MaCh3  2.4.2
Reference Guide
HistogramUtils.h
Go to the documentation of this file.
1 #pragma once
2 
3 // MaCh3 includes
6 
8 // ROOT include
9 #include "TGraphAsymmErrors.h"
10 #include "TObjString.h"
11 #include "TRandom3.h"
13 
17 
19 double OverflowIntegral(TH2Poly* poly);
20 
22 double NoOverflowIntegral(TH2Poly* poly);
23 
25 TH1D* PolyProjectionX(TObject* poly, std::string TempName, const std::vector<double>& xbins, const bool computeErrors = false);
27 TH1D* PolyProjectionY(TObject* poly, std::string TempName, const std::vector<double>& ybins, const bool computeErrors = false);
28 
30 TH2D* ConvertTH2PolyToTH2D(TH2Poly *poly, TH2D *TH2Dhist);
32 TH2Poly* ConvertTH2DtoTH2Poly(TH2D *TH2Dhist);
33 
35 TH2Poly* NormalisePoly(TH2Poly* Histogram);
36 
39 void NormaliseTH2Poly(TH2Poly* Histogram);
40 
42 TH2Poly* PolyScaleWidth(TH2Poly *Histogram, double scale);
43 
45 double PolyIntegralWidth(TH2Poly *Histogram);
46 
48 template<class HistType> HistType* RatioHists(HistType* NumHist, HistType* DenomHist);
49 
51 TH2Poly* RatioPolys(TH2Poly* NumPoly, TH2Poly* DenomPoly);
52 
57 TH2Poly* MakePolyHist(const std::string& name, const std::vector<double>& BinArray_x, const std::vector<double>& BinArray_y);
58 
62 void CheckTH2PolyFileVersion(TFile *file);
63 
65 void RemoveFitter(TH1D* hist, const std::string& name);
66 
68 void MakeFluctuatedHistogramStandard(TH1D *FluctHist, TH1D* PolyHist, TRandom3* rand);
70 void MakeFluctuatedHistogramAlternative(TH1D *FluctHist, TH1D* PolyHist, TRandom3* rand);
71 
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 
89 std::unique_ptr<TGraphAsymmErrors> MakeAsymGraph(TH1* sigmaArrayLeft, TH1* sigmaArrayCentr, TH1* sigmaArrayRight, const std::string& title);
90 
94 void FastViolinFill(TH2D* violin, TH1D* hist_1d);
95 
101 template <typename Derived, typename Base>
102 std::vector<Base*> CastVector(const std::vector<Derived*>& inputVec) {
103  std::vector<Base*> outputVec;
104  // Reserve space for efficiency
105  outputVec.reserve(inputVec.size());
106  for (auto* ptr : inputVec) {
107  outputVec.push_back(static_cast<Base*>(ptr));
108  }
109  return outputVec;
110 }
111 
113 inline std::string file_exists(std::string filename) {
114  std::ifstream infile(filename.c_str());
115  if (!infile.good()) {
116  MACH3LOG_ERROR("*** ERROR ***");
117  MACH3LOG_ERROR("File {} does not exist", filename);
118  MACH3LOG_ERROR("Please try again");
119  MACH3LOG_ERROR("*************");
120  throw MaCh3Exception(__FILE__ , __LINE__ );
121  }
122  return filename;
123 }
124 
126 double returnCherenkovThresholdMomentum(const int PDG);
127 
129 double CalculateQ2(double PLep, double PUpd, double EnuTrue, double InitialQ2 = 0.0);
130 
132 double CalculateEnu(double PLep, double cosTheta, double EB, bool neutrino);
133 
134 namespace M3 {
141 template <typename ObjectType>
142 std::unique_ptr<ObjectType> Clone(const ObjectType* obj, const std::string& name = "") {
143  std::string cloneName = name.empty() ? obj->GetName() : name;
144 
145  std::unique_ptr<ObjectType> Hist(static_cast<ObjectType*>(obj->Clone(cloneName.c_str())));
146  // Disable ROOT memory management because it causes lot of headache especially as smart pointers are much smarter
147  Hist->SetDirectory(nullptr);
148 
149  return Hist;
150 }
151 
152 
163 TFile* Open(const std::string& Name, const std::string& Type, const std::string& File, const int Line);
164 
166 void ScaleHistogram(TH1* Sample_Hist, const double scale);
167 
168 } //end M3
#define _MaCh3_Safe_Include_Start_
KS: Avoiding warning checking for headers.
Definition: Core.h:126
#define _MaCh3_Safe_Include_End_
KS: Restore warning checking after including external headers.
Definition: Core.h:140
std::unique_ptr< TGraphAsymmErrors > MakeAsymGraph(TH1 *sigmaArrayLeft, TH1 *sigmaArrayCentr, TH1 *sigmaArrayRight, const std::string &title)
Used by sigma variation, check how 1 sigma changes spectra.
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...
std::vector< Base * > CastVector(const std::vector< Derived * > &inputVec)
Converts a vector of pointers from a derived type to a base type.
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".
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.
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.
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 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.