Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 09:40:00

0001 //  Copyright (c) 2013 Christopher Kormanyos
0002 //  Use, modification and distribution are subject to the
0003 //  Boost Software License, Version 1.0. (See accompanying file
0004 //  LICENSE_1_0.txt or copy at http://www.boost.org/LICENSE_1_0.txt)
0005 //
0006 // This work is based on an earlier work:
0007 // "Algorithm 910: A Portable C++ Multiple-Precision System for Special-Function Calculations",
0008 // in ACM TOMS, {VOL 37, ISSUE 4, (February 2011)} (C) ACM, 2011. http://doi.acm.org/10.1145/1916461.1916469
0009 //
0010 // This header contains implementation details for estimating the zeros
0011 // of cylindrical Bessel and Neumann functions on the positive real axis.
0012 // Support is included for both positive as well as negative order.
0013 // Various methods are used to estimate the roots. These include
0014 // empirical curve fitting and McMahon's asymptotic approximation
0015 // for small order, uniform asymptotic expansion for large order,
0016 // and iteration and root interlacing for negative order.
0017 //
0018 #ifndef BOOST_MATH_BESSEL_JY_ZERO_2013_01_18_HPP_
0019   #define BOOST_MATH_BESSEL_JY_ZERO_2013_01_18_HPP_
0020 
0021   #include <algorithm>
0022   #include <boost/math/constants/constants.hpp>
0023   #include <boost/math/special_functions/math_fwd.hpp>
0024   #include <boost/math/special_functions/cbrt.hpp>
0025   #include <boost/math/special_functions/detail/airy_ai_bi_zero.hpp>
0026 
0027   namespace boost { namespace math {
0028   namespace detail
0029   {
0030     namespace bessel_zero
0031     {
0032       template<class T>
0033       T equation_nist_10_21_19(const T& v, const T& a)
0034       {
0035         // Get the initial estimate of the m'th root of Jv or Yv.
0036         // This subroutine is used for the order m with m > 1.
0037         // The order m has been used to create the input parameter a.
0038 
0039         // This is Eq. 10.21.19 in the NIST Handbook.
0040         const T mu                  = (v * v) * 4U;
0041         const T mu_minus_one        = mu - T(1);
0042         const T eight_a_inv         = T(1) / (a * 8U);
0043         const T eight_a_inv_squared = eight_a_inv * eight_a_inv;
0044 
0045         const T term3 = ((mu_minus_one *  4U) *     ((mu *    7U) -     T(31U) )) / 3U;
0046         const T term5 = ((mu_minus_one * 32U) *   ((((mu *   83U) -    T(982U) ) * mu) +    T(3779U) )) / 15U;
0047         const T term7 = ((mu_minus_one * 64U) * ((((((mu * 6949U) - T(153855UL)) * mu) + T(1585743UL)) * mu) - T(6277237UL))) / 105U;
0048 
0049         return a + ((((                      - term7
0050                        * eight_a_inv_squared - term5)
0051                        * eight_a_inv_squared - term3)
0052                        * eight_a_inv_squared - mu_minus_one)
0053                        * eight_a_inv);
0054       }
0055 
0056       template<typename T>
0057       class equation_as_9_3_39_and_its_derivative
0058       {
0059       public:
0060         explicit equation_as_9_3_39_and_its_derivative(const T& zt) : zeta(zt) { }
0061 
0062         equation_as_9_3_39_and_its_derivative(const equation_as_9_3_39_and_its_derivative&) = default;
0063 
0064         boost::math::tuple<T, T> operator()(const T& z) const
0065         {
0066           BOOST_MATH_STD_USING // ADL of std names, needed for acos, sqrt.
0067 
0068           // Return the function of zeta that is implicitly defined
0069           // in A&S Eq. 9.3.39 as a function of z. The function is
0070           // returned along with its derivative with respect to z.
0071 
0072           const T zsq_minus_one_sqrt = sqrt((z * z) - T(1));
0073 
0074           const T the_function(
0075               zsq_minus_one_sqrt
0076             - (  acos(T(1) / z) + ((T(2) / 3U) * (zeta * sqrt(zeta)))));
0077 
0078           const T its_derivative(zsq_minus_one_sqrt / z);
0079 
0080           return boost::math::tuple<T, T>(the_function, its_derivative);
0081         }
0082 
0083       private:
0084         const equation_as_9_3_39_and_its_derivative& operator=(const equation_as_9_3_39_and_its_derivative&) = delete;
0085         const T zeta;
0086       };
0087 
0088       template<class T, class Policy>
0089       static T equation_as_9_5_26(const T& v, const T& ai_bi_root, const Policy& pol)
0090       {
0091         BOOST_MATH_STD_USING // ADL of std names, needed for log, sqrt.
0092 
0093         // Obtain the estimate of the m'th zero of Jv or Yv.
0094         // The order m has been used to create the input parameter ai_bi_root.
0095         // Here, v is larger than about 2.2. The estimate is computed
0096         // from Abramowitz and Stegun Eqs. 9.5.22 and 9.5.26, page 371.
0097         //
0098         // The inversion of z as a function of zeta is mentioned in the text
0099         // following A&S Eq. 9.5.26. Here, we accomplish the inversion by
0100         // performing a Taylor expansion of Eq. 9.3.39 for large z to order 2
0101         // and solving the resulting quadratic equation, thereby taking
0102         // the positive root of the quadratic.
0103         // In other words: (2/3)(-zeta)^(3/2) approx = z + 1/(2z) - pi/2.
0104         // This leads to: z^2 - [(2/3)(-zeta)^(3/2) + pi/2]z + 1/2 = 0.
0105         //
0106         // With this initial estimate, Newton-Raphson iteration is used
0107         // to refine the value of the estimate of the root of z
0108         // as a function of zeta.
0109 
0110         const T v_pow_third(boost::math::cbrt(v, pol));
0111         const T v_pow_minus_two_thirds(T(1) / (v_pow_third * v_pow_third));
0112 
0113         // Obtain zeta using the order v combined with the m'th root of
0114         // an airy function, as shown in  A&S Eq. 9.5.22.
0115         const T zeta = v_pow_minus_two_thirds * (-ai_bi_root);
0116 
0117         const T zeta_sqrt = sqrt(zeta);
0118 
0119         // Set up a quadratic equation based on the Taylor series
0120         // expansion mentioned above.
0121         const T b = -((((zeta * zeta_sqrt) * 2U) / 3U) + boost::math::constants::half_pi<T>());
0122 
0123         // Solve the quadratic equation, taking the positive root.
0124         const T z_estimate = (-b + sqrt((b * b) - T(2))) / 2U;
0125 
0126         // Establish the range, the digits, and the iteration limit
0127         // for the upcoming root-finding.
0128         const T range_zmin = (std::max<T>)(z_estimate - T(1), T(1));
0129         const T range_zmax = z_estimate + T(1);
0130 
0131         const auto my_digits10 = static_cast<int>(static_cast<float>(boost::math::tools::digits<T>() * 0.301F));
0132 
0133         // Select the maximum allowed iterations based on the number
0134         // of decimal digits in the numeric type T, being at least 12.
0135         const auto iterations_allowed = static_cast<std::uintmax_t>((std::max)(12, my_digits10 * 2));
0136 
0137         std::uintmax_t iterations_used = iterations_allowed;
0138 
0139         // Calculate the root of z as a function of zeta.
0140         const T z = boost::math::tools::newton_raphson_iterate(
0141           boost::math::detail::bessel_zero::equation_as_9_3_39_and_its_derivative<T>(zeta),
0142           z_estimate,
0143           range_zmin,
0144           range_zmax,
0145           (std::min)(boost::math::tools::digits<T>(), boost::math::tools::digits<float>()),
0146           iterations_used);
0147 
0148         static_cast<void>(iterations_used);
0149 
0150         // Continue with the implementation of A&S Eq. 9.3.39.
0151         const T zsq_minus_one      = (z * z) - T(1);
0152         const T zsq_minus_one_sqrt = sqrt(zsq_minus_one);
0153 
0154         // This is A&S Eq. 9.3.42.
0155         const T b0_term_5_24 = T(5) / ((zsq_minus_one * zsq_minus_one_sqrt) * 24U);
0156         const T b0_term_1_8  = T(1) / ( zsq_minus_one_sqrt * 8U);
0157         const T b0_term_5_48 = T(5) / ((zeta * zeta) * 48U);
0158 
0159         const T b0 = -b0_term_5_48 + ((b0_term_5_24 + b0_term_1_8) / zeta_sqrt);
0160 
0161         // This is the second line of A&S Eq. 9.5.26 for f_k with k = 1.
0162         const T f1 = ((z * zeta_sqrt) * b0) / zsq_minus_one_sqrt;
0163 
0164         // This is A&S Eq. 9.5.22 expanded to k = 1 (i.e., one term in the series).
0165         return (v * z) + (f1 / v);
0166       }
0167 
0168       namespace cyl_bessel_j_zero_detail
0169       {
0170         template<class T, class Policy>
0171         T equation_nist_10_21_40_a(const T& v, const Policy& pol)
0172         {
0173           const T v_pow_third(boost::math::cbrt(v, pol));
0174           const T v_pow_minus_two_thirds(T(1) / (v_pow_third * v_pow_third));
0175 
0176           return v * (((((                         + T(0.043)
0177                           * v_pow_minus_two_thirds - T(0.0908))
0178                           * v_pow_minus_two_thirds - T(0.00397))
0179                           * v_pow_minus_two_thirds + T(1.033150))
0180                           * v_pow_minus_two_thirds + T(1.8557571))
0181                           * v_pow_minus_two_thirds + T(1));
0182         }
0183 
0184         template<class T, class Policy>
0185         class function_object_jv
0186         {
0187         public:
0188           function_object_jv(const T& v,
0189                              const Policy& pol) : my_v(v),
0190                                                   my_pol(pol) { }
0191 
0192           function_object_jv(const function_object_jv&) = default;
0193 
0194           T operator()(const T& x) const
0195           {
0196             return boost::math::cyl_bessel_j(my_v, x, my_pol);
0197           }
0198 
0199         private:
0200           const T my_v;
0201           const Policy& my_pol;
0202           const function_object_jv& operator=(const function_object_jv&) = delete;
0203         };
0204 
0205         template<class T, class Policy>
0206         class function_object_jv_and_jv_prime
0207         {
0208         public:
0209           function_object_jv_and_jv_prime(const T& v,
0210                                           const bool order_is_zero,
0211                                           const Policy& pol) : my_v(v),
0212                                                                my_order_is_zero(order_is_zero),
0213                                                                my_pol(pol) { }
0214 
0215           function_object_jv_and_jv_prime(const function_object_jv_and_jv_prime&) = default;
0216 
0217           boost::math::tuple<T, T> operator()(const T& x) const
0218           {
0219             // Obtain Jv(x) and Jv'(x).
0220             // Chris's original code called the Bessel function implementation layer direct, 
0221             // but that circumvented optimizations for integer-orders.  Call the documented
0222             // top level functions instead, and let them sort out which implementation to use.
0223             T j_v;
0224             T j_v_prime;
0225 
0226             if(my_order_is_zero)
0227             {
0228               j_v       =  boost::math::cyl_bessel_j(0, x, my_pol);
0229               j_v_prime = -boost::math::cyl_bessel_j(1, x, my_pol);
0230             }
0231             else
0232             {
0233                       j_v       = boost::math::cyl_bessel_j(  my_v,      x, my_pol);
0234               const T j_v_m1     (boost::math::cyl_bessel_j(T(my_v - 1), x, my_pol));
0235                       j_v_prime = j_v_m1 - ((my_v * j_v) / x);
0236             }
0237 
0238             // Return a tuple containing both Jv(x) and Jv'(x).
0239             return boost::math::make_tuple(j_v, j_v_prime);
0240           }
0241 
0242         private:
0243           const T my_v;
0244           const bool my_order_is_zero;
0245           const Policy& my_pol;
0246           const function_object_jv_and_jv_prime& operator=(const function_object_jv_and_jv_prime&) = delete;
0247         };
0248 
0249         template<class T> bool my_bisection_unreachable_tolerance(const T&, const T&) { return false; }
0250 
0251         template<class T, class Policy>
0252         T initial_guess(const T& v, const int m, const Policy& pol)
0253         {
0254           BOOST_MATH_STD_USING // ADL of std names, needed for floor.
0255 
0256           // Compute an estimate of the m'th root of cyl_bessel_j.
0257 
0258           T guess;
0259 
0260           // There is special handling for negative order.
0261           if(v < 0)
0262           {
0263             if((m == 1) && (v > -0.5F))
0264             {
0265               // For small, negative v, use the results of empirical curve fitting.
0266               // Mathematica(R) session for the coefficients:
0267               //  Table[{n, BesselJZero[n, 1]}, {n, -(1/2), 0, 1/10}]
0268               //  N[%, 20]
0269               //  Fit[%, {n^0, n^1, n^2, n^3, n^4, n^5, n^6}, n]
0270               guess = (((((    - T(0.2321156900729)
0271                            * v - T(0.1493247777488))
0272                            * v - T(0.15205419167239))
0273                            * v + T(0.07814930561249))
0274                            * v - T(0.17757573537688))
0275                            * v + T(1.542805677045663))
0276                            * v + T(2.40482555769577277);
0277 
0278               return guess;
0279             }
0280 
0281             // Create the positive order and extract its positive floor integer part.
0282             const T vv(-v);
0283             const T vv_floor(floor(vv));
0284 
0285             // The to-be-found root is bracketed by the roots of the
0286             // Bessel function whose reflected, positive integer order
0287             // is less than, but nearest to vv.
0288 
0289             T root_hi = boost::math::detail::bessel_zero::cyl_bessel_j_zero_detail::initial_guess(vv_floor, m, pol);
0290             T root_lo;
0291 
0292             if(m == 1)
0293             {
0294               // The estimate of the first root for negative order is found using
0295               // an adaptive range-searching algorithm.
0296               root_lo = T(root_hi - 0.1F);
0297 
0298               const bool hi_end_of_bracket_is_negative = (boost::math::cyl_bessel_j(v, root_hi, pol) < 0);
0299 
0300               while((root_lo > boost::math::tools::epsilon<T>()))
0301               {
0302                 const bool lo_end_of_bracket_is_negative = (boost::math::cyl_bessel_j(v, root_lo, pol) < 0);
0303 
0304                 if(hi_end_of_bracket_is_negative != lo_end_of_bracket_is_negative)
0305                 {
0306                   break;
0307                 }
0308 
0309                 root_hi = root_lo;
0310 
0311                 // Decrease the lower end of the bracket using an adaptive algorithm.
0312                 if(root_lo > 0.5F)
0313                 {
0314                   root_lo -= 0.5F;
0315                 }
0316                 else
0317                 {
0318                   root_lo *= 0.75F;
0319                 }
0320               }
0321             }
0322             else
0323             {
0324               root_lo = boost::math::detail::bessel_zero::cyl_bessel_j_zero_detail::initial_guess(vv_floor, m - 1, pol);
0325             }
0326 
0327             // Perform several steps of bisection iteration to refine the guess.
0328             std::uintmax_t number_of_iterations(12U);
0329 
0330             // Do the bisection iteration.
0331             const boost::math::tuple<T, T> guess_pair =
0332                boost::math::tools::bisect(
0333                   boost::math::detail::bessel_zero::cyl_bessel_j_zero_detail::function_object_jv<T, Policy>(v, pol),
0334                   root_lo,
0335                   root_hi,
0336                   boost::math::detail::bessel_zero::cyl_bessel_j_zero_detail::my_bisection_unreachable_tolerance<T>,
0337                   number_of_iterations);
0338 
0339             return (boost::math::get<0>(guess_pair) + boost::math::get<1>(guess_pair)) / 2U;
0340           }
0341 
0342           if(m == 1U)
0343           {
0344             // Get the initial estimate of the first root.
0345 
0346             if(v < 2.2F)
0347             {
0348               // For small v, use the results of empirical curve fitting.
0349               // Mathematica(R) session for the coefficients:
0350               //  Table[{n, BesselJZero[n, 1]}, {n, 0, 22/10, 1/10}]
0351               //  N[%, 20]
0352               //  Fit[%, {n^0, n^1, n^2, n^3, n^4, n^5, n^6}, n]
0353               guess = (((((    - T(0.0008342379046010)
0354                            * v + T(0.007590035637410))
0355                            * v - T(0.030640914772013))
0356                            * v + T(0.078232088020106))
0357                            * v - T(0.169668712590620))
0358                            * v + T(1.542187960073750))
0359                            * v + T(2.4048359915254634);
0360             }
0361             else
0362             {
0363               // For larger v, use the first line of Eqs. 10.21.40 in the NIST Handbook.
0364               guess = boost::math::detail::bessel_zero::cyl_bessel_j_zero_detail::equation_nist_10_21_40_a(v, pol);
0365             }
0366           }
0367           else
0368           {
0369             if(v < 2.2F)
0370             {
0371               // Use Eq. 10.21.19 in the NIST Handbook.
0372               const T a(((v + T(m * 2U)) - T(0.5)) * boost::math::constants::half_pi<T>());
0373 
0374               guess = boost::math::detail::bessel_zero::equation_nist_10_21_19(v, a);
0375             }
0376             else
0377             {
0378               // Get an estimate of the m'th root of airy_ai.
0379               const T airy_ai_root(boost::math::detail::airy_zero::airy_ai_zero_detail::initial_guess<T>(m, pol));
0380 
0381               // Use Eq. 9.5.26 in the A&S Handbook.
0382               guess = boost::math::detail::bessel_zero::equation_as_9_5_26(v, airy_ai_root, pol);
0383             }
0384           }
0385 
0386           return guess;
0387         }
0388       } // namespace cyl_bessel_j_zero_detail
0389 
0390       namespace cyl_neumann_zero_detail
0391       {
0392         template<class T, class Policy>
0393         T equation_nist_10_21_40_b(const T& v, const Policy& pol)
0394         {
0395           const T v_pow_third(boost::math::cbrt(v, pol));
0396           const T v_pow_minus_two_thirds(T(1) / (v_pow_third * v_pow_third));
0397 
0398           return v * (((((                         - T(0.001)
0399                           * v_pow_minus_two_thirds - T(0.0060))
0400                           * v_pow_minus_two_thirds + T(0.01198))
0401                           * v_pow_minus_two_thirds + T(0.260351))
0402                           * v_pow_minus_two_thirds + T(0.9315768))
0403                           * v_pow_minus_two_thirds + T(1));
0404         }
0405 
0406         template<class T, class Policy>
0407         class function_object_yv
0408         {
0409         public:
0410           function_object_yv(const T& v,
0411                              const Policy& pol) : my_v(v),
0412                                                   my_pol(pol) { }
0413 
0414           function_object_yv(const function_object_yv&) = default;
0415 
0416           T operator()(const T& x) const
0417           {
0418             return boost::math::cyl_neumann(my_v, x, my_pol);
0419           }
0420 
0421         private:
0422           const T my_v;
0423           const Policy& my_pol;
0424           const function_object_yv& operator=(const function_object_yv&) = delete;
0425         };
0426 
0427         template<class T, class Policy>
0428         class function_object_yv_and_yv_prime
0429         {
0430         public:
0431           function_object_yv_and_yv_prime(const T& v,
0432                                           const Policy& pol) : my_v(v),
0433                                                                my_pol(pol) { }
0434 
0435           function_object_yv_and_yv_prime(const function_object_yv_and_yv_prime&) = default;
0436 
0437           boost::math::tuple<T, T> operator()(const T& x) const
0438           {
0439             const T half_epsilon(boost::math::tools::epsilon<T>() / 2U);
0440 
0441             const bool order_is_zero = ((my_v > -half_epsilon) && (my_v < +half_epsilon));
0442 
0443             // Obtain Yv(x) and Yv'(x).
0444             // Chris's original code called the Bessel function implementation layer direct, 
0445             // but that circumvented optimizations for integer-orders.  Call the documented
0446             // top level functions instead, and let them sort out which implementation to use.
0447             T y_v;
0448             T y_v_prime;
0449 
0450             if(order_is_zero)
0451             {
0452               y_v       =  boost::math::cyl_neumann(0, x, my_pol);
0453               y_v_prime = -boost::math::cyl_neumann(1, x, my_pol);
0454             }
0455             else
0456             {
0457                       y_v       = boost::math::cyl_neumann(  my_v,      x, my_pol);
0458               const T y_v_m1     (boost::math::cyl_neumann(T(my_v - 1), x, my_pol));
0459                       y_v_prime = y_v_m1 - ((my_v * y_v) / x);
0460             }
0461 
0462             // Return a tuple containing both Yv(x) and Yv'(x).
0463             return boost::math::make_tuple(y_v, y_v_prime);
0464           }
0465 
0466         private:
0467           const T my_v;
0468           const Policy& my_pol;
0469           const function_object_yv_and_yv_prime& operator=(const function_object_yv_and_yv_prime&) = delete;
0470         };
0471 
0472         template<class T> bool my_bisection_unreachable_tolerance(const T&, const T&) { return false; }
0473 
0474         template<class T, class Policy>
0475         T initial_guess(const T& v, const int m, const Policy& pol)
0476         {
0477           BOOST_MATH_STD_USING // ADL of std names, needed for floor.
0478 
0479           // Compute an estimate of the m'th root of cyl_neumann.
0480 
0481           T guess;
0482 
0483           // There is special handling for negative order.
0484           if(v < 0)
0485           {
0486             // Create the positive order and extract its positive floor and ceiling integer parts.
0487             const T vv(-v);
0488             const T vv_floor(floor(vv));
0489 
0490             // The to-be-found root is bracketed by the roots of the
0491             // Bessel function whose reflected, positive integer order
0492             // is less than, but nearest to vv.
0493 
0494             // The special case of negative, half-integer order uses
0495             // the relation between Yv and spherical Bessel functions
0496             // in order to obtain the bracket for the root.
0497             // In these special cases, cyl_neumann(-n/2, x) = sph_bessel_j(+n/2, x)
0498             // for v = -n/2.
0499 
0500             T root_hi;
0501             T root_lo;
0502 
0503             if(m == 1)
0504             {
0505               // The estimate of the first root for negative order is found using
0506               // an adaptive range-searching algorithm.
0507               // Take special precautions for the discontinuity at negative,
0508               // half-integer orders and use different brackets above and below these.
0509               if(T(vv - vv_floor) < 0.5F)
0510               {
0511                 root_hi = boost::math::detail::bessel_zero::cyl_neumann_zero_detail::initial_guess(vv_floor, m, pol);
0512               }
0513               else
0514               {
0515                 root_hi = boost::math::detail::bessel_zero::cyl_bessel_j_zero_detail::initial_guess(T(vv_floor + 0.5F), m, pol);
0516               }
0517 
0518               root_lo = T(root_hi - 0.1F);
0519 
0520               const bool hi_end_of_bracket_is_negative = (boost::math::cyl_neumann(v, root_hi, pol) < 0);
0521 
0522               while((root_lo > boost::math::tools::epsilon<T>()))
0523               {
0524                 const bool lo_end_of_bracket_is_negative = (boost::math::cyl_neumann(v, root_lo, pol) < 0);
0525 
0526                 if(hi_end_of_bracket_is_negative != lo_end_of_bracket_is_negative)
0527                 {
0528                   break;
0529                 }
0530 
0531                 root_hi = root_lo;
0532 
0533                 // Decrease the lower end of the bracket using an adaptive algorithm.
0534                 if(root_lo > 0.5F)
0535                 {
0536                   root_lo -= 0.5F;
0537                 }
0538                 else
0539                 {
0540                   root_lo *= 0.75F;
0541                 }
0542               }
0543             }
0544             else
0545             {
0546               if(T(vv - vv_floor) < 0.5F)
0547               {
0548                 root_lo  = boost::math::detail::bessel_zero::cyl_neumann_zero_detail::initial_guess(vv_floor, m - 1, pol);
0549                 root_hi = boost::math::detail::bessel_zero::cyl_neumann_zero_detail::initial_guess(vv_floor, m, pol);
0550                 root_lo += 0.01F;
0551                 root_hi += 0.01F;
0552               }
0553               else
0554               {
0555                 root_lo = boost::math::detail::bessel_zero::cyl_bessel_j_zero_detail::initial_guess(T(vv_floor + 0.5F), m - 1, pol);
0556                 root_hi = boost::math::detail::bessel_zero::cyl_bessel_j_zero_detail::initial_guess(T(vv_floor + 0.5F), m, pol);
0557                 root_lo += 0.01F;
0558                 root_hi += 0.01F;
0559               }
0560             }
0561 
0562             // Perform several steps of bisection iteration to refine the guess.
0563             std::uintmax_t number_of_iterations(12U);
0564 
0565             // Do the bisection iteration.
0566             const boost::math::tuple<T, T> guess_pair =
0567                boost::math::tools::bisect(
0568                   boost::math::detail::bessel_zero::cyl_neumann_zero_detail::function_object_yv<T, Policy>(v, pol),
0569                   root_lo,
0570                   root_hi,
0571                   boost::math::detail::bessel_zero::cyl_neumann_zero_detail::my_bisection_unreachable_tolerance<T>,
0572                   number_of_iterations);
0573 
0574             return (boost::math::get<0>(guess_pair) + boost::math::get<1>(guess_pair)) / 2U;
0575           }
0576 
0577           if(m == 1U)
0578           {
0579             // Get the initial estimate of the first root.
0580 
0581             if(v < 2.2F)
0582             {
0583               // For small v, use the results of empirical curve fitting.
0584               // Mathematica(R) session for the coefficients:
0585               //  Table[{n, BesselYZero[n, 1]}, {n, 0, 22/10, 1/10}]
0586               //  N[%, 20]
0587               //  Fit[%, {n^0, n^1, n^2, n^3, n^4, n^5, n^6}, n]
0588               guess = (((((    - T(0.0025095909235652)
0589                            * v + T(0.021291887049053))
0590                            * v - T(0.076487785486526))
0591                            * v + T(0.159110268115362))
0592                            * v - T(0.241681668765196))
0593                            * v + T(1.4437846310885244))
0594                            * v + T(0.89362115190200490);
0595             }
0596             else
0597             {
0598               // For larger v, use the second line of Eqs. 10.21.40 in the NIST Handbook.
0599               guess = boost::math::detail::bessel_zero::cyl_neumann_zero_detail::equation_nist_10_21_40_b(v, pol);
0600             }
0601           }
0602           else
0603           {
0604             if(v < 2.2F)
0605             {
0606               // Use Eq. 10.21.19 in the NIST Handbook.
0607               const T a(((v + T(m * 2U)) - T(1.5)) * boost::math::constants::half_pi<T>());
0608 
0609               guess = boost::math::detail::bessel_zero::equation_nist_10_21_19(v, a);
0610             }
0611             else
0612             {
0613               // Get an estimate of the m'th root of airy_bi.
0614               const T airy_bi_root(boost::math::detail::airy_zero::airy_bi_zero_detail::initial_guess<T>(m, pol));
0615 
0616               // Use Eq. 9.5.26 in the A&S Handbook.
0617               guess = boost::math::detail::bessel_zero::equation_as_9_5_26(v, airy_bi_root, pol);
0618             }
0619           }
0620 
0621           return guess;
0622         }
0623       } // namespace cyl_neumann_zero_detail
0624     } // namespace bessel_zero
0625   } } } // namespace boost::math::detail
0626 
0627 #endif // BOOST_MATH_BESSEL_JY_ZERO_2013_01_18_HPP_