BCIntegrate Class Referenceabstract

A class for handling numerical operations for models. More...

#include <BCIntegrate.h>

Inheritance diagram for BCIntegrate:
[legend]
Collaboration diagram for BCIntegrate:
[legend]

Public Types

Enumerators
enum  BCOptimizationMethod {
  kOptEmpty, kOptSimAnn, kOptMetropolis, kOptMinuit,
  kOptDefault, NOptMethods
}
 An enumerator for the mode finding algorithm. More...
 
enum  BCIntegrationMethod {
  kIntEmpty, kIntMonteCarlo, kIntCuba, kIntGrid,
  kIntLaplace, kIntDefault, NIntMethods
}
 An enumerator for integration algorithms. More...
 
enum  BCMarginalizationMethod {
  kMargEmpty, kMargMetropolis, kMargMonteCarlo, kMargGrid,
  kMargDefault, NMargMethods
}
 An enumerator for marginalization algorithms. More...
 
enum  BCSASchedule { kSACauchy, kSABoltzmann, kSACustom, NSAMethods }
 An enumerator for the Simulated Annealing schedule. More...
 
enum  BCCubaMethod {
  kCubaVegas, kCubaSuave, kCubaDivonne, kCubaCuhre,
  kCubaDefault, NCubaMethods
}
 An enumerator for Cuba integration methods. More...
 
Function pointer types
typedef void(BCIntegrate::* tRandomizer) (std::vector< double > &) const
 A pointer for a function that chooses a next random point.
 
typedef double(BCIntegrate::* tEvaluator) (std::vector< double > &, const std::vector< double > &, bool &)
 A pointer for a function that evaluates at a point.
 
typedef void(* tIntegralUpdater) (const std::vector< double > &, const int &, double &, double &)
 A pointer for a function that updates the integral and absolute precision.
 
- Public Types inherited from BCEngineMCMC
enum  Precision {
  kLow, kQuick, kMedium, kHigh,
  kVeryHigh
}
 An enumerator for the status of a test. More...
 
enum  Phase { kPreRun = -1, kUnsetPhase = 0, kMainRun = +1 }
 An enumerator for the phase of the Markov chain. More...
 
enum  InitialPositionScheme { kInitCenter = 0, kInitRandomUniform = 1, kInitUserDefined = 2, kInitRandomPrior = 3 }
 An enumerator for markov-chain position initialization. More...
 

Public Member Functions

Constructors and destructors
 BCIntegrate (const std::string &name="model")
 Default constructor.
 
 BCIntegrate (const std::string &filename, const std::string &name, bool loadObservables=true)
 Read in MCMC constructor. More...
 
 BCIntegrate (const BCIntegrate &other)
 Copy constructor.
 
BCIntegrateoperator= (const BCIntegrate &)
 Copy-assignment operator.
 
virtual ~BCIntegrate ()
 Destructor.
 
Member functions (get)
double GetIntegral () const
 
BCIntegrate::BCOptimizationMethod GetOptimizationMethod () const
 
BCIntegrate::BCIntegrationMethod GetIntegrationMethod () const
 
BCIntegrate::BCMarginalizationMethod GetMarginalizationMethod () const
 
BCIntegrate::BCSASchedule GetSASchedule () const
 
void GetRandomVectorInParameterSpace (std::vector< double > &x) const
 Fills a vector of random numbers x[i] between fMin[i] and fMax[i] into a vector. More...
 
double GetRandomPoint (std::vector< double > &x)
 Fills a vector of (flat) random numbers in the limits of the parameters and returns the probability at that point. More...
 
int GetNIterationsMin () const
 
int GetNIterationsMax () const
 
int GetNIterationsPrecisionCheck () const
 
int GetNIterations () const
 
double GetRelativePrecision () const
 
double GetAbsolutePrecision () const
 
BCCubaMethod GetCubaIntegrationMethod () const
 
const BCCubaOptions::VegasGetCubaVegasOptions () const
 
const BCCubaOptions::SuaveGetCubaSuaveOptions () const
 
const BCCubaOptions::DivonneGetCubaDivonneOptions () const
 
const BCCubaOptions::CuhreGetCubaCuhreOptions () const
 
TH1 * GetSlice (std::vector< unsigned > indices, unsigned &nIterations, double &log_max_val, const std::vector< double > parameters=std::vector< double >(0), int nbins=0, bool normalize=true)
 Returns a one-dimensional slice of the pdf at the point and along a specific direction. More...
 
TH1 * GetSlice (const std::string &name, unsigned &nIterations, double &log_max_val, const std::vector< double > parameters=std::vector< double >(0), int nbins=0, bool normalize=true)
 Returns a one-dimensional slice of the pdf at the point and along a specific direction. More...
 
TH1 * GetSlice (unsigned index, unsigned &nIterations, double &log_max_val, const std::vector< double > parameters=std::vector< double >(0), int nbins=0, bool normalize=true)
 Returns a one-dimensional slice of the pdf at the point and along a specific direction. More...
 
TH2 * GetSlice (const std::string &name1, const std::string &name2, unsigned &nIterations, double &log_max_val, const std::vector< double > parameters=std::vector< double >(0), int nbins=0, bool normalize=true)
 Returns a two-dimensional slice of the pdf at the point and along two specified directions. More...
 
TH2 * GetSlice (unsigned index1, unsigned index2, unsigned &nIterations, double &log_max_val, const std::vector< double > parameters=std::vector< double >(0), int nbins=0, bool normalize=true)
 Returns a two-dimensional slice of the pdf at the point and along two specified directions. More...
 
double GetError () const
 
TMinuitMinimizer & GetMinuit ()
 
double GetSAT0 () const
 Returns the Simulated Annealing starting temperature. More...
 
double GetSATmin () const
 Returns the Simulated Annealing threshhold temperature. More...
 
virtual const std::vector< double > & GetBestFitParameters () const
 
const std::vector< double > & GetBestFitParameterErrors () const
 Returns the set of errors on the values of the parameters at the mode.
 
double GetLogMaximum () const
 Returns the posterior at the mode. More...
 
Member functions (set)
void SetFlagIgnorePrevOptimization (bool flag)
 
void SetOptimizationMethod (BCIntegrate::BCOptimizationMethod method)
 
void SetIntegrationMethod (BCIntegrate::BCIntegrationMethod method)
 
void SetMarginalizationMethod (BCIntegrate::BCMarginalizationMethod method)
 
void SetSASchedule (BCIntegrate::BCSASchedule schedule)
 
void SetNIterationsMin (int niterations)
 
void SetNIterationsMax (int niterations)
 
void SetNIterationsPrecisionCheck (int niterations)
 
void SetRelativePrecision (double relprecision)
 
void SetAbsolutePrecision (double absprecision)
 Set absolute precision of the numerical integation.
 
void SetCubaIntegrationMethod (BCCubaMethod type)
 Set Cuba integration method.
 
void SetCubaOptions (const BCCubaOptions::Vegas &options)
 Set options for CUBA's Vegas. More...
 
void SetCubaOptions (const BCCubaOptions::Suave &options)
 Set options for CUBA's Suave. More...
 
void SetCubaOptions (const BCCubaOptions::Divonne &options)
 Set options for CUBA's Divonne. More...
 
void SetCubaOptions (const BCCubaOptions::Cuhre &options)
 Set options for CUBA's Cuhre. More...
 
void SetSAT0 (double T0)
 Set starting temperature for Simulated Annealing. More...
 
void SetSATmin (double Tmin)
 Set threshold temperature for Simulated Annealing. More...
 
- Public Member Functions inherited from BCEngineMCMC
 BCEngineMCMC (const std::string &name="model")
 Default constructor. More...
 
 BCEngineMCMC (const BCEngineMCMC &enginemcmc)
 Copy constructor. More...
 
 BCEngineMCMC (const std::string &filename, const std::string &name, bool loadObservables=true)
 Read in MCMC constructor. More...
 
BCEngineMCMCoperator= (const BCEngineMCMC &)
 Copy-assignment operator.
 
virtual ~BCEngineMCMC ()
 Destructor. More...
 
const std::string & GetName () const
 
const std::string & GetSafeName () const
 
unsigned GetNChains () const
 
unsigned GetNLag () const
 
int GetCurrentIteration () const
 
unsigned GetCurrentChain () const
 
int GetNIterationsConvergenceGlobal () const
 
unsigned GetNIterationsPreRun () const
 
unsigned GetNIterationsPreRunMin () const
 
