libstdc++
Mathematical Special Functions
Collaboration diagram for Mathematical Special Functions:

Functions

template<typename _Tp >
__gnu_cxx::__promote< _Tp >::__type __gnu_cxx::airy_ai (_Tp __x)
 
float __gnu_cxx::airy_aif (float __x)
 
long double __gnu_cxx::airy_ail (long double __x)
 
template<typename _Tp >
__gnu_cxx::__promote< _Tp >::__type __gnu_cxx::airy_bi (_Tp __x)
 
float __gnu_cxx::airy_bif (float __x)
 
long double __gnu_cxx::airy_bil (long double __x)
 
template<typename _Tp >
__gnu_cxx::__promote< _Tp >::__type std::assoc_laguerre (unsigned int __n, unsigned int __m, _Tp __x)
 
float std::assoc_laguerref (unsigned int __n, unsigned int __m, float __x)
 
long double std::assoc_laguerrel (unsigned int __n, unsigned int __m, long double __x)
 
template<typename _Tp >
__gnu_cxx::__promote< _Tp >::__type std::assoc_legendre (unsigned int __l, unsigned int __m, _Tp __x)
 
float std::assoc_legendref (unsigned int __l, unsigned int __m, float __x)
 
long double std::assoc_legendrel (unsigned int __l, unsigned int __m, long double __x)
 
template<typename _Tpa , typename _Tpb >
__gnu_cxx::__promote_2< _Tpa, _Tpb >::__type std::beta (_Tpa __a, _Tpb __b)
 
float std::betaf (float __a, float __b)
 
long double std::betal (long double __a, long double __b)
 
template<typename _Tp >
__gnu_cxx::__promote< _Tp >::__type std::comp_ellint_1 (_Tp __k)
 
float std::comp_ellint_1f (float __k)
 
long double std::comp_ellint_1l (long double __k)
 
template<typename _Tp >
__gnu_cxx::__promote< _Tp >::__type std::comp_ellint_2 (_Tp __k)
 
float std::comp_ellint_2f (float __k)
 
long double std::comp_ellint_2l (long double __k)
 
template<typename _Tp , typename _Tpn >
__gnu_cxx::__promote_2< _Tp, _Tpn >::__type std::comp_ellint_3 (_Tp __k, _Tpn __nu)
 
float std::comp_ellint_3f (float __k, float __nu)
 
long double std::comp_ellint_3l (long double __k, long double __nu)
 
template<typename _Tpa , typename _Tpc , typename _Tp >
__gnu_cxx::__promote_3< _Tpa, _Tpc, _Tp >::__type __gnu_cxx::conf_hyperg (_Tpa __a, _Tpc __c, _Tp __x)
 
float __gnu_cxx::conf_hypergf (float __a, float __c, float __x)
 
long double __gnu_cxx::conf_hypergl (long double __a, long double __c, long double __x)
 
template<typename _Tpnu , typename _Tp >
__gnu_cxx::__promote_2< _Tpnu, _Tp >::__type std::cyl_bessel_i (_Tpnu __nu, _Tp __x)
 
float std::cyl_bessel_if (float __nu, float __x)
 
long double std::cyl_bessel_il (long double __nu, long double __x)
 
template<typename _Tpnu , typename _Tp >
__gnu_cxx::__promote_2< _Tpnu, _Tp >::__type std::cyl_bessel_j (_Tpnu __nu, _Tp __x)
 
float std::cyl_bessel_jf (float __nu, float __x)
 
long double std::cyl_bessel_jl (long double __nu, long double __x)
 
template<typename _Tpnu , typename _Tp >
__gnu_cxx::__promote_2< _Tpnu, _Tp >::__type std::cyl_bessel_k (_Tpnu __nu, _Tp __x)
 
float std::cyl_bessel_kf (float __nu, float __x)
 
long double std::cyl_bessel_kl (long double __nu, long double __x)
 
template<typename _Tpnu , typename _Tp >
__gnu_cxx::__promote_2< _Tpnu, _Tp >::__type std::cyl_neumann (_Tpnu __nu, _Tp __x)
 
float std::cyl_neumannf (float __nu, float __x)
 
long double std::cyl_neumannl (long double __nu, long double __x)
 
template<typename _Tp , typename _Tpp >
__gnu_cxx::__promote_2< _Tp, _Tpp >::__type std::ellint_1 (_Tp __k, _Tpp __phi)
 
float std::ellint_1f (float __k, float __phi)
 
long double std::ellint_1l (long double __k, long double __phi)
 
template<typename _Tp , typename _Tpp >
__gnu_cxx::__promote_2< _Tp, _Tpp >::__type std::ellint_2 (_Tp __k, _Tpp __phi)
 
float std::ellint_2f (float __k, float __phi)
 
long double std::ellint_2l (long double __k, long double __phi)
 
template<typename _Tp , typename _Tpn , typename _Tpp >
__gnu_cxx::__promote_3< _Tp, _Tpn, _Tpp >::__type std::ellint_3 (_Tp __k, _Tpn __nu, _Tpp __phi)
 
float std::ellint_3f (float __k, float __nu, float __phi)
 
long double std::ellint_3l (long double __k, long double __nu, long double __phi)
 
template<typename _Tp >
__gnu_cxx::__promote< _Tp >::__type std::expint (_Tp __x)
 
float std::expintf (float __x)
 
long double std::expintl (long double __x)
 
template<typename _Tp >
__gnu_cxx::__promote< _Tp >::__type std::hermite (unsigned int __n, _Tp __x)
 
float std::hermitef (unsigned int __n, float __x)
 
long double std::hermitel (unsigned int __n, long double __x)
 
template<typename _Tpa , typename _Tpb , typename _Tpc , typename _Tp >
__gnu_cxx::__promote_4< _Tpa, _Tpb, _Tpc, _Tp >::__type __gnu_cxx::hyperg (_Tpa __a, _Tpb __b, _Tpc __c, _Tp __x)
 
