MaCh3  2.5.1
Reference Guide
FitterBase.cpp
Go to the documentation of this file.
1 #include "FitterBase.h"
2 
4 #include "TRandom.h"
5 #include "TStopwatch.h"
6 #include "TTree.h"
7 #include "TGraphAsymmErrors.h"
9 
10 #pragma GCC diagnostic ignored "-Wuseless-cast"
11 
12 // *************************
13 // Initialise the Manager and make it an object of FitterBase class
14 // Now we can dump Manager settings to the output file
15 FitterBase::FitterBase(Manager * const man) : 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 }
64 
65 // *************************
66 // Destructor: close the logger and output file
68 // *************************
69  SaveOutput();
70  if(outputFile != nullptr) delete outputFile;
71  outputFile = nullptr;
72  MACH3LOG_DEBUG("Closing MaCh3 Fitter Engine");
73 }
74 
75 // *******************
76 // Prepare the output tree
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 }
148 
149 // *******************
150 // Prepare the output tree
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 }
219 
220 // *******************
222 // *******************
223  for (size_t i = 0; i < samples.size(); ++i) {
224  samples[i]->CleanMemoryBeforeFit();
225  }
226 }
227 
228 // *******************
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 }
257 
258 // *************************
259 // Add SampleHandler object to the Markov Chain
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_ERROR("SampleHandler with name '{}' already exists!", sample->GetName());
280  MACH3LOG_ERROR("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 }
293 
294 // *************************
295 // Add flux systematics, cross-section systematics, ND systematics to the chain
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 }
344 
345 // *******************
346 void FitterBase::StartFromPreviousFit(const std::string& FitName) {
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  delete ConfigCov;
375  }
376 
377  CovarianceFolder->Close();
378  delete CovarianceFolder;
379 
380  std::vector<double> branch_vals;
381  std::vector<std::string> branch_name;
382  systematics[s]->MatchMaCh3OutputBranches(posts, branch_vals, branch_name);
383  posts->GetEntry(posts->GetEntries()-1);
384 
385  systematics[s]->SetParameters(branch_vals);
386  systematics[s]->AcceptStep();
387 
388  MACH3LOG_INFO("Printing new starting values for: {}", systematics[s]->GetName());
389  systematics[s]->PrintPreFitCurrPropValues();
390 
391  // Resetting branch addressed to nullptr as we don't want to write into a delected vector out of scope...
392  for (int i = 0; i < systematics[s]->GetNumParams(); ++i) {
393  posts->SetBranchAddress(systematics[s]->GetParName(i).c_str(), nullptr);
394  }
395  }
396  logLCurr = log_val;
397  logLProp = log_val;
398  delete posts;
399  infile->Close();
400  delete infile;
401 
402  for (size_t s = 0; s < systematics.size(); ++s) {
403  if(systematics[s]->GetDoAdaption()){ //Use separate throw matrix for xsec
404  systematics[s]->SetNumberOfSteps(step_val);
405  }
406  }
407 }
408 
409 // *******************
410 // Process the MCMC output to get postfit etc
412 // *******************
413  if (fitMan == nullptr) return;
414 
415  // Process the MCMC
416  if (GetFromManager<bool>(fitMan->raw()["General"]["ProcessMCMC"], false, __FILE__ , __LINE__)){
417  // Make the processor
418  MCMCProcessor Processor(std::string(outputFile->GetName()));
419 
420  Processor.Initialise();
421  // Make the TVectorD pointers which hold the processed output
422  TVectorD *Central = nullptr;
423  TVectorD *Errors = nullptr;
424  TVectorD *Central_Gauss = nullptr;
425  TVectorD *Errors_Gauss = nullptr;
426  TVectorD *Peaks = nullptr;
427 
428  // Make the postfit
429  Processor.GetPostfit(Central, Errors, Central_Gauss, Errors_Gauss, Peaks);
430  Processor.DrawPostfit();
431 
432  // Make the TMatrix pointers which hold the processed output
433  TMatrixDSym *Covariance = nullptr;
434  TMatrixDSym *Correlation = nullptr;
435 
436  // Make the covariance matrix
437  Processor.GetCovariance(Covariance, Correlation);
438  Processor.DrawCovariance();
439 
440  std::vector<TString> BranchNames = Processor.GetBranchNames();
441 
442  // Re-open the TFile
443  if (!outputFile->IsOpen()) {
444  MACH3LOG_INFO("Opening output again to update with means..");
445  outputFile = new TFile(Get<std::string>(fitMan->raw()["General"]["OutputFile"], __FILE__, __LINE__).c_str(), "UPDATE");
446  }
447  Central->Write("PDF_Means");
448  Errors->Write("PDF_Errors");
449  Central_Gauss->Write("Gauss_Means");
450  Errors_Gauss->Write("Errors_Gauss");
451  Covariance->Write("Covariance");
452  Correlation->Write("Correlation");
453  }
454 }
455 
456 // *************************
457 // Run Drag Race
458 void FitterBase::DragRace(const int NLaps) {
459 // *************************
460  MACH3LOG_INFO("Let the Race Begin!");
461  MACH3LOG_INFO("All tests will be performed with {} threads", M3::GetNThreads());
462 
463  // Reweight the MC
464  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs)
465  {
466  TStopwatch clockRace;
467  clockRace.Start();
468  for(int Lap = 0; Lap < NLaps; ++Lap) {
469  samples[ivs]->Reweight();
470  }
471  clockRace.Stop();
472  MACH3LOG_INFO("It took {:.4f} s to reweights {} times sample: {}", clockRace.RealTime(), NLaps, samples[ivs]->GetName());
473  MACH3LOG_INFO("On average {:.6f}", clockRace.RealTime()/NLaps);
474  }
475 
476  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs)
477  {
478  TStopwatch clockRace;
479  clockRace.Start();
480  for(int Lap = 0; Lap < NLaps; ++Lap) {
481  samples[ivs]->GetLikelihood();
482  }
483  clockRace.Stop();
484  MACH3LOG_INFO("It took {:.4f} s to calculate GetLikelihood {} times sample: {}", clockRace.RealTime(), NLaps, samples[ivs]->GetName());
485  MACH3LOG_INFO("On average {:.6f}", clockRace.RealTime()/NLaps);
486  }
487  // 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
488  std::vector<std::vector<double>> StepsValuesBefore(systematics.size());
489  for (size_t s = 0; s < systematics.size(); ++s) {
490  StepsValuesBefore[s] = systematics[s]->GetProposed();
491  }
492  for (size_t s = 0; s < systematics.size(); ++s) {
493  TStopwatch clockRace;
494  clockRace.Start();
495  for(int Lap = 0; Lap < NLaps; ++Lap) {
496  systematics[s]->ProposeStep();
497  }
498  clockRace.Stop();
499  MACH3LOG_INFO("It took {:.4f} s to propose step {} times cov: {}", clockRace.RealTime(), NLaps, systematics[s]->GetName());
500  MACH3LOG_INFO("On average {:.6f}", clockRace.RealTime()/NLaps);
501  }
502  for (size_t s = 0; s < systematics.size(); ++s) {
503  systematics[s]->SetParameters(StepsValuesBefore[s]);
504  }
505 
506  for (size_t s = 0; s < systematics.size(); ++s) {
507  TStopwatch clockRace;
508  clockRace.Start();
509  for(int Lap = 0; Lap < NLaps; ++Lap) {
510  systematics[s]->GetLikelihood();
511  }
512  clockRace.Stop();
513  MACH3LOG_INFO("It took {:.4f} s to calculate get likelihood {} times cov: {}", clockRace.RealTime(), NLaps, systematics[s]->GetName());
514  MACH3LOG_INFO("On average {:.6f}", clockRace.RealTime()/NLaps);
515  }
516  MACH3LOG_INFO("End of race");
517 }
518 
519 // *************************
520 bool FitterBase::GetScanRange(std::map<std::string, std::vector<double>>& scanRanges) const {
521 // *************************
522  bool isScanRanges = false;
523  // YSP: Set up a mapping to store parameters with user-specified ranges, suggested by D. Barrow
524  if(fitMan->raw()["LLHScan"]["ScanRanges"]){
525  YAML::Node scanRangesList = fitMan->raw()["LLHScan"]["ScanRanges"];
526  for (auto it = scanRangesList.begin(); it != scanRangesList.end(); ++it) {
527  std::string itname = it->first.as<std::string>();
528  std::vector<double> itrange = it->second.as<std::vector<double>>();
529  // Set the mapping as param_name:param_range
530  scanRanges[itname] = itrange;
531  }
532  isScanRanges = true;
533  } else {
534  MACH3LOG_INFO("There are no user-defined parameter ranges, so I'll use default param bounds for LLH Scans");
535  }
536  return isScanRanges;
537 }
538 
539 // *************************
540 bool FitterBase::CheckSkipParameter(const std::vector<std::string>& SkipVector, const std::string& ParamName) const {
541 // *************************
542  bool skip = false;
543  for(unsigned int is = 0; is < SkipVector.size(); ++is)
544  {
545  if(ParamName.substr(0, SkipVector[is].length()) == SkipVector[is])
546  {
547  skip = true;
548  break;
549  }
550  }
551  return skip;
552 }
553 
554 
555 // *************************
556 void FitterBase::GetParameterScanRange(const ParameterHandlerBase* cov, const int i, double& CentralValue,
557  double& lower, double& upper, const int n_points, const std::string& suffix) const {
558 // *************************
559  // YSP: Set up a mapping to store parameters with user-specified ranges, suggested by D. Barrow
560  std::map<std::string, std::vector<double>> scanRanges;
561  const bool isScanRanges = GetScanRange(scanRanges);
562 
563  double nSigma = GetFromManager<int>(fitMan->raw()["LLHScan"]["LLHScanSigma"], 1., __FILE__, __LINE__);
564  bool IsPCA = cov->IsPCA();
565 
566  // Get the parameter name
567  std::string name = cov->GetParFancyName(i);
568  if (IsPCA) name += "_PCA";
569 
570  // Get the parameter priors and bounds
571  CentralValue = cov->GetParProp(i);
572  if (IsPCA) CentralValue = cov->GetPCAHandler()->GetParPropPCA(i);
573 
574  double prior = cov->GetParInit(i);
575  if (IsPCA) prior = cov->GetPCAHandler()->GetPreFitValuePCA(i);
576 
577  if (std::abs(CentralValue - prior) > 1e-10) {
578  MACH3LOG_INFO("For {} scanning around value {} rather than prior {}", name, CentralValue, prior);
579  }
580 
581  // Get the covariance matrix and do the +/- nSigma
582  // Set lower and upper bounds relative the CentralValue
583  // Set the parameter ranges between which LLH points are scanned
584  lower = CentralValue - nSigma*cov->GetDiagonalError(i);
585  upper = CentralValue + nSigma*cov->GetDiagonalError(i);
586  // If PCA, transform these parameter values to the PCA basis
587  if (IsPCA) {
588  lower = CentralValue - nSigma*std::sqrt((cov->GetPCAHandler()->GetEigenValues())(i));
589  upper = CentralValue + nSigma*std::sqrt((cov->GetPCAHandler()->GetEigenValues())(i));
590  MACH3LOG_INFO("eval {} = {:.2f}", i, cov->GetPCAHandler()->GetEigenValues()(i));
591  MACH3LOG_INFO("CV {} = {:.2f}", i, CentralValue);
592  MACH3LOG_INFO("lower {} = {:.2f}", i, lower);
593  MACH3LOG_INFO("upper {} = {:.2f}", i, upper);
594  MACH3LOG_INFO("nSigma = {:.2f}", nSigma);
595  }
596 
597  // Implementation suggested by D. Barrow
598  // If param ranges are specified in scanRanges node, extract it from there
599  if(isScanRanges){
600  // Find matching entries through std::maps
601  auto it = scanRanges.find(name);
602  if (it != scanRanges.end() && it->second.size() == 2) { //Making sure the range is has only two entries
603  lower = it->second[0];
604  upper = it->second[1];
605  MACH3LOG_INFO("Found matching param name for setting specified range for {}", name);
606  MACH3LOG_INFO("Range for {} = [{:.2f}, {:.2f}]", name, lower, upper);
607  }
608  }
609 
610  // Cross-section and flux parameters have boundaries that we scan between, check that these are respected in setting lower and upper variables
611  // This also applies for other parameters like osc, etc.
612  lower = std::max(lower, cov->GetLowerBound(i));
613  upper = std::min(upper, cov->GetUpperBound(i));
614  MACH3LOG_INFO("Scanning {} {} with {} steps, from [{:.2f} , {:.2f}], CV = {:.2f}", suffix, name, n_points, lower, upper, CentralValue);
615 }
616 
617 // *************************
618 // Run LLH scan
620 // *************************
621  // Save the settings into the output file
622  SaveSettings();
623 
624  MACH3LOG_INFO("Starting {}", __func__);
625 
626  //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.
627  bool PlotLLHScanBySample = GetFromManager<bool>(fitMan->raw()["LLHScan"]["LLHScanBySample"], false, __FILE__ , __LINE__);
628  auto SkipVector = GetFromManager<std::vector<std::string>>(fitMan->raw()["LLHScan"]["LLHScanSkipVector"], {}, __FILE__ , __LINE__);
629 
630  // Now finally get onto the LLH scan stuff
631  // Very similar code to MCMC but never start MCMC; just scan over the parameter space
632  std::vector<TDirectory *> Cov_LLH(systematics.size());
633  for(unsigned int ivc = 0; ivc < systematics.size(); ++ivc )
634  {
635  std::string NameTemp = systematics[ivc]->GetName();
636  NameTemp = NameTemp.substr(0, NameTemp.find("_cov")) + "_LLH";
637  Cov_LLH[ivc] = outputFile->mkdir(NameTemp.c_str());
638  }
639 
640  std::vector<TDirectory *> SampleClass_LLH(samples.size());
641  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs )
642  {
643  std::string NameTemp = samples[ivs]->GetName();
644  SampleClass_LLH[ivs] = outputFile->mkdir(NameTemp.c_str());
645  }
646 
647  TDirectory *Sample_LLH = outputFile->mkdir("Sample_LLH");
648  TDirectory *Total_LLH = outputFile->mkdir("Total_LLH");
649 
650  std::vector<TDirectory *>SampleSplit_LLH;
651  if(PlotLLHScanBySample)
652  {
653  SampleSplit_LLH.resize(TotalNSamples);
654  int SampleIterator = 0;
655  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs )
656  {
657  for(int is = 0; is < samples[ivs]->GetNSamples(); ++is )
658  {
659  SampleSplit_LLH[SampleIterator] = outputFile->mkdir((samples[ivs]->GetSampleTitle(is)+ "_LLH").c_str());
660  SampleIterator++;
661  }
662  }
663  }
664  // Number of points we do for each LLH scan
665  const int n_points = GetFromManager<int>(fitMan->raw()["LLHScan"]["LLHScanPoints"], 100, __FILE__ , __LINE__);
666 
667  // We print 5 reweights
668  const int countwidth = int(double(n_points)/double(5));
669 
670  // Loop over the covariance classes
671  for (ParameterHandlerBase *cov : systematics)
672  {
673  // Scan over all the parameters
674  // Get the number of parameters
675  int npars = cov->GetNumParams();
676  bool IsPCA = cov->IsPCA();
677  if (IsPCA) npars = cov->GetNParameters();
678  for (int i = 0; i < npars; ++i)
679  {
680  // Get the parameter name
681  std::string name = cov->GetParFancyName(i);
682  if (IsPCA) name += "_PCA";
683  // KS: Check if we want to skip this parameter
684  if(CheckSkipParameter(SkipVector, name)) continue;
685  // Get the parameter central and bounds
686  double CentralValue, lower, upper;
687  GetParameterScanRange(cov, i, CentralValue, lower, upper, n_points);
688  // Make the TH1D
689  auto hScan = std::make_unique<TH1D>((name + "_full").c_str(), (name + "_full").c_str(), n_points, lower, upper);
690  hScan->SetTitle((std::string("2LLH_full, ") + name + ";" + name + "; -2(ln L_{sample} + ln L_{xsec+flux} + ln L_{det})").c_str());
691  hScan->SetDirectory(nullptr);
692 
693  auto hScanSam = std::make_unique<TH1D>((name + "_sam").c_str(), (name + "_sam").c_str(), n_points, lower, upper);
694  hScanSam->SetTitle((std::string("2LLH_sam, ") + name + ";" + name + "; -2(ln L_{sample})").c_str());
695  hScanSam->SetDirectory(nullptr);
696 
697  std::vector<std::unique_ptr<TH1D>> hScanSample(samples.size());
698  std::vector<double> nSamLLH(samples.size(), 0.0);
699  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs )
700  {
701  std::string NameTemp = samples[ivs]->GetName();
702  hScanSample[ivs] = std::make_unique<TH1D>((name+"_"+NameTemp).c_str(), (name+"_" + NameTemp).c_str(), n_points, lower, upper);
703  hScanSample[ivs]->SetDirectory(nullptr);
704  hScanSample[ivs]->SetTitle(("2LLH_" + NameTemp + ", " + name + ";" + name + "; -2(ln L_{" + NameTemp +"})").c_str());
705  }
706 
707  std::vector<std::unique_ptr<TH1D>> hScanCov(systematics.size());
708  std::vector<double> nCovLLH(systematics.size(), 0.0);
709  for(unsigned int ivc = 0; ivc < systematics.size(); ++ivc )
710  {
711  std::string NameTemp = systematics[ivc]->GetName();
712  NameTemp = NameTemp.substr(0, NameTemp.find("_cov"));
713  hScanCov[ivc] = std::make_unique<TH1D>((name + "_" + NameTemp).c_str(), (name + "_" + NameTemp).c_str(), n_points, lower, upper);
714  hScanCov[ivc]->SetDirectory(nullptr);
715  hScanCov[ivc]->SetTitle(("2LLH_" + NameTemp + ", " + name + ";" + name + "; -2(ln L_{" + NameTemp +"})").c_str());
716  }
717 
718  std::vector<std::unique_ptr<TH1D>> hScanSamSplit;
719  std::vector<double> sampleSplitllh;
720  if(PlotLLHScanBySample)
721  {
722  int SampleIterator = 0;
723  hScanSamSplit.resize(TotalNSamples);
724  sampleSplitllh.resize(TotalNSamples);
725  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs )
726  {
727  for(int is = 0; is < samples[ivs]->GetNSamples(); ++is )
728  {
729  auto histName = name + samples[ivs]->GetSampleTitle(is);
730  auto histTitle = std::string("2LLH_sam, ") + name + ";" + name + "; -2(ln L_{sample})";
731  hScanSamSplit[SampleIterator] = std::make_unique<TH1D>(histName.c_str(), histTitle.c_str(), n_points, lower, upper);
732  hScanSamSplit[SampleIterator]->SetDirectory(nullptr);
733  SampleIterator++;
734  }
735  }
736  }
737 
738  // Scan over the parameter space
739  for (int j = 0; j < n_points; ++j)
740  {
741  if (j % countwidth == 0)
742  M3::Utils::PrintProgressBar(j, n_points);
743 
744  // For PCA we have to do it differently
745  if (IsPCA) {
746  cov->GetPCAHandler()->SetParPropPCA(i, hScan->GetBinCenter(j+1));
747  } else {
748  // Set the parameter
749  cov->SetParProp(i, hScan->GetBinCenter(j+1));
750  }
751 
752  // Reweight the MC
753  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs ){
754  samples[ivs]->Reweight();
755  }
756  //Total LLH
757  double totalllh = 0.;
758 
759  // Get the -log L likelihoods
760  double samplellh = 0.;
761 
762  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs ) {
763  nSamLLH[ivs] = samples[ivs]->GetLikelihood();
764  samplellh += nSamLLH[ivs];
765  }
766 
767  for(unsigned int ivc = 0; ivc < systematics.size(); ++ivc ) {
768  nCovLLH[ivc] = systematics[ivc]->GetLikelihood();
769  totalllh += nCovLLH[ivc];
770  }
771 
772  totalllh += samplellh;
773 
774  if(PlotLLHScanBySample)
775  {
776  int SampleIterator = 0;
777  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs )
778  {
779  for(int is = 0; is < samples[ivs]->GetNSamples(); ++is)
780  {
781  sampleSplitllh[SampleIterator] = samples[ivs]->GetSampleLikelihood(is);
782  SampleIterator++;
783  }
784  }
785  }
786 
787  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs ) {
788  hScanSample[ivs]->SetBinContent(j+1, 2*nSamLLH[ivs]);
789  }
790  for(unsigned int ivc = 0; ivc < systematics.size(); ++ivc ) {
791  hScanCov[ivc]->SetBinContent(j+1, 2*nCovLLH[ivc]);
792  }
793 
794  hScanSam->SetBinContent(j+1, 2*samplellh);
795  hScan->SetBinContent(j+1, 2*totalllh);
796 
797  if(PlotLLHScanBySample)
798  {
799  int SampleIterator = 0;
800  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs )
801  {
802  for(int is = 0; is < samples[ivs]->GetNSamples(); ++is)
803  {
804  hScanSamSplit[SampleIterator]->SetBinContent(j+1, 2*sampleSplitllh[SampleIterator]);
805  SampleIterator++;
806  }
807  }
808  }
809  }
810  for(unsigned int ivc = 0; ivc < systematics.size(); ++ivc )
811  {
812  Cov_LLH[ivc]->cd();
813  hScanCov[ivc]->Write();
814  }
815 
816  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs )
817  {
818  SampleClass_LLH[ivs]->cd();
819  hScanSample[ivs]->Write();
820  }
821  Sample_LLH->cd();
822  hScanSam->Write();
823  Total_LLH->cd();
824  hScan->Write();
825 
826  if(PlotLLHScanBySample)
827  {
828  int SampleIterator = 0;
829  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs )
830  {
831  for(int is = 0; is < samples[ivs]->GetNSamples(); ++is)
832  {
833  SampleSplit_LLH[SampleIterator]->cd();
834  hScanSamSplit[SampleIterator]->Write();
835  SampleIterator++;
836  }
837  }
838  }
839 
840  // Reset the parameters to their CentralValue central values
841  if (IsPCA) {
842  cov->GetPCAHandler()->SetParPropPCA(i, CentralValue);
843  } else {
844  cov->SetParProp(i, CentralValue);
845  }
846  }//end loop over systematics
847  }//end loop covariance classes
848 
849  for(unsigned int ivc = 0; ivc < systematics.size(); ++ivc )
850  {
851  Cov_LLH[ivc]->Write();
852  delete Cov_LLH[ivc];
853  }
854 
855  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs )
856  {
857  SampleClass_LLH[ivs]->Write();
858  delete SampleClass_LLH[ivs];
859  }
860 
861  Sample_LLH->Write();
862  delete Sample_LLH;
863 
864  Total_LLH->Write();
865  delete Total_LLH;
866 
867  if(PlotLLHScanBySample)
868  {
869  int SampleIterator = 0;
870  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs )
871  {
872  for(int is = 0; is < samples[ivs]->GetNSamples(); ++is )
873  {
874  SampleSplit_LLH[SampleIterator]->Write();
875  delete SampleSplit_LLH[SampleIterator];
876  SampleIterator++;
877  }
878  }
879  }
880 }
881 
882 // *************************
883 //LLH scan is good first estimate of step scale
884 void FitterBase::GetStepScaleBasedOnLLHScan(const std::string& outputFileName) {
885 // *************************
886  TFile* outputFileLLH = nullptr;
887  bool ownsfile = false;
888  if(outputFileName != ""){
889  outputFileLLH = M3::Open(outputFileName, "READ", __FILE__, __LINE__);
890  ownsfile = true;
891  } else {
892  outputFileLLH = outputFile;
893  }
894  TDirectory *Sample_LLH = outputFileLLH->Get<TDirectory>("Sample_LLH");
895  MACH3LOG_INFO("Starting Get Step Scale Based On LLHScan");
896 
897  if(!Sample_LLH || Sample_LLH->IsZombie())
898  {
899  MACH3LOG_WARN("Couldn't find Sample_LLH, it looks like LLH scan wasn't run, will do this now");
900  RunLLHScan();
901  Sample_LLH = outputFileLLH->Get<TDirectory>("Sample_LLH");
902  }
903 
904  for (ParameterHandlerBase *cov : systematics)
905  {
906  const int npars = cov->GetNumParams();
907  std::vector<double> StepScale(npars);
908  for (int i = 0; i < npars; ++i)
909  {
910  std::string name = cov->GetParFancyName(i);
911  StepScale[i] = cov->GetIndivStepScale(i);
912  TH1D* LLHScan = Sample_LLH->Get<TH1D>((name+"_sam").c_str());
913  if(LLHScan == nullptr)
914  {
915  MACH3LOG_WARN("Couldn't find LLH scan, for {}, skipping", name);
916  continue;
917  }
918  const double LLH_val = std::max(LLHScan->GetBinContent(1), LLHScan->GetBinContent(LLHScan->GetNbinsX()));
919  //If there is no sensitivity leave it
920  if(LLH_val < 0.001) continue;
921 
922  // EM: assuming that the likelihood is gaussian, approximate sigma value is given by variation/sqrt(-2LLH)
923  // can evaluate this at any point, simple to evaluate it in the first bin of the LLH scan
924  // KS: We assume variation is 1 sigma, each dial has different scale so it becomes faff...
925  const double Var = 1.;
926  const double approxSigma = std::abs(Var)/std::sqrt(LLH_val);
927  const double GlobalScale = cov->GetGlobalStepScale();
928  // 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)
929  const double TargetStep = approxSigma * 2.38 / std::sqrt(npars);
930  // KS: Need to divide by currently used gloalStepScale
931  const double NewStepScale = TargetStep / GlobalScale;
932 
933  StepScale[i] = NewStepScale;
934  MACH3LOG_DEBUG("Sigma: {}", approxSigma);
935  MACH3LOG_DEBUG("Target Step Size (before accounting for global step size): {}", TargetStep);
936  MACH3LOG_DEBUG("Optimal Step Size: {}", NewStepScale);
937  }
938  cov->SetIndivStepScale(StepScale);
939  cov->SaveUpdatedMatrixConfig();
940  }
941  if(ownsfile && outputFileLLH != nullptr) delete outputFileLLH;
942 }
943 
944 // *************************
945 // Run 2D LLH scan
947 // *************************
948  // Save the settings into the output file
949  SaveSettings();
950 
951  MACH3LOG_INFO("Starting 2D LLH Scan");
952 
953  TDirectory *Sample_2DLLH = outputFile->mkdir("Sample_2DLLH");
954  auto SkipVector = GetFromManager<std::vector<std::string>>(fitMan->raw()["LLHScan"]["LLHScanSkipVector"], {}, __FILE__ , __LINE__);;
955 
956  // Number of points we do for each LLH scan
957  const int n_points = GetFromManager<int>(fitMan->raw()["LLHScan"]["2DLLHScanPoints"], 20, __FILE__ , __LINE__);
958  // We print 5 reweights
959  const int countwidth = int(double(n_points)/double(5));
960 
961  // Loop over the covariance classes
962  for (ParameterHandlerBase *cov : systematics)
963  {
964  // Scan over all the parameters
965  // Get the number of parameters
966  int npars = cov->GetNumParams();
967  bool IsPCA = cov->IsPCA();
968  if (IsPCA) npars = cov->GetNParameters();
969 
970  for (int i = 0; i < npars; ++i)
971  {
972  std::string name_x = cov->GetParFancyName(i);
973  if (IsPCA) name_x += "_PCA";
974  // Get the parameter central and bounds
975  double central_x, lower_x, upper_x;
976  GetParameterScanRange(cov, i, central_x, lower_x, upper_x, n_points, "X");
977 
978  // KS: Check if we want to skip this parameter
979  if(CheckSkipParameter(SkipVector, name_x)) continue;
980 
981  for (int j = 0; j < i; ++j)
982  {
983  std::string name_y = cov->GetParFancyName(j);
984  if (IsPCA) name_y += "_PCA";
985  // KS: Check if we want to skip this parameter
986  if(CheckSkipParameter(SkipVector, name_y)) continue;
987 
988  // Get the parameter central and bounds
989  double central_y, lower_y, upper_y;
990  GetParameterScanRange(cov, j, central_y, lower_y, upper_y, n_points, "Y");
991 
992  auto hScanSam = std::make_unique<TH2D>((name_x + "_" + name_y + "_sam").c_str(), (name_x + "_" + name_y + "_sam").c_str(),
993  n_points, lower_x, upper_x, n_points, lower_y, upper_y);
994  hScanSam->SetDirectory(nullptr);
995  hScanSam->GetXaxis()->SetTitle(name_x.c_str());
996  hScanSam->GetYaxis()->SetTitle(name_y.c_str());
997  hScanSam->GetZaxis()->SetTitle("2LLH_sam");
998 
999  // Scan over the parameter space
1000  for (int x = 0; x < n_points; ++x)
1001  {
1002  if (x % countwidth == 0)
1003  M3::Utils::PrintProgressBar(x, n_points);
1004 
1005  for (int y = 0; y < n_points; ++y)
1006  {
1007  // For PCA we have to do it differently
1008  if (IsPCA) {
1009  cov->GetPCAHandler()->SetParPropPCA(i, hScanSam->GetXaxis()->GetBinCenter(x+1));
1010  cov->GetPCAHandler()->SetParPropPCA(j, hScanSam->GetYaxis()->GetBinCenter(y+1));
1011  } else {
1012  // Set the parameter
1013  cov->SetParProp(i, hScanSam->GetXaxis()->GetBinCenter(x+1));
1014  cov->SetParProp(j, hScanSam->GetYaxis()->GetBinCenter(y+1));
1015  }
1016  // Reweight the MC
1017  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs) {
1018  samples[ivs]->Reweight();
1019  }
1020 
1021  // Get the -log L likelihoods
1022  double samplellh = 0;
1023  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs) {
1024  samplellh += samples[ivs]->GetLikelihood();
1025  }
1026  hScanSam->SetBinContent(x+1, y+1, 2*samplellh);
1027  }// end loop over y points
1028  } // end loop over x points
1029 
1030  Sample_2DLLH->cd();
1031  hScanSam->Write();
1032  // Reset the parameters to their central central values
1033  if (IsPCA) {
1034  cov->GetPCAHandler()->SetParPropPCA(i, central_x);
1035  cov->GetPCAHandler()->SetParPropPCA(j, central_y);
1036  } else {
1037  cov->SetParProp(i, central_x);
1038  cov->SetParProp(j, central_y);
1039  }
1040  } //end loop over systematics y
1041  }//end loop over systematics X
1042  }//end loop covariance classes
1043  Sample_2DLLH->Write();
1044  delete Sample_2DLLH;
1045 }
1046 
1047 // *************************
1048 // Run a general multi-dimensional LLH scan
1050 // *************************
1051  // Save the settings into the output file
1052  SaveSettings();
1053 
1054  MACH3LOG_INFO("Starting {}", __func__);
1055 
1056  //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.
1057  bool PlotLLHScanBySample = GetFromManager<bool>(fitMan->raw()["LLHScan"]["LLHScanBySample"], false, __FILE__ , __LINE__);
1058  auto ParamsOfInterest = GetFromManager<std::vector<std::string>>(fitMan->raw()["LLHScan"]["LLHParameters"], {}, __FILE__, __LINE__);
1059 
1060  if(ParamsOfInterest.empty()) {
1061  MACH3LOG_WARN("There were no LLH parameters of interest specified to run the LLHMap! LLHMap will not run at all ...");
1062  return;
1063  }
1064 
1065  // Parameters IDs within the covariance objects
1066  // ParamsCovIDs = {Name, CovObj, IDinCovObj}
1067  std::vector<std::tuple<std::string, ParameterHandlerBase*, int>> ParamsCovIDs;
1068  for(auto& p : ParamsOfInterest) {
1069  bool found = false;
1070  for(auto cov : systematics) {
1071  for(int c = 0; c < cov->GetNumParams(); ++c) {
1072  if(cov->GetParName(c) == p || cov->GetParFancyName(c) == p) {
1073  bool add = true;
1074  for(auto& pc : ParamsCovIDs) {
1075  if(std::get<1>(pc) == cov && std::get<2>(pc) == c)
1076  {
1077  MACH3LOG_WARN("Parameter {} as {}({}) listed multiple times for LLHMap, omitting and using only once!", p, cov->GetName(), c);
1078  add = false;
1079  break;
1080  }
1081  }
1082 
1083  if(add)
1084  ParamsCovIDs.push_back(std::make_tuple(p, cov, c));
1085 
1086  found = true;
1087  break;
1088  }
1089  }
1090  if(found)
1091  break;
1092  }
1093  if(found)
1094  MACH3LOG_INFO("Parameter {} found in {} at an index {}.", p, std::get<1>(ParamsCovIDs.back())->GetName(), std::get<2>(ParamsCovIDs.back()));
1095  else
1096  MACH3LOG_WARN("Parameter {} not found in any of the systematic covariance objects. Will not scan over this one!", p);
1097  }
1098 
1099  // ParamsRanges["parameter"] = {nPoints, {low, high}}
1100  std::map<std::string, std::pair<int, std::pair<double, double>>> ParamsRanges;
1101 
1102  MACH3LOG_INFO("======================================================================================");
1103  MACH3LOG_INFO("Performing a general multi-dimensional LogL map scan over following parameters ranges:");
1104  MACH3LOG_INFO("======================================================================================");
1105  unsigned long TotalPoints = 1;
1106 
1107  double nSigma = GetFromManager<int>(fitMan->raw()["LLHScan"]["LLHScanSigma"], 1., __FILE__, __LINE__);
1108 
1109  // TN: Setting up the scan ranges might look like a re-implementation of the
1110  // FitterBase::GetScanRange, but I guess the use-case here is a bit different.
1111  // Anyway, just in case, we can discuss and rewrite to everyone's liking!
1112  for(auto& p : ParamsCovIDs) {
1113  // Auxiliary vars to help readability
1114  std::string name = std::get<0>(p);
1115  int i = std::get<2>(p);
1116  ParameterHandlerBase* cov = std::get<1>(p);
1117 
1118  ParamsRanges[name].first = GetFromManager<int>(fitMan->raw()["LLHScan"]["LLHScanPoints"], 20, __FILE__, __LINE__);
1119  if(CheckNodeExists(fitMan->raw(),"LLHScan","ScanPoints"))
1120  ParamsRanges[name].first = GetFromManager<int>(fitMan->raw()["LLHScan"]["ScanPoints"][name], ParamsRanges[name].first, __FILE__, __LINE__);
1121 
1122  // Get the parameter priors and bounds
1123  double CentralValue = cov->GetParProp(i);
1124 
1125  bool IsPCA = cov->IsPCA();
1126  if (IsPCA)
1127  CentralValue = cov->GetPCAHandler()->GetParPropPCA(i);
1128 
1129  double prior = cov->GetParInit(i);
1130  if (IsPCA)
1131  prior = cov->GetPCAHandler()->GetPreFitValuePCA(i);
1132 
1133  if (std::abs(CentralValue - prior) > 1e-10) {
1134  MACH3LOG_INFO("For {} scanning around value {} rather than prior {}", name, CentralValue, prior);
1135  }
1136  // Get the covariance matrix and do the +/- nSigma
1137  // Set lower and upper bounds relative the CentralValue
1138  // Set the parameter ranges between which LLH points are scanned
1139  double lower = CentralValue - nSigma*cov->GetDiagonalError(i);
1140  double upper = CentralValue + nSigma*cov->GetDiagonalError(i);
1141  // If PCA, transform these parameter values to the PCA basis
1142  if (IsPCA) {
1143  lower = CentralValue - nSigma*std::sqrt((cov->GetPCAHandler()->GetEigenValues())(i));
1144  upper = CentralValue + nSigma*std::sqrt((cov->GetPCAHandler()->GetEigenValues())(i));
1145  MACH3LOG_INFO("eval {} = {:.2f}", i, cov->GetPCAHandler()->GetEigenValues()(i));
1146  MACH3LOG_INFO("CV {} = {:.2f}", i, CentralValue);
1147  MACH3LOG_INFO("lower {} = {:.2f}", i, lower);
1148  MACH3LOG_INFO("upper {} = {:.2f}", i, upper);
1149  MACH3LOG_INFO("nSigma = {:.2f}", nSigma);
1150  }
1151 
1152  ParamsRanges[name].second = {lower,upper};
1153 
1154  if(CheckNodeExists(fitMan->raw(),"LLHScan","ScanRanges"))
1155  ParamsRanges[name].second = GetFromManager<std::pair<double,double>>(fitMan->raw()["LLHScan"]["ScanRanges"][name], ParamsRanges[name].second, __FILE__, __LINE__);
1156 
1157  MACH3LOG_INFO("{} from {:.4f} to {:.4f} with a {:.5f} step ({} points total)",
1158  name, ParamsRanges[name].second.first, ParamsRanges[name].second.second,
1159  (ParamsRanges[name].second.second - ParamsRanges[name].second.first)/(ParamsRanges[name].first - 1.),
1160  ParamsRanges[name].first);
1161 
1162  TotalPoints *= ParamsRanges[name].first;
1163  }
1164 
1165  // TN: Waiting for C++ 20 std::format() function
1166  MACH3LOG_INFO("In total, looping over {} points, from {} parameters. Estimates for run time:", TotalPoints, ParamsCovIDs.size());
1167  MACH3LOG_INFO(" 1 s per point = {} hours", double(TotalPoints)/3600.);
1168  MACH3LOG_INFO(" 0.1 s per point = {} hours", double(TotalPoints)/36000.);
1169  MACH3LOG_INFO("0.01 s per point = {} hours", double(TotalPoints)/360000.);
1170  MACH3LOG_INFO("==================================================================================");
1171 
1172  const int countwidth = int(double(TotalPoints)/double(20));
1173 
1174  // Tree to store LogL values
1175  auto LLHMap = new TTree("llhmap", "LLH Map");
1176 
1177  std::vector<double> CovLogL(systematics.size());
1178  for(unsigned int ivc = 0; ivc < systematics.size(); ++ivc )
1179  {
1180  std::string NameTemp = systematics[ivc]->GetName();
1181  NameTemp = NameTemp.substr(0, NameTemp.find("_cov")) + "_LLH";
1182  LLHMap->Branch(NameTemp.c_str(), &CovLogL[ivc]);
1183  }
1184 
1185  std::vector<double> SampleClassLogL(samples.size());
1186  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs )
1187  {
1188  std::string NameTemp = samples[ivs]->GetName()+"_LLH";
1189  LLHMap->Branch(NameTemp.c_str(), &SampleClassLogL[ivs]);
1190  }
1191 
1192  double SampleLogL, TotalLogL;
1193  LLHMap->Branch("Sample_LLH", &SampleLogL);
1194  LLHMap->Branch("Total_LLH", &TotalLogL);
1195 
1196  std::vector<double>SampleSplitLogL;
1197  if(PlotLLHScanBySample)
1198  {
1199  SampleSplitLogL.resize(TotalNSamples);
1200  int SampleIterator = 0;
1201  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs )
1202  {
1203  for(int is = 0; is < samples[ivs]->GetNSamples(); ++is )
1204  {
1205  std::string NameTemp =samples[ivs]->GetSampleTitle(is)+"_LLH";
1206  LLHMap->Branch(NameTemp.c_str(), &SampleSplitLogL[SampleIterator]);
1207  SampleIterator++;
1208  }
1209  }
1210  }
1211 
1212  std::vector<double> ParamsValues(ParamsCovIDs.size());
1213  for(unsigned int i=0; i < ParamsCovIDs.size(); ++i)
1214  LLHMap->Branch(std::get<0>(ParamsCovIDs[i]).c_str(), &ParamsValues[i]);
1215 
1216  // Setting up the scan
1217  // Starting at index {0,0,0,...}
1218  std::vector<unsigned long> idx(ParamsCovIDs.size(), 0);
1219 
1220  // loop over scanned points sp
1221  for(unsigned long sp = 0; sp < TotalPoints; ++sp)
1222  {
1223  // At each point need to find the indices and test values to calculate LogL
1224  for(unsigned int n = 0; n < ParamsCovIDs.size(); ++n)
1225  {
1226  // Auxiliaries
1227  std::string name = std::get<0>(ParamsCovIDs[n]);
1228  int points = ParamsRanges[name].first;
1229  double low = ParamsRanges[name].second.first;
1230  double high = ParamsRanges[name].second.second;
1231 
1232  // Find the n-th index of the sp-th scanned point
1233  unsigned long dev = 1;
1234  for(unsigned int m = 0; m <= n; ++m)
1235  dev *= ParamsRanges[std::get<0>(ParamsCovIDs[m])].first;
1236 
1237  idx[n] = sp % dev;
1238  if (n > 0)
1239  idx[n] = idx[n] / ( dev / points );
1240 
1241  // Parameter test value = low + ( high - low ) * idx / ( #points - 1 )
1242  ParamsValues[n] = low + (high-low) * double(idx[n])/double(points-1);
1243 
1244  // Now set the covariance objects
1245  // Auxiliary
1246  ParameterHandlerBase* cov = std::get<1>(ParamsCovIDs[n]);
1247  int i = std::get<2>(ParamsCovIDs[n]);
1248 
1249  if(cov->IsPCA())
1250  cov->GetPCAHandler()->SetParPropPCA(i, ParamsValues[n]);
1251  else
1252  cov->SetParProp(i, ParamsValues[n]);
1253  }
1254 
1255  // Reweight samples and calculate LogL
1256  TotalLogL = .0;
1257  SampleLogL = .0;
1258 
1259  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs)
1260  samples[ivs]->Reweight();
1261 
1262  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs)
1263  {
1264  SampleClassLogL[ivs] = 2.*samples[ivs]->GetLikelihood();
1265  SampleLogL += SampleClassLogL[ivs];
1266  }
1267  TotalLogL += SampleLogL;
1268 
1269  // CovObjs LogL
1270  for(unsigned int ivc = 0; ivc < systematics.size(); ++ivc )
1271  {
1272  CovLogL[ivc] = 2.*systematics[ivc]->GetLikelihood();
1273  TotalLogL += CovLogL[ivc];
1274  }
1275 
1276  if(PlotLLHScanBySample)
1277  {
1278  int SampleIterator = 0;
1279  for(unsigned int ivs = 0; ivs < samples.size(); ++ivs )
1280  {
1281  for(int is = 0; is < samples[ivs]->GetNSamples(); ++is)
1282  {
1283  SampleSplitLogL[SampleIterator] = 2.*samples[ivs]->GetSampleLikelihood(is);
1284  SampleIterator++;
1285  }
1286  }
1287  }
1288 
1289  LLHMap->Fill();
1290 
1291  if (sp % countwidth == 0)
1292  M3::Utils::PrintProgressBar(sp, TotalPoints);
1293  }
1294 
1295  outputFile->cd();
1296  LLHMap->Write();
1297 }
1298 
1299 // *************************
1300 // For comparison with P-Theta we usually have to apply different parameter values then usual 1, 3 sigma
1301 void FitterBase::CustomRange(const std::string& ParName, const double sigma, double& ParamShiftValue) const {
1302 // *************************
1303  if(!fitMan->raw()["SigmaVar"]["CustomRange"]) return;
1304 
1305  auto Config = fitMan->raw()["SigmaVar"]["CustomRange"];
1306 
1307  const auto sigmaStr = std::to_string(static_cast<int>(std::round(sigma)));
1308 
1309  if (Config[ParName] && Config[ParName][sigmaStr]) {
1310  ParamShiftValue = Config[ParName][sigmaStr].as<double>();
1311  MACH3LOG_INFO(" ::: setting custom range from config ::: {} -> {}", ParName, ParamShiftValue);
1312  }
1313 }
1314 
1315 // *************************
1317 void WriteHistograms(TH1 *hist, const std::string& baseName) {
1318 // *************************
1319  if (!hist) return;
1320  hist->SetTitle(baseName.c_str());
1321  // Get the class name of the histogram
1322  TString className = hist->ClassName();
1323 
1324  // Set the appropriate axis title based on the histogram type
1325  if (className.Contains("TH1")) {
1326  hist->GetYaxis()->SetTitle("Events");
1327  } else if (className.Contains("TH2")) {
1328  hist->GetZaxis()->SetTitle("Events");
1329  }
1330  hist->Write(baseName.c_str());
1331 }
1332 
1333 // *************************
1336  const std::string& suffix,
1337  const bool by_mode,
1338  const bool by_channel,
1339  const std::vector<TDirectory*>& SampleDir) {
1340 // *************************
1341  MaCh3Modes *modes = sample->GetMaCh3Modes();
1342  for (int iSample = 0; iSample < sample->GetNSamples(); ++iSample) {
1343  SampleDir[iSample]->cd();
1344  const std::string sampleName = sample->GetSampleTitle(iSample);
1345  for(int iDim1 = 0; iDim1 < sample->GetNDim(iSample); iDim1++) {
1346  std::string ProjectionName = sample->GetKinVarName(iSample, iDim1);
1347  std::string ProjectionSuffix = "_1DProj" + std::to_string(iDim1);
1348 
1349  // Probably a better way of handling this logic
1350  if (by_mode) {
1351  for (int iMode = 0; iMode < modes->GetNModes(); ++iMode) {
1352  auto modeHist = sample->Get1DVarHistByModeAndChannel(iSample, ProjectionName, iMode);
1353  WriteHistograms(modeHist.get(), sampleName + "_" + modes->GetMaCh3ModeName(iMode) + ProjectionSuffix + suffix);
1354  }
1355  }
1356 
1357  if (by_channel) {
1358  for (int iChan = 0; iChan < sample->GetNOscChannels(iSample); ++iChan) {
1359  auto chanHist = sample->Get1DVarHistByModeAndChannel(iSample, ProjectionName, -1, iChan); // -1 skips over mode plotting
1360  WriteHistograms(chanHist.get(), sampleName + "_" + sample->GetFlavourName(iSample, iChan) + ProjectionSuffix + suffix);
1361  }
1362  }
1363 
1364  if (by_mode && by_channel) {
1365  for (int iMode = 0; iMode < modes->GetNModes(); ++iMode) {
1366  for (int iChan = 0; iChan < sample->GetNOscChannels(iSample); ++iChan) {
1367  auto hist = sample->Get1DVarHistByModeAndChannel(iSample, ProjectionName, iMode, iChan);
1368  WriteHistograms(hist.get(), sampleName + "_" + modes->GetMaCh3ModeName(iMode) + "_" + sample->GetFlavourName(iSample, iChan) + ProjectionSuffix + suffix);
1369  }
1370  }
1371  }
1372 
1373  if (!by_mode && !by_channel) {
1374  auto hist = sample->Get1DVarHist(iSample, ProjectionName);
1375  WriteHistograms(hist.get(), sampleName + ProjectionSuffix + suffix);
1376  // Only for 2D and Beyond
1377  for (int iDim2 = iDim1 + 1; iDim2 < sample->GetNDim(iSample); ++iDim2) {
1378  // Get the names for the two dimensions
1379  std::string XVarName = sample->GetKinVarName(iSample, iDim1);
1380  std::string YVarName = sample->GetKinVarName(iSample, iDim2);
1381 
1382  // Get the 2D histogram for this pair
1383  auto hist2D = sample->Get2DVarHist(iSample, XVarName, YVarName);
1384 
1385  // Write the histogram
1386  std::string suffix2D = "_2DProj_" + std::to_string(iDim1) + "_vs_" + std::to_string(iDim2) + suffix;
1387  WriteHistograms(hist2D.get(), sampleName + suffix2D);
1388  }
1389  }
1390  }
1391  }
1392 }
1393 
1394 // *************************
1396 // *************************
1397  // Save the settings into the output file
1398  SaveSettings();
1399 
1400  bool plot_by_mode = GetFromManager<bool>(fitMan->raw()["SigmaVar"]["PlotByMode"], false);
1401  bool plot_by_channel = GetFromManager<bool>(fitMan->raw()["SigmaVar"]["PlotByChannel"], false);
1402  auto SkipVector = GetFromManager<std::vector<std::string>>(fitMan->raw()["SigmaVar"]["SkipVector"], {}, __FILE__ , __LINE__);
1403 
1404  if (plot_by_mode) MACH3LOG_INFO("Plotting by sample and mode");
1405  if (plot_by_channel) MACH3LOG_INFO("Plotting by sample and channel");
1406  if (!plot_by_mode && !plot_by_channel) MACH3LOG_INFO("Plotting by sample only");
1407  if (plot_by_mode && plot_by_channel) MACH3LOG_INFO("Plotting by sample, mode and channel");
1408 
1409  auto SigmaArray = GetFromManager<std::vector<double>>(fitMan->raw()["SigmaVar"]["SigmaArray"], {-3, -1, 0, 1, 3}, __FILE__ , __LINE__);
1410  if (std::find(SigmaArray.begin(), SigmaArray.end(), 0.0) == SigmaArray.end()) {
1411  MACH3LOG_ERROR(":: SigmaArray does not contain 0! Current contents: {} ::", fmt::join(SigmaArray, ", "));
1412  throw MaCh3Exception(__FILE__, __LINE__);
1413  }
1414 
1415  TDirectory* SigmaDir = outputFile->mkdir("SigmaVar");
1416  outputFile->cd();
1417 
1418  for (size_t s = 0; s < systematics.size(); ++s)
1419  {
1420  for(int i = 0; i < systematics[s]->GetNumParams(); i++)
1421  {
1422  std::string ParName = systematics[s]->GetParFancyName(i);
1423  // KS: Check if we want to skip this parameter
1424  if(CheckSkipParameter(SkipVector, ParName)) continue;
1425 
1426  MACH3LOG_INFO(":: Param {} ::", systematics[s]->GetParFancyName(i));
1427 
1428  TDirectory* ParamDir = SigmaDir->mkdir(ParName.c_str());
1429  ParamDir->cd();
1430 
1431  const double ParamCentralValue = systematics[s]->GetParProp(i);
1432  const double Prior = systematics[s]->GetParInit(i);
1433  const double ParamLower = systematics[s]->GetLowerBound(i);
1434  const double ParamUpper = systematics[s]->GetUpperBound(i);
1435 
1436  if (std::abs(ParamCentralValue - Prior) > 1e-10) {
1437  MACH3LOG_INFO("For {} scanning around value {} rather than prior {}", ParName, ParamCentralValue, Prior);
1438  }
1439 
1440  for(unsigned int iSample = 0; iSample < samples.size(); ++iSample)
1441  {
1442  auto* MaCh3Sample = samples[iSample];
1443  std::vector<TDirectory*> SampleDir(MaCh3Sample->GetNSamples());
1444  for (int SampleIndex = 0; SampleIndex < MaCh3Sample->GetNSamples(); ++SampleIndex) {
1445  SampleDir[SampleIndex] = ParamDir->mkdir(MaCh3Sample->GetSampleTitle(SampleIndex).c_str());
1446  }
1447 
1448  for (size_t j = 0; j < SigmaArray.size(); ++j) {
1449  double sigma = SigmaArray[j];
1450 
1451  double ParamShiftValue = ParamCentralValue + sigma * std::sqrt((*systematics[s]->GetCovMatrix())(i,i));
1452  ParamShiftValue = std::max(std::min(ParamShiftValue, ParamUpper), ParamLower);
1453 
1455  CustomRange(ParName, sigma, ParamShiftValue);
1456 
1457  MACH3LOG_INFO(" - set to {:<5.2f} ({:<2} sigma shift)", ParamShiftValue, sigma);
1458  systematics[s]->SetParProp(i, ParamShiftValue);
1459 
1460  std::ostringstream valStream;
1461  valStream << std::fixed << std::setprecision(2) << ParamShiftValue;
1462  std::string valueStr = valStream.str();
1463 
1464  std::ostringstream sigmaStream;
1465  sigmaStream << std::fixed << std::setprecision(2) << std::abs(sigma);
1466  std::string sigmaStr = sigmaStream.str();
1467 
1468  std::string suffix;
1469  if (sigma == 0) {
1470  suffix = "_" + ParName + "_nom_val_" + valueStr;
1471  } else {
1472  std::string sign = (sigma > 0) ? "p" : "n";
1473  suffix = "_" + ParName + "_sig_" + sign + sigmaStr + "_val_" + valueStr;
1474  }
1475 
1476  systematics[s]->SetParProp(i, ParamShiftValue);
1477  MaCh3Sample->Reweight();
1478 
1479  WriteHistogramsByMode(MaCh3Sample, suffix, plot_by_mode, plot_by_channel, SampleDir);
1480  }
1481  for (int subSampleIndex = 0; subSampleIndex < MaCh3Sample->GetNSamples(); ++subSampleIndex) {
1482  SampleDir[subSampleIndex]->Close();
1483  delete SampleDir[subSampleIndex];
1484  }
1485  ParamDir->cd();
1486  }
1487 
1488  systematics[s]->SetParProp(i, ParamCentralValue);
1489  MACH3LOG_INFO(" - set back to CV {:<5.2f}", ParamCentralValue);
1490  MACH3LOG_INFO("");
1491  ParamDir->Close();
1492  delete ParamDir;
1493  SigmaDir->cd();
1494  } // end loop over params
1495  } // end loop over systemics
1496  SigmaDir->Close();
1497  delete SigmaDir;
1498 
1499  outputFile->cd();
1500 }
#define _MaCh3_Safe_Include_Start_
KS: Avoiding warning checking for headers.
Definition: Core.h:126
#define _MaCh3_Safe_Include_End_
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 WriteHistograms(TH1 *hist, const std::string &baseName)
Helper to write histograms.
#define MACH3LOG_CRITICAL
Definition: MaCh3Logger.h:38
#define MACH3LOG_DEBUG
Definition: MaCh3Logger.h:34
#define MACH3LOG_ERROR
Definition: MaCh3Logger.h:37
#define MACH3LOG_INFO
Definition: MaCh3Logger.h:35
#define MACH3LOG_WARN
Definition: MaCh3Logger.h:36
std::vector< TString > BranchNames
Definition: RHat.cpp:54
TMacro YAMLtoTMacro(const YAML::Node &yaml_node, const std::string &name)
Convert a YAML node to a ROOT TMacro object.
Definition: YamlHelper.h:167
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
bool CheckNodeExists(const YAML::Node &node, Args... args)
KS: Wrapper function to call the recursive helper.
Definition: YamlHelper.h:60
void RunLLHScan()
Perform a 1D likelihood scan.
Definition: FitterBase.cpp:619
FitterBase(Manager *const fitMan)
Constructor.
Definition: FitterBase.cpp:15
void AddSystObj(ParameterHandlerBase *cov)
This function adds a Covariance object to the analysis framework. The Covariance object will be utili...
Definition: FitterBase.cpp:296
std::string GetName() const
Get name of class.
Definition: FitterBase.h:72
std::unique_ptr< TRandom3 > random
Random number.
Definition: FitterBase.h:146
double logLProp
proposed likelihood
Definition: FitterBase.h:117
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:540
void ProcessMCMC()
Process MCMC output.
Definition: FitterBase.cpp:411
int accCount
counts accepted steps
Definition: FitterBase.h:121
bool OutputPrepared
Checks if output prepared not repeat some operations.
Definition: FitterBase.h:167
void SaveOutput()
Save output and close files.
Definition: FitterBase.cpp:229
TFile * outputFile
Output.
Definition: FitterBase.h:149
void SaveSettings()
Save the settings that the MCMC was run with.
Definition: FitterBase.cpp:77
unsigned int step
current state
Definition: FitterBase.h:113
void PrepareOutput()
Prepare the output file.
Definition: FitterBase.cpp:151
bool SettingsSaved
Checks if setting saved not repeat some operations.
Definition: FitterBase.h:165
double accProb
current acceptance prob
Definition: FitterBase.h:119
virtual void StartFromPreviousFit(const std::string &FitName)
Allow to start from previous fit/chain.
Definition: FitterBase.cpp:346
bool FileSaved
Checks if file saved not repeat some operations.
Definition: FitterBase.h:163
std::string AlgorithmName
Name of fitting algorithm that is being used.
Definition: FitterBase.h:170
std::vector< double > sample_llh
store the llh breakdowns
Definition: FitterBase.h:126
void RunSigmaVar()
Perform a 1D/2D sigma var for all samples.
std::vector< SampleHandlerInterface * > samples
Sample holder.
Definition: FitterBase.h:131
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:556
double stepTime
Time of single step.
Definition: FitterBase.h:143
std::unique_ptr< TStopwatch > clock
tells global time how long fit took
Definition: FitterBase.h:139
Manager * fitMan
The manager for configuration handling.
Definition: FitterBase.h:110
unsigned int stepStart
step start, by default 0 if we start from previous chain then it will be different
Definition: FitterBase.h:123
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:520
std::unique_ptr< TStopwatch > stepClock
tells how long single step/fit iteration took
Definition: FitterBase.h:141
TDirectory * CovFolder
Output cov folder.
Definition: FitterBase.h:151
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...
TDirectory * SampleFolder
Output sample folder.
Definition: FitterBase.h:153
void DragRace(const int NLaps=100)
Calculates the required time for each sample or covariance object in a drag race simulation....
Definition: FitterBase.cpp:458
void Run2DLLHScan()
Perform a 2D likelihood scan.
Definition: FitterBase.cpp:946
unsigned int TotalNSamples
Total number of samples used, single SampleHandler can store more than one analysis sample!
Definition: FitterBase.h:133
double logLCurr
current likelihood
Definition: FitterBase.h:115
void RunLLHMap()
Perform a general multi-dimensional likelihood scan.
std::vector< double > syst_llh
systematic llh breakdowns
Definition: FitterBase.h:128
int auto_save
auto save every N steps
Definition: FitterBase.h:157
bool fTestLikelihood
Necessary for some fitting algorithms like PSO.
Definition: FitterBase.h:160
void GetStepScaleBasedOnLLHScan(const std::string &filename="")
LLH scan is good first estimate of step scale.
Definition: FitterBase.cpp:884
virtual ~FitterBase()
Destructor for the FitterBase class.
Definition: FitterBase.cpp:67
TTree * outTree
Output tree with posteriors.
Definition: FitterBase.h:155
void SanitiseInputs()
Remove obsolete memory and make other checks before fit starts.
Definition: FitterBase.cpp:221
void AddSampleHandler(SampleHandlerInterface *sample)
This function adds a sample PDF object to the analysis framework. The sample PDF object will be utili...
Definition: FitterBase.cpp:260
std::vector< ParameterHandlerBase * > systematics
Systematic holder.
Definition: FitterBase.h:136
Class responsible for processing MCMC chains, performing diagnostics, generating plots,...
Definition: MCMCProcessor.h:61
const std::vector< TString > & GetBranchNames() const
Get the vector of branch names from root file.
void Initialise()
Scan chain, what parameters we have and load information from covariance matrices.
void GetPostfit(TVectorD *&Central, TVectorD *&Errors, TVectorD *&Central_Gauss, TVectorD *&Errors_Gauss, TVectorD *&Peaks)
Get the post-fit results (arithmetic and Gaussian)
void DrawCovariance()
Draw the post-fit covariances.
void DrawPostfit()
Draw the post-fit comparisons.
void GetCovariance(TMatrixDSym *&Cov, TMatrixDSym *&Corr)
Get the post-fit covariances and correlations.
Custom exception class used throughout MaCh3.
KS: Class describing MaCh3 modes used in the analysis, it is being initialised from config.
Definition: MaCh3Modes.h:142
int GetNModes() const
KS: Get number of modes, keep in mind actual number is +1 greater due to unknown category.
Definition: MaCh3Modes.h:155
std::string GetMaCh3ModeName(const int Index) const
KS: Get normal name of mode, if mode not known you will get UNKNOWN_BAD.
Definition: MaCh3Modes.cpp:156
The manager class is responsible for managing configurations and settings.
Definition: Manager.h:16
YAML::Node const & raw() const
Return config.
Definition: Manager.h:47
void SaveSettings(TFile *const OutputFile) const
Add manager useful information's to TFile, in most cases to Fitter.
Definition: Manager.cpp:40
double GetParPropPCA(const int i) const
Get current parameter value using PCA.
Definition: PCAHandler.h:154
void SetParPropPCA(const int i, const double value)
Set proposed value for parameter in PCA base.
Definition: PCAHandler.h:136
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
Base class responsible for handling of systematic error parameters. Capable of using PCA or using ada...
int GetNumParams() const
Get total number of parameters.
TH2D * GetCorrelationMatrix() const
KS: Convert covariance matrix to correlation matrix and return TH2D which can be used for fancy plott...
std::string GetParName(const int i) const
Get name of parameter.
void SetParProp(const int i, const double val)
Set proposed parameter value.
double GetUpperBound(const int i) const
Get upper parameter bound in which it is physically valid.
TMatrixDSym * GetCovMatrix() const
Return covariance matrix.
std::string GetParFancyName(const int i) const
Get fancy name of the Parameter.
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
std::string GetName() const
Get name of covariance.
M3::float_t GetParProp(const int i) const
Get proposed parameter value.
double GetParInit(const int i) const
Get prior parameter value.
double GetDiagonalError(const int i) const
Get diagonal error for ith parameter.
YAML::Node GetConfig() const
Getter to return a copy of the YAML node.
Class responsible for handling implementation of samples used in analysis, reweighting and returning ...
virtual std::unique_ptr< TH1 > Get1DVarHistByModeAndChannel(const int iSample, const std::string &ProjectionVar_Str, const int kModeToFill=-1, const int kChannelToFill=-1, const int WeightStyle=0)=0
Build a 1D histogram for a given variable, optionally filtered by mode and channel.
virtual std::string GetName() const =0
Get name for Sample Handler.
virtual void SaveAdditionalInfo(TDirectory *Dir)
Store additional info in a chain.
virtual std::unique_ptr< TH1 > Get1DVarHist(const int iSample, const std::string &ProjectionVar, const std::vector< KinematicCut > &EventSelectionVec={}, int WeightStyle=0, const std::vector< KinematicCut > &SubEventSelectionVec={})=0
Return 1D projection of MC into given 1D variable (doesn't have to be variable used in the fit)
virtual std::string GetFlavourName(const int iSample, const int iChannel) const =0
Get the flavour name for a given sample and oscillation channel.
MaCh3Modes * GetMaCh3Modes() const
Return pointer to MaCh3 modes.
virtual int GetNOscChannels(const int iSample) const =0
Get number of oscillation channels for a single sample.
virtual M3::int_t GetNSamples()
returns total number of samples
virtual std::string GetSampleTitle(const int iSample) const =0
Get fancy title for specified samples.
virtual std::string GetKinVarName(const int iSample, const int Dimension) const =0
Return Kinematic Variable name for specified sample and dimension for example "Reconstructed_Neutrino...
virtual std::unique_ptr< TH2 > Get2DVarHist(const int iSample, const std::string &ProjectionVarX, const std::string &ProjectionVarY, const std::vector< KinematicCut > &EventSelectionVec={}, const int WeightStyle=0, const std::vector< KinematicCut > &SubEventSelectionVec={})=0
Build a 2D projection of MC events into specified variables.
virtual int GetNDim(const int Sample) const =0
DB Get what dimensionality binning for given sample has.
int getValue(const std::string &Type)
CW: Get info like RAM.
Definition: Monitor.cpp:252
void PrintProgressBar(const Long64_t Done, const Long64_t All)
KS: Simply print progress bar.
Definition: Monitor.cpp:229
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
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.
int GetNThreads()
number of threads which we need for example for TRandom3
Definition: Monitor.cpp:372