MaCh3  2.4.2
Reference Guide
Functions | Variables
ProcessMCMC.cpp File Reference

Main exectable responsible for different types of MCMC processing like drawing posteriors, triangle plots etc. Actual implantation of methods is in MCMCProcessor. More...

#include "Fitters/OscProcessor.h"
#include "Manager/Manager.h"
Include dependency graph for ProcessMCMC.cpp:

Go to the source code of this file.

Functions

void ProcessMCMC (const std::string &inputFile)
 Main function processing MCMC and Producing plots. More...
 
void MultipleProcessMCMC ()
 Function producing comparison of posterior and more betwen a few MCMC chains. More...
 
void CalcBayesFactor (MCMCProcessor *Processor)
 KS: Calculate Bayes factor for a given hypothesis, most informative are those related to osc params. However, it make relative easy interpretation for switch dials. More...
 
void CalcSavageDickey (MCMCProcessor *Processor)
 
void CalcBipolarPlot (MCMCProcessor *Processor)
 
void CalcParameterEvolution (MCMCProcessor *Processor)
 
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 calc again and compare. By performing such operation several hundred times we can check when matrix becomes stable. More...
 
TH2D * TMatrixIntoTH2D (TMatrixDSym *Matrix, const std::string &title)
 KS: Convert TMatrix to TH2D, mostly useful for making fancy plots. More...
 
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. More...
 
int main (int argc, char *argv[])
 
std::map< std::string, std::pair< double, double > > GetCustomBinning (const YAML::Node &Settings)
 Parse custom binning edges from a YAML configuration node. More...
 

Variables

int nFiles
 
std::vector< std::string > FileNames
 
std::vector< std::string > TitleNames
 
std::string config
 

Detailed Description

Main exectable responsible for different types of MCMC processing like drawing posteriors, triangle plots etc. Actual implantation of methods is in MCMCProcessor.

Author
Kamil Skwarczynski

Definition in file ProcessMCMC.cpp.

Function Documentation

◆ CalcBayesFactor()

void CalcBayesFactor ( MCMCProcessor Processor)
inline

KS: Calculate Bayes factor for a given hypothesis, most informative are those related to osc params. However, it make relative easy interpretation for switch dials.

Definition at line 421 of file ProcessMCMC.cpp.

422 {
423  YAML::Node card_yaml = M3OpenConfig(config.c_str());
424  YAML::Node Settings = card_yaml["ProcessMCMC"];
425 
426  std::vector<std::string> ParNames;
427  std::vector<std::vector<double>> Model1Bounds;
428  std::vector<std::vector<double>> Model2Bounds;
429  std::vector<std::vector<std::string>> ModelNames;
430  for (const auto& dg : Settings["BayesFactor"])
431  {
432  ParNames.push_back(dg[0].as<std::string>());
433  ModelNames.push_back(dg[1].as<std::vector<std::string>>());
434  Model1Bounds.push_back(dg[2].as<std::vector<double>>());
435  Model2Bounds.push_back(dg[3].as<std::vector<double>>());
436  }
437 
438  Processor->GetBayesFactor(ParNames, Model1Bounds, Model2Bounds, ModelNames);
439 }
std::string config
Definition: ProcessMCMC.cpp:30
#define M3OpenConfig(filename)
Macro to simplify calling LoadYaml with file and line info.
Definition: YamlHelper.h:589
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.

◆ CalcBipolarPlot()

void CalcBipolarPlot ( MCMCProcessor Processor)
inline

Definition at line 474 of file ProcessMCMC.cpp.

475 {
476  YAML::Node card_yaml = M3OpenConfig(config.c_str());
477  YAML::Node Settings = card_yaml["ProcessMCMC"];
478 
479  std::vector<std::string> ParNames;
480  for (const auto& d : Settings["BipolarPlot"])
481  {
482  ParNames.push_back(d[0].as<std::string>());
483  }
484  Processor->GetPolarPlot(ParNames);
485 }
void GetPolarPlot(const std::vector< std::string > &ParNames)
Make funny polar plot.

◆ CalcParameterEvolution()

void CalcParameterEvolution ( MCMCProcessor Processor)
inline

Definition at line 459 of file ProcessMCMC.cpp.

460 {
461  YAML::Node card_yaml = M3OpenConfig(config.c_str());
462  YAML::Node Settings = card_yaml["ProcessMCMC"];
463 
464  std::vector<std::string> ParNames;
465  std::vector<int> Intervals;
466  for (const auto& d : Settings["ParameterEvolution"])
467  {
468  ParNames.push_back(d[0].as<std::string>());
469  Intervals.push_back(d[1].as<int>());
470  }
471  Processor->ParameterEvolution(ParNames, Intervals);
472 }
void ParameterEvolution(const std::vector< std::string > &Names, const std::vector< int > &NIntervals)
Make .gif of parameter evolution.

