Special functions


Special Functions from MathCore

double ROOT::Math::erf (double x)
double ROOT::Math::erfc (double x)
double ROOT::Math::tgamma (double x)
double ROOT::Math::lgamma (double x)
double ROOT::Math::inc_gamma (double a, double x)
double ROOT::Math::inc_gamma_c (double a, double x)
double ROOT::Math::beta (double x, double y)
double ROOT::Math::inc_beta (double x, double a, double b)
double ROOT::Math::sinint (double x)
double ROOT::Math::cosint (double x)

Special Functions from MathMore

double ROOT::Math::assoc_laguerre (unsigned n, double m, double x)
double ROOT::Math::assoc_legendre (unsigned l, unsigned m, double x)
double ROOT::Math::comp_ellint_1 (double k)
double ROOT::Math::comp_ellint_2 (double k)
double ROOT::Math::comp_ellint_3 (double n, double k)
double ROOT::Math::conf_hyperg (double a, double b, double z)
double ROOT::Math::conf_hypergU (double a, double b, double z)
double ROOT::Math::cyl_bessel_i (double nu, double x)
double ROOT::Math::cyl_bessel_j (double nu, double x)
double ROOT::Math::cyl_bessel_k (double nu, double x)
double ROOT::Math::cyl_neumann (double nu, double x)
double ROOT::Math::ellint_1 (double k, double phi)
double ROOT::Math::ellint_2 (double k, double phi)
double ROOT::Math::ellint_3 (double n, double k, double phi)
double ROOT::Math::expint (double x)
double ROOT::Math::hyperg (double a, double b, double c, double x)
double ROOT::Math::laguerre (unsigned n, double x)
double ROOT::Math::legendre (unsigned l, double x)
double ROOT::Math::riemann_zeta (double x)
double ROOT::Math::sph_bessel (unsigned n, double x)
double ROOT::Math::sph_legendre (unsigned l, unsigned m, double theta)
double ROOT::Math::sph_neumann (unsigned n, double x)
double ROOT::Math::airy_Ai (double x)
double ROOT::Math::airy_Bi (double x)
double ROOT::Math::airy_Ai_deriv (double x)
double ROOT::Math::airy_Bi_deriv (double x)
double ROOT::Math::airy_zero_Ai (unsigned int s)
double ROOT::Math::airy_zero_Bi (unsigned int s)
double ROOT::Math::airy_zero_Ai_deriv (unsigned int s)
double ROOT::Math::airy_zero_Bi_deriv (unsigned int s)
double ROOT::Math::wigner_3j (int ja, int jb, int jc, int ma, int mb, int mc)
double ROOT::Math::wigner_6j (int ja, int jb, int jc, int jd, int je, int jf)
double ROOT::Math::wigner_9j (int ja, int jb, int jc, int jd, int je, int jf, int jg, int jh, int ji)

Detailed Description

Special mathematical functions. The naming and numbering of the functions is taken from Matt Austern, (Draft) Technical Report on Standard Library Extensions, N1687=04-0127, September 10, 2004

Author:
Created by Andras Zsenei on Mon Nov 8 2004

Function Documentation

double ROOT::Math::airy_Ai ( double  x  ) 

Calculates the Airy function Ai

\[ Ai(x) = \frac{1}{\pi} \int\limits_{0}^{\infty} \cos(xt + t^3/3) dt \]

For detailed description see Mathworld and Abramowitz&Stegun, Sect. 10.4. The implementation used is that of GSL.

Definition at line 364 of file SpecFuncMathMore.cxx.

Referenced by G__G__MathMore_99_0_26(), G__setup_memfuncROOTcLcLMath(), and testSpecFunc().

double ROOT::Math::airy_Ai_deriv ( double  x  ) 

Calculates the derivative of the Airy function Ai

\[ Ai(x) = \frac{1}{\pi} \int\limits_{0}^{\infty} \cos(xt + t^3/3) dt \]

For detailed description see Mathworld and Abramowitz&Stegun, Sect. 10.4. The implementation used is that of GSL.

Definition at line 380 of file SpecFuncMathMore.cxx.

Referenced by G__G__MathMore_99_0_28(), G__setup_memfuncROOTcLcLMath(), and testSpecFunc().

double ROOT::Math::airy_Bi ( double  x  ) 

Calculates the Airy function Bi

\[ Bi(x) = \frac{1}{\pi} \int\limits_{0}^{\infty} [\exp(xt-t^3/3) + \cos(xt + t^3/3)] dt \]

For detailed description see Mathworld and Abramowitz&Stegun, Sect. 10.4. The implementation used is that of GSL.

Definition at line 372 of file SpecFuncMathMore.cxx.

