hypergeometric_1F1.hpp 31 KB

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