1 // Special functions -*- C++ -*-
3 // Copyright (C) 2006-2013 Free Software Foundation, Inc.
5 // This file is part of the GNU ISO C++ Library. This library is free
6 // software; you can redistribute it and/or modify it under the
7 // terms of the GNU General Public License as published by the
8 // Free Software Foundation; either version 3, or (at your option)
11 // This library is distributed in the hope that it will be useful,
12 // but WITHOUT ANY WARRANTY; without even the implied warranty of
13 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 // GNU General Public License for more details.
16 // Under Section 7 of GPL version 3, you are granted additional
17 // permissions described in the GCC Runtime Library Exception, version
18 // 3.1, as published by the Free Software Foundation.
20 // You should have received a copy of the GNU General Public License and
21 // a copy of the GCC Runtime Library Exception along with this program;
22 // see the files COPYING3 and COPYING.RUNTIME respectively. If not, see
23 // <http://www.gnu.org/licenses/>.
25 /** @file tr1/modified_bessel_func.tcc
26 * This is an internal header file, included by other library headers.
27 * Do not attempt to use it directly. @headername{tr1/cmath}
31 // ISO C++ 14882 TR1: 5.2 Special functions
34 // Written by Edward Smith-Rowland.
37 // (1) Handbook of Mathematical Functions,
38 // Ed. Milton Abramowitz and Irene A. Stegun,
39 // Dover Publications,
40 // Section 9, pp. 355-434, Section 10 pp. 435-478
41 // (2) The Gnu Scientific Library, http://www.gnu.org/software/gsl
42 // (3) Numerical Recipes in C, by W. H. Press, S. A. Teukolsky,
43 // W. T. Vetterling, B. P. Flannery, Cambridge University Press (1992),
44 // 2nd ed, pp. 246-249.
46 #ifndef _GLIBCXX_TR1_MODIFIED_BESSEL_FUNC_TCC
47 #define _GLIBCXX_TR1_MODIFIED_BESSEL_FUNC_TCC 1
49 #include "special_function_util.h"
51 namespace std _GLIBCXX_VISIBILITY(default)
55 // [5.2] Special functions
57 // Implementation-space details.
60 _GLIBCXX_BEGIN_NAMESPACE_VERSION
63 * @brief Compute the modified Bessel functions @f$ I_\nu(x) @f$ and
64 * @f$ K_\nu(x) @f$ and their first derivatives
65 * @f$ I'_\nu(x) @f$ and @f$ K'_\nu(x) @f$ respectively.
66 * These four functions are computed together for numerical
69 * @param __nu The order of the Bessel functions.
70 * @param __x The argument of the Bessel functions.
71 * @param __Inu The output regular modified Bessel function.
72 * @param __Knu The output irregular modified Bessel function.
73 * @param __Ipnu The output derivative of the regular
74 * modified Bessel function.
75 * @param __Kpnu The output derivative of the irregular
76 * modified Bessel function.
78 template <typename _Tp>
80 __bessel_ik(_Tp __nu, _Tp __x,
81 _Tp & __Inu, _Tp & __Knu, _Tp & __Ipnu, _Tp & __Kpnu)
90 else if (__nu == _Tp(1))
100 __Knu = std::numeric_limits<_Tp>::infinity();
101 __Kpnu = -std::numeric_limits<_Tp>::infinity();
105 const _Tp __eps = std::numeric_limits<_Tp>::epsilon();
106 const _Tp __fp_min = _Tp(10) * std::numeric_limits<_Tp>::epsilon();
107 const int __max_iter = 15000;
108 const _Tp __x_min = _Tp(2);
110 const int __nl = static_cast<int>(__nu + _Tp(0.5L));
112 const _Tp __mu = __nu - __nl;
113 const _Tp __mu2 = __mu * __mu;
114 const _Tp __xi = _Tp(1) / __x;
115 const _Tp __xi2 = _Tp(2) * __xi;
116 _Tp __h = __nu * __xi;
117 if ( __h < __fp_min )
119 _Tp __b = __xi2 * __nu;
123 for ( __i = 1; __i <= __max_iter; ++__i )
126 __d = _Tp(1) / (__b + __d);
127 __c = __b + _Tp(1) / __c;
128 const _Tp __del = __c * __d;
130 if (std::abs(__del - _Tp(1)) < __eps)
133 if (__i > __max_iter)
134 std::__throw_runtime_error(__N("Argument x too large "
136 "try asymptotic expansion."));
137 _Tp __Inul = __fp_min;
138 _Tp __Ipnul = __h * __Inul;
139 _Tp __Inul1 = __Inul;
140 _Tp __Ipnu1 = __Ipnul;
141 _Tp __fact = __nu * __xi;
142 for (int __l = __nl; __l >= 1; --__l)
144 const _Tp __Inutemp = __fact * __Inul + __Ipnul;
146 __Ipnul = __fact * __Inutemp + __Inul;
149 _Tp __f = __Ipnul / __Inul;
153 const _Tp __x2 = __x / _Tp(2);
154 const _Tp __pimu = __numeric_constants<_Tp>::__pi() * __mu;
155 const _Tp __fact = (std::abs(__pimu) < __eps
156 ? _Tp(1) : __pimu / std::sin(__pimu));
157 _Tp __d = -std::log(__x2);
158 _Tp __e = __mu * __d;
159 const _Tp __fact2 = (std::abs(__e) < __eps
160 ? _Tp(1) : std::sinh(__e) / __e);
161 _Tp __gam1, __gam2, __gampl, __gammi;
162 __gamma_temme(__mu, __gam1, __gam2, __gampl, __gammi);
164 * (__gam1 * std::cosh(__e) + __gam2 * __fact2 * __d);
167 _Tp __p = __e / (_Tp(2) * __gampl);
168 _Tp __q = _Tp(1) / (_Tp(2) * __e * __gammi);
173 for (__i = 1; __i <= __max_iter; ++__i)
175 __ff = (__i * __ff + __p + __q) / (__i * __i - __mu2);
179 const _Tp __del = __c * __ff;
181 const _Tp __del1 = __c * (__p - __i * __ff);
183 if (std::abs(__del) < __eps * std::abs(__sum))
186 if (__i > __max_iter)
187 std::__throw_runtime_error(__N("Bessel k series failed to converge "
190 __Knu1 = __sum1 * __xi2;
194 _Tp __b = _Tp(2) * (_Tp(1) + __x);
195 _Tp __d = _Tp(1) / __b;
200 _Tp __a1 = _Tp(0.25L) - __mu2;
201 _Tp __q = __c = __a1;
203 _Tp __s = _Tp(1) + __q * __delh;
205 for (__i = 2; __i <= __max_iter; ++__i)
207 __a -= 2 * (__i - 1);
208 __c = -__a * __c / __i;
209 const _Tp __qnew = (__q1 - __b * __q2) / __a;
214 __d = _Tp(1) / (__b + __a * __d);
215 __delh = (__b * __d - _Tp(1)) * __delh;
217 const _Tp __dels = __q * __delh;
219 if ( std::abs(__dels / __s) < __eps )
222 if (__i > __max_iter)
223 std::__throw_runtime_error(__N("Steed's method failed "
226 __Kmu = std::sqrt(__numeric_constants<_Tp>::__pi() / (_Tp(2) * __x))
227 * std::exp(-__x) / __s;
228 __Knu1 = __Kmu * (__mu + __x + _Tp(0.5L) - __h) * __xi;
231 _Tp __Kpmu = __mu * __xi * __Kmu - __Knu1;
232 _Tp __Inumu = __xi / (__f * __Kmu - __Kpmu);
233 __Inu = __Inumu * __Inul1 / __Inul;
234 __Ipnu = __Inumu * __Ipnu1 / __Inul;
235 for ( __i = 1; __i <= __nl; ++__i )
237 const _Tp __Knutemp = (__mu + __i) * __xi2 * __Knu1 + __Kmu;
242 __Kpnu = __nu * __xi * __Kmu - __Knu1;
249 * @brief Return the regular modified Bessel function of order
250 * \f$ \nu \f$: \f$ I_{\nu}(x) \f$.
252 * The regular modified cylindrical Bessel function is:
254 * I_{\nu}(x) = \sum_{k=0}^{\infty}
255 * \frac{(x/2)^{\nu + 2k}}{k!\Gamma(\nu+k+1)}
258 * @param __nu The order of the regular modified Bessel function.
259 * @param __x The argument of the regular modified Bessel function.
260 * @return The output regular modified Bessel function.
262 template<typename _Tp>
264 __cyl_bessel_i(_Tp __nu, _Tp __x)
266 if (__nu < _Tp(0) || __x < _Tp(0))
267 std::__throw_domain_error(__N("Bad argument "
268 "in __cyl_bessel_i."));
269 else if (__isnan(__nu) || __isnan(__x))
270 return std::numeric_limits<_Tp>::quiet_NaN();
271 else if (__x * __x < _Tp(10) * (__nu + _Tp(1)))
272 return __cyl_bessel_ij_series(__nu, __x, +_Tp(1), 200);
275 _Tp __I_nu, __K_nu, __Ip_nu, __Kp_nu;
276 __bessel_ik(__nu, __x, __I_nu, __K_nu, __Ip_nu, __Kp_nu);
283 * @brief Return the irregular modified Bessel function
284 * \f$ K_{\nu}(x) \f$ of order \f$ \nu \f$.
286 * The irregular modified Bessel function is defined by:
288 * K_{\nu}(x) = \frac{\pi}{2}
289 * \frac{I_{-\nu}(x) - I_{\nu}(x)}{\sin \nu\pi}
291 * where for integral \f$ \nu = n \f$ a limit is taken:
292 * \f$ lim_{\nu \to n} \f$.
294 * @param __nu The order of the irregular modified Bessel function.
295 * @param __x The argument of the irregular modified Bessel function.
296 * @return The output irregular modified Bessel function.
298 template<typename _Tp>
300 __cyl_bessel_k(_Tp __nu, _Tp __x)
302 if (__nu < _Tp(0) || __x < _Tp(0))
303 std::__throw_domain_error(__N("Bad argument "
304 "in __cyl_bessel_k."));
305 else if (__isnan(__nu) || __isnan(__x))
306 return std::numeric_limits<_Tp>::quiet_NaN();
309 _Tp __I_nu, __K_nu, __Ip_nu, __Kp_nu;
310 __bessel_ik(__nu, __x, __I_nu, __K_nu, __Ip_nu, __Kp_nu);
317 * @brief Compute the spherical modified Bessel functions
318 * @f$ i_n(x) @f$ and @f$ k_n(x) @f$ and their first
319 * derivatives @f$ i'_n(x) @f$ and @f$ k'_n(x) @f$
322 * @param __n The order of the modified spherical Bessel function.
323 * @param __x The argument of the modified spherical Bessel function.
324 * @param __i_n The output regular modified spherical Bessel function.
325 * @param __k_n The output irregular modified spherical
327 * @param __ip_n The output derivative of the regular modified
328 * spherical Bessel function.
329 * @param __kp_n The output derivative of the irregular modified
330 * spherical Bessel function.
332 template <typename _Tp>
334 __sph_bessel_ik(unsigned int __n, _Tp __x,
335 _Tp & __i_n, _Tp & __k_n, _Tp & __ip_n, _Tp & __kp_n)
337 const _Tp __nu = _Tp(__n) + _Tp(0.5L);
339 _Tp __I_nu, __Ip_nu, __K_nu, __Kp_nu;
340 __bessel_ik(__nu, __x, __I_nu, __K_nu, __Ip_nu, __Kp_nu);
342 const _Tp __factor = __numeric_constants<_Tp>::__sqrtpio2()
345 __i_n = __factor * __I_nu;
346 __k_n = __factor * __K_nu;
347 __ip_n = __factor * __Ip_nu - __i_n / (_Tp(2) * __x);
348 __kp_n = __factor * __Kp_nu - __k_n / (_Tp(2) * __x);
355 * @brief Compute the Airy functions
356 * @f$ Ai(x) @f$ and @f$ Bi(x) @f$ and their first
357 * derivatives @f$ Ai'(x) @f$ and @f$ Bi(x) @f$
360 * @param __n The order of the Airy functions.
361 * @param __x The argument of the Airy functions.
362 * @param __i_n The output Airy function.
363 * @param __k_n The output Airy function.
364 * @param __ip_n The output derivative of the Airy function.
365 * @param __kp_n The output derivative of the Airy function.
367 template <typename _Tp>
369 __airy(_Tp __x, _Tp & __Ai, _Tp & __Bi, _Tp & __Aip, _Tp & __Bip)
371 const _Tp __absx = std::abs(__x);
372 const _Tp __rootx = std::sqrt(__absx);
373 const _Tp __z = _Tp(2) * __absx * __rootx / _Tp(3);
376 return std::numeric_limits<_Tp>::quiet_NaN();
377 else if (__x > _Tp(0))
379 _Tp __I_nu, __Ip_nu, __K_nu, __Kp_nu;
381 __bessel_ik(_Tp(1) / _Tp(3), __z, __I_nu, __K_nu, __Ip_nu, __Kp_nu);
382 __Ai = __rootx * __K_nu
383 / (__numeric_constants<_Tp>::__sqrt3()
384 * __numeric_constants<_Tp>::__pi());
385 __Bi = __rootx * (__K_nu / __numeric_constants<_Tp>::__pi()
386 + _Tp(2) * __I_nu / __numeric_constants<_Tp>::__sqrt3());
388 __bessel_ik(_Tp(2) / _Tp(3), __z, __I_nu, __K_nu, __Ip_nu, __Kp_nu);
389 __Aip = -__x * __K_nu
390 / (__numeric_constants<_Tp>::__sqrt3()
391 * __numeric_constants<_Tp>::__pi());
392 __Bip = __x * (__K_nu / __numeric_constants<_Tp>::__pi()
394 / __numeric_constants<_Tp>::__sqrt3());
396 else if (__x < _Tp(0))
398 _Tp __J_nu, __Jp_nu, __N_nu, __Np_nu;
400 __bessel_jn(_Tp(1) / _Tp(3), __z, __J_nu, __N_nu, __Jp_nu, __Np_nu);
401 __Ai = __rootx * (__J_nu
402 - __N_nu / __numeric_constants<_Tp>::__sqrt3()) / _Tp(2);
403 __Bi = -__rootx * (__N_nu
404 + __J_nu / __numeric_constants<_Tp>::__sqrt3()) / _Tp(2);
406 __bessel_jn(_Tp(2) / _Tp(3), __z, __J_nu, __N_nu, __Jp_nu, __Np_nu);
407 __Aip = __absx * (__N_nu / __numeric_constants<_Tp>::__sqrt3()
409 __Bip = __absx * (__J_nu / __numeric_constants<_Tp>::__sqrt3()
415 // Abramowitz & Stegun, page 446 section 10.4.4 on Airy functions.
416 // The number is Ai(0) = 3^{-2/3}/\Gamma(2/3).
417 __Ai = _Tp(0.35502805388781723926L);
418 __Bi = __Ai * __numeric_constants<_Tp>::__sqrt3();
421 // Abramowitz & Stegun, page 446 section 10.4.5 on Airy functions.
422 // The number is Ai'(0) = -3^{-1/3}/\Gamma(1/3).
423 __Aip = -_Tp(0.25881940379280679840L);
424 __Bip = -__Aip * __numeric_constants<_Tp>::__sqrt3();
430 _GLIBCXX_END_NAMESPACE_VERSION
431 } // namespace std::tr1::__detail
435 #endif // _GLIBCXX_TR1_MODIFIED_BESSEL_FUNC_TCC