MaCh3  2.4.2
Reference Guide
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
FitterBase Class Referenceabstract

Base class for implementing fitting algorithms. More...

#include <Fitters/FitterBase.h>

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

Public Member Functions

 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...
 
virtual void RunMCMC ()=0
 The specific fitting algorithm implemented in this function depends on the derived class. It could be Markov Chain Monte Carlo (MCMC), MinuitFit, or another algorithm. 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 RunLLHMap ()
 Perform a general multi-dimensional 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 1D/2D sigma var for all samples. More...
 
virtual void StartFromPreviousFit (const std::string &FitName)
 Allow to start from previous fit/chain. More...
 
std::string GetName () const
 Get name of class. More...
 

Protected Member Functions

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 GetScanRange (std::map< std::string, std::vector< double >> &scanRanges) const
 YSP: Set up a mapping to store parameters with user-specified ranges, suggested by D. Barrow. More...
 
void GetParameterScanRange (const ParameterHandlerBase *cov, const int i, double &CentralValue, double &lower, double &upper, const int n_points, const std::string &suffix="") const
 Helper function to get parameter scan range, central value. 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) const
 For comparison with other fitting frameworks (like P-Theta) we usually have to apply different parameter values then usual 1, 3 sigma. More...
 

Protected Attributes

ManagerfitMan
 The manager for configuration handling. 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, single SampleHandler can store more than one analysis sample! 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

Base class for implementing fitting algorithms.

This class wraps MaCh3 classes like SampleHandlerBase and ParameterHandlerBase. It serves as a base for different fitting algorithms and for validation techniques such as LLH scans.

Definition at line 26 of file FitterBase.h.

Constructor & Destructor Documentation

◆ FitterBase()

_MaCh3_Safe_Include_Start_ _MaCh3_Safe_Include_End_ FitterBase::FitterBase ( Manager *const  fitMan)

Constructor.

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

Definition at line 16 of file FitterBase.cpp.

16  : fitMan(man) {
17 // *************************
18  AlgorithmName = "";
19  //Get mach3 modes from Manager
20  random = std::make_unique<TRandom3>(Get<int>(fitMan->raw()["General"]["Seed"], __FILE__, __LINE__));
21 
22  // Counter of the accepted # of steps
23  accCount = 0;
24  step = 0;
25  stepStart = 0;
26 
27  clock = std::make_unique<TStopwatch>();
28  stepClock = std::make_unique<TStopwatch>();
29  #ifdef DEBUG
30  // Fit summary and debug info
31  debug = GetFromManager<bool>(fitMan->raw()["General"]["Debug"], false, __FILE__ , __LINE__);
32  #endif
33 
34  auto outfile = Get<std::string>(fitMan->raw()["General"]["OutputFile"], __FILE__ , __LINE__);
35  // Save output every auto_save steps
36  //you don't want this too often https://root.cern/root/html606/TTree_8cxx_source.html#l01229
37  auto_save = Get<int>(fitMan->raw()["General"]["MCMC"]["AutoSave"], __FILE__ , __LINE__);
38 
39  // Set the output file
40  outputFile = M3::Open(outfile, "RECREATE", __FILE__, __LINE__);
41  outputFile->cd();
42  // Set output tree
43  outTree = new TTree("posteriors", "Posterior_Distributions");
44  // Auto-save every 200MB, the bigger the better https://root.cern/root/html606/TTree_8cxx_source.html#l01229
45  outTree->SetAutoSave(-200E6);
46 
47  FileSaved = false;
48  SettingsSaved = false;
49  OutputPrepared = false;
50 
51  //Create TDirectory
52  CovFolder = outputFile->mkdir("CovarianceFolder");
53  outputFile->cd();
54  SampleFolder = outputFile->mkdir("SampleFolder");
55  outputFile->cd();
56 
57  #ifdef DEBUG
58  // Prepare the output log file
59  if (debug) debugFile.open((outfile+".log").c_str());
60  #endif
61 
62  TotalNSamples = 0;
63  fTestLikelihood = GetFromManager<bool>(fitMan->raw()["General"]["Fitter"]["FitTestLikelihood"], false, __FILE__ , __LINE__);
64 }
std::unique_ptr< TRandom3 > random
Random number.
Definition: FitterBase.h:144
int accCount
counts accepted steps
Definition: FitterBase.h:119
bool OutputPrepared
Checks if output prepared not repeat some operations.
Definition: FitterBase.h:165
TFile * outputFile
Output.
Definition: FitterBase.h:147
unsigned int step
current state
Definition: FitterBase.h:111
bool SettingsSaved
Checks if setting saved not repeat some operations.
Definition: FitterBase.h:163
bool FileSaved
Checks if file saved not repeat some operations.
Definition: FitterBase.h:161
std::string AlgorithmName
Name of fitting algorithm that is being used.
Definition: FitterBase.h:168
std::unique_ptr< TStopwatch > clock
tells global time how long fit took
Definition: FitterBase.h:137
Manager * fitMan
The manager for configuration handling.
Definition: FitterBase.h:108
unsigned int stepStart
step start, by default 0 if we start from previous chain then it will be different
Definition: FitterBase.h:121
std::unique_ptr< TStopwatch > stepClock
tells how long single step/fit iteration took
Definition: FitterBase.h:139
TDirectory * CovFolder
Output cov folder.
Definition: FitterBase.h:149
TDirectory * SampleFolder
Output sample folder.
Definition: FitterBase.h:151
unsigned int TotalNSamples
Total number of samples used, single SampleHandler can store more than one analysis sample!
Definition: FitterBase.h:131
int auto_save
auto save every N steps
Definition: FitterBase.h:155
bool fTestLikelihood
Necessary for some fitting algorithms like PSO.
Definition: FitterBase.h:158
TTree * outTree
Output tree with posteriors.
Definition: FitterBase.h:153
YAML::Node const & raw() const
Return config.
Definition: Manager.h:41
TFile * Open(const std::string &Name, const std::string &Type, const std::string &File, const int Line)
Opens a ROOT file with the given name and mode.

◆ ~FitterBase()

FitterBase::~FitterBase ( )
virtual

Destructor for the FitterBase class.

Definition at line 68 of file FitterBase.cpp.

68  {
69 // *************************
70  SaveOutput();
71 
72  if(outputFile != nullptr) delete outputFile;
73  outputFile = nullptr;
74  MACH3LOG_DEBUG("Closing MaCh3 Fitter Engine");
75 }
#define MACH3LOG_DEBUG
Definition: MaCh3Logger.h:34
void SaveOutput()
Save output and close files.
Definition: FitterBase.cpp:231

Member Function Documentation

◆ AddSampleHandler()

void FitterBase::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.

Parameters
sampleA pointer to a sample PDF object derived from ParameterHandlerBase.

Definition at line 262 of file FitterBase.cpp.

262  {
263 // *************************
264  // Check if any subsample name collides with already-registered subsamples
265  for (const auto &s : samples) {
266  for (int iExisting = 0; iExisting < s->GetNsamples(); ++iExisting) {
267  for (int iNew = 0; iNew < sample->GetNsamples(); ++iNew) {
268  if (s->GetSampleTitle(iExisting) == sample->GetSampleTitle(iNew)) {
270  "Duplicate sample title '{}' in handler {} detected: "
271  "same title exist in handler ", sample->GetSampleTitle(iNew),
272  sample->GetName(), s->GetName());
273  throw MaCh3Exception(__FILE__, __LINE__);
274  }
275  }
276  }
277  }
278 
279  for (const auto &s : samples) {
280  if (s->GetName() == sample->GetName()) {
281  MACH3LOG_WARN("SampleHandler with name '{}' already exists!", sample->GetName());
282  MACH3LOG_WARN("Is it intended?");
283  //throw MaCh3Exception(__FILE__ , __LINE__ );
284  }
285  }
286  // Save additional info from samples
287  SampleFolder->cd();
288 
290  TotalNSamples += sample->GetNsamples();
291  MACH3LOG_INFO("Adding {} object, with {} samples", sample->GetName(), sample->GetNsamples());
292  samples.push_back(sample);
293  outputFile->cd();
294 }
#define MACH3LOG_ERROR
Definition: MaCh3Logger.h:37
#define MACH3LOG_INFO
Definition: MaCh3Logger.h:35
#define MACH3LOG_WARN
Definition: MaCh3Logger.h:36
std::vector< SampleHandlerBase * > samples
Sample holder.
Definition: FitterBase.h:129
Custom exception class used throughout MaCh3.
virtual void SaveAdditionalInfo(TDirectory *Dir)
Store additional info in a chan.
virtual M3::int_t GetNsamples()
virtual std::string GetSampleTitle(const int Sample) const =0
virtual std::string GetName() const =0

◆ AddSystObj()

void FitterBase::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.

Parameters
covA pointer to a Covariance object derived from ParameterHandlerBase.

Definition at line 298 of file FitterBase.cpp.

