|
libstdc++
|
Functions | |
| 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 _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 _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) |
A collection of advanced mathematical special functions, defined by ISO/IEC IS 29124.
|
inline |
Return the associated Laguerre polynomial of nonnegative order n, nonnegative degree m and real argument x:
.
The associated Laguerre function of real degree
,
, is defined by
where
is the Pochhammer symbol and
is the confluent hypergeometric function.
The associated Laguerre polynomial is defined for integral degree
by:
where the Laguerre polynomial is defined by:
and
.
n | _Tp | The floating-point type of the argument __x. |
| __n | The order of the Laguerre function, __n >= 0. |
| __m | The degree of the Laguerre function, __m >= 0. |
| __x | The argument of the Laguerre function, __x >= 0. |
| std::domain_error | if __x < 0. |
|
inline |
|
inline |
|
inline |
Return the associated Legendre function of degree l and order m.
The associated Legendre function is derived from the Legendre function
by the Rodrigues formula:
l | _Tp | The floating-point type of the argument __x. |
| __l | The degree __l >= 0. |
| __m | The order __m <= l. |
| __x | The argument, abs(__x) <= 1. |
| std::domain_error | if abs(__x) > 1. |
|
inline |
|
inline |
|
inline |
Return the beta function,
, for real parameters a, b.
The beta function is defined by
where
and 
| _Tpa | The floating-point type of the parameter __a. |
| _Tpb | The floating-point type of the parameter __b. |
| __a | The first argument of the beta function, __a > 0 . |
| __b | The second argument of the beta function, __b > 0 . |
| std::domain_error | if __a < 0 or __b < 0 . |
|
inline |
|
inline |
|
inline |
Return the complete elliptic integral of the first kind
for real modulus k.
The complete elliptic integral of the first kind is defined as
where
is the incomplete elliptic integral of the first kind and the modulus
.
| _Tp | The floating-point type of the modulus __k. |
| __k | The modulus, abs(__k) <= 1 |
| std::domain_error | if abs(__k) > 1 . |
|
inline |
|
inline |
|
inline |
Return the complete elliptic integral of the second kind
for real modulus k.
The complete elliptic integral of the second kind is defined as
where
is the incomplete elliptic integral of the second kind and the modulus
.
| _Tp | The floating-point type of the modulus __k. |
| __k | The modulus, abs(__k) <= 1 |
| std::domain_error | if abs(__k) > 1. |
|
inline |
|
inline |
|
inline |
Return the complete elliptic integral of the third kind
for real modulus k.
The complete elliptic integral of the third kind is defined as
where
is the incomplete elliptic integral of the second kind and the modulus
.
| _Tp | The floating-point type of the modulus __k. |
| _Tpn | The floating-point type of the argument __nu. |
| __k | The modulus, abs(__k) <= 1 |
| __nu | The argument |
| std::domain_error | if abs(__k) > 1. |
|
inline |
|
inline |
|
inline |
Return the regular modified Bessel function
for real order
and argument
.
The regular modified cylindrical Bessel function is:
| _Tpnu | The floating-point type of the order __nu. |
| _Tp | The floating-point type of the argument __x. |
| __nu | The order |
| __x | The argument, __x >= 0 |
| std::domain_error | if __x < 0 . |
|
inline |
|
inline |
|
inline |
Return the Bessel function
of real order
and argument
.
The cylindrical Bessel function is:
| _Tpnu | The floating-point type of the order __nu. |
| _Tp | The floating-point type of the argument __x. |
| __nu | The order |
| __x | The argument, __x >= 0 |
| std::domain_error | if __x < 0 . |
|
inline |
|
inline |
|
inline |
Return the irregular modified Bessel function
of real order
and argument
.
The irregular modified Bessel function is defined by:
where for integral
a limit is taken:
. For negative argument we have simply:
| _Tpnu | The floating-point type of the order __nu. |
| _Tp | The floating-point type of the argument __x. |
| __nu | The order |
| __x | The argument, __x >= 0 |
| std::domain_error | if __x < 0 . |
|
inline |
|
inline |
|
inline |
Return the Neumann function
of real order
and argument
.
The Neumann function is defined by:
where
and for integral order
a limit is taken:
.
| _Tpnu | The floating-point type of the order __nu. |
| _Tp | The floating-point type of the argument __x. |
| __nu | The order |
| __x | The argument, __x >= 0 |
| std::domain_error | if __x < 0 . |
|
inline |
|
inline |
|
inline |
Return the incomplete elliptic integral of the first kind
for real modulus
and angle
.
The incomplete elliptic integral of the first kind is defined as
For
this becomes the complete elliptic integral of the first kind,
.
| _Tp | The floating-point type of the modulus __k. |
| _Tpp | The floating-point type of the angle __phi. |
| __k | The modulus, abs(__k) <= 1 |
| __phi | The integral limit argument in radians |
| std::domain_error | if abs(__k) > 1 . |
|
inline |
|
inline |
|
inline |
Return the incomplete elliptic integral of the second kind
.
The incomplete elliptic integral of the second kind is defined as
For
this becomes the complete elliptic integral of the second kind,
.
| _Tp | The floating-point type of the modulus __k. |
| _Tpp | The floating-point type of the angle __phi. |
| __k | The modulus, abs(__k) <= 1 |
| __phi | The integral limit argument in radians |
| std::domain_error | if abs(__k) > 1 . |
|
inline |
|
inline |
|
inline |
Return the incomplete elliptic integral of the third kind
.
The incomplete elliptic integral of the third kind is defined by:
For
this becomes the complete elliptic integral of the third kind,
.
| _Tp | The floating-point type of the modulus __k. |
| _Tpn | The floating-point type of the argument __nu. |
| _Tpp | The floating-point type of the angle __phi. |
| __k | The modulus, abs(__k) <= 1 |
| __nu | The second argument |
| __phi | The integral limit argument in radians |
| std::domain_error | if abs(__k) > 1 . |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
Return the Hermite polynomial
of order n and real argument x.
The Hermite polynomial is defined by:
The Hermite polynomial obeys a reflection formula:
| _Tp | The floating-point type of the argument __x. |
| __n | The order |
| __x | The argument |
|
inline |
|
inline |
|
inline |
Returns the Laguerre polynomial
of nonnegative degree n and real argument
.
The Laguerre polynomial is defined by:
| _Tp | The floating-point type of the argument __x. |
| __n | The nonnegative order |
| __x | The argument __x >= 0 |
| std::domain_error | if __x < 0 . |
|
inline |
|
inline |
|
inline |
Return the Legendre polynomial
of nonnegative degree
and real argument
.
The Legendre function of order
and argument
,
, is defined by:
| _Tp | The floating-point type of the argument __x. |
| __l | The degree |
| __x | The argument abs(__x) <= 1 |
| std::domain_error | if abs(__x) > 1 |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
Return the spherical Bessel function
of nonnegative order n and real argument
.
The spherical Bessel function is defined by:
| _Tp | The floating-point type of the argument __x. |
| __n | The integral order n >= 0 |
| __x | The real argument x >= 0 |
| std::domain_error | if __x < 0 . |
|
inline |
|
inline |
|
inline |
Return the spherical Legendre function of nonnegative integral degree l and order m and real angle
in radians.
The spherical Legendre function is defined by
| _Tp | The floating-point type of the angle __theta. |
| __l | The order __l >= 0 |
| __m | The degree __m >= 0 and __m <= __l |
| __theta | The radian polar angle argument |
|
inline |
|
inline |
|
inline |
Return the spherical Neumann function of integral order
and real argument
.
The spherical Neumann function is defined by
| _Tp | The floating-point type of the argument __x. |
| __n | The integral order n >= 0 |
| __x | The real argument __x >= 0 |
| std::domain_error | if __x < 0 . |
|
inline |