◆ CalcSavageDickey()

void CalcSavageDickey ( MCMCProcessor Processor)
inline

Definition at line 441 of file ProcessMCMC.cpp.

442 {
443  YAML::Node card_yaml = M3OpenConfig(config.c_str());
444  YAML::Node Settings = card_yaml["ProcessMCMC"];
445 
446  std::vector<std::string> ParNames;
447  std::vector<double> EvaluationPoint;
448  std::vector<std::vector<double>> Bounds;
449 
450  for (const auto& d : Settings["SavageDickey"])
451  {
452  ParNames.push_back(d[0].as<std::string>());
453  EvaluationPoint.push_back(d[1].as<double>());
454  Bounds.push_back(d[2].as<std::vector<double>>());
455  }
456  Processor->GetSavageDickey(ParNames, EvaluationPoint, Bounds);
457 }
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.

◆ DiagnoseCovarianceMatrix()

void DiagnoseCovarianceMatrix ( MCMCProcessor Processor,
const std::string &  inputFile 
)
inline

KS: You validate stability of posterior covariance matrix, you set burn calc cov matrix increase burn calc again and compare. By performing such operation several hundred times we can check when matrix becomes stable.

Definition at line 507 of file ProcessMCMC.cpp.

508 {
509  //Turn of plots from Processor
510  Processor->SetPrintToPDF(false);
511  // Open a TCanvas to write the posterior onto
512  auto Canvas = std::make_unique<TCanvas>("Canvas", "Canvas", 0, 0, 1024, 1024);
513  Canvas->SetGrid();
514  gStyle->SetOptStat(0);
515  gStyle->SetOptTitle(0);
516  Canvas->SetTickx();
517  Canvas->SetTicky();
518  Canvas->SetBottomMargin(0.1f);
519  Canvas->SetTopMargin(0.05f);
520  Canvas->SetRightMargin(0.15f);
521  Canvas->SetLeftMargin(0.10f);
522 
523  //KS: Fancy colours
524  const int NRGBs = 10;
525  TColor::InitializeColors();
526  Double_t stops[NRGBs] = { 0.00, 0.10, 0.25, 0.35, 0.50, 0.60, 0.65, 0.75, 0.90, 1.00 };
527  Double_t red[NRGBs] = { 0.50, 1.00, 1.00, 0.25, 0.00, 0.10, 0.50, 1.00, 0.75, 0.55 };
528  Double_t green[NRGBs] = { 0.00, 0.25, 1.00, 0.25, 0.00, 0.60, 0.90, 1.00, 0.75, 0.75 };
529  Double_t blue[NRGBs] = { 0.00, 0.25, 1.00, 1.00, 0.50, 0.60, 0.90, 1.00, 0.05, 0.05 };
530  TColor::CreateGradientColorTable(NRGBs, stops, red, green, blue, 255);
531  gStyle->SetNumberContours(255);
532 
533  std::string OutName = inputFile;
534  OutName = OutName.substr(0, OutName.find(".root"));
535  Canvas->Print(Form("Correlation_%s.pdf[", OutName.c_str()), "pdf");
536  Canvas->Print(Form("Covariance_%s.pdf[", OutName.c_str()), "pdf");
537 
538  YAML::Node card_yaml = M3OpenConfig(config.c_str());
539  YAML::Node Settings = card_yaml["ProcessMCMC"];
540 
541  const int entries = int(Processor->GetnSteps());
542  const int NIntervals = GetFromManager<int>(Settings["NIntervals"], 5);
543  const int IntervalsSize = entries/NIntervals;
544  //We start with burn from 0 (no burn in at all)
545  int BurnIn = 0;
546  MACH3LOG_INFO("Diagnosing matrices with entries={}, NIntervals={} and IntervalsSize={}", entries, NIntervals, IntervalsSize);
547 
548  TMatrixDSym *Covariance = nullptr;
549  TMatrixDSym *Correlation = nullptr;
550 
551  TH2D *CovariancePreviousHist = nullptr;
552  TH2D *CorrelationPreviousHist = nullptr;
553 
554  TH2D *CovarianceHist = nullptr;
555  TH2D *CorrelationHist = nullptr;
556 
557  //KS: Get first covariances, we need two for comparison...
558  Processor->SetStepCut(BurnIn);
559  Processor->GetCovariance(Covariance, Correlation);
560 
561  CovariancePreviousHist = TMatrixIntoTH2D(Covariance, "Covariance");
562  CorrelationPreviousHist = TMatrixIntoTH2D(Correlation, "Correlation");
563 
564  delete Covariance;
565  Covariance = nullptr;
566  delete Correlation;
567  Correlation = nullptr;
568 
569  //KS: Loop over all desired cuts
570  for(int k = 1; k < NIntervals; ++k)
571  {
572  BurnIn = k*IntervalsSize;
573  Processor->SetStepCut(BurnIn);
574  Processor->GetCovariance(Covariance, Correlation);
575  Processor->Reset2DPosteriors();
576 
577  CovarianceHist = TMatrixIntoTH2D(Covariance, "Covariance");
578  CorrelationHist = TMatrixIntoTH2D(Correlation, "Correlation");
579 
580  TH2D *CovarianceDiff = static_cast<TH2D*>(CovarianceHist->Clone("Covariance_Ratio"));
581  TH2D *CorrelationDiff = static_cast<TH2D*>(CorrelationHist->Clone("Correlation_Ratio"));
582 
583  //KS: Bit messy but quite often covariance is 0 is divided by 0 is problematic so
584  #ifdef MULTITHREAD
585  #pragma omp parallel for
586  #endif
587  for (int j = 1; j < CovarianceDiff->GetXaxis()->GetNbins()+1; ++j)
588  {
589  for (int i = 1; i < CovarianceDiff->GetYaxis()->GetNbins()+1; ++i)
590  {
591  if( std::fabs (CovarianceDiff->GetBinContent(j, i)) < 1.e-5 && std::fabs (CovariancePreviousHist->GetBinContent(j, i)) < 1.e-5)
592  {
593  CovarianceDiff->SetBinContent(j, i, M3::_BAD_DOUBLE_);
594  CovariancePreviousHist->SetBinContent(j, i, M3::_BAD_DOUBLE_);
595  }
596  if( std::fabs (CorrelationDiff->GetBinContent(j, i)) < 1.e-5 && std::fabs (CorrelationPreviousHist->GetBinContent(j, i)) < 1.e-5)
597  {
598  CorrelationDiff->SetBinContent(j, i, M3::_BAD_DOUBLE_);
599  CorrelationPreviousHist->SetBinContent(j, i, M3::_BAD_DOUBLE_);
600  }
601  }
602  }
603  //Divide matrices
604  CovarianceDiff->Divide(CovariancePreviousHist);
605  CorrelationDiff->Divide(CorrelationPreviousHist);
606 
607  //Now it is time for fancy names etc.
608  for (int j = 0; j < CovarianceDiff->GetXaxis()->GetNbins(); ++j)
609  {
610  TString Title = "";
611  double Prior = 1.0;
612  double PriorError = 1.0;
613 
614  Processor->GetNthParameter(j, Prior, PriorError, Title);
615 
616  CovarianceDiff->GetXaxis()->SetBinLabel(j+1, Title);
617  CovarianceDiff->GetYaxis()->SetBinLabel(j+1, Title);
618  CorrelationDiff->GetXaxis()->SetBinLabel(j+1, Title);
619  CorrelationDiff->GetYaxis()->SetBinLabel(j+1, Title);
620  }
621  CovarianceDiff->GetXaxis()->SetLabelSize(0.015f);
622  CovarianceDiff->GetYaxis()->SetLabelSize(0.015f);
623  CorrelationDiff->GetXaxis()->SetLabelSize(0.015f);
624  CorrelationDiff->GetYaxis()->SetLabelSize(0.015f);
625 
626  std::stringstream ss;
627  ss << "BCut_";
628  ss << BurnIn;
629  ss << "/";
630  ss << "BCut_";
631  ss << (k-1)*IntervalsSize;
632  std::string str = ss.str();
633 
634  TString Title = "Cov " + str;
635  CovarianceDiff->GetZaxis()->SetTitle( Title );
636  Title = "Corr " + str;
637  CorrelationDiff->GetZaxis()->SetTitle(Title);
638 
639  CovarianceDiff->SetMinimum(-2);
640  CovarianceDiff->SetMaximum(2);
641  CorrelationDiff->SetMinimum(-2);
642  CorrelationDiff->SetMaximum(2);
643 
644  Canvas->cd();
645  CovarianceDiff->Draw("colz");
646  Canvas->Print(Form("Covariance_%s.pdf", OutName.c_str()), "pdf");
647 
648  CorrelationDiff->Draw("colz");
649  Canvas->Print(Form("Correlation_%s.pdf", OutName.c_str()), "pdf");
650 
651  //KS: Current hist become previous as we need it for further comparison
652  delete CovariancePreviousHist;
653  CovariancePreviousHist = static_cast<TH2D*>(CovarianceHist->Clone());
654  delete CorrelationPreviousHist;
655  CorrelationPreviousHist = static_cast<TH2D*>(CorrelationHist->Clone());
656 
657  delete CovarianceHist;
658  CovarianceHist = nullptr;
659  delete CorrelationHist;
660  CorrelationHist = nullptr;
661 
662  delete CovarianceDiff;
663  delete CorrelationDiff;
664  delete Covariance;
665  Covariance = nullptr;
666  delete Correlation;
667  Correlation = nullptr;
668  }
669  Canvas->cd();
670  Canvas->Print(Form("Covariance_%s.pdf]", OutName.c_str()), "pdf");
671  Canvas->Print(Form("Correlation_%s.pdf]", OutName.c_str()), "pdf");
672 
673  Processor->SetPrintToPDF(true);
674  if(Covariance != nullptr) delete Covariance;
675  if(Correlation != nullptr) delete Correlation;
676  if(CovariancePreviousHist != nullptr) delete CovariancePreviousHist;
677  if(CorrelationPreviousHist != nullptr) delete CorrelationPreviousHist;
678  if(CovarianceHist != nullptr) delete CovarianceHist;
679  if(CorrelationHist != nullptr) delete CorrelationHist;
680 }
#define MACH3LOG_INFO
Definition: MaCh3Logger.h:35
TH2D * TMatrixIntoTH2D(TMatrixDSym *Matrix, const std::string &title)
KS: Convert TMatrix to TH2D, mostly useful for making fancy plots.
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 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 GetCovariance(TMatrixDSym *&Cov, TMatrixDSym *&Corr)
Get the post-fit covariances and correlations.
void SetStepCut(const std::string &Cuts)
Set the step cutting by string.
constexpr static const double _BAD_DOUBLE_
Default value used for double initialisation.
Definition: Core.h:53

