MaCh3  2.5.1
Reference Guide
Functions | Variables
GetPostfitParamPlots.cpp File Reference
#include <iostream>
#include <sstream>
#include <iomanip>
#include <algorithm>
#include "PlottingUtils/PlottingUtils.h"
#include "PlottingUtils/PlottingManager.h"
#include "TROOT.h"
#include "TGaxis.h"
#include "TString.h"
#include "TStyle.h"
#include "TH1.h"
#include "TH2.h"
#include "TF1.h"
#include "TLegend.h"
#include "TPad.h"
#include "TCanvas.h"
#include "TTree.h"
#include "TFile.h"
#include "TVectorD.h"
#include "TCandle.h"
#include "TFrame.h"
#include "TGraphAsymmErrors.h"
Include dependency graph for GetPostfitParamPlots.cpp:

Go to the source code of this file.

Functions

void copyParToBlockHist (const int localBin, const std::string &paramName, TH1D *blockHist, const std::string &type, const int fileId, const bool setLabels=true)
 
void InitializePads (const TCanvas *canvas, TPad *&pad3, TPad *&pad4)
 
void CopyViolinToBlock (TH2D *FullViolin, TH2D *ReducedViolin, const std::vector< std::string > &ParamNames)
 
template<typename HistType >
void PrettifyTitles (HistType *hist)
 
bool ReadSettings (const std::shared_ptr< TFile > &File1)
 
std::unique_ptr< TH1D > makeRatio (TH1D *PrefitCopy, TH1D *PostfitCopy, bool setAxes)
 
bool IsInvalidHist (const TH1D *hist, double invalid=M3::_BAD_DOUBLE_)
 
void DrawPlots (TCanvas *plotCanv, TH1D *PrefitCopy, const std::vector< std::unique_ptr< TH1D >> &PostfitVec, TPad *mainPad, TPad *ratioPad, const std::string &OutName)
 
void MakeParameterPlots (TCanvas *canv)
 
void MakeFluxPlots (TCanvas *canv)
 
void MakeNDDetPlots ()
 
void MakeRidgePlots ()
 
void GetPostfitParamPlots ()
 
std::unique_ptr< TGraphAsymmErrors > MakeTGraphAsymmErrors (const std::shared_ptr< TFile > &File, std::vector< int > Index={})
 
void GetViolinPlots ()
 KS: Make fancy violin plots. More...
 
void Get2DComparison (const std::string &FileName1, const std::string &FileName2)
 KS: Make comparison of 2D Posteriors. More...
 
int main (int argc, char *argv[])
 

Variables

M3::Plotting::PlottingManager * PlotMan = nullptr
 
int NDParameters
 
int NDParametersStartingPos
 
std::vector< int > NDSamplesBins
 
std::vector< std::string > NDSamplesNames
 
std::string SaveName
 
std::string plotType
 

Detailed Description

This script generates post-fit parameter plots. The central postfit value is taken as the Highest Posterior Density (HPD), but can be easily changed to another method such as Gaussian. Be cautious as parameter names and the number of parameters per plot are currently hardcoded.

Usage:

./GetPostfitParamPlots ProcessMCMC_Output1.root <ProcessMCMC_Output2.root> <ProcessMCMC_Output3.root>
void GetPostfitParamPlots()
Author
Clarence Wret
Will Parker
Kamil Skwarczynski
Ewan Miller

Definition in file GetPostfitParamPlots.cpp.

Function Documentation

◆ copyParToBlockHist()

void copyParToBlockHist ( const int  localBin,
const std::string &  paramName,
TH1D *  blockHist,
const std::string &  type,
const int  fileId,
const bool  setLabels = true 
)

Definition at line 62 of file GetPostfitParamPlots.cpp.

63  {
64  // Set the values in the sub-histograms
65  MACH3LOG_DEBUG("copying data from at local bin {}: for parameter {}", localBin, paramName);
66  MACH3LOG_DEBUG(" Fitter specific name: {}", PlotMan->input().translateName(fileId, M3::Plotting::kPostFit, paramName));
67  MACH3LOG_DEBUG(" value: {:.4f}", PlotMan->input().getPostFitValue(fileId, paramName, type));
68  MACH3LOG_DEBUG(" error: {:.4f}", PlotMan->input().getPostFitError(fileId, paramName, type));
69 
70  blockHist->SetBinContent(localBin +1, PlotMan->input().getPostFitValue(fileId, paramName, type));
71  blockHist->SetBinError(localBin +1, PlotMan->input().getPostFitError(fileId, paramName, type));
72 
73  if(setLabels){
74  blockHist->GetXaxis()->SetBinLabel(localBin +1, paramName.c_str());
75  blockHist->GetXaxis()->LabelsOption("v");
76  }
77 }
M3::Plotting::PlottingManager * PlotMan
#define MACH3LOG_DEBUG
Definition: MaCh3Logger.h:34

◆ CopyViolinToBlock()

void CopyViolinToBlock ( TH2D *  FullViolin,
TH2D *  ReducedViolin,
const std::vector< std::string > &  ParamNames 
)

Definition at line 97 of file GetPostfitParamPlots.cpp.

97  {
98  for(unsigned int i = 0; i < ParamNames.size(); i++)
99  {
100  int ParamBinId = M3::_BAD_INT_;
101  for (int ix = 0; ix < FullViolin->GetXaxis()->GetNbins(); ++ix) {
102  if(FullViolin->GetXaxis()->GetBinLabel(ix+1) == ParamNames[i])
103  {
104  ParamBinId = ix+1;
105  break;
106  }
107  }
108  if(ParamBinId == M3::_BAD_INT_) {
109  MACH3LOG_WARN("Didn't find param {}", ParamNames[i]);
110  continue;
111  }
112  //KS Fill content of reduced violin
113  for (int iy = 0; iy < FullViolin->GetYaxis()->GetNbins(); ++iy) {
114  ReducedViolin->SetBinContent(i+1, iy+1, FullViolin->GetBinContent(ParamBinId, iy+1));
115  ReducedViolin->GetXaxis()->SetBinLabel(i+1, ParamNames[i].c_str());
116  }
117  }
118  ReducedViolin->SetFillColor(FullViolin->GetFillColor());
119  ReducedViolin->SetFillColorAlpha(FullViolin->GetMarkerColor(), 0.35);
120  ReducedViolin->SetLineColor(FullViolin->GetMarkerColor());
121 
122  ReducedViolin->SetMarkerColor(FullViolin->GetMarkerColor());
123  ReducedViolin->SetMarkerStyle(FullViolin->GetMarkerStyle());
124  ReducedViolin->SetMarkerSize(FullViolin->GetMarkerSize());
125 
126  ReducedViolin->GetYaxis()->SetTitleOffset(FullViolin->GetTitleOffset());
127  ReducedViolin->GetYaxis()->SetTitle(FullViolin->GetYaxis()->GetTitle());
128  ReducedViolin->GetXaxis()->LabelsOption("v");
129 }
#define MACH3LOG_WARN
Definition: MaCh3Logger.h:36
constexpr static const int _BAD_INT_
Default value used for int initialisation.
Definition: Core.h:55

◆ DrawPlots()

void DrawPlots ( TCanvas *  plotCanv,
TH1D *  PrefitCopy,
const std::vector< std::unique_ptr< TH1D >> &  PostfitVec,
TPad *  mainPad,
TPad *  ratioPad,
const std::string &  OutName 
)

Definition at line 259 of file GetPostfitParamPlots.cpp.

