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);
61 return gsl_sf_legendre_Plm(
l,
m,
x);
71 return gsl_sf_ellint_Kcomp(k, GSL_PREC_DOUBLE);
82 return gsl_sf_ellint_Ecomp(k, GSL_PREC_DOUBLE);
134 return gsl_sf_ellint_P(
PI/2.0, k, -
n, GSL_PREC_DOUBLE);
145 return gsl_sf_hyperg_1F1(
a,
b, z);
153 return gsl_sf_hyperg_U(
a,
b, z);
164 return gsl_sf_bessel_Inu(nu,
x);
175 return gsl_sf_bessel_Jnu(nu,
x);
186 return gsl_sf_bessel_Knu(nu,
x);
198 return gsl_sf_bessel_Ynu(nu,
x);
210 return gsl_sf_ellint_F(phi, k, GSL_PREC_DOUBLE);
222 return gsl_sf_ellint_E(phi, k, GSL_PREC_DOUBLE);
276 return gsl_sf_ellint_P(phi, k, -
n, GSL_PREC_DOUBLE);
287 return gsl_sf_expint_Ei(
x);
296 return gsl_sf_expint_En(
n,
x);
316 return gsl_sf_hyperg_2F1(
a,
b,
c,
x);
326 return gsl_sf_laguerre_n(
n, 0,
x);
337 return gsl_sf_legendre_Pl(
l,
x);
348 return gsl_sf_zeta(
x);
359 return gsl_sf_bessel_jl(
n,
x);
370 return gsl_sf_legendre_sphPlm(
l,
m,
std::cos(theta));
382 return gsl_sf_bessel_yl(
n,
x);
390 return gsl_sf_airy_Ai(
x, GSL_PREC_DOUBLE);
398 return gsl_sf_airy_Bi(
x, GSL_PREC_DOUBLE);
406 return gsl_sf_airy_Ai_deriv(
x, GSL_PREC_DOUBLE);
414 return gsl_sf_airy_Bi_deriv(
x, GSL_PREC_DOUBLE);
422 return gsl_sf_airy_zero_Ai(
s);
430 return gsl_sf_airy_zero_Bi(
s);
438 return gsl_sf_airy_zero_Ai_deriv(
s);
446 return gsl_sf_airy_zero_Bi_deriv(
s);
452double wigner_3j(
int ja,
int jb,
int jc,
int ma,
int mb,
int mc) {
453 return gsl_sf_coupling_3j(ja,jb,jc,ma,mb,mc);
456double wigner_6j(
int ja,
int jb,
int jc,
int jd,
int je,
int jf) {
457 return gsl_sf_coupling_6j(ja,jb,jc,jd,je,jf);
460double wigner_9j(
int ja,
int jb,
int jc,
int jd,
int je,
int jf,
int jg,
int jh,
int ji) {
461 return gsl_sf_coupling_9j(ja,jb,jc,jd,je,jf,jg,jh,ji);
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 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 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.
Namespace for new Math classes and functions.
double legendre(unsigned l, unsigned m, double x)
static constexpr double s