298  {
299 // *************************
300  MACH3LOG_INFO("Adding systematic object {}, with {} params", cov->GetName(), cov->GetNumParams());
301  // KS: Need to make sure we don't have params with same name, otherwise ROOT I/O and parts of MaCh3 will be terribly confused...
302  for (size_t s = 0; s < systematics.size(); ++s)
303  {
304  for (int iPar = 0; iPar < systematics[s]->GetNumParams(); ++iPar)
305  {
306  for (int i = 0; i < cov->GetNumParams(); ++i)
307  {
308  if(systematics[s]->GetParName(iPar) == cov->GetParName(i)){
309  MACH3LOG_ERROR("ParameterHandler {} has param '{}' which already exists in in {}, with name {}",
310  cov->GetName(), cov->GetParName(i), systematics[s]->GetName(), systematics[s]->GetParName(iPar));
311  throw MaCh3Exception(__FILE__ , __LINE__ );
312  }
313  // Same for fancy name
314  if(systematics[s]->GetParFancyName(iPar) == cov->GetParFancyName(i)){
315  MACH3LOG_ERROR("ParameterHandler {} has param '{}' which already exists in {}, with name {}",
316  cov->GetName(), cov->GetParFancyName(i), systematics[s]->GetName(), systematics[s]->GetParFancyName(iPar));
317  throw MaCh3Exception(__FILE__ , __LINE__ );
318  }
319  }
320  }
321  }
322 
323  systematics.push_back(cov);
324 
325  CovFolder->cd();
326  std::vector<double> n_vec(cov->GetNumParams());
327  for (int i = 0; i < cov->GetNumParams(); ++i) {
328  n_vec[i] = cov->GetParInit(i);
329  }
330  cov->GetCovMatrix()->Write(cov->GetName().c_str());
331 
332  TH2D* CorrMatrix = cov->GetCorrelationMatrix();
333  CorrMatrix->Write((cov->GetName() + std::string("_Corr")).c_str());
334  delete CorrMatrix;
335 
336  // If we have yaml config file for covariance let's save it
337  YAML::Node Config = cov->GetConfig();
338  if(!Config.IsNull())
339  {
340  TMacro ConfigSave = YAMLtoTMacro(Config, (std::string("Config_") + cov->GetName()));
341  ConfigSave.Write();
342  }
343 
344  outputFile->cd();
345 }
TMacro YAMLtoTMacro(const YAML::Node &yaml_node, const std::string &name)
Convert a YAML node to a ROOT TMacro object.
Definition: YamlHelper.h:167
std::vector< ParameterHandlerBase * > systematics
Systematic holder.
Definition: FitterBase.h:134
int GetNumParams() const
Get total number of parameters.
std::string GetParName(const int i) const
Get name of parameter.
TMatrixDSym * GetCovMatrix() const
Return covariance matrix.
std::string GetParFancyName(const int i) const
Get fancy name of the Parameter.
TH2D * GetCorrelationMatrix()
KS: Convert covariance matrix to correlation matrix and return TH2D which can be used for fancy plott...
std::string GetName() const
Get name of covariance.
double GetParInit(const int i) const
Get prior parameter value.
YAML::Node GetConfig() const
Getter to return a copy of the YAML node.

◆ CheckSkipParameter()

bool FitterBase::CheckSkipParameter ( const std::vector< std::string > &  SkipVector,
const std::string &  ParamName 
) const
protected

KS: Check whether we want to skip parameter using skip vector.

Definition at line 543 of file FitterBase.cpp.

543  {
544 // *************************
545  bool skip = false;
546  for(unsigned int is = 0; is < SkipVector.size(); ++is)
547  {
548  if(ParamName.substr(0, SkipVector[is].length()) == SkipVector[is])
549  {
550  skip = true;
551  break;
552  }
553  }
554  return skip;
555 }

◆ CustomRange()

void FitterBase::CustomRange ( const std::string &  ParName,
const double  sigma,
double &  ParamShiftValue 
) const
protected

For comparison with other fitting frameworks (like P-Theta) we usually have to apply different parameter values then usual 1, 3 sigma.

Example YAML format:

SigmaVar:
Q2_norm_7: { "3": 2.0 }
SRC_Norm_O: { "-1": 0.5, "1": 1.5, "3": 2.0 }

Definition at line 1291 of file FitterBase.cpp.

1291  {
1292 // *************************
1293  if(!fitMan->raw()["SigmaVar"]["CustomRange"]) return;
1294 
1295  auto Config = fitMan->raw()["SigmaVar"]["CustomRange"];
1296 
1297  const auto sigmaStr = std::to_string(static_cast<int>(std::round(sigma)));
1298 
1299  if (Config[ParName] && Config[ParName][sigmaStr]) {
1300  ParamShiftValue = Config[ParName][sigmaStr].as<double>();
1301  MACH3LOG_INFO(" ::: setting custom range from config ::: {} -> {}", ParName, ParamShiftValue);
1302  }
1303 }

◆ DragRace()

void FitterBase::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.

Parameters
NLapsnumber of laps, every part of Fitter will be tested with given number of laps and you will get total and average time

Definition at line 461 of file FitterBase.cpp.

461  {
462 // *************************
463  MACH3LOG_INFO("Let the Race Begin!");
464  MACH3LOG_INFO("All tests will be performed with {} threads", M3::GetNThreads());
465 
466  // Reweight the MC
467  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs)
468  {
469  TStopwatch clockRace;
470  clockRace.Start();
471  for(int Lap = 0; Lap < NLaps; ++Lap) {
472  samples[ivs]->Reweight();
473  }
474  clockRace.Stop();
475  MACH3LOG_INFO("It took {:.4f} s to reweights {} times sample: {}", clockRace.RealTime(), NLaps, samples[ivs]->GetName());
476  MACH3LOG_INFO("On average {:.6f}", clockRace.RealTime()/NLaps);
477  }
478 
479  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs)
480  {
481  TStopwatch clockRace;
482  clockRace.Start();
483  for(int Lap = 0; Lap < NLaps; ++Lap) {
484  samples[ivs]->GetLikelihood();
485  }
486  clockRace.Stop();
487  MACH3LOG_INFO("It took {:.4f} s to calculate GetLikelihood {} times sample: {}", clockRace.RealTime(), NLaps, samples[ivs]->GetName());
488  MACH3LOG_INFO("On average {:.6f}", clockRace.RealTime()/NLaps);
489  }
490  // Get vector of proposed steps. If we want to run LLH scan or something else after we need to revert changes after proposing steps multiple times
491  std::vector<std::vector<double>> StepsValuesBefore(systematics.size());
492  for (size_t s = 0; s < systematics.size(); ++s) {
493  StepsValuesBefore[s] = systematics[s]->GetProposed();
494  }
495  for (size_t s = 0; s < systematics.size(); ++s) {
496  TStopwatch clockRace;
497  clockRace.Start();
498  for(int Lap = 0; Lap < NLaps; ++Lap) {
499  systematics[s]->ProposeStep();
500  }
501  clockRace.Stop();
502  MACH3LOG_INFO("It took {:.4f} s to propose step {} times cov: {}", clockRace.RealTime(), NLaps, systematics[s]->GetName());
503  MACH3LOG_INFO("On average {:.6f}", clockRace.RealTime()/NLaps);
504  }
505  for (size_t s = 0; s < systematics.size(); ++s) {
506  systematics[s]->SetParameters(StepsValuesBefore[s]);
507  }
508 
509  for (size_t s = 0; s < systematics.size(); ++s) {
510  TStopwatch clockRace;
511  clockRace.Start();
512  for(int Lap = 0; Lap < NLaps; ++Lap) {
513  systematics[s]->GetLikelihood();
514  }
515  clockRace.Stop();
516  MACH3LOG_INFO("It took {:.4f} s to calculate get likelihood {} times cov: {}", clockRace.RealTime(), NLaps, systematics[s]->GetName());
517  MACH3LOG_INFO("On average {:.6f}", clockRace.RealTime()/NLaps);
518  }
519  MACH3LOG_INFO("End of race");
520 }
int GetNThreads()
number of threads which we need for example for TRandom3
Definition: Monitor.cpp:372

◆ GetName()

std::string FitterBase::GetName ( ) const
inline

Get name of class.

Definition at line 70 of file FitterBase.h.

70 {return AlgorithmName;};

◆ GetParameterScanRange()

void FitterBase::GetParameterScanRange ( const ParameterHandlerBase cov,
const int  i,
double &  CentralValue,
double &  lower,
double &  upper,
const int  n_points,
const std::string &  suffix = "" 
) const
protected

Helper function to get parameter scan range, central value.

Definition at line 559 of file FitterBase.cpp.

560  {
561 // *************************
562  // YSP: Set up a mapping to store parameters with user-specified ranges, suggested by D. Barrow
563  std::map<std::string, std::vector<double>> scanRanges;
564  const bool isScanRanges = GetScanRange(scanRanges);
565 
566  double nSigma = GetFromManager<int>(fitMan->raw()["LLHScan"]["LLHScanSigma"], 1., __FILE__, __LINE__);
567  bool IsPCA = cov->IsPCA();
568 
569  // Get the parameter name
570  std::string name = cov->GetParFancyName(i);
571  if (IsPCA) name += "_PCA";
572 
573  // Get the parameter priors and bounds
574  CentralValue = cov->GetParProp(i);
575  if (IsPCA) CentralValue = cov->GetPCAHandler()->GetParPropPCA(i);
576 
577  double prior = cov->GetParInit(i);
578  if (IsPCA) prior = cov->GetPCAHandler()->GetPreFitValuePCA(i);
579 
580  if (std::abs(CentralValue - prior) > 1e-10) {
581  MACH3LOG_INFO("For {} scanning around value {} rather than prior {}", name, CentralValue, prior);
582  }
583 
584  // Get the covariance matrix and do the +/- nSigma
585  // Set lower and upper bounds relative the CentralValue
586  // Set the parameter ranges between which LLH points are scanned
587  lower = CentralValue - nSigma*cov->GetDiagonalError(i);
588  upper = CentralValue + nSigma*cov->GetDiagonalError(i);
589  // If PCA, transform these parameter values to the PCA basis
590  if (IsPCA) {
591  lower = CentralValue - nSigma*std::sqrt((cov->GetPCAHandler()->GetEigenValues())(i));
592  upper = CentralValue + nSigma*std::sqrt((cov->GetPCAHandler()->GetEigenValues())(i));
593  MACH3LOG_INFO("eval {} = {:.2f}", i, cov->GetPCAHandler()->GetEigenValues()(i));
594  MACH3LOG_INFO("CV {} = {:.2f}", i, CentralValue);
595  MACH3LOG_INFO("lower {} = {:.2f}", i, lower);
596  MACH3LOG_INFO("upper {} = {:.2f}", i, upper);
597  MACH3LOG_INFO("nSigma = {:.2f}", nSigma);
598  }
599 
600  // Implementation suggested by D. Barrow
601  // If param ranges are specified in scanRanges node, extract it from there
602  if(isScanRanges){
603  // Find matching entries through std::maps
604  auto it = scanRanges.find(name);
605  if (it != scanRanges.end() && it->second.size() == 2) { //Making sure the range is has only two entries
606  lower = it->second[0];
607  upper = it->second[1];
608  MACH3LOG_INFO("Found matching param name for setting specified range for {}", name);
609  MACH3LOG_INFO("Range for {} = [{:.2f}, {:.2f}]", name, lower, upper);
610  }
611  }
612 
613  // Cross-section and flux parameters have boundaries that we scan between, check that these are respected in setting lower and upper variables
614  // This also applies for other parameters like osc, etc.
615  lower = std::max(lower, cov->GetLowerBound(i));
616  upper = std::min(upper, cov->GetUpperBound(i));
617  MACH3LOG_INFO("Scanning {} {} with {} steps, from [{:.2f} , {:.2f}], CV = {:.2f}", suffix, name, n_points, lower, upper, CentralValue);
618 }
bool GetScanRange(std::map< std::string, std::vector< double >> &scanRanges) const
YSP: Set up a mapping to store parameters with user-specified ranges, suggested by D....
Definition: FitterBase.cpp:523
double GetParPropPCA(const int i) const
Get current parameter value using PCA.
Definition: PCAHandler.h:174
double GetPreFitValuePCA(const int i) const
Get current parameter value using PCA.
Definition: PCAHandler.h:181
const TVectorD GetEigenValues() const
Get eigen values for all parameters, if you want for decomposed only parameters use GetEigenValuesMas...
Definition: PCAHandler.h:200
double GetUpperBound(const int i) const
Get upper parameter bound in which it is physically valid.
PCAHandler * GetPCAHandler() const
Get pointer for PCAHandler.
double GetLowerBound(const int i) const
Get lower parameter bound in which it is physically valid.
bool IsPCA() const
is PCA, can use to query e.g. LLH scans
double GetDiagonalError(const int i) const
Get diagonal error for ith parameter.
double GetParProp(const int i) const
Get proposed parameter value.