260  {
261  // KS: If plot is invalid for some reason simply do not plot
262  if(IsInvalidHist(PrefitCopy)) return;
263 
264  // Draw!
265  plotCanv->cd();
266  mainPad->Draw();
267  mainPad->cd();
268  PrefitCopy->GetYaxis()->SetTitle("Parameter Value");
269 
270  PrefitCopy->GetYaxis()->SetLabelSize(0.);
271  PrefitCopy->GetYaxis()->SetTitleSize(0.05);
272  PrefitCopy->GetYaxis()->SetTitleOffset(1.3);
273  PrefitCopy->Draw("e2");
274 
275  for (int fileId = 0; fileId < static_cast<int>(PostfitVec.size()); fileId++) {
276  TH1D *postFitHist = PostfitVec[fileId].get();
277 
278  postFitHist->SetMarkerColor(TColor::GetColorPalette(fileId));
279  postFitHist->SetLineColor(TColor::GetColorPalette(fileId));
280  postFitHist->SetMarkerStyle(7);
281  postFitHist->SetLineStyle(1+fileId);
282  postFitHist->SetLineWidth(PlotMan->getOption<int>("plotLineWidth"));
283 
284  postFitHist->Draw("e1, same");
285  }
286 
287  plotCanv->Update();
288  auto axis = std::make_unique<TGaxis>(PrefitCopy->GetXaxis()->GetBinLowEdge(PrefitCopy->GetXaxis()->GetFirst()), gPad->GetUymin()+0.01,
289  PrefitCopy->GetXaxis()->GetBinLowEdge(PrefitCopy->GetXaxis()->GetFirst()), gPad->GetUymax(),
290  gPad->GetUymin()+0.01, gPad->GetUymax(), 510, "");
291  axis->SetLabelFont(43);
292  axis->SetLabelSize(25);
293  axis->Draw();
294 
295  plotCanv->cd();
296  ratioPad->Draw();
297  ratioPad->cd();
298 
299  std::vector<std::unique_ptr<TH1D>> ratioHists;
300  // save pointers to these so we can delete them once we are done
301  ratioHists.push_back(makeRatio(PrefitCopy, PostfitVec[0].get(), true));
302 
303  ratioHists[0]->Draw("p");
304  for(int postFitIdx = 1; postFitIdx < static_cast<int>(PostfitVec.size()); postFitIdx++){
305  ratioHists.push_back(makeRatio(PrefitCopy, PostfitVec[postFitIdx].get(), true));
306 
307  ratioHists[postFitIdx]->SetMarkerColor(TColor::GetColorPalette(postFitIdx));
308  ratioHists[postFitIdx]->SetLineColor(TColor::GetColorPalette(postFitIdx));
309  ratioHists[postFitIdx]->SetMarkerStyle(7);
310  ratioHists[postFitIdx]->SetLineStyle(1+postFitIdx);
311  ratioHists[postFitIdx]->SetLineWidth(PlotMan->getOption<int>("plotLineWidth"));
312 
313  ratioHists[postFitIdx]->Draw("p same");
314  }
315 
316  // draw lines across the plot at +-1 and 0
317  TLine line(ratioHists[0]->GetXaxis()->GetBinLowEdge(ratioHists[0]->GetXaxis()->GetFirst()), 0.0, ratioHists[0]->GetXaxis()->GetBinLowEdge(ratioHists[0]->GetXaxis()->GetLast()+1), 0.0);
318  TLine line2(ratioHists[0]->GetXaxis()->GetBinLowEdge(ratioHists[0]->GetXaxis()->GetFirst()), 1.0, ratioHists[0]->GetXaxis()->GetBinLowEdge(ratioHists[0]->GetXaxis()->GetLast()+1), 1.0);
319  TLine line3(ratioHists[0]->GetXaxis()->GetBinLowEdge(ratioHists[0]->GetXaxis()->GetFirst()), -1.0, ratioHists[0]->GetXaxis()->GetBinLowEdge(ratioHists[0]->GetXaxis()->GetLast()+1), -1.0);
320 
321  line.SetLineColor(kRed);
322  line.SetLineStyle(kDashed);
323  line.SetLineWidth(PlotMan->getOption<int>("refLineWidth"));
324  line2.SetLineColor(kRed);
325  line2.SetLineStyle(kDashed);
326  line2.SetLineWidth(PlotMan->getOption<int>("refLineWidth"));
327  line3.SetLineColor(kRed);
328  line3.SetLineStyle(kDashed);
329  line3.SetLineWidth(PlotMan->getOption<int>("refLineWidth"));
330 
331  line.Draw("same");
332  line2.Draw("same");
333  line3.Draw("same");
334 
335  plotCanv->Print((OutName).c_str());
336 }
std::unique_ptr< TH1D > makeRatio(TH1D *PrefitCopy, TH1D *PostfitCopy, bool setAxes)
bool IsInvalidHist(const TH1D *hist, double invalid=M3::_BAD_DOUBLE_)

◆ Get2DComparison()

void Get2DComparison ( const std::string &  FileName1,
const std::string &  FileName2 
)

KS: Make comparison of 2D Posteriors.

Definition at line 995 of file GetPostfitParamPlots.cpp.

996 {
997  auto canvas = std::make_unique<TCanvas>("canvas", "canvas", 0, 0, 1024, 1024);
998  canvas->SetBottomMargin(0.1f);
999  canvas->SetTopMargin(0.05f);
1000  canvas->SetRightMargin(0.03f);
1001  canvas->SetLeftMargin(0.15f);
1002 
1003  // Open the two ROOT files
1004  TFile* File1 = M3::Open(FileName1, "READ", __FILE__, __LINE__);
1005  TFile* File2 = M3::Open(FileName2, "READ", __FILE__, __LINE__);
1006 
1007  // Get the Post_2d_hists directory from both files
1008  TDirectory* Dir1 = File1->Get<TDirectory>("Post_2d_hists");
1009  TDirectory* Dir2 = File2->Get<TDirectory>("Post_2d_hists");
1010 
1011  if (!Dir1 || !Dir2) {
1012  MACH3LOG_WARN("Post_2d_hists directory not found in one or both files while running {}.", __func__);
1013  File1->Close();
1014  delete File1;
1015  File2->Close();
1016  delete File2;
1017  return;
1018  }
1019 
1020  // Get the list of keys in the first directory
1021  TIter next1(Dir1->GetListOfKeys());
1022  TKey* key1 = nullptr;
1023 
1024  // Prepare the output PDF filename
1025  std::string SaveName2D = "2DComparison_" + FileName1 + "_" + FileName2;
1026  SaveName2D = SaveName2D.substr(0, SaveName2D.find(".root"));
1027  SaveName2D = SaveName2D + ".pdf";
1028 
1029  canvas->Print((SaveName2D+"[").c_str());
1030  // Loop over keys in the first directory
1031  while ((key1 = static_cast<TKey*>(next1()))) {
1032  TString histName = key1->GetName();
1033 
1034  // Check if the key is a TH2D
1035  if (TString(key1->GetClassName()) == "TH2D") {
1036  TH2D* hist1 = static_cast<TH2D*>(key1->ReadObj());
1037 
1038  // Try to get the histogram with the same name from the second directory
1039  TH2D* hist2 = static_cast<TH2D*>(Dir2->Get(histName));
1040 
1041  if (hist2) {
1042  hist1->SetTitle("");
1043  hist1->SetTitle("");
1044 
1045  // Prettify axis titles
1046  std::string Xtitle = PlotMan->style().prettifyParamName(hist1->GetXaxis()->GetTitle());
1047  std::string Ytitle = PlotMan->style().prettifyParamName(hist1->GetYaxis()->GetTitle());
1048 
1049  // Adjust the axis ranges of hist1 to include both histograms
1050  double xmin = std::min(hist1->GetXaxis()->GetXmin(), hist2->GetXaxis()->GetXmin());
1051  double xmax = std::max(hist1->GetXaxis()->GetXmax(), hist2->GetXaxis()->GetXmax());
1052  double ymin = std::min(hist1->GetYaxis()->GetXmin(), hist2->GetYaxis()->GetXmin());
1053  double ymax = std::max(hist1->GetYaxis()->GetXmax(), hist2->GetYaxis()->GetXmax());
1054 
1055  hist1->GetXaxis()->SetRangeUser(xmin, xmax);
1056  hist1->GetYaxis()->SetRangeUser(ymin, ymax);
1057 
1058  hist1->GetXaxis()->SetTitle(Xtitle.c_str());
1059  hist1->GetYaxis()->SetTitle(Ytitle.c_str());
1060 
1061  hist1->SetLineColor(kBlue);
1062  hist1->SetLineStyle(kSolid);
1063  hist1->SetLineWidth(2);
1064 
1065  hist2->SetLineColor(kRed);
1066  hist2->SetLineStyle(kDashed);
1067  hist2->SetLineWidth(2);
1068 
1069  hist1->Draw("CONT3");
1070  hist2->Draw("CONT3 SAME");
1071 
1072  auto Legend = std::make_unique<TLegend>(0.20, 0.7, 0.4, 0.92);
1073  Legend->AddEntry(hist1, PlotMan->getFileLabel(0).c_str(), "l");
1074  Legend->AddEntry(hist2, PlotMan->getFileLabel(1).c_str(), "l");
1075  Legend->SetTextSize(0.03);
1076  Legend->SetLineColor(0);
1077  Legend->SetLineStyle(0);
1078  Legend->SetFillColor(0);
1079  Legend->SetFillStyle(0);
1080  Legend->SetBorderSize(0);
1081  Legend->Draw("SAME");
1082  canvas->Print((SaveName2D).c_str());
1083  }
1084  }
1085  }
1086  canvas->Print((SaveName2D+"]").c_str());
1087 
1088  File1->Close();
1089  delete File1;
1090  File2->Close();
1091  delete File2;
1092 }
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.

