NCFS-Pack
A generic (astro)particle physics analysis framework
 
Loading...
Searching...
No Matches
NcMath.h
Go to the documentation of this file.
1#ifndef NcMath_h
2#define NcMath_h
3// Copyright(c) 1998 NCFS/IIHE, All Rights Reserved.
4// See cxx source for full Copyright notice.
5
6#include <math.h>
7
8#include "TObject.h"
9#include "TF1.h"
10#include "TString.h"
11#include "TMath.h"
12#include "TH1I.h"
13#include "TH1F.h"
14#include "TH1D.h"
15#include "TList.h"
16#include "TLine.h"
17#include "TLegend.h"
18#include "TFeldmanCousins.h"
19#include "TGraph.h"
20
21#include "NcRandom.h"
22
24
25class NcMath : public TObject
26{
27 public:
28 NcMath(); // Default constructor
29 virtual ~NcMath(); // Destructor
30 NcMath(const NcMath& m); // Copy constructor
31 Double_t Zeta(Double_t x,Int_t nterms=100000) const; // The Riemann Zeta function for x>1.
32 Double_t Gamma(Double_t z) const; // Standard gamma function Gamma(z)
33 Double_t Gamma(Double_t a,Double_t x,Int_t mode=0) const; // Incomplete gamma function P(a,x) or gamma(a,x)
34 Double_t LnGamma(Double_t z) const; // Compute ln[Gamma(z)]
35 Double_t LnGamma(Double_t a,Double_t x,Int_t mode=0) const;// Compute ln of the incomplete gamma function P(a,x) or gamma(a,x)
36 Double_t Erf(Double_t x) const; // Error function erf(x)
37 Double_t Erfc(Double_t x) const; // Complementary error function erfc(x)
38 Double_t Prob(Double_t chi2,Int_t ndf,Int_t mode=1) const; // Compute Chi-squared probability
39 Double_t BesselI(Int_t n,Double_t x) const; // Compute integer order mod. Bessel function I_n(x)
40 Double_t BesselK(Int_t n,Double_t x) const; // Compute integer order mod. Bessel function K_n(x)
41 TF1 Chi2Dist(Int_t ndf) const; // Provide the Chi-squared distribution function
42 TF1 Chi2CDF(Int_t ndf) const; // Provide the Chi-squared cumulative distribution function
43 TF1 StudentDist(Double_t ndf) const; // Provide the Student's T distribution function
44 TF1 StudentCDF(Double_t ndf) const; // Provide the Student's T cumulative distribution function
45 TF1 FratioDist(Int_t ndf1,Int_t ndf2) const; // Provide the Fratio distribution function
46 TF1 FratioCDF(Int_t ndf1,Int_t ndf2) const; // Provide the Fratio cumulative distribution function
47 TF1 BinomialDist(Int_t n,Double_t p) const; // Provide the Binomial PDF p(k|n,p) for k successes
48 TF1 BinomialCDF(Int_t n,Double_t p) const; // Provide the Binomial CDF for p(k|n,p)
49 TF1 NegBinomialnDist(Int_t k,Double_t p) const; // Provide the Negative Binomial PDF p(n|k,p) for n trials
50 TF1 NegBinomialnCDF(Int_t k,Double_t p) const; // Provide the Negative Binomial CDF for p(n|k,p)
51 TF1 NegBinomialxDist(Int_t k,Double_t p) const; // Provide the Negative Binomial PDF p(x|k,p) for x failures
52 TF1 NegBinomialxCDF(Int_t k,Double_t p) const; // Provide the Negative Binomial CDF for p(x|k,p)
53 TF1 PoissonDist(Double_t mu) const; // Provide the Poisson p(n|mu) pdf for mean mu
54 TF1 PoissonCDF(Double_t mu) const; // Provide the Poisson p(n|mu) CDF for mean mu
55 TF1 PoissonDist(Double_t r,Double_t dt) const; // Provide the Poisson p(n|r,dt) pdf for rate r and time interval dt
56 TF1 PoissonCDF(Double_t r,Double_t dt) const; // Provide the Poisson p(n|r,dt) CDF for rate r and time interval dt
57 TF1 PoissonDtDist(Double_t r,Int_t n) const; // Provide the Poisson related p(dt|r,n) pdf for rate r and number of events n
58 TF1 PoissonDtCDF(Double_t r,Int_t n) const; // Provide the Poisson related p(dt|r,n) CDF for rate r and number of events n
59 TF1 GammaDtDist(Double_t r,Double_t z) const; // Provide the p(dt|r,z) pdf for given rate r and observed occurrences z
60 TF1 GaussDist(Double_t mu,Double_t sigma) const; // Provide the Gauss p(x|mu,sigma) pdf for mean mu and std dev. sigma
61 TF1 GaussCDF(Double_t mu,Double_t sigma) const; // Provide the Gauss p(x|mu,sigma) CDF for mean mu and std dev. sigma
62 TF1 RayleighDist(Double_t sigma) const; // Provide the Rayleigh p(r|sigma) pdf for a scale factor sigma
63 TF1 RayleighCDF(Double_t sigma) const; // Provide the Rayleigh p(r|sigma) CDF for a scale factor sigma
64 TF1 Rayleigh3Dist(Double_t sigma) const; // Provide the 3-dimensional Rayleigh p(r|sigma) pdf for a scale factor sigma
65 Double_t GetStatistic(TF1 f,TString name,Int_t n=0,Double_t vref=0,Int_t npx=1000) const; // Provide statistic specified by "name" for the function "f"
66 TGraph GetCDF(TF1 f,Int_t npx=1000) const; // Provide the CDF for the specified function "f" as a TGraph
67 Double_t GaussProb(Double_t q,Double_t mean=0,Double_t sigma=1,Int_t isig=0) const; // P(|x-mean|<=|q-mean|) for Gauss pdf
68 Double_t GaussPvalue(Double_t q,Double_t mean=0,Double_t sigma=1,Int_t sides=2,Int_t isig=0) const; // P-value of q for Gauss pdf
69 Double_t Chi2Pvalue(Double_t chi2,Int_t ndf,Int_t sides=0,Int_t sigma=0,Int_t mode=1) const; // Chi-squared P-value
70 Double_t StudentPvalue(Double_t t,Double_t ndf,Int_t sides=0,Int_t sigma=0) const; // Student's P-value
71 Double_t FratioPvalue(Double_t f,Int_t ndf1,Int_t ndf2,Int_t sides=0,Int_t sigma=0) const; // F ratio P-value
72 Double_t BinomialPvalue(Int_t k,Int_t n,Double_t p,Int_t sides=0,Int_t sigma=0,Int_t mode=0) const; // Bin. P-value
73 Double_t PoissonPvalue(Int_t k,Double_t mu,Int_t sides=0,Int_t sigma=0) const; // Poisson P-value
74 Double_t PoissonPvalue(Int_t k,Double_t r,Double_t dt,Int_t sides=0,Int_t sigma=0) const; // Poisson P-value
75 Double_t PoissonDtPvalue(Double_t dt,Double_t r,Int_t n,Int_t sides=0,Int_t sigma=0) const; // Poisson dt P-value
76 Double_t NegBinomialnPvalue(Int_t n,Int_t k,Double_t p,Int_t sides=0,Int_t sigma=0,Int_t mode=0) const; // NegBin. P-value for n trials
77 Double_t NegBinomialxPvalue(Int_t x,Int_t k,Double_t p,Int_t sides=0,Int_t sigma=0,Int_t mode=0) const; // NegBin. P-value for x failures
78 Double_t Log(Double_t B,Double_t x) const; // Compute log_B(x) with base B
79 Double_t Nfac(Int_t n,Int_t mode=0) const; // Compute n!
80 Double_t LnNfac(Int_t n,Int_t mode=2) const; // Compute ln(n!)
81 Double_t LogNfac(Int_t n,Int_t mode=2) const; // Compute log_10(n!)
82 Double_t Rfac(Double_t r) const; // Compute r! for a fractional value r
83 Double_t LnRfac(Double_t r) const; // Compute ln(r!) for a fractional value r
84 Double_t LogRfac(Double_t r) const; // Compute log_10(r!) for a fractional value r
85 Double_t PsiValue(Int_t m,Int_t* n,Double_t* p=0,Int_t f=0) const; // Bayesian Psi value of a counting exp. w.r.t. hypothesis
86 Double_t PsiValue(Int_t m,Double_t* n,Double_t* p=0,Int_t f=0) const; // Bayesian Psi value of a counting exp. w.r.t. hypothesis
87 Double_t PsiValue(TH1* his,TH1* hyp=0,TF1* pdf=0,Int_t f=0) const; // Bayesian Psi value of a counting exp. w.r.t. hypothesis
88 Double_t PsiExtreme(Double_t n,Int_t m,Double_t* p=0,Int_t k=0) const; // Extreme Bayesian Psi values w.r.t. some hypothesis
89 Double_t PsiExtreme(TH1* his,TH1* hyp=0,TF1* pdf=0,Int_t k=0) const; // Extreme Bayesian Psi values w.r.t. some hypothesis
90 Double_t PsiPvalue(Double_t psi0,Double_t nr,Double_t n,Int_t m,Double_t* p=0,Int_t f=0,Double_t* na=0,TH1F* psih=0,Int_t ncut=0,Double_t* nrx=0,Int_t mark=1); // P-value of psi0 for specified B_m exp.
91 Double_t PsiPvalue(Double_t psi0,Double_t nr,TH1* his,TH1* hyp=0,TF1* pdf=0,Int_t f=0,Double_t* na=0,TH1F* psih=0,Int_t ncut=0,Double_t* nrx=0,Int_t mark=1); // P-value of psi0 for specified B_m exp.
92 Double_t Chi2Value(Int_t m,Int_t* n,Double_t* p=0,Int_t* ndf=0) const; // Frequentist Chi2 value of a counting exp. w.r.t. hypothesis
93 Double_t Chi2Value(Int_t m,Double_t* n,Double_t* p=0,Int_t* ndf=0) const; // Frequentist Chi2 value of a counting exp. w.r.t. hypothesis
94 Double_t Chi2Value(TH1* his,TH1* hyp=0,TF1* pdf=0,Int_t* ndf=0) const; // Frequentist Chi2 value of a counting exp. w.r.t. hypothesis
95 Double_t MeanMu(Double_t cl,Double_t nbkg,Int_t mode,TF1* w=0,TFeldmanCousins* f=0,Int_t nmax=0) const; // Average Feldman-Cousins upper/lower limit
96 Double_t LiMaSignificance(Double_t Non,Double_t Ton,Double_t Noff,Double_t Toff,Double_t Ra=1,Double_t Re=1) const; // Li and Ma significance
97
98 protected:
99 Double_t GamSer(Double_t a,Double_t x) const; // Compute P(a,x) via serial representation
100 Double_t GamCf(Double_t a,Double_t x) const; // Compute P(a,x) via continued fractions
101 Double_t BesselI0(Double_t x) const; // Compute modified Bessel function I_0(x)
102 Double_t BesselK0(Double_t x) const; // Compute modified Bessel function K_0(x)
103 Double_t BesselI1(Double_t x) const; // Compute modified Bessel function I_1(x)
104 Double_t BesselK1(Double_t x) const; // Compute modified Bessel function K_1(x)
105
106 ClassDef(NcMath,14) // Various mathematical tools for scientific analysis.
107
108};
109#endif
Double_t StudentPvalue(Double_t t, Double_t ndf, Int_t sides=0, Int_t sigma=0) const
Definition NcMath.cxx:1875
Double_t PsiValue(Int_t m, Int_t *n, Double_t *p=0, Int_t f=0) const
Definition NcMath.cxx:2829
TF1 PoissonCDF(Double_t mu) const
Definition NcMath.cxx:1177
Double_t GetStatistic(TF1 f, TString name, Int_t n=0, Double_t vref=0, Int_t npx=1000) const
Definition NcMath.cxx:1541
Double_t BesselK(Int_t n, Double_t x) const
Definition NcMath.cxx:678
Double_t Log(Double_t B, Double_t x) const
Definition NcMath.cxx:2613
TF1 Chi2CDF(Int_t ndf) const
Definition NcMath.cxx:801
Double_t NegBinomialnPvalue(Int_t n, Int_t k, Double_t p, Int_t sides=0, Int_t sigma=0, Int_t mode=0) const
Definition NcMath.cxx:2385
TF1 StudentDist(Double_t ndf) const
Definition NcMath.cxx:824
Double_t GaussPvalue(Double_t q, Double_t mean=0, Double_t sigma=1, Int_t sides=2, Int_t isig=0) const
Definition NcMath.cxx:1732
TF1 PoissonDist(Double_t mu) const
Definition NcMath.cxx:1150
TF1 PoissonDtCDF(Double_t r, Int_t n) const
Definition NcMath.cxx:1293
Double_t LnGamma(Double_t z) const
Definition NcMath.cxx:193
Double_t BinomialPvalue(Int_t k, Int_t n, Double_t p, Int_t sides=0, Int_t sigma=0, Int_t mode=0) const
Definition NcMath.cxx:2056
Double_t Chi2Pvalue(Double_t chi2, Int_t ndf, Int_t sides=0, Int_t sigma=0, Int_t mode=1) const
Definition NcMath.cxx:1784
TF1 Chi2Dist(Int_t ndf) const
Definition NcMath.cxx:777
TF1 NegBinomialxCDF(Int_t k, Double_t p) const
Definition NcMath.cxx:1118
Double_t PoissonDtPvalue(Double_t dt, Double_t r, Int_t n, Int_t sides=0, Int_t sigma=0) const
Definition NcMath.cxx:2295
TF1 BinomialCDF(Int_t n, Double_t p) const
Definition NcMath.cxx:992
Double_t Nfac(Int_t n, Int_t mode=0) const
Definition NcMath.cxx:2631
Double_t MeanMu(Double_t cl, Double_t nbkg, Int_t mode, TF1 *w=0, TFeldmanCousins *f=0, Int_t nmax=0) const
Definition NcMath.cxx:4109
Double_t BesselI0(Double_t x) const
Definition NcMath.cxx:497
Double_t PoissonPvalue(Int_t k, Double_t mu, Int_t sides=0, Int_t sigma=0) const
Definition NcMath.cxx:2160
Double_t LiMaSignificance(Double_t Non, Double_t Ton, Double_t Noff, Double_t Toff, Double_t Ra=1, Double_t Re=1) const
Definition NcMath.cxx:4172
Double_t Gamma(Double_t z) const
Definition NcMath.cxx:121
Double_t LnRfac(Double_t r) const
Definition NcMath.cxx:2786
TF1 GammaDtDist(Double_t r, Double_t z) const
Definition NcMath.cxx:1324
Double_t GaussProb(Double_t q, Double_t mean=0, Double_t sigma=1, Int_t isig=0) const
Definition NcMath.cxx:1700
Double_t Erfc(Double_t x) const
Definition NcMath.cxx:383
TF1 BinomialDist(Int_t n, Double_t p) const
Definition NcMath.cxx:960
NcMath()
Definition NcMath.cxx:61
Double_t BesselI(Int_t n, Double_t x) const
Definition NcMath.cxx:720
Double_t FratioPvalue(Double_t f, Int_t ndf1, Int_t ndf2, Int_t sides=0, Int_t sigma=0) const
Definition NcMath.cxx:1963
Double_t LnNfac(Int_t n, Int_t mode=2) const
Definition NcMath.cxx:2688
Double_t Rfac(Double_t r) const
Definition NcMath.cxx:2764
Double_t Zeta(Double_t x, Int_t nterms=100000) const
Definition NcMath.cxx:88
TF1 Rayleigh3Dist(Double_t sigma) const
Definition NcMath.cxx:1490
Double_t PsiExtreme(Double_t n, Int_t m, Double_t *p=0, Int_t k=0) const
Definition NcMath.cxx:3138
TF1 NegBinomialxDist(Int_t k, Double_t p) const
Definition NcMath.cxx:1084
Double_t BesselK0(Double_t x) const
Definition NcMath.cxx:540
TF1 RayleighCDF(Double_t sigma) const
Definition NcMath.cxx:1468
virtual ~NcMath()
Definition NcMath.cxx:70
TGraph GetCDF(TF1 f, Int_t npx=1000) const
Definition NcMath.cxx:1649
Double_t LogRfac(Double_t r) const
Definition NcMath.cxx:2805
Double_t Prob(Double_t chi2, Int_t ndf, Int_t mode=1) const
Definition NcMath.cxx:421
TF1 GaussCDF(Double_t mu, Double_t sigma) const
Definition NcMath.cxx:1387
Double_t NegBinomialxPvalue(Int_t x, Int_t k, Double_t p, Int_t sides=0, Int_t sigma=0, Int_t mode=0) const
Definition NcMath.cxx:2498
Double_t BesselK1(Double_t x) const
Definition NcMath.cxx:631
Double_t Erf(Double_t x) const
Definition NcMath.cxx:370
TF1 RayleighDist(Double_t sigma) const
Definition NcMath.cxx:1415
TF1 GaussDist(Double_t mu, Double_t sigma) const
Definition NcMath.cxx:1357
TF1 NegBinomialnCDF(Int_t k, Double_t p) const
Definition NcMath.cxx:1054
TF1 NegBinomialnDist(Int_t k, Double_t p) const
Definition NcMath.cxx:1022
Double_t PsiPvalue(Double_t psi0, Double_t nr, Double_t n, Int_t m, Double_t *p=0, Int_t f=0, Double_t *na=0, TH1F *psih=0, Int_t ncut=0, Double_t *nrx=0, Int_t mark=1)
Definition NcMath.cxx:3528
TF1 FratioCDF(Int_t ndf1, Int_t ndf2) const
Definition NcMath.cxx:926
TF1 PoissonDtDist(Double_t r, Int_t n) const
Definition NcMath.cxx:1260
Double_t LogNfac(Int_t n, Int_t mode=2) const
Definition NcMath.cxx:2740
TF1 FratioDist(Int_t ndf1, Int_t ndf2) const
Definition NcMath.cxx:889
Double_t GamSer(Double_t a, Double_t x) const
Definition NcMath.cxx:271
TF1 StudentCDF(Double_t ndf) const
Definition NcMath.cxx:857
Double_t Chi2Value(Int_t m, Int_t *n, Double_t *p=0, Int_t *ndf=0) const
Definition NcMath.cxx:3867
Double_t GamCf(Double_t a, Double_t x) const
Definition NcMath.cxx:316
Double_t BesselI1(Double_t x) const
Definition NcMath.cxx:587