◆ GetScanRange()

bool FitterBase::GetScanRange ( std::map< std::string, std::vector< double >> &  scanRanges) const
protected

YSP: Set up a mapping to store parameters with user-specified ranges, suggested by D. Barrow.

Parameters
scanRangesA map with user specified parameter ranges

Definition at line 523 of file FitterBase.cpp.

523  {
524 // *************************
525  bool isScanRanges = false;
526  // YSP: Set up a mapping to store parameters with user-specified ranges, suggested by D. Barrow
527  if(fitMan->raw()["LLHScan"]["ScanRanges"]){
528  YAML::Node scanRangesList = fitMan->raw()["LLHScan"]["ScanRanges"];
529  for (auto it = scanRangesList.begin(); it != scanRangesList.end(); ++it) {
530  std::string itname = it->first.as<std::string>();
531  std::vector<double> itrange = it->second.as<std::vector<double>>();
532  // Set the mapping as param_name:param_range
533  scanRanges[itname] = itrange;
534  }
535  isScanRanges = true;
536  } else {
537  MACH3LOG_INFO("There are no user-defined parameter ranges, so I'll use default param bounds for LLH Scans");
538  }
539  return isScanRanges;
540 }

◆ GetStepScaleBasedOnLLHScan()

void FitterBase::GetStepScaleBasedOnLLHScan ( )

LLH scan is good first estimate of step scale.

Definition at line 887 of file FitterBase.cpp.

887  {
888 // *************************
889  TDirectory *Sample_LLH = outputFile->Get<TDirectory>("Sample_LLH");
890  MACH3LOG_INFO("Starting Get Step Scale Based On LLHScan");
891 
892  if(Sample_LLH->IsZombie())
893  {
894  MACH3LOG_WARN("Couldn't find Sample_LLH, it looks like LLH scan wasn't run, will do this now");
895  RunLLHScan();
896  }
897 
898  for (ParameterHandlerBase *cov : systematics)
899  {
900  const int npars = cov->GetNumParams();
901  std::vector<double> StepScale(npars);
902  for (int i = 0; i < npars; ++i)
903  {
904  std::string name = cov->GetParFancyName(i);
905 
906  StepScale[i] = cov->GetIndivStepScale(i);
907  TH1D* LLHScan = Sample_LLH->Get<TH1D>((name+"_sam").c_str());
908  if(LLHScan == nullptr)
909  {
910  MACH3LOG_WARN("Couldn't find LLH scan, for {}, skipping", name);
911  continue;
912  }
913  const double LLH_val = std::max(LLHScan->GetBinContent(1), LLHScan->GetBinContent(LLHScan->GetNbinsX()));
914  //If there is no sensitivity leave it
915  if(LLH_val < 0.001) continue;
916 
917  // EM: assuming that the likelihood is gaussian, approximate sigma value is given by variation/sqrt(-2LLH)
918  // can evaluate this at any point, simple to evaluate it in the first bin of the LLH scan
919  // KS: We assume variation is 1 sigma, each dial has different scale so it becomes faff...
920  const double Var = 1.;
921  const double approxSigma = TMath::Abs(Var)/std::sqrt(LLH_val);
922 
923  // Based on Ewan comment I just took the 1sigma width from the LLH, assuming it was Gaussian, but then had to also scale by 2.38/sqrt(N_params)
924  const double NewStepScale = approxSigma * 2.38/std::sqrt(npars);
925  StepScale[i] = NewStepScale;
926  MACH3LOG_DEBUG("Sigma: {}", approxSigma);
927  MACH3LOG_DEBUG("optimal Step Size: {}", NewStepScale);
928  }
929  cov->SetIndivStepScale(StepScale);
930  cov->SaveUpdatedMatrixConfig();
931  }
932 }
void RunLLHScan()
Perform a 1D likelihood scan.
Definition: FitterBase.cpp:622
Base class responsible for handling of systematic error parameters. Capable of using PCA or using ada...

◆ PrepareOutput()

void FitterBase::PrepareOutput ( )
protected

Prepare the output file.

Definition at line 153 of file FitterBase.cpp.

153  {
154 // *******************
155  if(OutputPrepared) return;
156  //MS: Check if we are fitting the test likelihood, rather than T2K likelihood, and only setup T2K output if not
157  if(!fTestLikelihood)
158  {
159  // Check that we have added samples
160  if (!samples.size()) {
161  MACH3LOG_CRITICAL("No samples Found! Is this really what you wanted?");
162  //throw MaCh3Exception(__FILE__ , __LINE__ );
163  }
164 
165  // Do we want to save proposal? This will break plotting scripts and is heave for disk space and step time. Only use when debugging
166  bool SaveProposal = GetFromManager<bool>(fitMan->raw()["General"]["SaveProposal"], false, __FILE__ , __LINE__);
167 
168  if(SaveProposal) MACH3LOG_INFO("Will save in the chain proposal parameters and LogL");
169  // Prepare the output trees
170  for (ParameterHandlerBase *cov : systematics) {
171  cov->SetBranches(*outTree, SaveProposal);
172  }
173 
174  outTree->Branch("LogL", &logLCurr, "LogL/D");
175  if(SaveProposal) outTree->Branch("LogLProp", &logLProp, "LogLProp/D");
176  outTree->Branch("accProb", &accProb, "accProb/D");
177  outTree->Branch("step", &step, "step/i");
178  outTree->Branch("stepTime", &stepTime, "stepTime/D");
179 
180  // Store individual likelihood components
181  // Using double means double as large output file!
182  sample_llh.resize(samples.size());
183  syst_llh.resize(systematics.size());
184 
185  for (size_t i = 0; i < samples.size(); ++i) {
186  std::stringstream oss, oss2;
187  oss << "LogL_sample_" << i;
188  oss2 << oss.str() << "/D";
189  outTree->Branch(oss.str().c_str(), &sample_llh[i], oss2.str().c_str());
190  }
191 
192  for (size_t i = 0; i < systematics.size(); ++i) {
193  std::stringstream oss, oss2;
194  oss << "LogL_systematic_" << systematics[i]->GetName();
195  oss2 << oss.str() << "/D";
196  outTree->Branch(oss.str().c_str(), &syst_llh[i], oss2.str().c_str());
197  }
198  }
199  else
200  {
201  outTree->Branch("LogL", &logLCurr, "LogL/D");
202  outTree->Branch("accProb", &accProb, "accProb/D");
203  outTree->Branch("step", &step, "step/i");
204  outTree->Branch("stepTime", &stepTime, "stepTime/D");
205  }
206 
207  MACH3LOG_INFO("-------------------- Starting MCMC --------------------");
208  #ifdef DEBUG
209  if (debug) {
210  debugFile << "----- Starting MCMC -----" << std::endl;
211  }
212  #endif
213  // Time the progress
214 
215 
216  clock->Start();
217 
218 
219  OutputPrepared = true;
220 }
#define MACH3LOG_CRITICAL
Definition: MaCh3Logger.h:38
double logLProp
proposed likelihood
Definition: FitterBase.h:115
double accProb
current acceptance prob
Definition: FitterBase.h:117
std::vector< double > sample_llh
store the llh breakdowns
Definition: FitterBase.h:124
double stepTime
Time of single step.
Definition: FitterBase.h:141
double logLCurr
current likelihood
Definition: FitterBase.h:113
std::vector< double > syst_llh
systematic llh breakdowns
Definition: FitterBase.h:126

◆ ProcessMCMC()

void FitterBase::ProcessMCMC ( )
protected

Process MCMC output.

Definition at line 414 of file FitterBase.cpp.