◆ GetPostfitParamPlots()

void GetPostfitParamPlots ( )

Definition at line 728 of file GetPostfitParamPlots.cpp.

729 {
730  SaveName = PlotMan->getOutputName();
731 
732  //KS: By default we take HPD values, but by changing "plotType" you can use for example Gauss
733  plotType = "HPD";
734  //plotType = "gaus";
735 
736  MACH3LOG_INFO("Plotting {} errors", plotType);
737 
738  // if we have one MaCh3 nd file then we can get settings from it
739  bool plotNDDet = false;
740  for (size_t fileId = 0; fileId < PlotMan->input().getNInputFiles(); fileId++) {
741  if(!ReadSettings(PlotMan->input().getFile(0).file)) {
742  MACH3LOG_INFO("at least one file provided does not have 'settings' tree indicating it is not MaCh3 ND file");
743  MACH3LOG_INFO(" sadly this means I cannot plot ND Det parameters as this is only supported for MaCh3 ND files for now... sorry :(");
744  plotNDDet = false;
745  }
746  }
747 
748  auto canv = std::make_unique<TCanvas>("canv", "canv", 1024, 1024);
749  //gStyle->SetPalette(51);
750  gStyle->SetOptStat(0); //Set 0 to disable statistic box
751  canv->SetLeftMargin(0.12);
752  canv->SetBottomMargin(0.12);
753  canv->SetTopMargin(0.08);
754  canv->SetRightMargin(0.04);
755  canv->Print((SaveName+"[").c_str());
756 
757  // Make a Legend page
758  auto leg = std::make_unique<TLegend>(0.0, 0.0, 1.0, 1.0);
759  // make a dummy TH1 to set out legend
760  auto Prefit = std::make_unique<TH1D>();
761  Prefit->SetDirectory(nullptr);
762  PlotMan->style().setTH1Style(Prefit.get(), PlotMan->getOption<std::string>("prefitHistStyle"));
763  leg->AddEntry(Prefit.get(), "Prior", "lpf");
764 
765  std::vector<std::unique_ptr<TH1D>> postFitHist_tmp(PlotMan->getNFiles());
766  for(unsigned int fileId = 0; fileId < PlotMan->getNFiles(); fileId++){
767  postFitHist_tmp[fileId] = std::make_unique<TH1D>();
768  postFitHist_tmp[fileId]->SetDirectory(nullptr);
769 
770  postFitHist_tmp[fileId]->SetMarkerColor(TColor::GetColorPalette(fileId));
771  postFitHist_tmp[fileId]->SetLineColor(TColor::GetColorPalette(fileId));
772  postFitHist_tmp[fileId]->SetMarkerStyle(7);
773  postFitHist_tmp[fileId]->SetLineStyle(1+fileId);
774  postFitHist_tmp[fileId]->SetLineWidth(PlotMan->getOption<int>("plotLineWidth"));
775  leg->AddEntry(postFitHist_tmp[fileId].get(), PlotMan->getFileLabel(fileId).c_str(), "lpf");
776  }
777 
778  canv->cd();
779  canv->Clear();
780  leg->Draw();
781  canv->Print((SaveName).c_str());
782 
783  MakeParameterPlots(canv.get());
784 
785  MakeFluxPlots(canv.get());
786 
787  //KS: By default we don't run ProcessMCMC with PlotDet as this take some time, in case we did let's make fancy plots
788  if(plotNDDet & (NDParameters > 0)) MakeNDDetPlots();
789 
790  canv->Print((SaveName+"]").c_str());
791 
792  MakeRidgePlots();
793 }
void MakeFluxPlots(TCanvas *canv)
void MakeParameterPlots(TCanvas *canv)
void MakeNDDetPlots()
bool ReadSettings(const std::shared_ptr< TFile > &File1)
std::string SaveName
void MakeRidgePlots()
int NDParameters
std::string plotType
#define MACH3LOG_INFO
Definition: MaCh3Logger.h:35

◆ GetViolinPlots()

void GetViolinPlots ( )

KS: Make fancy violin plots.

Definition at line 830 of file GetPostfitParamPlots.cpp.

