ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
SpecFuncMathMore.h
Go to the documentation of this file.
1 // @(#)root/mathmore:$Id$
2 // Authors: L. Moneta, A. Zsenei 08/2005
3 
4 // Authors: Andras Zsenei & Lorenzo Moneta 06/2005
5 
6  /**********************************************************************
7  * *
8  * Copyright (c) 2004 ROOT Foundation, CERN/PH-SFT *
9  * *
10  * This library is free software; you can redistribute it and/or *
11  * modify it under the terms of the GNU General Public License *
12  * as published by the Free Software Foundation; either version 2 *
13  * of the License, or (at your option) any later version. *
14  * *
15  * This library is distributed in the hope that it will be useful, *
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of *
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU *
18  * General Public License for more details. *
19  * *
20  * You should have received a copy of the GNU General Public License *
21  * along with this library (see file COPYING); if not, write *
22  * to the Free Software Foundation, Inc., 59 Temple Place, Suite *
23  * 330, Boston, MA 02111-1307 USA, or contact the author. *
24  * *
25  **********************************************************************/
26 
27 #if defined(__CINT__) && !defined(__MAKECINT__)
28 // avoid to include header file when using CINT
29 #ifndef _WIN32
30 #include "../lib/libMathMore.so"
31 #else
32 #include "../bin/libMathMore.dll"
33 #endif
34 
35 #else
36 
37 
38 /**
39 
40 Special mathematical functions.
41 The naming and numbering of the functions is taken from
42 <A HREF="http://www.open-std.org/jtc1/sc22/wg21/docs/papers/2004/n1687.pdf">
43 Matt Austern,
44 (Draft) Technical Report on Standard Library Extensions,
45 N1687=04-0127, September 10, 2004</A>
46 
47 @author Created by Andras Zsenei on Mon Nov 8 2004
48 
49 @defgroup SpecFunc Special functions
50 
51 */
52 
53 
54 
55 
56 
57 #ifndef ROOT_Math_SpecFuncMathMore
58 #define ROOT_Math_SpecFuncMathMore
59 
60 
61 
62 
63 namespace ROOT {
64 namespace Math {
65 
66  /** @name Special Functions from MathMore */
67 
68 
69  /**
70 
71 
72  Computes the generalized Laguerre polynomials for
73  \f$ n \geq 0, m > -1 \f$.
74  They are defined in terms of the confluent hypergeometric function.
75  For integer values of m they can be defined in terms of the Laguerre polynomials \f$L_n(x)\f$:
76 
77  \f[ L_{n}^{m}(x) = (-1)^{m} \frac{d^m}{dx^m} L_{n+m}(x) \f]
78 
79 
80  For detailed description see
81  <A HREF="http://mathworld.wolfram.com/LaguerrePolynomial.html">Mathworld</A>.
82  The implementation used is that of
83  <A HREF="http://www.gnu.org/software/gsl/manual/html_node/Laguerre-Functions.html">GSL</A>.
84 
85  This function is an extension of C++0x, also consistent in GSL,
86  Abramowitz and Stegun 1972 and MatheMathica that uses non-integer values for m.
87  C++0x calls for 'int m', more restrictive than necessary.
88  The definition for was incorrect in 'n1687.pdf', but fixed in
89  <A HREF="http://www.open-std.org/jtc1/sc22/wg21/docs/papers/2005/n1836.pdf">n1836.pdf</A>,
90  the most recent draft as of 2007-07-01
91 
92 
93  @ingroup SpecFunc
94 
95  */
96  // [5.2.1.1] associated Laguerre polynomials
97 
98  double assoc_laguerre(unsigned n, double m, double x);
99 
100 
101 
102 
103  /**
104 
105  Computes the associated Legendre polynomials.
106 
107  \f[ P_{l}^{m}(x) = (1-x^2)^{m/2} \frac{d^m}{dx^m} P_{l}(x) \f]
108 
109  with \f$m \geq 0\f$, \f$ l \geq m \f$ and \f$ |x|<1 \f$.
110  There are two sign conventions for associated Legendre polynomials.
111  As is the case with the above formula, some authors (e.g., Arfken
112  1985, pp. 668-669) omit the Condon-Shortley phase \f$(-1)^m\f$,
113  while others include it (e.g., Abramowitz and Stegun 1972).
114  One possible way to distinguish the two conventions is due to
115  Abramowitz and Stegun (1972, p. 332), who use the notation
116 
117  \f[ P_{lm} (x) = (-1)^m P_{l}^{m} (x)\f]
118 
119  to distinguish the two. For detailed description see
120  <A HREF="http://mathworld.wolfram.com/LegendrePolynomial.html">
121  Mathworld</A>. The implementation used is that of
122  <A HREF="http://www.gnu.org/software/gsl/manual/html_node/Associated-Legendre-Polynomials-and-Spherical-Harmonics.html">GSL</A>.
123 
124  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,
125  sph_legendre(l,m,theta).
126 
127  @ingroup SpecFunc
128 
129  */
130  // [5.2.1.2] associated Legendre functions
131 
132  double assoc_legendre(unsigned l, unsigned m, double x);
133 
134 
135 
136 
137 
138  /**
139 
140  Calculates the complete elliptic integral of the first kind.
141 
142  \f[ K(k) = F(k, \pi / 2) = \int_{0}^{\pi /2} \frac{d \theta}{\sqrt{1 - k^2 \sin^2{\theta}}} \f]
143 
144  with \f$0 \leq k^2 \leq 1\f$. For detailed description see
145  <A HREF="http://mathworld.wolfram.com/CompleteEllipticIntegraloftheFirstKind.html">
146  Mathworld</A>. The implementation used is that of
147  <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC100">GSL</A>.
148 
149  @ingroup SpecFunc
150 
151  */
152  // [5.2.1.4] (complete) elliptic integral of the first kind
153 
154  double comp_ellint_1(double k);
155 
156 
157 
158 
159  /**
160 
161  Calculates the complete elliptic integral of the second kind.
162 
163  \f[ E(k) = E(k , \pi / 2) = \int_{0}^{\pi /2} \sqrt{1 - k^2 \sin^2{\theta}} d \theta \f]
164 
165  with \f$0 \leq k^2 \leq 1\f$. For detailed description see
166  <A HREF="http://mathworld.wolfram.com/CompleteEllipticIntegraloftheSecondKind.html">
167  Mathworld</A>. The implementation used is that of
168  <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC100">GSL</A>.
169 
170  @ingroup SpecFunc
171 
172  */
173  // [5.2.1.5] (complete) elliptic integral of the second kind
174 
175  double comp_ellint_2(double k);
176 
177 
178 
179 
180  /**
181 
182  Calculates the complete elliptic integral of the third kind.
183 
184  \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]
185 
186  with \f$0 \leq k^2 \leq 1\f$. There are two sign conventions
187  for elliptic integrals of the third kind. Some authors (Abramowitz
188  and Stegun,
189  <A HREF="http://mathworld.wolfram.com/EllipticIntegraloftheThirdKind.html">
190  Mathworld</A>,
191  <A HREF="http://www.open-std.org/jtc1/sc22/wg21/docs/papers/2004/n1687.pdf">
192  C++ standard proposal</A>) use the above formula, while others
193  (<A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC95">
194  GSL</A>, <A HREF="http://planetmath.org/encyclopedia/EllipticIntegralsAndJacobiEllipticFunctions.html">
195  Planetmath</A> and
196  <A HREF="http://wwwasdoc.web.cern.ch/wwwasdoc/shortwrupsdir/c346/top.html">
197  CERNLIB</A>) use the + sign in front of n in the denominator. In
198  order to be C++ compliant, the present library uses the former
199  convention. The implementation used is that of
200  <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC101">GSL</A>.
201 
202  @ingroup SpecFunc
203 
204  */
205  // [5.2.1.6] (complete) elliptic integral of the third kind
206  double comp_ellint_3(double n, double k);
207 
208 
209 
210 
211  /**
212 
213  Calculates the confluent hypergeometric functions of the first kind.
214 
215  \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]
216 
217  For detailed description see
218  <A HREF="http://mathworld.wolfram.com/ConfluentHypergeometricFunctionoftheFirstKind.html">
219  Mathworld</A>. The implementation used is that of
220  <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC125">GSL</A>.
221 
222  @ingroup SpecFunc
223 
224  */
225  // [5.2.1.7] confluent hypergeometric functions
226 
227  double conf_hyperg(double a, double b, double z);
228 
229 
230  /**
231 
232  Calculates the confluent hypergeometric functions of the second kind, known also as Kummer function of the second kind,
233  it is related to the confluent hypergeometric functions of the first kind.
234 
235  \f[ U(a,b,z) = \frac{ \pi}{ \sin{\pi b } } \left[ \frac{ _{1}F_{1}(a,b,z) } {\Gamma(a-b+1) }
236  - \frac{ z^{1-b} { _{1}F_{1}}(a-b+1,2-b,z)}{\Gamma(a)} \right] \f]
237 
238  For detailed description see
239  <A HREF="http://mathworld.wolfram.com/ConfluentHypergeometricFunctionoftheSecondKind.html">
240  Mathworld</A>. The implementation used is that of
241  <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC125">GSL</A>.
242  This function is not part of the C++ standard proposal
243 
244  @ingroup SpecFunc
245 
246  */
247  // confluent hypergeometric functions of second type
248 
249  double conf_hypergU(double a, double b, double z);
250 
251 
252 
253  /**
254 
255  Calculates the modified Bessel function of the first kind
256  (also called regular modified (cylindrical) Bessel function).
257 
258  \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]
259 
260  for \f$x>0, \nu > 0\f$. For detailed description see
261  <A HREF="http://mathworld.wolfram.com/ModifiedBesselFunctionoftheFirstKind.html">
262  Mathworld</A>. The implementation used is that of
263  <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC71">GSL</A>.
264 
265  @ingroup SpecFunc
266 
267  */
268  // [5.2.1.8] regular modified cylindrical Bessel functions
269 
270  double cyl_bessel_i(double nu, double x);
271 
272 
273 
274 
275  /**
276 
277  Calculates the (cylindrical) Bessel functions of the first kind (also
278  called regular (cylindrical) Bessel functions).
279 
280  \f[ J_{\nu} (x) = \sum_{k=0}^{\infty} \frac{(-1)^k(\frac{1}{2}x)^{\nu + 2k}}{k! \Gamma(\nu + k + 1)} \f]
281 
282  For detailed description see
283  <A HREF="http://mathworld.wolfram.com/BesselFunctionoftheFirstKind.html">
284  Mathworld</A>. The implementation used is that of
285  <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC69">GSL</A>.
286 
287  @ingroup SpecFunc
288 
289  */
290  // [5.2.1.9] cylindrical Bessel functions (of the first kind)
291 
292  double cyl_bessel_j(double nu, double x);
293 
294 
295 
296 
297 
298  /**
299 
300  Calculates the modified Bessel functions of the second kind
301  (also called irregular modified (cylindrical) Bessel functions).
302 
303  \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}}
304 & \mbox{for integral $\nu$} \end{array} \right. \f]
305 
306  for \f$x>0, \nu > 0\f$. For detailed description see
307  <A HREF="http://mathworld.wolfram.com/ModifiedBesselFunctionoftheSecondKind.html">
308  Mathworld</A>. The implementation used is that of
309  <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC72">GSL</A>.
310 
311  @ingroup SpecFunc
312 
313  */
314  // [5.2.1.10] irregular modified cylindrical Bessel functions
315 
316  double cyl_bessel_k(double nu, double x);
317 
318 
319 
320 
321  /**
322 
323  Calculates the (cylindrical) Bessel functions of the second kind
324  (also called irregular (cylindrical) Bessel functions or
325  (cylindrical) Neumann functions).
326 
327  \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]
328 
329  For detailed description see
330  <A HREF="http://mathworld.wolfram.com/BesselFunctionoftheSecondKind.html">
331  Mathworld</A>. The implementation used is that of
332  <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC70">GSL</A>.
333 
334  @ingroup SpecFunc
335 
336  */
337  // [5.2.1.11] cylindrical Neumann functions;
338  // cylindrical Bessel functions (of the second kind)
339 
340  double cyl_neumann(double nu, double x);
341 
342 
343 
344 
345  /**
346 
347  Calculates the incomplete elliptic integral of the first kind.
348 
349  \f[ F(k, \phi) = \int_{0}^{\phi} \frac{d \theta}{\sqrt{1 - k^2 \sin^2{\theta}}} \f]
350 
351  with \f$0 \leq k^2 \leq 1\f$. For detailed description see
352  <A HREF="http://mathworld.wolfram.com/EllipticIntegraloftheFirstKind.html">
353  Mathworld</A>. The implementation used is that of
354  <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC101">GSL</A>.
355 
356  @param k
357  @param phi angle in radians
358 
359  @ingroup SpecFunc
360 
361  */
362  // [5.2.1.12] (incomplete) elliptic integral of the first kind
363  // phi in radians
364 
365  double ellint_1(double k, double phi);
366 
367 
368 
369 
370  /**
371 
372  Calculates the complete elliptic integral of the second kind.
373 
374  \f[ E(k , \phi) = \int_{0}^{\phi} \sqrt{1 - k^2 \sin^2{\theta}} d \theta \f]
375 
376  with \f$0 \leq k^2 \leq 1\f$. For detailed description see
377  <A HREF="http://mathworld.wolfram.com/EllipticIntegraloftheSecondKind.html">
378  Mathworld</A>. The implementation used is that of
379  <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC101">GSL</A>.
380 
381  @param k
382  @param phi angle in radians
383 
384  @ingroup SpecFunc
385 
386  */
387  // [5.2.1.13] (incomplete) elliptic integral of the second kind
388  // phi in radians
389 
390  double ellint_2(double k, double phi);
391 
392 
393 
394 
395  /**
396 
397  Calculates the complete elliptic integral of the third kind.
398 
399  \f[ \Pi (n, k, \phi) = \int_{0}^{\phi} \frac{d \theta}{(1 - n \sin^2{\theta})\sqrt{1 - k^2 \sin^2{\theta}}} \f]
400 
401  with \f$0 \leq k^2 \leq 1\f$. There are two sign conventions
402  for elliptic integrals of the third kind. Some authors (Abramowitz
403  and Stegun,
404  <A HREF="http://mathworld.wolfram.com/EllipticIntegraloftheThirdKind.html">
405  Mathworld</A>,
406  <A HREF="http://www.open-std.org/jtc1/sc22/wg21/docs/papers/2004/n1687.pdf">
407  C++ standard proposal</A>) use the above formula, while others
408  (<A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC95">
409  GSL</A>, <A HREF="http://planetmath.org/encyclopedia/EllipticIntegralsAndJacobiEllipticFunctions.html">
410  Planetmath</A> and
411  <A HREF="http://wwwasdoc.web.cern.ch/wwwasdoc/shortwrupsdir/c346/top.html">
412  CERNLIB</A>) use the + sign in front of n in the denominator. In
413  order to be C++ compliant, the present library uses the former
414  convention. The implementation used is that of
415  <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC101">GSL</A>.
416 
417  @param n
418  @param k
419  @param phi angle in radians
420 
421  @ingroup SpecFunc
422 
423  */
424  // [5.2.1.14] (incomplete) elliptic integral of the third kind
425  // phi in radians
426 
427  double ellint_3(double n, double k, double phi);
428 
429 
430 
431 
432  /**
433 
434  Calculates the exponential integral.
435 
436  \f[ Ei(x) = - \int_{-x}^{\infty} \frac{e^{-t}}{t} dt \f]
437 
438  For detailed description see
439  <A HREF="http://mathworld.wolfram.com/ExponentialIntegral.html">
440  Mathworld</A>. The implementation used is that of
441  <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC115">GSL</A>.
442 
443  @ingroup SpecFunc
444 
445  */
446  // [5.2.1.15] exponential integral
447 
448  double expint(double x);
449 
450 
451 
452  // [5.2.1.16] Hermite polynomials
453 
454  //double hermite(unsigned n, double x);
455 
456 
457 
458 
459 
460  /**
461 
462  Calculates Gauss' hypergeometric function.
463 
464  \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]
465 
466  For detailed description see
467  <A HREF="http://mathworld.wolfram.com/HypergeometricFunction.html">
468  Mathworld</A>. The implementation used is that of
469  <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC125">GSL</A>.
470 
471  @ingroup SpecFunc
472 
473  */
474  // [5.2.1.17] hypergeometric functions
475 
476  double hyperg(double a, double b, double c, double x);
477 
478 
479 
480  /**
481 
482  Calculates the Laguerre polynomials
483 
484  \f[ P_{l}(x) = \frac{ e^x}{n!} \frac{d^n}{dx^n} (x^n - e^{-x}) \f]
485 
486  for \f$x \geq 0 \f$ in the Rodrigues representation.
487  They corresponds to the associated Laguerre polynomial of order m=0.
488  See Abramowitz and Stegun, (22.5.16)
489  For detailed description see
490  <A HREF="http://mathworld.wolfram.com/LaguerrePolynomial.html">
491  Mathworld</A>.
492  The are implemented using the associated Laguerre polynomial of order m=0.
493 
494  @ingroup SpecFunc
495 
496  */
497  // [5.2.1.18] Laguerre polynomials
498 
499  double laguerre(unsigned n, double x);
500 
501 
502  /**
503 
504  Calculates the Legendre polynomials.
505 
506  \f[ P_{l}(x) = \frac{1}{2^l l!} \frac{d^l}{dx^l} (x^2 - 1)^l \f]
507 
508  for \f$l \geq 0, |x|\leq1\f$ in the Rodrigues representation.
509  For detailed description see
510  <A HREF="http://mathworld.wolfram.com/LegendrePolynomial.html">
511  Mathworld</A>. The implementation used is that of
512  <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC129">GSL</A>.
513 
514  @ingroup SpecFunc
515 
516  */
517  // [5.2.1.19] Legendre polynomials
518 
519  double legendre(unsigned l, double x);
520 
521 
522 
523 
524  /**
525 
526  Calculates the Riemann zeta function.
527 
528  \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]
529 
530  For detailed description see
531  <A HREF="http://mathworld.wolfram.com/RiemannZetaFunction.html">
532  Mathworld</A>. The implementation used is that of
533  <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC149">GSL</A>.
534 
535  CHECK WHETHER THE IMPLEMENTATION CALCULATES X<1
536 
537  @ingroup SpecFunc
538 
539  */
540  // [5.2.1.20] Riemann zeta function
541 
542  double riemann_zeta(double x);
543 
544 
545  /**
546 
547  Calculates the spherical Bessel functions of the first kind
548  (also called regular spherical Bessel functions).
549 
550  \f[ j_{n}(x) = \sqrt{\frac{\pi}{2x}} J_{n+1/2}(x) \f]
551 
552  For detailed description see
553  <A HREF="http://mathworld.wolfram.com/SphericalBesselFunctionoftheFirstKind.html">
554  Mathworld</A>. The implementation used is that of
555  <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC73">GSL</A>.
556 
557  @ingroup SpecFunc
558 
559  */
560  // [5.2.1.21] spherical Bessel functions of the first kind
561 
562  double sph_bessel(unsigned n, double x);
563 
564 
565  /**
566 
567  Computes the spherical (normalized) associated Legendre polynomials,
568  or spherical harmonic without azimuthal dependence (\f$e^(im\phi)\f$).
569 
570  \f[ Y_l^m(theta,0) = \sqrt{(2l+1)/(4\pi)} \sqrt{(l-m)!/(l+m)!} P_l^m(cos \theta) \f]
571 
572  for \f$m \geq 0, l \geq m\f$,
573  where the Condon-Shortley phase \f$(-1)^m\f$ is included in P_l^m(x)
574  This function is consistent with both C++0x and GSL,
575  even though there is a discrepancy in where to include the phase.
576  There is no reference in Abramowitz and Stegun.
577 
578 
579  @ingroup SpecFunc
580 
581  */
582 
583  // [5.2.1.22] spherical associated Legendre functions
584 
585  double sph_legendre(unsigned l, unsigned m, double theta);
586 
587 
588  /**
589 
590  Calculates the spherical Bessel functions of the second kind
591  (also called irregular spherical Bessel functions or
592  spherical Neumann functions).
593 
594  \f[ n_n(x) = y_n(x) = \sqrt{\frac{\pi}{2x}} N_{n+1/2}(x) \f]
595 
596  For detailed description see
597  <A HREF="http://mathworld.wolfram.com/SphericalBesselFunctionoftheSecondKind.html">
598  Mathworld</A>. The implementation used is that of
599  <A HREF="http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC74">GSL</A>.
600 
601  @ingroup SpecFunc
602 
603  */
604  // [5.2.1.23] spherical Neumann functions
605 
606  double sph_neumann(unsigned n, double x);
607 
608  /**
609 
610  Calculates the Airy function Ai
611 
612  \f[ Ai(x) = \frac{1}{\pi} \int\limits_{0}^{\infty} \cos(xt + t^3/3) dt \f]
613 
614  For detailed description see
615  <A HREF="http://mathworld.wolfram.com/AiryFunctions.html">
616  Mathworld</A>
617  and <A HREF="http://www.nrbook.com/abramowitz_and_stegun/page_446.htm">Abramowitz&Stegun, Sect. 10.4</A>.
618  The implementation used is that of
619  <A HREF="http://www.gnu.org/software/gsl/manual/html_node/Airy-Functions.html">GSL</A>.
620 
621  @ingroup SpecFunc
622 
623  */
624  // Airy function Ai
625 
626  double airy_Ai(double x);
627 
628  /**
629 
630  Calculates the Airy function Bi
631 
632  \f[ Bi(x) = \frac{1}{\pi} \int\limits_{0}^{\infty} [\exp(xt-t^3/3) + \cos(xt + t^3/3)] dt \f]
633 
634  For detailed description see
635  <A HREF="http://mathworld.wolfram.com/AiryFunctions.html">
636  Mathworld</A>
637  and <A HREF="http://www.nrbook.com/abramowitz_and_stegun/page_446.htm">Abramowitz&Stegun, Sect. 10.4</A>.
638  The implementation used is that of
639  <A HREF="http://www.gnu.org/software/gsl/manual/html_node/Airy-Functions.html">GSL</A>.
640 
641  @ingroup SpecFunc
642 
643  */
644  // Airy function Bi
645 
646  double airy_Bi(double x);
647 
648  /**
649 
650  Calculates the derivative of the Airy function Ai
651 
652  \f[ Ai(x) = \frac{1}{\pi} \int\limits_{0}^{\infty} \cos(xt + t^3/3) dt \f]
653 
654  For detailed description see
655  <A HREF="http://mathworld.wolfram.com/AiryFunctions.html">
656  Mathworld</A>
657  and <A HREF="http://www.nrbook.com/abramowitz_and_stegun/page_446.htm">Abramowitz&Stegun, Sect. 10.4</A>.
658  The implementation used is that of
659  <A HREF="http://www.gnu.org/software/gsl/manual/html_node/Derivatives-of-Airy-Functions.html">GSL</A>.
660 
661  @ingroup SpecFunc
662 
663  */
664  // Derivative of the Airy function Ai
665 
666  double airy_Ai_deriv(double x);
667 
668  /**
669 
670  Calculates the derivative of the Airy function Bi
671 
672  \f[ Bi(x) = \frac{1}{\pi} \int\limits_{0}^{\infty} [\exp(xt-t^3/3) + \cos(xt + t^3/3)] dt \f]
673 
674  For detailed description see
675  <A HREF="http://mathworld.wolfram.com/AiryFunctions.html">
676  Mathworld</A>
677  and <A HREF="http://www.nrbook.com/abramowitz_and_stegun/page_446.htm">Abramowitz&Stegun, Sect. 10.4</A>.
678  The implementation used is that of
679  <A HREF="http://www.gnu.org/software/gsl/manual/html_node/Derivatives-of-Airy-Functions.html">GSL</A>.
680 
681  @ingroup SpecFunc
682 
683  */
684  // Derivative of the Airy function Bi
685 
686  double airy_Bi_deriv(double x);
687 
688  /**
689 
690  Calculates the zeroes of the Airy function Ai
691 
692  \f[ Ai(x) = \frac{1}{\pi} \int\limits_{0}^{\infty} \cos(xt + t^3/3) dt \f]
693 
694  For detailed description see
695  <A HREF="http://mathworld.wolfram.com/AiryFunctionZeros.html">
696  Mathworld</A>
697  and <A HREF="http://www.nrbook.com/abramowitz_and_stegun/page_446.htm">Abramowitz&Stegun, Sect. 10.4</A>.
698  The implementation used is that of
699  <A HREF="http://www.gnu.org/software/gsl/manual/html_node/Zeros-of-Airy-Functions.html">GSL</A>.
700 
701  @ingroup SpecFunc
702 
703  */
704  // s-th zero of the Airy function Ai
705 
706  double airy_zero_Ai(unsigned int s);
707 
708  /**
709 
710  Calculates the zeroes of the Airy function Bi
711 
712  \f[ Bi(x) = \frac{1}{\pi} \int\limits_{0}^{\infty} [\exp(xt-t^3/3) + \cos(xt + t^3/3)] dt \f]
713 
714  For detailed description see
715  <A HREF="http://mathworld.wolfram.com/AiryFunctionZeros.html">
716  Mathworld</A>
717  and <A HREF="http://www.nrbook.com/abramowitz_and_stegun/page_446.htm">Abramowitz&Stegun, Sect. 10.4</A>.
718  The implementation used is that of
719  <A HREF="http://www.gnu.org/software/gsl/manual/html_node/Zeros-of-Airy-Functions.html">GSL</A>.
720 
721  @ingroup SpecFunc
722 
723  */
724  // s-th zero of the Airy function Bi
725 
726  double airy_zero_Bi(unsigned int s);
727 
728  /**
729 
730  Calculates the zeroes of the derivative of the Airy function Ai
731 
732  \f[ Ai(x) = \frac{1}{\pi} \int\limits_{0}^{\infty} \cos(xt + t^3/3) dt \f]
733 
734  For detailed description see
735  <A HREF="http://mathworld.wolfram.com/AiryFunctionZeros.html">
736  Mathworld</A>
737  and <A HREF="http://www.nrbook.com/abramowitz_and_stegun/page_446.htm">Abramowitz&Stegun, Sect. 10.4</A>.
738  The implementation used is that of
739  <A HREF="http://www.gnu.org/software/gsl/manual/html_node/Zeros-of-Derivatives-of-Airy-Functions.html">GSL</A>.
740 
741  @ingroup SpecFunc
742 
743  */
744  // s-th zero of the derivative of the Airy function Ai
745 
746  double airy_zero_Ai_deriv(unsigned int s);
747 
748  /**
749 
750  Calculates the zeroes of the derivative of the Airy function Bi
751 
752  \f[ Bi(x) = \frac{1}{\pi} \int\limits_{0}^{\infty} [\exp(xt-t^3/3) + \cos(xt + t^3/3)] dt \f]
753 
754  For detailed description see
755  <A HREF="http://mathworld.wolfram.com/AiryFunctionZeros.html">
756  Mathworld</A>
757  and <A HREF="http://www.nrbook.com/abramowitz_and_stegun/page_446.htm">Abramowitz&Stegun, Sect. 10.4</A>.
758  The implementation used is that of
759  <A HREF="http://www.gnu.org/software/gsl/manual/html_node/Zeros-of-Derivatives-of-Airy-Functions.html">GSL</A>.
760 
761  @ingroup SpecFunc
762 
763  */
764  // s-th zero of the derivative of the Airy function Bi
765 
766  double airy_zero_Bi_deriv(unsigned int s);
767 
768  /**
769 
770  Calculates the Wigner 3j coupling coefficients
771 
772  (ja jb jc
773  ma mb mc)
774 
775  where ja,ma,...etc are integers or half integers.
776  The function takes as input arguments only integers which corresponds
777  to half integer units, e.g two_ja = 2 * ja
778 
779  For detailed description see
780  <A HREF="http://mathworld.wolfram.com/Wigner3j-Symbol.html.html">
781  Mathworld</A>.
782  The implementation used is that of
783  <A HREF="http://www.gnu.org/software/gsl/manual/html_node/3_002dj-Symbols.html#g_t3_002dj-Symbols">GSL</A>.
784 
785  @ingroup SpecFunc
786 
787  */
788 
789  double wigner_3j(int two_ja, int two_jb, int two_jc, int two_ma, int two_mb, int two_mc);
790 
791  /**
792 
793  Calculates the Wigner 6j coupling coefficients
794 
795  (ja jb jc
796  jd je jf)
797 
798  where ja,jb,...etc are integers or half integers.
799  The function takes as input arguments only integers which corresponds
800  to half integer units, e.g two_ja = 2 * ja
801 
802  For detailed description see
803  <A HREF="http://mathworld.wolfram.com/Wigner6j-Symbol.html">
804  Mathworld</A>.
805  The implementation used is that of
806  <A HREF="http://www.gnu.org/software/gsl/manual/html_node/6_002dj-Symbols.html#g_t6_002dj-Symbols">GSL</A>.
807 
808  @ingroup SpecFunc
809 
810  */
811 
812  double wigner_6j(int two_ja, int two_jb, int two_jc, int two_jd, int two_je, int two_jf);
813 
814  /**
815 
816  Calculates the Wigner 9j coupling coefficients
817 
818  (ja jb jc
819  jd je jf
820  jg jh ji)
821 
822  where ja,jb...etc are integers or half integers.
823  The function takes as input arguments only integers which corresponds
824  to half integer units, e.g two_ja = 2 * ja
825 
826 
827  For detailed description see
828  <A HREF="http://mathworld.wolfram.com/Wigner9j-Symbol.html">
829  Mathworld</A>.
830  The implementation used is that of
831  <A HREF="http://www.gnu.org/software/gsl/manual/html_node/9_002dj-Symbols.html#g_t9_002dj-Symbols">GSL</A>.
832 
833  @ingroup SpecFunc
834 
835  */
836 
837  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);
838 
839 
840 
841 } // namespace Math
842 } // namespace ROOT
843 
844 
845 #endif //ROOT_Math_SpecFuncMathMore
846 
847 #endif // if defined (__CINT__) && !defined(__MAKECINT__)
double airy_zero_Bi(unsigned int s)
Calculates the zeroes of the Airy function Bi.
double sph_legendre(unsigned l, unsigned m, double theta)
Computes the spherical (normalized) associated Legendre polynomials, or spherical harmonic without az...
Float_t theta
Definition: shapesAnim.C:5
return c
double cyl_bessel_j(double nu, double x)
Calculates the (cylindrical) Bessel functions of the first kind (also called regular (cylindrical) Be...
double airy_zero_Ai_deriv(unsigned int s)
Calculates the zeroes of the derivative of the Airy function Ai.
TArc * a
Definition: textangle.C:12
double ellint_2(double k, double phi)
Calculates the complete elliptic integral of the second kind.
double wigner_6j(int two_ja, int two_jb, int two_jc, int two_jd, int two_je, int two_jf)
Calculates the Wigner 6j coupling coefficients.
double legendre(unsigned l, double x)
Calculates the Legendre polynomials.
double comp_ellint_1(double k)
Calculates the complete elliptic integral of the first kind.
double sph_bessel(unsigned n, double x)
Calculates the spherical Bessel functions of the first kind (also called regular spherical Bessel fun...
Double_t x[n]
Definition: legend1.C:17
double conf_hyperg(double a, double b, double z)
Calculates the confluent hypergeometric functions of the first kind.
double hyperg(double a, double b, double c, double x)
Calculates Gauss' hypergeometric function.
Float_t z[5]
Definition: Ifit.C:16
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)
Calculates the Wigner 9j coupling coefficients.
double airy_zero_Bi_deriv(unsigned int s)
Calculates the zeroes of the derivative of the Airy function Bi.
double expint(double x)
Calculates the exponential integral.
double ellint_3(double n, double k, double phi)
Calculates the complete elliptic integral of the third kind.
double riemann_zeta(double x)
Calculates the Riemann zeta function.
double comp_ellint_2(double k)
Calculates the complete elliptic integral of the second kind.
TMarker * m
Definition: textangle.C:8
double assoc_laguerre(unsigned n, double m, double x)
Computes the generalized Laguerre polynomials for .
double assoc_legendre(unsigned l, unsigned m, double x)
Computes the associated Legendre polynomials.
TLine * l
Definition: textangle.C:4
double airy_zero_Ai(unsigned int s)
Calculates the zeroes of the Airy function Ai.
double airy_Bi(double x)
Calculates the Airy function Bi.
double laguerre(unsigned n, double x)
Calculates the Laguerre polynomials.
double airy_Bi_deriv(double x)
Calculates the derivative of the Airy function Bi.
double comp_ellint_3(double n, double k)
Calculates the complete elliptic integral of the third kind.
Float_t phi
Definition: shapesAnim.C:6
double airy_Ai_deriv(double x)
Calculates the derivative of the Airy function Ai.
double ellint_1(double k, double phi)
Calculates the incomplete elliptic integral of the first kind.
double cyl_neumann(double nu, double x)
Calculates the (cylindrical) Bessel functions of the second kind (also called irregular (cylindrical)...
double wigner_3j(int two_ja, int two_jb, int two_jc, int two_ma, int two_mb, int two_mc)
Calculates the Wigner 3j coupling coefficients.
double sph_neumann(unsigned n, double x)
Calculates the spherical Bessel functions of the second kind (also called irregular spherical Bessel ...
double conf_hypergU(double a, double b, double z)
Calculates the confluent hypergeometric functions of the second kind, known also as Kummer function o...
double cyl_bessel_i(double nu, double x)
Calculates the modified Bessel function of the first kind (also called regular modified (cylindrical)...
double airy_Ai(double x)
Calculates the Airy function Ai.
const Int_t n
Definition: legend1.C:16
double cyl_bessel_k(double nu, double x)
Calculates the modified Bessel functions of the second kind (also called irregular modified (cylindri...