float __gnu_cxx::hypergf (float __a, float __b, float __c, float __x)
 
long double __gnu_cxx::hypergl (long double __a, long double __b, long double __c, long double __x)
 
template<typename _Tp >
__gnu_cxx::__promote< _Tp >::__type std::laguerre (unsigned int __n, _Tp __x)
 
float std::laguerref (unsigned int __n, float __x)
 
long double std::laguerrel (unsigned int __n, long double __x)
 
template<typename _Tp >
__gnu_cxx::__promote< _Tp >::__type std::legendre (unsigned int __l, _Tp __x)
 
float std::legendref (unsigned int __l, float __x)
 
long double std::legendrel (unsigned int __l, long double __x)
 
template<typename _Tp >
__gnu_cxx::__promote< _Tp >::__type std::riemann_zeta (_Tp __s)
 
float std::riemann_zetaf (float __s)
 
long double std::riemann_zetal (long double __s)
 
template<typename _Tp >
__gnu_cxx::__promote< _Tp >::__type std::sph_bessel (unsigned int __n, _Tp __x)
 
float std::sph_besself (unsigned int __n, float __x)
 
long double std::sph_bessell (unsigned int __n, long double __x)
 
template<typename _Tp >
__gnu_cxx::__promote< _Tp >::__type std::sph_legendre (unsigned int __l, unsigned int __m, _Tp __theta)
 
float std::sph_legendref (unsigned int __l, unsigned int __m, float __theta)
 
long double std::sph_legendrel (unsigned int __l, unsigned int __m, long double __theta)
 
template<typename _Tp >
__gnu_cxx::__promote< _Tp >::__type std::sph_neumann (unsigned int __n, _Tp __x)
 
float std::sph_neumannf (unsigned int __n, float __x)
 
long double std::sph_neumannl (unsigned int __n, long double __x)
 

Detailed Description

Mathematical Special Functions

A collection of advanced mathematical special functions, defined by ISO/IEC IS 29124 and then added to ISO C++ 2017.

Introduction and History

The first significant library upgrade on the road to C++2011, TR1, included a set of 23 mathematical functions that significantly extended the standard transcendental functions inherited from C and declared in <cmath>.

Although most components from TR1 were eventually adopted for C++11 these math functions were left behind out of concern for implementability. The math functions were published as a separate international standard IS 29124 - Extensions to the C++ Library to Support Mathematical Special Functions.

For C++17 these functions were incorporated into the main standard.

Contents

The following functions are implemented in namespace std:

The hypergeometric functions were stricken from the TR29124 and C++17 versions of this math library because of implementation concerns. However, since they were in the TR1 version and since they are popular we kept them as an extension in namespace __gnu_cxx:

Argument Promotion

The arguments suppled to the non-suffixed functions will be promoted according to the following rules:

  1. If any argument intended to be floating point is given an integral value That integral value is promoted to double.
  2. All floating point arguments are promoted up to the largest floating point precision among them.

NaN Arguments

If any of the floating point arguments supplied to these functions is invalid or NaN (std::numeric_limits<Tp>::quiet_NaN), the value NaN is returned.

Implementation

We strive to implement the underlying math with type generic algorithms to the greatest extent possible. In practice, the functions are thin wrappers that dispatch to function templates. Type dependence is controlled with std::numeric_limits and functions thereof.

We don't promote float to double or double to long double reflexively. The goal is for float functions to operate more quickly, at the cost of float accuracy and possibly a smaller domain of validity. Similaryly, long double should give you more dynamic range and slightly more pecision than double on many systems.

Testing

These functions have been tested against equivalent implementations from the Gnu Scientific Library, GSL and Boost and the ratio

\[
  \frac{|f - f_{test}|}{|f_{test}|}
\]

is generally found to be within 10-15 for 64-bit double on linux-x86_64 systems over most of the ranges of validity.

Todo:
Provide accuracy comparisons on a per-function basis for a small number of targets.

General Bibliography

See also
Abramowitz and Stegun: Handbook of Mathematical Functions, with Formulas, Graphs, and Mathematical Tables Edited by Milton Abramowitz and Irene A. Stegun, National Bureau of Standards Applied Mathematics Series - 55 Issued June 1964, Tenth Printing, December 1972, with corrections Electronic versions of A&S abound including both pdf and navigable html.
for example http://people.math.sfu.ca/~cbm/aands/
The old A&S has been redone as the NIST Digital Library of Mathematical Functions: http://dlmf.nist.gov/ This version is far more navigable and includes more recent work.
An Atlas of Functions: with Equator, the Atlas Function Calculator 2nd Edition, by Oldham, Keith B., Myland, Jan, Spanier, Jerome
Asymptotics and Special Functions by Frank W. J. Olver, Academic Press, 1974
Numerical Recipes in C, The Art of Scientific Computing, by William H. Press, Second Ed., Saul A. Teukolsky, William T. Vetterling, and Brian P. Flannery, Cambridge University Press, 1992
The Special Functions and Their Approximations: Volumes 1 and 2, by Yudell L. Luke, Academic Press, 1969

Function Documentation

◆ airy_ai()

template<typename _Tp >
__gnu_cxx::__promote< _Tp >::__type __gnu_cxx::airy_ai ( _Tp  __x)
inline

Return the Airy function $ Ai(x) $ of real argument x.

Definition at line 1240 of file specfun.h.

◆ airy_aif()

float __gnu_cxx::airy_aif ( float  __x)
inline

Return the Airy function $ Ai(x) $ of float argument x.

Definition at line 1217 of file specfun.h.

◆ airy_ail()

long double __gnu_cxx::airy_ail ( long double  __x)
inline

Return the Airy function $ Ai(x) $ of long double argument x.

Definition at line 1228 of file specfun.h.

◆ airy_bi()

