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