![]() |
MaCh3
2.2.3
Reference Guide
|
MCMC algorithm implementing the Metropolis–Rosenbluth–Rosenbluth–Teller–Teller (MR \(^2\)T \(^2\)) method. More...
#include <Fitters/MR2T2.h>
Public Member Functions | |
MR2T2 (manager *const manager) | |
Constructor. More... | |
virtual | ~MR2T2 ()=default |
Destructor. More... | |
![]() | |
MCMCBase (manager *const fitMan) | |
Constructor. More... | |
virtual | ~MCMCBase ()=default |
Destructor. More... | |
void | RunMCMC () override |
Actual implementation of MCMC fitting algorithm. More... | |
void | StartFromPreviousFit (const std::string &FitName) override |
Allow to start from previous fit/chain. More... | |
void | setChainLength (unsigned int L) |
Set how long chain should be. More... | |
![]() | |
FitterBase (manager *const fitMan) | |
Constructor. More... | |
virtual | ~FitterBase () |
Destructor for the FitterBase class. More... | |
void | AddSampleHandler (SampleHandlerBase *sample) |
This function adds a sample PDF object to the analysis framework. The sample PDF object will be utilized in fitting procedures or likelihood scans. More... | |
void | AddSystObj (ParameterHandlerBase *cov) |
This function adds a Covariance object to the analysis framework. The Covariance object will be utilized in fitting procedures or likelihood scans. More... | |
void | DragRace (const int NLaps=100) |
Calculates the required time for each sample or covariance object in a drag race simulation. Inspired by Dan's feature. More... | |
void | RunLLHScan () |
Perform a 1D likelihood scan. More... | |
void | GetStepScaleBasedOnLLHScan () |
LLH scan is good first estimate of step scale. More... | |
void | Run2DLLHScan () |
Perform a 2D likelihood scan. More... | |
void | RunSigmaVar () |
Perform a 2D and 1D sigma var for all samples. More... | |
void | RunSigmaVarFD () |
Perform a 1D sigma var for all samples. More... | |
std::string | GetName () const |
Get name of class. More... | |
Protected Member Functions | |
void | DoStep () override |
The MCMC step proposal and acceptance. More... | |
void | ProposeStep () override |
Propose a step. More... | |
double | AcceptanceProbability () override |
Step acceptance probability. More... | |
![]() | |
void | DoMCMCStep () |
The full StartStep->DoStep->EndStep chain. More... | |
void | PreStepProcess () |
Actions before step proposal [start stopwatch]. More... | |
void | PostStepProcess () |
Actions after step proposal [end stopwatch, fill tree]. More... | |
bool | IsStepAccepted (const double acc_prob) |
Is step accepted? More... | |
void | AcceptStep () |
Accept a step. More... | |
void | AdaptiveStep () |
Adaptive MCMC step. More... | |
void | PrintProgress () |
Print the progress. More... | |
![]() | |
void | ProcessMCMC () |
Process MCMC output. More... | |
void | PrepareOutput () |
Prepare the output file. More... | |
void | SaveOutput () |
Save output and close files. More... | |
void | SanitiseInputs () |
Remove obsolete memory and make other checks before fit starts. More... | |
void | SaveSettings () |
Save the settings that the MCMC was run with. More... | |
bool | GetScaneRange (std::map< std::string, std::vector< double >> &scanRanges) |
YSP: Set up a mapping to store parameters with user-specified ranges, suggested by D. Barrow. More... | |
bool | CheckSkipParameter (const std::vector< std::string > &SkipVector, const std::string &ParamName) const |
KS: Check whether we want to skip parameter using skip vector. More... | |
void | CustomRange (const std::string &ParName, const double sigma, double &ParamShiftValue) |
For comparison with P-Theta we usually have to apply different parameter values then usual 1, 3 sigma. More... | |
Additional Inherited Members | |
![]() | |
bool | out_of_bounds |
Do we reject based on hitting boundaries in systs. More... | |
bool | accept |
Accept. More... | |
unsigned int | chainLength |
number of steps in chain More... | |
bool | anneal |
simulated annealing More... | |
double | AnnealTemp |
simulated annealing temperature More... | |
![]() | |
manager * | fitMan |
The manager. More... | |
unsigned int | step |
current state More... | |
double | logLCurr |
current likelihood More... | |
double | logLProp |
proposed likelihood More... | |
double | accProb |
current acceptance prob More... | |
int | accCount |
counts accepted steps More... | |
unsigned int | stepStart |
step start, by default 0 if we start from previous chain then it will be different More... | |
std::vector< double > | sample_llh |
store the llh breakdowns More... | |
std::vector< double > | syst_llh |
systematic llh breakdowns More... | |
std::vector< SampleHandlerBase * > | samples |
Sample holder. More... | |
unsigned int | TotalNSamples |
Total number of samples used. More... | |
std::vector< ParameterHandlerBase * > | systematics |
Systematic holder. More... | |
std::unique_ptr< TStopwatch > | clock |
tells global time how long fit took More... | |
std::unique_ptr< TStopwatch > | stepClock |
tells how long single step/fit iteration took More... | |
double | stepTime |
Time of single step. More... | |
std::unique_ptr< TRandom3 > | random |
Random number. More... | |
TFile * | outputFile |
Output. More... | |
TDirectory * | CovFolder |
Output cov folder. More... | |
TDirectory * | SampleFolder |
Output sample folder. More... | |
TTree * | outTree |
Output tree with posteriors. More... | |
int | auto_save |
auto save every N steps More... | |
bool | fTestLikelihood |
Necessary for some fitting algorithms like PSO. More... | |
bool | FileSaved |
Checks if file saved not repeat some operations. More... | |
bool | SettingsSaved |
Checks if setting saved not repeat some operations. More... | |
bool | OutputPrepared |
Checks if output prepared not repeat some operations. More... | |
std::string | AlgorithmName |
Name of fitting algorithm that is being used. More... | |
MCMC algorithm implementing the Metropolis–Rosenbluth–Rosenbluth–Teller–Teller (MR \(^2\)T \(^2\)) method.
The acceptance probability for a proposed step is given by
\[ \alpha = \min \Biggl( 1, \exp\bigl( \log \mathcal{L}_\text{prop} - \log \mathcal{L}_\text{curr} \bigr) \Biggr), \]
where \(\log \mathcal{L}_\text{prop}\) is the log-likelihood of the proposed state and \(\log \mathcal{L}_\text{curr}\) is the log-likelihood of the current state.
If annealing is enabled, the acceptance probability is modified as
\[ \alpha_\text{anneal} = \min \Biggl( 1, \exp\Bigl( -\frac{\log \mathcal{L}_\text{prop} - \log \mathcal{L}_\text{curr}} {\exp(-\text{step}/T_\text{anneal})} \Bigr) \Biggr), \]
where \(T_\text{anneal}\) is the current annealing temperature and \(\text{step}\) is the iteration index.
MR2T2::MR2T2 | ( | manager *const | manager | ) |
Constructor.
fitMan | A pointer to a manager object, which will handle all settings. |
Definition at line 5 of file MR2T2.cpp.
|
virtualdefault |
Destructor.
|
overrideprotectedvirtual |
|
overrideprotectedvirtual |
The MCMC step proposal and acceptance.
Implements MCMCBase.
Definition at line 11 of file MR2T2.cpp.
|
overrideprotectedvirtual |
Propose a step.
Implements MCMCBase.
Definition at line 25 of file MR2T2.cpp.