414  {
415 // *******************
416  if (fitMan == nullptr) return;
417 
418  // Process the MCMC
419  if (GetFromManager<bool>(fitMan->raw()["General"]["ProcessMCMC"], false, __FILE__ , __LINE__)){
420  // Make the processor
421  MCMCProcessor Processor(std::string(outputFile->GetName()));
422 
423  Processor.Initialise();
424  // Make the TVectorD pointers which hold the processed output
425  TVectorD *Central = nullptr;
426  TVectorD *Errors = nullptr;
427  TVectorD *Central_Gauss = nullptr;
428  TVectorD *Errors_Gauss = nullptr;
429  TVectorD *Peaks = nullptr;
430 
431  // Make the postfit
432  Processor.GetPostfit(Central, Errors, Central_Gauss, Errors_Gauss, Peaks);
433  Processor.DrawPostfit();
434 
435  // Make the TMatrix pointers which hold the processed output
436  TMatrixDSym *Covariance = nullptr;
437  TMatrixDSym *Correlation = nullptr;
438 
439  // Make the covariance matrix
440  Processor.GetCovariance(Covariance, Correlation);
441  Processor.DrawCovariance();
442 
443  std::vector<TString> BranchNames = Processor.GetBranchNames();
444 
445  // Re-open the TFile
446  if (!outputFile->IsOpen()) {
447  MACH3LOG_INFO("Opening output again to update with means..");
448  outputFile = new TFile(Get<std::string>(fitMan->raw()["General"]["OutputFile"], __FILE__, __LINE__).c_str(), "UPDATE");
449  }
450  Central->Write("PDF_Means");
451  Errors->Write("PDF_Errors");
452  Central_Gauss->Write("Gauss_Means");
453  Errors_Gauss->Write("Errors_Gauss");
454  Covariance->Write("Covariance");
455  Correlation->Write("Correlation");
456  }
457 }
std::vector< TString > BranchNames
Definition: RHat.cpp:54
Class responsible for processing MCMC chains, performing diagnostics, generating plots,...
Definition: MCMCProcessor.h:58

◆ Run2DLLHScan()

void FitterBase::Run2DLLHScan ( )

Perform a 2D likelihood scan.

Warning
This operation may take a significant amount of time, especially for complex models.

Definition at line 936 of file FitterBase.cpp.

936  {
937 // *************************
938  // Save the settings into the output file
939  SaveSettings();
940 
941  MACH3LOG_INFO("Starting 2D LLH Scan");
942 
943  TDirectory *Sample_2DLLH = outputFile->mkdir("Sample_2DLLH");
944  auto SkipVector = GetFromManager<std::vector<std::string>>(fitMan->raw()["LLHScan"]["LLHScanSkipVector"], {}, __FILE__ , __LINE__);;
945 
946  // Number of points we do for each LLH scan
947  const int n_points = GetFromManager<int>(fitMan->raw()["LLHScan"]["2DLLHScanPoints"], 20, __FILE__ , __LINE__);
948  // We print 5 reweights
949  const int countwidth = int(double(n_points)/double(5));
950 
951  // Loop over the covariance classes
952  for (ParameterHandlerBase *cov : systematics)
953  {
954  // Scan over all the parameters
955  // Get the number of parameters
956  int npars = cov->GetNumParams();
957  bool IsPCA = cov->IsPCA();
958  if (IsPCA) npars = cov->GetNParameters();
959 
960  for (int i = 0; i < npars; ++i)
961  {
962  std::string name_x = cov->GetParFancyName(i);
963  if (IsPCA) name_x += "_PCA";
964  // Get the parameter central and bounds
965  double central_x, lower_x, upper_x;
966  GetParameterScanRange(cov, i, central_x, lower_x, upper_x, n_points, "X");
967 
968  // KS: Check if we want to skip this parameter
969  if(CheckSkipParameter(SkipVector, name_x)) continue;
970 
971  for (int j = 0; j < i; ++j)
972  {
973  std::string name_y = cov->GetParFancyName(j);
974  if (IsPCA) name_y += "_PCA";
975  // KS: Check if we want to skip this parameter
976  if(CheckSkipParameter(SkipVector, name_y)) continue;
977 
978  // Get the parameter central and bounds
979  double central_y, lower_y, upper_y;
980  GetParameterScanRange(cov, j, central_y, lower_y, upper_y, n_points, "Y");
981 
982  auto hScanSam = std::make_unique<TH2D>((name_x + "_" + name_y + "_sam").c_str(), (name_x + "_" + name_y + "_sam").c_str(),
983  n_points, lower_x, upper_x, n_points, lower_y, upper_y);
984  hScanSam->SetDirectory(nullptr);
985  hScanSam->GetXaxis()->SetTitle(name_x.c_str());
986  hScanSam->GetYaxis()->SetTitle(name_y.c_str());
987  hScanSam->GetZaxis()->SetTitle("2LLH_sam");
988 
989  // Scan over the parameter space
990  for (int x = 0; x < n_points; ++x)
991  {
992  if (x % countwidth == 0)
993  MaCh3Utils::PrintProgressBar(x, n_points);
994 
995  for (int y = 0; y < n_points; ++y)
996  {
997  // For PCA we have to do it differently
998  if (IsPCA) {
999  cov->GetPCAHandler()->SetParPropPCA(i, hScanSam->GetXaxis()->GetBinCenter(x+1));
1000  cov->GetPCAHandler()->SetParPropPCA(j, hScanSam->GetYaxis()->GetBinCenter(y+1));
1001  } else {
1002  // Set the parameter
1003  cov->SetParProp(i, hScanSam->GetXaxis()->GetBinCenter(x+1));
1004  cov->SetParProp(j, hScanSam->GetYaxis()->GetBinCenter(y+1));
1005  }
1006  // Reweight the MC
1007  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs) {
1008  samples[ivs]->Reweight();
1009  }
1010 
1011  // Get the -log L likelihoods
1012  double samplellh = 0;
1013  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs) {
1014  samplellh += samples[ivs]->GetLikelihood();
1015  }
1016  hScanSam->SetBinContent(x+1, y+1, 2*samplellh);
1017  }// end loop over y points
1018  } // end loop over x points
1019 
1020  Sample_2DLLH->cd();
1021  hScanSam->Write();
1022  // Reset the parameters to their central central values
1023  if (IsPCA) {
1024  cov->GetPCAHandler()->SetParPropPCA(i, central_x);
1025  cov->GetPCAHandler()->SetParPropPCA(j, central_y);
1026  } else {
1027  cov->SetParProp(i, central_x);
1028  cov->SetParProp(j, central_y);
1029  }
1030  } //end loop over systematics y
1031  }//end loop over systematics X
1032  }//end loop covariance classes
1033  Sample_2DLLH->Write();
1034  delete Sample_2DLLH;
1035 }
bool CheckSkipParameter(const std::vector< std::string > &SkipVector, const std::string &ParamName) const
KS: Check whether we want to skip parameter using skip vector.
Definition: FitterBase.cpp:543
void SaveSettings()
Save the settings that the MCMC was run with.
Definition: FitterBase.cpp:79
void GetParameterScanRange(const ParameterHandlerBase *cov, const int i, double &CentralValue, double &lower, double &upper, const int n_points, const std::string &suffix="") const
Helper function to get parameter scan range, central value.
Definition: FitterBase.cpp:559
void PrintProgressBar(const Long64_t Done, const Long64_t All)
KS: Simply print progress bar.
Definition: Monitor.cpp:228

◆ RunLLHMap()

void FitterBase::RunLLHMap ( )

Perform a general multi-dimensional likelihood scan.

Definition at line 1039 of file FitterBase.cpp.

