1*38fd1498Szrj // Special functions -*- C++ -*- 2*38fd1498Szrj 3*38fd1498Szrj // Copyright (C) 2006-2018 Free Software Foundation, Inc. 4*38fd1498Szrj // 5*38fd1498Szrj // This file is part of the GNU ISO C++ Library. This library is free 6*38fd1498Szrj // software; you can redistribute it and/or modify it under the 7*38fd1498Szrj // terms of the GNU General Public License as published by the 8*38fd1498Szrj // Free Software Foundation; either version 3, or (at your option) 9*38fd1498Szrj // any later version. 10*38fd1498Szrj // 11*38fd1498Szrj // This library is distributed in the hope that it will be useful, 12*38fd1498Szrj // but WITHOUT ANY WARRANTY; without even the implied warranty of 13*38fd1498Szrj // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 14*38fd1498Szrj // GNU General Public License for more details. 15*38fd1498Szrj // 16*38fd1498Szrj // Under Section 7 of GPL version 3, you are granted additional 17*38fd1498Szrj // permissions described in the GCC Runtime Library Exception, version 18*38fd1498Szrj // 3.1, as published by the Free Software Foundation. 19*38fd1498Szrj 20*38fd1498Szrj // You should have received a copy of the GNU General Public License and 21*38fd1498Szrj // a copy of the GCC Runtime Library Exception along with this program; 22*38fd1498Szrj // see the files COPYING3 and COPYING.RUNTIME respectively. If not, see 23*38fd1498Szrj // <http://www.gnu.org/licenses/>. 24*38fd1498Szrj 25*38fd1498Szrj /** @file tr1/modified_bessel_func.tcc 26*38fd1498Szrj * This is an internal header file, included by other library headers. 27*38fd1498Szrj * Do not attempt to use it directly. @headername{tr1/cmath} 28*38fd1498Szrj */ 29*38fd1498Szrj 30*38fd1498Szrj // 31*38fd1498Szrj // ISO C++ 14882 TR1: 5.2 Special functions 32*38fd1498Szrj // 33*38fd1498Szrj 34*38fd1498Szrj // Written by Edward Smith-Rowland. 35*38fd1498Szrj // 36*38fd1498Szrj // References: 37*38fd1498Szrj // (1) Handbook of Mathematical Functions, 38*38fd1498Szrj // Ed. Milton Abramowitz and Irene A. Stegun, 39*38fd1498Szrj // Dover Publications, 40*38fd1498Szrj // Section 9, pp. 355-434, Section 10 pp. 435-478 41*38fd1498Szrj // (2) The Gnu Scientific Library, http://www.gnu.org/software/gsl 42*38fd1498Szrj // (3) Numerical Recipes in C, by W. H. Press, S. A. Teukolsky, 43*38fd1498Szrj // W. T. Vetterling, B. P. Flannery, Cambridge University Press (1992), 44*38fd1498Szrj // 2nd ed, pp. 246-249. 45*38fd1498Szrj 46*38fd1498Szrj #ifndef _GLIBCXX_TR1_MODIFIED_BESSEL_FUNC_TCC 47*38fd1498Szrj #define _GLIBCXX_TR1_MODIFIED_BESSEL_FUNC_TCC 1 48*38fd1498Szrj 49*38fd1498Szrj #include "special_function_util.h" 50*38fd1498Szrj 51*38fd1498Szrj namespace std _GLIBCXX_VISIBILITY(default) 52*38fd1498Szrj { 53*38fd1498Szrj _GLIBCXX_BEGIN_NAMESPACE_VERSION 54*38fd1498Szrj 55*38fd1498Szrj #if _GLIBCXX_USE_STD_SPEC_FUNCS 56*38fd1498Szrj #elif defined(_GLIBCXX_TR1_CMATH) 57*38fd1498Szrj namespace tr1 58*38fd1498Szrj { 59*38fd1498Szrj #else 60*38fd1498Szrj # error do not include this header directly, use <cmath> or <tr1/cmath> 61*38fd1498Szrj #endif 62*38fd1498Szrj // [5.2] Special functions 63*38fd1498Szrj 64*38fd1498Szrj // Implementation-space details. 65*38fd1498Szrj namespace __detail 66*38fd1498Szrj { 67*38fd1498Szrj /** 68*38fd1498Szrj * @brief Compute the modified Bessel functions @f$ I_\nu(x) @f$ and 69*38fd1498Szrj * @f$ K_\nu(x) @f$ and their first derivatives 70*38fd1498Szrj * @f$ I'_\nu(x) @f$ and @f$ K'_\nu(x) @f$ respectively. 71*38fd1498Szrj * These four functions are computed together for numerical 72*38fd1498Szrj * stability. 73*38fd1498Szrj * 74*38fd1498Szrj * @param __nu The order of the Bessel functions. 75*38fd1498Szrj * @param __x The argument of the Bessel functions. 76*38fd1498Szrj * @param __Inu The output regular modified Bessel function. 77*38fd1498Szrj * @param __Knu The output irregular modified Bessel function. 78*38fd1498Szrj * @param __Ipnu The output derivative of the regular 79*38fd1498Szrj * modified Bessel function. 80*38fd1498Szrj * @param __Kpnu The output derivative of the irregular 81*38fd1498Szrj * modified Bessel function. 82*38fd1498Szrj */ 83*38fd1498Szrj template <typename _Tp> 84*38fd1498Szrj void __bessel_ik(_Tp __nu,_Tp __x,_Tp & __Inu,_Tp & __Knu,_Tp & __Ipnu,_Tp & __Kpnu)85*38fd1498Szrj __bessel_ik(_Tp __nu, _Tp __x, 86*38fd1498Szrj _Tp & __Inu, _Tp & __Knu, _Tp & __Ipnu, _Tp & __Kpnu) 87*38fd1498Szrj { 88*38fd1498Szrj if (__x == _Tp(0)) 89*38fd1498Szrj { 90*38fd1498Szrj if (__nu == _Tp(0)) 91*38fd1498Szrj { 92*38fd1498Szrj __Inu = _Tp(1); 93*38fd1498Szrj __Ipnu = _Tp(0); 94*38fd1498Szrj } 95*38fd1498Szrj else if (__nu == _Tp(1)) 96*38fd1498Szrj { 97*38fd1498Szrj __Inu = _Tp(0); 98*38fd1498Szrj __Ipnu = _Tp(0.5L); 99*38fd1498Szrj } 100*38fd1498Szrj else 101*38fd1498Szrj { 102*38fd1498Szrj __Inu = _Tp(0); 103*38fd1498Szrj __Ipnu = _Tp(0); 104*38fd1498Szrj } 105*38fd1498Szrj __Knu = std::numeric_limits<_Tp>::infinity(); 106*38fd1498Szrj __Kpnu = -std::numeric_limits<_Tp>::infinity(); 107*38fd1498Szrj return; 108*38fd1498Szrj } 109*38fd1498Szrj 110*38fd1498Szrj const _Tp __eps = std::numeric_limits<_Tp>::epsilon(); 111*38fd1498Szrj const _Tp __fp_min = _Tp(10) * std::numeric_limits<_Tp>::epsilon(); 112*38fd1498Szrj const int __max_iter = 15000; 113*38fd1498Szrj const _Tp __x_min = _Tp(2); 114*38fd1498Szrj 115*38fd1498Szrj const int __nl = static_cast<int>(__nu + _Tp(0.5L)); 116*38fd1498Szrj 117*38fd1498Szrj const _Tp __mu = __nu - __nl; 118*38fd1498Szrj const _Tp __mu2 = __mu * __mu; 119*38fd1498Szrj const _Tp __xi = _Tp(1) / __x; 120*38fd1498Szrj const _Tp __xi2 = _Tp(2) * __xi; 121*38fd1498Szrj _Tp __h = __nu * __xi; 122*38fd1498Szrj if ( __h < __fp_min ) 123*38fd1498Szrj __h = __fp_min; 124*38fd1498Szrj _Tp __b = __xi2 * __nu; 125*38fd1498Szrj _Tp __d = _Tp(0); 126*38fd1498Szrj _Tp __c = __h; 127*38fd1498Szrj int __i; 128*38fd1498Szrj for ( __i = 1; __i <= __max_iter; ++__i ) 129*38fd1498Szrj { 130*38fd1498Szrj __b += __xi2; 131*38fd1498Szrj __d = _Tp(1) / (__b + __d); 132*38fd1498Szrj __c = __b + _Tp(1) / __c; 133*38fd1498Szrj const _Tp __del = __c * __d; 134*38fd1498Szrj __h *= __del; 135*38fd1498Szrj if (std::abs(__del - _Tp(1)) < __eps) 136*38fd1498Szrj break; 137*38fd1498Szrj } 138*38fd1498Szrj if (__i > __max_iter) 139*38fd1498Szrj std::__throw_runtime_error(__N("Argument x too large " 140*38fd1498Szrj "in __bessel_ik; " 141*38fd1498Szrj "try asymptotic expansion.")); 142*38fd1498Szrj _Tp __Inul = __fp_min; 143*38fd1498Szrj _Tp __Ipnul = __h * __Inul; 144*38fd1498Szrj _Tp __Inul1 = __Inul; 145*38fd1498Szrj _Tp __Ipnu1 = __Ipnul; 146*38fd1498Szrj _Tp __fact = __nu * __xi; 147*38fd1498Szrj for (int __l = __nl; __l >= 1; --__l) 148*38fd1498Szrj { 149*38fd1498Szrj const _Tp __Inutemp = __fact * __Inul + __Ipnul; 150*38fd1498Szrj __fact -= __xi; 151*38fd1498Szrj __Ipnul = __fact * __Inutemp + __Inul; 152*38fd1498Szrj __Inul = __Inutemp; 153*38fd1498Szrj } 154*38fd1498Szrj _Tp __f = __Ipnul / __Inul; 155*38fd1498Szrj _Tp __Kmu, __Knu1; 156*38fd1498Szrj if (__x < __x_min) 157*38fd1498Szrj { 158*38fd1498Szrj const _Tp __x2 = __x / _Tp(2); 159*38fd1498Szrj const _Tp __pimu = __numeric_constants<_Tp>::__pi() * __mu; 160*38fd1498Szrj const _Tp __fact = (std::abs(__pimu) < __eps 161*38fd1498Szrj ? _Tp(1) : __pimu / std::sin(__pimu)); 162*38fd1498Szrj _Tp __d = -std::log(__x2); 163*38fd1498Szrj _Tp __e = __mu * __d; 164*38fd1498Szrj const _Tp __fact2 = (std::abs(__e) < __eps 165*38fd1498Szrj ? _Tp(1) : std::sinh(__e) / __e); 166*38fd1498Szrj _Tp __gam1, __gam2, __gampl, __gammi; 167*38fd1498Szrj __gamma_temme(__mu, __gam1, __gam2, __gampl, __gammi); 168*38fd1498Szrj _Tp __ff = __fact 169*38fd1498Szrj * (__gam1 * std::cosh(__e) + __gam2 * __fact2 * __d); 170*38fd1498Szrj _Tp __sum = __ff; 171*38fd1498Szrj __e = std::exp(__e); 172*38fd1498Szrj _Tp __p = __e / (_Tp(2) * __gampl); 173*38fd1498Szrj _Tp __q = _Tp(1) / (_Tp(2) * __e * __gammi); 174*38fd1498Szrj _Tp __c = _Tp(1); 175*38fd1498Szrj __d = __x2 * __x2; 176*38fd1498Szrj _Tp __sum1 = __p; 177*38fd1498Szrj int __i; 178*38fd1498Szrj for (__i = 1; __i <= __max_iter; ++__i) 179*38fd1498Szrj { 180*38fd1498Szrj __ff = (__i * __ff + __p + __q) / (__i * __i - __mu2); 181*38fd1498Szrj __c *= __d / __i; 182*38fd1498Szrj __p /= __i - __mu; 183*38fd1498Szrj __q /= __i + __mu; 184*38fd1498Szrj const _Tp __del = __c * __ff; 185*38fd1498Szrj __sum += __del; 186*38fd1498Szrj const _Tp __del1 = __c * (__p - __i * __ff); 187*38fd1498Szrj __sum1 += __del1; 188*38fd1498Szrj if (std::abs(__del) < __eps * std::abs(__sum)) 189*38fd1498Szrj break; 190*38fd1498Szrj } 191*38fd1498Szrj if (__i > __max_iter) 192*38fd1498Szrj std::__throw_runtime_error(__N("Bessel k series failed to converge " 193*38fd1498Szrj "in __bessel_ik.")); 194*38fd1498Szrj __Kmu = __sum; 195*38fd1498Szrj __Knu1 = __sum1 * __xi2; 196*38fd1498Szrj } 197*38fd1498Szrj else 198*38fd1498Szrj { 199*38fd1498Szrj _Tp __b = _Tp(2) * (_Tp(1) + __x); 200*38fd1498Szrj _Tp __d = _Tp(1) / __b; 201*38fd1498Szrj _Tp __delh = __d; 202*38fd1498Szrj _Tp __h = __delh; 203*38fd1498Szrj _Tp __q1 = _Tp(0); 204*38fd1498Szrj _Tp __q2 = _Tp(1); 205*38fd1498Szrj _Tp __a1 = _Tp(0.25L) - __mu2; 206*38fd1498Szrj _Tp __q = __c = __a1; 207*38fd1498Szrj _Tp __a = -__a1; 208*38fd1498Szrj _Tp __s = _Tp(1) + __q * __delh; 209*38fd1498Szrj int __i; 210*38fd1498Szrj for (__i = 2; __i <= __max_iter; ++__i) 211*38fd1498Szrj { 212*38fd1498Szrj __a -= 2 * (__i - 1); 213*38fd1498Szrj __c = -__a * __c / __i; 214*38fd1498Szrj const _Tp __qnew = (__q1 - __b * __q2) / __a; 215*38fd1498Szrj __q1 = __q2; 216*38fd1498Szrj __q2 = __qnew; 217*38fd1498Szrj __q += __c * __qnew; 218*38fd1498Szrj __b += _Tp(2); 219*38fd1498Szrj __d = _Tp(1) / (__b + __a * __d); 220*38fd1498Szrj __delh = (__b * __d - _Tp(1)) * __delh; 221*38fd1498Szrj __h += __delh; 222*38fd1498Szrj const _Tp __dels = __q * __delh; 223*38fd1498Szrj __s += __dels; 224*38fd1498Szrj if ( std::abs(__dels / __s) < __eps ) 225*38fd1498Szrj break; 226*38fd1498Szrj } 227*38fd1498Szrj if (__i > __max_iter) 228*38fd1498Szrj std::__throw_runtime_error(__N("Steed's method failed " 229*38fd1498Szrj "in __bessel_ik.")); 230*38fd1498Szrj __h = __a1 * __h; 231*38fd1498Szrj __Kmu = std::sqrt(__numeric_constants<_Tp>::__pi() / (_Tp(2) * __x)) 232*38fd1498Szrj * std::exp(-__x) / __s; 233*38fd1498Szrj __Knu1 = __Kmu * (__mu + __x + _Tp(0.5L) - __h) * __xi; 234*38fd1498Szrj } 235*38fd1498Szrj 236*38fd1498Szrj _Tp __Kpmu = __mu * __xi * __Kmu - __Knu1; 237*38fd1498Szrj _Tp __Inumu = __xi / (__f * __Kmu - __Kpmu); 238*38fd1498Szrj __Inu = __Inumu * __Inul1 / __Inul; 239*38fd1498Szrj __Ipnu = __Inumu * __Ipnu1 / __Inul; 240*38fd1498Szrj for ( __i = 1; __i <= __nl; ++__i ) 241*38fd1498Szrj { 242*38fd1498Szrj const _Tp __Knutemp = (__mu + __i) * __xi2 * __Knu1 + __Kmu; 243*38fd1498Szrj __Kmu = __Knu1; 244*38fd1498Szrj __Knu1 = __Knutemp; 245*38fd1498Szrj } 246*38fd1498Szrj __Knu = __Kmu; 247*38fd1498Szrj __Kpnu = __nu * __xi * __Kmu - __Knu1; 248*38fd1498Szrj 249*38fd1498Szrj return; 250*38fd1498Szrj } 251*38fd1498Szrj 252*38fd1498Szrj 253*38fd1498Szrj /** 254*38fd1498Szrj * @brief Return the regular modified Bessel function of order 255*38fd1498Szrj * \f$ \nu \f$: \f$ I_{\nu}(x) \f$. 256*38fd1498Szrj * 257*38fd1498Szrj * The regular modified cylindrical Bessel function is: 258*38fd1498Szrj * @f[ 259*38fd1498Szrj * I_{\nu}(x) = \sum_{k=0}^{\infty} 260*38fd1498Szrj * \frac{(x/2)^{\nu + 2k}}{k!\Gamma(\nu+k+1)} 261*38fd1498Szrj * @f] 262*38fd1498Szrj * 263*38fd1498Szrj * @param __nu The order of the regular modified Bessel function. 264*38fd1498Szrj * @param __x The argument of the regular modified Bessel function. 265*38fd1498Szrj * @return The output regular modified Bessel function. 266*38fd1498Szrj */ 267*38fd1498Szrj template<typename _Tp> 268*38fd1498Szrj _Tp __cyl_bessel_i(_Tp __nu,_Tp __x)269*38fd1498Szrj __cyl_bessel_i(_Tp __nu, _Tp __x) 270*38fd1498Szrj { 271*38fd1498Szrj if (__nu < _Tp(0) || __x < _Tp(0)) 272*38fd1498Szrj std::__throw_domain_error(__N("Bad argument " 273*38fd1498Szrj "in __cyl_bessel_i.")); 274*38fd1498Szrj else if (__isnan(__nu) || __isnan(__x)) 275*38fd1498Szrj return std::numeric_limits<_Tp>::quiet_NaN(); 276*38fd1498Szrj else if (__x * __x < _Tp(10) * (__nu + _Tp(1))) 277*38fd1498Szrj return __cyl_bessel_ij_series(__nu, __x, +_Tp(1), 200); 278*38fd1498Szrj else 279*38fd1498Szrj { 280*38fd1498Szrj _Tp __I_nu, __K_nu, __Ip_nu, __Kp_nu; 281*38fd1498Szrj __bessel_ik(__nu, __x, __I_nu, __K_nu, __Ip_nu, __Kp_nu); 282*38fd1498Szrj return __I_nu; 283*38fd1498Szrj } 284*38fd1498Szrj } 285*38fd1498Szrj 286*38fd1498Szrj 287*38fd1498Szrj /** 288*38fd1498Szrj * @brief Return the irregular modified Bessel function 289*38fd1498Szrj * \f$ K_{\nu}(x) \f$ of order \f$ \nu \f$. 290*38fd1498Szrj * 291*38fd1498Szrj * The irregular modified Bessel function is defined by: 292*38fd1498Szrj * @f[ 293*38fd1498Szrj * K_{\nu}(x) = \frac{\pi}{2} 294*38fd1498Szrj * \frac{I_{-\nu}(x) - I_{\nu}(x)}{\sin \nu\pi} 295*38fd1498Szrj * @f] 296*38fd1498Szrj * where for integral \f$ \nu = n \f$ a limit is taken: 297*38fd1498Szrj * \f$ lim_{\nu \to n} \f$. 298*38fd1498Szrj * 299*38fd1498Szrj * @param __nu The order of the irregular modified Bessel function. 300*38fd1498Szrj * @param __x The argument of the irregular modified Bessel function. 301*38fd1498Szrj * @return The output irregular modified Bessel function. 302*38fd1498Szrj */ 303*38fd1498Szrj template<typename _Tp> 304*38fd1498Szrj _Tp __cyl_bessel_k(_Tp __nu,_Tp __x)305*38fd1498Szrj __cyl_bessel_k(_Tp __nu, _Tp __x) 306*38fd1498Szrj { 307*38fd1498Szrj if (__nu < _Tp(0) || __x < _Tp(0)) 308*38fd1498Szrj std::__throw_domain_error(__N("Bad argument " 309*38fd1498Szrj "in __cyl_bessel_k.")); 310*38fd1498Szrj else if (__isnan(__nu) || __isnan(__x)) 311*38fd1498Szrj return std::numeric_limits<_Tp>::quiet_NaN(); 312*38fd1498Szrj else 313*38fd1498Szrj { 314*38fd1498Szrj _Tp __I_nu, __K_nu, __Ip_nu, __Kp_nu; 315*38fd1498Szrj __bessel_ik(__nu, __x, __I_nu, __K_nu, __Ip_nu, __Kp_nu); 316*38fd1498Szrj return __K_nu; 317*38fd1498Szrj } 318*38fd1498Szrj } 319*38fd1498Szrj 320*38fd1498Szrj 321*38fd1498Szrj /** 322*38fd1498Szrj * @brief Compute the spherical modified Bessel functions 323*38fd1498Szrj * @f$ i_n(x) @f$ and @f$ k_n(x) @f$ and their first 324*38fd1498Szrj * derivatives @f$ i'_n(x) @f$ and @f$ k'_n(x) @f$ 325*38fd1498Szrj * respectively. 326*38fd1498Szrj * 327*38fd1498Szrj * @param __n The order of the modified spherical Bessel function. 328*38fd1498Szrj * @param __x The argument of the modified spherical Bessel function. 329*38fd1498Szrj * @param __i_n The output regular modified spherical Bessel function. 330*38fd1498Szrj * @param __k_n The output irregular modified spherical 331*38fd1498Szrj * Bessel function. 332*38fd1498Szrj * @param __ip_n The output derivative of the regular modified 333*38fd1498Szrj * spherical Bessel function. 334*38fd1498Szrj * @param __kp_n The output derivative of the irregular modified 335*38fd1498Szrj * spherical Bessel function. 336*38fd1498Szrj */ 337*38fd1498Szrj template <typename _Tp> 338*38fd1498Szrj void __sph_bessel_ik(unsigned int __n,_Tp __x,_Tp & __i_n,_Tp & __k_n,_Tp & __ip_n,_Tp & __kp_n)339*38fd1498Szrj __sph_bessel_ik(unsigned int __n, _Tp __x, 340*38fd1498Szrj _Tp & __i_n, _Tp & __k_n, _Tp & __ip_n, _Tp & __kp_n) 341*38fd1498Szrj { 342*38fd1498Szrj const _Tp __nu = _Tp(__n) + _Tp(0.5L); 343*38fd1498Szrj 344*38fd1498Szrj _Tp __I_nu, __Ip_nu, __K_nu, __Kp_nu; 345*38fd1498Szrj __bessel_ik(__nu, __x, __I_nu, __K_nu, __Ip_nu, __Kp_nu); 346*38fd1498Szrj 347*38fd1498Szrj const _Tp __factor = __numeric_constants<_Tp>::__sqrtpio2() 348*38fd1498Szrj / std::sqrt(__x); 349*38fd1498Szrj 350*38fd1498Szrj __i_n = __factor * __I_nu; 351*38fd1498Szrj __k_n = __factor * __K_nu; 352*38fd1498Szrj __ip_n = __factor * __Ip_nu - __i_n / (_Tp(2) * __x); 353*38fd1498Szrj __kp_n = __factor * __Kp_nu - __k_n / (_Tp(2) * __x); 354*38fd1498Szrj 355*38fd1498Szrj return; 356*38fd1498Szrj } 357*38fd1498Szrj 358*38fd1498Szrj 359*38fd1498Szrj /** 360*38fd1498Szrj * @brief Compute the Airy functions 361*38fd1498Szrj * @f$ Ai(x) @f$ and @f$ Bi(x) @f$ and their first 362*38fd1498Szrj * derivatives @f$ Ai'(x) @f$ and @f$ Bi(x) @f$ 363*38fd1498Szrj * respectively. 364*38fd1498Szrj * 365*38fd1498Szrj * @param __x The argument of the Airy functions. 366*38fd1498Szrj * @param __Ai The output Airy function of the first kind. 367*38fd1498Szrj * @param __Bi The output Airy function of the second kind. 368*38fd1498Szrj * @param __Aip The output derivative of the Airy function 369*38fd1498Szrj * of the first kind. 370*38fd1498Szrj * @param __Bip The output derivative of the Airy function 371*38fd1498Szrj * of the second kind. 372*38fd1498Szrj */ 373*38fd1498Szrj template <typename _Tp> 374*38fd1498Szrj void __airy(_Tp __x,_Tp & __Ai,_Tp & __Bi,_Tp & __Aip,_Tp & __Bip)375*38fd1498Szrj __airy(_Tp __x, _Tp & __Ai, _Tp & __Bi, _Tp & __Aip, _Tp & __Bip) 376*38fd1498Szrj { 377*38fd1498Szrj const _Tp __absx = std::abs(__x); 378*38fd1498Szrj const _Tp __rootx = std::sqrt(__absx); 379*38fd1498Szrj const _Tp __z = _Tp(2) * __absx * __rootx / _Tp(3); 380*38fd1498Szrj const _Tp _S_NaN = std::numeric_limits<_Tp>::quiet_NaN(); 381*38fd1498Szrj const _Tp _S_inf = std::numeric_limits<_Tp>::infinity(); 382*38fd1498Szrj 383*38fd1498Szrj if (__isnan(__x)) 384*38fd1498Szrj __Bip = __Aip = __Bi = __Ai = std::numeric_limits<_Tp>::quiet_NaN(); 385*38fd1498Szrj else if (__z == _S_inf) 386*38fd1498Szrj { 387*38fd1498Szrj __Aip = __Ai = _Tp(0); 388*38fd1498Szrj __Bip = __Bi = _S_inf; 389*38fd1498Szrj } 390*38fd1498Szrj else if (__z == -_S_inf) 391*38fd1498Szrj __Bip = __Aip = __Bi = __Ai = _Tp(0); 392*38fd1498Szrj else if (__x > _Tp(0)) 393*38fd1498Szrj { 394*38fd1498Szrj _Tp __I_nu, __Ip_nu, __K_nu, __Kp_nu; 395*38fd1498Szrj 396*38fd1498Szrj __bessel_ik(_Tp(1) / _Tp(3), __z, __I_nu, __K_nu, __Ip_nu, __Kp_nu); 397*38fd1498Szrj __Ai = __rootx * __K_nu 398*38fd1498Szrj / (__numeric_constants<_Tp>::__sqrt3() 399*38fd1498Szrj * __numeric_constants<_Tp>::__pi()); 400*38fd1498Szrj __Bi = __rootx * (__K_nu / __numeric_constants<_Tp>::__pi() 401*38fd1498Szrj + _Tp(2) * __I_nu / __numeric_constants<_Tp>::__sqrt3()); 402*38fd1498Szrj 403*38fd1498Szrj __bessel_ik(_Tp(2) / _Tp(3), __z, __I_nu, __K_nu, __Ip_nu, __Kp_nu); 404*38fd1498Szrj __Aip = -__x * __K_nu 405*38fd1498Szrj / (__numeric_constants<_Tp>::__sqrt3() 406*38fd1498Szrj * __numeric_constants<_Tp>::__pi()); 407*38fd1498Szrj __Bip = __x * (__K_nu / __numeric_constants<_Tp>::__pi() 408*38fd1498Szrj + _Tp(2) * __I_nu 409*38fd1498Szrj / __numeric_constants<_Tp>::__sqrt3()); 410*38fd1498Szrj } 411*38fd1498Szrj else if (__x < _Tp(0)) 412*38fd1498Szrj { 413*38fd1498Szrj _Tp __J_nu, __Jp_nu, __N_nu, __Np_nu; 414*38fd1498Szrj 415*38fd1498Szrj __bessel_jn(_Tp(1) / _Tp(3), __z, __J_nu, __N_nu, __Jp_nu, __Np_nu); 416*38fd1498Szrj __Ai = __rootx * (__J_nu 417*38fd1498Szrj - __N_nu / __numeric_constants<_Tp>::__sqrt3()) / _Tp(2); 418*38fd1498Szrj __Bi = -__rootx * (__N_nu 419*38fd1498Szrj + __J_nu / __numeric_constants<_Tp>::__sqrt3()) / _Tp(2); 420*38fd1498Szrj 421*38fd1498Szrj __bessel_jn(_Tp(2) / _Tp(3), __z, __J_nu, __N_nu, __Jp_nu, __Np_nu); 422*38fd1498Szrj __Aip = __absx * (__N_nu / __numeric_constants<_Tp>::__sqrt3() 423*38fd1498Szrj + __J_nu) / _Tp(2); 424*38fd1498Szrj __Bip = __absx * (__J_nu / __numeric_constants<_Tp>::__sqrt3() 425*38fd1498Szrj - __N_nu) / _Tp(2); 426*38fd1498Szrj } 427*38fd1498Szrj else 428*38fd1498Szrj { 429*38fd1498Szrj // Reference: 430*38fd1498Szrj // Abramowitz & Stegun, page 446 section 10.4.4 on Airy functions. 431*38fd1498Szrj // The number is Ai(0) = 3^{-2/3}/\Gamma(2/3). 432*38fd1498Szrj __Ai = _Tp(0.35502805388781723926L); 433*38fd1498Szrj __Bi = __Ai * __numeric_constants<_Tp>::__sqrt3(); 434*38fd1498Szrj 435*38fd1498Szrj // Reference: 436*38fd1498Szrj // Abramowitz & Stegun, page 446 section 10.4.5 on Airy functions. 437*38fd1498Szrj // The number is Ai'(0) = -3^{-1/3}/\Gamma(1/3). 438*38fd1498Szrj __Aip = -_Tp(0.25881940379280679840L); 439*38fd1498Szrj __Bip = -__Aip * __numeric_constants<_Tp>::__sqrt3(); 440*38fd1498Szrj } 441*38fd1498Szrj 442*38fd1498Szrj return; 443*38fd1498Szrj } 444*38fd1498Szrj } // namespace __detail 445*38fd1498Szrj #if ! _GLIBCXX_USE_STD_SPEC_FUNCS && defined(_GLIBCXX_TR1_CMATH) 446*38fd1498Szrj } // namespace tr1 447*38fd1498Szrj #endif 448*38fd1498Szrj 449*38fd1498Szrj _GLIBCXX_END_NAMESPACE_VERSION 450*38fd1498Szrj } 451*38fd1498Szrj 452*38fd1498Szrj #endif // _GLIBCXX_TR1_MODIFIED_BESSEL_FUNC_TCC 453