bessel_ik.hpp 14 KB

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