1039  {
1040 // *************************
1041  // Save the settings into the output file
1042  SaveSettings();
1043 
1044  MACH3LOG_INFO("Starting {}", __func__);
1045 
1046  //KS: Turn it on if you want LLH scan for each ND sample separately, which increase time significantly but can be useful for validating new samples or dials.
1047  bool PlotLLHScanBySample = GetFromManager<bool>(fitMan->raw()["LLHScan"]["LLHScanBySample"], false, __FILE__ , __LINE__);
1048  auto ParamsOfInterest = GetFromManager<std::vector<std::string>>(fitMan->raw()["LLHScan"]["LLHParameters"], {}, __FILE__, __LINE__);
1049 
1050  if(ParamsOfInterest.empty()) {
1051  MACH3LOG_WARN("There were no LLH parameters of interest specified to run the LLHMap! LLHMap will not run at all ...");
1052  return;
1053  }
1054 
1055  // Parameters IDs within the covariance objects
1056  // ParamsCovIDs = {Name, CovObj, IDinCovObj}
1057  std::vector<std::tuple<std::string, ParameterHandlerBase*, int>> ParamsCovIDs;
1058  for(auto& p : ParamsOfInterest) {
1059  bool found = false;
1060  for(auto cov : systematics) {
1061  for(int c = 0; c < cov->GetNumParams(); ++c) {
1062  if(cov->GetParName(c) == p || cov->GetParFancyName(c) == p) {
1063  bool add = true;
1064  for(auto& pc : ParamsCovIDs) {
1065  if(std::get<1>(pc) == cov && std::get<2>(pc) == c)
1066  {
1067  MACH3LOG_WARN("Parameter {} as {}({}) listed multiple times for LLHMap, omitting and using only once!", p, cov->GetName(), c);
1068  add = false;
1069  break;
1070  }
1071  }
1072 
1073  if(add)
1074  ParamsCovIDs.push_back(std::make_tuple(p, cov, c));
1075 
1076  found = true;
1077  break;
1078  }
1079  }
1080  if(found)
1081  break;
1082  }
1083  if(found)
1084  MACH3LOG_INFO("Parameter {} found in {} at an index {}.", p, std::get<1>(ParamsCovIDs.back())->GetName(), std::get<2>(ParamsCovIDs.back()));
1085  else
1086  MACH3LOG_WARN("Parameter {} not found in any of the systematic covariance objects. Will not scan over this one!", p);
1087  }
1088 
1089  // ParamsRanges["parameter"] = {nPoints, {low, high}}
1090  std::map<std::string, std::pair<int, std::pair<double, double>>> ParamsRanges;
1091 
1092  MACH3LOG_INFO("======================================================================================");
1093  MACH3LOG_INFO("Performing a general multi-dimensional LogL map scan over following parameters ranges:");
1094  MACH3LOG_INFO("======================================================================================");
1095  unsigned long TotalPoints = 1;
1096 
1097  double nSigma = GetFromManager<int>(fitMan->raw()["LLHScan"]["LLHScanSigma"], 1., __FILE__, __LINE__);
1098 
1099  // TN: Setting up the scan ranges might look like a re-implementation of the
1100  // FitterBase::GetScanRange, but I guess the use-case here is a bit different.
1101  // Anyway, just in case, we can discuss and rewrite to everyone's liking!
1102  for(auto& p : ParamsCovIDs) {
1103  // Auxiliary vars to help readability
1104  std::string name = std::get<0>(p);
1105  int i = std::get<2>(p);
1106  ParameterHandlerBase* cov = std::get<1>(p);
1107 
1108  ParamsRanges[name].first = GetFromManager<int>(fitMan->raw()["LLHScan"]["LLHScanPoints"], 20, __FILE__, __LINE__);
1109  if(CheckNodeExists(fitMan->raw(),"LLHScan","ScanPoints"))
1110  ParamsRanges[name].first = GetFromManager<int>(fitMan->raw()["LLHScan"]["ScanPoints"][name], ParamsRanges[name].first, __FILE__, __LINE__);
1111 
1112  // Get the parameter priors and bounds
1113  double CentralValue = cov->GetParProp(i);
1114 
1115  bool IsPCA = cov->IsPCA();
1116  if (IsPCA)
1117  CentralValue = cov->GetPCAHandler()->GetParPropPCA(i);
1118 
1119  double prior = cov->GetParInit(i);
1120  if (IsPCA)
1121  prior = cov->GetPCAHandler()->GetPreFitValuePCA(i);
1122 
1123  if (std::abs(CentralValue - prior) > 1e-10) {
1124  MACH3LOG_INFO("For {} scanning around value {} rather than prior {}", name, CentralValue, prior);
1125  }
1126  // Get the covariance matrix and do the +/- nSigma
1127  // Set lower and upper bounds relative the CentralValue
1128  // Set the parameter ranges between which LLH points are scanned
1129  double lower = CentralValue - nSigma*cov->GetDiagonalError(i);
1130  double upper = CentralValue + nSigma*cov->GetDiagonalError(i);
1131  // If PCA, transform these parameter values to the PCA basis
1132  if (IsPCA) {
1133  lower = CentralValue - nSigma*std::sqrt((cov->GetPCAHandler()->GetEigenValues())(i));
1134  upper = CentralValue + nSigma*std::sqrt((cov->GetPCAHandler()->GetEigenValues())(i));
1135  MACH3LOG_INFO("eval {} = {:.2f}", i, cov->GetPCAHandler()->GetEigenValues()(i));
1136  MACH3LOG_INFO("CV {} = {:.2f}", i, CentralValue);
1137  MACH3LOG_INFO("lower {} = {:.2f}", i, lower);
1138  MACH3LOG_INFO("upper {} = {:.2f}", i, upper);
1139  MACH3LOG_INFO("nSigma = {:.2f}", nSigma);
1140  }
1141 
1142  ParamsRanges[name].second = {lower,upper};
1143 
1144  if(CheckNodeExists(fitMan->raw(),"LLHScan","ScanRanges"))
1145  ParamsRanges[name].second = GetFromManager<std::pair<double,double>>(fitMan->raw()["LLHScan"]["ScanRanges"][name], ParamsRanges[name].second, __FILE__, __LINE__);
1146 
1147  MACH3LOG_INFO("{} from {:.4f} to {:.4f} with a {:.5f} step ({} points total)",
1148  name, ParamsRanges[name].second.first, ParamsRanges[name].second.second,
1149  (ParamsRanges[name].second.second - ParamsRanges[name].second.first)/(ParamsRanges[name].first - 1.),
1150  ParamsRanges[name].first);
1151 
1152  TotalPoints *= ParamsRanges[name].first;
1153  }
1154 
1155  // TN: Waiting for C++ 20 std::format() function
1156  MACH3LOG_INFO("In total, looping over {} points, from {} parameters. Estimates for run time:", TotalPoints, ParamsCovIDs.size());
1157  MACH3LOG_INFO(" 1 s per point = {} hours", double(TotalPoints)/3600.);
1158  MACH3LOG_INFO(" 0.1 s per point = {} hours", double(TotalPoints)/36000.);
1159  MACH3LOG_INFO("0.01 s per point = {} hours", double(TotalPoints)/360000.);
1160  MACH3LOG_INFO("==================================================================================");
1161 
1162  const int countwidth = int(double(TotalPoints)/double(20));
1163 
1164  // Tree to store LogL values
1165  auto LLHMap = new TTree("llhmap", "LLH Map");
1166 
1167  std::vector<double> CovLogL(systematics.size());
1168  for(unsigned int ivc = 0; ivc < systematics.size(); ++ivc )
1169  {
1170  std::string NameTemp = systematics[ivc]->GetName();
1171  NameTemp = NameTemp.substr(0, NameTemp.find("_cov")) + "_LLH";
1172  LLHMap->Branch(NameTemp.c_str(), &CovLogL[ivc]);
1173  }
1174 
1175  std::vector<double> SampleClassLogL(samples.size());
1176  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs )
1177  {
1178  std::string NameTemp = samples[ivs]->GetName()+"_LLH";
1179  LLHMap->Branch(NameTemp.c_str(), &SampleClassLogL[ivs]);
1180  }
1181 
1182  double SampleLogL, TotalLogL;
1183  LLHMap->Branch("Sample_LLH", &SampleLogL);
1184  LLHMap->Branch("Total_LLH", &TotalLogL);
1185 
1186  std::vector<double>SampleSplitLogL;
1187  if(PlotLLHScanBySample)
1188  {
1189  SampleSplitLogL.resize(TotalNSamples);
1190  int SampleIterator = 0;
1191  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs )
1192  {
1193  for(int is = 0; is < samples[ivs]->GetNsamples(); ++is )
1194  {
1195  std::string NameTemp =samples[ivs]->GetSampleTitle(is)+"_LLH";
1196  LLHMap->Branch(NameTemp.c_str(), &SampleSplitLogL[SampleIterator]);
1197  SampleIterator++;
1198  }
1199  }
1200  }
1201 
1202  std::vector<double> ParamsValues(ParamsCovIDs.size());
1203  for(unsigned int i=0; i < ParamsCovIDs.size(); ++i)
1204  LLHMap->Branch(std::get<0>(ParamsCovIDs[i]).c_str(), &ParamsValues[i]);
1205 
1206  // Setting up the scan
1207  // Starting at index {0,0,0,...}
1208  std::vector<unsigned long> idx(ParamsCovIDs.size(), 0);
1209 
1210  // loop over scanned points sp
1211  for(unsigned long sp = 0; sp < TotalPoints; ++sp)
1212  {
1213  // At each point need to find the indices and test values to calculate LogL
1214  for(unsigned int n = 0; n < ParamsCovIDs.size(); ++n)
1215  {
1216  // Auxiliaries
1217  std::string name = std::get<0>(ParamsCovIDs[n]);
1218  int points = ParamsRanges[name].first;
1219  double low = ParamsRanges[name].second.first;
1220  double high = ParamsRanges[name].second.second;
1221 
1222  // Find the n-th index of the sp-th scanned point
1223  unsigned long dev = 1;
1224  for(unsigned int m = 0; m <= n; ++m)
1225  dev *= ParamsRanges[std::get<0>(ParamsCovIDs[m])].first;
1226 
1227  idx[n] = sp % dev;
1228  if (n > 0)
1229  idx[n] = idx[n] / ( dev / points );
1230 
1231  // Parameter test value = low + ( high - low ) * idx / ( #points - 1 )
1232  ParamsValues[n] = low + (high-low) * double(idx[n])/double(points-1);
1233 
1234  // Now set the covariance objects
1235  // Auxiliary
1236  ParameterHandlerBase* cov = std::get<1>(ParamsCovIDs[n]);
1237  int i = std::get<2>(ParamsCovIDs[n]);
1238 
1239  if(cov->IsPCA())
1240  cov->GetPCAHandler()->SetParPropPCA(i, ParamsValues[n]);
1241  else
1242  cov->SetParProp(i, ParamsValues[n]);
1243  }
1244 
1245  // Reweight samples and calculate LogL
1246  TotalLogL = .0;
1247  SampleLogL = .0;
1248 
1249  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs)
1250  samples[ivs]->Reweight();
1251 
1252  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs)
1253  {
1254  SampleClassLogL[ivs] = 2.*samples[ivs]->GetLikelihood();
1255  SampleLogL += SampleClassLogL[ivs];
1256  }
1257  TotalLogL += SampleLogL;
1258 
1259  // CovObjs LogL
1260  for(unsigned int ivc = 0; ivc < systematics.size(); ++ivc )
1261  {
1262  CovLogL[ivc] = 2.*systematics[ivc]->GetLikelihood();
1263  TotalLogL += CovLogL[ivc];
1264  }
1265 
1266  if(PlotLLHScanBySample)
1267  {
1268  int SampleIterator = 0;
1269  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs )
1270  {
1271  for(int is = 0; is < samples[ivs]->GetNsamples(); ++is)
1272  {
1273  SampleSplitLogL[SampleIterator] = 2.*samples[ivs]->GetSampleLikelihood(is);
1274  SampleIterator++;
1275  }
1276  }
1277  }
1278 
1279  LLHMap->Fill();
1280 
1281  if (sp % countwidth == 0)
1282  MaCh3Utils::PrintProgressBar(sp, TotalPoints);
1283  }
1284 
1285  outputFile->cd();
1286  LLHMap->Write();
1287 }
bool CheckNodeExists(const YAML::Node &node, Args... args)
KS: Wrapper function to call the recursive helper.
Definition: YamlHelper.h:60
void SetParPropPCA(const int i, const double value)
Set proposed value for parameter in PCA base.
Definition: PCAHandler.h:156
void SetParProp(const int i, const double val)
Set proposed parameter value.

◆ RunLLHScan()

void FitterBase::RunLLHScan ( )

Perform a 1D likelihood scan.

Definition at line 622 of file FitterBase.cpp.

