MaCh3  2.5.1
Reference Guide
Functions | Variables
PlotSigmaVariation.cpp File Reference
#include "PlottingUtils/PlottingUtils.h"
#include "PlottingUtils/PlottingManager.h"
Include dependency graph for PlotSigmaVariation.cpp:

Go to the source code of this file.

Functions

void FindKnot (std::vector< double > &SigmaValues, const std::string &dirname, const std::string &subdirname, const std::string &ProjName, std::string histname)
 Histograms have name like ND_CC0pi_1DProj0_Norm_Param_0_sig_n3.00_val_0.25. This code is trying to extract sigma names. More...
 
std::unique_ptr< TLegend > MakeLegend (double x1, double y1, double x2, double y2, double textSize=0.04)
 
void ScanInput (std::vector< std::string > &DialNameVecr, std::vector< std::string > &SampleNameVec, std::vector< int > &SampleDimVec, std::vector< double > &SigmaValues, const std::string &filename)
 Scan inputs to figure out dial name and used sample names. More...
 
bool SkipDirectory (const std::vector< std::string > &ExcludeString, const std::vector< std::string > &IncludeString, const std::string &dirname)
 Check whether to skip directory or not based on defined strings. More...
 
std::vector< double > GetDialValues (const std::vector< std::unique_ptr< TH1D >> &Poly)
 Extracts dial value for from histogram title. More...
 
void InitializePads (TCanvas *canv, TPad *&pad1, TPad *&pad2, double Pad1Bottom=0.25, double Pad2Top=0.25)
 
void MakeRatio (const std::vector< std::unique_ptr< TH1D >> &Poly, std::vector< std::unique_ptr< TH1D >> &Ratio)
 
void PlotRatio (const std::vector< std::unique_ptr< TH1D >> &Poly, const std::unique_ptr< TCanvas > &canv, const std::string &Title, const std::string &outfilename)
 
void CompareSigVar1D (const std::string &filename, const YAML::Node &Settings)
 
void PlotRatio2D (const std::vector< std::unique_ptr< TH2 >> &Poly, const std::unique_ptr< TCanvas > &canv, const std::string &Title, const std::string &outfilename)
 
void CompareSigVar2D (const std::string &filename, const YAML::Node &Settings)
 
void PlotEventRate (const std::vector< std::vector< std::unique_ptr< TH1D >>> &Poly, const std::unique_ptr< TCanvas > &canv, const std::string &Title, const std::string &outfilename)
 
void MakeEventRatePlot (const std::string &filename, const YAML::Node &Settings)
 
void PlotSigVar1D (const std::vector< std::vector< std::unique_ptr< TH1D >>> &Projection, const std::unique_ptr< TCanvas > &canv, const std::string &Title, const std::string &outfilename, const std::vector< std::string > &ParamNames, const std::vector< int > &ParamColour)
 
void OverlaySigVar1D (const std::string &filename, const YAML::Node &Settings)
 
int main (int argc, char **argv)
 

Variables

std::vector< std::string > DialNameVector
 
std::vector< std::string > SampleNameVector
 
std::vector< int > SampleMaxDim
 
std::vector< double > sigmaArray
 
int PriorKnot = M3::_BAD_INT_
 
constexpr const int NVars = 5
 
constexpr Color_t Colours [NVars] = {kRed, kGreen+1, kBlack, kBlue+1, kOrange+1}
 
constexpr ELineStyle Style [NVars] = {kDotted, kDashed, kSolid, kDashDotted, kDashDotted}
 
M3::Plotting::PlottingManager * PlotMan = nullptr
 

Detailed Description

Author
Kamil Skwarczynski
Todo:
add maybe ratio for PlotSigVar1D

Definition in file PlotSigmaVariation.cpp.

Function Documentation

◆ CompareSigVar1D()

void CompareSigVar1D ( const std::string &  filename,
const YAML::Node &  Settings 
)

Definition at line 358 of file PlotSigmaVariation.cpp.

359 {
360  //Get input file, make canvas and output file
361  auto canvas = std::make_unique<TCanvas>("canv", "canv", 1080, 1080);
362  TFile *infile = M3::Open(filename, "OPEN", __FILE__, __LINE__);
363  TDirectoryFile *SigmaDir = infile->Get<TDirectoryFile>("SigmaVar");
364 
365  std::string outfilename = filename.substr(0, filename.find(".root"));
366  outfilename = outfilename + "_RatioPlots1d.pdf";
367  gErrorIgnoreLevel = kWarning;
368  canvas->Print((outfilename+"[").c_str());
369 
370  auto IncludeString = GetFromManager<std::vector<std::string>>(Settings["IncludeString"], {});
371  auto ExcludeString = GetFromManager<std::vector<std::string>>(Settings["ExcludeString"], {});
372  TDirectory *dir = nullptr;
373 
374  for(size_t id = 0; id < DialNameVector.size(); id++)
375  {
376  for(size_t is = 0; is < SampleNameVector.size(); is++)
377  {
378  if(SkipDirectory(ExcludeString, IncludeString, (DialNameVector[id] + "/" + SampleNameVector[is]).c_str())) continue;
379  MACH3LOG_INFO("{} Entering {}/{}", __func__, DialNameVector[id], SampleNameVector[is]);
380  SigmaDir->cd((DialNameVector[id] + "/" + SampleNameVector[is]).c_str());
381 
382  //set dir to current directory
383  dir = gDirectory;
384 
385  // Loop over dimensions
386  for(int iDim = 0; iDim <= SampleMaxDim[is]; iDim++)
387  {
388  MACH3LOG_DEBUG("Starting loop over dimension {}", iDim);
389  //make -3,-1,0,1,3 polys
390  std::vector<std::unique_ptr<TH1D>> Projection;
391  TIter nextsub(dir->GetListOfKeys());
392  TKey *subsubkey = nullptr;
393 
394  //loop over items in directory, hard code which th2poly we want
395  while ((subsubkey = static_cast<TKey*>(nextsub())))
396  {
397  auto name = std::string(subsubkey->GetName());
398  auto classname = std::string(subsubkey->GetClassName());
399  // Looking
400  const std::string ProjectionName = "_1DProj" + std::to_string(iDim);
401  const bool IsProjection = (name.find(ProjectionName) != std::string::npos);
402  if (classname == "TH1D" && IsProjection)
403  {
404  name = DialNameVector[id] + "/" + SampleNameVector[is] + "/" + name;
405  Projection.emplace_back(M3::Clone(SigmaDir->Get<TH1D>(name.c_str())));
406  MACH3LOG_DEBUG("Adding hist {}", name);
407  }
408  }
409  std::string Title = PlotMan->style().prettifyParamName(DialNameVector[id]) + " "
410  + PlotMan->style().prettifySampleName(SampleNameVector[is]);
411  PlotRatio(Projection, canvas, Title, outfilename);
412  gDirectory->cd("..");
413  }
414  }
415  }
416 
417  canvas->Print((outfilename+"]").c_str());
418  infile->Close();
419  delete infile;
420 }
#define MACH3LOG_DEBUG
Definition: MaCh3Logger.h:34
#define MACH3LOG_INFO
Definition: MaCh3Logger.h:35
void PlotRatio(const std::vector< std::unique_ptr< TH1D >> &Poly, const std::unique_ptr< TCanvas > &canv, const std::string &Title, const std::string &outfilename)
std::vector< int > SampleMaxDim
std::vector< std::string > DialNameVector
std::vector< std::string > SampleNameVector
M3::Plotting::PlottingManager * PlotMan
bool SkipDirectory(const std::vector< std::string > &ExcludeString, const std::vector< std::string > &IncludeString, const std::string &dirname)
Check whether to skip directory or not based on defined strings.
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.
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.

