MaCh3  2.5.0
Reference Guide
StatisticalUtils.h
Go to the documentation of this file.
1 #pragma once
2 
3 // C++ includes
4 #include <iostream>
5 #include <sstream>
6 #include <iomanip>
7 #include <vector>
8 #include <algorithm>
9 #include <cmath>
10 
11 //MaCh3 includes
12 #include "Manager/Manager.h"
13 #include "Samples/HistogramUtils.h"
14 
16 // ROOT includes
17 #include "TCanvas.h"
18 #include "TLine.h"
19 #include "TROOT.h"
20 #include "TStyle.h"
22 
26 
27 
28 namespace M3 {
30  enum kInfCrit {
31  kBIC,
32  kDIC,
34  kInfCrits
35  };
36 }
37 
40 std::string GetJeffreysScale(const double BayesFactor);
41 
67 std::string GetDunneKaboth(const double BayesFactor);
68 
70 double GetSigmaValue(const int sigma);
71 
73 double GetBIC(const double llh, const int data, const int nPars);
74 
77 double GetNeffective(const int N1, const int N2);
78 
83 void CheckBonferoniCorrectedpValue(const std::vector<std::string>& SampleNameVec,
84  const std::vector<double>& PValVec,
85  const double Threshold = 0.05);
86 
90 double GetAndersonDarlingTestStat(const double CumulativeData, const double CumulativeMC, const double CumulativeJoint);
91 
93 int GetNumberOfRuns(const std::vector<int>& GroupClasifier);
94 
100 double GetBetaParameter(const double data, const double mc, const double w2, TestStatistic TestStat);
101 
105 double GetSubOptimality(const std::vector<double>& EigenValues, const int TotalTarameters);
106 
111 void GetArithmetic(TH1D * const hist, double& Mean, double& Error);
112 
118 void GetGaussian(TH1D*& hist, TF1* gauss, double& Mean, double& Error);
119 
127 void GetHPD(TH1D* const hist, double& Mean, double& Error, double& Error_p, double& Error_m, const double coverage = 0.6827);
128 
133 void GetCredibleInterval(const std::unique_ptr<TH1D>& hist, std::unique_ptr<TH1D>& hpost_copy, const double coverage = 0.6827);
134 
140 void GetCredibleIntervalSig(const std::unique_ptr<TH1D>& hist, std::unique_ptr<TH1D>& hpost_copy, const bool CredibleInSigmas, const double coverage = 0.6827);
141 
145 void GetCredibleRegion(std::unique_ptr<TH2D>& hist2D, const double coverage = 0.6827);
146 
151 void GetCredibleRegionSig(std::unique_ptr<TH2D>& hist2D, const bool CredibleInSigmas, const double coverage = 0.6827);
152 
155 double GetIQR(TH1D *Hist);
156 
162 double ComputeKLDivergence(TH2Poly* DataPoly, TH2Poly* PolyMC);
163 
168 double FisherCombinedPValue(const std::vector<double>& pvalues);
169 
176 void ThinningMCMC(const std::string& FilePath, const int ThinningCut);
177 
189 double GetZScore(const double value, const double mean, const double stddev);
190 
198 double GetPValueFromZScore(const double zScore);
199 
202 double GetModeError(TH1D* hpost);
203 
207 void Get2DBayesianpValue(TH2D *Histogram);
208 
209 
210 void PassErrorToRatioPlot(TH1D* RatioHist, TH1D* Hist1, TH1D* DataHist);
#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:141
double ** Mean
Definition: RHat.cpp:63
TestStatistic
Make an enum of the test statistic that we're using.
double GetSigmaValue(const int sigma)
KS: Convert sigma from normal distribution into percentage.
double ComputeKLDivergence(TH2Poly *DataPoly, TH2Poly *PolyMC)
Compute the Kullback-Leibler divergence between two TH2Poly histograms.
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,...
double GetSubOptimality(const std::vector< double > &EigenValues, const int TotalTarameters)
Based on .
void GetCredibleRegion(std::unique_ptr< TH2D > &hist2D, const double coverage=0.6827)
KS: Set 2D contour within some coverage.
double GetBIC(const double llh, const int data, const int nPars)
Get the Bayesian Information Criterion (BIC) or Schwarz information criterion (also SIC,...
double GetPValueFromZScore(const double zScore)
Compute the P-value from a given Z-score.
void Get2DBayesianpValue(TH2D *Histogram)
Calculates the 2D Bayesian p-value and generates a visualization.
double GetZScore(const double value, const double mean, const double stddev)
Compute the Z-score for a given value.
void GetGaussian(TH1D *&hist, TF1 *gauss, double &Mean, double &Error)
CW: Fit Gaussian to posterior.
double GetModeError(TH1D *hpost)
Get the mode error from a TH1D.
double GetAndersonDarlingTestStat(const double CumulativeData, const double CumulativeMC, const double CumulativeJoint)
double FisherCombinedPValue(const std::vector< double > &pvalues)
KS: Combine p-values using Fisher's method.
void PassErrorToRatioPlot(TH1D *RatioHist, TH1D *Hist1, TH1D *DataHist)
void GetCredibleRegionSig(std::unique_ptr< TH2D > &hist2D, const bool CredibleInSigmas, const double coverage=0.6827)
KS: Set 2D contour within some coverage.
double GetIQR(TH1D *Hist)
Interquartile Range (IQR)
void ThinningMCMC(const std::string &FilePath, const int ThinningCut)
Thin MCMC Chain, to save space and maintain low autocorrelations.
std::string GetJeffreysScale(const double BayesFactor)
KS: Following H. Jeffreys .
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)
double GetNeffective(const int N1, const int N2)
KS: See 14.3.10 in Numerical Recipes in C .
void GetArithmetic(TH1D *const hist, double &Mean, double &Error)
CW: Get Arithmetic mean from posterior.
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,...
std::string GetDunneKaboth(const double BayesFactor)
Convert a Bayes factor into an approximate particle-physics significance level using the Dunne–Kaboth...
int GetNumberOfRuns(const std::vector< int > &GroupClasifier)
KS: https://esjeevanand.uccollege.edu.in/wp-content/uploads/sites/114/2020/08/NON-PARAMTERIC-TEST-6....
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.
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.
Main namespace for MaCh3 software.
kInfCrit
KS: Different Information Criterion tests mostly based Gelman paper.
@ kWAIC
Watanabe-Akaike information criterion.
@ kInfCrits
This only enumerates.
@ kBIC
Bayesian Information Criterion.
@ kDIC
Deviance Information Criterion.