◆ GetCustomBinning()

std::map<std::string, std::pair<double, double> > GetCustomBinning ( const YAML::Node &  Settings)

Parse custom binning edges from a YAML configuration node.

This function reads the CustomBinEdges section of the YAML config and returns a mapping of parameter names to their lower and upper edges.

The expected YAML syntax is:

CustomBinEdges:
delta_cp: [-3.141592, 3.141592]
another_param: [min, max]
Parameters
SettingsYAML configuration node containing the optional CustomBinEdges section.

Definition at line 93 of file ProcessMCMC.cpp.

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 }
#define MACH3LOG_DEBUG
Definition: MaCh3Logger.h:34
#define MACH3LOG_ERROR
Definition: MaCh3Logger.h:37
Custom exception class used throughout MaCh3.

◆ GetTrianglePlot()

void GetTrianglePlot ( MCMCProcessor Processor)
inline

Definition at line 487 of file ProcessMCMC.cpp.

487  {
488  YAML::Node card_yaml = M3OpenConfig(config.c_str());
489  YAML::Node Settings = card_yaml["ProcessMCMC"];
490 
491  for (const auto& dg : Settings["TrianglePlot"])
492  {
493  std::string ParName = dg[0].as<std::string>();
494 
495  std::vector<std::string> NameVec = dg[1].as<std::vector<std::string>>();
496  Processor->MakeTrianglePlot(NameVec,
497  GetFromManager<std::vector<double>>(Settings["CredibleIntervals"], {0.99, 0.90, 0.68}),
498  GetFromManager<std::vector<short int>>(Settings["CredibleIntervalsColours"], {436, 430, 422}),
499  GetFromManager<std::vector<double>>(Settings["CredibleRegions"], {0.99, 0.90, 0.68}),
500  GetFromManager<std::vector<short int>>(Settings["CredibleRegionStyle"], {2, 1, 3}),
501  GetFromManager<std::vector<short int>>(Settings["CredibleRegionColor"], {413, 406, 416}),
502  GetFromManager<bool>(Settings["CredibleInSigmas"], false));
503  }
504 }
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
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.

