MaCh3  2.2.3
Reference Guide
Public Member Functions | Protected Member Functions | List of all members
MR2T2 Class Reference

MCMC algorithm implementing the Metropolis–Rosenbluth–Rosenbluth–Teller–Teller (MR \(^2\)T \(^2\)) method. More...

#include <Fitters/MR2T2.h>

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

Public Member Functions

 MR2T2 (manager *const manager)
 Constructor. More...
 
virtual ~MR2T2 ()=default
 Destructor. More...
 
- Public Member Functions inherited from MCMCBase
 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...
 
- Public Member Functions inherited from FitterBase
 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...
 
- Protected Member Functions inherited from MCMCBase
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...
 
- Protected Member Functions inherited from FitterBase
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

- Protected Attributes inherited from MCMCBase
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...
 
- Protected Attributes inherited from FitterBase
managerfitMan
 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...
 

Detailed Description

MCMC algorithm implementing the Metropolis–Rosenbluth–Rosenbluth–Teller–Teller (MR \(^2\)T \(^2\)) method.

Note
MR \(^2\)T \(^2\) is also known as Metropolis-Hastings [16]

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.

Author
Asher Kaboth

Definition at line 25 of file MR2T2.h.

Constructor & Destructor Documentation

◆ MR2T2()

MR2T2::MR2T2 ( manager *const  manager)

Constructor.

Parameters
fitManA pointer to a manager object, which will handle all settings.

Definition at line 5 of file MR2T2.cpp.

5  : MCMCBase(man) {
6 // ***************
7  AlgorithmName = "MR2T2";
8 }
MaCh3Plotting::PlottingManager * man
std::string AlgorithmName
Name of fitting algorithm that is being used.
Definition: FitterBase.h:170
MCMCBase(manager *const fitMan)
Constructor.
Definition: MCMCBase.cpp:6

◆ ~MR2T2()

virtual MR2T2::~MR2T2 ( )
virtualdefault

Destructor.

Member Function Documentation

◆ AcceptanceProbability()

double MR2T2::AcceptanceProbability ( )
overrideprotectedvirtual

Step acceptance probability.

Implements MCMCBase.

Definition at line 103 of file MR2T2.cpp.

103  {
104 // **********************
105  // Set the acceptance probability to zero
106  double acc_prob = 0.0;
107 
108  // Calculate acceptance probability
109  if (anneal)
110  acc_prob = std::min(1., std::exp(-(logLProp - logLCurr) / (std::exp(-step / AnnealTemp))));
111  else
112  acc_prob = std::min(1., std::exp(logLCurr - logLProp));
113 
114  return acc_prob;
115 }
double logLProp
proposed likelihood
Definition: FitterBase.h:117
unsigned int step
current state
Definition: FitterBase.h:113
double logLCurr
current likelihood
Definition: FitterBase.h:115
bool anneal
simulated annealing
Definition: MCMCBase.h:70
double AnnealTemp
simulated annealing temperature
Definition: MCMCBase.h:72

◆ DoStep()

void MR2T2::DoStep ( )
overrideprotectedvirtual

The MCMC step proposal and acceptance.

Implements MCMCBase.

Definition at line 11 of file MR2T2.cpp.

11  {
12 // *******************
13  ProposeStep();
14  // Does the MCMC accept this step?
16 
18  {
19  AcceptStep();
20  }
21 }
double accProb
current acceptance prob
Definition: FitterBase.h:119
void AcceptStep()
Accept a step.
Definition: MCMCBase.cpp:178
bool out_of_bounds
Do we reject based on hitting boundaries in systs.
Definition: MCMCBase.h:61
bool IsStepAccepted(const double acc_prob)
Is step accepted?
Definition: MCMCBase.cpp:160
void ProposeStep() override
Propose a step.
Definition: MR2T2.cpp:25
double AcceptanceProbability() override
Step acceptance probability.
Definition: MR2T2.cpp:103

◆ ProposeStep()

void MR2T2::ProposeStep ( )
overrideprotectedvirtual

Propose a step.

Implements MCMCBase.

Definition at line 25 of file MR2T2.cpp.

25  {
26 // *******************
27  // Initial likelihood
28  out_of_bounds = false;
29 
30  double llh = 0.0;
31 
32  // Loop over the systematics and propose the initial step
33  for (size_t s = 0; s < systematics.size(); ++s)
34  {
35  // Could throw the initial value here to do MCMC stability studies
36  // Propose the steps for the systematics
37  systematics[s]->ProposeStep();
38 
39  // Get the likelihood from the systematics
40  syst_llh[s] = systematics[s]->GetLikelihood();
41  llh += syst_llh[s];
42 
43 #ifdef DEBUG
44  if (debug)
45  debugFile << "LLH after " << systematics[s]->GetName() << " " << llh << std::endl;
46 #endif
47  }
48 
49  // Check if we've hit a boundary in the systematics
50  // In this case we can save time by not having to reconfigure the simulation
51  if (llh >= M3::_LARGE_LOGL_)
52  {
53  out_of_bounds = true;
54 #ifdef DEBUG
55  if (debug)
56  debugFile << "Rejecting based on boundary" << std::endl;
57 #endif
58  }
59 
60  // Only reweight when we have a good parameter configuration
61  // This speeds things up considerably because for every bad parameter configuration we don't have to reweight the MC
62  if (!out_of_bounds)
63  {
64  // Could multi-thread this
65  // But since sample reweight is multi-threaded it's probably better to do that
66  for (size_t i = 0; i < samples.size(); ++i)
67  {
68  samples[i]->Reweight();
69  }
70 
71  // DB for atmospheric event by event sample migration, need to fully reweight all samples to allow event passing prior to likelihood evaluation
72  for (size_t i = 0; i < samples.size(); ++i)
73  {
74  // Get the sample likelihoods and add them
75  sample_llh[i] = samples[i]->GetLikelihood();
76  llh += sample_llh[i];
77 #ifdef DEBUG
78  if (debug)
79  debugFile << "LLH after sample " << i << " " << llh << std::endl;
80 #endif
81  }
82 
83  // For when we don't have to reweight, set sample to madness
84  }
85  else
86  {
87  for (size_t i = 0; i < samples.size(); ++i)
88  {
89  // Set the sample_llh[i] to be madly high also to signify a step out of bounds
91 #ifdef DEBUG
92  if (debug)
93  debugFile << "LLH after REJECT sample " << i << " " << llh << std::endl;
94 #endif
95  }
96  }
97  // Save the proposed likelihood (class member)
98  logLProp = llh;
99 }
std::vector< SampleHandlerBase * > samples
Sample holder.
Definition: FitterBase.h:131
std::vector< double > sample_llh
store the llh breakdowns
Definition: FitterBase.h:126
std::vector< double > syst_llh
systematic llh breakdowns
Definition: FitterBase.h:128
std::vector< ParameterHandlerBase * > systematics
Systematic holder.
Definition: FitterBase.h:136
constexpr static const double _LARGE_LOGL_
Large Likelihood is used it parameter go out of physical boundary, this indicates in MCMC that such s...
Definition: Core.h:73

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