622  {
623 // *************************
624  // Save the settings into the output file
625  SaveSettings();
626 
627  MACH3LOG_INFO("Starting {}", __func__);
628 
629  //KS: Turn it on if you want LLH scan for each ND sample separately, which increase time significantly but can be useful for validating new samples or dials.
630  bool PlotLLHScanBySample = GetFromManager<bool>(fitMan->raw()["LLHScan"]["LLHScanBySample"], false, __FILE__ , __LINE__);
631  auto SkipVector = GetFromManager<std::vector<std::string>>(fitMan->raw()["LLHScan"]["LLHScanSkipVector"], {}, __FILE__ , __LINE__);
632 
633  // Now finally get onto the LLH scan stuff
634  // Very similar code to MCMC but never start MCMC; just scan over the parameter space
635  std::vector<TDirectory *> Cov_LLH(systematics.size());
636  for(unsigned int ivc = 0; ivc < systematics.size(); ++ivc )
637  {
638  std::string NameTemp = systematics[ivc]->GetName();
639  NameTemp = NameTemp.substr(0, NameTemp.find("_cov")) + "_LLH";
640  Cov_LLH[ivc] = outputFile->mkdir(NameTemp.c_str());
641  }
642 
643  std::vector<TDirectory *> SampleClass_LLH(samples.size());
644  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs )
645  {
646  std::string NameTemp = samples[ivs]->GetName();
647  SampleClass_LLH[ivs] = outputFile->mkdir(NameTemp.c_str());
648  }
649 
650  TDirectory *Sample_LLH = outputFile->mkdir("Sample_LLH");
651  TDirectory *Total_LLH = outputFile->mkdir("Total_LLH");
652 
653  std::vector<TDirectory *>SampleSplit_LLH;
654  if(PlotLLHScanBySample)
655  {
656  SampleSplit_LLH.resize(TotalNSamples);
657  int SampleIterator = 0;
658  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs )
659  {
660  for(int is = 0; is < samples[ivs]->GetNsamples(); ++is )
661  {
662  SampleSplit_LLH[SampleIterator] = outputFile->mkdir((samples[ivs]->GetSampleTitle(is)+ "_LLH").c_str());
663  SampleIterator++;
664  }
665  }
666  }
667  // Number of points we do for each LLH scan
668  const int n_points = GetFromManager<int>(fitMan->raw()["LLHScan"]["LLHScanPoints"], 100, __FILE__ , __LINE__);
669 
670  // We print 5 reweights
671  const int countwidth = int(double(n_points)/double(5));
672 
673  // Loop over the covariance classes
674  for (ParameterHandlerBase *cov : systematics)
675  {
676  // Scan over all the parameters
677  // Get the number of parameters
678  int npars = cov->GetNumParams();
679  bool IsPCA = cov->IsPCA();
680  if (IsPCA) npars = cov->GetNParameters();
681  for (int i = 0; i < npars; ++i)
682  {
683  // Get the parameter name
684  std::string name = cov->GetParFancyName(i);
685  if (IsPCA) name += "_PCA";
686  // KS: Check if we want to skip this parameter
687  if(CheckSkipParameter(SkipVector, name)) continue;
688  // Get the parameter central and bounds
689  double CentralValue, lower, upper;
690  GetParameterScanRange(cov, i, CentralValue, lower, upper, n_points);
691  // Make the TH1D
692  auto hScan = std::make_unique<TH1D>((name + "_full").c_str(), (name + "_full").c_str(), n_points, lower, upper);
693  hScan->SetTitle((std::string("2LLH_full, ") + name + ";" + name + "; -2(ln L_{sample} + ln L_{xsec+flux} + ln L_{det})").c_str());
694  hScan->SetDirectory(nullptr);
695 
696  auto hScanSam = std::make_unique<TH1D>((name + "_sam").c_str(), (name + "_sam").c_str(), n_points, lower, upper);
697  hScanSam->SetTitle((std::string("2LLH_sam, ") + name + ";" + name + "; -2(ln L_{sample})").c_str());
698  hScanSam->SetDirectory(nullptr);
699 
700  std::vector<std::unique_ptr<TH1D>> hScanSample(samples.size());
701  std::vector<double> nSamLLH(samples.size(), 0.0);
702  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs )
703  {
704  std::string NameTemp = samples[ivs]->GetName();
705  hScanSample[ivs] = std::make_unique<TH1D>((name+"_"+NameTemp).c_str(), (name+"_" + NameTemp).c_str(), n_points, lower, upper);
706  hScanSample[ivs]->SetDirectory(nullptr);
707  hScanSample[ivs]->SetTitle(("2LLH_" + NameTemp + ", " + name + ";" + name + "; -2(ln L_{" + NameTemp +"})").c_str());
708  }
709 
710  std::vector<std::unique_ptr<TH1D>> hScanCov(systematics.size());
711  std::vector<double> nCovLLH(systematics.size(), 0.0);
712  for(unsigned int ivc = 0; ivc < systematics.size(); ++ivc )
713  {
714  std::string NameTemp = systematics[ivc]->GetName();
715  NameTemp = NameTemp.substr(0, NameTemp.find("_cov"));
716  hScanCov[ivc] = std::make_unique<TH1D>((name + "_" + NameTemp).c_str(), (name + "_" + NameTemp).c_str(), n_points, lower, upper);
717  hScanCov[ivc]->SetDirectory(nullptr);
718  hScanCov[ivc]->SetTitle(("2LLH_" + NameTemp + ", " + name + ";" + name + "; -2(ln L_{" + NameTemp +"})").c_str());
719  }
720 
721  std::vector<std::unique_ptr<TH1D>> hScanSamSplit;
722  std::vector<double> sampleSplitllh;
723  if(PlotLLHScanBySample)
724  {
725  int SampleIterator = 0;
726  hScanSamSplit.resize(TotalNSamples);
727  sampleSplitllh.resize(TotalNSamples);
728  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs )
729  {
730  for(int is = 0; is < samples[ivs]->GetNsamples(); ++is )
731  {
732  auto histName = name + samples[ivs]->GetSampleTitle(is);
733  auto histTitle = std::string("2LLH_sam, ") + name + ";" + name + "; -2(ln L_{sample})";
734  hScanSamSplit[SampleIterator] = std::make_unique<TH1D>(histName.c_str(), histTitle.c_str(), n_points, lower, upper);
735  hScanSamSplit[SampleIterator]->SetDirectory(nullptr);
736  SampleIterator++;
737  }
738  }
739  }
740 
741  // Scan over the parameter space
742  for (int j = 0; j < n_points; ++j)
743  {
744  if (j % countwidth == 0)
745  MaCh3Utils::PrintProgressBar(j, n_points);
746 
747  // For PCA we have to do it differently
748  if (IsPCA) {
749  cov->GetPCAHandler()->SetParPropPCA(i, hScan->GetBinCenter(j+1));
750  } else {
751  // Set the parameter
752  cov->SetParProp(i, hScan->GetBinCenter(j+1));
753  }
754 
755  // Reweight the MC
756  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs ){
757  samples[ivs]->Reweight();
758  }
759  //Total LLH
760  double totalllh = 0.;
761 
762  // Get the -log L likelihoods
763  double samplellh = 0.;
764 
765  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs ) {
766  nSamLLH[ivs] = samples[ivs]->GetLikelihood();
767  samplellh += nSamLLH[ivs];
768  }
769 
770  for(unsigned int ivc = 0; ivc < systematics.size(); ++ivc ) {
771  nCovLLH[ivc] = systematics[ivc]->GetLikelihood();
772  totalllh += nCovLLH[ivc];
773  }
774 
775  totalllh += samplellh;
776 
777  if(PlotLLHScanBySample)
778  {
779  int SampleIterator = 0;
780  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs )
781  {
782  for(int is = 0; is < samples[ivs]->GetNsamples(); ++is)
783  {
784  sampleSplitllh[SampleIterator] = samples[ivs]->GetSampleLikelihood(is);
785  SampleIterator++;
786  }
787  }
788  }
789 
790  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs ) {
791  hScanSample[ivs]->SetBinContent(j+1, 2*nSamLLH[ivs]);
792  }
793  for(unsigned int ivc = 0; ivc < systematics.size(); ++ivc ) {
794  hScanCov[ivc]->SetBinContent(j+1, 2*nCovLLH[ivc]);
795  }
796 
797  hScanSam->SetBinContent(j+1, 2*samplellh);
798  hScan->SetBinContent(j+1, 2*totalllh);
799 
800  if(PlotLLHScanBySample)
801  {
802  int SampleIterator = 0;
803  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs )
804  {
805  for(int is = 0; is < samples[ivs]->GetNsamples(); ++is)
806  {
807  hScanSamSplit[SampleIterator]->SetBinContent(j+1, 2*sampleSplitllh[SampleIterator]);
808  SampleIterator++;
809  }
810  }
811  }
812  }
813  for(unsigned int ivc = 0; ivc < systematics.size(); ++ivc )
814  {
815  Cov_LLH[ivc]->cd();
816  hScanCov[ivc]->Write();
817  }
818 
819  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs )
820  {
821  SampleClass_LLH[ivs]->cd();
822  hScanSample[ivs]->Write();
823  }
824  Sample_LLH->cd();
825  hScanSam->Write();
826  Total_LLH->cd();
827  hScan->Write();
828 
829  if(PlotLLHScanBySample)
830  {
831  int SampleIterator = 0;
832  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs )
833  {
834  for(int is = 0; is < samples[ivs]->GetNsamples(); ++is)
835  {
836  SampleSplit_LLH[SampleIterator]->cd();
837  hScanSamSplit[SampleIterator]->Write();
838  SampleIterator++;
839  }
840  }
841  }
842 
843  // Reset the parameters to their CentralValue central values
844  if (IsPCA) {
845  cov->GetPCAHandler()->SetParPropPCA(i, CentralValue);
846  } else {
847  cov->SetParProp(i, CentralValue);
848  }
849  }//end loop over systematics
850  }//end loop covariance classes
851 
852  for(unsigned int ivc = 0; ivc < systematics.size(); ++ivc )
853  {
854  Cov_LLH[ivc]->Write();
855  delete Cov_LLH[ivc];
856  }
857 
858  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs )
859  {
860  SampleClass_LLH[ivs]->Write();
861  delete SampleClass_LLH[ivs];
862  }
863 
864  Sample_LLH->Write();
865  delete Sample_LLH;
866 
867  Total_LLH->Write();
868  delete Total_LLH;
869 
870  if(PlotLLHScanBySample)
871  {
872  int SampleIterator = 0;
873  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs )
874  {
875  for(int is = 0; is < samples[ivs]->GetNsamples(); ++is )
876  {
877  SampleSplit_LLH[SampleIterator]->Write();
878  delete SampleSplit_LLH[SampleIterator];
879  SampleIterator++;
880  }
881  }
882  }
883 }

