hypergeometric_1F1_large_abz.hpp 22 KB

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