MaCh3  2.5.0
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 (SampleHandlerInterface *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 (const std::string &filename="")
 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< SampleHandlerInterface * > 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 SampleHandlerInterface 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 15 of file FitterBase.cpp.

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

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

Member Function Documentation

◆ AddSampleHandler()

void FitterBase::AddSampleHandler ( SampleHandlerInterface 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 260 of file FitterBase.cpp.

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

◆ 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 296 of file FitterBase.cpp.

296  {
297 // *************************
298  MACH3LOG_INFO("Adding systematic object {}, with {} params", cov->GetName(), cov->GetNumParams());
299  // 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...
300  for (size_t s = 0; s < systematics.size(); ++s)
301  {
302  for (int iPar = 0; iPar < systematics[s]->GetNumParams(); ++iPar)
303  {
304  for (int i = 0; i < cov->GetNumParams(); ++i)
305  {
306  if(systematics[s]->GetParName(iPar) == cov->GetParName(i)){
307  MACH3LOG_ERROR("ParameterHandler {} has param '{}' which already exists in in {}, with name {}",
308  cov->GetName(), cov->GetParName(i), systematics[s]->GetName(), systematics[s]->GetParName(iPar));
309  throw MaCh3Exception(__FILE__ , __LINE__ );
310  }
311  // Same for fancy name
312  if(systematics[s]->GetParFancyName(iPar) == cov->GetParFancyName(i)){
313  MACH3LOG_ERROR("ParameterHandler {} has param '{}' which already exists in {}, with name {}",
314  cov->GetName(), cov->GetParFancyName(i), systematics[s]->GetName(), systematics[s]->GetParFancyName(iPar));
315  throw MaCh3Exception(__FILE__ , __LINE__ );
316  }
317  }
318  }
319  }
320 
321  systematics.push_back(cov);
322 
323  CovFolder->cd();
324  std::vector<double> n_vec(cov->GetNumParams());
325  for (int i = 0; i < cov->GetNumParams(); ++i) {
326  n_vec[i] = cov->GetParInit(i);
327  }
328  cov->GetCovMatrix()->Write(cov->GetName().c_str());
329 
330  TH2D* CorrMatrix = cov->GetCorrelationMatrix();
331  CorrMatrix->Write((cov->GetName() + std::string("_Corr")).c_str());
332  delete CorrMatrix;
333 
334  // If we have yaml config file for covariance let's save it
335  YAML::Node Config = cov->GetConfig();
336  if(!Config.IsNull())
337  {
338  TMacro ConfigSave = YAMLtoTMacro(Config, (std::string("Config_") + cov->GetName()));
339  ConfigSave.Write();
340  }
341 
342  outputFile->cd();
343 }
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:135
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 541 of file FitterBase.cpp.

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

◆ 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 1299 of file FitterBase.cpp.

1299  {
1300 // *************************
1301  if(!fitMan->raw()["SigmaVar"]["CustomRange"]) return;
1302 
1303  auto Config = fitMan->raw()["SigmaVar"]["CustomRange"];
1304 
1305  const auto sigmaStr = std::to_string(static_cast<int>(std::round(sigma)));
1306 
1307  if (Config[ParName] && Config[ParName][sigmaStr]) {
1308  ParamShiftValue = Config[ParName][sigmaStr].as<double>();
1309  MACH3LOG_INFO(" ::: setting custom range from config ::: {} -> {}", ParName, ParamShiftValue);
1310  }
1311 }

◆ 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 459 of file FitterBase.cpp.

459  {
460 // *************************
461  MACH3LOG_INFO("Let the Race Begin!");
462  MACH3LOG_INFO("All tests will be performed with {} threads", M3::GetNThreads());
463 
464  // Reweight the MC
465  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs)
466  {
467  TStopwatch clockRace;
468  clockRace.Start();
469  for(int Lap = 0; Lap < NLaps; ++Lap) {
470  samples[ivs]->Reweight();
471  }
472  clockRace.Stop();
473  MACH3LOG_INFO("It took {:.4f} s to reweights {} times sample: {}", clockRace.RealTime(), NLaps, samples[ivs]->GetName());
474  MACH3LOG_INFO("On average {:.6f}", clockRace.RealTime()/NLaps);
475  }
476 
477  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs)
478  {
479  TStopwatch clockRace;
480  clockRace.Start();
481  for(int Lap = 0; Lap < NLaps; ++Lap) {
482  samples[ivs]->GetLikelihood();
483  }
484  clockRace.Stop();
485  MACH3LOG_INFO("It took {:.4f} s to calculate GetLikelihood {} times sample: {}", clockRace.RealTime(), NLaps, samples[ivs]->GetName());
486  MACH3LOG_INFO("On average {:.6f}", clockRace.RealTime()/NLaps);
487  }
488  // 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
489  std::vector<std::vector<double>> StepsValuesBefore(systematics.size());
490  for (size_t s = 0; s < systematics.size(); ++s) {
491  StepsValuesBefore[s] = systematics[s]->GetProposed();
492  }
493  for (size_t s = 0; s < systematics.size(); ++s) {
494  TStopwatch clockRace;
495  clockRace.Start();
496  for(int Lap = 0; Lap < NLaps; ++Lap) {
497  systematics[s]->ProposeStep();
498  }
499  clockRace.Stop();
500  MACH3LOG_INFO("It took {:.4f} s to propose step {} times cov: {}", clockRace.RealTime(), NLaps, systematics[s]->GetName());
501  MACH3LOG_INFO("On average {:.6f}", clockRace.RealTime()/NLaps);
502  }
503  for (size_t s = 0; s < systematics.size(); ++s) {
504  systematics[s]->SetParameters(StepsValuesBefore[s]);
505  }
506 
507  for (size_t s = 0; s < systematics.size(); ++s) {
508  TStopwatch clockRace;
509  clockRace.Start();
510  for(int Lap = 0; Lap < NLaps; ++Lap) {
511  systematics[s]->GetLikelihood();
512  }
513  clockRace.Stop();
514  MACH3LOG_INFO("It took {:.4f} s to calculate get likelihood {} times cov: {}", clockRace.RealTime(), NLaps, systematics[s]->GetName());
515  MACH3LOG_INFO("On average {:.6f}", clockRace.RealTime()/NLaps);
516  }
517  MACH3LOG_INFO("End of race");
518 }
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 71 of file FitterBase.h.