◆ KolmogorovSmirnovTest()

void KolmogorovSmirnovTest ( const std::vector< std::unique_ptr< MCMCProcessor >> &  Processor,
const std::unique_ptr< TCanvas > &  Posterior,
const TString &  canvasname 
)
inline

KS: Perform KS test to check if two posteriors for the same parameter came from the same distribution.

Definition at line 698 of file ProcessMCMC.cpp.

701 {
702  constexpr Color_t CumulativeColor[] = {kBlue-1, kRed, kGreen+2};
703  constexpr Style_t CumulativeStyle[] = {kSolid, kDashed, kDotted};
704 
705  for(int i = 0; i < Processor[0]->GetNParams(); ++i)
706  {
707  // This holds the posterior density
708  std::vector<std::unique_ptr<TH1D>> hpost(nFiles);
709  std::vector<std::unique_ptr<TH1D>> CumulativeDistribution(nFiles);
710 
711  TString Title;
712  double Prior = 1.0;
713  double PriorError = 1.0;
714 
715  Processor[0]->GetNthParameter(i, Prior, PriorError, Title);
716  bool Skip = false;
717  for (int ik = 0 ; ik < nFiles; ik++)
718  {
719  int Index = 0;
720  if(ik == 0 ) Index = i;
721  else
722  {
723  // KS: If somehow this chain doesn't given params we skip it
724  Index = Processor[ik]->GetParamIndexFromName(hpost[0]->GetTitle());
725  if(Index == M3::_BAD_INT_)
726  {
727  Skip = true;
728  break;
729  }
730  }
731  hpost[ik] = M3::Clone(Processor[ik]->GetHpost(Index));
732  CumulativeDistribution[ik] = M3::Clone(Processor[ik]->GetHpost(Index));
733  CumulativeDistribution[ik]->Fill(0., 0.);
734  CumulativeDistribution[ik]->Reset();
735  CumulativeDistribution[ik]->SetMaximum(1.);
736  TString TempTitle = Title+" Kolmogorov Smirnov";
737  CumulativeDistribution[ik]->SetTitle(TempTitle);
738 
739  TempTitle = Title+" Value";
740  CumulativeDistribution[ik]->GetXaxis()->SetTitle(TempTitle);
741  CumulativeDistribution[ik]->GetYaxis()->SetTitle("Cumulative Probability");
742 
743  CumulativeDistribution[ik]->SetLineWidth(2);
744  CumulativeDistribution[ik]->SetLineColor(CumulativeColor[ik]);
745  CumulativeDistribution[ik]->SetLineStyle(CumulativeStyle[ik]);
746  }
747 
748  // Don't plot if this is a fixed histogram (i.e. the peak is the whole integral)
749  if(hpost[0]->GetMaximum() == hpost[0]->Integral()*1.5 || Skip) {
750  continue;
751  }
752 
753  for (int ik = 0 ; ik < nFiles; ik++)
754  {
755  const int NumberOfBins = hpost[ik]->GetXaxis()->GetNbins();
756  double Cumulative = 0;
757  const double Integral = hpost[ik]->Integral();
758  for (int j = 1; j < NumberOfBins+1; ++j)
759  {
760  Cumulative += hpost[ik]->GetBinContent(j)/Integral;
761  CumulativeDistribution[ik]->SetBinContent(j, Cumulative);
762  }
763  //KS: Set overflow to 1 just in case
764  CumulativeDistribution[ik]->SetBinContent(NumberOfBins+1, 1.);
765  }
766 
767  std::vector<int> TestStatBin(nFiles, 0);
768  std::vector<double> TestStatD(nFiles, -999);
769  std::vector<std::unique_ptr<TLine>> LineD(nFiles);
770  //Find KS statistic
771  for (int ik = 1 ; ik < nFiles; ik++)
772  {
773  const int NumberOfBins = CumulativeDistribution[0]->GetXaxis()->GetNbins();
774  for (int j = 1; j < NumberOfBins+1; ++j)
775  {
776  const double BinValue = CumulativeDistribution[0]->GetBinCenter(j);
777  const int BinNumber = CumulativeDistribution[ik]->FindBin(BinValue);
778  //KS: Calculate D statistic for this bin, only save it if it's bigger than previously found value
779  double TempDstat = std::fabs(CumulativeDistribution[0]->GetBinContent(j) - CumulativeDistribution[ik]->GetBinContent(BinNumber));
780  if(TempDstat > TestStatD[ik])
781  {
782  TestStatD[ik] = TempDstat;
783  TestStatBin[ik] = j;
784  }
785  }
786  }
787 
788  for (int ik = 0 ; ik < nFiles; ik++)
789  {
790  LineD[ik] = std::make_unique<TLine>(CumulativeDistribution[0]->GetBinCenter(TestStatBin[ik]), 0, CumulativeDistribution[0]->GetBinCenter(TestStatBin[ik]), CumulativeDistribution[0]->GetBinContent(TestStatBin[ik]));
791  LineD[ik]->SetLineColor(CumulativeColor[ik]);
792  LineD[ik]->SetLineWidth(2.0);
793  }
794  CumulativeDistribution[0]->Draw();
795  for (int ik = 0 ; ik < nFiles; ik++)
796  CumulativeDistribution[ik]->Draw("SAME");
797 
798  auto leg = std::make_unique<TLegend>(0.15, 0.7, 0.5, 0.90);
799  leg->SetTextSize(0.04f);
800  for (int ik = 0; ik < nFiles; ik++)
801  leg->AddEntry(CumulativeDistribution[ik].get(), TitleNames[ik].c_str(), "l");
802  for (int ik = 1; ik < nFiles; ik++)
803  leg->AddEntry(LineD[ik].get(), Form("#Delta D = %.4f", TestStatD[ik]), "l");
804 
805  leg->SetLineColor(0);
806  leg->SetLineStyle(0);
807  leg->SetFillColor(0);
808  leg->SetFillStyle(0);
809  leg->Draw("SAME");
810 
811  for (int ik = 1; ik < nFiles; ik++)
812  LineD[ik]->Draw("sam");
813 
814  Posterior->cd();
815  Posterior->Print(canvasname);
816  } //End loop over parameter
817 }
std::vector< std::string > TitleNames
Definition: ProcessMCMC.cpp:29
int nFiles
Definition: ProcessMCMC.cpp:27
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 int _BAD_INT_
Default value used for int initialisation.
Definition: Core.h:55

