hypergeometric_1F1_large_abz.hpp 22 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484
  1. ///////////////////////////////////////////////////////////////////////////////
  2. // Copyright 2018 John Maddock
  3. // Distributed under the Boost
  4. // 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_HYPERGEOMETRIC_1F1_LARGE_ABZ_HPP_
  7. #define BOOST_HYPERGEOMETRIC_1F1_LARGE_ABZ_HPP_
  8. #include <boost/math/special_functions/detail/hypergeometric_1F1_bessel.hpp>
  9. #include <boost/math/special_functions/detail/hypergeometric_series.hpp>
  10. #include <boost/math/special_functions/gamma.hpp>
  11. namespace boost { namespace math { namespace detail {
  12. template <class T>
  13. inline bool is_negative_integer(const T& x)
  14. {
  15. using std::floor;
  16. return (x <= 0) && (floor(x) == x);
  17. }
  18. template <class T, class Policy>
  19. struct hypergeometric_1F1_igamma_series
  20. {
  21. enum{ cache_size = 64 };
  22. typedef T result_type;
  23. hypergeometric_1F1_igamma_series(const T& alpha, const T& delta, const T& x, const Policy& pol)
  24. : delta_poch(-delta), alpha_poch(alpha), x(x), k(0), cache_offset(0), pol(pol)
  25. {
  26. BOOST_MATH_STD_USING
  27. T log_term = log(x) * -alpha;
  28. log_scaling = itrunc(log_term - 3 - boost::math::tools::log_min_value<T>() / 50);
  29. term = exp(log_term - log_scaling);
  30. refill_cache();
  31. }
  32. T operator()()
  33. {
  34. if (k - cache_offset >= cache_size)
  35. {
  36. cache_offset += cache_size;
  37. refill_cache();
  38. }
  39. T result = term * gamma_cache[k - cache_offset];
  40. term *= delta_poch * alpha_poch / (++k * x);
  41. delta_poch += 1;
  42. alpha_poch += 1;
  43. return result;
  44. }
  45. void refill_cache()
  46. {
  47. typedef typename lanczos::lanczos<T, Policy>::type lanczos_type;
  48. gamma_cache[cache_size - 1] = boost::math::gamma_p(alpha_poch + ((int)cache_size - 1), x, pol);
  49. for (int i = cache_size - 1; i > 0; --i)
  50. {
  51. gamma_cache[i - 1] = gamma_cache[i] >= 1 ? T(1) : T(gamma_cache[i] + regularised_gamma_prefix(T(alpha_poch + (i - 1)), x, pol, lanczos_type()) / (alpha_poch + (i - 1)));
  52. }
  53. }
  54. T delta_poch, alpha_poch, x, term;
  55. T gamma_cache[cache_size];
  56. int k;
  57. int log_scaling;
  58. int cache_offset;
  59. Policy pol;
  60. };
  61. template <class T, class Policy>
  62. T hypergeometric_1F1_igamma(const T& a, const T& b, const T& x, const T& b_minus_a, const Policy& pol, int& log_scaling)
  63. {
  64. BOOST_MATH_STD_USING
  65. if (b_minus_a == 0)
  66. {
  67. // special case: M(a,a,z) == exp(z)
  68. int scale = itrunc(x, pol);
  69. log_scaling += scale;
  70. return exp(x - scale);
  71. }
  72. hypergeometric_1F1_igamma_series<T, Policy> s(b_minus_a, a - 1, x, pol);
  73. log_scaling += s.log_scaling;
  74. boost::uintmax_t max_iter = boost::math::policies::get_max_series_iterations<Policy>();
  75. T result = boost::math::tools::sum_series(s, boost::math::policies::get_epsilon<T, Policy>(), max_iter);
  76. boost::math::policies::check_series_iterations<T>("boost::math::tgamma<%1%>(%1%,%1%)", max_iter, pol);
  77. T log_prefix = x + boost::math::lgamma(b, pol) - boost::math::lgamma(a, pol);
  78. int scale = itrunc(log_prefix);
  79. log_scaling += scale;
  80. return result * exp(log_prefix - scale);
  81. }
  82. template <class T, class Policy>
  83. T hypergeometric_1F1_shift_on_a(T h, const T& a_local, const T& b_local, const T& x, int a_shift, const Policy& pol, int& log_scaling)
  84. {
  85. BOOST_MATH_STD_USING
  86. T a = a_local + a_shift;
  87. if (a_shift == 0)
  88. return h;
  89. else if (a_shift > 0)
  90. {
  91. //
  92. // Forward recursion on a is stable as long as 2a-b+z > 0.
  93. // If 2a-b+z < 0 then backwards recursion is stable even though
  94. // the function may be strictly increasing with a. Potentially
  95. // we may need to split the recurrnce in 2 sections - one using
  96. // forward recursion, and one backwards.
  97. //
  98. // We will get the next seed value from the ratio
  99. // on b as that's stable and quick to compute.
  100. //
  101. T crossover_a = (b_local - x) / 2;
  102. int crossover_shift = itrunc(crossover_a - a_local);
  103. if (crossover_shift > 1)
  104. {
  105. //
  106. // Forwards recursion will start off unstable, but may switch to the stable direction later.
  107. // Start in the middle and go in both directions:
  108. //
  109. if (crossover_shift > a_shift)
  110. crossover_shift = a_shift;
  111. crossover_a = a_local + crossover_shift;
  112. boost::math::detail::hypergeometric_1F1_recurrence_b_coefficients<T> b_coef(crossover_a, b_local, x);
  113. boost::uintmax_t max_iter = boost::math::policies::get_max_series_iterations<Policy>();
  114. T b_ratio = boost::math::tools::function_ratio_from_backwards_recurrence(b_coef, boost::math::policies::get_epsilon<T, Policy>(), max_iter);
  115. boost::math::policies::check_series_iterations<T>("boost::math::hypergeometric_1F1_large_abz<%1%>(%1%,%1%,%1%)", max_iter, pol);
  116. //
  117. // Convert to a ratio:
  118. // (1+a-b)M(a, b, z) - aM(a+1, b, z) + (b-1)M(a, b-1, z) = 0
  119. //
  120. // hence: M(a+1,b,z) = ((1+a-b) / a) M(a,b,z) + ((b-1) / a) M(a,b,z)/b_ratio
  121. //
  122. T first = 1;
  123. T second = ((1 + crossover_a - b_local) / crossover_a) + ((b_local - 1) / crossover_a) / b_ratio;
  124. //
  125. // Recurse down to a_local, compare values and re-nomralise first and second:
  126. //
  127. boost::math::detail::hypergeometric_1F1_recurrence_a_coefficients<T> a_coef(crossover_a, b_local, x);
  128. int backwards_scale = 0;
  129. T comparitor = boost::math::tools::apply_recurrence_relation_backward(a_coef, crossover_shift, second, first, &backwards_scale);
  130. log_scaling -= backwards_scale;
  131. if ((h < 1) && (tools::max_value<T>() * h > comparitor))
  132. {
  133. // Need to rescale!
  134. int scale = itrunc(log(h), pol) + 1;
  135. h *= exp(T(-scale));
  136. log_scaling += scale;
  137. }
  138. comparitor /= h;
  139. first /= comparitor;
  140. second /= comparitor;
  141. //
  142. // Now we can recurse forwards for the rest of the range:
  143. //
  144. if (crossover_shift < a_shift)
  145. {
  146. boost::math::detail::hypergeometric_1F1_recurrence_a_coefficients<T> a_coef_2(crossover_a + 1, b_local, x);
  147. h = boost::math::tools::apply_recurrence_relation_forward(a_coef_2, a_shift - crossover_shift - 1, first, second, &log_scaling);
  148. }
  149. else
  150. h = first;
  151. }
  152. else
  153. {
  154. //
  155. // Regular case where forwards iteration is stable right from the start:
  156. //
  157. boost::math::detail::hypergeometric_1F1_recurrence_b_coefficients<T> b_coef(a_local, b_local, x);
  158. boost::uintmax_t max_iter = boost::math::policies::get_max_series_iterations<Policy>();
  159. T b_ratio = boost::math::tools::function_ratio_from_backwards_recurrence(b_coef, boost::math::policies::get_epsilon<T, Policy>(), max_iter);
  160. boost::math::policies::check_series_iterations<T>("boost::math::hypergeometric_1F1_large_abz<%1%>(%1%,%1%,%1%)", max_iter, pol);
  161. //
  162. // Convert to a ratio:
  163. // (1+a-b)M(a, b, z) - aM(a+1, b, z) + (b-1)M(a, b-1, z) = 0
  164. //
  165. // hence: M(a+1,b,z) = ((1+a-b) / a) M(a,b,z) + ((b-1) / a) M(a,b,z)/b_ratio
  166. //
  167. T second = ((1 + a_local - b_local) / a_local) * h + ((b_local - 1) / a_local) * h / b_ratio;
  168. boost::math::detail::hypergeometric_1F1_recurrence_a_coefficients<T> a_coef(a_local + 1, b_local, x);
  169. h = boost::math::tools::apply_recurrence_relation_forward(a_coef, --a_shift, h, second, &log_scaling);
  170. }
  171. }
  172. else
  173. {
  174. //
  175. // We've calculated h for a larger value of a than we want, and need to recurse down.
  176. // However, only forward iteration is stable, so calculate the ratio, compare values,
  177. // and normalise. Note that we calculate the ratio on b and convert to a since the
  178. // direction is the minimal solution for N->+INF.
  179. //
  180. // IMPORTANT: this is only currently called for a > b and therefore forwards iteration
  181. // is the only stable direction as we will only iterate down until a ~ b, but we
  182. // will check this with an assert:
  183. //
  184. BOOST_ASSERT(2 * a - b_local + x > 0);
  185. boost::math::detail::hypergeometric_1F1_recurrence_b_coefficients<T> b_coef(a, b_local, x);
  186. boost::uintmax_t max_iter = boost::math::policies::get_max_series_iterations<Policy>();
  187. T b_ratio = boost::math::tools::function_ratio_from_backwards_recurrence(b_coef, boost::math::policies::get_epsilon<T, Policy>(), max_iter);
  188. boost::math::policies::check_series_iterations<T>("boost::math::hypergeometric_1F1_large_abz<%1%>(%1%,%1%,%1%)", max_iter, pol);
  189. //
  190. // Convert to a ratio:
  191. // (1+a-b)M(a, b, z) - aM(a+1, b, z) + (b-1)M(a, b-1, z) = 0
  192. //
  193. // hence: M(a+1,b,z) = (1+a-b) / a M(a,b,z) + (b-1) / a M(a,b,z)/ (M(a,b,z)/M(a,b-1,z))
  194. //
  195. T first = 1; // arbitrary value;
  196. T second = ((1 + a - b_local) / a) + ((b_local - 1) / a) * (1 / b_ratio);
  197. if (a_shift == -1)
  198. h = h / second;
  199. else
  200. {
  201. boost::math::detail::hypergeometric_1F1_recurrence_a_coefficients<T> a_coef(a + 1, b_local, x);
  202. T comparitor = boost::math::tools::apply_recurrence_relation_forward(a_coef, -(a_shift + 1), first, second);
  203. if (boost::math::tools::min_value<T>() * comparitor > h)
  204. {
  205. // Ooops, need to rescale h:
  206. int rescale = itrunc(log(fabs(h)));
  207. T scale = exp(T(-rescale));
  208. h *= scale;
  209. log_scaling += rescale;
  210. }
  211. h = h / comparitor;
  212. }
  213. }
  214. return h;
  215. }
  216. template <class T, class Policy>
  217. T hypergeometric_1F1_shift_on_b(T h, const T& a, const T& b_local, const T& x, int b_shift, const Policy& pol, int& log_scaling)
  218. {
  219. BOOST_MATH_STD_USING
  220. T b = b_local + b_shift;
  221. if (b_shift == 0)
  222. return h;
  223. else if (b_shift > 0)
  224. {
  225. //
  226. // We get here for b_shift > 0 when b > z. We can't use forward recursion on b - it's unstable,
  227. // so grab the ratio and work backwards to b - b_shift and normalise.
  228. //
  229. boost::math::detail::hypergeometric_1F1_recurrence_b_coefficients<T> b_coef(a, b, x);
  230. boost::uintmax_t max_iter = boost::math::policies::get_max_series_iterations<Policy>();
  231. T first = 1; // arbitrary value;
  232. T second = 1 / boost::math::tools::function_ratio_from_backwards_recurrence(b_coef, boost::math::policies::get_epsilon<T, Policy>(), max_iter);
  233. boost::math::policies::check_series_iterations<T>("boost::math::hypergeometric_1F1_large_abz<%1%>(%1%,%1%,%1%)", max_iter, pol);
  234. if (b_shift == 1)
  235. h = h / second;
  236. else
  237. {
  238. //
  239. // Reset coefficients and recurse:
  240. //
  241. boost::math::detail::hypergeometric_1F1_recurrence_b_coefficients<T> b_coef_2(a, b - 1, x);
  242. int local_scale = 0;
  243. T comparitor = boost::math::tools::apply_recurrence_relation_backward(b_coef_2, --b_shift, first, second, &local_scale);
  244. log_scaling -= local_scale;
  245. if (boost::math::tools::min_value<T>() * comparitor > h)
  246. {
  247. // Ooops, need to rescale h:
  248. int rescale = itrunc(log(fabs(h)));
  249. T scale = exp(T(-rescale));
  250. h *= scale;
  251. log_scaling += rescale;
  252. }
  253. h = h / comparitor;
  254. }
  255. }
  256. else
  257. {
  258. T second;
  259. if (a == b_local)
  260. {
  261. // recurrence is trivial for a == b and method of ratios fails as the c-term goes to zero:
  262. second = -b_local * (1 - b_local - x) * h / (b_local * (b_local - 1));
  263. }
  264. else
  265. {
  266. BOOST_ASSERT(!is_negative_integer(b - a));
  267. boost::math::detail::hypergeometric_1F1_recurrence_b_coefficients<T> b_coef(a, b_local, x);
  268. boost::uintmax_t max_iter = boost::math::policies::get_max_series_iterations<Policy>();
  269. second = h / boost::math::tools::function_ratio_from_backwards_recurrence(b_coef, boost::math::policies::get_epsilon<T, Policy>(), max_iter);
  270. boost::math::policies::check_series_iterations<T>("boost::math::hypergeometric_1F1_large_abz<%1%>(%1%,%1%,%1%)", max_iter, pol);
  271. }
  272. if (b_shift == -1)
  273. h = second;
  274. else
  275. {
  276. boost::math::detail::hypergeometric_1F1_recurrence_b_coefficients<T> b_coef_2(a, b_local - 1, x);
  277. h = boost::math::tools::apply_recurrence_relation_backward(b_coef_2, -(++b_shift), h, second, &log_scaling);
  278. }
  279. }
  280. return h;
  281. }
  282. template <class T, class Policy>
  283. T hypergeometric_1F1_large_igamma(const T& a, const T& b, const T& x, const T& b_minus_a, const Policy& pol, int& log_scaling)
  284. {
  285. BOOST_MATH_STD_USING
  286. //
  287. // We need a < b < z in order to ensure there's at least a chance of convergence,
  288. // we can use recurrence relations to shift forwards on a+b or just a to achieve this,
  289. // for decent accuracy, try to keep 2b - 1 < a < 2b < z
  290. //
  291. int b_shift = b * 2 < x ? 0 : itrunc(b - x / 2);
  292. int a_shift = a > b - b_shift ? -itrunc(b - b_shift - a - 1) : -itrunc(b - b_shift - a);
  293. if (a_shift < 0)
  294. {
  295. // Might as well do all the shifting on b as scale a downwards:
  296. b_shift -= a_shift;
  297. a_shift = 0;
  298. }
  299. T a_local = a - a_shift;
  300. T b_local = b - b_shift;
  301. T b_minus_a_local = (a_shift == 0) && (b_shift == 0) ? b_minus_a : b_local - a_local;
  302. int local_scaling = 0;
  303. T h = hypergeometric_1F1_igamma(a_local, b_local, x, b_minus_a_local, pol, local_scaling);
  304. log_scaling += local_scaling;
  305. //
  306. // Apply shifts on a and b as required:
  307. //
  308. h = hypergeometric_1F1_shift_on_a(h, a_local, b_local, x, a_shift, pol, log_scaling);
  309. h = hypergeometric_1F1_shift_on_b(h, a, b_local, x, b_shift, pol, log_scaling);
  310. return h;
  311. }
  312. template <class T, class Policy>
  313. T hypergeometric_1F1_large_series(const T& a, const T& b, const T& z, const Policy& pol, int& log_scaling)
  314. {
  315. BOOST_MATH_STD_USING
  316. //
  317. // We make a small, and b > z:
  318. //
  319. int a_shift(0), b_shift(0);
  320. if (a * z > b)
  321. {
  322. a_shift = itrunc(a) - 5;
  323. b_shift = b < z ? itrunc(b - z - 1) : 0;
  324. }
  325. //
  326. // If a_shift is trivially small, there's really not much point in losing
  327. // accuracy to save a couple of iterations:
  328. //
  329. if (a_shift < 5)
  330. a_shift = 0;
  331. T a_local = a - a_shift;
  332. T b_local = b - b_shift;
  333. T h = boost::math::detail::hypergeometric_1F1_generic_series(a_local, b_local, z, pol, log_scaling, "hypergeometric_1F1_large_series<%1%>(a,b,z)");
  334. //
  335. // Apply shifts on a and b as required:
  336. //
  337. if (a_shift && (a_local == 0))
  338. {
  339. //
  340. // Shifting on a via method of ratios in hypergeometric_1F1_shift_on_a fails when
  341. // a_local == 0. However, the value of h calculated was trivial (unity), so
  342. // calculate a second 1F1 for a == 1 and recurse as normal:
  343. //
  344. int scale = 0;
  345. T h2 = boost::math::detail::hypergeometric_1F1_generic_series(T(a_local + 1), b_local, z, pol, scale, "hypergeometric_1F1_large_series<%1%>(a,b,z)");
  346. if (scale != log_scaling)
  347. {
  348. h2 *= exp(T(scale - log_scaling));
  349. }
  350. boost::math::detail::hypergeometric_1F1_recurrence_a_coefficients<T> coef(a_local + 1, b_local, z);
  351. h = boost::math::tools::apply_recurrence_relation_forward(coef, a_shift - 1, h, h2, &log_scaling);
  352. h = hypergeometric_1F1_shift_on_b(h, a, b_local, z, b_shift, pol, log_scaling);
  353. }
  354. else
  355. {
  356. h = hypergeometric_1F1_shift_on_a(h, a_local, b_local, z, a_shift, pol, log_scaling);
  357. h = hypergeometric_1F1_shift_on_b(h, a, b_local, z, b_shift, pol, log_scaling);
  358. }
  359. return h;
  360. }
  361. template <class T, class Policy>
  362. T hypergeometric_1F1_large_13_3_6_series(const T& a, const T& b, const T& z, const Policy& pol, int& log_scaling)
  363. {
  364. BOOST_MATH_STD_USING
  365. //
  366. // A&S 13.3.6 is good only when a ~ b, but isn't too fussy on the size of z.
  367. // So shift b to match a (b shifting seems to be more stable via method of ratios).
  368. //
  369. int b_shift = itrunc(b - a);
  370. T b_local = b - b_shift;
  371. T h = boost::math::detail::hypergeometric_1F1_AS_13_3_6(a, b_local, z, T(b_local - a), pol, log_scaling);
  372. return hypergeometric_1F1_shift_on_b(h, a, b_local, z, b_shift, pol, log_scaling);
  373. }
  374. template <class T, class Policy>
  375. T hypergeometric_1F1_large_abz(const T& a, const T& b, const T& z, const Policy& pol, int& log_scaling)
  376. {
  377. BOOST_MATH_STD_USING
  378. //
  379. // This is the selection logic to pick the "best" method.
  380. // We have a,b,z >> 0 and need to comute the approximate cost of each method
  381. // and then select whichever wins out.
  382. //
  383. enum method
  384. {
  385. method_series = 0,
  386. method_shifted_series,
  387. method_gamma,
  388. method_bessel
  389. };
  390. //
  391. // Cost of direct series, is the approx number of terms required until we hit convergence:
  392. //
  393. T current_cost = (sqrt(16 * z * (3 * a + z) + 9 * b * b - 24 * b * z) - 3 * b + 4 * z) / 6;
  394. method current_method = method_series;
  395. //
  396. // Cost of shifted series, is the number of recurrences required to move to a zone where
  397. // the series is convergent right from the start.
  398. // Note that recurrence relations fail for very small b, and too many recurrences on a
  399. // will completely destroy all our digits.
  400. // Also note that the method fails when b-a is a negative integer unless b is already
  401. // larger than z and thus does not need shifting.
  402. //
  403. T cost = a + ((b < z) ? T(z - b) : T(0));
  404. if((b > 1) && (cost < current_cost) && ((b > z) || !is_negative_integer(b-a)))
  405. {
  406. current_method = method_shifted_series;
  407. current_cost = cost;
  408. }
  409. //
  410. // Cost for gamma function method is the number of recurrences required to move it
  411. // into a convergent zone, note that recurrence relations fail for very small b.
  412. // Also add on a fudge factor to account for the fact that this method is both
  413. // more expensive to compute (requires gamma functions), and less accurate than the
  414. // methods above:
  415. //
  416. T b_shift = fabs(b * 2 < z ? T(0) : T(b - z / 2));
  417. T a_shift = fabs(a > b - b_shift ? T(-(b - b_shift - a - 1)) : T(-(b - b_shift - a)));
  418. cost = 1000 + b_shift + a_shift;
  419. if((b > 1) && (cost <= current_cost))
  420. {
  421. current_method = method_gamma;
  422. current_cost = cost;
  423. }
  424. //
  425. // Cost for bessel approximation is the number of recurrences required to make a ~ b,
  426. // Note that recurrence relations fail for very small b. We also have issue with large
  427. // z: either overflow/numeric instability or else the series goes divergent. We seem to be
  428. // OK for z smaller than log_max_value<Quad> though, maybe we can stretch a little further
  429. // but that's not clear...
  430. // Also need to add on a fudge factor to the cost to account for the fact that we need
  431. // to calculate the Bessel functions, this is not quite as high as the gamma function
  432. // method above as this is generally more accurate and so prefered if the methods are close:
  433. //
  434. cost = 50 + fabs(b - a);
  435. if((b > 1) && (cost <= current_cost) && (z < tools::log_max_value<T>()) && (z < 11356) && (b - a != 0.5f))
  436. {
  437. current_method = method_bessel;
  438. current_cost = cost;
  439. }
  440. switch (current_method)
  441. {
  442. case method_series:
  443. return detail::hypergeometric_1F1_generic_series(a, b, z, pol, log_scaling, "hypergeometric_1f1_large_abz<%1%>(a,b,z)");
  444. case method_shifted_series:
  445. return detail::hypergeometric_1F1_large_series(a, b, z, pol, log_scaling);
  446. case method_gamma:
  447. return detail::hypergeometric_1F1_large_igamma(a, b, z, T(b - a), pol, log_scaling);
  448. case method_bessel:
  449. return detail::hypergeometric_1F1_large_13_3_6_series(a, b, z, pol, log_scaling);
  450. }
  451. return 0; // We don't get here!
  452. }
  453. } } } // namespaces
  454. #endif // BOOST_HYPERGEOMETRIC_1F1_LARGE_ABZ_HPP_