bessel_jy_zero.hpp 24 KB

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