template<typename _Tp >
__gnu_cxx::__promote< _Tp >::__type __gnu_cxx::airy_bi ( _Tp  __x)
inline

Return the Airy function $ Bi(x) $ of real argument x.

Definition at line 1275 of file specfun.h.

◆ airy_bif()

float __gnu_cxx::airy_bif ( float  __x)
inline

Return the Airy function $ Bi(x) $ of float argument x.

Definition at line 1252 of file specfun.h.

◆ airy_bil()

long double __gnu_cxx::airy_bil ( long double  __x)
inline

Return the Airy function $ Bi(x) $ of long double argument x.

Definition at line 1263 of file specfun.h.

◆ assoc_laguerre()

template<typename _Tp >
__gnu_cxx::__promote< _Tp >::__type std::assoc_laguerre ( unsigned int  __n,
unsigned int  __m,
_Tp  __x 
)
inline

Return the associated Laguerre polynomial of nonnegative order n, nonnegative degree m and real argument x: $ L_n^m(x) $.

The associated Laguerre function of real degree $ \alpha $, $ L_n^\alpha(x) $, is defined by

\[
         L_n^\alpha(x) = \frac{(\alpha + 1)_n}{n!}
                         {}_1F_1(-n; \alpha + 1; x)
\]

where $ (\alpha)_n $ is the Pochhammer symbol and $ {}_1F_1(a; c; x) $ is the confluent hypergeometric function.

The associated Laguerre polynomial is defined for integral degree $ \alpha = m $ by:

\[
         L_n^m(x) = (-1)^m \frac{d^m}{dx^m} L_{n + m}(x)
\]

where the Laguerre polynomial is defined by:

\[
         L_n(x) = \frac{e^x}{n!} \frac{d^n}{dx^n} (x^ne^{-x})
\]

and $ x >= 0 $.

See also
laguerre for details of the Laguerre function of degree n
Template Parameters
_TpThe floating-point type of the argument __x.
Parameters
__nThe order of the Laguerre function, __n >= 0.
__mThe degree of the Laguerre function, __m >= 0.
__xThe argument of the Laguerre function, __x >= 0.
Exceptions
std::domain_errorif __x < 0.

Definition at line 250 of file specfun.h.

◆ assoc_laguerref()

float std::assoc_laguerref ( unsigned int  __n,
unsigned int  __m,
float  __x 
)
inline

Return the associated Laguerre polynomial of order n, degree m: $ L_n^m(x) $ for float argument.

See also
assoc_laguerre for more details.

Definition at line 204 of file specfun.h.

◆ assoc_laguerrel()

long double std::assoc_laguerrel ( unsigned int  __n,
unsigned int  __m,
long double  __x 
)
inline

Return the associated Laguerre polynomial of order n, degree m: $ L_n^m(x) $.

See also
assoc_laguerre for more details.

Definition at line 214 of file specfun.h.

◆ assoc_legendre()

template<typename _Tp >
__gnu_cxx::__promote< _Tp >::__type std::assoc_legendre ( unsigned int  __l,
unsigned int  __m,
_Tp  __x 
)
inline

Return the associated Legendre function of degree l and order m.

The associated Legendre function is derived from the Legendre function $ P_l(x) $ by the Rodrigues formula:

\[
  P_l^m(x) = (1 - x^2)^{m/2}\frac{d^m}{dx^m}P_l(x)
\]

See also
legendre for details of the Legendre function of degree l
Template Parameters
_TpThe floating-point type of the argument __x.
Parameters
__lThe degree __l >= 0.
__mThe order __m <= l.
__xThe argument, abs(__x) <= 1.
Exceptions
std::domain_errorif abs(__x) > 1.

Definition at line 296 of file specfun.h.

◆ assoc_legendref()

float std::assoc_legendref ( unsigned int  __l,
unsigned int  __m,
float  __x 
)
inline

Return the associated Legendre function of degree l and order m for float argument.

See also
assoc_legendre for more details.

Definition at line 265 of file specfun.h.

◆ assoc_legendrel()

long double std::assoc_legendrel ( unsigned int  __l,
unsigned int  __m,
long double  __x 
)
inline

Return the associated Legendre function of degree l and order m.

See also
assoc_legendre for more details.

Definition at line 274 of file specfun.h.

◆ beta()

template<typename _Tpa , typename _Tpb >
__gnu_cxx::__promote_2< _Tpa, _Tpb >::__type std::beta ( _Tpa  __a,
_Tpb  __b 
)
inline

Return the beta function, $B(a,b)$, for real parameters a, b.

The beta function is defined by

\[
  B(a,b) = \int_0^1 t^{a - 1} (1 - t)^{b - 1} dt
         = \frac{\Gamma(a)\Gamma(b)}{\Gamma(a+b)}
\]

where $ a > 0 $ and $ b > 0 $

Template Parameters
_TpaThe floating-point type of the parameter __a.
_TpbThe floating-point type of the parameter __b.
Parameters
__aThe first argument of the beta function, __a > 0 .
__bThe second argument of the beta function, __b > 0 .
Exceptions
std::domain_errorif __a < 0 or __b < 0 .

Definition at line 341 of file specfun.h.

◆ betaf()

float std::betaf ( float  __a,
float  __b 
)
inline

Return the beta function, $ B(a,b) $, for float parameters a, b.

See also
beta for more details.

Definition at line 310 of file specfun.h.

◆ betal()

long double std::betal ( long double  __a,
long double  __b 
)
inline

Return the beta function, $B(a,b)$, for long double parameters a, b.

See also
beta for more details.

Definition at line 320 of file specfun.h.

◆ comp_ellint_1()

template<typename _Tp >
__gnu_cxx::__promote< _Tp >::__type std::comp_ellint_1 ( _Tp  __k)
inline

Return the complete elliptic integral of the first kind $ K(k) $ for real modulus k.

The complete elliptic integral of the first kind is defined as