831 {
832  //KS: Should be in some config... either way it control whether you plot symmetric or asymmetric error bars
833  bool PlotAssym = true;
834 
835  //KS: No idea why but ROOT changed treatment of violin in R6. If you have non uniform binning this will results in very hard to see violin plots.
836  TCandle::SetScaledViolin(false);
837 
838  std::string OutputName = "";
839  for(unsigned int fileId = 0; fileId < PlotMan->getNFiles(); fileId++){
840  MACH3LOG_INFO("File {}: {} ", fileId, PlotMan->getFileName(fileId));
841  OutputName += PlotMan->getFileName(fileId);
842  OutputName = OutputName.substr(0, OutputName.find(".root"));
843  }
844  MACH3LOG_INFO("Making Violin Plot");
845  OutputName += "_Violin";
846  if(PlotAssym) OutputName += "_Assym";
847 
848  auto canvas = std::make_unique<TCanvas>("canv", "canv", 1024, 1024);
849  canvas->SetGrid();
850  gStyle->SetOptStat(0);
851  //KS: Remove errors on X axis as they are confusing in violin type of plot
852  if(!PlotAssym) gStyle->SetErrorX(0.0001);
853  canvas->SetTickx();
854  canvas->SetTicky();
855  canvas->SetBottomMargin(0.25);
856  canvas->SetTopMargin(0.08);
857  canvas->SetRightMargin(0.03);
858  canvas->SetLeftMargin(0.10);
859  canvas->Print((OutputName+".pdf[").c_str());
860  canvas->SetGrid();
861 
862  if(PlotMan->input().getFile(0).file->Get<TH2D>( "param_violin_prior" ) == nullptr)
863  {
864  MACH3LOG_WARN("Couldn't find violin plot, make sure method from MCMCProcessor is being called");
865  return;
866  }
867  std::unique_ptr<TH2D> ViolinPre = M3::Clone(PlotMan->input().getFile(0).file->Get<TH2D>( "param_violin_prior" ));
868  // Do some fancy replacements
869  ViolinPre->SetFillColor(kRed);
870  ViolinPre->SetFillColorAlpha(kRed, 0.35);
871  ViolinPre->SetMarkerColor(kRed);
872  ViolinPre->SetMarkerStyle(20);
873  ViolinPre->SetMarkerSize(0.5);
874 
875  ViolinPre->GetYaxis()->SetTitleOffset(1.3);
876  ViolinPre->GetYaxis()->SetTitle("Parameter Value");
877  ViolinPre->GetXaxis()->LabelsOption("v");
878 
879  std::unique_ptr<TH1D> Postfit = M3::Clone(PlotMan->input().getFile(0).file->Get<TH1D>( ("param_xsec_"+plotType).c_str() ));
880  Postfit->SetMarkerColor(kRed);
881  Postfit->SetLineColor(kRed);
882  Postfit->SetMarkerStyle(7);
883 
884  std::vector<std::unique_ptr<TH2D>> Violin(PlotMan->getNFiles());
885  //KS: I know hardcoded but we can figure out later...
886  const std::vector<Color_t> FillColors = { kBlue, kGreen, kMagenta };
887  for(unsigned int fileId = 0; fileId < PlotMan->getNFiles(); fileId++) {
888  Violin[fileId] = M3::Clone(PlotMan->input().getFile(fileId).file->Get<TH2D>( "param_violin" ));
889  if(Violin[fileId] == nullptr)
890  {
891  MACH3LOG_ERROR("Couldn't find violin plot, make sure method from MCMCProcessor is being called");
892  return;
893  }
894  Violin[fileId]->SetMarkerColor(FillColors.at(fileId));
895  Violin[fileId]->SetLineColor(FillColors.at(fileId));
896  Violin[fileId]->SetFillColor(FillColors.at(fileId));
897  Violin[fileId]->SetFillColorAlpha(FillColors.at(fileId), 0.35);
898 
899  if(fileId == 0){
900  Violin[fileId]->SetMarkerStyle(20);
901  Violin[fileId]->SetMarkerSize(0.5);
902  }
903  }
904 
905  std::unique_ptr<TGraphAsymmErrors> PostGraphAll = MakeTGraphAsymmErrors(PlotMan->input().getFile(0).file);
906  // Make a Legend page
907  auto leg = std::make_unique<TLegend>(0.0, 0.0, 1.0, 1.0);
908  if (ViolinPre != nullptr) leg->AddEntry(ViolinPre.get(), "Prior", "lpf");
909  for(unsigned int fileId = 0; fileId < PlotMan->getNFiles(); fileId++) {
910  leg->AddEntry(Violin[fileId].get(), PlotMan->getFileLabel(fileId).c_str(), "lpf");
911  }
912  if(PlotAssym) leg->AddEntry(PostGraphAll.get(), "HPD Assym", "lp");
913  else leg->AddEntry(Postfit.get(), "HPD", "lpf");
914 
915  canvas->cd();
916  canvas->Clear();
917  leg->Draw();
918  canvas->Print((OutputName+".pdf").c_str());
919 
920  // get the names of the blocks of parameters to group together
921  std::vector<std::string> const blockNames = PlotMan->getOption<std::vector<std::string>>("paramGroups");
922  const int nPlots = static_cast<int>(blockNames.size());
923 
924  for (int i = 0; i < nPlots; i++)
925  {
926  // get the configuration for this parameter
927  std::string blockName = blockNames[i];
928  YAML::Node paramBlock = PlotMan->getOption(blockName);
929  std::string blockTitle = paramBlock[0].as<std::string>();
930  std::vector<double> blockLimits = paramBlock[1].as<std::vector<double>>();
931  std::vector<std::string> blockContents = paramBlock[2].as<std::vector<std::string>>();
932 
933  // get num of params in the block
934  const int nParams = static_cast<int>(blockContents.size());
935 
936  // set some plot things
937  auto blockHist_prefit = std::make_unique<TH2D>((blockName + "_Prefit").c_str(), blockTitle.c_str(), nParams, 0.0, static_cast<double>(nParams),
938  ViolinPre->GetYaxis()->GetNbins(), ViolinPre->GetYaxis()->GetXmin(), ViolinPre->GetYaxis()->GetXmax());
939  CopyViolinToBlock(ViolinPre.get(), blockHist_prefit.get(), blockContents);
940  // set the y axis limits we got from config
941  blockHist_prefit->GetYaxis()->SetRangeUser(blockLimits[0], blockLimits[1]);
942 
943  std::vector<std::unique_ptr<TH2D>> blockHist(PlotMan->getNFiles());
944  for(unsigned int fileId = 0; fileId < PlotMan->getNFiles(); fileId++) {
945  blockHist[fileId] = std::make_unique<TH2D>((blockTitle + "Violin" + fileId).Data(), (blockTitle + "Violin" + fileId).Data(),
946  nParams, 0.0, static_cast<double>(nParams), Violin[fileId]->GetYaxis()->GetNbins(),
947  Violin[fileId]->GetYaxis()->GetXmin(), Violin[fileId]->GetYaxis()->GetXmax());
948  CopyViolinToBlock(Violin[fileId].get(), blockHist[fileId].get(), blockContents);
949  blockHist[fileId]->GetYaxis()->SetRangeUser(blockLimits[fileId], blockLimits[1]);
950  }
951  // Do some fancy replacements
952  PrettifyTitles(blockHist_prefit.get());
953 
954  // set some plot things
955  auto blockHist_Best = std::make_unique<TH1D>(blockName.c_str(), blockTitle.c_str(),
956  nParams, 0.0, static_cast<double>(nParams));
957  // set the errors for the prefit block hist
958  for(int localBin=0; localBin < nParams; localBin ++){
959  // the "local" bin is the params index within the group of parameters
960  std::string paramName = blockContents[localBin];
961  copyParToBlockHist(localBin, paramName, blockHist_Best.get(), "", 0);
962  }
963 
964  std::vector<int> Index;
965  for(unsigned int is = 0; is < blockContents.size(); is++) {
966  int ParamBinId = M3::_BAD_INT_;
967  for (int ix = 0; ix < ViolinPre->GetXaxis()->GetNbins(); ++ix) {
968  if(ViolinPre->GetXaxis()->GetBinLabel(ix+1) == blockContents[is]) {
969  ParamBinId = ix;
970  break;
971  }
972  }
973  Index.push_back(ParamBinId);
974  }
975  std::unique_ptr<TGraphAsymmErrors> PostGraph = MakeTGraphAsymmErrors(PlotMan->input().getFile(0).file, Index);
976  PostGraph->SetMarkerColor(kBlack);
977  PostGraph->SetLineColor(kBlack);
978  PostGraph->SetMarkerStyle(7);
979  PostGraph->SetLineWidth(2);
980  PostGraph->SetLineStyle(kSolid);
981 
982  blockHist_prefit->Draw("violinX(03100300)");
983  for(unsigned int fileId = 0; fileId < PlotMan->getNFiles(); fileId++) {
984  blockHist[fileId]->Draw("violinX(03100300) SAME");
985  }
986 
987  if(PlotAssym) PostGraph->Draw("P SAME");
988  else Postfit->Draw("SAME");
989  canvas->Print((OutputName+".pdf").c_str());
990  }
991  canvas->Print((OutputName+".pdf]").c_str());
992 }
void copyParToBlockHist(const int localBin, const std::string &paramName, TH1D *blockHist, const std::string &type, const int fileId, const bool setLabels=true)
void PrettifyTitles(HistType *hist)
std::unique_ptr< TGraphAsymmErrors > MakeTGraphAsymmErrors(const std::shared_ptr< TFile > &File, std::vector< int > Index={})
void CopyViolinToBlock(TH2D *FullViolin, TH2D *ReducedViolin, const std::vector< std::string > &ParamNames)
#define MACH3LOG_ERROR
Definition: MaCh3Logger.h:37
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.

◆ InitializePads()

void InitializePads ( const TCanvas *  canvas,
TPad *&  pad3,
TPad *&  pad4 
)

Definition at line 79 of file GetPostfitParamPlots.cpp.

79  {
80  // Initialize TPad p3
81  pad3 = new TPad("Top", "Top", 0.0, 0.4, 1.0, 1.0);
82  pad3->SetLeftMargin(canvas->GetLeftMargin());
83  pad3->SetRightMargin(canvas->GetRightMargin());
84  pad3->SetTopMargin(canvas->GetTopMargin());
85  pad3->SetBottomMargin(0);
86  pad3->SetGrid();
87 
88  // Initialize TPad p4
89  pad4 = new TPad("Bottom", "Bottom", 0.0, 0.0, 1.0, 0.4);
90  pad4->SetLeftMargin(canvas->GetLeftMargin());
91  pad4->SetRightMargin(canvas->GetRightMargin());
92  pad4->SetTopMargin(0);
93  pad4->SetBottomMargin(0.75);
94  pad4->SetGrid();
95 }