◆ RunMCMC()

virtual void FitterBase::RunMCMC ( )
pure virtual

The specific fitting algorithm implemented in this function depends on the derived class. It could be Markov Chain Monte Carlo (MCMC), MinuitFit, or another algorithm.

Implemented in PyLikelihoodFit, PyFitterBase, PSO, PredictiveThrower, MinuitFit, and MCMCBase.

◆ RunSigmaVar()

void FitterBase::RunSigmaVar ( )

Perform a 1D/2D sigma var for all samples.

Apply custom range to make easier comparison with p-theta

Definition at line 1393 of file FitterBase.cpp.

1393  {
1394 // *************************
1395  // Save the settings into the output file
1396  SaveSettings();
1397 
1398  bool plot_by_mode = GetFromManager<bool>(fitMan->raw()["SigmaVar"]["PlotByMode"], false);
1399  bool plot_by_channel = GetFromManager<bool>(fitMan->raw()["SigmaVar"]["PlotByChannel"], false);
1400  auto SkipVector = GetFromManager<std::vector<std::string>>(fitMan->raw()["SigmaVar"]["SkipVector"], {}, __FILE__ , __LINE__);
1401 
1402  if (plot_by_mode) MACH3LOG_INFO("Plotting by sample and mode");
1403  if (plot_by_channel) MACH3LOG_INFO("Plotting by sample and channel");
1404  if (!plot_by_mode && !plot_by_channel) MACH3LOG_INFO("Plotting by sample only");
1405  if (plot_by_mode && plot_by_channel) MACH3LOG_INFO("Plotting by sample, mode and channel");
1406 
1407  auto SigmaArray = GetFromManager<std::vector<double>>(fitMan->raw()["SigmaVar"]["SigmaArray"], {-3, -1, 0, 1, 3}, __FILE__ , __LINE__);
1408  if (std::find(SigmaArray.begin(), SigmaArray.end(), 0.0) == SigmaArray.end()) {
1409  MACH3LOG_ERROR(":: SigmaArray does not contain 0! Current contents: {} ::", fmt::join(SigmaArray, ", "));
1410  throw MaCh3Exception(__FILE__, __LINE__);
1411  }
1412 
1413  TDirectory* SigmaDir = outputFile->mkdir("SigmaVar");
1414  outputFile->cd();
1415 
1416  for (size_t s = 0; s < systematics.size(); ++s)
1417  {
1418  for(int i = 0; i < systematics[s]->GetNumParams(); i++)
1419  {
1420  std::string ParName = systematics[s]->GetParFancyName(i);
1421  // KS: Check if we want to skip this parameter
1422  if(CheckSkipParameter(SkipVector, ParName)) continue;
1423 
1424  MACH3LOG_INFO(":: Param {} ::", systematics[s]->GetParFancyName(i));
1425 
1426  TDirectory* ParamDir = SigmaDir->mkdir(ParName.c_str());
1427  ParamDir->cd();
1428 
1429  const double ParamCentralValue = systematics[s]->GetParProp(i);
1430  const double Prior = systematics[s]->GetParInit(i);
1431  const double ParamLower = systematics[s]->GetLowerBound(i);
1432  const double ParamUpper = systematics[s]->GetUpperBound(i);
1433 
1434  if (std::abs(ParamCentralValue - Prior) > 1e-10) {
1435  MACH3LOG_INFO("For {} scanning around value {} rather than prior {}", ParName, ParamCentralValue, Prior);
1436  }
1437 
1438  for(unsigned int iSample = 0; iSample < samples.size(); ++iSample)
1439  {
1440  auto* MaCh3Sample = samples[iSample];
1441  std::vector<TDirectory*> SampleDir(MaCh3Sample->GetNsamples());
1442  for (int SampleIndex = 0; SampleIndex < MaCh3Sample->GetNsamples(); ++SampleIndex) {
1443  SampleDir[SampleIndex] = ParamDir->mkdir(MaCh3Sample->GetSampleTitle(SampleIndex).c_str());
1444  }
1445 
1446  for (size_t j = 0; j < SigmaArray.size(); ++j) {
1447  double sigma = SigmaArray[j];
1448 
1449  double ParamShiftValue = ParamCentralValue + sigma * std::sqrt((*systematics[s]->GetCovMatrix())(i,i));
1450  ParamShiftValue = std::max(std::min(ParamShiftValue, ParamUpper), ParamLower);
1451 
1453  CustomRange(ParName, sigma, ParamShiftValue);
1454 
1455  MACH3LOG_INFO(" - set to {:<5.2f} ({:<2} sigma shift)", ParamShiftValue, sigma);
1456  systematics[s]->SetParProp(i, ParamShiftValue);
1457 
1458  std::ostringstream valStream;
1459  valStream << std::fixed << std::setprecision(2) << ParamShiftValue;
1460  std::string valueStr = valStream.str();
1461 
1462  std::ostringstream sigmaStream;
1463  sigmaStream << std::fixed << std::setprecision(2) << std::abs(sigma);
1464  std::string sigmaStr = sigmaStream.str();
1465 
1466  std::string suffix;
1467  if (sigma == 0) {
1468  suffix = "_" + ParName + "_nom_val_" + valueStr;
1469  } else {
1470  std::string sign = (sigma > 0) ? "p" : "n";
1471  suffix = "_" + ParName + "_sig_" + sign + sigmaStr + "_val_" + valueStr;
1472  }
1473 
1474  systematics[s]->SetParProp(i, ParamShiftValue);
1475  MaCh3Sample->Reweight();
1476 
1477  WriteHistogramsByMode(MaCh3Sample, suffix, plot_by_mode, plot_by_channel, SampleDir);
1478  }
1479  for (int subSampleIndex = 0; subSampleIndex < MaCh3Sample->GetNsamples(); ++subSampleIndex) {
1480  SampleDir[subSampleIndex]->Close();
1481  delete SampleDir[subSampleIndex];
1482  }
1483  ParamDir->cd();
1484  }
1485 
1486  systematics[s]->SetParProp(i, ParamCentralValue);
1487  MACH3LOG_INFO(" - set back to CV {:<5.2f}", ParamCentralValue);
1488  MACH3LOG_INFO("");
1489  ParamDir->Close();
1490  delete ParamDir;
1491  SigmaDir->cd();
1492  } // end loop over params
1493  } // end loop over systemics
1494  SigmaDir->Close();
1495  delete SigmaDir;
1496 
1497  outputFile->cd();
1498 }
void WriteHistogramsByMode(SampleHandlerBase *sample, const std::string &suffix, const bool by_mode, const bool by_channel, const std::vector< TDirectory * > &SampleDir)
Generic histogram writer - should make main code more palatable.
void CustomRange(const std::string &ParName, const double sigma, double &ParamShiftValue) const
For comparison with other fitting frameworks (like P-Theta) we usually have to apply different parame...

◆ SanitiseInputs()

void FitterBase::SanitiseInputs ( )
protected

Remove obsolete memory and make other checks before fit starts.

Todo:
consider expanding into ParmaterHandler and add more sanitisers

Definition at line 223 of file FitterBase.cpp.

223  {
224 // *******************
225  for (size_t i = 0; i < samples.size(); ++i) {
226  samples[i]->CleanMemoryBeforeFit();
227  }
228 }

◆ SaveOutput()

void FitterBase::SaveOutput ( )
protected

Save output and close files.

Definition at line 231 of file FitterBase.cpp.

231  {
232 // *******************
233  if(FileSaved) return;
234  //Stop Clock
235  clock->Stop();
236 
237  //KS: Some version of ROOT keep spamming about accessing already deleted object which is wrong and not helpful...
238  int originalErrorLevel = gErrorIgnoreLevel;
239  gErrorIgnoreLevel = kFatal;
240 
241  outputFile->cd();
242  outTree->Write();
243 
244  MACH3LOG_INFO("{} steps took {:.2e} seconds to complete. ({:.2e}s / step).", step - stepStart, clock->RealTime(), clock->RealTime() / static_cast<double>(step - stepStart));
245  MACH3LOG_INFO("{} steps were accepted.", accCount);
246  #ifdef DEBUG
247  if (debug)
248  {
249  debugFile << "\n\n" << step << " steps took " << clock->RealTime() << " seconds to complete. (" << clock->RealTime() / step << "s / step).\n" << accCount<< " steps were accepted." << std::endl;
250  debugFile.close();
251  }
252  #endif
253 
254  outputFile->Close();
255  FileSaved = true;
256 
257  gErrorIgnoreLevel = originalErrorLevel;
258 }

◆ SaveSettings()

void FitterBase::SaveSettings ( )
protected

Save the settings that the MCMC was run with.

Definition at line 79 of file FitterBase.cpp.

