Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-09-14 08:36:51

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