bessel_jy.hpp 22 KB

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