79  {
80 // *******************
81  if(SettingsSaved) return;
82 
83  outputFile->cd();
84 
85  TDirectory* MaCh3Version = outputFile->mkdir("MaCh3Engine");
86  MaCh3Version->cd();
87 
88  if (std::getenv("MaCh3_ROOT") == nullptr) {
89  MACH3LOG_ERROR("Need MaCh3_ROOT environment variable");
90  MACH3LOG_ERROR("Please remember about source bin/setup.MaCh3.sh");
91  throw MaCh3Exception(__FILE__ , __LINE__ );
92  }
93 
94  if (std::getenv("MACH3") == nullptr) {
95  MACH3LOG_ERROR("Need MACH3 environment variable");
96  throw MaCh3Exception(__FILE__ , __LINE__ );
97  }
98 
99  std::string header_path = std::string(std::getenv("MACH3"));
100  header_path += "/version.h";
101  FILE* file = fopen(header_path.c_str(), "r");
102  //KS: It is better to use experiment specific header file. If given experiment didn't provide it we gonna use one given by Core MaCh3.
103  if (!file) {
104  header_path = std::string(std::getenv("MaCh3_ROOT"));
105  header_path += "/version.h";
106  } else {
107  fclose(file);
108  }
109 
110  // EM: embed the cmake generated version.h file
111  TMacro versionHeader("version_header", "version_header");
112  versionHeader.ReadFile(header_path.c_str());
113  versionHeader.Write();
114 
115  if(GetName() == ""){
116  MACH3LOG_ERROR("Name of currently used algorithm is {}", GetName());
117  MACH3LOG_ERROR("Have you forgotten to modify AlgorithmName?");
118  throw MaCh3Exception(__FILE__ , __LINE__ );
119  }
120  TNamed Engine(GetName(), GetName());
121  Engine.Write(GetName().c_str());
122 
123  MaCh3Version->Write();
124  delete MaCh3Version;
125 
126  outputFile->cd();
127 
129 
130  MACH3LOG_WARN("\033[0;31mCurrent Total RAM usage is {:.2f} GB\033[0m", MaCh3Utils::getValue("VmRSS") / 1048576.0);
131  MACH3LOG_WARN("\033[0;31mOut of Total available RAM {:.2f} GB\033[0m", MaCh3Utils::getValue("MemTotal") / 1048576.0);
132 
133  MACH3LOG_INFO("#####Current Setup#####");
134  MACH3LOG_INFO("Number of covariances: {}", systematics.size());
135  for(unsigned int i = 0; i < systematics.size(); ++i)
136  MACH3LOG_INFO("{}: Cov name: {}, it has {} params", i, systematics[i]->GetName(), systematics[i]->GetNumParams());
137  MACH3LOG_INFO("Number of SampleHandlers: {}", samples.size());
138  for(unsigned int i = 0; i < samples.size(); ++i) {
139  MACH3LOG_INFO("{}: SampleHandler name: {}, it has {} samples",i , samples[i]->GetName(), samples[i]->GetNsamples());
140  for(int iSam = 0; iSam < samples[i]->GetNsamples(); ++iSam) {
141  MACH3LOG_INFO(" {}: Sample name: {}, with {} osc channels",iSam , samples[i]->GetSampleTitle(iSam), samples[i]->GetNOscChannels(iSam));
142  }
143  }
144  //TN: Have to close the folder in order to write it to disk before SaveOutput is called in the destructor
145  CovFolder->Close();
146  SampleFolder->Close();
147 
148  SettingsSaved = true;
149 }
std::string GetName() const
Get name of class.
Definition: FitterBase.h:70
void SaveSettings(TFile *const OutputFile) const
Add manager useful information's to TFile, in most cases to Fitter.
Definition: Manager.cpp:49
int getValue(const std::string &Type)
CW: Get info like RAM.
Definition: Monitor.cpp:251

◆ StartFromPreviousFit()

void FitterBase::StartFromPreviousFit ( const std::string &  FitName)
virtual

Allow to start from previous fit/chain.

Parameters
FitNameName of previous chain

Reimplemented in MCMCBase.

Definition at line 348 of file FitterBase.cpp.

348  {
349 // *******************
350  MACH3LOG_INFO("Getting starting position from {}", FitName);
351  TFile *infile = M3::Open(FitName, "READ", __FILE__, __LINE__);
352  TTree *posts = infile->Get<TTree>("posteriors");
353  unsigned int step_val = 0;
354  double log_val = M3::_LARGE_LOGL_;
355  posts->SetBranchAddress("step",&step_val);
356  posts->SetBranchAddress("LogL",&log_val);
357 
358  for (size_t s = 0; s < systematics.size(); ++s)
359  {
360  TDirectory* CovarianceFolder = infile->Get<TDirectory>("CovarianceFolder");
361 
362  std::string ConfigName = "Config_" + systematics[s]->GetName();
363  TMacro *ConfigCov = CovarianceFolder->Get<TMacro>(ConfigName.c_str());
364  // KS: Not every covariance uses yaml, if it uses yaml make sure they are identical
365  if (ConfigCov != nullptr) {
366  // Config which was in MCMC from which we are starting
367  YAML::Node CovSettings = TMacroToYAML(*ConfigCov);
368  // Config from currently used cov object
369  YAML::Node ConfigCurrent = systematics[s]->GetConfig();
370 
371  if (!compareYAMLNodes(CovSettings, ConfigCurrent))
372  {
373  MACH3LOG_ERROR("Yaml configs in previous chain (from path {}) and current one are different", FitName);
374  throw MaCh3Exception(__FILE__ , __LINE__ );
375  }
376 
377  delete ConfigCov;
378  }
379 
380  CovarianceFolder->Close();
381  delete CovarianceFolder;
382 
383  std::vector<double> branch_vals;
384  std::vector<std::string> branch_name;
385  systematics[s]->MatchMaCh3OutputBranches(posts, branch_vals, branch_name);
386  posts->GetEntry(posts->GetEntries()-1);
387 
388  systematics[s]->SetParameters(branch_vals);
389  systematics[s]->AcceptStep();
390 
391  MACH3LOG_INFO("Printing new starting values for: {}", systematics[s]->GetName());
392  systematics[s]->PrintNominalCurrProp();
393 
394  // Resetting branch addressed to nullptr as we don't want to write into a delected vector out of scope...
395  for (int i = 0; i < systematics[s]->GetNumParams(); ++i) {
396  posts->SetBranchAddress(systematics[s]->GetParName(i).c_str(), nullptr);
397  }
398  }
399  logLCurr = log_val;
400  logLProp = log_val;
401  delete posts;
402  infile->Close();
403  delete infile;
404 
405  for (size_t s = 0; s < systematics.size(); ++s) {
406  if(systematics[s]->GetDoAdaption()){ //Use separate throw matrix for xsec
407  systematics[s]->SetNumberOfSteps(step_val);
408  }
409  }
410 }
YAML::Node TMacroToYAML(const TMacro &macro)
KS: Convert a ROOT TMacro object to a YAML node.
Definition: YamlHelper.h:152
bool compareYAMLNodes(const YAML::Node &node1, const YAML::Node &node2, bool Mute=false)
Compare if yaml nodes are identical.
Definition: YamlHelper.h:186
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:80

Member Data Documentation

◆ accCount

int FitterBase::accCount
protected

counts accepted steps

Definition at line 119 of file FitterBase.h.

◆ accProb

double FitterBase::accProb
protected

current acceptance prob

Definition at line 117 of file FitterBase.h.

◆ AlgorithmName

std::string FitterBase::AlgorithmName
protected

Name of fitting algorithm that is being used.

Definition at line 168 of file FitterBase.h.

◆ auto_save

int FitterBase::auto_save
protected

auto save every N steps

Definition at line 155 of file FitterBase.h.

◆ clock

std::unique_ptr<TStopwatch> FitterBase::clock
protected

tells global time how long fit took

Definition at line 137 of file FitterBase.h.

◆ CovFolder

TDirectory* FitterBase::CovFolder
protected

Output cov folder.

Definition at line 149 of file FitterBase.h.

◆ FileSaved

bool FitterBase::FileSaved
protected

Checks if file saved not repeat some operations.

Definition at line 161 of file FitterBase.h.

◆ fitMan

Manager* FitterBase::fitMan
protected

The manager for configuration handling.

Definition at line 108 of file FitterBase.h.

◆ fTestLikelihood

bool FitterBase::fTestLikelihood
protected

Necessary for some fitting algorithms like PSO.

Definition at line 158 of file FitterBase.h.

◆ logLCurr

double FitterBase::logLCurr
protected

current likelihood

Definition at line 113 of file FitterBase.h.

◆ logLProp

double FitterBase::logLProp
protected

proposed likelihood

Definition at line 115 of file FitterBase.h.

◆ outputFile

TFile* FitterBase::outputFile
protected

Output.

Definition at line 147 of file FitterBase.h.

◆ OutputPrepared

bool FitterBase::OutputPrepared
protected

Checks if output prepared not repeat some operations.

Definition at line 165 of file FitterBase.h.

◆ outTree

TTree* FitterBase::outTree
protected

Output tree with posteriors.

Definition at line 153 of file FitterBase.h.

◆ random

std::unique_ptr<TRandom3> FitterBase::random
protected

Random number.

Definition at line 144 of file FitterBase.h.

◆ sample_llh

std::vector<double> FitterBase::sample_llh
protected

store the llh breakdowns

Definition at line 124 of file FitterBase.h.

◆ SampleFolder

TDirectory* FitterBase::SampleFolder
protected

Output sample folder.

Definition at line 151 of file FitterBase.h.

◆ samples

std::vector<SampleHandlerBase*> FitterBase::samples
protected

Sample holder.

Definition at line 129 of file FitterBase.h.

◆ SettingsSaved

bool FitterBase::SettingsSaved
protected

Checks if setting saved not repeat some operations.

Definition at line 163 of file FitterBase.h.

◆ step

unsigned int FitterBase::step
protected

current state

Definition at line 111 of file FitterBase.h.

◆ stepClock

std::unique_ptr<TStopwatch> FitterBase::stepClock
protected

tells how long single step/fit iteration took

Definition at line 139 of file FitterBase.h.

◆ stepStart

unsigned int FitterBase::stepStart
protected

step start, by default 0 if we start from previous chain then it will be different

Definition at line 121 of file FitterBase.h.

◆ stepTime

double FitterBase::stepTime
protected

Time of single step.

Definition at line 141 of file FitterBase.h.

◆ syst_llh

std::vector<double> FitterBase::syst_llh
protected

systematic llh breakdowns

Definition at line 126 of file FitterBase.h.

◆ systematics

std::vector<ParameterHandlerBase*> FitterBase::systematics
protected

Systematic holder.

Definition at line 134 of file FitterBase.h.

◆ TotalNSamples

unsigned int FitterBase::TotalNSamples
protected

Total number of samples used, single SampleHandler can store more than one analysis sample!

Definition at line 131 of file FitterBase.h.


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