71 {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 557 of file FitterBase.cpp.

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

◆ 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 521 of file FitterBase.cpp.

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

◆ GetStepScaleBasedOnLLHScan()

void FitterBase::GetStepScaleBasedOnLLHScan ( const std::string &  filename = "")

LLH scan is good first estimate of step scale.

Definition at line 885 of file FitterBase.cpp.

885  {
886 // *************************
887  TFile* outputFileLLH = nullptr;
888  bool ownsfile = false;
889  if(outputFileName != ""){
890  outputFileLLH = M3::Open(outputFileName, "READ", __FILE__, __LINE__);
891  ownsfile = true;
892  }
893  else
894  outputFileLLH = outputFile;
895 
896  TDirectory *Sample_LLH = outputFileLLH->Get<TDirectory>("Sample_LLH");
897  MACH3LOG_INFO("Starting Get Step Scale Based On LLHScan");
898 
899  if(!Sample_LLH || Sample_LLH->IsZombie())
900  {
901  MACH3LOG_WARN("Couldn't find Sample_LLH, it looks like LLH scan wasn't run, will do this now");
902  RunLLHScan();
903  Sample_LLH = outputFileLLH->Get<TDirectory>("Sample_LLH");
904  }
905 
906  for (ParameterHandlerBase *cov : systematics)
907  {
908  const int npars = cov->GetNumParams();
909  std::vector<double> StepScale(npars);
910  for (int i = 0; i < npars; ++i)
911  {
912  std::string name = cov->GetParFancyName(i);
913  StepScale[i] = cov->GetIndivStepScale(i);
914  TH1D* LLHScan = Sample_LLH->Get<TH1D>((name+"_sam").c_str());
915  if(LLHScan == nullptr)
916  {
917  MACH3LOG_WARN("Couldn't find LLH scan, for {}, skipping", name);
918  continue;
919  }
920  const double LLH_val = std::max(LLHScan->GetBinContent(1), LLHScan->GetBinContent(LLHScan->GetNbinsX()));
921  //If there is no sensitivity leave it
922  if(LLH_val < 0.001) continue;
923 
924  // EM: assuming that the likelihood is gaussian, approximate sigma value is given by variation/sqrt(-2LLH)
925  // can evaluate this at any point, simple to evaluate it in the first bin of the LLH scan
926  // KS: We assume variation is 1 sigma, each dial has different scale so it becomes faff...
927  const double Var = 1.;
928  const double approxSigma = TMath::Abs(Var)/std::sqrt(LLH_val);
929 
930  // 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)
931  const double NewStepScale = approxSigma * 2.38/std::sqrt(npars);
932  StepScale[i] = NewStepScale;
933  MACH3LOG_DEBUG("Sigma: {}", approxSigma);
934  MACH3LOG_DEBUG("optimal Step Size: {}", NewStepScale);
935  }
936  cov->SetIndivStepScale(StepScale);
937  cov->SaveUpdatedMatrixConfig();
938  }
939  if(ownsfile && outputFileLLH != nullptr) delete outputFileLLH;
940 }
void RunLLHScan()
Perform a 1D likelihood scan.
Definition: FitterBase.cpp:620
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 151 of file FitterBase.cpp.

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

◆ ProcessMCMC()

void FitterBase::ProcessMCMC ( )
protected

Process MCMC output.

Definition at line 412 of file FitterBase.cpp.

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

◆ 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 944 of file FitterBase.cpp.

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

◆ RunLLHMap()

void FitterBase::RunLLHMap ( )

Perform a general multi-dimensional likelihood scan.

Author
Tomas Nosek

Definition at line 1047 of file FitterBase.cpp.

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

◆ RunLLHScan()

void FitterBase::RunLLHScan ( )

Perform a 1D likelihood scan.

Definition at line 620 of file FitterBase.cpp.

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

◆ 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 MCMCBase, PyLikelihoodFit, PyFitterBase, PredictiveThrower, PSO, and MinuitFit.

◆ 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(SampleHandlerInterface *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 221 of file FitterBase.cpp.

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

◆ SaveOutput()

void FitterBase::SaveOutput ( )
protected

Save output and close files.

Definition at line 229 of file FitterBase.cpp.

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

◆ SaveSettings()

void FitterBase::SaveSettings ( )
protected

Save the settings that the MCMC was run with.

Definition at line 77 of file FitterBase.cpp.

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

◆ 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 346 of file FitterBase.cpp.

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

◆ accProb

double FitterBase::accProb
protected

current acceptance prob

Definition at line 118 of file FitterBase.h.

◆ AlgorithmName

std::string FitterBase::AlgorithmName
protected

Name of fitting algorithm that is being used.

Definition at line 169 of file FitterBase.h.

◆ auto_save

int FitterBase::auto_save
protected

auto save every N steps

Definition at line 156 of file FitterBase.h.

◆ clock

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

tells global time how long fit took

Definition at line 138 of file FitterBase.h.

◆ CovFolder

TDirectory* FitterBase::CovFolder
protected

Output cov folder.

Definition at line 150 of file FitterBase.h.

◆ FileSaved

bool FitterBase::FileSaved
protected

Checks if file saved not repeat some operations.

Definition at line 162 of file FitterBase.h.

◆ fitMan

Manager* FitterBase::fitMan
protected

The manager for configuration handling.

Definition at line 109 of file FitterBase.h.

◆ fTestLikelihood

bool FitterBase::fTestLikelihood
protected

Necessary for some fitting algorithms like PSO.

Definition at line 159 of file FitterBase.h.

◆ logLCurr

double FitterBase::logLCurr
protected

current likelihood

Definition at line 114 of file FitterBase.h.

◆ logLProp

double FitterBase::logLProp
protected

proposed likelihood

Definition at line 116 of file FitterBase.h.

◆ outputFile

TFile* FitterBase::outputFile
protected

Output.

Definition at line 148 of file FitterBase.h.

◆ OutputPrepared

bool FitterBase::OutputPrepared
protected

Checks if output prepared not repeat some operations.

Definition at line 166 of file FitterBase.h.

◆ outTree

TTree* FitterBase::outTree
protected

Output tree with posteriors.

Definition at line 154 of file FitterBase.h.

◆ random

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

Random number.

Definition at line 145 of file FitterBase.h.

◆ sample_llh

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

store the llh breakdowns

Definition at line 125 of file FitterBase.h.

◆ SampleFolder

TDirectory* FitterBase::SampleFolder
protected

Output sample folder.

Definition at line 152 of file FitterBase.h.

◆ samples

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

Sample holder.

Definition at line 130 of file FitterBase.h.

◆ SettingsSaved

bool FitterBase::SettingsSaved
protected

Checks if setting saved not repeat some operations.

Definition at line 164 of file FitterBase.h.

◆ step

unsigned int FitterBase::step
protected

current state

Definition at line 112 of file FitterBase.h.

◆ stepClock

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

tells how long single step/fit iteration took

Definition at line 140 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 122 of file FitterBase.h.

◆ stepTime

double FitterBase::stepTime
protected

Time of single step.

Definition at line 142 of file FitterBase.h.

◆ syst_llh

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

systematic llh breakdowns

Definition at line 127 of file FitterBase.h.

◆ systematics

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

Systematic holder.

Definition at line 135 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 132 of file FitterBase.h.


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