ROOT  6.06/08
Reference Guide
SpecFuncMathMore.cxx
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) 2005 , LCG ROOT MathLib Team *
9  * *
10  * *
11  **********************************************************************/
12 
13 #include <cmath>
14 
15 #ifndef PI
16 #define PI 3.14159265358979323846264338328 /* pi */
17 #endif
18 
19 
20 #include "gsl/gsl_sf_bessel.h"
21 #include "gsl/gsl_sf_legendre.h"
22 #include "gsl/gsl_sf_laguerre.h"
23 #include "gsl/gsl_sf_hyperg.h"
24 #include "gsl/gsl_sf_ellint.h"
25 //#include "gsl/gsl_sf_gamma.h"
26 #include "gsl/gsl_sf_expint.h"
27 #include "gsl/gsl_sf_zeta.h"
28 #include "gsl/gsl_sf_airy.h"
29 #include "gsl/gsl_sf_coupling.h"
30 
31 
32 namespace ROOT {
33 namespace Math {
34 
35 
36 
37 
38 // [5.2.1.1] associated Laguerre polynomials
39 // (26.x.12)
40 
41 double assoc_laguerre(unsigned n, double m, double x) {
42 
43  return gsl_sf_laguerre_n(n, m, x);
44 
45 }
46 
47 
48 
49 // [5.2.1.2] associated Legendre functions
50 // (26.x.8)
51 
52 double assoc_legendre(unsigned l, unsigned m, double x) {
53  // add (-1)^m
54  return (m%2 == 0) ? gsl_sf_legendre_Plm(l, m, x) : -gsl_sf_legendre_Plm(l, m, x);
55 
56 }
57 
58 
59 
60 
61 
62 // [5.2.1.4] (complete) elliptic integral of the first kind
63 // (26.x.15.2)
64 
65 double comp_ellint_1(double k) {
66 
67  return gsl_sf_ellint_Kcomp(k, GSL_PREC_DOUBLE);
68 
69 }
70 
71 
72 
73 // [5.2.1.5] (complete) elliptic integral of the second kind
74 // (26.x.16.2)
75 
76 double comp_ellint_2(double k) {
77 
78  return gsl_sf_ellint_Ecomp(k, GSL_PREC_DOUBLE);
79 
80 }
81 
82 
83 
84 // [5.2.1.6] (complete) elliptic integral of the third kind
85 // (26.x.17.2)
86 /**
87 Complete elliptic integral of the third kind.
88 
89 There are two different definitions used for the elliptic
90 integral of the third kind:
91 \f[
92 P(\phi,k,n) = \int_0^\phi \frac{dt}{(1 + n \sin^2{t})\sqrt{1 - k^2 \sin^2{t}}}
93 \f]
94 and
95 \f[
96 P(\phi,k,n) = \int_0^\phi \frac{dt}{(1 - n \sin^2{t})\sqrt{1 - k^2 \sin^2{t}}}
97 \f]
98 the former is adopted by
99 
100 - GSL
101  http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC95
102 
103 - Planetmath
104  http://planetmath.org/encyclopedia/EllipticIntegralsAndJacobiEllipticFunctions.html
105 
106 - CERNLIB
107  https://cern-tex.web.cern.ch/cern-tex/shortwrupsdir/c346/top.html
108 
109  while the latter is used by
110 
111 - Abramowitz and Stegun
112 
113 - Mathematica
114  http://mathworld.wolfram.com/EllipticIntegraloftheThirdKind.html
115 
116 - C++ standard
117  http://www.open-std.org/jtc1/sc22/wg21/docs/papers/2004/n1687.pdf
118 
119  in order to be C++ compliant, we decided to use the latter, hence the
120  change of the sign in the function call to GSL.
121 
122  */
123 
124 double comp_ellint_3(double n, double k) {
125 
126  return gsl_sf_ellint_P(PI/2.0, k, -n, GSL_PREC_DOUBLE);
127 
128 }
129 
130 
131 
132 // [5.2.1.7] confluent hypergeometric functions
133 // (26.x.14)
134 
135 double conf_hyperg(double a, double b, double z) {
136 
137  return gsl_sf_hyperg_1F1(a, b, z);
138 
139 }
140 
141 // confluent hypergeometric functions of second type
142 
143 double conf_hypergU(double a, double b, double z) {
144 
145  return gsl_sf_hyperg_U(a, b, z);
146 
147 }
148 
149 
150 
151 // [5.2.1.8] regular modified cylindrical Bessel functions
152 // (26.x.4.1)
153 
154 double cyl_bessel_i(double nu, double x) {
155 
156  return gsl_sf_bessel_Inu(nu, x);
157 
158 }
159 
160 
161 
162 // [5.2.1.9] cylindrical Bessel functions (of the first kind)
163 // (26.x.2)
164 
165 double cyl_bessel_j(double nu, double x) {
166 
167  return gsl_sf_bessel_Jnu(nu, x);
168 
169 }
170 
171 
172 
173 // [5.2.1.10] irregular modified cylindrical Bessel functions
174 // (26.x.4.2)
175 
176 double cyl_bessel_k(double nu, double x) {
177 
178  return gsl_sf_bessel_Knu(nu, x);
179 
180 }
181 
182 
183 
184 // [5.2.1.11] cylindrical Neumann functions;
185 // cylindrical Bessel functions (of the second kind)
186 // (26.x.3)
187 
188 double cyl_neumann(double nu, double x) {
189 
190  return gsl_sf_bessel_Ynu(nu, x);
191 
192 }
193 
194 
195 
196 // [5.2.1.12] (incomplete) elliptic integral of the first kind
197 // phi in radians
198 // (26.x.15.1)
199 
200 double ellint_1(double k, double phi) {
201 
202  return gsl_sf_ellint_F(phi, k, GSL_PREC_DOUBLE);
203 
204 }
205 
206 
207 
208 // [5.2.1.13] (incomplete) elliptic integral of the second kind
209 // phi in radians
210 // (26.x.16.1)
211 
212 double ellint_2(double k, double phi) {
213 
214  return gsl_sf_ellint_E(phi, k, GSL_PREC_DOUBLE);
215 
216 }
217 
218 
219 
220 // [5.2.1.14] (incomplete) elliptic integral of the third kind
221 // phi in radians
222 // (26.x.17.1)
223 /**
224 Incomplete elliptic integral of the third kind.
225 
226 There are two different definitions used for the elliptic
227 integral of the third kind:
228 \f[
229 P(\phi,k,n) = \int_0^\phi \frac{dt}{(1 + n \sin^2{t})\sqrt{1 - k^2 \sin^2{t}}}
230 \f]
231 and
232 \f[
233 P(\phi,k,n) = \int_0^\phi \frac{dt}{(1 - n \sin^2{t})\sqrt{1 - k^2 \sin^2{t}}}
234 \f]
235 the former is adopted by
236 
237 - GSL
238  http://www.gnu.org/software/gsl/manual/gsl-ref_7.html#SEC95
239 
240 - Planetmath
241  http://planetmath.org/encyclopedia/EllipticIntegralsAndJacobiEllipticFunctions.html
242 
243 - CERNLIB
244  https://cern-tex.web.cern.ch/cern-tex/shortwrupsdir/c346/top.html
245 
246  while the latter is used by
247 
248 - Abramowitz and Stegun
249 
250 - Mathematica
251  http://mathworld.wolfram.com/EllipticIntegraloftheThirdKind.html
252 
253 - C++ standard
254  http://www.open-std.org/jtc1/sc22/wg21/docs/papers/2004/n1687.pdf
255 
256  in order to be C++ compliant, we decided to use the latter, hence the
257  change of the sign in the function call to GSL.
258 
259  */
260 
261 double ellint_3(double n, double k, double phi) {
262 
263  return gsl_sf_ellint_P(phi, k, -n, GSL_PREC_DOUBLE);
264 
265 }
266 
267 
268 
269 // [5.2.1.15] exponential integral
270 // (26.x.20)
271 
272 double expint(double x) {
273 
274  return gsl_sf_expint_Ei(x);
275 
276 }
277 
278 
279 
280 // [5.2.1.16] Hermite polynomials
281 // (26.x.10)
282 
283 //double hermite(unsigned n, double x) {
284 //}
285 
286 
287 
288 
289 // [5.2.1.17] hypergeometric functions
290 // (26.x.13)
291 
292 double hyperg(double a, double b, double c, double x) {
293 
294  return gsl_sf_hyperg_2F1(a, b, c, x);
295 
296 }
297 
298 
299 
300 // [5.2.1.18] Laguerre polynomials
301 // (26.x.11)
302 
303 double laguerre(unsigned n, double x) {
304  return gsl_sf_laguerre_n(n, 0, x);
305 }
306 
307 
308 
309 
310 // [5.2.1.19] Legendre polynomials
311 // (26.x.7)
312 
313 double legendre(unsigned l, double x) {
314 
315  return gsl_sf_legendre_Pl(l, x);
316 
317 }
318 
319 
320 
321 // [5.2.1.20] Riemann zeta function
322 // (26.x.22)
323 
324 double riemann_zeta(double x) {
325 
326  return gsl_sf_zeta(x);
327 
328 }
329 
330 
331 
332 // [5.2.1.21] spherical Bessel functions of the first kind
333 // (26.x.5)
334 
335 double sph_bessel(unsigned n, double x) {
336 
337  return gsl_sf_bessel_jl(n, x);
338 
339 }
340 
341 
342 
343 // [5.2.1.22] spherical associated Legendre functions
344 // (26.x.9) ??????????
345 
346 double sph_legendre(unsigned l, unsigned m, double theta) {
347 
348  return gsl_sf_legendre_sphPlm(l, m, std::cos(theta));
349 
350 }
351 
352 
353 
354 
355 // [5.2.1.23] spherical Neumann functions
356 // (26.x.6)
357 
358 double sph_neumann(unsigned n, double x) {
359 
360  return gsl_sf_bessel_yl(n, x);
361 
362 }
363 
364 // Airy function Ai
365 
366 double airy_Ai(double x) {
367 
368  return gsl_sf_airy_Ai(x, GSL_PREC_DOUBLE);
369 
370 }
371 
372 // Airy function Bi
373 
374 double airy_Bi(double x) {
375 
376  return gsl_sf_airy_Bi(x, GSL_PREC_DOUBLE);
377 
378 }
379 
380 // Derivative of the Airy function Ai
381 
382 double airy_Ai_deriv(double x) {
383 
384  return gsl_sf_airy_Ai_deriv(x, GSL_PREC_DOUBLE);
385 
386 }
387 
388 // Derivative of the Airy function Bi
389 
390 double airy_Bi_deriv(double x) {
391 
392  return gsl_sf_airy_Bi_deriv(x, GSL_PREC_DOUBLE);
393 
394 }
395 
396 // s-th zero of the Airy function Ai
397 
398 double airy_zero_Ai(unsigned int s) {
399 
400  return gsl_sf_airy_zero_Ai(s);
401 
402 }
403 
404 // s-th zero of the Airy function Bi
405 
406 double airy_zero_Bi(unsigned int s) {
407 
408  return gsl_sf_airy_zero_Bi(s);
409 
410 }
411 
412 // s-th zero of the derivative of the Airy function Ai
413 
414 double airy_zero_Ai_deriv(unsigned int s) {
415 
416  return gsl_sf_airy_zero_Ai_deriv(s);
417 
418 }
419 
420 // s-th zero of the derivative of the Airy function Bi
421 
422 double airy_zero_Bi_deriv(unsigned int s) {
423 
424  return gsl_sf_airy_zero_Bi_deriv(s);
425 
426 }
427 
428 // wigner coefficient
429 
430 double wigner_3j(int ja, int jb, int jc, int ma, int mb, int mc) {
431  return gsl_sf_coupling_3j(ja,jb,jc,ma,mb,mc);
432 }
433 
434 double wigner_6j(int ja, int jb, int jc, int jd, int je, int jf) {
435  return gsl_sf_coupling_6j(ja,jb,jc,jd,je,jf);
436 }
437 
438 double wigner_9j(int ja, int jb, int jc, int jd, int je, int jf, int jg, int jh, int ji) {
439  return gsl_sf_coupling_9j(ja,jb,jc,jd,je,jf,jg,jh,ji);
440 }
441 
442 } // namespace Math
443 } // namespace ROOT
double airy_zero_Bi(unsigned int s)
Calculates the zeroes of the Airy function Bi.
Namespace for new ROOT classes and functions.
Definition: ROOT.py:1
double sph_legendre(unsigned l, unsigned m, double theta)
Computes the spherical (normalized) associated Legendre polynomials, or spherical harmonic without az...
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 cos(double)
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&#39; hypergeometric function.
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.
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)...
Namespace for new Math classes and functions.
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
#define PI
double cyl_bessel_k(double nu, double x)
Calculates the modified Bessel functions of the second kind (also called irregular modified (cylindri...