9#include "TGraphAsymmErrors.h"
10#include "TObjString.h"
25TH1D*
PolyProjectionX(TObject* poly, std::string TempName,
const std::vector<double>& xbins,
const bool computeErrors =
false);
27TH1D*
PolyProjectionY(TObject* poly, std::string TempName,
const std::vector<double>& ybins,
const bool computeErrors =
false);
48template<
class HistType> HistType*
RatioHists(HistType* NumHist, HistType* DenomHist);
51TH2Poly*
RatioPolys(TH2Poly* NumPoly, TH2Poly* DenomPoly);
57TH2Poly*
MakePolyHist(
const std::string& name,
const std::vector<double>& BinArray_x,
const std::vector<double>& BinArray_y);
86std::unique_ptr<TGraphAsymmErrors>
MakeAsymGraph(TH1* sigmaArrayLeft, TH1* sigmaArrayCentr, TH1* sigmaArrayRight,
const std::string& title);
98template <
typename Derived,
typename Base>
99std::vector<Base*>
CastVector(
const std::vector<Derived*>& inputVec) {
100 std::vector<Base*> outputVec;
102 outputVec.reserve(inputVec.size());
103 for (
auto* ptr : inputVec) {
104 outputVec.push_back(
static_cast<Base*
>(ptr));
111 std::ifstream infile(filename.c_str());
112 if (!infile.good()) {
126double CalculateQ2(
double PLep,
double PUpd,
double EnuTrue,
double InitialQ2 = 0.0);
129double CalculateEnu(
double PLep,
double cosTheta,
double EB,
bool neutrino);
137template <
typename ObjectType>
138std::unique_ptr<ObjectType>
Clone(
const ObjectType* obj,
const std::string& name =
"") {
139 std::string cloneName = name.empty() ? obj->GetName() : name;
141 std::unique_ptr<ObjectType> Hist(
static_cast<ObjectType*
>(obj->Clone(cloneName.c_str())));
143 Hist->SetDirectory(
nullptr);
157TFile*
Open(
const std::string& Name,
const std::string& Type,
const std::string& File,
const int Line);
#define _MaCh3_Safe_Include_Start_
KS: Avoiding warning checking for headers.
#define _MaCh3_Safe_Include_End_
KS: Restore warning checking after including external headers.
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...
TH2Poly * RatioPolys(TH2Poly *NumPoly, TH2Poly *DenomPoly)
Helper to make ratio of TH2Polys.
void CheckTH2PolyFileVersion(TFile *file)
KS: ROOT changes something with binning when moving from ROOT 5 to ROOT 6. If you open ROOT5 produced...
TH1D * PolyProjectionX(TObject *poly, std::string TempName, const std::vector< double > &xbins, const bool computeErrors=false)
WP: Poly Projectors.
void MakeFluctuatedHistogramAlternative(TH1D *FluctHist, TH1D *PolyHist, TRandom3 *rand)
Make Poisson fluctuation of TH1D hist using slow method which is only for cross-check.
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 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.
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 * PolyScaleWidth(TH2Poly *Histogram, double scale)
WP: Helper to scale th2poly analogous to th2d scale with option "width".
std::vector< Base * > CastVector(const std::vector< Derived * > &inputVec)
Converts a vector of pointers from a derived type to a base type.
double returnCherenkovThresholdMomentum(int PDG)
DB Get the Cherenkov momentum threshold in MeV.
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.
HistType * RatioHists(HistType *NumHist, HistType *DenomHist)
Helper to make ratio histograms.
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 * 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.
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 * NormalisePoly(TH2Poly *Histogram)
WP: Helper to Normalise histograms.
TH2D * ConvertTH2PolyToTH2D(TH2Poly *poly, TH2D *TH2Dhist)
KS: Convert TH2D to TH2Poly.
std::string file_exists(std::string filename)
Helper to check if files exist or not.
TH1D * PolyProjectionY(TObject *poly, std::string TempName, const std::vector< double > &ybins, const bool computeErrors=false)
WP: Poly Projectors.
Custom exception class for MaCh3 errors.
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.
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.