\[
  K(k) = F(k,\pi/2) = \int_0^{\pi/2}\frac{d\theta}
                                             {\sqrt{1 - k^2 sin^2\theta}}
\]

where $ F(k,\phi) $ is the incomplete elliptic integral of the first kind and the modulus $ |k| <= 1 $.

See also
ellint_1 for details of the incomplete elliptic function of the first kind.
Template Parameters
_TpThe floating-point type of the modulus __k.
Parameters
__kThe modulus, abs(__k) <= 1
Exceptions
std::domain_errorif abs(__k) > 1 .

Definition at line 389 of file specfun.h.

◆ comp_ellint_1f()

float std::comp_ellint_1f ( float  __k)
inline

Return the complete elliptic integral of the first kind $ E(k) $ for float modulus k.

See also
comp_ellint_1 for details.

Definition at line 356 of file specfun.h.

◆ comp_ellint_1l()

long double std::comp_ellint_1l ( long double  __k)
inline

Return the complete elliptic integral of the first kind $ E(k) $ for long double modulus k.

See also
comp_ellint_1 for details.

Definition at line 366 of file specfun.h.

◆ comp_ellint_2()

template<typename _Tp >
__gnu_cxx::__promote< _Tp >::__type std::comp_ellint_2 ( _Tp  __k)
inline

Return the complete elliptic integral of the second kind $ E(k) $ for real modulus k.

The complete elliptic integral of the second kind is defined as

\[
  E(k) = E(k,\pi/2) = \int_0^{\pi/2}\sqrt{1 - k^2 sin^2\theta}
\]

where $ E(k,\phi) $ is the incomplete elliptic integral of the second kind and the modulus $ |k| <= 1 $.

See also
ellint_2 for details of the incomplete elliptic function of the second kind.
Template Parameters
_TpThe floating-point type of the modulus __k.
Parameters
__kThe modulus, abs(__k) <= 1
Exceptions
std::domain_errorif abs(__k) > 1.

Definition at line 436 of file specfun.h.

◆ comp_ellint_2f()

float std::comp_ellint_2f ( float  __k)
inline

Return the complete elliptic integral of the second kind $ E(k) $ for float modulus k.

See also
comp_ellint_2 for details.

Definition at line 404 of file specfun.h.

◆ comp_ellint_2l()

long double std::comp_ellint_2l ( long double  __k)
inline

Return the complete elliptic integral of the second kind $ E(k) $ for long double modulus k.

See also
comp_ellint_2 for details.

Definition at line 414 of file specfun.h.

◆ comp_ellint_3()

template<typename _Tp , typename _Tpn >
__gnu_cxx::__promote_2< _Tp, _Tpn >::__type std::comp_ellint_3 ( _Tp  __k,
_Tpn  __nu 
)
inline

Return the complete elliptic integral of the third kind $ \Pi(k,\nu) = \Pi(k,\nu,\pi/2) $ for real modulus k.

The complete elliptic integral of the third kind is defined as

\[
  \Pi(k,\nu) = \Pi(k,\nu,\pi/2) = \int_0^{\pi/2}
                     \frac{d\theta}
                   {(1 - \nu \sin^2\theta)\sqrt{1 - k^2 \sin^2\theta}}
\]

where $ \Pi(k,\nu,\phi) $ is the incomplete elliptic integral of the second kind and the modulus $ |k| <= 1 $.

See also
ellint_3 for details of the incomplete elliptic function of the third kind.
Template Parameters
_TpThe floating-point type of the modulus __k.
_TpnThe floating-point type of the argument __nu.
Parameters
__kThe modulus, abs(__k) <= 1
__nuThe argument
Exceptions
std::domain_errorif abs(__k) > 1.

Definition at line 487 of file specfun.h.

◆ comp_ellint_3f()

float std::comp_ellint_3f ( float  __k,
float  __nu 
)
inline

Return the complete elliptic integral of the third kind $ \Pi(k,\nu) $ for float modulus k.

See also
comp_ellint_3 for details.

Definition at line 451 of file specfun.h.

◆ comp_ellint_3l()

long double std::comp_ellint_3l ( long double  __k,
long double  __nu 
)
inline

Return the complete elliptic integral of the third kind $ \Pi(k,\nu) $ for long double modulus k.

See also
comp_ellint_3 for details.

Definition at line 461 of file specfun.h.

◆ conf_hyperg()

template<typename _Tpa , typename _Tpc , typename _Tp >
__gnu_cxx::__promote_3< _Tpa, _Tpc, _Tp >::__type __gnu_cxx::conf_hyperg ( _Tpa  __a,
_Tpc  __c,
_Tp  __x 
)
inline

Return the confluent hypergeometric function $ {}_1F_1(a;c;x) $ of real numeratorial parameter a, denominatorial parameter c, and argument x.

The confluent hypergeometric function is defined by

\[
   {}_1F_1(a;c;x) = \sum_{n=0}^{\infty} \frac{(a)_n x^n}{(c)_n n!}
\]

where the Pochhammer symbol is $ (x)_k = (x)(x+1)...(x+k-1) $, $ (x)_0 = 1 $

Parameters
__aThe numeratorial parameter
__cThe denominatorial parameter
__xThe argument

Definition at line 1325 of file specfun.h.

◆ conf_hypergf()

float __gnu_cxx::conf_hypergf ( float  __a,
float  __c,
float  __x 
)
inline

Return the confluent hypergeometric function $ {}_1F_1(a;c;x) $ of float numeratorial parameter a, denominatorial parameter c, and argument x.

See also
conf_hyperg for details.

Definition at line 1293 of file specfun.h.

◆ conf_hypergl()

long double __gnu_cxx::conf_hypergl ( long double  __a,
long double  __c,
long double  __x 
)
inline

Return the confluent hypergeometric function $ {}_1F_1(a;c;x) $ of long double numeratorial parameter a, denominatorial parameter c, and argument x.

See also
conf_hyperg for details.

