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