◆ main()

int main ( int  argc,
char *  argv[] 
)

Definition at line 32 of file ProcessMCMC.cpp.

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 }
void SetMaCh3LoggerFormat()
Set messaging format of the logger.
Definition: MaCh3Logger.h:61
void ProcessMCMC(const std::string &inputFile)
Main function processing MCMC and Producing plots.
void MultipleProcessMCMC()
Function producing comparison of posterior and more betwen a few MCMC chains.
std::vector< std::string > FileNames
Definition: ProcessMCMC.cpp:28
bool CheckNodeExists(const YAML::Node &node, Args... args)
KS: Wrapper function to call the recursive helper.
Definition: YamlHelper.h:60

◆ MultipleProcessMCMC()

void MultipleProcessMCMC ( )
inline

Function producing comparison of posterior and more betwen a few MCMC chains.

Definition at line 212 of file ProcessMCMC.cpp.

213 {
214  YAML::Node card_yaml = M3OpenConfig(config.c_str());
215  YAML::Node Settings = card_yaml["ProcessMCMC"];
216 
217  constexpr Color_t PosteriorColor[] = {kBlue-1, kRed, kGreen+2};
218  //constexpr Style_t PosteriorStyle[] = {kSolid, kDashed, kDotted};
219  nFiles = int(FileNames.size());
220  std::vector<std::unique_ptr<MCMCProcessor>> Processor(nFiles);
221 
222  if(!Settings["BurnInSteps"]) {
223  MACH3LOG_WARN("BurnInSteps not set, defaulting to 20%");
224  }
225 
226  for (int ik = 0; ik < nFiles; ik++)
227  {
228  MACH3LOG_INFO("File for study: {}", FileNames[ik]);
229  // Make the processor
230  Processor[ik] = std::make_unique<MCMCProcessor>(FileNames[ik]);
231  Processor[ik]->SetOutputSuffix(("_" + std::to_string(ik)).c_str());
232 
233  Processor[ik]->SetExcludedTypes(GetFromManager<std::vector<std::string>>(Settings["ExcludedTypes"], {}));
234  Processor[ik]->SetExcludedNames(GetFromManager<std::vector<std::string>>(Settings["ExcludedNames"], {}));
235  Processor[ik]->SetExcludedGroups(GetFromManager<std::vector<std::string>>(Settings["ExcludedGroups"], {}));
236 
237  //Apply additional cuts to 1D posterior
238  Processor[ik]->SetPosterior1DCut(GetFromManager<std::string>(Settings["Posterior1DCut"], ""));
239 
240  Processor[ik]->SetPlotRelativeToPrior(GetFromManager<bool>(Settings["PlotRelativeToPrior"], false));
241  Processor[ik]->SetFancyNames(GetFromManager<bool>(Settings["FancyNames"], true));
242  Processor[ik]->Initialise();
243 
244  if(Settings["BurnInSteps"]) {
245  Processor[ik]->SetStepCut(Settings["BurnInSteps"].as<int>());
246  }else {
247  Processor[ik]->SetStepCut(static_cast<int>(Processor[ik]->GetnSteps()/5));
248  }
249 
250  if(Settings["MaxEntries"]) {
251  Processor[ik]->SetEntries(Get<int>(Settings["MaxEntries"], __FILE__, __LINE__));
252  }
253  if(Settings["NBins"]) {
254  Processor[ik]->SetNBins(Get<int>(Settings["NBins"], __FILE__, __LINE__));
255  }
256  }
257 
258  Processor[0]->MakePostfit(GetCustomBinning(Settings));
259  Processor[0]->DrawPostfit();
260  // Get edges from first histogram to ensure all params use same binning
261  std::map<std::string, std::pair<double, double>> ParamEdges;
262  for(int i = 0; i < Processor[0]->GetNParams(); ++i) {
263  // Get the histogram for the i-th parameter
264  TH1D* hist = Processor[0]->GetHpost(i);
265  if (!hist) {
266  MACH3LOG_DEBUG("Histogram for parameter {} is null.", i);
267  continue;
268  }
269 
270  // Get the parameter name (title of the histogram)
271  std::string paramName = hist->GetTitle();
272 
273  // Get the axis limits (edges)
274  TAxis* axis = hist->GetXaxis();
275  double xmin = axis->GetXmin();
276  double xmax = axis->GetXmax();
277 
278  MACH3LOG_DEBUG("Adding bin edges for {} equal to {:.4f}, {:.4f}",paramName, xmin, xmax);
279  // Insert into the map
280  ParamEdges[paramName] = std::make_pair(xmin, xmax);
281  }
282 
283  //KS: Multithreading here is very tempting but there are some issues with root that need to be resovled :(
284  for (int ik = 1; ik < nFiles; ik++)
285  {
286  // Make the postfit
287  Processor[ik]->MakePostfit(ParamEdges);
288  Processor[ik]->DrawPostfit();
289  }
290 
291  // Open a TCanvas to write the posterior onto
292  auto Posterior = std::make_unique<TCanvas>("PosteriorMulti", "PosteriorMulti", 0, 0, 1024, 1024);
293  gStyle->SetOptStat(0);
294  gStyle->SetOptTitle(0);
295  Posterior->SetGrid();
296  Posterior->SetBottomMargin(0.1f);
297  Posterior->SetTopMargin(0.05f);
298  Posterior->SetRightMargin(0.03f);
299  Posterior->SetLeftMargin(0.15f);
300 
301  // First filename: keep path, just remove ".root"
302  // Would be nice to specify outpath in a later update
303  size_t pos = FileNames[0].rfind(".root");
304  std::string base = (pos == std::string::npos) ? FileNames[0] : FileNames[0].substr(0, pos);
305  TString canvasname = base;
306 
307  // Remaining filenames: strip path and ".root"
308  // 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
309  for (int ik = 1; ik < nFiles; ik++) {
310  pos = FileNames[ik].find_last_of('/');
311  base = (pos == std::string::npos) ? FileNames[ik] : FileNames[ik].substr(pos + 1);
312  pos = base.rfind(".root");
313  if (pos != std::string::npos) base = base.substr(0, pos);
314  canvasname += "_" + TString(base);
315  }
316 
317  canvasname = canvasname +".pdf[";
318 
319  Posterior->Print(canvasname);
320  // Once the pdf file is open no longer need to bracket
321  canvasname.ReplaceAll("[","");
322 
323  for(int i = 0; i < Processor[0]->GetNParams(); ++i)
324  {
325  // This holds the posterior density
326  std::vector<std::unique_ptr<TH1D>> hpost(nFiles);
327  std::vector<std::unique_ptr<TLine>> hpd(nFiles);
328  hpost[0] = M3::Clone(Processor[0]->GetHpost(i));
329  hpost[0]->GetYaxis()->SetTitle("Posterior Density");
330  bool Skip = false;
331  for (int ik = 1 ; ik < nFiles; ik++)
332  {
333  // KS: If somehow this chain doesn't given params we skip it
334  const int Index = Processor[ik]->GetParamIndexFromName(hpost[0]->GetTitle());
335  if(Index == M3::_BAD_INT_)
336  {
337  Skip = true;
338  break;
339  }
340  hpost[ik] = M3::Clone(Processor[ik]->GetHpost(Index));
341  }
342 
343  // Don't plot if this is a fixed histogram (i.e. the peak is the whole integral)
344  if(hpost[0]->GetMaximum() == hpost[0]->Integral()*1.5 || Skip) {
345  continue;
346  }
347  for (int ik = 0; ik < nFiles; ik++)
348  {
349  RemoveFitter(hpost[ik].get(), "Gauss");
350 
351  // Set some nice colours
352  hpost[ik]->SetLineColor(PosteriorColor[ik]);
353  //hpost[ik]->SetLineStyle(PosteriorStyle[ik]);
354  hpost[ik]->SetLineWidth(2);
355 
356  // Area normalise the distributions
357  hpost[ik]->Scale(1./hpost[ik]->Integral());
358  }
359  TString Title;
360  double Prior = 1.0;
361  double PriorError = 1.0;
362 
363  Processor[0]->GetNthParameter(i, Prior, PriorError, Title);
364 
365  // Now make the TLine for the Asimov
366  auto Asimov = std::make_unique<TLine>(Prior, hpost[0]->GetMinimum(), Prior, hpost[0]->GetMaximum());
367  Asimov->SetLineColor(kRed-3);
368  Asimov->SetLineWidth(2);
369  Asimov->SetLineStyle(kDashed);
370 
371  // Make a nice little TLegend
372  auto leg = std::make_unique<TLegend>(0.20, 0.7, 0.6, 0.97);
373  leg->SetTextSize(0.03f);
374  leg->SetFillColor(0);
375  leg->SetFillStyle(0);
376  leg->SetLineColor(0);
377  leg->SetLineStyle(0);
378  TString asimovLeg = Form("#splitline{Prior}{x = %.2f , #sigma = %.2f}", Prior, PriorError);
379  leg->AddEntry(Asimov.get(), asimovLeg, "l");
380 
381  for (int ik = 0; ik < nFiles; ik++)
382  {
383  TString rebinLeg = Form("#splitline{%s}{#mu = %.2f, #sigma = %.2f}", TitleNames[ik].c_str(), hpost[ik]->GetMean(), hpost[ik]->GetRMS());
384  leg->AddEntry(hpost[ik].get(), rebinLeg, "l");
385 
386  hpd[ik] = std::make_unique<TLine>(hpost[ik]->GetBinCenter(hpost[ik]->GetMaximumBin()), hpost[ik]->GetMinimum(),
387  hpost[ik]->GetBinCenter(hpost[ik]->GetMaximumBin()), hpost[ik]->GetMaximum());
388  hpd[ik]->SetLineColor(hpost[ik]->GetLineColor());
389  hpd[ik]->SetLineWidth(2);
390  hpd[ik]->SetLineStyle(kSolid);
391  }
392 
393  // Find the maximum value to nicely resize hist
394  double maximum = 0;
395  for (int ik = 0; ik < nFiles; ik++) maximum = std::max(maximum, hpost[ik]->GetMaximum());
396  for (int ik = 0; ik < nFiles; ik++) hpost[ik]->SetMaximum(1.3*maximum);
397 
398  hpost[0]->Draw("hist");
399  for (int ik = 1; ik < nFiles; ik++) hpost[ik]->Draw("hist same");
400  Asimov->Draw("same");
401  for (int ik = 0; ik < nFiles; ik++) hpd[ik]->Draw("same");
402  leg->Draw("same");
403  Posterior->cd();
404  Posterior->Print(canvasname);
405  }//End loop over parameters
406 
407  // Finally draw the parameter plot onto the PDF
408  // Close the .pdf file with all the posteriors
409  Posterior->cd();
410  Posterior->Clear();
411 
412  if(GetFromManager<bool>(Settings["PerformKStest"], true)) KolmogorovSmirnovTest(Processor, Posterior, canvasname);
413 
414  // Close the pdf file
415  MACH3LOG_INFO("Closing pdf {}", canvasname);
416  canvasname+="]";
417  Posterior->Print(canvasname);
418 }
void RemoveFitter(TH1D *hist, const std::string &name)
KS: Remove fitted TF1 from hist to make comparison easier.
#define MACH3LOG_WARN
Definition: MaCh3Logger.h:36
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...
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