Referenced by G__G__MathMore_99_0_27(), G__setup_memfuncROOTcLcLMath(), and testSpecFunc().

double ROOT::Math::airy_Bi_deriv ( double  x  ) 

Calculates the derivative of the Airy function Bi

\[ Bi(x) = \frac{1}{\pi} \int\limits_{0}^{\infty} [\exp(xt-t^3/3) + \cos(xt + t^3/3)] dt \]

For detailed description see Mathworld and Abramowitz&Stegun, Sect. 10.4. The implementation used is that of GSL.

Definition at line 388 of file SpecFuncMathMore.cxx.

Referenced by G__G__MathMore_99_0_29(), G__setup_memfuncROOTcLcLMath(), and testSpecFunc().

double ROOT::Math::airy_zero_Ai ( unsigned int  s  ) 

Calculates the zeroes of the Airy function Ai

\[ Ai(x) = \frac{1}{\pi} \int\limits_{0}^{\infty} \cos(xt + t^3/3) dt \]

For detailed description see Mathworld and Abramowitz&Stegun, Sect. 10.4. The implementation used is that of GSL.

Definition at line 396 of file SpecFuncMathMore.cxx.

Referenced by G__G__MathMore_99_0_30(), G__setup_memfuncROOTcLcLMath(), and testSpecFunc().

double ROOT::Math::airy_zero_Ai_deriv ( unsigned int  s  ) 

Calculates the zeroes of the derivative of the Airy function Ai

\[ Ai(x) = \frac{1}{\pi} \int\limits_{0}^{\infty} \cos(xt + t^3/3) dt \]

For detailed description see Mathworld and Abramowitz&Stegun, Sect. 10.4. The implementation used is that of GSL.

Definition at line 412 of file SpecFuncMathMore.cxx.

Referenced by G__G__MathMore_99_0_32(), G__setup_memfuncROOTcLcLMath(), and testSpecFunc().

double ROOT::Math::airy_zero_Bi ( unsigned int  s  ) 

Calculates the zeroes of the Airy function Bi

\[ Bi(x) = \frac{1}{\pi} \int\limits_{0}^{\infty} [\exp(xt-t^3/3) + \cos(xt + t^3/3)] dt \]

For detailed description see Mathworld and Abramowitz&Stegun, Sect. 10.4. The implementation used is that of GSL.

Definition at line 404 of file SpecFuncMathMore.cxx.

Referenced by G__G__MathMore_99_0_31(), G__setup_memfuncROOTcLcLMath(), and testSpecFunc().

double ROOT::Math::airy_zero_Bi_deriv ( unsigned int  s  ) 

Calculates the zeroes of the derivative of the Airy function Bi

\[ Bi(x) = \frac{1}{\pi} \int\limits_{0}^{\infty} [\exp(xt-t^3/3) + \cos(xt + t^3/3)] dt \]

For detailed description see Mathworld and Abramowitz&Stegun, Sect. 10.4. The implementation used is that of GSL.

Definition at line 420 of file SpecFuncMathMore.cxx.

Referenced by G__G__MathMore_99_0_33(), G__setup_memfuncROOTcLcLMath(), and testSpecFunc().

double ROOT::Math::assoc_laguerre ( unsigned  n,
double  m,
double  x 
)

Computes the generalized Laguerre polynomials for $ n \geq 0, m > -1 $. They are defined in terms of the confluent hypergeometric function. For integer values of m they can be defined in terms of the Laguerre polynomials $L_n(x)$:

\[ L_{n}^{m}(x) = (-1)^{m} \frac{d^m}{dx^m} L_{n+m}(x) \]

For detailed description see Mathworld. The implementation used is that of GSL.

This function is an extension of C++0x, also consistent in GSL, Abramowitz and Stegun 1972 and MatheMathica that uses non-integer values for m. C++0x calls for 'int m', more restrictive than necessary. The definition for was incorrect in 'n1687.pdf', but fixed in n1836.pdf, the most recent draft as of 2007-07-01

Definition at line 41 of file SpecFuncMathMore.cxx.

Referenced by G__G__MathMore_99_0_4(), G__setup_memfuncROOTcLcLMath(), RooMathMoreReg::RooMathMoreReg(), and testSpecFunc().

double ROOT::Math::assoc_legendre ( unsigned  l,
unsigned  m,
double  x 
)

Computes the associated Legendre polynomials.

\[ P_{l}^{m}(x) = (1-x^2)^{m/2} \frac{d^m}{dx^m} P_{l}(x) \]

with $m \geq 0$, $ l \geq m $ and $ |x|<1 $. There are two sign conventions for associated Legendre polynomials. As is the case with the above formula, some authors (e.g., Arfken 1985, pp. 668-669) omit the Condon-Shortley phase $(-1)^m$, while others include it (e.g., Abramowitz and Stegun 1972). One possible way to distinguish the two conventions is due to Abramowitz and Stegun (1972, p. 332), who use the notation

