]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RALICE/AliMath.cxx
Configs for DC05 added (Test version !).
[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 // According to the value of the parameter "mode" various algorithms
292 // can be selected.
293 //
294 // mode = 0 : Calculations are based on the incomplete gamma function P(a,x),
295 //            where a=ndf/2 and x=chi2/2.
296 //
297 //        1 : Same as for mode=0. However, in case ndf=1 an exact expression
298 //            based on the error function Erf() is used.
299 //
300 //        2 : Same as for mode=0. However, in case ndf>30 a Gaussian approximation
301 //            is used instead of the gamma function.
302 //
303 // When invoked as Prob(chi2,ndf) the default mode=1 is used.
304 //
305 // P(a,x) represents the probability that the observed Chi-squared
306 // for a correct model should be less than the value chi2.
307 //
308 // The returned probability corresponds to 1-P(a,x),
309 // which denotes the probability that an observed Chi-squared exceeds
310 // the value chi2 by chance, even for a correct model.
311 //
312 //--- NvE 14-nov-1998 UU-SAP Utrecht
313  
314  if (ndf <= 0) return 0; // Set CL to zero in case ndf<=0
315  
316  if (chi2 <= 0.)
317  {
318   if (chi2 < 0.)
319   {
320     return 0;
321   }
322   else
323   {
324    return 1;
325   }
326  }
327
328  Double_t v=-1.;
329
330  switch (mode)
331  {
332   case 1: // Exact expression for ndf=1 as alternative for the gamma function
333    if (ndf==1) v=1.-Erf(sqrt(chi2)/sqrt(2.));
334    break;
335  
336   case 2: // Gaussian approximation for large ndf (i.e. ndf>30) as alternative for the gamma function
337    if (ndf>30)
338    {
339     Double_t q=sqrt(2.*chi2)-sqrt(double(2*ndf-1));
340     if (q>0.) v=0.5*(1.-Erf(q/sqrt(2.)));
341    }
342    break;
343  }
344  
345  if (v<0.)
346  {
347   // Evaluate the incomplete gamma function
348   Double_t a=double(ndf)/2.;
349   Double_t x=chi2/2.;
350   v=1.-Gamma(a,x);
351  }
352
353  return v;
354 }
355 ///////////////////////////////////////////////////////////////////////////
356 Double_t AliMath::BesselI0(Double_t x) const
357 {
358 // Computation of the modified Bessel function I_0(x) for any real x.
359 //
360 // The algorithm is based on the article by Abramowitz and Stegun [1]
361 // as denoted in Numerical Recipes 2nd ed. on p. 230 (W.H.Press et al.).
362 //
363 // [1] M.Abramowitz and I.A.Stegun, Handbook of Mathematical Functions,
364 //     Applied Mathematics Series vol. 55 (1964), Washington.  
365 //
366 //--- NvE 12-mar-2000 UU-SAP Utrecht
367
368  // Parameters of the polynomial approximation  
369  const Double_t kp1=1.0,          kp2=3.5156229,    kp3=3.0899424,
370                 kp4=1.2067492,    kp5=0.2659732,    kp6=3.60768e-2,  kp7=4.5813e-3;
371
372  const Double_t kq1= 0.39894228,  kq2= 1.328592e-2, kq3= 2.25319e-3,
373                 kq4=-1.57565e-3,  kq5= 9.16281e-3,  kq6=-2.057706e-2,
374                 kq7= 2.635537e-2, kq8=-1.647633e-2, kq9= 3.92377e-3; 
375
376  Double_t ax=fabs(x);
377
378  Double_t y=0,result=0;
379
380  if (ax < 3.75)
381  {
382   y=pow(x/3.75,2);
383   result=kp1+y*(kp2+y*(kp3+y*(kp4+y*(kp5+y*(kp6+y*kp7)))));
384  }
385  else
386  {
387   y=3.75/ax;
388   result=(exp(ax)/sqrt(ax))
389          *(kq1+y*(kq2+y*(kq3+y*(kq4+y*(kq5+y*(kq6+y*(kq7+y*(kq8+y*kq9))))))));
390  }
391
392  return result;
393 }
394 ///////////////////////////////////////////////////////////////////////////
395 Double_t AliMath::BesselK0(Double_t x) const
396 {
397 // Computation of the modified Bessel function K_0(x) for positive real x.
398 //
399 // The algorithm is based on the article by Abramowitz and Stegun [1]
400 // as denoted in Numerical Recipes 2nd ed. on p. 230 (W.H.Press et al.).
401 //
402 // [1] M.Abramowitz and I.A.Stegun, Handbook of Mathematical Functions,
403 //     Applied Mathematics Series vol. 55 (1964), Washington.  
404 //
405 //--- NvE 12-mar-2000 UU-SAP Utrecht
406
407  // Parameters of the polynomial approximation  
408  const Double_t kp1=-0.57721566,  kp2=0.42278420,   kp3=0.23069756,
409                 kp4= 3.488590e-2, kp5=2.62698e-3,   kp6=1.0750e-4,    kp7=7.4e-6;
410
411  const Double_t kq1= 1.25331414,  kq2=-7.832358e-2, kq3= 2.189568e-2,
412                 kq4=-1.062446e-2, kq5= 5.87872e-3,  kq6=-2.51540e-3,  kq7=5.3208e-4;
413
414  if (x <= 0)
415  {
416   cout << " *BesselK0* Invalid argument x = " << x << endl;
417   return 0;
418  }
419
420  Double_t y=0,result=0;
421
422  if (x <= 2)
423  {
424   y=x*x/4.;
425   result=(-log(x/2.)*BesselI0(x))
426          +(kp1+y*(kp2+y*(kp3+y*(kp4+y*(kp5+y*(kp6+y*kp7))))));
427  }
428  else
429  {
430   y=2./x;
431   result=(exp(-x)/sqrt(x))
432          *(kq1+y*(kq2+y*(kq3+y*(kq4+y*(kq5+y*(kq6+y*kq7))))));
433  }
434
435  return result;
436 }
437 ///////////////////////////////////////////////////////////////////////////
438 Double_t AliMath::BesselI1(Double_t x) const
439 {
440 // Computation of the modified Bessel function I_1(x) for any real x.
441 //
442 // The algorithm is based on the article by Abramowitz and Stegun [1]
443 // as denoted in Numerical Recipes 2nd ed. on p. 230 (W.H.Press et al.).
444 //
445 // [1] M.Abramowitz and I.A.Stegun, Handbook of Mathematical Functions,
446 //     Applied Mathematics Series vol. 55 (1964), Washington.  
447 //
448 //--- NvE 12-mar-2000 UU-SAP Utrecht
449
450  // Parameters of the polynomial approximation  
451  const Double_t kp1=0.5,          kp2=0.87890594,   kp3=0.51498869,
452                 kp4=0.15084934,   kp5=2.658733e-2,  kp6=3.01532e-3,  kp7=3.2411e-4;
453
454  const Double_t kq1= 0.39894228,  kq2=-3.988024e-2, kq3=-3.62018e-3,
455                 kq4= 1.63801e-3,  kq5=-1.031555e-2, kq6= 2.282967e-2,
456                 kq7=-2.895312e-2, kq8= 1.787654e-2, kq9=-4.20059e-3; 
457
458  Double_t ax=fabs(x);
459
460  Double_t y=0,result=0;
461
462  if (ax < 3.75)
463  {
464   y=pow(x/3.75,2);
465   result=x*(kp1+y*(kp2+y*(kp3+y*(kp4+y*(kp5+y*(kp6+y*kp7))))));
466  }
467  else
468  {
469   y=3.75/ax;
470   result=(exp(ax)/sqrt(ax))
471          *(kq1+y*(kq2+y*(kq3+y*(kq4+y*(kq5+y*(kq6+y*(kq7+y*(kq8+y*kq9))))))));
472   if (x < 0) result=-result;
473  }
474
475  return result;
476 }
477 ///////////////////////////////////////////////////////////////////////////
478 Double_t AliMath::BesselK1(Double_t x) const
479 {
480 // Computation of the modified Bessel function K_1(x) for positive real x.
481 //
482 // The algorithm is based on the article by Abramowitz and Stegun [1]
483 // as denoted in Numerical Recipes 2nd ed. on p. 230 (W.H.Press et al.).
484 //
485 // [1] M.Abramowitz and I.A.Stegun, Handbook of Mathematical Functions,
486 //     Applied Mathematics Series vol. 55 (1964), Washington.  
487 //
488 //--- NvE 12-mar-2000 UU-SAP Utrecht
489
490  // Parameters of the polynomial approximation  
491  const Double_t kp1= 1.,          kp2= 0.15443144,  kp3=-0.67278579,
492                 kp4=-0.18156897,  kp5=-1.919402e-2, kp6=-1.10404e-3,  kp7=-4.686e-5;
493
494  const Double_t kq1= 1.25331414,  kq2= 0.23498619,  kq3=-3.655620e-2,
495                 kq4= 1.504268e-2, kq5=-7.80353e-3,  kq6= 3.25614e-3,  kq7=-6.8245e-4;
496
497  if (x <= 0)
498  {
499   cout << " *BesselK1* Invalid argument x = " << x << endl;
500   return 0;
501  }
502
503  Double_t y=0,result=0;
504
505  if (x <= 2)
506  {
507   y=x*x/4.;
508   result=(log(x/2.)*BesselI1(x))
509          +(1./x)*(kp1+y*(kp2+y*(kp3+y*(kp4+y*(kp5+y*(kp6+y*kp7))))));
510  }
511  else
512  {
513   y=2./x;
514   result=(exp(-x)/sqrt(x))
515          *(kq1+y*(kq2+y*(kq3+y*(kq4+y*(kq5+y*(kq6+y*kq7))))));
516  }
517
518  return result;
519 }
520 ///////////////////////////////////////////////////////////////////////////
521 Double_t AliMath::BesselK(Int_t n,Double_t x) const
522 {
523 // Computation of the Integer Order Modified Bessel function K_n(x)
524 // for n=0,1,2,... and positive real x.
525 //
526 // The algorithm uses the recurrence relation
527 //
528 //               K_n+1(x) = (2n/x)*K_n(x) + K_n-1(x) 
529 //
530 // as denoted in Numerical Recipes 2nd ed. on p. 232 (W.H.Press et al.).
531 //
532 //--- NvE 12-mar-2000 UU-SAP Utrecht
533
534  if (x <= 0 || n < 0)
535  {
536   cout << " *BesselK* Invalid argument(s) (n,x) = (" << n << " , " << x << ")" << endl;
537   return 0;
538  }
539
540  if (n==0) return BesselK0(x);
541
542  if (n==1) return BesselK1(x);
543
544  // Perform upward recurrence for all x
545  Double_t tox=2./x;
546  Double_t bkm=BesselK0(x);
547  Double_t bk=BesselK1(x);
548  Double_t bkp=0;
549  for (Int_t j=1; j<n; j++)
550  {
551   bkp=bkm+double(j)*tox*bk;
552   bkm=bk;
553   bk=bkp;
554  }
555
556  return bk;
557 }
558 ///////////////////////////////////////////////////////////////////////////
559 Double_t AliMath::BesselI(Int_t n,Double_t x) const
560 {
561 // Computation of the Integer Order Modified Bessel function I_n(x)
562 // for n=0,1,2,... and any real x.
563 //
564 // The algorithm uses the recurrence relation
565 //
566 //               I_n+1(x) = (-2n/x)*I_n(x) + I_n-1(x) 
567 //
568 // as denoted in Numerical Recipes 2nd ed. on p. 232 (W.H.Press et al.).
569 //
570 //--- NvE 12-mar-2000 UU-SAP Utrecht
571
572  Int_t iacc=40; // Increase to enhance accuracy
573  Double_t bigno=1.e10, bigni=1.e-10;
574
575  if (n < 0)
576  {
577   cout << " *BesselI* Invalid argument (n,x) = (" << n << " , " << x << ")" << endl;
578   return 0;
579  }
580
581  if (n==0) return BesselI0(x);
582
583  if (n==1) return BesselI1(x);
584
585  if (fabs(x) < 1.e-10) return 0;
586
587  Double_t tox=2./fabs(x);
588  Double_t bip=0,bim=0;
589  Double_t bi=1;
590  Double_t result=0;
591  Int_t m=2*((n+int(sqrt(float(iacc*n))))); // Downward recurrence from even m
592  for (Int_t j=m; j<=1; j--)
593  {
594   bim=bip+double(j)*tox*bi;
595   bip=bi;
596   bi=bim;
597   if (fabs(bi) > bigno) // Renormalise to prevent overflows
598   {
599    result*=bigni;
600    bi*=bigni;
601    bip*=bigni;
602   }
603   if (j==n) result=bip;
604  }
605
606  result*=BesselI0(x)/bi; // Normalise with I0(x)
607  if ((x < 0) && (n%2 == 1)) result=-result;
608
609  return result;
610 }
611 ///////////////////////////////////////////////////////////////////////////