unsigned GetNIterationsPreRunMax () const
 
unsigned GetNIterationsRun () const
 
unsigned GetNIterationsPreRunCheck () const
 
unsigned GetPreRunCheckClear ()
 
double GetMinimumEfficiency () const
 
double GetMaximumEfficiency () const
 
double GetScaleFactorLowerLimit () const
 
double GetScaleFactorUpperLimit () const
 
const std::vector< std::vector< double > > & GetScaleFactors () const
 
const ChainStateGetChainState (unsigned c) const
 
const std::vector< double > & Getx (unsigned c) const
 
double Getx (unsigned c, unsigned p) const
 
double GetLogProbx (unsigned c) const
 
BCEngineMCMC::Phase GetPhase () const
 
BCEngineMCMC::InitialPositionScheme GetInitialPositionScheme () const
 
unsigned GetInitialPositionAttemptLimit () const
 
bool GetProposeMultivariate () const
 
double GetProposalFunctionDof () const
 
unsigned GetMultivariateCovarianceUpdates () const
 
double GetMultivariateCovarianceUpdateLambda () const
 
double GetMultivariateEpsilon () const
 
double GetMultivariateScaleMultiplier () const
 
double GetRValueParametersCriterion () const
 
const std::vector< double > & GetRValueParameters () const
 
double GetRValueParameters (unsigned index) const
 
bool GetCorrectRValueForSamplingVariability () const
 Flag for correcting convergence checking for initial sampling variability. More...
 
bool GetFlagRun () const
 
TTree * GetMarkovChainTree () const
 Retrieve the tree containing the Markov chain. More...
 
TTree * GetParameterTree () const
 Retrieve the tree containing the parameter information. More...
 
TFile * GetOutputFile () const
 Retrieve output file for MCMC. More...
 
const BCEngineMCMC::StatisticsGetStatistics () const
 Get combined statistics for all chains. More...
 
const BCEngineMCMC::StatisticsGetStatistics (unsigned c) const
 Get MCMC statistics for one chain. More...
 
const std::vector< BCEngineMCMC::Statistics > & GetStatisticsVector () const
 Get vector of MCMC statistics for each chain separately. More...
 
bool GetRescaleHistogramRangesAfterPreRun () const
 
double GetHistogramRescalePadding () const
 
virtual std::vector< unsigned > GetH1DPrintOrder () const
 
virtual std::vector< std::pair< unsigned, unsigned > > GetH2DPrintOrder () const
 
bool MarginalizedHistogramExists (unsigned index) const
 
bool MarginalizedHistogramExists (unsigned index1, unsigned index2) const
 
TH1 * GetMarginalizedHistogram (const std::string &name) const
 Obtain the individual marginalized distributions with respect to one parameter as a ROOT TH1. More...
 
TH1 * GetMarginalizedHistogram (unsigned index) const
 Obtain the individual marginalized distributions with respect to one parameter as a ROOT TH1. More...
 
TH2 * GetMarginalizedHistogram (const std::string &name1, const std::string &name2) const
 Obtain the individual marginalized distributions with respect to two parameters as a ROOT TH2. More...
 
TH2 * GetMarginalizedHistogram (unsigned index1, unsigned index2) const
 Obtain the individual marginalized distributions with respect to two parameters as a ROOT TH2. More...
 
BCH1D GetMarginalized (const std::string &name) const
 Obtain the individual marginalized distributions with respect to one parameter. More...
 
BCH1D GetMarginalized (unsigned index) const
 Obtain the individual marginalized distributions with respect to one parameter. More...
 
BCH2D GetMarginalized (const std::string &name1, const std::string &name2) const
 Obtain the individual marginalized distributions with respect to two parameters. More...
 
BCH2D GetMarginalized (unsigned index1, unsigned index2) const
 Obtain the individual marginalized distributions with respect to two parameters. More...
 
unsigned GetMaximumParameterNameLength (bool observables=true) const
 
BCVariableGetVariable (unsigned index)
 
const BCVariableGetVariable (unsigned index) const
 
unsigned GetNVariables () const
 
BCParameterSetGetParameters ()
 
const BCParameterSetGetParameters () const
 
BCParameterGetParameter (unsigned index)
 
const BCParameterGetParameter (unsigned index) const
 
BCParameterGetParameter (const std::string &name)
 
const BCParameterGetParameter (const std::string &name) const
 
unsigned GetNParameters () const
 
unsigned GetNFixedParameters () const
 
unsigned GetNFreeParameters () const
 
BCObservableSetGetObservables ()
 
const BCObservableSetGetObservables () const
 
BCObservableGetObservable (unsigned index)
 
const BCObservableGetObservable (unsigned index) const
 
BCObservableGetObservable (const std::string &name)
 
const BCObservableGetObservable (const std::string &name) const
 
unsigned GetNObservables () const
 
const std::vector< double > & GetLocalModes (bool force_recalculation=false)
 
bool GetReuseObservables () const
 
BCH1DGetBCH1DdrawingOptions ()
 
BCH2DGetBCH2DdrawingOptions ()
 
void SetName (const std::string &name)
 Sets the name of the engine. More...
 
void SetScaleFactorLowerLimit (double l)
 Set scale factor lower limit.
 
void SetScaleFactorUpperLimit (double l)
 Set scale factor upper limit.
 
void SetInitialScaleFactors (const std::vector< double > &scale)
 Set the initial scale factors for the factorized proposal function. More...
 
void SetNChains (unsigned n)
 Sets the number of Markov chains which are run in parallel. More...
 
void SetNLag (unsigned n)
 Sets the lag of the Markov chains.
 
void SetNIterationsPreRunMax (unsigned n)
 Sets the maximum number of iterations in the pre-run. More...
 
void SetNIterationsRun (unsigned n)
 Sets the number of iterations. More...
 
void SetNIterationsPreRunMin (unsigned n)
 Sets the minimum number of iterations in the pre-run.
 
void SetNIterationsPreRunCheck (unsigned n)
 Sets the number of iterations between scale adjustments and convergence checks in the pre-run. More...
 
void SetPreRunCheckClear (unsigned n)
 Sets the number of prerun checks to make inbetween statistics clearing. More...
 
void SetMinimumEfficiency (double efficiency)
 Sets the minimum efficiency required for a chain. More...
 
void SetMaximumEfficiency (double efficiency)
 Sets the maximum efficiency required for a chain. More...
 
void SetRandomSeed (unsigned seed)
 Set the random number seed.
 
void SetInitialPositions (const std::vector< double > &x0s)
 Sets the initial positions for all chains. More...
 
void SetInitialPositions (const std::vector< std::vector< double > > &x0s)
 Sets the initial positions for all chains. More...
 
void SetInitialPositionScheme (BCEngineMCMC::InitialPositionScheme scheme)
 Sets flag which defines initial position. More...
 
void SetInitialPositionAttemptLimit (unsigned n)
 Sets maximum number of attempts to find a valid initial position. More...
 
void SetProposeMultivariate (bool flag)
 Set flag to true to turn on the multivariate proposal for MCMC based on (Haario et al., 2001) where the covariance is learned from the prerun. More...
 
void SetProposalFunctionDof (double dof=1)
 Set the degree of freedom of the proposal function for MCMC. More...
 
void SetMultivariateCovarianceUpdateLambda (double l)
 Set weighting for multivariate proposal function covariance update. More...
 
void SetMultivariateEpsilon (double epsilon)
 Sets multivariate-proposal-function cholesky-decomposition nudge. More...
 
void SetMultivariateScaleMultiplier (double s)
 Sets multivariate-proposal-function scale multiplier. More...
 
void SetFlagFillHistograms (bool flag)
 Sets whether to fill histograms. More...
 
void SetFlagFillHistograms (bool flag_1d, bool flag_2d)
 Sets the whether to fill histograms. More...
 
void SetFillHistogramParPar (unsigned x, unsigned y, bool flag=true)
 Sets whether to fill particular H2 histogram: par(y) vs. More...
 
void SetFillHistogramParPar (const std::string &x, const std::string &y, bool flag=true)
 Sets whether to fill particular H2 histogram: par(y) vs. More...
 
void SetFillHistogramParObs (unsigned x, unsigned y, bool flag=true)
 Sets whether to fill particular H2 histogram: obs(y) vs. More...
 
void SetFillHistogramParObs (const std::string &x, const std::string &y, bool flag=true)
 Sets whether to fill particular H2 histogram: obs(y) vs. More...
 