\[ P_{lm} (x) = (-1)^m P_{l}^{m} (x)\]

to distinguish the two. For detailed description see Mathworld. The implementation used is that of GSL.

The definition uses is the one of C++0x, $ P_{lm}$, while GSL and MatheMatica use the $P_{l}^{m}$ definition. Note that C++0x and GSL definitions agree instead for the normalized associated Legendre polynomial, sph_legendre(l,m,theta).

Definition at line 52 of file SpecFuncMathMore.cxx.

Referenced by RooLegendre::evaluate(), G__G__MathMore_99_0_5(), G__setup_memfuncROOTcLcLMath(), RooMathMoreReg::RooMathMoreReg(), and testSpecFunc().

double ROOT::Math::beta ( double  x,
double  y 
)

Calculates the beta function.

\[ B(x,y) = \frac{\Gamma(x) \Gamma(y)}{\Gamma(x+y)} \]

for x>0 and y>0. For detailed description see Mathworld.

Definition at line 111 of file SpecFuncMathCore.cxx.

References exp(), ROOT::Math::Cephes::lgam(), and ROOT::Math::lgamma().

Referenced by TMath::Beta(), TMath::BetaDist(), TDecompSVD::Bidiagonalize(), TEfficiency::Combine(), TGraphAsymmErrors::Divide(), G__G__MathCore_170_0_7(), G__setup_memfuncROOTcLcLMath(), TGenPhaseSpace::Generate(), TEfficiency::GetEfficiency(), TEfficiency::GetEfficiencyErrorLow(), TEfficiency::GetEfficiencyErrorUp(), RooHistError::getInterval(), TF1::GetQuantiles(), ROOT::Math::KelvinFunctions::Kei(), ROOT::Math::KelvinFunctions::Ker(), main(), ROOT::Minuit2::SimplexBuilder::Minimum(), TMinuit::mnsimp(), TSpectrum2Painter::Paint(), TGaxis::PaintAxis(), ROOT::Math::Boost::Rectify(), ROOT::Math::BoostX::Rectify(), ROOT::Math::BoostY::Rectify(), ROOT::Math::BoostZ::Rectify(), rf105_funcbinding(), RooMathCoreReg::RooMathCoreReg(), MyEvent::ScatterAngle(), testBetaFunction(), TestBasic105::testCode(), testLorentzVector(), testSpecFunc(), testSpecFuncBeta(), and testStatFunctions().

double ROOT::Math::comp_ellint_1 ( double  k  ) 

Calculates the complete elliptic integral of the first kind.

\[ K(k) = F(k, \pi / 2) = \int_{0}^{\pi /2} \frac{d \theta}{\sqrt{1 - k^2 \sin^2{\theta}}} \]

with $0 \leq k^2 \leq 1$. For detailed description see Mathworld. The implementation used is that of GSL.

Definition at line 65 of file SpecFuncMathMore.cxx.

Referenced by G__G__MathMore_99_0_6(), G__setup_memfuncROOTcLcLMath(), RooMathMoreReg::RooMathMoreReg(), and testSpecFunc().

double ROOT::Math::comp_ellint_2 ( double  k  ) 

Calculates the complete elliptic integral of the second kind.

\[ E(k) = E(k , \pi / 2) = \int_{0}^{\pi /2} \sqrt{1 - k^2 \sin^2{\theta}} d \theta \]

with $0 \leq k^2 \leq 1$. For detailed description see Mathworld. The implementation used is that of GSL.

Definition at line 76 of file SpecFuncMathMore.cxx.

Referenced by G__G__MathMore_99_0_7(), G__setup_memfuncROOTcLcLMath(), RooMathMoreReg::RooMathMoreReg(), and testSpecFunc().

double ROOT::Math::comp_ellint_3 ( double  n,
double  k 
)

Calculates the complete elliptic integral of the third kind.

\[ \Pi (n, k, \pi / 2) = \int_{0}^{\pi /2} \frac{d \theta}{(1 - n \sin^2{\theta})\sqrt{1 - k^2 \sin^2{\theta}}} \]

with $0 \leq k^2 \leq 1$. There are two sign conventions for elliptic integrals of the third kind. Some authors (Abramowitz and Stegun, Mathworld, C++ standard proposal) use the above formula, while others (GSL, Planetmath and CERNLIB) use the + sign in front of n in the denominator. In order to be C++ compliant, the present library uses the former convention. The implementation used is that of GSL.

Definition at line 123 of file SpecFuncMathMore.cxx.

References PI.

