An engine class for Markov Chain Monte Carlo. More...
#include <BCEngineMCMC.h>
Classes | |
struct | ChainState |
A struct for holding a state in a Markov chain. More... | |
struct | Statistics |
A struct for holding statistical information about samples. More... | |
Public Types | |
Enumerators | |
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 destructor | |
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... | |
BCEngineMCMC & | operator= (const BCEngineMCMC &) |
Copy-assignment operator. | |
virtual | ~BCEngineMCMC () |
Destructor. More... | |
Getters | |
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 ChainState & | GetChainState (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::Statistics & | GetStatistics () const |
Get combined statistics for all chains. More... | |
const BCEngineMCMC::Statistics & | GetStatistics (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 |
BCVariable & | GetVariable (unsigned index) |
const BCVariable & | GetVariable (unsigned index) const |
unsigned | GetNVariables () const |
BCParameterSet & | GetParameters () |
const BCParameterSet & | GetParameters () const |
BCParameter & | GetParameter (unsigned index) |
const BCParameter & | GetParameter (unsigned index) const |
BCParameter & | GetParameter (const std::string &name) |
const BCParameter & | GetParameter (const std::string &name) const |
unsigned | GetNParameters () const |
unsigned | GetNFixedParameters () const |
unsigned | GetNFreeParameters () const |
BCObservableSet & | GetObservables () |
const BCObservableSet & | GetObservables () const |
BCObservable & | GetObservable (unsigned index) |
const BCObservable & | GetObservable (unsigned index) const |
BCObservable & | GetObservable (const std::string &name) |
const BCObservable & | GetObservable (const std::string &name) const |
unsigned | GetNObservables () const |
virtual const std::vector< double > & | GetBestFitParameters () const |
virtual const std::vector< double > & | GetBestFitParameterErrors () const |
const std::vector< double > & | GetLocalModes (bool force_recalculation=false) |
virtual double | GetLogMaximum () const |
bool | GetReuseObservables () const |
BCH1D & | GetBCH1DdrawingOptions () |
BCH2D & | GetBCH2DdrawingOptions () |
Setters | |
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... | |
Prior setting functions (all deprecated). | |
The priors are not used for the likelihood calculation by BCEngineMCMC, but are used for initializing the positions of the chains. | |
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 () |
Output functions | |
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... | |
Protected Member Functions | |
virtual std::string | GetBestFitSummary (unsigned i) const |
Get string summarizing best fit for single variable. More... | |
virtual void | PrintBestFitSummary () const |
Print best fit to log. More... | |
virtual void | PrintMarginalizationSummary () const |
Print marginalization to log. More... | |
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 | |
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< ChainState > | fMCMCStates |
The current states of each Markov chain. More... | |
std::vector< BCEngineMCMC::Statistics > | fMCMCStatistics |
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... | |
Friends | |
Structs | |
class | BCModel |
Miscellaneous methods | |
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 ¶meter) |
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. | |
virtual double | LogEval (const std::vector< double > ¶meters)=0 |
Needs to be overloaded in the derived class. More... | |
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 | ResetResults () |
Reset the MCMC variables. 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 ¶meterTreeName="", 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... | |
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
An engine class for Markov Chain Monte Carlo.
- Version
- 1.0
- Date
- 08.2008
This class represents an engine class for Markov Chain Monte Carlo (MCMC). One or more chains can be defined simultaneously.
Definition at line 55 of file BCEngineMCMC.h.
Member Enumeration Documentation
An enumerator for markov-chain position initialization.
Definition at line 82 of file BCEngineMCMC.h.
enum BCEngineMCMC::Phase |
An enumerator for the phase of the Markov chain.
Negative values are in pre-run, 0 is unset, positive is main run.
Enumerator | |
---|---|
kPreRun |
In pre-run. |
kUnsetPhase |
Unset. |
kMainRun |
In main run. |
Definition at line 75 of file BCEngineMCMC.h.
An enumerator for the status of a test.
See BCEngineMCMC::SetPrecision for full details
Definition at line 65 of file BCEngineMCMC.h.
Constructor & Destructor Documentation
BCEngineMCMC::BCEngineMCMC | ( | const std::string & | name = "model" | ) |
Default constructor.
Definition at line 55 of file BCEngineMCMC.cxx.
BCEngineMCMC::BCEngineMCMC | ( | const BCEngineMCMC & | enginemcmc | ) |
Copy constructor.
Definition at line 126 of file BCEngineMCMC.cxx.
BCEngineMCMC::BCEngineMCMC | ( | const std::string & | filename, |
const std::string & | name, | ||
bool | loadObservables = true |
||
) |
Read in MCMC constructor.
- Parameters
-
filename Path of file holding model. name Name 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.loadObservables Flag for whether to load observables for file (true; default) or to let user respecify observables.
Definition at line 90 of file BCEngineMCMC.cxx.
|
virtual |
Destructor.
Definition at line 301 of file BCEngineMCMC.cxx.
Member Function Documentation
bool BCEngineMCMC::AcceptOrRejectPoint | ( | unsigned | chain, |
unsigned | parameter | ||
) |
Accept or rejects a point for a chain and updates efficiency.
- Parameters
-
chain chain index parameter index of parameter to update efficiency for
- Returns
- Whether proposed point was accepted (true) or previous point was kept (false).
Definition at line 1580 of file BCEngineMCMC.cxx.
|
inlinevirtual |
- Deprecated:
- Instead use GetObservables().Add(...) Adds a user-calculated observable.
- Parameters
-
name name of observable min minimum value of the observable max maximum value of the observable latexname Optional latexname used for plotting unitstring Unit string to be printed for observable.
- Returns
- Success of action.
Definition at line 1266 of file BCEngineMCMC.h.
|
inlinevirtual |
- Deprecated:
- Instead use GetObservables().Add(obs) Adds a user-calculated observable to the model.
- Parameters
-
obs A user-calculated observable
- Returns
- Success of action.
Definition at line 1274 of file BCEngineMCMC.h.
|
inlinevirtual |
- Deprecated:
- Instead use GetParameters().Add(...) Adds a parameter.
- Parameters
-
name Name of parameter min minimum value of the parameter max maximum value of the parameter latexname Optional latexname used for plotting unitstring Unit string to be printed for parameter.
- Returns
- Success of action.
Definition at line 1246 of file BCEngineMCMC.h.
|
inlinevirtual |
- Deprecated:
- Instead use GetParameters().Add(parameter) Adds a parameter to the model.
- Parameters
-
parameter A model parameter
- Returns
- Success of action.
Definition at line 1254 of file BCEngineMCMC.h.
void BCEngineMCMC::CalculateCholeskyDecompositions | ( | ) |
Calculate Cholesky decompositions needed for multivariate proposal function.
Definition at line 1471 of file BCEngineMCMC.cxx.
|
inlinevirtual |
Evaluates user-defined observables.
To be overloaded by user to calculate user observables.
- Parameters
-
pars Parameter set to evaluate observables for.
Reimplemented in BCPriorModel.
Definition at line 1292 of file BCEngineMCMC.h.
void BCEngineMCMC::CloseOutputFile | ( | ) |
Close the root output file.
Definition at line 1744 of file BCEngineMCMC.cxx.
|
virtual |
Create histograms from parameter and observable sets.
- Parameters
-
rescale_ranges Rescale axis ranges to range reached by MCMC tree if true
Definition at line 2596 of file BCEngineMCMC.cxx.
bool BCEngineMCMC::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.
- Parameters
-
i0 Index of first parameter to print. npar Number of parameters to print, set to 0 to print all. interval_content Probability mass to display in smallest X interval band quantile_values Vector of quantile values to draw rescale_ranges Flag for rescaling to range surveyed by MCMC chains
- Returns
- Success of action.
Definition at line 3026 of file BCEngineMCMC.cxx.
|
virtual |
Evaluates user-defined observables at current state of chain and stores results in fMCMCState.
- Parameters
-
chain Chain to evaluate observables for.
Definition at line 2352 of file BCEngineMCMC.cxx.
|
inline |
Definition at line 739 of file BCEngineMCMC.h.
|
inline |
Definition at line 744 of file BCEngineMCMC.h.
|
virtual |
- Returns
- vector of standard errors (sqrt(variance)) of 1D marginals
Reimplemented in BCIntegrate.
Definition at line 639 of file BCEngineMCMC.cxx.
|
virtual |
- Returns
- vector of parameter values at global mode.
Reimplemented in BCIntegrate.
Definition at line 633 of file BCEngineMCMC.cxx.
|
protectedvirtual |
Get string summarizing best fit for single variable.
- Parameters
-
i Index of variable to summarize.
Reimplemented in BCIntegrate.
Definition at line 2792 of file BCEngineMCMC.cxx.
|
inline |
- Parameters
-
c index of the Markov chain
- Returns
- state of chains
Definition at line 364 of file BCEngineMCMC.h.
|
inline |
Flag for correcting convergence checking for initial sampling variability.
Definition at line 451 of file BCEngineMCMC.h.
unsigned BCEngineMCMC::GetCurrentChain | ( | ) | const |
- Returns
- current chain index. If no chain is currently selected (for example, not within mcmc), return 0.
Definition at line 606 of file BCEngineMCMC.cxx.
|
inline |
- Returns
- current iterations
Definition at line 289 of file BCEngineMCMC.h.
|
inline |
- Returns
- if MCMC has been performed or not
Definition at line 456 of file BCEngineMCMC.h.
|
virtual |
- Returns
- vector indices for order of printing 1D histograms.
Definition at line 2891 of file BCEngineMCMC.cxx.
|
virtual |
- Returns
- vector of index pairs for order of printing 2D histograms.
Definition at line 2902 of file BCEngineMCMC.cxx.
|
inline |
- Returns
- Factor by which to enlarge histogram ranges when rescaling to add padding beyond range.
Definition at line 497 of file BCEngineMCMC.h.
|
inline |
- Returns
- maximum number of allowed attempts to set the initial position
Definition at line 398 of file BCEngineMCMC.h.
|
inline |
- Returns
- flag which defined initial position
Definition at line 393 of file BCEngineMCMC.h.
const std::vector< double > & BCEngineMCMC::GetLocalModes | ( | bool | force_recalculation = false | ) |
- Returns
- vector of the local modes of parameters and observables
- Parameters
-
force_recalculation flag for forcing recalculation of local modes from histograms.
Definition at line 648 of file BCEngineMCMC.cxx.
|
inlinevirtual |
- Returns
- The log of the value at the mode.
Reimplemented in BCIntegrate.
Definition at line 729 of file BCEngineMCMC.h.
|
inline |
- Returns
- log of the probability of the current points of the Markov chain.
- Parameters
-
c chain index
Definition at line 383 of file BCEngineMCMC.h.
|
inline |
Obtain the individual marginalized distributions with respect to one parameter.
- Note
- The most efficient method is to access by index.
- Parameters
-
name The parameter's name
- Returns
- 1D marginalized probability
Definition at line 561 of file BCEngineMCMC.h.
BCH1D BCEngineMCMC::GetMarginalized | ( | unsigned | index | ) | const |
Obtain the individual marginalized distributions with respect to one parameter.
- Parameters
-
index The parameter index
- Returns
- 1D marginalized probability
Definition at line 573 of file BCEngineMCMC.cxx.
|
inline |
Obtain the individual marginalized distributions with respect to two parameters.
- Note
- The most efficient method is to access by indices.
- Parameters
-
name1 Name of first parameter name2 Name of second parameter
- Returns
- 2D marginalized probability
Definition at line 578 of file BCEngineMCMC.h.
BCH2D BCEngineMCMC::GetMarginalized | ( | unsigned | index1, |
unsigned | index2 | ||
) | const |
Obtain the individual marginalized distributions with respect to two parameters.
- Parameters
-
index1 Index of first parameter index2 Index of second parameter
- Returns
- 2D marginalized probability
Definition at line 587 of file BCEngineMCMC.cxx.
|
inline |
Obtain the individual marginalized distributions with respect to one parameter as a ROOT TH1.
- Note
- The most efficient method is to access by index.
- Parameters
-
name The parameter's name
- Returns
- 1D marginalized probability
Definition at line 527 of file BCEngineMCMC.h.
TH1 * BCEngineMCMC::GetMarginalizedHistogram | ( | unsigned | index | ) | const |
Obtain the individual marginalized distributions with respect to one parameter as a ROOT TH1.
- Parameters
-
index The parameter index
- Returns
- 1D marginalized probability
Definition at line 528 of file BCEngineMCMC.cxx.
|
inline |
Obtain the individual marginalized distributions with respect to two parameters as a ROOT TH2.
- Note
- The most efficient method is to access by indices.
- Parameters
-
name1 Name of first parameter name2 Name of second parameter
- Returns
- 2D marginalized probability
Definition at line 544 of file BCEngineMCMC.h.
TH2 * BCEngineMCMC::GetMarginalizedHistogram | ( | unsigned | index1, |
unsigned | index2 | ||
) | const |
Obtain the individual marginalized distributions with respect to two parameters as a ROOT TH2.
- Parameters
-
index1 Index of first parameter index2 Index of second parameter
- Returns
- 2D marginalized probability
Definition at line 546 of file BCEngineMCMC.cxx.
|
inline |
Retrieve the tree containing the Markov chain.
Definition at line 461 of file BCEngineMCMC.h.
|
inline |
- Returns
- maximum efficiency required for a chain.
Definition at line 343 of file BCEngineMCMC.h.
|
inline |
- Parameters
-
observables Whether to check max length of user-defined observable names.
- Returns
- length of longest parameter name.
Definition at line 592 of file BCEngineMCMC.h.
|
inline |
- Returns
- minimum efficiency required for a chain.
Definition at line 338 of file BCEngineMCMC.h.
|
inline |
- Returns
- number of updates to multivariate-proposal-function covariance performed.
Definition at line 416 of file BCEngineMCMC.h.
|
inline |
- Returns
- multivariate-proposal-function Cholesky decomposition nudge size.
Definition at line 426 of file BCEngineMCMC.h.
|
inline |
- Returns
- multivariate-proposal-function scale multiplier.
Definition at line 431 of file BCEngineMCMC.h.
|
inline |
- Returns
- The name of the engine.
Definition at line 269 of file BCEngineMCMC.h.
|
inline |
- Returns
- number of Markov chains
Definition at line 279 of file BCEngineMCMC.h.
bool BCEngineMCMC::GetNewPointMetropolis | ( | ) |
Generate a new point using the Metropolis algorithm for all chains.
- Returns
- Whether proposed point was accepted (true) or previous point was kept (false).
Definition at line 1647 of file BCEngineMCMC.cxx.
bool BCEngineMCMC::GetNewPointMetropolis | ( | unsigned | chain | ) |
Generate a new point using the Metropolis algorithm for one chain.
- Parameters
-
chain chain index
- Returns
- Whether proposed point was accepted (true) or previous point was kept (false).
Definition at line 1632 of file BCEngineMCMC.cxx.
bool BCEngineMCMC::GetNewPointMetropolis | ( | unsigned | chain, |
unsigned | parameter | ||
) |
Generate a new point using the Metropolis algorithm for one chain, varying only one parameter's value.
- Parameters
-
chain Chain index parameter Index of single parameter to update.
- Returns
- Whether proposed point was accepted (true) or previous point was kept (false).
Definition at line 1617 of file BCEngineMCMC.cxx.
|
inline |
- Deprecated:
- Instead call GetParameters().GetNFixedParameters()
- Returns
- The number of fixed parameters.
Definition at line 662 of file BCEngineMCMC.h.
|
inline |
- Deprecated:
- Instead call GetParameters().GetNFreeParameters()
- Returns
- The number of free parameters.
Definition at line 668 of file BCEngineMCMC.h.
|
inline |
- Returns
- number of iterations needed for all chains to converge simultaneously. A value of -1 indicates that the chains did not converge during a prerun. A value of 0 indicates that chain criteria were set by hand and convergence was not checked.
Definition at line 301 of file BCEngineMCMC.h.
unsigned BCEngineMCMC::GetNIterationsPreRun | ( | ) | const |
- Returns
- the number of iterations the prerun actually took. If the prerun wasn't run, that's 0. Else it is bounded by GetNIterationsPreRunMin() and GetNIterationsPreRunMax().
Definition at line 622 of file BCEngineMCMC.cxx.
|
inline |
- Returns
- number of iterations between scale adjustments and convergence checking during pre-run.
Definition at line 328 of file BCEngineMCMC.h.
|
inline |
- Returns
- maximum number of pre-run iterations for a Markov chain
Definition at line 318 of file BCEngineMCMC.h.
|
inline |
- Returns
- minimum number of pre-run iterations for a Markov chain
Definition at line 313 of file BCEngineMCMC.h.
|
inline |
- Returns
- number of iterations for a Markov chain
Definition at line 323 of file BCEngineMCMC.h.
|
inline |
- Returns
- lag of the Markov chains
Definition at line 284 of file BCEngineMCMC.h.
|
inline |
- Returns
- The number of user-defined observables.
Definition at line 711 of file BCEngineMCMC.h.
|
inline |
- Returns
- The number of parameters of the model.
Definition at line 656 of file BCEngineMCMC.h.
|
inline |
- Returns
- The number of parameters of the model.
Definition at line 613 of file BCEngineMCMC.h.
|
inline |
- Deprecated:
- Instead call GetObservables().At(index)
- Parameters
-
index The index of the observable in the observable set.
- Returns
- The user-defined observable.
Definition at line 685 of file BCEngineMCMC.h.
|
inline |
- Deprecated:
- Instead call GetObservbles().At(index)
- Parameters
-
index The index of the observable in the observable set.
- Returns
- The user-defined observable.
Definition at line 692 of file BCEngineMCMC.h.
|
inline |
- Deprecated:
- Instead call GetObservables().Get(name)
- Parameters
-
name The name of the observable in the observable set.
- Returns
- The user-defined observable.
Definition at line 699 of file BCEngineMCMC.h.
|
inline |
- Deprecated:
- Instead call GetObservables().Get(name)
- Parameters
-
name The name of the observable in the observable set.
- Returns
- The user-defined observable.
Definition at line 706 of file BCEngineMCMC.h.
|
inline |
- Returns
- Observable set.
Definition at line 673 of file BCEngineMCMC.h.
|
inline |
- Returns
- Observable set.
Definition at line 678 of file BCEngineMCMC.h.
|
inline |
Retrieve output file for MCMC.
Definition at line 471 of file BCEngineMCMC.h.
|
inline |
- Deprecated:
- Instead call GetParameters().At(index)
- Parameters
-
index The index of the parameter in the parameter set.
- Returns
- The parameter.
Definition at line 630 of file BCEngineMCMC.h.
|
inline |
- Deprecated:
- Instead call GetParameters().At(index)
- Parameters
-
index The index of the parameter in the parameter set.
- Returns
- The parameter.
Definition at line 637 of file BCEngineMCMC.h.
|
inline |
- Deprecated:
- Instead call GetParameters().Get(name)
- Parameters
-
name The name of the parameter in the parameter set.
- Returns
- The parameter.
Definition at line 644 of file BCEngineMCMC.h.
|
inline |
- Deprecated:
- Instead call GetParameters().Get(name)
- Parameters
-
name The name of the parameter in the parameter set.
- Returns
- The parameter.
Definition at line 651 of file BCEngineMCMC.h.
|
inline |
- Returns
- Parameter set.
Definition at line 618 of file BCEngineMCMC.h.
|
inline |
- Returns
- Parameter set.
Definition at line 623 of file BCEngineMCMC.h.
|
inline |
Retrieve the tree containing the parameter information.
Definition at line 466 of file BCEngineMCMC.h.
|
inline |
- Returns
- pointer to the phase of a run.
Definition at line 388 of file BCEngineMCMC.h.
|
inline |
- Returns
- number of iterations between clearings of statistics for convergence checking.
Definition at line 333 of file BCEngineMCMC.h.
|
inline |
- Returns
- Degree of freedom of multivariate proposal function. Anything >0 is the degree of freedom of a multivariate Student's t distribution, <= 0 is a multivariate Gaussian
Definition at line 411 of file BCEngineMCMC.h.
bool BCEngineMCMC::GetProposalPointMetropolis | ( | unsigned | chain, |
std::vector< double > & | x | ||
) |
Return a proposal point for the Metropolis algorithm.
- Parameters
-
chain chain index x proposal point
- Returns
- flag indicating whether the new point lies within the allowed range
Definition at line 1528 of file BCEngineMCMC.cxx.
bool BCEngineMCMC::GetProposalPointMetropolis | ( | unsigned | chain, |
unsigned | parameter, | ||
std::vector< double > & | x | ||
) |
Return a proposal point for the Metropolis algorithm.
- Parameters
-
chain chain index parameter Index of single parameter to update. x proposal point
- Returns
- flag indicating whether the new point lies within the allowed range
Definition at line 1557 of file BCEngineMCMC.cxx.
|
inline |
- Returns
- whether to use a multivariate proposal function.
Definition at line 403 of file BCEngineMCMC.h.
|
inline |
- Returns
- Flag for whether to rescale histogram ranges to fit MCMC reach after pre-run.
Definition at line 492 of file BCEngineMCMC.h.
|
inline |
- Returns
- Flag whether to reuse user-defined observables from MCMC tree when looping through it.
Definition at line 734 of file BCEngineMCMC.h.
|
inline |
- Returns
- vector of R values for parameters
Definition at line 441 of file BCEngineMCMC.h.
|
inline |
- Returns
- R-value for a parameter
- Parameters
-
index parameter index
Definition at line 447 of file BCEngineMCMC.h.
|
inline |
- Returns
- R-value criterion for parameters
Definition at line 436 of file BCEngineMCMC.h.
|
inline |
- Returns
- The name of the engine with spaces removed.
Definition at line 274 of file BCEngineMCMC.h.
|
inline |
- Returns
- proposal function scale factor lower limit
Definition at line 348 of file BCEngineMCMC.h.
|
inline |
- Returns
- scale factor for all parameters and chains
Definition at line 358 of file BCEngineMCMC.h.
|
inline |
- Returns
- proposal function scale factor upper limit
Definition at line 353 of file BCEngineMCMC.h.
|
inline |
Get combined statistics for all chains.
Definition at line 476 of file BCEngineMCMC.h.
|
inline |
Get MCMC statistics for one chain.
- Parameters
-
c Chain to get statistics of.
Definition at line 482 of file BCEngineMCMC.h.
|
inline |
Get vector of MCMC statistics for each chain separately.
Definition at line 487 of file BCEngineMCMC.h.
|
inline |
- Parameters
-
index The index of the variable running first over 0,...,N_parameters in the ParameterSet, and then over N_parameters,...,N_parameters+N_observables in the ObservableSet
- Returns
- The variable.
Definition at line 600 of file BCEngineMCMC.h.
|
inline |
- Parameters
-
index The index of the observable running first over 0,...,N_parameters in the ParameterSet, and then over N_parameters,...,N_parameters+N_observables in the ObservableSet
- Returns
- The observable.
Definition at line 608 of file BCEngineMCMC.h.
|
inline |
- Parameters
-
c index of the Markov chain
- Returns
- current point of the Markov chain
Definition at line 370 of file BCEngineMCMC.h.
|
inline |
- Parameters
-
c chain index p parameter index
- Returns
- parameter of the Markov chain
Definition at line 377 of file BCEngineMCMC.h.
|
virtual |
Initialize the trees containing the Markov chains and parameter info.
- Parameters
-
replacetree Flag to delete and recreate tree object if already existing. replacefile Flag to delete and recreate file object if already existing.
Reimplemented in BCModel.
Definition at line 703 of file BCEngineMCMC.cxx.
void BCEngineMCMC::LoadMCMC | ( | const std::string & | filename, |
std::string | mcmcTreeName = "" , |
||
std::string | parameterTreeName = "" , |
||
bool | loadObservables = true |
||
) |
Load previous MCMC run.
- Parameters
-
filename Pathname of file containing model trees. mcmcTreeName Name of tree inside file containing MCMC, empty string (default) loads [modelname]_mcmc. parameterTreeName Name of tree inside file containing parameter list, empty string (default) loads [modelname]_parameters. loadObservables Flag for whether to load observables from parameter list and MCMC trees.
Definition at line 1138 of file BCEngineMCMC.cxx.
void BCEngineMCMC::LoadMCMC | ( | TTree * | mcmcTree, |
TTree * | parTree, | ||
bool | loadObservables = true |
||
) |
Load previous MCMC run.
- Parameters
-
mcmcTree Tree containing MCMC samples. parTree Tree containing definition of parameters. loadObservables Flag for whether to load observables.
Definition at line 1218 of file BCEngineMCMC.cxx.
void BCEngineMCMC::LoadMCMCParameters | ( | TTree & | partree | ) |
Load MCMC parameters from parameter tree: nchains, proposal function type, scales.
- Parameters
-
partree Tree holding parameter information.
Definition at line 999 of file BCEngineMCMC.cxx.
void BCEngineMCMC::LoadParametersFromTree | ( | TTree * | partree, |
bool | loadObservables = true |
||
) |
Load parameters and observables from tree.
- Parameters
-
partree Tree holding parameter information. loadObservables Flag for whether to also load observables.
Definition at line 918 of file BCEngineMCMC.cxx.
|
pure virtual |
Needs to be overloaded in the derived class.
- Parameters
-
parameters Parameter set to evaluate at.
- Returns
- natural logarithm of the function to map with MCMC
Implemented in BCIntegrate, and BCModel.
|
inline |
- Parameters
-
index Index of histogram of which to check existence
- Returns
- Whether the marginalized histogram exists.
Definition at line 511 of file BCEngineMCMC.h.
|
inline |
- Parameters
-
index1 X Index of histogram of which to check existence. index2 Y Index of histogram of which to check existence.
- Returns
- Whether the marginalized histogram exists.
Definition at line 518 of file BCEngineMCMC.h.
|
inlinevirtual |
Interface allowing to execute arbitrary code for each new point of the MCMC whether it is accepted or not.
This method needs to be overloaded in the derived class
- Note
- This method is called for every call to the likelihood.
- Parameters
-
point point that was generated and checked ichain index of the chain accepted flag whether or not the point was accepted for the chain
Definition at line 1441 of file BCEngineMCMC.h.
void BCEngineMCMC::MCMCInitialize | ( | ) |
Resets all containers used in MCMC and initializes starting points.
Definition at line 2393 of file BCEngineMCMC.cxx.
|
inlinevirtual |
User hook called from MCMCInitialize().
MCMCUserInitialize() is called after all settings for the upcoming MCMC run are fixed but before the initial point is chosen and before any call to user methods such as LogLikelihood() or LogAPrioriProbability() are issued. MCMCUserInitialize() is useful for example to allocate separate copies of objects within the user model, one per chain, for thread safety.
- Note
- Any error inside MCMCUserInitialize() should be signaled via an exception.
Reimplemented in BCFitter.
Definition at line 1404 of file BCEngineMCMC.h.
|
inlinevirtual |
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.
With Remarginalize(), there is no lag, and MCMCUserIterationInterface() is called on every iteration before filling histograms.
In both instances of use, CalculateObservables() has already been called.
This function is not necessarily called for each call to the likelihood. Its frequency is affected by the lag; and for the factorized proposal function, it is called only after all parameters have been updated.
This method needs to be overloaded by the user to do anything.
- Note
- For uncertainty propagation, it is much more convenient BCObservable and CalculateObservables()
Reimplemented in BCMTF, and BCFitter.
Definition at line 1428 of file BCEngineMCMC.h.
bool BCEngineMCMC::Metropolis | ( | ) |
Runs Metropolis algorithm.
- Returns
- Success of action.
Definition at line 2183 of file BCEngineMCMC.cxx.
bool BCEngineMCMC::MetropolisPreRun | ( | ) |
Runs a pre run for the Metropolis algorithm.
- Returns
- Success of action.
Definition at line 1762 of file BCEngineMCMC.cxx.
|
virtual |
Check parameter tree against model.
- Parameters
-
partree Tree of parameters to check against model. checkObservables Flag for whether to check observables.
- Returns
- Whether tree's parameters match those of model.
Definition at line 1065 of file BCEngineMCMC.cxx.
void BCEngineMCMC::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.
Definition at line 1419 of file BCEngineMCMC.cxx.
unsigned BCEngineMCMC::PrintAllMarginalized | ( | const std::string & | filename, |
unsigned | hdiv = 1 , |
||
unsigned | vdiv = 1 |
||
) | const |
Print all marginalizations.
- Parameters
-
filename Path to file to print to hdiv Number of columns of plots per page vdiv Number of rows of plots per page
- Returns
- Number of plots printed
Definition at line 2932 of file BCEngineMCMC.cxx.
|
protectedvirtual |
Print best fit to log.
Reimplemented in BCIntegrate.
Definition at line 2774 of file BCEngineMCMC.cxx.
bool BCEngineMCMC::PrintCorrelationMatrix | ( | const std::string & | filename = "matrix.pdf" | ) | const |
Print a correlation matrix for the parameters.
- Parameters
-
filename Path to file to print correlation matrix to.
- Returns
- Success of action.
Definition at line 3278 of file BCEngineMCMC.cxx.
bool BCEngineMCMC::PrintCorrelationPlot | ( | const std::string & | filename = "correlation.pdf" , |
bool | include_observables = true |
||
) | const |
Print a correlation plot for the parameters.
- Parameters
-
filename Path to file to print correlation plot to. include_observables Flag for including observables (default: true)
- Returns
- Success of action.
Definition at line 3397 of file BCEngineMCMC.cxx.
|
protectedvirtual |
Print marginalization to log.
Reimplemented in BCIntegrate.
Definition at line 2809 of file BCEngineMCMC.cxx.
|
protectedvirtual |
Print model summary to log.
Definition at line 2746 of file BCEngineMCMC.cxx.
bool BCEngineMCMC::PrintParameterLatex | ( | const std::string & | filename | ) | const |
Print a LaTeX table of the parameters.
- Parameters
-
filename Path to file tp print LaTeX table of parameters to.
- Returns
- Success of action.
Definition at line 3533 of file BCEngineMCMC.cxx.
unsigned BCEngineMCMC::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.
- Parameters
-
filename Path to filename to print to. npar Number of parameters per page, print all on one page if set to zero or negative interval_content Probability mass to display in smallest X interval band quantile_values Vector of quantile values to draw rescale_ranges Flag for rescaling to range surveyed by MCMC chains
- Returns
- Number of pages printed.
Definition at line 2986 of file BCEngineMCMC.cxx.
void BCEngineMCMC::PrintParameters | ( | const std::vector< double > & | P, |
void(*)(const std::string &) | output = BCLog::OutSummary |
||
) | const |
Print parameters.
- Parameters
-
P vector of the parameter values to be printed output pointer to the output function to be used, which defaults to BCLog::OutSummary
Definition at line 2881 of file BCEngineMCMC.cxx.
|
virtual |
Prints a summary to the logs.
Definition at line 2736 of file BCEngineMCMC.cxx.
|
virtual |
The default proposal function is a Breit-Wigner random walk.
It can be overloaded by the user to set the proposal function.
- Parameters
-
ichain the chain index ipar the parameter index
- Returns
- the unscaled proposal point
Definition at line 1519 of file BCEngineMCMC.cxx.
|
virtual |
Marginalize from TTree.
- Parameters
-
autorange Flag for automatically choosing variable ranges to.
Definition at line 1253 of file BCEngineMCMC.cxx.
|
virtual |
Reset the MCMC variables.
Reimplemented in BCIntegrate.
Definition at line 2364 of file BCEngineMCMC.cxx.
|
static |
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)
- Parameters
-
means Vector of means of sample batches. variances Vector of variances of samples batches. n Number of samples in each batch. correctForSamplingVariability Flag to control correcting R value for initial sampling variability.
- Returns
- R value for set of batches of samples.
Definition at line 2106 of file BCEngineMCMC.cxx.
|
inline |
Set flag to correct convergence checking for initial sampling variability.
Definition at line 975 of file BCEngineMCMC.h.
|
protected |
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).
- Parameters
-
x Index of variable for horizontal axis. y Index of variable for vertical axis. flag Whether to fill 2D histogram.
Definition at line 2574 of file BCEngineMCMC.cxx.
|
inline |
Sets whether to fill particular H2 histogram: obs(y) vs.
obs(x)
- Parameters
-
x Index of observable for horizontal axis. y Index of observable for vertical axis. flag Whether to fill histogram.
Definition at line 941 of file BCEngineMCMC.h.
|
inline |
Sets whether to fill particular H2 histogram: obs(y) vs.
obs(x)
- Parameters
-
x Name of observable for horizontal axis. y Name of observable for vertical axis. flag Whether to fill histogram.
Definition at line 948 of file BCEngineMCMC.h.
|
inline |
Sets whether to fill particular H2 histogram: par(y) vs.
obs(x)
- Parameters
-
x Index of observable for horizontal axis. y Index of parameter for vertical axis. flag Whether to fill histogram.
Definition at line 955 of file BCEngineMCMC.h.
|
inline |
Sets whether to fill particular H2 histogram: par(y) vs.
obs(x)
- Parameters
-
x Name of observable for horizontal axis. y Name of parameter for vertical axis. flag Whether to fill histogram.
Definition at line 962 of file BCEngineMCMC.h.
|
inline |
Sets whether to fill particular H2 histogram: obs(y) vs.
par(x)
- Parameters
-
x Index of parameter for horizontal axis. y Index of observable for vertical axis. flag Whether to fill histogram.
Definition at line 927 of file BCEngineMCMC.h.
|
inline |
Sets whether to fill particular H2 histogram: obs(y) vs.
par(x)
- Parameters
-
x Name of parameter for horizontal axis. y Name of observable for vertical axis. flag Whether to fill histogram.
Definition at line 934 of file BCEngineMCMC.h.
|
inline |
Sets whether to fill particular H2 histogram: par(y) vs.
par(x)
- Parameters
-
x Index of parameter for horizontal axis. y Index of parameter for vertical axis. flag Whether to fill histogram.
Definition at line 913 of file BCEngineMCMC.h.
|
inline |
Sets whether to fill particular H2 histogram: par(y) vs.
par(x)
- Parameters
-
x Name of parameter for horizontal axis. y Name of parameter for vertical axis. flag Whether to fill histogram.
Definition at line 920 of file BCEngineMCMC.h.
|
inline |
Sets whether to fill histograms.
Applies to 1D and 2D histogams of all parameters and observables so far added.
Definition at line 901 of file BCEngineMCMC.h.
|
inline |
Sets the whether to fill histograms.
Applies to all parameters and observables so far added.
Definition at line 906 of file BCEngineMCMC.h.
|
inline |
Set if a (new) prerun should be performed.
Definition at line 966 of file BCEngineMCMC.h.
|
inline |
Set enlargement factor of range for when rescaling.
Definition at line 1009 of file BCEngineMCMC.h.
|
inline |
Sets maximum number of attempts to find a valid initial position.
Definition at line 840 of file BCEngineMCMC.h.
void BCEngineMCMC::SetInitialPositions | ( | const std::vector< double > & | x0s | ) |
Sets the initial positions for all chains.
- Parameters
-
x0s initial positions for all chains.
Definition at line 673 of file BCEngineMCMC.cxx.
|
inline |
Sets the initial positions for all chains.
- Parameters
-
x0s initial positions for all chains.
Definition at line 830 of file BCEngineMCMC.h.
|
inline |
Sets flag which defines initial position.
Definition at line 835 of file BCEngineMCMC.h.
|
inline |
Set the initial scale factors for the factorized proposal function.
The same factors are used for every chain and get updated during the pre-run.
- Parameters
-
scale a vector of doubles containing the scale factors for each parameter.
Definition at line 770 of file BCEngineMCMC.h.
|
inline |
Sets the maximum efficiency required for a chain.
Definition at line 815 of file BCEngineMCMC.h.
|
inline |
Sets the minimum efficiency required for a chain.
Definition at line 810 of file BCEngineMCMC.h.
|
inline |
Set weighting for multivariate proposal function covariance update.
value forced into [0, 1]
Definition at line 886 of file BCEngineMCMC.h.
|
inline |
Sets multivariate-proposal-function cholesky-decomposition nudge.
Definition at line 891 of file BCEngineMCMC.h.
|
inline |
Sets multivariate-proposal-function scale multiplier.
Definition at line 896 of file BCEngineMCMC.h.
void BCEngineMCMC::SetName | ( | const std::string & | name | ) |
Sets the name of the engine.
- Parameters
-
name Name of the engine
Definition at line 308 of file BCEngineMCMC.cxx.
|
inline |
Set the number of bins for the marginalized distribution of all parameters.
- Parameters
-
nbins Number of bins
Definition at line 994 of file BCEngineMCMC.h.
|
inline |
Sets the number of Markov chains which are run in parallel.
Definition at line 775 of file BCEngineMCMC.h.
|
inline |
Sets the number of iterations between scale adjustments and convergence checks in the pre-run.
Definition at line 800 of file BCEngineMCMC.h.
|
inline |
Sets the maximum number of iterations in the pre-run.
Definition at line 784 of file BCEngineMCMC.h.
|
inline |
Sets the number of iterations.
Definition at line 789 of file BCEngineMCMC.h.
void BCEngineMCMC::SetPrecision | ( | BCEngineMCMC::Precision | precision | ) |
Set the precision for the MCMC run.
Definition at line 353 of file BCEngineMCMC.cxx.
|
inline |
Copy precision for the MCMC run from other model.
Definition at line 984 of file BCEngineMCMC.h.
void BCEngineMCMC::SetPrecision | ( | const BCEngineMCMC & | other | ) |
Copy precision for the MCMC run from other model.
Definition at line 410 of file BCEngineMCMC.cxx.
|
inline |
Sets the number of prerun checks to make inbetween statistics clearing.
Definition at line 805 of file BCEngineMCMC.h.
void BCEngineMCMC::SetPrior | ( | unsigned | index, |
TF1 & | f, | ||
bool | logL = true |
||
) |
- Deprecated:
- Instead call: GetPrarameter(index)->SetPrior(new BCTF1Prior(f)) Set prior for a parameter.
- Parameters
-
index The parameter index f A function describing the prior logL Whether function is of log of prior (true) or just prior (false)
Definition at line 315 of file BCEngineMCMC.cxx.
|
inline |
- Deprecated:
- Instead call: GetParameter(name)->SetPrior(new BCTF1Prior(f)) Set prior for a parameter.
- Parameters
-
name The parameter name f A function describing the prior logL Whether function is of log of prior (true) or just prior (false)
Definition at line 1075 of file BCEngineMCMC.h.
void BCEngineMCMC::SetPrior | ( | unsigned | index, |
TH1 & | h, | ||
bool | interpolate = false |
||
) |
- Deprecated:
- Instead call: GetParameter(index)->SetPrior(new BCTH1Prior(h,interpolate)) Set prior for a parameter.
- Parameters
-
index parameter index h histogram describing the prior interpolate whether or not to use linear interpolation
Definition at line 324 of file BCEngineMCMC.cxx.
|
inline |
- Deprecated:
- Instead call: GetParameter(name)->SetPrior(new BCTH1Prior(h,interpolate)) Set prior for a parameter.
- Parameters
-
name parameter name h histogram describing the prior interpolate whether or not to use linear interpolation
- Returns
- success of action.
Definition at line 1145 of file BCEngineMCMC.h.
|
inline |
- Deprecated:
- Instead call: GetParameter(index)->SetPriorConstant() Set constant prior for this parameter
- Parameters
-
index the index of the parameter
Definition at line 1051 of file BCEngineMCMC.h.
|
inline |
- Deprecated:
- Instead call: GetParameter(name)->SetPriorConstant() Set constant prior for this parameter
- Parameters
-
name the name of the parameter
Definition at line 1058 of file BCEngineMCMC.h.
|
inline |
- Deprecated:
- Instead call: GetParameters().SetPriorConstantAll()
Definition at line 1150 of file BCEngineMCMC.h.
|
inline |
- Deprecated:
- Instead call: GetParameter(index)->Fix(value) Fixes parameter to value.
- Parameters
-
index The parameter index value The position of the delta function.
Definition at line 1083 of file BCEngineMCMC.h.
|
inline |
- Deprecated:
- Instead call: GetParameter(name)->Fix(value) Fixes parameter to value.
- Parameters
-
name The parameter name value The position of the delta function.
Definition at line 1091 of file BCEngineMCMC.h.
void BCEngineMCMC::SetPriorGauss | ( | unsigned | index, |
double | mean, | ||
double | sigma | ||
) |
- Deprecated:
- Instead call: GetParameter(index)->SetPrior(new BCGaussianPrior(mean,sigma)) Set Gaussian prior for a parameter.
- Parameters
-
index The parameter index mean The mean of the Gaussian sigma The sigma of the Gaussian
Definition at line 330 of file BCEngineMCMC.cxx.
|
inline |
- Deprecated:
- Instead call: GetParameter(name)->SetPrior(new BCGaussianPrior(mean,sigma)) Set Gaussian prior for a parameter.
- Parameters
-
name The parameter name mean The mean of the Gaussian sigma The sigma of the Gaussian
Definition at line 1108 of file BCEngineMCMC.h.
void BCEngineMCMC::SetPriorGauss | ( | unsigned | index, |
double | mode, | ||
double | sigma_below, | ||
double | sigma_above | ||
) |
- Deprecated:
- Instead call: GetParameter(index)->SetPrior(new BCSplitGaussianPrior(mode,sigma_below,sigma_above)) Set Gaussian prior for a parameter with two different widths.
- Parameters
-
index The parameter index mode The mode of the Gaussian sigma_below Standard deviation below mode. sigma_above Standard deviation above mode.
Definition at line 336 of file BCEngineMCMC.cxx.
|
inline |
- Deprecated:
- Instead call: GetParameter(name)->SetPrior(new BCSplitGaussianPrior(mode,sigma_below,sigma_above)) Set Gaussian prior for a parameter with two different widths.
- Parameters
-
name The parameter name mode The mode of the Gaussian sigma_below Standard deviation below mode. sigma_above Standard deviation above mode.
Definition at line 1127 of file BCEngineMCMC.h.
|
inline |
Set the degree of freedom of the proposal function for MCMC.
The default dof == 1
represents a Cauchy distribution. For any positive value of dof
, a Student's t distribution with the corresponding degree of freedom is used. For dof <= 0
, a Gaussian is used.
A small positive degree of freedom leads to fat tails in the proposal. This makes it easier to make a large jump in a single iteration but generally leads to a lower acceptance rate being optimal.
Definition at line 878 of file BCEngineMCMC.h.
|
inline |
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.
If flag == false
, use a factorized proposal in which every parameter is varied individually, one after the other. This means for N parameters, N calls to the target density occur until every parameter has been (attempted to be) varied exactly once. In contrast, with the multivariate proposal all parameters are varied simultaneously and a move can occur for a single call to the target density.
For both multivariate and factorized proposal, the acceptance rate of the Markov chain is tuned to lie within the limits given by SetMinimumEfficiency() and SetMaximumEfficiency().
Definition at line 859 of file BCEngineMCMC.h.
|
inline |
Set flag for rescaling histogram ranges after pre-run.
Definition at line 1004 of file BCEngineMCMC.h.
|
inline |
- Parameters
-
flag Flag whether to reuse user-defined observables from MCMC tree when looping through it.
Definition at line 999 of file BCEngineMCMC.h.
void BCEngineMCMC::UpdateChainIndex | ( | int | chain | ) |
Keep track of which chain is currently computed (within a thread).
- Warning
- Call this method only if you know what you are doing!
Definition at line 3750 of file BCEngineMCMC.cxx.
|
protected |
return appropriate update interval
- Parameters
-
N total number of iterations to be accomplished.
- Returns
- Appropriate interval to output after.
Definition at line 3639 of file BCEngineMCMC.cxx.
|
virtual |
Update multivariate proposal function covariances.
- Parameters
-
a update control factor.
Definition at line 1431 of file BCEngineMCMC.cxx.
|
virtual |
Update multivariate proposal function covariances.
Definition at line 1459 of file BCEngineMCMC.cxx.
|
protected |
Update Paramater TTree with scales and efficiencies.
Definition at line 817 of file BCEngineMCMC.cxx.
bool BCEngineMCMC::ValidMCMCTree | ( | TTree * | tree, |
bool | checkObservables = true |
||
) | const |
Check tree structure for MCMC tree.
- Parameters
-
tree MCMC Tree to check validity of. checkObservables Flag for whether to check observables.
- Returns
- Validity of tree.
Definition at line 874 of file BCEngineMCMC.cxx.
bool BCEngineMCMC::ValidParameterTree | ( | TTree * | tree | ) | const |
Check tree structure for parameter tree.
- Parameters
-
tree Parameter tree to check validity of.
- Returns
- Validty of tree.
Definition at line 894 of file BCEngineMCMC.cxx.
void BCEngineMCMC::WriteMarginalizedDistributions | ( | const std::string & | filename, |
const std::string & | option, | ||
bool | closeExistingFile = false |
||
) |
Write marginalization histograms to file.
- Parameters
-
filename Path to write file to. option Options passed to ROOT's TFile::Open. closeExistingFile if file is already open, whether to close it after writing.
Definition at line 464 of file BCEngineMCMC.cxx.
|
inline |
Turn on/off writing of Markov chain to root file.
If setting true, you must first set filename with function with filename arguments.
- Parameters
-
flag Flag for writing Markov chain (run and prerun) to ROOT file.
Definition at line 1016 of file BCEngineMCMC.h.
void BCEngineMCMC::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.
- Parameters
-
filename Name of file to write chain to. option file-open options (TFile), must be "NEW", "CREATE", "RECREATE", or "UPDATE" (i.e. writeable). flag_run Flag for writing run Markov chain to ROOT file (true) or not (false). flag_prerun Flag for writing prerun Markov chain to ROOT file (true) or not (false).
Definition at line 446 of file BCEngineMCMC.cxx.
void BCEngineMCMC::WriteMarkovChainPreRun | ( | bool | flag | ) |
Turn on/off writing of Markov chain to root file during prerun.
If setting either true, you must first set filename with function with filename arguments.
- Parameters
-
flag Flag for writing prerun Markov chain to ROOT file (true) or not (false).
Definition at line 438 of file BCEngineMCMC.cxx.
void BCEngineMCMC::WriteMarkovChainRun | ( | bool | flag | ) |
Turn on/off writing of Markov chain to root file during run.
If setting either true, you must first set filename with function with filename arguments.
- Parameters
-
flag Flag for writing run Markov chain to ROOT file (true) or not (false).
Definition at line 430 of file BCEngineMCMC.cxx.
Member Data Documentation
|
protected |
A BCH1D (with no histogram) for storing BCH1D drawing options.
Definition at line 1817 of file BCEngineMCMC.h.
|
protected |
A BCH2D (with no histogram) for storing BCH2D drawing options.
Definition at line 1821 of file BCEngineMCMC.h.
|
protected |
flag for correcting R value for initial sampling variability.
Definition at line 1760 of file BCEngineMCMC.h.
|
protected |
Vector of 2D marginalized distributions.
Definition at line 1779 of file BCEngineMCMC.h.
|
protected |
factor for enlarging range of histograms when rescaling.
Definition at line 1829 of file BCEngineMCMC.h.
|
protected |
Maximum number of attempts to make to set the initial position.
Definition at line 1731 of file BCEngineMCMC.h.
|
protected |
Variable which defines the initial position.
See enum MCMCInitialPosition for possible values.
Definition at line 1727 of file BCEngineMCMC.h.
|
protected |
Vector of local modes.
Definition at line 1813 of file BCEngineMCMC.h.
|
protected |
The current iteration number.
If not called within the running of the algorithm, return -1.
Definition at line 1618 of file BCEngineMCMC.h.
|
protected |
The intial position of each Markov chain.
The length of the vectors is equal to fMCMCNChains * fMCMCNParameters. First, the values of the first Markov chain are saved, then those of the second and so on
Definition at line 1714 of file BCEngineMCMC.h.
|
protected |
User-provided initial values of the scale factors of the factorized proposal function.
Definition at line 1682 of file BCEngineMCMC.h.
|
protected |
Number of iterations between scale adjustments and convergence checks in pre-run.
Definition at line 1623 of file BCEngineMCMC.h.
|
protected |
Output file for for writing MCMC Tree.
Definition at line 1656 of file BCEngineMCMC.h.
|
protected |
Output filename for for writing MCMC Tree.
Definition at line 1660 of file BCEngineMCMC.h.
|
protected |
Output file open option for for writing MCMC Tree.
Definition at line 1664 of file BCEngineMCMC.h.
|
protected |
The phase of the run.
Definition at line 1744 of file BCEngineMCMC.h.
|
protected |
Number of iterations between clearing of convergence stats in pre-run.
Definition at line 1627 of file BCEngineMCMC.h.
|
protected |
Degree of freedom of Student's t proposal.
If <= 0, use Gaussian proposal.
Definition at line 1740 of file BCEngineMCMC.h.
|
protected |
Scale factors for proposal functions.
First index is for the chain, second index is for the parameters (if not using Multivariate proposal function).
Definition at line 1678 of file BCEngineMCMC.h.
|
protected |
Flag for using multivariate proposal function.
Definition at line 1735 of file BCEngineMCMC.h.
|
protected |
The current states of each Markov chain.
Definition at line 1748 of file BCEngineMCMC.h.
|
protected |
Statistics for each Markov chain.
Definition at line 1752 of file BCEngineMCMC.h.
|
protected |
Statistics across all Markov chains.
Definition at line 1756 of file BCEngineMCMC.h.
|
protected |
The tree containing the Markov chains.
Definition at line 1789 of file BCEngineMCMC.h.
|
protected |
flag for whether MCMC Tree successfully loaded.
Definition at line 1793 of file BCEngineMCMC.h.
|
protected |
flag for whether to reuse MCMC Tree's observables.
Definition at line 1797 of file BCEngineMCMC.h.
|
protected |
weighting parameter for multivariate-proposal-function covariance update.
Definition at line 1699 of file BCEngineMCMC.h.
|
protected |
Number of multivariate-proposal-function covariance updates performed.
Definition at line 1695 of file BCEngineMCMC.h.
|
protected |
multivariate-proposal-function cholesky-decomposition nudge.
Definition at line 1703 of file BCEngineMCMC.h.
|
protected |
Cholesky decompositions for multivariate proposal function.
Index is over chain.
Definition at line 1691 of file BCEngineMCMC.h.
|
protected |
Covariance matrices used in multivariate proposal functions.
Definition at line 1686 of file BCEngineMCMC.h.
|
protected |
factor to multiply or divide scale factors by in adjusting multivariate-proposal-function scales.
Definition at line 1707 of file BCEngineMCMC.h.
|
protected |
Name of the engine.
Definition at line 1593 of file BCEngineMCMC.h.
|
protected |
The tree containing the parameter information.
Definition at line 1809 of file BCEngineMCMC.h.
|
protected |
Vector of pairs of indices for which 2D histograms should be stored.
a negative index indicates an observable, with observable zero as -1, observable one as -2, etc.
Definition at line 1785 of file BCEngineMCMC.h.
|
protected |
flag for rescaling of histograms after pre-run.
Definition at line 1825 of file BCEngineMCMC.h.
|
protected |
Safe name of the engine for use in naming ROOT objects.
Definition at line 1597 of file BCEngineMCMC.h.
The documentation for this class was generated from the following files:
- /root/bat/BAT/BCEngineMCMC.h
- /root/bat/src/BCEngineMCMC.cxx