specfun.h 46 KB

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