Definition at line 1304 of file specfun.h.

◆ cyl_bessel_i()

template<typename _Tpnu , typename _Tp >
__gnu_cxx::__promote_2< _Tpnu, _Tp >::__type std::cyl_bessel_i ( _Tpnu  __nu,
_Tp  __x 
)
inline

Return the regular modified Bessel function $ I_{\nu}(x) $ for real order $ \nu $ and argument $ x >= 0 $.

The regular modified cylindrical Bessel function is:

\[
 I_{\nu}(x) = i^{-\nu}J_\nu(ix) = \sum_{k=0}^{\infty}
                \frac{(x/2)^{\nu + 2k}}{k!\Gamma(\nu+k+1)}
\]

Template Parameters
_TpnuThe floating-point type of the order __nu.
_TpThe floating-point type of the argument __x.
Parameters
__nuThe order
__xThe argument, __x >= 0
Exceptions
std::domain_errorif __x < 0 .

Definition at line 533 of file specfun.h.

◆ cyl_bessel_if()

float std::cyl_bessel_if ( float  __nu,
float  __x 
)
inline

Return the regular modified Bessel function $ I_{\nu}(x) $ for float order $ \nu $ and argument $ x >= 0 $.

See also
cyl_bessel_i for setails.

Definition at line 502 of file specfun.h.

◆ cyl_bessel_il()

long double std::cyl_bessel_il ( long double  __nu,
long double  __x 
)
inline

Return the regular modified Bessel function $ I_{\nu}(x) $ for long double order $ \nu $ and argument $ x >= 0 $.

See also
cyl_bessel_i for setails.

Definition at line 512 of file specfun.h.

◆ cyl_bessel_j()

template<typename _Tpnu , typename _Tp >
__gnu_cxx::__promote_2< _Tpnu, _Tp >::__type std::cyl_bessel_j ( _Tpnu  __nu,
_Tp  __x 
)
inline

Return the Bessel function $ J_{\nu}(x) $ of real order $ \nu $ and argument $ x >= 0 $.

The cylindrical Bessel function is:

\[
   J_{\nu}(x) = \sum_{k=0}^{\infty}
             \frac{(-1)^k (x/2)^{\nu + 2k}}{k!\Gamma(\nu+k+1)}
\]

Template Parameters
_TpnuThe floating-point type of the order __nu.
_TpThe floating-point type of the argument __x.
Parameters
__nuThe order
__xThe argument, __x >= 0
Exceptions
std::domain_errorif __x < 0 .

Definition at line 579 of file specfun.h.

◆ cyl_bessel_jf()

float std::cyl_bessel_jf ( float  __nu,
float  __x 
)
inline

Return the Bessel function of the first kind $ J_{\nu}(x) $ for float order $ \nu $ and argument $ x >= 0 $.

See also
cyl_bessel_j for setails.

Definition at line 548 of file specfun.h.

◆ cyl_bessel_jl()

long double std::cyl_bessel_jl ( long double  __nu,
long double  __x 
)
inline

Return the Bessel function of the first kind $ J_{\nu}(x) $ for long double order $ \nu $ and argument $ x >= 0 $.

See also
cyl_bessel_j for setails.

Definition at line 558 of file specfun.h.

◆ cyl_bessel_k()

template<typename _Tpnu , typename _Tp >
__gnu_cxx::__promote_2< _Tpnu, _Tp >::__type std::cyl_bessel_k ( _Tpnu  __nu,
_Tp  __x 
)
inline

Return the irregular modified Bessel function $ K_{\nu}(x) $ of real order $ \nu $ and argument $ x $.

The irregular modified Bessel function is defined by:

\[
        K_{\nu}(x) = \frac{\pi}{2}
                     \frac{I_{-\nu}(x) - I_{\nu}(x)}{\sin \nu\pi}
\]

where for integral $ \nu = n $ a limit is taken: $ lim_{\nu \to n} $. For negative argument we have simply:

\[
        K_{-\nu}(x) = K_{\nu}(x)
\]

Template Parameters
_TpnuThe floating-point type of the order __nu.
_TpThe floating-point type of the argument __x.
Parameters
__nuThe order
__xThe argument, __x >= 0
Exceptions
std::domain_errorif __x < 0 .

Definition at line 631 of file specfun.h.

◆ cyl_bessel_kf()

float std::cyl_bessel_kf ( float  __nu,
float  __x 
)
inline

Return the irregular modified Bessel function $ K_{\nu}(x) $ for float order $ \nu $ and argument $ x >= 0 $.

See also
cyl_bessel_k for setails.

Definition at line 594 of file specfun.h.

◆ cyl_bessel_kl()

long double std::cyl_bessel_kl ( long double  __nu,
long double  __x 
)
inline

Return the irregular modified Bessel function $ K_{\nu}(x) $ for long double order $ \nu $ and argument $ x >= 0 $.

See also
cyl_bessel_k for setails.

Definition at line 604 of file specfun.h.

◆ cyl_neumann()

template<typename _Tpnu , typename _Tp >
__gnu_cxx::__promote_2< _Tpnu, _Tp >::__type std::cyl_neumann ( _Tpnu  __nu,
_Tp  __x 
)
inline

Return the Neumann function $ N_{\nu}(x) $ of real order $ \nu $ and argument $ x >= 0 $.

The Neumann function is defined by:

\[
   N_{\nu}(x) = \frac{J_{\nu}(x) \cos \nu\pi - J_{-\nu}(x)}
                     {\sin \nu\pi}
\]

where $ x >= 0 $ and for integral order $ \nu = n $ a limit is taken: $ lim_{\nu \to n} $.

Template Parameters
_TpnuThe floating-point type of the order __nu.
_TpThe floating-point type of the argument __x.
Parameters
__nuThe order
__xThe argument, __x >= 0
Exceptions
std::domain_errorif __x < 0 .

Definition at line 679 of file specfun.h.

