bessel_jy_zero.hpp 25 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627
  1. // Copyright (c) 2013 Christopher Kormanyos
  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. //
  6. // This work is based on an earlier work:
  7. // "Algorithm 910: A Portable C++ Multiple-Precision System for Special-Function Calculations",
  8. // in ACM TOMS, {VOL 37, ISSUE 4, (February 2011)} (C) ACM, 2011. http://doi.acm.org/10.1145/1916461.1916469
  9. //
  10. // This header contains implementation details for estimating the zeros
  11. // of cylindrical Bessel and Neumann functions on the positive real axis.
  12. // Support is included for both positive as well as negative order.
  13. // Various methods are used to estimate the roots. These include
  14. // empirical curve fitting and McMahon's asymptotic approximation
  15. // for small order, uniform asymptotic expansion for large order,
  16. // and iteration and root interlacing for negative order.
  17. //
  18. #ifndef BOOST_MATH_BESSEL_JY_ZERO_2013_01_18_HPP_
  19. #define BOOST_MATH_BESSEL_JY_ZERO_2013_01_18_HPP_
  20. #include <algorithm>
  21. #include <boost/math/constants/constants.hpp>
  22. #include <boost/math/special_functions/math_fwd.hpp>
  23. #include <boost/math/special_functions/cbrt.hpp>
  24. #include <boost/math/special_functions/detail/airy_ai_bi_zero.hpp>
  25. namespace boost { namespace math {
  26. namespace detail
  27. {
  28. namespace bessel_zero
  29. {
  30. template<class T>
  31. T equation_nist_10_21_19(const T& v, const T& a)
  32. {
  33. // Get the initial estimate of the m'th root of Jv or Yv.
  34. // This subroutine is used for the order m with m > 1.
  35. // The order m has been used to create the input parameter a.
  36. // This is Eq. 10.21.19 in the NIST Handbook.
  37. const T mu = (v * v) * 4U;
  38. const T mu_minus_one = mu - T(1);
  39. const T eight_a_inv = T(1) / (a * 8U);
  40. const T eight_a_inv_squared = eight_a_inv * eight_a_inv;
  41. const T term3 = ((mu_minus_one * 4U) * ((mu * 7U) - T(31U) )) / 3U;
  42. const T term5 = ((mu_minus_one * 32U) * ((((mu * 83U) - T(982U) ) * mu) + T(3779U) )) / 15U;
  43. const T term7 = ((mu_minus_one * 64U) * ((((((mu * 6949U) - T(153855UL)) * mu) + T(1585743UL)) * mu) - T(6277237UL))) / 105U;
  44. return a + (((( - term7
  45. * eight_a_inv_squared - term5)
  46. * eight_a_inv_squared - term3)
  47. * eight_a_inv_squared - mu_minus_one)
  48. * eight_a_inv);
  49. }
  50. template<typename T>
  51. class equation_as_9_3_39_and_its_derivative
  52. {
  53. public:
  54. explicit equation_as_9_3_39_and_its_derivative(const T& zt) : zeta(zt) { }
  55. equation_as_9_3_39_and_its_derivative(const equation_as_9_3_39_and_its_derivative&) = default;
  56. boost::math::tuple<T, T> operator()(const T& z) const
  57. {
  58. BOOST_MATH_STD_USING // ADL of std names, needed for acos, sqrt.
  59. // Return the function of zeta that is implicitly defined
  60. // in A&S Eq. 9.3.39 as a function of z. The function is
  61. // returned along with its derivative with respect to z.
  62. const T zsq_minus_one_sqrt = sqrt((z * z) - T(1));
  63. const T the_function(
  64. zsq_minus_one_sqrt
  65. - ( acos(T(1) / z) + ((T(2) / 3U) * (zeta * sqrt(zeta)))));
  66. const T its_derivative(zsq_minus_one_sqrt / z);
  67. return boost::math::tuple<T, T>(the_function, its_derivative);
  68. }
  69. private:
  70. const equation_as_9_3_39_and_its_derivative& operator=(const equation_as_9_3_39_and_its_derivative&) = delete;
  71. const T zeta;
  72. };
  73. template<class T, class Policy>
  74. static T equation_as_9_5_26(const T& v, const T& ai_bi_root, const Policy& pol)
  75. {
  76. BOOST_MATH_STD_USING // ADL of std names, needed for log, sqrt.
  77. // Obtain the estimate of the m'th zero of Jv or Yv.
  78. // The order m has been used to create the input parameter ai_bi_root.
  79. // Here, v is larger than about 2.2. The estimate is computed
  80. // from Abramowitz and Stegun Eqs. 9.5.22 and 9.5.26, page 371.
  81. //
  82. // The inversion of z as a function of zeta is mentioned in the text
  83. // following A&S Eq. 9.5.26. Here, we accomplish the inversion by
  84. // performing a Taylor expansion of Eq. 9.3.39 for large z to order 2
  85. // and solving the resulting quadratic equation, thereby taking
  86. // the positive root of the quadratic.
  87. // In other words: (2/3)(-zeta)^(3/2) approx = z + 1/(2z) - pi/2.
  88. // This leads to: z^2 - [(2/3)(-zeta)^(3/2) + pi/2]z + 1/2 = 0.
  89. //
  90. // With this initial estimate, Newton-Raphson iteration is used
  91. // to refine the value of the estimate of the root of z
  92. // as a function of zeta.
  93. const T v_pow_third(boost::math::cbrt(v, pol));
  94. const T v_pow_minus_two_thirds(T(1) / (v_pow_third * v_pow_third));
  95. // Obtain zeta using the order v combined with the m'th root of
  96. // an airy function, as shown in A&S Eq. 9.5.22.
  97. const T zeta = v_pow_minus_two_thirds * (-ai_bi_root);
  98. const T zeta_sqrt = sqrt(zeta);
  99. // Set up a quadratic equation based on the Taylor series
  100. // expansion mentioned above.
  101. const T b = -((((zeta * zeta_sqrt) * 2U) / 3U) + boost::math::constants::half_pi<T>());
  102. // Solve the quadratic equation, taking the positive root.
  103. const T z_estimate = (-b + sqrt((b * b) - T(2))) / 2U;
  104. // Establish the range, the digits, and the iteration limit
  105. // for the upcoming root-finding.
  106. const T range_zmin = (std::max<T>)(z_estimate - T(1), T(1));
  107. const T range_zmax = z_estimate + T(1);
  108. const auto my_digits10 = static_cast<int>(static_cast<float>(boost::math::tools::digits<T>() * 0.301F));
  109. // Select the maximum allowed iterations based on the number
  110. // of decimal digits in the numeric type T, being at least 12.
  111. const auto iterations_allowed = static_cast<std::uintmax_t>((std::max)(12, my_digits10 * 2));
  112. std::uintmax_t iterations_used = iterations_allowed;
  113. // Calculate the root of z as a function of zeta.
  114. const T z = boost::math::tools::newton_raphson_iterate(
  115. boost::math::detail::bessel_zero::equation_as_9_3_39_and_its_derivative<T>(zeta),
  116. z_estimate,
  117. range_zmin,
  118. range_zmax,
  119. (std::min)(boost::math::tools::digits<T>(), boost::math::tools::digits<float>()),
  120. iterations_used);
  121. static_cast<void>(iterations_used);
  122. // Continue with the implementation of A&S Eq. 9.3.39.
  123. const T zsq_minus_one = (z * z) - T(1);
  124. const T zsq_minus_one_sqrt = sqrt(zsq_minus_one);
  125. // This is A&S Eq. 9.3.42.
  126. const T b0_term_5_24 = T(5) / ((zsq_minus_one * zsq_minus_one_sqrt) * 24U);
  127. const T b0_term_1_8 = T(1) / ( zsq_minus_one_sqrt * 8U);
  128. const T b0_term_5_48 = T(5) / ((zeta * zeta) * 48U);
  129. const T b0 = -b0_term_5_48 + ((b0_term_5_24 + b0_term_1_8) / zeta_sqrt);
  130. // This is the second line of A&S Eq. 9.5.26 for f_k with k = 1.
  131. const T f1 = ((z * zeta_sqrt) * b0) / zsq_minus_one_sqrt;
  132. // This is A&S Eq. 9.5.22 expanded to k = 1 (i.e., one term in the series).
  133. return (v * z) + (f1 / v);
  134. }
  135. namespace cyl_bessel_j_zero_detail
  136. {
  137. template<class T, class Policy>
  138. T equation_nist_10_21_40_a(const T& v, const Policy& pol)
  139. {
  140. const T v_pow_third(boost::math::cbrt(v, pol));
  141. const T v_pow_minus_two_thirds(T(1) / (v_pow_third * v_pow_third));
  142. return v * ((((( + T(0.043)
  143. * v_pow_minus_two_thirds - T(0.0908))
  144. * v_pow_minus_two_thirds - T(0.00397))
  145. * v_pow_minus_two_thirds + T(1.033150))
  146. * v_pow_minus_two_thirds + T(1.8557571))
  147. * v_pow_minus_two_thirds + T(1));
  148. }
  149. template<class T, class Policy>
  150. class function_object_jv
  151. {
  152. public:
  153. function_object_jv(const T& v,
  154. const Policy& pol) : my_v(v),
  155. my_pol(pol) { }
  156. function_object_jv(const function_object_jv&) = default;
  157. T operator()(const T& x) const
  158. {
  159. return boost::math::cyl_bessel_j(my_v, x, my_pol);
  160. }
  161. private:
  162. const T my_v;
  163. const Policy& my_pol;
  164. const function_object_jv& operator=(const function_object_jv&) = delete;
  165. };
  166. template<class T, class Policy>
  167. class function_object_jv_and_jv_prime
  168. {
  169. public:
  170. function_object_jv_and_jv_prime(const T& v,
  171. const bool order_is_zero,
  172. const Policy& pol) : my_v(v),
  173. my_order_is_zero(order_is_zero),
  174. my_pol(pol) { }
  175. function_object_jv_and_jv_prime(const function_object_jv_and_jv_prime&) = default;
  176. boost::math::tuple<T, T> operator()(const T& x) const
  177. {
  178. // Obtain Jv(x) and Jv'(x).
  179. // Chris's original code called the Bessel function implementation layer direct,
  180. // but that circumvented optimizations for integer-orders. Call the documented
  181. // top level functions instead, and let them sort out which implementation to use.
  182. T j_v;
  183. T j_v_prime;
  184. if(my_order_is_zero)
  185. {
  186. j_v = boost::math::cyl_bessel_j(0, x, my_pol);
  187. j_v_prime = -boost::math::cyl_bessel_j(1, x, my_pol);
  188. }
  189. else
  190. {
  191. j_v = boost::math::cyl_bessel_j( my_v, x, my_pol);
  192. const T j_v_m1 (boost::math::cyl_bessel_j(T(my_v - 1), x, my_pol));
  193. j_v_prime = j_v_m1 - ((my_v * j_v) / x);
  194. }
  195. // Return a tuple containing both Jv(x) and Jv'(x).
  196. return boost::math::make_tuple(j_v, j_v_prime);
  197. }
  198. private:
  199. const T my_v;
  200. const bool my_order_is_zero;
  201. const Policy& my_pol;
  202. const function_object_jv_and_jv_prime& operator=(const function_object_jv_and_jv_prime&) = delete;
  203. };
  204. template<class T> bool my_bisection_unreachable_tolerance(const T&, const T&) { return false; }
  205. template<class T, class Policy>
  206. T initial_guess(const T& v, const int m, const Policy& pol)
  207. {
  208. BOOST_MATH_STD_USING // ADL of std names, needed for floor.
  209. // Compute an estimate of the m'th root of cyl_bessel_j.
  210. T guess;
  211. // There is special handling for negative order.
  212. if(v < 0)
  213. {
  214. if((m == 1) && (v > -0.5F))
  215. {
  216. // For small, negative v, use the results of empirical curve fitting.
  217. // Mathematica(R) session for the coefficients:
  218. // Table[{n, BesselJZero[n, 1]}, {n, -(1/2), 0, 1/10}]
  219. // N[%, 20]
  220. // Fit[%, {n^0, n^1, n^2, n^3, n^4, n^5, n^6}, n]
  221. guess = ((((( - T(0.2321156900729)
  222. * v - T(0.1493247777488))
  223. * v - T(0.15205419167239))
  224. * v + T(0.07814930561249))
  225. * v - T(0.17757573537688))
  226. * v + T(1.542805677045663))
  227. * v + T(2.40482555769577277);
  228. return guess;
  229. }
  230. // Create the positive order and extract its positive floor integer part.
  231. const T vv(-v);
  232. const T vv_floor(floor(vv));
  233. // The to-be-found root is bracketed by the roots of the
  234. // Bessel function whose reflected, positive integer order
  235. // is less than, but nearest to vv.
  236. T root_hi = boost::math::detail::bessel_zero::cyl_bessel_j_zero_detail::initial_guess(vv_floor, m, pol);
  237. T root_lo;
  238. if(m == 1)
  239. {
  240. // The estimate of the first root for negative order is found using
  241. // an adaptive range-searching algorithm.
  242. root_lo = T(root_hi - 0.1F);
  243. const bool hi_end_of_bracket_is_negative = (boost::math::cyl_bessel_j(v, root_hi, pol) < 0);
  244. while((root_lo > boost::math::tools::epsilon<T>()))
  245. {
  246. const bool lo_end_of_bracket_is_negative = (boost::math::cyl_bessel_j(v, root_lo, pol) < 0);
  247. if(hi_end_of_bracket_is_negative != lo_end_of_bracket_is_negative)
  248. {
  249. break;
  250. }
  251. root_hi = root_lo;
  252. // Decrease the lower end of the bracket using an adaptive algorithm.
  253. if(root_lo > 0.5F)
  254. {
  255. root_lo -= 0.5F;
  256. }
  257. else
  258. {
  259. root_lo *= 0.75F; // LCOV_EXCL_LINE probably unreachable, but hard to prove?
  260. }
  261. }
  262. }
  263. else
  264. {
  265. root_lo = boost::math::detail::bessel_zero::cyl_bessel_j_zero_detail::initial_guess(vv_floor, m - 1, pol);
  266. }
  267. // Perform several steps of bisection iteration to refine the guess.
  268. std::uintmax_t number_of_iterations(12U);
  269. // Do the bisection iteration.
  270. const boost::math::tuple<T, T> guess_pair =
  271. boost::math::tools::bisect(
  272. boost::math::detail::bessel_zero::cyl_bessel_j_zero_detail::function_object_jv<T, Policy>(v, pol),
  273. root_lo,
  274. root_hi,
  275. boost::math::detail::bessel_zero::cyl_bessel_j_zero_detail::my_bisection_unreachable_tolerance<T>,
  276. number_of_iterations);
  277. return (boost::math::get<0>(guess_pair) + boost::math::get<1>(guess_pair)) / 2U;
  278. }
  279. if(m == 1U)
  280. {
  281. // Get the initial estimate of the first root.
  282. if(v < 2.2F)
  283. {
  284. // For small v, use the results of empirical curve fitting.
  285. // Mathematica(R) session for the coefficients:
  286. // Table[{n, BesselJZero[n, 1]}, {n, 0, 22/10, 1/10}]
  287. // N[%, 20]
  288. // Fit[%, {n^0, n^1, n^2, n^3, n^4, n^5, n^6}, n]
  289. guess = ((((( - T(0.0008342379046010)
  290. * v + T(0.007590035637410))
  291. * v - T(0.030640914772013))
  292. * v + T(0.078232088020106))
  293. * v - T(0.169668712590620))
  294. * v + T(1.542187960073750))
  295. * v + T(2.4048359915254634);
  296. }
  297. else
  298. {
  299. // For larger v, use the first line of Eqs. 10.21.40 in the NIST Handbook.
  300. guess = boost::math::detail::bessel_zero::cyl_bessel_j_zero_detail::equation_nist_10_21_40_a(v, pol);
  301. }
  302. }
  303. else
  304. {
  305. if(v < 2.2F)
  306. {
  307. // Use Eq. 10.21.19 in the NIST Handbook.
  308. const T a(((v + T(m * 2U)) - T(0.5)) * boost::math::constants::half_pi<T>());
  309. guess = boost::math::detail::bessel_zero::equation_nist_10_21_19(v, a);
  310. }
  311. else
  312. {
  313. // Get an estimate of the m'th root of airy_ai.
  314. const T airy_ai_root(boost::math::detail::airy_zero::airy_ai_zero_detail::initial_guess<T>(m, pol));
  315. // Use Eq. 9.5.26 in the A&S Handbook.
  316. guess = boost::math::detail::bessel_zero::equation_as_9_5_26(v, airy_ai_root, pol);
  317. }
  318. }
  319. return guess;
  320. }
  321. } // namespace cyl_bessel_j_zero_detail
  322. namespace cyl_neumann_zero_detail
  323. {
  324. template<class T, class Policy>
  325. T equation_nist_10_21_40_b(const T& v, const Policy& pol)
  326. {
  327. const T v_pow_third(boost::math::cbrt(v, pol));
  328. const T v_pow_minus_two_thirds(T(1) / (v_pow_third * v_pow_third));
  329. return v * ((((( - T(0.001)
  330. * v_pow_minus_two_thirds - T(0.0060))
  331. * v_pow_minus_two_thirds + T(0.01198))
  332. * v_pow_minus_two_thirds + T(0.260351))
  333. * v_pow_minus_two_thirds + T(0.9315768))
  334. * v_pow_minus_two_thirds + T(1));
  335. }
  336. template<class T, class Policy>
  337. class function_object_yv
  338. {
  339. public:
  340. function_object_yv(const T& v,
  341. const Policy& pol) : my_v(v),
  342. my_pol(pol) { }
  343. function_object_yv(const function_object_yv&) = default;
  344. T operator()(const T& x) const
  345. {
  346. return boost::math::cyl_neumann(my_v, x, my_pol);
  347. }
  348. private:
  349. const T my_v;
  350. const Policy& my_pol;
  351. const function_object_yv& operator=(const function_object_yv&) = delete;
  352. };
  353. template<class T, class Policy>
  354. class function_object_yv_and_yv_prime
  355. {
  356. public:
  357. function_object_yv_and_yv_prime(const T& v,
  358. const Policy& pol) : my_v(v),
  359. my_pol(pol) { }
  360. function_object_yv_and_yv_prime(const function_object_yv_and_yv_prime&) = default;
  361. boost::math::tuple<T, T> operator()(const T& x) const
  362. {
  363. const T half_epsilon(boost::math::tools::epsilon<T>() / 2U);
  364. const bool order_is_zero = ((my_v > -half_epsilon) && (my_v < +half_epsilon));
  365. // Obtain Yv(x) and Yv'(x).
  366. // Chris's original code called the Bessel function implementation layer direct,
  367. // but that circumvented optimizations for integer-orders. Call the documented
  368. // top level functions instead, and let them sort out which implementation to use.
  369. T y_v;
  370. T y_v_prime;
  371. if(order_is_zero)
  372. {
  373. y_v = boost::math::cyl_neumann(0, x, my_pol);
  374. y_v_prime = -boost::math::cyl_neumann(1, x, my_pol);
  375. }
  376. else
  377. {
  378. y_v = boost::math::cyl_neumann( my_v, x, my_pol);
  379. const T y_v_m1 (boost::math::cyl_neumann(T(my_v - 1), x, my_pol));
  380. y_v_prime = y_v_m1 - ((my_v * y_v) / x);
  381. }
  382. // Return a tuple containing both Yv(x) and Yv'(x).
  383. return boost::math::make_tuple(y_v, y_v_prime);
  384. }
  385. private:
  386. const T my_v;
  387. const Policy& my_pol;
  388. const function_object_yv_and_yv_prime& operator=(const function_object_yv_and_yv_prime&) = delete;
  389. };
  390. template<class T> bool my_bisection_unreachable_tolerance(const T&, const T&) { return false; }
  391. template<class T, class Policy>
  392. T initial_guess(const T& v, const int m, const Policy& pol)
  393. {
  394. BOOST_MATH_STD_USING // ADL of std names, needed for floor.
  395. // Compute an estimate of the m'th root of cyl_neumann.
  396. T guess;
  397. // There is special handling for negative order.
  398. if(v < 0)
  399. {
  400. // Create the positive order and extract its positive floor and ceiling integer parts.
  401. const T vv(-v);
  402. const T vv_floor(floor(vv));
  403. // The to-be-found root is bracketed by the roots of the
  404. // Bessel function whose reflected, positive integer order
  405. // is less than, but nearest to vv.
  406. // The special case of negative, half-integer order uses
  407. // the relation between Yv and spherical Bessel functions
  408. // in order to obtain the bracket for the root.
  409. // In these special cases, cyl_neumann(-n/2, x) = sph_bessel_j(+n/2, x)
  410. // for v = -n/2.
  411. T root_hi;
  412. T root_lo;
  413. if(m == 1)
  414. {
  415. // The estimate of the first root for negative order is found using
  416. // an adaptive range-searching algorithm.
  417. // Take special precautions for the discontinuity at negative,
  418. // half-integer orders and use different brackets above and below these.
  419. if(T(vv - vv_floor) < 0.5F)
  420. {
  421. root_hi = boost::math::detail::bessel_zero::cyl_neumann_zero_detail::initial_guess(vv_floor, m, pol);
  422. }
  423. else
  424. {
  425. root_hi = boost::math::detail::bessel_zero::cyl_bessel_j_zero_detail::initial_guess(T(vv_floor + 0.5F), m, pol);
  426. }
  427. root_lo = T(root_hi - 0.1F);
  428. const bool hi_end_of_bracket_is_negative = (boost::math::cyl_neumann(v, root_hi, pol) < 0);
  429. while((root_lo > boost::math::tools::epsilon<T>()))
  430. {
  431. const bool lo_end_of_bracket_is_negative = (boost::math::cyl_neumann(v, root_lo, pol) < 0);
  432. if(hi_end_of_bracket_is_negative != lo_end_of_bracket_is_negative)
  433. {
  434. break;
  435. }
  436. root_hi = root_lo;
  437. // Decrease the lower end of the bracket using an adaptive algorithm.
  438. if(root_lo > 0.5F)
  439. {
  440. root_lo -= 0.5F;
  441. }
  442. else
  443. {
  444. root_lo *= 0.75F; // LCOV_EXCL_LINE probably unreachable, but hard to prove?
  445. }
  446. }
  447. }
  448. else
  449. {
  450. if(T(vv - vv_floor) < 0.5F)
  451. {
  452. root_lo = boost::math::detail::bessel_zero::cyl_neumann_zero_detail::initial_guess(vv_floor, m - 1, pol);
  453. root_hi = boost::math::detail::bessel_zero::cyl_neumann_zero_detail::initial_guess(vv_floor, m, pol);
  454. root_lo += 0.01F;
  455. root_hi += 0.01F;
  456. }
  457. else
  458. {
  459. root_lo = boost::math::detail::bessel_zero::cyl_bessel_j_zero_detail::initial_guess(T(vv_floor + 0.5F), m - 1, pol);
  460. root_hi = boost::math::detail::bessel_zero::cyl_bessel_j_zero_detail::initial_guess(T(vv_floor + 0.5F), m, pol);
  461. root_lo += 0.01F;
  462. root_hi += 0.01F;
  463. }
  464. }
  465. // Perform several steps of bisection iteration to refine the guess.
  466. std::uintmax_t number_of_iterations(12U);
  467. // Do the bisection iteration.
  468. const boost::math::tuple<T, T> guess_pair =
  469. boost::math::tools::bisect(
  470. boost::math::detail::bessel_zero::cyl_neumann_zero_detail::function_object_yv<T, Policy>(v, pol),
  471. root_lo,
  472. root_hi,
  473. boost::math::detail::bessel_zero::cyl_neumann_zero_detail::my_bisection_unreachable_tolerance<T>,
  474. number_of_iterations);
  475. return (boost::math::get<0>(guess_pair) + boost::math::get<1>(guess_pair)) / 2U;
  476. }
  477. if(m == 1U)
  478. {
  479. // Get the initial estimate of the first root.
  480. if(v < 2.2F)
  481. {
  482. // For small v, use the results of empirical curve fitting.
  483. // Mathematica(R) session for the coefficients:
  484. // Table[{n, BesselYZero[n, 1]}, {n, 0, 22/10, 1/10}]
  485. // N[%, 20]
  486. // Fit[%, {n^0, n^1, n^2, n^3, n^4, n^5, n^6}, n]
  487. guess = ((((( - T(0.0025095909235652)
  488. * v + T(0.021291887049053))
  489. * v - T(0.076487785486526))
  490. * v + T(0.159110268115362))
  491. * v - T(0.241681668765196))
  492. * v + T(1.4437846310885244))
  493. * v + T(0.89362115190200490);
  494. }
  495. else
  496. {
  497. // For larger v, use the second line of Eqs. 10.21.40 in the NIST Handbook.
  498. guess = boost::math::detail::bessel_zero::cyl_neumann_zero_detail::equation_nist_10_21_40_b(v, pol);
  499. }
  500. }
  501. else
  502. {
  503. if(v < 2.2F)
  504. {
  505. // Use Eq. 10.21.19 in the NIST Handbook.
  506. const T a(((v + T(m * 2U)) - T(1.5)) * boost::math::constants::half_pi<T>());
  507. guess = boost::math::detail::bessel_zero::equation_nist_10_21_19(v, a);
  508. }
  509. else
  510. {
  511. // Get an estimate of the m'th root of airy_bi.
  512. const T airy_bi_root(boost::math::detail::airy_zero::airy_bi_zero_detail::initial_guess<T>(m, pol));
  513. // Use Eq. 9.5.26 in the A&S Handbook.
  514. guess = boost::math::detail::bessel_zero::equation_as_9_5_26(v, airy_bi_root, pol);
  515. }
  516. }
  517. return guess;
  518. }
  519. } // namespace cyl_neumann_zero_detail
  520. } // namespace bessel_zero
  521. } } } // namespace boost::math::detail
  522. #endif // BOOST_MATH_BESSEL_JY_ZERO_2013_01_18_HPP_