Referenced by G__G__MathMore_99_0_8(), G__setup_memfuncROOTcLcLMath(), and testSpecFunc().

double ROOT::Math::conf_hyperg ( double  a,
double  b,
double  z 
)

Calculates the confluent hypergeometric functions of the first kind.

\[ _{1}F_{1}(a;b;z) = \frac{\Gamma(b)}{\Gamma(a)} \sum_{n=0}^{\infty} \frac{\Gamma(a+n)}{\Gamma(b+n)} \frac{z^n}{n!} \]

For detailed description see Mathworld. The implementation used is that of GSL.

Definition at line 134 of file SpecFuncMathMore.cxx.

Referenced by G__G__MathMore_99_0_9(), G__setup_memfuncROOTcLcLMath(), RooMathMoreReg::RooMathMoreReg(), and testSpecFunc().

double ROOT::Math::conf_hypergU ( double  a,
double  b,
double  z 
)

Calculates the confluent hypergeometric functions of the second kind, known also as Kummer function of the second kind, it is related to the confluent hypergeometric functions of the first kind.

\[ U(a,b,z) = \frac{ \pi}{ \sin{\pi b } } \left[ \frac{ _{1}F_{1}(a,b,z) } {\Gamma(a-b+1) } - \frac{ z^{1-b} { _{1}F_{1}}(a-b+1,2-b,z)}{\Gamma(a)} \right] \]

For detailed description see Mathworld. The implementation used is that of GSL. This function is not part of the C++ standard proposal

Definition at line 142 of file SpecFuncMathMore.cxx.

Referenced by G__G__MathMore_99_0_10(), G__setup_memfuncROOTcLcLMath(), and RooMathMoreReg::RooMathMoreReg().

double ROOT::Math::cosint ( double  x  ) 

Calculates the real part of the cosine integral (Ci).

For x<0, the imaginary part is i and has to be added by the user, for x>0 the imaginary part of Ci(x) is 0.

\[ Ci(x) = - \int_{x}^{\infty} \frac{\cos t}{t} dt = \gamma + \ln x + \int_{0}^{x} \frac{\cos t - 1}{t} dt\]

For detailed description see Mathworld. The implementation used is that of CERNLIB, based on Y.L. Luke, The special functions and their approximations, v.II, (Academic Press, New York l969) 325-326.

Definition at line 212 of file SpecFuncMathCore.cxx.

References c, cos(), h, i, RootCsg::infinity, log(), p, and sin().

Referenced by G__G__MathCore_170_0_10(), G__setup_memfuncROOTcLcLMath(), ROOT::Math::VavilovAccurate::Set(), and testSiCi().

double ROOT::Math::cyl_bessel_i ( double  nu,
double  x 
)

Calculates the modified Bessel function of the first kind (also called regular modified (cylindrical) Bessel function).

\[ I_{\nu} (x) = i^{-\nu} J_{\nu}(ix) = \sum_{k=0}^{\infty} \frac{(\frac{1}{2}x)^{\nu + 2k}}{k! \Gamma(\nu + k + 1)} \]

for $x>0, \nu > 0$. For detailed description see Mathworld. The implementation used is that of GSL.

Definition at line 153 of file SpecFuncMathMore.cxx.

Referenced by G__G__MathMore_99_0_11(), G__setup_memfuncROOTcLcLMath(), ROOT::Math::noncentral_chisquared_pdf(), RooMathMoreReg::RooMathMoreReg(), and testSpecFunc().

double ROOT::Math::cyl_bessel_j ( double  nu,
double  x 
)

Calculates the (cylindrical) Bessel functions of the first kind (also called regular (cylindrical) Bessel functions).

\[ J_{\nu} (x) = \sum_{k=0}^{\infty} \frac{(-1)^k(\frac{1}{2}x)^{\nu + 2k}}{k! \Gamma(\nu + k + 1)} \]

For detailed description see Mathworld. The implementation used is that of GSL.

Definition at line 164 of file SpecFuncMathMore.cxx.

Referenced by G__G__MathMore_99_0_12(), G__setup_memfuncROOTcLcLMath(), RooMathMoreReg::RooMathMoreReg(), and testSpecFunc().

double ROOT::Math::cyl_bessel_k ( double  nu,
double  x 
)

Calculates the modified Bessel functions of the second kind (also called irregular modified (cylindrical) Bessel functions).

