hypergeometric_1F1.hpp 32 KB

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