xref: /dflybsd-src/contrib/gcc-8.0/libstdc++-v3/include/bits/specfun.h (revision 38fd149817dfbff97799f62fcb70be98c4e32523)
1*38fd1498Szrj // Mathematical Special Functions for -*- C++ -*-
2*38fd1498Szrj 
3*38fd1498Szrj // Copyright (C) 2006-2018 Free Software Foundation, Inc.
4*38fd1498Szrj //
5*38fd1498Szrj // This file is part of the GNU ISO C++ Library.  This library is free
6*38fd1498Szrj // software; you can redistribute it and/or modify it under the
7*38fd1498Szrj // terms of the GNU General Public License as published by the
8*38fd1498Szrj // Free Software Foundation; either version 3, or (at your option)
9*38fd1498Szrj // any later version.
10*38fd1498Szrj 
11*38fd1498Szrj // This library is distributed in the hope that it will be useful,
12*38fd1498Szrj // but WITHOUT ANY WARRANTY; without even the implied warranty of
13*38fd1498Szrj // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
14*38fd1498Szrj // GNU General Public License for more details.
15*38fd1498Szrj 
16*38fd1498Szrj // Under Section 7 of GPL version 3, you are granted additional
17*38fd1498Szrj // permissions described in the GCC Runtime Library Exception, version
18*38fd1498Szrj // 3.1, as published by the Free Software Foundation.
19*38fd1498Szrj 
20*38fd1498Szrj // You should have received a copy of the GNU General Public License and
21*38fd1498Szrj // a copy of the GCC Runtime Library Exception along with this program;
22*38fd1498Szrj // see the files COPYING3 and COPYING.RUNTIME respectively.  If not, see
23*38fd1498Szrj // <http://www.gnu.org/licenses/>.
24*38fd1498Szrj 
25*38fd1498Szrj /** @file bits/specfun.h
26*38fd1498Szrj  *  This is an internal header file, included by other library headers.
27*38fd1498Szrj  *  Do not attempt to use it directly. @headername{cmath}
28*38fd1498Szrj  */
29*38fd1498Szrj 
30*38fd1498Szrj #ifndef _GLIBCXX_BITS_SPECFUN_H
31*38fd1498Szrj #define _GLIBCXX_BITS_SPECFUN_H 1
32*38fd1498Szrj 
33*38fd1498Szrj #pragma GCC visibility push(default)
34*38fd1498Szrj 
35*38fd1498Szrj #include <bits/c++config.h>
36*38fd1498Szrj 
37*38fd1498Szrj #define __STDCPP_MATH_SPEC_FUNCS__ 201003L
38*38fd1498Szrj 
39*38fd1498Szrj #define __cpp_lib_math_special_functions 201603L
40*38fd1498Szrj 
41*38fd1498Szrj #if __cplusplus <= 201403L && __STDCPP_WANT_MATH_SPEC_FUNCS__ == 0
42*38fd1498Szrj # error include <cmath> and define __STDCPP_WANT_MATH_SPEC_FUNCS__
43*38fd1498Szrj #endif
44*38fd1498Szrj 
45*38fd1498Szrj #include <bits/stl_algobase.h>
46*38fd1498Szrj #include <limits>
47*38fd1498Szrj #include <type_traits>
48*38fd1498Szrj 
49*38fd1498Szrj #include <tr1/gamma.tcc>
50*38fd1498Szrj #include <tr1/bessel_function.tcc>
51*38fd1498Szrj #include <tr1/beta_function.tcc>
52*38fd1498Szrj #include <tr1/ell_integral.tcc>
53*38fd1498Szrj #include <tr1/exp_integral.tcc>
54*38fd1498Szrj #include <tr1/hypergeometric.tcc>
55*38fd1498Szrj #include <tr1/legendre_function.tcc>
56*38fd1498Szrj #include <tr1/modified_bessel_func.tcc>
57*38fd1498Szrj #include <tr1/poly_hermite.tcc>
58*38fd1498Szrj #include <tr1/poly_laguerre.tcc>
59*38fd1498Szrj #include <tr1/riemann_zeta.tcc>
60*38fd1498Szrj 
_GLIBCXX_VISIBILITY(default)61*38fd1498Szrj namespace std _GLIBCXX_VISIBILITY(default)
62*38fd1498Szrj {
63*38fd1498Szrj _GLIBCXX_BEGIN_NAMESPACE_VERSION
64*38fd1498Szrj 
65*38fd1498Szrj   /**
66*38fd1498Szrj    * @defgroup mathsf Mathematical Special Functions
67*38fd1498Szrj    * @ingroup numerics
68*38fd1498Szrj    *
69*38fd1498Szrj    * A collection of advanced mathematical special functions,
70*38fd1498Szrj    * defined by ISO/IEC IS 29124.
71*38fd1498Szrj    * @{
72*38fd1498Szrj    */
73*38fd1498Szrj 
74*38fd1498Szrj   /**
75*38fd1498Szrj    * @mainpage Mathematical Special Functions
76*38fd1498Szrj    *
77*38fd1498Szrj    * @section intro Introduction and History
78*38fd1498Szrj    * The first significant library upgrade on the road to C++2011,
79*38fd1498Szrj    * <a href="http://www.open-std.org/JTC1/SC22/WG21/docs/papers/2005/n1836.pdf">
80*38fd1498Szrj    * TR1</a>, included a set of 23 mathematical functions that significantly
81*38fd1498Szrj    * extended the standard transcendental functions inherited from C and declared
82*38fd1498Szrj    * in @<cmath@>.
83*38fd1498Szrj    *
84*38fd1498Szrj    * Although most components from TR1 were eventually adopted for C++11 these
85*38fd1498Szrj    * math functions were left behind out of concern for implementability.
86*38fd1498Szrj    * The math functions were published as a separate international standard
87*38fd1498Szrj    * <a href="http://www.open-std.org/JTC1/SC22/WG21/docs/papers/2010/n3060.pdf">
88*38fd1498Szrj    * IS 29124 - Extensions to the C++ Library to Support Mathematical Special
89*38fd1498Szrj    * Functions</a>.
90*38fd1498Szrj    *
91*38fd1498Szrj    * For C++17 these functions were incorporated into the main standard.
92*38fd1498Szrj    *
93*38fd1498Szrj    * @section contents Contents
94*38fd1498Szrj    * The following functions are implemented in namespace @c std:
95*38fd1498Szrj    * - @ref assoc_laguerre "assoc_laguerre - Associated Laguerre functions"
96*38fd1498Szrj    * - @ref assoc_legendre "assoc_legendre - Associated Legendre functions"
97*38fd1498Szrj    * - @ref beta "beta - Beta functions"
98*38fd1498Szrj    * - @ref comp_ellint_1 "comp_ellint_1 - Complete elliptic functions of the first kind"
99*38fd1498Szrj    * - @ref comp_ellint_2 "comp_ellint_2 - Complete elliptic functions of the second kind"
100*38fd1498Szrj    * - @ref comp_ellint_3 "comp_ellint_3 - Complete elliptic functions of the third kind"
101*38fd1498Szrj    * - @ref cyl_bessel_i "cyl_bessel_i - Regular modified cylindrical Bessel functions"
102*38fd1498Szrj    * - @ref cyl_bessel_j "cyl_bessel_j - Cylindrical Bessel functions of the first kind"
103*38fd1498Szrj    * - @ref cyl_bessel_k "cyl_bessel_k - Irregular modified cylindrical Bessel functions"
104*38fd1498Szrj    * - @ref cyl_neumann "cyl_neumann - Cylindrical Neumann functions or Cylindrical Bessel functions of the second kind"
105*38fd1498Szrj    * - @ref ellint_1 "ellint_1 - Incomplete elliptic functions of the first kind"
106*38fd1498Szrj    * - @ref ellint_2 "ellint_2 - Incomplete elliptic functions of the second kind"
107*38fd1498Szrj    * - @ref ellint_3 "ellint_3 - Incomplete elliptic functions of the third kind"
108*38fd1498Szrj    * - @ref expint "expint - The exponential integral"
109*38fd1498Szrj    * - @ref hermite "hermite - Hermite polynomials"
110*38fd1498Szrj    * - @ref laguerre "laguerre - Laguerre functions"
111*38fd1498Szrj    * - @ref legendre "legendre - Legendre polynomials"
112*38fd1498Szrj    * - @ref riemann_zeta "riemann_zeta - The Riemann zeta function"
113*38fd1498Szrj    * - @ref sph_bessel "sph_bessel - Spherical Bessel functions"
114*38fd1498Szrj    * - @ref sph_legendre "sph_legendre - Spherical Legendre functions"
115*38fd1498Szrj    * - @ref sph_neumann "sph_neumann - Spherical Neumann functions"
116*38fd1498Szrj    *
117*38fd1498Szrj    * The hypergeometric functions were stricken from the TR29124 and C++17
118*38fd1498Szrj    * versions of this math library because of implementation concerns.
119*38fd1498Szrj    * However, since they were in the TR1 version and since they are popular
120*38fd1498Szrj    * we kept them as an extension in namespace @c __gnu_cxx:
121*38fd1498Szrj    * - @ref __gnu_cxx::conf_hyperg "conf_hyperg - Confluent hypergeometric functions"
122*38fd1498Szrj    * - @ref __gnu_cxx::hyperg "hyperg - Hypergeometric functions"
123*38fd1498Szrj    *
124*38fd1498Szrj    * @section general General Features
125*38fd1498Szrj    *
126*38fd1498Szrj    * @subsection promotion Argument Promotion
127*38fd1498Szrj    * The arguments suppled to the non-suffixed functions will be promoted
128*38fd1498Szrj    * according to the following rules:
129*38fd1498Szrj    * 1. If any argument intended to be floating point is given an integral value
130*38fd1498Szrj    * That integral value is promoted to double.
131*38fd1498Szrj    * 2. All floating point arguments are promoted up to the largest floating
132*38fd1498Szrj    *    point precision among them.
133*38fd1498Szrj    *
134*38fd1498Szrj    * @subsection NaN NaN Arguments
135*38fd1498Szrj    * If any of the floating point arguments supplied to these functions is
136*38fd1498Szrj    * invalid or NaN (std::numeric_limits<Tp>::quiet_NaN),
137*38fd1498Szrj    * the value NaN is returned.
138*38fd1498Szrj    *
139*38fd1498Szrj    * @section impl Implementation
140*38fd1498Szrj    *
141*38fd1498Szrj    * We strive to implement the underlying math with type generic algorithms
142*38fd1498Szrj    * to the greatest extent possible.  In practice, the functions are thin
143*38fd1498Szrj    * wrappers that dispatch to function templates. Type dependence is
144*38fd1498Szrj    * controlled with std::numeric_limits and functions thereof.
145*38fd1498Szrj    *
146*38fd1498Szrj    * We don't promote @c float to @c double or @c double to <tt>long double</tt>
147*38fd1498Szrj    * reflexively.  The goal is for @c float functions to operate more quickly,
148*38fd1498Szrj    * at the cost of @c float accuracy and possibly a smaller domain of validity.
149*38fd1498Szrj    * Similaryly, <tt>long double</tt> should give you more dynamic range
150*38fd1498Szrj    * and slightly more pecision than @c double on many systems.
151*38fd1498Szrj    *
152*38fd1498Szrj    * @section testing Testing
153*38fd1498Szrj    *
154*38fd1498Szrj    * These functions have been tested against equivalent implementations
155*38fd1498Szrj    * from the <a href="http://www.gnu.org/software/gsl">
156*38fd1498Szrj    * Gnu Scientific Library, GSL</a> and
157*38fd1498Szrj    * <a href="http://www.boost.org/doc/libs/1_60_0/libs/math/doc/html/index.html>Boost</a>
158*38fd1498Szrj    * and the ratio
159*38fd1498Szrj    * @f[
160*38fd1498Szrj    *   \frac{|f - f_{test}|}{|f_{test}|}
161*38fd1498Szrj    * @f]
162*38fd1498Szrj    * is generally found to be within 10^-15 for 64-bit double on linux-x86_64 systems
163*38fd1498Szrj    * over most of the ranges of validity.
164*38fd1498Szrj    *
165*38fd1498Szrj    * @todo Provide accuracy comparisons on a per-function basis for a small
166*38fd1498Szrj    *       number of targets.
167*38fd1498Szrj    *
168*38fd1498Szrj    * @section bibliography General Bibliography
169*38fd1498Szrj    *
170*38fd1498Szrj    * @see Abramowitz and Stegun: Handbook of Mathematical Functions,
171*38fd1498Szrj    * with Formulas, Graphs, and Mathematical Tables
172*38fd1498Szrj    * Edited by Milton Abramowitz and Irene A. Stegun,
173*38fd1498Szrj    * National Bureau of Standards  Applied Mathematics Series - 55
174*38fd1498Szrj    * Issued June 1964, Tenth Printing, December 1972, with corrections
175*38fd1498Szrj    * Electronic versions of A&S abound including both pdf and navigable html.
176*38fd1498Szrj    * @see for example  http://people.math.sfu.ca/~cbm/aands/
177*38fd1498Szrj    *
178*38fd1498Szrj    * @see The old A&S has been redone as the
179*38fd1498Szrj    * NIST Digital Library of Mathematical Functions: http://dlmf.nist.gov/
180*38fd1498Szrj    * This version is far more navigable and includes more recent work.
181*38fd1498Szrj    *
182*38fd1498Szrj    * @see An Atlas of Functions: with Equator, the Atlas Function Calculator
183*38fd1498Szrj    * 2nd Edition, by Oldham, Keith B., Myland, Jan, Spanier, Jerome
184*38fd1498Szrj    *
185*38fd1498Szrj    * @see Asymptotics and Special Functions by Frank W. J. Olver,
186*38fd1498Szrj    * Academic Press, 1974
187*38fd1498Szrj    *
188*38fd1498Szrj    * @see Numerical Recipes in C, The Art of Scientific Computing,
189*38fd1498Szrj    * by William H. Press, Second Ed., Saul A. Teukolsky,
190*38fd1498Szrj    * William T. Vetterling, and Brian P. Flannery,
191*38fd1498Szrj    * Cambridge University Press, 1992
192*38fd1498Szrj    *
193*38fd1498Szrj    * @see The Special Functions and Their Approximations: Volumes 1 and 2,
194*38fd1498Szrj    * by Yudell L. Luke, Academic Press, 1969
195*38fd1498Szrj    */
196*38fd1498Szrj 
197*38fd1498Szrj   // Associated Laguerre polynomials
198*38fd1498Szrj 
199*38fd1498Szrj   /**
200*38fd1498Szrj    * Return the associated Laguerre polynomial of order @c n,
201*38fd1498Szrj    * degree @c m: @f$ L_n^m(x) @f$ for @c float argument.
202*38fd1498Szrj    *
203*38fd1498Szrj    * @see assoc_laguerre for more details.
204*38fd1498Szrj    */
205*38fd1498Szrj   inline float
206*38fd1498Szrj   assoc_laguerref(unsigned int __n, unsigned int __m, float __x)
207*38fd1498Szrj   { return __detail::__assoc_laguerre<float>(__n, __m, __x); }
208*38fd1498Szrj 
209*38fd1498Szrj   /**
210*38fd1498Szrj    * Return the associated Laguerre polynomial of order @c n,
211*38fd1498Szrj    * degree @c m: @f$ L_n^m(x) @f$.
212*38fd1498Szrj    *
213*38fd1498Szrj    * @see assoc_laguerre for more details.
214*38fd1498Szrj    */
215*38fd1498Szrj   inline long double
216*38fd1498Szrj   assoc_laguerrel(unsigned int __n, unsigned int __m, long double __x)
217*38fd1498Szrj   { return __detail::__assoc_laguerre<long double>(__n, __m, __x); }
218*38fd1498Szrj 
219*38fd1498Szrj   /**
220*38fd1498Szrj    * Return the associated Laguerre polynomial of nonnegative order @c n,
221*38fd1498Szrj    * nonnegative degree @c m and real argument @c x: @f$ L_n^m(x) @f$.
222*38fd1498Szrj    *
223*38fd1498Szrj    * The associated Laguerre function of real degree @f$ \alpha @f$,
224*38fd1498Szrj    * @f$ L_n^\alpha(x) @f$, is defined by
225*38fd1498Szrj    * @f[
226*38fd1498Szrj    * 	 L_n^\alpha(x) = \frac{(\alpha + 1)_n}{n!}
227*38fd1498Szrj    * 			 {}_1F_1(-n; \alpha + 1; x)
228*38fd1498Szrj    * @f]
229*38fd1498Szrj    * where @f$ (\alpha)_n @f$ is the Pochhammer symbol and
230*38fd1498Szrj    * @f$ {}_1F_1(a; c; x) @f$ is the confluent hypergeometric function.
231*38fd1498Szrj    *
232*38fd1498Szrj    * The associated Laguerre polynomial is defined for integral
233*38fd1498Szrj    * degree @f$ \alpha = m @f$ by:
234*38fd1498Szrj    * @f[
235*38fd1498Szrj    * 	 L_n^m(x) = (-1)^m \frac{d^m}{dx^m} L_{n + m}(x)
236*38fd1498Szrj    * @f]
237*38fd1498Szrj    * where the Laguerre polynomial is defined by:
238*38fd1498Szrj    * @f[
239*38fd1498Szrj    * 	 L_n(x) = \frac{e^x}{n!} \frac{d^n}{dx^n} (x^ne^{-x})
240*38fd1498Szrj    * @f]
241*38fd1498Szrj    * and @f$ x >= 0 @f$.
242*38fd1498Szrj    * @see laguerre for details of the Laguerre function of degree @c n
243*38fd1498Szrj    *
244*38fd1498Szrj    * @tparam _Tp The floating-point type of the argument @c __x.
245*38fd1498Szrj    * @param __n The order of the Laguerre function, <tt>__n >= 0</tt>.
246*38fd1498Szrj    * @param __m The degree of the Laguerre function, <tt>__m >= 0</tt>.
247*38fd1498Szrj    * @param __x The argument of the Laguerre function, <tt>__x >= 0</tt>.
248*38fd1498Szrj    * @throw std::domain_error if <tt>__x < 0</tt>.
249*38fd1498Szrj    */
250*38fd1498Szrj   template<typename _Tp>
251*38fd1498Szrj     inline typename __gnu_cxx::__promote<_Tp>::__type
252*38fd1498Szrj     assoc_laguerre(unsigned int __n, unsigned int __m, _Tp __x)
253*38fd1498Szrj     {
254*38fd1498Szrj       typedef typename __gnu_cxx::__promote<_Tp>::__type __type;
255*38fd1498Szrj       return __detail::__assoc_laguerre<__type>(__n, __m, __x);
256*38fd1498Szrj     }
257*38fd1498Szrj 
258*38fd1498Szrj   // Associated Legendre functions
259*38fd1498Szrj 
260*38fd1498Szrj   /**
261*38fd1498Szrj    * Return the associated Legendre function of degree @c l and order @c m
262*38fd1498Szrj    * for @c float argument.
263*38fd1498Szrj    *
264*38fd1498Szrj    * @see assoc_legendre for more details.
265*38fd1498Szrj    */
266*38fd1498Szrj   inline float
267*38fd1498Szrj   assoc_legendref(unsigned int __l, unsigned int __m, float __x)
268*38fd1498Szrj   { return __detail::__assoc_legendre_p<float>(__l, __m, __x); }
269*38fd1498Szrj 
270*38fd1498Szrj   /**
271*38fd1498Szrj    * Return the associated Legendre function of degree @c l and order @c m.
272*38fd1498Szrj    *
273*38fd1498Szrj    * @see assoc_legendre for more details.
274*38fd1498Szrj    */
275*38fd1498Szrj   inline long double
276*38fd1498Szrj   assoc_legendrel(unsigned int __l, unsigned int __m, long double __x)
277*38fd1498Szrj   { return __detail::__assoc_legendre_p<long double>(__l, __m, __x); }
278*38fd1498Szrj 
279*38fd1498Szrj 
280*38fd1498Szrj   /**
281*38fd1498Szrj    * Return the associated Legendre function of degree @c l and order @c m.
282*38fd1498Szrj    *
283*38fd1498Szrj    * The associated Legendre function is derived from the Legendre function
284*38fd1498Szrj    * @f$ P_l(x) @f$ by the Rodrigues formula:
285*38fd1498Szrj    * @f[
286*38fd1498Szrj    *   P_l^m(x) = (1 - x^2)^{m/2}\frac{d^m}{dx^m}P_l(x)
287*38fd1498Szrj    * @f]
288*38fd1498Szrj    * @see legendre for details of the Legendre function of degree @c l
289*38fd1498Szrj    *
290*38fd1498Szrj    * @tparam _Tp The floating-point type of the argument @c __x.
291*38fd1498Szrj    * @param  __l  The degree <tt>__l >= 0</tt>.
292*38fd1498Szrj    * @param  __m  The order <tt>__m <= l</tt>.
293*38fd1498Szrj    * @param  __x  The argument, <tt>abs(__x) <= 1</tt>.
294*38fd1498Szrj    * @throw std::domain_error if <tt>abs(__x) > 1</tt>.
295*38fd1498Szrj    */
296*38fd1498Szrj   template<typename _Tp>
297*38fd1498Szrj     inline typename __gnu_cxx::__promote<_Tp>::__type
298*38fd1498Szrj     assoc_legendre(unsigned int __l, unsigned int __m, _Tp __x)
299*38fd1498Szrj     {
300*38fd1498Szrj       typedef typename __gnu_cxx::__promote<_Tp>::__type __type;
301*38fd1498Szrj       return __detail::__assoc_legendre_p<__type>(__l, __m, __x);
302*38fd1498Szrj     }
303*38fd1498Szrj 
304*38fd1498Szrj   // Beta functions
305*38fd1498Szrj 
306*38fd1498Szrj   /**
307*38fd1498Szrj    * Return the beta function, @f$ B(a,b) @f$, for @c float parameters @c a, @c b.
308*38fd1498Szrj    *
309*38fd1498Szrj    * @see beta for more details.
310*38fd1498Szrj    */
311*38fd1498Szrj   inline float
312*38fd1498Szrj   betaf(float __a, float __b)
313*38fd1498Szrj   { return __detail::__beta<float>(__a, __b); }
314*38fd1498Szrj 
315*38fd1498Szrj   /**
316*38fd1498Szrj    * Return the beta function, @f$B(a,b)@f$, for long double
317*38fd1498Szrj    * parameters @c a, @c b.
318*38fd1498Szrj    *
319*38fd1498Szrj    * @see beta for more details.
320*38fd1498Szrj    */
321*38fd1498Szrj   inline long double
322*38fd1498Szrj   betal(long double __a, long double __b)
323*38fd1498Szrj   { return __detail::__beta<long double>(__a, __b); }
324*38fd1498Szrj 
325*38fd1498Szrj   /**
326*38fd1498Szrj    * Return the beta function, @f$B(a,b)@f$, for real parameters @c a, @c b.
327*38fd1498Szrj    *
328*38fd1498Szrj    * The beta function is defined by
329*38fd1498Szrj    * @f[
330*38fd1498Szrj    *   B(a,b) = \int_0^1 t^{a - 1} (1 - t)^{b - 1} dt
331*38fd1498Szrj    *          = \frac{\Gamma(a)\Gamma(b)}{\Gamma(a+b)}
332*38fd1498Szrj    * @f]
333*38fd1498Szrj    * where @f$ a > 0 @f$ and @f$ b > 0 @f$
334*38fd1498Szrj    *
335*38fd1498Szrj    * @tparam _Tpa The floating-point type of the parameter @c __a.
336*38fd1498Szrj    * @tparam _Tpb The floating-point type of the parameter @c __b.
337*38fd1498Szrj    * @param __a The first argument of the beta function, <tt> __a > 0 </tt>.
338*38fd1498Szrj    * @param __b The second argument of the beta function, <tt> __b > 0 </tt>.
339*38fd1498Szrj    * @throw std::domain_error if <tt> __a < 0 </tt> or <tt> __b < 0 </tt>.
340*38fd1498Szrj    */
341*38fd1498Szrj   template<typename _Tpa, typename _Tpb>
342*38fd1498Szrj     inline typename __gnu_cxx::__promote_2<_Tpa, _Tpb>::__type
343*38fd1498Szrj     beta(_Tpa __a, _Tpb __b)
344*38fd1498Szrj     {
345*38fd1498Szrj       typedef typename __gnu_cxx::__promote_2<_Tpa, _Tpb>::__type __type;
346*38fd1498Szrj       return __detail::__beta<__type>(__a, __b);
347*38fd1498Szrj     }
348*38fd1498Szrj 
349*38fd1498Szrj   // Complete elliptic integrals of the first kind
350*38fd1498Szrj 
351*38fd1498Szrj   /**
352*38fd1498Szrj    * Return the complete elliptic integral of the first kind @f$ E(k) @f$
353*38fd1498Szrj    * for @c float modulus @c k.
354*38fd1498Szrj    *
355*38fd1498Szrj    * @see comp_ellint_1 for details.
356*38fd1498Szrj    */
357*38fd1498Szrj   inline float
358*38fd1498Szrj   comp_ellint_1f(float __k)
359*38fd1498Szrj   { return __detail::__comp_ellint_1<float>(__k); }
360*38fd1498Szrj 
361*38fd1498Szrj   /**
362*38fd1498Szrj    * Return the complete elliptic integral of the first kind @f$ E(k) @f$
363*38fd1498Szrj    * for long double modulus @c k.
364*38fd1498Szrj    *
365*38fd1498Szrj    * @see comp_ellint_1 for details.
366*38fd1498Szrj    */
367*38fd1498Szrj   inline long double
368*38fd1498Szrj   comp_ellint_1l(long double __k)
369*38fd1498Szrj   { return __detail::__comp_ellint_1<long double>(__k); }
370*38fd1498Szrj 
371*38fd1498Szrj   /**
372*38fd1498Szrj    * Return the complete elliptic integral of the first kind
373*38fd1498Szrj    * @f$ K(k) @f$ for real modulus @c k.
374*38fd1498Szrj    *
375*38fd1498Szrj    * The complete elliptic integral of the first kind is defined as
376*38fd1498Szrj    * @f[
377*38fd1498Szrj    *   K(k) = F(k,\pi/2) = \int_0^{\pi/2}\frac{d\theta}
378*38fd1498Szrj    * 					     {\sqrt{1 - k^2 sin^2\theta}}
379*38fd1498Szrj    * @f]
380*38fd1498Szrj    * where @f$ F(k,\phi) @f$ is the incomplete elliptic integral of the
381*38fd1498Szrj    * first kind and the modulus @f$ |k| <= 1 @f$.
382*38fd1498Szrj    * @see ellint_1 for details of the incomplete elliptic function
383*38fd1498Szrj    * of the first kind.
384*38fd1498Szrj    *
385*38fd1498Szrj    * @tparam _Tp The floating-point type of the modulus @c __k.
386*38fd1498Szrj    * @param  __k  The modulus, <tt> abs(__k) <= 1 </tt>
387*38fd1498Szrj    * @throw std::domain_error if <tt> abs(__k) > 1 </tt>.
388*38fd1498Szrj    */
389*38fd1498Szrj   template<typename _Tp>
390*38fd1498Szrj     inline typename __gnu_cxx::__promote<_Tp>::__type
391*38fd1498Szrj     comp_ellint_1(_Tp __k)
392*38fd1498Szrj     {
393*38fd1498Szrj       typedef typename __gnu_cxx::__promote<_Tp>::__type __type;
394*38fd1498Szrj       return __detail::__comp_ellint_1<__type>(__k);
395*38fd1498Szrj     }
396*38fd1498Szrj 
397*38fd1498Szrj   // Complete elliptic integrals of the second kind
398*38fd1498Szrj 
399*38fd1498Szrj   /**
400*38fd1498Szrj    * Return the complete elliptic integral of the second kind @f$ E(k) @f$
401*38fd1498Szrj    * for @c float modulus @c k.
402*38fd1498Szrj    *
403*38fd1498Szrj    * @see comp_ellint_2 for details.
404*38fd1498Szrj    */
405*38fd1498Szrj   inline float
406*38fd1498Szrj   comp_ellint_2f(float __k)
407*38fd1498Szrj   { return __detail::__comp_ellint_2<float>(__k); }
408*38fd1498Szrj 
409*38fd1498Szrj   /**
410*38fd1498Szrj    * Return the complete elliptic integral of the second kind @f$ E(k) @f$
411*38fd1498Szrj    * for long double modulus @c k.
412*38fd1498Szrj    *
413*38fd1498Szrj    * @see comp_ellint_2 for details.
414*38fd1498Szrj    */
415*38fd1498Szrj   inline long double
416*38fd1498Szrj   comp_ellint_2l(long double __k)
417*38fd1498Szrj   { return __detail::__comp_ellint_2<long double>(__k); }
418*38fd1498Szrj 
419*38fd1498Szrj   /**
420*38fd1498Szrj    * Return the complete elliptic integral of the second kind @f$ E(k) @f$
421*38fd1498Szrj    * for real modulus @c k.
422*38fd1498Szrj    *
423*38fd1498Szrj    * The complete elliptic integral of the second kind is defined as
424*38fd1498Szrj    * @f[
425*38fd1498Szrj    *   E(k) = E(k,\pi/2) = \int_0^{\pi/2}\sqrt{1 - k^2 sin^2\theta}
426*38fd1498Szrj    * @f]
427*38fd1498Szrj    * where @f$ E(k,\phi) @f$ is the incomplete elliptic integral of the
428*38fd1498Szrj    * second kind and the modulus @f$ |k| <= 1 @f$.
429*38fd1498Szrj    * @see ellint_2 for details of the incomplete elliptic function
430*38fd1498Szrj    * of the second kind.
431*38fd1498Szrj    *
432*38fd1498Szrj    * @tparam _Tp The floating-point type of the modulus @c __k.
433*38fd1498Szrj    * @param  __k  The modulus, @c abs(__k) <= 1
434*38fd1498Szrj    * @throw std::domain_error if @c abs(__k) > 1.
435*38fd1498Szrj    */
436*38fd1498Szrj   template<typename _Tp>
437*38fd1498Szrj     inline typename __gnu_cxx::__promote<_Tp>::__type
438*38fd1498Szrj     comp_ellint_2(_Tp __k)
439*38fd1498Szrj     {
440*38fd1498Szrj       typedef typename __gnu_cxx::__promote<_Tp>::__type __type;
441*38fd1498Szrj       return __detail::__comp_ellint_2<__type>(__k);
442*38fd1498Szrj     }
443*38fd1498Szrj 
444*38fd1498Szrj   // Complete elliptic integrals of the third kind
445*38fd1498Szrj 
446*38fd1498Szrj   /**
447*38fd1498Szrj    * @brief Return the complete elliptic integral of the third kind
448*38fd1498Szrj    * @f$ \Pi(k,\nu) @f$ for @c float modulus @c k.
449*38fd1498Szrj    *
450*38fd1498Szrj    * @see comp_ellint_3 for details.
451*38fd1498Szrj    */
452*38fd1498Szrj   inline float
453*38fd1498Szrj   comp_ellint_3f(float __k, float __nu)
454*38fd1498Szrj   { return __detail::__comp_ellint_3<float>(__k, __nu); }
455*38fd1498Szrj 
456*38fd1498Szrj   /**
457*38fd1498Szrj    * @brief Return the complete elliptic integral of the third kind
458*38fd1498Szrj    * @f$ \Pi(k,\nu) @f$ for <tt>long double</tt> modulus @c k.
459*38fd1498Szrj    *
460*38fd1498Szrj    * @see comp_ellint_3 for details.
461*38fd1498Szrj    */
462*38fd1498Szrj   inline long double
463*38fd1498Szrj   comp_ellint_3l(long double __k, long double __nu)
464*38fd1498Szrj   { return __detail::__comp_ellint_3<long double>(__k, __nu); }
465*38fd1498Szrj 
466*38fd1498Szrj   /**
467*38fd1498Szrj    * Return the complete elliptic integral of the third kind
468*38fd1498Szrj    * @f$ \Pi(k,\nu) = \Pi(k,\nu,\pi/2) @f$ for real modulus @c k.
469*38fd1498Szrj    *
470*38fd1498Szrj    * The complete elliptic integral of the third kind is defined as
471*38fd1498Szrj    * @f[
472*38fd1498Szrj    *   \Pi(k,\nu) = \Pi(k,\nu,\pi/2) = \int_0^{\pi/2}
473*38fd1498Szrj    * 		     \frac{d\theta}
474*38fd1498Szrj    * 		   {(1 - \nu \sin^2\theta)\sqrt{1 - k^2 \sin^2\theta}}
475*38fd1498Szrj    * @f]
476*38fd1498Szrj    * where @f$ \Pi(k,\nu,\phi) @f$ is the incomplete elliptic integral of the
477*38fd1498Szrj    * second kind and the modulus @f$ |k| <= 1 @f$.
478*38fd1498Szrj    * @see ellint_3 for details of the incomplete elliptic function
479*38fd1498Szrj    * of the third kind.
480*38fd1498Szrj    *
481*38fd1498Szrj    * @tparam _Tp The floating-point type of the modulus @c __k.
482*38fd1498Szrj    * @tparam _Tpn The floating-point type of the argument @c __nu.
483*38fd1498Szrj    * @param  __k  The modulus, @c abs(__k) <= 1
484*38fd1498Szrj    * @param  __nu  The argument
485*38fd1498Szrj    * @throw std::domain_error if @c abs(__k) > 1.
486*38fd1498Szrj    */
487*38fd1498Szrj   template<typename _Tp, typename _Tpn>
488*38fd1498Szrj     inline typename __gnu_cxx::__promote_2<_Tp, _Tpn>::__type
489*38fd1498Szrj     comp_ellint_3(_Tp __k, _Tpn __nu)
490*38fd1498Szrj     {
491*38fd1498Szrj       typedef typename __gnu_cxx::__promote_2<_Tp, _Tpn>::__type __type;
492*38fd1498Szrj       return __detail::__comp_ellint_3<__type>(__k, __nu);
493*38fd1498Szrj     }
494*38fd1498Szrj 
495*38fd1498Szrj   // Regular modified cylindrical Bessel functions
496*38fd1498Szrj 
497*38fd1498Szrj   /**
498*38fd1498Szrj    * Return the regular modified Bessel function @f$ I_{\nu}(x) @f$
499*38fd1498Szrj    * for @c float order @f$ \nu @f$ and argument @f$ x >= 0 @f$.
500*38fd1498Szrj    *
501*38fd1498Szrj    * @see cyl_bessel_i for setails.
502*38fd1498Szrj    */
503*38fd1498Szrj   inline float
504*38fd1498Szrj   cyl_bessel_if(float __nu, float __x)
505*38fd1498Szrj   { return __detail::__cyl_bessel_i<float>(__nu, __x); }
506*38fd1498Szrj 
507*38fd1498Szrj   /**
508*38fd1498Szrj    * Return the regular modified Bessel function @f$ I_{\nu}(x) @f$
509*38fd1498Szrj    * for <tt>long double</tt> order @f$ \nu @f$ and argument @f$ x >= 0 @f$.
510*38fd1498Szrj    *
511*38fd1498Szrj    * @see cyl_bessel_i for setails.
512*38fd1498Szrj    */
513*38fd1498Szrj   inline long double
514*38fd1498Szrj   cyl_bessel_il(long double __nu, long double __x)
515*38fd1498Szrj   { return __detail::__cyl_bessel_i<long double>(__nu, __x); }
516*38fd1498Szrj 
517*38fd1498Szrj   /**
518*38fd1498Szrj    * Return the regular modified Bessel function @f$ I_{\nu}(x) @f$
519*38fd1498Szrj    * for real order @f$ \nu @f$ and argument @f$ x >= 0 @f$.
520*38fd1498Szrj    *
521*38fd1498Szrj    * The regular modified cylindrical Bessel function is:
522*38fd1498Szrj    * @f[
523*38fd1498Szrj    *  I_{\nu}(x) = i^{-\nu}J_\nu(ix) = \sum_{k=0}^{\infty}
524*38fd1498Szrj    * 		\frac{(x/2)^{\nu + 2k}}{k!\Gamma(\nu+k+1)}
525*38fd1498Szrj    * @f]
526*38fd1498Szrj    *
527*38fd1498Szrj    * @tparam _Tpnu The floating-point type of the order @c __nu.
528*38fd1498Szrj    * @tparam _Tp The floating-point type of the argument @c __x.
529*38fd1498Szrj    * @param  __nu  The order
530*38fd1498Szrj    * @param  __x   The argument, <tt> __x >= 0 </tt>
531*38fd1498Szrj    * @throw std::domain_error if <tt> __x < 0 </tt>.
532*38fd1498Szrj    */
533*38fd1498Szrj   template<typename _Tpnu, typename _Tp>
534*38fd1498Szrj     inline typename __gnu_cxx::__promote_2<_Tpnu, _Tp>::__type
535*38fd1498Szrj     cyl_bessel_i(_Tpnu __nu, _Tp __x)
536*38fd1498Szrj     {
537*38fd1498Szrj       typedef typename __gnu_cxx::__promote_2<_Tpnu, _Tp>::__type __type;
538*38fd1498Szrj       return __detail::__cyl_bessel_i<__type>(__nu, __x);
539*38fd1498Szrj     }
540*38fd1498Szrj 
541*38fd1498Szrj   // Cylindrical Bessel functions (of the first kind)
542*38fd1498Szrj 
543*38fd1498Szrj   /**
544*38fd1498Szrj    * Return the Bessel function of the first kind @f$ J_{\nu}(x) @f$
545*38fd1498Szrj    * for @c float order @f$ \nu @f$ and argument @f$ x >= 0 @f$.
546*38fd1498Szrj    *
547*38fd1498Szrj    * @see cyl_bessel_j for setails.
548*38fd1498Szrj    */
549*38fd1498Szrj   inline float
550*38fd1498Szrj   cyl_bessel_jf(float __nu, float __x)
551*38fd1498Szrj   { return __detail::__cyl_bessel_j<float>(__nu, __x); }
552*38fd1498Szrj 
553*38fd1498Szrj   /**
554*38fd1498Szrj    * Return the Bessel function of the first kind @f$ J_{\nu}(x) @f$
555*38fd1498Szrj    * for <tt>long double</tt> order @f$ \nu @f$ and argument @f$ x >= 0 @f$.
556*38fd1498Szrj    *
557*38fd1498Szrj    * @see cyl_bessel_j for setails.
558*38fd1498Szrj    */
559*38fd1498Szrj   inline long double
560*38fd1498Szrj   cyl_bessel_jl(long double __nu, long double __x)
561*38fd1498Szrj   { return __detail::__cyl_bessel_j<long double>(__nu, __x); }
562*38fd1498Szrj 
563*38fd1498Szrj   /**
564*38fd1498Szrj    * Return the Bessel function @f$ J_{\nu}(x) @f$ of real order @f$ \nu @f$
565*38fd1498Szrj    * and argument @f$ x >= 0 @f$.
566*38fd1498Szrj    *
567*38fd1498Szrj    * The cylindrical Bessel function is:
568*38fd1498Szrj    * @f[
569*38fd1498Szrj    *    J_{\nu}(x) = \sum_{k=0}^{\infty}
570*38fd1498Szrj    *              \frac{(-1)^k (x/2)^{\nu + 2k}}{k!\Gamma(\nu+k+1)}
571*38fd1498Szrj    * @f]
572*38fd1498Szrj    *
573*38fd1498Szrj    * @tparam _Tpnu The floating-point type of the order @c __nu.
574*38fd1498Szrj    * @tparam _Tp The floating-point type of the argument @c __x.
575*38fd1498Szrj    * @param  __nu  The order
576*38fd1498Szrj    * @param  __x   The argument, <tt> __x >= 0 </tt>
577*38fd1498Szrj    * @throw std::domain_error if <tt> __x < 0 </tt>.
578*38fd1498Szrj    */
579*38fd1498Szrj   template<typename _Tpnu, typename _Tp>
580*38fd1498Szrj     inline typename __gnu_cxx::__promote_2<_Tpnu, _Tp>::__type
581*38fd1498Szrj     cyl_bessel_j(_Tpnu __nu, _Tp __x)
582*38fd1498Szrj     {
583*38fd1498Szrj       typedef typename __gnu_cxx::__promote_2<_Tpnu, _Tp>::__type __type;
584*38fd1498Szrj       return __detail::__cyl_bessel_j<__type>(__nu, __x);
585*38fd1498Szrj     }
586*38fd1498Szrj 
587*38fd1498Szrj   // Irregular modified cylindrical Bessel functions
588*38fd1498Szrj 
589*38fd1498Szrj   /**
590*38fd1498Szrj    * Return the irregular modified Bessel function @f$ K_{\nu}(x) @f$
591*38fd1498Szrj    * for @c float order @f$ \nu @f$ and argument @f$ x >= 0 @f$.
592*38fd1498Szrj    *
593*38fd1498Szrj    * @see cyl_bessel_k for setails.
594*38fd1498Szrj    */
595*38fd1498Szrj   inline float
596*38fd1498Szrj   cyl_bessel_kf(float __nu, float __x)
597*38fd1498Szrj   { return __detail::__cyl_bessel_k<float>(__nu, __x); }
598*38fd1498Szrj 
599*38fd1498Szrj   /**
600*38fd1498Szrj    * Return the irregular modified Bessel function @f$ K_{\nu}(x) @f$
601*38fd1498Szrj    * for <tt>long double</tt> order @f$ \nu @f$ and argument @f$ x >= 0 @f$.
602*38fd1498Szrj    *
603*38fd1498Szrj    * @see cyl_bessel_k for setails.
604*38fd1498Szrj    */
605*38fd1498Szrj   inline long double
606*38fd1498Szrj   cyl_bessel_kl(long double __nu, long double __x)
607*38fd1498Szrj   { return __detail::__cyl_bessel_k<long double>(__nu, __x); }
608*38fd1498Szrj 
609*38fd1498Szrj   /**
610*38fd1498Szrj    * Return the irregular modified Bessel function @f$ K_{\nu}(x) @f$
611*38fd1498Szrj    * of real order @f$ \nu @f$ and argument @f$ x @f$.
612*38fd1498Szrj    *
613*38fd1498Szrj    * The irregular modified Bessel function is defined by:
614*38fd1498Szrj    * @f[
615*38fd1498Szrj    * 	K_{\nu}(x) = \frac{\pi}{2}
616*38fd1498Szrj    * 		     \frac{I_{-\nu}(x) - I_{\nu}(x)}{\sin \nu\pi}
617*38fd1498Szrj    * @f]
618*38fd1498Szrj    * where for integral @f$ \nu = n @f$ a limit is taken:
619*38fd1498Szrj    * @f$ lim_{\nu \to n} @f$.
620*38fd1498Szrj    * For negative argument we have simply:
621*38fd1498Szrj    * @f[
622*38fd1498Szrj    * 	K_{-\nu}(x) = K_{\nu}(x)
623*38fd1498Szrj    * @f]
624*38fd1498Szrj    *
625*38fd1498Szrj    * @tparam _Tpnu The floating-point type of the order @c __nu.
626*38fd1498Szrj    * @tparam _Tp The floating-point type of the argument @c __x.
627*38fd1498Szrj    * @param  __nu  The order
628*38fd1498Szrj    * @param  __x   The argument, <tt> __x >= 0 </tt>
629*38fd1498Szrj    * @throw std::domain_error if <tt> __x < 0 </tt>.
630*38fd1498Szrj    */
631*38fd1498Szrj   template<typename _Tpnu, typename _Tp>
632*38fd1498Szrj     inline typename __gnu_cxx::__promote_2<_Tpnu, _Tp>::__type
633*38fd1498Szrj     cyl_bessel_k(_Tpnu __nu, _Tp __x)
634*38fd1498Szrj     {
635*38fd1498Szrj       typedef typename __gnu_cxx::__promote_2<_Tpnu, _Tp>::__type __type;
636*38fd1498Szrj       return __detail::__cyl_bessel_k<__type>(__nu, __x);
637*38fd1498Szrj     }
638*38fd1498Szrj 
639*38fd1498Szrj   // Cylindrical Neumann functions
640*38fd1498Szrj 
641*38fd1498Szrj   /**
642*38fd1498Szrj    * Return the Neumann function @f$ N_{\nu}(x) @f$
643*38fd1498Szrj    * of @c float order @f$ \nu @f$ and argument @f$ x @f$.
644*38fd1498Szrj    *
645*38fd1498Szrj    * @see cyl_neumann for setails.
646*38fd1498Szrj    */
647*38fd1498Szrj   inline float
648*38fd1498Szrj   cyl_neumannf(float __nu, float __x)
649*38fd1498Szrj   { return __detail::__cyl_neumann_n<float>(__nu, __x); }
650*38fd1498Szrj 
651*38fd1498Szrj   /**
652*38fd1498Szrj    * Return the Neumann function @f$ N_{\nu}(x) @f$
653*38fd1498Szrj    * of <tt>long double</tt> order @f$ \nu @f$ and argument @f$ x @f$.
654*38fd1498Szrj    *
655*38fd1498Szrj    * @see cyl_neumann for setails.
656*38fd1498Szrj    */
657*38fd1498Szrj   inline long double
658*38fd1498Szrj   cyl_neumannl(long double __nu, long double __x)
659*38fd1498Szrj   { return __detail::__cyl_neumann_n<long double>(__nu, __x); }
660*38fd1498Szrj 
661*38fd1498Szrj   /**
662*38fd1498Szrj    * Return the Neumann function @f$ N_{\nu}(x) @f$
663*38fd1498Szrj    * of real order @f$ \nu @f$ and argument @f$ x >= 0 @f$.
664*38fd1498Szrj    *
665*38fd1498Szrj    * The Neumann function is defined by:
666*38fd1498Szrj    * @f[
667*38fd1498Szrj    *    N_{\nu}(x) = \frac{J_{\nu}(x) \cos \nu\pi - J_{-\nu}(x)}
668*38fd1498Szrj    *                      {\sin \nu\pi}
669*38fd1498Szrj    * @f]
670*38fd1498Szrj    * where @f$ x >= 0 @f$ and for integral order @f$ \nu = n @f$
671*38fd1498Szrj    * a limit is taken: @f$ lim_{\nu \to n} @f$.
672*38fd1498Szrj    *
673*38fd1498Szrj    * @tparam _Tpnu The floating-point type of the order @c __nu.
674*38fd1498Szrj    * @tparam _Tp The floating-point type of the argument @c __x.
675*38fd1498Szrj    * @param  __nu  The order
676*38fd1498Szrj    * @param  __x   The argument, <tt> __x >= 0 </tt>
677*38fd1498Szrj    * @throw std::domain_error if <tt> __x < 0 </tt>.
678*38fd1498Szrj    */
679*38fd1498Szrj   template<typename _Tpnu, typename _Tp>
680*38fd1498Szrj     inline typename __gnu_cxx::__promote_2<_Tpnu, _Tp>::__type
681*38fd1498Szrj     cyl_neumann(_Tpnu __nu, _Tp __x)
682*38fd1498Szrj     {
683*38fd1498Szrj       typedef typename __gnu_cxx::__promote_2<_Tpnu, _Tp>::__type __type;
684*38fd1498Szrj       return __detail::__cyl_neumann_n<__type>(__nu, __x);
685*38fd1498Szrj     }
686*38fd1498Szrj 
687*38fd1498Szrj   // Incomplete elliptic integrals of the first kind
688*38fd1498Szrj 
689*38fd1498Szrj   /**
690*38fd1498Szrj    * Return the incomplete elliptic integral of the first kind @f$ E(k,\phi) @f$
691*38fd1498Szrj    * for @c float modulus @f$ k @f$ and angle @f$ \phi @f$.
692*38fd1498Szrj    *
693*38fd1498Szrj    * @see ellint_1 for details.
694*38fd1498Szrj    */
695*38fd1498Szrj   inline float
696*38fd1498Szrj   ellint_1f(float __k, float __phi)
697*38fd1498Szrj   { return __detail::__ellint_1<float>(__k, __phi); }
698*38fd1498Szrj 
699*38fd1498Szrj   /**
700*38fd1498Szrj    * Return the incomplete elliptic integral of the first kind @f$ E(k,\phi) @f$
701*38fd1498Szrj    * for <tt>long double</tt> modulus @f$ k @f$ and angle @f$ \phi @f$.
702*38fd1498Szrj    *
703*38fd1498Szrj    * @see ellint_1 for details.
704*38fd1498Szrj    */
705*38fd1498Szrj   inline long double
706*38fd1498Szrj   ellint_1l(long double __k, long double __phi)
707*38fd1498Szrj   { return __detail::__ellint_1<long double>(__k, __phi); }
708*38fd1498Szrj 
709*38fd1498Szrj   /**
710*38fd1498Szrj    * Return the incomplete elliptic integral of the first kind @f$ F(k,\phi) @f$
711*38fd1498Szrj    * for @c real modulus @f$ k @f$ and angle @f$ \phi @f$.
712*38fd1498Szrj    *
713*38fd1498Szrj    * The incomplete elliptic integral of the first kind is defined as
714*38fd1498Szrj    * @f[
715*38fd1498Szrj    *   F(k,\phi) = \int_0^{\phi}\frac{d\theta}
716*38fd1498Szrj    * 				     {\sqrt{1 - k^2 sin^2\theta}}
717*38fd1498Szrj    * @f]
718*38fd1498Szrj    * For  @f$ \phi= \pi/2 @f$ this becomes the complete elliptic integral of
719*38fd1498Szrj    * the first kind, @f$ K(k) @f$.  @see comp_ellint_1.
720*38fd1498Szrj    *
721*38fd1498Szrj    * @tparam _Tp The floating-point type of the modulus @c __k.
722*38fd1498Szrj    * @tparam _Tpp The floating-point type of the angle @c __phi.
723*38fd1498Szrj    * @param  __k  The modulus, <tt> abs(__k) <= 1 </tt>
724*38fd1498Szrj    * @param  __phi  The integral limit argument in radians
725*38fd1498Szrj    * @throw std::domain_error if <tt> abs(__k) > 1 </tt>.
726*38fd1498Szrj    */
727*38fd1498Szrj   template<typename _Tp, typename _Tpp>
728*38fd1498Szrj     inline typename __gnu_cxx::__promote_2<_Tp, _Tpp>::__type
729*38fd1498Szrj     ellint_1(_Tp __k, _Tpp __phi)
730*38fd1498Szrj     {
731*38fd1498Szrj       typedef typename __gnu_cxx::__promote_2<_Tp, _Tpp>::__type __type;
732*38fd1498Szrj       return __detail::__ellint_1<__type>(__k, __phi);
733*38fd1498Szrj     }
734*38fd1498Szrj 
735*38fd1498Szrj   // Incomplete elliptic integrals of the second kind
736*38fd1498Szrj 
737*38fd1498Szrj   /**
738*38fd1498Szrj    * @brief Return the incomplete elliptic integral of the second kind
739*38fd1498Szrj    * @f$ E(k,\phi) @f$ for @c float argument.
740*38fd1498Szrj    *
741*38fd1498Szrj    * @see ellint_2 for details.
742*38fd1498Szrj    */
743*38fd1498Szrj   inline float
744*38fd1498Szrj   ellint_2f(float __k, float __phi)
745*38fd1498Szrj   { return __detail::__ellint_2<float>(__k, __phi); }
746*38fd1498Szrj 
747*38fd1498Szrj   /**
748*38fd1498Szrj    * @brief Return the incomplete elliptic integral of the second kind
749*38fd1498Szrj    * @f$ E(k,\phi) @f$.
750*38fd1498Szrj    *
751*38fd1498Szrj    * @see ellint_2 for details.
752*38fd1498Szrj    */
753*38fd1498Szrj   inline long double
754*38fd1498Szrj   ellint_2l(long double __k, long double __phi)
755*38fd1498Szrj   { return __detail::__ellint_2<long double>(__k, __phi); }
756*38fd1498Szrj 
757*38fd1498Szrj   /**
758*38fd1498Szrj    * Return the incomplete elliptic integral of the second kind
759*38fd1498Szrj    * @f$ E(k,\phi) @f$.
760*38fd1498Szrj    *
761*38fd1498Szrj    * The incomplete elliptic integral of the second kind is defined as
762*38fd1498Szrj    * @f[
763*38fd1498Szrj    *   E(k,\phi) = \int_0^{\phi} \sqrt{1 - k^2 sin^2\theta}
764*38fd1498Szrj    * @f]
765*38fd1498Szrj    * For  @f$ \phi= \pi/2 @f$ this becomes the complete elliptic integral of
766*38fd1498Szrj    * the second kind, @f$ E(k) @f$.  @see comp_ellint_2.
767*38fd1498Szrj    *
768*38fd1498Szrj    * @tparam _Tp The floating-point type of the modulus @c __k.
769*38fd1498Szrj    * @tparam _Tpp The floating-point type of the angle @c __phi.
770*38fd1498Szrj    * @param  __k  The modulus, <tt> abs(__k) <= 1 </tt>
771*38fd1498Szrj    * @param  __phi  The integral limit argument in radians
772*38fd1498Szrj    * @return  The elliptic function of the second kind.
773*38fd1498Szrj    * @throw std::domain_error if <tt> abs(__k) > 1 </tt>.
774*38fd1498Szrj    */
775*38fd1498Szrj   template<typename _Tp, typename _Tpp>
776*38fd1498Szrj     inline typename __gnu_cxx::__promote_2<_Tp, _Tpp>::__type
777*38fd1498Szrj     ellint_2(_Tp __k, _Tpp __phi)
778*38fd1498Szrj     {
779*38fd1498Szrj       typedef typename __gnu_cxx::__promote_2<_Tp, _Tpp>::__type __type;
780*38fd1498Szrj       return __detail::__ellint_2<__type>(__k, __phi);
781*38fd1498Szrj     }
782*38fd1498Szrj 
783*38fd1498Szrj   // Incomplete elliptic integrals of the third kind
784*38fd1498Szrj 
785*38fd1498Szrj   /**
786*38fd1498Szrj    * @brief Return the incomplete elliptic integral of the third kind
787*38fd1498Szrj    * @f$ \Pi(k,\nu,\phi) @f$ for @c float argument.
788*38fd1498Szrj    *
789*38fd1498Szrj    * @see ellint_3 for details.
790*38fd1498Szrj    */
791*38fd1498Szrj   inline float
792*38fd1498Szrj   ellint_3f(float __k, float __nu, float __phi)
793*38fd1498Szrj   { return __detail::__ellint_3<float>(__k, __nu, __phi); }
794*38fd1498Szrj 
795*38fd1498Szrj   /**
796*38fd1498Szrj    * @brief Return the incomplete elliptic integral of the third kind
797*38fd1498Szrj    * @f$ \Pi(k,\nu,\phi) @f$.
798*38fd1498Szrj    *
799*38fd1498Szrj    * @see ellint_3 for details.
800*38fd1498Szrj    */
801*38fd1498Szrj   inline long double
802*38fd1498Szrj   ellint_3l(long double __k, long double __nu, long double __phi)
803*38fd1498Szrj   { return __detail::__ellint_3<long double>(__k, __nu, __phi); }
804*38fd1498Szrj 
805*38fd1498Szrj   /**
806*38fd1498Szrj    * @brief Return the incomplete elliptic integral of the third kind
807*38fd1498Szrj    * @f$ \Pi(k,\nu,\phi) @f$.
808*38fd1498Szrj    *
809*38fd1498Szrj    * The incomplete elliptic integral of the third kind is defined by:
810*38fd1498Szrj    * @f[
811*38fd1498Szrj    *   \Pi(k,\nu,\phi) = \int_0^{\phi}
812*38fd1498Szrj    * 			 \frac{d\theta}
813*38fd1498Szrj    * 			 {(1 - \nu \sin^2\theta)
814*38fd1498Szrj    * 			  \sqrt{1 - k^2 \sin^2\theta}}
815*38fd1498Szrj    * @f]
816*38fd1498Szrj    * For  @f$ \phi= \pi/2 @f$ this becomes the complete elliptic integral of
817*38fd1498Szrj    * the third kind, @f$ \Pi(k,\nu) @f$.  @see comp_ellint_3.
818*38fd1498Szrj    *
819*38fd1498Szrj    * @tparam _Tp The floating-point type of the modulus @c __k.
820*38fd1498Szrj    * @tparam _Tpn The floating-point type of the argument @c __nu.
821*38fd1498Szrj    * @tparam _Tpp The floating-point type of the angle @c __phi.
822*38fd1498Szrj    * @param  __k  The modulus, <tt> abs(__k) <= 1 </tt>
823*38fd1498Szrj    * @param  __nu  The second argument
824*38fd1498Szrj    * @param  __phi  The integral limit argument in radians
825*38fd1498Szrj    * @return  The elliptic function of the third kind.
826*38fd1498Szrj    * @throw std::domain_error if <tt> abs(__k) > 1 </tt>.
827*38fd1498Szrj    */
828*38fd1498Szrj   template<typename _Tp, typename _Tpn, typename _Tpp>
829*38fd1498Szrj     inline typename __gnu_cxx::__promote_3<_Tp, _Tpn, _Tpp>::__type
830*38fd1498Szrj     ellint_3(_Tp __k, _Tpn __nu, _Tpp __phi)
831*38fd1498Szrj     {
832*38fd1498Szrj       typedef typename __gnu_cxx::__promote_3<_Tp, _Tpn, _Tpp>::__type __type;
833*38fd1498Szrj       return __detail::__ellint_3<__type>(__k, __nu, __phi);
834*38fd1498Szrj     }
835*38fd1498Szrj 
836*38fd1498Szrj   // Exponential integrals
837*38fd1498Szrj 
838*38fd1498Szrj   /**
839*38fd1498Szrj    * Return the exponential integral @f$ Ei(x) @f$ for @c float argument @c x.
840*38fd1498Szrj    *
841*38fd1498Szrj    * @see expint for details.
842*38fd1498Szrj    */
843*38fd1498Szrj   inline float
844*38fd1498Szrj   expintf(float __x)
845*38fd1498Szrj   { return __detail::__expint<float>(__x); }
846*38fd1498Szrj 
847*38fd1498Szrj   /**
848*38fd1498Szrj    * Return the exponential integral @f$ Ei(x) @f$
849*38fd1498Szrj    * for <tt>long double</tt> argument @c x.
850*38fd1498Szrj    *
851*38fd1498Szrj    * @see expint for details.
852*38fd1498Szrj    */
853*38fd1498Szrj   inline long double
854*38fd1498Szrj   expintl(long double __x)
855*38fd1498Szrj   { return __detail::__expint<long double>(__x); }
856*38fd1498Szrj 
857*38fd1498Szrj   /**
858*38fd1498Szrj    * Return the exponential integral @f$ Ei(x) @f$ for @c real argument @c x.
859*38fd1498Szrj    *
860*38fd1498Szrj    * The exponential integral is given by
861*38fd1498Szrj    * \f[
862*38fd1498Szrj    *   Ei(x) = -\int_{-x}^\infty \frac{e^t}{t} dt
863*38fd1498Szrj    * \f]
864*38fd1498Szrj    *
865*38fd1498Szrj    * @tparam _Tp The floating-point type of the argument @c __x.
866*38fd1498Szrj    * @param  __x  The argument of the exponential integral function.
867*38fd1498Szrj    */
868*38fd1498Szrj   template<typename _Tp>
869*38fd1498Szrj     inline typename __gnu_cxx::__promote<_Tp>::__type
870*38fd1498Szrj     expint(_Tp __x)
871*38fd1498Szrj     {
872*38fd1498Szrj       typedef typename __gnu_cxx::__promote<_Tp>::__type __type;
873*38fd1498Szrj       return __detail::__expint<__type>(__x);
874*38fd1498Szrj     }
875*38fd1498Szrj 
876*38fd1498Szrj   // Hermite polynomials
877*38fd1498Szrj 
878*38fd1498Szrj   /**
879*38fd1498Szrj    * Return the Hermite polynomial @f$ H_n(x) @f$ of nonnegative order n
880*38fd1498Szrj    * and float argument @c x.
881*38fd1498Szrj    *
882*38fd1498Szrj    * @see hermite for details.
883*38fd1498Szrj    */
884*38fd1498Szrj   inline float
885*38fd1498Szrj   hermitef(unsigned int __n, float __x)
886*38fd1498Szrj   { return __detail::__poly_hermite<float>(__n, __x); }
887*38fd1498Szrj 
888*38fd1498Szrj   /**
889*38fd1498Szrj    * Return the Hermite polynomial @f$ H_n(x) @f$ of nonnegative order n
890*38fd1498Szrj    * and <tt>long double</tt> argument @c x.
891*38fd1498Szrj    *
892*38fd1498Szrj    * @see hermite for details.
893*38fd1498Szrj    */
894*38fd1498Szrj   inline long double
895*38fd1498Szrj   hermitel(unsigned int __n, long double __x)
896*38fd1498Szrj   { return __detail::__poly_hermite<long double>(__n, __x); }
897*38fd1498Szrj 
898*38fd1498Szrj   /**
899*38fd1498Szrj    * Return the Hermite polynomial @f$ H_n(x) @f$ of order n
900*38fd1498Szrj    * and @c real argument @c x.
901*38fd1498Szrj    *
902*38fd1498Szrj    * The Hermite polynomial is defined by:
903*38fd1498Szrj    * @f[
904*38fd1498Szrj    *   H_n(x) = (-1)^n e^{x^2} \frac{d^n}{dx^n} e^{-x^2}
905*38fd1498Szrj    * @f]
906*38fd1498Szrj    *
907*38fd1498Szrj    * The Hermite polynomial obeys a reflection formula:
908*38fd1498Szrj    * @f[
909*38fd1498Szrj    *   H_n(-x) = (-1)^n H_n(x)
910*38fd1498Szrj    * @f]
911*38fd1498Szrj    *
912*38fd1498Szrj    * @tparam _Tp The floating-point type of the argument @c __x.
913*38fd1498Szrj    * @param __n The order
914*38fd1498Szrj    * @param __x The argument
915*38fd1498Szrj    */
916*38fd1498Szrj   template<typename _Tp>
917*38fd1498Szrj     inline typename __gnu_cxx::__promote<_Tp>::__type
918*38fd1498Szrj     hermite(unsigned int __n, _Tp __x)
919*38fd1498Szrj     {
920*38fd1498Szrj       typedef typename __gnu_cxx::__promote<_Tp>::__type __type;
921*38fd1498Szrj       return __detail::__poly_hermite<__type>(__n, __x);
922*38fd1498Szrj     }
923*38fd1498Szrj 
924*38fd1498Szrj   // Laguerre polynomials
925*38fd1498Szrj 
926*38fd1498Szrj   /**
927*38fd1498Szrj    * Returns the Laguerre polynomial @f$ L_n(x) @f$ of nonnegative degree @c n
928*38fd1498Szrj    * and @c float argument  @f$ x >= 0 @f$.
929*38fd1498Szrj    *
930*38fd1498Szrj    * @see laguerre for more details.
931*38fd1498Szrj    */
932*38fd1498Szrj   inline float
933*38fd1498Szrj   laguerref(unsigned int __n, float __x)
934*38fd1498Szrj   { return __detail::__laguerre<float>(__n, __x); }
935*38fd1498Szrj 
936*38fd1498Szrj   /**
937*38fd1498Szrj    * Returns the Laguerre polynomial @f$ L_n(x) @f$ of nonnegative degree @c n
938*38fd1498Szrj    * and <tt>long double</tt> argument @f$ x >= 0 @f$.
939*38fd1498Szrj    *
940*38fd1498Szrj    * @see laguerre for more details.
941*38fd1498Szrj    */
942*38fd1498Szrj   inline long double
943*38fd1498Szrj   laguerrel(unsigned int __n, long double __x)
944*38fd1498Szrj   { return __detail::__laguerre<long double>(__n, __x); }
945*38fd1498Szrj 
946*38fd1498Szrj   /**
947*38fd1498Szrj    * Returns the Laguerre polynomial @f$ L_n(x) @f$
948*38fd1498Szrj    * of nonnegative degree @c n and real argument @f$ x >= 0 @f$.
949*38fd1498Szrj    *
950*38fd1498Szrj    * The Laguerre polynomial is defined by:
951*38fd1498Szrj    * @f[
952*38fd1498Szrj    * 	 L_n(x) = \frac{e^x}{n!} \frac{d^n}{dx^n} (x^ne^{-x})
953*38fd1498Szrj    * @f]
954*38fd1498Szrj    *
955*38fd1498Szrj    * @tparam _Tp The floating-point type of the argument @c __x.
956*38fd1498Szrj    * @param __n The nonnegative order
957*38fd1498Szrj    * @param __x The argument <tt> __x >= 0 </tt>
958*38fd1498Szrj    * @throw std::domain_error if <tt> __x < 0 </tt>.
959*38fd1498Szrj    */
960*38fd1498Szrj   template<typename _Tp>
961*38fd1498Szrj     inline typename __gnu_cxx::__promote<_Tp>::__type
962*38fd1498Szrj     laguerre(unsigned int __n, _Tp __x)
963*38fd1498Szrj     {
964*38fd1498Szrj       typedef typename __gnu_cxx::__promote<_Tp>::__type __type;
965*38fd1498Szrj       return __detail::__laguerre<__type>(__n, __x);
966*38fd1498Szrj     }
967*38fd1498Szrj 
968*38fd1498Szrj   // Legendre polynomials
969*38fd1498Szrj 
970*38fd1498Szrj   /**
971*38fd1498Szrj    * Return the Legendre polynomial @f$ P_l(x) @f$ of nonnegative
972*38fd1498Szrj    * degree @f$ l @f$ and @c float argument @f$ |x| <= 0 @f$.
973*38fd1498Szrj    *
974*38fd1498Szrj    * @see legendre for more details.
975*38fd1498Szrj    */
976*38fd1498Szrj   inline float
977*38fd1498Szrj   legendref(unsigned int __l, float __x)
978*38fd1498Szrj   { return __detail::__poly_legendre_p<float>(__l, __x); }
979*38fd1498Szrj 
980*38fd1498Szrj   /**
981*38fd1498Szrj    * Return the Legendre polynomial @f$ P_l(x) @f$ of nonnegative
982*38fd1498Szrj    * degree @f$ l @f$ and <tt>long double</tt> argument @f$ |x| <= 0 @f$.
983*38fd1498Szrj    *
984*38fd1498Szrj    * @see legendre for more details.
985*38fd1498Szrj    */
986*38fd1498Szrj   inline long double
987*38fd1498Szrj   legendrel(unsigned int __l, long double __x)
988*38fd1498Szrj   { return __detail::__poly_legendre_p<long double>(__l, __x); }
989*38fd1498Szrj 
990*38fd1498Szrj   /**
991*38fd1498Szrj    * Return the Legendre polynomial @f$ P_l(x) @f$ of nonnegative
992*38fd1498Szrj    * degree @f$ l @f$ and real argument @f$ |x| <= 0 @f$.
993*38fd1498Szrj    *
994*38fd1498Szrj    * The Legendre function of order @f$ l @f$ and argument @f$ x @f$,
995*38fd1498Szrj    * @f$ P_l(x) @f$, is defined by:
996*38fd1498Szrj    * @f[
997*38fd1498Szrj    *   P_l(x) = \frac{1}{2^l l!}\frac{d^l}{dx^l}(x^2 - 1)^{l}
998*38fd1498Szrj    * @f]
999*38fd1498Szrj    *
1000*38fd1498Szrj    * @tparam _Tp The floating-point type of the argument @c __x.
1001*38fd1498Szrj    * @param __l The degree @f$ l >= 0 @f$
1002*38fd1498Szrj    * @param __x The argument @c abs(__x) <= 1
1003*38fd1498Szrj    * @throw std::domain_error if @c abs(__x) > 1
1004*38fd1498Szrj    */
1005*38fd1498Szrj   template<typename _Tp>
1006*38fd1498Szrj     inline typename __gnu_cxx::__promote<_Tp>::__type
1007*38fd1498Szrj     legendre(unsigned int __l, _Tp __x)
1008*38fd1498Szrj     {
1009*38fd1498Szrj       typedef typename __gnu_cxx::__promote<_Tp>::__type __type;
1010*38fd1498Szrj       return __detail::__poly_legendre_p<__type>(__l, __x);
1011*38fd1498Szrj     }
1012*38fd1498Szrj 
1013*38fd1498Szrj   // Riemann zeta functions
1014*38fd1498Szrj 
1015*38fd1498Szrj   /**
1016*38fd1498Szrj    * Return the Riemann zeta function @f$ \zeta(s) @f$
1017*38fd1498Szrj    * for @c float argument @f$ s @f$.
1018*38fd1498Szrj    *
1019*38fd1498Szrj    * @see riemann_zeta for more details.
1020*38fd1498Szrj    */
1021*38fd1498Szrj   inline float
1022*38fd1498Szrj   riemann_zetaf(float __s)
1023*38fd1498Szrj   { return __detail::__riemann_zeta<float>(__s); }
1024*38fd1498Szrj 
1025*38fd1498Szrj   /**
1026*38fd1498Szrj    * Return the Riemann zeta function @f$ \zeta(s) @f$
1027*38fd1498Szrj    * for <tt>long double</tt> argument @f$ s @f$.
1028*38fd1498Szrj    *
1029*38fd1498Szrj    * @see riemann_zeta for more details.
1030*38fd1498Szrj    */
1031*38fd1498Szrj   inline long double
1032*38fd1498Szrj   riemann_zetal(long double __s)
1033*38fd1498Szrj   { return __detail::__riemann_zeta<long double>(__s); }
1034*38fd1498Szrj 
1035*38fd1498Szrj   /**
1036*38fd1498Szrj    * Return the Riemann zeta function @f$ \zeta(s) @f$
1037*38fd1498Szrj    * for real argument @f$ s @f$.
1038*38fd1498Szrj    *
1039*38fd1498Szrj    * The Riemann zeta function is defined by:
1040*38fd1498Szrj    * @f[
1041*38fd1498Szrj    * 	\zeta(s) = \sum_{k=1}^{\infty} k^{-s} \hbox{ for } s > 1
1042*38fd1498Szrj    * @f]
1043*38fd1498Szrj    * and
1044*38fd1498Szrj    * @f[
1045*38fd1498Szrj    * 	\zeta(s) = \frac{1}{1-2^{1-s}}\sum_{k=1}^{\infty}(-1)^{k-1}k^{-s}
1046*38fd1498Szrj    *              \hbox{ for } 0 <= s <= 1
1047*38fd1498Szrj    * @f]
1048*38fd1498Szrj    * For s < 1 use the reflection formula:
1049*38fd1498Szrj    * @f[
1050*38fd1498Szrj    * 	\zeta(s) = 2^s \pi^{s-1} \sin(\frac{\pi s}{2}) \Gamma(1-s) \zeta(1-s)
1051*38fd1498Szrj    * @f]
1052*38fd1498Szrj    *
1053*38fd1498Szrj    * @tparam _Tp The floating-point type of the argument @c __s.
1054*38fd1498Szrj    * @param __s The argument <tt> s != 1 </tt>
1055*38fd1498Szrj    */
1056*38fd1498Szrj   template<typename _Tp>
1057*38fd1498Szrj     inline typename __gnu_cxx::__promote<_Tp>::__type
1058*38fd1498Szrj     riemann_zeta(_Tp __s)
1059*38fd1498Szrj     {
1060*38fd1498Szrj       typedef typename __gnu_cxx::__promote<_Tp>::__type __type;
1061*38fd1498Szrj       return __detail::__riemann_zeta<__type>(__s);
1062*38fd1498Szrj     }
1063*38fd1498Szrj 
1064*38fd1498Szrj   // Spherical Bessel functions
1065*38fd1498Szrj 
1066*38fd1498Szrj   /**
1067*38fd1498Szrj    * Return the spherical Bessel function @f$ j_n(x) @f$ of nonnegative order n
1068*38fd1498Szrj    * and @c float argument @f$ x >= 0 @f$.
1069*38fd1498Szrj    *
1070*38fd1498Szrj    * @see sph_bessel for more details.
1071*38fd1498Szrj    */
1072*38fd1498Szrj   inline float
1073*38fd1498Szrj   sph_besself(unsigned int __n, float __x)
1074*38fd1498Szrj   { return __detail::__sph_bessel<float>(__n, __x); }
1075*38fd1498Szrj 
1076*38fd1498Szrj   /**
1077*38fd1498Szrj    * Return the spherical Bessel function @f$ j_n(x) @f$ of nonnegative order n
1078*38fd1498Szrj    * and <tt>long double</tt> argument @f$ x >= 0 @f$.
1079*38fd1498Szrj    *
1080*38fd1498Szrj    * @see sph_bessel for more details.
1081*38fd1498Szrj    */
1082*38fd1498Szrj   inline long double
1083*38fd1498Szrj   sph_bessell(unsigned int __n, long double __x)
1084*38fd1498Szrj   { return __detail::__sph_bessel<long double>(__n, __x); }
1085*38fd1498Szrj 
1086*38fd1498Szrj   /**
1087*38fd1498Szrj    * Return the spherical Bessel function @f$ j_n(x) @f$ of nonnegative order n
1088*38fd1498Szrj    * and real argument @f$ x >= 0 @f$.
1089*38fd1498Szrj    *
1090*38fd1498Szrj    * The spherical Bessel function is defined by:
1091*38fd1498Szrj    * @f[
1092*38fd1498Szrj    *  j_n(x) = \left(\frac{\pi}{2x} \right) ^{1/2} J_{n+1/2}(x)
1093*38fd1498Szrj    * @f]
1094*38fd1498Szrj    *
1095*38fd1498Szrj    * @tparam _Tp The floating-point type of the argument @c __x.
1096*38fd1498Szrj    * @param  __n  The integral order <tt> n >= 0 </tt>
1097*38fd1498Szrj    * @param  __x  The real argument <tt> x >= 0 </tt>
1098*38fd1498Szrj    * @throw std::domain_error if <tt> __x < 0 </tt>.
1099*38fd1498Szrj    */
1100*38fd1498Szrj   template<typename _Tp>
1101*38fd1498Szrj     inline typename __gnu_cxx::__promote<_Tp>::__type
1102*38fd1498Szrj     sph_bessel(unsigned int __n, _Tp __x)
1103*38fd1498Szrj     {
1104*38fd1498Szrj       typedef typename __gnu_cxx::__promote<_Tp>::__type __type;
1105*38fd1498Szrj       return __detail::__sph_bessel<__type>(__n, __x);
1106*38fd1498Szrj     }
1107*38fd1498Szrj 
1108*38fd1498Szrj   // Spherical associated Legendre functions
1109*38fd1498Szrj 
1110*38fd1498Szrj   /**
1111*38fd1498Szrj    * Return the spherical Legendre function of nonnegative integral
1112*38fd1498Szrj    * degree @c l and order @c m and float angle @f$ \theta @f$ in radians.
1113*38fd1498Szrj    *
1114*38fd1498Szrj    * @see sph_legendre for details.
1115*38fd1498Szrj    */
1116*38fd1498Szrj   inline float
1117*38fd1498Szrj   sph_legendref(unsigned int __l, unsigned int __m, float __theta)
1118*38fd1498Szrj   { return __detail::__sph_legendre<float>(__l, __m, __theta); }
1119*38fd1498Szrj 
1120*38fd1498Szrj   /**
1121*38fd1498Szrj    * Return the spherical Legendre function of nonnegative integral
1122*38fd1498Szrj    * degree @c l and order @c m and <tt>long double</tt> angle @f$ \theta @f$
1123*38fd1498Szrj    * in radians.
1124*38fd1498Szrj    *
1125*38fd1498Szrj    * @see sph_legendre for details.
1126*38fd1498Szrj    */
1127*38fd1498Szrj   inline long double
1128*38fd1498Szrj   sph_legendrel(unsigned int __l, unsigned int __m, long double __theta)
1129*38fd1498Szrj   { return __detail::__sph_legendre<long double>(__l, __m, __theta); }
1130*38fd1498Szrj 
1131*38fd1498Szrj   /**
1132*38fd1498Szrj    * Return the spherical Legendre function of nonnegative integral
1133*38fd1498Szrj    * degree @c l and order @c m and real angle @f$ \theta @f$ in radians.
1134*38fd1498Szrj    *
1135*38fd1498Szrj    * The spherical Legendre function is defined by
1136*38fd1498Szrj    * @f[
1137*38fd1498Szrj    *  Y_l^m(\theta,\phi) = (-1)^m[\frac{(2l+1)}{4\pi}
1138*38fd1498Szrj    *                              \frac{(l-m)!}{(l+m)!}]
1139*38fd1498Szrj    *                   P_l^m(\cos\theta) \exp^{im\phi}
1140*38fd1498Szrj    * @f]
1141*38fd1498Szrj    *
1142*38fd1498Szrj    * @tparam _Tp The floating-point type of the angle @c __theta.
1143*38fd1498Szrj    * @param __l The order <tt> __l >= 0 </tt>
1144*38fd1498Szrj    * @param __m The degree <tt> __m >= 0 </tt> and <tt> __m <= __l </tt>
1145*38fd1498Szrj    * @param __theta The radian polar angle argument
1146*38fd1498Szrj    */
1147*38fd1498Szrj   template<typename _Tp>
1148*38fd1498Szrj     inline typename __gnu_cxx::__promote<_Tp>::__type
1149*38fd1498Szrj     sph_legendre(unsigned int __l, unsigned int __m, _Tp __theta)
1150*38fd1498Szrj     {
1151*38fd1498Szrj       typedef typename __gnu_cxx::__promote<_Tp>::__type __type;
1152*38fd1498Szrj       return __detail::__sph_legendre<__type>(__l, __m, __theta);
1153*38fd1498Szrj     }
1154*38fd1498Szrj 
1155*38fd1498Szrj   // Spherical Neumann functions
1156*38fd1498Szrj 
1157*38fd1498Szrj   /**
1158*38fd1498Szrj    * Return the spherical Neumann function of integral order @f$ n >= 0 @f$
1159*38fd1498Szrj    * and @c float argument @f$ x >= 0 @f$.
1160*38fd1498Szrj    *
1161*38fd1498Szrj    * @see sph_neumann for details.
1162*38fd1498Szrj    */
1163*38fd1498Szrj   inline float
1164*38fd1498Szrj   sph_neumannf(unsigned int __n, float __x)
1165*38fd1498Szrj   { return __detail::__sph_neumann<float>(__n, __x); }
1166*38fd1498Szrj 
1167*38fd1498Szrj   /**
1168*38fd1498Szrj    * Return the spherical Neumann function of integral order @f$ n >= 0 @f$
1169*38fd1498Szrj    * and <tt>long double</tt> @f$ x >= 0 @f$.
1170*38fd1498Szrj    *
1171*38fd1498Szrj    * @see sph_neumann for details.
1172*38fd1498Szrj    */
1173*38fd1498Szrj   inline long double
1174*38fd1498Szrj   sph_neumannl(unsigned int __n, long double __x)
1175*38fd1498Szrj   { return __detail::__sph_neumann<long double>(__n, __x); }
1176*38fd1498Szrj 
1177*38fd1498Szrj   /**
1178*38fd1498Szrj    * Return the spherical Neumann function of integral order @f$ n >= 0 @f$
1179*38fd1498Szrj    * and real argument @f$ x >= 0 @f$.
1180*38fd1498Szrj    *
1181*38fd1498Szrj    * The spherical Neumann function is defined by
1182*38fd1498Szrj    * @f[
1183*38fd1498Szrj    *    n_n(x) = \left(\frac{\pi}{2x} \right) ^{1/2} N_{n+1/2}(x)
1184*38fd1498Szrj    * @f]
1185*38fd1498Szrj    *
1186*38fd1498Szrj    * @tparam _Tp The floating-point type of the argument @c __x.
1187*38fd1498Szrj    * @param  __n  The integral order <tt> n >= 0 </tt>
1188*38fd1498Szrj    * @param  __x  The real argument <tt> __x >= 0 </tt>
1189*38fd1498Szrj    * @throw std::domain_error if <tt> __x < 0 </tt>.
1190*38fd1498Szrj    */
1191*38fd1498Szrj   template<typename _Tp>
1192*38fd1498Szrj     inline typename __gnu_cxx::__promote<_Tp>::__type
1193*38fd1498Szrj     sph_neumann(unsigned int __n, _Tp __x)
1194*38fd1498Szrj     {
1195*38fd1498Szrj       typedef typename __gnu_cxx::__promote<_Tp>::__type __type;
1196*38fd1498Szrj       return __detail::__sph_neumann<__type>(__n, __x);
1197*38fd1498Szrj     }
1198*38fd1498Szrj 
1199*38fd1498Szrj   // @} group mathsf
1200*38fd1498Szrj 
1201*38fd1498Szrj _GLIBCXX_END_NAMESPACE_VERSION
1202*38fd1498Szrj } // namespace std
1203*38fd1498Szrj 
1204*38fd1498Szrj #ifndef __STRICT_ANSI__
_GLIBCXX_VISIBILITY(default)1205*38fd1498Szrj namespace __gnu_cxx _GLIBCXX_VISIBILITY(default)
1206*38fd1498Szrj {
1207*38fd1498Szrj _GLIBCXX_BEGIN_NAMESPACE_VERSION
1208*38fd1498Szrj 
1209*38fd1498Szrj   // Airy functions
1210*38fd1498Szrj 
1211*38fd1498Szrj   /**
1212*38fd1498Szrj    * Return the Airy function @f$ Ai(x) @f$ of @c float argument x.
1213*38fd1498Szrj    */
1214*38fd1498Szrj   inline float
1215*38fd1498Szrj   airy_aif(float __x)
1216*38fd1498Szrj   {
1217*38fd1498Szrj     float __Ai, __Bi, __Aip, __Bip;
1218*38fd1498Szrj     std::__detail::__airy<float>(__x, __Ai, __Bi, __Aip, __Bip);
1219*38fd1498Szrj     return __Ai;
1220*38fd1498Szrj   }
1221*38fd1498Szrj 
1222*38fd1498Szrj   /**
1223*38fd1498Szrj    * Return the Airy function @f$ Ai(x) @f$ of <tt>long double</tt> argument x.
1224*38fd1498Szrj    */
1225*38fd1498Szrj   inline long double
1226*38fd1498Szrj   airy_ail(long double __x)
1227*38fd1498Szrj   {
1228*38fd1498Szrj     long double __Ai, __Bi, __Aip, __Bip;
1229*38fd1498Szrj     std::__detail::__airy<long double>(__x, __Ai, __Bi, __Aip, __Bip);
1230*38fd1498Szrj     return __Ai;
1231*38fd1498Szrj   }
1232*38fd1498Szrj 
1233*38fd1498Szrj   /**
1234*38fd1498Szrj    * Return the Airy function @f$ Ai(x) @f$ of real argument x.
1235*38fd1498Szrj    */
1236*38fd1498Szrj   template<typename _Tp>
1237*38fd1498Szrj     inline typename __gnu_cxx::__promote<_Tp>::__type
1238*38fd1498Szrj     airy_ai(_Tp __x)
1239*38fd1498Szrj     {
1240*38fd1498Szrj       typedef typename __gnu_cxx::__promote<_Tp>::__type __type;
1241*38fd1498Szrj       __type __Ai, __Bi, __Aip, __Bip;
1242*38fd1498Szrj       std::__detail::__airy<__type>(__x, __Ai, __Bi, __Aip, __Bip);
1243*38fd1498Szrj       return __Ai;
1244*38fd1498Szrj     }
1245*38fd1498Szrj 
1246*38fd1498Szrj   /**
1247*38fd1498Szrj    * Return the Airy function @f$ Bi(x) @f$ of @c float argument x.
1248*38fd1498Szrj    */
1249*38fd1498Szrj   inline float
1250*38fd1498Szrj   airy_bif(float __x)
1251*38fd1498Szrj   {
1252*38fd1498Szrj     float __Ai, __Bi, __Aip, __Bip;
1253*38fd1498Szrj     std::__detail::__airy<float>(__x, __Ai, __Bi, __Aip, __Bip);
1254*38fd1498Szrj     return __Bi;
1255*38fd1498Szrj   }
1256*38fd1498Szrj 
1257*38fd1498Szrj   /**
1258*38fd1498Szrj    * Return the Airy function @f$ Bi(x) @f$ of <tt>long double</tt> argument x.
1259*38fd1498Szrj    */
1260*38fd1498Szrj   inline long double
1261*38fd1498Szrj   airy_bil(long double __x)
1262*38fd1498Szrj   {
1263*38fd1498Szrj     long double __Ai, __Bi, __Aip, __Bip;
1264*38fd1498Szrj     std::__detail::__airy<long double>(__x, __Ai, __Bi, __Aip, __Bip);
1265*38fd1498Szrj     return __Bi;
1266*38fd1498Szrj   }
1267*38fd1498Szrj 
1268*38fd1498Szrj   /**
1269*38fd1498Szrj    * Return the Airy function @f$ Bi(x) @f$ of real argument x.
1270*38fd1498Szrj    */
1271*38fd1498Szrj   template<typename _Tp>
1272*38fd1498Szrj     inline typename __gnu_cxx::__promote<_Tp>::__type
1273*38fd1498Szrj     airy_bi(_Tp __x)
1274*38fd1498Szrj     {
1275*38fd1498Szrj       typedef typename __gnu_cxx::__promote<_Tp>::__type __type;
1276*38fd1498Szrj       __type __Ai, __Bi, __Aip, __Bip;
1277*38fd1498Szrj       std::__detail::__airy<__type>(__x, __Ai, __Bi, __Aip, __Bip);
1278*38fd1498Szrj       return __Bi;
1279*38fd1498Szrj     }
1280*38fd1498Szrj 
1281*38fd1498Szrj   // Confluent hypergeometric functions
1282*38fd1498Szrj 
1283*38fd1498Szrj   /**
1284*38fd1498Szrj    * Return the confluent hypergeometric function @f$ {}_1F_1(a;c;x) @f$
1285*38fd1498Szrj    * of @c float numeratorial parameter @c a, denominatorial parameter @c c,
1286*38fd1498Szrj    * and argument @c x.
1287*38fd1498Szrj    *
1288*38fd1498Szrj    * @see conf_hyperg for details.
1289*38fd1498Szrj    */
1290*38fd1498Szrj   inline float
1291*38fd1498Szrj   conf_hypergf(float __a, float __c, float __x)
1292*38fd1498Szrj   { return std::__detail::__conf_hyperg<float>(__a, __c, __x); }
1293*38fd1498Szrj 
1294*38fd1498Szrj   /**
1295*38fd1498Szrj    * Return the confluent hypergeometric function @f$ {}_1F_1(a;c;x) @f$
1296*38fd1498Szrj    * of <tt>long double</tt> numeratorial parameter @c a,
1297*38fd1498Szrj    * denominatorial parameter @c c, and argument @c x.
1298*38fd1498Szrj    *
1299*38fd1498Szrj    * @see conf_hyperg for details.
1300*38fd1498Szrj    */
1301*38fd1498Szrj   inline long double
1302*38fd1498Szrj   conf_hypergl(long double __a, long double __c, long double __x)
1303*38fd1498Szrj   { return std::__detail::__conf_hyperg<long double>(__a, __c, __x); }
1304*38fd1498Szrj 
1305*38fd1498Szrj   /**
1306*38fd1498Szrj    * Return the confluent hypergeometric function @f$ {}_1F_1(a;c;x) @f$
1307*38fd1498Szrj    * of real numeratorial parameter @c a, denominatorial parameter @c c,
1308*38fd1498Szrj    * and argument @c x.
1309*38fd1498Szrj    *
1310*38fd1498Szrj    * The confluent hypergeometric function is defined by
1311*38fd1498Szrj    * @f[
1312*38fd1498Szrj    *    {}_1F_1(a;c;x) = \sum_{n=0}^{\infty} \frac{(a)_n x^n}{(c)_n n!}
1313*38fd1498Szrj    * @f]
1314*38fd1498Szrj    * where the Pochhammer symbol is @f$ (x)_k = (x)(x+1)...(x+k-1) @f$,
1315*38fd1498Szrj    * @f$ (x)_0 = 1 @f$
1316*38fd1498Szrj    *
1317*38fd1498Szrj    * @param __a The numeratorial parameter
1318*38fd1498Szrj    * @param __c The denominatorial parameter
1319*38fd1498Szrj    * @param __x The argument
1320*38fd1498Szrj    */
1321*38fd1498Szrj   template<typename _Tpa, typename _Tpc, typename _Tp>
1322*38fd1498Szrj     inline typename __gnu_cxx::__promote_3<_Tpa, _Tpc, _Tp>::__type
1323*38fd1498Szrj     conf_hyperg(_Tpa __a, _Tpc __c, _Tp __x)
1324*38fd1498Szrj     {
1325*38fd1498Szrj       typedef typename __gnu_cxx::__promote_3<_Tpa, _Tpc, _Tp>::__type __type;
1326*38fd1498Szrj       return std::__detail::__conf_hyperg<__type>(__a, __c, __x);
1327*38fd1498Szrj     }
1328*38fd1498Szrj 
1329*38fd1498Szrj   // Hypergeometric functions
1330*38fd1498Szrj 
1331*38fd1498Szrj   /**
1332*38fd1498Szrj    * Return the hypergeometric function @f$ {}_2F_1(a,b;c;x) @f$
1333*38fd1498Szrj    * of @ float numeratorial parameters @c a and @c b,
1334*38fd1498Szrj    * denominatorial parameter @c c, and argument @c x.
1335*38fd1498Szrj    *
1336*38fd1498Szrj    * @see hyperg for details.
1337*38fd1498Szrj    */
1338*38fd1498Szrj   inline float
1339*38fd1498Szrj   hypergf(float __a, float __b, float __c, float __x)
1340*38fd1498Szrj   { return std::__detail::__hyperg<float>(__a, __b, __c, __x); }
1341*38fd1498Szrj 
1342*38fd1498Szrj   /**
1343*38fd1498Szrj    * Return the hypergeometric function @f$ {}_2F_1(a,b;c;x) @f$
1344*38fd1498Szrj    * of <tt>long double</tt> numeratorial parameters @c a and @c b,
1345*38fd1498Szrj    * denominatorial parameter @c c, and argument @c x.
1346*38fd1498Szrj    *
1347*38fd1498Szrj    * @see hyperg for details.
1348*38fd1498Szrj    */
1349*38fd1498Szrj   inline long double
1350*38fd1498Szrj   hypergl(long double __a, long double __b, long double __c, long double __x)
1351*38fd1498Szrj   { return std::__detail::__hyperg<long double>(__a, __b, __c, __x); }
1352*38fd1498Szrj 
1353*38fd1498Szrj   /**
1354*38fd1498Szrj    * Return the hypergeometric function @f$ {}_2F_1(a,b;c;x) @f$
1355*38fd1498Szrj    * of real numeratorial parameters @c a and @c b,
1356*38fd1498Szrj    * denominatorial parameter @c c, and argument @c x.
1357*38fd1498Szrj    *
1358*38fd1498Szrj    * The hypergeometric function is defined by
1359*38fd1498Szrj    * @f[
1360*38fd1498Szrj    *    {}_2F_1(a;c;x) = \sum_{n=0}^{\infty} \frac{(a)_n (b)_n x^n}{(c)_n n!}
1361*38fd1498Szrj    * @f]
1362*38fd1498Szrj    * where the Pochhammer symbol is @f$ (x)_k = (x)(x+1)...(x+k-1) @f$,
1363*38fd1498Szrj    * @f$ (x)_0 = 1 @f$
1364*38fd1498Szrj    *
1365*38fd1498Szrj    * @param __a The first numeratorial parameter
1366*38fd1498Szrj    * @param __b The second numeratorial parameter
1367*38fd1498Szrj    * @param __c The denominatorial parameter
1368*38fd1498Szrj    * @param __x The argument
1369*38fd1498Szrj    */
1370*38fd1498Szrj   template<typename _Tpa, typename _Tpb, typename _Tpc, typename _Tp>
1371*38fd1498Szrj     inline typename __gnu_cxx::__promote_4<_Tpa, _Tpb, _Tpc, _Tp>::__type
1372*38fd1498Szrj     hyperg(_Tpa __a, _Tpb __b, _Tpc __c, _Tp __x)
1373*38fd1498Szrj     {
1374*38fd1498Szrj       typedef typename __gnu_cxx::__promote_4<_Tpa, _Tpb, _Tpc, _Tp>
1375*38fd1498Szrj 		::__type __type;
1376*38fd1498Szrj       return std::__detail::__hyperg<__type>(__a, __b, __c, __x);
1377*38fd1498Szrj     }
1378*38fd1498Szrj 
1379*38fd1498Szrj _GLIBCXX_END_NAMESPACE_VERSION
1380*38fd1498Szrj } // namespace __gnu_cxx
1381*38fd1498Szrj #endif // __STRICT_ANSI__
1382*38fd1498Szrj 
1383*38fd1498Szrj #pragma GCC visibility pop
1384*38fd1498Szrj 
1385*38fd1498Szrj #endif // _GLIBCXX_BITS_SPECFUN_H
1386