◆ cyl_neumannf()

float std::cyl_neumannf ( float  __nu,
float  __x 
)
inline

Return the Neumann function $ N_{\nu}(x) $ of float order $ \nu $ and argument $ x $.

See also
cyl_neumann for setails.

Definition at line 646 of file specfun.h.

◆ cyl_neumannl()

long double std::cyl_neumannl ( long double  __nu,
long double  __x 
)
inline

Return the Neumann function $ N_{\nu}(x) $ of long double order $ \nu $ and argument $ x $.

See also
cyl_neumann for setails.

Definition at line 656 of file specfun.h.

◆ ellint_1()

template<typename _Tp , typename _Tpp >
__gnu_cxx::__promote_2< _Tp, _Tpp >::__type std::ellint_1 ( _Tp  __k,
_Tpp  __phi 
)
inline

Return the incomplete elliptic integral of the first kind $ F(k,\phi) $ for real modulus $ k $ and angle $ \phi $.

The incomplete elliptic integral of the first kind is defined as

\[
  F(k,\phi) = \int_0^{\phi}\frac{d\theta}
                                     {\sqrt{1 - k^2 sin^2\theta}}
\]

For $ \phi= \pi/2 $ this becomes the complete elliptic integral of the first kind, $ K(k) $.

See also
comp_ellint_1.
Template Parameters
_TpThe floating-point type of the modulus __k.
_TppThe floating-point type of the angle __phi.
Parameters
__kThe modulus, abs(__k) <= 1
__phiThe integral limit argument in radians
Exceptions
std::domain_errorif abs(__k) > 1 .

Definition at line 727 of file specfun.h.

◆ ellint_1f()

float std::ellint_1f ( float  __k,
float  __phi 
)
inline

Return the incomplete elliptic integral of the first kind $ E(k,\phi) $ for float modulus $ k $ and angle $ \phi $.

See also
ellint_1 for details.

Definition at line 694 of file specfun.h.

◆ ellint_1l()

long double std::ellint_1l ( long double  __k,
long double  __phi 
)
inline

Return the incomplete elliptic integral of the first kind $ E(k,\phi) $ for long double modulus $ k $ and angle $ \phi $.

See also
ellint_1 for details.

Definition at line 704 of file specfun.h.

◆ ellint_2()

template<typename _Tp , typename _Tpp >
__gnu_cxx::__promote_2< _Tp, _Tpp >::__type std::ellint_2 ( _Tp  __k,
_Tpp  __phi 
)
inline

Return the incomplete elliptic integral of the second kind $ E(k,\phi) $.

The incomplete elliptic integral of the second kind is defined as

\[
  E(k,\phi) = \int_0^{\phi} \sqrt{1 - k^2 sin^2\theta}
\]

For $ \phi= \pi/2 $ this becomes the complete elliptic integral of the second kind, $ E(k) $.

See also
comp_ellint_2.
Template Parameters
_TpThe floating-point type of the modulus __k.
_TppThe floating-point type of the angle __phi.
Parameters
__kThe modulus, abs(__k) <= 1
__phiThe integral limit argument in radians
Returns
The elliptic function of the second kind.
Exceptions
std::domain_errorif abs(__k) > 1 .

Definition at line 775 of file specfun.h.

◆ ellint_2f()

float std::ellint_2f ( float  __k,
float  __phi 
)
inline

Return the incomplete elliptic integral of the second kind $ E(k,\phi) $ for float argument.

See also
ellint_2 for details.

Definition at line 742 of file specfun.h.

◆ ellint_2l()

long double std::ellint_2l ( long double  __k,
long double  __phi 
)
inline

Return the incomplete elliptic integral of the second kind $ E(k,\phi) $.

See also
ellint_2 for details.

Definition at line 752 of file specfun.h.

◆ ellint_3()

template<typename _Tp , typename _Tpn , typename _Tpp >
__gnu_cxx::__promote_3< _Tp, _Tpn, _Tpp >::__type std::ellint_3 ( _Tp  __k,
_Tpn  __nu,
_Tpp  __phi 
)
inline

Return the incomplete elliptic integral of the third kind $ \Pi(k,\nu,\phi) $.

The incomplete elliptic integral of the third kind is defined by:

\[
  \Pi(k,\nu,\phi) = \int_0^{\phi}
                         \frac{d\theta}
                         {(1 - \nu \sin^2\theta)
                          \sqrt{1 - k^2 \sin^2\theta}}
\]

For $ \phi= \pi/2 $ this becomes the complete elliptic integral of the third kind, $ \Pi(k,\nu) $.

See also
comp_ellint_3.
Template Parameters
_TpThe floating-point type of the modulus __k.
_TpnThe floating-point type of the argument __nu.
_TppThe floating-point type of the angle __phi.
Parameters
__kThe modulus, abs(__k) <= 1
__nuThe second argument
__phiThe integral limit argument in radians
Returns
The elliptic function of the third kind.
Exceptions
std::domain_errorif abs(__k) > 1 .

Definition at line 828 of file specfun.h.

◆ ellint_3f()

float std::ellint_3f ( float  __k,
float  __nu,
float  __phi 
)
inline

Return the incomplete elliptic integral of the third kind $ \Pi(k,\nu,\phi) $ for float argument.

See also
ellint_3 for details.

Definition at line 790 of file specfun.h.

◆ ellint_3l()

long double std::ellint_3l ( long double  __k,
long double  __nu,
long double  __phi 
)
inline

Return the incomplete elliptic integral of the third kind $ \Pi(k,\nu,\phi) $.

See also
ellint_3 for details.

Definition at line 800 of file specfun.h.

◆ expint()

template<typename _Tp >
__gnu_cxx::__promote< _Tp >::__type std::expint ( _Tp  __x)
inline

Return the exponential integral $ Ei(x) $ for real argument x.

The exponential integral is given by