\[ K_{\nu} (x) = \frac{\pi}{2} i^{\nu + 1} (J_{\nu} (ix) + iN(ix)) = \left\{ \begin{array}{cl} \frac{\pi}{2} \frac{I_{-\nu}(x) - I_{\nu}(x)}{\sin{\nu \pi}} & \mbox{for non-integral $\nu$} \\ \frac{\pi}{2} \lim{\mu \to \nu} \frac{I_{-\mu}(x) - I_{\mu}(x)}{\sin{\mu \pi}} & \mbox{for integral $\nu$} \end{array} \right. \]

for $x>0, \nu > 0$. For detailed description see Mathworld. The implementation used is that of GSL.

Definition at line 175 of file SpecFuncMathMore.cxx.

Referenced by G__G__MathMore_99_0_13(), G__setup_memfuncROOTcLcLMath(), RooMathMoreReg::RooMathMoreReg(), and testSpecFunc().

double ROOT::Math::cyl_neumann ( double  nu,
double  x 
)

Calculates the (cylindrical) Bessel functions of the second kind (also called irregular (cylindrical) Bessel functions or (cylindrical) Neumann functions).

\[ N_{\nu} (x) = Y_{\nu} (x) = \left\{ \begin{array}{cl} \frac{J_{\nu} \cos{\nu \pi}-J_{-\nu}(x)}{\sin{\nu \pi}} & \mbox{for non-integral $\nu$} \\ \lim{\mu \to \nu} \frac{J_{\mu} \cos{\mu \pi}-J_{-\mu}(x)}{\sin{\mu \pi}} & \mbox{for integral $\nu$} \end{array} \right. \]

For detailed description see Mathworld. The implementation used is that of GSL.

Definition at line 187 of file SpecFuncMathMore.cxx.

Referenced by G__G__MathMore_99_0_14(), G__setup_memfuncROOTcLcLMath(), RooMathMoreReg::RooMathMoreReg(), and testSpecFunc().

double ROOT::Math::ellint_1 ( double  k,
double  phi 
)

Calculates the incomplete elliptic integral of the first kind.

\[ F(k, \phi) = \int_{0}^{\phi} \frac{d \theta}{\sqrt{1 - k^2 \sin^2{\theta}}} \]

with $0 \leq k^2 \leq 1$. For detailed description see Mathworld. The implementation used is that of GSL.

Parameters:
k 
phi angle in radians

Definition at line 199 of file SpecFuncMathMore.cxx.

Referenced by G__G__MathMore_99_0_15(), G__setup_memfuncROOTcLcLMath(), RooMathMoreReg::RooMathMoreReg(), and testSpecFunc().

double ROOT::Math::ellint_2 ( double  k,
double  phi 
)

Calculates the complete elliptic integral of the second kind.

\[ E(k , \phi) = \int_{0}^{\phi} \sqrt{1 - k^2 \sin^2{\theta}} d \theta \]

with $0 \leq k^2 \leq 1$. For detailed description see Mathworld. The implementation used is that of GSL.

Parameters:
k 
phi angle in radians

Definition at line 211 of file SpecFuncMathMore.cxx.

Referenced by G__G__MathMore_99_0_16(), G__setup_memfuncROOTcLcLMath(), RooMathMoreReg::RooMathMoreReg(), and testSpecFunc().

double ROOT::Math::ellint_3 ( double  n,
double  k,
double  phi 
)

Calculates the complete elliptic integral of the third kind.

\[ \Pi (n, k, \phi) = \int_{0}^{\phi} \frac{d \theta}{(1 - n \sin^2{\theta})\sqrt{1 - k^2 \sin^2{\theta}}} \]

with $0 \leq k^2 \leq 1$. There are two sign conventions for elliptic integrals of the third kind. Some authors (Abramowitz and Stegun, Mathworld, C++ standard proposal) use the above formula, while others (GSL, Planetmath and CERNLIB) use the + sign in front of n in the denominator. In order to be C++ compliant, the present library uses the former convention. The implementation used is that of GSL.

Parameters:
n 
k 
phi angle in radians

Definition at line 259 of file SpecFuncMathMore.cxx.

Referenced by G__G__MathMore_99_0_17(), G__setup_memfuncROOTcLcLMath(), RooMathMoreReg::RooMathMoreReg(), and testSpecFunc().

double ROOT::Math::erf ( double  x  ) 

Error function encountered in integrating the normal distribution.

\[ erf(x) = \frac{2}{\sqrt{\pi}} \int_{0}^{x} e^{-t^2} dt \]

For detailed description see Mathworld. The implementation used is that of GSL. This function is provided only for convenience, in case your standard C++ implementation does not support it. If it does, please use these standard version!

Definition at line 59 of file SpecFuncMathCore.cxx.

References ROOT::Math::Cephes::erf(), ROOT::Math::Cephes::erfc(), ROOT::Math::Cephes::erfT, ROOT::Math::Cephes::erfU, ROOT::Math::Polynomial1eval(), ROOT::Math::Polynomialeval(), and y.

Referenced by TMath::Erf(), ROOT::Math::erf(), ROOT::Math::Cephes::erfc(), G__G__MathCore_170_0_1(), G__setup_memfuncROOTcLcLMath(), ROOT::Math::lognormal_cdf(), ROOT::Math::lognormal_cdf_c(), ROOT::Math::normal_cdf(), ROOT::Math::normal_cdf_c(), rf105_funcbinding(), RooMathCoreReg::RooMathCoreReg(), TestBasic105::testCode(), testSpecFunc(), and testSpecFuncErf().

double ROOT::Math::erfc ( double  x  ) 

Complementary error function.

\[ erfc(x) = 1 - erf(x) = \frac{2}{\sqrt{\pi}} \int_{x}^{\infty} e^{-t^2} dt \]

For detailed description see Mathworld. The implementation used is that of Cephes from S. Moshier.

Definition at line 44 of file SpecFuncMathCore.cxx.

References ROOT::Math::Cephes::erf(), ROOT::Math::Cephes::erfc(), ROOT::Math::Cephes::erfP, ROOT::Math::Cephes::erfQ, ROOT::Math::Cephes::erfR, ROOT::Math::Cephes::erfS, exp(), kMAXLOG, p, ROOT::Math::Polynomial1eval(), ROOT::Math::Polynomialeval(), x, and y.

Referenced by ROOT::Math::Cephes::erf(), TMath::Erfc(), ROOT::Math::erfc(), G__G__MathCore_170_0_2(), G__setup_memfuncROOTcLcLMath(), ROOT::Math::lognormal_cdf(), ROOT::Math::lognormal_cdf_c(), ROOT::Math::normal_cdf(), ROOT::Math::normal_cdf_c(), RooMathCoreReg::RooMathCoreReg(), testSpecFunc(), and testSpecFuncErf().

double ROOT::Math::expint ( double  x  ) 

Calculates the exponential integral.

\[ Ei(x) = - \int_{-x}^{\infty} \frac{e^{-t}}{t} dt \]

For detailed description see Mathworld. The implementation used is that of GSL.

Definition at line 270 of file SpecFuncMathMore.cxx.

Referenced by ROOT::Math::VavilovAccurate::E1plLog(), G__G__MathMore_99_0_18(), G__setup_memfuncROOTcLcLMath(), RooMathMoreReg::RooMathMoreReg(), and testSpecFunc().

double ROOT::Math::hyperg ( double  a,
double  b,
double  c,
double  x 
)

Calculates Gauss' hypergeometric function.

\[ _{2}F_{1}(a,b;c;x) = \frac{\Gamma(c)}{\Gamma(a) \Gamma(b)} \sum_{n=0}^{\infty} \frac{\Gamma(a+n)\Gamma(b+n)}{\Gamma(c+n)} \frac{x^n}{n!} \]

For detailed description see Mathworld. The implementation used is that of GSL.

Definition at line 290 of file SpecFuncMathMore.cxx.

Referenced by G__G__MathMore_99_0_19(), G__setup_memfuncROOTcLcLMath(), and testSpecFunc().

double ROOT::Math::inc_beta ( double  x,
double  a,
double  b 
)

Calculates the normalized (regularized) incomplete beta function.

\[ B(x, a, b ) = \frac{ \int_{0}^{x} u^{a-1} (1-u)^{b-1} du } { B(a,b) } \]

for 0<=x<=1, a>0, and b>0. For detailed description see Mathworld. The implementation used is that of Cephes from S. Moshier.

Definition at line 115 of file SpecFuncMathCore.cxx.

References ROOT::Math::Cephes::incbet().

Referenced by ROOT::Math::beta_cdf(), ROOT::Math::beta_cdf_c(), TMath::BetaIncomplete(), ROOT::Math::fdistribution_cdf(), ROOT::Math::fdistribution_cdf_c(), G__G__MathCore_170_0_8(), G__setup_memfuncROOTcLcLMath(), RooMathCoreReg::RooMathCoreReg(), ROOT::Math::tdistribution_cdf(), ROOT::Math::tdistribution_cdf_c(), testSpecFunc(), and testSpecFuncBetaI().

double ROOT::Math::inc_gamma ( double  a,
double  x 
)

Calculates the normalized (regularized) lower incomplete gamma function (lower integral)

\[ P(a, x) = \frac{ 1} {\Gamma(a) } \int_{0}^{x} t^{a-1} e^{-t} dt \]

For a detailed description see Mathworld. The implementation used is that of Cephes from S. Moshier. In this implementation both a and x must be positive. If a is negative 1.0 is returned for every x. This is correct only if a is negative integer. For a>0 and x<0 0 is returned (this is correct only for a>0 and x=0).

Definition at line 99 of file SpecFuncMathCore.cxx.

References ROOT::Math::Cephes::igam().

Referenced by ROOT::Math::chisquared_cdf(), G__G__MathCore_170_0_5(), G__setup_memfuncROOTcLcLMath(), TMath::Gamma(), ROOT::Math::gamma_cdf(), RooMathCoreReg::RooMathCoreReg(), testSpecFunc(), and testSpecFuncGamma().

double ROOT::Math::inc_gamma_c ( double  a,
double  x 
)

Calculates the normalized (regularized) upper incomplete gamma function (upper integral)

\[ Q(a, x) = \frac{ 1} {\Gamma(a) } \int_{x}^{\infty} t^{a-1} e^{-t} dt \]

For a detailed description see Mathworld. The implementation used is that of Cephes from S. Moshier. In this implementation both a and x must be positive. If a is negative, 0 is returned for every x. This is correct only if a is negative integer. For a>0 and x<0 1 is returned (this is correct only for a>0 and x=0).

Definition at line 103 of file SpecFuncMathCore.cxx.

References ROOT::Math::Cephes::igamc().

Referenced by ROOT::Math::chisquared_cdf_c(), G__G__MathCore_170_0_6(), G__setup_memfuncROOTcLcLMath(), ROOT::Math::gamma_cdf_c(), RooMathCoreReg::RooMathCoreReg(), and testSpecFunc().

double ROOT::Math::laguerre ( unsigned  n,
double  x 
)

Calculates the Laguerre polynomials

\[ P_{l}(x) = \frac{ e^x}{n!} \frac{d^n}{dx^n} (x^n - e^{-x}) \]

for $x \geq 0 $ in the Rodrigues representation. They corresponds to the associated Laguerre polynomial of order m=0. See Abramowitz and Stegun, (22.5.16) For detailed description see Mathworld. The are implemented using the associated Laguerre polynomial of order m=0.

Definition at line 301 of file SpecFuncMathMore.cxx.

Referenced by G__G__MathMore_99_0_20(), G__setup_memfuncROOTcLcLMath(), RooMathMoreReg::RooMathMoreReg(), and testSpecFunc().

double ROOT::Math::legendre ( unsigned  l,
double  x 
)

Calculates the Legendre polynomials.

\[ P_{l}(x) = \frac{1}{2^l l!} \frac{d^l}{dx^l} (x^2 - 1)^l \]

for $l \geq 0, |x|\leq1$ in the Rodrigues representation. For detailed description see Mathworld. The implementation used is that of GSL.

Definition at line 311 of file SpecFuncMathMore.cxx.

Referenced by G__G__MathMore_99_0_21(), G__setup_memfuncROOTcLcLMath(), RooMathMoreReg::RooMathMoreReg(), and testSpecFunc().

double ROOT::Math::lgamma ( double  x  ) 

Calculates the logarithm of the gamma function

The implementation used is that of Cephes from S. Moshier.

Definition at line 74 of file SpecFuncMathCore.cxx.

References ROOT::Math::Cephes::lgam().

Referenced by ROOT::Math::beta(), ROOT::Math::beta_pdf(), ROOT::Math::binomial_pdf(), ROOT::Math::chisquared_pdf(), ROOT::Math::fdistribution_pdf(), G__G__MathCore_170_0_4(), G__setup_memfuncROOTcLcLMath(), ROOT::Math::gamma_pdf(), TMath::LnGamma(), ROOT::Math::negative_binomial_pdf(), ROOT::Math::poisson_pdf(), RooMathCoreReg::RooMathCoreReg(), ROOT::Math::tdistribution_pdf(), testSpecFunc(), and testSpecFuncGamma().

double ROOT::Math::riemann_zeta ( double  x  ) 

Calculates the Riemann zeta function.

\[ \zeta (x) = \left\{ \begin{array}{cl} \sum_{k=1}^{\infty}k^{-x} & \mbox{for $x > 1$} \\ 2^x \pi^{x-1} \sin{(\frac{1}{2}\pi x)} \Gamma(1-x) \zeta (1-x) & \mbox{for $x < 1$} \end{array} \right. \]

For detailed description see Mathworld. The implementation used is that of GSL.

CHECK WHETHER THE IMPLEMENTATION CALCULATES X<1

Definition at line 322 of file SpecFuncMathMore.cxx.

Referenced by G__G__MathMore_99_0_22(), G__setup_memfuncROOTcLcLMath(), RooMathMoreReg::RooMathMoreReg(), and testSpecFunc().

double ROOT::Math::sinint ( double  x  ) 

Calculates the sine integral.

\[ Si(x) = - \int_{0}^{x} \frac{\sin t}{t} dt \]

For detailed description see Mathworld. The implementation used is that of CERNLIB, based on Y.L. Luke, The special functions and their approximations, v.II, (Academic Press, New York l969) 325-326.

Definition at line 122 of file SpecFuncMathCore.cxx.

References cos(), h, i, p, PI, r8, s, sin(), and y.

Referenced by G__G__MathCore_170_0_9(), G__setup_memfuncROOTcLcLMath(), ROOT::Math::VavilovAccurate::Set(), and testSiCi().

double ROOT::Math::sph_bessel ( unsigned  n,
double  x 
)

Calculates the spherical Bessel functions of the first kind (also called regular spherical Bessel functions).

\[ j_{n}(x) = \sqrt{\frac{\pi}{2x}} J_{n+1/2}(x) \]

For detailed description see Mathworld. The implementation used is that of GSL.

Definition at line 333 of file SpecFuncMathMore.cxx.

Referenced by G__G__MathMore_99_0_23(), G__setup_memfuncROOTcLcLMath(), RooMathMoreReg::RooMathMoreReg(), and testSpecFunc().

double ROOT::Math::sph_legendre ( unsigned  l,
unsigned  m,
double  theta 
)

Computes the spherical (normalized) associated Legendre polynomials, or spherical harmonic without azimuthal dependence ($e^(im\phi)$).

\[ Y_l^m(theta,0) = \sqrt{(2l+1)/(4\pi)} \sqrt{(l-m)!/(l+m)!} P_l^m(cos \theta) \]

for $m \geq 0, l \geq m$, where the Condon-Shortley phase $(-1)^m$ is included in P_l^m(x) This function is consistent with both C++0x and GSL, even though there is a discrepancy in where to include the phase. There is no reference in Abramowitz and Stegun.

Definition at line 344 of file SpecFuncMathMore.cxx.

References cos().

Referenced by G__G__MathMore_99_0_24(), G__setup_memfuncROOTcLcLMath(), RooMathMoreReg::RooMathMoreReg(), and testSpecFunc().

double ROOT::Math::sph_neumann ( unsigned  n,
double  x 
)

Calculates the spherical Bessel functions of the second kind (also called irregular spherical Bessel functions or spherical Neumann functions).

\[ n_n(x) = y_n(x) = \sqrt{\frac{\pi}{2x}} N_{n+1/2}(x) \]

For detailed description see Mathworld. The implementation used is that of GSL.

Definition at line 356 of file SpecFuncMathMore.cxx.

Referenced by G__G__MathMore_99_0_25(), G__setup_memfuncROOTcLcLMath(), RooMathMoreReg::RooMathMoreReg(), and testSpecFunc().

double ROOT::Math::tgamma ( double  x  ) 

The gamma function is defined to be the extension of the factorial to real numbers.

\[ \Gamma(x) = \int_{0}^{\infty} t^{x-1} e^{-t} dt \]

For detailed description see Mathworld. The implementation used is that of Cephes from S. Moshier.

Definition at line 89 of file SpecFuncMathCore.cxx.

References ROOT::Math::Cephes::gamma().

Referenced by GammaFunction::DoEval(), G__G__MathCore_170_0_3(), G__setup_memfuncROOTcLcLMath(), TMath::Gamma(), gamma_func(), ROOT::Math::noncentral_chisquared_pdf(), RooMathCoreReg::RooMathCoreReg(), testSpecFunc(), and testSpecFuncGamma().

double ROOT::Math::wigner_3j ( int  ja,
int  jb,
int  jc,
int  ma,
int  mb,
int  mc 
)

Calculates the Wigner 3j coupling coefficients

For detailed description see Mathworld. The implementation used is that of GSL.

Definition at line 428 of file SpecFuncMathMore.cxx.

Referenced by G__G__MathMore_99_0_34(), and G__setup_memfuncROOTcLcLMath().

double ROOT::Math::wigner_6j ( int  ja,
int  jb,
int  jc,
int  jd,
int  je,
int  jf 
)

Calculates the Wigner 6j coupling coefficients

For detailed description see Mathworld. The implementation used is that of GSL.

Definition at line 432 of file SpecFuncMathMore.cxx.

Referenced by G__G__MathMore_99_0_35(), and G__setup_memfuncROOTcLcLMath().

double ROOT::Math::wigner_9j ( int  ja,
int  jb,
int  jc,
int  jd,
int  je,
int  jf,
int  jg,
int  jh,
int  ji 
)

Calculates the Wigner 9j coupling coefficients

For detailed description see Mathworld. The implementation used is that of GSL.

Definition at line 436 of file SpecFuncMathMore.cxx.

Referenced by G__G__MathMore_99_0_36(), and G__setup_memfuncROOTcLcLMath().


Generated on Tue Jul 5 16:09:42 2011 for ROOT_528-00b_version by  doxygen 1.5.1