void std::tr1::__detail::__airy | ( | const _Tp | __x, | |
_Tp & | __Ai, | |||
_Tp & | __Bi, | |||
_Tp & | __Aip, | |||
_Tp & | __Bip | |||
) | [inline] |
Compute the Airy functions and
and their first derivatives
and
respectively.
__n | The order of the Airy functions. | |
__x | The argument of the Airy functions. | |
__i_n | The output Airy function. | |
__k_n | The output Airy function. | |
__ip_n | The output derivative of the Airy function. | |
__kp_n | The output derivative of the Airy function. |
Definition at line 370 of file modified_bessel_func.tcc.
References __bessel_ik(), __bessel_jn(), std::tr1::__detail::__numeric_constants< _Tp >::__pi(), std::tr1::__detail::__numeric_constants< _Tp >::__sqrt3(), std::abs(), and std::sqrt().
_Tp std::tr1::__detail::__assoc_laguerre | ( | const unsigned int | __n, | |
const unsigned int | __m, | |||
const _Tp | __x | |||
) | [inline] |
This routine returns the associated Laguerre polynomial of order n, degree m: .
The associated Laguerre polynomial is defined for integral by:
where the Laguerre polynomial is defined by:
__n | The order of the Laguerre polynomial. | |
__m | The degree of the Laguerre polynomial. | |
__x | The argument of the Laguerre polynomial. |
Definition at line 297 of file poly_laguerre.tcc.
_Tp std::tr1::__detail::__assoc_legendre_p | ( | const unsigned int | __l, | |
const unsigned int | __m, | |||
const _Tp | __x | |||
) | [inline] |
Return the associated Legendre function by recursion on .
The associated Legendre function is derived from the Legendre function by the Rodrigues formula:
l | The order of the associated Legendre function. ![]() | |
m | The order of the associated Legendre function. ![]() | |
x | The argument of the associated Legendre function. ![]() |
Definition at line 133 of file legendre_function.tcc.
References __poly_legendre_p(), and std::sqrt().
_Tp std::tr1::__detail::__bernoulli | ( | const int | __n | ) | [inline] |
_Tp std::tr1::__detail::__bernoulli_series | ( | unsigned int | __n | ) | [inline] |
This returns Bernoulli numbers from a table or by summation for larger values.
Recursion is unstable.
__n | the order n of the Bernoulli number. |
Definition at line 70 of file gamma.tcc.
References std::pow().
void std::tr1::__detail::__bessel_ik | ( | const _Tp | __nu, | |
const _Tp | __x, | |||
_Tp & | __Inu, | |||
_Tp & | __Knu, | |||
_Tp & | __Ipnu, | |||
_Tp & | __Kpnu | |||
) | [inline] |
Compute the modified Bessel functions and
and their first derivatives
and
respectively. These four functions are computed together for numerical stability.
__nu | The order of the Bessel functions. | |
__x | The argument of the Bessel functions. | |
__Inu | The output regular modified Bessel function. | |
__Knu | The output irregular modified Bessel function. | |
__Ipnu | The output derivative of the regular modified Bessel function. | |
__Kpnu | The output derivative of the irregular modified Bessel function. |
Definition at line 81 of file modified_bessel_func.tcc.
References __gamma_temme(), std::tr1::__detail::__numeric_constants< _Tp >::__pi(), std::abs(), std::cosh(), std::exp(), std::log(), std::sin(), std::sinh(), and std::sqrt().
Referenced by __airy(), __cyl_bessel_i(), __cyl_bessel_k(), and __sph_bessel_ik().
void std::tr1::__detail::__bessel_jn | ( | const _Tp | __nu, | |
const _Tp | __x, | |||
_Tp & | __Jnu, | |||
_Tp & | __Nnu, | |||
_Tp & | __Jpnu, | |||
_Tp & | __Npnu | |||
) | [inline] |
Compute the Bessel and Neumann
functions and their first derivatives
and
respectively. These four functions are computed together for numerical stability.
__nu | The order of the Bessel functions. | |
__x | The argument of the Bessel functions. | |
__Jnu | The output Bessel function of the first kind. | |
__Nnu | The output Neumann function (Bessel function of the second kind). | |
__Jpnu | The output derivative of the Bessel function of the first kind. | |
__Npnu | The output derivative of the Neumann function. |
Definition at line 128 of file bessel_function.tcc.
References __gamma_temme(), std::tr1::__detail::__numeric_constants< _Tp >::__pi(), std::abs(), std::cosh(), std::exp(), std::log(), std::max(), std::sin(), std::sinh(), and std::sqrt().
Referenced by __airy(), __cyl_bessel_j(), __cyl_neumann_n(), and __sph_bessel_jn().
_Tp std::tr1::__detail::__beta | ( | _Tp | __x, | |
_Tp | __y | |||
) | [inline] |
Return the beta function .
The beta function is defined by
__x | The first argument of the beta function. | |
__y | The second argument of the beta function. |
Definition at line 185 of file beta_function.tcc.
References __beta_lgamma().
Referenced by __ellint_rj().
_Tp std::tr1::__detail::__beta_gamma | ( | _Tp | __x, | |
_Tp | __y | |||
) | [inline] |
Return the beta function: .
The beta function is defined by
__x | The first argument of the beta function. | |
__y | The second argument of the beta function. |
Definition at line 75 of file beta_function.tcc.
References __gamma().
_Tp std::tr1::__detail::__beta_lgamma | ( | _Tp | __x, | |
_Tp | __y | |||
) | [inline] |
Return the beta function using the log gamma functions.
The beta function is defined by
__x | The first argument of the beta function. | |
__y | The second argument of the beta function. |
Definition at line 123 of file beta_function.tcc.
References __log_gamma(), and std::exp().
Referenced by __beta().
_Tp std::tr1::__detail::__beta_product | ( | _Tp | __x, | |
_Tp | __y | |||
) | [inline] |
Return the beta function using the product form.
The beta function is defined by
__x | The first argument of the beta function. | |
__y | The second argument of the beta function. |
Definition at line 154 of file beta_function.tcc.
_Tp std::tr1::__detail::__bincoef | ( | const unsigned int | __n, | |
const unsigned int | __k | |||
) | [inline] |
Return the binomial coefficient. The binomial coefficient is given by:
.
__n | The first argument of the binomial coefficient. | |
__k | The second argument of the binomial coefficient. |
Definition at line 310 of file gamma.tcc.
References std::exp(), and std::log().
_Tp std::tr1::__detail::__comp_ellint_1 | ( | const _Tp | __k | ) | [inline] |
Return the complete elliptic integral of the first kind using the Carlson formulation.
The complete elliptic integral of the first kind is defined as
where is the incomplete elliptic integral of the first kind.
__k | The argument of the complete elliptic function. |
Definition at line 191 of file ell_integral.tcc.
References __ellint_rf(), and std::abs().
Referenced by __ellint_1().
_Tp std::tr1::__detail::__comp_ellint_1_series | ( | const _Tp | __k | ) | [inline] |
Return the complete elliptic integral of the first kind by series expansion.
The complete elliptic integral of the first kind is defined as
This routine is not bad as long as |k| is somewhat smaller than 1 but is not is good as the Carlson elliptic integral formulation.
__k | The argument of the complete elliptic function. |
Definition at line 153 of file ell_integral.tcc.
References std::tr1::__detail::__numeric_constants< _Tp >::__pi_2().
_Tp std::tr1::__detail::__comp_ellint_2 | ( | const _Tp | __k | ) | [inline] |
Return the complete elliptic integral of the second kind using the Carlson formulation.
The complete elliptic integral of the second kind is defined as
__k | The argument of the complete elliptic function. |
Definition at line 402 of file ell_integral.tcc.
References __ellint_rd(), __ellint_rf(), and std::abs().
Referenced by __ellint_2().
_Tp std::tr1::__detail::__comp_ellint_2_series | ( | const _Tp | __k | ) | [inline] |
Return the complete elliptic integral of the second kind by series expansion.
The complete elliptic integral of the second kind is defined as
This routine is not bad as long as |k| is somewhat smaller than 1 but is not is good as the Carlson elliptic integral formulation.
__k | The argument of the complete elliptic function. |
Definition at line 266 of file ell_integral.tcc.
References std::tr1::__detail::__numeric_constants< _Tp >::__pi_2().
_Tp std::tr1::__detail::__comp_ellint_3 | ( | const _Tp | __k, | |
const _Tp | __nu | |||
) | [inline] |
Return the complete elliptic integral of the third kind using the Carlson formulation.
The complete elliptic integral of the third kind is defined as
__k | The argument of the elliptic function. | |
__nu | The second argument of the elliptic function. |
Definition at line 670 of file ell_integral.tcc.
References __ellint_rf(), __ellint_rj(), and std::abs().
Referenced by __ellint_3().
_Tp std::tr1::__detail::__conf_hyperg | ( | const _Tp | __a, | |
const _Tp | __c, | |||
const _Tp | __x | |||
) | [inline] |
Return the confluent hypogeometric function .
__a | The "numerator" parameter. | |
__c | The "denominator" parameter. | |
__x | The argument of the confluent hypergeometric function. |
Definition at line 222 of file hypergeometric.tcc.
References __conf_hyperg_luke(), __conf_hyperg_series(), and std::exp().
_Tp std::tr1::__detail::__conf_hyperg_luke | ( | const _Tp | __a, | |
const _Tp | __c, | |||
const _Tp | __xin | |||
) | [inline] |
Return the hypogeometric function by an iterative procedure described in Luke, Algorithms for the Computation of Mathematical Functions.
Like the case of the 2F1 rational approximations, these are probably guaranteed to converge for x < 0, barring gross numerical instability in the pre-asymptotic regime.
Definition at line 115 of file hypergeometric.tcc.
References std::abs(), and std::pow().
Referenced by __conf_hyperg().
_Tp std::tr1::__detail::__conf_hyperg_series | ( | const _Tp | __a, | |
const _Tp | __c, | |||
const _Tp | __x | |||
) | [inline] |
This routine returns the confluent hypergeometric function by series expansion.
If a and b are integers and a < 0 and either b > 0 or b < a then the series is a polynomial with a finite number of terms. If b is an integer and b <= 0 the confluent hypergeometric function is undefined.
__a | The "numerator" parameter. | |
__c | The "denominator" parameter. | |
__x | The argument of the confluent hypergeometric function. |
Definition at line 78 of file hypergeometric.tcc.
References std::abs().
Referenced by __conf_hyperg().
_Tp std::tr1::__detail::__cyl_bessel_i | ( | const _Tp | __nu, | |
const _Tp | __x | |||
) | [inline] |
Return the regular modified Bessel function of order :
.
The regular modified cylindrical Bessel function is:
__nu | The order of the regular modified Bessel function. | |
__x | The argument of the regular modified Bessel function. |
Definition at line 265 of file modified_bessel_func.tcc.
References __bessel_ik(), and __cyl_bessel_ij_series().
_Tp std::tr1::__detail::__cyl_bessel_ij_series | ( | const _Tp | __nu, | |
const _Tp | __x, | |||
const _Tp | __sgn, | |||
const unsigned int | __max_iter | |||
) | [inline] |
This routine returns the cylindrical Bessel functions of order :
or
by series expansion.
The modified cylindrical Bessel function is:
where or
for
or
respectively.
See Abramowitz & Stegun, 9.1.10 Abramowitz & Stegun, 9.6.7 (1) Handbook of Mathematical Functions, ed. Milton Abramowitz and Irene A. Stegun, Dover Publications, Equation 9.1.10 p. 360 and Equation 9.6.10 p. 375
__nu | The order of the Bessel function. | |
__x | The argument of the Bessel function. | |
__sgn | The sign of the alternate terms -1 for the Bessel function of the first kind. +1 for the modified Bessel function of the first kind. |
Definition at line 410 of file bessel_function.tcc.
References __log_gamma(), std::abs(), std::exp(), and std::log().
Referenced by __cyl_bessel_i(), and __cyl_bessel_j().
_Tp std::tr1::__detail::__cyl_bessel_j | ( | const _Tp | __nu, | |
const _Tp | __x | |||
) | [inline] |
Return the Bessel function of order :
.
The cylindrical Bessel function is:
__nu | The order of the Bessel function. | |
__x | The argument of the Bessel function. |
Definition at line 454 of file bessel_function.tcc.
References __bessel_jn(), __cyl_bessel_ij_series(), and __cyl_bessel_jn_asymp().
void std::tr1::__detail::__cyl_bessel_jn_asymp | ( | const _Tp | __nu, | |
const _Tp | __x, | |||
_Tp & | __Jnu, | |||
_Tp & | __Nnu | |||
) | [inline] |
This routine computes the asymptotic cylindrical Bessel and Neumann functions of order nu: ,
.
References: (1) Handbook of Mathematical Functions, ed. Milton Abramowitz and Irene A. Stegun, Dover Publications, Section 9 p. 364, Equations 9.2.5-9.2.10
__nu | The order of the Bessel functions. | |
__x | The argument of the Bessel functions. | |
__Jnu | The output Bessel function of the first kind. | |
__Nnu | The output Neumann function (Bessel function of the second kind). |
Definition at line 353 of file bessel_function.tcc.
References std::cos(), std::sin(), and std::sqrt().
Referenced by __cyl_bessel_j(), and __cyl_neumann_n().
_Tp std::tr1::__detail::__cyl_bessel_k | ( | const _Tp | __nu, | |
const _Tp | __x | |||
) | [inline] |
Return the irregular modified Bessel function of order
.
The irregular modified Bessel function is defined by:
where for integral a limit is taken:
.
__nu | The order of the irregular modified Bessel function. | |
__x | The argument of the irregular modified Bessel function. |
Definition at line 301 of file modified_bessel_func.tcc.
References __bessel_ik().
_Tp std::tr1::__detail::__cyl_neumann_n | ( | const _Tp | __nu, | |
const _Tp | __x | |||
) | [inline] |
Return the Neumann function of order :
.
The Neumann function is defined by:
where for integral a limit is taken:
.
__nu | The order of the Neumann function. | |
__x | The argument of the Neumann function. |
Definition at line 496 of file bessel_function.tcc.
References __bessel_jn(), and __cyl_bessel_jn_asymp().
_Tp std::tr1::__detail::__ellint_1 | ( | const _Tp | __k, | |
const _Tp | __phi | |||
) | [inline] |
Return the incomplete elliptic integral of the first kind using the Carlson formulation.
The incomplete elliptic integral of the first kind is defined as
__k | The argument of the elliptic function. | |
__phi | The integral limit argument of the elliptic function. |
Definition at line 219 of file ell_integral.tcc.
References __comp_ellint_1(), __ellint_rf(), std::tr1::__detail::__numeric_constants< _Tp >::__pi(), std::abs(), std::cos(), and std::sin().
_Tp std::tr1::__detail::__ellint_2 | ( | const _Tp | __k, | |
const _Tp | __phi | |||
) | [inline] |
Return the incomplete elliptic integral of the second kind using the Carlson formulation.
The incomplete elliptic integral of the second kind is defined as
__k | The argument of the elliptic function. | |
__phi | The integral limit argument of the elliptic function. |
Definition at line 436 of file ell_integral.tcc.
References __comp_ellint_2(), __ellint_rd(), __ellint_rf(), std::tr1::__detail::__numeric_constants< _Tp >::__pi(), std::abs(), std::cos(), and std::sin().
_Tp std::tr1::__detail::__ellint_3 | ( | const _Tp | __k, | |
const _Tp | __nu, | |||
const _Tp | __phi | |||
) | [inline] |
Return the incomplete elliptic integral of the third kind using the Carlson formulation.
The incomplete elliptic integral of the third kind is defined as
__k | The argument of the elliptic function. | |
__nu | The second argument of the elliptic function. | |
__phi | The integral limit argument of the elliptic function. |
Definition at line 710 of file ell_integral.tcc.
References __comp_ellint_3(), __ellint_rf(), __ellint_rj(), std::tr1::__detail::__numeric_constants< _Tp >::__pi(), std::abs(), std::cos(), and std::sin().
_Tp std::tr1::__detail::__ellint_rc | ( | const _Tp | __x, | |
const _Tp | __y | |||
) | [inline] |
Return the Carlson elliptic function where
is the Carlson elliptic function of the first kind.
The Carlson elliptic function is defined by:
Based on Carlson's algorithms:
__x | The first argument. | |
__y | The second argument. |
Definition at line 495 of file ell_integral.tcc.
References std::abs(), std::max(), std::min(), std::pow(), and std::sqrt().
Referenced by __ellint_rj().
_Tp std::tr1::__detail::__ellint_rd | ( | const _Tp | __x, | |
const _Tp | __y, | |||
const _Tp | __z | |||
) | [inline] |
Return the Carlson elliptic function of the second kind where
is the Carlson elliptic function of the third kind.
The Carlson elliptic function of the second kind is defined by:
Based on Carlson's algorithms:
__x | The first of two symmetric arguments. | |
__y | The second of two symmetric arguments. | |
__z | The third argument. |
Definition at line 314 of file ell_integral.tcc.
References std::abs(), std::max(), std::min(), std::pow(), and std::sqrt().
Referenced by __comp_ellint_2(), and __ellint_2().
_Tp std::tr1::__detail::__ellint_rf | ( | const _Tp | __x, | |
const _Tp | __y, | |||
const _Tp | __z | |||
) | [inline] |
Return the Carlson elliptic function of the first kind.
The Carlson elliptic function of the first kind is defined by:
__x | The first of three symmetric arguments. | |
__y | The second of three symmetric arguments. | |
__z | The third of three symmetric arguments. |
Definition at line 74 of file ell_integral.tcc.
References std::abs(), std::max(), std::min(), std::pow(), and std::sqrt().
Referenced by __comp_ellint_1(), __comp_ellint_2(), __comp_ellint_3(), __ellint_1(), __ellint_2(), and __ellint_3().
_Tp std::tr1::__detail::__ellint_rj | ( | const _Tp | __x, | |
const _Tp | __y, | |||
const _Tp | __z, | |||
const _Tp | __p | |||
) | [inline] |
Return the Carlson elliptic function of the third kind.
The Carlson elliptic function of the third kind is defined by:
Based on Carlson's algorithms:
__x | The first of three symmetric arguments. | |
__y | The second of three symmetric arguments. | |
__z | The third of three symmetric arguments. | |
__p | The fourth argument. |
Definition at line 566 of file ell_integral.tcc.
References __beta(), __ellint_rc(), std::abs(), std::max(), std::min(), std::pow(), and std::sqrt().
Referenced by __comp_ellint_3(), and __ellint_3().
_Tp std::tr1::__detail::__expint | ( | const _Tp | __x | ) | [inline] |
Return the exponential integral .
The exponential integral is given by
__x | The argument of the exponential integral function. |
Definition at line 512 of file exp_integral.tcc.
References __expint_Ei().
_Tp std::tr1::__detail::__expint | ( | const unsigned int | __n, | |
const _Tp | __x | |||
) | [inline] |
Return the exponential integral .
The exponential integral is given by
This is something of an extension.
__n | The order of the exponential integral function. | |
__x | The argument of the exponential integral function. |
Definition at line 472 of file exp_integral.tcc.
References __expint_E1(), __expint_En_recursion(), and std::exp().
_Tp std::tr1::__detail::__expint_asymp | ( | const unsigned int | __n, | |
const _Tp | __x | |||
) | [inline] |
Return the exponential integral for large argument.
The exponential integral is given by
This is something of an extension.
__n | The order of the exponential integral function. | |
__x | The argument of the exponential integral function. |
Definition at line 404 of file exp_integral.tcc.
References std::abs(), and std::exp().
_Tp std::tr1::__detail::__expint_E1 | ( | const _Tp | __x | ) | [inline] |
Return the exponential integral .
The exponential integral is given by
__x | The argument of the exponential integral function. |
Definition at line 374 of file exp_integral.tcc.
References __expint_E1_asymp(), __expint_E1_series(), __expint_Ei(), and __expint_En_cont_frac().
Referenced by __expint(), __expint_Ei(), and __expint_En_recursion().
_Tp std::tr1::__detail::__expint_E1_asymp | ( | const _Tp | __x | ) | [inline] |
Return the exponential integral by asymptotic expansion.
The exponential integral is given by
__x | The argument of the exponential integral function. |
Definition at line 114 of file exp_integral.tcc.
References std::abs(), and std::exp().
Referenced by __expint_E1().
_Tp std::tr1::__detail::__expint_E1_series | ( | const _Tp | __x | ) | [inline] |
Return the exponential integral by series summation. This should be good for
.
The exponential integral is given by
__x | The argument of the exponential integral function. |
Definition at line 77 of file exp_integral.tcc.
References std::tr1::__detail::__numeric_constants< _Tp >::__gamma_e(), std::abs(), and std::log().
Referenced by __expint_E1().
_Tp std::tr1::__detail::__expint_Ei | ( | const _Tp | __x | ) | [inline] |
Return the exponential integral .
The exponential integral is given by
__x | The argument of the exponential integral function. |
Definition at line 350 of file exp_integral.tcc.
References __expint_E1(), __expint_Ei_asymp(), __expint_Ei_series(), and std::log().
Referenced by __expint(), and __expint_E1().
_Tp std::tr1::__detail::__expint_Ei_asymp | ( | const _Tp | __x | ) | [inline] |
Return the exponential integral by asymptotic expansion.
The exponential integral is given by
__x | The argument of the exponential integral function. |
Definition at line 317 of file exp_integral.tcc.
References std::exp().
Referenced by __expint_Ei().
_Tp std::tr1::__detail::__expint_Ei_series | ( | const _Tp | __x | ) | [inline] |
Return the exponential integral by series summation.
The exponential integral is given by
__x | The argument of the exponential integral function. |
Definition at line 286 of file exp_integral.tcc.
References std::tr1::__detail::__numeric_constants< _Tp >::__gamma_e(), and std::log().
Referenced by __expint_Ei().
_Tp std::tr1::__detail::__expint_En_cont_frac | ( | const unsigned int | __n, | |
const _Tp | __x | |||
) | [inline] |
Return the exponential integral by continued fractions.
The exponential integral is given by
__n | The order of the exponential integral function. | |
__x | The argument of the exponential integral function. |
Definition at line 197 of file exp_integral.tcc.
References std::abs(), std::exp(), and std::min().
Referenced by __expint_E1().
_Tp std::tr1::__detail::__expint_En_recursion | ( | const unsigned int | __n, | |
const _Tp | __x | |||
) | [inline] |
Return the exponential integral by recursion. Use upward recursion for
and downward recursion (Miller's algorithm) otherwise.
The exponential integral is given by
__n | The order of the exponential integral function. | |
__x | The argument of the exponential integral function. |
Definition at line 242 of file exp_integral.tcc.
References __expint_E1(), and std::exp().
Referenced by __expint().
_Tp std::tr1::__detail::__expint_En_series | ( | const unsigned int | __n, | |
const _Tp | __x | |||
) | [inline] |
Return the exponential integral by series summation.
The exponential integral is given by
__n | The order of the exponential integral function. | |
__x | The argument of the exponential integral function. |
Definition at line 151 of file exp_integral.tcc.
References std::tr1::__detail::__numeric_constants< _Tp >::__gamma_e(), __psi(), std::abs(), and std::log().
_Tp std::tr1::__detail::__expint_large_n | ( | const unsigned int | __n, | |
const _Tp | __x | |||
) | [inline] |
Return the exponential integral for large order.
The exponential integral is given by
This is something of an extension.
__n | The order of the exponential integral function. | |
__x | The argument of the exponential integral function. |
Definition at line 438 of file exp_integral.tcc.
References std::abs(), and std::exp().
_Tp std::tr1::__detail::__gamma | ( | const _Tp | __x | ) | [inline] |
Return .
__x | The argument of the gamma function. |
Definition at line 333 of file gamma.tcc.
References __log_gamma(), and std::exp().
Referenced by __beta_gamma(), and __gamma_temme().
void std::tr1::__detail::__gamma_temme | ( | const _Tp | __mu, | |
_Tp & | __gam1, | |||
_Tp & | __gam2, | |||
_Tp & | __gampl, | |||
_Tp & | __gammi | |||
) | [inline] |
Compute the gamma functions required by the Temme series expansions of and
.
and
where is
and
. is the nearest integer to
. The values of
and
are returned as well.
The accuracy requirements on this are exquisite.
__mu | The input parameter of the gamma functions. | |
__gam1 | The output function ![]() | |
__gam2 | The output function ![]() | |
__gampl | The output function ![]() | |
__gammi | The output function ![]() |
Definition at line 90 of file bessel_function.tcc.
References __gamma(), and std::abs().
Referenced by __bessel_ik(), and __bessel_jn().
_Tp std::tr1::__detail::__hurwitz_zeta | ( | const _Tp | __a, | |
const _Tp | __s | |||
) | [inline] |
Return the Hurwitz zeta function for all s != 1 and x > -1.
The Hurwitz zeta function is defined by:
The Riemann zeta function is a special case:
Definition at line 426 of file riemann_zeta.tcc.
References __hurwitz_zeta_glob().
Referenced by __psi().
_Tp std::tr1::__detail::__hurwitz_zeta_glob | ( | const _Tp | __a, | |
const _Tp | __s | |||
) | [inline] |
Return the Hurwitz zeta function for all s != 1 and x > -1.
The Hurwitz zeta function is defined by:
The Riemann zeta function is a special case:
This functions uses the double sum that converges for s != 1 and x > -1:
Definition at line 361 of file riemann_zeta.tcc.
References __log_gamma(), std::abs(), std::exp(), std::log(), and std::pow().
Referenced by __hurwitz_zeta().
_Tp std::tr1::__detail::__hyperg | ( | const _Tp | __a, | |
const _Tp | __b, | |||
const _Tp | __c, | |||
const _Tp | __x | |||
) | [inline] |
Return the hypogeometric function .
The hypogeometric function is defined by
__a | The first "numerator" parameter. | |
__a | The second "numerator" parameter. | |
__c | The "denominator" parameter. | |
__x | The argument of the confluent hypergeometric function. |
Definition at line 724 of file hypergeometric.tcc.
References __hyperg_luke(), __hyperg_reflect(), __hyperg_series(), std::abs(), and std::pow().
_Tp std::tr1::__detail::__hyperg_luke | ( | const _Tp | __a, | |
const _Tp | __b, | |||
const _Tp | __c, | |||
const _Tp | __xin | |||
) | [inline] |
Return the hypogeometric function by an iterative procedure described in Luke, Algorithms for the Computation of Mathematical Functions.
Definition at line 300 of file hypergeometric.tcc.
References std::abs(), and std::pow().
Referenced by __hyperg().
_Tp std::tr1::__detail::__hyperg_reflect | ( | const _Tp | __a, | |
const _Tp | __b, | |||
const _Tp | __c, | |||
const _Tp | __x | |||
) | [inline] |
Return the hypogeometric function by the reflection formulae in Abramowitz & Stegun formula 15.3.6 for d = c - a - b not integral and formula 15.3.11 for d = c - a - b integral. This assumes a, b, c != negative integer.
The hypogeometric function is defined by
The reflection formula for nonintegral is:
The reflection formula for integral is:
Definition at line 434 of file hypergeometric.tcc.
References std::tr1::__detail::__numeric_constants< _Tp >::__gamma_e(), __hyperg_series(), __log_gamma(), __log_gamma_sign(), __psi(), std::abs(), std::exp(), and std::log().
Referenced by __hyperg().
_Tp std::tr1::__detail::__hyperg_series | ( | const _Tp | __a, | |
const _Tp | __b, | |||
const _Tp | __c, | |||
const _Tp | __x | |||
) | [inline] |
Return the hypogeometric function by series expansion.
The hypogeometric function is defined by
This works and it's pretty fast.
__a | The first "numerator" parameter. | |
__a | The second "numerator" parameter. | |
__c | The "denominator" parameter. | |
__x | The argument of the confluent hypergeometric function. |
Definition at line 266 of file hypergeometric.tcc.
References std::abs().
Referenced by __hyperg(), and __hyperg_reflect().
_Tp std::tr1::__detail::__laguerre | ( | const unsigned int | __n, | |
const _Tp | __x | |||
) | [inline] |
This routine returns the Laguerre polynomial of order n: .
The Laguerre polynomial is defined by:
__n | The order of the Laguerre polynomial. | |
__x | The argument of the Laguerre polynomial. |
Definition at line 320 of file poly_laguerre.tcc.
_Tp std::tr1::__detail::__log_bincoef | ( | const unsigned int | __n, | |
const unsigned int | __k | |||
) | [inline] |
Return the logarithm of the binomial coefficient. The binomial coefficient is given by:
.
__n | The first argument of the binomial coefficient. | |
__k | The second argument of the binomial coefficient. |
Definition at line 279 of file gamma.tcc.
References __log_gamma(), and std::log().
_Tp std::tr1::__detail::__log_gamma | ( | const _Tp | __x | ) | [inline] |
Return . This will return values even for
. To recover the sign of
for any argument use __log_gamma_sign.
__x | The argument of the log of the gamma function. |
Definition at line 221 of file gamma.tcc.
References std::tr1::__detail::__numeric_constants< _Tp >::__lnpi(), __log_gamma_lanczos(), std::abs(), std::log(), and std::sin().
Referenced by __beta_lgamma(), __cyl_bessel_ij_series(), __gamma(), __hurwitz_zeta_glob(), __hyperg_reflect(), __log_bincoef(), __poly_laguerre_large_n(), __psi(), __riemann_zeta(), __riemann_zeta_glob(), and __sph_legendre().
_Tp std::tr1::__detail::__log_gamma_bernoulli | ( | const _Tp | __x | ) | [inline] |
Return by asymptotic expansion with Bernoulli number coefficients. This is like Sterling's approximation.
__x | The argument of the log of the gamma function. |
Definition at line 149 of file gamma.tcc.
References std::__lg(), and std::log().
_Tp std::tr1::__detail::__log_gamma_lanczos | ( | const _Tp | __x | ) | [inline] |
Return by the Lanczos method. This method dominates all others on the positive axis I think.
__x | The argument of the log of the gamma function. |
Definition at line 177 of file gamma.tcc.
References std::tr1::__detail::__numeric_constants< _Tp >::__euler(), and std::log().
Referenced by __log_gamma().
_Tp std::tr1::__detail::__log_gamma_sign | ( | const _Tp | __x | ) | [inline] |
Return the sign of . At nonpositive integers zero is returned.
__x | The argument of the gamma function. |
Definition at line 248 of file gamma.tcc.
References std::sin().
Referenced by __hyperg_reflect().
_Tp std::tr1::__detail::__poly_hermite | ( | const unsigned int | __n, | |
const _Tp | __x | |||
) | [inline] |
This routine returns the Hermite polynomial of order n: .
The Hermite polynomial is defined by:
__n | The order of the Hermite polynomial. | |
__x | The argument of the Hermite polynomial. |
Definition at line 112 of file poly_hermite.tcc.
References __poly_hermite_recursion().
_Tp std::tr1::__detail::__poly_hermite_recursion | ( | const unsigned int | __n, | |
const _Tp | __x | |||
) | [inline] |
This routine returns the Hermite polynomial of order n: by recursion on n.
The Hermite polynomial is defined by:
__n | The order of the Hermite polynomial. | |
__x | The argument of the Hermite polynomial. |
Definition at line 70 of file poly_hermite.tcc.
Referenced by __poly_hermite().
_Tp std::tr1::__detail::__poly_laguerre | ( | const unsigned int | __n, | |
const _Tpa | __alpha1, | |||
const _Tp | __x | |||
) | [inline] |
This routine returns the associated Laguerre polynomial of order n, degree :
.
The associated Laguerre function is defined by
where is the Pochhammer symbol and
is the confluent hypergeometric function.
The associated Laguerre polynomial is defined for integral by:
where the Laguerre polynomial is defined by:
__n | The order of the Laguerre function. | |
__alpha | The degree of the Laguerre function. | |
__x | The argument of the Laguerre function. |
Definition at line 244 of file poly_laguerre.tcc.
References __poly_laguerre_hyperg(), __poly_laguerre_large_n(), and __poly_laguerre_recursion().
_Tp std::tr1::__detail::__poly_laguerre_hyperg | ( | const unsigned int | __n, | |
const _Tpa | __alpha1, | |||
const _Tp | __x | |||
) | [inline] |
Evaluate the polynomial based on the confluent hypergeometric function in a safe way, with no restriction on the arguments.
The associated Laguerre function is defined by
where is the Pochhammer symbol and
is the confluent hypergeometric function.
This function assumes x != 0.
This is from the GNU Scientific Library.
Definition at line 127 of file poly_laguerre.tcc.
References std::abs().
Referenced by __poly_laguerre().
_Tp std::tr1::__detail::__poly_laguerre_large_n | ( | const unsigned | __n, | |
const _Tpa | __alpha1, | |||
const _Tp | __x | |||
) | [inline] |
This routine returns the associated Laguerre polynomial of order , degree
for large n. Abramowitz & Stegun, 13.5.21.
__n | The order of the Laguerre function. | |
__alpha | The degree of the Laguerre function. | |
__x | The argument of the Laguerre function. |
Definition at line 72 of file poly_laguerre.tcc.
References __log_gamma(), std::tr1::__detail::__numeric_constants< _Tp >::__pi_2(), std::exp(), std::log(), std::sin(), and std::sqrt().
Referenced by __poly_laguerre().
_Tp std::tr1::__detail::__poly_laguerre_recursion | ( | const unsigned int | __n, | |
const _Tpa | __alpha1, | |||
const _Tp | __x | |||
) | [inline] |
This routine returns the associated Laguerre polynomial of order , degree
:
by recursion.
The associated Laguerre function is defined by
where is the Pochhammer symbol and
is the confluent hypergeometric function.
The associated Laguerre polynomial is defined for integral by:
where the Laguerre polynomial is defined by:
__n | The order of the Laguerre function. | |
__alpha | The degree of the Laguerre function. | |
__x | The argument of the Laguerre function. |
Definition at line 184 of file poly_laguerre.tcc.
Referenced by __poly_laguerre().
_Tp std::tr1::__detail::__poly_legendre_p | ( | const unsigned int | __l, | |
const _Tp | __x | |||
) | [inline] |
Return the Legendre polynomial by recursion on order .
The Legendre function of and
,
, is defined by:
l | The order of the Legendre polynomial. ![]() | |
x | The argument of the Legendre polynomial. ![]() |
Definition at line 76 of file legendre_function.tcc.
Referenced by __assoc_legendre_p(), and __sph_legendre().
_Tp std::tr1::__detail::__psi | ( | const unsigned int | __n, | |
const _Tp | __x | |||
) | [inline] |
Return the polygamma function .
The polygamma function is related to the Hurwitz zeta function:
Definition at line 444 of file gamma.tcc.
References __hurwitz_zeta(), __log_gamma(), __psi(), and std::exp().
_Tp std::tr1::__detail::__psi | ( | const _Tp | __x | ) | [inline] |
Return the digamma function. The digamma or function is defined by
For negative argument the reflection formula is used:
.
Definition at line 415 of file gamma.tcc.
References std::tr1::__detail::__numeric_constants< _Tp >::__pi(), __psi_asymp(), __psi_series(), std::abs(), std::cos(), and std::sin().
Referenced by __expint_En_series(), __hyperg_reflect(), and __psi().
_Tp std::tr1::__detail::__psi_asymp | ( | const _Tp | __x | ) | [inline] |
Return the digamma function for large argument. The digamma or function is defined by
.
The asymptotic series is given by:
Definition at line 384 of file gamma.tcc.
References std::abs(), and std::log().
Referenced by __psi().
_Tp std::tr1::__detail::__psi_series | ( | const _Tp | __x | ) | [inline] |
Return the digamma function by series expansion. The digamma or function is defined by
.
The series is given by:
Definition at line 354 of file gamma.tcc.
References std::tr1::__detail::__numeric_constants< _Tp >::__gamma_e(), and std::abs().
Referenced by __psi().
_Tp std::tr1::__detail::__riemann_zeta | ( | const _Tp | __s | ) | [inline] |
Return the Riemann zeta function .
The Riemann zeta function is defined by:
For s < 1 use the reflection formula:
Definition at line 289 of file riemann_zeta.tcc.
References __log_gamma(), std::tr1::__detail::__numeric_constants< _Tp >::__pi(), __riemann_zeta_glob(), __riemann_zeta_product(), __riemann_zeta_sum(), std::exp(), std::pow(), and std::sin().
_Tp std::tr1::__detail::__riemann_zeta_alt | ( | const _Tp | __s | ) | [inline] |
Evaluate the Riemann zeta function by an alternate series for s > 0.
The Riemann zeta function is defined by:
For s < 1 use the reflection formula:
Definition at line 111 of file riemann_zeta.tcc.
References std::abs(), and std::pow().
_Tp std::tr1::__detail::__riemann_zeta_glob | ( | const _Tp | __s | ) | [inline] |
Evaluate the Riemann zeta function by series for all s != 1. Convergence is great until largish negative numbers. Then the convergence of the > 0 sum gets better.
The series is:
Havil 2003, p. 206.
The Riemann zeta function is defined by:
For s < 1 use the reflection formula:
Definition at line 153 of file riemann_zeta.tcc.
References __log_gamma(), std::tr1::__detail::__numeric_constants< _Tp >::__pi(), std::abs(), std::exp(), std::log(), std::pow(), and std::sin().
Referenced by __riemann_zeta().
_Tp std::tr1::__detail::__riemann_zeta_product | ( | const _Tp | __s | ) | [inline] |
Compute the Riemann zeta function using the product over prime factors.
where are the prime numbers.
The Riemann zeta function is defined by:
For s < 1 use the reflection formula:
Definition at line 248 of file riemann_zeta.tcc.
References std::pow().
Referenced by __riemann_zeta().
_Tp std::tr1::__detail::__riemann_zeta_sum | ( | const _Tp | __s | ) | [inline] |
Compute the Riemann zeta function by summation for s > 1.
The Riemann zeta function is defined by:
For s < 1 use the reflection formula:
Definition at line 74 of file riemann_zeta.tcc.
References std::pow().
Referenced by __riemann_zeta().
_Tp std::tr1::__detail::__sph_bessel | ( | const unsigned int | __n, | |
const _Tp | __x | |||
) | [inline] |
Return the spherical Bessel function of order n.
The spherical Bessel function is defined by:
__n | The order of the spherical Bessel function. | |
__x | The argument of the spherical Bessel function. |
Definition at line 568 of file bessel_function.tcc.
References __sph_bessel_jn().
void std::tr1::__detail::__sph_bessel_ik | ( | const unsigned int | __n, | |
const _Tp | __x, | |||
_Tp & | __i_n, | |||
_Tp & | __k_n, | |||
_Tp & | __ip_n, | |||
_Tp & | __kp_n | |||
) | [inline] |
Compute the spherical modified Bessel functions and
and their first derivatives
and
respectively.
__n | The order of the modified spherical Bessel function. | |
__x | The argument of the modified spherical Bessel function. | |
__i_n | The output regular modified spherical Bessel function. | |
__k_n | The output irregular modified spherical Bessel function. | |
__ip_n | The output derivative of the regular modified spherical Bessel function. | |
__kp_n | The output derivative of the irregular modified spherical Bessel function. |
Definition at line 335 of file modified_bessel_func.tcc.
References __bessel_ik(), std::tr1::__detail::__numeric_constants< _Tp >::__sqrtpio2(), and std::sqrt().
void std::tr1::__detail::__sph_bessel_jn | ( | const unsigned int | __n, | |
const _Tp | __x, | |||
_Tp & | __j_n, | |||
_Tp & | __n_n, | |||
_Tp & | __jp_n, | |||
_Tp & | __np_n | |||
) | [inline] |
Compute the spherical Bessel and Neumann
functions and their first derivatives
and
respectively.
__n | The order of the spherical Bessel function. | |
__x | The argument of the spherical Bessel function. | |
__j_n | The output spherical Bessel function. | |
__n_n | The output spherical Neumann function. | |
__jp_n | The output derivative of the spherical Bessel function. | |
__np_n | The output derivative of the spherical Neumann function. |
Definition at line 533 of file bessel_function.tcc.
References __bessel_jn(), std::tr1::__detail::__numeric_constants< _Tp >::__sqrtpio2(), and std::sqrt().
Referenced by __sph_bessel(), and __sph_neumann().
_Tp std::tr1::__detail::__sph_legendre | ( | const unsigned int | __l, | |
const unsigned int | __m, | |||
const _Tp | __theta | |||
) | [inline] |
Return the spherical associated Legendre function.
The spherical associated Legendre function of ,
, and
is defined as
where
is the spherical harmonic function and is the associated Legendre function.
This function differs from the associated Legendre function by argument () and by a normalization factor but this factor is rather large for large
and
and so this function is stable for larger differences of
and
.
l | The order of the spherical associated Legendre function. ![]() | |
m | The order of the spherical associated Legendre function. ![]() | |
theta | The radian angle argument of the spherical associated Legendre function. |
Definition at line 213 of file legendre_function.tcc.
References std::tr1::__detail::__numeric_constants< _Tp >::__lnpi(), __log_gamma(), __poly_legendre_p(), std::cos(), std::exp(), std::log(), and std::sqrt().
_Tp std::tr1::__detail::__sph_neumann | ( | const unsigned int | __n, | |
const _Tp | __x | |||
) | [inline] |
Return the spherical Neumann function .
The spherical Neumann function is defined by:
__n | The order of the spherical Neumann function. | |
__x | The argument of the spherical Neumann function. |
Definition at line 606 of file bessel_function.tcc.
References __sph_bessel_jn().