hypergeometric_1F1.hpp 33 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823
  1. ///////////////////////////////////////////////////////////////////////////////
  2. // Copyright 2014 Anton Bikineev
  3. // Copyright 2014 Christopher Kormanyos
  4. // Copyright 2014 John Maddock
  5. // Copyright 2014 Paul Bristow
  6. // Distributed under the Boost
  7. // Software License, Version 1.0. (See accompanying file
  8. // LICENSE_1_0.txt or copy at http://www.boost.org/LICENSE_1_0.txt)
  9. #ifndef BOOST_MATH_HYPERGEOMETRIC_1F1_HPP
  10. #define BOOST_MATH_HYPERGEOMETRIC_1F1_HPP
  11. #include <boost/math/tools/config.hpp>
  12. #include <boost/math/policies/policy.hpp>
  13. #include <boost/math/policies/error_handling.hpp>
  14. #include <boost/math/special_functions/detail/hypergeometric_series.hpp>
  15. #include <boost/math/special_functions/detail/hypergeometric_asym.hpp>
  16. #include <boost/math/special_functions/detail/hypergeometric_rational.hpp>
  17. #include <boost/math/special_functions/detail/hypergeometric_1F1_recurrence.hpp>
  18. #include <boost/math/special_functions/detail/hypergeometric_1F1_by_ratios.hpp>
  19. #include <boost/math/special_functions/detail/hypergeometric_pade.hpp>
  20. #include <boost/math/special_functions/detail/hypergeometric_1F1_bessel.hpp>
  21. #include <boost/math/special_functions/detail/hypergeometric_1F1_scaled_series.hpp>
  22. #include <boost/math/special_functions/detail/hypergeometric_pFq_checked_series.hpp>
  23. #include <boost/math/special_functions/detail/hypergeometric_1F1_addition_theorems_on_z.hpp>
  24. #include <boost/math/special_functions/detail/hypergeometric_1F1_large_abz.hpp>
  25. #include <boost/math/special_functions/detail/hypergeometric_1F1_small_a_negative_b_by_ratio.hpp>
  26. #include <boost/math/special_functions/detail/hypergeometric_1F1_negative_b_regions.hpp>
  27. namespace boost { namespace math { namespace detail {
  28. // check when 1F1 series can't decay to polynom
  29. template <class T>
  30. inline bool check_hypergeometric_1F1_parameters(const T& a, const T& b)
  31. {
  32. BOOST_MATH_STD_USING
  33. if ((b <= 0) && (b == floor(b)))
  34. {
  35. if ((a >= 0) || (a < b) || (a != floor(a)))
  36. return false;
  37. }
  38. return true;
  39. }
  40. template <class T, class Policy>
  41. T hypergeometric_1F1_divergent_fallback(const T& a, const T& b, const T& z, const Policy& pol, long long& log_scaling)
  42. {
  43. BOOST_MATH_STD_USING
  44. const char* function = "hypergeometric_1F1_divergent_fallback<%1%>(%1%,%1%,%1%)";
  45. //
  46. // We get here if either:
  47. // 1) We decide up front that Tricomi's method won't work, or:
  48. // 2) We've called Tricomi's method and it's failed.
  49. //
  50. if (b > 0)
  51. {
  52. // Commented out since recurrence seems to always be better?
  53. #if 0
  54. if ((z < b) && (a > -50))
  55. // Might as well use a recurrence in preference to z-recurrence:
  56. return hypergeometric_1F1_backward_recurrence_for_negative_a(a, b, z, pol, function, log_scaling);
  57. T z_limit = fabs((2 * a - b) / (sqrt(fabs(a))));
  58. int k = 1 + itrunc(z - z_limit);
  59. // If k is too large we destroy all the digits in the result:
  60. T convergence_at_50 = (b - a + 50) * k / (z * 50);
  61. if ((k > 0) && (k < 50) && (fabs(convergence_at_50) < 1) && (z > z_limit))
  62. {
  63. return boost::math::detail::hypergeometric_1f1_recurrence_on_z_minus_zero(a, b, T(z - k), k, pol, log_scaling);
  64. }
  65. #endif
  66. if (z < b)
  67. return hypergeometric_1F1_backward_recurrence_for_negative_a(a, b, z, pol, function, log_scaling);
  68. else
  69. return hypergeometric_1F1_backwards_recursion_on_b_for_negative_a(a, b, z, pol, function, log_scaling);
  70. }
  71. else // b < 0
  72. {
  73. if (a < 0)
  74. {
  75. if ((b < a) && (z < -b / 4))
  76. // Defensive programming: it is *almost* certain that we can never get here, proving that is hard though...
  77. return hypergeometric_1F1_from_function_ratio_negative_ab(a, b, z, pol, log_scaling); // LCOV_EXCL_LINE
  78. else
  79. {
  80. //
  81. // Solve (a+n)z/((b+n)n) == 1 for n, the number of iterations till the series starts to converge.
  82. // If this is well away from the origin then it's probably better to use the series to evaluate this.
  83. // Note that if sqr is negative then we have no solution, so assign an arbitrarily large value to the
  84. // number of iterations.
  85. //
  86. bool can_use_recursion = (z - b + 100 < boost::math::policies::get_max_series_iterations<Policy>()) && (100 - a < boost::math::policies::get_max_series_iterations<Policy>());
  87. T sqr = 4 * a * z + b * b - 2 * b * z + z * z;
  88. T iterations_to_convergence = sqr > 0 ? T(0.5f * (-sqrt(sqr) - b + z)) : T(-a - b);
  89. if(can_use_recursion && ((std::max)(a, b) + iterations_to_convergence > -300))
  90. return hypergeometric_1F1_backwards_recursion_on_b_for_negative_a(a, b, z, pol, function, log_scaling);
  91. //
  92. // When a < b and if we fall through to the series, then we get divergent behaviour when b crosses the origin
  93. // so ideally we would pick another method. Otherwise the terms immediately after b crosses the origin may
  94. // suffer catastrophic cancellation....
  95. //
  96. if((a < b) && can_use_recursion)
  97. return hypergeometric_1F1_backwards_recursion_on_b_for_negative_a(a, b, z, pol, function, log_scaling);
  98. }
  99. }
  100. else
  101. {
  102. //
  103. // Start by getting the domain of the recurrence relations, we get either:
  104. // -1 Backwards recursion is stable and the CF will converge to double precision.
  105. // +1 Forwards recursion is stable and the CF will converge to double precision.
  106. // 0 No man's land, we're not far enough away from the crossover point to get double precision from either CF.
  107. //
  108. // At higher than double precision we need to be further away from the crossover location to
  109. // get full converge, but it's not clear how much further - indeed at quad precision it's
  110. // basically impossible to ever get forwards iteration to work. Backwards seems to work
  111. // OK as long as a > 1 whatever the precision though.
  112. //
  113. int domain = hypergeometric_1F1_negative_b_recurrence_region(a, b, z);
  114. if ((domain < 0) && ((a > 1) || (boost::math::policies::digits<T, Policy>() <= 64)))
  115. return hypergeometric_1F1_from_function_ratio_negative_b(a, b, z, pol, log_scaling);
  116. else if (domain > 0)
  117. {
  118. if (boost::math::policies::digits<T, Policy>() <= 64)
  119. return hypergeometric_1F1_from_function_ratio_negative_b_forwards(a, b, z, pol, log_scaling);
  120. // LCOV_EXCL_START, what follows is multiprecision only
  121. #ifndef BOOST_MATH_NO_EXCEPTIONS
  122. try
  123. #endif
  124. {
  125. return hypergeometric_1F1_checked_series_impl(a, b, z, pol, log_scaling);
  126. }
  127. #ifndef BOOST_MATH_NO_EXCEPTIONS
  128. catch (const evaluation_error&)
  129. {
  130. //
  131. // The series failed, try the recursions instead and hope we get at least double precision:
  132. //
  133. return hypergeometric_1F1_from_function_ratio_negative_b_forwards(a, b, z, pol, log_scaling);
  134. }
  135. #endif
  136. // LCOV_EXCL_STOP
  137. }
  138. //
  139. // We could fall back to Tricomi's approximation if we're in the transition zone
  140. // between the above two regions. However, I've been unable to find any examples
  141. // where this is better than the series, and there are many cases where it leads to
  142. // quite grievous errors.
  143. /*
  144. else if (allow_tricomi)
  145. {
  146. T aa = a < 1 ? T(1) : a;
  147. if (z < fabs((2 * aa - b) / (sqrt(fabs(aa * b)))))
  148. return hypergeometric_1F1_AS_13_3_7_tricomi(a, b, z, pol, log_scaling);
  149. }
  150. */
  151. }
  152. }
  153. // If we get here, then we've run out of methods to try, use the checked series which will
  154. // raise an error if the result is garbage:
  155. return hypergeometric_1F1_checked_series_impl(a, b, z, pol, log_scaling);
  156. }
  157. #if 0
  158. // Archived, not used, see comments at call site.
  159. template <class T>
  160. bool is_convergent_negative_z_series(const T& a, const T& b, const T& z, const T& b_minus_a)
  161. {
  162. BOOST_MATH_STD_USING
  163. //
  164. // Filter out some cases we don't want first:
  165. //
  166. if((b_minus_a > 0) && (b > 0))
  167. {
  168. if (a < 0)
  169. return false;
  170. }
  171. //
  172. // Generic check: we have small initial divergence and are convergent after 10 terms:
  173. //
  174. if ((fabs(z * a / b) < 2) && (fabs(z * (a + 10) / ((b + 10) * 10)) < 1))
  175. {
  176. // Double check for divergence when we cross the origin on a and b:
  177. if (a < 0)
  178. {
  179. T n = 3 - floor(a);
  180. if (fabs((a + n) * z / ((b + n) * n)) < 1)
  181. {
  182. if (b < 0)
  183. {
  184. T m = 3 - floor(b);
  185. if (fabs((a + m) * z / ((b + m) * m)) < 1)
  186. return true;
  187. }
  188. else
  189. return true;
  190. }
  191. }
  192. else if (b < 0)
  193. {
  194. T n = 3 - floor(b);
  195. if (fabs((a + n) * z / ((b + n) * n)) < 1)
  196. return true;
  197. }
  198. }
  199. if ((b > 0) && (a < 0))
  200. {
  201. //
  202. // For a and z both negative, we're OK with some initial divergence as long as
  203. // it occurs before we hit the origin, as to start with all the terms have the
  204. // same sign.
  205. //
  206. // https://www.wolframalpha.com/input/?i=solve+(a%2Bn)z+%2F+((b%2Bn)n)+%3D%3D+1+for+n
  207. //
  208. T sqr = 4 * a * z + b * b - 2 * b * z + z * z;
  209. T iterations_to_convergence = sqr > 0 ? T(0.5f * (-sqrt(sqr) - b + z)) : T(-a + b);
  210. if (iterations_to_convergence < 0)
  211. iterations_to_convergence = 0.5f * (sqrt(sqr) - b + z);
  212. if (a + iterations_to_convergence < -50)
  213. {
  214. // Need to check for divergence when we cross the origin on a:
  215. if (a > -1)
  216. return true;
  217. T n = 300 - floor(a);
  218. if(fabs((a + n) * z / ((b + n) * n)) < 1)
  219. return true;
  220. }
  221. }
  222. return false;
  223. }
  224. #endif
  225. template <class T>
  226. inline T cyl_bessel_i_shrinkage_rate(const T& z)
  227. {
  228. // Approximately the ratio I_10.5(z/2) / I_9.5(z/2), this gives us an idea of how quickly
  229. // the Bessel terms in A&S 13.6.4 are converging:
  230. if (z < -160)
  231. return 1;
  232. if (z < -40)
  233. return 0.75f;
  234. if (z < -20)
  235. return 0.5f;
  236. if (z < -7)
  237. return 0.25f;
  238. if (z < -2)
  239. return 0.1f;
  240. return 0.05f;
  241. }
  242. template <class T>
  243. inline bool hypergeometric_1F1_is_13_3_6_region(const T& a, const T& b, const T& z)
  244. {
  245. BOOST_MATH_STD_USING
  246. if(fabs(a) == 0.5)
  247. return false;
  248. if ((z < 0) && (fabs(10 * a / b) < 1) && (fabs(a) < 50))
  249. {
  250. T shrinkage = cyl_bessel_i_shrinkage_rate(z);
  251. // We want the first term not too divergent, and convergence by term 10:
  252. if ((fabs((2 * a - 1) * (2 * a - b) / b) < 2) && (fabs(shrinkage * (2 * a + 9) * (2 * a - b + 10) / (10 * (b + 10))) < 0.75))
  253. return true;
  254. }
  255. return false;
  256. }
  257. template <class T>
  258. inline bool hypergeometric_1F1_need_kummer_reflection(const T& a, const T& b, const T& z)
  259. {
  260. BOOST_MATH_STD_USING
  261. //
  262. // Check to see if we should apply Kummer's relation or not:
  263. //
  264. if (z > 0)
  265. return false;
  266. if (z < -1)
  267. return true;
  268. //
  269. // When z is small and negative, things get more complex.
  270. // More often than not we do not need apply Kummer's relation and the
  271. // series is convergent as is, but we do need to check:
  272. //
  273. if (a > 0)
  274. {
  275. if (b > 0)
  276. {
  277. return fabs((a + 10) * z / (10 * (b + 10))) < 1; // Is the 10'th term convergent?
  278. }
  279. else
  280. {
  281. return true; // Likely to be divergent as b crosses the origin
  282. }
  283. }
  284. else // a < 0
  285. {
  286. if (b > 0)
  287. {
  288. return false; // Terms start off all positive and then by the time a crosses the origin we *must* be convergent.
  289. }
  290. else
  291. {
  292. return true; // Likely to be divergent as b crosses the origin, but hard to rationalise about!
  293. }
  294. }
  295. }
  296. template <class T, class Policy>
  297. T hypergeometric_1F1_imp(const T& a, const T& b, const T& z, const Policy& pol, long long& log_scaling)
  298. {
  299. BOOST_MATH_STD_USING // exp, fabs, sqrt
  300. static const char* const function = "boost::math::hypergeometric_1F1<%1%,%1%,%1%>(%1%,%1%,%1%)";
  301. if ((z == 0) || (a == 0))
  302. return T(1);
  303. // undefined result:
  304. if (!detail::check_hypergeometric_1F1_parameters(a, b))
  305. return policies::raise_domain_error<T>(function, "Function is indeterminate for negative integer b = %1%.", b, pol);
  306. // other checks:
  307. if (a == -1)
  308. {
  309. T r = 1 - (z / b);
  310. if (fabs(r) < 0.5)
  311. r = (b - z) / b;
  312. return r;
  313. }
  314. const T b_minus_a = b - a;
  315. // 0f0 a == b case;
  316. if (b_minus_a == 0)
  317. {
  318. if ((a < 0) && (floor(a) == a))
  319. {
  320. // Special case, use the truncated series to match what Mathematica does.
  321. if ((a < -20) && (z > 0) && (z < 1))
  322. {
  323. // https://functions.wolfram.com/HypergeometricFunctions/Hypergeometric1F1/03/01/04/02/0002/
  324. return exp(z) * boost::math::gamma_q(1 - a, z, pol);
  325. }
  326. // https://functions.wolfram.com/HypergeometricFunctions/Hypergeometric1F1/03/01/04/02/0003/
  327. return hypergeometric_1F1_checked_series_impl(a, b, z, pol, log_scaling);
  328. }
  329. long long scale = lltrunc(z, pol);
  330. log_scaling += scale;
  331. return exp(z - scale);
  332. }
  333. // Special case for b-a = -1, we don't use for small a as it throws the digits of a away and leads to large errors:
  334. if ((b_minus_a == -1) && (fabs(a) > 0.5))
  335. {
  336. // for negative small integer a it is reasonable to use truncated series - polynomial
  337. if ((a < 0) && (a == ceil(a)) && (a > -50))
  338. return detail::hypergeometric_1F1_generic_series(a, b, z, pol, log_scaling, function);
  339. log_scaling = lltrunc(floor(z));
  340. T local_z = z - log_scaling;
  341. return (b + z) * exp(local_z) / b;
  342. }
  343. if ((a == 1) && (b == 2))
  344. return boost::math::expm1(z, pol) / z;
  345. if ((b - a == b) && (fabs(z / b) < policies::get_epsilon<T, Policy>()))
  346. return 1;
  347. //
  348. // Special case for A&S 13.3.6:
  349. //
  350. if (z < 0)
  351. {
  352. if (hypergeometric_1F1_is_13_3_6_region(a, b, z))
  353. {
  354. // a is tiny compared to b, and z < 0
  355. // 13.3.6 appears to be the most efficient and often the most accurate method.
  356. T r = boost::math::detail::hypergeometric_1F1_AS_13_3_6(b_minus_a, b, T(-z), a, pol, log_scaling);
  357. long long scale = lltrunc(z, pol);
  358. log_scaling += scale;
  359. return r * exp(z - scale);
  360. }
  361. if ((b < 0) && (fabs(a) < 1e-2))
  362. {
  363. //
  364. // This is a tricky area, potentially we have no good method at all:
  365. //
  366. if (b - ceil(b) == a)
  367. {
  368. // Fractional parts of a and b are genuinely equal, we might as well
  369. // apply Kummer's relation and get a truncated series:
  370. long long scaling = lltrunc(z);
  371. T r = exp(z - scaling) * detail::hypergeometric_1F1_imp<T>(b_minus_a, b, -z, pol, log_scaling);
  372. log_scaling += scaling;
  373. return r;
  374. }
  375. if ((b < -1) && (max_b_for_1F1_small_a_negative_b_by_ratio(z) < b))
  376. return hypergeometric_1F1_small_a_negative_b_by_ratio(a, b, z, pol, log_scaling);
  377. if ((b > -1) && (b < -0.5f))
  378. {
  379. // Recursion is meta-stable:
  380. T first = hypergeometric_1F1_imp(a, T(b + 2), z, pol);
  381. T second = hypergeometric_1F1_imp(a, T(b + 1), z, pol);
  382. return tools::apply_recurrence_relation_backward(hypergeometric_1F1_recurrence_small_b_coefficients<T>(a, b, z, 1), 1, first, second);
  383. }
  384. //
  385. // We've got nothing left but 13.3.6, even though it may be initially divergent:
  386. //
  387. T r = boost::math::detail::hypergeometric_1F1_AS_13_3_6(b_minus_a, b, T(-z), a, pol, log_scaling);
  388. long long scale = lltrunc(z, pol);
  389. log_scaling += scale;
  390. return r * exp(z - scale);
  391. }
  392. }
  393. //
  394. // Asymptotic expansion for large z
  395. // TODO: check region for higher precision types.
  396. // Use recurrence relations to move to this region when a and b are also large.
  397. //
  398. if (detail::hypergeometric_1F1_asym_region(a, b, z, pol))
  399. {
  400. long long saved_scale = log_scaling;
  401. #ifndef BOOST_MATH_NO_EXCEPTIONS
  402. try
  403. #endif
  404. {
  405. return hypergeometric_1F1_asym_large_z_series(a, b, z, pol, log_scaling);
  406. }
  407. #ifndef BOOST_MATH_NO_EXCEPTIONS
  408. catch (const evaluation_error&)
  409. {
  410. }
  411. #endif
  412. //
  413. // Very occasionally our convergence criteria don't quite go to full precision
  414. // and we have to try another method:
  415. //
  416. log_scaling = saved_scale;
  417. }
  418. if ((fabs(a * z / b) < 3.5) && (fabs(z * 100) < fabs(b)) && ((fabs(a) > 1e-2) || (b < -5)))
  419. return detail::hypergeometric_1F1_rational(a, b, z, pol);
  420. if (hypergeometric_1F1_need_kummer_reflection(a, b, z))
  421. {
  422. if (a == 1)
  423. return detail::hypergeometric_1F1_pade(b, z, pol);
  424. #if 0
  425. //
  426. // Commented out: is_convergent_negative_z_series is fine so far as it goes
  427. // but there appear to be no cases that use it, and in extremis, we will
  428. // fall through to the series evaluation anyway.
  429. //
  430. if (is_convergent_negative_z_series(a, b, z, b_minus_a))
  431. {
  432. if ((boost::math::sign(b_minus_a) == boost::math::sign(b)) && ((b > 0) || (b < -200)))
  433. {
  434. // Series is close enough to convergent that we should be OK,
  435. // In this domain b - a ~ b and since 1F1[a, a, z] = e^z 1F1[b-a, b, -z]
  436. // and 1F1[a, a, -z] = e^-z the result must necessarily be somewhere near unity.
  437. // We have to rule out b small and negative because if b crosses the origin early
  438. // in the series (before we're pretty much converged) then all bets are off.
  439. // Note that this can go badly wrong when b and z are both large and negative,
  440. // in that situation the series goes in waves of large and small values which
  441. // may or may not cancel out. Likewise the initial part of the series may or may
  442. // not converge, and even if it does may or may not give a correct answer!
  443. // For example 1F1[-small, -1252.5, -1043.7] can loose up to ~800 digits due to
  444. // cancellation and is basically incalculable via this method.
  445. return hypergeometric_1F1_checked_series_impl(a, b, z, pol, log_scaling);
  446. }
  447. }
  448. #endif
  449. if ((b < 0) && (floor(b) == b))
  450. {
  451. // Negative integer b, so a must be a negative integer too.
  452. // Kummer's transformation fails here!
  453. if(a > -50)
  454. return detail::hypergeometric_1F1_generic_series(a, b, z, pol, log_scaling, function);
  455. // Is there anything better than this??
  456. return hypergeometric_1F1_imp(a, float_next(b), z, pol, log_scaling);
  457. }
  458. else
  459. {
  460. // Let's otherwise make z positive (almost always)
  461. // by Kummer's transformation
  462. // (we also don't transform if z belongs to [-1,0])
  463. // Also note that Kummer's transformation fails when b is
  464. // a negative integer, although this seems to be unmentioned
  465. // in the literature...
  466. long long scaling = lltrunc(z);
  467. T r = exp(z - scaling) * detail::hypergeometric_1F1_imp<T>(b_minus_a, b, -z, pol, log_scaling);
  468. log_scaling += scaling;
  469. return r;
  470. }
  471. }
  472. //
  473. // Check for initial divergence:
  474. //
  475. bool series_is_divergent = (a + 1) * z / (b + 1) < -1;
  476. if (series_is_divergent && (a < 0) && (b < 0) && (a > -1))
  477. series_is_divergent = false; // Best off taking the series in this situation
  478. //
  479. // If series starts off non-divergent, and becomes divergent later
  480. // then it's because both a and b are negative, so check for later
  481. // divergence as well:
  482. //
  483. if (!series_is_divergent && (a < 0) && (b < 0) && (b > a))
  484. {
  485. //
  486. // We need to exclude situations where we're over the initial "hump"
  487. // in the series terms (ie series has already converged by the time
  488. // b crosses the origin:
  489. //
  490. //T fa = fabs(a);
  491. //T fb = fabs(b);
  492. T convergence_point = sqrt((a - 1) * (a - b)) - a;
  493. if (-b < convergence_point)
  494. {
  495. T n = -floor(b);
  496. series_is_divergent = (a + n) * z / ((b + n) * n) < -1;
  497. }
  498. }
  499. else if (!series_is_divergent && (b < 0) && (a > 0))
  500. {
  501. // Series almost always become divergent as b crosses the origin:
  502. series_is_divergent = true;
  503. }
  504. if (series_is_divergent && (b < -1) && (b > -5) && (a > b))
  505. series_is_divergent = false; // don't bother with divergence, series will be OK
  506. //
  507. // Test for alternating series due to negative a,
  508. // in particular, see if the series is initially divergent
  509. // If so use the recurrence relation on a:
  510. //
  511. if (series_is_divergent)
  512. {
  513. if((a < 0) && (floor(a) == a) && (-a < policies::get_max_series_iterations<Policy>()))
  514. // This works amazingly well for negative integer a:
  515. return hypergeometric_1F1_backward_recurrence_for_negative_a(a, b, z, pol, function, log_scaling);
  516. //
  517. // In what follows we have to set limits on how large z can be otherwise
  518. // the Bessel series become large and divergent and all the digits cancel out.
  519. // The criteria are distinctly empiracle rather than based on a firm analysis
  520. // of the terms in the series.
  521. //
  522. if (b > 0)
  523. {
  524. T z_limit = fabs((2 * a - b) / (sqrt(fabs(a))));
  525. if ((z < z_limit) && hypergeometric_1F1_is_tricomi_viable_positive_b(a, b, z))
  526. return detail::hypergeometric_1F1_AS_13_3_7_tricomi(a, b, z, pol, log_scaling);
  527. }
  528. else // b < 0
  529. {
  530. if (a < 0)
  531. {
  532. T z_limit = fabs((2 * a - b) / (sqrt(fabs(a))));
  533. //
  534. // I hate these hard limits, but they're about the best we can do to try and avoid
  535. // Bessel function internal failures: these will be caught and handled
  536. // but up the expense of this function call:
  537. //
  538. if (((z < z_limit) || (a > -500)) && ((b > -500) || (b - 2 * a > 0)) && (z < -a))
  539. {
  540. //
  541. // Outside this domain we will probably get better accuracy from the recursive methods.
  542. //
  543. if(!(((a < b) && (z > -b)) || (z > z_limit)))
  544. return detail::hypergeometric_1F1_AS_13_3_7_tricomi(a, b, z, pol, log_scaling);
  545. //
  546. // When b and z are both very small, we get large errors from the recurrence methods
  547. // in the fallbacks. Tricomi seems to work well here, as does direct series evaluation
  548. // at least some of the time. Picking the right method is not easy, and sometimes this
  549. // is much worse than the fallback. Overall though, it's a reasonable choice that keeps
  550. // the very worst errors under control.
  551. //
  552. if(b > -1)
  553. return detail::hypergeometric_1F1_AS_13_3_7_tricomi(a, b, z, pol, log_scaling);
  554. }
  555. }
  556. //
  557. // We previously used Tricomi here, but it appears to be worse than
  558. // the recurrence-based algorithms in hypergeometric_1F1_divergent_fallback.
  559. /*
  560. else
  561. {
  562. T aa = a < 1 ? T(1) : a;
  563. if (z < fabs((2 * aa - b) / (sqrt(fabs(aa * b)))))
  564. return detail::hypergeometric_1F1_AS_13_3_7_tricomi(a, b, z, pol, log_scaling);
  565. }*/
  566. }
  567. return hypergeometric_1F1_divergent_fallback(a, b, z, pol, log_scaling);
  568. }
  569. if (hypergeometric_1F1_is_13_3_6_region(b_minus_a, b, T(-z)))
  570. {
  571. // b_minus_a is tiny compared to b, and -z < 0
  572. // 13.3.6 appears to be the most efficient and often the most accurate method.
  573. return boost::math::detail::hypergeometric_1F1_AS_13_3_6(a, b, z, b_minus_a, pol, log_scaling);
  574. }
  575. #if 0
  576. if ((a > 0) && (b > 0) && (a * z / b > 2))
  577. {
  578. //
  579. // Series is initially divergent and slow to converge, see if applying
  580. // Kummer's relation can improve things:
  581. //
  582. if (is_convergent_negative_z_series(b_minus_a, b, T(-z), b_minus_a))
  583. {
  584. long long scaling = lltrunc(z);
  585. T r = exp(z - scaling) * detail::hypergeometric_1F1_checked_series_impl(b_minus_a, b, T(-z), pol, log_scaling);
  586. log_scaling += scaling;
  587. return r;
  588. }
  589. }
  590. #endif
  591. if ((a > 0) && (b > 0) && (a * z > 50))
  592. return detail::hypergeometric_1F1_large_abz(a, b, z, pol, log_scaling);
  593. if (b < 0)
  594. return detail::hypergeometric_1F1_checked_series_impl(a, b, z, pol, log_scaling);
  595. return detail::hypergeometric_1F1_generic_series(a, b, z, pol, log_scaling, function);
  596. }
  597. template <class T, class Policy>
  598. inline T hypergeometric_1F1_imp(const T& a, const T& b, const T& z, const Policy& pol)
  599. {
  600. BOOST_MATH_STD_USING // exp, fabs, sqrt
  601. long long log_scaling = 0;
  602. T result = hypergeometric_1F1_imp(a, b, z, pol, log_scaling);
  603. //
  604. // Actual result will be result * e^log_scaling.
  605. //
  606. static const thread_local long long max_scaling = lltrunc(boost::math::tools::log_max_value<T>()) - 2;
  607. static const thread_local T max_scale_factor = exp(T(max_scaling));
  608. while (log_scaling > max_scaling)
  609. {
  610. result *= max_scale_factor;
  611. log_scaling -= max_scaling;
  612. }
  613. while (log_scaling < -max_scaling)
  614. {
  615. result /= max_scale_factor;
  616. log_scaling += max_scaling;
  617. }
  618. if (log_scaling)
  619. result *= exp(T(log_scaling));
  620. return result;
  621. }
  622. template <class T, class Policy>
  623. inline T log_hypergeometric_1F1_imp(const T& a, const T& b, const T& z, int* sign, const Policy& pol)
  624. {
  625. BOOST_MATH_STD_USING // exp, fabs, sqrt
  626. long long log_scaling = 0;
  627. T result = hypergeometric_1F1_imp(a, b, z, pol, log_scaling);
  628. if (sign)
  629. *sign = result < 0 ? -1 : 1;
  630. result = log(fabs(result)) + log_scaling;
  631. return result;
  632. }
  633. template <class T, class Policy>
  634. inline T hypergeometric_1F1_regularized_imp(const T& a, const T& b, const T& z, const Policy& pol)
  635. {
  636. BOOST_MATH_STD_USING // exp, fabs, sqrt
  637. long long log_scaling = 0;
  638. T result = hypergeometric_1F1_imp(a, b, z, pol, log_scaling);
  639. //
  640. // Actual result will be result * e^log_scaling / tgamma(b).
  641. //
  642. int result_sign = 1;
  643. T scale = log_scaling - boost::math::lgamma(b, &result_sign, pol);
  644. static const thread_local T max_scaling = boost::math::tools::log_max_value<T>() - 2;
  645. static const thread_local T max_scale_factor = exp(max_scaling);
  646. while (scale > max_scaling)
  647. {
  648. if((fabs(result) > 1) && (fabs(tools::max_value<T>()) / result <= max_scale_factor))
  649. return policies::raise_overflow_error<T>("hypergeometric_1F1_regularized", nullptr, pol);
  650. // This is *probably* unreachable:
  651. // LCOV_EXCL_START
  652. result *= max_scale_factor;
  653. scale -= max_scaling;
  654. // LCOV_EXCL_STOP
  655. }
  656. while (scale < -max_scaling)
  657. {
  658. result /= max_scale_factor;
  659. scale += max_scaling;
  660. }
  661. if (scale != 0)
  662. {
  663. scale = exp(scale);
  664. if ((scale > 1) && (fabs(result) > 1) && (fabs(tools::max_value<T>() / result) <= scale))
  665. return policies::raise_overflow_error<T>("hypergeometric_1F1_regularized", nullptr, pol);
  666. result *= scale;
  667. }
  668. return result * result_sign;
  669. }
  670. } // namespace detail
  671. template <class T1, class T2, class T3, class Policy>
  672. inline typename tools::promote_args<T1, T2, T3>::type hypergeometric_1F1(T1 a, T2 b, T3 z, const Policy& /* pol */)
  673. {
  674. BOOST_FPU_EXCEPTION_GUARD
  675. typedef typename tools::promote_args<T1, T2, T3>::type result_type;
  676. typedef typename policies::evaluation<result_type, Policy>::type value_type;
  677. typedef typename policies::normalise<
  678. Policy,
  679. policies::promote_float<false>,
  680. policies::promote_double<false>,
  681. policies::discrete_quantile<>,
  682. policies::assert_undefined<> >::type forwarding_policy;
  683. return policies::checked_narrowing_cast<result_type, Policy>(
  684. detail::hypergeometric_1F1_imp<value_type>(
  685. static_cast<value_type>(a),
  686. static_cast<value_type>(b),
  687. static_cast<value_type>(z),
  688. forwarding_policy()),
  689. "boost::math::hypergeometric_1F1<%1%>(%1%,%1%,%1%)");
  690. }
  691. template <class T1, class T2, class T3>
  692. inline typename tools::promote_args<T1, T2, T3>::type hypergeometric_1F1(T1 a, T2 b, T3 z)
  693. {
  694. return hypergeometric_1F1(a, b, z, policies::policy<>());
  695. }
  696. template <class T1, class T2, class T3, class Policy>
  697. inline typename tools::promote_args<T1, T2, T3>::type hypergeometric_1F1_regularized(T1 a, T2 b, T3 z, const Policy& /* pol */)
  698. {
  699. BOOST_FPU_EXCEPTION_GUARD
  700. typedef typename tools::promote_args<T1, T2, T3>::type result_type;
  701. typedef typename policies::evaluation<result_type, Policy>::type value_type;
  702. typedef typename policies::normalise<
  703. Policy,
  704. policies::promote_float<false>,
  705. policies::promote_double<false>,
  706. policies::discrete_quantile<>,
  707. policies::assert_undefined<> >::type forwarding_policy;
  708. return policies::checked_narrowing_cast<result_type, Policy>(
  709. detail::hypergeometric_1F1_regularized_imp<value_type>(
  710. static_cast<value_type>(a),
  711. static_cast<value_type>(b),
  712. static_cast<value_type>(z),
  713. forwarding_policy()),
  714. "boost::math::hypergeometric_1F1<%1%>(%1%,%1%,%1%)");
  715. }
  716. template <class T1, class T2, class T3>
  717. inline typename tools::promote_args<T1, T2, T3>::type hypergeometric_1F1_regularized(T1 a, T2 b, T3 z)
  718. {
  719. return hypergeometric_1F1_regularized(a, b, z, policies::policy<>());
  720. }
  721. template <class T1, class T2, class T3, class Policy>
  722. inline typename tools::promote_args<T1, T2, T3>::type log_hypergeometric_1F1(T1 a, T2 b, T3 z, const Policy& /* pol */)
  723. {
  724. BOOST_FPU_EXCEPTION_GUARD
  725. typedef typename tools::promote_args<T1, T2, T3>::type result_type;
  726. typedef typename policies::evaluation<result_type, Policy>::type value_type;
  727. typedef typename policies::normalise<
  728. Policy,
  729. policies::promote_float<false>,
  730. policies::promote_double<false>,
  731. policies::discrete_quantile<>,
  732. policies::assert_undefined<> >::type forwarding_policy;
  733. return policies::checked_narrowing_cast<result_type, Policy>(
  734. detail::log_hypergeometric_1F1_imp<value_type>(
  735. static_cast<value_type>(a),
  736. static_cast<value_type>(b),
  737. static_cast<value_type>(z),
  738. 0,
  739. forwarding_policy()),
  740. "boost::math::hypergeometric_1F1<%1%>(%1%,%1%,%1%)");
  741. }
  742. template <class T1, class T2, class T3>
  743. inline typename tools::promote_args<T1, T2, T3>::type log_hypergeometric_1F1(T1 a, T2 b, T3 z)
  744. {
  745. return log_hypergeometric_1F1(a, b, z, policies::policy<>());
  746. }
  747. template <class T1, class T2, class T3, class Policy>
  748. inline typename tools::promote_args<T1, T2, T3>::type log_hypergeometric_1F1(T1 a, T2 b, T3 z, int* sign, const Policy& /* pol */)
  749. {
  750. BOOST_FPU_EXCEPTION_GUARD
  751. typedef typename tools::promote_args<T1, T2, T3>::type result_type;
  752. typedef typename policies::evaluation<result_type, Policy>::type value_type;
  753. typedef typename policies::normalise<
  754. Policy,
  755. policies::promote_float<false>,
  756. policies::promote_double<false>,
  757. policies::discrete_quantile<>,
  758. policies::assert_undefined<> >::type forwarding_policy;
  759. return policies::checked_narrowing_cast<result_type, Policy>(
  760. detail::log_hypergeometric_1F1_imp<value_type>(
  761. static_cast<value_type>(a),
  762. static_cast<value_type>(b),
  763. static_cast<value_type>(z),
  764. sign,
  765. forwarding_policy()),
  766. "boost::math::hypergeometric_1F1<%1%>(%1%,%1%,%1%)");
  767. }
  768. template <class T1, class T2, class T3>
  769. inline typename tools::promote_args<T1, T2, T3>::type log_hypergeometric_1F1(T1 a, T2 b, T3 z, int* sign)
  770. {
  771. return log_hypergeometric_1F1(a, b, z, sign, policies::policy<>());
  772. }
  773. } } // namespace boost::math
  774. #endif // BOOST_MATH_HYPERGEOMETRIC_HPP