KelvinFunctions.cxx

Go to the documentation of this file.
00001 // @(#)root/mathmore:$Id: KelvinFunctions.cxx 20882 2007-11-19 11:31:26Z rdm $
00002 // The functions in this class have been imported by Jason Detwiler (jasondet@gmail.com) from
00003 // CodeCogs GNU General Public License Agreement
00004 // Copyright (C) 2004-2005 CodeCogs, Zyba Ltd, Broadwood, Holford, TA5 1DU,
00005 // England.
00006 //
00007 // This program is free software; you can redistribute it and/or modify it
00008 // under
00009 // the terms of the GNU General Public License as published by CodeCogs. 
00010 // You must retain a copy of this licence in all copies. 
00011 //
00012 // This program is distributed in the hope that it will be useful, but
00013 // WITHOUT ANY
00014 // WARRANTY; without even the implied warranty of MERCHANTABILITY or
00015 // FITNESS FOR A
00016 // PARTICULAR PURPOSE. See the Adapted GNU General Public License for more
00017 // details.
00018 //
00019 // *** THIS SOFTWARE CAN NOT BE USED FOR COMMERCIAL GAIN. ***
00020 // ---------------------------------------------------------------------------------
00021 
00022 #include "Math/KelvinFunctions.h"
00023 #include <math.h>
00024 
00025 
00026 namespace ROOT {
00027 namespace Math {
00028 
00029 double KelvinFunctions::fgMin = 20;
00030 double KelvinFunctions::fgEpsilon = 1.e-20;
00031 
00032 double kSqrt2 = 1.4142135623730950488016887242097;
00033 double kPi    = 3.14159265358979323846;
00034 double kEulerGamma = 0.577215664901532860606512090082402431042;
00035 
00036 
00037 /* Begin_Html
00038 <center><h2>KelvinFunctions</h2></center>
00039 
00040 <p>
00041 This class calculates the Kelvin functions Ber(x), Bei(x), Ker(x),
00042 Kei(x), and their first derivatives. 
00043 </p>
00044 
00045 End_Html */
00046 
00047 //______________________________________________________________________________
00048 double KelvinFunctions::Ber(double x) 
00049 {
00050    // Begin_Latex
00051    // Ber(x) = Ber_{0}(x) = Re#left[J_{0}#left(x e^{3#pii/4}#right)#right]
00052    // End_Latex
00053    // where x is real, and Begin_Latex J_{0}(z) End_Latex is the zeroth-order Bessel
00054    // function of the first kind.
00055    // 
00056    // If x < fgMin (=20), Ber(x) is computed according to its polynomial
00057    // approximation
00058    // Begin_Latex
00059    // Ber(x) = 1 + #sum_{n #geq 1}#frac{(-1)^{n}(x/2)^{4n}}{[(2n)!]^{2}}
00060    // End_Latex
00061    // For x > fgMin, Ber(x) is computed according to its asymptotic
00062    // expansion:
00063    // Begin_Latex
00064    // Ber(x) = #frac{e^{x/#sqrt{2}}}{#sqrt{2#pix}} [F1(x) cos#alpha + G1(x) sin#alpha] - #frac{1}{#pi}Kei(x)
00065    // End_Latex
00066    // where Begin_Latex #alpha = #frac{x}{#sqrt{2}} - #frac{#pi}{8} End_Latex. 
00067    // See also F1(x) and G1(x).
00068    //
00069    // Begin_Macro
00070    // {
00071    //   TCanvas *c = new TCanvas("c","c",0,0,500,300);
00072    //   TF1 *fBer = new TF1("fBer","ROOT::Math::KelvinFunctions::Ber(x)",-10,10);
00073    //   fBer->Draw();
00074    //   return c;
00075    // }
00076    // End_Macro
00077 
00078    if (fabs(x) < fgEpsilon) return 1;
00079  
00080    if (fabs(x) < fgMin) {
00081       double sum, factorial = 1, n = 1;
00082       double term = 1, x_factor = x * x * x * x * 0.0625;
00083  
00084       sum = 1;
00085  
00086       do {
00087          factorial = 4 * n * n * (2 * n - 1) * (2 * n - 1);
00088          term *= (-1) / factorial * x_factor;
00089          sum += term;
00090          n += 1;
00091          if (n > 1000) break;
00092       } while (fabs(term) > fgEpsilon * sum);
00093  
00094       return sum;
00095    } else {
00096       double alpha = x / kSqrt2 - kPi / 8;
00097       double value = F1(x) * cos(alpha) + G1(x) * sin(alpha);
00098       value *= exp(x / kSqrt2) / sqrt(2 * kPi * x);
00099       value -= Kei(x) / kPi;
00100       return value;
00101    }
00102 }
00103  
00104 //______________________________________________________________________________
00105 double KelvinFunctions::Bei(double x) 
00106 {
00107    // Begin_Latex
00108    // Bei(x) = Bei_{0}(x) = Im#left[J_{0}#left(x e^{3#pii/4}#right)#right]
00109    // End_Latex
00110    // where x is real, and Begin_Latex J_{0}(z) End_Latex is the zeroth-order Bessel
00111    // function of the first kind.
00112    // 
00113    // If x < fgMin (=20), Bei(x) is computed according to its polynomial
00114    // approximation
00115    // Begin_Latex
00116    // Bei(x) = #sum_{n #geq 0}#frac{(-1)^{n}(x/2)^{4n+2}}{[(2n+1)!]^{2}}
00117    // End_Latex
00118    // For x > fgMin, Bei(x) is computed according to its asymptotic
00119    // expansion:
00120    // Begin_Latex
00121    // Bei(x) = #frac{e^{x/#sqrt{2}}}{#sqrt{2#pix}} [F1(x) sin#alpha + G1(x) cos#alpha] - #frac{1}{#pi}Ker(x)
00122    // End_Latex
00123    // where Begin_Latex #alpha = #frac{x}{#sqrt{2}} - #frac{#pi}{8} End_Latex
00124    // See also F1(x) and G1(x).
00125    //
00126    // Begin_Macro
00127    // {
00128    //   TCanvas *c = new TCanvas("c","c",0,0,500,300);
00129    //   TF1 *fBei = new TF1("fBei","ROOT::Math::KelvinFunctions::Bei(x)",-10,10);
00130    //   fBei->Draw();
00131    //   return c;
00132    // }
00133    // End_Macro
00134    
00135 
00136    if (fabs(x) < fgEpsilon) return 0;
00137  
00138    if (fabs(x) < fgMin) {
00139       double sum, factorial = 1, n = 1;
00140       double term = x * x * 0.25, x_factor = term * term;
00141  
00142       sum = term;
00143  
00144       do {
00145          factorial = 4 * n * n * (2 * n + 1) * (2 * n + 1);
00146          term *= (-1) / factorial * x_factor;
00147          sum += term;
00148          n += 1;
00149          if (n > 1000) break;
00150       } while (fabs(term) > fgEpsilon * sum);
00151  
00152       return sum;
00153    } else {
00154       double alpha = x / kSqrt2 - kPi / 8;
00155       double value = F1(x) * sin(alpha) + G1(x) * cos(alpha);
00156       value *= exp(x / kSqrt2) / sqrt(2 * kPi * x);
00157       value += Ker(x) / kPi;
00158       return value;
00159    }
00160 }
00161  
00162  
00163 //______________________________________________________________________________
00164 double KelvinFunctions::Ker(double x) 
00165 {
00166    // Begin_Latex
00167    // Ker(x) = Ker_{0}(x) = Re#left[K_{0}#left(x e^{3#pii/4}#right)#right]
00168    // End_Latex
00169    // where x is real, and Begin_Latex K_{0}(z) End_Latex is the zeroth-order modified
00170    // Bessel function of the second kind.
00171    //
00172    // If x < fgMin (=20), Ker(x) is computed according to its polynomial
00173    // approximation
00174    // Begin_Latex
00175    // 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}
00176    // End_Latex
00177    // where Begin_Latex #gamma = 0.577215664... End_Latex is the Euler-Mascheroni constant, 
00178    // Begin_Latex #delta = #pi End_Latex for x < 0 and is otherwise zero, and 
00179    // Begin_Latex 
00180    // H_{n} = #sum_{k = 1}^{n} #frac{1}{k} 
00181    // End_Latex
00182    // For x > fgMin, Ker(x) is computed according to its asymptotic
00183    // expansion:
00184    // Begin_Latex
00185    // Ker(x) = #sqrt{#frac{#pi}{2x}} e^{-x/#sqrt{2}} [F2(x) cos#beta + G2(x) sin#beta]
00186    // End_Latex
00187    // where Begin_Latex #beta = #frac{x}{#sqrt{2}} + #frac{#pi}{8} End_Latex
00188    // See also F2(x) and G2(x).
00189    //
00190    // Begin_Macro
00191    // {
00192    //   TCanvas *c = new TCanvas("c","c",0,0,500,300);
00193    //   TF1 *fKer = new TF1("fKer","ROOT::Math::KelvinFunctions::Ker(x)",-10,10);
00194    //   fKer->Draw();
00195    //   return c;
00196    // }
00197    // End_Macro
00198 
00199    if (fabs(x) < fgEpsilon) return 1E+100;
00200  
00201    if (fabs(x) < fgMin) {
00202       double term = 1, x_factor = x * x * x * x * 0.0625;
00203       double factorial = 1, harmonic = 0, n = 1, sum;
00204       double delta = 0;
00205       if(x < 0) delta = kPi;
00206 
00207       sum  = - (log(fabs(x) * 0.5) + kEulerGamma) * Ber(x) + (kPi * 0.25 - delta) * Bei(x);
00208  
00209       do {
00210          factorial = 4 * n * n * (2 * n - 1) * (2 * n - 1);
00211          term *= (-1) / factorial * x_factor;
00212          harmonic += 1 / (2 * n - 1 ) + 1 / (2 * n);
00213          sum += term * harmonic;
00214          n += 1;
00215          if (n > 1000) break;
00216       } while (fabs(term * harmonic) > fgEpsilon * sum);
00217  
00218       return sum;
00219    } else {
00220       double beta = x / kSqrt2 + kPi / 8;
00221       double value = F2(x) * cos(beta) - G2(x) * sin(beta);
00222       value *= sqrt(kPi / (2 * x)) * exp(- x / kSqrt2);
00223       return value;
00224    }
00225 }
00226  
00227  
00228  
00229 //______________________________________________________________________________
00230 double KelvinFunctions::Kei(double x) 
00231 {
00232    // Begin_Latex
00233    // Kei(x) = Kei_{0}(x) = Im#left[K_{0}#left(x e^{3#pii/4}#right)#right]
00234    // End_Latex
00235    // where x is real, and Begin_Latex K_{0}(z) End_Latex is the zeroth-order modified
00236    // Bessel function of the second kind.
00237    //
00238    // If x < fgMin (=20), Kei(x) is computed according to its polynomial
00239    // approximation
00240    // Begin_Latex
00241    // 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}
00242    // End_Latex
00243    // where Begin_Latex #gamma = 0.577215664... End_Latex is the Euler-Mascheroni constant,
00244    // Begin_Latex #delta = #pi End_Latex for x < 0 and is otherwise zero, and 
00245    // Begin_Latex 
00246    // H_{n} = #sum_{k = 1}^{n} #frac{1}{k} 
00247    // End_Latex
00248    // For x > fgMin, Kei(x) is computed according to its asymptotic
00249    // expansion:
00250    // Begin_Latex
00251    // Kei(x) = - #sqrt{#frac{#pi}{2x}} e^{-x/#sqrt{2}} [F2(x) sin#beta + G2(x) cos#beta]
00252    // End_Latex
00253    // where Begin_Latex #beta = #frac{x}{#sqrt{2}} + #frac{#pi}{8} End_Latex
00254    // See also F2(x) and G2(x).
00255    //
00256    // Begin_Macro
00257    // {
00258    //   TCanvas *c = new TCanvas("c","c",0,0,500,300);
00259    //   TF1 *fKei = new TF1("fKei","ROOT::Math::KelvinFunctions::Kei(x)",-10,10);
00260    //   fKei->Draw();
00261    //   return c;
00262    // }
00263    // End_Macro
00264 
00265    if (fabs(x) < fgEpsilon) return (-0.25 * kPi);
00266  
00267    if (fabs(x) < fgMin) {
00268       double term = x * x * 0.25, x_factor = term * term;
00269       double factorial = 1, harmonic = 1, n = 1, sum;
00270       double delta = 0;
00271       if(x < 0) delta = kPi;
00272 
00273       sum  = term - (log(fabs(x) * 0.5) + kEulerGamma) * Bei(x) - (kPi * 0.25 - delta) * Ber(x);
00274  
00275       do {
00276          factorial = 4 * n * n * (2 * n + 1) * (2 * n + 1);
00277          term *= (-1) / factorial * x_factor;
00278          harmonic += 1 / (2 * n) + 1 / (2 * n + 1);
00279          sum += term * harmonic;
00280          n += 1;
00281          if (n > 1000) break;
00282       } while (fabs(term * harmonic) > fgEpsilon * sum);
00283  
00284       return sum;
00285    } else {
00286       double beta = x / kSqrt2 + kPi / 8;
00287       double value = - F2(x) * sin(beta) - G2(x) * cos(beta);
00288       value *= sqrt(kPi / (2 * x)) * exp(- x / kSqrt2);
00289       return value;
00290    }
00291 }
00292  
00293  
00294  
00295 //______________________________________________________________________________
00296 double KelvinFunctions::DBer(double x) 
00297 {
00298    // Calculates the first derivative of Ber(x).
00299    //
00300    // If x < fgMin (=20), DBer(x) is computed according to the derivative of
00301    // the polynomial approximation of Ber(x). Otherwise it is computed
00302    // according to its asymptotic expansion
00303    // Begin_Latex
00304    // #frac{d}{dx} Ber(x) = M cos#left(#theta - #frac{#pi}{4}#right)
00305    // End_Latex
00306    // See also M(x) and Theta(x).
00307    //
00308    // Begin_Macro
00309    // {
00310    //   TCanvas *c = new TCanvas("c","c",0,0,500,300);
00311    //   TF1 *fDBer = new TF1("fDBer","ROOT::Math::KelvinFunctions::DBer(x)",-10,10);
00312    //   fDBer->Draw();
00313    //   return c;
00314    // }
00315    // End_Macro
00316    if (fabs(x) < fgEpsilon) return 0;
00317  
00318    if (fabs(x) < fgMin) {
00319       double sum, factorial = 1, n = 1;
00320       double term = - x * x * x * 0.0625, x_factor = - term * x;
00321  
00322       sum = term;
00323  
00324       do {
00325          factorial = 4 * n * (n + 1) * (2 * n + 1) * (2 * n + 1);
00326          term *= (-1) / factorial * x_factor;
00327          sum += term;
00328          n += 1;
00329          if (n > 1000) break;
00330       } while (fabs(term) > fgEpsilon * sum);
00331  
00332       return sum;
00333    }
00334    else return (M(x) * sin(Theta(x) - kPi / 4));
00335 }
00336  
00337  
00338  
00339 //______________________________________________________________________________
00340 double KelvinFunctions::DBei(double x) 
00341 {
00342    // Calculates the first derivative of Bei(x).
00343    //
00344    // If x < fgMin (=20), DBei(x) is computed according to the derivative of
00345    // the polynomial approximation of Bei(x). Otherwise it is computed
00346    // according to its asymptotic expansion
00347    // Begin_Latex
00348    // #frac{d}{dx} Bei(x) = M sin#left(#theta - #frac{#pi}{4}#right)
00349    // End_Latex
00350    // See also M(x) and Theta(x).
00351    //
00352    // Begin_Macro
00353    // {
00354    //   TCanvas *c = new TCanvas("c","c",0,0,500,300);
00355    //   TF1 *fDBei = new TF1("fDBei","ROOT::Math::KelvinFunctions::DBei(x)",-10,10);
00356    //   fDBei->Draw();
00357    //   return c;
00358    // }
00359    // End_Macro
00360    if (fabs(x) < fgEpsilon) return 0;
00361   
00362    if (fabs(x) < fgMin) {
00363       double sum, factorial = 1, n = 1;
00364       double term = x * 0.5, x_factor = x * x * x * x * 0.0625;
00365  
00366       sum = term;
00367  
00368       do {
00369          factorial = 4 * n * n * (2 * n - 1) * (2 * n + 1);
00370          term *= (-1) * x_factor / factorial;
00371          sum += term;
00372          n += 1;
00373          if (n > 1000) break;
00374       } while (fabs(term) > fgEpsilon * sum);
00375  
00376       return sum;
00377    }
00378    else return (M(x) * cos(Theta(x) - kPi / 4));
00379 }
00380  
00381  
00382  
00383 //______________________________________________________________________________
00384 double KelvinFunctions::DKer(double x) 
00385 {
00386    // Calculates the first derivative of Ker(x).
00387    //
00388    // If x < fgMin (=20), DKer(x) is computed according to the derivative of
00389    // the polynomial approximation of Ker(x). Otherwise it is computed
00390    // according to its asymptotic expansion
00391    // Begin_Latex
00392    // #frac{d}{dx} Ker(x) = N cos#left(#phi - #frac{#pi}{4}#right)
00393    // End_Latex
00394    // See also N(x) and Phi(x).
00395    //
00396    // Begin_Macro
00397    // {
00398    //   TCanvas *c = new TCanvas("c","c",0,0,500,300);
00399    //   TF1 *fDKer = new TF1("fDKer","ROOT::Math::KelvinFunctions::DKer(x)",-10,10);
00400    //   fDKer->Draw();
00401    //   return c;
00402    // }
00403    // End_Macro
00404    if (fabs(x) < fgEpsilon) return -1E+100;
00405  
00406    if (fabs(x) < fgMin) {
00407       double term = - x * x * x * 0.0625, x_factor = - term * x;
00408       double factorial = 1, harmonic = 1.5, n = 1, sum;
00409       double delta = 0;
00410       if(x < 0) delta = kPi;
00411  
00412       sum = 1.5 * term - Ber(x) / x - (log(fabs(x) * 0.5) + kEulerGamma) * DBer(x) + (0.25 * kPi - delta) * DBei(x);
00413  
00414       do {
00415          factorial = 4 * n * (n + 1) * (2 * n + 1) * (2 * n + 1);
00416          term *= (-1) / factorial * x_factor;
00417          harmonic += 1 / (2 * n + 1 ) + 1 / (2 * n + 2);
00418          sum += term * harmonic;
00419          n += 1;
00420          if (n > 1000) break;
00421       } while (fabs(term * harmonic) > fgEpsilon * sum);
00422  
00423       return sum;
00424    }
00425    else return N(x) * sin(Phi(x) - kPi / 4);
00426 }
00427  
00428  
00429  
00430 //______________________________________________________________________________
00431 double KelvinFunctions::DKei(double x) 
00432 {
00433    // Calculates the first derivative of Kei(x).
00434    //
00435    // If x < fgMin (=20), DKei(x) is computed according to the derivative of
00436    // the polynomial approximation of Kei(x). Otherwise it is computed
00437    // according to its asymptotic expansion
00438    // Begin_Latex
00439    // #frac{d}{dx} Kei(x) = N sin#left(#phi - #frac{#pi}{4}#right)
00440    // End_Latex
00441    // See also N(x) and Phi(x).
00442    //
00443    // Begin_Macro
00444    // {
00445    //   TCanvas *c = new TCanvas("c","c",0,0,500,300);
00446    //   TF1 *fDKei = new TF1("fDKei","ROOT::Math::KelvinFunctions::DKei(x)",-10,10);
00447    //   fDKei->Draw();
00448    //   return c;
00449    // }
00450    // End_Macro
00451    if (fabs(x) < fgEpsilon) return 0;
00452  
00453    if (fabs(x) < fgMin) {
00454       double term = 0.5 * x, x_factor = x * x * x * x * 0.0625;
00455       double factorial = 1, harmonic = 1, n = 1, sum;
00456       double delta = 0;
00457       if(x < 0) delta = kPi;
00458  
00459       sum  = term - Bei(x) / x - (log(fabs(x) * 0.5) + kEulerGamma) * DBei(x) - (kPi * 0.25 - delta) * DBer(x);
00460  
00461       do {
00462          factorial = 4 * n * n * (2 * n - 1) * (2 * n + 1);
00463          term *= (-1) / factorial * x_factor;
00464          harmonic += 1 / (2 * n) + 1 / (2 * n + 1);
00465          sum += term * harmonic;
00466          n += 1;
00467          if (n > 1000) break;
00468       } while (fabs(term * harmonic) > fgEpsilon * sum);
00469  
00470       return sum;
00471    }
00472    else return N(x) * cos(Phi(x) - kPi / 4);
00473 }
00474  
00475  
00476  
00477 //______________________________________________________________________________
00478 double KelvinFunctions::F1(double x) 
00479 {
00480    // Utility function appearing in the calculations of the Kelvin
00481    // functions Bei(x) and Ber(x) (and their derivatives). F1(x) is given by
00482    // Begin_Latex
00483    // 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)
00484    // End_Latex
00485    double sum, term;
00486    double prod = 1, x_factor = 8 * x, factorial = 1, n = 2;
00487   
00488    sum = kSqrt2 / (16 * x);
00489  
00490    do {
00491       factorial *= n;
00492       prod *= (2 * n - 1) * (2 * n - 1);
00493       x_factor *= 8 * x;
00494       term = prod / (factorial * x_factor) * cos(0.25 * n * kPi);
00495       sum += term;
00496       n += 1;
00497       if (n > 1000) break;
00498    } while (fabs(term) > fgEpsilon * sum);
00499 
00500    sum += 1;
00501  
00502    return sum;
00503 }
00504  
00505 //______________________________________________________________________________
00506 double KelvinFunctions::F2(double x) 
00507 {
00508    // Utility function appearing in the calculations of the Kelvin
00509    // functions Kei(x) and Ker(x) (and their derivatives). F2(x) is given by
00510    // Begin_Latex
00511    // 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)
00512    // End_Latex
00513    double sum, term;
00514    double prod = 1, x_factor = 8 * x, factorial = 1, n = 2;
00515  
00516    sum = kSqrt2 / (16 * x);
00517  
00518    do {
00519       factorial *= - n;
00520       prod *= (2 * n - 1) * (2 * n - 1);
00521       x_factor *= 8 * x;
00522       term = (prod / (factorial * x_factor)) * cos(0.25 * n * kPi);
00523       sum += term;
00524       n += 1;
00525       if (n > 1000) break;
00526    } while (fabs(term) > fgEpsilon * sum);
00527 
00528    sum += 1;
00529  
00530    return sum;
00531 }
00532  
00533  
00534  
00535 //______________________________________________________________________________
00536 double KelvinFunctions::G1(double x) 
00537 {
00538    // Utility function appearing in the calculations of the Kelvin
00539    // functions Bei(x) and Ber(x) (and their derivatives). G1(x) is given by
00540    // Begin_Latex
00541    // G1(x) = #sum_{n #geq 1} #frac{#prod_{m=1}^{n}(2m - 1)^{2}}{n! (8x)^{n}} sin#left(#frac{n#pi}{4}#right)
00542    // End_Latex
00543    double sum, term;
00544    double prod = 1, x_factor = 8 * x, factorial = 1, n = 2;
00545  
00546    sum = kSqrt2 / (16 * x);
00547  
00548    do {
00549       factorial *= n;
00550       prod *= (2 * n - 1) * (2 * n - 1);
00551       x_factor *= 8 * x;
00552       term = prod / (factorial * x_factor) * sin(0.25 * n * kPi);
00553       sum += term;
00554       n += 1;
00555       if (n > 1000) break;
00556    } while (fabs(term) > fgEpsilon * sum);
00557  
00558    return sum;
00559 }
00560  
00561 //______________________________________________________________________________
00562 double KelvinFunctions::G2(double x) 
00563 {
00564    // Utility function appearing in the calculations of the Kelvin
00565    // functions Kei(x) and Ker(x) (and their derivatives). G2(x) is given by
00566    // Begin_Latex
00567    // 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)
00568    // End_Latex
00569    double sum, term;
00570    double prod = 1, x_factor = 8 * x, factorial = 1, n = 2;
00571  
00572    sum = kSqrt2 / (16 * x);
00573  
00574    do {
00575       factorial *= - n;
00576       prod *= (2 * n - 1) * (2 * n - 1);
00577       x_factor *= 8 * x;
00578       term = prod / (factorial * x_factor) * sin(0.25 * n * kPi);
00579       sum += term;
00580       n += 1;
00581       if (n > 1000) break;
00582    } while (fabs(term) > fgEpsilon * sum);
00583  
00584    return sum;
00585 }
00586  
00587  
00588  
00589 //______________________________________________________________________________
00590 double KelvinFunctions::M(double x) 
00591 {
00592    // Utility function appearing in the asymptotic expansions of DBer(x) and
00593    // DBei(x). M(x) is given by
00594    // Begin_Latex
00595    // M(x) = #frac{e^{x/#sqrt{2}}}{#sqrt{2#pix}}#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)
00596    // End_Latex
00597    double value = 1 + 1 / (8 * kSqrt2 * x) + 1 / (256 * x * x) - 399 / (6144 * kSqrt2 * x * x * x);
00598    value *= exp(x / kSqrt2) / sqrt(2 * kPi * x);
00599    return value;
00600 }
00601  
00602  
00603  
00604 //______________________________________________________________________________
00605 double KelvinFunctions::Theta(double x) 
00606 {
00607    // Utility function appearing in the asymptotic expansions of DBer(x) and
00608    // DBei(x). Begin_Latex #theta(x) #End_Latex is given by
00609    // Begin_Latex
00610    // #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)
00611    // End_Latex
00612    double value = x / kSqrt2 - kPi / 8;
00613    value -= 1 / (8 * kSqrt2 * x) + 1 / (16 * x * x) + 25 / (384 * kSqrt2 * x * x * x);
00614    return value;
00615 }
00616  
00617  
00618  
00619 //______________________________________________________________________________
00620 double KelvinFunctions::N(double x) 
00621 {
00622    // Utility function appearing in the asymptotic expansions of DKer(x) and
00623    // DKei(x). (x) is given by
00624    // Begin_Latex
00625    // 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)
00626    // End_Latex
00627    double value = 1 - 1 / (8 * kSqrt2 * x) + 1 / (256 * x * x) + 399 / (6144 * kSqrt2 * x * x * x);
00628    value *= exp(- x / kSqrt2) * sqrt(kPi / (2 * x));
00629    return value;
00630 }
00631  
00632  
00633  
00634 //______________________________________________________________________________
00635 double KelvinFunctions::Phi(double x) 
00636 {
00637    // Utility function appearing in the asymptotic expansions of DKer(x) and
00638    // DKei(x). Begin_Latex #phi(x) #End_Latex is given by
00639    // Begin_Latex
00640    // #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)
00641    // End_Latex
00642    double value = - x / kSqrt2 - kPi / 8;
00643    value += 1 / (8 * kSqrt2 * x) - 1 / (16 * x * x) + 25 / (384 * kSqrt2 * x * x * x);
00644    return value;
00645 }
00646 
00647 
00648 } // namespace Math
00649 } // namespace ROOT
00650  
00651  
00652  

Generated on Tue Jul 5 14:34:57 2011 for ROOT_528-00b_version by  doxygen 1.5.1