specfun.h 46 KB

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