]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RALICE/AliMath.cxx
TPCNoiseMapComponent included into build (Kelly)
[u/mrichter/AliRoot.git] / RALICE / AliMath.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 // $Id$
17
18 ///////////////////////////////////////////////////////////////////////////
19 // Class AliMath
20 // Various mathematical tools which may be very convenient while
21 // performing physics analysis.
22 //
23 // Example : Probability of a Chi-squared value
24 // =========
25 //
26 // AliMath M;
27 // Float_t chi2=20;            // The chi-squared value
28 // Int_t ndf=12;               // The number of degrees of freedom
29 // Float_t p=M.Prob(chi2,ndf); // The probability that at least a Chi-squared
30 //                             // value of chi2 will be observed, even for a
31 //                             // correct model
32 //
33 //--- Author: Nick van Eijndhoven 14-nov-1998 UU-SAP Utrecht
34 //- Modified: NvE $Date$ UU-SAP Utrecht
35 ///////////////////////////////////////////////////////////////////////////
36
37 #include "AliMath.h"
38 #include "Riostream.h"
39  
40 ClassImp(AliMath) // Class implementation to enable ROOT I/O
41  
42 AliMath::AliMath() : TObject()
43 {
44 // Default constructor
45 }
46 ///////////////////////////////////////////////////////////////////////////
47 AliMath::~AliMath()
48 {
49 // Destructor
50 }
51 ///////////////////////////////////////////////////////////////////////////
52 AliMath::AliMath(const AliMath& m) : TObject(m)
53 {
54 // Copy constructor
55 }
56 ///////////////////////////////////////////////////////////////////////////
57 Double_t AliMath::Gamma(Double_t z) const
58 {
59 // Computation of gamma(z) for all z>0.
60 //
61 // The algorithm is based on the article by C.Lanczos [1] as denoted in
62 // Numerical Recipes 2nd ed. on p. 207 (W.H.Press et al.).
63 //
64 // [1] C.Lanczos, SIAM Journal of Numerical Analysis B1 (1964), 86.
65 //
66 //--- Nve 14-nov-1998 UU-SAP Utrecht
67  
68  if (z<=0.)
69  {
70   cout << "*Gamma(z)* Wrong argument z = " << z << endl;
71   return 0;
72  }
73  
74  Double_t v=LnGamma(z);
75  return exp(v);
76 }
77 ///////////////////////////////////////////////////////////////////////////
78 Double_t AliMath::Gamma(Double_t a,Double_t x) const
79 {
80 // Computation of the incomplete gamma function P(a,x)
81 //
82 // The algorithm is based on the formulas and code as denoted in
83 // Numerical Recipes 2nd ed. on p. 210-212 (W.H.Press et al.).
84 //
85 //--- Nve 14-nov-1998 UU-SAP Utrecht
86  
87  if (a<=0.)
88  {
89   cout << "*Gamma(a,x)* Invalid argument a = " << a << endl;
90   return 0;
91  }
92  
93  if (x<=0.)
94  {
95   if (x<0) cout << "*Gamma(a,x)* Invalid argument x = " << x << endl;
96   return 0;
97  }
98  
99  if (x<(a+1.))
100  {
101   return GamSer(a,x);
102  }
103  else
104  {
105   return GamCf(a,x);
106  }
107 }
108 ///////////////////////////////////////////////////////////////////////////
109 Double_t AliMath::LnGamma(Double_t z) const
110 {
111 // Computation of ln[gamma(z)] for all z>0.
112 //
113 // The algorithm is based on the article by C.Lanczos [1] as denoted in
114 // Numerical Recipes 2nd ed. on p. 207 (W.H.Press et al.).
115 //
116 // [1] C.Lanczos, SIAM Journal of Numerical Analysis B1 (1964), 86.
117 //
118 // The accuracy of the result is better than 2e-10.
119 //
120 //--- Nve 14-nov-1998 UU-SAP Utrecht
121  
122  if (z<=0.)
123  {
124   cout << "*LnGamma(z)* Wrong argument z = " << z << endl;
125   return 0;
126  }
127  
128  // Coefficients for the series expansion
129  Double_t c[7];
130  c[0]=  2.5066282746310005;
131  c[1]= 76.18009172947146;
132  c[2]=-86.50532032941677;
133  c[3]= 24.01409824083091;
134  c[4]= -1.231739572450155;
135  c[5]=  0.1208650973866179e-2;
136  c[6]= -0.5395239384953e-5;
137  
138  Double_t x=z;
139  Double_t y=x;
140  Double_t tmp=x+5.5;
141  tmp=(x+0.5)*log(tmp)-tmp;
142  Double_t ser=1.000000000190015;
143  for (Int_t i=1; i<7; i++)
144  {
145   y+=1.;
146   ser+=c[i]/y;
147  }
148  Double_t v=tmp+log(c[0]*ser/x);
149  return v;
150 }
151 ///////////////////////////////////////////////////////////////////////////
152 Double_t AliMath::GamSer(Double_t a,Double_t x) const
153 {
154 // Computation of the incomplete gamma function P(a,x)
155 // via its series representation.
156 //
157 // The algorithm is based on the formulas and code as denoted in
158 // Numerical Recipes 2nd ed. on p. 210-212 (W.H.Press et al.).
159 //
160 //--- Nve 14-nov-1998 UU-SAP Utrecht
161  
162  Int_t itmax=100;   // Maximum number of iterations
163  Double_t eps=3.e-7; // Relative accuracy
164  
165  if (a<=0.)
166  {
167   cout << "*GamSer(a,x)* Invalid argument a = " << a << endl;
168   return 0;
169  }
170  
171  if (x<=0.)
172  {
173   if (x<0) cout << "*GamSer(a,x)* Invalid argument x = " << x << endl;
174   return 0;
175  }
176  
177  Double_t gln=LnGamma(a);
178  Double_t ap=a;
179  Double_t sum=1./a;
180  Double_t del=sum;
181  for (Int_t n=1; n<=itmax; n++)
182  {
183   ap+=1.;
184   del=del*x/ap;
185   sum+=del;
186   if (fabs(del)<fabs(sum*eps)) break;
187   if (n==itmax) cout << "*GamSer(a,x)* a too large or itmax too small" << endl;
188  }
189  Double_t v=sum*exp(-x+a*log(x)-gln);
190  return v;
191 }
192 ///////////////////////////////////////////////////////////////////////////
193 Double_t AliMath::GamCf(Double_t a,Double_t x) const
194 {
195 // Computation of the incomplete gamma function P(a,x)
196 // via its continued fraction representation.
197 //
198 // The algorithm is based on the formulas and code as denoted in
199 // Numerical Recipes 2nd ed. on p. 210-212 (W.H.Press et al.).
200 //
201 //--- Nve 14-nov-1998 UU-SAP Utrecht
202  
203  Int_t itmax=100;      // Maximum number of iterations
204  Double_t eps=3.e-7;    // Relative accuracy
205  Double_t fpmin=1.e-30; // Smallest Double_t value allowed here
206  
207  if (a<=0.)
208  {
209   cout << "*GamCf(a,x)* Invalid argument a = " << a << endl;
210   return 0;
211  }
212  
213  if (x<=0.)
214  {
215   if (x<0) cout << "*GamCf(a,x)* Invalid argument x = " << x << endl;
216   return 0;
217  }
218  
219  Double_t gln=LnGamma(a);
220  Double_t b=x+1.-a;
221  Double_t c=1./fpmin;
222  Double_t d=1./b;
223  Double_t h=d;
224  Double_t an,del;
225  for (Int_t i=1; i<=itmax; i++)
226  {
227   an=double(-i)*(double(i)-a);
228   b+=2.;
229   d=an*d+b;
230   if (fabs(d)<fpmin) d=fpmin;
231   c=b+an/c;
232   if (fabs(c)<fpmin) c=fpmin;
233   d=1./d;
234   del=d*c;
235   h=h*del;
236   if (fabs(del-1.)<eps) break;
237   if (i==itmax) cout << "*GamCf(a,x)* a too large or itmax too small" << endl;
238  }
239  Double_t v=exp(-x+a*log(x)-gln)*h;
240  return (1.-v);
241 }
242 ///////////////////////////////////////////////////////////////////////////
243 Double_t AliMath::Erf(Double_t x) const
244 {
245 // Computation of the error function erf(x).
246 //
247 //--- NvE 14-nov-1998 UU-SAP Utrecht
248  
249  return (1.-Erfc(x));
250 }
251 ///////////////////////////////////////////////////////////////////////////
252 Double_t AliMath::Erfc(Double_t x) const
253 {
254 // Computation of the complementary error function erfc(x).
255 //
256 // The algorithm is based on a Chebyshev fit as denoted in
257 // Numerical Recipes 2nd ed. on p. 214 (W.H.Press et al.).
258 //
259 // The fractional error is always less than 1.2e-7.
260 //
261 //--- Nve 14-nov-1998 UU-SAP Utrecht
262  
263  // The parameters of the Chebyshev fit
264  const Double_t ka1=-1.26551223,  ka2=1.00002368,
265                 ka3= 0.37409196,  ka4=0.09678418,
266                 ka5=-0.18628806,  ka6=0.27886807,
267                 ka7=-1.13520398,  ka8=1.48851587,
268                 ka9=-0.82215223, ka10=0.17087277;
269  
270  Double_t v=1.; // The return value
271  
272  Double_t z=fabs(x);
273  
274  if (z <= 0.) return v; // erfc(0)=1
275  
276  Double_t t=1./(1.+0.5*z);
277  
278  v=t*exp((-z*z)
279    +ka1+t*(ka2+t*(ka3+t*(ka4+t*(ka5+t*(ka6+t*(ka7+t*(ka8+t*(ka9+t*ka10)))))))));
280  
281  if (x < 0.) v=2.-v; // erfc(-x)=2-erfc(x)
282  
283  return v;
284 }
285 ///////////////////////////////////////////////////////////////////////////
286 Double_t AliMath::Prob(Double_t chi2,Int_t ndf,Int_t mode) const
287 {
288 // Computation of the probability for a certain Chi-squared (chi2)
289 // and number of degrees of freedom (ndf).
290 //
291 // A more clear and flexible facility is offered by Chi2Pvalue. 
292 //
293 // According to the value of the parameter "mode" various algorithms
294 // can be selected.
295 //
296 // mode = 0 : Calculations are based on the incomplete gamma function P(a,x),
297 //            where a=ndf/2 and x=chi2/2.
298 //
299 //        1 : Same as for mode=0. However, in case ndf=1 an exact expression
300 //            based on the error function Erf() is used.
301 //
302 //        2 : Same as for mode=0. However, in case ndf>30 a Gaussian approximation
303 //            is used instead of the gamma function.
304 //
305 // When invoked as Prob(chi2,ndf) the default mode=1 is used.
306 //
307 // P(a,x) represents the probability that the observed Chi-squared
308 // for a correct model should be less than the value chi2.
309 //
310 // The returned probability corresponds to 1-P(a,x),
311 // which denotes the probability that an observed Chi-squared exceeds
312 // the value chi2 by chance, even for a correct model.
313 //
314 //--- NvE 14-nov-1998 UU-SAP Utrecht
315  
316  if (ndf <= 0) return 0; // Set CL to zero in case ndf<=0
317  
318  if (chi2 <= 0.)
319  {
320   if (chi2 < 0.)
321   {
322     return 0;
323   }
324   else
325   {
326    return 1;
327   }
328  }
329
330  Double_t v=-1.;
331
332  switch (mode)
333  {
334   case 1: // Exact expression for ndf=1 as alternative for the gamma function
335    if (ndf==1) v=1.-Erf(sqrt(chi2)/sqrt(2.));
336    break;
337  
338   case 2: // Gaussian approximation for large ndf (i.e. ndf>30) as alternative for the gamma function
339    if (ndf>30)
340    {
341     Double_t q=sqrt(2.*chi2)-sqrt(double(2*ndf-1));
342     if (q>0.) v=0.5*(1.-Erf(q/sqrt(2.)));
343    }
344    break;
345  }
346  
347  if (v<0.)
348  {
349   // Evaluate the incomplete gamma function
350   Double_t a=double(ndf)/2.;
351   Double_t x=chi2/2.;
352   v=1.-Gamma(a,x);
353  }
354
355  return v;
356 }
357 ///////////////////////////////////////////////////////////////////////////
358 Double_t AliMath::BesselI0(Double_t x) const
359 {
360 // Computation of the modified Bessel function I_0(x) for any real x.
361 //
362 // The algorithm is based on the article by Abramowitz and Stegun [1]
363 // as denoted in Numerical Recipes 2nd ed. on p. 230 (W.H.Press et al.).
364 //
365 // [1] M.Abramowitz and I.A.Stegun, Handbook of Mathematical Functions,
366 //     Applied Mathematics Series vol. 55 (1964), Washington.  
367 //
368 //--- NvE 12-mar-2000 UU-SAP Utrecht
369
370  // Parameters of the polynomial approximation  
371  const Double_t kp1=1.0,          kp2=3.5156229,    kp3=3.0899424,
372                 kp4=1.2067492,    kp5=0.2659732,    kp6=3.60768e-2,  kp7=4.5813e-3;
373
374  const Double_t kq1= 0.39894228,  kq2= 1.328592e-2, kq3= 2.25319e-3,
375                 kq4=-1.57565e-3,  kq5= 9.16281e-3,  kq6=-2.057706e-2,
376                 kq7= 2.635537e-2, kq8=-1.647633e-2, kq9= 3.92377e-3; 
377
378  Double_t ax=fabs(x);
379
380  Double_t y=0,result=0;
381
382  if (ax < 3.75)
383  {
384   y=pow(x/3.75,2);
385   result=kp1+y*(kp2+y*(kp3+y*(kp4+y*(kp5+y*(kp6+y*kp7)))));
386  }
387  else
388  {
389   y=3.75/ax;
390   result=(exp(ax)/sqrt(ax))
391          *(kq1+y*(kq2+y*(kq3+y*(kq4+y*(kq5+y*(kq6+y*(kq7+y*(kq8+y*kq9))))))));
392  }
393
394  return result;
395 }
396 ///////////////////////////////////////////////////////////////////////////
397 Double_t AliMath::BesselK0(Double_t x) const
398 {
399 // Computation of the modified Bessel function K_0(x) for positive real x.
400 //
401 // The algorithm is based on the article by Abramowitz and Stegun [1]
402 // as denoted in Numerical Recipes 2nd ed. on p. 230 (W.H.Press et al.).
403 //
404 // [1] M.Abramowitz and I.A.Stegun, Handbook of Mathematical Functions,
405 //     Applied Mathematics Series vol. 55 (1964), Washington.  
406 //
407 //--- NvE 12-mar-2000 UU-SAP Utrecht
408
409  // Parameters of the polynomial approximation  
410  const Double_t kp1=-0.57721566,  kp2=0.42278420,   kp3=0.23069756,
411                 kp4= 3.488590e-2, kp5=2.62698e-3,   kp6=1.0750e-4,    kp7=7.4e-6;
412
413  const Double_t kq1= 1.25331414,  kq2=-7.832358e-2, kq3= 2.189568e-2,
414                 kq4=-1.062446e-2, kq5= 5.87872e-3,  kq6=-2.51540e-3,  kq7=5.3208e-4;
415
416  if (x <= 0)
417  {
418   cout << " *BesselK0* Invalid argument x = " << x << endl;
419   return 0;
420  }
421
422  Double_t y=0,result=0;
423
424  if (x <= 2)
425  {
426   y=x*x/4.;
427   result=(-log(x/2.)*BesselI0(x))
428          +(kp1+y*(kp2+y*(kp3+y*(kp4+y*(kp5+y*(kp6+y*kp7))))));
429  }
430  else
431  {
432   y=2./x;
433   result=(exp(-x)/sqrt(x))
434          *(kq1+y*(kq2+y*(kq3+y*(kq4+y*(kq5+y*(kq6+y*kq7))))));
435  }
436
437  return result;
438 }
439 ///////////////////////////////////////////////////////////////////////////
440 Double_t AliMath::BesselI1(Double_t x) const
441 {
442 // Computation of the modified Bessel function I_1(x) for any real x.
443 //
444 // The algorithm is based on the article by Abramowitz and Stegun [1]
445 // as denoted in Numerical Recipes 2nd ed. on p. 230 (W.H.Press et al.).
446 //
447 // [1] M.Abramowitz and I.A.Stegun, Handbook of Mathematical Functions,
448 //     Applied Mathematics Series vol. 55 (1964), Washington.  
449 //
450 //--- NvE 12-mar-2000 UU-SAP Utrecht
451
452  // Parameters of the polynomial approximation  
453  const Double_t kp1=0.5,          kp2=0.87890594,   kp3=0.51498869,
454                 kp4=0.15084934,   kp5=2.658733e-2,  kp6=3.01532e-3,  kp7=3.2411e-4;
455
456  const Double_t kq1= 0.39894228,  kq2=-3.988024e-2, kq3=-3.62018e-3,
457                 kq4= 1.63801e-3,  kq5=-1.031555e-2, kq6= 2.282967e-2,
458                 kq7=-2.895312e-2, kq8= 1.787654e-2, kq9=-4.20059e-3; 
459
460  Double_t ax=fabs(x);
461
462  Double_t y=0,result=0;
463
464  if (ax < 3.75)
465  {
466   y=pow(x/3.75,2);
467   result=x*(kp1+y*(kp2+y*(kp3+y*(kp4+y*(kp5+y*(kp6+y*kp7))))));
468  }
469  else
470  {
471   y=3.75/ax;
472   result=(exp(ax)/sqrt(ax))
473          *(kq1+y*(kq2+y*(kq3+y*(kq4+y*(kq5+y*(kq6+y*(kq7+y*(kq8+y*kq9))))))));
474   if (x < 0) result=-result;
475  }
476
477  return result;
478 }
479 ///////////////////////////////////////////////////////////////////////////
480 Double_t AliMath::BesselK1(Double_t x) const
481 {
482 // Computation of the modified Bessel function K_1(x) for positive real x.
483 //
484 // The algorithm is based on the article by Abramowitz and Stegun [1]
485 // as denoted in Numerical Recipes 2nd ed. on p. 230 (W.H.Press et al.).
486 //
487 // [1] M.Abramowitz and I.A.Stegun, Handbook of Mathematical Functions,
488 //     Applied Mathematics Series vol. 55 (1964), Washington.  
489 //
490 //--- NvE 12-mar-2000 UU-SAP Utrecht
491
492  // Parameters of the polynomial approximation  
493  const Double_t kp1= 1.,          kp2= 0.15443144,  kp3=-0.67278579,
494                 kp4=-0.18156897,  kp5=-1.919402e-2, kp6=-1.10404e-3,  kp7=-4.686e-5;
495
496  const Double_t kq1= 1.25331414,  kq2= 0.23498619,  kq3=-3.655620e-2,
497                 kq4= 1.504268e-2, kq5=-7.80353e-3,  kq6= 3.25614e-3,  kq7=-6.8245e-4;
498
499  if (x <= 0)
500  {
501   cout << " *BesselK1* Invalid argument x = " << x << endl;
502   return 0;
503  }
504
505  Double_t y=0,result=0;
506
507  if (x <= 2)
508  {
509   y=x*x/4.;
510   result=(log(x/2.)*BesselI1(x))
511          +(1./x)*(kp1+y*(kp2+y*(kp3+y*(kp4+y*(kp5+y*(kp6+y*kp7))))));
512  }
513  else
514  {
515   y=2./x;
516   result=(exp(-x)/sqrt(x))
517          *(kq1+y*(kq2+y*(kq3+y*(kq4+y*(kq5+y*(kq6+y*kq7))))));
518  }
519
520  return result;
521 }
522 ///////////////////////////////////////////////////////////////////////////
523 Double_t AliMath::BesselK(Int_t n,Double_t x) const
524 {
525 // Computation of the Integer Order Modified Bessel function K_n(x)
526 // for n=0,1,2,... and positive real x.
527 //
528 // The algorithm uses the recurrence relation
529 //
530 //               K_n+1(x) = (2n/x)*K_n(x) + K_n-1(x) 
531 //
532 // as denoted in Numerical Recipes 2nd ed. on p. 232 (W.H.Press et al.).
533 //
534 //--- NvE 12-mar-2000 UU-SAP Utrecht
535
536  if (x <= 0 || n < 0)
537  {
538   cout << " *BesselK* Invalid argument(s) (n,x) = (" << n << " , " << x << ")" << endl;
539   return 0;
540  }
541
542  if (n==0) return BesselK0(x);
543
544  if (n==1) return BesselK1(x);
545
546  // Perform upward recurrence for all x
547  Double_t tox=2./x;
548  Double_t bkm=BesselK0(x);
549  Double_t bk=BesselK1(x);
550  Double_t bkp=0;
551  for (Int_t j=1; j<n; j++)
552  {
553   bkp=bkm+double(j)*tox*bk;
554   bkm=bk;
555   bk=bkp;
556  }
557
558  return bk;
559 }
560 ///////////////////////////////////////////////////////////////////////////
561 Double_t AliMath::BesselI(Int_t n,Double_t x) const
562 {
563 // Computation of the Integer Order Modified Bessel function I_n(x)
564 // for n=0,1,2,... and any real x.
565 //
566 // The algorithm uses the recurrence relation
567 //
568 //               I_n+1(x) = (-2n/x)*I_n(x) + I_n-1(x) 
569 //
570 // as denoted in Numerical Recipes 2nd ed. on p. 232 (W.H.Press et al.).
571 //
572 //--- NvE 12-mar-2000 UU-SAP Utrecht
573
574  Int_t iacc=40; // Increase to enhance accuracy
575  Double_t bigno=1.e10, bigni=1.e-10;
576
577  if (n < 0)
578  {
579   cout << " *BesselI* Invalid argument (n,x) = (" << n << " , " << x << ")" << endl;
580   return 0;
581  }
582
583  if (n==0) return BesselI0(x);
584
585  if (n==1) return BesselI1(x);
586
587  if (fabs(x) < 1.e-10) return 0;
588
589  Double_t tox=2./fabs(x);
590  Double_t bip=0,bim=0;
591  Double_t bi=1;
592  Double_t result=0;
593  Int_t m=2*((n+int(sqrt(float(iacc*n))))); // Downward recurrence from even m
594  for (Int_t j=m; j<=1; j--)
595  {
596   bim=bip+double(j)*tox*bi;
597   bip=bi;
598   bi=bim;
599   if (fabs(bi) > bigno) // Renormalise to prevent overflows
600   {
601    result*=bigni;
602    bi*=bigni;
603    bip*=bigni;
604   }
605   if (j==n) result=bip;
606  }
607
608  result*=BesselI0(x)/bi; // Normalise with I0(x)
609  if ((x < 0) && (n%2 == 1)) result=-result;
610
611  return result;
612 }
613 ///////////////////////////////////////////////////////////////////////////
614 TF1 AliMath::Chi2Dist(Int_t ndf) const
615 {
616 // Provide the Chi-squared distribution function corresponding to the
617 // specified ndf degrees of freedom.
618 //
619 // Details can be found in the excellent textbook of Phil Gregory
620 // "Bayesian Logical Data Analysis for the Physical Sciences".
621 //
622 // Note : <chi2>=ndf  Var(chi2)=2*ndf
623  
624  TF1 chi2dist("chi2dist","1./(TMath::Gamma([0]/2.)*pow(2,[0]/2.))*pow(x,[0]/2.-1.)*exp(-x/2.)");
625  TString title="#chi^{2} distribution for ndf = ";
626  title+=ndf; 
627  chi2dist.SetTitle(title.Data());
628  chi2dist.SetParName(0,"ndf");
629  chi2dist.SetParameter(0,ndf);
630
631  return chi2dist;
632 }
633 ///////////////////////////////////////////////////////////////////////////
634 TF1 AliMath::StudentDist(Double_t ndf) const
635 {
636 // Provide the Student's T distribution function corresponding to the
637 // specified ndf degrees of freedom.
638 //
639 // In a frequentist approach, the Student's T distribution is particularly
640 // useful in making inferences about the mean of an underlying population
641 // based on the data from a random sample.
642 //
643 // In a Bayesian context it is used to characterise the posterior PDF
644 // for a particular state of information. 
645 //
646 // Note : ndf is not restricted to integer values
647 //
648 // Details can be found in the excellent textbook of Phil Gregory
649 // "Bayesian Logical Data Analysis for the Physical Sciences".
650 //
651 // Note : <T>=0  Var(T)=ndf/(ndf-2)
652  
653  TF1 tdist("tdist","(TMath::Gamma(([0]+1.)/2.)/(sqrt(pi*[0])*TMath::Gamma([0]/2.)))*pow(1.+x*x/[0],-([0]+1.)/2.)");
654  TString title="Student's t distribution for ndf = ";
655  title+=ndf; 
656  tdist.SetTitle(title.Data());
657  tdist.SetParName(0,"ndf");
658  tdist.SetParameter(0,ndf);
659
660  return tdist;
661 }
662 ///////////////////////////////////////////////////////////////////////////
663 TF1 AliMath::FratioDist(Int_t ndf1,Int_t ndf2) const
664 {
665 // Provide the F (ratio) distribution function corresponding to the
666 // specified ndf1 and ndf2 degrees of freedom of the two samples.
667 //
668 // In a frequentist approach, the F (ratio) distribution is particularly useful
669 // in making inferences about the ratio of the variances of two underlying
670 // populations based on a measurement of the variance of a random sample taken
671 // from each one of the two populations.
672 // So the F test provides a means to investigate the degree of equality of
673 // two populations.
674 //
675 // Details can be found in the excellent textbook of Phil Gregory
676 // "Bayesian Logical Data Analysis for the Physical Sciences".
677 //
678 // Note : <F>=ndf2/(ndf2-2)  Var(F)=2*ndf2*ndf2*(ndf2+ndf1-2)/(ndf1*(ndf2-1)*(ndf2-1)*(ndf2-4))
679  
680  TF1 fdist("fdist",
681  "(TMath::Gamma(([0]+[1])/2.)/(TMath::Gamma([0]/2.)*TMath::Gamma([1]/2.)))*pow([0]/[1],[0]/2.)*pow(x,([0]-2.)/2.)/pow(1.+x*[0]/[1],([0]+[1])/2.)");
682  TString title="F (ratio) distribution for ndf1 = ";
683  title+=ndf1;
684  title+=" ndf2 = "; 
685  title+=ndf2;
686  fdist.SetTitle(title.Data());
687  fdist.SetParName(0,"ndf1");
688  fdist.SetParameter(0,ndf1);
689  fdist.SetParName(1,"ndf2");
690  fdist.SetParameter(1,ndf2);
691
692  return fdist;
693 }
694 ///////////////////////////////////////////////////////////////////////////
695 TF1 AliMath::BinomialDist(Int_t n,Double_t p) const
696 {
697 // Provide the Binomial distribution function corresponding to the
698 // specified number of trials n and probability p of success.
699 //
700 // p(k|n,p) = probability to obtain exactly k successes in n trials. 
701 //
702 // Details can be found in the excelent textbook of Phil Gregory
703 // "Bayesian Logical Data Analysis for the Physical Sciences".
704 //
705 // Note : <k>=n*p  Var(k)=n*p*(1-p)
706  
707  TF1 bindist("bindist","TMath::Binomial(int([0]),int(x))*pow([1],int(x))*pow(1.-[1],int([0])-int(x))");
708  TString title="Binomial distribution for n = ";
709  title+=n;
710  title+=" p = "; 
711  title+=p;
712  bindist.SetTitle(title.Data());
713  bindist.SetParName(0,"n");
714  bindist.SetParameter(0,n);
715  bindist.SetParName(1,"p");
716  bindist.SetParameter(1,p);
717
718  return bindist;
719 }
720 ///////////////////////////////////////////////////////////////////////////
721 TF1 AliMath::NegBinomialDist(Int_t k,Double_t p) const
722 {
723 // Provide the Negative Binomial distribution function corresponding to the
724 // specified number of successes k and probability p of success.
725 //
726 // p(n|k,p) = probability to have reached k successes after n trials.
727 //
728 // Details can be found in the excelent textbook of Phil Gregory
729 // "Bayesian Logical Data Analysis for the Physical Sciences".
730 //
731 // Note : <n>=k*(1-p)/p  Var(n)=k*(1-p)/(p*p) 
732  
733  TF1 negbindist("negbindist","TMath::Binomial(int(x)-1,int([0])-1)*pow([1],int([0]))*pow(1.-[1],int(x)-int([0]))");
734  TString title="Negative Binomial distribution for k = ";
735  title+=k;
736  title+=" p = "; 
737  title+=p;
738  negbindist.SetTitle(title.Data());
739  negbindist.SetParName(0,"k");
740  negbindist.SetParameter(0,k);
741  negbindist.SetParName(1,"p");
742  negbindist.SetParameter(1,p);
743
744  return negbindist;
745 }
746 ///////////////////////////////////////////////////////////////////////////
747 TF1 AliMath::PoissonDist(Double_t mu) const
748 {
749 // Provide the Poisson distribution function corresponding to the
750 // specified average rate (in time or space) mu of occurrences.
751 //
752 // p(n|mu) = probability for n occurrences in a certain time or space interval.
753 //
754 // Details can be found in the excelent textbook of Phil Gregory
755 // "Bayesian Logical Data Analysis for the Physical Sciences".
756 //
757 // Note : <n>=mu  Var(n)=mu
758  
759  TF1 poissdist("poissdist","exp(-[0])*pow([0],int(x))/TMath::Factorial(int(x))");
760  TString title="Poisson distribution for mu = ";
761  title+=mu;
762  poissdist.SetTitle(title.Data());
763  poissdist.SetParName(0,"mu");
764  poissdist.SetParameter(0,mu);
765
766  return poissdist;
767 }
768 ///////////////////////////////////////////////////////////////////////////
769 Double_t AliMath::Chi2Pvalue(Double_t chi2,Int_t ndf,Int_t sides,Int_t sigma,Int_t mode) const
770 {
771 // Computation of the P-value for a certain specified Chi-squared (chi2) value 
772 // for a Chi-squared distribution with ndf degrees of freedom.
773 //
774 // The P-value for a certain Chi-squared value chi2 corresponds to the fraction of
775 // repeatedly drawn equivalent samples from a certain population, which is expected
776 // to yield a Chi-squared value exceeding (less than) the value chi2 for an
777 // upper (lower) tail test in case a certain hypothesis is true.
778 //
779 // Further details can be found in the excellent textbook of Phil Gregory
780 // "Bayesian Logical Data Analysis for the Physical Sciences".
781 //
782 // Note : <Chi2>=ndf  Var(Chi2)=2*ndf
783 // 
784 // With the "sides" parameter a one-sided or two-sided test can be selected
785 // using either the upper or lower tail contents.
786 // In case of automatic upper/lower selection the decision is made on basis
787 // of the location of the input chi2 value w.r.t. <Chi2> of the distribution. 
788 //
789 // sides =  1 : One-sided test using the upper tail contents
790 //          2 : Two-sided test using the upper tail contents
791 //         -1 : One-sided test using the lower tail contents
792 //         -2 : Two-sided test using the lower tail contents
793 //          0 : One-sided test using the auto-selected upper or lower tail contents
794 //          3 : Two-sided test using the auto-selected upper or lower tail contents
795 //
796 // The argument "sigma" allows for the following return values :
797 //
798 // sigma = 0 : P-value is returned as the above specified fraction
799 //         1 : The difference chi2-<Chi2> expressed in units of sigma
800 //             Note : This difference may be negative.
801 //  
802 // According to the value of the parameter "mode" various algorithms
803 // can be selected.
804 //
805 // mode = 0 : Calculations are based on the incomplete gamma function.
806 //
807 //        1 : Same as for mode=0. However, in case ndf=1 an exact expression
808 //            based on the error function Erf() is used.
809 //
810 //        2 : Same as for mode=0. However, in case ndf>30 a Gaussian approximation
811 //            is used instead of the gamma function.
812 //
813 // The default values are sides=0, sigma=0 and mode=1.
814 //
815 //--- NvE 21-may-2005 Utrecht University
816
817  if (ndf<=0) return 0;
818
819  Double_t mean=ndf;
820
821  if (!sides) // Automatic one-sided test
822  {
823   sides=1;
824   if (chi2<mean) sides=-1;
825  }
826
827  if (sides==3) // Automatic two-sided test
828  {
829   sides=2;
830   if (chi2<mean) sides=-2;
831  }
832
833  Double_t val=0;
834  if (sigma) // P-value in units of sigma
835  {
836   Double_t s=sqrt(double(2*ndf));
837   val=(chi2-mean)/s;
838  }
839  else // P-value from tail contents
840  {
841   if (sides>0) // Upper tail
842   { 
843    val=Prob(chi2,ndf,mode);
844   }
845   else // Lower tail
846   {
847    val=1.-Prob(chi2,ndf,mode);
848   }
849  }
850
851  if (abs(sides)==2) val=val*2.;
852
853  return val;
854 }
855 ///////////////////////////////////////////////////////////////////////////
856 Double_t AliMath::StudentPvalue(Double_t t,Double_t ndf,Int_t sides,Int_t sigma) const
857 {
858 // Computation of the P-value for a certain specified Student's t value 
859 // for a Student's T distribution with ndf degrees of freedom.
860 //
861 // In a frequentist approach, the Student's T distribution is particularly
862 // useful in making inferences about the mean of an underlying population
863 // based on the data from a random sample.
864 //
865 // The P-value for a certain t value corresponds to the fraction of
866 // repeatedly drawn equivalent samples from a certain population, which is expected
867 // to yield a T value exceeding (less than) the value t for an upper (lower)
868 // tail test in case a certain hypothesis is true.
869 //
870 // Further details can be found in the excellent textbook of Phil Gregory
871 // "Bayesian Logical Data Analysis for the Physical Sciences".
872 //
873 // Note : <T>=0  Var(T)=ndf/(ndf-2)
874 // 
875 // With the "sides" parameter a one-sided or two-sided test can be selected
876 // using either the upper or lower tail contents.
877 // In case of automatic upper/lower selection the decision is made on basis
878 // of the location of the input t value w.r.t. <T> of the distribution. 
879 //
880 // sides =  1 : One-sided test using the upper tail contents
881 //          2 : Two-sided test using the upper tail contents
882 //         -1 : One-sided test using the lower tail contents
883 //         -2 : Two-sided test using the lower tail contents
884 //          0 : One-sided test using the auto-selected upper or lower tail contents
885 //          3 : Two-sided test using the auto-selected upper or lower tail contents
886 //
887 // The argument "sigma" allows for the following return values :
888 //
889 // sigma = 0 : P-value is returned as the above specified fraction
890 //         1 : The difference t-<T> expressed in units of sigma
891 //             Note : This difference may be negative and sigma
892 //                    is only defined for ndf>2.
893 //  
894 // The default values are sides=0 and sigma=0.
895 //
896 //--- NvE 21-may-2005 Utrecht University
897
898  if (ndf<=0) return 0;
899
900  Double_t mean=0;
901
902  if (!sides) // Automatic one-sided test
903  {
904   sides=1;
905   if (t<mean) sides=-1;
906  }
907
908  if (sides==3) // Automatic two-sided test
909  {
910   sides=2;
911   if (t<mean) sides=-2;
912  }
913
914  Double_t val=0;
915  if (sigma) // P-value in units of sigma
916  { 
917   if (ndf>2) // Sigma is only defined for ndf>2
918   {
919    Double_t s=sqrt(ndf/(ndf-2.));
920    val=t/s;
921   }
922  }
923  else // P-value from tail contents
924  {
925   if (sides>0) // Upper tail
926   {
927    val=1.-TMath::StudentI(t,ndf);
928   }
929   else // Lower tail
930   {
931    val=TMath::StudentI(t,ndf);
932   }
933  }
934
935  if (abs(sides)==2) val=val*2.;
936
937  return val;
938 }
939 ///////////////////////////////////////////////////////////////////////////
940 Double_t AliMath::FratioPvalue(Double_t f,Int_t ndf1,Int_t ndf2,Int_t sides,Int_t sigma) const
941 {
942 // Computation of the P-value for a certain specified F ratio f value 
943 // for an F (ratio) distribution with ndf1 and ndf2 degrees of freedom
944 // for the two samples X,Y respectively to be compared in the ratio X/Y.
945 //
946 // In a frequentist approach, the F (ratio) distribution is particularly useful
947 // in making inferences about the ratio of the variances of two underlying
948 // populations based on a measurement of the variance of a random sample taken
949 // from each one of the two populations.
950 // So the F test provides a means to investigate the degree of equality of
951 // two populations.
952 //
953 // The P-value for a certain f value corresponds to the fraction of
954 // repeatedly drawn equivalent samples from a certain population, which is expected
955 // to yield an F value exceeding (less than) the value f for an upper (lower)
956 // tail test in case a certain hypothesis is true.
957 //
958 // Further details can be found in the excellent textbook of Phil Gregory
959 // "Bayesian Logical Data Analysis for the Physical Sciences".
960 //
961 // Note : <F>=ndf2/(ndf2-2)  Var(F)=2*ndf2*ndf2*(ndf2+ndf1-2)/(ndf1*(ndf2-1)*(ndf2-1)*(ndf2-4))
962 // 
963 // With the "sides" parameter a one-sided or two-sided test can be selected
964 // using either the upper or lower tail contents.
965 // In case of automatic upper/lower selection the decision is made on basis
966 // of the location of the input f value w.r.t. <F> of the distribution. 
967 //
968 // sides =  1 : One-sided test using the upper tail contents
969 //          2 : Two-sided test using the upper tail contents
970 //         -1 : One-sided test using the lower tail contents
971 //         -2 : Two-sided test using the lower tail contents
972 //          0 : One-sided test using the auto-selected upper or lower tail contents
973 //          3 : Two-sided test using the auto-selected upper or lower tail contents
974 //
975 // The argument "sigma" allows for the following return values :
976 //
977 // sigma = 0 : P-value is returned as the above specified fraction
978 //         1 : The difference f-<F> expressed in units of sigma
979 //             Note : This difference may be negative and sigma
980 //                    is only defined for ndf2>4.
981 //  
982 // The default values are sides=0 and sigma=0.
983 //
984 //--- NvE 21-may-2005 Utrecht University
985
986  if (ndf1<=0 || ndf2<=0 || f<=0) return 0;
987
988  Double_t mean=double(ndf2/(ndf2-2));
989
990  if (!sides) // Automatic one-sided test
991  {
992   sides=1;
993   if (f<mean) sides=-1;
994  }
995
996  if (sides==3) // Automatic two-sided test
997  {
998   sides=2;
999   if (f<mean) sides=-2;
1000  }
1001
1002  Double_t val=0;
1003  if (sigma) // P-value in units of sigma
1004  { 
1005   // Sigma is only defined for ndf2>4
1006   if (ndf2>4)
1007   {
1008    Double_t s=sqrt(double(ndf2*ndf2*(2*ndf2+2*ndf1-4))/double(ndf1*pow(double(ndf2-1),2)*(ndf2-4)));
1009    val=(f-mean)/s;
1010   }
1011  }
1012  else // P-value from tail contents 
1013  {
1014   if (sides>0) // Upper tail
1015   {
1016    val=1.-TMath::FDistI(f,ndf1,ndf2);
1017   }
1018   else // Lower tail
1019   {
1020    val=TMath::FDistI(f,ndf1,ndf2);
1021   }
1022  }
1023
1024  if (abs(sides)==2) val=val*2.;
1025
1026  return val;
1027 }
1028 ///////////////////////////////////////////////////////////////////////////
1029 Double_t AliMath::BinomialPvalue(Int_t k,Int_t n,Double_t p,Int_t sides,Int_t sigma,Int_t mode) const
1030 {
1031 // Computation of the P-value for a certain specified number of successes k
1032 // for a Binomial distribution with n trials and success probability p.
1033 //
1034 // The P-value for a certain number of successes k corresponds to the fraction of
1035 // repeatedly drawn equivalent samples from a certain population, which is expected
1036 // to yield a number of successes exceeding (less than) the value k for an
1037 // upper (lower) tail test in case a certain hypothesis is true.
1038 //
1039 // Further details can be found in the excellent textbook of Phil Gregory
1040 // "Bayesian Logical Data Analysis for the Physical Sciences".
1041 //
1042 // Note : <K>=n*p  Var(K)=n*p*(1-p)
1043 // 
1044 // With the "sides" parameter a one-sided or two-sided test can be selected
1045 // using either the upper or lower tail contents.
1046 // In case of automatic upper/lower selection the decision is made on basis
1047 // of the location of the input k value w.r.t. <K> of the distribution. 
1048 //
1049 // sides =  1 : One-sided test using the upper tail contents
1050 //          2 : Two-sided test using the upper tail contents
1051 //         -1 : One-sided test using the lower tail contents
1052 //         -2 : Two-sided test using the lower tail contents
1053 //          0 : One-sided test using the auto-selected upper or lower tail contents
1054 //          3 : Two-sided test using the auto-selected upper or lower tail contents
1055 //
1056 // The argument "sigma" allows for the following return values :
1057 //
1058 // sigma = 0 : P-value is returned as the above specified fraction
1059 //         1 : The difference k-<K> expressed in units of sigma
1060 //             Note : This difference may be negative.
1061 //
1062 // mode = 0 : Incomplete Beta function will be used to calculate the tail content.
1063 //        1 : Straightforward summation of the Binomial terms is used.
1064 //
1065 // The Incomplete Beta function in general provides the most accurate values.
1066 //
1067 // The default values are sides=0, sigma=0 and mode=0.
1068 //
1069 //--- NvE 24-may-2005 Utrecht University
1070
1071  Double_t mean=double(n)*p;
1072
1073  if (!sides) // Automatic one-sided test
1074  {
1075   sides=1;
1076   if (k<mean) sides=-1;
1077  }
1078
1079  if (sides==3) // Automatic two-sided test
1080  {
1081   sides=2;
1082   if (k<mean) sides=-2;
1083  }
1084
1085  Double_t val=0;
1086
1087  if (sigma) // P-value in units of sigma
1088  {
1089   Double_t s=sqrt(double(n)*p*(1.-p));
1090   val=(double(k)-mean)/s;
1091  }
1092  else // P-value from tail contents
1093  {
1094   if (sides>0)
1095   {
1096    if (!mode) // Use Incomplete Beta function
1097    {
1098     val=TMath::BetaIncomplete(p,k+1,n-k);
1099    }
1100    else // Use straightforward summation
1101    {
1102     for (Int_t i=k+1; i<=n; i++)
1103     {
1104      val+=TMath::Binomial(n,i)*pow(p,i)*pow(1.-p,n-i);
1105     }
1106    }
1107   }
1108   else
1109   {
1110    if (!mode) // Use Incomplete Beta function
1111    {
1112     val=1.-TMath::BetaIncomplete(p,k+1,n-k);
1113    }
1114    else
1115    {
1116     for (Int_t j=0; j<=k; j++)
1117     {
1118      val+=TMath::Binomial(n,j)*pow(p,j)*pow(1.-p,n-j);
1119     }
1120    }
1121   }
1122  }
1123
1124  if (abs(sides)==2) val=val*2.;
1125
1126  return val;
1127 }
1128 ///////////////////////////////////////////////////////////////////////////
1129 Double_t AliMath::PoissonPvalue(Int_t k,Double_t mu,Int_t sides,Int_t sigma) const
1130 {
1131 // Computation of the P-value for a certain specified number of occurrences k
1132 // for a Poisson distribution with a specified average rate (in time or space)
1133 // mu of occurrences.
1134 //
1135 // The P-value for a certain number of occurrences k corresponds to the fraction of
1136 // repeatedly drawn equivalent samples from a certain population, which is expected
1137 // to yield a number of occurrences exceeding (less than) the value k for an
1138 // upper (lower) tail test in case a certain hypothesis is true.
1139 //
1140 // Further details can be found in the excellent textbook of Phil Gregory
1141 // "Bayesian Logical Data Analysis for the Physical Sciences".
1142 //
1143 // Note : <K>=mu  Var(K)=mu
1144 // 
1145 // With the "sides" parameter a one-sided or two-sided test can be selected
1146 // using either the upper or lower tail contents.
1147 // In case of automatic upper/lower selection the decision is made on basis
1148 // of the location of the input k value w.r.t. <K> of the distribution. 
1149 //
1150 // sides =  1 : One-sided test using the upper tail contents
1151 //          2 : Two-sided test using the upper tail contents
1152 //         -1 : One-sided test using the lower tail contents
1153 //         -2 : Two-sided test using the lower tail contents
1154 //          0 : One-sided test using the auto-selected upper or lower tail contents
1155 //          3 : Two-sided test using the auto-selected upper or lower tail contents
1156 //
1157 // The argument "sigma" allows for the following return values :
1158 //
1159 // sigma = 0 : P-value is returned as the above specified fraction
1160 //         1 : The difference k-<K> expressed in units of sigma
1161 //             Note : This difference may be negative.
1162 //
1163 // The default values are sides=0 and sigma=0.
1164 //
1165 // Note : The tail contents are given by the incomplete Gamma function P(a,x).
1166 //        Lower tail contents = 1-P(k,mu)
1167 //        Upper tail contents = P(k,mu)
1168 //
1169 //--- NvE 24-may-2005 Utrecht University
1170
1171  Double_t mean=mu;
1172
1173  if (!sides) // Automatic one-sided test
1174  {
1175   sides=1;
1176   if (k<mean) sides=-1;
1177  }
1178
1179  if (sides==3) // Automatic two-sided test
1180  {
1181   sides=2;
1182   if (k<mean) sides=-2;
1183  }
1184
1185  Double_t val=0;
1186
1187  if (sigma) // P-value in units of sigma
1188  {
1189   Double_t s=sqrt(mu);
1190   val=(double(k)-mean)/s;
1191  }
1192  else // P-value from tail contents
1193  {
1194   if (sides>0) // Upper tail
1195   {
1196    val=Gamma(k,mu);
1197   }
1198   else // Lower tail
1199   {
1200    val=1.-Gamma(k,mu);
1201   }
1202  }
1203
1204  if (abs(sides)==2) val=val*2.;
1205
1206  return val;
1207 }
1208 ///////////////////////////////////////////////////////////////////////////
1209 Double_t AliMath::NegBinomialPvalue(Int_t n,Int_t k,Double_t p,Int_t sides,Int_t sigma,Int_t mode) const
1210 {
1211 // Computation of the P-value for a certain specified number of trials n
1212 // for a Negative Binomial distribution where exactly k successes are to
1213 // be reached which have each a probability p.
1214 //
1215 // The P-value for a certain number of trials n corresponds to the fraction of
1216 // repeatedly drawn equivalent samples from a certain population, which is expected
1217 // to yield a number of trials exceeding (less than) the value n for an
1218 // upper (lower) tail test in case a certain hypothesis is true.
1219 //
1220 // Further details can be found in the excellent textbook of Phil Gregory
1221 // "Bayesian Logical Data Analysis for the Physical Sciences".
1222 //
1223 // Note : <N>=k*(1-p)/p  Var(N)=k*(1-p)/(p*p) 
1224 // 
1225 // With the "sides" parameter a one-sided or two-sided test can be selected
1226 // using either the upper or lower tail contents.
1227 // In case of automatic upper/lower selection the decision is made on basis
1228 // of the location of the input n value w.r.t. <N> of the distribution. 
1229 //
1230 // sides =  1 : One-sided test using the upper tail contents
1231 //          2 : Two-sided test using the upper tail contents
1232 //         -1 : One-sided test using the lower tail contents
1233 //         -2 : Two-sided test using the lower tail contents
1234 //          0 : One-sided test using the auto-selected upper or lower tail contents
1235 //          3 : Two-sided test using the auto-selected upper or lower tail contents
1236 //
1237 // The argument "sigma" allows for the following return values :
1238 //
1239 // sigma = 0 : P-value is returned as the above specified fraction
1240 //         1 : The difference n-<N> expressed in units of sigma
1241 //             Note : This difference may be negative.
1242 //
1243 // mode = 0 : Incomplete Beta function will be used to calculate the tail content.
1244 //        1 : Straightforward summation of the Negative Binomial terms is used.
1245 //
1246 // The Incomplete Beta function in general provides the most accurate values.
1247 //
1248 // The default values are sides=0, sigma=0 and mode=0.
1249 //
1250 //--- NvE 24-may-2005 Utrecht University
1251
1252  Double_t mean=double(k)*(1.-p)/p;
1253
1254  if (!sides) // Automatic one-sided test
1255  {
1256   sides=1;
1257   if (n<mean) sides=-1;
1258  }
1259
1260  if (sides==3) // Automatic two-sided test
1261  {
1262   sides=2;
1263   if (n<mean) sides=-2;
1264  }
1265
1266  Double_t val=0;
1267
1268  if (sigma) // P-value in units of sigma
1269  {
1270   Double_t s=sqrt(double(k)*(1.-p)/(p*p));
1271   val=(double(n)-mean)/s;
1272  }
1273  else // P-value from tail contents
1274  {
1275   if (sides>0) // Upper tail
1276   {
1277    if (!mode) // Use Incomplete Beta function
1278    {
1279     val=1.-TMath::BetaIncomplete(p,k,n-k);
1280    }
1281    else // Use straightforward summation
1282    {
1283     for (Int_t i=1; i<n; i++)
1284     {
1285      val+=TMath::Binomial(i-1,k-1)*pow(p,k)*pow(1.-p,i-k);
1286     }
1287     val=1.-val;
1288    }
1289   }
1290   else // Lower tail
1291   {
1292    if (!mode) // Use Incomplete Beta function
1293    {
1294     val=TMath::BetaIncomplete(p,k,n-k);
1295    }
1296    else
1297    {
1298     for (Int_t j=1; j<n; j++)
1299     {
1300      val+=TMath::Binomial(j-1,k-1)*pow(p,k)*pow(1.-p,j-k);
1301     }
1302    }
1303   }
1304  }
1305
1306  if (abs(sides)==2) val=val*2.;
1307
1308  return val;
1309 }
1310 ///////////////////////////////////////////////////////////////////////////
1311 Double_t AliMath::Nfac(Int_t n,Int_t mode) const
1312 {
1313 // Compute n!.
1314 // The algorithm can be selected by the "mode" input argument.
1315 //
1316 // mode : 0 ==> Calculation by means of straightforward multiplication
1317 //      : 1 ==> Calculation by means of Stirling's approximation
1318 //      : 2 ==> Calculation by means of n!=Gamma(n+1)
1319 //
1320 // For large n the calculation modes 1 and 2 will in general be faster.
1321 // By default mode=0 is used.
1322 // For n<0 the value 0 will be returned.
1323 //
1324 // Note : Because of Double_t value overflow the maximum value is n=170.
1325 //
1326 //--- NvE 20-jan-2007 Utrecht University
1327
1328  if (n<0) return 0;
1329  if (n==0) return 1;
1330
1331  Double_t twopi=2.*acos(-1.);
1332  Double_t z=0;
1333  Double_t nfac=1;
1334  Int_t i=n;
1335  
1336  switch (mode)
1337  {
1338   case 0: // Straightforward multiplication
1339    while (i>1)
1340    {
1341     nfac*=Double_t(i);
1342     i--;
1343    }
1344    break;
1345
1346   case 1: // Stirling's approximation 
1347    z=n;
1348    nfac=sqrt(twopi)*pow(z,z+0.5)*exp(-z)*(1.+1./(12.*z));
1349    break;
1350
1351   case 2: // Use of Gamma(n+1)
1352    z=n+1;
1353    nfac=Gamma(z);
1354    break;
1355
1356   default:
1357    nfac=0;
1358    break;
1359  }
1360
1361  return nfac;
1362 }
1363 ///////////////////////////////////////////////////////////////////////////
1364 Double_t AliMath::LnNfac(Int_t n,Int_t mode) const
1365 {
1366 // Compute ln(n!).
1367 // The algorithm can be selected by the "mode" input argument.
1368 //
1369 // mode : 0 ==>  Calculation via evaluation of n! followed by taking ln(n!)
1370 //      : 1 ==>  Calculation via Stirling's approximation ln(n!)=0.5*ln(2*pi)+(n+0.5)*ln(n)-n+1/(12*n)
1371 //      : 2 ==>  Calculation by means of ln(n!)=LnGamma(n+1)
1372 //
1373 // Note : Because of Double_t value overflow the maximum value is n=170 for mode=0.
1374 //
1375 // For mode=2 rather accurate results are obtained for both small and large n.
1376 // By default mode=2 is used.
1377 // For n<1 the value 0 will be returned.
1378 //
1379 //--- NvE 20-jan-2007 Utrecht University
1380
1381  if (n<=1) return 0;
1382
1383  Double_t twopi=2.*acos(-1.);
1384  Double_t z=0;
1385  Double_t lognfac=0;
1386  
1387  switch (mode)
1388  {
1389   case 0: // Straightforward ln(n!)
1390    z=Nfac(n);
1391    lognfac=log(z);
1392    break;
1393
1394   case 1: // Stirling's approximation 
1395    z=n;
1396    lognfac=0.5*log(twopi)+(z+0.5)*log(z)-z+1./(12.*z);
1397    break;
1398
1399   case 2: // Use of LnGamma(n+1)
1400    z=n+1;
1401    lognfac=LnGamma(z);
1402    break;
1403
1404   default:
1405    lognfac=0;
1406    break;
1407  }
1408
1409  return lognfac;
1410 }
1411 ///////////////////////////////////////////////////////////////////////////
1412 Double_t AliMath::LogNfac(Int_t n,Int_t mode) const
1413 {
1414 // Compute log_10(n!).
1415 // First ln(n!) is evaluated via invokation of LnNfac(n,mode).
1416 // Then the algorithm log_10(z)=ln(z)*log_10(e) is used.
1417 //
1418 //--- NvE 20-jan-2007 Utrecht University
1419
1420  Double_t e=exp(1.);
1421
1422  Double_t val=LnNfac(n,mode);
1423  val*=log10(e);
1424
1425  return val;
1426 }
1427 ///////////////////////////////////////////////////////////////////////////
1428 Double_t AliMath::PsiValue(Int_t m,Int_t* n,Double_t* p,Int_t f) const
1429 {
1430 // Provide the Bayesian Psi value of observations of a counting experiment
1431 // w.r.t. a Bernoulli class hypothesis B_m.
1432 // The hypothesis B_m represents a counting experiment with m different
1433 // possible outcomes and is completely defined by the probabilities
1434 // of the various outcomes (and the requirement that the sum of all these
1435 // probabilities equals 1).
1436 //
1437 // The Psi value provides (in dB scale) the amount of support that the
1438 // data can maximally give to any Bernoulli class hypothesis different
1439 // from the currently specified B_m.
1440 //
1441 // To be specific : Psi=-10*log[p(D|B_m I)]
1442 //
1443 // where p(D|B_m I) represents the likelihood of the data D under the condition
1444 // that B_m and some prior I are true.
1445 //
1446 // A Psi value of zero indicates a perfect match between the observations
1447 // and the specified hypothesis.
1448 // Further mathematical details can be found in astro-ph/0702029.
1449 //
1450 // m : The number of different possible outcomes of the counting experiment
1451 // n : The observed number of different outcomes
1452 // p : The probabilities of the different outcomes according to the hypothesis
1453 // f : Flag to indicate the use of a frequentist (Stirling) approximation (f=1)
1454 //     or the exact Bayesian expression (f=0).
1455 //
1456 // In case no probabilities are given (i.e. pk=0), a uniform distribution
1457 // is assumed.
1458 //
1459 // The default values are pk=0 and f=0.
1460 //
1461 // In the case of inconsistent input, a Psi value of -1 is returned.
1462 //
1463 //--- NvE 03-oct-2007 Utrecht University
1464
1465  Double_t psi=-1;
1466
1467  if (m<=0 || !n) return psi;
1468
1469  Int_t ntot=0;
1470  for (Int_t j=0; j<m; j++)
1471  {
1472   if (n[j]>0) ntot+=n[j];
1473  }
1474
1475  psi=0;
1476  Double_t pk=1./float(m); // Prob. of getting outcome k for a uniform distr.
1477  for (Int_t i=0; i<m; i++)
1478  {
1479   if (p) pk=p[i]; // Probabilities from a specific B_m hypothesis
1480   if (n[i]>0 && pk>0)
1481   {
1482    if (!f) // Exact Bayesian expression
1483    {
1484     psi+=double(n[i])*log10(pk)-LogNfac(n[i]);
1485    }
1486    else // Frequentist (Stirling) approximation
1487    {
1488     if (ntot>0) psi+=double(n[i])*log10(n[i]/(ntot*pk));
1489    }
1490   }
1491  }
1492
1493  if (!f) // Exact Bayesian expression
1494  {
1495   psi+=LogNfac(ntot);
1496   psi*=-10.;
1497  }
1498  else // Frequentist (Stirling) approximation
1499  {
1500   psi*=10.;
1501  }
1502
1503  return psi;
1504 }
1505 ///////////////////////////////////////////////////////////////////////////
1506 Double_t AliMath::PsiValue(TH1F* his,TH1F* hyp,TF1* pdf,Int_t f) const
1507 {
1508 // Provide the Bayesian Psi value of observations of a counting experiment
1509 // (in histogram format) w.r.t. a Bernoulli class hypothesis B_m.
1510 // The hypothesis B_m represents a counting experiment with m different
1511 // possible outcomes and is completely defined by the probabilities
1512 // of the various outcomes (and the requirement that the sum of all these
1513 // probabilities equals 1).
1514 // The specification of a hypothesis B_m can be provided either in
1515 // histogram format (hyp) or via a probability distribution function (pdf),
1516 // as outlined below.
1517 // Note : The pdf does not need to be normalised.
1518 //
1519 // The Psi value provides (in dB scale) the amount of support that the
1520 // data can maximally give to any Bernoulli class hypothesis different
1521 // from the currently specified B_m.
1522 //
1523 // To be specific : Psi=-10*log[p(D|B_m I)]
1524 //
1525 // where p(D|B_m I) represents the likelihood of the data D under the condition
1526 // that B_m and some prior I are true.
1527 //
1528 // A Psi value of zero indicates a perfect match between the observations
1529 // and the specified hypothesis.
1530 // Further mathematical details can be found in astro-ph/0702029.
1531 //
1532 // his : The experimental observations in histogram format
1533 // hyp : Hypothetical observations according to some hypothesis
1534 // pdf : Probability distribution function for the hypothesis
1535 // f   : Flag to indicate the use of a frequentist (Stirling) approximation (f=1)
1536 //       or the exact Bayesian expression (f=0).
1537 //
1538 // In case no hypothesis is specified (i.e. hyp=0 and pdf=0), a uniform
1539 // background distribution is assumed.
1540 //
1541 // Default values are : hyp=0, pdf=0 and f=0.
1542 //
1543 // In the case of inconsistent input, a Psi value of -1 is returned.
1544 //
1545 //--- NvE 03-oct-2007 Utrecht University
1546
1547  Double_t psi=-1;
1548
1549  if (!his) return psi;
1550
1551  TAxis* xaxis=his->GetXaxis();
1552  Double_t xmin=xaxis->GetXmin();
1553  Double_t xmax=xaxis->GetXmax();
1554  Double_t range=xmax-xmin;
1555  Int_t nbins=his->GetNbinsX();
1556  Double_t nensig=his->GetEntries();
1557
1558  if (nbins<=0 || nensig<=0 || range<=0) return psi;
1559
1560  Int_t* n=new Int_t[nbins];
1561  Double_t* p=new Double_t[nbins]; 
1562  Int_t nk;
1563  Double_t pk;
1564
1565  // Uniform hypothesis distribution
1566  if (!hyp && !pdf)
1567  {
1568   for (Int_t i=1; i<=nbins; i++)
1569   {
1570    nk=int(his->GetBinContent(i));
1571    pk=his->GetBinWidth(i)/range;
1572    n[i-1]=0;
1573    p[i-1]=0;
1574    if (nk>0) n[i-1]=nk;
1575    if (pk>0) p[i-1]=pk;
1576   }
1577   psi=PsiValue(nbins,n,p,f);
1578  }
1579
1580  // Hypothesis specified via a pdf
1581  if (pdf && !hyp)
1582  {
1583   pdf->SetRange(xmin,xmax);
1584   Double_t ftot=pdf->Integral(xmin,xmax);
1585   if (ftot>0)
1586   {
1587    Double_t x1,x2;
1588    for (Int_t ipdf=1; ipdf<=nbins; ipdf++)
1589    {
1590     nk=int(his->GetBinContent(ipdf));
1591     x1=his->GetBinLowEdge(ipdf);
1592     x2=x1+his->GetBinWidth(ipdf);
1593     pk=pdf->Integral(x1,x2)/ftot;
1594     n[ipdf-1]=0;
1595     p[ipdf-1]=0;
1596     if (nk>0) n[ipdf-1]=nk;
1597     if (pk>0) p[ipdf-1]=pk;
1598    }
1599    psi=PsiValue(nbins,n,p,f);
1600   }
1601  }
1602
1603  // Hypothesis specified via a histogram
1604  if (hyp && !pdf)
1605  {
1606   TH1F* href=new TH1F(*his);
1607   href->Reset();
1608   Double_t nenhyp=hyp->GetEntries();
1609   Float_t x,y;
1610   for (Int_t ihyp=1; ihyp<=hyp->GetNbinsX(); ihyp++)
1611   {
1612    x=hyp->GetBinCenter(ihyp);
1613    y=hyp->GetBinContent(ihyp);
1614    href->Fill(x,y);
1615   }
1616   for (Int_t j=1; j<=nbins; j++)
1617   {
1618    nk=int(his->GetBinContent(j));
1619    pk=href->GetBinContent(j)/nenhyp;
1620    n[j-1]=0;
1621    p[j-1]=0;
1622    if (nk>0) n[j-1]=nk;
1623    if (pk>0) p[j-1]=pk;
1624   }
1625   psi=PsiValue(nbins,n,p,f);
1626   delete href;
1627  }
1628
1629  delete [] n;
1630  delete [] p;
1631
1632  return psi;
1633 }
1634 ///////////////////////////////////////////////////////////////////////////
1635 Double_t AliMath::Chi2Value(Int_t m,Int_t* n,Double_t* p,Int_t* ndf) const
1636 {
1637 // Provide the frequentist chi-squared value of observations of a counting
1638 // experiment w.r.t. a Bernoulli class hypothesis B_m.
1639 // The hypothesis B_m represents a counting experiment with m different
1640 // possible outcomes and is completely defined by the probabilities
1641 // of the various outcomes (and the requirement that the sum of all these
1642 // probabilities equals 1).
1643 //
1644 // Further mathematical details can be found in astro-ph/0702029.
1645 //
1646 // m   : The number of different possible outcomes of the counting experiment
1647 // n   : The observed number of different outcomes
1648 // p   : The probabilities of the different outcomes according to the hypothesis
1649 // ndf : The returned number of degrees of freedom
1650 //
1651 // Note : Both the arrays "n" and (when provided) "p" should be of dimension "m".
1652 //
1653 // In case no probabilities are given (i.e. pk=0), a uniform distribution
1654 // is assumed.
1655 //
1656 // The default values are pk=0 and ndf=0.
1657 //
1658 // In the case of inconsistent input, a chi-squared and ndf value of -1 is returned.
1659 //
1660 //--- NvE 03-oct-2007 Utrecht University
1661
1662  Double_t chi=-1;
1663  if (ndf) *ndf=-1;
1664
1665  if (m<=0 || !n) return chi;
1666
1667  Int_t ntot=0;
1668  for (Int_t j=0; j<m; j++)
1669  {
1670   if (n[j]>0) ntot+=n[j];
1671  }
1672
1673  chi=0;
1674  Double_t pk=1./float(m); // Prob. of getting outcome k for a uniform distr.
1675  for (Int_t i=0; i<m; i++)
1676  {
1677   if (p) pk=p[i]; // Probabilities from a specific B_m hypothesis
1678   if (n[i]>0 && pk>0 && ntot>0)
1679   {
1680    chi+=pow(double(n[i])-double(ntot)*pk,2)/(ntot*pk);
1681    if (ndf) (*ndf)=(*ndf)+1;
1682   }
1683  }
1684
1685  return chi;
1686 }
1687 ///////////////////////////////////////////////////////////////////////////
1688 Double_t AliMath::Chi2Value(TH1F* his,TH1F* hyp,TF1* pdf,Int_t* ndf) const
1689 {
1690 // Provide the frequentist chi-squared value of observations of a counting
1691 // experiment (in histogram format) w.r.t. a Bernoulli class hypothesis B_m.
1692 // The hypothesis B_m represents a counting experiment with m different
1693 // possible outcomes and is completely defined by the probabilities
1694 // of the various outcomes (and the requirement that the sum of all these
1695 // probabilities equals 1).
1696 // The specification of a hypothesis B_m can be provided either in
1697 // histogram format (hyp) or via a probability distribution function (pdf),
1698 // as outlined below.
1699 // Note : The pdf does not need to be normalised.
1700 //
1701 // Further mathematical details can be found in astro-ph/0702029.
1702 //
1703 // his : The experimental observations in histogram format
1704 // hyp : Hypothetical observations according to some hypothesis
1705 // pdf : Probability distribution function for the hypothesis
1706 // ndf : The returned number of degrees of freedom
1707 //
1708 // In case no hypothesis is specified (i.e. hyp=0 and pdf=0), a uniform
1709 // background distribution is assumed.
1710 //
1711 // Default values are : hyp=0, pdf=0 and ndf=0.
1712 //
1713 // In the case of inconsistent input, a chi-squared and ndf value of -1 is returned.
1714 //
1715 //--- NvE 03-oct-2007 Utrecht University
1716
1717  Double_t chi=-1;
1718
1719  if (!his) return chi;
1720
1721  TAxis* xaxis=his->GetXaxis();
1722  Double_t xmin=xaxis->GetXmin();
1723  Double_t xmax=xaxis->GetXmax();
1724  Double_t range=xmax-xmin;
1725  Int_t nbins=his->GetNbinsX();
1726  Double_t nensig=his->GetEntries();
1727
1728  if (nbins<=0 || nensig<=0 || range<=0) return chi;
1729
1730  Int_t* n=new Int_t[nbins];
1731  Double_t* p=new Double_t[nbins]; 
1732  Int_t nk;
1733  Double_t pk;
1734
1735  // Uniform hypothesis distribution
1736  if (!hyp && !pdf)
1737  {
1738   for (Int_t i=1; i<=nbins; i++)
1739   {
1740    nk=int(his->GetBinContent(i));
1741    pk=his->GetBinWidth(i)/range;
1742    n[i-1]=0;
1743    p[i-1]=0;
1744    if (nk>0) n[i-1]=nk;
1745    if (pk>0) p[i-1]=pk;
1746   }
1747   chi=Chi2Value(nbins,n,p,ndf);
1748  }
1749
1750  // Hypothesis specified via a pdf
1751  if (pdf && !hyp)
1752  {
1753   pdf->SetRange(xmin,xmax);
1754   Double_t ftot=pdf->Integral(xmin,xmax);
1755   if (ftot>0)
1756   {
1757    Double_t x1,x2;
1758    for (Int_t ipdf=1; ipdf<=nbins; ipdf++)
1759    {
1760     nk=int(his->GetBinContent(ipdf));
1761     x1=his->GetBinLowEdge(ipdf);
1762     x2=x1+his->GetBinWidth(ipdf);
1763     pk=pdf->Integral(x1,x2)/ftot;
1764     n[ipdf-1]=0;
1765     p[ipdf-1]=0;
1766     if (nk>0) n[ipdf-1]=nk;
1767     if (pk>0) p[ipdf-1]=pk;
1768    }
1769    chi=Chi2Value(nbins,n,p,ndf);
1770   }
1771  }
1772
1773  // Hypothesis specified via a histogram
1774  if (hyp && !pdf)
1775  {
1776   TH1F* href=new TH1F(*his);
1777   href->Reset();
1778   Double_t nenhyp=hyp->GetEntries();
1779   Float_t x,y;
1780   for (Int_t ihyp=1; ihyp<=hyp->GetNbinsX(); ihyp++)
1781   {
1782    x=hyp->GetBinCenter(ihyp);
1783    y=hyp->GetBinContent(ihyp);
1784    href->Fill(x,y);
1785   }
1786   for (Int_t j=1; j<=nbins; j++)
1787   {
1788    nk=int(his->GetBinContent(j));
1789    pk=href->GetBinContent(j)/nenhyp;
1790    n[j-1]=0;
1791    p[j-1]=0;
1792    if (nk>0) n[j-1]=nk;
1793    if (pk>0) p[j-1]=pk;
1794   }
1795   chi=Chi2Value(nbins,n,p,ndf);
1796   delete href;
1797  }
1798
1799  delete [] n;
1800  delete [] p;
1801
1802  return chi;
1803 }
1804 ///////////////////////////////////////////////////////////////////////////