void SetFillHistogramObsObs (unsigned x, unsigned y, bool flag=true)
 Sets whether to fill particular H2 histogram: obs(y) vs. More...
 
void SetFillHistogramObsObs (const std::string &x, const std::string &y, bool flag=true)
 Sets whether to fill particular H2 histogram: obs(y) vs. More...
 
void SetFillHistogramObsPar (unsigned x, unsigned y, bool flag=true)
 Sets whether to fill particular H2 histogram: par(y) vs. More...
 
void SetFillHistogramObsPar (const std::string &x, const std::string &y, bool flag=true)
 Sets whether to fill particular H2 histogram: par(y) vs. More...
 
void SetFlagPreRun (bool flag)
 Set if a (new) prerun should be performed. More...
 
void SetRValueParametersCriterion (double r)
 Sets the parameter R-value criterion for convergence of all chains.
 
void SetCorrectRValueForSamplingVariability (bool flag=true)
 Set flag to correct convergence checking for initial sampling variability. More...
 
void SetPrecision (BCEngineMCMC::Precision precision)
 Set the precision for the MCMC run. More...
 
void SetPrecision (const BCEngineMCMC *other)
 Copy precision for the MCMC run from other model. More...
 
void SetPrecision (const BCEngineMCMC &other)
 Copy precision for the MCMC run from other model. More...
 
void SetNbins (unsigned int nbins)
 Set the number of bins for the marginalized distribution of all parameters. More...
 
void SetReuseObservables (bool flag)
 
void SetRescaleHistogramRangesAfterPreRun (bool flag=true)
 Set flag for rescaling histogram ranges after pre-run. More...
 
void SetHistogramRescalingPadding (double factor)
 Set enlargement factor of range for when rescaling. More...
 
void WriteMarkovChain (bool flag)
 Turn on/off writing of Markov chain to root file. More...
 
void WriteMarkovChainRun (bool flag)
 Turn on/off writing of Markov chain to root file during run. More...
 
void WriteMarkovChainPreRun (bool flag)
 Turn on/off writing of Markov chain to root file during prerun. More...
 
void WriteMarkovChain (const std::string &filename, const std::string &option, bool flag_run=true, bool flag_prerun=true)
 Turn on writing of Markov chain to root file. More...
 
void SetPriorConstant (unsigned index)
 
void SetPriorConstant (const std::string &name)
 
void SetPrior (unsigned index, TF1 &f, bool logL=true)
 
void SetPrior (const std::string &name, TF1 &f, bool logL=true)
 
void SetPriorDelta (unsigned index, double value)
 
void SetPriorDelta (const std::string &name, double value)
 
void SetPriorGauss (unsigned index, double mean, double sigma)
 
void SetPriorGauss (const std::string &name, double mean, double sigma)
 
void SetPriorGauss (unsigned index, double mode, double sigma_below, double sigma_above)
 
void SetPriorGauss (const std::string &name, double mode, double sigma_below, double sigma_above)
 
void SetPrior (unsigned index, TH1 &h, bool interpolate=false)
 
void SetPrior (const std::string &name, TH1 &h, bool interpolate=false)
 
void SetPriorConstantAll ()
 
void WriteMarginalizedDistributions (const std::string &filename, const std::string &option, bool closeExistingFile=false)
 Write marginalization histograms to file. More...
 
virtual void PrintSummary () const
 Prints a summary to the logs. More...
 
void PrintParameters (const std::vector< double > &P, void(*output)(const std::string &)=BCLog::OutSummary) const
 Print parameters. More...
 
unsigned PrintAllMarginalized (const std::string &filename, unsigned hdiv=1, unsigned vdiv=1) const
 Print all marginalizations. More...
 
unsigned PrintParameterPlot (const std::string &filename, int npar=10, double interval_content=68e-2, std::vector< double > quantile_values=std::vector< double >(0), bool rescale_ranges=true) const
 Print a summary plot for the parameters and user-defined observables. More...
 
bool DrawParameterPlot (unsigned i0, unsigned npar=0, double interval_content=68e-2, std::vector< double > quantile_values=std::vector< double >(0), bool rescale_ranges=true) const
 Draw a summary plot for the parameters in the range provided to current pad. More...
 
bool PrintCorrelationMatrix (const std::string &filename="matrix.pdf") const
 Print a correlation matrix for the parameters. More...
 
bool PrintCorrelationPlot (const std::string &filename="correlation.pdf", bool include_observables=true) const
 Print a correlation plot for the parameters. More...
 
bool PrintParameterLatex (const std::string &filename) const
 Print a LaTeX table of the parameters. More...
 
virtual void CreateHistograms (bool rescale_ranges=false)
 Create histograms from parameter and observable sets. More...
 
virtual void InitializeMarkovChainTree (bool replacetree=false, bool replacefile=false)
 Initialize the trees containing the Markov chains and parameter info. More...
 
virtual bool AddParameter (const std::string &name, double min, double max, const std::string &latexname="", const std::string &unitstring="")
 
virtual bool AddParameter (BCParameter &parameter)
 
virtual bool AddObservable (const std::string &name, double min, double max, const std::string &latexname="", const std::string &unitstring="")
 
virtual bool AddObservable (BCObservable &obs)
 
virtual void EvaluateObservables ()
 Evaluates user-defined observables at current state of all chains and stores results in fMCMCState.
 
virtual void EvaluateObservables (unsigned chain)
 Evaluates user-defined observables at current state of chain and stores results in fMCMCState. More...
 
virtual void CalculateObservables (const std::vector< double > &pars)
 Evaluates user-defined observables. More...
 
virtual double ProposalFunction (unsigned ichain, unsigned ipar)
 The default proposal function is a Breit-Wigner random walk. More...
 
bool GetProposalPointMetropolis (unsigned chain, std::vector< double > &x)
 Return a proposal point for the Metropolis algorithm. More...
 
bool GetProposalPointMetropolis (unsigned chain, unsigned parameter, std::vector< double > &x)
 Return a proposal point for the Metropolis algorithm. More...
 
bool GetNewPointMetropolis ()
 Generate a new point using the Metropolis algorithm for all chains. More...
 
bool GetNewPointMetropolis (unsigned chain)
 Generate a new point using the Metropolis algorithm for one chain. More...
 
bool GetNewPointMetropolis (unsigned chain, unsigned parameter)
 Generate a new point using the Metropolis algorithm for one chain, varying only one parameter's value. More...
 
bool AcceptOrRejectPoint (unsigned chain, unsigned parameter)
 Accept or rejects a point for a chain and updates efficiency. More...
 
void InChainFillHistograms (const ChainState &cs)
 Fill marginalized distributions from a chain state.
 
void InChainFillHistograms ()
 Fill marginalized distributions from all chain states.
 
void InChainFillTree (const ChainState &cs, unsigned chain_number)
 Write a chain state to the tree.
 
void InChainFillTree ()
 Write all chain states to the tree.
 
bool Metropolis ()
 Runs Metropolis algorithm. More...
 
bool MetropolisPreRun ()
 Runs a pre run for the Metropolis algorithm. More...
 
void MCMCInitialize ()
 Resets all containers used in MCMC and initializes starting points. More...
 
virtual void MCMCUserInitialize ()
 User hook called from MCMCInitialize(). More...
 
virtual void MCMCUserIterationInterface ()
 Interface allowing to execute arbitrary code for each iteration of the MCMC while running the chains after applying the lag but before filling histograms or the output tree. More...
 
virtual void MCMCCurrentPointInterface (const std::vector< double > &point, int ichain, bool accepted)
 Interface allowing to execute arbitrary code for each new point of the MCMC whether it is accepted or not. More...
 
void LoadParametersFromTree (TTree *partree, bool loadObservables=true)
 Load parameters and observables from tree. More...
 
void LoadMCMCParameters (TTree &partree)
 Load MCMC parameters from parameter tree: nchains, proposal function type, scales. More...
 
virtual bool ParameterTreeMatchesModel (TTree *partree, bool checkObservables=true)
 Check parameter tree against model. More...
 
void LoadMCMC (const std::string &filename, std::string mcmcTreeName="", std::string parameterTreeName="", bool loadObservables=true)
 Load previous MCMC run. More...
 
void LoadMCMC (TTree *mcmcTree, TTree *parTree, bool loadObservables=true)
 Load previous MCMC run. More...
 
bool ValidMCMCTree (TTree *tree, bool checkObservables=true) const
 Check tree structure for MCMC tree. More...
 