◆ CompareSigVar2D()

void CompareSigVar2D ( const std::string &  filename,
const YAML::Node &  Settings 
)

Definition at line 473 of file PlotSigmaVariation.cpp.

474 {
475  //Get input file, make canvas and output file
476  auto canvas = std::make_unique<TCanvas>("canv", "canv", 1080, 1080);
477  TFile *infile = M3::Open(filename, "OPEN", __FILE__, __LINE__);
478  TDirectoryFile *SigmaDir = infile->Get<TDirectoryFile>("SigmaVar");
479 
480  std::string outfilename = filename.substr(0, filename.find(".root"));
481  outfilename = outfilename + "_RatioPlots2d.pdf";
482  gErrorIgnoreLevel = kWarning;
483  canvas->Print((outfilename+"[").c_str());
484 
485  auto IncludeString = GetFromManager<std::vector<std::string>>(Settings["IncludeString"], {});
486  auto ExcludeString = GetFromManager<std::vector<std::string>>(Settings["ExcludeString"], {});
487  TDirectory *dir = nullptr;
488 
489  for(size_t id = 0; id < DialNameVector.size(); id++)
490  {
491  for(size_t is = 0; is < SampleNameVector.size(); is++)
492  {
493  if(SkipDirectory(ExcludeString, IncludeString, (DialNameVector[id] + "/" + SampleNameVector[is]).c_str())) continue;
494  MACH3LOG_INFO("{} Entering {}/{}", __func__, DialNameVector[id], SampleNameVector[is]);
495  SigmaDir->cd((DialNameVector[id] + "/" + SampleNameVector[is]).c_str());
496 
497  //set dir to current directory
498  dir = gDirectory;
499 
500  int nDim = SampleMaxDim[is];
501 
502  TIter nextsub(dir->GetListOfKeys());
503  TKey *subsubkey = nullptr;
504 
505  // Loop over all unique dimension pairs
506  for (int iDim1 = 0; iDim1 <= nDim; ++iDim1) {
507  for (int iDim2 = iDim1 + 1; iDim2 <= nDim; ++iDim2) {
508  // Reset iterator for each dimension pair
509  nextsub.Reset();
510  //make -3,-1,0,1,3 polys
511  std::vector<std::unique_ptr<TH2>> Projection;
512 
513  //loop over items in directory, hard code which th2poly we want
514  while ((subsubkey = static_cast<TKey*>(nextsub())))
515  {
516  auto name = std::string(subsubkey->GetName());
517  auto classname = std::string(subsubkey->GetClassName());
518  // Looking
519  const std::string ProjectionName = "_2DProj_" + std::to_string(iDim1) + "_vs_" + std::to_string(iDim2);
520  const bool IsProjection = (name.find(ProjectionName) != std::string::npos);
521  if ((classname == "TH2D" || classname == "TH2Poly")&& IsProjection)
522  {
523  name = DialNameVector[id] + "/" + SampleNameVector[is] + "/" + name;
524  Projection.emplace_back(M3::Clone(SigmaDir->Get<TH2>(name.c_str())));
525  MACH3LOG_DEBUG("Adding hist {}", name);
526  }
527  }
528  std::string Title = PlotMan->style().prettifyParamName(DialNameVector[id]) + " "
529  + PlotMan->style().prettifySampleName(SampleNameVector[is]);
530  if(Projection.size() == sigmaArray.size()) PlotRatio2D(Projection, canvas, Title, outfilename);
531  }
532  }
533  gDirectory->cd("..");
534  }
535  }
536  canvas->Print((outfilename+"]").c_str());
537  infile->Close();
538  delete infile;
539 }
void PlotRatio2D(const std::vector< std::unique_ptr< TH2 >> &Poly, const std::unique_ptr< TCanvas > &canv, const std::string &Title, const std::string &outfilename)
std::vector< double > sigmaArray

◆ FindKnot()

void FindKnot ( std::vector< double > &  SigmaValues,
const std::string &  dirname,
const std::string &  subdirname,
const std::string &  ProjName,
std::string  histname 
)

Histograms have name like ND_CC0pi_1DProj0_Norm_Param_0_sig_n3.00_val_0.25. This code is trying to extract sigma names.

Definition at line 28 of file PlotSigmaVariation.cpp.

32  {
33  auto StripPrefix = [](std::string& str, const std::string& prefix) {
34  if (str.find(prefix + "_") == 0) {
35  str.erase(0, prefix.length() + 1);
36  if (str.find(prefix) == 0) {
37  MACH3LOG_ERROR("Failed to strip prefix '{}' from string '{}'", prefix, str);
38  throw MaCh3Exception(__FILE__, __LINE__);
39  }
40  } else {
41  MACH3LOG_ERROR("String '{}' does not start with expected prefix '{}'", str, prefix);
42  throw MaCh3Exception(__FILE__, __LINE__);
43  }
44  };
45 
46  // Remove sample and dial name to avoid potential issues
47  StripPrefix(histname, subdirname);
48  StripPrefix(histname, ProjName);
49  StripPrefix(histname, dirname);
50 
51  MACH3LOG_DEBUG("Name afters striping {}", histname);
52 
53  double sigma = 0.0;
54  // Find the "_sig_" part in the name
55  size_t sig_pos = histname.find("_sig_");
56 
57  // Extract the part after "_sig_"
58  std::string sigma_part = histname.substr(sig_pos + 5);
59 
60  // Check if it starts with 'n' (negative) or 'p' (positive) or 'nom' (0 or prior)
61  if (histname.find("nom_") != std::string::npos) {
62  sigma = 0.0;
63  PriorKnot = static_cast<int>(SigmaValues.size());
64  MACH3LOG_DEBUG("Found prior knot {}", PriorKnot);
65  } else if (sigma_part.size() > 0 && sigma_part[0] == 'n') {
66  sigma = -std::stod(sigma_part.substr(1));
67  } else if (sigma_part.size() > 0 && sigma_part[0] == 'p') {
68  sigma = std::stod(sigma_part.substr(1));
69  } else {
70  MACH3LOG_ERROR("Weirdly formatted string {}", sigma_part);
71  throw MaCh3Exception(__FILE__ , __LINE__ );
72  }
73  MACH3LOG_DEBUG("Adding sigma {}", sigma);
74 
75  SigmaValues.push_back(sigma);
76 }
#define MACH3LOG_ERROR
Definition: MaCh3Logger.h:37
int PriorKnot
Custom exception class used throughout MaCh3.

◆ GetDialValues()

std::vector<double> GetDialValues ( const std::vector< std::unique_ptr< TH1D >> &  Poly)

Extracts dial value for from histogram title.

Definition at line 203 of file PlotSigmaVariation.cpp.

203  {
204  std::vector<double> values;
205  for (const auto& hist : Poly) {
206  std::string title = hist->GetTitle();
207  auto pos = title.rfind("_val_");
208  if (pos != std::string::npos) {
209  std::string val_str = title.substr(pos + 5); // skip "_val_"
210  double val = std::stod(val_str);
211  values.push_back(val);
212  MACH3LOG_DEBUG("Extracted dial value {} from title '{}'", val, title);
213  } else {
214  MACH3LOG_DEBUG("Failed to extract dial value from title '{}'", title);
215  throw MaCh3Exception(__FILE__, __LINE__);
216  }
217  }
218  return values;
219 }

◆ InitializePads()

void InitializePads ( TCanvas *  canv,
TPad *&  pad1,
TPad *&  pad2,
double  Pad1Bottom = 0.25,
double  Pad2Top = 0.25 
)

Definition at line 222 of file PlotSigmaVariation.cpp.

223 {
224  // Delete existing pads if they exist
225  if (pad1) delete pad1;
226  if (pad2) delete pad2;
227 
228  // Allocate new pads
229  pad1 = new TPad("pad1", "pad1", 0., Pad2Top, 1., 1.0);
230  pad2 = new TPad("pad2", "pad2", 0., 0., 1., Pad1Bottom);
231 
232  // Append pads to canvas
233  pad1->AppendPad();
234  pad2->AppendPad();
235 
236  // Set margins for pad1
237  pad1->SetLeftMargin(canv->GetLeftMargin());
238  pad1->SetRightMargin(canv->GetRightMargin());
239  pad1->SetTopMargin(canv->GetTopMargin());
240  pad1->SetBottomMargin(0);
241 
242  // Set margins for pad2
243  pad2->SetLeftMargin(canv->GetLeftMargin());
244  pad2->SetRightMargin(canv->GetRightMargin());
245  pad2->SetTopMargin(0);
246  pad2->SetBottomMargin(0.30);
247 
248  // Enable grid for both pads
249  pad1->SetGrid();
250  pad2->SetGrid();
251 }

◆ main()

int main ( int  argc,
char **  argv 
)

Definition at line 809 of file PlotSigmaVariation.cpp.

810 {
812  if (argc != 3)
813  {
814  MACH3LOG_ERROR("Need two inputs: output of sigma var and config");
815  throw MaCh3Exception(__FILE__, __LINE__);
816  }
817  std::string filename = argv[1];
818  std::string ConfigName = argv[2];
819  MACH3LOG_INFO("Running {} with {} {}", argv[0], filename, ConfigName);
820  // Load the YAML file
821  YAML::Node Config = M3OpenConfig(ConfigName);
822 
823  // Access the "MatrixPlotter" section
824  YAML::Node settings = Config["PlotSigmaVariation"];
825 
827 
828  PlotMan = new M3::Plotting::PlottingManager();
829  PlotMan->initialise();
830 
831  CompareSigVar1D(filename, settings);
832  CompareSigVar2D(filename, settings);
833  MakeEventRatePlot(filename, settings);
834  OverlaySigVar1D(filename, settings);
835 
836  if(PlotMan) delete PlotMan;
837  return 0;
838 }
void SetMaCh3LoggerFormat()
Set messaging format of the logger.
Definition: MaCh3Logger.h:61
void CompareSigVar2D(const std::string &filename, const YAML::Node &Settings)
void OverlaySigVar1D(const std::string &filename, const YAML::Node &Settings)
void CompareSigVar1D(const std::string &filename, const YAML::Node &Settings)
void ScanInput(std::vector< std::string > &DialNameVecr, std::vector< std::string > &SampleNameVec, std::vector< int > &SampleDimVec, std::vector< double > &SigmaValues, const std::string &filename)
Scan inputs to figure out dial name and used sample names.
void MakeEventRatePlot(const std::string &filename, const YAML::Node &Settings)
#define M3OpenConfig(filename)
Macro to simplify calling LoadYaml with file and line info.
Definition: YamlHelper.h:589

◆ MakeEventRatePlot()

void MakeEventRatePlot ( const std::string &  filename,
const YAML::Node &  Settings 
)

Definition at line 625 of file PlotSigmaVariation.cpp.

626 {
627  (void) Settings;
628  //Get input file, make canvas and output file
629  auto canvas = std::make_unique<TCanvas>("canv", "canv", 1080, 1080);
630  canvas->SetRightMargin(0.05);
631  TFile *infile = M3::Open(filename, "OPEN", __FILE__, __LINE__);
632  TDirectoryFile *SigmaDir = infile->Get<TDirectoryFile>("SigmaVar");
633 
634  std::string outfilename = filename.substr(0, filename.find(".root"));
635  outfilename = outfilename + "_EventRate.pdf";
636  gErrorIgnoreLevel = kWarning;
637  canvas->Print((outfilename+"[").c_str());
638 
639  TDirectory *dir = nullptr;
640  for(size_t id = 0; id < DialNameVector.size(); id++)
641  {
642  std::vector<std::vector<std::unique_ptr<TH1D>>> Projection(SampleNameVector.size());
643  for(size_t is = 0; is < SampleNameVector.size(); is++)
644  {
645  MACH3LOG_INFO("{} Entering {}/{}", __func__, DialNameVector[id], SampleNameVector[is]);
646  SigmaDir->cd((DialNameVector[id] + "/" + SampleNameVector[is]).c_str());
647 
648  //set dir to current directory
649  dir = gDirectory;
650 
651  //make -3,-1,0,1,3 polys
652  TIter nextsub(dir->GetListOfKeys());
653  TKey *subsubkey = nullptr;
654 
655  //loop over items in directory, hard code which th2poly we want
656  while ((subsubkey = static_cast<TKey*>(nextsub())))
657  {
658  auto name = std::string(subsubkey->GetName());
659  auto classname = std::string(subsubkey->GetClassName());
660  // Looking
661  const std::string ProjectionName = "_1DProj0";
662  const bool IsProjection = (name.find(ProjectionName) != std::string::npos);
663  if (classname == "TH1D" && IsProjection)
664  {
665  name = DialNameVector[id] + "/" + SampleNameVector[is] + "/" + name;
666  Projection[is].emplace_back(M3::Clone(SigmaDir->Get<TH1D>(name.c_str())));
667  MACH3LOG_DEBUG("Adding hist {}", name);
668  }
669  }
670  gDirectory->cd("..");
671  }
672  std::string Title = PlotMan->style().prettifyParamName(DialNameVector[id]);
673  PlotEventRate(Projection, canvas, Title, outfilename);
674  }
675 
676  canvas->Print((outfilename+"]").c_str());
677  infile->Close();
678  delete infile;
679 }
void PlotEventRate(const std::vector< std::vector< std::unique_ptr< TH1D >>> &Poly, const std::unique_ptr< TCanvas > &canv, const std::string &Title, const std::string &outfilename)

◆ MakeLegend()

std::unique_ptr<TLegend> MakeLegend ( double  x1,
double  y1,
double  x2,
double  y2,
double  textSize = 0.04 
)

Definition at line 78 of file PlotSigmaVariation.cpp.

80 {
81  auto leg = std::make_unique<TLegend>(x1, y1, x2, y2);
82  leg->SetTextSize(textSize);
83  leg->SetLineColor(0);
84  leg->SetLineStyle(0);
85  leg->SetFillColor(0);
86  leg->SetFillStyle(0);
87  leg->SetBorderSize(0);
88 
89  return leg;
90 }

◆ MakeRatio()

void MakeRatio ( const std::vector< std::unique_ptr< TH1D >> &  Poly,
std::vector< std::unique_ptr< TH1D >> &  Ratio 
)

Definition at line 253 of file PlotSigmaVariation.cpp.

254  {
255  size_t ratio_index = 0; // Track position in Ratio vector
256  for (int i = 0; i < static_cast<int>(Poly.size()); ++i) {
257  if (i == PriorKnot) continue; // Skip PriorKnot
258  Ratio[ratio_index] = M3::Clone(Poly[i].get());
259  Ratio[ratio_index]->Divide(Poly[PriorKnot].get());
260  ratio_index++;
261  }
262 
263  Ratio[0]->GetYaxis()->SetTitle("Ratio to Prior");
264  Ratio[0]->SetBit(TH1D::kNoTitle);
265 
266  M3::Plotting::SetSymmetricRatioRange(Ratio);
267 }

◆ OverlaySigVar1D()

void OverlaySigVar1D ( const std::string &  filename,
const YAML::Node &  Settings 
)

Definition at line 741 of file PlotSigmaVariation.cpp.

742 {
743  auto ParamNames = GetFromManager<std::vector<std::string>>(Settings["ParamNames"], {});
744  auto SigColours = GetFromManager<std::vector<int>>(Settings["Coulour"], {});
745 
746  if(ParamNames.size() == 0) return;
747 
748  if(ParamNames.size() != SigColours.size()){
749  MACH3LOG_ERROR("Wrong Size");
750  throw MaCh3Exception(__FILE__, __LINE__);
751  }
752 
753  //Get input file, make canvas and output file
754  auto canvas = std::make_unique<TCanvas>("canv", "canv", 1080, 1080);
755  TFile *infile = M3::Open(filename, "OPEN", __FILE__, __LINE__);
756  TDirectoryFile *SigmaDir = infile->Get<TDirectoryFile>("SigmaVar");
757 
758  std::string outfilename = filename.substr(0, filename.find(".root"));
759  outfilename = outfilename + "_Overlay1d.pdf";
760  gErrorIgnoreLevel = kWarning;
761  canvas->Print((outfilename+"[").c_str());
762 
763  TDirectory *dir = nullptr;
764 
765  for(size_t is = 0; is < SampleNameVector.size(); is++)
766  {
767  // Loop over dimensions
768  for(int iDim = 0; iDim <= SampleMaxDim[is]; iDim++)
769  {
770  std::vector<std::vector<std::unique_ptr<TH1D>>> Projection(ParamNames.size());
771  for(size_t id = 0; id < ParamNames.size(); id++)
772  {
773  MACH3LOG_INFO("{} Entering {}/{}", __func__, ParamNames[id], SampleNameVector[is]);
774  SigmaDir->cd((ParamNames[id] + "/" + SampleNameVector[is]).c_str());
775 
776  //set dir to current directory
777  dir = gDirectory;
778  //make -3,-1,0,1,3 polys
779  TIter nextsub(dir->GetListOfKeys());
780  TKey *subsubkey = nullptr;
781 
782  //loop over items in directory, hard code which th2poly we want
783  while ((subsubkey = static_cast<TKey*>(nextsub())))
784  {
785  auto name = std::string(subsubkey->GetName());
786  auto classname = std::string(subsubkey->GetClassName());
787  // Looking
788  const std::string ProjectionName = "_1DProj" + std::to_string(iDim);
789  const bool IsProjection = (name.find(ProjectionName) != std::string::npos);
790  if (classname == "TH1D" && IsProjection)
791  {
792  name = ParamNames[id] + "/" + SampleNameVector[is] + "/" + name;
793  Projection[id].emplace_back(M3::Clone(SigmaDir->Get<TH1D>(name.c_str())));
794  MACH3LOG_DEBUG("Adding hist {}", name);
795  }
796  }
797  }
798  std::string Title = PlotMan->style().prettifySampleName(SampleNameVector[is]);
799  PlotSigVar1D(Projection, canvas, Title, outfilename, ParamNames, SigColours);
800  gDirectory->cd("..");
801  }
802  }
803 
804  canvas->Print((outfilename+"]").c_str());
805  infile->Close();
806  delete infile;
807 }
void PlotSigVar1D(const std::vector< std::vector< std::unique_ptr< TH1D >>> &Projection, const std::unique_ptr< TCanvas > &canv, const std::string &Title, const std::string &outfilename, const std::vector< std::string > &ParamNames, const std::vector< int > &ParamColour)

◆ PlotEventRate()

void PlotEventRate ( const std::vector< std::vector< std::unique_ptr< TH1D >>> &  Poly,
const std::unique_ptr< TCanvas > &  canv,
const std::string &  Title,
const std::string &  outfilename 
)

Definition at line 542 of file PlotSigmaVariation.cpp.

546 {
547  std::vector<std::unique_ptr<TH1D>> EvenRates(sigmaArray.size());
548  for(int ih = 0; ih < static_cast<int>(sigmaArray.size()); ih++)
549  {
550  EvenRates[ih] = std::make_unique<TH1D>(Title.c_str(), Title.c_str(), SampleNameVector.size(), 0, SampleNameVector.size());
551  EvenRates[ih]->SetDirectory(nullptr);
552  EvenRates[ih]->GetYaxis()->SetTitleOffset(1.4);
553  EvenRates[ih]->GetYaxis()->SetTitle("Events");
554  EvenRates[ih]->SetLineWidth(2.);
555  EvenRates[ih]->SetLineStyle(Style[ih]);
556  EvenRates[ih]->SetLineColor(Colours[ih]);
557 
558  for(size_t iSample = 0; iSample < SampleNameVector.size(); iSample ++) {
559  EvenRates[ih]->SetBinContent(iSample+1, Poly[iSample][ih]->Integral());
560 
561  std::string SamName = PlotMan->style().prettifySampleName(SampleNameVector[iSample]);
562  EvenRates[ih]->GetXaxis()->SetBinLabel(iSample+1, SamName.c_str());
563  }
564  }
565 
566  TPad* pad1 = nullptr;
567  TPad* pad2 = nullptr;
568  InitializePads(canv.get(), pad1, pad2, 0.50, 0.50);
569  pad2->SetBottomMargin(0.60);
570 
571  pad1->cd();
572  double max = 0;
573  for(int ik = 0; ik < static_cast<int>(sigmaArray.size()); ++ik) {
574  max = std::max(max, EvenRates[ik]->GetMaximum());
575  }
576  EvenRates[0]->SetTitle(Title.c_str());
577  EvenRates[0]->SetMaximum(max*1.2);
578  EvenRates[0]->Draw("HIST");
579  for(int ik = 1; ik < static_cast<int>(sigmaArray.size()); ++ik)
580  {
581  EvenRates[ik]->Draw("HIST SAME");
582  }
583 
584  auto DialValues = GetDialValues(Poly[0]);
585  auto leg = MakeLegend(0.55, 0.55, 0.8, 0.88, 0.04);
586  for (int j = 0; j < static_cast<int>(sigmaArray.size()); j++)
587  {
588  if(j == PriorKnot) {
589  leg->AddEntry(EvenRates[j].get(), Form("Prior (%.2f)", DialValues[j]), "l");
590  } else {
591  leg->AddEntry(EvenRates[j].get(), Form("%.0f#sigma (%.2f)", sigmaArray[j], DialValues[j]), "l");
592  }
593  }
594  leg->Draw("SAME");
595 
596  pad2->cd();
597  TLine line(EvenRates[0]->GetXaxis()->GetBinLowEdge(EvenRates[0]->GetXaxis()->GetFirst()),
598  1.0, EvenRates[0]->GetXaxis()->GetBinUpEdge(EvenRates[0]->GetXaxis()->GetLast()), 1.0);
599  std::vector<std::unique_ptr<TH1D>> Ratio(sigmaArray.size()-1);
600  MakeRatio(EvenRates, Ratio);
601  Ratio[0]->GetXaxis()->SetTitleSize(0.08);
602  Ratio[0]->GetYaxis()->SetTitleOffset(0.4);
603  Ratio[0]->GetYaxis()->SetTitleSize(0.06);
604 
605  Ratio[0]->GetXaxis()->SetLabelSize(0.08);
606  Ratio[0]->GetYaxis()->SetLabelSize(0.04);
607  Ratio[0]->GetXaxis()->LabelsOption("v");
608  Ratio[0]->Draw("HIST");
609 
610  for(int ik = 1; ik < static_cast<int>(sigmaArray.size())-1; ++ik)
611  {
612  Ratio[ik]->Draw("HIST SAME");
613  }
614 
615  line.SetLineWidth(2);
616  line.SetLineColor(kBlack);
617  line.Draw("SAME");
618 
619  canv->Print((outfilename).c_str());
620 
621  delete pad1;
622  delete pad2;
623 }
constexpr ELineStyle Style[NVars]
constexpr Color_t Colours[NVars]
std::vector< double > GetDialValues(const std::vector< std::unique_ptr< TH1D >> &Poly)
Extracts dial value for from histogram title.
std::unique_ptr< TLegend > MakeLegend(double x1, double y1, double x2, double y2, double textSize=0.04)
void InitializePads(TCanvas *canv, TPad *&pad1, TPad *&pad2, double Pad1Bottom=0.25, double Pad2Top=0.25)
void MakeRatio(const std::vector< std::unique_ptr< TH1D >> &Poly, std::vector< std::unique_ptr< TH1D >> &Ratio)

◆ PlotRatio()

void PlotRatio ( const std::vector< std::unique_ptr< TH1D >> &  Poly,
const std::unique_ptr< TCanvas > &  canv,
const std::string &  Title,
const std::string &  outfilename 
)

Definition at line 269 of file PlotSigmaVariation.cpp.

273 {
274  canv->Clear();
275  gStyle->SetDrawBorder(0);
276  gStyle->SetTitleBorderSize(2);
277  gStyle->SetOptStat(0); //Set 0 to disable statistic box
278  canv->SetGrid();
279  canv->SetTopMargin(0.10);
280  canv->SetBottomMargin(0.08);
281  canv->SetRightMargin(0.05);
282  canv->SetLeftMargin(0.12);
283 
284  TPad* pad1 = nullptr;
285  TPad* pad2 = nullptr;
286  InitializePads(canv.get(), pad1, pad2);
287 
288  pad1->cd();
289 
290  auto DialValues = GetDialValues(Poly);
291  double max = 0;
292  for(int ik = 0; ik < static_cast<int>(sigmaArray.size()); ++ik)
293  {
294  Poly[ik]->SetLineWidth(2.);
295  Poly[ik]->SetLineColor(Colours[ik]);
296  Poly[ik]->SetLineStyle(Style[ik]);
297  auto BinWidthScale = PlotMan->style().getBinWidthScale(Poly[ik]->GetXaxis()->GetTitle());
298  Poly[ik]->GetYaxis()->SetTitle(fmt::format("Events/{:.0f}", BinWidthScale).c_str());
299  M3::ScaleHistogram(Poly[ik].get(), BinWidthScale);
300  max = std::max(max, Poly[ik]->GetMaximum());
301  }
302  Poly[0]->SetTitle(Title.c_str());
303  Poly[0]->SetMaximum(max*1.2);
304  Poly[0]->Draw("HIST");
305  for(int ik = 1; ik < static_cast<int>(sigmaArray.size()); ++ik)
306  {
307  Poly[ik]->Draw("HIST SAME");
308  }
309 
310  std::vector<double> Integral(sigmaArray.size());
311  for(int ik = 0; ik < static_cast<int>(sigmaArray.size()); ++ik)
312  Integral[ik] = Poly[ik]->Integral();
313 
314  auto leg = MakeLegend(0.55, 0.55, 0.8, 0.88, 0.04);
315  leg->SetTextSize(0.04);
316  for (int j = 0; j < static_cast<int>(sigmaArray.size()); j++)
317  {
318  if(j == PriorKnot) {
319  leg->AddEntry(Poly[j].get(), Form("Prior (%.2f), #int=%.2f", DialValues[j], Integral[j]), "l");
320  } else {
321  leg->AddEntry(Poly[j].get(), Form("%.0f#sigma (%.2f), #int=%.2f", sigmaArray[j], DialValues[j], Integral[j]), "l");
322  }
323  }
324  leg->Draw("SAME");
325 
326  pad2->cd();
327 
328  TLine line(Poly[0]->GetXaxis()->GetBinLowEdge(Poly[0]->GetXaxis()->GetFirst()),
329  1.0, Poly[0]->GetXaxis()->GetBinUpEdge(Poly[0]->GetXaxis()->GetLast()), 1.0);
330  std::vector<std::unique_ptr<TH1D>> Ratio(sigmaArray.size()-1);
331 
332  MakeRatio(Poly, Ratio);
333 
334  auto PrettyX = PlotMan->style().prettifyKinematicName(Ratio[0]->GetXaxis()->GetTitle());
335  Ratio[0]->GetXaxis()->SetTitle(PrettyX.c_str());
336  Ratio[0]->GetXaxis()->SetTitleSize(0.12);
337  Ratio[0]->GetYaxis()->SetTitleOffset(0.4);
338  Ratio[0]->GetYaxis()->SetTitleSize(0.10);
339 
340  Ratio[0]->GetXaxis()->SetLabelSize(0.10);
341  Ratio[0]->GetYaxis()->SetLabelSize(0.10);
342  Ratio[0]->Draw("HIST");
343  for(int ik = 1; ik < static_cast<int>(sigmaArray.size())-1; ++ik)
344  {
345  Ratio[ik]->Draw("HIST SAME");
346  }
347 
348  line.SetLineWidth(2);
349  line.SetLineColor(kBlack);
350  line.Draw("SAME");
351 
352  canv->Print((outfilename).c_str());
353 
354  delete pad1;
355  delete pad2;
356 }
void ScaleHistogram(TH1 *Sample_Hist, const double scale)
Scale histogram to get divided by bin width.

◆ PlotRatio2D()

void PlotRatio2D ( const std::vector< std::unique_ptr< TH2 >> &  Poly,
const std::unique_ptr< TCanvas > &  canv,
const std::string &  Title,
const std::string &  outfilename 
)

Definition at line 422 of file PlotSigmaVariation.cpp.

426 {
427  canv->Clear();
428  gStyle->SetDrawBorder(0);
429  gStyle->SetTitleBorderSize(2);
430  gStyle->SetOptStat(0); //Set 0 to disable statistic box
431  canv->SetGrid();
432  canv->SetTopMargin(0.10);
433  canv->SetBottomMargin(0.10);
434  canv->SetLeftMargin(0.12);
435  canv->SetRightMargin(0.20);
436 
437  constexpr int NRGBs = 5;
438  TColor::InitializeColors();
439  Double_t stops[NRGBs] = { 0.00, 0.25, 0.50, 0.75, 1.00 };
440  Double_t red[NRGBs] = { 0.00, 0.25, 1.00, 1.00, 0.50 };
441  Double_t green[NRGBs] = { 0.00, 0.25, 1.00, 0.25, 0.00 };
442  Double_t blue[NRGBs] = { 0.50, 1.00, 1.00, 0.25, 0.00 };
443  TColor::CreateGradientColorTable(5, stops, red, green, blue, 255);
444  gStyle->SetNumberContours(255);
445 
446  for (int i = 0; i < static_cast<int>(Poly.size()); ++i) {
447  if (i == PriorKnot) continue; // Skip PriorKnot
448  std::unique_ptr<TH2> Ratio = M3::Clone(Poly[i].get());
449  Ratio->Divide(Poly[PriorKnot].get());
450  Ratio->SetTitle((Title + " " + std::to_string(static_cast<int>(sigmaArray[i])) + "sigma").c_str());
451 
452  const double maxz = Ratio->GetMaximum();
453  const double minz = Ratio->GetMinimum();
454  if (std::fabs(1-maxz) > std::fabs(1-minz))
455  Ratio->GetZaxis()->SetRangeUser(1-std::fabs(1-maxz),1+std::fabs(1-maxz));
456  else
457  Ratio->GetZaxis()->SetRangeUser(1-std::fabs(1-minz),1+std::fabs(1-minz));
458  Ratio->GetXaxis()->SetTitleOffset(1.1);
459  Ratio->GetYaxis()->SetTitleOffset(1.1);
460  Ratio->GetZaxis()->SetTitleOffset(1.5);
461  Ratio->GetZaxis()->SetTitle("Ratio to Prior");
462 
463  auto PrettyX = PlotMan->style().prettifyKinematicName(Ratio->GetXaxis()->GetTitle());
464  Ratio->GetXaxis()->SetTitle(PrettyX.c_str());
465  auto PrettyY = PlotMan->style().prettifyKinematicName(Ratio->GetYaxis()->GetTitle());
466  Ratio->GetYaxis()->SetTitle(PrettyY.c_str());
467 
468  Ratio->Draw("COLZ");
469  canv->Print((outfilename).c_str());
470  }
471 }

◆ PlotSigVar1D()

void PlotSigVar1D ( const std::vector< std::vector< std::unique_ptr< TH1D >>> &  Projection,
const std::unique_ptr< TCanvas > &  canv,
const std::string &  Title,
const std::string &  outfilename,
const std::vector< std::string > &  ParamNames,
const std::vector< int > &  ParamColour 
)

Definition at line 681 of file PlotSigmaVariation.cpp.

687 {
688  canv->Clear();
689  gStyle->SetDrawBorder(0);
690  gStyle->SetTitleBorderSize(2);
691  gStyle->SetOptStat(0); //Set 0 to disable statistic box
692  canv->SetGrid();
693  canv->SetTopMargin(0.10);
694  canv->SetBottomMargin(0.08);
695  canv->SetRightMargin(0.05);
696  canv->SetLeftMargin(0.12);
697 
698  auto PriorHist = Projection[0][PriorKnot].get();
699  PriorHist->SetTitle(Title.c_str());
700 
701  auto BinWidthScale = PlotMan->style().getBinWidthScale(PriorHist->GetXaxis()->GetTitle());
702  PriorHist->GetYaxis()->SetTitle(fmt::format("Events/{:.0f}", BinWidthScale).c_str());
703  PriorHist->Draw("HIST");
704  PriorHist->SetLineWidth(2.);
705  PriorHist->SetLineColor(kBlack);
706  M3::ScaleHistogram(PriorHist, BinWidthScale);
707 
708  auto PrettyX = PlotMan->style().prettifyKinematicName(PriorHist->GetXaxis()->GetTitle());
709  PriorHist->GetXaxis()->SetTitle(PrettyX.c_str());
710  for(int ik = 0; ik < static_cast<int>(sigmaArray.size()); ++ik)
711  {
712  if(ik == PriorKnot) continue;
713 
714  double max = 0;
715  for(size_t nParam = 0; nParam < Projection.size(); nParam++) {
716  Projection[nParam][ik]->SetLineWidth(2.);
717  Projection[nParam][ik]->SetLineColor(ParamColour[nParam]);
718  Projection[nParam][ik]->SetLineStyle(kDotted);
719  Projection[nParam][ik]->GetYaxis()->SetTitle(fmt::format("Events/{:.0f}", BinWidthScale).c_str());
720  M3::ScaleHistogram(Projection[nParam][ik].get(), BinWidthScale);
721  max = std::max(max, Projection[nParam][ik]->GetMaximum());
722  }
723  PriorHist->SetMaximum(max*1.2);
724 
725  PriorHist->Draw("HIST");
726  for(size_t nParam = 0; nParam < Projection.size(); nParam++) {
727  Projection[nParam][ik]->Draw("HIST SAME");
728  }
729  auto leg = MakeLegend(0.50, 0.55, 0.70, 0.75, 0.035);
730  leg->AddEntry(PriorHist, "Prior", "l");
731  for(size_t nParam = 0; nParam < Projection.size(); nParam++) {
732  leg->AddEntry(Projection[nParam][ik].get(), Form("%s (%.0f#sigma)", ParamNames[nParam].c_str(), sigmaArray[ik]), "l");
733  }
734  leg->Draw("SAME");
735 
736  canv->Print((outfilename).c_str());
737  canv->cd();
738  }
739 }

◆ ScanInput()

void ScanInput ( std::vector< std::string > &  DialNameVecr,
std::vector< std::string > &  SampleNameVec,
std::vector< int > &  SampleDimVec,
std::vector< double > &  SigmaValues,
const std::string &  filename 
)

Scan inputs to figure out dial name and used sample names.

Definition at line 93 of file PlotSigmaVariation.cpp.

98 {
99  MACH3LOG_DEBUG("Starting {}", __func__);
100  TFile *infile = M3::Open(filename, "OPEN", __FILE__, __LINE__);
101  TDirectoryFile *SigmaDir = infile->Get<TDirectoryFile>("SigmaVar");
102 
103  //Get all entries in input file
104  TIter next(SigmaDir->GetListOfKeys());
105  TKey *key = nullptr;
106 
107  // Loop over directory with dial names
108  while ((key = static_cast<TKey*>(next()))) {
109  // get directory names, ignore flux
110  auto classname = std::string(key->GetClassName());
111  auto dirname = std::string(key->GetName());
112 
113  if (classname != "TDirectoryFile") continue;
114  dirname = std::string(key->GetName());
115 
116  SigmaDir->cd(dirname.c_str());
117  TIter nextsub(gDirectory->GetListOfKeys());
118  TKey *subkey = nullptr;
119 
120  DialNameVecr.push_back(dirname);
121  MACH3LOG_DEBUG("Entering Dial {}", dirname);
122 
123  if(SampleNameVec.size() != 0) continue;
124  //loop over directories with sample names
125  while ((subkey = static_cast<TKey*>(nextsub())))
126  {
127  auto subdirname = std::string(subkey->GetName());
128  SampleNameVec.push_back(subdirname);
129  SampleDimVec.push_back(0);
130  MACH3LOG_DEBUG("Entering Sample {}", subdirname);
131  SigmaDir->cd((dirname + "/" + subdirname).c_str());
132 
133  TKey *subsubkey = nullptr;
134  TIter nextsubsub(gDirectory->GetListOfKeys());
135 
136  // Check if we already filled sigma vector
137  bool FillSigma = false;
138  if(SigmaValues.size() == 0) FillSigma = true;
139  // loop over histograms
140  while ((subsubkey = static_cast<TKey*>(nextsubsub())))
141  {
142  auto subsubdirname = std::string(subsubkey->GetTitle());
143  MACH3LOG_DEBUG("Entering Hist {}", subsubdirname);
144  std::string histname = subsubdirname;
145 
146  classname = std::string(subsubkey->GetClassName());
147 
148  if (classname != "TH1D") continue;
149  // Find if there is more dimensions
150  size_t proj_pos = histname.find("_1DProj");
151  int proj_number = -1;
152 
153  if (proj_pos != std::string::npos) {
154  size_t number_start = proj_pos + 7; // skip "_1DProj"
155  size_t number_end = histname.find_first_not_of("0123456789", number_start);
156  proj_number = std::stoi(histname.substr(number_start, number_end - number_start));
157  }
158  SampleDimVec.back() = std::max(proj_number, SampleDimVec.back());
159  MACH3LOG_DEBUG("Found dimension {} with dimension", proj_number);
160 
161  // KS: Extract knot position from hist only we haven't done this before and for projection X
162  // Sigma are same for projection Y and Z and beyond
163  if(FillSigma && proj_number == 0) FindKnot(SigmaValues, dirname, subdirname, "1DProj" + std::to_string(proj_number), histname);
164  }
165  }
166  }
167 
168  if(PriorKnot == M3::_BAD_INT_){
169  MACH3LOG_ERROR("Didn't find prior knot, something is not right...");
170  throw MaCh3Exception(__FILE__ , __LINE__ );
171  }
172 
173  if(sigmaArray.size() < NVars){
174  MACH3LOG_ERROR("Found sigma {}, while I have some hardcoding for {}",sigmaArray.size(), NVars);
175  throw MaCh3Exception(__FILE__ , __LINE__ );
176  }
177 
178  if(SampleNameVec.size() != SampleDimVec.size()) {
179  MACH3LOG_ERROR("Sample name vec ({}) and sample dimension vec ({}) have different sizes, something is not right");
180  throw MaCh3Exception(__FILE__ , __LINE__ );
181  }
182 
183  infile->Close();
184  delete infile;
185 }
void FindKnot(std::vector< double > &SigmaValues, const std::string &dirname, const std::string &subdirname, const std::string &ProjName, std::string histname)
Histograms have name like ND_CC0pi_1DProj0_Norm_Param_0_sig_n3.00_val_0.25. This code is trying to ex...
constexpr const int NVars
constexpr static const int _BAD_INT_
Default value used for int initialisation.
Definition: Core.h:55

◆ SkipDirectory()

bool SkipDirectory ( const std::vector< std::string > &  ExcludeString,
const std::vector< std::string > &  IncludeString,
const std::string &  dirname 
)

Check whether to skip directory or not based on defined strings.

Definition at line 188 of file PlotSigmaVariation.cpp.

189 {
190  bool Skip = false;
191  for(unsigned int i = 0; i < ExcludeString.size(); i++)
192  {
193  if (dirname.find(ExcludeString[i]) != std::string::npos){ Skip = true; break; }
194  }
195  for(unsigned int i = 0; i < IncludeString.size(); i++)
196  {
197  if (!(dirname.find(IncludeString[i]) != std::string::npos)){ Skip = true; break; }
198  }
199  return Skip;
200 }

Variable Documentation

◆ Colours

constexpr Color_t Colours[NVars] = {kRed, kGreen+1, kBlack, kBlue+1, kOrange+1}
constexpr

Definition at line 21 of file PlotSigmaVariation.cpp.

◆ DialNameVector

std::vector<std::string> DialNameVector

Definition at line 13 of file PlotSigmaVariation.cpp.

◆ NVars

constexpr const int NVars = 5
constexpr

Definition at line 20 of file PlotSigmaVariation.cpp.

◆ PlotMan

M3::Plotting::PlottingManager* PlotMan = nullptr
Warning
KS: keep raw pointer or ensure manual delete of PlotMan. If spdlog in automatically deleted before PlotMan then destructor has some spdlog and this could cause segfault

Definition at line 25 of file PlotSigmaVariation.cpp.

◆ PriorKnot

int PriorKnot = M3::_BAD_INT_

Definition at line 18 of file PlotSigmaVariation.cpp.

◆ SampleMaxDim

std::vector<int> SampleMaxDim

Definition at line 15 of file PlotSigmaVariation.cpp.

◆ SampleNameVector

std::vector<std::string> SampleNameVector

Definition at line 14 of file PlotSigmaVariation.cpp.

◆ sigmaArray

std::vector<double> sigmaArray

Definition at line 16 of file PlotSigmaVariation.cpp.

◆ Style

constexpr ELineStyle Style[NVars] = {kDotted, kDashed, kSolid, kDashDotted, kDashDotted}
constexpr

Definition at line 22 of file PlotSigmaVariation.cpp.