![]() |
MaCh3 2.2.1
Reference Guide
|
#include "Fitters/StatisticalUtils.h"
Go to the source code of this file.
Functions | |
std::string | GetJeffreysScale (const double BayesFactor) |
KS: Following H. Jeffreys [17]. | |
std::string | GetDunneKaboth (const double BayesFactor) |
KS: Based on Table 1 in https://www.t2k.org/docs/technotes/435. | |
double | GetSigmaValue (const int sigma) |
KS: Convert sigma from normal distribution into percentage. | |
double | GetBIC (const double llh, const int data, const int nPars) |
Get the Bayesian Information Criterion (BIC) or Schwarz information criterion (also SIC, SBC, SBIC) | |
double | GetNeffective (const int N1, const int N2) |
KS: See 14.3.10 in Numerical Recipes in C [23]. | |
void | CheckBonferoniCorrectedpValue (const std::vector< std::string > &SampleNameVec, const std::vector< double > &PValVec, const double Threshold) |
KS: For more see https://www.t2k.org/docs/technotes/429/TN429_v8#page=63. | |
double | GetAndersonDarlingTestStat (const double CumulativeData, const double CumulativeMC, const double CumulativeJoint) |
int | GetNumberOfRuns (const std::vector< int > &GroupClasifier) |
KS: https://esjeevanand.uccollege.edu.in/wp-content/uploads/sites/114/2020/08/NON-PARAMTERIC-TEST-6.pdf. | |
double | GetBetaParameter (const double data, const double mc, const double w2, TestStatistic TestStat) |
KS: Calculate Beta parameter which will be different based on specified test statistic. | |
double | GetSubOptimality (const std::vector< double > &EigenValues, const int TotalTarameters) |
Based on [24]. | |
void | GetArithmetic (TH1D *const hist, double &Mean, double &Error) |
CW: Get Arithmetic mean from posterior. | |
void | GetGaussian (TH1D *&hist, TF1 *gauss, double &Mean, double &Error) |
CW: Fit Gaussian to posterior. | |
void | GetHPD (TH1D *const hist, double &Mean, double &Error, double &Error_p, double &Error_m, const double coverage) |
Get Highest Posterior Density (HPD) | |
void | GetCredibleInterval (const std::unique_ptr< TH1D > &hist, std::unique_ptr< TH1D > &hpost_copy, const double coverage) |
KS: Get 1D histogram within credible interval, hpost_copy has to have the same binning, I don't do Copy() as this will lead to problems if this is used under multithreading. | |
void | GetCredibleIntervalSig (const std::unique_ptr< TH1D > &hist, std::unique_ptr< TH1D > &hpost_copy, const bool CredibleInSigmas, const double coverage) |
KS: Get 1D histogram within credible interval, hpost_copy has to have the same binning, I don't do Copy() as this will lead to problems if this is used under multithreading. | |
void | GetCredibleRegion (std::unique_ptr< TH2D > &hist2D, const double coverage) |
KS: Set 2D contour within some coverage. | |
void | GetCredibleRegionSig (std::unique_ptr< TH2D > &hist2D, const bool CredibleInSigmas, const double coverage) |
KS: Set 2D contour within some coverage. | |
double | GetIQR (TH1D *Hist) |
Interquartile Range (IQR) | |
double | ComputeKLDivergence (TH2Poly *DataPoly, TH2Poly *PolyMC) |
Compute the Kullback-Leibler divergence between two TH2Poly histograms. | |
double | FisherCombinedPValue (const std::vector< double > &pvalues) |
KS: Combine p-values using Fisher's method. | |
void | ThinningMCMC (const std::string &FilePath, const int ThinningCut) |
Thin MCMC Chain, to save space and maintain low autocorrelations. | |
double | GetZScore (const double value, const double mean, const double stddev) |
Compute the Z-score for a given value. | |
double | GetPValueFromZScore (const double zScore) |
Compute the P-value from a given Z-score. | |
double | GetModeError (TH1D *hpost) |
Get the mode error from a TH1D. | |
void | Get2DBayesianpValue (TH2D *Histogram) |
Calculates the 2D Bayesian p-value and generates a visualization. | |
void CheckBonferoniCorrectedpValue | ( | const std::vector< std::string > & | SampleNameVec, |
const std::vector< double > & | PValVec, | ||
const double | Threshold = 0.05 |
||
) |
KS: For more see https://www.t2k.org/docs/technotes/429/TN429_v8#page=63.
SampleNameVec | vector of sample names |
PValVec | pvalue for each sample |
Threshold | pvalue accepted threshold, usually 5% |
Definition at line 90 of file StatisticalUtils.cpp.
double ComputeKLDivergence | ( | TH2Poly * | DataPoly, |
TH2Poly * | PolyMC | ||
) |
Compute the Kullback-Leibler divergence between two TH2Poly histograms.
DataPoly | Pointer to the data histogram (TH2Poly). |
PolyMC | Pointer to the Monte Carlo histogram (TH2Poly). |
Definition at line 456 of file StatisticalUtils.cpp.
double FisherCombinedPValue | ( | const std::vector< double > & | pvalues | ) |
KS: Combine p-values using Fisher's method.
pvalues | A vector of individual p-values to combine. |
Definition at line 471 of file StatisticalUtils.cpp.
void Get2DBayesianpValue | ( | TH2D * | Histogram | ) |
Calculates the 2D Bayesian p-value and generates a visualization.
Histogram | A pointer to a TH2D histogram object. The function modifies the histogram's title to include the p-value information. |
TempCanvas->Write()
. Definition at line 596 of file StatisticalUtils.cpp.
double GetAndersonDarlingTestStat | ( | const double | CumulativeData, |
const double | CumulativeMC, | ||
const double | CumulativeJoint | ||
) |
CumulativeData | Value of CDF for data |
CumulativeMC | Value of CDF for MC |
CumulativeJoint | Value of CDF for joint data and MC distribution |
Definition at line 128 of file StatisticalUtils.cpp.
void GetArithmetic | ( | TH1D *const | hist, |
double & | Mean, | ||
double & | Error | ||
) |
CW: Get Arithmetic mean from posterior.
hist | histograms from which we extract arithmetic mean |
Mean | Arithmetic Mean value |
Error | Arithmetic Error value |
Definition at line 201 of file StatisticalUtils.cpp.
double GetBetaParameter | ( | const double | data, |
const double | mc, | ||
const double | w2, | ||
TestStatistic | TestStat | ||
) |
KS: Calculate Beta parameter which will be different based on specified test statistic.
data | Number of data events in a bin |
mc | Number of MC events in a bin |
w2 | Value of weight squared in a bin |
TestStat | Test statistic based on which we calculate beta |
Definition at line 154 of file StatisticalUtils.cpp.
double GetBIC | ( | const double | llh, |
const int | data, | ||
const int | nPars | ||
) |
Get the Bayesian Information Criterion (BIC) or Schwarz information criterion (also SIC, SBC, SBIC)
Definition at line 68 of file StatisticalUtils.cpp.
void GetCredibleInterval | ( | const std::unique_ptr< TH1D > & | hist, |
std::unique_ptr< TH1D > & | hpost_copy, | ||
const double | coverage = 0.6827 |
||
) |
KS: Get 1D histogram within credible interval, hpost_copy has to have the same binning, I don't do Copy() as this will lead to problems if this is used under multithreading.
hist | histograms based on which we calculate credible interval |
coverage | What is defined coverage, by default 0.6827 (1 sigma) |
Loop over histogram bins with highest number of entries until covered 90 or 68.3%
Get bin of highest content and save the number of entries reached so far
Replace bin value by -1 so it is not looped over as being maximum bin again
Definition at line 316 of file StatisticalUtils.cpp.
void GetCredibleIntervalSig | ( | const std::unique_ptr< TH1D > & | hist, |
std::unique_ptr< TH1D > & | hpost_copy, | ||
const bool | CredibleInSigmas, | ||
const double | coverage = 0.6827 |
||
) |
KS: Get 1D histogram within credible interval, hpost_copy has to have the same binning, I don't do Copy() as this will lead to problems if this is used under multithreading.
hist | histograms based on which we calculate credible interval |
CredibleInSigmas | Whether interval is in sigmas or percentage |
coverage | What is defined coverage, by default 0.6827 (1 sigma) |
Definition at line 367 of file StatisticalUtils.cpp.
void GetCredibleRegion | ( | std::unique_ptr< TH2D > & | hist2D, |
const double | coverage = 0.6827 |
||
) |
KS: Set 2D contour within some coverage.
hist2D | histograms based on which we calculate credible regions |
coverage | What is defined coverage, by default 0.6827 (1 sigma) |
Loop over histogram bins with highest number of entries until covered 90 or 68.3%
Get bin of highest content and save the number of entries reached so far
Replace bin value by -1 so it is not looped over as being maximum bin again
Definition at line 380 of file StatisticalUtils.cpp.
void GetCredibleRegionSig | ( | std::unique_ptr< TH2D > & | hist2D, |
const bool | CredibleInSigmas, | ||
const double | coverage = 0.6827 |
||
) |
KS: Set 2D contour within some coverage.
hist2D | histograms based on which we calculate credible regions |
CredibleInSigmas | Whether interval is in sigmas or percentage |
coverage | What is defined coverage, by default 0.6827 (1 sigma) |
Definition at line 431 of file StatisticalUtils.cpp.
std::string GetDunneKaboth | ( | const double | BayesFactor | ) |
KS: Based on Table 1 in https://www.t2k.org/docs/technotes/435.
BayesFactor | Obtained value of Bayes factor |
Definition at line 20 of file StatisticalUtils.cpp.
void GetGaussian | ( | TH1D *& | hist, |
TF1 * | gauss, | ||
double & | Mean, | ||
double & | Error | ||
) |
CW: Fit Gaussian to posterior.
hist | histograms to which we fit gaussian |
gauss | tf1 with gaussian, we pass pointer to make things faster |
Mean | Gaussian Mean value |
Error | Gaussian Error value |
Definition at line 208 of file StatisticalUtils.cpp.
void GetHPD | ( | TH1D *const | hist, |
double & | Mean, | ||
double & | Error, | ||
double & | Error_p, | ||
double & | Error_m, | ||
const double | coverage = 0.6827 |
||
) |
Get Highest Posterior Density (HPD)
hist | histograms from which we HPD |
Mean | HPD Mean value |
Error | HPD Error value |
Error_p | HPD Negative (left hand side) Error value |
Error_m | HPD Positive (right hand side) Error value |
coverage | What is defined coverage, by default 0.6827 (1 sigma) |
Definition at line 228 of file StatisticalUtils.cpp.
double GetIQR | ( | TH1D * | Hist | ) |
Interquartile Range (IQR)
hist | histograms from which we IQR |
Definition at line 443 of file StatisticalUtils.cpp.
std::string GetJeffreysScale | ( | const double | BayesFactor | ) |
KS: Following H. Jeffreys [17].
BayesFactor | Obtained value of Bayes factor |
Definition at line 5 of file StatisticalUtils.cpp.
double GetModeError | ( | TH1D * | hpost | ) |
Get the mode error from a TH1D.
hpost | hist from which we extract mode error |
Definition at line 542 of file StatisticalUtils.cpp.
double GetNeffective | ( | const int | N1, |
const int | N2 | ||
) |
KS: See 14.3.10 in Numerical Recipes in C [23].
Definition at line 81 of file StatisticalUtils.cpp.
int GetNumberOfRuns | ( | const std::vector< int > & | GroupClasifier | ) |
KS: https://esjeevanand.uccollege.edu.in/wp-content/uploads/sites/114/2020/08/NON-PARAMTERIC-TEST-6.pdf.
Definition at line 137 of file StatisticalUtils.cpp.
double GetPValueFromZScore | ( | const double | zScore | ) |
Compute the P-value from a given Z-score.
The P-value represents the probability of observing a value as extreme as the given Z-score under the null hypothesis of a standard normal distribution.
zScore | The Z-score for which to compute the P-value. |
Definition at line 535 of file StatisticalUtils.cpp.
double GetSigmaValue | ( | const int | sigma | ) |
KS: Convert sigma from normal distribution into percentage.
Definition at line 36 of file StatisticalUtils.cpp.
double GetSubOptimality | ( | const std::vector< double > & | EigenValues, |
const int | TotalTarameters | ||
) |
Based on [24].
EigenValues | Eigen values of covariance matrix |
Definition at line 184 of file StatisticalUtils.cpp.
double GetZScore | ( | const double | value, |
const double | mean, | ||
const double | stddev | ||
) |
Compute the Z-score for a given value.
The Z-score indicates how many standard deviations a value is from the mean. A positive Z-score means the value is above the mean, while a negative Z-score means it is below the mean.
value | The data point for which to compute the Z-score. |
mean | The mean of the data set. |
stddev | The standard deviation of the data set. Must be non-zero. |
Definition at line 529 of file StatisticalUtils.cpp.
void ThinningMCMC | ( | const std::string & | FilePath, |
const int | ThinningCut | ||
) |
Thin MCMC Chain, to save space and maintain low autocorrelations.
FilePath | Path to MCMC chain you want to thin |
ThinningCut | every which entry you want to thin [19] |
Definition at line 487 of file StatisticalUtils.cpp.