MaCh3  2.5.0
Reference Guide
HistogramUtils.h
Go to the documentation of this file.
1 #pragma once
2 
3 // MaCh3 includes
6 #include "Manager/YamlHelper.h"
7 
8 //KS: Joy of forward declaration https://gieseanw.wordpress.com/2018/02/25/the-joys-of-forward-declarations-results-from-the-real-world/
9 class TRandom3;
10 class TObject;
11 
15 
17 double OverflowIntegral(TH2Poly* poly);
18 
20 double NoOverflowIntegral(TH2Poly* poly);
21 
23 TH1D* PolyProjectionX(TObject* poly, std::string TempName, const std::vector<double>& xbins, const bool computeErrors = false);
25 TH1D* PolyProjectionY(TObject* poly, std::string TempName, const std::vector<double>& ybins, const bool computeErrors = false);
26 
28 TH2D* ConvertTH2PolyToTH2D(TH2Poly *poly, TH2D *TH2Dhist);
30 TH2Poly* ConvertTH2DtoTH2Poly(TH2D *TH2Dhist);
31 
33 TH2Poly* NormalisePoly(TH2Poly* Histogram);
34 
37 void NormaliseTH2Poly(TH2Poly* Histogram);
38 
40 TH2Poly* PolyScaleWidth(TH2Poly *Histogram, double scale);
41 
43 double PolyIntegralWidth(TH2Poly *Histogram);
44 
46 template<class HistType> HistType* RatioHists(HistType* NumHist, HistType* DenomHist);
47 
49 TH2Poly* RatioPolys(TH2Poly* NumPoly, TH2Poly* DenomPoly);
50 
55 TH2Poly* MakePolyHist(const std::string& name, const std::vector<double>& BinArray_x, const std::vector<double>& BinArray_y);
56 
60 void CheckTH2PolyFileVersion(TFile *file);
61 
63 void RemoveFitter(TH1D* hist, const std::string& name);
64 
66 void MakeFluctuatedHistogramStandard(TH1D *FluctHist, TH1D* PolyHist, TRandom3* rand);
68 void MakeFluctuatedHistogramAlternative(TH1D *FluctHist, TH1D* PolyHist, TRandom3* rand);
69 
71 void MakeFluctuatedHistogramStandard(TH2D *FluctHist, TH2D* PolyHist, TRandom3* rand);
73 void MakeFluctuatedHistogramAlternative(TH2D *FluctHist, TH2D* PolyHist, TRandom3* rand);
74 
76 void MakeFluctuatedHistogramStandard(TH2Poly *FluctHist, TH2Poly* PolyHist, TRandom3* rand);
78 void MakeFluctuatedHistogramAlternative(TH2Poly *FluctHist, TH2Poly* PolyHist, TRandom3* rand);
79 
81 int GetRandomPoly2(const TH2Poly* PolyHist, TRandom3* rand);
82 
86 void FastViolinFill(TH2D* violin, TH1D* hist_1d);
87 
93 template <typename Derived, typename Base>
94 std::vector<std::unique_ptr<Base>>
95 CastVector(std::vector<std::unique_ptr<Derived>>&& inputVec) {
96  std::vector<std::unique_ptr<Base>> outputVec;
97  outputVec.reserve(inputVec.size());
98  for (auto& ptr : inputVec) {
99  outputVec.push_back(std::unique_ptr<Base>(std::move(ptr)));
100  }
101  return outputVec;
102 }
103 
105 inline std::string file_exists(std::string filename) {
106  std::ifstream infile(filename.c_str());
107  if (!infile.good()) {
108  MACH3LOG_ERROR("File {} does not exist", filename);
109  MACH3LOG_ERROR("Please try again");
110  throw MaCh3Exception(__FILE__ , __LINE__ );
111  }
112  return filename;
113 }
114 
117 double returnCherenkovThresholdMomentum(const int PDG);
118 
124 double CalculateQ2(double PLep, double PUpd, double EnuTrue, double InitialQ2 = 0.0);
125 
131 double CalculateEnu(double PLep, double cosTheta, double EB, bool neutrino);
132 
136 std::unique_ptr<TH1D> MakeSummaryFromSpectra(const TH2D* Spectra, const std::string& name);
137 
138 namespace M3 {
145 template <typename ObjectType>
146 std::unique_ptr<ObjectType> Clone(const ObjectType* obj, const std::string& name = "") {
147  std::string cloneName = name.empty() ? obj->GetName() : name;
148 
149  std::unique_ptr<ObjectType> Hist(static_cast<ObjectType*>(obj->Clone(cloneName.c_str())));
150  // Disable ROOT memory management because it causes lot of headache especially as smart pointers are much smarter
151  Hist->SetDirectory(nullptr);
152 
153  return Hist;
154 }
155 
166 TFile* Open(const std::string& Name, const std::string& Type, const std::string& File, const int Line);
167 
169 void ScaleHistogram(TH1* Sample_Hist, const double scale);
170 
172 void CheckBinningMatch(const TH1D* Hist1, const TH1D* Hist2, const std::string& File, const int Line);
174 void CheckBinningMatch(const TH2D* Hist1, const TH2D* Hist2, const std::string& File, const int Line);
176 void CheckBinningMatch(TH2Poly* Hist1, TH2Poly* Hist2, const std::string& File, const int Line);
178 YAML::Node PolyToYaml(TH2Poly* Hist, const std::string& YamlName, const std::string& File, const int Line);
179 
180 } //end M3
double CalculateQ2(double PLep, double PUpd, double EnuTrue, double InitialQ2=0.0)
Recalculate Q^2 after Eb shift. Takes in shifted lepton momentum, lepton angle, and true neutrino ene...
void CheckTH2PolyFileVersion(TFile *file)
KS: ROOT changes something with binning when moving from ROOT 5 to ROOT 6. If you open ROOT5 produced...
void MakeFluctuatedHistogramAlternative(TH1D *FluctHist, TH1D *PolyHist, TRandom3 *rand)
Make Poisson fluctuation of TH1D hist using slow method which is only for cross-check.
TH2Poly * RatioPolys(TH2Poly *NumPoly, TH2Poly *DenomPoly)
Helper to make ratio of TH2Polys.
double CalculateEnu(double PLep, double cosTheta, double EB, bool neutrino)
Recalculate Enu after Eb shift. Takes in shifted lepton momentum, lepton angle, and binding energy ch...
void NormaliseTH2Poly(TH2Poly *Histogram)
Helper to Normalise histograms.
double returnCherenkovThresholdMomentum(const int PDG)
DB Get the Cherenkov momentum threshold in MeV.
double OverflowIntegral(TH2Poly *poly)
WP: Helper function for calculating unbinned Integral of TH2Poly i.e including overflow.
void MakeFluctuatedHistogramStandard(TH1D *FluctHist, TH1D *PolyHist, TRandom3 *rand)
Make Poisson fluctuation of TH1D hist using default fast method.
TH2D * ConvertTH2PolyToTH2D(TH2Poly *poly, TH2D *TH2Dhist)
KS: Convert TH2D to TH2Poly.
int GetRandomPoly2(const TH2Poly *PolyHist, TRandom3 *rand)
KS: ROOT developers were too lazy do develop getRanom2 for TH2Poly, this implementation is based on l...
void FastViolinFill(TH2D *violin, TH1D *hist_1d)
KS: Fill Violin histogram with entry from a toy.
TH2Poly * NormalisePoly(TH2Poly *Histogram)
WP: Helper to Normalise histograms.
HistType * RatioHists(HistType *NumHist, HistType *DenomHist)
Helper to make ratio histograms.
TH2Poly * PolyScaleWidth(TH2Poly *Histogram, double scale)
WP: Helper to scale th2poly analogous to th2d scale with option "width".
void RemoveFitter(TH1D *hist, const std::string &name)
KS: Remove fitted TF1 from hist to make comparison easier.
double PolyIntegralWidth(TH2Poly *Histogram)
WP: Helper to calc integral of th2poly analogous to th2d integra; with option "width".
std::vector< std::unique_ptr< Base > > CastVector(std::vector< std::unique_ptr< Derived >> &&inputVec)
Converts a vector of pointers from a derived type to a base type.
TH2Poly * ConvertTH2DtoTH2Poly(TH2D *TH2Dhist)
KS: Convert TH2Poly to TH2D.
double NoOverflowIntegral(TH2Poly *poly)
WP: Helper function for calculating binned Integral of TH2Poly i.e not including overflow.
TH2Poly * MakePolyHist(const std::string &name, const std::vector< double > &BinArray_x, const std::vector< double > &BinArray_y)
WP: Helper function to create TH2Poly histogram with uniform binning.
TH1D * PolyProjectionY(TObject *poly, std::string TempName, const std::vector< double > &ybins, const bool computeErrors=false)
WP: Poly Projectors.
std::unique_ptr< TH1D > MakeSummaryFromSpectra(const TH2D *Spectra, const std::string &name)
Build a 1D posterior-predictive summary from a violin spectrum.
TH1D * PolyProjectionX(TObject *poly, std::string TempName, const std::vector< double > &xbins, const bool computeErrors=false)
WP: Poly Projectors.
std::string file_exists(std::string filename)
Helper to check if files exist or not.
#define MACH3LOG_ERROR
Definition: MaCh3Logger.h:37
Definitions of generic parameter structs and utility templates for MaCh3.
Utility functions for handling YAML nodes.
Custom exception class used throughout MaCh3.
Main namespace for MaCh3 software.
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.
void CheckBinningMatch(const TH1D *Hist1, const TH1D *Hist2, const std::string &File, const int Line)
KS: Helper function check if data and MC binning matches.
void ScaleHistogram(TH1 *Sample_Hist, const double scale)
Scale histogram to get divided by bin width.
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.
YAML::Node PolyToYaml(TH2Poly *Hist, const std::string &YamlName, const std::string &File, const int Line)
KS: Convert TH2Poly into yaml config accepted by MaCh3.