Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 10:10:21

0001 // @(#)root/mathmore:$Id$
0002 // Authors: L. Moneta, A. Zsenei   08/2005
0003 
0004 // Authors: Andras Zsenei & Lorenzo Moneta   06/2005
0005 
0006  /**********************************************************************
0007   *                                                                    *
0008   * Copyright (c) 2004 ROOT Foundation,  CERN/PH-SFT                   *
0009   *                                                                    *
0010   * This library is free software; you can redistribute it and/or      *
0011   * modify it under the terms of the GNU General Public License        *
0012   * as published by the Free Software Foundation; either version 2     *
0013   * of the License, or (at your option) any later version.             *
0014   *                                                                    *
0015   * This library is distributed in the hope that it will be useful,    *
0016   * but WITHOUT ANY WARRANTY; without even the implied warranty of     *
0017   * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU   *
0018   * General Public License for more details.                           *
0019   *                                                                    *
0020   * You should have received a copy of the GNU General Public License  *
0021   * along with this library (see file COPYING); if not, write          *
0022   * to the Free Software Foundation, Inc., 59 Temple Place, Suite      *
0023   * 330, Boston, MA 02111-1307 USA, or contact the author.             *
0024   *                                                                    *
0025   **********************************************************************/
0026 
0027 /**
0028 
0029 Special mathematical functions.
0030 The naming and numbering of the functions is taken from
0031 <A HREF="http://www.open-std.org/jtc1/sc22/wg21/docs/papers/2004/n1687.pdf">
0032 Matt Austern,
0033 (Draft) Technical Report on Standard Library Extensions,
0034 N1687=04-0127, September 10, 2004</A>
0035 
0036 @author Created by Andras Zsenei on Mon Nov 8 2004
0037 
0038 @defgroup SpecFunc Special functions
0039 
0040 */
0041 
0042 
0043 
0044 
0045 
0046 #ifndef ROOT_Math_SpecFuncMathMore
0047 #define ROOT_Math_SpecFuncMathMore
0048 
0049 
0050 
0051 
0052 namespace ROOT {
0053 namespace Math {
0054 
0055    /** @name Special Functions from MathMore */
0056 
0057 
0058   /**
0059 
0060 
0061   Computes the generalized Laguerre polynomials for
0062   \f$ n \geq 0, m > -1 \f$.
0063   They are defined in terms of the confluent hypergeometric function.
0064   For integer values of m they can be defined in terms of the Laguerre polynomials \f$L_n(x)\f$:
0065 
0066   \f[ L_{n}^{m}(x) = (-1)^{m} \frac{d^m}{dx^m} L_{n+m}(x) \f]
0067 
0068 
0069   For detailed description see
0070   <A HREF="http://mathworld.wolfram.com/LaguerrePolynomial.html">Mathworld</A>.
0071   The implementation used is that of
0072   <A HREF="http://www.gnu.org/software/gsl/manual/html_node/Laguerre-Functions.html">GSL</A>.
0073 
0074   This function is an extension of C++0x, also consistent in GSL,
0075   Abramowitz and Stegun 1972 and MatheMathica that uses non-integer values for m.
0076   C++0x calls for 'int m', more restrictive than necessary.
0077   The definition for was incorrect in 'n1687.pdf', but fixed in
0078   <A HREF="http://www.open-std.org/jtc1/sc22/wg21/docs/papers/2005/n1836.pdf">n1836.pdf</A>,
0079   the most recent draft as of 2007-07-01
0080 
0081 
0082   @ingroup SpecFunc
0083 
0084   */
0085   // [5.2.1.1] associated Laguerre polynomials
0086 
0087   double assoc_laguerre(unsigned n, double m, double x);
0088 
0089 
0090 
0091 
0092   /**
0093 
0094   Computes the associated Legendre polynomials.
0095 
0096   \f[ P_{l}^{m}(x) = (1-x^2)^{m/2} \frac{d^m}{dx^m} P_{l}(x) \f]
0097 
0098   with \f$m \geq 0\f$, \f$ l \geq m \f$ and \f$ |x|<1 \f$.
0099   There are two sign conventions for associated Legendre polynomials.
0100   As is the case with the above formula, some authors (e.g., Arfken
0101   1985, pp. 668-669) omit the Condon-Shortley phase \f$(-1)^m\f$,
0102   while others include it (e.g., Abramowitz and Stegun 1972).
0103   One possible way to distinguish the two conventions is due to
0104   Abramowitz and Stegun (1972, p. 332), who use the notation
0105 
0106   \f[ P_{lm} (x) = (-1)^m P_{l}^{m} (x)\f]
0107 
0108   to distinguish the two. For detailed description see
0109   <A HREF="http://mathworld.wolfram.com/LegendrePolynomial.html">
0110   Mathworld</A>. The implementation used is that of
0111   <A HREF="http://www.gnu.org/software/gsl/manual/html_node/Associated-Legendre-Polynomials-and-Spherical-Harmonics.html">GSL</A>.
0112 
0113   The definition uses is the one of C++0x, \f$ P_{lm}\f$, while GSL and MatheMatica use the \f$P_{l}^{m}\f$ definition. Note that C++0x and GSL definitions agree instead for the normalized associated Legendre polynomial,
0114   sph_legendre(l,m,theta).
0115 
0116   @ingroup SpecFunc
0117 
0118   */
0119   // [5.2.1.2] associated Legendre functions
0120 
0121   double assoc_legendre(unsigned l, unsigned m, double x);
0122 
0123   // Shortcut for RooFit to call the gsl legendre functions without the branches in the above implementation.
0124   namespace internal{
0125     double legendre(unsigned l, unsigned m, double x);
0126   }
0127 
0128 
0129   /**
0130 
0131   Calculates the complete elliptic integral of the first kind.
0132 
0133   \f[ K(k) = F(k, \pi / 2) = \int_{0}^{\pi /2} \frac{d \theta}{\sqrt{1 - k^2 \sin^2{\theta}}} \f]
0134 
0135   with \f$0 \leq k^2 \leq 1\f$. For detailed description see
0136   <A HREF="http://mathworld.wolfram.com/CompleteEllipticIntegraloftheFirstKind.html">
0137   Mathworld</A>. The implementation used is that of
0138   <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC100">GSL</A>.
0139 
0140   @ingroup SpecFunc
0141 
0142   */
0143   // [5.2.1.4] (complete) elliptic integral of the first kind
0144 
0145   double comp_ellint_1(double k);
0146 
0147 
0148 
0149 
0150   /**
0151 
0152   Calculates the complete elliptic integral of the second kind.
0153 
0154   \f[ E(k) = E(k , \pi / 2) = \int_{0}^{\pi /2} \sqrt{1 - k^2 \sin^2{\theta}} d \theta  \f]
0155 
0156   with \f$0 \leq k^2 \leq 1\f$. For detailed description see
0157   <A HREF="http://mathworld.wolfram.com/CompleteEllipticIntegraloftheSecondKind.html">
0158   Mathworld</A>. The implementation used is that of
0159   <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC100">GSL</A>.
0160 
0161   @ingroup SpecFunc
0162 
0163   */
0164   // [5.2.1.5] (complete) elliptic integral of the second kind
0165 
0166   double comp_ellint_2(double k);
0167 
0168 
0169 
0170 
0171   /**
0172 
0173   Calculates the complete elliptic integral of the third kind.
0174 
0175   \f[ \Pi (n, k, \pi / 2) = \int_{0}^{\pi /2} \frac{d \theta}{(1 - n \sin^2{\theta})\sqrt{1 - k^2 \sin^2{\theta}}}  \f]
0176 
0177   with \f$0 \leq k^2 \leq 1\f$. There are two sign conventions
0178   for elliptic integrals of the third kind. Some authors (Abramowitz
0179   and Stegun,
0180   <A HREF="http://mathworld.wolfram.com/EllipticIntegraloftheThirdKind.html">
0181   Mathworld</A>,
0182   <A HREF="http://www.open-std.org/jtc1/sc22/wg21/docs/papers/2004/n1687.pdf">
0183   C++ standard proposal</A>) use the above formula, while others
0184   (<A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC95">
0185   GSL</A>, <A HREF="http://planetmath.org/encyclopedia/EllipticIntegralsAndJacobiEllipticFunctions.html">
0186   Planetmath</A> and
0187   <A HREF="https://cern-tex.web.cern.ch/cern-tex/shortwrupsdir/c346/top.html">
0188   CERNLIB</A>) use the + sign in front of n in the denominator. In
0189   order to be C++ compliant, the present library uses the former
0190   convention. The implementation used is that of
0191   <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC101">GSL</A>.
0192 
0193   @ingroup SpecFunc
0194 
0195   */
0196   // [5.2.1.6] (complete) elliptic integral of the third kind
0197   double comp_ellint_3(double n, double k);
0198 
0199 
0200 
0201 
0202   /**
0203 
0204   Calculates the confluent hypergeometric functions of the first kind.
0205 
0206   \f[ _{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!}  \f]
0207 
0208   For detailed description see
0209   <A HREF="http://mathworld.wolfram.com/ConfluentHypergeometricFunctionoftheFirstKind.html">
0210   Mathworld</A>. The implementation used is that of
0211   <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC125">GSL</A>.
0212 
0213   @ingroup SpecFunc
0214 
0215   */
0216   // [5.2.1.7] confluent hypergeometric functions
0217 
0218   double conf_hyperg(double a, double b, double z);
0219 
0220 
0221   /**
0222 
0223   Calculates the confluent hypergeometric functions of the second kind, known also as Kummer function of the second kind,
0224   it is related to the confluent hypergeometric functions of the first kind.
0225 
0226   \f[ U(a,b,z)  = \frac{ \pi}{ \sin{\pi b } } \left[ \frac{ _{1}F_{1}(a,b,z) } {\Gamma(a-b+1) }
0227             - \frac{ z^{1-b} { _{1}F_{1}}(a-b+1,2-b,z)}{\Gamma(a)} \right]  \f]
0228 
0229   For detailed description see
0230   <A HREF="http://mathworld.wolfram.com/ConfluentHypergeometricFunctionoftheSecondKind.html">
0231   Mathworld</A>. The implementation used is that of
0232   <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC125">GSL</A>.
0233   This function is not part of the C++ standard proposal
0234 
0235   @ingroup SpecFunc
0236 
0237   */
0238   // confluent hypergeometric functions of second type
0239 
0240   double conf_hypergU(double a, double b, double z);
0241 
0242 
0243 
0244   /**
0245 
0246   Calculates the modified Bessel function of the first kind
0247   (also called regular modified (cylindrical) Bessel function).
0248 
0249   \f[ I_{\nu} (x) = i^{-\nu} J_{\nu}(ix) = \sum_{k=0}^{\infty} \frac{(\frac{1}{2}x)^{\nu + 2k}}{k! \Gamma(\nu + k + 1)} \f]
0250 
0251   for \f$x>0, \nu > 0\f$. For detailed description see
0252   <A HREF="http://mathworld.wolfram.com/ModifiedBesselFunctionoftheFirstKind.html">
0253   Mathworld</A>. The implementation used is that of
0254   <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC71">GSL</A>.
0255 
0256   @ingroup SpecFunc
0257 
0258   */
0259   // [5.2.1.8] regular modified cylindrical Bessel functions
0260 
0261   double cyl_bessel_i(double nu, double x);
0262 
0263 
0264 
0265 
0266   /**
0267 
0268   Calculates the (cylindrical) Bessel functions of the first kind (also
0269   called regular (cylindrical) Bessel functions).
0270 
0271   \f[ J_{\nu} (x) = \sum_{k=0}^{\infty} \frac{(-1)^k(\frac{1}{2}x)^{\nu + 2k}}{k! \Gamma(\nu + k + 1)} \f]
0272 
0273   For detailed description see
0274   <A HREF="http://mathworld.wolfram.com/BesselFunctionoftheFirstKind.html">
0275   Mathworld</A>. The implementation used is that of
0276   <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC69">GSL</A>.
0277 
0278   @ingroup SpecFunc
0279 
0280   */
0281   // [5.2.1.9] cylindrical Bessel functions (of the first kind)
0282 
0283   double cyl_bessel_j(double nu, double x);
0284 
0285 
0286 
0287 
0288 
0289   /**
0290 
0291   Calculates the modified Bessel functions of the second kind
0292   (also called irregular modified (cylindrical) Bessel functions).
0293 
0294   \f[ 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}}
0295 & \mbox{for integral $\nu$} \end{array}  \right.  \f]
0296 
0297   for \f$x>0, \nu > 0\f$. For detailed description see
0298   <A HREF="http://mathworld.wolfram.com/ModifiedBesselFunctionoftheSecondKind.html">
0299   Mathworld</A>. The implementation used is that of
0300   <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC72">GSL</A>.
0301 
0302   @ingroup SpecFunc
0303 
0304   */
0305   // [5.2.1.10] irregular modified cylindrical Bessel functions
0306 
0307   double cyl_bessel_k(double nu, double x);
0308 
0309 
0310 
0311 
0312   /**
0313 
0314   Calculates the (cylindrical) Bessel functions of the second kind
0315   (also called irregular (cylindrical) Bessel functions or
0316   (cylindrical) Neumann functions).
0317 
0318   \f[ 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.  \f]
0319 
0320    For detailed description see
0321   <A HREF="http://mathworld.wolfram.com/BesselFunctionoftheSecondKind.html">
0322   Mathworld</A>. The implementation used is that of
0323   <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC70">GSL</A>.
0324 
0325   @ingroup SpecFunc
0326 
0327   */
0328   // [5.2.1.11] cylindrical Neumann functions;
0329   // cylindrical Bessel functions (of the second kind)
0330 
0331   double cyl_neumann(double nu, double x);
0332 
0333 
0334 
0335 
0336   /**
0337 
0338   Calculates the incomplete elliptic integral of the first kind.
0339 
0340   \f[ F(k, \phi) = \int_{0}^{\phi} \frac{d \theta}{\sqrt{1 - k^2 \sin^2{\theta}}} \f]
0341 
0342   with \f$0 \leq k^2 \leq 1\f$. For detailed description see
0343   <A HREF="http://mathworld.wolfram.com/EllipticIntegraloftheFirstKind.html">
0344   Mathworld</A>. The implementation used is that of
0345   <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC101">GSL</A>.
0346 
0347   @param k
0348   @param phi angle in radians
0349 
0350   @ingroup SpecFunc
0351 
0352   */
0353   // [5.2.1.12] (incomplete) elliptic integral of the first kind
0354   // phi in radians
0355 
0356   double ellint_1(double k, double phi);
0357 
0358 
0359 
0360 
0361   /**
0362 
0363   Calculates the complete elliptic integral of the second kind.
0364 
0365   \f[ E(k , \phi) = \int_{0}^{\phi} \sqrt{1 - k^2 \sin^2{\theta}} d \theta  \f]
0366 
0367   with \f$0 \leq k^2 \leq 1\f$. For detailed description see
0368   <A HREF="http://mathworld.wolfram.com/EllipticIntegraloftheSecondKind.html">
0369   Mathworld</A>. The implementation used is that of
0370   <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC101">GSL</A>.
0371 
0372   @param k
0373   @param phi angle in radians
0374 
0375   @ingroup SpecFunc
0376 
0377   */
0378   // [5.2.1.13] (incomplete) elliptic integral of the second kind
0379   // phi in radians
0380 
0381   double ellint_2(double k, double phi);
0382 
0383 
0384 
0385 
0386   /**
0387 
0388   Calculates the complete elliptic integral of the third kind.
0389 
0390   \f[ \Pi (n, k, \phi) = \int_{0}^{\phi} \frac{d \theta}{(1 - n \sin^2{\theta})\sqrt{1 - k^2 \sin^2{\theta}}}  \f]
0391 
0392   with \f$0 \leq k^2 \leq 1\f$. There are two sign conventions
0393   for elliptic integrals of the third kind. Some authors (Abramowitz
0394   and Stegun,
0395   <A HREF="http://mathworld.wolfram.com/EllipticIntegraloftheThirdKind.html">
0396   Mathworld</A>,
0397   <A HREF="http://www.open-std.org/jtc1/sc22/wg21/docs/papers/2004/n1687.pdf">
0398   C++ standard proposal</A>) use the above formula, while others
0399   (<A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC95">
0400   GSL</A>, <A HREF="http://planetmath.org/encyclopedia/EllipticIntegralsAndJacobiEllipticFunctions.html">
0401   Planetmath</A> and
0402   <A HREF="https://cern-tex.web.cern.ch/cern-tex/shortwrupsdir/c346/top.html">
0403   CERNLIB</A>) use the + sign in front of n in the denominator. In
0404   order to be C++ compliant, the present library uses the former
0405   convention. The implementation used is that of
0406   <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC101">GSL</A>.
0407 
0408   @param n
0409   @param k
0410   @param phi angle in radians
0411 
0412   @ingroup SpecFunc
0413 
0414   */
0415   // [5.2.1.14] (incomplete) elliptic integral of the third kind
0416   // phi in radians
0417 
0418   double ellint_3(double n, double k, double phi);
0419 
0420 
0421 
0422 
0423   /**
0424 
0425   Calculates the exponential integral.
0426 
0427   \f[ Ei(x) = - \int_{-x}^{\infty} \frac{e^{-t}}{t} dt \f]
0428 
0429   For detailed description see
0430   <A HREF="http://mathworld.wolfram.com/ExponentialIntegral.html">
0431   Mathworld</A>. The implementation used is that of
0432   <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC115">GSL</A>.
0433 
0434   @ingroup SpecFunc
0435 
0436   */
0437   // [5.2.1.15] exponential integral
0438 
0439   double expint(double x);
0440   double expint_n(int n, double x);
0441 
0442 
0443 
0444   // [5.2.1.16] Hermite polynomials
0445 
0446   //double hermite(unsigned n, double x);
0447 
0448 
0449 
0450 
0451 
0452   /**
0453 
0454   Calculates Gauss' hypergeometric function.
0455 
0456   \f[ _{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!} \f]
0457 
0458   For detailed description see
0459   <A HREF="http://mathworld.wolfram.com/HypergeometricFunction.html">
0460   Mathworld</A>. The implementation used is that of
0461   <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC125">GSL</A>.
0462 
0463   @ingroup SpecFunc
0464 
0465   */
0466   // [5.2.1.17] hypergeometric functions
0467 
0468   double hyperg(double a, double b, double c, double x);
0469 
0470 
0471 
0472   /**
0473 
0474   Calculates the Laguerre polynomials
0475 
0476   \f[ P_{l}(x) = \frac{ e^x}{n!} \frac{d^n}{dx^n}  (x^n - e^{-x}) \f]
0477 
0478   for \f$x \geq 0 \f$ in the Rodrigues representation.
0479   They corresponds to the associated Laguerre polynomial of order m=0.
0480   See Abramowitz and Stegun, (22.5.16)
0481   For detailed description see
0482   <A HREF="http://mathworld.wolfram.com/LaguerrePolynomial.html">
0483   Mathworld</A>.
0484   The are implemented using the associated Laguerre polynomial of order m=0.
0485 
0486   @ingroup SpecFunc
0487 
0488   */
0489   // [5.2.1.18] Laguerre polynomials
0490 
0491   double laguerre(unsigned n, double x);
0492 
0493 
0494   /**
0495 
0496   Calculates the Lambert W function on branch 0
0497 
0498   The Lambert W functions are defined to be the solution of the equation
0499 
0500   \f[ W(x) \exp(W(x)) = x \f]
0501 
0502   For detailed description see
0503   <A HREF="https://mathworld.wolfram.com/LambertW-Function.html">Mathworld</A>
0504   or <A HREF="https://en.wikipedia.org/wiki/Lambert_W_function">Wikipedia</A>.
0505 
0506   This function implements the Lambert W function on branch 0, which is real
0507   valued and defined for \f$ x \geq -1/e \f$ with \f$ W_0(x) \geq -1 \f$.
0508 
0509   @ingroup SpecFunc
0510 
0511   */
0512 
0513   double lambert_W0(double x);
0514 
0515 
0516   /**
0517 
0518   Calculates the Lambert W function on branch -1
0519 
0520   The Lambert W functions are defined to be the solution of the equation
0521 
0522   \f[ W(x) \exp(W(x)) = x \f]
0523 
0524   For detailed description see
0525   <A HREF="https://mathworld.wolfram.com/LambertW-Function.html">Mathworld</A>
0526   or <A HREF="https://en.wikipedia.org/wiki/Lambert_W_function">Wikipedia</A>.
0527 
0528   This function implements the Lambert W function on branch -1, which is real
0529   valued and defined for \f$ -1/e \seq x < 0 \f$ with
0530   \f$ W_{-1}(x) \seq -1 \f$.
0531 
0532   @ingroup SpecFunc
0533 
0534   */
0535 
0536   double lambert_Wm1(double x);
0537 
0538 
0539   /**
0540 
0541   Calculates the Legendre polynomials.
0542 
0543   \f[ P_{l}(x) = \frac{1}{2^l l!} \frac{d^l}{dx^l}  (x^2 - 1)^l \f]
0544 
0545   for \f$l \geq 0, |x|\leq1\f$ in the Rodrigues representation.
0546   For detailed description see
0547   <A HREF="http://mathworld.wolfram.com/LegendrePolynomial.html">
0548   Mathworld</A>. The implementation used is that of
0549   <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC129">GSL</A>.
0550 
0551   @ingroup SpecFunc
0552 
0553   */
0554   // [5.2.1.19] Legendre polynomials
0555 
0556   double legendre(unsigned l, double x);
0557 
0558 
0559 
0560 
0561   /**
0562 
0563   Calculates the Riemann zeta function.
0564 
0565   \f[ \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. \f]
0566 
0567   For detailed description see
0568   <A HREF="http://mathworld.wolfram.com/RiemannZetaFunction.html">
0569   Mathworld</A>. The implementation used is that of
0570   <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC149">GSL</A>.
0571 
0572   CHECK WHETHER THE IMPLEMENTATION CALCULATES X<1
0573 
0574   @ingroup SpecFunc
0575 
0576   */
0577   // [5.2.1.20] Riemann zeta function
0578 
0579   double riemann_zeta(double x);
0580 
0581 
0582   /**
0583 
0584   Calculates the spherical Bessel functions of the first kind
0585   (also called regular spherical Bessel functions).
0586 
0587   \f[ j_{n}(x) = \sqrt{\frac{\pi}{2x}} J_{n+1/2}(x) \f]
0588 
0589   For detailed description see
0590   <A HREF="http://mathworld.wolfram.com/SphericalBesselFunctionoftheFirstKind.html">
0591   Mathworld</A>. The implementation used is that of
0592   <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC73">GSL</A>.
0593 
0594   @ingroup SpecFunc
0595 
0596   */
0597   // [5.2.1.21] spherical Bessel functions of the first kind
0598 
0599   double sph_bessel(unsigned n, double x);
0600 
0601 
0602   /**
0603 
0604   Computes the spherical (normalized)  associated Legendre polynomials,
0605   or spherical harmonic without azimuthal dependence (\f$e^(im\phi)\f$).
0606 
0607   \f[ Y_l^m(theta,0) = \sqrt{(2l+1)/(4\pi)} \sqrt{(l-m)!/(l+m)!} P_l^m(cos \theta) \f]
0608 
0609   for \f$m \geq 0, l \geq m\f$,
0610   where the Condon-Shortley phase  \f$(-1)^m\f$ is included in P_l^m(x)
0611   This function is consistent with both C++0x and GSL,
0612   even though there is a discrepancy in where to include the phase.
0613   There is no reference in Abramowitz and Stegun.
0614 
0615 
0616   @ingroup SpecFunc
0617 
0618   */
0619 
0620   // [5.2.1.22] spherical associated Legendre functions
0621 
0622   double sph_legendre(unsigned l, unsigned m, double theta);
0623 
0624 
0625   /**
0626 
0627   Calculates the spherical Bessel functions of the second kind
0628   (also called irregular spherical Bessel functions or
0629   spherical Neumann functions).
0630 
0631   \f[ n_n(x) = y_n(x) = \sqrt{\frac{\pi}{2x}} N_{n+1/2}(x) \f]
0632 
0633   For detailed description see
0634   <A HREF="http://mathworld.wolfram.com/SphericalBesselFunctionoftheSecondKind.html">
0635   Mathworld</A>. The implementation used is that of
0636   <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC74">GSL</A>.
0637 
0638   @ingroup SpecFunc
0639 
0640   */
0641   // [5.2.1.23] spherical Neumann functions
0642 
0643   double sph_neumann(unsigned n, double x);
0644 
0645   /**
0646 
0647   Calculates the Airy function Ai
0648 
0649   \f[ Ai(x) = \frac{1}{\pi} \int\limits_{0}^{\infty} \cos(xt + t^3/3) dt \f]
0650 
0651   For detailed description see
0652   <A HREF="http://mathworld.wolfram.com/AiryFunctions.html">
0653   Mathworld</A>
0654   and <A HREF="http://www.nrbook.com/abramowitz_and_stegun/page_446.htm">Abramowitz&Stegun, Sect. 10.4</A>.
0655   The implementation used is that of
0656   <A HREF="http://www.gnu.org/software/gsl/manual/html_node/Airy-Functions.html">GSL</A>.
0657 
0658   @ingroup SpecFunc
0659 
0660   */
0661   // Airy function Ai
0662 
0663   double airy_Ai(double x);
0664 
0665   /**
0666 
0667   Calculates the Airy function Bi
0668 
0669   \f[ Bi(x) = \frac{1}{\pi} \int\limits_{0}^{\infty} [\exp(xt-t^3/3) + \cos(xt + t^3/3)] dt \f]
0670 
0671   For detailed description see
0672   <A HREF="http://mathworld.wolfram.com/AiryFunctions.html">
0673   Mathworld</A>
0674   and <A HREF="http://www.nrbook.com/abramowitz_and_stegun/page_446.htm">Abramowitz&Stegun, Sect. 10.4</A>.
0675   The implementation used is that of
0676   <A HREF="http://www.gnu.org/software/gsl/manual/html_node/Airy-Functions.html">GSL</A>.
0677 
0678   @ingroup SpecFunc
0679 
0680   */
0681   // Airy function Bi
0682 
0683   double airy_Bi(double x);
0684 
0685   /**
0686 
0687   Calculates the derivative of the Airy function Ai
0688 
0689   \f[ Ai(x) = \frac{1}{\pi} \int\limits_{0}^{\infty} \cos(xt + t^3/3) dt \f]
0690 
0691   For detailed description see
0692   <A HREF="http://mathworld.wolfram.com/AiryFunctions.html">
0693   Mathworld</A>
0694   and <A HREF="http://www.nrbook.com/abramowitz_and_stegun/page_446.htm">Abramowitz&Stegun, Sect. 10.4</A>.
0695   The implementation used is that of
0696   <A HREF="http://www.gnu.org/software/gsl/manual/html_node/Derivatives-of-Airy-Functions.html">GSL</A>.
0697 
0698   @ingroup SpecFunc
0699 
0700   */
0701   // Derivative of the Airy function Ai
0702 
0703   double airy_Ai_deriv(double x);
0704 
0705   /**
0706 
0707   Calculates the derivative of the Airy function Bi
0708 
0709   \f[ Bi(x) = \frac{1}{\pi} \int\limits_{0}^{\infty} [\exp(xt-t^3/3) + \cos(xt + t^3/3)] dt \f]
0710 
0711   For detailed description see
0712   <A HREF="http://mathworld.wolfram.com/AiryFunctions.html">
0713   Mathworld</A>
0714   and <A HREF="http://www.nrbook.com/abramowitz_and_stegun/page_446.htm">Abramowitz&Stegun, Sect. 10.4</A>.
0715   The implementation used is that of
0716   <A HREF="http://www.gnu.org/software/gsl/manual/html_node/Derivatives-of-Airy-Functions.html">GSL</A>.
0717 
0718   @ingroup SpecFunc
0719 
0720   */
0721   // Derivative of the Airy function Bi
0722 
0723   double airy_Bi_deriv(double x);
0724 
0725   /**
0726 
0727   Calculates the zeroes of the Airy function Ai
0728 
0729   \f[ Ai(x) = \frac{1}{\pi} \int\limits_{0}^{\infty} \cos(xt + t^3/3) dt \f]
0730 
0731   For detailed description see
0732   <A HREF="http://mathworld.wolfram.com/AiryFunctionZeros.html">
0733   Mathworld</A>
0734   and <A HREF="http://www.nrbook.com/abramowitz_and_stegun/page_446.htm">Abramowitz&Stegun, Sect. 10.4</A>.
0735   The implementation used is that of
0736   <A HREF="http://www.gnu.org/software/gsl/manual/html_node/Zeros-of-Airy-Functions.html">GSL</A>.
0737 
0738   @ingroup SpecFunc
0739 
0740   */
0741   // s-th zero of the Airy function Ai
0742 
0743   double airy_zero_Ai(unsigned int s);
0744 
0745   /**
0746 
0747   Calculates the zeroes of the Airy function Bi
0748 
0749   \f[ Bi(x) = \frac{1}{\pi} \int\limits_{0}^{\infty} [\exp(xt-t^3/3) + \cos(xt + t^3/3)] dt \f]
0750 
0751   For detailed description see
0752   <A HREF="http://mathworld.wolfram.com/AiryFunctionZeros.html">
0753   Mathworld</A>
0754   and <A HREF="http://www.nrbook.com/abramowitz_and_stegun/page_446.htm">Abramowitz&Stegun, Sect. 10.4</A>.
0755   The implementation used is that of
0756   <A HREF="http://www.gnu.org/software/gsl/manual/html_node/Zeros-of-Airy-Functions.html">GSL</A>.
0757 
0758   @ingroup SpecFunc
0759 
0760   */
0761   // s-th zero of the Airy function Bi
0762 
0763   double airy_zero_Bi(unsigned int s);
0764 
0765   /**
0766 
0767   Calculates the zeroes of the derivative of the Airy function Ai
0768 
0769   \f[ Ai(x) = \frac{1}{\pi} \int\limits_{0}^{\infty} \cos(xt + t^3/3) dt \f]
0770 
0771   For detailed description see
0772   <A HREF="http://mathworld.wolfram.com/AiryFunctionZeros.html">
0773   Mathworld</A>
0774   and <A HREF="http://www.nrbook.com/abramowitz_and_stegun/page_446.htm">Abramowitz&Stegun, Sect. 10.4</A>.
0775   The implementation used is that of
0776   <A HREF="http://www.gnu.org/software/gsl/manual/html_node/Zeros-of-Derivatives-of-Airy-Functions.html">GSL</A>.
0777 
0778   @ingroup SpecFunc
0779 
0780   */
0781   // s-th zero of the derivative of the Airy function Ai
0782 
0783   double airy_zero_Ai_deriv(unsigned int s);
0784 
0785   /**
0786 
0787   Calculates the zeroes of the derivative of the Airy function Bi
0788 
0789   \f[ Bi(x) = \frac{1}{\pi} \int\limits_{0}^{\infty} [\exp(xt-t^3/3) + \cos(xt + t^3/3)] dt \f]
0790 
0791   For detailed description see
0792   <A HREF="http://mathworld.wolfram.com/AiryFunctionZeros.html">
0793   Mathworld</A>
0794   and <A HREF="http://www.nrbook.com/abramowitz_and_stegun/page_446.htm">Abramowitz&Stegun, Sect. 10.4</A>.
0795   The implementation used is that of
0796   <A HREF="http://www.gnu.org/software/gsl/manual/html_node/Zeros-of-Derivatives-of-Airy-Functions.html">GSL</A>.
0797 
0798   @ingroup SpecFunc
0799 
0800   */
0801   // s-th zero of the derivative of the Airy function Bi
0802 
0803   double airy_zero_Bi_deriv(unsigned int s);
0804 
0805   /**
0806 
0807    Calculates the Wigner 3j coupling coefficients
0808 
0809            (ja jb jc
0810            ma mb mc)
0811 
0812    where ja,ma,...etc are integers or half integers.
0813    The function takes as input arguments only integers which corresponds
0814    to half integer units, e.g two_ja = 2 * ja
0815 
0816    For detailed description see
0817    <A HREF="http://mathworld.wolfram.com/Wigner3j-Symbol.html.html">
0818    Mathworld</A>.
0819    The implementation used is that of
0820    <A HREF="http://www.gnu.org/software/gsl/manual/html_node/3_002dj-Symbols.html#g_t3_002dj-Symbols">GSL</A>.
0821 
0822    @ingroup SpecFunc
0823 
0824   */
0825 
0826   double wigner_3j(int two_ja, int two_jb, int two_jc, int two_ma, int two_mb, int two_mc);
0827 
0828   /**
0829 
0830   Calculates the Wigner 6j coupling coefficients
0831 
0832           (ja jb jc
0833            jd je jf)
0834 
0835    where ja,jb,...etc are integers or half integers.
0836    The function takes as input arguments only integers which corresponds
0837    to half integer units, e.g two_ja = 2 * ja
0838 
0839    For detailed description see
0840    <A HREF="http://mathworld.wolfram.com/Wigner6j-Symbol.html">
0841    Mathworld</A>.
0842    The implementation used is that of
0843    <A HREF="http://www.gnu.org/software/gsl/manual/html_node/6_002dj-Symbols.html#g_t6_002dj-Symbols">GSL</A>.
0844 
0845    @ingroup SpecFunc
0846 
0847   */
0848 
0849   double wigner_6j(int two_ja, int two_jb, int two_jc, int two_jd, int two_je, int two_jf);
0850 
0851   /**
0852 
0853   Calculates the Wigner 9j coupling coefficients
0854 
0855           (ja jb jc
0856            jd je jf
0857            jg jh ji)
0858 
0859    where ja,jb...etc are integers or half integers.
0860    The function takes as input arguments only integers which corresponds
0861    to half integer units, e.g two_ja = 2 * ja
0862 
0863 
0864    For detailed description see
0865    <A HREF="http://mathworld.wolfram.com/Wigner9j-Symbol.html">
0866    Mathworld</A>.
0867    The implementation used is that of
0868    <A HREF="http://www.gnu.org/software/gsl/manual/html_node/9_002dj-Symbols.html#g_t9_002dj-Symbols">GSL</A>.
0869 
0870    @ingroup SpecFunc
0871 
0872   */
0873 
0874    double wigner_9j(int two_ja, int two_jb, int two_jc, int two_jd, int two_je, int two_jf, int two_jg, int two_jh, int two_ji);
0875 
0876 
0877 
0878 } // namespace Math
0879 } // namespace ROOT
0880 
0881 
0882 #endif //ROOT_Math_SpecFuncMathMore