]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RALICE/AliMath.cxx
29-nov-2002 NvE Typo corrected in AliCollider.cxx.
[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  
39 ClassImp(AliMath) // Class implementation to enable ROOT I/O
40  
41 AliMath::AliMath()
42 {
43 // Default constructor
44 }
45 ///////////////////////////////////////////////////////////////////////////
46 AliMath::~AliMath()
47 {
48 // Destructor
49 }
50 ///////////////////////////////////////////////////////////////////////////
51 Double_t AliMath::Gamma(Double_t z)
52 {
53 // Computation of gamma(z) for all z>0.
54 //
55 // The algorithm is based on the article by C.Lanczos [1] as denoted in
56 // Numerical Recipes 2nd ed. on p. 207 (W.H.Press et al.).
57 //
58 // [1] C.Lanczos, SIAM Journal of Numerical Analysis B1 (1964), 86.
59 //
60 //--- Nve 14-nov-1998 UU-SAP Utrecht
61  
62  if (z<=0.)
63  {
64   cout << "*Gamma(z)* Wrong argument z = " << z << endl;
65   return 0;
66  }
67  
68  Double_t v=LnGamma(z);
69  return exp(v);
70 }
71 ///////////////////////////////////////////////////////////////////////////
72 Double_t AliMath::Gamma(Double_t a,Double_t x)
73 {
74 // Computation of the incomplete gamma function P(a,x)
75 //
76 // The algorithm is based on the formulas and code as denoted in
77 // Numerical Recipes 2nd ed. on p. 210-212 (W.H.Press et al.).
78 //
79 //--- Nve 14-nov-1998 UU-SAP Utrecht
80  
81  if (a<=0.)
82  {
83   cout << "*Gamma(a,x)* Invalid argument a = " << a << endl;
84   return 0;
85  }
86  
87  if (x<=0.)
88  {
89   if (x<0) cout << "*Gamma(a,x)* Invalid argument x = " << x << endl;
90   return 0;
91  }
92  
93  if (x<(a+1.))
94  {
95   return GamSer(a,x);
96  }
97  else
98  {
99   return GamCf(a,x);
100  }
101 }
102 ///////////////////////////////////////////////////////////////////////////
103 Double_t AliMath::LnGamma(Double_t z)
104 {
105 // Computation of ln[gamma(z)] for all z>0.
106 //
107 // The algorithm is based on the article by C.Lanczos [1] as denoted in
108 // Numerical Recipes 2nd ed. on p. 207 (W.H.Press et al.).
109 //
110 // [1] C.Lanczos, SIAM Journal of Numerical Analysis B1 (1964), 86.
111 //
112 // The accuracy of the result is better than 2e-10.
113 //
114 //--- Nve 14-nov-1998 UU-SAP Utrecht
115  
116  if (z<=0.)
117  {
118   cout << "*LnGamma(z)* Wrong argument z = " << z << endl;
119   return 0;
120  }
121  
122  // Coefficients for the series expansion
123  Double_t c[7];
124  c[0]=  2.5066282746310005;
125  c[1]= 76.18009172947146;
126  c[2]=-86.50532032941677;
127  c[3]= 24.01409824083091;
128  c[4]= -1.231739572450155;
129  c[5]=  0.1208650973866179e-2;
130  c[6]= -0.5395239384953e-5;
131  
132  Double_t x=z;
133  Double_t y=x;
134  Double_t tmp=x+5.5;
135  tmp=(x+0.5)*log(tmp)-tmp;
136  Double_t ser=1.000000000190015;
137  for (Int_t i=1; i<7; i++)
138  {
139   y+=1.;
140   ser+=c[i]/y;
141  }
142  Double_t v=tmp+log(c[0]*ser/x);
143  return v;
144 }
145 ///////////////////////////////////////////////////////////////////////////
146 Double_t AliMath::GamSer(Double_t a,Double_t x)
147 {
148 // Computation of the incomplete gamma function P(a,x)
149 // via its series representation.
150 //
151 // The algorithm is based on the formulas and code as denoted in
152 // Numerical Recipes 2nd ed. on p. 210-212 (W.H.Press et al.).
153 //
154 //--- Nve 14-nov-1998 UU-SAP Utrecht
155  
156  Int_t itmax=100;   // Maximum number of iterations
157  Double_t eps=3.e-7; // Relative accuracy
158  
159  if (a<=0.)
160  {
161   cout << "*GamSer(a,x)* Invalid argument a = " << a << endl;
162   return 0;
163  }
164  
165  if (x<=0.)
166  {
167   if (x<0) cout << "*GamSer(a,x)* Invalid argument x = " << x << endl;
168   return 0;
169  }
170  
171  Double_t gln=LnGamma(a);
172  Double_t ap=a;
173  Double_t sum=1./a;
174  Double_t del=sum;
175  for (Int_t n=1; n<=itmax; n++)
176  {
177   ap+=1.;
178   del=del*x/ap;
179   sum+=del;
180   if (fabs(del)<fabs(sum*eps)) break;
181   if (n==itmax) cout << "*GamSer(a,x)* a too large or itmax too small" << endl;
182  }
183  Double_t v=sum*exp(-x+a*log(x)-gln);
184  return v;
185 }
186 ///////////////////////////////////////////////////////////////////////////
187 Double_t AliMath::GamCf(Double_t a,Double_t x)
188 {
189 // Computation of the incomplete gamma function P(a,x)
190 // via its continued fraction representation.
191 //
192 // The algorithm is based on the formulas and code as denoted in
193 // Numerical Recipes 2nd ed. on p. 210-212 (W.H.Press et al.).
194 //
195 //--- Nve 14-nov-1998 UU-SAP Utrecht
196  
197  Int_t itmax=100;      // Maximum number of iterations
198  Double_t eps=3.e-7;    // Relative accuracy
199  Double_t fpmin=1.e-30; // Smallest Double_t value allowed here
200  
201  if (a<=0.)
202  {
203   cout << "*GamCf(a,x)* Invalid argument a = " << a << endl;
204   return 0;
205  }
206  
207  if (x<=0.)
208  {
209   if (x<0) cout << "*GamCf(a,x)* Invalid argument x = " << x << endl;
210   return 0;
211  }
212  
213  Double_t gln=LnGamma(a);
214  Double_t b=x+1.-a;
215  Double_t c=1./fpmin;
216  Double_t d=1./b;
217  Double_t h=d;
218  Double_t an,del;
219  for (Int_t i=1; i<=itmax; i++)
220  {
221   an=double(-i)*(double(i)-a);
222   b+=2.;
223   d=an*d+b;
224   if (fabs(d)<fpmin) d=fpmin;
225   c=b+an/c;
226   if (fabs(c)<fpmin) c=fpmin;
227   d=1./d;
228   del=d*c;
229   h=h*del;
230   if (fabs(del-1.)<eps) break;
231   if (i==itmax) cout << "*GamCf(a,x)* a too large or itmax too small" << endl;
232  }
233  Double_t v=exp(-x+a*log(x)-gln)*h;
234  return (1.-v);
235 }
236 ///////////////////////////////////////////////////////////////////////////
237 Double_t AliMath::Erf(Double_t x)
238 {
239 // Computation of the error function erf(x).
240 //
241 //--- NvE 14-nov-1998 UU-SAP Utrecht
242  
243  return (1.-Erfc(x));
244 }
245 ///////////////////////////////////////////////////////////////////////////
246 Double_t AliMath::Erfc(Double_t x)
247 {
248 // Computation of the complementary error function erfc(x).
249 //
250 // The algorithm is based on a Chebyshev fit as denoted in
251 // Numerical Recipes 2nd ed. on p. 214 (W.H.Press et al.).
252 //
253 // The fractional error is always less than 1.2e-7.
254 //
255 //--- Nve 14-nov-1998 UU-SAP Utrecht
256  
257  // The parameters of the Chebyshev fit
258  const Double_t a1=-1.26551223,  a2=1.00002368,
259                 a3= 0.37409196,  a4=0.09678418,
260                 a5=-0.18628806,  a6=0.27886807,
261                 a7=-1.13520398,  a8=1.48851587,
262                 a9=-0.82215223, a10=0.17087277;
263  
264  Double_t v=1.; // The return value
265  
266  Double_t z=fabs(x);
267  
268  if (z <= 0.) return v; // erfc(0)=1
269  
270  Double_t t=1./(1.+0.5*z);
271  
272  v=t*exp((-z*z)
273          +a1+t*(a2+t*(a3+t*(a4+t*(a5+t*(a6+t*(a7+t*(a8+t*(a9+t*a10)))))))));
274  
275  if (x < 0.) v=2.-v; // erfc(-x)=2-erfc(x)
276  
277  return v;
278 }
279 ///////////////////////////////////////////////////////////////////////////
280 Double_t AliMath::Prob(Double_t chi2,Int_t ndf)
281 {
282 // Computation of the probability for a certain Chi-squared (chi2)
283 // and number of degrees of freedom (ndf).
284 //
285 // Calculations are based on the incomplete gamma function P(a,x),
286 // where a=ndf/2 and x=chi2/2.
287 //
288 // P(a,x) represents the probability that the observed Chi-squared
289 // for a correct model should be less than the value chi2.
290 //
291 // The returned probability corresponds to 1-P(a,x),
292 // which denotes the probability that an observed Chi-squared exceeds
293 // the value chi2 by chance, even for a correct model.
294 //
295 //--- NvE 14-nov-1998 UU-SAP Utrecht
296  
297  if (ndf <= 0) return 0; // Set CL to zero in case ndf<=0
298  
299  if (chi2 <= 0.)
300  {
301   if (chi2 < 0.)
302   {
303     return 0;
304   }
305   else
306   {
307    return 1;
308   }
309  }
310  
311 // Alternative which is exact
312 // This code may be activated in case the gamma function gives problems
313 // if (ndf==1)
314 // {
315 //  Double_t v=1.-Erf(sqrt(chi2)/sqrt(2.));
316 //  return v;
317 // }
318  
319 // Gaussian approximation for large ndf
320 // This code may be activated in case the gamma function shows a problem
321 // Double_t q=sqrt(2.*chi2)-sqrt(double(2*ndf-1));
322 // if (n>30 && q>0.)
323 // {
324 //  Double_t v=0.5*(1.-Erf(q/sqrt(2.)));
325 //  return v;
326 // }
327  
328  // Evaluate the incomplete gamma function
329  Double_t a=double(ndf)/2.;
330  Double_t x=chi2/2.;
331  return (1.-Gamma(a,x));
332 }
333 ///////////////////////////////////////////////////////////////////////////
334 Double_t AliMath::BesselI0(Double_t x)
335 {
336 // Computation of the modified Bessel function I_0(x) for any real x.
337 //
338 // The algorithm is based on the article by Abramowitz and Stegun [1]
339 // as denoted in Numerical Recipes 2nd ed. on p. 230 (W.H.Press et al.).
340 //
341 // [1] M.Abramowitz and I.A.Stegun, Handbook of Mathematical Functions,
342 //     Applied Mathematics Series vol. 55 (1964), Washington.  
343 //
344 //--- NvE 12-mar-2000 UU-SAP Utrecht
345
346  // Parameters of the polynomial approximation  
347  const Double_t p1=1.0,          p2=3.5156229,    p3=3.0899424,
348                 p4=1.2067492,    p5=0.2659732,    p6=3.60768e-2,  p7=4.5813e-3;
349
350  const Double_t q1= 0.39894228,  q2= 1.328592e-2, q3= 2.25319e-3,
351                 q4=-1.57565e-3,  q5= 9.16281e-3,  q6=-2.057706e-2,
352                 q7= 2.635537e-2, q8=-1.647633e-2, q9= 3.92377e-3; 
353
354  Double_t ax=fabs(x);
355
356  Double_t y=0,result=0;
357
358  if (ax < 3.75)
359  {
360   y=pow(x/3.75,2);
361   result=p1+y*(p2+y*(p3+y*(p4+y*(p5+y*(p6+y*p7)))));
362  }
363  else
364  {
365   y=3.75/ax;
366   result=(exp(ax)/sqrt(ax))*(q1+y*(q2+y*(q3+y*(q4+y*(q5+y*(q6+y*(q7+y*(q8+y*q9))))))));
367  }
368
369  return result;
370 }
371 ///////////////////////////////////////////////////////////////////////////
372 Double_t AliMath::BesselK0(Double_t x)
373 {
374 // Computation of the modified Bessel function K_0(x) for positive real x.
375 //
376 // The algorithm is based on the article by Abramowitz and Stegun [1]
377 // as denoted in Numerical Recipes 2nd ed. on p. 230 (W.H.Press et al.).
378 //
379 // [1] M.Abramowitz and I.A.Stegun, Handbook of Mathematical Functions,
380 //     Applied Mathematics Series vol. 55 (1964), Washington.  
381 //
382 //--- NvE 12-mar-2000 UU-SAP Utrecht
383
384  // Parameters of the polynomial approximation  
385  const Double_t p1=-0.57721566,  p2=0.42278420,   p3=0.23069756,
386                 p4= 3.488590e-2, p5=2.62698e-3,   p6=1.0750e-4,    p7=7.4e-5;
387
388  const Double_t q1= 1.25331414,  q2=-7.832358e-2, q3= 2.189568e-2,
389                 q4=-1.062446e-2, q5= 5.87872e-3,  q6=-2.51540e-3,  q7=5.3208e-4;
390
391  if (x <= 0)
392  {
393   cout << " *BesselK0* Invalid argument x = " << x << endl;
394   return 0;
395  }
396
397  Double_t y=0,result=0;
398
399  if (x <= 2)
400  {
401   y=x*x/4.;
402   result=(-log(x/2.)*BesselI0(x))+(p1+y*(p2+y*(p3+y*(p4+y*(p5+y*(p6+y*p7))))));
403  }
404  else
405  {
406   y=2./x;
407   result=(exp(-x)/sqrt(x))*(q1+y*(q2+y*(q3+y*(q4+y*(q5+y*(q6+y*q7))))));
408  }
409
410  return result;
411 }
412 ///////////////////////////////////////////////////////////////////////////
413 Double_t AliMath::BesselI1(Double_t x)
414 {
415 // Computation of the modified Bessel function I_1(x) for any real x.
416 //
417 // The algorithm is based on the article by Abramowitz and Stegun [1]
418 // as denoted in Numerical Recipes 2nd ed. on p. 230 (W.H.Press et al.).
419 //
420 // [1] M.Abramowitz and I.A.Stegun, Handbook of Mathematical Functions,
421 //     Applied Mathematics Series vol. 55 (1964), Washington.  
422 //
423 //--- NvE 12-mar-2000 UU-SAP Utrecht
424
425  // Parameters of the polynomial approximation  
426  const Double_t p1=0.5,          p2=0.87890594,   p3=0.51498869,
427                 p4=0.15084934,   p5=2.658733e-2,  p6=3.01532e-3,  p7=3.2411e-4;
428
429  const Double_t q1= 0.39894228,  q2=-3.988024e-2, q3=-3.62018e-3,
430                 q4= 1.63801e-3,  q5=-1.031555e-2, q6= 2.282967e-2,
431                 q7=-2.895312e-2, q8= 1.787654e-2, q9=-4.20059e-3; 
432
433  Double_t ax=fabs(x);
434
435  Double_t y=0,result=0;
436
437  if (ax < 3.75)
438  {
439   y=pow(x/3.75,2);
440   result=x*(p1+y*(p2+y*(p3+y*(p4+y*(p5+y*(p6+y*p7))))));
441  }
442  else
443  {
444   y=3.75/ax;
445   result=(exp(ax)/sqrt(ax))*(q1+y*(q2+y*(q3+y*(q4+y*(q5+y*(q6+y*(q7+y*(q8+y*q9))))))));
446   if (x < 0) result=-result;
447  }
448
449  return result;
450 }
451 ///////////////////////////////////////////////////////////////////////////
452 Double_t AliMath::BesselK1(Double_t x)
453 {
454 // Computation of the modified Bessel function K_1(x) for positive real x.
455 //
456 // The algorithm is based on the article by Abramowitz and Stegun [1]
457 // as denoted in Numerical Recipes 2nd ed. on p. 230 (W.H.Press et al.).
458 //
459 // [1] M.Abramowitz and I.A.Stegun, Handbook of Mathematical Functions,
460 //     Applied Mathematics Series vol. 55 (1964), Washington.  
461 //
462 //--- NvE 12-mar-2000 UU-SAP Utrecht
463
464  // Parameters of the polynomial approximation  
465  const Double_t p1= 1.,          p2= 0.15443144,  p3=-0.67278579,
466                 p4=-0.18156897,  p5=-1.919402e-2, p6=-1.10404e-3,  p7=-4.686e-5;
467
468  const Double_t q1= 1.25331414,  q2= 0.23498619,  q3=-3.655620e-2,
469                 q4= 1.504268e-2, q5=-7.80353e-3,  q6= 3.25614e-3,  q7=-6.8245e-4;
470
471  if (x <= 0)
472  {
473   cout << " *BesselK1* Invalid argument x = " << x << endl;
474   return 0;
475  }
476
477  Double_t y=0,result=0;
478
479  if (x <= 2)
480  {
481   y=x*x/4.;
482   result=(log(x/2.)*BesselI1(x))+(1./x)*(p1+y*(p2+y*(p3+y*(p4+y*(p5+y*(p6+y*p7))))));
483  }
484  else
485  {
486   y=2./x;
487   result=(exp(-x)/sqrt(x))*(q1+y*(q2+y*(q3+y*(q4+y*(q5+y*(q6+y*q7))))));
488  }
489
490  return result;
491 }
492 ///////////////////////////////////////////////////////////////////////////
493 Double_t AliMath::BesselK(Int_t n,Double_t x)
494 {
495 // Computation of the Integer Order Modified Bessel function K_n(x)
496 // for n=0,1,2,... and positive real x.
497 //
498 // The algorithm uses the recurrence relation
499 //
500 //               K_n+1(x) = (2n/x)*K_n(x) + K_n-1(x) 
501 //
502 // as denoted in Numerical Recipes 2nd ed. on p. 232 (W.H.Press et al.).
503 //
504 //--- NvE 12-mar-2000 UU-SAP Utrecht
505
506  if (x <= 0 || n < 0)
507  {
508   cout << " *BesselK* Invalid argument(s) (n,x) = (" << n << " , " << x << ")" << endl;
509   return 0;
510  }
511
512  if (n==0) return BesselK0(x);
513
514  if (n==1) return BesselK1(x);
515
516  // Perform upward recurrence for all x
517  Double_t tox=2./x;
518  Double_t bkm=BesselK0(x);
519  Double_t bk=BesselK1(x);
520  Double_t bkp=0;
521  for (Int_t j=1; j<n; j++)
522  {
523   bkp=bkm+double(j)*tox*bk;
524   bkm=bk;
525   bk=bkp;
526  }
527
528  return bk;
529 }
530 ///////////////////////////////////////////////////////////////////////////
531 Double_t AliMath::BesselI(Int_t n,Double_t x)
532 {
533 // Computation of the Integer Order Modified Bessel function I_n(x)
534 // for n=0,1,2,... and any real x.
535 //
536 // The algorithm uses the recurrence relation
537 //
538 //               I_n+1(x) = (-2n/x)*I_n(x) + I_n-1(x) 
539 //
540 // as denoted in Numerical Recipes 2nd ed. on p. 232 (W.H.Press et al.).
541 //
542 //--- NvE 12-mar-2000 UU-SAP Utrecht
543
544  Int_t iacc=40; // Increase to enhance accuracy
545  Double_t bigno=1.e10, bigni=1.e-10;
546
547  if (n < 0)
548  {
549   cout << " *BesselI* Invalid argument (n,x) = (" << n << " , " << x << ")" << endl;
550   return 0;
551  }
552
553  if (n==0) return BesselI0(x);
554
555  if (n==1) return BesselI1(x);
556
557  if (fabs(x) < 1.e-10) return 0;
558
559  Double_t tox=2./fabs(x);
560  Double_t bip=0,bim=0;
561  Double_t bi=1;
562  Double_t result=0;
563  Int_t m=2*((n+int(sqrt(float(iacc*n))))); // Downward recurrence from even m
564  for (Int_t j=m; j<=1; j--)
565  {
566   bim=bip+double(j)*tox*bi;
567   bip=bi;
568   bi=bim;
569   if (fabs(bi) > bigno) // Renormalise to prevent overflows
570   {
571    result*=bigni;
572    bi*=bigni;
573    bip*=bigni;
574   }
575   if (j==n) result=bip;
576  }
577
578  result*=BesselI0(x)/bi; // Normalise with I0(x)
579  if ((x < 0) && (n%2 == 1)) result=-result;
580
581  return result;
582 }
583 ///////////////////////////////////////////////////////////////////////////