◆ IsInvalidHist()

bool IsInvalidHist ( const TH1D *  hist,
double  invalid = M3::_BAD_DOUBLE_ 
)

Definition at line 240 of file GetPostfitParamPlots.cpp.

241 {
242  if (!hist) return true;
243 
244  const int nBins = hist->GetNbinsX();
245 
246  for (int i = 1; i <= nBins; i++) {
247  double val = hist->GetBinContent(i);
248 
249  // If ANY bin is valid → histogram is usable
250  if (val != invalid) {
251  return false;
252  }
253  }
254 
255  // All bins are invalid
256  return true;
257 }

◆ main()

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

Definition at line 1094 of file GetPostfitParamPlots.cpp.

1095 {
1097  // Avoid Info in <TCanvas::Print>
1098  gErrorIgnoreLevel = kWarning;
1099 
1100  PlotMan = new M3::Plotting::PlottingManager();
1101  PlotMan->parseInputs(argc, argv);
1102  #ifdef MACH3_DEBUG
1103  PlotMan->input().getFile(0).file->ls();
1104  #endif
1105  PlotMan->setExec("GetPostfitParamPlots");
1106 
1107  PlotMan->style().setPalette(PlotMan->getOption<std::string>("colorPalette"));
1108 
1110  GetViolinPlots();
1111 
1112  if (PlotMan->input().getNInputFiles() == 2)
1113  {
1114  std::string filename1 = PlotMan->getFileName(0);
1115  std::string filename2 = PlotMan->getFileName(1);
1116  Get2DComparison(filename1, filename2);
1117  }
1118 
1119  if(PlotMan) delete PlotMan;
1120  return 0;
1121 }
void Get2DComparison(const std::string &FileName1, const std::string &FileName2)
KS: Make comparison of 2D Posteriors.
void GetViolinPlots()
KS: Make fancy violin plots.
void SetMaCh3LoggerFormat()
Set messaging format of the logger.
Definition: MaCh3Logger.h:61

◆ MakeFluxPlots()

void MakeFluxPlots ( TCanvas *  canv)

Definition at line 400 of file GetPostfitParamPlots.cpp.

401 {
402  // these for non named params where we don't need as much space
403  auto p1 = std::make_unique<TPad>("p1", "p1", 0.0, 0.3, 1.0, 1.0);
404  auto p2 = std::make_unique<TPad>("p2", "p2", 0.0, 0.0, 1.0, 0.3);
405  p1->SetLeftMargin(canv->GetLeftMargin());
406  p1->SetRightMargin(canv->GetRightMargin());
407  p1->SetTopMargin(canv->GetTopMargin());
408  p1->SetBottomMargin(0);
409  p1->SetGrid();
410  p2->SetLeftMargin(canv->GetLeftMargin());
411  p2->SetRightMargin(canv->GetRightMargin());
412  p2->SetTopMargin(0);
413  p2->SetBottomMargin(0.25);
414  p2->SetGrid();
415 
416  p1->SetLogx(true);
417  p2->SetLogx(true);
418 
419  // get the names of the blocks of parameters to group together
420  std::vector<std::string> const fluxBlockNames = PlotMan->getOption<std::vector<std::string>>("fluxGroups");
421  auto const fluxBinningTable = PlotMan->getOption("FluxBinning");
422 
423  const int FluxPlots = static_cast<int>(fluxBlockNames.size());
424 
425  for (int i = 0; i < FluxPlots; i++)
426  {
427  // get the configuration for this block
428  std::string fluxBlockName = fluxBlockNames[i];
429  YAML::Node paramBlock = PlotMan->getOption(fluxBlockName);
430  std::string blockTitle = paramBlock[0].as<std::string>();
431  std::vector<double> blockLimits = paramBlock[1].as<std::vector<double>>();
432  std::string blockBinningName = paramBlock[2].as<std::string>();
433  std::vector<int> blockContents = paramBlock[3].as<std::vector<int>>();
434 
435  // get the binning for this block of flux params
436  std::vector<double> binning = fluxBinningTable[blockBinningName].as<std::vector<double>>();
437 
438  // get num of params in the block
439  int nParams = blockContents[1] - blockContents[0] +1;
440  // check for sanity
441  if(nParams <= 0 || blockContents.size() > 2){
442  MACH3LOG_CRITICAL("Invalid flux parameter block endpoints specified for {}", fluxBlockName);
443  MACH3LOG_CRITICAL(" Should have the form [<low index>, <up index>]");
444  throw MaCh3Exception(__FILE__ , __LINE__ );
445  }
446  if (nParams != static_cast<int>(binning.size()) - 1) {
447  MACH3LOG_CRITICAL("Binning provided for flux param block {} does not match the number of parameters specified for the block", fluxBlockName);
448  MACH3LOG_CRITICAL(" Provided {} parameters but {} bins", nParams, binning.size() -1);
449  throw MaCh3Exception(__FILE__ , __LINE__ );
450  }
451 
452  auto blockHist_prefit = std::make_unique<TH1D>(fluxBlockName.c_str(), blockTitle.c_str(), nParams, binning.data());
453  blockHist_prefit->GetYaxis()->SetTitle("Parameter Variation");
454  blockHist_prefit->GetXaxis()->SetTitle("E_{#nu} (GeV)");
455  blockHist_prefit->GetXaxis()->SetTitleOffset(blockHist_prefit->GetXaxis()->GetTitleOffset()*1.2);
456  PlotMan->style().setTH1Style(blockHist_prefit.get(), PlotMan->getOption<std::string>("prefitHistStyle"));
457  // set the errors for the prefit block hist
458  for(int fluxParId = blockContents[0]; fluxParId <= blockContents[1]; fluxParId++){
459  int localBin = fluxParId - blockContents[0];
460  std::string paramName = "b_" + std::to_string(fluxParId);
461  copyParToBlockHist(localBin, paramName, blockHist_prefit.get(), "Prior", 0, false);
462  }
463 
464  // now set for the postfit blocks for all files
465  std::vector <std::unique_ptr<TH1D>> blockHist_postfit_Vec;
466  for(unsigned int fileId = 0; fileId < PlotMan->getNFiles(); fileId++){
467  auto blockHist_postfit = std::make_unique<TH1D>(fluxBlockName.c_str(), blockTitle.c_str(), nParams, binning.data());
468 
469  for(int fluxParId = blockContents[0]; fluxParId <= blockContents[1]; fluxParId++){
470  int localBin = fluxParId - blockContents[0];
471  std::string paramName = "b_" + std::to_string(fluxParId);
472 
473  copyParToBlockHist(localBin, paramName, blockHist_postfit.get(), "", fileId, false);
474  }
475  blockHist_postfit_Vec.push_back(std::move(blockHist_postfit));
476  }
477  // set the y axis limits we got from config
478  blockHist_prefit->GetYaxis()->SetRangeUser(blockLimits[0], blockLimits[1]);
479  DrawPlots(canv, blockHist_prefit.get(), blockHist_postfit_Vec, p1.get(), p2.get(), SaveName);
480  }
481 
482  canv->cd();
483  canv->SetLogx(false);
484  canv->SetBottomMargin(canv->GetBottomMargin()*1.7);
485 }
void DrawPlots(TCanvas *plotCanv, TH1D *PrefitCopy, const std::vector< std::unique_ptr< TH1D >> &PostfitVec, TPad *mainPad, TPad *ratioPad, const std::string &OutName)
#define MACH3LOG_CRITICAL
Definition: MaCh3Logger.h:38
Custom exception class used throughout MaCh3.

◆ MakeNDDetPlots()

void MakeNDDetPlots ( )
Warning
This is legacy functions and will become deprecated

Definition at line 488 of file GetPostfitParamPlots.cpp.

