11debfc3dSmrg // Mathematical Special Functions for -*- C++ -*-
21debfc3dSmrg
3*8feb0f0bSmrg // Copyright (C) 2006-2020 Free Software Foundation, Inc.
41debfc3dSmrg //
51debfc3dSmrg // This file is part of the GNU ISO C++ Library. This library is free
61debfc3dSmrg // software; you can redistribute it and/or modify it under the
71debfc3dSmrg // terms of the GNU General Public License as published by the
81debfc3dSmrg // Free Software Foundation; either version 3, or (at your option)
91debfc3dSmrg // any later version.
101debfc3dSmrg
111debfc3dSmrg // This library is distributed in the hope that it will be useful,
121debfc3dSmrg // but WITHOUT ANY WARRANTY; without even the implied warranty of
131debfc3dSmrg // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
141debfc3dSmrg // GNU General Public License for more details.
151debfc3dSmrg
161debfc3dSmrg // Under Section 7 of GPL version 3, you are granted additional
171debfc3dSmrg // permissions described in the GCC Runtime Library Exception, version
181debfc3dSmrg // 3.1, as published by the Free Software Foundation.
191debfc3dSmrg
201debfc3dSmrg // You should have received a copy of the GNU General Public License and
211debfc3dSmrg // a copy of the GCC Runtime Library Exception along with this program;
221debfc3dSmrg // see the files COPYING3 and COPYING.RUNTIME respectively. If not, see
231debfc3dSmrg // <http://www.gnu.org/licenses/>.
241debfc3dSmrg
251debfc3dSmrg /** @file bits/specfun.h
261debfc3dSmrg * This is an internal header file, included by other library headers.
271debfc3dSmrg * Do not attempt to use it directly. @headername{cmath}
281debfc3dSmrg */
291debfc3dSmrg
301debfc3dSmrg #ifndef _GLIBCXX_BITS_SPECFUN_H
311debfc3dSmrg #define _GLIBCXX_BITS_SPECFUN_H 1
321debfc3dSmrg
331debfc3dSmrg #pragma GCC visibility push(default)
341debfc3dSmrg
351debfc3dSmrg #include <bits/c++config.h>
361debfc3dSmrg
371debfc3dSmrg #define __STDCPP_MATH_SPEC_FUNCS__ 201003L
381debfc3dSmrg
391debfc3dSmrg #define __cpp_lib_math_special_functions 201603L
401debfc3dSmrg
411debfc3dSmrg #if __cplusplus <= 201403L && __STDCPP_WANT_MATH_SPEC_FUNCS__ == 0
421debfc3dSmrg # error include <cmath> and define __STDCPP_WANT_MATH_SPEC_FUNCS__
431debfc3dSmrg #endif
441debfc3dSmrg
451debfc3dSmrg #include <bits/stl_algobase.h>
461debfc3dSmrg #include <limits>
471debfc3dSmrg #include <type_traits>
481debfc3dSmrg
491debfc3dSmrg #include <tr1/gamma.tcc>
501debfc3dSmrg #include <tr1/bessel_function.tcc>
511debfc3dSmrg #include <tr1/beta_function.tcc>
521debfc3dSmrg #include <tr1/ell_integral.tcc>
531debfc3dSmrg #include <tr1/exp_integral.tcc>
541debfc3dSmrg #include <tr1/hypergeometric.tcc>
551debfc3dSmrg #include <tr1/legendre_function.tcc>
561debfc3dSmrg #include <tr1/modified_bessel_func.tcc>
571debfc3dSmrg #include <tr1/poly_hermite.tcc>
581debfc3dSmrg #include <tr1/poly_laguerre.tcc>
591debfc3dSmrg #include <tr1/riemann_zeta.tcc>
601debfc3dSmrg
_GLIBCXX_VISIBILITY(default)611debfc3dSmrg namespace std _GLIBCXX_VISIBILITY(default)
621debfc3dSmrg {
631debfc3dSmrg _GLIBCXX_BEGIN_NAMESPACE_VERSION
641debfc3dSmrg
651debfc3dSmrg /**
661debfc3dSmrg * @defgroup mathsf Mathematical Special Functions
671debfc3dSmrg * @ingroup numerics
681debfc3dSmrg *
69*8feb0f0bSmrg * @section mathsf_desc Mathematical Special Functions
701debfc3dSmrg *
71*8feb0f0bSmrg * A collection of advanced mathematical special functions,
72*8feb0f0bSmrg * defined by ISO/IEC IS 29124 and then added to ISO C++ 2017.
73*8feb0f0bSmrg *
74*8feb0f0bSmrg *
75*8feb0f0bSmrg * @subsection mathsf_intro Introduction and History
761debfc3dSmrg * The first significant library upgrade on the road to C++2011,
771debfc3dSmrg * <a href="http://www.open-std.org/JTC1/SC22/WG21/docs/papers/2005/n1836.pdf">
781debfc3dSmrg * TR1</a>, included a set of 23 mathematical functions that significantly
791debfc3dSmrg * extended the standard transcendental functions inherited from C and declared
801debfc3dSmrg * in @<cmath@>.
811debfc3dSmrg *
821debfc3dSmrg * Although most components from TR1 were eventually adopted for C++11 these
831debfc3dSmrg * math functions were left behind out of concern for implementability.
841debfc3dSmrg * The math functions were published as a separate international standard
851debfc3dSmrg * <a href="http://www.open-std.org/JTC1/SC22/WG21/docs/papers/2010/n3060.pdf">
861debfc3dSmrg * IS 29124 - Extensions to the C++ Library to Support Mathematical Special
871debfc3dSmrg * Functions</a>.
881debfc3dSmrg *
891debfc3dSmrg * For C++17 these functions were incorporated into the main standard.
901debfc3dSmrg *
91*8feb0f0bSmrg * @subsection mathsf_contents Contents
921debfc3dSmrg * The following functions are implemented in namespace @c std:
931debfc3dSmrg * - @ref assoc_laguerre "assoc_laguerre - Associated Laguerre functions"
941debfc3dSmrg * - @ref assoc_legendre "assoc_legendre - Associated Legendre functions"
951debfc3dSmrg * - @ref beta "beta - Beta functions"
961debfc3dSmrg * - @ref comp_ellint_1 "comp_ellint_1 - Complete elliptic functions of the first kind"
971debfc3dSmrg * - @ref comp_ellint_2 "comp_ellint_2 - Complete elliptic functions of the second kind"
981debfc3dSmrg * - @ref comp_ellint_3 "comp_ellint_3 - Complete elliptic functions of the third kind"
991debfc3dSmrg * - @ref cyl_bessel_i "cyl_bessel_i - Regular modified cylindrical Bessel functions"
1001debfc3dSmrg * - @ref cyl_bessel_j "cyl_bessel_j - Cylindrical Bessel functions of the first kind"
1011debfc3dSmrg * - @ref cyl_bessel_k "cyl_bessel_k - Irregular modified cylindrical Bessel functions"
1021debfc3dSmrg * - @ref cyl_neumann "cyl_neumann - Cylindrical Neumann functions or Cylindrical Bessel functions of the second kind"
1031debfc3dSmrg * - @ref ellint_1 "ellint_1 - Incomplete elliptic functions of the first kind"
1041debfc3dSmrg * - @ref ellint_2 "ellint_2 - Incomplete elliptic functions of the second kind"
1051debfc3dSmrg * - @ref ellint_3 "ellint_3 - Incomplete elliptic functions of the third kind"
1061debfc3dSmrg * - @ref expint "expint - The exponential integral"
1071debfc3dSmrg * - @ref hermite "hermite - Hermite polynomials"
1081debfc3dSmrg * - @ref laguerre "laguerre - Laguerre functions"
1091debfc3dSmrg * - @ref legendre "legendre - Legendre polynomials"
1101debfc3dSmrg * - @ref riemann_zeta "riemann_zeta - The Riemann zeta function"
1111debfc3dSmrg * - @ref sph_bessel "sph_bessel - Spherical Bessel functions"
1121debfc3dSmrg * - @ref sph_legendre "sph_legendre - Spherical Legendre functions"
1131debfc3dSmrg * - @ref sph_neumann "sph_neumann - Spherical Neumann functions"
1141debfc3dSmrg *
1151debfc3dSmrg * The hypergeometric functions were stricken from the TR29124 and C++17
1161debfc3dSmrg * versions of this math library because of implementation concerns.
1171debfc3dSmrg * However, since they were in the TR1 version and since they are popular
1181debfc3dSmrg * we kept them as an extension in namespace @c __gnu_cxx:
119a2dc1f3fSmrg * - @ref __gnu_cxx::conf_hyperg "conf_hyperg - Confluent hypergeometric functions"
120a2dc1f3fSmrg * - @ref __gnu_cxx::hyperg "hyperg - Hypergeometric functions"
1211debfc3dSmrg *
122*8feb0f0bSmrg * <!-- @subsection mathsf_general General Features -->
1231debfc3dSmrg *
124*8feb0f0bSmrg * @subsection mathsf_promotion Argument Promotion
1251debfc3dSmrg * The arguments suppled to the non-suffixed functions will be promoted
1261debfc3dSmrg * according to the following rules:
1271debfc3dSmrg * 1. If any argument intended to be floating point is given an integral value
1281debfc3dSmrg * That integral value is promoted to double.
1291debfc3dSmrg * 2. All floating point arguments are promoted up to the largest floating
1301debfc3dSmrg * point precision among them.
1311debfc3dSmrg *
132*8feb0f0bSmrg * @subsection mathsf_NaN NaN Arguments
1331debfc3dSmrg * If any of the floating point arguments supplied to these functions is
1341debfc3dSmrg * invalid or NaN (std::numeric_limits<Tp>::quiet_NaN),
1351debfc3dSmrg * the value NaN is returned.
1361debfc3dSmrg *
137*8feb0f0bSmrg * @subsection mathsf_impl Implementation
1381debfc3dSmrg *
1391debfc3dSmrg * We strive to implement the underlying math with type generic algorithms
1401debfc3dSmrg * to the greatest extent possible. In practice, the functions are thin
1411debfc3dSmrg * wrappers that dispatch to function templates. Type dependence is
1421debfc3dSmrg * controlled with std::numeric_limits and functions thereof.
1431debfc3dSmrg *
1441debfc3dSmrg * We don't promote @c float to @c double or @c double to <tt>long double</tt>
1451debfc3dSmrg * reflexively. The goal is for @c float functions to operate more quickly,
1461debfc3dSmrg * at the cost of @c float accuracy and possibly a smaller domain of validity.
1471debfc3dSmrg * Similaryly, <tt>long double</tt> should give you more dynamic range
1481debfc3dSmrg * and slightly more pecision than @c double on many systems.
1491debfc3dSmrg *
150*8feb0f0bSmrg * @subsection mathsf_testing Testing
1511debfc3dSmrg *
1521debfc3dSmrg * These functions have been tested against equivalent implementations
1531debfc3dSmrg * from the <a href="http://www.gnu.org/software/gsl">
1541debfc3dSmrg * Gnu Scientific Library, GSL</a> and
155*8feb0f0bSmrg * <a href="http://www.boost.org/doc/libs/1_60_0/libs/math/doc/html/index.html">Boost</a>
1561debfc3dSmrg * and the ratio
1571debfc3dSmrg * @f[
1581debfc3dSmrg * \frac{|f - f_{test}|}{|f_{test}|}
1591debfc3dSmrg * @f]
160*8feb0f0bSmrg * is generally found to be within 10<sup>-15</sup> for 64-bit double on
161*8feb0f0bSmrg * linux-x86_64 systems over most of the ranges of validity.
1621debfc3dSmrg *
1631debfc3dSmrg * @todo Provide accuracy comparisons on a per-function basis for a small
1641debfc3dSmrg * number of targets.
1651debfc3dSmrg *
166*8feb0f0bSmrg * @subsection mathsf_bibliography General Bibliography
1671debfc3dSmrg *
1681debfc3dSmrg * @see Abramowitz and Stegun: Handbook of Mathematical Functions,
1691debfc3dSmrg * with Formulas, Graphs, and Mathematical Tables
1701debfc3dSmrg * Edited by Milton Abramowitz and Irene A. Stegun,
1711debfc3dSmrg * National Bureau of Standards Applied Mathematics Series - 55
1721debfc3dSmrg * Issued June 1964, Tenth Printing, December 1972, with corrections
1731debfc3dSmrg * Electronic versions of A&S abound including both pdf and navigable html.
1741debfc3dSmrg * @see for example http://people.math.sfu.ca/~cbm/aands/
1751debfc3dSmrg *
1761debfc3dSmrg * @see The old A&S has been redone as the
1771debfc3dSmrg * NIST Digital Library of Mathematical Functions: http://dlmf.nist.gov/
1781debfc3dSmrg * This version is far more navigable and includes more recent work.
1791debfc3dSmrg *
1801debfc3dSmrg * @see An Atlas of Functions: with Equator, the Atlas Function Calculator
1811debfc3dSmrg * 2nd Edition, by Oldham, Keith B., Myland, Jan, Spanier, Jerome
1821debfc3dSmrg *
1831debfc3dSmrg * @see Asymptotics and Special Functions by Frank W. J. Olver,
1841debfc3dSmrg * Academic Press, 1974
1851debfc3dSmrg *
1861debfc3dSmrg * @see Numerical Recipes in C, The Art of Scientific Computing,
1871debfc3dSmrg * by William H. Press, Second Ed., Saul A. Teukolsky,
1881debfc3dSmrg * William T. Vetterling, and Brian P. Flannery,
1891debfc3dSmrg * Cambridge University Press, 1992
1901debfc3dSmrg *
1911debfc3dSmrg * @see The Special Functions and Their Approximations: Volumes 1 and 2,
1921debfc3dSmrg * by Yudell L. Luke, Academic Press, 1969
193*8feb0f0bSmrg *
194*8feb0f0bSmrg * @{
1951debfc3dSmrg */
1961debfc3dSmrg
1971debfc3dSmrg // Associated Laguerre polynomials
1981debfc3dSmrg
1991debfc3dSmrg /**
2001debfc3dSmrg * Return the associated Laguerre polynomial of order @c n,
2011debfc3dSmrg * degree @c m: @f$ L_n^m(x) @f$ for @c float argument.
2021debfc3dSmrg *
2031debfc3dSmrg * @see assoc_laguerre for more details.
2041debfc3dSmrg */
2051debfc3dSmrg inline float
2061debfc3dSmrg assoc_laguerref(unsigned int __n, unsigned int __m, float __x)
2071debfc3dSmrg { return __detail::__assoc_laguerre<float>(__n, __m, __x); }
2081debfc3dSmrg
2091debfc3dSmrg /**
2101debfc3dSmrg * Return the associated Laguerre polynomial of order @c n,
2111debfc3dSmrg * degree @c m: @f$ L_n^m(x) @f$.
2121debfc3dSmrg *
2131debfc3dSmrg * @see assoc_laguerre for more details.
2141debfc3dSmrg */
2151debfc3dSmrg inline long double
2161debfc3dSmrg assoc_laguerrel(unsigned int __n, unsigned int __m, long double __x)
2171debfc3dSmrg { return __detail::__assoc_laguerre<long double>(__n, __m, __x); }
2181debfc3dSmrg
2191debfc3dSmrg /**
2201debfc3dSmrg * Return the associated Laguerre polynomial of nonnegative order @c n,
2211debfc3dSmrg * nonnegative degree @c m and real argument @c x: @f$ L_n^m(x) @f$.
2221debfc3dSmrg *
2231debfc3dSmrg * The associated Laguerre function of real degree @f$ \alpha @f$,
2241debfc3dSmrg * @f$ L_n^\alpha(x) @f$, is defined by
2251debfc3dSmrg * @f[
2261debfc3dSmrg * L_n^\alpha(x) = \frac{(\alpha + 1)_n}{n!}
2271debfc3dSmrg * {}_1F_1(-n; \alpha + 1; x)
2281debfc3dSmrg * @f]
2291debfc3dSmrg * where @f$ (\alpha)_n @f$ is the Pochhammer symbol and
2301debfc3dSmrg * @f$ {}_1F_1(a; c; x) @f$ is the confluent hypergeometric function.
2311debfc3dSmrg *
2321debfc3dSmrg * The associated Laguerre polynomial is defined for integral
2331debfc3dSmrg * degree @f$ \alpha = m @f$ by:
2341debfc3dSmrg * @f[
2351debfc3dSmrg * L_n^m(x) = (-1)^m \frac{d^m}{dx^m} L_{n + m}(x)
2361debfc3dSmrg * @f]
2371debfc3dSmrg * where the Laguerre polynomial is defined by:
2381debfc3dSmrg * @f[
2391debfc3dSmrg * L_n(x) = \frac{e^x}{n!} \frac{d^n}{dx^n} (x^ne^{-x})
2401debfc3dSmrg * @f]
2411debfc3dSmrg * and @f$ x >= 0 @f$.
2421debfc3dSmrg * @see laguerre for details of the Laguerre function of degree @c n
2431debfc3dSmrg *
2441debfc3dSmrg * @tparam _Tp The floating-point type of the argument @c __x.
2451debfc3dSmrg * @param __n The order of the Laguerre function, <tt>__n >= 0</tt>.
2461debfc3dSmrg * @param __m The degree of the Laguerre function, <tt>__m >= 0</tt>.
2471debfc3dSmrg * @param __x The argument of the Laguerre function, <tt>__x >= 0</tt>.
2481debfc3dSmrg * @throw std::domain_error if <tt>__x < 0</tt>.
2491debfc3dSmrg */
2501debfc3dSmrg template<typename _Tp>
2511debfc3dSmrg inline typename __gnu_cxx::__promote<_Tp>::__type
2521debfc3dSmrg assoc_laguerre(unsigned int __n, unsigned int __m, _Tp __x)
2531debfc3dSmrg {
2541debfc3dSmrg typedef typename __gnu_cxx::__promote<_Tp>::__type __type;
2551debfc3dSmrg return __detail::__assoc_laguerre<__type>(__n, __m, __x);
2561debfc3dSmrg }
2571debfc3dSmrg
2581debfc3dSmrg // Associated Legendre functions
2591debfc3dSmrg
2601debfc3dSmrg /**
2611debfc3dSmrg * Return the associated Legendre function of degree @c l and order @c m
2621debfc3dSmrg * for @c float argument.
2631debfc3dSmrg *
2641debfc3dSmrg * @see assoc_legendre for more details.
2651debfc3dSmrg */
2661debfc3dSmrg inline float
2671debfc3dSmrg assoc_legendref(unsigned int __l, unsigned int __m, float __x)
2681debfc3dSmrg { return __detail::__assoc_legendre_p<float>(__l, __m, __x); }
2691debfc3dSmrg
2701debfc3dSmrg /**
2711debfc3dSmrg * Return the associated Legendre function of degree @c l and order @c m.
2721debfc3dSmrg *
2731debfc3dSmrg * @see assoc_legendre for more details.
2741debfc3dSmrg */
2751debfc3dSmrg inline long double
2761debfc3dSmrg assoc_legendrel(unsigned int __l, unsigned int __m, long double __x)
2771debfc3dSmrg { return __detail::__assoc_legendre_p<long double>(__l, __m, __x); }
2781debfc3dSmrg
2791debfc3dSmrg
2801debfc3dSmrg /**
2811debfc3dSmrg * Return the associated Legendre function of degree @c l and order @c m.
2821debfc3dSmrg *
2831debfc3dSmrg * The associated Legendre function is derived from the Legendre function
2841debfc3dSmrg * @f$ P_l(x) @f$ by the Rodrigues formula:
2851debfc3dSmrg * @f[
2861debfc3dSmrg * P_l^m(x) = (1 - x^2)^{m/2}\frac{d^m}{dx^m}P_l(x)
2871debfc3dSmrg * @f]
2881debfc3dSmrg * @see legendre for details of the Legendre function of degree @c l
2891debfc3dSmrg *
2901debfc3dSmrg * @tparam _Tp The floating-point type of the argument @c __x.
2911debfc3dSmrg * @param __l The degree <tt>__l >= 0</tt>.
2921debfc3dSmrg * @param __m The order <tt>__m <= l</tt>.
2931debfc3dSmrg * @param __x The argument, <tt>abs(__x) <= 1</tt>.
2941debfc3dSmrg * @throw std::domain_error if <tt>abs(__x) > 1</tt>.
2951debfc3dSmrg */
2961debfc3dSmrg template<typename _Tp>
2971debfc3dSmrg inline typename __gnu_cxx::__promote<_Tp>::__type
2981debfc3dSmrg assoc_legendre(unsigned int __l, unsigned int __m, _Tp __x)
2991debfc3dSmrg {
3001debfc3dSmrg typedef typename __gnu_cxx::__promote<_Tp>::__type __type;
3011debfc3dSmrg return __detail::__assoc_legendre_p<__type>(__l, __m, __x);
3021debfc3dSmrg }
3031debfc3dSmrg
3041debfc3dSmrg // Beta functions
3051debfc3dSmrg
3061debfc3dSmrg /**
3071debfc3dSmrg * Return the beta function, @f$ B(a,b) @f$, for @c float parameters @c a, @c b.
3081debfc3dSmrg *
3091debfc3dSmrg * @see beta for more details.
3101debfc3dSmrg */
3111debfc3dSmrg inline float
3121debfc3dSmrg betaf(float __a, float __b)
3131debfc3dSmrg { return __detail::__beta<float>(__a, __b); }
3141debfc3dSmrg
3151debfc3dSmrg /**
3161debfc3dSmrg * Return the beta function, @f$B(a,b)@f$, for long double
3171debfc3dSmrg * parameters @c a, @c b.
3181debfc3dSmrg *
3191debfc3dSmrg * @see beta for more details.
3201debfc3dSmrg */
3211debfc3dSmrg inline long double
3221debfc3dSmrg betal(long double __a, long double __b)
3231debfc3dSmrg { return __detail::__beta<long double>(__a, __b); }
3241debfc3dSmrg
3251debfc3dSmrg /**
3261debfc3dSmrg * Return the beta function, @f$B(a,b)@f$, for real parameters @c a, @c b.
3271debfc3dSmrg *
3281debfc3dSmrg * The beta function is defined by
3291debfc3dSmrg * @f[
3301debfc3dSmrg * B(a,b) = \int_0^1 t^{a - 1} (1 - t)^{b - 1} dt
3311debfc3dSmrg * = \frac{\Gamma(a)\Gamma(b)}{\Gamma(a+b)}
3321debfc3dSmrg * @f]
3331debfc3dSmrg * where @f$ a > 0 @f$ and @f$ b > 0 @f$
3341debfc3dSmrg *
3351debfc3dSmrg * @tparam _Tpa The floating-point type of the parameter @c __a.
3361debfc3dSmrg * @tparam _Tpb The floating-point type of the parameter @c __b.
3371debfc3dSmrg * @param __a The first argument of the beta function, <tt> __a > 0 </tt>.
3381debfc3dSmrg * @param __b The second argument of the beta function, <tt> __b > 0 </tt>.
3391debfc3dSmrg * @throw std::domain_error if <tt> __a < 0 </tt> or <tt> __b < 0 </tt>.
3401debfc3dSmrg */
3411debfc3dSmrg template<typename _Tpa, typename _Tpb>
3421debfc3dSmrg inline typename __gnu_cxx::__promote_2<_Tpa, _Tpb>::__type
3431debfc3dSmrg beta(_Tpa __a, _Tpb __b)
3441debfc3dSmrg {
3451debfc3dSmrg typedef typename __gnu_cxx::__promote_2<_Tpa, _Tpb>::__type __type;
3461debfc3dSmrg return __detail::__beta<__type>(__a, __b);
3471debfc3dSmrg }
3481debfc3dSmrg
3491debfc3dSmrg // Complete elliptic integrals of the first kind
3501debfc3dSmrg
3511debfc3dSmrg /**
3521debfc3dSmrg * Return the complete elliptic integral of the first kind @f$ E(k) @f$
3531debfc3dSmrg * for @c float modulus @c k.
3541debfc3dSmrg *
3551debfc3dSmrg * @see comp_ellint_1 for details.
3561debfc3dSmrg */
3571debfc3dSmrg inline float
3581debfc3dSmrg comp_ellint_1f(float __k)
3591debfc3dSmrg { return __detail::__comp_ellint_1<float>(__k); }
3601debfc3dSmrg
3611debfc3dSmrg /**
3621debfc3dSmrg * Return the complete elliptic integral of the first kind @f$ E(k) @f$
3631debfc3dSmrg * for long double modulus @c k.
3641debfc3dSmrg *
3651debfc3dSmrg * @see comp_ellint_1 for details.
3661debfc3dSmrg */
3671debfc3dSmrg inline long double
3681debfc3dSmrg comp_ellint_1l(long double __k)
3691debfc3dSmrg { return __detail::__comp_ellint_1<long double>(__k); }
3701debfc3dSmrg
3711debfc3dSmrg /**
3721debfc3dSmrg * Return the complete elliptic integral of the first kind
3731debfc3dSmrg * @f$ K(k) @f$ for real modulus @c k.
3741debfc3dSmrg *
3751debfc3dSmrg * The complete elliptic integral of the first kind is defined as
3761debfc3dSmrg * @f[
3771debfc3dSmrg * K(k) = F(k,\pi/2) = \int_0^{\pi/2}\frac{d\theta}
3781debfc3dSmrg * {\sqrt{1 - k^2 sin^2\theta}}
3791debfc3dSmrg * @f]
3801debfc3dSmrg * where @f$ F(k,\phi) @f$ is the incomplete elliptic integral of the
3811debfc3dSmrg * first kind and the modulus @f$ |k| <= 1 @f$.
3821debfc3dSmrg * @see ellint_1 for details of the incomplete elliptic function
3831debfc3dSmrg * of the first kind.
3841debfc3dSmrg *
3851debfc3dSmrg * @tparam _Tp The floating-point type of the modulus @c __k.
3861debfc3dSmrg * @param __k The modulus, <tt> abs(__k) <= 1 </tt>
3871debfc3dSmrg * @throw std::domain_error if <tt> abs(__k) > 1 </tt>.
3881debfc3dSmrg */
3891debfc3dSmrg template<typename _Tp>
3901debfc3dSmrg inline typename __gnu_cxx::__promote<_Tp>::__type
3911debfc3dSmrg comp_ellint_1(_Tp __k)
3921debfc3dSmrg {
3931debfc3dSmrg typedef typename __gnu_cxx::__promote<_Tp>::__type __type;
3941debfc3dSmrg return __detail::__comp_ellint_1<__type>(__k);
3951debfc3dSmrg }
3961debfc3dSmrg
3971debfc3dSmrg // Complete elliptic integrals of the second kind
3981debfc3dSmrg
3991debfc3dSmrg /**
4001debfc3dSmrg * Return the complete elliptic integral of the second kind @f$ E(k) @f$
4011debfc3dSmrg * for @c float modulus @c k.
4021debfc3dSmrg *
4031debfc3dSmrg * @see comp_ellint_2 for details.
4041debfc3dSmrg */
4051debfc3dSmrg inline float
4061debfc3dSmrg comp_ellint_2f(float __k)
4071debfc3dSmrg { return __detail::__comp_ellint_2<float>(__k); }
4081debfc3dSmrg
4091debfc3dSmrg /**
4101debfc3dSmrg * Return the complete elliptic integral of the second kind @f$ E(k) @f$
4111debfc3dSmrg * for long double modulus @c k.
4121debfc3dSmrg *
4131debfc3dSmrg * @see comp_ellint_2 for details.
4141debfc3dSmrg */
4151debfc3dSmrg inline long double
4161debfc3dSmrg comp_ellint_2l(long double __k)
4171debfc3dSmrg { return __detail::__comp_ellint_2<long double>(__k); }
4181debfc3dSmrg
4191debfc3dSmrg /**
4201debfc3dSmrg * Return the complete elliptic integral of the second kind @f$ E(k) @f$
4211debfc3dSmrg * for real modulus @c k.
4221debfc3dSmrg *
4231debfc3dSmrg * The complete elliptic integral of the second kind is defined as
4241debfc3dSmrg * @f[
4251debfc3dSmrg * E(k) = E(k,\pi/2) = \int_0^{\pi/2}\sqrt{1 - k^2 sin^2\theta}
4261debfc3dSmrg * @f]
4271debfc3dSmrg * where @f$ E(k,\phi) @f$ is the incomplete elliptic integral of the
4281debfc3dSmrg * second kind and the modulus @f$ |k| <= 1 @f$.
4291debfc3dSmrg * @see ellint_2 for details of the incomplete elliptic function
4301debfc3dSmrg * of the second kind.
4311debfc3dSmrg *
4321debfc3dSmrg * @tparam _Tp The floating-point type of the modulus @c __k.
4331debfc3dSmrg * @param __k The modulus, @c abs(__k) <= 1
4341debfc3dSmrg * @throw std::domain_error if @c abs(__k) > 1.
4351debfc3dSmrg */
4361debfc3dSmrg template<typename _Tp>
4371debfc3dSmrg inline typename __gnu_cxx::__promote<_Tp>::__type
4381debfc3dSmrg comp_ellint_2(_Tp __k)
4391debfc3dSmrg {
4401debfc3dSmrg typedef typename __gnu_cxx::__promote<_Tp>::__type __type;
4411debfc3dSmrg return __detail::__comp_ellint_2<__type>(__k);
4421debfc3dSmrg }
4431debfc3dSmrg
4441debfc3dSmrg // Complete elliptic integrals of the third kind
4451debfc3dSmrg
4461debfc3dSmrg /**
4471debfc3dSmrg * @brief Return the complete elliptic integral of the third kind
4481debfc3dSmrg * @f$ \Pi(k,\nu) @f$ for @c float modulus @c k.
4491debfc3dSmrg *
4501debfc3dSmrg * @see comp_ellint_3 for details.
4511debfc3dSmrg */
4521debfc3dSmrg inline float
4531debfc3dSmrg comp_ellint_3f(float __k, float __nu)
4541debfc3dSmrg { return __detail::__comp_ellint_3<float>(__k, __nu); }
4551debfc3dSmrg
4561debfc3dSmrg /**
4571debfc3dSmrg * @brief Return the complete elliptic integral of the third kind
4581debfc3dSmrg * @f$ \Pi(k,\nu) @f$ for <tt>long double</tt> modulus @c k.
4591debfc3dSmrg *
4601debfc3dSmrg * @see comp_ellint_3 for details.
4611debfc3dSmrg */
4621debfc3dSmrg inline long double
4631debfc3dSmrg comp_ellint_3l(long double __k, long double __nu)
4641debfc3dSmrg { return __detail::__comp_ellint_3<long double>(__k, __nu); }
4651debfc3dSmrg
4661debfc3dSmrg /**
4671debfc3dSmrg * Return the complete elliptic integral of the third kind
4681debfc3dSmrg * @f$ \Pi(k,\nu) = \Pi(k,\nu,\pi/2) @f$ for real modulus @c k.
4691debfc3dSmrg *
4701debfc3dSmrg * The complete elliptic integral of the third kind is defined as
4711debfc3dSmrg * @f[
4721debfc3dSmrg * \Pi(k,\nu) = \Pi(k,\nu,\pi/2) = \int_0^{\pi/2}
4731debfc3dSmrg * \frac{d\theta}
4741debfc3dSmrg * {(1 - \nu \sin^2\theta)\sqrt{1 - k^2 \sin^2\theta}}
4751debfc3dSmrg * @f]
4761debfc3dSmrg * where @f$ \Pi(k,\nu,\phi) @f$ is the incomplete elliptic integral of the
4771debfc3dSmrg * second kind and the modulus @f$ |k| <= 1 @f$.
4781debfc3dSmrg * @see ellint_3 for details of the incomplete elliptic function
4791debfc3dSmrg * of the third kind.
4801debfc3dSmrg *
4811debfc3dSmrg * @tparam _Tp The floating-point type of the modulus @c __k.
4821debfc3dSmrg * @tparam _Tpn The floating-point type of the argument @c __nu.
4831debfc3dSmrg * @param __k The modulus, @c abs(__k) <= 1
4841debfc3dSmrg * @param __nu The argument
4851debfc3dSmrg * @throw std::domain_error if @c abs(__k) > 1.
4861debfc3dSmrg */
4871debfc3dSmrg template<typename _Tp, typename _Tpn>
4881debfc3dSmrg inline typename __gnu_cxx::__promote_2<_Tp, _Tpn>::__type
4891debfc3dSmrg comp_ellint_3(_Tp __k, _Tpn __nu)
4901debfc3dSmrg {
4911debfc3dSmrg typedef typename __gnu_cxx::__promote_2<_Tp, _Tpn>::__type __type;
4921debfc3dSmrg return __detail::__comp_ellint_3<__type>(__k, __nu);
4931debfc3dSmrg }
4941debfc3dSmrg
4951debfc3dSmrg // Regular modified cylindrical Bessel functions
4961debfc3dSmrg
4971debfc3dSmrg /**
4981debfc3dSmrg * Return the regular modified Bessel function @f$ I_{\nu}(x) @f$
4991debfc3dSmrg * for @c float order @f$ \nu @f$ and argument @f$ x >= 0 @f$.
5001debfc3dSmrg *
5011debfc3dSmrg * @see cyl_bessel_i for setails.
5021debfc3dSmrg */
5031debfc3dSmrg inline float
5041debfc3dSmrg cyl_bessel_if(float __nu, float __x)
5051debfc3dSmrg { return __detail::__cyl_bessel_i<float>(__nu, __x); }
5061debfc3dSmrg
5071debfc3dSmrg /**
5081debfc3dSmrg * Return the regular modified Bessel function @f$ I_{\nu}(x) @f$
5091debfc3dSmrg * for <tt>long double</tt> order @f$ \nu @f$ and argument @f$ x >= 0 @f$.
5101debfc3dSmrg *
5111debfc3dSmrg * @see cyl_bessel_i for setails.
5121debfc3dSmrg */
5131debfc3dSmrg inline long double
5141debfc3dSmrg cyl_bessel_il(long double __nu, long double __x)
5151debfc3dSmrg { return __detail::__cyl_bessel_i<long double>(__nu, __x); }
5161debfc3dSmrg
5171debfc3dSmrg /**
5181debfc3dSmrg * Return the regular modified Bessel function @f$ I_{\nu}(x) @f$
5191debfc3dSmrg * for real order @f$ \nu @f$ and argument @f$ x >= 0 @f$.
5201debfc3dSmrg *
5211debfc3dSmrg * The regular modified cylindrical Bessel function is:
5221debfc3dSmrg * @f[
5231debfc3dSmrg * I_{\nu}(x) = i^{-\nu}J_\nu(ix) = \sum_{k=0}^{\infty}
5241debfc3dSmrg * \frac{(x/2)^{\nu + 2k}}{k!\Gamma(\nu+k+1)}
5251debfc3dSmrg * @f]
5261debfc3dSmrg *
5271debfc3dSmrg * @tparam _Tpnu The floating-point type of the order @c __nu.
5281debfc3dSmrg * @tparam _Tp The floating-point type of the argument @c __x.
5291debfc3dSmrg * @param __nu The order
5301debfc3dSmrg * @param __x The argument, <tt> __x >= 0 </tt>
5311debfc3dSmrg * @throw std::domain_error if <tt> __x < 0 </tt>.
5321debfc3dSmrg */
5331debfc3dSmrg template<typename _Tpnu, typename _Tp>
5341debfc3dSmrg inline typename __gnu_cxx::__promote_2<_Tpnu, _Tp>::__type
5351debfc3dSmrg cyl_bessel_i(_Tpnu __nu, _Tp __x)
5361debfc3dSmrg {
5371debfc3dSmrg typedef typename __gnu_cxx::__promote_2<_Tpnu, _Tp>::__type __type;
5381debfc3dSmrg return __detail::__cyl_bessel_i<__type>(__nu, __x);
5391debfc3dSmrg }
5401debfc3dSmrg
5411debfc3dSmrg // Cylindrical Bessel functions (of the first kind)
5421debfc3dSmrg
5431debfc3dSmrg /**
5441debfc3dSmrg * Return the Bessel function of the first kind @f$ J_{\nu}(x) @f$
5451debfc3dSmrg * for @c float order @f$ \nu @f$ and argument @f$ x >= 0 @f$.
5461debfc3dSmrg *
5471debfc3dSmrg * @see cyl_bessel_j for setails.
5481debfc3dSmrg */
5491debfc3dSmrg inline float
5501debfc3dSmrg cyl_bessel_jf(float __nu, float __x)
5511debfc3dSmrg { return __detail::__cyl_bessel_j<float>(__nu, __x); }
5521debfc3dSmrg
5531debfc3dSmrg /**
5541debfc3dSmrg * Return the Bessel function of the first kind @f$ J_{\nu}(x) @f$
5551debfc3dSmrg * for <tt>long double</tt> order @f$ \nu @f$ and argument @f$ x >= 0 @f$.
5561debfc3dSmrg *
5571debfc3dSmrg * @see cyl_bessel_j for setails.
5581debfc3dSmrg */
5591debfc3dSmrg inline long double
5601debfc3dSmrg cyl_bessel_jl(long double __nu, long double __x)
5611debfc3dSmrg { return __detail::__cyl_bessel_j<long double>(__nu, __x); }
5621debfc3dSmrg
5631debfc3dSmrg /**
5641debfc3dSmrg * Return the Bessel function @f$ J_{\nu}(x) @f$ of real order @f$ \nu @f$
5651debfc3dSmrg * and argument @f$ x >= 0 @f$.
5661debfc3dSmrg *
5671debfc3dSmrg * The cylindrical Bessel function is:
5681debfc3dSmrg * @f[
5691debfc3dSmrg * J_{\nu}(x) = \sum_{k=0}^{\infty}
5701debfc3dSmrg * \frac{(-1)^k (x/2)^{\nu + 2k}}{k!\Gamma(\nu+k+1)}
5711debfc3dSmrg * @f]
5721debfc3dSmrg *
5731debfc3dSmrg * @tparam _Tpnu The floating-point type of the order @c __nu.
5741debfc3dSmrg * @tparam _Tp The floating-point type of the argument @c __x.
5751debfc3dSmrg * @param __nu The order
5761debfc3dSmrg * @param __x The argument, <tt> __x >= 0 </tt>
5771debfc3dSmrg * @throw std::domain_error if <tt> __x < 0 </tt>.
5781debfc3dSmrg */
5791debfc3dSmrg template<typename _Tpnu, typename _Tp>
5801debfc3dSmrg inline typename __gnu_cxx::__promote_2<_Tpnu, _Tp>::__type
5811debfc3dSmrg cyl_bessel_j(_Tpnu __nu, _Tp __x)
5821debfc3dSmrg {
5831debfc3dSmrg typedef typename __gnu_cxx::__promote_2<_Tpnu, _Tp>::__type __type;
5841debfc3dSmrg return __detail::__cyl_bessel_j<__type>(__nu, __x);
5851debfc3dSmrg }
5861debfc3dSmrg
5871debfc3dSmrg // Irregular modified cylindrical Bessel functions
5881debfc3dSmrg
5891debfc3dSmrg /**
5901debfc3dSmrg * Return the irregular modified Bessel function @f$ K_{\nu}(x) @f$
5911debfc3dSmrg * for @c float order @f$ \nu @f$ and argument @f$ x >= 0 @f$.
5921debfc3dSmrg *
5931debfc3dSmrg * @see cyl_bessel_k for setails.
5941debfc3dSmrg */
5951debfc3dSmrg inline float
5961debfc3dSmrg cyl_bessel_kf(float __nu, float __x)
5971debfc3dSmrg { return __detail::__cyl_bessel_k<float>(__nu, __x); }
5981debfc3dSmrg
5991debfc3dSmrg /**
6001debfc3dSmrg * Return the irregular modified Bessel function @f$ K_{\nu}(x) @f$
6011debfc3dSmrg * for <tt>long double</tt> order @f$ \nu @f$ and argument @f$ x >= 0 @f$.
6021debfc3dSmrg *
6031debfc3dSmrg * @see cyl_bessel_k for setails.
6041debfc3dSmrg */
6051debfc3dSmrg inline long double
6061debfc3dSmrg cyl_bessel_kl(long double __nu, long double __x)
6071debfc3dSmrg { return __detail::__cyl_bessel_k<long double>(__nu, __x); }
6081debfc3dSmrg
6091debfc3dSmrg /**
6101debfc3dSmrg * Return the irregular modified Bessel function @f$ K_{\nu}(x) @f$
6111debfc3dSmrg * of real order @f$ \nu @f$ and argument @f$ x @f$.
6121debfc3dSmrg *
6131debfc3dSmrg * The irregular modified Bessel function is defined by:
6141debfc3dSmrg * @f[
6151debfc3dSmrg * K_{\nu}(x) = \frac{\pi}{2}
6161debfc3dSmrg * \frac{I_{-\nu}(x) - I_{\nu}(x)}{\sin \nu\pi}
6171debfc3dSmrg * @f]
6181debfc3dSmrg * where for integral @f$ \nu = n @f$ a limit is taken:
6191debfc3dSmrg * @f$ lim_{\nu \to n} @f$.
6201debfc3dSmrg * For negative argument we have simply:
6211debfc3dSmrg * @f[
6221debfc3dSmrg * K_{-\nu}(x) = K_{\nu}(x)
6231debfc3dSmrg * @f]
6241debfc3dSmrg *
6251debfc3dSmrg * @tparam _Tpnu The floating-point type of the order @c __nu.
6261debfc3dSmrg * @tparam _Tp The floating-point type of the argument @c __x.
6271debfc3dSmrg * @param __nu The order
6281debfc3dSmrg * @param __x The argument, <tt> __x >= 0 </tt>
6291debfc3dSmrg * @throw std::domain_error if <tt> __x < 0 </tt>.
6301debfc3dSmrg */
6311debfc3dSmrg template<typename _Tpnu, typename _Tp>
6321debfc3dSmrg inline typename __gnu_cxx::__promote_2<_Tpnu, _Tp>::__type
6331debfc3dSmrg cyl_bessel_k(_Tpnu __nu, _Tp __x)
6341debfc3dSmrg {
6351debfc3dSmrg typedef typename __gnu_cxx::__promote_2<_Tpnu, _Tp>::__type __type;
6361debfc3dSmrg return __detail::__cyl_bessel_k<__type>(__nu, __x);
6371debfc3dSmrg }
6381debfc3dSmrg
6391debfc3dSmrg // Cylindrical Neumann functions
6401debfc3dSmrg
6411debfc3dSmrg /**
6421debfc3dSmrg * Return the Neumann function @f$ N_{\nu}(x) @f$
6431debfc3dSmrg * of @c float order @f$ \nu @f$ and argument @f$ x @f$.
6441debfc3dSmrg *
6451debfc3dSmrg * @see cyl_neumann for setails.
6461debfc3dSmrg */
6471debfc3dSmrg inline float
6481debfc3dSmrg cyl_neumannf(float __nu, float __x)
6491debfc3dSmrg { return __detail::__cyl_neumann_n<float>(__nu, __x); }
6501debfc3dSmrg
6511debfc3dSmrg /**
6521debfc3dSmrg * Return the Neumann function @f$ N_{\nu}(x) @f$
6531debfc3dSmrg * of <tt>long double</tt> order @f$ \nu @f$ and argument @f$ x @f$.
6541debfc3dSmrg *
6551debfc3dSmrg * @see cyl_neumann for setails.
6561debfc3dSmrg */
6571debfc3dSmrg inline long double
6581debfc3dSmrg cyl_neumannl(long double __nu, long double __x)
6591debfc3dSmrg { return __detail::__cyl_neumann_n<long double>(__nu, __x); }
6601debfc3dSmrg
6611debfc3dSmrg /**
6621debfc3dSmrg * Return the Neumann function @f$ N_{\nu}(x) @f$
6631debfc3dSmrg * of real order @f$ \nu @f$ and argument @f$ x >= 0 @f$.
6641debfc3dSmrg *
6651debfc3dSmrg * The Neumann function is defined by:
6661debfc3dSmrg * @f[
6671debfc3dSmrg * N_{\nu}(x) = \frac{J_{\nu}(x) \cos \nu\pi - J_{-\nu}(x)}
6681debfc3dSmrg * {\sin \nu\pi}
6691debfc3dSmrg * @f]
6701debfc3dSmrg * where @f$ x >= 0 @f$ and for integral order @f$ \nu = n @f$
6711debfc3dSmrg * a limit is taken: @f$ lim_{\nu \to n} @f$.
6721debfc3dSmrg *
6731debfc3dSmrg * @tparam _Tpnu The floating-point type of the order @c __nu.
6741debfc3dSmrg * @tparam _Tp The floating-point type of the argument @c __x.
6751debfc3dSmrg * @param __nu The order
6761debfc3dSmrg * @param __x The argument, <tt> __x >= 0 </tt>
6771debfc3dSmrg * @throw std::domain_error if <tt> __x < 0 </tt>.
6781debfc3dSmrg */
6791debfc3dSmrg template<typename _Tpnu, typename _Tp>
6801debfc3dSmrg inline typename __gnu_cxx::__promote_2<_Tpnu, _Tp>::__type
6811debfc3dSmrg cyl_neumann(_Tpnu __nu, _Tp __x)
6821debfc3dSmrg {
6831debfc3dSmrg typedef typename __gnu_cxx::__promote_2<_Tpnu, _Tp>::__type __type;
6841debfc3dSmrg return __detail::__cyl_neumann_n<__type>(__nu, __x);
6851debfc3dSmrg }
6861debfc3dSmrg
6871debfc3dSmrg // Incomplete elliptic integrals of the first kind
6881debfc3dSmrg
6891debfc3dSmrg /**
6901debfc3dSmrg * Return the incomplete elliptic integral of the first kind @f$ E(k,\phi) @f$
6911debfc3dSmrg * for @c float modulus @f$ k @f$ and angle @f$ \phi @f$.
6921debfc3dSmrg *
6931debfc3dSmrg * @see ellint_1 for details.
6941debfc3dSmrg */
6951debfc3dSmrg inline float
6961debfc3dSmrg ellint_1f(float __k, float __phi)
6971debfc3dSmrg { return __detail::__ellint_1<float>(__k, __phi); }
6981debfc3dSmrg
6991debfc3dSmrg /**
7001debfc3dSmrg * Return the incomplete elliptic integral of the first kind @f$ E(k,\phi) @f$
7011debfc3dSmrg * for <tt>long double</tt> modulus @f$ k @f$ and angle @f$ \phi @f$.
7021debfc3dSmrg *
7031debfc3dSmrg * @see ellint_1 for details.
7041debfc3dSmrg */
7051debfc3dSmrg inline long double
7061debfc3dSmrg ellint_1l(long double __k, long double __phi)
7071debfc3dSmrg { return __detail::__ellint_1<long double>(__k, __phi); }
7081debfc3dSmrg
7091debfc3dSmrg /**
7101debfc3dSmrg * Return the incomplete elliptic integral of the first kind @f$ F(k,\phi) @f$
7111debfc3dSmrg * for @c real modulus @f$ k @f$ and angle @f$ \phi @f$.
7121debfc3dSmrg *
7131debfc3dSmrg * The incomplete elliptic integral of the first kind is defined as
7141debfc3dSmrg * @f[
7151debfc3dSmrg * F(k,\phi) = \int_0^{\phi}\frac{d\theta}
7161debfc3dSmrg * {\sqrt{1 - k^2 sin^2\theta}}
7171debfc3dSmrg * @f]
7181debfc3dSmrg * For @f$ \phi= \pi/2 @f$ this becomes the complete elliptic integral of
7191debfc3dSmrg * the first kind, @f$ K(k) @f$. @see comp_ellint_1.
7201debfc3dSmrg *
7211debfc3dSmrg * @tparam _Tp The floating-point type of the modulus @c __k.
7221debfc3dSmrg * @tparam _Tpp The floating-point type of the angle @c __phi.
7231debfc3dSmrg * @param __k The modulus, <tt> abs(__k) <= 1 </tt>
7241debfc3dSmrg * @param __phi The integral limit argument in radians
7251debfc3dSmrg * @throw std::domain_error if <tt> abs(__k) > 1 </tt>.
7261debfc3dSmrg */
7271debfc3dSmrg template<typename _Tp, typename _Tpp>
7281debfc3dSmrg inline typename __gnu_cxx::__promote_2<_Tp, _Tpp>::__type
7291debfc3dSmrg ellint_1(_Tp __k, _Tpp __phi)
7301debfc3dSmrg {
7311debfc3dSmrg typedef typename __gnu_cxx::__promote_2<_Tp, _Tpp>::__type __type;
7321debfc3dSmrg return __detail::__ellint_1<__type>(__k, __phi);
7331debfc3dSmrg }
7341debfc3dSmrg
7351debfc3dSmrg // Incomplete elliptic integrals of the second kind
7361debfc3dSmrg
7371debfc3dSmrg /**
7381debfc3dSmrg * @brief Return the incomplete elliptic integral of the second kind
7391debfc3dSmrg * @f$ E(k,\phi) @f$ for @c float argument.
7401debfc3dSmrg *
7411debfc3dSmrg * @see ellint_2 for details.
7421debfc3dSmrg */
7431debfc3dSmrg inline float
7441debfc3dSmrg ellint_2f(float __k, float __phi)
7451debfc3dSmrg { return __detail::__ellint_2<float>(__k, __phi); }
7461debfc3dSmrg
7471debfc3dSmrg /**
7481debfc3dSmrg * @brief Return the incomplete elliptic integral of the second kind
7491debfc3dSmrg * @f$ E(k,\phi) @f$.
7501debfc3dSmrg *
7511debfc3dSmrg * @see ellint_2 for details.
7521debfc3dSmrg */
7531debfc3dSmrg inline long double
7541debfc3dSmrg ellint_2l(long double __k, long double __phi)
7551debfc3dSmrg { return __detail::__ellint_2<long double>(__k, __phi); }
7561debfc3dSmrg
7571debfc3dSmrg /**
7581debfc3dSmrg * Return the incomplete elliptic integral of the second kind
7591debfc3dSmrg * @f$ E(k,\phi) @f$.
7601debfc3dSmrg *
7611debfc3dSmrg * The incomplete elliptic integral of the second kind is defined as
7621debfc3dSmrg * @f[
7631debfc3dSmrg * E(k,\phi) = \int_0^{\phi} \sqrt{1 - k^2 sin^2\theta}
7641debfc3dSmrg * @f]
7651debfc3dSmrg * For @f$ \phi= \pi/2 @f$ this becomes the complete elliptic integral of
7661debfc3dSmrg * the second kind, @f$ E(k) @f$. @see comp_ellint_2.
7671debfc3dSmrg *
7681debfc3dSmrg * @tparam _Tp The floating-point type of the modulus @c __k.
7691debfc3dSmrg * @tparam _Tpp The floating-point type of the angle @c __phi.
7701debfc3dSmrg * @param __k The modulus, <tt> abs(__k) <= 1 </tt>
7711debfc3dSmrg * @param __phi The integral limit argument in radians
7721debfc3dSmrg * @return The elliptic function of the second kind.
7731debfc3dSmrg * @throw std::domain_error if <tt> abs(__k) > 1 </tt>.
7741debfc3dSmrg */
7751debfc3dSmrg template<typename _Tp, typename _Tpp>
7761debfc3dSmrg inline typename __gnu_cxx::__promote_2<_Tp, _Tpp>::__type
7771debfc3dSmrg ellint_2(_Tp __k, _Tpp __phi)
7781debfc3dSmrg {
7791debfc3dSmrg typedef typename __gnu_cxx::__promote_2<_Tp, _Tpp>::__type __type;
7801debfc3dSmrg return __detail::__ellint_2<__type>(__k, __phi);
7811debfc3dSmrg }
7821debfc3dSmrg
7831debfc3dSmrg // Incomplete elliptic integrals of the third kind
7841debfc3dSmrg
7851debfc3dSmrg /**
7861debfc3dSmrg * @brief Return the incomplete elliptic integral of the third kind
7871debfc3dSmrg * @f$ \Pi(k,\nu,\phi) @f$ for @c float argument.
7881debfc3dSmrg *
7891debfc3dSmrg * @see ellint_3 for details.
7901debfc3dSmrg */
7911debfc3dSmrg inline float
7921debfc3dSmrg ellint_3f(float __k, float __nu, float __phi)
7931debfc3dSmrg { return __detail::__ellint_3<float>(__k, __nu, __phi); }
7941debfc3dSmrg
7951debfc3dSmrg /**
7961debfc3dSmrg * @brief Return the incomplete elliptic integral of the third kind
7971debfc3dSmrg * @f$ \Pi(k,\nu,\phi) @f$.
7981debfc3dSmrg *
7991debfc3dSmrg * @see ellint_3 for details.
8001debfc3dSmrg */
8011debfc3dSmrg inline long double
8021debfc3dSmrg ellint_3l(long double __k, long double __nu, long double __phi)
8031debfc3dSmrg { return __detail::__ellint_3<long double>(__k, __nu, __phi); }
8041debfc3dSmrg
8051debfc3dSmrg /**
8061debfc3dSmrg * @brief Return the incomplete elliptic integral of the third kind
8071debfc3dSmrg * @f$ \Pi(k,\nu,\phi) @f$.
8081debfc3dSmrg *
8091debfc3dSmrg * The incomplete elliptic integral of the third kind is defined by:
8101debfc3dSmrg * @f[
8111debfc3dSmrg * \Pi(k,\nu,\phi) = \int_0^{\phi}
8121debfc3dSmrg * \frac{d\theta}
8131debfc3dSmrg * {(1 - \nu \sin^2\theta)
8141debfc3dSmrg * \sqrt{1 - k^2 \sin^2\theta}}
8151debfc3dSmrg * @f]
8161debfc3dSmrg * For @f$ \phi= \pi/2 @f$ this becomes the complete elliptic integral of
8171debfc3dSmrg * the third kind, @f$ \Pi(k,\nu) @f$. @see comp_ellint_3.
8181debfc3dSmrg *
8191debfc3dSmrg * @tparam _Tp The floating-point type of the modulus @c __k.
8201debfc3dSmrg * @tparam _Tpn The floating-point type of the argument @c __nu.
8211debfc3dSmrg * @tparam _Tpp The floating-point type of the angle @c __phi.
8221debfc3dSmrg * @param __k The modulus, <tt> abs(__k) <= 1 </tt>
8231debfc3dSmrg * @param __nu The second argument
8241debfc3dSmrg * @param __phi The integral limit argument in radians
8251debfc3dSmrg * @return The elliptic function of the third kind.
8261debfc3dSmrg * @throw std::domain_error if <tt> abs(__k) > 1 </tt>.
8271debfc3dSmrg */
8281debfc3dSmrg template<typename _Tp, typename _Tpn, typename _Tpp>
8291debfc3dSmrg inline typename __gnu_cxx::__promote_3<_Tp, _Tpn, _Tpp>::__type
8301debfc3dSmrg ellint_3(_Tp __k, _Tpn __nu, _Tpp __phi)
8311debfc3dSmrg {
8321debfc3dSmrg typedef typename __gnu_cxx::__promote_3<_Tp, _Tpn, _Tpp>::__type __type;
8331debfc3dSmrg return __detail::__ellint_3<__type>(__k, __nu, __phi);
8341debfc3dSmrg }
8351debfc3dSmrg
8361debfc3dSmrg // Exponential integrals
8371debfc3dSmrg
8381debfc3dSmrg /**
8391debfc3dSmrg * Return the exponential integral @f$ Ei(x) @f$ for @c float argument @c x.
8401debfc3dSmrg *
8411debfc3dSmrg * @see expint for details.
8421debfc3dSmrg */
8431debfc3dSmrg inline float
8441debfc3dSmrg expintf(float __x)
8451debfc3dSmrg { return __detail::__expint<float>(__x); }
8461debfc3dSmrg
8471debfc3dSmrg /**
8481debfc3dSmrg * Return the exponential integral @f$ Ei(x) @f$
8491debfc3dSmrg * for <tt>long double</tt> argument @c x.
8501debfc3dSmrg *
8511debfc3dSmrg * @see expint for details.
8521debfc3dSmrg */
8531debfc3dSmrg inline long double
8541debfc3dSmrg expintl(long double __x)
8551debfc3dSmrg { return __detail::__expint<long double>(__x); }
8561debfc3dSmrg
8571debfc3dSmrg /**
8581debfc3dSmrg * Return the exponential integral @f$ Ei(x) @f$ for @c real argument @c x.
8591debfc3dSmrg *
8601debfc3dSmrg * The exponential integral is given by
8611debfc3dSmrg * \f[
8621debfc3dSmrg * Ei(x) = -\int_{-x}^\infty \frac{e^t}{t} dt
8631debfc3dSmrg * \f]
8641debfc3dSmrg *
8651debfc3dSmrg * @tparam _Tp The floating-point type of the argument @c __x.
8661debfc3dSmrg * @param __x The argument of the exponential integral function.
8671debfc3dSmrg */
8681debfc3dSmrg template<typename _Tp>
8691debfc3dSmrg inline typename __gnu_cxx::__promote<_Tp>::__type
8701debfc3dSmrg expint(_Tp __x)
8711debfc3dSmrg {
8721debfc3dSmrg typedef typename __gnu_cxx::__promote<_Tp>::__type __type;
8731debfc3dSmrg return __detail::__expint<__type>(__x);
8741debfc3dSmrg }
8751debfc3dSmrg
8761debfc3dSmrg // Hermite polynomials
8771debfc3dSmrg
8781debfc3dSmrg /**
8791debfc3dSmrg * Return the Hermite polynomial @f$ H_n(x) @f$ of nonnegative order n
8801debfc3dSmrg * and float argument @c x.
8811debfc3dSmrg *
8821debfc3dSmrg * @see hermite for details.
8831debfc3dSmrg */
8841debfc3dSmrg inline float
8851debfc3dSmrg hermitef(unsigned int __n, float __x)
8861debfc3dSmrg { return __detail::__poly_hermite<float>(__n, __x); }
8871debfc3dSmrg
8881debfc3dSmrg /**
8891debfc3dSmrg * Return the Hermite polynomial @f$ H_n(x) @f$ of nonnegative order n
8901debfc3dSmrg * and <tt>long double</tt> argument @c x.
8911debfc3dSmrg *
8921debfc3dSmrg * @see hermite for details.
8931debfc3dSmrg */
8941debfc3dSmrg inline long double
8951debfc3dSmrg hermitel(unsigned int __n, long double __x)
8961debfc3dSmrg { return __detail::__poly_hermite<long double>(__n, __x); }
8971debfc3dSmrg
8981debfc3dSmrg /**
8991debfc3dSmrg * Return the Hermite polynomial @f$ H_n(x) @f$ of order n
9001debfc3dSmrg * and @c real argument @c x.
9011debfc3dSmrg *
9021debfc3dSmrg * The Hermite polynomial is defined by:
9031debfc3dSmrg * @f[
9041debfc3dSmrg * H_n(x) = (-1)^n e^{x^2} \frac{d^n}{dx^n} e^{-x^2}
9051debfc3dSmrg * @f]
9061debfc3dSmrg *
9071debfc3dSmrg * The Hermite polynomial obeys a reflection formula:
9081debfc3dSmrg * @f[
9091debfc3dSmrg * H_n(-x) = (-1)^n H_n(x)
9101debfc3dSmrg * @f]
9111debfc3dSmrg *
9121debfc3dSmrg * @tparam _Tp The floating-point type of the argument @c __x.
9131debfc3dSmrg * @param __n The order
9141debfc3dSmrg * @param __x The argument
9151debfc3dSmrg */
9161debfc3dSmrg template<typename _Tp>
9171debfc3dSmrg inline typename __gnu_cxx::__promote<_Tp>::__type
9181debfc3dSmrg hermite(unsigned int __n, _Tp __x)
9191debfc3dSmrg {
9201debfc3dSmrg typedef typename __gnu_cxx::__promote<_Tp>::__type __type;
9211debfc3dSmrg return __detail::__poly_hermite<__type>(__n, __x);
9221debfc3dSmrg }
9231debfc3dSmrg
9241debfc3dSmrg // Laguerre polynomials
9251debfc3dSmrg
9261debfc3dSmrg /**
9271debfc3dSmrg * Returns the Laguerre polynomial @f$ L_n(x) @f$ of nonnegative degree @c n
9281debfc3dSmrg * and @c float argument @f$ x >= 0 @f$.
9291debfc3dSmrg *
9301debfc3dSmrg * @see laguerre for more details.
9311debfc3dSmrg */
9321debfc3dSmrg inline float
9331debfc3dSmrg laguerref(unsigned int __n, float __x)
9341debfc3dSmrg { return __detail::__laguerre<float>(__n, __x); }
9351debfc3dSmrg
9361debfc3dSmrg /**
9371debfc3dSmrg * Returns the Laguerre polynomial @f$ L_n(x) @f$ of nonnegative degree @c n
9381debfc3dSmrg * and <tt>long double</tt> argument @f$ x >= 0 @f$.
9391debfc3dSmrg *
9401debfc3dSmrg * @see laguerre for more details.
9411debfc3dSmrg */
9421debfc3dSmrg inline long double
9431debfc3dSmrg laguerrel(unsigned int __n, long double __x)
9441debfc3dSmrg { return __detail::__laguerre<long double>(__n, __x); }
9451debfc3dSmrg
9461debfc3dSmrg /**
9471debfc3dSmrg * Returns the Laguerre polynomial @f$ L_n(x) @f$
9481debfc3dSmrg * of nonnegative degree @c n and real argument @f$ x >= 0 @f$.
9491debfc3dSmrg *
9501debfc3dSmrg * The Laguerre polynomial is defined by:
9511debfc3dSmrg * @f[
9521debfc3dSmrg * L_n(x) = \frac{e^x}{n!} \frac{d^n}{dx^n} (x^ne^{-x})
9531debfc3dSmrg * @f]
9541debfc3dSmrg *
9551debfc3dSmrg * @tparam _Tp The floating-point type of the argument @c __x.
9561debfc3dSmrg * @param __n The nonnegative order
9571debfc3dSmrg * @param __x The argument <tt> __x >= 0 </tt>
9581debfc3dSmrg * @throw std::domain_error if <tt> __x < 0 </tt>.
9591debfc3dSmrg */
9601debfc3dSmrg template<typename _Tp>
9611debfc3dSmrg inline typename __gnu_cxx::__promote<_Tp>::__type
9621debfc3dSmrg laguerre(unsigned int __n, _Tp __x)
9631debfc3dSmrg {
9641debfc3dSmrg typedef typename __gnu_cxx::__promote<_Tp>::__type __type;
9651debfc3dSmrg return __detail::__laguerre<__type>(__n, __x);
9661debfc3dSmrg }
9671debfc3dSmrg
9681debfc3dSmrg // Legendre polynomials
9691debfc3dSmrg
9701debfc3dSmrg /**
9711debfc3dSmrg * Return the Legendre polynomial @f$ P_l(x) @f$ of nonnegative
9721debfc3dSmrg * degree @f$ l @f$ and @c float argument @f$ |x| <= 0 @f$.
9731debfc3dSmrg *
9741debfc3dSmrg * @see legendre for more details.
9751debfc3dSmrg */
9761debfc3dSmrg inline float
9771debfc3dSmrg legendref(unsigned int __l, float __x)
9781debfc3dSmrg { return __detail::__poly_legendre_p<float>(__l, __x); }
9791debfc3dSmrg
9801debfc3dSmrg /**
9811debfc3dSmrg * Return the Legendre polynomial @f$ P_l(x) @f$ of nonnegative
9821debfc3dSmrg * degree @f$ l @f$ and <tt>long double</tt> argument @f$ |x| <= 0 @f$.
9831debfc3dSmrg *
9841debfc3dSmrg * @see legendre for more details.
9851debfc3dSmrg */
9861debfc3dSmrg inline long double
9871debfc3dSmrg legendrel(unsigned int __l, long double __x)
9881debfc3dSmrg { return __detail::__poly_legendre_p<long double>(__l, __x); }
9891debfc3dSmrg
9901debfc3dSmrg /**
9911debfc3dSmrg * Return the Legendre polynomial @f$ P_l(x) @f$ of nonnegative
9921debfc3dSmrg * degree @f$ l @f$ and real argument @f$ |x| <= 0 @f$.
9931debfc3dSmrg *
9941debfc3dSmrg * The Legendre function of order @f$ l @f$ and argument @f$ x @f$,
9951debfc3dSmrg * @f$ P_l(x) @f$, is defined by:
9961debfc3dSmrg * @f[
9971debfc3dSmrg * P_l(x) = \frac{1}{2^l l!}\frac{d^l}{dx^l}(x^2 - 1)^{l}
9981debfc3dSmrg * @f]
9991debfc3dSmrg *
10001debfc3dSmrg * @tparam _Tp The floating-point type of the argument @c __x.
10011debfc3dSmrg * @param __l The degree @f$ l >= 0 @f$
10021debfc3dSmrg * @param __x The argument @c abs(__x) <= 1
10031debfc3dSmrg * @throw std::domain_error if @c abs(__x) > 1
10041debfc3dSmrg */
10051debfc3dSmrg template<typename _Tp>
10061debfc3dSmrg inline typename __gnu_cxx::__promote<_Tp>::__type
10071debfc3dSmrg legendre(unsigned int __l, _Tp __x)
10081debfc3dSmrg {
10091debfc3dSmrg typedef typename __gnu_cxx::__promote<_Tp>::__type __type;
10101debfc3dSmrg return __detail::__poly_legendre_p<__type>(__l, __x);
10111debfc3dSmrg }
10121debfc3dSmrg
10131debfc3dSmrg // Riemann zeta functions
10141debfc3dSmrg
10151debfc3dSmrg /**
10161debfc3dSmrg * Return the Riemann zeta function @f$ \zeta(s) @f$
10171debfc3dSmrg * for @c float argument @f$ s @f$.
10181debfc3dSmrg *
10191debfc3dSmrg * @see riemann_zeta for more details.
10201debfc3dSmrg */
10211debfc3dSmrg inline float
10221debfc3dSmrg riemann_zetaf(float __s)
10231debfc3dSmrg { return __detail::__riemann_zeta<float>(__s); }
10241debfc3dSmrg
10251debfc3dSmrg /**
10261debfc3dSmrg * Return the Riemann zeta function @f$ \zeta(s) @f$
10271debfc3dSmrg * for <tt>long double</tt> argument @f$ s @f$.
10281debfc3dSmrg *
10291debfc3dSmrg * @see riemann_zeta for more details.
10301debfc3dSmrg */
10311debfc3dSmrg inline long double
10321debfc3dSmrg riemann_zetal(long double __s)
10331debfc3dSmrg { return __detail::__riemann_zeta<long double>(__s); }
10341debfc3dSmrg
10351debfc3dSmrg /**
10361debfc3dSmrg * Return the Riemann zeta function @f$ \zeta(s) @f$
10371debfc3dSmrg * for real argument @f$ s @f$.
10381debfc3dSmrg *
10391debfc3dSmrg * The Riemann zeta function is defined by:
10401debfc3dSmrg * @f[
10411debfc3dSmrg * \zeta(s) = \sum_{k=1}^{\infty} k^{-s} \hbox{ for } s > 1
10421debfc3dSmrg * @f]
10431debfc3dSmrg * and
10441debfc3dSmrg * @f[
10451debfc3dSmrg * \zeta(s) = \frac{1}{1-2^{1-s}}\sum_{k=1}^{\infty}(-1)^{k-1}k^{-s}
10461debfc3dSmrg * \hbox{ for } 0 <= s <= 1
10471debfc3dSmrg * @f]
10481debfc3dSmrg * For s < 1 use the reflection formula:
10491debfc3dSmrg * @f[
10501debfc3dSmrg * \zeta(s) = 2^s \pi^{s-1} \sin(\frac{\pi s}{2}) \Gamma(1-s) \zeta(1-s)
10511debfc3dSmrg * @f]
10521debfc3dSmrg *
10531debfc3dSmrg * @tparam _Tp The floating-point type of the argument @c __s.
10541debfc3dSmrg * @param __s The argument <tt> s != 1 </tt>
10551debfc3dSmrg */
10561debfc3dSmrg template<typename _Tp>
10571debfc3dSmrg inline typename __gnu_cxx::__promote<_Tp>::__type
10581debfc3dSmrg riemann_zeta(_Tp __s)
10591debfc3dSmrg {
10601debfc3dSmrg typedef typename __gnu_cxx::__promote<_Tp>::__type __type;
10611debfc3dSmrg return __detail::__riemann_zeta<__type>(__s);
10621debfc3dSmrg }
10631debfc3dSmrg
10641debfc3dSmrg // Spherical Bessel functions
10651debfc3dSmrg
10661debfc3dSmrg /**
10671debfc3dSmrg * Return the spherical Bessel function @f$ j_n(x) @f$ of nonnegative order n
10681debfc3dSmrg * and @c float argument @f$ x >= 0 @f$.
10691debfc3dSmrg *
10701debfc3dSmrg * @see sph_bessel for more details.
10711debfc3dSmrg */
10721debfc3dSmrg inline float
10731debfc3dSmrg sph_besself(unsigned int __n, float __x)
10741debfc3dSmrg { return __detail::__sph_bessel<float>(__n, __x); }
10751debfc3dSmrg
10761debfc3dSmrg /**
10771debfc3dSmrg * Return the spherical Bessel function @f$ j_n(x) @f$ of nonnegative order n
10781debfc3dSmrg * and <tt>long double</tt> argument @f$ x >= 0 @f$.
10791debfc3dSmrg *
10801debfc3dSmrg * @see sph_bessel for more details.
10811debfc3dSmrg */
10821debfc3dSmrg inline long double
10831debfc3dSmrg sph_bessell(unsigned int __n, long double __x)
10841debfc3dSmrg { return __detail::__sph_bessel<long double>(__n, __x); }
10851debfc3dSmrg
10861debfc3dSmrg /**
10871debfc3dSmrg * Return the spherical Bessel function @f$ j_n(x) @f$ of nonnegative order n
10881debfc3dSmrg * and real argument @f$ x >= 0 @f$.
10891debfc3dSmrg *
10901debfc3dSmrg * The spherical Bessel function is defined by:
10911debfc3dSmrg * @f[
10921debfc3dSmrg * j_n(x) = \left(\frac{\pi}{2x} \right) ^{1/2} J_{n+1/2}(x)
10931debfc3dSmrg * @f]
10941debfc3dSmrg *
10951debfc3dSmrg * @tparam _Tp The floating-point type of the argument @c __x.
10961debfc3dSmrg * @param __n The integral order <tt> n >= 0 </tt>
10971debfc3dSmrg * @param __x The real argument <tt> x >= 0 </tt>
10981debfc3dSmrg * @throw std::domain_error if <tt> __x < 0 </tt>.
10991debfc3dSmrg */
11001debfc3dSmrg template<typename _Tp>
11011debfc3dSmrg inline typename __gnu_cxx::__promote<_Tp>::__type
11021debfc3dSmrg sph_bessel(unsigned int __n, _Tp __x)
11031debfc3dSmrg {
11041debfc3dSmrg typedef typename __gnu_cxx::__promote<_Tp>::__type __type;
11051debfc3dSmrg return __detail::__sph_bessel<__type>(__n, __x);
11061debfc3dSmrg }
11071debfc3dSmrg
11081debfc3dSmrg // Spherical associated Legendre functions
11091debfc3dSmrg
11101debfc3dSmrg /**
11111debfc3dSmrg * Return the spherical Legendre function of nonnegative integral
11121debfc3dSmrg * degree @c l and order @c m and float angle @f$ \theta @f$ in radians.
11131debfc3dSmrg *
11141debfc3dSmrg * @see sph_legendre for details.
11151debfc3dSmrg */
11161debfc3dSmrg inline float
11171debfc3dSmrg sph_legendref(unsigned int __l, unsigned int __m, float __theta)
11181debfc3dSmrg { return __detail::__sph_legendre<float>(__l, __m, __theta); }
11191debfc3dSmrg
11201debfc3dSmrg /**
11211debfc3dSmrg * Return the spherical Legendre function of nonnegative integral
11221debfc3dSmrg * degree @c l and order @c m and <tt>long double</tt> angle @f$ \theta @f$
11231debfc3dSmrg * in radians.
11241debfc3dSmrg *
11251debfc3dSmrg * @see sph_legendre for details.
11261debfc3dSmrg */
11271debfc3dSmrg inline long double
11281debfc3dSmrg sph_legendrel(unsigned int __l, unsigned int __m, long double __theta)
11291debfc3dSmrg { return __detail::__sph_legendre<long double>(__l, __m, __theta); }
11301debfc3dSmrg
11311debfc3dSmrg /**
11321debfc3dSmrg * Return the spherical Legendre function of nonnegative integral
11331debfc3dSmrg * degree @c l and order @c m and real angle @f$ \theta @f$ in radians.
11341debfc3dSmrg *
11351debfc3dSmrg * The spherical Legendre function is defined by
11361debfc3dSmrg * @f[
11371debfc3dSmrg * Y_l^m(\theta,\phi) = (-1)^m[\frac{(2l+1)}{4\pi}
11381debfc3dSmrg * \frac{(l-m)!}{(l+m)!}]
11391debfc3dSmrg * P_l^m(\cos\theta) \exp^{im\phi}
11401debfc3dSmrg * @f]
11411debfc3dSmrg *
11421debfc3dSmrg * @tparam _Tp The floating-point type of the angle @c __theta.
11431debfc3dSmrg * @param __l The order <tt> __l >= 0 </tt>
11441debfc3dSmrg * @param __m The degree <tt> __m >= 0 </tt> and <tt> __m <= __l </tt>
11451debfc3dSmrg * @param __theta The radian polar angle argument
11461debfc3dSmrg */
11471debfc3dSmrg template<typename _Tp>
11481debfc3dSmrg inline typename __gnu_cxx::__promote<_Tp>::__type
11491debfc3dSmrg sph_legendre(unsigned int __l, unsigned int __m, _Tp __theta)
11501debfc3dSmrg {
11511debfc3dSmrg typedef typename __gnu_cxx::__promote<_Tp>::__type __type;
11521debfc3dSmrg return __detail::__sph_legendre<__type>(__l, __m, __theta);
11531debfc3dSmrg }
11541debfc3dSmrg
11551debfc3dSmrg // Spherical Neumann functions
11561debfc3dSmrg
11571debfc3dSmrg /**
11581debfc3dSmrg * Return the spherical Neumann function of integral order @f$ n >= 0 @f$
11591debfc3dSmrg * and @c float argument @f$ x >= 0 @f$.
11601debfc3dSmrg *
11611debfc3dSmrg * @see sph_neumann for details.
11621debfc3dSmrg */
11631debfc3dSmrg inline float
11641debfc3dSmrg sph_neumannf(unsigned int __n, float __x)
11651debfc3dSmrg { return __detail::__sph_neumann<float>(__n, __x); }
11661debfc3dSmrg
11671debfc3dSmrg /**
11681debfc3dSmrg * Return the spherical Neumann function of integral order @f$ n >= 0 @f$
11691debfc3dSmrg * and <tt>long double</tt> @f$ x >= 0 @f$.
11701debfc3dSmrg *
11711debfc3dSmrg * @see sph_neumann for details.
11721debfc3dSmrg */
11731debfc3dSmrg inline long double
11741debfc3dSmrg sph_neumannl(unsigned int __n, long double __x)
11751debfc3dSmrg { return __detail::__sph_neumann<long double>(__n, __x); }
11761debfc3dSmrg
11771debfc3dSmrg /**
11781debfc3dSmrg * Return the spherical Neumann function of integral order @f$ n >= 0 @f$
11791debfc3dSmrg * and real argument @f$ x >= 0 @f$.
11801debfc3dSmrg *
11811debfc3dSmrg * The spherical Neumann function is defined by
11821debfc3dSmrg * @f[
11831debfc3dSmrg * n_n(x) = \left(\frac{\pi}{2x} \right) ^{1/2} N_{n+1/2}(x)
11841debfc3dSmrg * @f]
11851debfc3dSmrg *
11861debfc3dSmrg * @tparam _Tp The floating-point type of the argument @c __x.
11871debfc3dSmrg * @param __n The integral order <tt> n >= 0 </tt>
11881debfc3dSmrg * @param __x The real argument <tt> __x >= 0 </tt>
11891debfc3dSmrg * @throw std::domain_error if <tt> __x < 0 </tt>.
11901debfc3dSmrg */
11911debfc3dSmrg template<typename _Tp>
11921debfc3dSmrg inline typename __gnu_cxx::__promote<_Tp>::__type
11931debfc3dSmrg sph_neumann(unsigned int __n, _Tp __x)
11941debfc3dSmrg {
11951debfc3dSmrg typedef typename __gnu_cxx::__promote<_Tp>::__type __type;
11961debfc3dSmrg return __detail::__sph_neumann<__type>(__n, __x);
11971debfc3dSmrg }
11981debfc3dSmrg
1199*8feb0f0bSmrg /// @} group mathsf
12001debfc3dSmrg
12011debfc3dSmrg _GLIBCXX_END_NAMESPACE_VERSION
12021debfc3dSmrg } // namespace std
12031debfc3dSmrg
1204a2dc1f3fSmrg #ifndef __STRICT_ANSI__
_GLIBCXX_VISIBILITY(default)12051debfc3dSmrg namespace __gnu_cxx _GLIBCXX_VISIBILITY(default)
12061debfc3dSmrg {
1207a2dc1f3fSmrg _GLIBCXX_BEGIN_NAMESPACE_VERSION
1208a2dc1f3fSmrg
1209*8feb0f0bSmrg /** @addtogroup mathsf
1210*8feb0f0bSmrg * @{
1211*8feb0f0bSmrg */
1212*8feb0f0bSmrg
1213a2dc1f3fSmrg // Airy functions
1214a2dc1f3fSmrg
1215a2dc1f3fSmrg /**
1216a2dc1f3fSmrg * Return the Airy function @f$ Ai(x) @f$ of @c float argument x.
1217a2dc1f3fSmrg */
1218a2dc1f3fSmrg inline float
1219a2dc1f3fSmrg airy_aif(float __x)
1220a2dc1f3fSmrg {
1221a2dc1f3fSmrg float __Ai, __Bi, __Aip, __Bip;
1222a2dc1f3fSmrg std::__detail::__airy<float>(__x, __Ai, __Bi, __Aip, __Bip);
1223a2dc1f3fSmrg return __Ai;
1224a2dc1f3fSmrg }
1225a2dc1f3fSmrg
1226a2dc1f3fSmrg /**
1227a2dc1f3fSmrg * Return the Airy function @f$ Ai(x) @f$ of <tt>long double</tt> argument x.
1228a2dc1f3fSmrg */
1229a2dc1f3fSmrg inline long double
1230a2dc1f3fSmrg airy_ail(long double __x)
1231a2dc1f3fSmrg {
1232a2dc1f3fSmrg long double __Ai, __Bi, __Aip, __Bip;
1233a2dc1f3fSmrg std::__detail::__airy<long double>(__x, __Ai, __Bi, __Aip, __Bip);
1234a2dc1f3fSmrg return __Ai;
1235a2dc1f3fSmrg }
1236a2dc1f3fSmrg
1237a2dc1f3fSmrg /**
1238a2dc1f3fSmrg * Return the Airy function @f$ Ai(x) @f$ of real argument x.
1239a2dc1f3fSmrg */
1240a2dc1f3fSmrg template<typename _Tp>
1241a2dc1f3fSmrg inline typename __gnu_cxx::__promote<_Tp>::__type
1242a2dc1f3fSmrg airy_ai(_Tp __x)
1243a2dc1f3fSmrg {
1244a2dc1f3fSmrg typedef typename __gnu_cxx::__promote<_Tp>::__type __type;
1245a2dc1f3fSmrg __type __Ai, __Bi, __Aip, __Bip;
1246a2dc1f3fSmrg std::__detail::__airy<__type>(__x, __Ai, __Bi, __Aip, __Bip);
1247a2dc1f3fSmrg return __Ai;
1248a2dc1f3fSmrg }
1249a2dc1f3fSmrg
1250a2dc1f3fSmrg /**
1251a2dc1f3fSmrg * Return the Airy function @f$ Bi(x) @f$ of @c float argument x.
1252a2dc1f3fSmrg */
1253a2dc1f3fSmrg inline float
1254a2dc1f3fSmrg airy_bif(float __x)
1255a2dc1f3fSmrg {
1256a2dc1f3fSmrg float __Ai, __Bi, __Aip, __Bip;
1257a2dc1f3fSmrg std::__detail::__airy<float>(__x, __Ai, __Bi, __Aip, __Bip);
1258a2dc1f3fSmrg return __Bi;
1259a2dc1f3fSmrg }
1260a2dc1f3fSmrg
1261a2dc1f3fSmrg /**
1262a2dc1f3fSmrg * Return the Airy function @f$ Bi(x) @f$ of <tt>long double</tt> argument x.
1263a2dc1f3fSmrg */
1264a2dc1f3fSmrg inline long double
1265a2dc1f3fSmrg airy_bil(long double __x)
1266a2dc1f3fSmrg {
1267a2dc1f3fSmrg long double __Ai, __Bi, __Aip, __Bip;
1268a2dc1f3fSmrg std::__detail::__airy<long double>(__x, __Ai, __Bi, __Aip, __Bip);
1269a2dc1f3fSmrg return __Bi;
1270a2dc1f3fSmrg }
1271a2dc1f3fSmrg
1272a2dc1f3fSmrg /**
1273a2dc1f3fSmrg * Return the Airy function @f$ Bi(x) @f$ of real argument x.
1274a2dc1f3fSmrg */
1275a2dc1f3fSmrg template<typename _Tp>
1276a2dc1f3fSmrg inline typename __gnu_cxx::__promote<_Tp>::__type
1277a2dc1f3fSmrg airy_bi(_Tp __x)
1278a2dc1f3fSmrg {
1279a2dc1f3fSmrg typedef typename __gnu_cxx::__promote<_Tp>::__type __type;
1280a2dc1f3fSmrg __type __Ai, __Bi, __Aip, __Bip;
1281a2dc1f3fSmrg std::__detail::__airy<__type>(__x, __Ai, __Bi, __Aip, __Bip);
1282a2dc1f3fSmrg return __Bi;
1283a2dc1f3fSmrg }
12841debfc3dSmrg
12851debfc3dSmrg // Confluent hypergeometric functions
12861debfc3dSmrg
12871debfc3dSmrg /**
12881debfc3dSmrg * Return the confluent hypergeometric function @f$ {}_1F_1(a;c;x) @f$
12891debfc3dSmrg * of @c float numeratorial parameter @c a, denominatorial parameter @c c,
12901debfc3dSmrg * and argument @c x.
12911debfc3dSmrg *
12921debfc3dSmrg * @see conf_hyperg for details.
12931debfc3dSmrg */
12941debfc3dSmrg inline float
12951debfc3dSmrg conf_hypergf(float __a, float __c, float __x)
12961debfc3dSmrg { return std::__detail::__conf_hyperg<float>(__a, __c, __x); }
12971debfc3dSmrg
12981debfc3dSmrg /**
12991debfc3dSmrg * Return the confluent hypergeometric function @f$ {}_1F_1(a;c;x) @f$
13001debfc3dSmrg * of <tt>long double</tt> numeratorial parameter @c a,
13011debfc3dSmrg * denominatorial parameter @c c, and argument @c x.
13021debfc3dSmrg *
13031debfc3dSmrg * @see conf_hyperg for details.
13041debfc3dSmrg */
13051debfc3dSmrg inline long double
13061debfc3dSmrg conf_hypergl(long double __a, long double __c, long double __x)
13071debfc3dSmrg { return std::__detail::__conf_hyperg<long double>(__a, __c, __x); }
13081debfc3dSmrg
13091debfc3dSmrg /**
13101debfc3dSmrg * Return the confluent hypergeometric function @f$ {}_1F_1(a;c;x) @f$
13111debfc3dSmrg * of real numeratorial parameter @c a, denominatorial parameter @c c,
13121debfc3dSmrg * and argument @c x.
13131debfc3dSmrg *
13141debfc3dSmrg * The confluent hypergeometric function is defined by
13151debfc3dSmrg * @f[
13161debfc3dSmrg * {}_1F_1(a;c;x) = \sum_{n=0}^{\infty} \frac{(a)_n x^n}{(c)_n n!}
13171debfc3dSmrg * @f]
13181debfc3dSmrg * where the Pochhammer symbol is @f$ (x)_k = (x)(x+1)...(x+k-1) @f$,
13191debfc3dSmrg * @f$ (x)_0 = 1 @f$
13201debfc3dSmrg *
13211debfc3dSmrg * @param __a The numeratorial parameter
13221debfc3dSmrg * @param __c The denominatorial parameter
13231debfc3dSmrg * @param __x The argument
13241debfc3dSmrg */
13251debfc3dSmrg template<typename _Tpa, typename _Tpc, typename _Tp>
13261debfc3dSmrg inline typename __gnu_cxx::__promote_3<_Tpa, _Tpc, _Tp>::__type
13271debfc3dSmrg conf_hyperg(_Tpa __a, _Tpc __c, _Tp __x)
13281debfc3dSmrg {
13291debfc3dSmrg typedef typename __gnu_cxx::__promote_3<_Tpa, _Tpc, _Tp>::__type __type;
13301debfc3dSmrg return std::__detail::__conf_hyperg<__type>(__a, __c, __x);
13311debfc3dSmrg }
13321debfc3dSmrg
13331debfc3dSmrg // Hypergeometric functions
13341debfc3dSmrg
13351debfc3dSmrg /**
13361debfc3dSmrg * Return the hypergeometric function @f$ {}_2F_1(a,b;c;x) @f$
13371debfc3dSmrg * of @ float numeratorial parameters @c a and @c b,
13381debfc3dSmrg * denominatorial parameter @c c, and argument @c x.
13391debfc3dSmrg *
13401debfc3dSmrg * @see hyperg for details.
13411debfc3dSmrg */
13421debfc3dSmrg inline float
13431debfc3dSmrg hypergf(float __a, float __b, float __c, float __x)
13441debfc3dSmrg { return std::__detail::__hyperg<float>(__a, __b, __c, __x); }
13451debfc3dSmrg
13461debfc3dSmrg /**
13471debfc3dSmrg * Return the hypergeometric function @f$ {}_2F_1(a,b;c;x) @f$
13481debfc3dSmrg * of <tt>long double</tt> numeratorial parameters @c a and @c b,
13491debfc3dSmrg * denominatorial parameter @c c, and argument @c x.
13501debfc3dSmrg *
13511debfc3dSmrg * @see hyperg for details.
13521debfc3dSmrg */
13531debfc3dSmrg inline long double
13541debfc3dSmrg hypergl(long double __a, long double __b, long double __c, long double __x)
13551debfc3dSmrg { return std::__detail::__hyperg<long double>(__a, __b, __c, __x); }
13561debfc3dSmrg
13571debfc3dSmrg /**
13581debfc3dSmrg * Return the hypergeometric function @f$ {}_2F_1(a,b;c;x) @f$
13591debfc3dSmrg * of real numeratorial parameters @c a and @c b,
13601debfc3dSmrg * denominatorial parameter @c c, and argument @c x.
13611debfc3dSmrg *
13621debfc3dSmrg * The hypergeometric function is defined by
13631debfc3dSmrg * @f[
13641debfc3dSmrg * {}_2F_1(a;c;x) = \sum_{n=0}^{\infty} \frac{(a)_n (b)_n x^n}{(c)_n n!}
13651debfc3dSmrg * @f]
13661debfc3dSmrg * where the Pochhammer symbol is @f$ (x)_k = (x)(x+1)...(x+k-1) @f$,
13671debfc3dSmrg * @f$ (x)_0 = 1 @f$
13681debfc3dSmrg *
13691debfc3dSmrg * @param __a The first numeratorial parameter
13701debfc3dSmrg * @param __b The second numeratorial parameter
13711debfc3dSmrg * @param __c The denominatorial parameter
13721debfc3dSmrg * @param __x The argument
13731debfc3dSmrg */
13741debfc3dSmrg template<typename _Tpa, typename _Tpb, typename _Tpc, typename _Tp>
13751debfc3dSmrg inline typename __gnu_cxx::__promote_4<_Tpa, _Tpb, _Tpc, _Tp>::__type
13761debfc3dSmrg hyperg(_Tpa __a, _Tpb __b, _Tpc __c, _Tp __x)
13771debfc3dSmrg {
13781debfc3dSmrg typedef typename __gnu_cxx::__promote_4<_Tpa, _Tpb, _Tpc, _Tp>
13791debfc3dSmrg ::__type __type;
13801debfc3dSmrg return std::__detail::__hyperg<__type>(__a, __b, __c, __x);
13811debfc3dSmrg }
13821debfc3dSmrg
1383*8feb0f0bSmrg /// @}
1384a2dc1f3fSmrg _GLIBCXX_END_NAMESPACE_VERSION
13851debfc3dSmrg } // namespace __gnu_cxx
1386a2dc1f3fSmrg #endif // __STRICT_ANSI__
13871debfc3dSmrg
13881debfc3dSmrg #pragma GCC visibility pop
13891debfc3dSmrg
13901debfc3dSmrg #endif // _GLIBCXX_BITS_SPECFUN_H
1391