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