489 {
490  MACH3LOG_INFO("ND detector parameters: {}", NDParameters);
491  auto canv = std::make_unique<TCanvas>("canv", "canv", 1024, 1024);
492  //gStyle->SetPalette(51);
493  gStyle->SetOptStat(0); //Set 0 to disable statistic box
494  canv->SetLeftMargin(0.12);
495  canv->SetBottomMargin(0.12);
496  canv->SetTopMargin(0.08);
497  canv->SetRightMargin(0.04);
498 
499  // make a dummy TH1 to set out legend
500  TH1D* Prefit = new TH1D();
501  PlotMan->style().setTH1Style(Prefit, PlotMan->getOption<std::string>("prefitHistStyle"));
502  Prefit->GetYaxis()->SetTitleOffset(Prefit->GetYaxis()->GetTitleOffset()*1.2);
503 
504  TPad* pTop = nullptr;
505  TPad* pDown = nullptr;
506  InitializePads(canv.get(), pTop, pDown);
507 
508  int NDbinCounter = NDParametersStartingPos;
509  int Start = NDbinCounter;
510 
511  MACH3LOG_INFO("Running on {} samples", NDSamplesNames.size());
512 
513  for (unsigned int i = 0; i < NDSamplesNames.size(); ++i)
514  {
515  MACH3LOG_DEBUG("--- On sample {}", NDSamplesNames[i]);
516  NDbinCounter += NDSamplesBins[i];
517 
518  std::vector<std::unique_ptr<TH1D>> PostfitNDDetHistVec(PlotMan->getNFiles());
519  TH1D *PreFitNDDetHist = PlotMan->input().getFile(0).file->Get<TH1D>(Form("param_%s_prefit", NDSamplesNames[i].c_str()));
520  PlotMan->style().setTH1Style(PreFitNDDetHist, PlotMan->getOption<std::string>("prefitHistStyle"));
521 
522  std::string temp = NDSamplesNames[i].c_str();
523  while (temp.find("_") != std::string::npos) {
524  temp.replace(temp.find("_"), 1, std::string(" "));
525  }
526  PreFitNDDetHist->SetTitle(temp.c_str());
527  PreFitNDDetHist->GetXaxis()->SetRangeUser(Start, NDbinCounter);
528 
529  MACH3LOG_DEBUG(" Start bin: {} :: End bin: {}", Start, NDbinCounter);
530  // set the x range for the postfits
531  for(unsigned int fileId = 0; fileId < PlotMan->getNFiles(); fileId++){
532  PostfitNDDetHistVec[fileId] = M3::Clone(PlotMan->input().getFile(fileId).file->Get<TH1D>(Form("param_%s_%s", NDSamplesNames[i].c_str(), plotType.c_str())));
533  }
534 
535  //KS: We don't' need name for every nd param
536  for(int j = 0; j < NDSamplesBins[i]; ++j)
537  {
538  bool ProductOfTen = false;
539  if(j % 10) ProductOfTen = true;
540  if(j != 0 && ProductOfTen) PreFitNDDetHist->GetXaxis()->SetBinLabel(Start+j+1, " ");
541  else{
542  PreFitNDDetHist->GetXaxis()->SetBinLabel(Start+j+1, Form("Det Variation Bin %i", Start+j));
543  }
544  }
545 
546  PreFitNDDetHist->GetYaxis()->SetRangeUser(PlotMan->getOption<double>("detParYRange_low"), PlotMan->getOption<double>("detParYRange_high"));
547 
548  Start += NDSamplesBins[i];
549 
550  DrawPlots(canv.get(), PreFitNDDetHist, PostfitNDDetHistVec, pTop, pDown, SaveName);
551  canv->Update();
552  }
553  delete pTop;
554  delete pDown;
555  delete Prefit;
556 }
void InitializePads(const TCanvas *canvas, TPad *&pad3, TPad *&pad4)
std::vector< int > NDSamplesBins
int NDParametersStartingPos
std::vector< std::string > NDSamplesNames

◆ MakeParameterPlots()

void MakeParameterPlots ( TCanvas *  canv)

Definition at line 338 of file GetPostfitParamPlots.cpp.

339 {
340  // get the names of the blocks of parameters to group together
341  std::vector<std::string> const blockNames = PlotMan->getOption<std::vector<std::string>>("paramGroups");
342  const int nPlots = static_cast<int>(blockNames.size());
343 
344  TPad *p3 = nullptr;
345  TPad *p4 = nullptr;
346  // these for named parameters where we need a nice big gap at the botto to fit the names
347  InitializePads(canv, p3, p4);
348 
349  for (int i = 0; i < nPlots; i++)
350  {
351  // get the configuration for this parameter
352  std::string blockName = blockNames[i];
353  YAML::Node paramBlock = PlotMan->getOption(blockName);
354  auto blockTitle = paramBlock[0].as<std::string>();
355  auto blockLimits = paramBlock[1].as<std::vector<double>>();
356  auto blockContents = paramBlock[2].as<std::vector<std::string>>();
357 
358  // get num of params in the block
359  const int nParams = static_cast<int>(blockContents.size());
360 
361  // set some plot things
362  auto blockHist_prefit = std::make_unique<TH1D>(blockName.c_str(), blockTitle.c_str(),
363  nParams, 0.0, static_cast<double>(nParams));
364 
365  PlotMan->style().setTH1Style(blockHist_prefit.get(), PlotMan->getOption<std::string>("prefitHistStyle"));
366 
367  // set the errors for the prefit block hist
368  for(int localBin=0; localBin < nParams; localBin ++){
369  // the "local" bin is the params index within the group of parameters
370  std::string paramName = blockContents[localBin];
371  copyParToBlockHist(localBin, paramName, blockHist_prefit.get(), "Prior", 0);
372  }
373 
374  // now set for the postfit blocks for all files
375  std::vector <std::unique_ptr<TH1D>> blockHist_postfit_Vec;
376  for(unsigned int fileId = 0; fileId < PlotMan->getNFiles(); fileId++) {
377  auto blockHist_postfit = std::make_unique<TH1D>((blockName + PlotMan->getFileName(fileId)).c_str(),
378  blockTitle.c_str(), nParams, 0.0, static_cast<double>(nParams));
379 
380  // loop through all the parameters in this block and set the contents in the blocks TH1
381  for(int localBin=0; localBin < nParams; localBin ++){
382  // the "local" bin is the params index within the group of parameters
383  std::string paramName = blockContents[localBin];
384  copyParToBlockHist(localBin, paramName, blockHist_postfit.get(), "", fileId);
385  }
386  blockHist_postfit_Vec.push_back(std::move(blockHist_postfit));
387  }
388  // set the y axis limits we got from config
389  blockHist_prefit->GetYaxis()->SetRangeUser(blockLimits[0], blockLimits[1]);
390 
391  // Do some fancy replacements
392  PrettifyTitles(blockHist_prefit.get());
393 
394  DrawPlots(canv, blockHist_prefit.get(), blockHist_postfit_Vec, p3, p4, SaveName);
395  }
396  delete p3;
397  delete p4;
398 }

◆ makeRatio()

std::unique_ptr<TH1D> makeRatio ( TH1D *  PrefitCopy,
TH1D *  PostfitCopy,
bool  setAxes 
)

Definition at line 172 of file GetPostfitParamPlots.cpp.