bool ValidParameterTree (TTree *tree) const
 Check tree structure for parameter tree. More...
 
void CloseOutputFile ()
 Close the root output file. More...
 
virtual void Remarginalize (bool autorange=true)
 Marginalize from TTree. More...
 
void PrepareToContinueMarginalization (const std::string &filename, const std::string &mcmcTreeName="", const std::string &parameterTreeName="", bool loadObservables=true, bool autorange=true)
 Continue the marginalization already stored in another file. More...
 
virtual bool UpdateMultivariateProposalFunctionCovariances (double a)
 Update multivariate proposal function covariances. More...
 
virtual bool UpdateMultivariateProposalFunctionCovariances ()
 Update multivariate proposal function covariances. More...
 
void CalculateCholeskyDecompositions ()
 Calculate Cholesky decompositions needed for multivariate proposal function. More...
 
void UpdateChainIndex (int chain)
 Keep track of which chain is currently computed (within a thread). More...
 

Protected Member Functions

virtual std::string GetBestFitSummary (unsigned i) const
 Get string summarizing best fit for single variable. More...
 
unsigned IntegrationOutputFrequency () const
 Determine frequency of output during integration.
 
void LogOutputAtEndOfIntegration (double integral, double absprecision, double relprecision, int nIterations)
 Helper method to output at end of integration. More...
 
void LogOutputAtIntegrationStatusUpdate (BCIntegrationMethod type, double integral, double absprecision, int nIterations)
 Helper method to output integration status. More...
 
void LogOutputAtStartOfIntegration (BCIntegrationMethod type, BCCubaMethod cubatype)
 Helper method to output at beginning of integration. More...
 
virtual void PrintBestFitSummary () const
 Print best fit to log.
 
virtual void PrintMarginalizationSummary () const
 Print marginalization to log. More...
 
- Protected Member Functions inherited from BCEngineMCMC
virtual void PrintModelSummary () const
 Print model summary to log. More...
 
void SetFillHistogram (int x, int y, bool flag)
 Set whether to fill 2D histogram y vs x: positive indices for parameters; negative for observables, starting at -1 and going more negative—observable index = -(index+1). More...
 
unsigned UpdateFrequency (unsigned N) const
 return appropriate update interval More...
 
void UpdateParameterTree ()
 Update Paramater TTree with scales and efficiencies. More...
 

Protected Attributes

bool fFlagIgnorePrevOptimization
 Flag for ignoring older results of optimization.
 
bool fFlagMarginalized
 flag indicating if the model was marginalized
 
double fSALogProb
 Log probability of current simulated annealing iteration. More...
 
int fSANIterations
 Number of iterations for simualted annealing. More...
 
double fSAT0
 Starting temperature for Simulated Annealing.
 
double fSATemperature
 Current temperature of simulated annealing algorithm. More...
 
double fSATmin
 Minimal/Threshold temperature for Simulated Annealing.
 
std::vector< double > fSAx
 Current simulated annealing parameter point. More...
 
- Protected Attributes inherited from BCEngineMCMC
BCH1D fBCH1DdrawingOptions
 A BCH1D (with no histogram) for storing BCH1D drawing options. More...
 
BCH2D fBCH2DdrawingOptions
 A BCH2D (with no histogram) for storing BCH2D drawing options. More...
 
bool fCorrectRValueForSamplingVariability
 flag for correcting R value for initial sampling variability. More...
 
std::vector< TH1 * > fH1Marginalized
 Vector of 1D marginalized distributions.
 
std::vector< std::vector< TH2 * > > fH2Marginalized
 Vector of 2D marginalized distributions. More...
 
double fHistogramRescalePadding
 factor for enlarging range of histograms when rescaling. More...
 
unsigned fInitialPositionAttemptLimit
 Maximum number of attempts to make to set the initial position. More...
 
BCEngineMCMC::InitialPositionScheme fInitialPositionScheme
 Variable which defines the initial position. More...
 
std::vector< double > fLocalModes
 Vector of local modes. More...
 
int fMCMCCurrentIteration
 The current iteration number. More...
 
double fMCMCEfficiencyMax
 The maximum allowed efficiency for MCMC.
 
double fMCMCEfficiencyMin
 The minimum required efficiency for MCMC.
 
bool fMCMCFlagWriteChainToFile
 Flag to write Markov chains to file.
 
bool fMCMCFlagWritePreRunToFile
 Flag to write pre run to file.
 
std::vector< std::vector< double > > fMCMCInitialPosition
 The intial position of each Markov chain. More...
 
std::vector< double > fMCMCInitialScaleFactors
 User-provided initial values of the scale factors of the factorized proposal function. More...
 
unsigned fMCMCNChains
 Number of Markov chains ran in parallel.
 
int fMCMCNIterationsConvergenceGlobal
 Number of iterations needed for all chains to convergence simultaneously.
 
unsigned fMCMCNIterationsPreRunCheck
 Number of iterations between scale adjustments and convergence checks in pre-run. More...
 
unsigned fMCMCNIterationsPreRunMax
 Maximum number of iterations for a Markov chain prerun.
 
unsigned fMCMCNIterationsPreRunMin
 Minimum number of iterations for the pre-run.
 
unsigned fMCMCNIterationsRun
 Number of iterations for a Markov chain run.
 
unsigned fMCMCNLag
 The lag for the Markov Chain.
 
TFile * fMCMCOutputFile
 Output file for for writing MCMC Tree. More...
 
std::string fMCMCOutputFilename
 Output filename for for writing MCMC Tree. More...
 
std::string fMCMCOutputFileOption
 Output file open option for for writing MCMC Tree. More...
 
BCEngineMCMC::Phase fMCMCPhase
 The phase of the run. More...
 
unsigned fMCMCPreRunCheckClear
 Number of iterations between clearing of convergence stats in pre-run. More...
 
double fMCMCProposalFunctionDof
 Degree of freedom of Student's t proposal. More...
 
std::vector< std::vector< double > > fMCMCProposalFunctionScaleFactor
 Scale factors for proposal functions. More...
 
bool fMCMCProposeMultivariate
 Flag for using multivariate proposal function. More...
 
std::vector< double > fMCMCRValueParameters
 The R-values for each parameter.
 
double fMCMCRValueParametersCriterion
 The R-value criterion for convergence of parameters.
 
double fMCMCScaleFactorLowerLimit
 Lower limit for scale factors.
 
double fMCMCScaleFactorUpperLimit
 Upper limit for scale factors.
 
std::vector< ChainStatefMCMCStates
 The current states of each Markov chain. More...
 
std::vector< BCEngineMCMC::StatisticsfMCMCStatistics
 Statistics for each Markov chain. More...
 
BCEngineMCMC::Statistics fMCMCStatistics_AllChains
 Statistics across all Markov chains. More...
 
TTree * fMCMCTree
 The tree containing the Markov chains. More...
 
unsigned int fMCMCTree_Chain
 Chain number for storing into tree.
 
ChainState fMCMCTree_State
 MC state object for storing into tree.
 
bool fMCMCTreeLoaded
 flag for whether MCMC Tree successfully loaded. More...
 
bool fMCMCTreeReuseObservables
 flag for whether to reuse MCMC Tree's observables. More...
 
double fMultivariateCovarianceUpdateLambda
 weighting parameter for multivariate-proposal-function covariance update. More...
 
unsigned fMultivariateCovarianceUpdates
 Number of multivariate-proposal-function covariance updates performed. More...
 
double fMultivariateEpsilon
 multivariate-proposal-function cholesky-decomposition nudge. More...
 
std::vector< TMatrixD > fMultivariateProposalFunctionCholeskyDecomposition
 Cholesky decompositions for multivariate proposal function. More...
 
std::vector< TMatrixDSym > fMultivariateProposalFunctionCovariance
 Covariance matrices used in multivariate proposal functions. More...
 
double fMultivariateScaleMultiplier
 factor to multiply or divide scale factors by in adjusting multivariate-proposal-function scales. More...
 
std::string fName
 Name of the engine. More...
 
BCAux::BCTrash< TObject > fObjectTrash
 Storage for plot objects with proper clean-up.
 
BCObservableSet fObservables
 User-calculated Observables Set.
 
BCParameterSet fParameters
 Parameter settings.
 
TTree * fParameterTree
 The tree containing the parameter information. More...
 
TRandom3 fRandom
 Random number generator.
 
std::vector< std::pair< int, int > > fRequestedH2
 Vector of pairs of indices for which 2D histograms should be stored. More...
 
