bessel_ik.hpp 14 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451
  1. // Copyright (c) 2006 Xiaogang Zhang
  2. // Use, modification and distribution are subject to the
  3. // Boost Software License, Version 1.0. (See accompanying file
  4. // LICENSE_1_0.txt or copy at http://www.boost.org/LICENSE_1_0.txt)
  5. #ifndef BOOST_MATH_BESSEL_IK_HPP
  6. #define BOOST_MATH_BESSEL_IK_HPP
  7. #ifdef _MSC_VER
  8. #pragma once
  9. #endif
  10. #include <boost/math/special_functions/round.hpp>
  11. #include <boost/math/special_functions/gamma.hpp>
  12. #include <boost/math/special_functions/sin_pi.hpp>
  13. #include <boost/math/constants/constants.hpp>
  14. #include <boost/math/policies/error_handling.hpp>
  15. #include <boost/math/tools/config.hpp>
  16. // Modified Bessel functions of the first and second kind of fractional order
  17. namespace boost { namespace math {
  18. namespace detail {
  19. template <class T, class Policy>
  20. struct cyl_bessel_i_small_z
  21. {
  22. typedef T result_type;
  23. cyl_bessel_i_small_z(T v_, T z_) : k(0), v(v_), mult(z_*z_/4)
  24. {
  25. BOOST_MATH_STD_USING
  26. term = 1;
  27. }
  28. T operator()()
  29. {
  30. T result = term;
  31. ++k;
  32. term *= mult / k;
  33. term /= k + v;
  34. return result;
  35. }
  36. private:
  37. unsigned k;
  38. T v;
  39. T term;
  40. T mult;
  41. };
  42. template <class T, class Policy>
  43. inline T bessel_i_small_z_series(T v, T x, const Policy& pol)
  44. {
  45. BOOST_MATH_STD_USING
  46. T prefix;
  47. if(v < max_factorial<T>::value)
  48. {
  49. prefix = pow(x / 2, v) / boost::math::tgamma(v + 1, pol);
  50. }
  51. else
  52. {
  53. prefix = v * log(x / 2) - boost::math::lgamma(v + 1, pol);
  54. prefix = exp(prefix);
  55. }
  56. if(prefix == 0)
  57. return prefix;
  58. cyl_bessel_i_small_z<T, Policy> s(v, x);
  59. boost::uintmax_t max_iter = policies::get_max_series_iterations<Policy>();
  60. #if BOOST_WORKAROUND(__BORLANDC__, BOOST_TESTED_AT(0x582))
  61. T zero = 0;
  62. T result = boost::math::tools::sum_series(s, boost::math::policies::get_epsilon<T, Policy>(), max_iter, zero);
  63. #else
  64. T result = boost::math::tools::sum_series(s, boost::math::policies::get_epsilon<T, Policy>(), max_iter);
  65. #endif
  66. policies::check_series_iterations<T>("boost::math::bessel_j_small_z_series<%1%>(%1%,%1%)", max_iter, pol);
  67. return prefix * result;
  68. }
  69. // Calculate K(v, x) and K(v+1, x) by method analogous to
  70. // Temme, Journal of Computational Physics, vol 21, 343 (1976)
  71. template <typename T, typename Policy>
  72. int temme_ik(T v, T x, T* K, T* K1, const Policy& pol)
  73. {
  74. T f, h, p, q, coef, sum, sum1, tolerance;
  75. T a, b, c, d, sigma, gamma1, gamma2;
  76. unsigned long k;
  77. BOOST_MATH_STD_USING
  78. using namespace boost::math::tools;
  79. using namespace boost::math::constants;
  80. // |x| <= 2, Temme series converge rapidly
  81. // |x| > 2, the larger the |x|, the slower the convergence
  82. BOOST_ASSERT(abs(x) <= 2);
  83. BOOST_ASSERT(abs(v) <= 0.5f);
  84. T gp = boost::math::tgamma1pm1(v, pol);
  85. T gm = boost::math::tgamma1pm1(-v, pol);
  86. a = log(x / 2);
  87. b = exp(v * a);
  88. sigma = -a * v;
  89. c = abs(v) < tools::epsilon<T>() ?
  90. T(1) : T(boost::math::sin_pi(v) / (v * pi<T>()));
  91. d = abs(sigma) < tools::epsilon<T>() ?
  92. T(1) : T(sinh(sigma) / sigma);
  93. gamma1 = abs(v) < tools::epsilon<T>() ?
  94. T(-euler<T>()) : T((0.5f / v) * (gp - gm) * c);
  95. gamma2 = (2 + gp + gm) * c / 2;
  96. // initial values
  97. p = (gp + 1) / (2 * b);
  98. q = (1 + gm) * b / 2;
  99. f = (cosh(sigma) * gamma1 + d * (-a) * gamma2) / c;
  100. h = p;
  101. coef = 1;
  102. sum = coef * f;
  103. sum1 = coef * h;
  104. BOOST_MATH_INSTRUMENT_VARIABLE(p);
  105. BOOST_MATH_INSTRUMENT_VARIABLE(q);
  106. BOOST_MATH_INSTRUMENT_VARIABLE(f);
  107. BOOST_MATH_INSTRUMENT_VARIABLE(sigma);
  108. BOOST_MATH_INSTRUMENT_CODE(sinh(sigma));
  109. BOOST_MATH_INSTRUMENT_VARIABLE(gamma1);
  110. BOOST_MATH_INSTRUMENT_VARIABLE(gamma2);
  111. BOOST_MATH_INSTRUMENT_VARIABLE(c);
  112. BOOST_MATH_INSTRUMENT_VARIABLE(d);
  113. BOOST_MATH_INSTRUMENT_VARIABLE(a);
  114. // series summation
  115. tolerance = tools::epsilon<T>();
  116. for (k = 1; k < policies::get_max_series_iterations<Policy>(); k++)
  117. {
  118. f = (k * f + p + q) / (k*k - v*v);
  119. p /= k - v;
  120. q /= k + v;
  121. h = p - k * f;
  122. coef *= x * x / (4 * k);
  123. sum += coef * f;
  124. sum1 += coef * h;
  125. if (abs(coef * f) < abs(sum) * tolerance)
  126. {
  127. break;
  128. }
  129. }
  130. policies::check_series_iterations<T>("boost::math::bessel_ik<%1%>(%1%,%1%) in temme_ik", k, pol);
  131. *K = sum;
  132. *K1 = 2 * sum1 / x;
  133. return 0;
  134. }
  135. // Evaluate continued fraction fv = I_(v+1) / I_v, derived from
  136. // Abramowitz and Stegun, Handbook of Mathematical Functions, 1972, 9.1.73
  137. template <typename T, typename Policy>
  138. int CF1_ik(T v, T x, T* fv, const Policy& pol)
  139. {
  140. T C, D, f, a, b, delta, tiny, tolerance;
  141. unsigned long k;
  142. BOOST_MATH_STD_USING
  143. // |x| <= |v|, CF1_ik converges rapidly
  144. // |x| > |v|, CF1_ik needs O(|x|) iterations to converge
  145. // modified Lentz's method, see
  146. // Lentz, Applied Optics, vol 15, 668 (1976)
  147. tolerance = 2 * tools::epsilon<T>();
  148. BOOST_MATH_INSTRUMENT_VARIABLE(tolerance);
  149. tiny = sqrt(tools::min_value<T>());
  150. BOOST_MATH_INSTRUMENT_VARIABLE(tiny);
  151. C = f = tiny; // b0 = 0, replace with tiny
  152. D = 0;
  153. for (k = 1; k < policies::get_max_series_iterations<Policy>(); k++)
  154. {
  155. a = 1;
  156. b = 2 * (v + k) / x;
  157. C = b + a / C;
  158. D = b + a * D;
  159. if (C == 0) { C = tiny; }
  160. if (D == 0) { D = tiny; }
  161. D = 1 / D;
  162. delta = C * D;
  163. f *= delta;
  164. BOOST_MATH_INSTRUMENT_VARIABLE(delta-1);
  165. if (abs(delta - 1) <= tolerance)
  166. {
  167. break;
  168. }
  169. }
  170. BOOST_MATH_INSTRUMENT_VARIABLE(k);
  171. policies::check_series_iterations<T>("boost::math::bessel_ik<%1%>(%1%,%1%) in CF1_ik", k, pol);
  172. *fv = f;
  173. return 0;
  174. }
  175. // Calculate K(v, x) and K(v+1, x) by evaluating continued fraction
  176. // z1 / z0 = U(v+1.5, 2v+1, 2x) / U(v+0.5, 2v+1, 2x), see
  177. // Thompson and Barnett, Computer Physics Communications, vol 47, 245 (1987)
  178. template <typename T, typename Policy>
  179. int CF2_ik(T v, T x, T* Kv, T* Kv1, const Policy& pol)
  180. {
  181. BOOST_MATH_STD_USING
  182. using namespace boost::math::constants;
  183. T S, C, Q, D, f, a, b, q, delta, tolerance, current, prev;
  184. unsigned long k;
  185. // |x| >= |v|, CF2_ik converges rapidly
  186. // |x| -> 0, CF2_ik fails to converge
  187. BOOST_ASSERT(abs(x) > 1);
  188. // Steed's algorithm, see Thompson and Barnett,
  189. // Journal of Computational Physics, vol 64, 490 (1986)
  190. tolerance = tools::epsilon<T>();
  191. a = v * v - 0.25f;
  192. b = 2 * (x + 1); // b1
  193. D = 1 / b; // D1 = 1 / b1
  194. f = delta = D; // f1 = delta1 = D1, coincidence
  195. prev = 0; // q0
  196. current = 1; // q1
  197. Q = C = -a; // Q1 = C1 because q1 = 1
  198. S = 1 + Q * delta; // S1
  199. BOOST_MATH_INSTRUMENT_VARIABLE(tolerance);
  200. BOOST_MATH_INSTRUMENT_VARIABLE(a);
  201. BOOST_MATH_INSTRUMENT_VARIABLE(b);
  202. BOOST_MATH_INSTRUMENT_VARIABLE(D);
  203. BOOST_MATH_INSTRUMENT_VARIABLE(f);
  204. for (k = 2; k < policies::get_max_series_iterations<Policy>(); k++) // starting from 2
  205. {
  206. // continued fraction f = z1 / z0
  207. a -= 2 * (k - 1);
  208. b += 2;
  209. D = 1 / (b + a * D);
  210. delta *= b * D - 1;
  211. f += delta;
  212. // series summation S = 1 + \sum_{n=1}^{\infty} C_n * z_n / z_0
  213. q = (prev - (b - 2) * current) / a;
  214. prev = current;
  215. current = q; // forward recurrence for q
  216. C *= -a / k;
  217. Q += C * q;
  218. S += Q * delta;
  219. //
  220. // Under some circumstances q can grow very small and C very
  221. // large, leading to under/overflow. This is particularly an
  222. // issue for types which have many digits precision but a narrow
  223. // exponent range. A typical example being a "double double" type.
  224. // To avoid this situation we can normalise q (and related prev/current)
  225. // and C. All other variables remain unchanged in value. A typical
  226. // test case occurs when x is close to 2, for example cyl_bessel_k(9.125, 2.125).
  227. //
  228. if(q < tools::epsilon<T>())
  229. {
  230. C *= q;
  231. prev /= q;
  232. current /= q;
  233. q = 1;
  234. }
  235. // S converges slower than f
  236. BOOST_MATH_INSTRUMENT_VARIABLE(Q * delta);
  237. BOOST_MATH_INSTRUMENT_VARIABLE(abs(S) * tolerance);
  238. BOOST_MATH_INSTRUMENT_VARIABLE(S);
  239. if (abs(Q * delta) < abs(S) * tolerance)
  240. {
  241. break;
  242. }
  243. }
  244. policies::check_series_iterations<T>("boost::math::bessel_ik<%1%>(%1%,%1%) in CF2_ik", k, pol);
  245. if(x >= tools::log_max_value<T>())
  246. *Kv = exp(0.5f * log(pi<T>() / (2 * x)) - x - log(S));
  247. else
  248. *Kv = sqrt(pi<T>() / (2 * x)) * exp(-x) / S;
  249. *Kv1 = *Kv * (0.5f + v + x + (v * v - 0.25f) * f) / x;
  250. BOOST_MATH_INSTRUMENT_VARIABLE(*Kv);
  251. BOOST_MATH_INSTRUMENT_VARIABLE(*Kv1);
  252. return 0;
  253. }
  254. enum{
  255. need_i = 1,
  256. need_k = 2
  257. };
  258. // Compute I(v, x) and K(v, x) simultaneously by Temme's method, see
  259. // Temme, Journal of Computational Physics, vol 19, 324 (1975)
  260. template <typename T, typename Policy>
  261. int bessel_ik(T v, T x, T* I, T* K, int kind, const Policy& pol)
  262. {
  263. // Kv1 = K_(v+1), fv = I_(v+1) / I_v
  264. // Ku1 = K_(u+1), fu = I_(u+1) / I_u
  265. T u, Iv, Kv, Kv1, Ku, Ku1, fv;
  266. T W, current, prev, next;
  267. bool reflect = false;
  268. unsigned n, k;
  269. int org_kind = kind;
  270. BOOST_MATH_INSTRUMENT_VARIABLE(v);
  271. BOOST_MATH_INSTRUMENT_VARIABLE(x);
  272. BOOST_MATH_INSTRUMENT_VARIABLE(kind);
  273. BOOST_MATH_STD_USING
  274. using namespace boost::math::tools;
  275. using namespace boost::math::constants;
  276. static const char* function = "boost::math::bessel_ik<%1%>(%1%,%1%)";
  277. if (v < 0)
  278. {
  279. reflect = true;
  280. v = -v; // v is non-negative from here
  281. kind |= need_k;
  282. }
  283. n = iround(v, pol);
  284. u = v - n; // -1/2 <= u < 1/2
  285. BOOST_MATH_INSTRUMENT_VARIABLE(n);
  286. BOOST_MATH_INSTRUMENT_VARIABLE(u);
  287. if (x < 0)
  288. {
  289. *I = *K = policies::raise_domain_error<T>(function,
  290. "Got x = %1% but real argument x must be non-negative, complex number result not supported.", x, pol);
  291. return 1;
  292. }
  293. if (x == 0)
  294. {
  295. Iv = (v == 0) ? static_cast<T>(1) : static_cast<T>(0);
  296. if(kind & need_k)
  297. {
  298. Kv = policies::raise_overflow_error<T>(function, 0, pol);
  299. }
  300. else
  301. {
  302. Kv = std::numeric_limits<T>::quiet_NaN(); // any value will do
  303. }
  304. if(reflect && (kind & need_i))
  305. {
  306. T z = (u + n % 2);
  307. Iv = boost::math::sin_pi(z, pol) == 0 ?
  308. Iv :
  309. policies::raise_overflow_error<T>(function, 0, pol); // reflection formula
  310. }
  311. *I = Iv;
  312. *K = Kv;
  313. return 0;
  314. }
  315. // x is positive until reflection
  316. W = 1 / x; // Wronskian
  317. if (x <= 2) // x in (0, 2]
  318. {
  319. temme_ik(u, x, &Ku, &Ku1, pol); // Temme series
  320. }
  321. else // x in (2, \infty)
  322. {
  323. CF2_ik(u, x, &Ku, &Ku1, pol); // continued fraction CF2_ik
  324. }
  325. BOOST_MATH_INSTRUMENT_VARIABLE(Ku);
  326. BOOST_MATH_INSTRUMENT_VARIABLE(Ku1);
  327. prev = Ku;
  328. current = Ku1;
  329. T scale = 1;
  330. T scale_sign = 1;
  331. for (k = 1; k <= n; k++) // forward recurrence for K
  332. {
  333. T fact = 2 * (u + k) / x;
  334. if((tools::max_value<T>() - fabs(prev)) / fact < fabs(current))
  335. {
  336. prev /= current;
  337. scale /= current;
  338. scale_sign *= boost::math::sign(current);
  339. current = 1;
  340. }
  341. next = fact * current + prev;
  342. prev = current;
  343. current = next;
  344. }
  345. Kv = prev;
  346. Kv1 = current;
  347. BOOST_MATH_INSTRUMENT_VARIABLE(Kv);
  348. BOOST_MATH_INSTRUMENT_VARIABLE(Kv1);
  349. if(kind & need_i)
  350. {
  351. T lim = (4 * v * v + 10) / (8 * x);
  352. lim *= lim;
  353. lim *= lim;
  354. lim /= 24;
  355. if((lim < tools::epsilon<T>() * 10) && (x > 100))
  356. {
  357. // x is huge compared to v, CF1 may be very slow
  358. // to converge so use asymptotic expansion for large
  359. // x case instead. Note that the asymptotic expansion
  360. // isn't very accurate - so it's deliberately very hard
  361. // to get here - probably we're going to overflow:
  362. Iv = asymptotic_bessel_i_large_x(v, x, pol);
  363. }
  364. else if((v > 0) && (x / v < 0.25))
  365. {
  366. Iv = bessel_i_small_z_series(v, x, pol);
  367. }
  368. else
  369. {
  370. CF1_ik(v, x, &fv, pol); // continued fraction CF1_ik
  371. Iv = scale * W / (Kv * fv + Kv1); // Wronskian relation
  372. }
  373. }
  374. else
  375. Iv = std::numeric_limits<T>::quiet_NaN(); // any value will do
  376. if (reflect)
  377. {
  378. T z = (u + n % 2);
  379. T fact = (2 / pi<T>()) * (boost::math::sin_pi(z) * Kv);
  380. if(fact == 0)
  381. *I = Iv;
  382. else if(tools::max_value<T>() * scale < fact)
  383. *I = (org_kind & need_i) ? T(sign(fact) * scale_sign * policies::raise_overflow_error<T>(function, 0, pol)) : T(0);
  384. else
  385. *I = Iv + fact / scale; // reflection formula
  386. }
  387. else
  388. {
  389. *I = Iv;
  390. }
  391. if(tools::max_value<T>() * scale < Kv)
  392. *K = (org_kind & need_k) ? T(sign(Kv) * scale_sign * policies::raise_overflow_error<T>(function, 0, pol)) : T(0);
  393. else
  394. *K = Kv / scale;
  395. BOOST_MATH_INSTRUMENT_VARIABLE(*I);
  396. BOOST_MATH_INSTRUMENT_VARIABLE(*K);
  397. return 0;
  398. }
  399. }}} // namespaces
  400. #endif // BOOST_MATH_BESSEL_IK_HPP