1 // Copyright (c) 2006 Xiaogang Zhang
2 // Use, modification and distribution are subject to the
3 // Boost Software License, Version 1.0. (See accompanying file
4 // LICENSE_1_0.txt or copy at http://www.boost.org/LICENSE_1_0.txt)
6 #ifndef BOOST_MATH_BESSEL_JY_HPP
7 #define BOOST_MATH_BESSEL_JY_HPP
13 #include <boost/math/tools/config.hpp>
14 #include <boost/math/special_functions/gamma.hpp>
15 #include <boost/math/special_functions/sign.hpp>
16 #include <boost/math/special_functions/hypot.hpp>
17 #include <boost/math/special_functions/sin_pi.hpp>
18 #include <boost/math/special_functions/cos_pi.hpp>
19 #include <boost/math/special_functions/detail/bessel_jy_asym.hpp>
20 #include <boost/math/special_functions/detail/bessel_jy_series.hpp>
21 #include <boost/math/constants/constants.hpp>
22 #include <boost/math/policies/error_handling.hpp>
23 #include <boost/mpl/if.hpp>
24 #include <boost/type_traits/is_floating_point.hpp>
27 // Bessel functions of the first and second kind of fractional order
29 namespace boost { namespace math {
34 // Simultaneous calculation of A&S 9.2.9 and 9.2.10
35 // for use in A&S 9.2.5 and 9.2.6.
36 // This series is quick to evaluate, but divergent unless
37 // x is very large, in fact it's pretty hard to figure out
38 // with any degree of precision when this series actually
39 // *will* converge!! Consequently, we may just have to
42 template <class T, class Policy>
43 bool hankel_PQ(T v, T x, T* p, T* q, const Policy& )
46 T tolerance = 2 * policies::get_epsilon<T, Policy>();
57 term *= (mu - sq * sq) / (k * z8);
61 T mult = (sq * sq - mu) / (k * z8);
62 ok = fabs(mult) < 0.5f;
68 while((fabs(term) > tolerance * *p) && ok);
72 // Calculate Y(v, x) and Y(v+1, x) by Temme's method, see
73 // Temme, Journal of Computational Physics, vol 21, 343 (1976)
74 template <typename T, typename Policy>
75 int temme_jy(T v, T x, T* Y, T* Y1, const Policy& pol)
77 T g, h, p, q, f, coef, sum, sum1, tolerance;
82 using namespace boost::math::tools;
83 using namespace boost::math::constants;
85 BOOST_ASSERT(fabs(v) <= 0.5f); // precondition for using this routine
87 T gp = boost::math::tgamma1pm1(v, pol);
88 T gm = boost::math::tgamma1pm1(-v, pol);
89 T spv = boost::math::sin_pi(v, pol);
90 T spv2 = boost::math::sin_pi(v/2, pol);
95 d = abs(sigma) < tools::epsilon<T>() ?
96 T(1) : sinh(sigma) / sigma;
97 e = abs(v) < tools::epsilon<T>() ? T(v*pi<T>()*pi<T>() / 2)
98 : T(2 * spv2 * spv2 / v);
100 T g1 = (v == 0) ? T(-euler<T>()) : T((gp - gm) / ((1 + gp) * (1 + gm) * 2 * v));
101 T g2 = (2 + gp + gm) / ((1 + gp) * (1 + gm) * 2);
102 T vspv = (fabs(v) < tools::epsilon<T>()) ? T(1/constants::pi<T>()) : T(v / spv);
103 f = (g1 * cosh(sigma) - g2 * a * d) * 2 * vspv;
105 p = vspv / (xp * (1 + gm));
106 q = vspv * xp / (1 + gp);
115 T coef_mult = -x * x / 4;
118 tolerance = policies::get_epsilon<T, Policy>();
119 for (k = 1; k < policies::get_max_series_iterations<Policy>(); k++)
121 f = (k * f + p + q) / (k*k - v2);
126 coef *= coef_mult / k;
129 if (abs(coef * g) < abs(sum) * tolerance)
134 policies::check_series_iterations<T>("boost::math::bessel_jy<%1%>(%1%,%1%) in temme_jy", k, pol);
141 // Evaluate continued fraction fv = J_(v+1) / J_v, see
142 // Abramowitz and Stegun, Handbook of Mathematical Functions, 1972, 9.1.73
143 template <typename T, typename Policy>
144 int CF1_jy(T v, T x, T* fv, int* sign, const Policy& pol)
146 T C, D, f, a, b, delta, tiny, tolerance;
152 // |x| <= |v|, CF1_jy converges rapidly
153 // |x| > |v|, CF1_jy needs O(|x|) iterations to converge
155 // modified Lentz's method, see
156 // Lentz, Applied Optics, vol 15, 668 (1976)
157 tolerance = 2 * policies::get_epsilon<T, Policy>();;
158 tiny = sqrt(tools::min_value<T>());
159 C = f = tiny; // b0 = 0, replace with tiny
161 for (k = 1; k < policies::get_max_series_iterations<Policy>() * 100; k++)
167 if (C == 0) { C = tiny; }
168 if (D == 0) { D = tiny; }
172 if (D < 0) { s = -s; }
173 if (abs(delta - 1) < tolerance)
176 policies::check_series_iterations<T>("boost::math::bessel_jy<%1%>(%1%,%1%) in CF1_jy", k / 100, pol);
178 *sign = s; // sign of denominator
183 // This algorithm was originally written by Xiaogang Zhang
184 // using std::complex to perform the complex arithmetic.
185 // However, that turns out to 10x or more slower than using
186 // all real-valued arithmetic, so it's been rewritten using
189 template <typename T, typename Policy>
190 int CF2_jy(T v, T x, T* p, T* q, const Policy& pol)
194 T Cr, Ci, Dr, Di, fr, fi, a, br, bi, delta_r, delta_i, temp;
198 // |x| >= |v|, CF2_jy converges rapidly
199 // |x| -> 0, CF2_jy fails to converge
200 BOOST_ASSERT(fabs(x) > 1);
202 // modified Lentz's method, complex numbers involved, see
203 // Lentz, Applied Optics, vol 15, 668 (1976)
204 T tolerance = 2 * policies::get_epsilon<T, Policy>();
205 tiny = sqrt(tools::min_value<T>());
210 a = (0.25f - v2) / x; // Note complex this one time only!
214 Ci = bi + a * Cr / temp;
218 if (fabs(Cr) + fabs(Ci) < tiny) { Cr = tiny; }
219 if (fabs(Dr) + fabs(Di) < tiny) { Dr = tiny; }
220 temp = Dr * Dr + Di * Di;
223 delta_r = Cr * Dr - Ci * Di;
224 delta_i = Ci * Dr + Cr * Di;
226 fr = temp * delta_r - fi * delta_i;
227 fi = temp * delta_i + fi * delta_r;
228 for (k = 2; k < policies::get_max_series_iterations<Policy>(); k++)
234 temp = Cr * Cr + Ci * Ci;
235 Cr = br + a * Cr / temp;
236 Ci = bi - a * Ci / temp;
239 if (fabs(Cr) + fabs(Ci) < tiny) { Cr = tiny; }
240 if (fabs(Dr) + fabs(Di) < tiny) { Dr = tiny; }
241 temp = Dr * Dr + Di * Di;
244 delta_r = Cr * Dr - Ci * Di;
245 delta_i = Ci * Dr + Cr * Di;
247 fr = temp * delta_r - fi * delta_i;
248 fi = temp * delta_i + fi * delta_r;
249 if (fabs(delta_r - 1) + fabs(delta_i) < tolerance)
252 policies::check_series_iterations<T>("boost::math::bessel_jy<%1%>(%1%,%1%) in CF2_jy", k, pol);
259 static const int need_j = 1;
260 static const int need_y = 2;
262 // Compute J(v, x) and Y(v, x) simultaneously by Steed's method, see
263 // Barnett et al, Computer Physics Communications, vol 8, 377 (1974)
264 template <typename T, typename Policy>
265 int bessel_jy(T v, T x, T* J, T* Y, int kind, const Policy& pol)
267 BOOST_ASSERT(x >= 0);
269 T u, Jv, Ju, Yv, Yv1, Yu, Yu1(0), fv, fu;
270 T W, p, q, gamma, current, prev, next;
271 bool reflect = false;
278 static const char* function = "boost::math::bessel_jy<%1%>(%1%,%1%)";
281 using namespace boost::math::tools;
282 using namespace boost::math::constants;
287 v = -v; // v is non-negative from here
289 if(v > static_cast<T>((std::numeric_limits<int>::max)()))
290 policies::raise_evaluation_error<T>(function, "Order of Bessel function is too large to evaluate: got %1%", v, pol);
292 u = v - n; // -1/2 <= u < 1/2
297 cp = boost::math::cos_pi(z, pol);
298 sp = boost::math::sin_pi(z, pol);
300 kind = need_j|need_y; // need both for reflection formula
307 else if((u == 0) || !reflect)
309 else if(kind & need_j)
310 *J = policies::raise_domain_error<T>(function, "Value of Bessel J_v(x) is complex-infinity at %1%", x, pol); // complex infinity
312 *J = std::numeric_limits<T>::quiet_NaN(); // any value will do, not using J.
314 if((kind & need_y) == 0)
315 *Y = std::numeric_limits<T>::quiet_NaN(); // any value will do, not using Y.
317 *Y = -policies::raise_overflow_error<T>(function, 0, pol);
319 *Y = policies::raise_domain_error<T>(function, "Value of Bessel Y_v(x) is complex-infinity at %1%", x, pol); // complex infinity
323 // x is positive until reflection
324 W = T(2) / (x * pi<T>()); // Wronskian
326 if(((kind & need_y) == 0) && ((x < 1) || (v > x * x / 4) || (x < 5)))
329 // This series will actually converge rapidly for all small
330 // x - say up to x < 20 - but the first few terms are large
331 // and divergent which leads to large errors :-(
333 Jv = bessel_j_small_z_series(v, x, pol);
334 Yv = std::numeric_limits<T>::quiet_NaN();
336 else if((x < 1) && (u != 0) && (log(policies::get_epsilon<T, Policy>() / 2) > v * log((x/2) * (x/2) / v)))
338 // Evaluate using series representations.
339 // This is particularly important for x << v as in this
340 // area temme_jy may be slow to converge, if it converges at all.
341 // Requires x is not an integer.
343 Jv = bessel_j_small_z_series(v, x, pol);
345 Jv = std::numeric_limits<T>::quiet_NaN();
346 if((org_kind&need_y && (!reflect || (cp != 0)))
347 || (org_kind & need_j && (reflect && (sp != 0))))
349 // Only calculate if we need it, and if the reflection formula will actually use it:
350 Yv = bessel_y_small_z_series(v, x, &Yv_scale, pol);
353 Yv = std::numeric_limits<T>::quiet_NaN();
355 else if((u == 0) && (x < policies::get_epsilon<T, Policy>()))
357 // Truncated series evaluation for small x and v an integer,
358 // much quicker in this area than temme_jy below.
360 Jv = bessel_j_small_z_series(v, x, pol);
362 Jv = std::numeric_limits<T>::quiet_NaN();
363 if((org_kind&need_y && (!reflect || (cp != 0)))
364 || (org_kind & need_j && (reflect && (sp != 0))))
366 // Only calculate if we need it, and if the reflection formula will actually use it:
367 Yv = bessel_yn_small_z(n, x, &Yv_scale, pol);
370 Yv = std::numeric_limits<T>::quiet_NaN();
372 else if(asymptotic_bessel_large_x_limit(v, x))
376 Yv = asymptotic_bessel_y_large_x_2(v, x);
379 Yv = std::numeric_limits<T>::quiet_NaN(); // any value will do, we're not using it.
382 Jv = asymptotic_bessel_j_large_x_2(v, x);
385 Jv = std::numeric_limits<T>::quiet_NaN(); // any value will do, we're not using it.
387 else if((x > 8) && hankel_PQ(v, x, &p, &q, pol))
390 // Hankel approximation: note that this method works best when x
391 // is large, but in that case we end up calculating sines and cosines
392 // of large values, with horrendous resulting accuracy. It is fast though
395 // Normally we calculate sin/cos(chi) where:
397 // chi = x - fmod(T(v / 2 + 0.25f), T(2)) * boost::math::constants::pi<T>();
399 // But this introduces large errors, so use sin/cos addition formulae to
402 T mod_v = fmod(T(v / 2 + 0.25f), T(2));
405 T sv = sin_pi(mod_v);
406 T cv = cos_pi(mod_v);
408 T sc = sx * cv - sv * cx; // == sin(chi);
409 T cc = cx * cv + sx * sv; // == cos(chi);
410 T chi = boost::math::constants::root_two<T>() / (boost::math::constants::root_pi<T>() * sqrt(x)); //sqrt(2 / (boost::math::constants::pi<T>() * x));
411 Yv = chi * (p * sc + q * cc);
412 Jv = chi * (p * cc - q * sc);
414 else if (x <= 2) // x in (0, 2]
416 if(temme_jy(u, x, &Yu, &Yu1, pol)) // Temme series
425 policies::check_series_iterations<T>(function, n, pol);
426 for (k = 1; k <= n; k++) // forward recurrence for Y
428 T fact = 2 * (u + k) / x;
429 if((tools::max_value<T>() - fabs(prev)) / fact < fabs(current))
435 next = fact * current - prev;
443 CF1_jy(v, x, &fv, &s, pol); // continued fraction CF1_jy
444 Jv = scale * W / (Yv * fv - Yv1); // Wronskian relation
447 Jv = std::numeric_limits<T>::quiet_NaN(); // any value will do, we're not using it.
450 else // x in (2, \infty)
455 CF1_jy(v, x, &fv, &s, pol);
456 // tiny initial value to prevent overflow
457 T init = sqrt(tools::min_value<T>());
458 prev = fv * s * init;
460 if(v < max_factorial<T>::value)
462 policies::check_series_iterations<T>(function, n, pol);
463 for (k = n; k > 0; k--) // backward recurrence for J
465 next = 2 * (u + k) * current / x - prev;
469 ratio = (s * init) / current; // scaling ratio
470 // can also call CF1_jy() to get fu, not much difference in precision
476 // When v is large we may get overflow in this calculation
477 // leading to NaN's and other nasty surprises:
479 policies::check_series_iterations<T>(function, n, pol);
481 for (k = n; k > 0; k--) // backward recurrence for J
483 T t = 2 * (u + k) / x;
484 if((t > 1) && (tools::max_value<T>() / t < current))
489 next = t * current - prev;
495 ratio = (s * init) / current; // scaling ratio
496 // can also call CF1_jy() to get fu, not much difference in precision
505 CF2_jy(u, x, &p, &q, pol); // continued fraction CF2_jy
506 T t = u / x - fu; // t = J'/J
509 // We can't allow gamma to cancel out to zero competely as it messes up
510 // the subsequent logic. So pretend that one bit didn't cancel out
511 // and set to a suitably small value. The only test case we've been able to
512 // find for this, is when v = 8.5 and x = 4*PI.
516 gamma = u * tools::epsilon<T>() / x;
518 Ju = sign(current) * sqrt(W / (q + gamma * (p - t)));
520 Jv = Ju * ratio; // normalization
523 Yu1 = Yu * (u/x - p - q/gamma);
530 policies::check_series_iterations<T>(function, n, pol);
531 for (k = 1; k <= n; k++) // forward recurrence for Y
533 T fact = 2 * (u + k) / x;
534 if((tools::max_value<T>() - fabs(prev)) / fact < fabs(current))
540 next = fact * current - prev;
547 Yv = std::numeric_limits<T>::quiet_NaN(); // any value will do, we're not using it.
552 if((sp != 0) && (tools::max_value<T>() * fabs(Yv_scale) < fabs(sp * Yv)))
553 *J = org_kind & need_j ? T(-sign(sp) * sign(Yv) * sign(Yv_scale) * policies::raise_overflow_error<T>(function, 0, pol)) : T(0);
555 *J = cp * Jv - (sp == 0 ? T(0) : T((sp * Yv) / Yv_scale)); // reflection formula
556 if((cp != 0) && (tools::max_value<T>() * fabs(Yv_scale) < fabs(cp * Yv)))
557 *Y = org_kind & need_y ? T(-sign(cp) * sign(Yv) * sign(Yv_scale) * policies::raise_overflow_error<T>(function, 0, pol)) : T(0);
559 *Y = (sp != 0 ? sp * Jv : T(0)) + (cp == 0 ? T(0) : T((cp * Yv) / Yv_scale));
564 if(tools::max_value<T>() * fabs(Yv_scale) < fabs(Yv))
565 *Y = org_kind & need_y ? T(sign(Yv) * sign(Yv_scale) * policies::raise_overflow_error<T>(function, 0, pol)) : T(0);
573 } // namespace detail
577 #endif // BOOST_MATH_BESSEL_JY_HPP