172  {
173  // set up the ratio hist
174  std::unique_ptr<TH1D> Ratio = M3::Clone(PrefitCopy);
175  Ratio->GetYaxis()->SetTitle("(x_{Post}-#mu_{Prior})/#sigma_{Prior}");
176  Ratio->SetMinimum(-3.7);
177  Ratio->SetMaximum(3.7);
178 
179  for (int j = 0; j < Ratio->GetXaxis()->GetNbins(); ++j)
180  {
181  if ( PrefitCopy->GetBinError(j+1) > 1.e-5 )
182  {
183  Ratio->SetBinContent(j+1, (PostfitCopy->GetBinContent(j+1)-PrefitCopy->GetBinContent(j+1))/PrefitCopy->GetBinError(j+1));
184 
185  double up = (PostfitCopy->GetBinContent(j+1)+PostfitCopy->GetBinError(j+1)-PrefitCopy->GetBinContent(j+1))/PrefitCopy->GetBinError(j+1);
186  double down = (PostfitCopy->GetBinContent(j+1)-PostfitCopy->GetBinError(j+1)-PrefitCopy->GetBinContent(j+1))/PrefitCopy->GetBinError(j+1);
187 
188  double maximum = up-Ratio->GetBinContent(j+1);
189  double minimum = Ratio->GetBinContent(j+1)-down;
190 
191  Ratio->SetBinError(j+1, std::max(maximum, minimum));
192  }
193  //KS: Most likely flat prior
194  else {
195  Ratio->SetBinContent(j+1, (PostfitCopy->GetBinContent(j+1)-PrefitCopy->GetBinContent(j+1)));
196 
197  double up = (PostfitCopy->GetBinContent(j+1)+PostfitCopy->GetBinError(j+1)-PrefitCopy->GetBinContent(j+1));
198  double down = (PostfitCopy->GetBinContent(j+1)-PostfitCopy->GetBinError(j+1)-PrefitCopy->GetBinContent(j+1));
199 
200  double maximum = up-Ratio->GetBinContent(j+1);
201  double minimum = Ratio->GetBinContent(j+1)-down;
202 
203  Ratio->SetBinError(j+1, std::max(maximum, minimum));
204  }
205  } //end loop over parameters
206 
207  if(setAxes){
208  Ratio->SetFillStyle(0);
209  Ratio->SetFillColor(0);
210 
211  Ratio->SetLineColor(PostfitCopy->GetLineColor());
212  if (Ratio->GetLineColor() == 0) Ratio->SetLineColor(kBlack);
213  Ratio->SetMarkerColor(PostfitCopy->GetMarkerColor());
214 
215  Ratio->SetLineWidth(PlotMan->getOption<int>("plotLineWidth"));
216  Ratio->SetTitle("");
217 
218  Ratio->SetMarkerSize(2);
219  Ratio->SetMarkerStyle(20);
220 
221  Ratio->GetYaxis()->SetTitleSize(25);
222  Ratio->GetYaxis()->SetTitleFont(43);
223  Ratio->GetYaxis()->SetTitleOffset(2.0);
224  Ratio->GetYaxis()->SetLabelFont(43);
225  Ratio->GetYaxis()->SetLabelSize(25);
226  Ratio->GetYaxis()->CenterTitle();
227  Ratio->GetYaxis()->SetNdivisions(5,2,0);
228 
229  Ratio->GetXaxis()->SetTitleSize(25);
230  Ratio->GetXaxis()->SetTitleFont(43);
231  Ratio->GetXaxis()->SetTitleOffset(4.0);
232  Ratio->GetXaxis()->SetLabelOffset(0.025);
233  Ratio->GetXaxis()->SetLabelFont(43);
234  Ratio->GetXaxis()->SetLabelSize(25);
235  }
236 
237  return Ratio;
238 }

◆ MakeRidgePlots()

void MakeRidgePlots ( )

Definition at line 558 of file GetPostfitParamPlots.cpp.

559 {
560  gStyle->SetPalette(51);
561 
562  auto blankCanv = std::make_unique<TCanvas>("blankCanv", "blankCanv", 2048, 2048);
563  blankCanv->SaveAs("RidgePlots.pdf[");
564 
565  // get the names of the blocks of parameters to group together
566  const auto blockNames = PlotMan->getOption<std::vector<std::string>>("paramGroups");
567  const int nPlots = static_cast<int>(blockNames.size());
568 
569  constexpr double padTopMargin = 0.9;
570  constexpr double padBottomMargin = 0.1;
571  constexpr double padOverlap = 0.9;
572  constexpr double ridgeLineWidth = 1.0;
573  for (int i = 0; i < nPlots; i++)
574  {
575  // get the configuration for this parameter
576  std::string blockName = blockNames[i];
577  auto const &paramBlock = PlotMan->getOption(blockName);
578  auto blockTitle = paramBlock[0].as<std::string>();
579  auto blockLimits = paramBlock[1].as<std::vector<double>>();
580  auto blockContents = paramBlock[2].as<std::vector<std::string>>();
581 
582  // the directory of histograms
583  TDirectoryFile *posteriorDir = PlotMan->input().getFile(0).file->Get<TDirectoryFile>("Post_1d_hists");
584 
585  // get num of params in the block
586  int nParams = static_cast<int>(blockContents.size());
587 
588  if (nParams == 1) {
589  MACH3LOG_WARN("{} doesn't work for single param", __func__);
590  continue;
591  }
592  auto ridgeCanv = std::make_unique<TCanvas>("RidgePlotCanv", "RidgePlotCanv", 2048, 2048);
593  ridgeCanv->Divide(1,1+nParams, 0.01, 0.0);
594 
595  auto title = std::make_unique<TLatex>();
596  title->SetTextAlign(21);
597  title->SetTextSize(0.03);
598  title->DrawLatex(0.5, 0.95, blockTitle.c_str());
599 
600  auto label = std::make_unique<TLatex>();
601  label->SetTextAlign(31);
602  label->SetTextSize(0.02);
603 
604  auto line = std::make_unique<TLine>();
605  line->SetLineColor(kBlack);
606  line->SetLineWidth(ridgeLineWidth);
607 
608  // use this to set the limits and also to plot the x axis and grid
609  auto axisPlot = std::make_unique<TH1D>("axis plot", "", 1, blockLimits[0], blockLimits[1]);
610 
611  std::vector<std::unique_ptr<TH1D>> axisPlot_holder(nParams);
612  std::vector<std::unique_ptr<TPad>> graph_holder(nParams);
613  for(int parId = 0; parId < nParams; parId++) {
614  std::string paramName = blockContents[parId];
615 
616  TH1D *posteriorDist = nullptr;
617  // get the list of objects in the directory
618  TIter next(posteriorDir->GetListOfKeys());
619  while (TKey* key = static_cast<TKey*>(next())) {
620  // check if the end of the param name matches with the MaCh3 name, do this so we exclude things like nds_ at the start of the name
621  std::string str(key->GetTitle());
622  std::string name = PlotMan->input().translateName(0, M3::Plotting::kPostFit, paramName);
623  uint pos = str.find(name);
624  bool foundPar = (pos == str.length() - name.length());
625 
626  MACH3LOG_TRACE("Looking for {} in {}", name, str);
627  if(foundPar){
628  MACH3LOG_TRACE("Found it");
629  posteriorDist = posteriorDir->Get<TH1D>(key->GetName());
630  }
631  }
632 
633  if(posteriorDist == nullptr){
634  MACH3LOG_WARN("Couldn't find parameter {} when making ridgeline plots", paramName);
635  MACH3LOG_WARN("It could be fixed param");
636  continue;
637  }
638 
639  // EM: do some funky scaling so that we always get evenly spaced pads in the range [bottomMargin, TopMargin] with the specified overlap
640  double padAnchor = padBottomMargin + (static_cast<double>(nParams - parId - 1) /
641  static_cast<double>(nParams - 1)) * (padTopMargin - padBottomMargin);
642  double padWidth = (padTopMargin - padBottomMargin) / static_cast<double>(nParams);
643  double norm = (padTopMargin - padBottomMargin);
644 
645  double padTop = padWidth * (1.0 + padOverlap) * (padTopMargin - padAnchor) / norm + padAnchor;
646  double padBottom = padAnchor - padWidth * (1.0 + padOverlap) * (padAnchor - padBottomMargin) / norm;
647 
648  auto pad = std::make_unique<TPad>(paramName.c_str(), "", 0.3, padBottom, 0.9, padTop, -1, 0, -1);
649  ridgeCanv->cd();
650 
651  pad->SetBottomMargin(0.0);
652  pad->SetTopMargin(0.0);
653  pad->SetLeftMargin(0.0);
654  pad->SetRightMargin(0.0);
655 
656  pad->Draw();
657  pad->cd();
658  pad->SetFillStyle(4000);
659 
660  gPad->SetFrameFillStyle(4000);
661  TF1* gausFunc = posteriorDist->GetFunction("Gauss");
662  if (posteriorDist->GetFunction("Gauss")) {
663  gausFunc->SetBit(TF1::kNotDraw);
664  }
665  posteriorDist->SetTitle("");
666  posteriorDist->SetLineWidth(ridgeLineWidth);
667 
668  auto axisPlot_tmp = M3::Clone(axisPlot.get(), Form("AxisPlot_%s", paramName.c_str()));
669  axisPlot_tmp->Draw("A");
670  posteriorDist->Draw("H SAME");
671 
672  axisPlot_tmp->GetYaxis()->SetRangeUser(0.0, 0.7 *posteriorDist->GetMaximum());
673  posteriorDist->SetLineColor(kWhite);
674  posteriorDist->SetFillColorAlpha(TColor::GetColorPalette(floor(static_cast<float>(parId) *
675  TColor::GetNumberOfColors() / static_cast<float>(nParams))), 0.85);
676 
677  posteriorDist->GetXaxis()->SetRangeUser(blockLimits[0], blockLimits[1]);
678  posteriorDist->GetYaxis()->SetTitle(paramName.c_str());
679 
680  //EM: Have to get the frame drawn by the histogram and then set it to be transparent... it took me an hour to get rid of this one line on a plot
681  gPad->Modified(); gPad->Update();
682  TFrame *frame = gPad->GetFrame();
683  frame->SetLineColorAlpha(0, 0.0);
684 
685  ridgeCanv->cd();
686  label->DrawLatexNDC(0.29, padBottom + 0.005, PlotMan->style().prettifyParamName(paramName).c_str());
687  line->DrawLine(0.1, padBottom, 0.9, padBottom);
688 
689  axisPlot_holder[parId] = std::move(axisPlot_tmp);
690  graph_holder[parId] = std::move(pad);
691  }
692 
693  ridgeCanv->cd();
694  ridgeCanv->SetGrid(1,1);
695  auto axisPad = std::make_unique<TPad>("AxisPad", "", 0.3, 0.0, 0.9, 1.0, -1, 0, -1);
696  axisPad->SetLeftMargin(0.0);
697  axisPad->SetRightMargin(0.0);
698  axisPad->Draw();
699  axisPad->cd();
700  axisPad->SetGrid(1,1);
701  axisPad->SetFrameFillStyle(4000);
702 
703  axisPlot->GetXaxis()->SetTickSize(0.01);
704  axisPlot->GetXaxis()->SetTitle("Parameter Variation");
705  axisPlot->GetYaxis()->SetLabelOffset(9999);
706  axisPlot->GetYaxis()->SetLabelSize(0);
707  axisPlot->GetYaxis()->SetTickSize(0);
708  axisPlot->GetYaxis()->SetAxisColor(0,0.0);
709  axisPlot->Draw("AXIS");
710  axisPlot->Draw("AXIG SAME");
711 
712  axisPlot->SetFillStyle(4000);
713  axisPad->SetFillStyle(4000);
714 
715  axisPad->SetGrid(1,1);
716  gPad->Modified(); gPad->Update();
717  gPad->SetFrameFillStyle(4000);
718 
719  gPad->Modified(); gPad->Update();
720  TFrame *frame = gPad->GetFrame();
721  frame->SetLineColorAlpha(0, 0.0);
722 
723  ridgeCanv->SaveAs("RidgePlots.pdf");
724  }
725  blankCanv->SaveAs("RidgePlots.pdf]");
726 }
#define MACH3LOG_TRACE
Definition: MaCh3Logger.h:33

