1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 ///////////////////////////////////////////////////////////////////////////
20 // Various mathematical tools which may be very convenient while
21 // performing physics analysis.
23 // Example : Probability of a Chi-squared value
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
33 //--- Author: Nick van Eijndhoven 14-nov-1998 UU-SAP Utrecht
34 //- Modified: NvE $Date$ UU-SAP Utrecht
35 ///////////////////////////////////////////////////////////////////////////
38 #include "Riostream.h"
40 ClassImp(AliMath) // Class implementation to enable ROOT I/O
42 AliMath::AliMath() : TObject()
44 // Default constructor
46 ///////////////////////////////////////////////////////////////////////////
51 ///////////////////////////////////////////////////////////////////////////
52 AliMath::AliMath(AliMath& m) : TObject(m)
56 ///////////////////////////////////////////////////////////////////////////
57 Double_t AliMath::Gamma(Double_t z)
59 // Computation of gamma(z) for all z>0.
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.).
64 // [1] C.Lanczos, SIAM Journal of Numerical Analysis B1 (1964), 86.
66 //--- Nve 14-nov-1998 UU-SAP Utrecht
70 cout << "*Gamma(z)* Wrong argument z = " << z << endl;
74 Double_t v=LnGamma(z);
77 ///////////////////////////////////////////////////////////////////////////
78 Double_t AliMath::Gamma(Double_t a,Double_t x)
80 // Computation of the incomplete gamma function P(a,x)
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.).
85 //--- Nve 14-nov-1998 UU-SAP Utrecht
89 cout << "*Gamma(a,x)* Invalid argument a = " << a << endl;
95 if (x<0) cout << "*Gamma(a,x)* Invalid argument x = " << x << endl;
108 ///////////////////////////////////////////////////////////////////////////
109 Double_t AliMath::LnGamma(Double_t z)
111 // Computation of ln[gamma(z)] for all z>0.
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.).
116 // [1] C.Lanczos, SIAM Journal of Numerical Analysis B1 (1964), 86.
118 // The accuracy of the result is better than 2e-10.
120 //--- Nve 14-nov-1998 UU-SAP Utrecht
124 cout << "*LnGamma(z)* Wrong argument z = " << z << endl;
128 // Coefficients for the series expansion
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;
141 tmp=(x+0.5)*log(tmp)-tmp;
142 Double_t ser=1.000000000190015;
143 for (Int_t i=1; i<7; i++)
148 Double_t v=tmp+log(c[0]*ser/x);
151 ///////////////////////////////////////////////////////////////////////////
152 Double_t AliMath::GamSer(Double_t a,Double_t x)
154 // Computation of the incomplete gamma function P(a,x)
155 // via its series representation.
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.).
160 //--- Nve 14-nov-1998 UU-SAP Utrecht
162 Int_t itmax=100; // Maximum number of iterations
163 Double_t eps=3.e-7; // Relative accuracy
167 cout << "*GamSer(a,x)* Invalid argument a = " << a << endl;
173 if (x<0) cout << "*GamSer(a,x)* Invalid argument x = " << x << endl;
177 Double_t gln=LnGamma(a);
181 for (Int_t n=1; n<=itmax; n++)
186 if (fabs(del)<fabs(sum*eps)) break;
187 if (n==itmax) cout << "*GamSer(a,x)* a too large or itmax too small" << endl;
189 Double_t v=sum*exp(-x+a*log(x)-gln);
192 ///////////////////////////////////////////////////////////////////////////
193 Double_t AliMath::GamCf(Double_t a,Double_t x)
195 // Computation of the incomplete gamma function P(a,x)
196 // via its continued fraction representation.
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.).
201 //--- Nve 14-nov-1998 UU-SAP Utrecht
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
209 cout << "*GamCf(a,x)* Invalid argument a = " << a << endl;
215 if (x<0) cout << "*GamCf(a,x)* Invalid argument x = " << x << endl;
219 Double_t gln=LnGamma(a);
225 for (Int_t i=1; i<=itmax; i++)
227 an=double(-i)*(double(i)-a);
230 if (fabs(d)<fpmin) d=fpmin;
232 if (fabs(c)<fpmin) c=fpmin;
236 if (fabs(del-1.)<eps) break;
237 if (i==itmax) cout << "*GamCf(a,x)* a too large or itmax too small" << endl;
239 Double_t v=exp(-x+a*log(x)-gln)*h;
242 ///////////////////////////////////////////////////////////////////////////
243 Double_t AliMath::Erf(Double_t x)
245 // Computation of the error function erf(x).
247 //--- NvE 14-nov-1998 UU-SAP Utrecht
251 ///////////////////////////////////////////////////////////////////////////
252 Double_t AliMath::Erfc(Double_t x)
254 // Computation of the complementary error function erfc(x).
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.).
259 // The fractional error is always less than 1.2e-7.
261 //--- Nve 14-nov-1998 UU-SAP Utrecht
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;
270 Double_t v=1.; // The return value
274 if (z <= 0.) return v; // erfc(0)=1
276 Double_t t=1./(1.+0.5*z);
279 +ka1+t*(ka2+t*(ka3+t*(ka4+t*(ka5+t*(ka6+t*(ka7+t*(ka8+t*(ka9+t*ka10)))))))));
281 if (x < 0.) v=2.-v; // erfc(-x)=2-erfc(x)
285 ///////////////////////////////////////////////////////////////////////////
286 Double_t AliMath::Prob(Double_t chi2,Int_t ndf)
288 // Computation of the probability for a certain Chi-squared (chi2)
289 // and number of degrees of freedom (ndf).
291 // Calculations are based on the incomplete gamma function P(a,x),
292 // where a=ndf/2 and x=chi2/2.
294 // P(a,x) represents the probability that the observed Chi-squared
295 // for a correct model should be less than the value chi2.
297 // The returned probability corresponds to 1-P(a,x),
298 // which denotes the probability that an observed Chi-squared exceeds
299 // the value chi2 by chance, even for a correct model.
301 //--- NvE 14-nov-1998 UU-SAP Utrecht
303 if (ndf <= 0) return 0; // Set CL to zero in case ndf<=0
317 // Alternative which is exact
318 // This code may be activated in case the gamma function gives problems
321 // Double_t v=1.-Erf(sqrt(chi2)/sqrt(2.));
325 // Gaussian approximation for large ndf
326 // This code may be activated in case the gamma function shows a problem
327 // Double_t q=sqrt(2.*chi2)-sqrt(double(2*ndf-1));
330 // Double_t v=0.5*(1.-Erf(q/sqrt(2.)));
334 // Evaluate the incomplete gamma function
335 Double_t a=double(ndf)/2.;
337 return (1.-Gamma(a,x));
339 ///////////////////////////////////////////////////////////////////////////
340 Double_t AliMath::BesselI0(Double_t x)
342 // Computation of the modified Bessel function I_0(x) for any real x.
344 // The algorithm is based on the article by Abramowitz and Stegun [1]
345 // as denoted in Numerical Recipes 2nd ed. on p. 230 (W.H.Press et al.).
347 // [1] M.Abramowitz and I.A.Stegun, Handbook of Mathematical Functions,
348 // Applied Mathematics Series vol. 55 (1964), Washington.
350 //--- NvE 12-mar-2000 UU-SAP Utrecht
352 // Parameters of the polynomial approximation
353 const Double_t kp1=1.0, kp2=3.5156229, kp3=3.0899424,
354 kp4=1.2067492, kp5=0.2659732, kp6=3.60768e-2, kp7=4.5813e-3;
356 const Double_t kq1= 0.39894228, kq2= 1.328592e-2, kq3= 2.25319e-3,
357 kq4=-1.57565e-3, kq5= 9.16281e-3, kq6=-2.057706e-2,
358 kq7= 2.635537e-2, kq8=-1.647633e-2, kq9= 3.92377e-3;
362 Double_t y=0,result=0;
367 result=kp1+y*(kp2+y*(kp3+y*(kp4+y*(kp5+y*(kp6+y*kp7)))));
372 result=(exp(ax)/sqrt(ax))
373 *(kq1+y*(kq2+y*(kq3+y*(kq4+y*(kq5+y*(kq6+y*(kq7+y*(kq8+y*kq9))))))));
378 ///////////////////////////////////////////////////////////////////////////
379 Double_t AliMath::BesselK0(Double_t x)
381 // Computation of the modified Bessel function K_0(x) for positive real x.
383 // The algorithm is based on the article by Abramowitz and Stegun [1]
384 // as denoted in Numerical Recipes 2nd ed. on p. 230 (W.H.Press et al.).
386 // [1] M.Abramowitz and I.A.Stegun, Handbook of Mathematical Functions,
387 // Applied Mathematics Series vol. 55 (1964), Washington.
389 //--- NvE 12-mar-2000 UU-SAP Utrecht
391 // Parameters of the polynomial approximation
392 const Double_t kp1=-0.57721566, kp2=0.42278420, kp3=0.23069756,
393 kp4= 3.488590e-2, kp5=2.62698e-3, kp6=1.0750e-4, kp7=7.4e-5;
395 const Double_t kq1= 1.25331414, kq2=-7.832358e-2, kq3= 2.189568e-2,
396 kq4=-1.062446e-2, kq5= 5.87872e-3, kq6=-2.51540e-3, kq7=5.3208e-4;
400 cout << " *BesselK0* Invalid argument x = " << x << endl;
404 Double_t y=0,result=0;
409 result=(-log(x/2.)*BesselI0(x))
410 +(kp1+y*(kp2+y*(kp3+y*(kp4+y*(kp5+y*(kp6+y*kp7))))));
415 result=(exp(-x)/sqrt(x))
416 *(kq1+y*(kq2+y*(kq3+y*(kq4+y*(kq5+y*(kq6+y*kq7))))));
421 ///////////////////////////////////////////////////////////////////////////
422 Double_t AliMath::BesselI1(Double_t x)
424 // Computation of the modified Bessel function I_1(x) for any real x.
426 // The algorithm is based on the article by Abramowitz and Stegun [1]
427 // as denoted in Numerical Recipes 2nd ed. on p. 230 (W.H.Press et al.).
429 // [1] M.Abramowitz and I.A.Stegun, Handbook of Mathematical Functions,
430 // Applied Mathematics Series vol. 55 (1964), Washington.
432 //--- NvE 12-mar-2000 UU-SAP Utrecht
434 // Parameters of the polynomial approximation
435 const Double_t kp1=0.5, kp2=0.87890594, kp3=0.51498869,
436 kp4=0.15084934, kp5=2.658733e-2, kp6=3.01532e-3, kp7=3.2411e-4;
438 const Double_t kq1= 0.39894228, kq2=-3.988024e-2, kq3=-3.62018e-3,
439 kq4= 1.63801e-3, kq5=-1.031555e-2, kq6= 2.282967e-2,
440 kq7=-2.895312e-2, kq8= 1.787654e-2, kq9=-4.20059e-3;
444 Double_t y=0,result=0;
449 result=x*(kp1+y*(kp2+y*(kp3+y*(kp4+y*(kp5+y*(kp6+y*kp7))))));
454 result=(exp(ax)/sqrt(ax))
455 *(kq1+y*(kq2+y*(kq3+y*(kq4+y*(kq5+y*(kq6+y*(kq7+y*(kq8+y*kq9))))))));
456 if (x < 0) result=-result;
461 ///////////////////////////////////////////////////////////////////////////
462 Double_t AliMath::BesselK1(Double_t x)
464 // Computation of the modified Bessel function K_1(x) for positive real x.
466 // The algorithm is based on the article by Abramowitz and Stegun [1]
467 // as denoted in Numerical Recipes 2nd ed. on p. 230 (W.H.Press et al.).
469 // [1] M.Abramowitz and I.A.Stegun, Handbook of Mathematical Functions,
470 // Applied Mathematics Series vol. 55 (1964), Washington.
472 //--- NvE 12-mar-2000 UU-SAP Utrecht
474 // Parameters of the polynomial approximation
475 const Double_t kp1= 1., kp2= 0.15443144, kp3=-0.67278579,
476 kp4=-0.18156897, kp5=-1.919402e-2, kp6=-1.10404e-3, kp7=-4.686e-5;
478 const Double_t kq1= 1.25331414, kq2= 0.23498619, kq3=-3.655620e-2,
479 kq4= 1.504268e-2, kq5=-7.80353e-3, kq6= 3.25614e-3, kq7=-6.8245e-4;
483 cout << " *BesselK1* Invalid argument x = " << x << endl;
487 Double_t y=0,result=0;
492 result=(log(x/2.)*BesselI1(x))
493 +(1./x)*(kp1+y*(kp2+y*(kp3+y*(kp4+y*(kp5+y*(kp6+y*kp7))))));
498 result=(exp(-x)/sqrt(x))
499 *(kq1+y*(kq2+y*(kq3+y*(kq4+y*(kq5+y*(kq6+y*kq7))))));
504 ///////////////////////////////////////////////////////////////////////////
505 Double_t AliMath::BesselK(Int_t n,Double_t x)
507 // Computation of the Integer Order Modified Bessel function K_n(x)
508 // for n=0,1,2,... and positive real x.
510 // The algorithm uses the recurrence relation
512 // K_n+1(x) = (2n/x)*K_n(x) + K_n-1(x)
514 // as denoted in Numerical Recipes 2nd ed. on p. 232 (W.H.Press et al.).
516 //--- NvE 12-mar-2000 UU-SAP Utrecht
520 cout << " *BesselK* Invalid argument(s) (n,x) = (" << n << " , " << x << ")" << endl;
524 if (n==0) return BesselK0(x);
526 if (n==1) return BesselK1(x);
528 // Perform upward recurrence for all x
530 Double_t bkm=BesselK0(x);
531 Double_t bk=BesselK1(x);
533 for (Int_t j=1; j<n; j++)
535 bkp=bkm+double(j)*tox*bk;
542 ///////////////////////////////////////////////////////////////////////////
543 Double_t AliMath::BesselI(Int_t n,Double_t x)
545 // Computation of the Integer Order Modified Bessel function I_n(x)
546 // for n=0,1,2,... and any real x.
548 // The algorithm uses the recurrence relation
550 // I_n+1(x) = (-2n/x)*I_n(x) + I_n-1(x)
552 // as denoted in Numerical Recipes 2nd ed. on p. 232 (W.H.Press et al.).
554 //--- NvE 12-mar-2000 UU-SAP Utrecht
556 Int_t iacc=40; // Increase to enhance accuracy
557 Double_t bigno=1.e10, bigni=1.e-10;
561 cout << " *BesselI* Invalid argument (n,x) = (" << n << " , " << x << ")" << endl;
565 if (n==0) return BesselI0(x);
567 if (n==1) return BesselI1(x);
569 if (fabs(x) < 1.e-10) return 0;
571 Double_t tox=2./fabs(x);
572 Double_t bip=0,bim=0;
575 Int_t m=2*((n+int(sqrt(float(iacc*n))))); // Downward recurrence from even m
576 for (Int_t j=m; j<=1; j--)
578 bim=bip+double(j)*tox*bi;
581 if (fabs(bi) > bigno) // Renormalise to prevent overflows
587 if (j==n) result=bip;
590 result*=BesselI0(x)/bi; // Normalise with I0(x)
591 if ((x < 0) && (n%2 == 1)) result=-result;
595 ///////////////////////////////////////////////////////////////////////////