]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RALICE/AliMath.cxx
4815c385393625515057ae103ed2370f1beacea7
[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  Double_t pi=acos(-1.);
654
655  TF1 tdist("tdist","(TMath::Gamma(([0]+1.)/2.)/(sqrt(pi*[0])*TMath::Gamma([0]/2.)))*pow(1.+x*x/[0],-([0]+1.)/2.)");
656  TString title="Student's t distribution for ndf = ";
657  title+=ndf; 
658  tdist.SetTitle(title.Data());
659  tdist.SetParName(0,"ndf");
660  tdist.SetParameter(0,ndf);
661
662  return tdist;
663 }
664 ///////////////////////////////////////////////////////////////////////////
665 TF1 AliMath::FratioDist(Int_t ndf1,Int_t ndf2) const
666 {
667 // Provide the F (ratio) distribution function corresponding to the
668 // specified ndf1 and ndf2 degrees of freedom of the two samples.
669 //
670 // In a frequentist approach, the F (ratio) distribution is particularly useful
671 // in making inferences about the ratio of the variances of two underlying
672 // populations based on a measurement of the variance of a random sample taken
673 // from each one of the two populations.
674 // So the F test provides a means to investigate the degree of equality of
675 // two populations.
676 //
677 // Details can be found in the excellent textbook of Phil Gregory
678 // "Bayesian Logical Data Analysis for the Physical Sciences".
679 //
680 // Note : <F>=ndf2/(ndf2-2)  Var(F)=2*ndf2*ndf2*(ndf2+ndf1-2)/(ndf1*(ndf2-1)*(ndf2-1)*(ndf2-4))
681  
682  TF1 fdist("fdist",
683  "(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.)");
684  TString title="F (ratio) distribution for ndf1 = ";
685  title+=ndf1;
686  title+=" ndf2 = "; 
687  title+=ndf2;
688  fdist.SetTitle(title.Data());
689  fdist.SetParName(0,"ndf1");
690  fdist.SetParameter(0,ndf1);
691  fdist.SetParName(1,"ndf2");
692  fdist.SetParameter(1,ndf2);
693
694  return fdist;
695 }
696 ///////////////////////////////////////////////////////////////////////////
697 TF1 AliMath::BinomialDist(Int_t n,Double_t p) const
698 {
699 // Provide the Binomial distribution function corresponding to the
700 // specified number of trials n and probability p of success.
701 //
702 // p(k|n,p) = probability to obtain exactly k successes in n trials. 
703 //
704 // Details can be found in the excelent textbook of Phil Gregory
705 // "Bayesian Logical Data Analysis for the Physical Sciences".
706 //
707 // Note : <k>=n*p  Var(k)=n*p*(1-p)
708  
709  TF1 bindist("bindist","TMath::Binomial(int([0]),int(x))*pow([1],int(x))*pow(1.-[1],int([0])-int(x))");
710  TString title="Binomial distribution for n = ";
711  title+=n;
712  title+=" p = "; 
713  title+=p;
714  bindist.SetTitle(title.Data());
715  bindist.SetParName(0,"n");
716  bindist.SetParameter(0,n);
717  bindist.SetParName(1,"p");
718  bindist.SetParameter(1,p);
719
720  return bindist;
721 }
722 ///////////////////////////////////////////////////////////////////////////
723 TF1 AliMath::NegBinomialDist(Int_t k,Double_t p) const
724 {
725 // Provide the Negative Binomial distribution function corresponding to the
726 // specified number of successes k and probability p of success.
727 //
728 // p(n|k,p) = probability to have reached k successes after n trials.
729 //
730 // Details can be found in the excelent textbook of Phil Gregory
731 // "Bayesian Logical Data Analysis for the Physical Sciences".
732 //
733 // Note : <n>=k*(1-p)/p  Var(n)=k*(1-p)/(p*p) 
734  
735  TF1 negbindist("negbindist","TMath::Binomial(int(x)-1,int([0])-1)*pow([1],int([0]))*pow(1.-[1],int(x)-int([0]))");
736  TString title="Negative Binomial distribution for k = ";
737  title+=k;
738  title+=" p = "; 
739  title+=p;
740  negbindist.SetTitle(title.Data());
741  negbindist.SetParName(0,"k");
742  negbindist.SetParameter(0,k);
743  negbindist.SetParName(1,"p");
744  negbindist.SetParameter(1,p);
745
746  return negbindist;
747 }
748 ///////////////////////////////////////////////////////////////////////////
749 TF1 AliMath::PoissonDist(Double_t mu) const
750 {
751 // Provide the Poisson distribution function corresponding to the
752 // specified average rate (in time or space) mu of occurrences.
753 //
754 // p(n|mu) = probability for n occurrences in a certain time or space interval.
755 //
756 // Details can be found in the excelent textbook of Phil Gregory
757 // "Bayesian Logical Data Analysis for the Physical Sciences".
758 //
759 // Note : <n>=mu  Var(n)=mu
760  
761  TF1 poissdist("poissdist","exp(-[0])*pow([0],int(x))/TMath::Factorial(int(x))");
762  TString title="Poisson distribution for mu = ";
763  title+=mu;
764  poissdist.SetTitle(title.Data());
765  poissdist.SetParName(0,"mu");
766  poissdist.SetParameter(0,mu);
767
768  return poissdist;
769 }
770 ///////////////////////////////////////////////////////////////////////////
771 Double_t AliMath::Chi2Pvalue(Double_t chi2,Int_t ndf,Int_t sides,Int_t sigma,Int_t mode) const
772 {
773 // Computation of the P-value for a certain specified Chi-squared (chi2) value 
774 // for a Chi-squared distribution with ndf degrees of freedom.
775 //
776 // The P-value for a certain Chi-squared value chi2 corresponds to the fraction of
777 // repeatedly drawn equivalent samples from a certain population, which is expected
778 // to yield a Chi-squared value exceeding (less than) the value chi2 for an
779 // upper (lower) tail test in case a certain hypothesis is true.
780 //
781 // Further details can be found in the excellent textbook of Phil Gregory
782 // "Bayesian Logical Data Analysis for the Physical Sciences".
783 //
784 // Note : <Chi2>=ndf  Var(Chi2)=2*ndf
785 // 
786 // With the "sides" parameter a one-sided or two-sided test can be selected
787 // using either the upper or lower tail contents.
788 // In case of automatic upper/lower selection the decision is made on basis
789 // of the location of the input chi2 value w.r.t. <Chi2> of the distribution. 
790 //
791 // sides =  1 : One-sided test using the upper tail contents
792 //          2 : Two-sided test using the upper tail contents
793 //         -1 : One-sided test using the lower tail contents
794 //         -2 : Two-sided test using the lower tail contents
795 //          0 : One-sided test using the auto-selected upper or lower tail contents
796 //          3 : Two-sided test using the auto-selected upper or lower tail contents
797 //
798 // The argument "sigma" allows for the following return values :
799 //
800 // sigma = 0 : P-value is returned as the above specified fraction
801 //         1 : The difference chi2-<Chi2> expressed in units of sigma
802 //             Note : This difference may be negative.
803 //  
804 // According to the value of the parameter "mode" various algorithms
805 // can be selected.
806 //
807 // mode = 0 : Calculations are based on the incomplete gamma function.
808 //
809 //        1 : Same as for mode=0. However, in case ndf=1 an exact expression
810 //            based on the error function Erf() is used.
811 //
812 //        2 : Same as for mode=0. However, in case ndf>30 a Gaussian approximation
813 //            is used instead of the gamma function.
814 //
815 // The default values are sides=0, sigma=0 and mode=1.
816 //
817 //--- NvE 21-may-2005 Utrecht University
818
819  if (ndf<=0) return 0;
820
821  Double_t mean=ndf;
822
823  if (!sides) // Automatic one-sided test
824  {
825   sides=1;
826   if (chi2<mean) sides=-1;
827  }
828
829  if (sides==3) // Automatic two-sided test
830  {
831   sides=2;
832   if (chi2<mean) sides=-2;
833  }
834
835  Double_t val=0;
836  if (sigma) // P-value in units of sigma
837  {
838   Double_t s=sqrt(double(2*ndf));
839   val=(chi2-mean)/s;
840  }
841  else // P-value from tail contents
842  {
843   if (sides>0) // Upper tail
844   { 
845    val=Prob(chi2,ndf,mode);
846   }
847   else // Lower tail
848   {
849    val=1.-Prob(chi2,ndf,mode);
850   }
851  }
852
853  if (abs(sides)==2) val=val*2.;
854
855  return val;
856 }
857 ///////////////////////////////////////////////////////////////////////////
858 Double_t AliMath::StudentPvalue(Double_t t,Double_t ndf,Int_t sides,Int_t sigma) const
859 {
860 // Computation of the P-value for a certain specified Student's t value 
861 // for a Student's T distribution with ndf degrees of freedom.
862 //
863 // In a frequentist approach, the Student's T distribution is particularly
864 // useful in making inferences about the mean of an underlying population
865 // based on the data from a random sample.
866 //
867 // The P-value for a certain t value corresponds to the fraction of
868 // repeatedly drawn equivalent samples from a certain population, which is expected
869 // to yield a T value exceeding (less than) the value t for an upper (lower)
870 // tail test in case a certain hypothesis is true.
871 //
872 // Further details can be found in the excellent textbook of Phil Gregory
873 // "Bayesian Logical Data Analysis for the Physical Sciences".
874 //
875 // Note : <T>=0  Var(T)=ndf/(ndf-2)
876 // 
877 // With the "sides" parameter a one-sided or two-sided test can be selected
878 // using either the upper or lower tail contents.
879 // In case of automatic upper/lower selection the decision is made on basis
880 // of the location of the input t value w.r.t. <T> of the distribution. 
881 //
882 // sides =  1 : One-sided test using the upper tail contents
883 //          2 : Two-sided test using the upper tail contents
884 //         -1 : One-sided test using the lower tail contents
885 //         -2 : Two-sided test using the lower tail contents
886 //          0 : One-sided test using the auto-selected upper or lower tail contents
887 //          3 : Two-sided test using the auto-selected upper or lower tail contents
888 //
889 // The argument "sigma" allows for the following return values :
890 //
891 // sigma = 0 : P-value is returned as the above specified fraction
892 //         1 : The difference t-<T> expressed in units of sigma
893 //             Note : This difference may be negative and sigma
894 //                    is only defined for ndf>2.
895 //  
896 // The default values are sides=0 and sigma=0.
897 //
898 //--- NvE 21-may-2005 Utrecht University
899
900  if (ndf<=0) return 0;
901
902  Double_t mean=0;
903
904  if (!sides) // Automatic one-sided test
905  {
906   sides=1;
907   if (t<mean) sides=-1;
908  }
909
910  if (sides==3) // Automatic two-sided test
911  {
912   sides=2;
913   if (t<mean) sides=-2;
914  }
915
916  Double_t val=0;
917  if (sigma) // P-value in units of sigma
918  { 
919   if (ndf>2) // Sigma is only defined for ndf>2
920   {
921    Double_t s=sqrt(ndf/(ndf-2.));
922    val=t/s;
923   }
924  }
925  else // P-value from tail contents
926  {
927   if (sides>0) // Upper tail
928   {
929    val=1.-TMath::StudentI(t,ndf);
930   }
931   else // Lower tail
932   {
933    val=TMath::StudentI(t,ndf);
934   }
935  }
936
937  if (abs(sides)==2) val=val*2.;
938
939  return val;
940 }
941 ///////////////////////////////////////////////////////////////////////////
942 Double_t AliMath::FratioPvalue(Double_t f,Int_t ndf1,Int_t ndf2,Int_t sides,Int_t sigma) const
943 {
944 // Computation of the P-value for a certain specified F ratio f value 
945 // for an F (ratio) distribution with ndf1 and ndf2 degrees of freedom
946 // for the two samples X,Y respectively to be compared in the ratio X/Y.
947 //
948 // In a frequentist approach, the F (ratio) distribution is particularly useful
949 // in making inferences about the ratio of the variances of two underlying
950 // populations based on a measurement of the variance of a random sample taken
951 // from each one of the two populations.
952 // So the F test provides a means to investigate the degree of equality of
953 // two populations.
954 //
955 // The P-value for a certain f value corresponds to the fraction of
956 // repeatedly drawn equivalent samples from a certain population, which is expected
957 // to yield an F value exceeding (less than) the value f for an upper (lower)
958 // tail test in case a certain hypothesis is true.
959 //
960 // Further details can be found in the excellent textbook of Phil Gregory
961 // "Bayesian Logical Data Analysis for the Physical Sciences".
962 //
963 // Note : <F>=ndf2/(ndf2-2)  Var(F)=2*ndf2*ndf2*(ndf2+ndf1-2)/(ndf1*(ndf2-1)*(ndf2-1)*(ndf2-4))
964 // 
965 // With the "sides" parameter a one-sided or two-sided test can be selected
966 // using either the upper or lower tail contents.
967 // In case of automatic upper/lower selection the decision is made on basis
968 // of the location of the input f value w.r.t. <F> of the distribution. 
969 //
970 // sides =  1 : One-sided test using the upper tail contents
971 //          2 : Two-sided test using the upper tail contents
972 //         -1 : One-sided test using the lower tail contents
973 //         -2 : Two-sided test using the lower tail contents
974 //          0 : One-sided test using the auto-selected upper or lower tail contents
975 //          3 : Two-sided test using the auto-selected upper or lower tail contents
976 //
977 // The argument "sigma" allows for the following return values :
978 //
979 // sigma = 0 : P-value is returned as the above specified fraction
980 //         1 : The difference f-<F> expressed in units of sigma
981 //             Note : This difference may be negative and sigma
982 //                    is only defined for ndf2>4.
983 //  
984 // The default values are sides=0 and sigma=0.
985 //
986 //--- NvE 21-may-2005 Utrecht University
987
988  if (ndf1<=0 || ndf2<=0 || f<=0) return 0;
989
990  Double_t mean=double(ndf2/(ndf2-2));
991
992  if (!sides) // Automatic one-sided test
993  {
994   sides=1;
995   if (f<mean) sides=-1;
996  }
997
998  if (sides==3) // Automatic two-sided test
999  {
1000   sides=2;
1001   if (f<mean) sides=-2;
1002  }
1003
1004  Double_t val=0;
1005  if (sigma) // P-value in units of sigma
1006  { 
1007   // Sigma is only defined for ndf2>4
1008   if (ndf2>4)
1009   {
1010    Double_t s=sqrt(double(ndf2*ndf2*(2*ndf2+2*ndf1-4))/double(ndf1*pow(ndf2-1,2)*(ndf2-4)));
1011    val=(f-mean)/s;
1012   }
1013  }
1014  else // P-value from tail contents 
1015  {
1016   if (sides>0) // Upper tail
1017   {
1018    val=1.-TMath::FDistI(f,ndf1,ndf2);
1019   }
1020   else // Lower tail
1021   {
1022    val=TMath::FDistI(f,ndf1,ndf2);
1023   }
1024  }
1025
1026  if (abs(sides)==2) val=val*2.;
1027
1028  return val;
1029 }
1030 ///////////////////////////////////////////////////////////////////////////
1031 Double_t AliMath::BinomialPvalue(Int_t k,Int_t n,Double_t p,Int_t sides,Int_t sigma,Int_t mode) const
1032 {
1033 // Computation of the P-value for a certain specified number of successes k
1034 // for a Binomial distribution with n trials and success probability p.
1035 //
1036 // The P-value for a certain number of successes k corresponds to the fraction of
1037 // repeatedly drawn equivalent samples from a certain population, which is expected
1038 // to yield a number of successes exceeding (less than) the value k for an
1039 // upper (lower) tail test in case a certain hypothesis is true.
1040 //
1041 // Further details can be found in the excellent textbook of Phil Gregory
1042 // "Bayesian Logical Data Analysis for the Physical Sciences".
1043 //
1044 // Note : <K>=n*p  Var(K)=n*p*(1-p)
1045 // 
1046 // With the "sides" parameter a one-sided or two-sided test can be selected
1047 // using either the upper or lower tail contents.
1048 // In case of automatic upper/lower selection the decision is made on basis
1049 // of the location of the input k value w.r.t. <K> of the distribution. 
1050 //
1051 // sides =  1 : One-sided test using the upper tail contents
1052 //          2 : Two-sided test using the upper tail contents
1053 //         -1 : One-sided test using the lower tail contents
1054 //         -2 : Two-sided test using the lower tail contents
1055 //          0 : One-sided test using the auto-selected upper or lower tail contents
1056 //          3 : Two-sided test using the auto-selected upper or lower tail contents
1057 //
1058 // The argument "sigma" allows for the following return values :
1059 //
1060 // sigma = 0 : P-value is returned as the above specified fraction
1061 //         1 : The difference k-<K> expressed in units of sigma
1062 //             Note : This difference may be negative.
1063 //
1064 // mode = 0 : Incomplete Beta function will be used to calculate the tail content.
1065 //        1 : Straightforward summation of the Binomial terms is used.
1066 //
1067 // The Incomplete Beta function in general provides the most accurate values.
1068 //
1069 // The default values are sides=0, sigma=0 and mode=0.
1070 //
1071 //--- NvE 24-may-2005 Utrecht University
1072
1073  Double_t mean=double(n)*p;
1074
1075  if (!sides) // Automatic one-sided test
1076  {
1077   sides=1;
1078   if (k<mean) sides=-1;
1079  }
1080
1081  if (sides==3) // Automatic two-sided test
1082  {
1083   sides=2;
1084   if (k<mean) sides=-2;
1085  }
1086
1087  Double_t val=0;
1088
1089  if (sigma) // P-value in units of sigma
1090  {
1091   Double_t s=sqrt(double(n)*p*(1.-p));
1092   val=(double(k)-mean)/s;
1093  }
1094  else // P-value from tail contents
1095  {
1096   if (sides>0)
1097   {
1098    if (!mode) // Use Incomplete Beta function
1099    {
1100     val=TMath::BetaIncomplete(p,k+1,n-k);
1101    }
1102    else // Use straightforward summation
1103    {
1104     for (Int_t i=k+1; i<=n; i++)
1105     {
1106      val+=TMath::Binomial(n,i)*pow(p,i)*pow(1.-p,n-i);
1107     }
1108    }
1109   }
1110   else
1111   {
1112    if (!mode) // Use Incomplete Beta function
1113    {
1114     val=1.-TMath::BetaIncomplete(p,k+1,n-k);
1115    }
1116    else
1117    {
1118     for (Int_t j=0; j<=k; j++)
1119     {
1120      val+=TMath::Binomial(n,j)*pow(p,j)*pow(1.-p,n-j);
1121     }
1122    }
1123   }
1124  }
1125
1126  if (abs(sides)==2) val=val*2.;
1127
1128  return val;
1129 }
1130 ///////////////////////////////////////////////////////////////////////////
1131 Double_t AliMath::PoissonPvalue(Int_t k,Double_t mu,Int_t sides,Int_t sigma) const
1132 {
1133 // Computation of the P-value for a certain specified number of occurrences k
1134 // for a Poisson distribution with a specified average rate (in time or space)
1135 // mu of occurrences.
1136 //
1137 // The P-value for a certain number of occurrences k corresponds to the fraction of
1138 // repeatedly drawn equivalent samples from a certain population, which is expected
1139 // to yield a number of occurrences exceeding (less than) the value k for an
1140 // upper (lower) tail test in case a certain hypothesis is true.
1141 //
1142 // Further details can be found in the excellent textbook of Phil Gregory
1143 // "Bayesian Logical Data Analysis for the Physical Sciences".
1144 //
1145 // Note : <K>=mu  Var(K)=mu
1146 // 
1147 // With the "sides" parameter a one-sided or two-sided test can be selected
1148 // using either the upper or lower tail contents.
1149 // In case of automatic upper/lower selection the decision is made on basis
1150 // of the location of the input k value w.r.t. <K> of the distribution. 
1151 //
1152 // sides =  1 : One-sided test using the upper tail contents
1153 //          2 : Two-sided test using the upper tail contents
1154 //         -1 : One-sided test using the lower tail contents
1155 //         -2 : Two-sided test using the lower tail contents
1156 //          0 : One-sided test using the auto-selected upper or lower tail contents
1157 //          3 : Two-sided test using the auto-selected upper or lower tail contents
1158 //
1159 // The argument "sigma" allows for the following return values :
1160 //
1161 // sigma = 0 : P-value is returned as the above specified fraction
1162 //         1 : The difference k-<K> expressed in units of sigma
1163 //             Note : This difference may be negative.
1164 //
1165 // The default values are sides=0 and sigma=0.
1166 //
1167 // Note : The tail contents are given by the incomplete Gamma function P(a,x).
1168 //        Lower tail contents = 1-P(k,mu)
1169 //        Upper tail contents = P(k,mu)
1170 //
1171 //--- NvE 24-may-2005 Utrecht University
1172
1173  Double_t mean=mu;
1174
1175  if (!sides) // Automatic one-sided test
1176  {
1177   sides=1;
1178   if (k<mean) sides=-1;
1179  }
1180
1181  if (sides==3) // Automatic two-sided test
1182  {
1183   sides=2;
1184   if (k<mean) sides=-2;
1185  }
1186
1187  Double_t val=0;
1188
1189  if (sigma) // P-value in units of sigma
1190  {
1191   Double_t s=sqrt(mu);
1192   val=(double(k)-mean)/s;
1193  }
1194  else // P-value from tail contents
1195  {
1196   if (sides>0) // Upper tail
1197   {
1198    val=Gamma(k,mu);
1199   }
1200   else // Lower tail
1201   {
1202    val=1.-Gamma(k,mu);
1203   }
1204  }
1205
1206  if (abs(sides)==2) val=val*2.;
1207
1208  return val;
1209 }
1210 ///////////////////////////////////////////////////////////////////////////
1211 Double_t AliMath::NegBinomialPvalue(Int_t n,Int_t k,Double_t p,Int_t sides,Int_t sigma,Int_t mode) const
1212 {
1213 // Computation of the P-value for a certain specified number of trials n
1214 // for a Negative Binomial distribution where exactly k successes are to
1215 // be reached which have each a probability p.
1216 //
1217 // The P-value for a certain number of trials n corresponds to the fraction of
1218 // repeatedly drawn equivalent samples from a certain population, which is expected
1219 // to yield a number of trials exceeding (less than) the value n for an
1220 // upper (lower) tail test in case a certain hypothesis is true.
1221 //
1222 // Further details can be found in the excellent textbook of Phil Gregory
1223 // "Bayesian Logical Data Analysis for the Physical Sciences".
1224 //
1225 // Note : <N>=k*(1-p)/p  Var(N)=k*(1-p)/(p*p) 
1226 // 
1227 // With the "sides" parameter a one-sided or two-sided test can be selected
1228 // using either the upper or lower tail contents.
1229 // In case of automatic upper/lower selection the decision is made on basis
1230 // of the location of the input n value w.r.t. <N> of the distribution. 
1231 //
1232 // sides =  1 : One-sided test using the upper tail contents
1233 //          2 : Two-sided test using the upper tail contents
1234 //         -1 : One-sided test using the lower tail contents
1235 //         -2 : Two-sided test using the lower tail contents
1236 //          0 : One-sided test using the auto-selected upper or lower tail contents
1237 //          3 : Two-sided test using the auto-selected upper or lower tail contents
1238 //
1239 // The argument "sigma" allows for the following return values :
1240 //
1241 // sigma = 0 : P-value is returned as the above specified fraction
1242 //         1 : The difference n-<N> expressed in units of sigma
1243 //             Note : This difference may be negative.
1244 //
1245 // mode = 0 : Incomplete Beta function will be used to calculate the tail content.
1246 //        1 : Straightforward summation of the Negative Binomial terms is used.
1247 //
1248 // The Incomplete Beta function in general provides the most accurate values.
1249 //
1250 // The default values are sides=0, sigma=0 and mode=0.
1251 //
1252 //--- NvE 24-may-2005 Utrecht University
1253
1254  Double_t mean=double(k)*(1.-p)/p;
1255
1256  if (!sides) // Automatic one-sided test
1257  {
1258   sides=1;
1259   if (n<mean) sides=-1;
1260  }
1261
1262  if (sides==3) // Automatic two-sided test
1263  {
1264   sides=2;
1265   if (n<mean) sides=-2;
1266  }
1267
1268  Double_t val=0;
1269
1270  if (sigma) // P-value in units of sigma
1271  {
1272   Double_t s=sqrt(double(k)*(1.-p)/(p*p));
1273   val=(double(n)-mean)/s;
1274  }
1275  else // P-value from tail contents
1276  {
1277   if (sides>0) // Upper tail
1278   {
1279    if (!mode) // Use Incomplete Beta function
1280    {
1281     val=1.-TMath::BetaIncomplete(p,k,n-k);
1282    }
1283    else // Use straightforward summation
1284    {
1285     for (Int_t i=1; i<n; i++)
1286     {
1287      val+=TMath::Binomial(i-1,k-1)*pow(p,k)*pow(1.-p,i-k);
1288     }
1289     val=1.-val;
1290    }
1291   }
1292   else // Lower tail
1293   {
1294    if (!mode) // Use Incomplete Beta function
1295    {
1296     val=TMath::BetaIncomplete(p,k,n-k);
1297    }
1298    else
1299    {
1300     for (Int_t j=1; j<n; j++)
1301     {
1302      val+=TMath::Binomial(j-1,k-1)*pow(p,k)*pow(1.-p,j-k);
1303     }
1304    }
1305   }
1306  }
1307
1308  if (abs(sides)==2) val=val*2.;
1309
1310  return val;
1311 }
1312 ///////////////////////////////////////////////////////////////////////////