]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - RALICE/AliMath.h
Coding convetion violations fixed.
[u/mrichter/AliRoot.git] / RALICE / AliMath.h
index 3ff93f4d4f50ddbf3f240c1ce27c3134f8b99895..e0c71b088590865af77b8c272b8dd54673360354 100644 (file)
@@ -8,6 +8,10 @@
 #include <math.h>
  
 #include "TObject.h"
+#include "TF1.h"
+#include "TString.h"
+#include "TMath.h"
+#include "TH1F.h"
  
 class AliMath : public TObject
 {
@@ -23,7 +27,26 @@ class AliMath : public TObject
   Double_t Prob(Double_t chi2,Int_t ndf,Int_t mode=1) const; // Compute Chi-squared probability
   Double_t BesselI(Int_t n,Double_t x) const;                // Compute integer order mod. Bessel function I_n(x)
   Double_t BesselK(Int_t n,Double_t x) const;                // Compute integer order mod. Bessel function K_n(x)
+  TF1 Chi2Dist(Int_t ndf) const;                             // Provide the Chi-squared distribution function
+  TF1 StudentDist(Double_t ndf) const;                       // Provide the Student's T distribution function
+  TF1 FratioDist(Int_t ndf1,Int_t ndf2) const;               // Provide the Fratio distribution function
+  TF1 BinomialDist(Int_t n,Double_t p) const;                // Provide the Binomial distribution function
+  TF1 NegBinomialDist(Int_t k,Double_t p) const;             // Provide the Negative Binomial distribution function
+  TF1 PoissonDist(Double_t mu) const;                        // Provide the Poisson distribution function
+  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
+  Double_t StudentPvalue(Double_t t,Double_t ndf,Int_t sides=0,Int_t sigma=0) const; // Student's P-value
+  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
+  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
+  Double_t PoissonPvalue(Int_t k,Double_t mu,Int_t sides=0,Int_t sigma=0) const; // Poisson P-value
+  Double_t NegBinomialPvalue(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
+  Double_t Nfac(Int_t n,Int_t mode=0) const;    // Compute n!
+  Double_t LnNfac(Int_t n,Int_t mode=2) const;  // Compute ln(n!) 
+  Double_t LogNfac(Int_t n,Int_t mode=2) const; // Compute log_10(n!) 
+  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
+  Double_t PsiValue(TH1F* his,TH1F* hyp=0,TF1* pdf=0,Int_t f=0) const; // Bayesian Psi value of a counting exp. w.r.t. hypothesis
+  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
+  Double_t Chi2Value(TH1F* his,TH1F* hyp=0,TF1* pdf=0,Int_t* ndf=0) const; // Frequentist Chi2 value of a counting exp. w.r.t. hypothesis
+
  protected:
   Double_t GamSer(Double_t a,Double_t x) const; // Compute P(a,x) via serial representation
   Double_t GamCf(Double_t a,Double_t x) const;  // Compute P(a,x) via continued fractions
@@ -32,7 +55,7 @@ class AliMath : public TObject
   Double_t BesselI1(Double_t x) const;          // Compute modified Bessel function I_1(x)
   Double_t BesselK1(Double_t x) const;          // Compute modified Bessel function K_1(x)
  
- ClassDef(AliMath,3) // Various mathematical tools for physics analysis.
+ ClassDef(AliMath,7) // Various mathematical tools for physics analysis.
  
 };
 #endif