bool fRescaleHistogramRangesAfterPreRun
 flag for rescaling of histograms after pre-run. More...
 
std::string fSafeName
 Safe name of the engine for use in naming ROOT objects. More...
 

Member functions (miscellaneous methods)

virtual double Eval (const std::vector< double > &x)=0
 Evaluate the unnormalized probability at a point in parameter space. More...
 
virtual double LogEval (const std::vector< double > &x)
 Evaluate the natural logarithm of the Eval function. More...
 
double Normalize ()
 Performs integration. More...
 
double Integrate (BCIntegrationMethod intmethod)
 Does the integration over the un-normalized probability. More...
 
double Integrate ()
 Perform the integration. More...
 
double Integrate (BCIntegrationMethod type, tRandomizer randomizer, tEvaluator evaluator, tIntegralUpdater updater, std::vector< double > &sums)
 Does the integration over the un-normalized probability. More...
 
double EvaluatorMC (std::vector< double > &sums, const std::vector< double > &point, bool &accepted)
 Evaluates integrator.
 
int MarginalizeAll ()
 Marginalize all probabilities wrt. More...
 
int MarginalizeAll (BCMarginalizationMethod margmethod)
 Marginalize all probabilities wrt. More...
 
virtual void MarginalizePreprocess ()
 Method executed for before marginalization. More...
 
virtual void MarginalizePostprocess ()
 Method executed after marginalization. More...
 
std::vector< double > FindMode (std::vector< double > start=std::vector< double >())
 Do the mode finding using a method set via SetOptimizationMethod. More...
 
std::vector< double > FindMode (BCIntegrate::BCOptimizationMethod optmethod, std::vector< double > start=std::vector< double >())
 Find mode using a specific method. More...
 
double SATemperature (double t) const
 Temperature annealing schedule for use with Simulated Annealing. More...
 
double SATemperatureBoltzmann (double t) const
 Temperature annealing schedule for use with Simulated Annealing. More...
 
double SATemperatureCauchy (double t) const
 Temperature annealing schedule for use with Simulated Annealing. More...
 
virtual double SATemperatureCustom (double t) const
 Temperature annealing schedule for use with Simulated Annealing. More...
 
std::vector< double > GetProposalPointSA (const std::vector< double > &x, int t) const
 Generates a new state in a neighbourhood around x that is to be accepted or rejected by the Simulated Annealing algorithm. More...
 
std::vector< double > GetProposalPointSABoltzmann (const std::vector< double > &x, int t) const
 Generates a new state in a neighbourhood around x that is to be accepted or rejected by the Simulated Annealing algorithm. More...
 
std::vector< double > GetProposalPointSACauchy (const std::vector< double > &x, int t) const
 Generates a new state in a neighbourhood around x that is to be accepted or rejected by the Simulated Annealing algorithm. More...
 
virtual std::vector< double > GetProposalPointSACustom (const std::vector< double > &x, int t) const
 Generates a new state in a neighbourhood around x that is to be accepted or rejected by the Simulated Annealing algorithm. More...
 
std::vector< double > SAHelperGetRandomPointOnHypersphere () const
 Generates a uniform distributed random point on the surface of a fNvar-dimensional Hypersphere. More...
 
double SAHelperGetRadialCauchy () const
 Generates the radial part of a n-dimensional Cauchy distribution. More...
 
double SAHelperSinusToNIntegral (int dim, double theta) const
 Returns the Integral of sin^dim from 0 to theta. More...
 
virtual void ResetResults ()
 Reset all information on the best-fit parameters. More...
 
std::string DumpIntegrationMethod (BCIntegrationMethod type) const
 Return string with the name for a given integration type. More...
 
std::string DumpCurrentIntegrationMethod () const
 Return string with the name for the currently set integration type. More...
 
std::string DumpUsedIntegrationMethod () const
 Return string with the name for the previously used integration type. More...
 
std::string DumpMarginalizationMethod (BCMarginalizationMethod type) const
 Return string with the name for a given marginalization type. More...
 
std::string DumpCurrentMarginalizationMethod () const
 Return string with the name for the currently set marginalization type. More...
 
std::string DumpUsedMarginalizationMethod () const
 Return string with the name for the marginalization type used. More...
 
std::string DumpOptimizationMethod (BCOptimizationMethod type) const
 Return string with the name for a given optimization type. More...
 
std::string DumpCurrentOptimizationMethod () const
 Return string with the name for the currently set optimization type. More...
 
std::string DumpUsedOptimizationMethod () const
 Return string with the name for the optimization type used to find the current mode. More...
 
std::string DumpCubaIntegrationMethod (BCCubaMethod type) const
 Return string with the name for a given Cuba integration type. More...
 
std::string DumpCubaIntegrationMethod () const
 Return string with the name for the currently set Cuba integration type. More...
 
void SetBestFitParameters (const std::vector< double > &x)
 Set best fit parameters values. More...
 
void SetBestFitParameters (const std::vector< double > &x, const double &new_value, double &old_value)
 
bool CheckMarginalizationAvailability (BCMarginalizationMethod type)
 Check availability of integration routine for marginalization. More...
 
bool CheckMarginalizationIndices (TH1 *hist, const std::vector< unsigned > &index)
 Check that indices of parameters to marginalize w/r/t are correct.
 
double IntegrateLaplace ()
 Integrate using the Laplace approximation. More...
 
static void IntegralUpdaterMC (const std::vector< double > &sums, const int &nIterations, double &integral, double &absprecision)
 Updates info about integrator.
 

Additional Inherited Members

- Static Public Member Functions inherited from BCEngineMCMC
static double RValue (const std::vector< double > &means, const std::vector< double > &variances, unsigned n, bool correctForSamplingVariability=true)
 Calculate R value of set of batches of samples—represented by their means and variances, all batches containing the same number of samples—according to Brooks & Gelman, "General Methods for Monitoring Convergence of Iterative Simulations," (1988) More...
 

Detailed Description

A class for handling numerical operations for models.

Author
Daniel Kollar
Kevin Kröninger
Version
1.0
Date
08.2008

This is a base class for a model class. It contains numerical methods to carry out the integration, marginalization, peak finding etc.

Definition at line 143 of file BCIntegrate.h.

Member Enumeration Documentation

An enumerator for Cuba integration methods.

Enumerator
kCubaVegas 

Vegas.

kCubaSuave 

Suave.

kCubaDivonne 

Divonne.

kCubaCuhre 

Cuhre.

kCubaDefault 

Default.

NCubaMethods 

number of available CUBA methods

Definition at line 196 of file BCIntegrate.h.

An enumerator for integration algorithms.

Enumerator
kIntEmpty 

No integration method set.

kIntMonteCarlo 

Sample mean method.

kIntCuba 

Use CUBA interface.

kIntGrid 

Integration by gridding of parameter space.

kIntLaplace 

Laplace approximation.

kIntDefault 

Default.

NIntMethods 

number of available integration methods

Definition at line 164 of file BCIntegrate.h.

An enumerator for marginalization algorithms.

Enumerator
kMargEmpty 

No marginalization method set.

kMargMetropolis 

Metropolis Hastings.

kMargMonteCarlo 

Sample mean Monte Carlo.

kMargGrid 

Marginalization by gridding of parameter space.

kMargDefault 

Default.

NMargMethods 

number of available marginalization methods

Definition at line 176 of file BCIntegrate.h.

An enumerator for the mode finding algorithm.

Enumerator
kOptEmpty 

No optimization method set.

kOptSimAnn 

Simulated annealing.

kOptMetropolis 

Metropolis Hastings.

kOptMinuit 

ROOT's Minuit.

kOptDefault 

Default.

NOptMethods 

number of available optimization methods

Definition at line 153 of file BCIntegrate.h.

An enumerator for the Simulated Annealing schedule.

Enumerator
kSACauchy 

Cauchy scheduler.

kSABoltzmann 

Boltzman scheduler.

kSACustom 

Custom scheduler.

NSAMethods 

number of available schedulers

Definition at line 187 of file BCIntegrate.h.

Constructor & Destructor Documentation

BCIntegrate::BCIntegrate ( const std::string &  filename,
const std::string &  name,
bool  loadObservables = true 
)

Read in MCMC constructor.

