SpecFuncMathMore.h 25 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848
  1. // @(#)root/mathmore:$Id$
  2. // Authors: L. Moneta, A. Zsenei 08/2005
  3. // Authors: Andras Zsenei & Lorenzo Moneta 06/2005
  4. /**********************************************************************
  5. * *
  6. * Copyright (c) 2004 ROOT Foundation, CERN/PH-SFT *
  7. * *
  8. * This library is free software; you can redistribute it and/or *
  9. * modify it under the terms of the GNU General Public License *
  10. * as published by the Free Software Foundation; either version 2 *
  11. * of the License, or (at your option) any later version. *
  12. * *
  13. * This library is distributed in the hope that it will be useful, *
  14. * but WITHOUT ANY WARRANTY; without even the implied warranty of *
  15. * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU *
  16. * General Public License for more details. *
  17. * *
  18. * You should have received a copy of the GNU General Public License *
  19. * along with this library (see file COPYING); if not, write *
  20. * to the Free Software Foundation, Inc., 59 Temple Place, Suite *
  21. * 330, Boston, MA 02111-1307 USA, or contact the author. *
  22. * *
  23. **********************************************************************/
  24. #if defined(__CINT__) && !defined(__MAKECINT__)
  25. // avoid to include header file when using CINT
  26. #ifndef _WIN32
  27. #include "../lib/libMathMore.so"
  28. #else
  29. #include "../bin/libMathMore.dll"
  30. #endif
  31. #else
  32. /**
  33. Special mathematical functions.
  34. The naming and numbering of the functions is taken from
  35. <A HREF="http://www.open-std.org/jtc1/sc22/wg21/docs/papers/2004/n1687.pdf">
  36. Matt Austern,
  37. (Draft) Technical Report on Standard Library Extensions,
  38. N1687=04-0127, September 10, 2004</A>
  39. @author Created by Andras Zsenei on Mon Nov 8 2004
  40. @defgroup SpecFunc Special functions
  41. */
  42. #ifndef ROOT_Math_SpecFuncMathMore
  43. #define ROOT_Math_SpecFuncMathMore
  44. namespace ROOT {
  45. namespace Math {
  46. /** @name Special Functions from MathMore */
  47. /**
  48. Computes the generalized Laguerre polynomials for
  49. \f$ n \geq 0, m > -1 \f$.
  50. They are defined in terms of the confluent hypergeometric function.
  51. For integer values of m they can be defined in terms of the Laguerre polynomials \f$L_n(x)\f$:
  52. \f[ L_{n}^{m}(x) = (-1)^{m} \frac{d^m}{dx^m} L_{n+m}(x) \f]
  53. For detailed description see
  54. <A HREF="http://mathworld.wolfram.com/LaguerrePolynomial.html">Mathworld</A>.
  55. The implementation used is that of
  56. <A HREF="http://www.gnu.org/software/gsl/manual/html_node/Laguerre-Functions.html">GSL</A>.
  57. This function is an extension of C++0x, also consistent in GSL,
  58. Abramowitz and Stegun 1972 and MatheMathica that uses non-integer values for m.
  59. C++0x calls for 'int m', more restrictive than necessary.
  60. The definition for was incorrect in 'n1687.pdf', but fixed in
  61. <A HREF="http://www.open-std.org/jtc1/sc22/wg21/docs/papers/2005/n1836.pdf">n1836.pdf</A>,
  62. the most recent draft as of 2007-07-01
  63. @ingroup SpecFunc
  64. */
  65. // [5.2.1.1] associated Laguerre polynomials
  66. double assoc_laguerre(unsigned n, double m, double x);
  67. /**
  68. Computes the associated Legendre polynomials.
  69. \f[ P_{l}^{m}(x) = (1-x^2)^{m/2} \frac{d^m}{dx^m} P_{l}(x) \f]
  70. with \f$m \geq 0\f$, \f$ l \geq m \f$ and \f$ |x|<1 \f$.
  71. There are two sign conventions for associated Legendre polynomials.
  72. As is the case with the above formula, some authors (e.g., Arfken
  73. 1985, pp. 668-669) omit the Condon-Shortley phase \f$(-1)^m\f$,
  74. while others include it (e.g., Abramowitz and Stegun 1972).
  75. One possible way to distinguish the two conventions is due to
  76. Abramowitz and Stegun (1972, p. 332), who use the notation
  77. \f[ P_{lm} (x) = (-1)^m P_{l}^{m} (x)\f]
  78. to distinguish the two. For detailed description see
  79. <A HREF="http://mathworld.wolfram.com/LegendrePolynomial.html">
  80. Mathworld</A>. The implementation used is that of
  81. <A HREF="http://www.gnu.org/software/gsl/manual/html_node/Associated-Legendre-Polynomials-and-Spherical-Harmonics.html">GSL</A>.
  82. 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,
  83. sph_legendre(l,m,theta).
  84. @ingroup SpecFunc
  85. */
  86. // [5.2.1.2] associated Legendre functions
  87. double assoc_legendre(unsigned l, unsigned m, double x);
  88. /**
  89. Calculates the complete elliptic integral of the first kind.
  90. \f[ K(k) = F(k, \pi / 2) = \int_{0}^{\pi /2} \frac{d \theta}{\sqrt{1 - k^2 \sin^2{\theta}}} \f]
  91. with \f$0 \leq k^2 \leq 1\f$. For detailed description see
  92. <A HREF="http://mathworld.wolfram.com/CompleteEllipticIntegraloftheFirstKind.html">
  93. Mathworld</A>. The implementation used is that of
  94. <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC100">GSL</A>.
  95. @ingroup SpecFunc
  96. */
  97. // [5.2.1.4] (complete) elliptic integral of the first kind
  98. double comp_ellint_1(double k);
  99. /**
  100. Calculates the complete elliptic integral of the second kind.
  101. \f[ E(k) = E(k , \pi / 2) = \int_{0}^{\pi /2} \sqrt{1 - k^2 \sin^2{\theta}} d \theta \f]
  102. with \f$0 \leq k^2 \leq 1\f$. For detailed description see
  103. <A HREF="http://mathworld.wolfram.com/CompleteEllipticIntegraloftheSecondKind.html">
  104. Mathworld</A>. The implementation used is that of
  105. <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC100">GSL</A>.
  106. @ingroup SpecFunc
  107. */
  108. // [5.2.1.5] (complete) elliptic integral of the second kind
  109. double comp_ellint_2(double k);
  110. /**
  111. Calculates the complete elliptic integral of the third kind.
  112. \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]
  113. with \f$0 \leq k^2 \leq 1\f$. There are two sign conventions
  114. for elliptic integrals of the third kind. Some authors (Abramowitz
  115. and Stegun,
  116. <A HREF="http://mathworld.wolfram.com/EllipticIntegraloftheThirdKind.html">
  117. Mathworld</A>,
  118. <A HREF="http://www.open-std.org/jtc1/sc22/wg21/docs/papers/2004/n1687.pdf">
  119. C++ standard proposal</A>) use the above formula, while others
  120. (<A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC95">
  121. GSL</A>, <A HREF="http://planetmath.org/encyclopedia/EllipticIntegralsAndJacobiEllipticFunctions.html">
  122. Planetmath</A> and
  123. <A HREF="http://wwwasdoc.web.cern.ch/wwwasdoc/shortwrupsdir/c346/top.html">
  124. CERNLIB</A>) use the + sign in front of n in the denominator. In
  125. order to be C++ compliant, the present library uses the former
  126. convention. The implementation used is that of
  127. <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC101">GSL</A>.
  128. @ingroup SpecFunc
  129. */
  130. // [5.2.1.6] (complete) elliptic integral of the third kind
  131. double comp_ellint_3(double n, double k);
  132. /**
  133. Calculates the confluent hypergeometric functions of the first kind.
  134. \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]
  135. For detailed description see
  136. <A HREF="http://mathworld.wolfram.com/ConfluentHypergeometricFunctionoftheFirstKind.html">
  137. Mathworld</A>. The implementation used is that of
  138. <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC125">GSL</A>.
  139. @ingroup SpecFunc
  140. */
  141. // [5.2.1.7] confluent hypergeometric functions
  142. double conf_hyperg(double a, double b, double z);
  143. /**
  144. Calculates the confluent hypergeometric functions of the second kind, known also as Kummer function of the second kind,
  145. it is related to the confluent hypergeometric functions of the first kind.
  146. \f[ U(a,b,z) = \frac{ \pi}{ \sin{\pi b } } \left[ \frac{ _{1}F_{1}(a,b,z) } {\Gamma(a-b+1) }
  147. - \frac{ z^{1-b} { _{1}F_{1}}(a-b+1,2-b,z)}{\Gamma(a)} \right] \f]
  148. For detailed description see
  149. <A HREF="http://mathworld.wolfram.com/ConfluentHypergeometricFunctionoftheSecondKind.html">
  150. Mathworld</A>. The implementation used is that of
  151. <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC125">GSL</A>.
  152. This function is not part of the C++ standard proposal
  153. @ingroup SpecFunc
  154. */
  155. // confluent hypergeometric functions of second type
  156. double conf_hypergU(double a, double b, double z);
  157. /**
  158. Calculates the modified Bessel function of the first kind
  159. (also called regular modified (cylindrical) Bessel function).
  160. \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]
  161. for \f$x>0, \nu > 0\f$. For detailed description see
  162. <A HREF="http://mathworld.wolfram.com/ModifiedBesselFunctionoftheFirstKind.html">
  163. Mathworld</A>. The implementation used is that of
  164. <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC71">GSL</A>.
  165. @ingroup SpecFunc
  166. */
  167. // [5.2.1.8] regular modified cylindrical Bessel functions
  168. double cyl_bessel_i(double nu, double x);
  169. /**
  170. Calculates the (cylindrical) Bessel functions of the first kind (also
  171. called regular (cylindrical) Bessel functions).
  172. \f[ J_{\nu} (x) = \sum_{k=0}^{\infty} \frac{(-1)^k(\frac{1}{2}x)^{\nu + 2k}}{k! \Gamma(\nu + k + 1)} \f]
  173. For detailed description see
  174. <A HREF="http://mathworld.wolfram.com/BesselFunctionoftheFirstKind.html">
  175. Mathworld</A>. The implementation used is that of
  176. <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC69">GSL</A>.
  177. @ingroup SpecFunc
  178. */
  179. // [5.2.1.9] cylindrical Bessel functions (of the first kind)
  180. double cyl_bessel_j(double nu, double x);
  181. /**
  182. Calculates the modified Bessel functions of the second kind
  183. (also called irregular modified (cylindrical) Bessel functions).
  184. \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}}
  185. & \mbox{for integral $\nu$} \end{array} \right. \f]
  186. for \f$x>0, \nu > 0\f$. For detailed description see
  187. <A HREF="http://mathworld.wolfram.com/ModifiedBesselFunctionoftheSecondKind.html">
  188. Mathworld</A>. The implementation used is that of
  189. <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC72">GSL</A>.
  190. @ingroup SpecFunc
  191. */
  192. // [5.2.1.10] irregular modified cylindrical Bessel functions
  193. double cyl_bessel_k(double nu, double x);
  194. /**
  195. Calculates the (cylindrical) Bessel functions of the second kind
  196. (also called irregular (cylindrical) Bessel functions or
  197. (cylindrical) Neumann functions).
  198. \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]
  199. For detailed description see
  200. <A HREF="http://mathworld.wolfram.com/BesselFunctionoftheSecondKind.html">
  201. Mathworld</A>. The implementation used is that of
  202. <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC70">GSL</A>.
  203. @ingroup SpecFunc
  204. */
  205. // [5.2.1.11] cylindrical Neumann functions;
  206. // cylindrical Bessel functions (of the second kind)
  207. double cyl_neumann(double nu, double x);
  208. /**
  209. Calculates the incomplete elliptic integral of the first kind.
  210. \f[ F(k, \phi) = \int_{0}^{\phi} \frac{d \theta}{\sqrt{1 - k^2 \sin^2{\theta}}} \f]
  211. with \f$0 \leq k^2 \leq 1\f$. For detailed description see
  212. <A HREF="http://mathworld.wolfram.com/EllipticIntegraloftheFirstKind.html">
  213. Mathworld</A>. The implementation used is that of
  214. <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC101">GSL</A>.
  215. @param k
  216. @param phi angle in radians
  217. @ingroup SpecFunc
  218. */
  219. // [5.2.1.12] (incomplete) elliptic integral of the first kind
  220. // phi in radians
  221. double ellint_1(double k, double phi);
  222. /**
  223. Calculates the complete elliptic integral of the second kind.
  224. \f[ E(k , \phi) = \int_{0}^{\phi} \sqrt{1 - k^2 \sin^2{\theta}} d \theta \f]
  225. with \f$0 \leq k^2 \leq 1\f$. For detailed description see
  226. <A HREF="http://mathworld.wolfram.com/EllipticIntegraloftheSecondKind.html">
  227. Mathworld</A>. The implementation used is that of
  228. <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC101">GSL</A>.
  229. @param k
  230. @param phi angle in radians
  231. @ingroup SpecFunc
  232. */
  233. // [5.2.1.13] (incomplete) elliptic integral of the second kind
  234. // phi in radians
  235. double ellint_2(double k, double phi);
  236. /**
  237. Calculates the complete elliptic integral of the third kind.
  238. \f[ \Pi (n, k, \phi) = \int_{0}^{\phi} \frac{d \theta}{(1 - n \sin^2{\theta})\sqrt{1 - k^2 \sin^2{\theta}}} \f]
  239. with \f$0 \leq k^2 \leq 1\f$. There are two sign conventions
  240. for elliptic integrals of the third kind. Some authors (Abramowitz
  241. and Stegun,
  242. <A HREF="http://mathworld.wolfram.com/EllipticIntegraloftheThirdKind.html">
  243. Mathworld</A>,
  244. <A HREF="http://www.open-std.org/jtc1/sc22/wg21/docs/papers/2004/n1687.pdf">
  245. C++ standard proposal</A>) use the above formula, while others
  246. (<A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC95">
  247. GSL</A>, <A HREF="http://planetmath.org/encyclopedia/EllipticIntegralsAndJacobiEllipticFunctions.html">
  248. Planetmath</A> and
  249. <A HREF="http://wwwasdoc.web.cern.ch/wwwasdoc/shortwrupsdir/c346/top.html">
  250. CERNLIB</A>) use the + sign in front of n in the denominator. In
  251. order to be C++ compliant, the present library uses the former
  252. convention. The implementation used is that of
  253. <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC101">GSL</A>.
  254. @param n
  255. @param k
  256. @param phi angle in radians
  257. @ingroup SpecFunc
  258. */
  259. // [5.2.1.14] (incomplete) elliptic integral of the third kind
  260. // phi in radians
  261. double ellint_3(double n, double k, double phi);
  262. /**
  263. Calculates the exponential integral.
  264. \f[ Ei(x) = - \int_{-x}^{\infty} \frac{e^{-t}}{t} dt \f]
  265. For detailed description see
  266. <A HREF="http://mathworld.wolfram.com/ExponentialIntegral.html">
  267. Mathworld</A>. The implementation used is that of
  268. <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC115">GSL</A>.
  269. @ingroup SpecFunc
  270. */
  271. // [5.2.1.15] exponential integral
  272. double expint(double x);
  273. // [5.2.1.16] Hermite polynomials
  274. //double hermite(unsigned n, double x);
  275. /**
  276. Calculates Gauss' hypergeometric function.
  277. \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]
  278. For detailed description see
  279. <A HREF="http://mathworld.wolfram.com/HypergeometricFunction.html">
  280. Mathworld</A>. The implementation used is that of
  281. <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC125">GSL</A>.
  282. @ingroup SpecFunc
  283. */
  284. // [5.2.1.17] hypergeometric functions
  285. double hyperg(double a, double b, double c, double x);
  286. /**
  287. Calculates the Laguerre polynomials
  288. \f[ P_{l}(x) = \frac{ e^x}{n!} \frac{d^n}{dx^n} (x^n - e^{-x}) \f]
  289. for \f$x \geq 0 \f$ in the Rodrigues representation.
  290. They corresponds to the associated Laguerre polynomial of order m=0.
  291. See Abramowitz and Stegun, (22.5.16)
  292. For detailed description see
  293. <A HREF="http://mathworld.wolfram.com/LaguerrePolynomial.html">
  294. Mathworld</A>.
  295. The are implemented using the associated Laguerre polynomial of order m=0.
  296. @ingroup SpecFunc
  297. */
  298. // [5.2.1.18] Laguerre polynomials
  299. double laguerre(unsigned n, double x);
  300. /**
  301. Calculates the Legendre polynomials.
  302. \f[ P_{l}(x) = \frac{1}{2^l l!} \frac{d^l}{dx^l} (x^2 - 1)^l \f]
  303. for \f$l \geq 0, |x|\leq1\f$ in the Rodrigues representation.
  304. For detailed description see
  305. <A HREF="http://mathworld.wolfram.com/LegendrePolynomial.html">
  306. Mathworld</A>. The implementation used is that of
  307. <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC129">GSL</A>.
  308. @ingroup SpecFunc
  309. */
  310. // [5.2.1.19] Legendre polynomials
  311. double legendre(unsigned l, double x);
  312. /**
  313. Calculates the Riemann zeta function.
  314. \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]
  315. For detailed description see
  316. <A HREF="http://mathworld.wolfram.com/RiemannZetaFunction.html">
  317. Mathworld</A>. The implementation used is that of
  318. <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC149">GSL</A>.
  319. CHECK WHETHER THE IMPLEMENTATION CALCULATES X<1
  320. @ingroup SpecFunc
  321. */
  322. // [5.2.1.20] Riemann zeta function
  323. double riemann_zeta(double x);
  324. /**
  325. Calculates the spherical Bessel functions of the first kind
  326. (also called regular spherical Bessel functions).
  327. \f[ j_{n}(x) = \sqrt{\frac{\pi}{2x}} J_{n+1/2}(x) \f]
  328. For detailed description see
  329. <A HREF="http://mathworld.wolfram.com/SphericalBesselFunctionoftheFirstKind.html">
  330. Mathworld</A>. The implementation used is that of
  331. <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC73">GSL</A>.
  332. @ingroup SpecFunc
  333. */
  334. // [5.2.1.21] spherical Bessel functions of the first kind
  335. double sph_bessel(unsigned n, double x);
  336. /**
  337. Computes the spherical (normalized) associated Legendre polynomials,
  338. or spherical harmonic without azimuthal dependence (\f$e^(im\phi)\f$).
  339. \f[ Y_l^m(theta,0) = \sqrt{(2l+1)/(4\pi)} \sqrt{(l-m)!/(l+m)!} P_l^m(cos \theta) \f]
  340. for \f$m \geq 0, l \geq m\f$,
  341. where the Condon-Shortley phase \f$(-1)^m\f$ is included in P_l^m(x)
  342. This function is consistent with both C++0x and GSL,
  343. even though there is a discrepancy in where to include the phase.
  344. There is no reference in Abramowitz and Stegun.
  345. @ingroup SpecFunc
  346. */
  347. // [5.2.1.22] spherical associated Legendre functions
  348. double sph_legendre(unsigned l, unsigned m, double theta);
  349. /**
  350. Calculates the spherical Bessel functions of the second kind
  351. (also called irregular spherical Bessel functions or
  352. spherical Neumann functions).
  353. \f[ n_n(x) = y_n(x) = \sqrt{\frac{\pi}{2x}} N_{n+1/2}(x) \f]
  354. For detailed description see
  355. <A HREF="http://mathworld.wolfram.com/SphericalBesselFunctionoftheSecondKind.html">
  356. Mathworld</A>. The implementation used is that of
  357. <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC74">GSL</A>.
  358. @ingroup SpecFunc
  359. */
  360. // [5.2.1.23] spherical Neumann functions
  361. double sph_neumann(unsigned n, double x);
  362. /**
  363. Calculates the Airy function Ai
  364. \f[ Ai(x) = \frac{1}{\pi} \int\limits_{0}^{\infty} \cos(xt + t^3/3) dt \f]
  365. For detailed description see
  366. <A HREF="http://mathworld.wolfram.com/AiryFunctions.html">
  367. Mathworld</A>
  368. and <A HREF="http://www.nrbook.com/abramowitz_and_stegun/page_446.htm">Abramowitz&Stegun, Sect. 10.4</A>.
  369. The implementation used is that of
  370. <A HREF="http://www.gnu.org/software/gsl/manual/html_node/Airy-Functions.html">GSL</A>.
  371. @ingroup SpecFunc
  372. */
  373. // Airy function Ai
  374. double airy_Ai(double x);
  375. /**
  376. Calculates the Airy function Bi
  377. \f[ Bi(x) = \frac{1}{\pi} \int\limits_{0}^{\infty} [\exp(xt-t^3/3) + \cos(xt + t^3/3)] dt \f]
  378. For detailed description see
  379. <A HREF="http://mathworld.wolfram.com/AiryFunctions.html">
  380. Mathworld</A>
  381. and <A HREF="http://www.nrbook.com/abramowitz_and_stegun/page_446.htm">Abramowitz&Stegun, Sect. 10.4</A>.
  382. The implementation used is that of
  383. <A HREF="http://www.gnu.org/software/gsl/manual/html_node/Airy-Functions.html">GSL</A>.
  384. @ingroup SpecFunc
  385. */
  386. // Airy function Bi
  387. double airy_Bi(double x);
  388. /**
  389. Calculates the derivative of the Airy function Ai
  390. \f[ Ai(x) = \frac{1}{\pi} \int\limits_{0}^{\infty} \cos(xt + t^3/3) dt \f]
  391. For detailed description see
  392. <A HREF="http://mathworld.wolfram.com/AiryFunctions.html">
  393. Mathworld</A>
  394. and <A HREF="http://www.nrbook.com/abramowitz_and_stegun/page_446.htm">Abramowitz&Stegun, Sect. 10.4</A>.
  395. The implementation used is that of
  396. <A HREF="http://www.gnu.org/software/gsl/manual/html_node/Derivatives-of-Airy-Functions.html">GSL</A>.
  397. @ingroup SpecFunc
  398. */
  399. // Derivative of the Airy function Ai
  400. double airy_Ai_deriv(double x);
  401. /**
  402. Calculates the derivative of the Airy function Bi
  403. \f[ Bi(x) = \frac{1}{\pi} \int\limits_{0}^{\infty} [\exp(xt-t^3/3) + \cos(xt + t^3/3)] dt \f]
  404. For detailed description see
  405. <A HREF="http://mathworld.wolfram.com/AiryFunctions.html">
  406. Mathworld</A>
  407. and <A HREF="http://www.nrbook.com/abramowitz_and_stegun/page_446.htm">Abramowitz&Stegun, Sect. 10.4</A>.
  408. The implementation used is that of
  409. <A HREF="http://www.gnu.org/software/gsl/manual/html_node/Derivatives-of-Airy-Functions.html">GSL</A>.
  410. @ingroup SpecFunc
  411. */
  412. // Derivative of the Airy function Bi
  413. double airy_Bi_deriv(double x);
  414. /**
  415. Calculates the zeroes of the Airy function Ai
  416. \f[ Ai(x) = \frac{1}{\pi} \int\limits_{0}^{\infty} \cos(xt + t^3/3) dt \f]
  417. For detailed description see
  418. <A HREF="http://mathworld.wolfram.com/AiryFunctionZeros.html">
  419. Mathworld</A>
  420. and <A HREF="http://www.nrbook.com/abramowitz_and_stegun/page_446.htm">Abramowitz&Stegun, Sect. 10.4</A>.
  421. The implementation used is that of
  422. <A HREF="http://www.gnu.org/software/gsl/manual/html_node/Zeros-of-Airy-Functions.html">GSL</A>.
  423. @ingroup SpecFunc
  424. */
  425. // s-th zero of the Airy function Ai
  426. double airy_zero_Ai(unsigned int s);
  427. /**
  428. Calculates the zeroes of the Airy function Bi
  429. \f[ Bi(x) = \frac{1}{\pi} \int\limits_{0}^{\infty} [\exp(xt-t^3/3) + \cos(xt + t^3/3)] dt \f]
  430. For detailed description see
  431. <A HREF="http://mathworld.wolfram.com/AiryFunctionZeros.html">
  432. Mathworld</A>
  433. and <A HREF="http://www.nrbook.com/abramowitz_and_stegun/page_446.htm">Abramowitz&Stegun, Sect. 10.4</A>.
  434. The implementation used is that of
  435. <A HREF="http://www.gnu.org/software/gsl/manual/html_node/Zeros-of-Airy-Functions.html">GSL</A>.
  436. @ingroup SpecFunc
  437. */
  438. // s-th zero of the Airy function Bi
  439. double airy_zero_Bi(unsigned int s);
  440. /**
  441. Calculates the zeroes of the derivative of the Airy function Ai
  442. \f[ Ai(x) = \frac{1}{\pi} \int\limits_{0}^{\infty} \cos(xt + t^3/3) dt \f]
  443. For detailed description see
  444. <A HREF="http://mathworld.wolfram.com/AiryFunctionZeros.html">
  445. Mathworld</A>
  446. and <A HREF="http://www.nrbook.com/abramowitz_and_stegun/page_446.htm">Abramowitz&Stegun, Sect. 10.4</A>.
  447. The implementation used is that of
  448. <A HREF="http://www.gnu.org/software/gsl/manual/html_node/Zeros-of-Derivatives-of-Airy-Functions.html">GSL</A>.
  449. @ingroup SpecFunc
  450. */
  451. // s-th zero of the derivative of the Airy function Ai
  452. double airy_zero_Ai_deriv(unsigned int s);
  453. /**
  454. Calculates the zeroes of the derivative of the Airy function Bi
  455. \f[ Bi(x) = \frac{1}{\pi} \int\limits_{0}^{\infty} [\exp(xt-t^3/3) + \cos(xt + t^3/3)] dt \f]
  456. For detailed description see
  457. <A HREF="http://mathworld.wolfram.com/AiryFunctionZeros.html">
  458. Mathworld</A>
  459. and <A HREF="http://www.nrbook.com/abramowitz_and_stegun/page_446.htm">Abramowitz&Stegun, Sect. 10.4</A>.
  460. The implementation used is that of
  461. <A HREF="http://www.gnu.org/software/gsl/manual/html_node/Zeros-of-Derivatives-of-Airy-Functions.html">GSL</A>.
  462. @ingroup SpecFunc
  463. */
  464. // s-th zero of the derivative of the Airy function Bi
  465. double airy_zero_Bi_deriv(unsigned int s);
  466. /**
  467. Calculates the Wigner 3j coupling coefficients
  468. (ja jb jc
  469. ma mb mc)
  470. where ja,ma,...etc are integers or half integers.
  471. The function takes as input arguments only integers which corresponds
  472. to half integer units, e.g two_ja = 2 * ja
  473. For detailed description see
  474. <A HREF="http://mathworld.wolfram.com/Wigner3j-Symbol.html.html">
  475. Mathworld</A>.
  476. The implementation used is that of
  477. <A HREF="http://www.gnu.org/software/gsl/manual/html_node/3_002dj-Symbols.html#g_t3_002dj-Symbols">GSL</A>.
  478. @ingroup SpecFunc
  479. */
  480. double wigner_3j(int two_ja, int two_jb, int two_jc, int two_ma, int two_mb, int two_mc);
  481. /**
  482. Calculates the Wigner 6j coupling coefficients
  483. (ja jb jc
  484. jd je jf)
  485. where ja,jb,...etc are integers or half integers.
  486. The function takes as input arguments only integers which corresponds
  487. to half integer units, e.g two_ja = 2 * ja
  488. For detailed description see
  489. <A HREF="http://mathworld.wolfram.com/Wigner6j-Symbol.html">
  490. Mathworld</A>.
  491. The implementation used is that of
  492. <A HREF="http://www.gnu.org/software/gsl/manual/html_node/6_002dj-Symbols.html#g_t6_002dj-Symbols">GSL</A>.
  493. @ingroup SpecFunc
  494. */
  495. double wigner_6j(int two_ja, int two_jb, int two_jc, int two_jd, int two_je, int two_jf);
  496. /**
  497. Calculates the Wigner 9j coupling coefficients
  498. (ja jb jc
  499. jd je jf
  500. jg jh ji)
  501. where ja,jb...etc are integers or half integers.
  502. The function takes as input arguments only integers which corresponds
  503. to half integer units, e.g two_ja = 2 * ja
  504. For detailed description see
  505. <A HREF="http://mathworld.wolfram.com/Wigner9j-Symbol.html">
  506. Mathworld</A>.
  507. The implementation used is that of
  508. <A HREF="http://www.gnu.org/software/gsl/manual/html_node/9_002dj-Symbols.html#g_t9_002dj-Symbols">GSL</A>.
  509. @ingroup SpecFunc
  510. */
  511. 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);
  512. } // namespace Math
  513. } // namespace ROOT
  514. #endif //ROOT_Math_SpecFuncMathMore
  515. #endif // if defined (__CINT__) && !defined(__MAKECINT__)