◆ MakeTGraphAsymmErrors()

std::unique_ptr<TGraphAsymmErrors> MakeTGraphAsymmErrors ( const std::shared_ptr< TFile > &  File,
std::vector< int >  Index = {} 
)

Definition at line 795 of file GetPostfitParamPlots.cpp.

795  {})
796 {
797  int GraphBins = Index.size();
798  std::vector<double> x(GraphBins);
799  std::vector<double> y(GraphBins);
800  std::vector<double> exl(GraphBins);
801  std::vector<double> eyl(GraphBins);
802  std::vector<double> exh(GraphBins);
803  std::vector<double> eyh(GraphBins);
804 
805  TH1D* PostHist = static_cast<TH1D*>(File->Get( ("param_xsec_"+plotType).c_str() ));
806 
807  auto Errors_HPD_Positive = static_cast<TVectorD*>(File->Get( "Errors_HPD_Positive" ));
808  auto Errors_HPD_Negative = static_cast<TVectorD*>(File->Get( "Errors_HPD_Negative" ));
809  //KS: I am tempted to multithread this...
810  for(int i = 0; i < GraphBins; ++i)
811  {
812  int Counter = Index.size() == 0 ? i : Index[i];
813  //KS: We are extracting value from three object each having different numbering scheme, I have checked carefully so this is correct please don't change all these +1 +0.5 etc. it just work...
814  x[i] = i + 0.5;
815  y[i] = PostHist->GetBinContent(Counter+1);
816 
817  //KS: We don't want x axis errors as they are confusing in Violin plot
818  exh[i] = 0.00001;
819  exl[i] = 0.00001;
820  eyh[i] = (*Errors_HPD_Positive)(Counter);
821  eyl[i] = (*Errors_HPD_Negative)(Counter);
822  }
823  auto PostGraph = std::make_unique<TGraphAsymmErrors>(GraphBins, x.data(), y.data(), exl.data(), exh.data(), eyl.data(), eyh.data());
824  PostGraph->SetTitle("");
825 
826  return PostGraph;
827 }

◆ PrettifyTitles()

template<typename HistType >
void PrettifyTitles ( HistType *  hist)

Definition at line 132 of file GetPostfitParamPlots.cpp.

132  {
133  int nBins = hist->GetXaxis()->GetNbins();
134  for(int i = 0; i < nBins; ++i) {
135  std::string title = hist->GetXaxis()->GetBinLabel(i+1);
136  title = PlotMan->style().prettifyParamName(title);
137  hist->GetXaxis()->SetBinLabel(i+1, title.c_str());
138  }
139 }

◆ ReadSettings()

bool ReadSettings ( const std::shared_ptr< TFile > &  File1)

Definition at line 141 of file GetPostfitParamPlots.cpp.

142 {
143  MACH3LOG_DEBUG("Reading settings for file {}", File1->GetName());
144  #ifdef MACH3_DEBUG
145  File1->ls();
146  #endif
147  TTree *Settings = (File1->Get<TTree>("Settings"));
148 
149  // can't find settings tree :(
150  if (!Settings) return false;
151 
152  MACH3LOG_DEBUG("Got settings tree");
153  M3::Utils::Print(Settings);
154 
155  Settings->SetBranchAddress("NDParameters", &NDParameters);
156  Settings->SetBranchAddress("NDParametersStartingPos", &NDParametersStartingPos);
157 
158  std::vector<int> *NDSamples_Bins = 0;
159  std::vector<std::string> *NDSamples_Names = 0;
160  Settings->SetBranchAddress("NDSamplesNames", &NDSamples_Names);
161  Settings->SetBranchAddress("NDSamplesBins", &NDSamples_Bins);
162 
163  Settings->GetEntry(0);
164 
165  NDSamplesNames = *NDSamples_Names;
166  NDSamplesBins = *NDSamples_Bins;
167 
168  MACH3LOG_DEBUG("Read settings tree successfully");
169  return true;
170 }
void Print(const TTree *tree)
Definition: Monitor.cpp:326

Variable Documentation

◆ NDParameters

int NDParameters

Definition at line 53 of file GetPostfitParamPlots.cpp.

◆ NDParametersStartingPos

int NDParametersStartingPos

Definition at line 54 of file GetPostfitParamPlots.cpp.

◆ NDSamplesBins

std::vector<int> NDSamplesBins

Definition at line 56 of file GetPostfitParamPlots.cpp.

◆ NDSamplesNames

std::vector<std::string> NDSamplesNames

Definition at line 57 of file GetPostfitParamPlots.cpp.

◆ PlotMan

M3::Plotting::PlottingManager* PlotMan = nullptr

Definition at line 51 of file GetPostfitParamPlots.cpp.

◆ plotType

std::string plotType

Definition at line 60 of file GetPostfitParamPlots.cpp.

◆ SaveName

std::string SaveName

Definition at line 59 of file GetPostfitParamPlots.cpp.