16 #define PI 3.14159265358979323846264338328 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" 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" 43 return gsl_sf_laguerre_n(n, m, x);
54 return (m%2 == 0) ? gsl_sf_legendre_Plm(l, m, x) : -gsl_sf_legendre_Plm(l, m, x);
67 return gsl_sf_ellint_Kcomp(k, GSL_PREC_DOUBLE);
78 return gsl_sf_ellint_Ecomp(k, GSL_PREC_DOUBLE);
130 return gsl_sf_ellint_P(
PI/2.0, k, -n, GSL_PREC_DOUBLE);
141 return gsl_sf_hyperg_1F1(a, b, z);
149 return gsl_sf_hyperg_U(a, b, z);
160 return gsl_sf_bessel_Inu(nu, x);
171 return gsl_sf_bessel_Jnu(nu, x);
182 return gsl_sf_bessel_Knu(nu, x);
194 return gsl_sf_bessel_Ynu(nu, x);
206 return gsl_sf_ellint_F(phi, k, GSL_PREC_DOUBLE);
218 return gsl_sf_ellint_E(phi, k, GSL_PREC_DOUBLE);
272 return gsl_sf_ellint_P(phi, k, -n, GSL_PREC_DOUBLE);
283 return gsl_sf_expint_Ei(x);
292 return gsl_sf_expint_En(n, x);
310 double hyperg(
double a,
double b,
double c,
double x) {
312 return gsl_sf_hyperg_2F1(a, b, c, x);
322 return gsl_sf_laguerre_n(n, 0, x);
333 return gsl_sf_legendre_Pl(l, x);
344 return gsl_sf_zeta(x);
355 return gsl_sf_bessel_jl(n, x);
366 return gsl_sf_legendre_sphPlm(l, m,
std::cos(theta));
378 return gsl_sf_bessel_yl(n, x);
386 return gsl_sf_airy_Ai(x, GSL_PREC_DOUBLE);
394 return gsl_sf_airy_Bi(x, GSL_PREC_DOUBLE);
402 return gsl_sf_airy_Ai_deriv(x, GSL_PREC_DOUBLE);
410 return gsl_sf_airy_Bi_deriv(x, GSL_PREC_DOUBLE);
418 return gsl_sf_airy_zero_Ai(s);
426 return gsl_sf_airy_zero_Bi(s);
434 return gsl_sf_airy_zero_Ai_deriv(s);
442 return gsl_sf_airy_zero_Bi_deriv(s);
448 double wigner_3j(
int ja,
int jb,
int jc,
int ma,
int mb,
int mc) {
449 return gsl_sf_coupling_3j(ja,jb,jc,ma,mb,mc);
452 double wigner_6j(
int ja,
int jb,
int jc,
int jd,
int je,
int jf) {
453 return gsl_sf_coupling_6j(ja,jb,jc,jd,je,jf);
456 double wigner_9j(
int ja,
int jb,
int jc,
int jd,
int je,
int jf,
int jg,
int jh,
int ji) {
457 return gsl_sf_coupling_9j(ja,jb,jc,jd,je,jf,jg,jh,ji);
double airy_zero_Bi(unsigned int s)
Calculates the zeroes of the Airy function Bi.
Namespace for new ROOT classes and functions.
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.
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 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 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.
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t Double_t Double_t Double_t Double_t Int_t m
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.
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.
double airy_zero_Ai(unsigned int s)
Calculates the zeroes of the Airy function Ai.
double expint_n(int n, double x)
double airy_Bi(double x)
Calculates the Airy function Bi.
double laguerre(unsigned n, double x)
Calculates the Laguerre polynomials.
* x
Deprecated and error prone model selection interface.
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.
you should not use this method at all Int_t Int_t z
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.
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t Double_t Double_t Double_t Double_t b
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.
double cyl_bessel_k(double nu, double x)
Calculates the modified Bessel functions of the second kind (also called irregular modified (cylindri...