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 Revision 1.3 1999/11/03 14:23:17 fca
19 New version of RALICE introduced
21 Revision 1.2 1999/09/29 09:24:28 fca
22 Introduction of the Copyright and cvs Log
26 ///////////////////////////////////////////////////////////////////////////
28 // Various mathematical tools which may be very convenient while
29 // performing physics analysis.
31 // Example : Probability of a Chi-squared value
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
41 //--- Author: Nick van Eijndhoven 14-nov-1998 UU-SAP Utrecht
42 //- Modified: NvE 12-mar-2000 UU-SAP Utrecht
43 ///////////////////////////////////////////////////////////////////////////
47 ClassImp(AliMath) // Class implementation to enable ROOT I/O
51 // Default constructor
53 ///////////////////////////////////////////////////////////////////////////
58 ///////////////////////////////////////////////////////////////////////////
59 Double_t AliMath::Gamma(Double_t z)
61 // Computation of gamma(z) for all z>0.
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.).
66 // [1] C.Lanczos, SIAM Journal of Numerical Analysis B1 (1964), 86.
68 //--- Nve 14-nov-1998 UU-SAP Utrecht
72 cout << "*Gamma(z)* Wrong argument z = " << z << endl;
76 Double_t v=LnGamma(z);
79 ///////////////////////////////////////////////////////////////////////////
80 Double_t AliMath::Gamma(Double_t a,Double_t x)
82 // Computation of the incomplete gamma function P(a,x)
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.).
87 //--- Nve 14-nov-1998 UU-SAP Utrecht
91 cout << "*Gamma(a,x)* Invalid argument a = " << a << endl;
97 if (x<0) cout << "*Gamma(a,x)* Invalid argument x = " << x << endl;
110 ///////////////////////////////////////////////////////////////////////////
111 Double_t AliMath::LnGamma(Double_t z)
113 // Computation of ln[gamma(z)] for all z>0.
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.).
118 // [1] C.Lanczos, SIAM Journal of Numerical Analysis B1 (1964), 86.
120 // The accuracy of the result is better than 2e-10.
122 //--- Nve 14-nov-1998 UU-SAP Utrecht
126 cout << "*LnGamma(z)* Wrong argument z = " << z << endl;
130 // Coefficients for the series expansion
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;
143 tmp=(x+0.5)*log(tmp)-tmp;
144 Double_t ser=1.000000000190015;
145 for (Int_t i=1; i<7; i++)
150 Double_t v=tmp+log(c[0]*ser/x);
153 ///////////////////////////////////////////////////////////////////////////
154 Double_t AliMath::GamSer(Double_t a,Double_t x)
156 // Computation of the incomplete gamma function P(a,x)
157 // via its series representation.
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.).
162 //--- Nve 14-nov-1998 UU-SAP Utrecht
164 Int_t itmax=100; // Maximum number of iterations
165 Double_t eps=3.e-7; // Relative accuracy
169 cout << "*GamSer(a,x)* Invalid argument a = " << a << endl;
175 if (x<0) cout << "*GamSer(a,x)* Invalid argument x = " << x << endl;
179 Double_t gln=LnGamma(a);
183 for (Int_t n=1; n<=itmax; n++)
188 if (fabs(del)<fabs(sum*eps)) break;
189 if (n==itmax) cout << "*GamSer(a,x)* a too large or itmax too small" << endl;
191 Double_t v=sum*exp(-x+a*log(x)-gln);
194 ///////////////////////////////////////////////////////////////////////////
195 Double_t AliMath::GamCf(Double_t a,Double_t x)
197 // Computation of the incomplete gamma function P(a,x)
198 // via its continued fraction representation.
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.).
203 //--- Nve 14-nov-1998 UU-SAP Utrecht
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
211 cout << "*GamCf(a,x)* Invalid argument a = " << a << endl;
217 if (x<0) cout << "*GamCf(a,x)* Invalid argument x = " << x << endl;
221 Double_t gln=LnGamma(a);
227 for (Int_t i=1; i<=itmax; i++)
229 an=double(-i)*(double(i)-a);
232 if (fabs(d)<fpmin) d=fpmin;
234 if (fabs(c)<fpmin) c=fpmin;
238 if (fabs(del-1.)<eps) break;
239 if (i==itmax) cout << "*GamCf(a,x)* a too large or itmax too small" << endl;
241 Double_t v=exp(-x+a*log(x)-gln)*h;
244 ///////////////////////////////////////////////////////////////////////////
245 Double_t AliMath::Erf(Double_t x)
247 // Computation of the error function erf(x).
249 //--- NvE 14-nov-1998 UU-SAP Utrecht
253 ///////////////////////////////////////////////////////////////////////////
254 Double_t AliMath::Erfc(Double_t x)
256 // Computation of the complementary error function erfc(x).
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.).
261 // The fractional error is always less than 1.2e-7.
263 //--- Nve 14-nov-1998 UU-SAP Utrecht
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;
272 Double_t v=1.; // The return value
276 if (z <= 0.) return v; // erfc(0)=1
278 Double_t t=1./(1.+0.5*z);
281 +a1+t*(a2+t*(a3+t*(a4+t*(a5+t*(a6+t*(a7+t*(a8+t*(a9+t*a10)))))))));
283 if (x < 0.) v=2.-v; // erfc(-x)=2-erfc(x)
287 ///////////////////////////////////////////////////////////////////////////
288 Double_t AliMath::Prob(Double_t chi2,Int_t ndf)
290 // Computation of the probability for a certain Chi-squared (chi2)
291 // and number of degrees of freedom (ndf).
293 // Calculations are based on the incomplete gamma function P(a,x),
294 // where a=ndf/2 and x=chi2/2.
296 // P(a,x) represents the probability that the observed Chi-squared
297 // for a correct model should be less than the value chi2.
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.
303 //--- NvE 14-nov-1998 UU-SAP Utrecht
305 if (ndf <= 0) return 0; // Set CL to zero in case ndf<=0
319 // Alternative which is exact
320 // This code may be activated in case the gamma function gives problems
323 // Double_t v=1.-Erf(sqrt(chi2)/sqrt(2.));
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));
332 // Double_t v=0.5*(1.-Erf(q/sqrt(2.)));
336 // Evaluate the incomplete gamma function
337 Double_t a=double(ndf)/2.;
339 return (1.-Gamma(a,x));
341 ///////////////////////////////////////////////////////////////////////////
342 Double_t AliMath::BesselI0(Double_t x)
344 // Computation of the modified Bessel function I_0(x) for any real x.
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.).
349 // [1] M.Abramowitz and I.A.Stegun, Handbook of Mathematical Functions,
350 // Applied Mathematics Series vol. 55 (1964), Washington.
352 //--- NvE 12-mar-2000 UU-SAP Utrecht
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;
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;
364 Double_t y=0,result=0;
369 result=p1+y*(p2+y*(p3+y*(p4+y*(p5+y*(p6+y*p7)))));
374 result=(exp(ax)/sqrt(ax))*(q1+y*(q2+y*(q3+y*(q4+y*(q5+y*(q6+y*(q7+y*(q8+y*q9))))))));
379 ///////////////////////////////////////////////////////////////////////////
380 Double_t AliMath::BesselK0(Double_t x)
382 // Computation of the modified Bessel function K_0(x) for positive real x.
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.).
387 // [1] M.Abramowitz and I.A.Stegun, Handbook of Mathematical Functions,
388 // Applied Mathematics Series vol. 55 (1964), Washington.
390 //--- NvE 12-mar-2000 UU-SAP Utrecht
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;
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;
401 cout << " *BesselK0* Invalid argument x = " << x << endl;
405 Double_t y=0,result=0;
410 result=(-log(x/2.)*BesselI0(x))+(p1+y*(p2+y*(p3+y*(p4+y*(p5+y*(p6+y*p7))))));
415 result=(exp(-x)/sqrt(x))*(q1+y*(q2+y*(q3+y*(q4+y*(q5+y*(q6+y*q7))))));
420 ///////////////////////////////////////////////////////////////////////////
421 Double_t AliMath::BesselI1(Double_t x)
423 // Computation of the modified Bessel function I_1(x) for any real x.
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.).
428 // [1] M.Abramowitz and I.A.Stegun, Handbook of Mathematical Functions,
429 // Applied Mathematics Series vol. 55 (1964), Washington.
431 //--- NvE 12-mar-2000 UU-SAP Utrecht
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;
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;
443 Double_t y=0,result=0;
448 result=x*(p1+y*(p2+y*(p3+y*(p4+y*(p5+y*(p6+y*p7))))));
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;
459 ///////////////////////////////////////////////////////////////////////////
460 Double_t AliMath::BesselK1(Double_t x)
462 // Computation of the modified Bessel function K_1(x) for positive real x.
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.).
467 // [1] M.Abramowitz and I.A.Stegun, Handbook of Mathematical Functions,
468 // Applied Mathematics Series vol. 55 (1964), Washington.
470 //--- NvE 12-mar-2000 UU-SAP Utrecht
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;
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;
481 cout << " *BesselK1* Invalid argument x = " << x << endl;
485 Double_t y=0,result=0;
490 result=(log(x/2.)*BesselI1(x))+(1./x)*(p1+y*(p2+y*(p3+y*(p4+y*(p5+y*(p6+y*p7))))));
495 result=(exp(-x)/sqrt(x))*(q1+y*(q2+y*(q3+y*(q4+y*(q5+y*(q6+y*q7))))));
500 ///////////////////////////////////////////////////////////////////////////
501 Double_t AliMath::BesselK(Int_t n,Double_t x)
503 // Computation of the Integer Order Modified Bessel function K_n(x)
504 // for n=0,1,2,... and positive real x.
506 // The algorithm uses the recurrence relation
508 // K_n+1(x) = (2n/x)*K_n(x) + K_n-1(x)
510 // as denoted in Numerical Recipes 2nd ed. on p. 232 (W.H.Press et al.).
512 //--- NvE 12-mar-2000 UU-SAP Utrecht
516 cout << " *BesselK* Invalid argument(s) (n,x) = (" << n << " , " << x << ")" << endl;
520 if (n==0) return BesselK0(x);
522 if (n==1) return BesselK1(x);
524 // Perform upward recurrence for all x
526 Double_t bkm=BesselK0(x);
527 Double_t bk=BesselK1(x);
529 for (Int_t j=1; j<n; j++)
531 bkp=bkm+double(j)*tox*bk;
538 ///////////////////////////////////////////////////////////////////////////
539 Double_t AliMath::BesselI(Int_t n,Double_t x)
541 // Computation of the Integer Order Modified Bessel function I_n(x)
542 // for n=0,1,2,... and any real x.
544 // The algorithm uses the recurrence relation
546 // I_n+1(x) = (-2n/x)*I_n(x) + I_n-1(x)
548 // as denoted in Numerical Recipes 2nd ed. on p. 232 (W.H.Press et al.).
550 //--- NvE 12-mar-2000 UU-SAP Utrecht
552 Int_t iacc=40; // Increase to enhance accuracy
553 Double_t bigno=1.e10, bigni=1.e-10;
557 cout << " *BesselI* Invalid argument (n,x) = (" << n << " , " << x << ")" << endl;
561 if (n==0) return BesselI0(x);
563 if (n==1) return BesselI1(x);
565 if (fabs(x) < 1.e-10) return 0;
567 Double_t tox=2./fabs(x);
568 Double_t bip=0,bim=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--)
574 bim=bip+double(j)*tox*bi;
577 if (fabs(bi) > bigno) // Renormalise to prevent overflows
583 if (j==n) result=bip;
586 result*=BesselI0(x)/bi; // Normalise with I0(x)
587 if ((x < 0) && (n%2 == 1)) result=-result;
591 ///////////////////////////////////////////////////////////////////////////