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