ROOT  6.06/08
Reference Guide
KelvinFunctions.cxx
Go to the documentation of this file.
1 // @(#)root/mathmore:$Id$
2 // The functions in this class have been imported by Jason Detwiler (jasondet@gmail.com) from
3 // CodeCogs GNU General Public License Agreement
4 // Copyright (C) 2004-2005 CodeCogs, Zyba Ltd, Broadwood, Holford, TA5 1DU,
5 // England.
6 //
7 // This program is free software; you can redistribute it and/or modify it
8 // under
9 // the terms of the GNU General Public License as published by CodeCogs.
10 // You must retain a copy of this licence in all copies.
11 //
12 // This program is distributed in the hope that it will be useful, but
13 // WITHOUT ANY
14 // WARRANTY; without even the implied warranty of MERCHANTABILITY or
15 // FITNESS FOR A
16 // PARTICULAR PURPOSE. See the Adapted GNU General Public License for more
17 // details.
18 //
19 // *** THIS SOFTWARE CAN NOT BE USED FOR COMMERCIAL GAIN. ***
20 // ---------------------------------------------------------------------------------
21 
22 #include "Math/KelvinFunctions.h"
23 #include <math.h>
24 
25 
26 namespace ROOT {
27 namespace Math {
28 
29 double KelvinFunctions::fgMin = 20;
30 double KelvinFunctions::fgEpsilon = 1.e-20;
31 
32 double kSqrt2 = 1.4142135623730950488016887242097;
33 double kPi = 3.14159265358979323846;
34 double kEulerGamma = 0.577215664901532860606512090082402431042;
35 
36 
37 /**
38 \class KelvinFunctions
39 
40 This class calculates the Kelvin functions Ber(x), Bei(x), Ker(x),
41 Kei(x), and their first derivatives.
42 */
43 
44 ////////////////////////////////////////////////////////////////////////////////
45 /// \f[
46 /// Ber(x) = Ber_{0}(x) = Re\left[J_{0}\left(x e^{3\pi i/4}\right)\right]
47 /// \f]
48 /// where x is real, and \f$J_{0}(z)\f$ is the zeroth-order Bessel
49 /// function of the first kind.
50 ///
51 /// If x < fgMin (=20), Ber(x) is computed according to its polynomial
52 /// approximation
53 /// \f[
54 /// Ber(x) = 1 + \sum_{n \geq 1}\frac{(-1)^{n}(x/2)^{4n}}{[(2n)!]^{2}}
55 /// \f]
56 /// For x > fgMin, Ber(x) is computed according to its asymptotic
57 /// expansion:
58 /// \f[
59 /// Ber(x) = \frac{e^{x/\sqrt{2}}}{\sqrt{2\pi x}} [F1(x) cos\alpha + G1(x) sin\alpha] - \frac{1}{\pi}Kei(x)
60 /// \f]
61 /// where \f$\alpha = \frac{x}{\sqrt{2}} - \frac{\pi}{8}\f$.
62 ///
63 /// See also F1() and G1().
64 ///
65 /// Begin_Macro
66 /// {
67 /// TCanvas *c = new TCanvas("c","c",0,0,500,300);
68 /// TF1 *fBer = new TF1("fBer","ROOT::Math::KelvinFunctions::Ber(x)",-10,10);
69 /// fBer->Draw();
70 /// return c;
71 /// }
72 /// End_Macro
73 
74 double KelvinFunctions::Ber(double x)
75 {
76  if (fabs(x) < fgEpsilon) return 1;
77 
78  if (fabs(x) < fgMin) {
79  double sum, factorial = 1, n = 1;
80  double term = 1, x_factor = x * x * x * x * 0.0625;
81 
82  sum = 1;
83 
84  do {
85  factorial = 4 * n * n * (2 * n - 1) * (2 * n - 1);
86  term *= (-1) / factorial * x_factor;
87  sum += term;
88  n += 1;
89  if (n > 1000) break;
90  } while (fabs(term) > fgEpsilon * sum);
91 
92  return sum;
93  } else {
94  double alpha = x / kSqrt2 - kPi / 8;
95  double value = F1(x) * cos(alpha) + G1(x) * sin(alpha);
96  value *= exp(x / kSqrt2) / sqrt(2 * kPi * x);
97  value -= Kei(x) / kPi;
98  return value;
99  }
100 }
101 
102 ////////////////////////////////////////////////////////////////////////////////
103 /// \f[
104 /// Bei(x) = Bei_{0}(x) = Im\left[J_{0}\left(x e^{3\pi i/4}\right)\right]
105 /// \f]
106 /// where x is real, and \f$J_{0}(z)\f$ is the zeroth-order Bessel
107 /// function of the first kind.
108 ///
109 /// If x < fgMin (=20), Bei(x) is computed according to its polynomial
110 /// approximation
111 /// \f[
112 /// Bei(x) = \sum_{n \geq 0}\frac{(-1)^{n}(x/2)^{4n+2}}{[(2n+1)!]^{2}}
113 /// \f]
114 /// For x > fgMin, Bei(x) is computed according to its asymptotic
115 /// expansion:
116 /// \f[
117 /// Bei(x) = \frac{e^{x/\sqrt{2}}}{\sqrt{2\pi x}} [F1(x) sin\alpha + G1(x) cos\alpha] - \frac{1}{\pi}Ker(x)
118 /// \f]
119 /// where \f$\alpha = \frac{x}{\sqrt{2}} - \frac{\pi}{8}\f$.
120 ///
121 /// See also F1() and G1().
122 ///
123 /// Begin_Macro
124 /// {
125 /// TCanvas *c = new TCanvas("c","c",0,0,500,300);
126 /// TF1 *fBei = new TF1("fBei","ROOT::Math::KelvinFunctions::Bei(x)",-10,10);
127 /// fBei->Draw();
128 /// return c;
129 /// }
130 /// End_Macro
131 
132 double KelvinFunctions::Bei(double x)
133 {
134 
135  if (fabs(x) < fgEpsilon) return 0;
136 
137  if (fabs(x) < fgMin) {
138  double sum, factorial = 1, n = 1;
139  double term = x * x * 0.25, x_factor = term * term;
140 
141  sum = term;
142 
143  do {
144  factorial = 4 * n * n * (2 * n + 1) * (2 * n + 1);
145  term *= (-1) / factorial * x_factor;
146  sum += term;
147  n += 1;
148  if (n > 1000) break;
149  } while (fabs(term) > fgEpsilon * sum);
150 
151  return sum;
152  } else {
153  double alpha = x / kSqrt2 - kPi / 8;
154  double value = F1(x) * sin(alpha) + G1(x) * cos(alpha);
155  value *= exp(x / kSqrt2) / sqrt(2 * kPi * x);
156  value += Ker(x) / kPi;
157  return value;
158  }
159 }
160 
161 
162 ////////////////////////////////////////////////////////////////////////////////
163 /// \f[
164 /// Ker(x) = Ker_{0}(x) = Re\left[K_{0}\left(x e^{3\pi i/4}\right)\right]
165 /// \f]
166 /// where x is real, and \f$K_{0}(z)\f$ is the zeroth-order modified
167 /// Bessel function of the second kind.
168 ///
169 /// If x < fgMin (=20), Ker(x) is computed according to its polynomial
170 /// approximation
171 /// \f[
172 /// Ker(x) = -\left(ln \frac{|x|}{2} + \gamma\right) Ber(x) + \left(\frac{\pi}{4} - \delta\right) Bei(x) + \sum_{n \geq 0} \frac{(-1)^{n}}{[(2n)!]^{2}} H_{2n} \left(\frac{x}{2}\right)^{4n}
173 /// \f]
174 /// where \f$\gamma = 0.577215664...\f$ is the Euler-Mascheroni constant,
175 /// \f$\delta = \pi\f$ for x < 0 and is otherwise zero, and
176 /// \f[
177 /// H_{n} = \sum_{k = 1}^{n} \frac{1}{k}
178 /// \f]
179 /// For x > fgMin, Ker(x) is computed according to its asymptotic
180 /// expansion:
181 /// \f[
182 /// Ker(x) = \sqrt{\frac{\pi}{2x}} e^{-x/\sqrt{2}} [F2(x) cos\beta + G2(x) sin\beta]
183 /// \f]
184 /// where \f$\beta = \frac{x}{\sqrt{2}} + \frac{\pi}{8}\f$.
185 ///
186 /// See also F2() and G2().
187 ///
188 /// Begin_Macro
189 /// {
190 /// TCanvas *c = new TCanvas("c","c",0,0,500,300);
191 /// TF1 *fKer = new TF1("fKer","ROOT::Math::KelvinFunctions::Ker(x)",-10,10);
192 /// fKer->Draw();
193 /// return c;
194 /// }
195 /// End_Macro
196 
197 double KelvinFunctions::Ker(double x)
198 {
199  if (fabs(x) < fgEpsilon) return 1E+100;
200 
201  if (fabs(x) < fgMin) {
202  double term = 1, x_factor = x * x * x * x * 0.0625;
203  double factorial = 1, harmonic = 0, n = 1, sum;
204  double delta = 0;
205  if(x < 0) delta = kPi;
206 
207  sum = - (log(fabs(x) * 0.5) + kEulerGamma) * Ber(x) + (kPi * 0.25 - delta) * Bei(x);
208 
209  do {
210  factorial = 4 * n * n * (2 * n - 1) * (2 * n - 1);
211  term *= (-1) / factorial * x_factor;
212  harmonic += 1 / (2 * n - 1 ) + 1 / (2 * n);
213  sum += term * harmonic;
214  n += 1;
215  if (n > 1000) break;
216  } while (fabs(term * harmonic) > fgEpsilon * sum);
217 
218  return sum;
219  } else {
220  double beta = x / kSqrt2 + kPi / 8;
221  double value = F2(x) * cos(beta) - G2(x) * sin(beta);
222  value *= sqrt(kPi / (2 * x)) * exp(- x / kSqrt2);
223  return value;
224  }
225 }
226 
227 
228 
229 ////////////////////////////////////////////////////////////////////////////////
230 /// \f[
231 /// Kei(x) = Kei_{0}(x) = Im\left[K_{0}\left(x e^{3\pi i/4}\right)\right]
232 /// \f]
233 /// where x is real, and \f$K_{0}(z)\f$ is the zeroth-order modified
234 /// Bessel function of the second kind.
235 ///
236 /// If x < fgMin (=20), Kei(x) is computed according to its polynomial
237 /// approximation
238 /// \f[
239 /// Kei(x) = -\left(ln \frac{x}{2} + \gamma\right) Bei(x) - \left(\frac{\pi}{4} - \delta\right) Ber(x) + \sum_{n \geq 0} \frac{(-1)^{n}}{[(2n)!]^{2}} H_{2n} \left(\frac{x}{2}\right)^{4n+2}
240 /// \f]
241 /// where \f$\gamma = 0.577215664...\f$ is the Euler-Mascheroni constant,
242 /// \f$\delta = \pi\f$ for x < 0 and is otherwise zero, and
243 /// \f[
244 /// H_{n} = \sum_{k = 1}^{n} \frac{1}{k}
245 /// \f]
246 /// For x > fgMin, Kei(x) is computed according to its asymptotic
247 /// expansion:
248 /// \f[
249 /// Kei(x) = - \sqrt{\frac{\pi}{2x}} e^{-x/\sqrt{2}} [F2(x) sin\beta + G2(x) cos\beta]
250 /// \f]
251 /// where \f$\beta = \frac{x}{\sqrt{2}} + \frac{\pi}{8}\f$.
252 ///
253 /// See also F2() and G2().
254 ///
255 /// Begin_Macro
256 /// {
257 /// TCanvas *c = new TCanvas("c","c",0,0,500,300);
258 /// TF1 *fKei = new TF1("fKei","ROOT::Math::KelvinFunctions::Kei(x)",-10,10);
259 /// fKei->Draw();
260 /// return c;
261 /// }
262 /// End_Macro
263 
264 double KelvinFunctions::Kei(double x)
265 {
266  if (fabs(x) < fgEpsilon) return (-0.25 * kPi);
267 
268  if (fabs(x) < fgMin) {
269  double term = x * x * 0.25, x_factor = term * term;
270  double factorial = 1, harmonic = 1, n = 1, sum;
271  double delta = 0;
272  if(x < 0) delta = kPi;
273 
274  sum = term - (log(fabs(x) * 0.5) + kEulerGamma) * Bei(x) - (kPi * 0.25 - delta) * Ber(x);
275 
276  do {
277  factorial = 4 * n * n * (2 * n + 1) * (2 * n + 1);
278  term *= (-1) / factorial * x_factor;
279  harmonic += 1 / (2 * n) + 1 / (2 * n + 1);
280  sum += term * harmonic;
281  n += 1;
282  if (n > 1000) break;
283  } while (fabs(term * harmonic) > fgEpsilon * sum);
284 
285  return sum;
286  } else {
287  double beta = x / kSqrt2 + kPi / 8;
288  double value = - F2(x) * sin(beta) - G2(x) * cos(beta);
289  value *= sqrt(kPi / (2 * x)) * exp(- x / kSqrt2);
290  return value;
291  }
292 }
293 
294 
295 
296 ////////////////////////////////////////////////////////////////////////////////
297 /// Calculates the first derivative of Ber(x).
298 ///
299 /// If x < fgMin (=20), DBer(x) is computed according to the derivative of
300 /// the polynomial approximation of Ber(x). Otherwise it is computed
301 /// according to its asymptotic expansion
302 /// \f[
303 /// \frac{d}{dx} Ber(x) = M cos\left(\theta - \frac{\pi}{4}\right)
304 /// \f]
305 /// See also M() and Theta().
306 ///
307 /// Begin_Macro
308 /// {
309 /// TCanvas *c = new TCanvas("c","c",0,0,500,300);
310 /// TF1 *fDBer = new TF1("fDBer","ROOT::Math::KelvinFunctions::DBer(x)",-10,10);
311 /// fDBer->Draw();
312 /// return c;
313 /// }
314 /// End_Macro
315 
316 double KelvinFunctions::DBer(double x)
317 {
318  if (fabs(x) < fgEpsilon) return 0;
319 
320  if (fabs(x) < fgMin) {
321  double sum, factorial = 1, n = 1;
322  double term = - x * x * x * 0.0625, x_factor = - term * x;
323 
324  sum = term;
325 
326  do {
327  factorial = 4 * n * (n + 1) * (2 * n + 1) * (2 * n + 1);
328  term *= (-1) / factorial * x_factor;
329  sum += term;
330  n += 1;
331  if (n > 1000) break;
332  } while (fabs(term) > fgEpsilon * sum);
333 
334  return sum;
335  }
336  else return (M(x) * sin(Theta(x) - kPi / 4));
337 }
338 
339 
340 
341 ////////////////////////////////////////////////////////////////////////////////
342 /// Calculates the first derivative of Bei(x).
343 ///
344 /// If x < fgMin (=20), DBei(x) is computed according to the derivative of
345 /// the polynomial approximation of Bei(x). Otherwise it is computed
346 /// according to its asymptotic expansion
347 /// \f[
348 /// \frac{d}{dx} Bei(x) = M sin\left(\theta - \frac{\pi}{4}\right)
349 /// \f]
350 /// See also M() and Theta().
351 ///
352 /// Begin_Macro
353 /// {
354 /// TCanvas *c = new TCanvas("c","c",0,0,500,300);
355 /// TF1 *fDBei = new TF1("fDBei","ROOT::Math::KelvinFunctions::DBei(x)",-10,10);
356 /// fDBei->Draw();
357 /// return c;
358 /// }
359 /// End_Macro
360 
361 double KelvinFunctions::DBei(double x)
362 {
363  if (fabs(x) < fgEpsilon) return 0;
364 
365  if (fabs(x) < fgMin) {
366  double sum, factorial = 1, n = 1;
367  double term = x * 0.5, x_factor = x * x * x * x * 0.0625;
368 
369  sum = term;
370 
371  do {
372  factorial = 4 * n * n * (2 * n - 1) * (2 * n + 1);
373  term *= (-1) * x_factor / factorial;
374  sum += term;
375  n += 1;
376  if (n > 1000) break;
377  } while (fabs(term) > fgEpsilon * sum);
378 
379  return sum;
380  }
381  else return (M(x) * cos(Theta(x) - kPi / 4));
382 }
383 
384 
385 
386 ////////////////////////////////////////////////////////////////////////////////
387 /// Calculates the first derivative of Ker(x).
388 ///
389 /// If x < fgMin (=20), DKer(x) is computed according to the derivative of
390 /// the polynomial approximation of Ker(x). Otherwise it is computed
391 /// according to its asymptotic expansion
392 /// \f[
393 /// \frac{d}{dx} Ker(x) = N cos\left(\phi - \frac{\pi}{4}\right)
394 /// \f]
395 /// See also N() and Phi().
396 ///
397 /// Begin_Macro
398 /// {
399 /// TCanvas *c = new TCanvas("c","c",0,0,500,300);
400 /// TF1 *fDKer = new TF1("fDKer","ROOT::Math::KelvinFunctions::DKer(x)",-10,10);
401 /// fDKer->Draw();
402 /// return c;
403 /// }
404 /// End_Macro
405 
406 double KelvinFunctions::DKer(double x)
407 {
408  if (fabs(x) < fgEpsilon) return -1E+100;
409 
410  if (fabs(x) < fgMin) {
411  double term = - x * x * x * 0.0625, x_factor = - term * x;
412  double factorial = 1, harmonic = 1.5, n = 1, sum;
413  double delta = 0;
414  if(x < 0) delta = kPi;
415 
416  sum = 1.5 * term - Ber(x) / x - (log(fabs(x) * 0.5) + kEulerGamma) * DBer(x) + (0.25 * kPi - delta) * DBei(x);
417 
418  do {
419  factorial = 4 * n * (n + 1) * (2 * n + 1) * (2 * n + 1);
420  term *= (-1) / factorial * x_factor;
421  harmonic += 1 / (2 * n + 1 ) + 1 / (2 * n + 2);
422  sum += term * harmonic;
423  n += 1;
424  if (n > 1000) break;
425  } while (fabs(term * harmonic) > fgEpsilon * sum);
426 
427  return sum;
428  }
429  else return N(x) * sin(Phi(x) - kPi / 4);
430 }
431 
432 
433 
434 ////////////////////////////////////////////////////////////////////////////////
435 /// Calculates the first derivative of Kei(x).
436 ///
437 /// If x < fgMin (=20), DKei(x) is computed according to the derivative of
438 /// the polynomial approximation of Kei(x). Otherwise it is computed
439 /// according to its asymptotic expansion
440 /// \f[
441 /// \frac{d}{dx} Kei(x) = N sin\left(\phi - \frac{\pi}{4}\right)
442 /// \f]
443 /// See also N() and Phi().
444 ///
445 /// Begin_Macro
446 /// {
447 /// TCanvas *c = new TCanvas("c","c",0,0,500,300);
448 /// TF1 *fDKei = new TF1("fDKei","ROOT::Math::KelvinFunctions::DKei(x)",-10,10);
449 /// fDKei->Draw();
450 /// return c;
451 /// }
452 /// End_Macro
453 
454 double KelvinFunctions::DKei(double x)
455 {
456  if (fabs(x) < fgEpsilon) return 0;
457 
458  if (fabs(x) < fgMin) {
459  double term = 0.5 * x, x_factor = x * x * x * x * 0.0625;
460  double factorial = 1, harmonic = 1, n = 1, sum;
461  double delta = 0;
462  if(x < 0) delta = kPi;
463 
464  sum = term - Bei(x) / x - (log(fabs(x) * 0.5) + kEulerGamma) * DBei(x) - (kPi * 0.25 - delta) * DBer(x);
465 
466  do {
467  factorial = 4 * n * n * (2 * n - 1) * (2 * n + 1);
468  term *= (-1) / factorial * x_factor;
469  harmonic += 1 / (2 * n) + 1 / (2 * n + 1);
470  sum += term * harmonic;
471  n += 1;
472  if (n > 1000) break;
473  } while (fabs(term * harmonic) > fgEpsilon * sum);
474 
475  return sum;
476  }
477  else return N(x) * cos(Phi(x) - kPi / 4);
478 }
479 
480 
481 
482 ////////////////////////////////////////////////////////////////////////////////
483 /// Utility function appearing in the calculations of the Kelvin
484 /// functions Bei(x) and Ber(x) (and their derivatives). F1(x) is given by
485 /// \f[
486 /// F1(x) = 1 + \sum_{n \geq 1} \frac{\prod_{m=1}^{n}(2m - 1)^{2}}{n! (8x)^{n}} cos\left(\frac{n\pi}{4}\right)
487 /// \f]
488 
489 double KelvinFunctions::F1(double x)
490 {
491  double sum, term;
492  double prod = 1, x_factor = 8 * x, factorial = 1, n = 2;
493 
494  sum = kSqrt2 / (16 * x);
495 
496  do {
497  factorial *= n;
498  prod *= (2 * n - 1) * (2 * n - 1);
499  x_factor *= 8 * x;
500  term = prod / (factorial * x_factor) * cos(0.25 * n * kPi);
501  sum += term;
502  n += 1;
503  if (n > 1000) break;
504  } while (fabs(term) > fgEpsilon * sum);
505 
506  sum += 1;
507 
508  return sum;
509 }
510 
511 ////////////////////////////////////////////////////////////////////////////////
512 /// Utility function appearing in the calculations of the Kelvin
513 /// functions Kei(x) and Ker(x) (and their derivatives). F2(x) is given by
514 /// \f[
515 /// F2(x) = 1 + \sum_{n \geq 1} (-1)^{n} \frac{\prod_{m=1}^{n}(2m - 1)^{2}}{n! (8x)^{n}} cos\left(\frac{n\pi}{4}\right)
516 /// \f]
517 
518 double KelvinFunctions::F2(double x)
519 {
520  double sum, term;
521  double prod = 1, x_factor = 8 * x, factorial = 1, n = 2;
522 
523  sum = kSqrt2 / (16 * x);
524 
525  do {
526  factorial *= - n;
527  prod *= (2 * n - 1) * (2 * n - 1);
528  x_factor *= 8 * x;
529  term = (prod / (factorial * x_factor)) * cos(0.25 * n * kPi);
530  sum += term;
531  n += 1;
532  if (n > 1000) break;
533  } while (fabs(term) > fgEpsilon * sum);
534 
535  sum += 1;
536 
537  return sum;
538 }
539 
540 
541 
542 ////////////////////////////////////////////////////////////////////////////////
543 /// Utility function appearing in the calculations of the Kelvin
544 /// functions Bei(x) and Ber(x) (and their derivatives). G1(x) is given by
545 /// \f[
546 /// G1(x) = \sum_{n \geq 1} \frac{\prod_{m=1}^{n}(2m - 1)^{2}}{n! (8x)^{n}} sin\left(\frac{n\pi}{4}\right)
547 /// \f]
548 
549 double KelvinFunctions::G1(double x)
550 {
551  double sum, term;
552  double prod = 1, x_factor = 8 * x, factorial = 1, n = 2;
553 
554  sum = kSqrt2 / (16 * x);
555 
556  do {
557  factorial *= n;
558  prod *= (2 * n - 1) * (2 * n - 1);
559  x_factor *= 8 * x;
560  term = prod / (factorial * x_factor) * sin(0.25 * n * kPi);
561  sum += term;
562  n += 1;
563  if (n > 1000) break;
564  } while (fabs(term) > fgEpsilon * sum);
565 
566  return sum;
567 }
568 
569 ////////////////////////////////////////////////////////////////////////////////
570 /// Utility function appearing in the calculations of the Kelvin
571 /// functions Kei(x) and Ker(x) (and their derivatives). G2(x) is given by
572 /// \f[
573 /// G2(x) = \sum_{n \geq 1} (-1)^{n} \frac{\prod_{m=1}^{n}(2m - 1)^{2}}{n! (8x)^{n}} sin\left(\frac{n\pi}{4}\right)
574 /// \f]
575 
576 double KelvinFunctions::G2(double x)
577 {
578  double sum, term;
579  double prod = 1, x_factor = 8 * x, factorial = 1, n = 2;
580 
581  sum = kSqrt2 / (16 * x);
582 
583  do {
584  factorial *= - n;
585  prod *= (2 * n - 1) * (2 * n - 1);
586  x_factor *= 8 * x;
587  term = prod / (factorial * x_factor) * sin(0.25 * n * kPi);
588  sum += term;
589  n += 1;
590  if (n > 1000) break;
591  } while (fabs(term) > fgEpsilon * sum);
592 
593  return sum;
594 }
595 
596 
597 
598 ////////////////////////////////////////////////////////////////////////////////
599 /// Utility function appearing in the asymptotic expansions of DBer(x) and
600 /// DBei(x). M(x) is given by
601 /// \f[
602 /// M(x) = \frac{e^{x/\sqrt{2}}}{\sqrt{2\pi x}}\left(1 + \frac{1}{8\sqrt{2} x} + \frac{1}{256 x^{2}} - \frac{399}{6144\sqrt{2} x^{3}} + O\left(\frac{1}{x^{4}}\right)\right)
603 /// \f]
604 
605 double KelvinFunctions::M(double x)
606 {
607  double value = 1 + 1 / (8 * kSqrt2 * x) + 1 / (256 * x * x) - 399 / (6144 * kSqrt2 * x * x * x);
608  value *= exp(x / kSqrt2) / sqrt(2 * kPi * x);
609  return value;
610 }
611 
612 
613 
614 ////////////////////////////////////////////////////////////////////////////////
615 /// Utility function appearing in the asymptotic expansions of DBer(x) and
616 /// DBei(x). \f$\theta(x)\f$ is given by
617 /// \f[
618 /// \theta(x) = \frac{x}{\sqrt{2}} - \frac{\pi}{8} - \frac{1}{8\sqrt{2} x} - \frac{1}{16 x^{2}} - \frac{25}{384\sqrt{2} x^{3}} + O\left(\frac{1}{x^{5}}\right)
619 /// \f]
620 
621 double KelvinFunctions::Theta(double x)
622 {
623  double value = x / kSqrt2 - kPi / 8;
624  value -= 1 / (8 * kSqrt2 * x) + 1 / (16 * x * x) + 25 / (384 * kSqrt2 * x * x * x);
625  return value;
626 }
627 
628 
629 
630 ////////////////////////////////////////////////////////////////////////////////
631 /// Utility function appearing in the asymptotic expansions of DKer(x) and
632 /// DKei(x). N(x) is given by
633 /// \f[
634 /// N(x) = \sqrt{\frac{\pi}{2x}} e^{-x/\sqrt{2}} \left(1 - \frac{1}{8\sqrt{2} x} + \frac{1}{256 x^{2}} + \frac{399}{6144\sqrt{2} x^{3}} + O\left(\frac{1}{x^{4}}\right)\right)
635 /// \f]
636 
637 double KelvinFunctions::N(double x)
638 {
639  double value = 1 - 1 / (8 * kSqrt2 * x) + 1 / (256 * x * x) + 399 / (6144 * kSqrt2 * x * x * x);
640  value *= exp(- x / kSqrt2) * sqrt(kPi / (2 * x));
641  return value;
642 }
643 
644 
645 
646 ////////////////////////////////////////////////////////////////////////////////
647 /// Utility function appearing in the asymptotic expansions of DKer(x) and
648 /// DKei(x). \f$\phi(x)\f$ is given by
649 /// \f[
650 /// \phi(x) = - \frac{x}{\sqrt{2}} - \frac{\pi}{8} + \frac{1}{8\sqrt{2} x} - \frac{1}{16 x^{2}} + \frac{25}{384\sqrt{2} x^{3}} + O\left(\frac{1}{x^{5}}\right)
651 /// \f]
652 
653 double KelvinFunctions::Phi(double x)
654 {
655  double value = - x / kSqrt2 - kPi / 8;
656  value += 1 / (8 * kSqrt2 * x) - 1 / (16 * x * x) + 25 / (384 * kSqrt2 * x * x * x);
657  return value;
658 }
659 
660 
661 } // namespace Math
662 } // namespace ROOT
663 
664 
665 
static double Ber(double x)
where x is real, and is the zeroth-order Bessel function of the first kind.
Namespace for new ROOT classes and functions.
Definition: ROOT.py:1
static double DBei(double x)
Calculates the first derivative of Bei(x).
static const double kSqrt2
double kEulerGamma
double cos(double)
VecExpr< UnaryOp< Sqrt< T >, VecExpr< A, T, D >, T >, T, D > sqrt(const VecExpr< A, T, D > &rhs)
double beta(double x, double y)
Calculates the beta function.
static double F1(double x)
Utility function appearing in the calculations of the Kelvin functions Bei(x) and Ber(x) (and their d...
Double_t x[n]
Definition: legend1.C:17
static double G1(double x)
Utility function appearing in the calculations of the Kelvin functions Bei(x) and Ber(x) (and their d...
double sin(double)
VecExpr< UnaryOp< Fabs< T >, VecExpr< A, T, D >, T >, T, D > fabs(const VecExpr< A, T, D > &rhs)
static double Kei(double x)
where x is real, and is the zeroth-order modified Bessel function of the second kind...
static double DKer(double x)
Calculates the first derivative of Ker(x).
Double_t E()
Definition: TMath.h:54
static double M(double x)
Utility function appearing in the asymptotic expansions of DBer(x) and DBei(x).
static double Bei(double x)
where x is real, and is the zeroth-order Bessel function of the first kind.
static double Theta(double x)
Utility function appearing in the asymptotic expansions of DBer(x) and DBei(x).
Namespace for new Math classes and functions.
static double G2(double x)
Utility function appearing in the calculations of the Kelvin functions Kei(x) and Ker(x) (and their d...
static double DBer(double x)
Calculates the first derivative of Ber(x).
static double Phi(double x)
Utility function appearing in the asymptotic expansions of DKer(x) and DKei(x).
static double Ker(double x)
where x is real, and is the zeroth-order modified Bessel function of the second kind...
static double DKei(double x)
Calculates the first derivative of Kei(x).
double exp(double)
static double F2(double x)
Utility function appearing in the calculations of the Kelvin functions Kei(x) and Ker(x) (and their d...
float value
Definition: math.cpp:443
const Int_t n
Definition: legend1.C:16
double log(double)
static double N(double x)
Utility function appearing in the asymptotic expansions of DKer(x) and DKei(x).