\[
  Ei(x) = -\int_{-x}^\infty \frac{e^t}{t} dt
\]

Template Parameters
_TpThe floating-point type of the argument __x.
Parameters
__xThe argument of the exponential integral function.

Definition at line 868 of file specfun.h.

◆ expintf()

float std::expintf ( float  __x)
inline

Return the exponential integral $ Ei(x) $ for float argument x.

See also
expint for details.

Definition at line 842 of file specfun.h.

◆ expintl()

long double std::expintl ( long double  __x)
inline

Return the exponential integral $ Ei(x) $ for long double argument x.

See also
expint for details.

Definition at line 852 of file specfun.h.

◆ hermite()

template<typename _Tp >
__gnu_cxx::__promote< _Tp >::__type std::hermite ( unsigned int  __n,
_Tp  __x 
)
inline

Return the Hermite polynomial $ H_n(x) $ of order n and real argument x.

The Hermite polynomial is defined by:

\[
  H_n(x) = (-1)^n e^{x^2} \frac{d^n}{dx^n} e^{-x^2}
\]

The Hermite polynomial obeys a reflection formula:

\[
  H_n(-x) = (-1)^n H_n(x)
\]

Template Parameters
_TpThe floating-point type of the argument __x.
Parameters
__nThe order
__xThe argument

Definition at line 916 of file specfun.h.

◆ hermitef()

float std::hermitef ( unsigned int  __n,
float  __x 
)
inline

Return the Hermite polynomial $ H_n(x) $ of nonnegative order n and float argument x.

See also
hermite for details.

Definition at line 883 of file specfun.h.

◆ hermitel()

long double std::hermitel ( unsigned int  __n,
long double  __x 
)
inline

Return the Hermite polynomial $ H_n(x) $ of nonnegative order n and long double argument x.

See also
hermite for details.

Definition at line 893 of file specfun.h.

◆ hyperg()

template<typename _Tpa , typename _Tpb , typename _Tpc , typename _Tp >
__gnu_cxx::__promote_4< _Tpa, _Tpb, _Tpc, _Tp >::__type __gnu_cxx::hyperg ( _Tpa  __a,
_Tpb  __b,
_Tpc  __c,
_Tp  __x 
)
inline

Return the hypergeometric function $ {}_2F_1(a,b;c;x) $ of real numeratorial parameters a and b, denominatorial parameter c, and argument x.

The hypergeometric function is defined by

\[
   {}_2F_1(a;c;x) = \sum_{n=0}^{\infty} \frac{(a)_n (b)_n x^n}{(c)_n n!}
\]

where the Pochhammer symbol is $ (x)_k = (x)(x+1)...(x+k-1) $, $ (x)_0 = 1 $

Parameters
__aThe first numeratorial parameter
__bThe second numeratorial parameter
__cThe denominatorial parameter
__xThe argument

Definition at line 1374 of file specfun.h.

◆ hypergf()

float __gnu_cxx::hypergf ( float  __a,
float  __b,
float  __c,
float  __x 
)
inline

Return the hypergeometric function $ {}_2F_1(a,b;c;x) $ of @ float numeratorial parameters a and b, denominatorial parameter c, and argument x.

See also
hyperg for details.

Definition at line 1341 of file specfun.h.

◆ hypergl()

long double __gnu_cxx::hypergl ( long double  __a,
long double  __b,
long double  __c,
long double  __x 
)
inline

Return the hypergeometric function $ {}_2F_1(a,b;c;x) $ of long double numeratorial parameters a and b, denominatorial parameter c, and argument x.

See also
hyperg for details.

Definition at line 1352 of file specfun.h.

◆ laguerre()

template<typename _Tp >
__gnu_cxx::__promote< _Tp >::__type std::laguerre ( unsigned int  __n,
_Tp  __x 
)
inline

Returns the Laguerre polynomial $ L_n(x) $ of nonnegative degree n and real argument $ x >= 0 $.

The Laguerre polynomial is defined by:

\[
         L_n(x) = \frac{e^x}{n!} \frac{d^n}{dx^n} (x^ne^{-x})
\]

Template Parameters
_TpThe floating-point type of the argument __x.
Parameters
__nThe nonnegative order
__xThe argument __x >= 0
Exceptions
std::domain_errorif __x < 0 .

Definition at line 960 of file specfun.h.

◆ laguerref()

float std::laguerref ( unsigned int  __n,
float  __x 
)
inline

Returns the Laguerre polynomial $ L_n(x) $ of nonnegative degree n and float argument $ x >= 0 $.

See also
laguerre for more details.

Definition at line 931 of file specfun.h.

◆ laguerrel()

long double std::laguerrel ( unsigned int  __n,
long double  __x 
)
inline

Returns the Laguerre polynomial $ L_n(x) $ of nonnegative degree n and long double argument $ x >= 0 $.

See also
laguerre for more details.

Definition at line 941 of file specfun.h.

◆ legendre()

template<typename _Tp >
__gnu_cxx::__promote< _Tp >::__type std::legendre ( unsigned int  __l,
_Tp  __x 
)
inline

Return the Legendre polynomial $ P_l(x) $ of nonnegative degree $ l $ and real argument $ |x| <= 0 $.

The Legendre function of order $ l $ and argument $ x $, $ P_l(x) $, is defined by:

\[
  P_l(x) = \frac{1}{2^l l!}\frac{d^l}{dx^l}(x^2 - 1)^{l}
\]

Template Parameters
_TpThe floating-point type of the argument __x.
Parameters
__lThe degree $ l >= 0 $
__xThe argument abs(__x) <= 1
Exceptions
std::domain_errorif abs(__x) > 1

Definition at line 1005 of file specfun.h.

◆ legendref()

float std::legendref ( unsigned int  __l,
float  __x 
)
inline

Return the Legendre polynomial $ P_l(x) $ of nonnegative degree $ l $ and float argument $ |x| <= 0 $.

See also
legendre for more details.

