MaCh3  2.5.1
Reference Guide
ProcessMCMC.cpp
Go to the documentation of this file.
1 //MaCh3 includes
2 #include "Fitters/OscProcessor.h"
3 #include "Manager/Manager.h"
4 
10 
12 inline void ProcessMCMC(const std::string& inputFile);
14 inline void MultipleProcessMCMC();
15 inline void CalcBayesFactor(MCMCProcessor* Processor);
16 inline void CalcSavageDickey(MCMCProcessor* Processor);
17 inline void CalcBipolarPlot(MCMCProcessor* Processor);
18 inline void CalcParameterEvolution(MCMCProcessor* Processor);
19 inline void GetTrianglePlot(MCMCProcessor* Processor);
20 inline void DiagnoseCovarianceMatrix(MCMCProcessor* Processor, const std::string& inputFile);
22 inline TH2D* TMatrixIntoTH2D(TMatrixDSym* Matrix, const std::string& title);
24 inline void KolmogorovSmirnovTest(const std::vector<std::unique_ptr<MCMCProcessor>>& Processor,
25  const std::unique_ptr<TCanvas>& Posterior,
26  const TString& canvasname);
27 int nFiles;
28 std::vector <std::string> FileNames;
29 std::vector <std::string> TitleNames;
30 std::string config;
31 
32 int main(int argc, char *argv[])
33 {
35  nFiles = 0;
36  if (argc != 3 && argc !=6 && argc != 8)
37  {
38  MACH3LOG_ERROR("How to use: ");
39  MACH3LOG_ERROR(" single chain: {} <Config> <MCMM_ND_Output.root>", argv[0]);
40  MACH3LOG_ERROR(" two chain: {} <Config> <MCMM_ND_Output_1.root> <Title 1> <MCMC_ND_Output_2.root> <Title 2>", argv[0]);
41  MACH3LOG_ERROR(" three chain: {} <Config> <MCMM_ND_Output_1.root> <Title 1> <MCMC_ND_Output_2.root> <Title 2> <MCMC_ND_Output_3.root> <Title 3>", argv[0]);
42  throw MaCh3Exception(__FILE__ , __LINE__ );
43  }
44 
45  config = argv[1];
46  YAML::Node card_yaml = M3OpenConfig(config);
47  if (!CheckNodeExists(card_yaml, "ProcessMCMC")) {
48  MACH3LOG_ERROR("The 'ProcessMCMC' node is not defined in the YAML configuration.");
49  throw MaCh3Exception(__FILE__ , __LINE__ );
50  }
51 
52  if (argc == 3)
53  {
54  MACH3LOG_INFO("Producing single fit output");
55  std::string filename = argv[2];
56  ProcessMCMC(filename);
57  }
58  // If we want to compare two or more fits (e.g. binning changes or introducing new params/priors)
59  else if (argc == 6 || argc == 8)
60  {
61  MACH3LOG_INFO("Producing two fit comparison");
62  FileNames.push_back(argv[2]);
63  TitleNames.push_back(argv[3]);
64 
65  FileNames.push_back(argv[4]);
66  TitleNames.push_back(argv[5]);
67  //KS: If there is third file add it
68  if(argc == 8)
69  {
70  FileNames.push_back(argv[6]);
71  TitleNames.push_back(argv[7]);
72  }
73 
75  }
76  return 0;
77 }
78 
93 std::map<std::string, std::pair<double, double>> GetCustomBinning(const YAML::Node& Settings)
94 {
95  std::map<std::string, std::pair<double, double>> CustomBinning;
96  if (Settings["CustomBinEdges"]) {
97  const YAML::Node& edges = Settings["CustomBinEdges"];
98 
99  for (const auto& node : edges) {
100  std::string key = node.first.as<std::string>();
101  auto values = node.second.as<std::vector<double>>();
102 
103  if (values.size() == 2) {
104  CustomBinning[key] = std::make_pair(values[0], values[1]);
105  MACH3LOG_DEBUG("Adding custom binning {} with {:.4f}, {:.4f}", key, values[0], values[1]);
106  } else {
107  MACH3LOG_ERROR("Invalid number of values for key: {}", key);
108  throw MaCh3Exception(__FILE__ , __LINE__ );
109  }
110  }
111  }
112  return CustomBinning;
113 }
114 
115 void ProcessMCMC(const std::string& inputFile)
116 {
117  MACH3LOG_INFO("File for study: {} with config {}", inputFile, config);
118  // Make the processor)
119  auto Processor = std::make_unique<OscProcessor>(inputFile);
120 
121  YAML::Node card_yaml = M3OpenConfig(config.c_str());
122  YAML::Node Settings = card_yaml["ProcessMCMC"];
123 
124  const bool PlotCorr = GetFromManager<bool>(Settings["PlotCorr"], false);
125 
126  Processor->SetExcludedTypes(GetFromManager<std::vector<std::string>>(Settings["ExcludedTypes"], {}));
127  Processor->SetExcludedNames(GetFromManager<std::vector<std::string>>(Settings["ExcludedNames"], {}));
128  Processor->SetExcludedGroups(GetFromManager<std::vector<std::string>>(Settings["ExcludedGroups"], {}));
129 
130  //Apply additional cuts to 1D posterior
131  Processor->SetPosterior1DCut(GetFromManager<std::string>(Settings["Posterior1DCut"], ""));
132 
133  if(PlotCorr) Processor->SetOutputSuffix("_drawCorr");
134  //KS:Turn off plotting detector and some other setting, should be via some config
135  Processor->SetPlotRelativeToPrior(GetFromManager<bool>(Settings["PlotRelativeToPrior"], false));
136  Processor->SetPrintToPDF(GetFromManager<bool>(Settings["PrintToPDF"], true));
137 
138  //KS: Whether you want prior error bands for parameters with flat prior or not
139  Processor->SetPlotErrorForFlatPrior(GetFromManager<bool>(Settings["PlotErrorForFlatPrior"], true));
140  Processor->SetFancyNames(GetFromManager<bool>(Settings["FancyNames"], true));
141  Processor->SetPlotBinValue(GetFromManager<bool>(Settings["PlotBinValue"], false));
142  //KS: Plot only 2D posteriors with correlations greater than 0.2
143  Processor->SetPost2DPlotThreshold(GetFromManager<double>(Settings["Post2DPlotThreshold"], 0.2));
144 
145  Processor->Initialise();
146 
147  if(Settings["BurnInSteps"]) {
148  Processor->SetStepCut(Settings["BurnInSteps"].as<int>());
149  } else {
150  MACH3LOG_WARN("BurnInSteps not set, defaulting to 20%");
151  Processor->SetStepCut(static_cast<int>(Processor->GetnSteps()/5));
152  }
153  if(Settings["MaxEntries"]) {
154  Processor->SetEntries(Get<int>(Settings["MaxEntries"], __FILE__, __LINE__));
155  }
156  if(Settings["NBins"]) {
157  Processor->SetNBins(Get<int>(Settings["NBins"], __FILE__, __LINE__));
158  }
159  if(Settings["Thinning"])
160  {
161  if(Settings["Thinning"][0].as<bool>()){
162  Processor->ThinMCMC(Settings["Thinning"][1].as<int>());
163  }
164  }
165  // Make the postfit
166  Processor->MakePostfit(GetCustomBinning(Settings));
167  Processor->DrawPostfit();
168  //KS: Should set via config whether you want below or not
169  if(GetFromManager<bool>(Settings["MakeCredibleIntervals"], true)) {
170  Processor->MakeCredibleIntervals(GetFromManager<std::vector<double>>(Settings["CredibleIntervals"], {0.99, 0.90, 0.68}),
171  GetFromManager<std::vector<short int>>(Settings["CredibleIntervalsColours"], {436, 430, 422}),
172  GetFromManager<bool>(Settings["CredibleInSigmas"], false));
173  }
174  if(GetFromManager<bool>(Settings["CalcBayesFactor"], true)) CalcBayesFactor(Processor.get());
175  if(GetFromManager<bool>(Settings["CalcSavageDickey"], true)) CalcSavageDickey(Processor.get());
176  if(GetFromManager<bool>(Settings["CalcBipolarPlot"], false)) CalcBipolarPlot(Processor.get());
177  if(GetFromManager<bool>(Settings["CalcParameterEvolution"], false)) CalcParameterEvolution(Processor.get());
178 
179  if(PlotCorr)
180  {
181  Processor->SetSmoothing(GetFromManager<bool>(Settings["Smoothing"], true));
182  // Make the covariance matrix
183  //We have different treatment for multithread
184  Processor->CacheSteps();
185  //KS: Since we cached let's make fancy violins :)
186  if(GetFromManager<bool>(Settings["MakeViolin"], true)) Processor->MakeViolin();
187  Processor->MakeCovariance_MP();
188 
189  Processor->DrawCovariance();
190  if(GetFromManager<bool>(Settings["MakeCovarianceYAML"], true)) Processor->MakeCovarianceYAML(GetFromManager<std::string>(Settings["CovarianceYAMLOutName"], "UpdatedCorrelationMatrix.yaml"), GetFromManager<std::string>(Settings["CovarianceYAMLMeansMethod"], "HPD"));
191 
192  auto const &MakeSubOptimality = Settings["MakeSubOptimality"];
193  if(MakeSubOptimality[0].as<bool>()) Processor->MakeSubOptimality(MakeSubOptimality[1].as<int>());
194 
195  if(GetFromManager<bool>(Settings["MakeCredibleRegions"], false)) {
196  Processor->MakeCredibleRegions(GetFromManager<std::vector<double>>(Settings["CredibleRegions"], {0.99, 0.90, 0.68}),
197  GetFromManager<std::vector<short int>>(Settings["CredibleRegionStyle"], {2, 1, 3}),
198  GetFromManager<std::vector<short int>>(Settings["CredibleRegionColor"], {413, 406, 416}),
199  GetFromManager<bool>(Settings["CredibleInSigmas"], false),
200  GetFromManager<bool>(Settings["Draw2DPosterior"], true),
201  GetFromManager<bool>(Settings["DrawBestFit"], true));
202  }
203  if(GetFromManager<bool>(Settings["GetTrianglePlot"], true)) GetTrianglePlot(Processor.get());
204 
205  //KS: When creating covariance matrix longest time is spend on caching every step, since we already cached we can run some fancy covariance stability diagnostic
206  if(GetFromManager<bool>(Settings["DiagnoseCovarianceMatrix"], false)) DiagnoseCovarianceMatrix(Processor.get(), inputFile);
207  }
208  Processor->ProduceChi2(GetFromManager<std::string>(Settings["Chi2Group"], "Osc"));
209  if(GetFromManager<bool>(Settings["JarlskogAnalysis"], true)) Processor->PerformJarlskogAnalysis();
210  if(GetFromManager<bool>(Settings["ProducePMNSElements"], true)) Processor->ProducePMNSElements();
211  if(GetFromManager<bool>(Settings["ProduceUnitarityTriangles"], true)) Processor->ProduceUnitarityTriangles();
212  if(GetFromManager<bool>(Settings["MakePiePlot"], true)) Processor->MakePiePlot();
213 }
214 
216 {
217  YAML::Node card_yaml = M3OpenConfig(config.c_str());
218  YAML::Node Settings = card_yaml["ProcessMCMC"];
219 
220  constexpr Color_t PosteriorColor[] = {kBlue-1, kRed, kGreen+2};
221  //constexpr Style_t PosteriorStyle[] = {kSolid, kDashed, kDotted};
222  nFiles = int(FileNames.size());
223  std::vector<std::unique_ptr<MCMCProcessor>> Processor(nFiles);
224 
225  if(!Settings["BurnInSteps"]) {
226  MACH3LOG_WARN("BurnInSteps not set, defaulting to 20%");
227  }
228 
229  for (int ik = 0; ik < nFiles; ik++)
230  {
231  MACH3LOG_INFO("File for study: {}", FileNames[ik]);
232  // Make the processor
233  Processor[ik] = std::make_unique<MCMCProcessor>(FileNames[ik]);
234  Processor[ik]->SetOutputSuffix(("_" + std::to_string(ik)).c_str());
235 
236  Processor[ik]->SetExcludedTypes(GetFromManager<std::vector<std::string>>(Settings["ExcludedTypes"], {}));
237  Processor[ik]->SetExcludedNames(GetFromManager<std::vector<std::string>>(Settings["ExcludedNames"], {}));
238  Processor[ik]->SetExcludedGroups(GetFromManager<std::vector<std::string>>(Settings["ExcludedGroups"], {}));
239 
240  //Apply additional cuts to 1D posterior
241  Processor[ik]->SetPosterior1DCut(GetFromManager<std::string>(Settings["Posterior1DCut"], ""));
242 
243  Processor[ik]->SetPlotRelativeToPrior(GetFromManager<bool>(Settings["PlotRelativeToPrior"], false));
244  Processor[ik]->SetFancyNames(GetFromManager<bool>(Settings["FancyNames"], true));
245  Processor[ik]->Initialise();
246 
247  if(Settings["BurnInSteps"]) {
248  Processor[ik]->SetStepCut(Settings["BurnInSteps"].as<int>());
249  }else {
250  Processor[ik]->SetStepCut(static_cast<int>(Processor[ik]->GetnSteps()/5));
251  }
252 
253  if(Settings["MaxEntries"]) {
254  Processor[ik]->SetEntries(Get<int>(Settings["MaxEntries"], __FILE__, __LINE__));
255  }
256  if(Settings["NBins"]) {
257  Processor[ik]->SetNBins(Get<int>(Settings["NBins"], __FILE__, __LINE__));
258  }
259  }
260 
261  Processor[0]->MakePostfit(GetCustomBinning(Settings));
262  Processor[0]->DrawPostfit();
263  // Get edges from first histogram to ensure all params use same binning
264  std::map<std::string, std::pair<double, double>> ParamEdges;
265  for(int i = 0; i < Processor[0]->GetNParams(); ++i) {
266  // Get the histogram for the i-th parameter
267  TH1D* hist = Processor[0]->GetHpost(i);
268  if (!hist) {
269  MACH3LOG_DEBUG("Histogram for parameter {} is null.", i);
270  continue;
271  }
272 
273  // Get the parameter name (title of the histogram)
274  std::string paramName = hist->GetTitle();
275 
276  // Get the axis limits (edges)
277  TAxis* axis = hist->GetXaxis();
278  double xmin = axis->GetXmin();
279  double xmax = axis->GetXmax();
280 
281  MACH3LOG_DEBUG("Adding bin edges for {} equal to {:.4f}, {:.4f}",paramName, xmin, xmax);
282  // Insert into the map
283  ParamEdges[paramName] = std::make_pair(xmin, xmax);
284  }
285 
286  //KS: Multithreading here is very tempting but there are some issues with root that need to be resovled :(
287  for (int ik = 1; ik < nFiles; ik++)
288  {
289  // Make the postfit
290  Processor[ik]->MakePostfit(ParamEdges);
291  Processor[ik]->DrawPostfit();
292  }
293 
294  // Open a TCanvas to write the posterior onto
295  auto Posterior = std::make_unique<TCanvas>("PosteriorMulti", "PosteriorMulti", 0, 0, 1024, 1024);
296  gStyle->SetOptStat(0);
297  gStyle->SetOptTitle(0);
298  Posterior->SetGrid();
299  Posterior->SetBottomMargin(0.1f);
300  Posterior->SetTopMargin(0.05f);
301  Posterior->SetRightMargin(0.03f);
302  Posterior->SetLeftMargin(0.15f);
303 
304  // First filename: keep path, just remove ".root"
305  // Would be nice to specify outpath in a later update
306  size_t pos = FileNames[0].rfind(".root");
307  std::string base = (pos == std::string::npos) ? FileNames[0] : FileNames[0].substr(0, pos);
308  TString canvasname = base;
309 
310  // Remaining filenames: strip path and ".root"
311  // So if you have /path/to/file1.root and /path/to/file2.root or /another/path/to/file2.root, canvasname = /path/to/file1_file2.root
312  for (int ik = 1; ik < nFiles; ik++) {
313  pos = FileNames[ik].find_last_of('/');
314  base = (pos == std::string::npos) ? FileNames[ik] : FileNames[ik].substr(pos + 1);
315  pos = base.rfind(".root");
316  if (pos != std::string::npos) base = base.substr(0, pos);
317  canvasname += "_" + TString(base);
318  }
319 
320  canvasname = canvasname +".pdf[";
321 
322  Posterior->Print(canvasname);
323  // Once the pdf file is open no longer need to bracket
324  canvasname.ReplaceAll("[","");
325 
326  for(int i = 0; i < Processor[0]->GetNParams(); ++i)
327  {
328  // This holds the posterior density
329  std::vector<std::unique_ptr<TH1D>> hpost(nFiles);
330  std::vector<std::unique_ptr<TLine>> hpd(nFiles);
331  hpost[0] = M3::Clone(Processor[0]->GetHpost(i));
332  hpost[0]->GetYaxis()->SetTitle("Posterior Density");
333  bool Skip = false;
334  for (int ik = 1 ; ik < nFiles; ik++)
335  {
336  // KS: If somehow this chain doesn't given params we skip it
337  const int Index = Processor[ik]->GetParamIndexFromName(hpost[0]->GetTitle());
338  if(Index == M3::_BAD_INT_)
339  {
340  Skip = true;
341  break;
342  }
343  hpost[ik] = M3::Clone(Processor[ik]->GetHpost(Index));
344  }
345 
346  // Don't plot if this is a fixed histogram (i.e. the peak is the whole integral)
347  if(hpost[0]->GetMaximum() == hpost[0]->Integral()*1.5 || Skip) {
348  continue;
349  }
350  for (int ik = 0; ik < nFiles; ik++)
351  {
352  RemoveFitter(hpost[ik].get(), "Gauss");
353 
354  // Set some nice colours
355  hpost[ik]->SetLineColor(PosteriorColor[ik]);
356  //hpost[ik]->SetLineStyle(PosteriorStyle[ik]);
357  hpost[ik]->SetLineWidth(2);
358 
359  // Area normalise the distributions
360  hpost[ik]->Scale(1./hpost[ik]->Integral());
361  }
362  TString Title;
363  double Prior = 1.0;
364  double PriorError = 1.0;
365 
366  Processor[0]->GetNthParameter(i, Prior, PriorError, Title);
367 
368  // Now make the TLine for the Asimov
369  auto Asimov = std::make_unique<TLine>(Prior, hpost[0]->GetMinimum(), Prior, hpost[0]->GetMaximum());
370  Asimov->SetLineColor(kRed-3);
371  Asimov->SetLineWidth(2);
372  Asimov->SetLineStyle(kDashed);
373 
374  // Make a nice little TLegend
375  auto leg = std::make_unique<TLegend>(0.20, 0.7, 0.6, 0.97);
376  leg->SetTextSize(0.03f);
377  leg->SetFillColor(0);
378  leg->SetFillStyle(0);
379  leg->SetLineColor(0);
380  leg->SetLineStyle(0);
381  TString asimovLeg = Form("#splitline{Prior}{x = %.2f , #sigma = %.2f}", Prior, PriorError);
382  leg->AddEntry(Asimov.get(), asimovLeg, "l");
383 
384  for (int ik = 0; ik < nFiles; ik++)
385  {
386  TString rebinLeg = Form("#splitline{%s}{#mu = %.2f, #sigma = %.2f}", TitleNames[ik].c_str(), hpost[ik]->GetMean(), hpost[ik]->GetRMS());
387  leg->AddEntry(hpost[ik].get(), rebinLeg, "l");
388 
389  hpd[ik] = std::make_unique<TLine>(hpost[ik]->GetBinCenter(hpost[ik]->GetMaximumBin()), hpost[ik]->GetMinimum(),
390  hpost[ik]->GetBinCenter(hpost[ik]->GetMaximumBin()), hpost[ik]->GetMaximum());
391  hpd[ik]->SetLineColor(hpost[ik]->GetLineColor());
392  hpd[ik]->SetLineWidth(2);
393  hpd[ik]->SetLineStyle(kSolid);
394  }
395 
396  // Find the maximum value to nicely resize hist
397  double maximum = 0;
398  for (int ik = 0; ik < nFiles; ik++) maximum = std::max(maximum, hpost[ik]->GetMaximum());
399  for (int ik = 0; ik < nFiles; ik++) hpost[ik]->SetMaximum(1.3*maximum);
400 
401  hpost[0]->Draw("hist");
402  for (int ik = 1; ik < nFiles; ik++) hpost[ik]->Draw("hist same");
403  Asimov->Draw("same");
404  for (int ik = 0; ik < nFiles; ik++) hpd[ik]->Draw("same");
405  leg->Draw("same");
406  Posterior->cd();
407  Posterior->Print(canvasname);
408  }//End loop over parameters
409 
410  // Finally draw the parameter plot onto the PDF
411  // Close the .pdf file with all the posteriors
412  Posterior->cd();
413  Posterior->Clear();
414 
415  if(GetFromManager<bool>(Settings["PerformKStest"], true)) KolmogorovSmirnovTest(Processor, Posterior, canvasname);
416 
417  // Close the pdf file
418  MACH3LOG_INFO("Closing pdf {}", canvasname);
419  canvasname+="]";
420  Posterior->Print(canvasname);
421 }
422 
425 {
426  YAML::Node card_yaml = M3OpenConfig(config.c_str());
427  YAML::Node Settings = card_yaml["ProcessMCMC"];
428 
429  std::vector<std::string> ParNames;
430  std::vector<std::vector<double>> Model1Bounds;
431  std::vector<std::vector<double>> Model2Bounds;
432  std::vector<std::vector<std::string>> ModelNames;
433  for (const auto& dg : Settings["BayesFactor"])
434  {
435  ParNames.push_back(dg[0].as<std::string>());
436  ModelNames.push_back(dg[1].as<std::vector<std::string>>());
437  Model1Bounds.push_back(dg[2].as<std::vector<double>>());
438  Model2Bounds.push_back(dg[3].as<std::vector<double>>());
439  }
440 
441  Processor->GetBayesFactor(ParNames, Model1Bounds, Model2Bounds, ModelNames);
442 }
443 
445 {
446  YAML::Node card_yaml = M3OpenConfig(config.c_str());
447  YAML::Node Settings = card_yaml["ProcessMCMC"];
448 
449  std::vector<std::string> ParNames;
450  std::vector<double> EvaluationPoint;
451  std::vector<std::vector<double>> Bounds;
452 
453  for (const auto& d : Settings["SavageDickey"])
454  {
455  ParNames.push_back(d[0].as<std::string>());
456  EvaluationPoint.push_back(d[1].as<double>());
457  Bounds.push_back(d[2].as<std::vector<double>>());
458  }
459  Processor->GetSavageDickey(ParNames, EvaluationPoint, Bounds);
460 }
461 
463 {
464  YAML::Node card_yaml = M3OpenConfig(config.c_str());
465  YAML::Node Settings = card_yaml["ProcessMCMC"];
466 
467  std::vector<std::string> ParNames;
468  std::vector<int> Intervals;
469  for (const auto& d : Settings["ParameterEvolution"])
470  {
471  ParNames.push_back(d[0].as<std::string>());
472  Intervals.push_back(d[1].as<int>());
473  }
474  Processor->ParameterEvolution(ParNames, Intervals);
475 }
476 
478 {
479  YAML::Node card_yaml = M3OpenConfig(config.c_str());
480  YAML::Node Settings = card_yaml["ProcessMCMC"];
481 
482  std::vector<std::string> ParNames;
483  for (const auto& d : Settings["BipolarPlot"])
484  {
485  ParNames.push_back(d[0].as<std::string>());
486  }
487  Processor->GetPolarPlot(ParNames);
488 }
489 
490 void GetTrianglePlot(MCMCProcessor* Processor) {
491  YAML::Node card_yaml = M3OpenConfig(config.c_str());
492  YAML::Node Settings = card_yaml["ProcessMCMC"];
493 
494  for (const auto& dg : Settings["TrianglePlot"])
495  {
496  std::string ParName = dg[0].as<std::string>();
497 
498  std::vector<std::string> NameVec = dg[1].as<std::vector<std::string>>();
499  Processor->MakeTrianglePlot(NameVec,
500  GetFromManager<std::vector<double>>(Settings["CredibleIntervals"], {0.99, 0.90, 0.68}),
501  GetFromManager<std::vector<short int>>(Settings["CredibleIntervalsColours"], {436, 430, 422}),
502  GetFromManager<std::vector<double>>(Settings["CredibleRegions"], {0.99, 0.90, 0.68}),
503  GetFromManager<std::vector<short int>>(Settings["CredibleRegionStyle"], {2, 1, 3}),
504  GetFromManager<std::vector<short int>>(Settings["CredibleRegionColor"], {413, 406, 416}),
505  GetFromManager<bool>(Settings["CredibleInSigmas"], false));
506  }
507 }
508 
510 void DiagnoseCovarianceMatrix(MCMCProcessor* Processor, const std::string& inputFile)
511 {
512  //Turn of plots from Processor
513  Processor->SetPrintToPDF(false);
514  // Open a TCanvas to write the posterior onto
515  auto Canvas = std::make_unique<TCanvas>("Canvas", "Canvas", 0, 0, 1024, 1024);
516  Canvas->SetGrid();
517  gStyle->SetOptStat(0);
518  gStyle->SetOptTitle(0);
519  Canvas->SetTickx();
520  Canvas->SetTicky();
521  Canvas->SetBottomMargin(0.1f);
522  Canvas->SetTopMargin(0.05f);
523  Canvas->SetRightMargin(0.15f);
524  Canvas->SetLeftMargin(0.10f);
525 
526  //KS: Fancy colours
527  const int NRGBs = 10;
528  TColor::InitializeColors();
529  Double_t stops[NRGBs] = { 0.00, 0.10, 0.25, 0.35, 0.50, 0.60, 0.65, 0.75, 0.90, 1.00 };
530  Double_t red[NRGBs] = { 0.50, 1.00, 1.00, 0.25, 0.00, 0.10, 0.50, 1.00, 0.75, 0.55 };
531  Double_t green[NRGBs] = { 0.00, 0.25, 1.00, 0.25, 0.00, 0.60, 0.90, 1.00, 0.75, 0.75 };
532  Double_t blue[NRGBs] = { 0.00, 0.25, 1.00, 1.00, 0.50, 0.60, 0.90, 1.00, 0.05, 0.05 };
533  TColor::CreateGradientColorTable(NRGBs, stops, red, green, blue, 255);
534  gStyle->SetNumberContours(255);
535 
536  std::string OutName = inputFile;
537  OutName = OutName.substr(0, OutName.find(".root"));
538  Canvas->Print(Form("Correlation_%s.pdf[", OutName.c_str()), "pdf");
539  Canvas->Print(Form("Covariance_%s.pdf[", OutName.c_str()), "pdf");
540 
541  YAML::Node card_yaml = M3OpenConfig(config.c_str());
542  YAML::Node Settings = card_yaml["ProcessMCMC"];
543 
544  const int entries = int(Processor->GetnSteps());
545  const int NIntervals = GetFromManager<int>(Settings["NIntervals"], 5);
546  const int IntervalsSize = entries/NIntervals;
547  //We start with burn from 0 (no burn in at all)
548  int BurnIn = 0;
549  MACH3LOG_INFO("Diagnosing matrices with entries={}, NIntervals={} and IntervalsSize={}", entries, NIntervals, IntervalsSize);
550 
551  TMatrixDSym *Covariance = nullptr;
552  TMatrixDSym *Correlation = nullptr;
553 
554  TH2D *CovariancePreviousHist = nullptr;
555  TH2D *CorrelationPreviousHist = nullptr;
556 
557  TH2D *CovarianceHist = nullptr;
558  TH2D *CorrelationHist = nullptr;
559 
560  //KS: Get first covariances, we need two for comparison...
561  Processor->SetStepCut(BurnIn);
562  Processor->GetCovariance(Covariance, Correlation);
563 
564  CovariancePreviousHist = TMatrixIntoTH2D(Covariance, "Covariance");
565  CorrelationPreviousHist = TMatrixIntoTH2D(Correlation, "Correlation");
566 
567  delete Covariance;
568  Covariance = nullptr;
569  delete Correlation;
570  Correlation = nullptr;
571 
572  //KS: Loop over all desired cuts
573  for(int k = 1; k < NIntervals; ++k)
574  {
575  BurnIn = k*IntervalsSize;
576  Processor->SetStepCut(BurnIn);
577  Processor->GetCovariance(Covariance, Correlation);
578  Processor->Reset2DPosteriors();
579 
580  CovarianceHist = TMatrixIntoTH2D(Covariance, "Covariance");
581  CorrelationHist = TMatrixIntoTH2D(Correlation, "Correlation");
582 
583  TH2D *CovarianceDiff = static_cast<TH2D*>(CovarianceHist->Clone("Covariance_Ratio"));
584  TH2D *CorrelationDiff = static_cast<TH2D*>(CorrelationHist->Clone("Correlation_Ratio"));
585 
586  //KS: Bit messy but quite often covariance is 0 is divided by 0 is problematic so
587  #ifdef MULTITHREAD
588  #pragma omp parallel for
589  #endif
590  for (int j = 1; j < CovarianceDiff->GetXaxis()->GetNbins()+1; ++j)
591  {
592  for (int i = 1; i < CovarianceDiff->GetYaxis()->GetNbins()+1; ++i)
593  {
594  if( std::fabs (CovarianceDiff->GetBinContent(j, i)) < 1.e-5 && std::fabs (CovariancePreviousHist->GetBinContent(j, i)) < 1.e-5)
595  {
596  CovarianceDiff->SetBinContent(j, i, M3::_BAD_DOUBLE_);
597  CovariancePreviousHist->SetBinContent(j, i, M3::_BAD_DOUBLE_);
598  }
599  if( std::fabs (CorrelationDiff->GetBinContent(j, i)) < 1.e-5 && std::fabs (CorrelationPreviousHist->GetBinContent(j, i)) < 1.e-5)
600  {
601  CorrelationDiff->SetBinContent(j, i, M3::_BAD_DOUBLE_);
602  CorrelationPreviousHist->SetBinContent(j, i, M3::_BAD_DOUBLE_);
603  }
604  }
605  }
606  //Divide matrices
607  CovarianceDiff->Divide(CovariancePreviousHist);
608  CorrelationDiff->Divide(CorrelationPreviousHist);
609 
610  //Now it is time for fancy names etc.
611  for (int j = 0; j < CovarianceDiff->GetXaxis()->GetNbins(); ++j)
612  {
613  TString Title = "";
614  double Prior = 1.0;
615  double PriorError = 1.0;
616 
617  Processor->GetNthParameter(j, Prior, PriorError, Title);
618 
619  CovarianceDiff->GetXaxis()->SetBinLabel(j+1, Title);
620  CovarianceDiff->GetYaxis()->SetBinLabel(j+1, Title);
621  CorrelationDiff->GetXaxis()->SetBinLabel(j+1, Title);
622  CorrelationDiff->GetYaxis()->SetBinLabel(j+1, Title);
623  }
624  CovarianceDiff->GetXaxis()->SetLabelSize(0.015f);
625  CovarianceDiff->GetYaxis()->SetLabelSize(0.015f);
626  CorrelationDiff->GetXaxis()->SetLabelSize(0.015f);
627  CorrelationDiff->GetYaxis()->SetLabelSize(0.015f);
628 
629  std::stringstream ss;
630  ss << "BCut_";
631  ss << BurnIn;
632  ss << "/";
633  ss << "BCut_";
634  ss << (k-1)*IntervalsSize;
635  std::string str = ss.str();
636 
637  TString Title = "Cov " + str;
638  CovarianceDiff->GetZaxis()->SetTitle( Title );
639  Title = "Corr " + str;
640  CorrelationDiff->GetZaxis()->SetTitle(Title);
641 
642  CovarianceDiff->SetMinimum(-2);
643  CovarianceDiff->SetMaximum(2);
644  CorrelationDiff->SetMinimum(-2);
645  CorrelationDiff->SetMaximum(2);
646 
647  Canvas->cd();
648  CovarianceDiff->Draw("colz");
649  Canvas->Print(Form("Covariance_%s.pdf", OutName.c_str()), "pdf");
650 
651  CorrelationDiff->Draw("colz");
652  Canvas->Print(Form("Correlation_%s.pdf", OutName.c_str()), "pdf");
653 
654  //KS: Current hist become previous as we need it for further comparison
655  delete CovariancePreviousHist;
656  CovariancePreviousHist = static_cast<TH2D*>(CovarianceHist->Clone());
657  delete CorrelationPreviousHist;
658  CorrelationPreviousHist = static_cast<TH2D*>(CorrelationHist->Clone());
659 
660  delete CovarianceHist;
661  CovarianceHist = nullptr;
662  delete CorrelationHist;
663  CorrelationHist = nullptr;
664 
665  delete CovarianceDiff;
666  delete CorrelationDiff;
667  delete Covariance;
668  Covariance = nullptr;
669  delete Correlation;
670  Correlation = nullptr;
671  }
672  Canvas->cd();
673  Canvas->Print(Form("Covariance_%s.pdf]", OutName.c_str()), "pdf");
674  Canvas->Print(Form("Correlation_%s.pdf]", OutName.c_str()), "pdf");
675 
676  Processor->SetPrintToPDF(true);
677  if(Covariance != nullptr) delete Covariance;
678  if(Correlation != nullptr) delete Correlation;
679  if(CovariancePreviousHist != nullptr) delete CovariancePreviousHist;
680  if(CorrelationPreviousHist != nullptr) delete CorrelationPreviousHist;
681  if(CovarianceHist != nullptr) delete CovarianceHist;
682  if(CorrelationHist != nullptr) delete CorrelationHist;
683 }
684 
685 //KS: Convert TMatrix to TH2D, mostly useful for making fancy plots
686 TH2D* TMatrixIntoTH2D(TMatrixDSym* Matrix, const std::string& title)
687 {
688  TH2D* hMatrix = new TH2D(title.c_str(), title.c_str(), Matrix->GetNrows(), 0.0, Matrix->GetNrows(), Matrix->GetNcols(), 0.0, Matrix->GetNcols());
689  for(int i = 0; i < Matrix->GetNrows(); i++)
690  {
691  for(int j = 0; j < Matrix->GetNcols(); j++)
692  {
693  //KS: +1 because there is offset in histogram relative to TMatrix
694  hMatrix->SetBinContent(i+1,j+1, (*Matrix)(i,j));
695  }
696  }
697  return hMatrix;
698 }
699 
700 // KS: Perform KS test to check if two posteriors for the same parameter came from the same distribution
701 void KolmogorovSmirnovTest(const std::vector<std::unique_ptr<MCMCProcessor>>& Processor,
702  const std::unique_ptr<TCanvas>& Posterior,
703  const TString& canvasname)
704 {
705  constexpr Color_t CumulativeColor[] = {kBlue-1, kRed, kGreen+2};
706  constexpr Style_t CumulativeStyle[] = {kSolid, kDashed, kDotted};
707 
708  for(int i = 0; i < Processor[0]->GetNParams(); ++i)
709  {
710  // This holds the posterior density
711  std::vector<std::unique_ptr<TH1D>> hpost(nFiles);
712  std::vector<std::unique_ptr<TH1D>> CumulativeDistribution(nFiles);
713 
714  TString Title;
715  double Prior = 1.0;
716  double PriorError = 1.0;
717 
718  Processor[0]->GetNthParameter(i, Prior, PriorError, Title);
719  bool Skip = false;
720  for (int ik = 0 ; ik < nFiles; ik++)
721  {
722  int Index = 0;
723  if(ik == 0 ) Index = i;
724  else
725  {
726  // KS: If somehow this chain doesn't given params we skip it
727  Index = Processor[ik]->GetParamIndexFromName(hpost[0]->GetTitle());
728  if(Index == M3::_BAD_INT_)
729  {
730  Skip = true;
731  break;
732  }
733  }
734  hpost[ik] = M3::Clone(Processor[ik]->GetHpost(Index));
735  CumulativeDistribution[ik] = M3::Clone(Processor[ik]->GetHpost(Index));
736  CumulativeDistribution[ik]->Fill(0., 0.);
737  CumulativeDistribution[ik]->Reset();
738  CumulativeDistribution[ik]->SetMaximum(1.);
739  TString TempTitle = Title+" Kolmogorov Smirnov";
740  CumulativeDistribution[ik]->SetTitle(TempTitle);
741 
742  TempTitle = Title+" Value";
743  CumulativeDistribution[ik]->GetXaxis()->SetTitle(TempTitle);
744  CumulativeDistribution[ik]->GetYaxis()->SetTitle("Cumulative Probability");
745 
746  CumulativeDistribution[ik]->SetLineWidth(2);
747  CumulativeDistribution[ik]->SetLineColor(CumulativeColor[ik]);
748  CumulativeDistribution[ik]->SetLineStyle(CumulativeStyle[ik]);
749  }
750 
751  // Don't plot if this is a fixed histogram (i.e. the peak is the whole integral)
752  if(hpost[0]->GetMaximum() == hpost[0]->Integral()*1.5 || Skip) {
753  continue;
754  }
755 
756  for (int ik = 0 ; ik < nFiles; ik++)
757  {
758  const int NumberOfBins = hpost[ik]->GetXaxis()->GetNbins();
759  double Cumulative = 0;
760  const double Integral = hpost[ik]->Integral();
761  for (int j = 1; j < NumberOfBins+1; ++j)
762  {
763  Cumulative += hpost[ik]->GetBinContent(j)/Integral;
764  CumulativeDistribution[ik]->SetBinContent(j, Cumulative);
765  }
766  //KS: Set overflow to 1 just in case
767  CumulativeDistribution[ik]->SetBinContent(NumberOfBins+1, 1.);
768  }
769 
770  std::vector<int> TestStatBin(nFiles, 0);
771  std::vector<double> TestStatD(nFiles, -999);
772  std::vector<std::unique_ptr<TLine>> LineD(nFiles);
773  //Find KS statistic
774  for (int ik = 1 ; ik < nFiles; ik++)
775  {
776  const int NumberOfBins = CumulativeDistribution[0]->GetXaxis()->GetNbins();
777  for (int j = 1; j < NumberOfBins+1; ++j)
778  {
779  const double BinValue = CumulativeDistribution[0]->GetBinCenter(j);
780  const int BinNumber = CumulativeDistribution[ik]->FindBin(BinValue);
781  //KS: Calculate D statistic for this bin, only save it if it's bigger than previously found value
782  double TempDstat = std::fabs(CumulativeDistribution[0]->GetBinContent(j) - CumulativeDistribution[ik]->GetBinContent(BinNumber));
783  if(TempDstat > TestStatD[ik])
784  {
785  TestStatD[ik] = TempDstat;
786  TestStatBin[ik] = j;
787  }
788  }
789  }
790 
791  for (int ik = 0 ; ik < nFiles; ik++)
792  {
793  LineD[ik] = std::make_unique<TLine>(CumulativeDistribution[0]->GetBinCenter(TestStatBin[ik]), 0, CumulativeDistribution[0]->GetBinCenter(TestStatBin[ik]), CumulativeDistribution[0]->GetBinContent(TestStatBin[ik]));
794  LineD[ik]->SetLineColor(CumulativeColor[ik]);
795  LineD[ik]->SetLineWidth(2.0);
796  }
797  CumulativeDistribution[0]->Draw();
798  for (int ik = 0 ; ik < nFiles; ik++)
799  CumulativeDistribution[ik]->Draw("SAME");
800 
801  auto leg = std::make_unique<TLegend>(0.15, 0.7, 0.5, 0.90);
802  leg->SetTextSize(0.04f);
803  for (int ik = 0; ik < nFiles; ik++)
804  leg->AddEntry(CumulativeDistribution[ik].get(), TitleNames[ik].c_str(), "l");
805  for (int ik = 1; ik < nFiles; ik++)
806  leg->AddEntry(LineD[ik].get(), Form("#Delta D = %.4f", TestStatD[ik]), "l");
807 
808  leg->SetLineColor(0);
809  leg->SetLineStyle(0);
810  leg->SetFillColor(0);
811  leg->SetFillStyle(0);
812  leg->Draw("SAME");
813 
814  for (int ik = 1; ik < nFiles; ik++)
815  LineD[ik]->Draw("sam");
816 
817  Posterior->cd();
818  Posterior->Print(canvasname);
819  } //End loop over parameter
820 }
void RemoveFitter(TH1D *hist, const std::string &name)
KS: Remove fitted TF1 from hist to make comparison easier.
#define MACH3LOG_DEBUG
Definition: MaCh3Logger.h:34
#define MACH3LOG_ERROR
Definition: MaCh3Logger.h:37
#define MACH3LOG_INFO
Definition: MaCh3Logger.h:35
void SetMaCh3LoggerFormat()
Set messaging format of the logger.
Definition: MaCh3Logger.h:61
#define MACH3LOG_WARN
Definition: MaCh3Logger.h:36
void GetTrianglePlot(MCMCProcessor *Processor)
void DiagnoseCovarianceMatrix(MCMCProcessor *Processor, const std::string &inputFile)
KS: You validate stability of posterior covariance matrix, you set burn calc cov matrix increase burn...
int main(int argc, char *argv[])
Definition: ProcessMCMC.cpp:32
std::string config
Definition: ProcessMCMC.cpp:30
TH2D * TMatrixIntoTH2D(TMatrixDSym *Matrix, const std::string &title)
KS: Convert TMatrix to TH2D, mostly useful for making fancy plots.
std::vector< std::string > TitleNames
Definition: ProcessMCMC.cpp:29
void CalcBipolarPlot(MCMCProcessor *Processor)
void CalcParameterEvolution(MCMCProcessor *Processor)
void ProcessMCMC(const std::string &inputFile)
Main function processing MCMC and Producing plots.
int nFiles
Definition: ProcessMCMC.cpp:27
void MultipleProcessMCMC()
Function producing comparison of posterior and more betwen a few MCMC chains.
std::vector< std::string > FileNames
Definition: ProcessMCMC.cpp:28
void KolmogorovSmirnovTest(const std::vector< std::unique_ptr< MCMCProcessor >> &Processor, const std::unique_ptr< TCanvas > &Posterior, const TString &canvasname)
KS: Perform KS test to check if two posteriors for the same parameter came from the same distribution...
void CalcBayesFactor(MCMCProcessor *Processor)
KS: Calculate Bayes factor for a given hypothesis, most informative are those related to osc params....
std::map< std::string, std::pair< double, double > > GetCustomBinning(const YAML::Node &Settings)
Parse custom binning edges from a YAML configuration node.
Definition: ProcessMCMC.cpp:93
void CalcSavageDickey(MCMCProcessor *Processor)
Type GetFromManager(const YAML::Node &node, const Type defval, const std::string File="", const int Line=1)
Get content of config file if node is not found take default value specified.
Definition: YamlHelper.h:329
bool CheckNodeExists(const YAML::Node &node, Args... args)
KS: Wrapper function to call the recursive helper.
Definition: YamlHelper.h:60
#define M3OpenConfig(filename)
Macro to simplify calling LoadYaml with file and line info.
Definition: YamlHelper.h:589
Class responsible for processing MCMC chains, performing diagnostics, generating plots,...
Definition: MCMCProcessor.h:61
void GetNthParameter(const int param, double &Prior, double &PriorError, TString &Title) const
Get properties of parameter by passing it number.
void Reset2DPosteriors()
Reset 2D posteriors, in case we would like to calculate in again with different BurnInCut.
void GetBayesFactor(const std::vector< std::string > &ParName, const std::vector< std::vector< double >> &Model1Bounds, const std::vector< std::vector< double >> &Model2Bounds, const std::vector< std::vector< std::string >> &ModelNames)
Calculate Bayes factor for vector of params, and model boundaries.
void MakeTrianglePlot(const std::vector< std::string > &ParNames, const std::vector< double > &CredibleIntervals={0.99, 0.90, 0.68 }, const std::vector< Color_t > &CredibleIntervalsColours={kCyan+4, kCyan-2, kCyan-10}, const std::vector< double > &CredibleRegions={0.99, 0.90, 0.68}, const std::vector< Style_t > &CredibleRegionStyle={kDashed, kSolid, kDotted}, const std::vector< Color_t > &CredibleRegionColor={kGreen-3, kGreen-10, kGreen}, const bool CredibleInSigmas=false)
Make fancy triangle plot for selected parameters.
void SetPrintToPDF(const bool PlotOrNot)
Whether to dump all plots into PDF.
Long64_t GetnSteps()
Get Number of Steps that Chain has, for merged chains will not be the same nEntries.
void GetPolarPlot(const std::vector< std::string > &ParNames)
Make funny polar plot.
void ParameterEvolution(const std::vector< std::string > &Names, const std::vector< int > &NIntervals)
Make .gif of parameter evolution.
void GetCovariance(TMatrixDSym *&Cov, TMatrixDSym *&Corr)
Get the post-fit covariances and correlations.
void SetStepCut(const std::string &Cuts)
Set the step cutting by string.
void GetSavageDickey(const std::vector< std::string > &ParName, const std::vector< double > &EvaluationPoint, const std::vector< std::vector< double >> &Bounds)
Calculate Bayes factor for point like hypothesis using SavageDickey.
Custom exception class used throughout MaCh3.
std::unique_ptr< ObjectType > Clone(const ObjectType *obj, const std::string &name="")
KS: Creates a copy of a ROOT-like object and wraps it in a smart pointer.
constexpr static const double _BAD_DOUBLE_
Default value used for double initialisation.
Definition: Core.h:53
constexpr static const int _BAD_INT_
Default value used for int initialisation.
Definition: Core.h:55