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