Definition at line 975 of file specfun.h.

◆ legendrel()

long double std::legendrel ( unsigned int  __l,
long double  __x 
)
inline

Return the Legendre polynomial $ P_l(x) $ of nonnegative degree $ l $ and long double argument $ |x| <= 0 $.

See also
legendre for more details.

Definition at line 985 of file specfun.h.

◆ riemann_zeta()

template<typename _Tp >
__gnu_cxx::__promote< _Tp >::__type std::riemann_zeta ( _Tp  __s)
inline

Return the Riemann zeta function $ \zeta(s) $ for real argument $ s $.

The Riemann zeta function is defined by:

\[
        \zeta(s) = \sum_{k=1}^{\infty} k^{-s} \hbox{ for } s > 1
\]

and

\[
        \zeta(s) = \frac{1}{1-2^{1-s}}\sum_{k=1}^{\infty}(-1)^{k-1}k^{-s}
             \hbox{ for } 0 <= s <= 1
\]

For s < 1 use the reflection formula:

\[
        \zeta(s) = 2^s \pi^{s-1} \sin(\frac{\pi s}{2}) \Gamma(1-s) \zeta(1-s)
\]

Template Parameters
_TpThe floating-point type of the argument __s.
Parameters
__sThe argument s != 1

Definition at line 1056 of file specfun.h.

◆ riemann_zetaf()

float std::riemann_zetaf ( float  __s)
inline

Return the Riemann zeta function $ \zeta(s) $ for float argument $ s $.

See also
riemann_zeta for more details.

Definition at line 1020 of file specfun.h.

◆ riemann_zetal()

long double std::riemann_zetal ( long double  __s)
inline

Return the Riemann zeta function $ \zeta(s) $ for long double argument $ s $.

See also
riemann_zeta for more details.

Definition at line 1030 of file specfun.h.

◆ sph_bessel()

template<typename _Tp >
__gnu_cxx::__promote< _Tp >::__type std::sph_bessel ( unsigned int  __n,
_Tp  __x 
)
inline

Return the spherical Bessel function $ j_n(x) $ of nonnegative order n and real argument $ x >= 0 $.

The spherical Bessel function is defined by:

\[
 j_n(x) = \left(\frac{\pi}{2x} \right) ^{1/2} J_{n+1/2}(x)
\]

Template Parameters
_TpThe floating-point type of the argument __x.
Parameters
__nThe integral order n >= 0
__xThe real argument x >= 0
Exceptions
std::domain_errorif __x < 0 .

Definition at line 1100 of file specfun.h.

◆ sph_besself()

float std::sph_besself ( unsigned int  __n,
float  __x 
)
inline

Return the spherical Bessel function $ j_n(x) $ of nonnegative order n and float argument $ x >= 0 $.

See also
sph_bessel for more details.

Definition at line 1071 of file specfun.h.

◆ sph_bessell()

long double std::sph_bessell ( unsigned int  __n,
long double  __x 
)
inline

Return the spherical Bessel function $ j_n(x) $ of nonnegative order n and long double argument $ x >= 0 $.

See also
sph_bessel for more details.

Definition at line 1081 of file specfun.h.

◆ sph_legendre()

template<typename _Tp >
__gnu_cxx::__promote< _Tp >::__type std::sph_legendre ( unsigned int  __l,
unsigned int  __m,
_Tp  __theta 
)
inline

Return the spherical Legendre function of nonnegative integral degree l and order m and real angle $ \theta $ in radians.

The spherical Legendre function is defined by

\[
 Y_l^m(\theta,\phi) = (-1)^m[\frac{(2l+1)}{4\pi}
                             \frac{(l-m)!}{(l+m)!}]
                  P_l^m(\cos\theta) \exp^{im\phi}
\]

Template Parameters
_TpThe floating-point type of the angle __theta.
Parameters
__lThe order __l >= 0
__mThe degree __m >= 0 and __m <= __l
__thetaThe radian polar angle argument

Definition at line 1147 of file specfun.h.

◆ sph_legendref()

float std::sph_legendref ( unsigned int  __l,
unsigned int  __m,
float  __theta 
)
inline

Return the spherical Legendre function of nonnegative integral degree l and order m and float angle $ \theta $ in radians.

See also
sph_legendre for details.

Definition at line 1115 of file specfun.h.

◆ sph_legendrel()

long double std::sph_legendrel ( unsigned int  __l,
unsigned int  __m,
long double  __theta 
)
inline

Return the spherical Legendre function of nonnegative integral degree l and order m and long double angle $ \theta $ in radians.

See also
sph_legendre for details.

Definition at line 1126 of file specfun.h.

◆ sph_neumann()

template<typename _Tp >
__gnu_cxx::__promote< _Tp >::__type std::sph_neumann ( unsigned int  __n,
_Tp  __x 
)
inline

Return the spherical Neumann function of integral order $ n >= 0 $ and real argument $ x >= 0 $.

The spherical Neumann function is defined by

\[
   n_n(x) = \left(\frac{\pi}{2x} \right) ^{1/2} N_{n+1/2}(x)
\]

Template Parameters
_TpThe floating-point type of the argument __x.
Parameters
__nThe integral order n >= 0
__xThe real argument __x >= 0
Exceptions
std::domain_errorif __x < 0 .

Definition at line 1191 of file specfun.h.

◆ sph_neumannf()

float std::sph_neumannf ( unsigned int  __n,
float  __x 
)
inline

Return the spherical Neumann function of integral order $ n >= 0 $ and float argument $ x >= 0 $.

See also
sph_neumann for details.

Definition at line 1162 of file specfun.h.

◆ sph_neumannl()

long double std::sph_neumannl ( unsigned int  __n,
long double  __x 
)
inline

Return the spherical Neumann function of integral order $ n >= 0 $ and long double $ x >= 0 $.

See also
sph_neumann for details.

Definition at line 1172 of file specfun.h.