◆ ProcessMCMC()

void ProcessMCMC ( const std::string &  inputFile)
inline

Main function processing MCMC and Producing plots.

Definition at line 115 of file ProcessMCMC.cpp.

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  if(GetFromManager<bool>(Settings["JarlskogAnalysis"], true)) Processor->PerformJarlskogAnalysis();
209  if(GetFromManager<bool>(Settings["MakePiePlot"], true)) Processor->MakePiePlot();
210 }
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...
void CalcBipolarPlot(MCMCProcessor *Processor)
void CalcParameterEvolution(MCMCProcessor *Processor)
void CalcBayesFactor(MCMCProcessor *Processor)
KS: Calculate Bayes factor for a given hypothesis, most informative are those related to osc params....
void CalcSavageDickey(MCMCProcessor *Processor)

◆ TMatrixIntoTH2D()

TH2D * TMatrixIntoTH2D ( TMatrixDSym *  Matrix,
const std::string &  title 
)
inline

KS: Convert TMatrix to TH2D, mostly useful for making fancy plots.

Definition at line 683 of file ProcessMCMC.cpp.

684 {
685  TH2D* hMatrix = new TH2D(title.c_str(), title.c_str(), Matrix->GetNrows(), 0.0, Matrix->GetNrows(), Matrix->GetNcols(), 0.0, Matrix->GetNcols());
686  for(int i = 0; i < Matrix->GetNrows(); i++)
687  {
688  for(int j = 0; j < Matrix->GetNcols(); j++)
689  {
690  //KS: +1 because there is offset in histogram relative to TMatrix
691  hMatrix->SetBinContent(i+1,j+1, (*Matrix)(i,j));
692  }
693  }
694  return hMatrix;
695 }

Variable Documentation

◆ config

std::string config

Definition at line 30 of file ProcessMCMC.cpp.

◆ FileNames

std::vector<std::string> FileNames

Definition at line 28 of file ProcessMCMC.cpp.

◆ nFiles

int nFiles

Definition at line 27 of file ProcessMCMC.cpp.

◆ TitleNames

std::vector<std::string> TitleNames

Definition at line 29 of file ProcessMCMC.cpp.