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