Parameters
filenamePath of file holding model.
nameName of model (file should contain TTree's [name]_mcmc and [name]_parameters.
if empty string is given, properly matching TTrees are searched for in the file.
loadObservablesFlag for whether to load observables for file (true; default) or to let user respecify observables.

Definition at line 166 of file BCIntegrate.cxx.

Member Function Documentation

bool BCIntegrate::CheckMarginalizationAvailability ( BCMarginalizationMethod  type)

Check availability of integration routine for marginalization.

Returns
availability of marginalization method

Definition at line 582 of file BCIntegrate.cxx.

std::string BCIntegrate::DumpCubaIntegrationMethod ( BCIntegrate::BCCubaMethod  type) const

Return string with the name for a given Cuba integration type.

Parameters
typecode for the Cuba integration type
Returns
string containing the name of the Cuba integration type

Definition at line 1928 of file BCIntegrate.cxx.

std::string BCIntegrate::DumpCubaIntegrationMethod ( ) const
inline

Return string with the name for the currently set Cuba integration type.

Returns
string containing the name of the Cuba integration type

Definition at line 783 of file BCIntegrate.h.

std::string BCIntegrate::DumpCurrentIntegrationMethod ( ) const
inline

Return string with the name for the currently set integration type.

Returns
string containing the name of the integration type

Definition at line 729 of file BCIntegrate.h.

std::string BCIntegrate::DumpCurrentMarginalizationMethod ( ) const
inline

Return string with the name for the currently set marginalization type.

Returns
string containing the name of the marginalization type

Definition at line 747 of file BCIntegrate.h.

std::string BCIntegrate::DumpCurrentOptimizationMethod ( ) const
inline

Return string with the name for the currently set optimization type.

Returns
string containing the name of the optimization type

Definition at line 765 of file BCIntegrate.h.

std::string BCIntegrate::DumpIntegrationMethod ( BCIntegrate::BCIntegrationMethod  type) const

Return string with the name for a given integration type.

Parameters
typecode for the integration type
Returns
string containing the name of the integration type

Definition at line 1871 of file BCIntegrate.cxx.

std::string BCIntegrate::DumpMarginalizationMethod ( BCIntegrate::BCMarginalizationMethod  type) const

Return string with the name for a given marginalization type.

Parameters
typecode for the marginalization type
Returns
string containing the name of the marginalization type

Definition at line 1890 of file BCIntegrate.cxx.

std::string BCIntegrate::DumpOptimizationMethod ( BCIntegrate::BCOptimizationMethod  type) const

Return string with the name for a given optimization type.

Parameters
typecode for the optimization type
Returns
string containing the name of the optimization type

Definition at line 1909 of file BCIntegrate.cxx.

std::string BCIntegrate::DumpUsedIntegrationMethod ( ) const
inline

Return string with the name for the previously used integration type.

Returns
string containing the name of the integration type

Definition at line 735 of file BCIntegrate.h.

std::string BCIntegrate::DumpUsedMarginalizationMethod ( ) const
inline

Return string with the name for the marginalization type used.

Returns
string containing the name of the marginalization type

Definition at line 753 of file BCIntegrate.h.

std::string BCIntegrate::DumpUsedOptimizationMethod ( ) const
inline

Return string with the name for the optimization type used to find the current mode.

Returns
string containing the name of the optimization type

Definition at line 771 of file BCIntegrate.h.

virtual double BCIntegrate::Eval ( const std::vector< double > &  x)
pure virtual

Evaluate the unnormalized probability at a point in parameter space.

Method needs to be overloaded by the user.

Parameters
xThe point in parameter space
Returns
The unnormalized probability

Implemented in BCModel.

std::vector< double > BCIntegrate::FindMode ( std::vector< double >  start = std::vector<double>())

Do the mode finding using a method set via SetOptimizationMethod.

Default is Minuit. The mode can be extracted using the GetBestFitParameters() method.

A starting point for the mode finding can be specified for Minuit. If not specified, the previously found maximum (typically from marginalization) is used as an initial point. If that is not available, then the Minuit default will be used (center of the parameter space).

Returns
The mode found.
Note
The result may not coincide with the result of GetBestFitParameters() if a previous optimization found a better value.

Definition at line 1022 of file BCIntegrate.cxx.

std::vector< double > BCIntegrate::FindMode ( BCIntegrate::BCOptimizationMethod  optmethod,
std::vector< double >  start = std::vector<double>() 
)

Find mode using a specific method.

The original method will be reset.

Parameters
optmethodthe optimization method
startthe starting point for the optimization algorithm
Returns
the mode
See also
std::vector<double> FindMode(std::vector<double> start = std::vector<double>(0));

Definition at line 1090 of file BCIntegrate.cxx.

double BCIntegrate::GetAbsolutePrecision ( ) const
inline
Returns
Requested absolute precision of the numerical integration

Definition at line 321 of file BCIntegrate.h.

const std::vector< double > & BCIntegrate::GetBestFitParameters ( ) const
virtual
Returns
vector of parameter and observable values at the mode.

Reimplemented from BCEngineMCMC.

Definition at line 274 of file BCIntegrate.cxx.

std::string BCIntegrate::GetBestFitSummary ( unsigned  i) const
protectedvirtual

Get string summarizing best fit for single variable.

Parameters
iIndex of variable to summarize.

Reimplemented from BCEngineMCMC.

Definition at line 2032 of file BCIntegrate.cxx.

const BCCubaOptions::Cuhre& BCIntegrate::GetCubaCuhreOptions ( ) const
inline
Returns
Options used for integration with CUBA's Cuhre

Definition at line 346 of file BCIntegrate.h.

const BCCubaOptions::Divonne& BCIntegrate::GetCubaDivonneOptions ( ) const
inline
Returns
Options used for integration with CUBA's Divonne

Definition at line 341 of file BCIntegrate.h.

BCCubaMethod BCIntegrate::GetCubaIntegrationMethod ( ) const
inline
Returns
Cuba Integration method

Definition at line 326 of file BCIntegrate.h.

const BCCubaOptions::Suave& BCIntegrate::GetCubaSuaveOptions ( ) const
inline
Returns
Options used for integration with CUBA's Suave

Definition at line 336 of file BCIntegrate.h.

const BCCubaOptions::Vegas& BCIntegrate::GetCubaVegasOptions ( ) const
inline
Returns
Options used for integration with CUBA's Vegas

Definition at line 331 of file BCIntegrate.h.

double BCIntegrate::GetError ( ) const
inline
Returns
The uncertainty in the most recent Monte Carlo integration

Definition at line 411 of file BCIntegrate.h.

double BCIntegrate::GetIntegral ( ) const
inline
Returns
The integral.

Definition at line 259 of file BCIntegrate.h.

BCIntegrate::BCIntegrationMethod BCIntegrate::GetIntegrationMethod ( ) const
inline
Returns
The current integration method

Definition at line 269 of file BCIntegrate.h.

double BCIntegrate::GetLogMaximum ( ) const
inlinevirtual

Returns the posterior at the mode.

Returns
the posterior.

Reimplemented from BCEngineMCMC.

Definition at line 442 of file BCIntegrate.h.

BCIntegrate::BCMarginalizationMethod BCIntegrate::GetMarginalizationMethod ( ) const
inline
Returns
The current marginalization method

Definition at line 274 of file BCIntegrate.h.

TMinuitMinimizer& BCIntegrate::GetMinuit ( )
inline
Returns
Minuit minimizer used for mode finding. The object is in a valid but unusable state if no mode finding has been performed.

Definition at line 418 of file BCIntegrate.h.

int BCIntegrate::GetNIterations ( ) const
inline
Returns
The number of iterations for the most recent Monte Carlo integration

Definition at line 311 of file BCIntegrate.h.

int BCIntegrate::GetNIterationsMax ( ) const
inline
Returns
The number of maximum iterations for integration

Definition at line 301 of file BCIntegrate.h.

int BCIntegrate::GetNIterationsMin ( ) const
inline
Returns
The number of minimum iterations for integration

Definition at line 296 of file BCIntegrate.h.

int BCIntegrate::GetNIterationsPrecisionCheck ( ) const
inline
Returns
The interval for checking precision in integration

Definition at line 306 of file BCIntegrate.h.

BCIntegrate::BCOptimizationMethod BCIntegrate::GetOptimizationMethod ( ) const
inline
Returns
The current optimization method

Definition at line 264 of file BCIntegrate.h.

std::vector< double > BCIntegrate::GetProposalPointSA ( const std::vector< double > &  x,
int  t 
) const

Generates a new state in a neighbourhood around x that is to be accepted or rejected by the Simulated Annealing algorithm.

Delegates the generation to the appropriate method according to fSASchedule.

Parameters
xlast state.
ttime iterator to determine current temperature.

Definition at line 1282 of file BCIntegrate.cxx.

std::vector< double > BCIntegrate::GetProposalPointSABoltzmann ( const std::vector< double > &  x,
int  t 
) const

Generates a new state in a neighbourhood around x that is to be accepted or rejected by the Simulated Annealing algorithm.

This method is used for Boltzmann annealing schedule.

Parameters
xlast state.
ttime iterator to determine current temperature.

Definition at line 1294 of file BCIntegrate.cxx.

std::vector< double > BCIntegrate::GetProposalPointSACauchy ( const std::vector< double > &  x,
int  t 
) const

Generates a new state in a neighbourhood around x that is to be accepted or rejected by the Simulated Annealing algorithm.

This method is used for Cauchy annealing schedule.

Parameters
xlast state.
ttime iterator to determine current temperature.

Definition at line 1314 of file BCIntegrate.cxx.

std::vector< double > BCIntegrate::GetProposalPointSACustom ( const std::vector< double > &  x,
int  t 
) const
virtual

Generates a new state in a neighbourhood around x that is to be accepted or rejected by the Simulated Annealing algorithm.

This is a virtual method to be overridden by a user-defined custom proposal function.

Parameters
xlast state.
ttime iterator to determine current temperature.

Definition at line 1356 of file BCIntegrate.cxx.

double BCIntegrate::GetRandomPoint ( std::vector< double > &  x)

Fills a vector of (flat) random numbers in the limits of the parameters and returns the probability at that point.

Parameters
xA vector of doubles to fill
Returns
The (unnormalized) probability at the random point

Definition at line 1006 of file BCIntegrate.cxx.

void BCIntegrate::GetRandomVectorInParameterSpace ( std::vector< double > &  x) const

Fills a vector of random numbers x[i] between fMin[i] and fMax[i] into a vector.

Parameters
xA vector of doubles to fill

Definition at line 1013 of file BCIntegrate.cxx.

double BCIntegrate::GetRelativePrecision ( ) const
inline
Returns
Requested relative precision of the numerical integation

Definition at line 316 of file BCIntegrate.h.

BCIntegrate::BCSASchedule BCIntegrate::GetSASchedule ( ) const
inline
Returns
The Simulated Annealing schedule

Definition at line 279 of file BCIntegrate.h.

double BCIntegrate::GetSAT0 ( ) const
inline

Returns the Simulated Annealing starting temperature.

Definition at line 423 of file BCIntegrate.h.

double BCIntegrate::GetSATmin ( ) const
inline

Returns the Simulated Annealing threshhold temperature.

Definition at line 428 of file BCIntegrate.h.

TH1 * BCIntegrate::GetSlice ( std::vector< unsigned >  indices,
unsigned &  nIterations,
double &  log_max_val,
const std::vector< double >  parameters = std::vector<double>(0),
int  nbins = 0,
bool  normalize = true 
)

Returns a one-dimensional slice of the pdf at the point and along a specific direction.

Parameters
indicesThe indices of the model parameters along which the slice is calculated.
nIterationsAdd the number of posterior evaluations performed.
log_max_valStores the log of the maximum value before normalizing
parametersThe point at which the other parameters are fixed.
nbinsThe number of bins of the 1D-histogram.
normalizeFlag for turning on normalization of histogram.
Returns
The slice histogram.

Definition at line 840 of file BCIntegrate.cxx.

TH1* BCIntegrate::GetSlice ( const std::string &  name,
unsigned &  nIterations,
double &  log_max_val,
const std::vector< double >  parameters = std::vector<double>(0),
int  nbins = 0,
bool  normalize = true 
)
inline

Returns a one-dimensional slice of the pdf at the point and along a specific direction.

Parameters
nameThe name of the model parameter along which the slice is calculated.
nIterationsAdd the number of posterior evaluations performed.
log_max_valStores the log of the maximum value before normalizing
parametersThe point at which the other parameters are fixed.
nbinsThe number of bins of the 1D-histogram.
normalizeFlag for turning on normalization of histogram.
Returns
The 1D slice.

Definition at line 369 of file BCIntegrate.h.

TH1* BCIntegrate::GetSlice ( unsigned  index,
unsigned &  nIterations,
double &  log_max_val,
const std::vector< double >  parameters = std::vector<double>(0),
int  nbins = 0,
bool  normalize = true 
)
inline

Returns a one-dimensional slice of the pdf at the point and along a specific direction.

Parameters
indexThe index of the model parameter along which the slice is calculated.
nIterationsAdd the number of posterior evaluations performed.
log_max_valStores the log of the maximum value before normalizing
parametersThe point at which the other parameters are fixed.
nbinsThe number of bins of the 1D-histogram.
normalizeFlag for turning on normalization of histogram.
Returns
The 1D slice.

Definition at line 381 of file BCIntegrate.h.

TH2* BCIntegrate::GetSlice ( const std::string &  name1,
const std::string &  name2,
unsigned &  nIterations,
double &  log_max_val,
const std::vector< double >  parameters = std::vector<double>(0),
int  nbins = 0,
bool  normalize = true 
)
inline

Returns a two-dimensional slice of the pdf at the point and along two specified directions.

Parameters
name1The name of the first model parameter along which the slice is calculated.
name2The name of the second model parameter along which the slice is calculated.
nIterationsAdd the number of posterior evaluations performed.
log_max_valStores the log of the maximum value before normalizing
parametersThe point at which the other parameters are fixed.
nbinsThe number of bins on each axis of the 2D-histogram.
normalizeFlag for turning on normalization of histogram.
Returns
The 2D slice.

Definition at line 394 of file BCIntegrate.h.

TH2 * BCIntegrate::GetSlice ( unsigned  index1,
unsigned  index2,
unsigned &  nIterations,
double &  log_max_val,
const std::vector< double >  parameters = std::vector<double>(0),
int  nbins = 0,
bool  normalize = true 
)

Returns a two-dimensional slice of the pdf at the point and along two specified directions.

Parameters
index1The index of the first model parameter along which the slice is calculated.
index2The index of the second model parameter along which the slice is calculated.
nIterationsAdd the number of posterior evaluations performed.
log_max_valStores the log of the maximum value before normalizing
parametersThe point at which the other parameters are fixed.
nbinsThe number of bins on each axis of the 2D-histogram.
normalizeFlag for turning on normalization of histogram.
Returns
The 2D slice.

Definition at line 998 of file BCIntegrate.cxx.

double BCIntegrate::Integrate ( BCIntegrationMethod  intmethod)

Does the integration over the un-normalized probability.

Parameters
intmethodThe integration method to used
Returns
The normalization value

Definition at line 385 of file BCIntegrate.cxx.

double BCIntegrate::Integrate ( )

Perform the integration.

Returns
the integral

Definition at line 306 of file BCIntegrate.cxx.

double BCIntegrate::Integrate ( BCIntegrationMethod  type,
tRandomizer  randomizer,
tEvaluator  evaluator,
tIntegralUpdater  updater,
std::vector< double > &  sums 
)

Does the integration over the un-normalized probability.

Parameters
typeThe integration method to used (for printing out status updates by name)
randomizerPointer to function to choose next random point
evaluatorPointer to function to evaluate point
updaterPointer to function to update integral and precision
sumsVector of doubles holding values used in integral calculation
Returns
The integral value

Definition at line 490 of file BCIntegrate.cxx.

double BCIntegrate::IntegrateLaplace ( )

Integrate using the Laplace approximation.

Take the result of a previous successful minuit run to estimate the covariance matrix. Else it runs minuit (again).

Returns
the integral on the log(!) scale. Note that is is different from the other integration methods that work on the linear scale.

Definition at line 1782 of file BCIntegrate.cxx.

virtual double BCIntegrate::LogEval ( const std::vector< double > &  x)
inlinevirtual

Evaluate the natural logarithm of the Eval function.

For better numerical stability, this method should also be overloaded by the user.

Parameters
xThe point in parameter space
Returns
log(Eval(x))

Implements BCEngineMCMC.

Reimplemented in BCModel.

Definition at line 555 of file BCIntegrate.h.

void BCIntegrate::LogOutputAtEndOfIntegration ( double  integral,
double  absprecision,
double  relprecision,
int  nIterations 
)
protected

Helper method to output at end of integration.

Definition at line 475 of file BCIntegrate.cxx.

void BCIntegrate::LogOutputAtIntegrationStatusUpdate ( BCIntegrationMethod  type,
double  integral,
double  absprecision,
int  nIterations 
)
protected

Helper method to output integration status.

Definition at line 484 of file BCIntegrate.cxx.

void BCIntegrate::LogOutputAtStartOfIntegration ( BCIntegrationMethod  type,
BCCubaMethod  cubatype 
)
protected

Helper method to output at beginning of integration.

Definition at line 420 of file BCIntegrate.cxx.

int BCIntegrate::MarginalizeAll ( )

Marginalize all probabilities wrt.

single parameters and all combinations of two parameters. The individual distributions can be retrieved using the GetMarginalized method.

Returns
Total number of marginalized distributions

Definition at line 635 of file BCIntegrate.cxx.

int BCIntegrate::MarginalizeAll ( BCIntegrate::BCMarginalizationMethod  margmethod)

Marginalize all probabilities wrt.

single parameters and all combinations of two parameters. The individual distributions can be retrieved using the GetMarginalized method.

Parameters
margmethodthe marginalization method.
Returns
Total number of marginalized distributions

Definition at line 821 of file BCIntegrate.cxx.

virtual void BCIntegrate::MarginalizePostprocess ( )
inlinevirtual

Method executed after marginalization.

User's code should be provided via overloading in the derived class

Definition at line 616 of file BCIntegrate.h.

virtual void BCIntegrate::MarginalizePreprocess ( )
inlinevirtual

Method executed for before marginalization.

User's code should be provided via overloading in the derived class

Reimplemented in BCFitter.

Definition at line 610 of file BCIntegrate.h.

double BCIntegrate::Normalize ( )
inline

Performs integration.

Definition at line 560 of file BCIntegrate.h.

void BCIntegrate::PrintMarginalizationSummary ( ) const
protectedvirtual

Print marginalization to log.

Reimplemented from BCEngineMCMC.

Definition at line 1997 of file BCIntegrate.cxx.

void BCIntegrate::ResetResults ( )
virtual

Reset all information on the best-fit parameters.

Reimplemented from BCEngineMCMC.

Definition at line 292 of file BCIntegrate.cxx.

double BCIntegrate::SAHelperGetRadialCauchy ( ) const

Generates the radial part of a n-dimensional Cauchy distribution.

Helper function for Cauchy annealing.

Definition at line 1404 of file BCIntegrate.cxx.

std::vector< double > BCIntegrate::SAHelperGetRandomPointOnHypersphere ( ) const

Generates a uniform distributed random point on the surface of a fNvar-dimensional Hypersphere.

Used as a helper to generate proposal points for Cauchy annealing.

Definition at line 1363 of file BCIntegrate.cxx.

double BCIntegrate::SAHelperSinusToNIntegral ( int  dim,
double  theta 
) const

Returns the Integral of sin^dim from 0 to theta.

Helper function needed for generating Cauchy distributions.

Definition at line 1460 of file BCIntegrate.cxx.

double BCIntegrate::SATemperature ( double  t) const

Temperature annealing schedule for use with Simulated Annealing.

Delegates to the appropriate method according to fSASchedule.

Parameters
titerator for lowering the temperature over time.

Definition at line 1251 of file BCIntegrate.cxx.

double BCIntegrate::SATemperatureBoltzmann ( double  t) const

Temperature annealing schedule for use with Simulated Annealing.

This method is used for Boltzmann annealing schedule.

Parameters
titerator for lowering the temperature over time.

Definition at line 1263 of file BCIntegrate.cxx.

double BCIntegrate::SATemperatureCauchy ( double  t) const

Temperature annealing schedule for use with Simulated Annealing.

This method is used for Cauchy annealing schedule.

Parameters
titerator for lowering the temperature over time.

Definition at line 1269 of file BCIntegrate.cxx.

double BCIntegrate::SATemperatureCustom ( double  t) const
virtual

Temperature annealing schedule for use with Simulated Annealing.

This is a virtual method to be overridden by a user-defined custom temperature schedule.

Parameters
titerator for lowering the temperature over time.

Definition at line 1275 of file BCIntegrate.cxx.

void BCIntegrate::SetBestFitParameters ( const std::vector< double > &  x)

Set best fit parameters values.

Parameters
xParameter set to designate as best fit parameters.

Definition at line 414 of file BCIntegrate.cxx.

void BCIntegrate::SetBestFitParameters ( const std::vector< double > &  x,
const double &  new_value,
double &  old_value 
)
Parameters
xSet best-fit parameters to x if new_value >= old_value
new_valueis the value of the function at x
old_valueis the old best fit value, updated to new_value, if it is larger

Definition at line 405 of file BCIntegrate.cxx.

void BCIntegrate::SetCubaOptions ( const BCCubaOptions::Vegas options)
inline

Set options for CUBA's Vegas.

Parameters
optionsOptions for CUBA

Definition at line 506 of file BCIntegrate.h.

void BCIntegrate::SetCubaOptions ( const BCCubaOptions::Suave options)
inline

Set options for CUBA's Suave.

Parameters
optionsOptions for CUBA

Definition at line 512 of file BCIntegrate.h.

void BCIntegrate::SetCubaOptions ( const BCCubaOptions::Divonne options)
inline

Set options for CUBA's Divonne.

Parameters
optionsOptions for CUBA

Definition at line 518 of file BCIntegrate.h.

void BCIntegrate::SetCubaOptions ( const BCCubaOptions::Cuhre options)
inline

Set options for CUBA's Cuhre.

Parameters
optionsOptions for CUBA

Definition at line 524 of file BCIntegrate.h.

void BCIntegrate::SetFlagIgnorePrevOptimization ( bool  flag)
inline
Parameters
flagFlag whether or not to ignore result of previous mode finding

Definition at line 451 of file BCIntegrate.h.

void BCIntegrate::SetIntegrationMethod ( BCIntegrate::BCIntegrationMethod  method)
Parameters
methodThe current integration method

Definition at line 1488 of file BCIntegrate.cxx.

void BCIntegrate::SetMarginalizationMethod ( BCIntegrate::BCMarginalizationMethod  method)
inline
Parameters
methodThe current marginalization method

Definition at line 465 of file BCIntegrate.h.

void BCIntegrate::SetNIterationsMax ( int  niterations)
inline
Parameters
niterationsThe maximum number of iterations for integration

Definition at line 480 of file BCIntegrate.h.

void BCIntegrate::SetNIterationsMin ( int  niterations)
inline
Parameters
niterationsThe maximum number of iterations for integration

Definition at line 475 of file BCIntegrate.h.

void BCIntegrate::SetNIterationsPrecisionCheck ( int  niterations)
inline
Parameters
niterationsinterval for checking precision in integration routines

Definition at line 485 of file BCIntegrate.h.

void BCIntegrate::SetOptimizationMethod ( BCIntegrate::BCOptimizationMethod  method)
inline
Parameters
methodThe current optimization method

Definition at line 456 of file BCIntegrate.h.

void BCIntegrate::SetRelativePrecision ( double  relprecision)
inline
Parameters
relprecisionThe relative precision envisioned for Monte Carlo integration

Definition at line 491 of file BCIntegrate.h.

void BCIntegrate::SetSASchedule ( BCIntegrate::BCSASchedule  schedule)
inline
Parameters
scheduleThe Simulated Annealing schedule

Definition at line 470 of file BCIntegrate.h.

void BCIntegrate::SetSAT0 ( double  T0)
inline

Set starting temperature for Simulated Annealing.

Parameters
T0starting temperature.

Definition at line 530 of file BCIntegrate.h.

void BCIntegrate::SetSATmin ( double  Tmin)
inline

Set threshold temperature for Simulated Annealing.

Parameters
Tminthreshold temperature.

Definition at line 536 of file BCIntegrate.h.

Member Data Documentation

double BCIntegrate::fSALogProb
protected

Log probability of current simulated annealing iteration.

Definition at line 858 of file BCIntegrate.h.

int BCIntegrate::fSANIterations
protected

Number of iterations for simualted annealing.

Definition at line 850 of file BCIntegrate.h.

double BCIntegrate::fSATemperature
protected

Current temperature of simulated annealing algorithm.

Definition at line 854 of file BCIntegrate.h.

std::vector<double> BCIntegrate::fSAx
protected

Current simulated annealing parameter point.

Definition at line 862 of file BCIntegrate.h.


The documentation for this class was generated from the following files: