MaCh3  2.4.2
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 const double ScalingFactor = 10
 
constexpr Color_t Colours [NVars] = {kRed, kGreen+1, kBlack, kBlue+1, kOrange+1}
 
constexpr ELineStyle Style [NVars] = {kDotted, kDashed, kSolid, kDashDotted, kDashDotted}
 
MaCh3Plotting::PlottingManager * PlotMan
 

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 376 of file PlotSigmaVariation.cpp.

377 {
378  //Get input file, make canvas and output file
379  auto canvas = std::make_unique<TCanvas>("canv", "canv", 1080, 1080);
380  TFile *infile = M3::Open(filename, "OPEN", __FILE__, __LINE__);
381  TDirectoryFile *SigmaDir = infile->Get<TDirectoryFile>("SigmaVar");
382 
383  std::string outfilename = filename.substr(0, filename.find(".root"));
384  outfilename = outfilename + "_RatioPlots1d.pdf";
385  gErrorIgnoreLevel = kWarning;
386  canvas->Print((outfilename+"[").c_str());
387 
388  auto IncludeString = GetFromManager<std::vector<std::string>>(Settings["IncludeString"], {});
389  auto ExcludeString = GetFromManager<std::vector<std::string>>(Settings["ExcludeString"], {});
390  TDirectory *dir = nullptr;
391 
392  for(size_t id = 0; id < DialNameVector.size(); id++)
393  {
394  for(size_t is = 0; is < SampleNameVector.size(); is++)
395  {
396  if(SkipDirectory(ExcludeString, IncludeString, (DialNameVector[id] + "/" + SampleNameVector[is]).c_str())) continue;
397  MACH3LOG_INFO("{} Entering {}/{}", __func__, DialNameVector[id], SampleNameVector[is]);
398  SigmaDir->cd((DialNameVector[id] + "/" + SampleNameVector[is]).c_str());
399 
400  //set dir to current directory
401  dir = gDirectory;
402 
403  // Loop over dimensions
404  for(int iDim = 0; iDim <= SampleMaxDim[is]; iDim++)
405  {
406  MACH3LOG_DEBUG("Starting loop over dimension {}", iDim);
407  //make -3,-1,0,1,3 polys
408  std::vector<std::unique_ptr<TH1D>> Projection;
409  TIter nextsub(dir->GetListOfKeys());
410  TKey *subsubkey = nullptr;
411 
412  //loop over items in directory, hard code which th2poly we want
413  while ((subsubkey = static_cast<TKey*>(nextsub())))
414  {
415  auto name = std::string(subsubkey->GetName());
416  auto classname = std::string(subsubkey->GetClassName());
417  // Looking
418  const std::string ProjectionName = "_1DProj" + std::to_string(iDim);
419  const bool IsProjection = (name.find(ProjectionName) != std::string::npos);
420  if (classname == "TH1D" && IsProjection)
421  {
422  name = DialNameVector[id] + "/" + SampleNameVector[is] + "/" + name;
423  Projection.emplace_back(M3::Clone(SigmaDir->Get<TH1D>(name.c_str())));
424  MACH3LOG_DEBUG("Adding hist {}", name);
425  }
426  }
427  std::string Title = PlotMan->style().prettifyParamName(DialNameVector[id]) + " "
428  + PlotMan->style().prettifySampleName(SampleNameVector[is]);
429  PlotRatio(Projection, canvas, Title, outfilename);
430  gDirectory->cd("..");
431  }
432  }
433  }
434 
435  canvas->Print((outfilename+"]").c_str());
436  infile->Close();
437  delete infile;
438 }
#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)
MaCh3Plotting::PlottingManager * PlotMan
std::vector< int > SampleMaxDim
std::vector< std::string > DialNameVector
std::vector< std::string > SampleNameVector
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 491 of file PlotSigmaVariation.cpp.

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

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

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

◆ InitializePads()

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

Definition at line 223 of file PlotSigmaVariation.cpp.

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

◆ main()

int main ( int  argc,
char **  argv 
)

Definition at line 824 of file PlotSigmaVariation.cpp.

825 {
827  if (argc != 3)
828  {
829  MACH3LOG_ERROR("Need two inputs: output of sigma var and config");
830  throw MaCh3Exception(__FILE__, __LINE__);
831  }
832  std::string filename = argv[1];
833  std::string ConfigName = argv[2];
834  MACH3LOG_INFO("Running {} with {} {}", argv[0], filename, ConfigName);
835  // Load the YAML file
836  YAML::Node Config = M3OpenConfig(ConfigName);
837 
838  // Access the "MatrixPlotter" section
839  YAML::Node settings = Config["PlotSigmaVariation"];
840 
842 
843  PlotMan = new MaCh3Plotting::PlottingManager();
844  PlotMan->initialise();
845 
846  CompareSigVar1D(filename, settings);
847  CompareSigVar2D(filename, settings);
848  MakeEventRatePlot(filename, settings);
849  OverlaySigVar1D(filename, settings);
850 
851  if(PlotMan) delete PlotMan;
852  return 0;
853 }
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 643 of file PlotSigmaVariation.cpp.

644 {
645  (void) Settings;
646  //Get input file, make canvas and output file
647  auto canvas = std::make_unique<TCanvas>("canv", "canv", 1080, 1080);
648  TFile *infile = M3::Open(filename, "OPEN", __FILE__, __LINE__);
649  TDirectoryFile *SigmaDir = infile->Get<TDirectoryFile>("SigmaVar");
650 
651  std::string outfilename = filename.substr(0, filename.find(".root"));
652  outfilename = outfilename + "_EventRate.pdf";
653  gErrorIgnoreLevel = kWarning;
654  canvas->Print((outfilename+"[").c_str());
655 
656  TDirectory *dir = nullptr;
657  for(size_t id = 0; id < DialNameVector.size(); id++)
658  {
659  std::vector<std::vector<std::unique_ptr<TH1D>>> Projection(SampleNameVector.size());
660  for(size_t is = 0; is < SampleNameVector.size(); is++)
661  {
662  MACH3LOG_INFO("{} Entering {}/{}", __func__, DialNameVector[id], SampleNameVector[is]);
663  SigmaDir->cd((DialNameVector[id] + "/" + SampleNameVector[is]).c_str());
664 
665  //set dir to current directory
666  dir = gDirectory;
667 
668  //make -3,-1,0,1,3 polys
669  TIter nextsub(dir->GetListOfKeys());
670  TKey *subsubkey = nullptr;
671 
672  //loop over items in directory, hard code which th2poly we want
673  while ((subsubkey = static_cast<TKey*>(nextsub())))
674  {
675  auto name = std::string(subsubkey->GetName());
676  auto classname = std::string(subsubkey->GetClassName());
677  // Looking
678  const std::string ProjectionName = "_1DProj0";
679  const bool IsProjection = (name.find(ProjectionName) != std::string::npos);
680  if (classname == "TH1D" && IsProjection)
681  {
682  name = DialNameVector[id] + "/" + SampleNameVector[is] + "/" + name;
683  Projection[is].emplace_back(M3::Clone(SigmaDir->Get<TH1D>(name.c_str())));
684  MACH3LOG_DEBUG("Adding hist {}", name);
685  }
686  }
687  gDirectory->cd("..");
688  }
689  std::string Title = PlotMan->style().prettifyParamName(DialNameVector[id]);
690  PlotEventRate(Projection, canvas, Title, outfilename);
691  }
692 
693  canvas->Print((outfilename+"]").c_str());
694  infile->Close();
695  delete infile;
696 }
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 79 of file PlotSigmaVariation.cpp.

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

◆ MakeRatio()

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

Definition at line 254 of file PlotSigmaVariation.cpp.

255  {
256  size_t ratio_index = 0; // Track position in Ratio vector
257  for (int i = 0; i < static_cast<int>(Poly.size()); ++i) {
258  if (i == PriorKnot) continue; // Skip PriorKnot
259  Ratio[ratio_index] = M3::Clone(Poly[i].get());
260  Ratio[ratio_index]->Divide(Poly[PriorKnot].get());
261  ratio_index++;
262  }
263 
264  Ratio[0]->GetYaxis()->SetTitle("Ratio to Prior");
265  Ratio[0]->SetBit(TH1D::kNoTitle);
266  Ratio[0]->SetBit(TH1D::kNoTitle);
267 
268  double maxz = -999;
269  double minz = +999;
270  for (int j = 0; j < static_cast<int>(sigmaArray.size())-1; j++)
271  {
272  for (int i = 1; i < Ratio[0]->GetXaxis()->GetNbins(); i++)
273  {
274  maxz = std::max(maxz, Ratio[j]->GetBinContent(i));
275  minz = std::min(minz, Ratio[j]->GetBinContent(i));
276  }
277  }
278  maxz = maxz*1.001;
279  minz = minz*1.001;
280 
281  if (std::fabs(1 - maxz) > std::fabs(1-minz))
282  Ratio[0]->GetYaxis()->SetRangeUser(1-std::fabs(1-maxz),1+std::fabs(1-maxz));
283  else
284  Ratio[0]->GetYaxis()->SetRangeUser(1-std::fabs(1-minz),1+std::fabs(1-minz));
285 }

◆ OverlaySigVar1D()

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

Definition at line 756 of file PlotSigmaVariation.cpp.

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

564 {
565  std::vector<std::unique_ptr<TH1D>> EvenRates(sigmaArray.size());
566  for(int ih = 0; ih < static_cast<int>(sigmaArray.size()); ih++)
567  {
568  EvenRates[ih] = std::make_unique<TH1D>(Title.c_str(), Title.c_str(), SampleNameVector.size(), 0, SampleNameVector.size());
569  EvenRates[ih]->SetDirectory(nullptr);
570  EvenRates[ih]->GetYaxis()->SetTitleOffset(1.4);
571  EvenRates[ih]->GetYaxis()->SetTitle("Events");
572  EvenRates[ih]->SetLineWidth(2.);
573  EvenRates[ih]->SetLineStyle(Style[ih]);
574  EvenRates[ih]->SetLineColor(Colours[ih]);
575 
576  for(size_t iSample = 0; iSample < SampleNameVector.size(); iSample ++) {
577  EvenRates[ih]->SetBinContent(iSample+1, Poly[iSample][ih]->Integral());
578 
579  std::string SamName = PlotMan->style().prettifySampleName(SampleNameVector[iSample]);
580  EvenRates[ih]->GetXaxis()->SetBinLabel(iSample+1, SamName.c_str());
581  }
582  }
583 
584  TPad* pad1 = nullptr;
585  TPad* pad2 = nullptr;
586  InitializePads(canv.get(), pad1, pad2, 0.50, 0.50);
587  pad2->SetBottomMargin(0.60);
588 
589  pad1->cd();
590  double max = 0;
591  for(int ik = 0; ik < static_cast<int>(sigmaArray.size()); ++ik) {
592  max = std::max(max, EvenRates[ik]->GetMaximum());
593  }
594  EvenRates[0]->SetTitle(Title.c_str());
595  EvenRates[0]->SetMaximum(max*1.2);
596  EvenRates[0]->Draw("HIST");
597  for(int ik = 1; ik < static_cast<int>(sigmaArray.size()); ++ik)
598  {
599  EvenRates[ik]->Draw("HIST SAME");
600  }
601 
602  auto DialValues = GetDialValues(Poly[0]);
603  auto leg = MakeLegend(0.55, 0.55, 0.8, 0.88, 0.04);
604  for (int j = 0; j < static_cast<int>(sigmaArray.size()); j++)
605  {
606  if(j == PriorKnot) {
607  leg->AddEntry(EvenRates[j].get(), Form("Prior (%.2f)", DialValues[j]), "l");
608  } else {
609  leg->AddEntry(EvenRates[j].get(), Form("%.0f#sigma (%.2f)", sigmaArray[j], DialValues[j]), "l");
610  }
611  }
612  leg->Draw("SAME");
613 
614  pad2->cd();
615  TLine line(EvenRates[0]->GetXaxis()->GetBinLowEdge(EvenRates[0]->GetXaxis()->GetFirst()),
616  1.0, EvenRates[0]->GetXaxis()->GetBinUpEdge(EvenRates[0]->GetXaxis()->GetLast()), 1.0);
617  std::vector<std::unique_ptr<TH1D>> Ratio(sigmaArray.size()-1);
618  MakeRatio(EvenRates, Ratio);
619  Ratio[0]->GetXaxis()->SetTitleSize(0.08);
620  Ratio[0]->GetYaxis()->SetTitleOffset(0.4);
621  Ratio[0]->GetYaxis()->SetTitleSize(0.06);
622 
623  Ratio[0]->GetXaxis()->SetLabelSize(0.08);
624  Ratio[0]->GetYaxis()->SetLabelSize(0.04);
625  Ratio[0]->GetXaxis()->LabelsOption("v");
626  Ratio[0]->Draw("HIST");
627 
628  for(int ik = 1; ik < static_cast<int>(sigmaArray.size())-1; ++ik)
629  {
630  Ratio[ik]->Draw("HIST SAME");
631  }
632 
633  line.SetLineWidth(2);
634  line.SetLineColor(kBlack);
635  line.Draw("SAME");
636 
637  canv->Print((outfilename).c_str());
638 
639  delete pad1;
640  delete pad2;
641 }
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 288 of file PlotSigmaVariation.cpp.

292 {
293  canv->Clear();
294  gStyle->SetDrawBorder(0);
295  gStyle->SetTitleBorderSize(2);
296  gStyle->SetOptStat(0); //Set 0 to disable statistic box
297  canv->SetGrid();
298  canv->SetTopMargin(0.10);
299  canv->SetBottomMargin(0.08);
300  canv->SetRightMargin(0.05);
301  canv->SetLeftMargin(0.12);
302 
303  TPad* pad1 = nullptr;
304  TPad* pad2 = nullptr;
305  InitializePads(canv.get(), pad1, pad2);
306 
307  pad1->cd();
308 
309  auto DialValues = GetDialValues(Poly);
310  double max = 0;
311  for(int ik = 0; ik < static_cast<int>(sigmaArray.size()); ++ik)
312  {
313  Poly[ik]->SetLineWidth(2.);
314  Poly[ik]->SetLineColor(Colours[ik]);
315  Poly[ik]->SetLineStyle(Style[ik]);
316  Poly[ik]->GetYaxis()->SetTitle(fmt::format("Events/{:.0f}", ScalingFactor).c_str());
317  M3::ScaleHistogram(Poly[ik].get(), ScalingFactor);
318  max = std::max(max, Poly[ik]->GetMaximum());
319  }
320  Poly[0]->SetTitle(Title.c_str());
321  Poly[0]->SetMaximum(max*1.2);
322  Poly[0]->Draw("HIST");
323  for(int ik = 1; ik < static_cast<int>(sigmaArray.size()); ++ik)
324  {
325  Poly[ik]->Draw("HIST SAME");
326  }
327 
328  std::vector<double> Integral(sigmaArray.size());
329  for(int ik = 0; ik < static_cast<int>(sigmaArray.size()); ++ik)
330  Integral[ik] = Poly[ik]->Integral();
331 
332  auto leg = MakeLegend(0.55, 0.55, 0.8, 0.88, 0.04);
333  leg->SetTextSize(0.04);
334  for (int j = 0; j < static_cast<int>(sigmaArray.size()); j++)
335  {
336  if(j == PriorKnot) {
337  leg->AddEntry(Poly[j].get(), Form("Prior (%.2f), #int=%.2f", DialValues[j], Integral[j]), "l");
338  } else {
339  leg->AddEntry(Poly[j].get(), Form("%.0f#sigma (%.2f), #int=%.2f", sigmaArray[j], DialValues[j], Integral[j]), "l");
340  }
341  }
342  leg->Draw("SAME");
343 
344  pad2->cd();
345 
346  TLine line(Poly[0]->GetXaxis()->GetBinLowEdge(Poly[0]->GetXaxis()->GetFirst()),
347  1.0, Poly[0]->GetXaxis()->GetBinUpEdge(Poly[0]->GetXaxis()->GetLast()), 1.0);
348  std::vector<std::unique_ptr<TH1D>> Ratio(sigmaArray.size()-1);
349 
350  MakeRatio(Poly, Ratio);
351 
352  auto PrettyX = PlotMan->style().prettifyKinematicName(Ratio[0]->GetXaxis()->GetTitle());
353  Ratio[0]->GetXaxis()->SetTitle(PrettyX.c_str());
354  Ratio[0]->GetXaxis()->SetTitleSize(0.12);
355  Ratio[0]->GetYaxis()->SetTitleOffset(0.4);
356  Ratio[0]->GetYaxis()->SetTitleSize(0.10);
357 
358  Ratio[0]->GetXaxis()->SetLabelSize(0.10);
359  Ratio[0]->GetYaxis()->SetLabelSize(0.10);
360  Ratio[0]->Draw("HIST");
361  for(int ik = 1; ik < static_cast<int>(sigmaArray.size())-1; ++ik)
362  {
363  Ratio[ik]->Draw("HIST SAME");
364  }
365 
366  line.SetLineWidth(2);
367  line.SetLineColor(kBlack);
368  line.Draw("SAME");
369 
370  canv->Print((outfilename).c_str());
371 
372  delete pad1;
373  delete pad2;
374 }
constexpr const double ScalingFactor
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 440 of file PlotSigmaVariation.cpp.

444 {
445  canv->Clear();
446  gStyle->SetDrawBorder(0);
447  gStyle->SetTitleBorderSize(2);
448  gStyle->SetOptStat(0); //Set 0 to disable statistic box
449  canv->SetGrid();
450  canv->SetTopMargin(0.10);
451  canv->SetBottomMargin(0.10);
452  canv->SetLeftMargin(0.12);
453  canv->SetRightMargin(0.20);
454 
455  constexpr int NRGBs = 5;
456  TColor::InitializeColors();
457  Double_t stops[NRGBs] = { 0.00, 0.25, 0.50, 0.75, 1.00 };
458  Double_t red[NRGBs] = { 0.00, 0.25, 1.00, 1.00, 0.50 };
459  Double_t green[NRGBs] = { 0.00, 0.25, 1.00, 0.25, 0.00 };
460  Double_t blue[NRGBs] = { 0.50, 1.00, 1.00, 0.25, 0.00 };
461  TColor::CreateGradientColorTable(5, stops, red, green, blue, 255);
462  gStyle->SetNumberContours(255);
463 
464  for (int i = 0; i < static_cast<int>(Poly.size()); ++i) {
465  if (i == PriorKnot) continue; // Skip PriorKnot
466  std::unique_ptr<TH2> Ratio = M3::Clone(Poly[i].get());
467  Ratio->Divide(Poly[PriorKnot].get());
468  Ratio->SetTitle((Title + " " + std::to_string(static_cast<int>(sigmaArray[i])) + "sigma").c_str());
469 
470  const double maxz = Ratio->GetMaximum();
471  const double minz = Ratio->GetMinimum();
472  if (std::fabs(1-maxz) > std::fabs(1-minz))
473  Ratio->GetZaxis()->SetRangeUser(1-std::fabs(1-maxz),1+std::fabs(1-maxz));
474  else
475  Ratio->GetZaxis()->SetRangeUser(1-std::fabs(1-minz),1+std::fabs(1-minz));
476  Ratio->GetXaxis()->SetTitleOffset(1.1);
477  Ratio->GetYaxis()->SetTitleOffset(1.1);
478  Ratio->GetZaxis()->SetTitleOffset(1.5);
479  Ratio->GetZaxis()->SetTitle("Ratio to Prior");
480 
481  auto PrettyX = PlotMan->style().prettifyKinematicName(Ratio->GetXaxis()->GetTitle());
482  Ratio->GetXaxis()->SetTitle(PrettyX.c_str());
483  auto PrettyY = PlotMan->style().prettifyKinematicName(Ratio->GetYaxis()->GetTitle());
484  Ratio->GetYaxis()->SetTitle(PrettyY.c_str());
485 
486  Ratio->Draw("COLZ");
487  canv->Print((outfilename).c_str());
488  }
489 }

◆ 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 698 of file PlotSigmaVariation.cpp.

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

◆ 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 94 of file PlotSigmaVariation.cpp.

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

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

Variable Documentation

◆ Colours

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

Definition at line 23 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

MaCh3Plotting::PlottingManager* PlotMan
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 27 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.

◆ ScalingFactor

constexpr const double ScalingFactor = 10
constexpr

Definition at line 21 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 24 of file PlotSigmaVariation.cpp.