bessel_jy.hpp 21 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589
  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)
  5. #ifndef BOOST_MATH_BESSEL_JY_HPP
  6. #define BOOST_MATH_BESSEL_JY_HPP
  7. #ifdef _MSC_VER
  8. #pragma once
  9. #endif
  10. #include <boost/math/tools/config.hpp>
  11. #include <boost/math/special_functions/gamma.hpp>
  12. #include <boost/math/special_functions/sign.hpp>
  13. #include <boost/math/special_functions/hypot.hpp>
  14. #include <boost/math/special_functions/sin_pi.hpp>
  15. #include <boost/math/special_functions/cos_pi.hpp>
  16. #include <boost/math/special_functions/detail/bessel_jy_asym.hpp>
  17. #include <boost/math/special_functions/detail/bessel_jy_series.hpp>
  18. #include <boost/math/constants/constants.hpp>
  19. #include <boost/math/policies/error_handling.hpp>
  20. #include <boost/mpl/if.hpp>
  21. #include <boost/type_traits/is_floating_point.hpp>
  22. #include <complex>
  23. // Bessel functions of the first and second kind of fractional order
  24. namespace boost { namespace math {
  25. namespace detail {
  26. //
  27. // Simultaneous calculation of A&S 9.2.9 and 9.2.10
  28. // for use in A&S 9.2.5 and 9.2.6.
  29. // This series is quick to evaluate, but divergent unless
  30. // x is very large, in fact it's pretty hard to figure out
  31. // with any degree of precision when this series actually
  32. // *will* converge!! Consequently, we may just have to
  33. // try it and see...
  34. //
  35. template <class T, class Policy>
  36. bool hankel_PQ(T v, T x, T* p, T* q, const Policy& )
  37. {
  38. BOOST_MATH_STD_USING
  39. T tolerance = 2 * policies::get_epsilon<T, Policy>();
  40. *p = 1;
  41. *q = 0;
  42. T k = 1;
  43. T z8 = 8 * x;
  44. T sq = 1;
  45. T mu = 4 * v * v;
  46. T term = 1;
  47. bool ok = true;
  48. do
  49. {
  50. term *= (mu - sq * sq) / (k * z8);
  51. *q += term;
  52. k += 1;
  53. sq += 2;
  54. T mult = (sq * sq - mu) / (k * z8);
  55. ok = fabs(mult) < 0.5f;
  56. term *= mult;
  57. *p += term;
  58. k += 1;
  59. sq += 2;
  60. }
  61. while((fabs(term) > tolerance * *p) && ok);
  62. return ok;
  63. }
  64. // Calculate Y(v, x) and Y(v+1, x) by Temme's method, see
  65. // Temme, Journal of Computational Physics, vol 21, 343 (1976)
  66. template <typename T, typename Policy>
  67. int temme_jy(T v, T x, T* Y, T* Y1, const Policy& pol)
  68. {
  69. T g, h, p, q, f, coef, sum, sum1, tolerance;
  70. T a, d, e, sigma;
  71. unsigned long k;
  72. BOOST_MATH_STD_USING
  73. using namespace boost::math::tools;
  74. using namespace boost::math::constants;
  75. BOOST_ASSERT(fabs(v) <= 0.5f); // precondition for using this routine
  76. T gp = boost::math::tgamma1pm1(v, pol);
  77. T gm = boost::math::tgamma1pm1(-v, pol);
  78. T spv = boost::math::sin_pi(v, pol);
  79. T spv2 = boost::math::sin_pi(v/2, pol);
  80. T xp = pow(x/2, v);
  81. a = log(x / 2);
  82. sigma = -a * v;
  83. d = abs(sigma) < tools::epsilon<T>() ?
  84. T(1) : sinh(sigma) / sigma;
  85. e = abs(v) < tools::epsilon<T>() ? T(v*pi<T>()*pi<T>() / 2)
  86. : T(2 * spv2 * spv2 / v);
  87. T g1 = (v == 0) ? T(-euler<T>()) : T((gp - gm) / ((1 + gp) * (1 + gm) * 2 * v));
  88. T g2 = (2 + gp + gm) / ((1 + gp) * (1 + gm) * 2);
  89. T vspv = (fabs(v) < tools::epsilon<T>()) ? T(1/constants::pi<T>()) : T(v / spv);
  90. f = (g1 * cosh(sigma) - g2 * a * d) * 2 * vspv;
  91. p = vspv / (xp * (1 + gm));
  92. q = vspv * xp / (1 + gp);
  93. g = f + e * q;
  94. h = p;
  95. coef = 1;
  96. sum = coef * g;
  97. sum1 = coef * h;
  98. T v2 = v * v;
  99. T coef_mult = -x * x / 4;
  100. // series summation
  101. tolerance = policies::get_epsilon<T, Policy>();
  102. for (k = 1; k < policies::get_max_series_iterations<Policy>(); k++)
  103. {
  104. f = (k * f + p + q) / (k*k - v2);
  105. p /= k - v;
  106. q /= k + v;
  107. g = f + e * q;
  108. h = p - k * g;
  109. coef *= coef_mult / k;
  110. sum += coef * g;
  111. sum1 += coef * h;
  112. if (abs(coef * g) < abs(sum) * tolerance)
  113. {
  114. break;
  115. }
  116. }
  117. policies::check_series_iterations<T>("boost::math::bessel_jy<%1%>(%1%,%1%) in temme_jy", k, pol);
  118. *Y = -sum;
  119. *Y1 = -2 * sum1 / x;
  120. return 0;
  121. }
  122. // Evaluate continued fraction fv = J_(v+1) / J_v, see
  123. // Abramowitz and Stegun, Handbook of Mathematical Functions, 1972, 9.1.73
  124. template <typename T, typename Policy>
  125. int CF1_jy(T v, T x, T* fv, int* sign, const Policy& pol)
  126. {
  127. T C, D, f, a, b, delta, tiny, tolerance;
  128. unsigned long k;
  129. int s = 1;
  130. BOOST_MATH_STD_USING
  131. // |x| <= |v|, CF1_jy converges rapidly
  132. // |x| > |v|, CF1_jy needs O(|x|) iterations to converge
  133. // modified Lentz's method, see
  134. // Lentz, Applied Optics, vol 15, 668 (1976)
  135. tolerance = 2 * policies::get_epsilon<T, Policy>();;
  136. tiny = sqrt(tools::min_value<T>());
  137. C = f = tiny; // b0 = 0, replace with tiny
  138. D = 0;
  139. for (k = 1; k < policies::get_max_series_iterations<Policy>() * 100; k++)
  140. {
  141. a = -1;
  142. b = 2 * (v + k) / x;
  143. C = b + a / C;
  144. D = b + a * D;
  145. if (C == 0) { C = tiny; }
  146. if (D == 0) { D = tiny; }
  147. D = 1 / D;
  148. delta = C * D;
  149. f *= delta;
  150. if (D < 0) { s = -s; }
  151. if (abs(delta - 1) < tolerance)
  152. { break; }
  153. }
  154. policies::check_series_iterations<T>("boost::math::bessel_jy<%1%>(%1%,%1%) in CF1_jy", k / 100, pol);
  155. *fv = -f;
  156. *sign = s; // sign of denominator
  157. return 0;
  158. }
  159. //
  160. // This algorithm was originally written by Xiaogang Zhang
  161. // using std::complex to perform the complex arithmetic.
  162. // However, that turns out to 10x or more slower than using
  163. // all real-valued arithmetic, so it's been rewritten using
  164. // real values only.
  165. //
  166. template <typename T, typename Policy>
  167. int CF2_jy(T v, T x, T* p, T* q, const Policy& pol)
  168. {
  169. BOOST_MATH_STD_USING
  170. T Cr, Ci, Dr, Di, fr, fi, a, br, bi, delta_r, delta_i, temp;
  171. T tiny;
  172. unsigned long k;
  173. // |x| >= |v|, CF2_jy converges rapidly
  174. // |x| -> 0, CF2_jy fails to converge
  175. BOOST_ASSERT(fabs(x) > 1);
  176. // modified Lentz's method, complex numbers involved, see
  177. // Lentz, Applied Optics, vol 15, 668 (1976)
  178. T tolerance = 2 * policies::get_epsilon<T, Policy>();
  179. tiny = sqrt(tools::min_value<T>());
  180. Cr = fr = -0.5f / x;
  181. Ci = fi = 1;
  182. //Dr = Di = 0;
  183. T v2 = v * v;
  184. a = (0.25f - v2) / x; // Note complex this one time only!
  185. br = 2 * x;
  186. bi = 2;
  187. temp = Cr * Cr + 1;
  188. Ci = bi + a * Cr / temp;
  189. Cr = br + a / temp;
  190. Dr = br;
  191. Di = bi;
  192. if (fabs(Cr) + fabs(Ci) < tiny) { Cr = tiny; }
  193. if (fabs(Dr) + fabs(Di) < tiny) { Dr = tiny; }
  194. temp = Dr * Dr + Di * Di;
  195. Dr = Dr / temp;
  196. Di = -Di / temp;
  197. delta_r = Cr * Dr - Ci * Di;
  198. delta_i = Ci * Dr + Cr * Di;
  199. temp = fr;
  200. fr = temp * delta_r - fi * delta_i;
  201. fi = temp * delta_i + fi * delta_r;
  202. for (k = 2; k < policies::get_max_series_iterations<Policy>(); k++)
  203. {
  204. a = k - 0.5f;
  205. a *= a;
  206. a -= v2;
  207. bi += 2;
  208. temp = Cr * Cr + Ci * Ci;
  209. Cr = br + a * Cr / temp;
  210. Ci = bi - a * Ci / temp;
  211. Dr = br + a * Dr;
  212. Di = bi + a * Di;
  213. if (fabs(Cr) + fabs(Ci) < tiny) { Cr = tiny; }
  214. if (fabs(Dr) + fabs(Di) < tiny) { Dr = tiny; }
  215. temp = Dr * Dr + Di * Di;
  216. Dr = Dr / temp;
  217. Di = -Di / temp;
  218. delta_r = Cr * Dr - Ci * Di;
  219. delta_i = Ci * Dr + Cr * Di;
  220. temp = fr;
  221. fr = temp * delta_r - fi * delta_i;
  222. fi = temp * delta_i + fi * delta_r;
  223. if (fabs(delta_r - 1) + fabs(delta_i) < tolerance)
  224. break;
  225. }
  226. policies::check_series_iterations<T>("boost::math::bessel_jy<%1%>(%1%,%1%) in CF2_jy", k, pol);
  227. *p = fr;
  228. *q = fi;
  229. return 0;
  230. }
  231. static const int need_j = 1;
  232. static const int need_y = 2;
  233. // Compute J(v, x) and Y(v, x) simultaneously by Steed's method, see
  234. // Barnett et al, Computer Physics Communications, vol 8, 377 (1974)
  235. template <typename T, typename Policy>
  236. int bessel_jy(T v, T x, T* J, T* Y, int kind, const Policy& pol)
  237. {
  238. BOOST_ASSERT(x >= 0);
  239. T u, Jv, Ju, Yv, Yv1, Yu, Yu1(0), fv, fu;
  240. T W, p, q, gamma, current, prev, next;
  241. bool reflect = false;
  242. unsigned n, k;
  243. int s;
  244. int org_kind = kind;
  245. T cp = 0;
  246. T sp = 0;
  247. static const char* function = "boost::math::bessel_jy<%1%>(%1%,%1%)";
  248. BOOST_MATH_STD_USING
  249. using namespace boost::math::tools;
  250. using namespace boost::math::constants;
  251. if (v < 0)
  252. {
  253. reflect = true;
  254. v = -v; // v is non-negative from here
  255. }
  256. if (v > static_cast<T>((std::numeric_limits<int>::max)()))
  257. {
  258. *J = *Y = policies::raise_evaluation_error<T>(function, "Order of Bessel function is too large to evaluate: got %1%", v, pol);
  259. return 1;
  260. }
  261. n = iround(v, pol);
  262. u = v - n; // -1/2 <= u < 1/2
  263. if(reflect)
  264. {
  265. T z = (u + n % 2);
  266. cp = boost::math::cos_pi(z, pol);
  267. sp = boost::math::sin_pi(z, pol);
  268. if(u != 0)
  269. kind = need_j|need_y; // need both for reflection formula
  270. }
  271. if(x == 0)
  272. {
  273. if(v == 0)
  274. *J = 1;
  275. else if((u == 0) || !reflect)
  276. *J = 0;
  277. else if(kind & need_j)
  278. *J = policies::raise_domain_error<T>(function, "Value of Bessel J_v(x) is complex-infinity at %1%", x, pol); // complex infinity
  279. else
  280. *J = std::numeric_limits<T>::quiet_NaN(); // any value will do, not using J.
  281. if((kind & need_y) == 0)
  282. *Y = std::numeric_limits<T>::quiet_NaN(); // any value will do, not using Y.
  283. else if(v == 0)
  284. *Y = -policies::raise_overflow_error<T>(function, 0, pol);
  285. else
  286. *Y = policies::raise_domain_error<T>(function, "Value of Bessel Y_v(x) is complex-infinity at %1%", x, pol); // complex infinity
  287. return 1;
  288. }
  289. // x is positive until reflection
  290. W = T(2) / (x * pi<T>()); // Wronskian
  291. T Yv_scale = 1;
  292. if(((kind & need_y) == 0) && ((x < 1) || (v > x * x / 4) || (x < 5)))
  293. {
  294. //
  295. // This series will actually converge rapidly for all small
  296. // x - say up to x < 20 - but the first few terms are large
  297. // and divergent which leads to large errors :-(
  298. //
  299. Jv = bessel_j_small_z_series(v, x, pol);
  300. Yv = std::numeric_limits<T>::quiet_NaN();
  301. }
  302. else if((x < 1) && (u != 0) && (log(policies::get_epsilon<T, Policy>() / 2) > v * log((x/2) * (x/2) / v)))
  303. {
  304. // Evaluate using series representations.
  305. // This is particularly important for x << v as in this
  306. // area temme_jy may be slow to converge, if it converges at all.
  307. // Requires x is not an integer.
  308. if(kind&need_j)
  309. Jv = bessel_j_small_z_series(v, x, pol);
  310. else
  311. Jv = std::numeric_limits<T>::quiet_NaN();
  312. if((org_kind&need_y && (!reflect || (cp != 0)))
  313. || (org_kind & need_j && (reflect && (sp != 0))))
  314. {
  315. // Only calculate if we need it, and if the reflection formula will actually use it:
  316. Yv = bessel_y_small_z_series(v, x, &Yv_scale, pol);
  317. }
  318. else
  319. Yv = std::numeric_limits<T>::quiet_NaN();
  320. }
  321. else if((u == 0) && (x < policies::get_epsilon<T, Policy>()))
  322. {
  323. // Truncated series evaluation for small x and v an integer,
  324. // much quicker in this area than temme_jy below.
  325. if(kind&need_j)
  326. Jv = bessel_j_small_z_series(v, x, pol);
  327. else
  328. Jv = std::numeric_limits<T>::quiet_NaN();
  329. if((org_kind&need_y && (!reflect || (cp != 0)))
  330. || (org_kind & need_j && (reflect && (sp != 0))))
  331. {
  332. // Only calculate if we need it, and if the reflection formula will actually use it:
  333. Yv = bessel_yn_small_z(n, x, &Yv_scale, pol);
  334. }
  335. else
  336. Yv = std::numeric_limits<T>::quiet_NaN();
  337. }
  338. else if(asymptotic_bessel_large_x_limit(v, x))
  339. {
  340. if(kind&need_y)
  341. {
  342. Yv = asymptotic_bessel_y_large_x_2(v, x);
  343. }
  344. else
  345. Yv = std::numeric_limits<T>::quiet_NaN(); // any value will do, we're not using it.
  346. if(kind&need_j)
  347. {
  348. Jv = asymptotic_bessel_j_large_x_2(v, x);
  349. }
  350. else
  351. Jv = std::numeric_limits<T>::quiet_NaN(); // any value will do, we're not using it.
  352. }
  353. else if((x > 8) && hankel_PQ(v, x, &p, &q, pol))
  354. {
  355. //
  356. // Hankel approximation: note that this method works best when x
  357. // is large, but in that case we end up calculating sines and cosines
  358. // of large values, with horrendous resulting accuracy. It is fast though
  359. // when it works....
  360. //
  361. // Normally we calculate sin/cos(chi) where:
  362. //
  363. // chi = x - fmod(T(v / 2 + 0.25f), T(2)) * boost::math::constants::pi<T>();
  364. //
  365. // But this introduces large errors, so use sin/cos addition formulae to
  366. // improve accuracy:
  367. //
  368. T mod_v = fmod(T(v / 2 + 0.25f), T(2));
  369. T sx = sin(x);
  370. T cx = cos(x);
  371. T sv = sin_pi(mod_v);
  372. T cv = cos_pi(mod_v);
  373. T sc = sx * cv - sv * cx; // == sin(chi);
  374. T cc = cx * cv + sx * sv; // == cos(chi);
  375. T chi = boost::math::constants::root_two<T>() / (boost::math::constants::root_pi<T>() * sqrt(x)); //sqrt(2 / (boost::math::constants::pi<T>() * x));
  376. Yv = chi * (p * sc + q * cc);
  377. Jv = chi * (p * cc - q * sc);
  378. }
  379. else if (x <= 2) // x in (0, 2]
  380. {
  381. if(temme_jy(u, x, &Yu, &Yu1, pol)) // Temme series
  382. {
  383. // domain error:
  384. *J = *Y = Yu;
  385. return 1;
  386. }
  387. prev = Yu;
  388. current = Yu1;
  389. T scale = 1;
  390. policies::check_series_iterations<T>(function, n, pol);
  391. for (k = 1; k <= n; k++) // forward recurrence for Y
  392. {
  393. T fact = 2 * (u + k) / x;
  394. if((tools::max_value<T>() - fabs(prev)) / fact < fabs(current))
  395. {
  396. scale /= current;
  397. prev /= current;
  398. current = 1;
  399. }
  400. next = fact * current - prev;
  401. prev = current;
  402. current = next;
  403. }
  404. Yv = prev;
  405. Yv1 = current;
  406. if(kind&need_j)
  407. {
  408. CF1_jy(v, x, &fv, &s, pol); // continued fraction CF1_jy
  409. Jv = scale * W / (Yv * fv - Yv1); // Wronskian relation
  410. }
  411. else
  412. Jv = std::numeric_limits<T>::quiet_NaN(); // any value will do, we're not using it.
  413. Yv_scale = scale;
  414. }
  415. else // x in (2, \infty)
  416. {
  417. // Get Y(u, x):
  418. T ratio;
  419. CF1_jy(v, x, &fv, &s, pol);
  420. // tiny initial value to prevent overflow
  421. T init = sqrt(tools::min_value<T>());
  422. BOOST_MATH_INSTRUMENT_VARIABLE(init);
  423. prev = fv * s * init;
  424. current = s * init;
  425. if(v < max_factorial<T>::value)
  426. {
  427. policies::check_series_iterations<T>(function, n, pol);
  428. for (k = n; k > 0; k--) // backward recurrence for J
  429. {
  430. next = 2 * (u + k) * current / x - prev;
  431. prev = current;
  432. current = next;
  433. }
  434. ratio = (s * init) / current; // scaling ratio
  435. // can also call CF1_jy() to get fu, not much difference in precision
  436. fu = prev / current;
  437. }
  438. else
  439. {
  440. //
  441. // When v is large we may get overflow in this calculation
  442. // leading to NaN's and other nasty surprises:
  443. //
  444. policies::check_series_iterations<T>(function, n, pol);
  445. bool over = false;
  446. for (k = n; k > 0; k--) // backward recurrence for J
  447. {
  448. T t = 2 * (u + k) / x;
  449. if((t > 1) && (tools::max_value<T>() / t < current))
  450. {
  451. over = true;
  452. break;
  453. }
  454. next = t * current - prev;
  455. prev = current;
  456. current = next;
  457. }
  458. if(!over)
  459. {
  460. ratio = (s * init) / current; // scaling ratio
  461. // can also call CF1_jy() to get fu, not much difference in precision
  462. fu = prev / current;
  463. }
  464. else
  465. {
  466. ratio = 0;
  467. fu = 1;
  468. }
  469. }
  470. CF2_jy(u, x, &p, &q, pol); // continued fraction CF2_jy
  471. T t = u / x - fu; // t = J'/J
  472. gamma = (p - t) / q;
  473. //
  474. // We can't allow gamma to cancel out to zero competely as it messes up
  475. // the subsequent logic. So pretend that one bit didn't cancel out
  476. // and set to a suitably small value. The only test case we've been able to
  477. // find for this, is when v = 8.5 and x = 4*PI.
  478. //
  479. if(gamma == 0)
  480. {
  481. gamma = u * tools::epsilon<T>() / x;
  482. }
  483. BOOST_MATH_INSTRUMENT_VARIABLE(current);
  484. BOOST_MATH_INSTRUMENT_VARIABLE(W);
  485. BOOST_MATH_INSTRUMENT_VARIABLE(q);
  486. BOOST_MATH_INSTRUMENT_VARIABLE(gamma);
  487. BOOST_MATH_INSTRUMENT_VARIABLE(p);
  488. BOOST_MATH_INSTRUMENT_VARIABLE(t);
  489. Ju = sign(current) * sqrt(W / (q + gamma * (p - t)));
  490. BOOST_MATH_INSTRUMENT_VARIABLE(Ju);
  491. Jv = Ju * ratio; // normalization
  492. Yu = gamma * Ju;
  493. Yu1 = Yu * (u/x - p - q/gamma);
  494. if(kind&need_y)
  495. {
  496. // compute Y:
  497. prev = Yu;
  498. current = Yu1;
  499. policies::check_series_iterations<T>(function, n, pol);
  500. for (k = 1; k <= n; k++) // forward recurrence for Y
  501. {
  502. T fact = 2 * (u + k) / x;
  503. if((tools::max_value<T>() - fabs(prev)) / fact < fabs(current))
  504. {
  505. prev /= current;
  506. Yv_scale /= current;
  507. current = 1;
  508. }
  509. next = fact * current - prev;
  510. prev = current;
  511. current = next;
  512. }
  513. Yv = prev;
  514. }
  515. else
  516. Yv = std::numeric_limits<T>::quiet_NaN(); // any value will do, we're not using it.
  517. }
  518. if (reflect)
  519. {
  520. if((sp != 0) && (tools::max_value<T>() * fabs(Yv_scale) < fabs(sp * Yv)))
  521. *J = org_kind & need_j ? T(-sign(sp) * sign(Yv) * (Yv_scale != 0 ? sign(Yv_scale) : 1) * policies::raise_overflow_error<T>(function, 0, pol)) : T(0);
  522. else
  523. *J = cp * Jv - (sp == 0 ? T(0) : T((sp * Yv) / Yv_scale)); // reflection formula
  524. if((cp != 0) && (tools::max_value<T>() * fabs(Yv_scale) < fabs(cp * Yv)))
  525. *Y = org_kind & need_y ? T(-sign(cp) * sign(Yv) * (Yv_scale != 0 ? sign(Yv_scale) : 1) * policies::raise_overflow_error<T>(function, 0, pol)) : T(0);
  526. else
  527. *Y = (sp != 0 ? sp * Jv : T(0)) + (cp == 0 ? T(0) : T((cp * Yv) / Yv_scale));
  528. }
  529. else
  530. {
  531. *J = Jv;
  532. if(tools::max_value<T>() * fabs(Yv_scale) < fabs(Yv))
  533. *Y = org_kind & need_y ? T(sign(Yv) * sign(Yv_scale) * policies::raise_overflow_error<T>(function, 0, pol)) : T(0);
  534. else
  535. *Y = Yv / Yv_scale;
  536. }
  537. return 0;
  538. }
  539. } // namespace detail
  540. }} // namespaces
  541. #endif // BOOST_MATH_BESSEL_JY_HPP