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 ///////////////////////////////////////////////////////////////////////////
39 ClassImp(AliMath) // Class implementation to enable ROOT I/O
43 // Default constructor
45 ///////////////////////////////////////////////////////////////////////////
50 ///////////////////////////////////////////////////////////////////////////
51 Double_t AliMath::Gamma(Double_t z)
53 // Computation of gamma(z) for all z>0.
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.).
58 // [1] C.Lanczos, SIAM Journal of Numerical Analysis B1 (1964), 86.
60 //--- Nve 14-nov-1998 UU-SAP Utrecht
64 cout << "*Gamma(z)* Wrong argument z = " << z << endl;
68 Double_t v=LnGamma(z);
71 ///////////////////////////////////////////////////////////////////////////
72 Double_t AliMath::Gamma(Double_t a,Double_t x)
74 // Computation of the incomplete gamma function P(a,x)
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.).
79 //--- Nve 14-nov-1998 UU-SAP Utrecht
83 cout << "*Gamma(a,x)* Invalid argument a = " << a << endl;
89 if (x<0) cout << "*Gamma(a,x)* Invalid argument x = " << x << endl;
102 ///////////////////////////////////////////////////////////////////////////
103 Double_t AliMath::LnGamma(Double_t z)
105 // Computation of ln[gamma(z)] for all z>0.
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.).
110 // [1] C.Lanczos, SIAM Journal of Numerical Analysis B1 (1964), 86.
112 // The accuracy of the result is better than 2e-10.
114 //--- Nve 14-nov-1998 UU-SAP Utrecht
118 cout << "*LnGamma(z)* Wrong argument z = " << z << endl;
122 // Coefficients for the series expansion
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;
135 tmp=(x+0.5)*log(tmp)-tmp;
136 Double_t ser=1.000000000190015;
137 for (Int_t i=1; i<7; i++)
142 Double_t v=tmp+log(c[0]*ser/x);
145 ///////////////////////////////////////////////////////////////////////////
146 Double_t AliMath::GamSer(Double_t a,Double_t x)
148 // Computation of the incomplete gamma function P(a,x)
149 // via its series representation.
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.).
154 //--- Nve 14-nov-1998 UU-SAP Utrecht
156 Int_t itmax=100; // Maximum number of iterations
157 Double_t eps=3.e-7; // Relative accuracy
161 cout << "*GamSer(a,x)* Invalid argument a = " << a << endl;
167 if (x<0) cout << "*GamSer(a,x)* Invalid argument x = " << x << endl;
171 Double_t gln=LnGamma(a);
175 for (Int_t n=1; n<=itmax; n++)
180 if (fabs(del)<fabs(sum*eps)) break;
181 if (n==itmax) cout << "*GamSer(a,x)* a too large or itmax too small" << endl;
183 Double_t v=sum*exp(-x+a*log(x)-gln);
186 ///////////////////////////////////////////////////////////////////////////
187 Double_t AliMath::GamCf(Double_t a,Double_t x)
189 // Computation of the incomplete gamma function P(a,x)
190 // via its continued fraction representation.
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.).
195 //--- Nve 14-nov-1998 UU-SAP Utrecht
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
203 cout << "*GamCf(a,x)* Invalid argument a = " << a << endl;
209 if (x<0) cout << "*GamCf(a,x)* Invalid argument x = " << x << endl;
213 Double_t gln=LnGamma(a);
219 for (Int_t i=1; i<=itmax; i++)
221 an=double(-i)*(double(i)-a);
224 if (fabs(d)<fpmin) d=fpmin;
226 if (fabs(c)<fpmin) c=fpmin;
230 if (fabs(del-1.)<eps) break;
231 if (i==itmax) cout << "*GamCf(a,x)* a too large or itmax too small" << endl;
233 Double_t v=exp(-x+a*log(x)-gln)*h;
236 ///////////////////////////////////////////////////////////////////////////
237 Double_t AliMath::Erf(Double_t x)
239 // Computation of the error function erf(x).
241 //--- NvE 14-nov-1998 UU-SAP Utrecht
245 ///////////////////////////////////////////////////////////////////////////
246 Double_t AliMath::Erfc(Double_t x)
248 // Computation of the complementary error function erfc(x).
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.).
253 // The fractional error is always less than 1.2e-7.
255 //--- Nve 14-nov-1998 UU-SAP Utrecht
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;
264 Double_t v=1.; // The return value
268 if (z <= 0.) return v; // erfc(0)=1
270 Double_t t=1./(1.+0.5*z);
273 +a1+t*(a2+t*(a3+t*(a4+t*(a5+t*(a6+t*(a7+t*(a8+t*(a9+t*a10)))))))));
275 if (x < 0.) v=2.-v; // erfc(-x)=2-erfc(x)
279 ///////////////////////////////////////////////////////////////////////////
280 Double_t AliMath::Prob(Double_t chi2,Int_t ndf)
282 // Computation of the probability for a certain Chi-squared (chi2)
283 // and number of degrees of freedom (ndf).
285 // Calculations are based on the incomplete gamma function P(a,x),
286 // where a=ndf/2 and x=chi2/2.
288 // P(a,x) represents the probability that the observed Chi-squared
289 // for a correct model should be less than the value chi2.
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.
295 //--- NvE 14-nov-1998 UU-SAP Utrecht
297 if (ndf <= 0) return 0; // Set CL to zero in case ndf<=0
311 // Alternative which is exact
312 // This code may be activated in case the gamma function gives problems
315 // Double_t v=1.-Erf(sqrt(chi2)/sqrt(2.));
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));
324 // Double_t v=0.5*(1.-Erf(q/sqrt(2.)));
328 // Evaluate the incomplete gamma function
329 Double_t a=double(ndf)/2.;
331 return (1.-Gamma(a,x));
333 ///////////////////////////////////////////////////////////////////////////
334 Double_t AliMath::BesselI0(Double_t x)
336 // Computation of the modified Bessel function I_0(x) for any real x.
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.).
341 // [1] M.Abramowitz and I.A.Stegun, Handbook of Mathematical Functions,
342 // Applied Mathematics Series vol. 55 (1964), Washington.
344 //--- NvE 12-mar-2000 UU-SAP Utrecht
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;
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;
356 Double_t y=0,result=0;
361 result=p1+y*(p2+y*(p3+y*(p4+y*(p5+y*(p6+y*p7)))));
366 result=(exp(ax)/sqrt(ax))*(q1+y*(q2+y*(q3+y*(q4+y*(q5+y*(q6+y*(q7+y*(q8+y*q9))))))));
371 ///////////////////////////////////////////////////////////////////////////
372 Double_t AliMath::BesselK0(Double_t x)
374 // Computation of the modified Bessel function K_0(x) for positive real x.
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.).
379 // [1] M.Abramowitz and I.A.Stegun, Handbook of Mathematical Functions,
380 // Applied Mathematics Series vol. 55 (1964), Washington.
382 //--- NvE 12-mar-2000 UU-SAP Utrecht
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;
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;
393 cout << " *BesselK0* Invalid argument x = " << x << endl;
397 Double_t y=0,result=0;
402 result=(-log(x/2.)*BesselI0(x))+(p1+y*(p2+y*(p3+y*(p4+y*(p5+y*(p6+y*p7))))));
407 result=(exp(-x)/sqrt(x))*(q1+y*(q2+y*(q3+y*(q4+y*(q5+y*(q6+y*q7))))));
412 ///////////////////////////////////////////////////////////////////////////
413 Double_t AliMath::BesselI1(Double_t x)
415 // Computation of the modified Bessel function I_1(x) for any real x.
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.).
420 // [1] M.Abramowitz and I.A.Stegun, Handbook of Mathematical Functions,
421 // Applied Mathematics Series vol. 55 (1964), Washington.
423 //--- NvE 12-mar-2000 UU-SAP Utrecht
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;
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;
435 Double_t y=0,result=0;
440 result=x*(p1+y*(p2+y*(p3+y*(p4+y*(p5+y*(p6+y*p7))))));
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;
451 ///////////////////////////////////////////////////////////////////////////
452 Double_t AliMath::BesselK1(Double_t x)
454 // Computation of the modified Bessel function K_1(x) for positive real x.
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.).
459 // [1] M.Abramowitz and I.A.Stegun, Handbook of Mathematical Functions,
460 // Applied Mathematics Series vol. 55 (1964), Washington.
462 //--- NvE 12-mar-2000 UU-SAP Utrecht
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;
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;
473 cout << " *BesselK1* Invalid argument x = " << x << endl;
477 Double_t y=0,result=0;
482 result=(log(x/2.)*BesselI1(x))+(1./x)*(p1+y*(p2+y*(p3+y*(p4+y*(p5+y*(p6+y*p7))))));
487 result=(exp(-x)/sqrt(x))*(q1+y*(q2+y*(q3+y*(q4+y*(q5+y*(q6+y*q7))))));
492 ///////////////////////////////////////////////////////////////////////////
493 Double_t AliMath::BesselK(Int_t n,Double_t x)
495 // Computation of the Integer Order Modified Bessel function K_n(x)
496 // for n=0,1,2,... and positive real x.
498 // The algorithm uses the recurrence relation
500 // K_n+1(x) = (2n/x)*K_n(x) + K_n-1(x)
502 // as denoted in Numerical Recipes 2nd ed. on p. 232 (W.H.Press et al.).
504 //--- NvE 12-mar-2000 UU-SAP Utrecht
508 cout << " *BesselK* Invalid argument(s) (n,x) = (" << n << " , " << x << ")" << endl;
512 if (n==0) return BesselK0(x);
514 if (n==1) return BesselK1(x);
516 // Perform upward recurrence for all x
518 Double_t bkm=BesselK0(x);
519 Double_t bk=BesselK1(x);
521 for (Int_t j=1; j<n; j++)
523 bkp=bkm+double(j)*tox*bk;
530 ///////////////////////////////////////////////////////////////////////////
531 Double_t AliMath::BesselI(Int_t n,Double_t x)
533 // Computation of the Integer Order Modified Bessel function I_n(x)
534 // for n=0,1,2,... and any real x.
536 // The algorithm uses the recurrence relation
538 // I_n+1(x) = (-2n/x)*I_n(x) + I_n-1(x)
540 // as denoted in Numerical Recipes 2nd ed. on p. 232 (W.H.Press et al.).
542 //--- NvE 12-mar-2000 UU-SAP Utrecht
544 Int_t iacc=40; // Increase to enhance accuracy
545 Double_t bigno=1.e10, bigni=1.e-10;
549 cout << " *BesselI* Invalid argument (n,x) = (" << n << " , " << x << ")" << endl;
553 if (n==0) return BesselI0(x);
555 if (n==1) return BesselI1(x);
557 if (fabs(x) < 1.e-10) return 0;
559 Double_t tox=2./fabs(x);
560 Double_t bip=0,bim=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--)
566 bim=bip+double(j)*tox*bi;
569 if (fabs(bi) > bigno) // Renormalise to prevent overflows
575 if (j==n) result=bip;
578 result*=BesselI0(x)/bi; // Normalise with I0(x)
579 if ((x < 0) && (n%2 == 1)) result=-result;
583 ///////////////////////////////////////////////////////////////////////////