![]() |
MaCh3 2.2.1
Reference Guide
|
#include "Samples/SampleStructs.h"
#include "Parameters/ParameterStructs.h"
#include "TGraphAsymmErrors.h"
#include "TObjString.h"
#include "TRandom3.h"
Go to the source code of this file.
Namespaces | |
namespace | M3 |
Functions | |
double | OverflowIntegral (TH2Poly *poly) |
WP: Helper function for calculating unbinned Integral of TH2Poly i.e including overflow. | |
double | NoOverflowIntegral (TH2Poly *poly) |
WP: Helper function for calculating binned Integral of TH2Poly i.e not including overflow. | |
TH1D * | PolyProjectionX (TObject *poly, std::string TempName, const std::vector< double > &xbins, const bool computeErrors=false) |
WP: Poly Projectors. | |
TH1D * | PolyProjectionY (TObject *poly, std::string TempName, const std::vector< double > &ybins, const bool computeErrors=false) |
WP: Poly Projectors. | |
TH2D * | ConvertTH2PolyToTH2D (TH2Poly *poly, TH2D *TH2Dhist) |
KS: Convert TH2D to TH2Poly. | |
TH2Poly * | ConvertTH2DtoTH2Poly (TH2D *TH2Dhist) |
KS: Convert TH2Poly to TH2D. | |
TH2Poly * | NormalisePoly (TH2Poly *Histogram) |
WP: Helper to Normalise histograms. | |
void | NormaliseTH2Poly (TH2Poly *Histogram) |
Helper to Normalise histograms. | |
TH2Poly * | PolyScaleWidth (TH2Poly *Histogram, double scale) |
WP: Helper to scale th2poly analogous to th2d scale with option "width". | |
double | PolyIntegralWidth (TH2Poly *Histogram) |
WP: Helper to calc integral of th2poly analogous to th2d integra; with option "width". | |
template<class HistType > | |
HistType * | RatioHists (HistType *NumHist, HistType *DenomHist) |
Helper to make ratio histograms. | |
TH2Poly * | RatioPolys (TH2Poly *NumPoly, TH2Poly *DenomPoly) |
Helper to make ratio of TH2Polys. | |
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. | |
void | CheckTH2PolyFileVersion (TFile *file) |
KS: ROOT changes something with binning when moving from ROOT 5 to ROOT 6. If you open ROOT5 produced file with ROOT6 you will be missing 9 last bins. | |
void | RemoveFitter (TH1D *hist, const std::string &name) |
KS: Remove fitted TF1 from hist to make comparison easier. | |
void | MakeFluctuatedHistogramStandard (TH1D *FluctHist, TH1D *PolyHist, TRandom3 *rand) |
Make Poisson fluctuation of TH1D hist using default fast method. | |
void | MakeFluctuatedHistogramAlternative (TH1D *FluctHist, TH1D *PolyHist, TRandom3 *rand) |
Make Poisson fluctuation of TH1D hist using slow method which is only for cross-check. | |
void | MakeFluctuatedHistogramStandard (TH2Poly *FluctHist, TH2Poly *PolyHist, TRandom3 *rand) |
Make Poisson fluctuation of TH2Poly hist using default fast method. | |
void | MakeFluctuatedHistogramAlternative (TH2Poly *FluctHist, TH2Poly *PolyHist, TRandom3 *rand) |
Make Poisson fluctuation of TH2Poly hist using slow method which is only for cross-check. | |
int | GetRandomPoly2 (const TH2Poly *PolyHist, TRandom3 *rand) |
KS: ROOT developers were too lazy do develop getRanom2 for TH2Poly, this implementation is based on link | |
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. | |
void | FastViolinFill (TH2D *violin, TH1D *hist_1d) |
KS: Fill Violin histogram with entry from a toy. | |
template<typename Derived , typename Base > | |
std::vector< Base * > | CastVector (const std::vector< Derived * > &inputVec) |
Converts a vector of pointers from a derived type to a base type. | |
std::string | file_exists (std::string filename) |
Helper to check if files exist or not. | |
double | returnCherenkovThresholdMomentum (int PDG) |
DB Get the Cherenkov momentum threshold in MeV. | |
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 energy. | |
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 change, and if nu/anu. | |
template<typename ObjectType > | |
std::unique_ptr< ObjectType > | M3::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 * | M3::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. | |
Definition in file HistogramUtils.h.
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 change, and if nu/anu.
Definition at line 589 of file HistogramUtils.cpp.
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 energy.
Definition at line 570 of file HistogramUtils.cpp.
std::vector< Base * > CastVector | ( | const std::vector< Derived * > & | inputVec | ) |
Converts a vector of pointers from a derived type to a base type.
Derived | The derived class type. |
Base | The base class type. |
inputVec | A std::vector of pointers to Derived objects. |
std::vector
of pointers to Base
objects. Definition at line 99 of file HistogramUtils.h.
void CheckTH2PolyFileVersion | ( | TFile * | file | ) |
KS: ROOT changes something with binning when moving from ROOT 5 to ROOT 6. If you open ROOT5 produced file with ROOT6 you will be missing 9 last bins.
However if you use ROOT6 and have ROOT6 file exactly the same code will work. Something have changed with how TH2Poly bins are stored in TFile
file | ROOT file that we will make version checks |
Definition at line 11 of file HistogramUtils.cpp.
TH2Poly * ConvertTH2DtoTH2Poly | ( | TH2D * | TH2Dhist | ) |
KS: Convert TH2Poly to TH2D.
Definition at line 282 of file HistogramUtils.cpp.
TH2D * ConvertTH2PolyToTH2D | ( | TH2Poly * | poly, |
TH2D * | TH2Dhist | ||
) |
KS: Convert TH2D to TH2Poly.
Definition at line 248 of file HistogramUtils.cpp.
void FastViolinFill | ( | TH2D * | violin, |
TH1D * | hist_1d | ||
) |
KS: Fill Violin histogram with entry from a toy.
violin | hist that will be filled |
hist_1d | refence hist from which we take entries to be filled |
Definition at line 548 of file HistogramUtils.cpp.
|
inline |
Helper to check if files exist or not.
Definition at line 110 of file HistogramUtils.h.
int GetRandomPoly2 | ( | const TH2Poly * | PolyHist, |
TRandom3 * | rand | ||
) |
KS: ROOT developers were too lazy do develop getRanom2 for TH2Poly, this implementation is based on link
Definition at line 466 of file HistogramUtils.cpp.
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.
sigmaArrayLeft | sigma var hist at -1 or -3 sigma shift |
sigmaArrayCentr | sigma var hist at prior values |
sigmaArrayRight | sigma var hist at +1 or +3 sigma shift |
title | A tittle for returned object |
TGraphAsymmErrors
object that visualizes the sigma variation of spectra, showing confidence intervals between different sigma shifts. Definition at line 516 of file HistogramUtils.cpp.
void MakeFluctuatedHistogramAlternative | ( | TH1D * | FluctHist, |
TH1D * | PolyHist, | ||
TRandom3 * | rand | ||
) |
Make Poisson fluctuation of TH1D hist using slow method which is only for cross-check.
Definition at line 443 of file HistogramUtils.cpp.
void MakeFluctuatedHistogramAlternative | ( | TH2Poly * | FluctHist, |
TH2Poly * | PolyHist, | ||
TRandom3 * | rand | ||
) |
Make Poisson fluctuation of TH2Poly hist using slow method which is only for cross-check.
Definition at line 495 of file HistogramUtils.cpp.
void MakeFluctuatedHistogramStandard | ( | TH1D * | FluctHist, |
TH1D * | PolyHist, | ||
TRandom3 * | rand | ||
) |
Make Poisson fluctuation of TH1D hist using default fast method.
Definition at line 405 of file HistogramUtils.cpp.
void MakeFluctuatedHistogramStandard | ( | TH2Poly * | FluctHist, |
TH2Poly * | PolyHist, | ||
TRandom3 * | rand | ||
) |
Make Poisson fluctuation of TH2Poly hist using default fast method.
Definition at line 424 of file HistogramUtils.cpp.
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.
name | This will be tittle of output histogram |
BinArray_x | Bin edges for X axis |
BinArray_y | Bin edges for Y axis |
Definition at line 370 of file HistogramUtils.cpp.
double NoOverflowIntegral | ( | TH2Poly * | poly | ) |
WP: Helper function for calculating binned Integral of TH2Poly i.e not including overflow.
Definition at line 50 of file HistogramUtils.cpp.
TH2Poly * NormalisePoly | ( | TH2Poly * | Histogram | ) |
WP: Helper to Normalise histograms.
Definition at line 198 of file HistogramUtils.cpp.
void NormaliseTH2Poly | ( | TH2Poly * | Histogram | ) |
Helper to Normalise histograms.
Histogram | hist which we normalise |
Definition at line 211 of file HistogramUtils.cpp.
double OverflowIntegral | ( | TH2Poly * | poly | ) |
WP: Helper function for calculating unbinned Integral of TH2Poly i.e including overflow.
Definition at line 35 of file HistogramUtils.cpp.
double PolyIntegralWidth | ( | TH2Poly * | Histogram | ) |
WP: Helper to calc integral of th2poly analogous to th2d integra; with option "width".
Definition at line 351 of file HistogramUtils.cpp.
TH1D * PolyProjectionX | ( | TObject * | poly, |
std::string | TempName, | ||
const std::vector< double > & | xbins, | ||
const bool | computeErrors = false |
||
) |
WP: Poly Projectors.
Definition at line 63 of file HistogramUtils.cpp.
TH1D * PolyProjectionY | ( | TObject * | poly, |
std::string | TempName, | ||
const std::vector< double > & | ybins, | ||
const bool | computeErrors = false |
||
) |
WP: Poly Projectors.
Definition at line 131 of file HistogramUtils.cpp.
TH2Poly * PolyScaleWidth | ( | TH2Poly * | Histogram, |
double | scale | ||
) |
WP: Helper to scale th2poly analogous to th2d scale with option "width".
Definition at line 331 of file HistogramUtils.cpp.
HistType * RatioHists | ( | HistType * | NumHist, |
HistType * | DenomHist | ||
) |
Helper to make ratio histograms.
Definition at line 223 of file HistogramUtils.cpp.
TH2Poly * RatioPolys | ( | TH2Poly * | NumPoly, |
TH2Poly * | DenomPoly | ||
) |
Helper to make ratio of TH2Polys.
Definition at line 235 of file HistogramUtils.cpp.
void RemoveFitter | ( | TH1D * | hist, |
const std::string & | name | ||
) |
KS: Remove fitted TF1 from hist to make comparison easier.
Definition at line 394 of file HistogramUtils.cpp.
double returnCherenkovThresholdMomentum | ( | int | PDG | ) |
DB Get the Cherenkov momentum threshold in MeV.
Definition at line 560 of file HistogramUtils.cpp.