Back to home page

EIC code displayed by LXR

 
 

    


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

0001 //  Copyright (c) 2006 Xiaogang Zhang
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 #ifndef BOOST_MATH_BESSEL_JY_HPP
0007 #define BOOST_MATH_BESSEL_JY_HPP
0008 
0009 #ifdef _MSC_VER
0010 #pragma once
0011 #endif
0012 
0013 #include <boost/math/tools/config.hpp>
0014 #include <boost/math/special_functions/gamma.hpp>
0015 #include <boost/math/special_functions/sign.hpp>
0016 #include <boost/math/special_functions/hypot.hpp>
0017 #include <boost/math/special_functions/sin_pi.hpp>
0018 #include <boost/math/special_functions/cos_pi.hpp>
0019 #include <boost/math/special_functions/detail/bessel_jy_asym.hpp>
0020 #include <boost/math/special_functions/detail/bessel_jy_series.hpp>
0021 #include <boost/math/constants/constants.hpp>
0022 #include <boost/math/policies/error_handling.hpp>
0023 #include <complex>
0024 
0025 // Bessel functions of the first and second kind of fractional order
0026 
0027 namespace boost { namespace math {
0028 
0029    namespace detail {
0030 
0031       //
0032       // Simultaneous calculation of A&S 9.2.9 and 9.2.10
0033       // for use in A&S 9.2.5 and 9.2.6.
0034       // This series is quick to evaluate, but divergent unless
0035       // x is very large, in fact it's pretty hard to figure out
0036       // with any degree of precision when this series actually
0037       // *will* converge!!  Consequently, we may just have to
0038       // try it and see...
0039       //
0040       template <class T, class Policy>
0041       bool hankel_PQ(T v, T x, T* p, T* q, const Policy& )
0042       {
0043          BOOST_MATH_STD_USING
0044             T tolerance = 2 * policies::get_epsilon<T, Policy>();
0045          *p = 1;
0046          *q = 0;
0047          T k = 1;
0048          T z8 = 8 * x;
0049          T sq = 1;
0050          T mu = 4 * v * v;
0051          T term = 1;
0052          bool ok = true;
0053          do
0054          {
0055             term *= (mu - sq * sq) / (k * z8);
0056             *q += term;
0057             k += 1;
0058             sq += 2;
0059             T mult = (sq * sq - mu) / (k * z8);
0060             ok = fabs(mult) < 0.5f;
0061             term *= mult;
0062             *p += term;
0063             k += 1;
0064             sq += 2;
0065          }
0066          while((fabs(term) > tolerance * *p) && ok);
0067          return ok;
0068       }
0069 
0070       // Calculate Y(v, x) and Y(v+1, x) by Temme's method, see
0071       // Temme, Journal of Computational Physics, vol 21, 343 (1976)
0072       template <typename T, typename Policy>
0073       int temme_jy(T v, T x, T* Y, T* Y1, const Policy& pol)
0074       {
0075          T g, h, p, q, f, coef, sum, sum1, tolerance;
0076          T a, d, e, sigma;
0077          unsigned long k;
0078 
0079          BOOST_MATH_STD_USING
0080             using namespace boost::math::tools;
0081          using namespace boost::math::constants;
0082 
0083          BOOST_MATH_ASSERT(fabs(v) <= 0.5f);  // precondition for using this routine
0084 
0085          T gp = boost::math::tgamma1pm1(v, pol);
0086          T gm = boost::math::tgamma1pm1(-v, pol);
0087          T spv = boost::math::sin_pi(v, pol);
0088          T spv2 = boost::math::sin_pi(v/2, pol);
0089          T xp = pow(x/2, v);
0090 
0091          a = log(x / 2);
0092          sigma = -a * v;
0093          d = abs(sigma) < tools::epsilon<T>() ?
0094             T(1) : sinh(sigma) / sigma;
0095          e = abs(v) < tools::epsilon<T>() ? T(v*pi<T>()*pi<T>() / 2)
0096             : T(2 * spv2 * spv2 / v);
0097 
0098          T g1 = (v == 0) ? T(-euler<T>()) : T((gp - gm) / ((1 + gp) * (1 + gm) * 2 * v));
0099          T g2 = (2 + gp + gm) / ((1 + gp) * (1 + gm) * 2);
0100          T vspv = (fabs(v) < tools::epsilon<T>()) ? T(1/constants::pi<T>()) : T(v / spv);
0101          f = (g1 * cosh(sigma) - g2 * a * d) * 2 * vspv;
0102 
0103          p = vspv / (xp * (1 + gm));
0104          q = vspv * xp / (1 + gp);
0105 
0106          g = f + e * q;
0107          h = p;
0108          coef = 1;
0109          sum = coef * g;
0110          sum1 = coef * h;
0111 
0112          T v2 = v * v;
0113          T coef_mult = -x * x / 4;
0114 
0115          // series summation
0116          tolerance = policies::get_epsilon<T, Policy>();
0117          for (k = 1; k < policies::get_max_series_iterations<Policy>(); k++)
0118          {
0119             f = (k * f + p + q) / (k*k - v2);
0120             p /= k - v;
0121             q /= k + v;
0122             g = f + e * q;
0123             h = p - k * g;
0124             coef *= coef_mult / k;
0125             sum += coef * g;
0126             sum1 += coef * h;
0127             if (abs(coef * g) < abs(sum) * tolerance)
0128             {
0129                break;
0130             }
0131          }
0132          policies::check_series_iterations<T>("boost::math::bessel_jy<%1%>(%1%,%1%) in temme_jy", k, pol);
0133          *Y = -sum;
0134          *Y1 = -2 * sum1 / x;
0135 
0136          return 0;
0137       }
0138 
0139       // Evaluate continued fraction fv = J_(v+1) / J_v, see
0140       // Abramowitz and Stegun, Handbook of Mathematical Functions, 1972, 9.1.73
0141       template <typename T, typename Policy>
0142       int CF1_jy(T v, T x, T* fv, int* sign, const Policy& pol)
0143       {
0144          T C, D, f, a, b, delta, tiny, tolerance;
0145          unsigned long k;
0146          int s = 1;
0147 
0148          BOOST_MATH_STD_USING
0149 
0150             // |x| <= |v|, CF1_jy converges rapidly
0151             // |x| > |v|, CF1_jy needs O(|x|) iterations to converge
0152 
0153             // modified Lentz's method, see
0154             // Lentz, Applied Optics, vol 15, 668 (1976)
0155             tolerance = 2 * policies::get_epsilon<T, Policy>();
0156          tiny = sqrt(tools::min_value<T>());
0157          C = f = tiny;                           // b0 = 0, replace with tiny
0158          D = 0;
0159          for (k = 1; k < policies::get_max_series_iterations<Policy>() * 100; k++)
0160          {
0161             a = -1;
0162             b = 2 * (v + k) / x;
0163             C = b + a / C;
0164             D = b + a * D;
0165             if (C == 0) { C = tiny; }
0166             if (D == 0) { D = tiny; }
0167             D = 1 / D;
0168             delta = C * D;
0169             f *= delta;
0170             if (D < 0) { s = -s; }
0171             if (abs(delta - 1) < tolerance)
0172             { break; }
0173          }
0174          policies::check_series_iterations<T>("boost::math::bessel_jy<%1%>(%1%,%1%) in CF1_jy", k / 100, pol);
0175          *fv = -f;
0176          *sign = s;                              // sign of denominator
0177 
0178          return 0;
0179       }
0180       //
0181       // This algorithm was originally written by Xiaogang Zhang
0182       // using std::complex to perform the complex arithmetic.
0183       // However, that turns out to 10x or more slower than using
0184       // all real-valued arithmetic, so it's been rewritten using
0185       // real values only.
0186       //
0187       template <typename T, typename Policy>
0188       int CF2_jy(T v, T x, T* p, T* q, const Policy& pol)
0189       {
0190          BOOST_MATH_STD_USING
0191 
0192             T Cr, Ci, Dr, Di, fr, fi, a, br, bi, delta_r, delta_i, temp;
0193          T tiny;
0194          unsigned long k;
0195 
0196          // |x| >= |v|, CF2_jy converges rapidly
0197          // |x| -> 0, CF2_jy fails to converge
0198          BOOST_MATH_ASSERT(fabs(x) > 1);
0199 
0200          // modified Lentz's method, complex numbers involved, see
0201          // Lentz, Applied Optics, vol 15, 668 (1976)
0202          T tolerance = 2 * policies::get_epsilon<T, Policy>();
0203          tiny = sqrt(tools::min_value<T>());
0204          Cr = fr = -0.5f / x;
0205          Ci = fi = 1;
0206          //Dr = Di = 0;
0207          T v2 = v * v;
0208          a = (0.25f - v2) / x; // Note complex this one time only!
0209          br = 2 * x;
0210          bi = 2;
0211          temp = Cr * Cr + 1;
0212          Ci = bi + a * Cr / temp;
0213          Cr = br + a / temp;
0214          Dr = br;
0215          Di = bi;
0216          if (fabs(Cr) + fabs(Ci) < tiny) { Cr = tiny; }
0217          if (fabs(Dr) + fabs(Di) < tiny) { Dr = tiny; }
0218          temp = Dr * Dr + Di * Di;
0219          Dr = Dr / temp;
0220          Di = -Di / temp;
0221          delta_r = Cr * Dr - Ci * Di;
0222          delta_i = Ci * Dr + Cr * Di;
0223          temp = fr;
0224          fr = temp * delta_r - fi * delta_i;
0225          fi = temp * delta_i + fi * delta_r;
0226          for (k = 2; k < policies::get_max_series_iterations<Policy>(); k++)
0227          {
0228             a = k - 0.5f;
0229             a *= a;
0230             a -= v2;
0231             bi += 2;
0232             temp = Cr * Cr + Ci * Ci;
0233             Cr = br + a * Cr / temp;
0234             Ci = bi - a * Ci / temp;
0235             Dr = br + a * Dr;
0236             Di = bi + a * Di;
0237             if (fabs(Cr) + fabs(Ci) < tiny) { Cr = tiny; }
0238             if (fabs(Dr) + fabs(Di) < tiny) { Dr = tiny; }
0239             temp = Dr * Dr + Di * Di;
0240             Dr = Dr / temp;
0241             Di = -Di / temp;
0242             delta_r = Cr * Dr - Ci * Di;
0243             delta_i = Ci * Dr + Cr * Di;
0244             temp = fr;
0245             fr = temp * delta_r - fi * delta_i;
0246             fi = temp * delta_i + fi * delta_r;
0247             if (fabs(delta_r - 1) + fabs(delta_i) < tolerance)
0248                break;
0249          }
0250          policies::check_series_iterations<T>("boost::math::bessel_jy<%1%>(%1%,%1%) in CF2_jy", k, pol);
0251          *p = fr;
0252          *q = fi;
0253 
0254          return 0;
0255       }
0256 
0257       static const int need_j = 1;
0258       static const int need_y = 2;
0259 
0260       // Compute J(v, x) and Y(v, x) simultaneously by Steed's method, see
0261       // Barnett et al, Computer Physics Communications, vol 8, 377 (1974)
0262       template <typename T, typename Policy>
0263       int bessel_jy(T v, T x, T* J, T* Y, int kind, const Policy& pol)
0264       {
0265          BOOST_MATH_ASSERT(x >= 0);
0266 
0267          T u, Jv, Ju, Yv, Yv1, Yu, Yu1(0), fv, fu;
0268          T W, p, q, gamma, current, prev, next;
0269          bool reflect = false;
0270          unsigned n, k;
0271          int s;
0272          int org_kind = kind;
0273          T cp = 0;
0274          T sp = 0;
0275 
0276          static const char* function = "boost::math::bessel_jy<%1%>(%1%,%1%)";
0277 
0278          BOOST_MATH_STD_USING
0279             using namespace boost::math::tools;
0280          using namespace boost::math::constants;
0281 
0282          if (v < 0)
0283          {
0284             reflect = true;
0285             v = -v;                             // v is non-negative from here
0286          }
0287          if (v > static_cast<T>((std::numeric_limits<int>::max)()))
0288          {
0289             *J = *Y = policies::raise_evaluation_error<T>(function, "Order of Bessel function is too large to evaluate: got %1%", v, pol);
0290             return 1;
0291          }
0292          n = iround(v, pol);
0293          u = v - n;                              // -1/2 <= u < 1/2
0294 
0295          if(reflect)
0296          {
0297             T z = (u + n % 2);
0298             cp = boost::math::cos_pi(z, pol);
0299             sp = boost::math::sin_pi(z, pol);
0300             if(u != 0)
0301                kind = need_j|need_y;               // need both for reflection formula
0302          }
0303 
0304          if(x == 0)
0305          {
0306             if(v == 0)
0307                *J = 1;
0308             else if((u == 0) || !reflect)
0309                *J = 0;
0310             else if(kind & need_j)
0311                *J = policies::raise_domain_error<T>(function, "Value of Bessel J_v(x) is complex-infinity at %1%", x, pol); // complex infinity
0312             else
0313                *J = std::numeric_limits<T>::quiet_NaN();  // any value will do, not using J.
0314 
0315             if((kind & need_y) == 0)
0316                *Y = std::numeric_limits<T>::quiet_NaN();  // any value will do, not using Y.
0317             else if(v == 0)
0318                *Y = -policies::raise_overflow_error<T>(function, nullptr, pol);
0319             else
0320                *Y = policies::raise_domain_error<T>(function, "Value of Bessel Y_v(x) is complex-infinity at %1%", x, pol); // complex infinity
0321             return 1;
0322          }
0323 
0324          // x is positive until reflection
0325          W = T(2) / (x * pi<T>());               // Wronskian
0326          T Yv_scale = 1;
0327          if(((kind & need_y) == 0) && ((x < 1) || (v > x * x / 4) || (x < 5)))
0328          {
0329             //
0330             // This series will actually converge rapidly for all small
0331             // x - say up to x < 20 - but the first few terms are large
0332             // and divergent which leads to large errors :-(
0333             //
0334             Jv = bessel_j_small_z_series(v, x, pol);
0335             Yv = std::numeric_limits<T>::quiet_NaN();
0336          }
0337          else if((x < 1) && (u != 0) && (log(policies::get_epsilon<T, Policy>() / 2) > v * log((x/2) * (x/2) / v)))
0338          {
0339             // Evaluate using series representations.
0340             // This is particularly important for x << v as in this
0341             // area temme_jy may be slow to converge, if it converges at all.
0342             // Requires x is not an integer.
0343             if(kind&need_j)
0344                Jv = bessel_j_small_z_series(v, x, pol);
0345             else
0346                Jv = std::numeric_limits<T>::quiet_NaN();
0347             if((org_kind&need_y && (!reflect || (cp != 0)))
0348                || (org_kind & need_j && (reflect && (sp != 0))))
0349             {
0350                // Only calculate if we need it, and if the reflection formula will actually use it:
0351                Yv = bessel_y_small_z_series(v, x, &Yv_scale, pol);
0352             }
0353             else
0354                Yv = std::numeric_limits<T>::quiet_NaN();
0355          }
0356          else if((u == 0) && (x < policies::get_epsilon<T, Policy>()))
0357          {
0358             // Truncated series evaluation for small x and v an integer,
0359             // much quicker in this area than temme_jy below.
0360             if(kind&need_j)
0361                Jv = bessel_j_small_z_series(v, x, pol);
0362             else
0363                Jv = std::numeric_limits<T>::quiet_NaN();
0364             if((org_kind&need_y && (!reflect || (cp != 0)))
0365                || (org_kind & need_j && (reflect && (sp != 0))))
0366             {
0367                // Only calculate if we need it, and if the reflection formula will actually use it:
0368                Yv = bessel_yn_small_z(n, x, &Yv_scale, pol);
0369             }
0370             else
0371                Yv = std::numeric_limits<T>::quiet_NaN();
0372          }
0373          else if(asymptotic_bessel_large_x_limit(v, x))
0374          {
0375             if(kind&need_y)
0376             {
0377                Yv = asymptotic_bessel_y_large_x_2(v, x, pol);
0378             }
0379             else
0380                Yv = std::numeric_limits<T>::quiet_NaN(); // any value will do, we're not using it.
0381             if(kind&need_j)
0382             {
0383                Jv = asymptotic_bessel_j_large_x_2(v, x, pol);
0384             }
0385             else
0386                Jv = std::numeric_limits<T>::quiet_NaN(); // any value will do, we're not using it.
0387          }
0388          else if((x > 8) && hankel_PQ(v, x, &p, &q, pol))
0389          {
0390             //
0391             // Hankel approximation: note that this method works best when x
0392             // is large, but in that case we end up calculating sines and cosines
0393             // of large values, with horrendous resulting accuracy.  It is fast though
0394             // when it works....
0395             //
0396             // Normally we calculate sin/cos(chi) where:
0397             //
0398             // chi = x - fmod(T(v / 2 + 0.25f), T(2)) * boost::math::constants::pi<T>();
0399             //
0400             // But this introduces large errors, so use sin/cos addition formulae to
0401             // improve accuracy:
0402             //
0403             T mod_v = fmod(T(v / 2 + 0.25f), T(2));
0404             T sx = sin(x);
0405             T cx = cos(x);
0406             T sv = boost::math::sin_pi(mod_v, pol);
0407             T cv = boost::math::cos_pi(mod_v, pol);
0408 
0409             T sc = sx * cv - sv * cx; // == sin(chi);
0410             T cc = cx * cv + sx * sv; // == cos(chi);
0411             T chi = boost::math::constants::root_two<T>() / (boost::math::constants::root_pi<T>() * sqrt(x)); //sqrt(2 / (boost::math::constants::pi<T>() * x));
0412             Yv = chi * (p * sc + q * cc);
0413             Jv = chi * (p * cc - q * sc);
0414          }
0415          else if (x <= 2)                           // x in (0, 2]
0416          {
0417             if(temme_jy(u, x, &Yu, &Yu1, pol))             // Temme series
0418             {
0419                // domain error:
0420                *J = *Y = Yu;
0421                return 1;
0422             }
0423             prev = Yu;
0424             current = Yu1;
0425             T scale = 1;
0426             policies::check_series_iterations<T>(function, n, pol);
0427             for (k = 1; k <= n; k++)            // forward recurrence for Y
0428             {
0429                T fact = 2 * (u + k) / x;
0430                if((tools::max_value<T>() - fabs(prev)) / fact < fabs(current))
0431                {
0432                   scale /= current;
0433                   prev /= current;
0434                   current = 1;
0435                }
0436                next = fact * current - prev;
0437                prev = current;
0438                current = next;
0439             }
0440             Yv = prev;
0441             Yv1 = current;
0442             if(kind&need_j)
0443             {
0444                CF1_jy(v, x, &fv, &s, pol);                 // continued fraction CF1_jy
0445                Jv = scale * W / (Yv * fv - Yv1);           // Wronskian relation
0446             }
0447             else
0448                Jv = std::numeric_limits<T>::quiet_NaN(); // any value will do, we're not using it.
0449             Yv_scale = scale;
0450          }
0451          else                                    // x in (2, \infty)
0452          {
0453             // Get Y(u, x):
0454 
0455             T ratio;
0456             CF1_jy(v, x, &fv, &s, pol);
0457             // tiny initial value to prevent overflow
0458             T init = sqrt(tools::min_value<T>());
0459             BOOST_MATH_INSTRUMENT_VARIABLE(init);
0460             prev = fv * s * init;
0461             current = s * init;
0462             if(v < max_factorial<T>::value)
0463             {
0464                policies::check_series_iterations<T>(function, n, pol);
0465                for (k = n; k > 0; k--)             // backward recurrence for J
0466                {
0467                   next = 2 * (u + k) * current / x - prev;
0468                   //
0469                   // We can't allow next to completely cancel out or the subsequent logic breaks.
0470                   // Pretend that one bit did not cancel:
0471                   if (next == 0)
0472                   {
0473                      next = prev * tools::epsilon<T>() / 2;
0474                   }
0475                   prev = current;
0476                   current = next;
0477                }
0478                ratio = (s * init) / current;     // scaling ratio
0479                // can also call CF1_jy() to get fu, not much difference in precision
0480                fu = prev / current;
0481             }
0482             else
0483             {
0484                //
0485                // When v is large we may get overflow in this calculation
0486                // leading to NaN's and other nasty surprises:
0487                //
0488                policies::check_series_iterations<T>(function, n, pol);
0489                bool over = false;
0490                for (k = n; k > 0; k--)             // backward recurrence for J
0491                {
0492                   T t = 2 * (u + k) / x;
0493                   if((t > 1) && (tools::max_value<T>() / t < current))
0494                   {
0495                      over = true;
0496                      break;
0497                   }
0498                   next = t * current - prev;
0499                   prev = current;
0500                   current = next;
0501                }
0502                if(!over)
0503                {
0504                   ratio = (s * init) / current;     // scaling ratio
0505                   // can also call CF1_jy() to get fu, not much difference in precision
0506                   fu = prev / current;
0507                }
0508                else
0509                {
0510                   ratio = 0;
0511                   fu = 1;
0512                }
0513             }
0514             CF2_jy(u, x, &p, &q, pol);                  // continued fraction CF2_jy
0515             T t = u / x - fu;                   // t = J'/J
0516             gamma = (p - t) / q;
0517             //
0518             // We can't allow gamma to cancel out to zero completely as it messes up
0519             // the subsequent logic.  So pretend that one bit didn't cancel out
0520             // and set to a suitably small value.  The only test case we've been able to
0521             // find for this, is when v = 8.5 and x = 4*PI.
0522             //
0523             if(gamma == 0)
0524             {
0525                gamma = u * tools::epsilon<T>() / x;
0526             }
0527             BOOST_MATH_INSTRUMENT_VARIABLE(current);
0528             BOOST_MATH_INSTRUMENT_VARIABLE(W);
0529             BOOST_MATH_INSTRUMENT_VARIABLE(q);
0530             BOOST_MATH_INSTRUMENT_VARIABLE(gamma);
0531             BOOST_MATH_INSTRUMENT_VARIABLE(p);
0532             BOOST_MATH_INSTRUMENT_VARIABLE(t);
0533             Ju = sign(current) * sqrt(W / (q + gamma * (p - t)));
0534             BOOST_MATH_INSTRUMENT_VARIABLE(Ju);
0535 
0536             Jv = Ju * ratio;                    // normalization
0537 
0538             Yu = gamma * Ju;
0539             Yu1 = Yu * (u/x - p - q/gamma);
0540 
0541             if(kind&need_y)
0542             {
0543                // compute Y:
0544                prev = Yu;
0545                current = Yu1;
0546                policies::check_series_iterations<T>(function, n, pol);
0547                for (k = 1; k <= n; k++)            // forward recurrence for Y
0548                {
0549                   T fact = 2 * (u + k) / x;
0550                   if((tools::max_value<T>() - fabs(prev)) / fact < fabs(current))
0551                   {
0552                      prev /= current;
0553                      Yv_scale /= current;
0554                      current = 1;
0555                   }
0556                   next = fact * current - prev;
0557                   prev = current;
0558                   current = next;
0559                }
0560                Yv = prev;
0561             }
0562             else
0563                Yv = std::numeric_limits<T>::quiet_NaN(); // any value will do, we're not using it.
0564          }
0565 
0566          if (reflect)
0567          {
0568             if((sp != 0) && (tools::max_value<T>() * fabs(Yv_scale) < fabs(sp * Yv)))
0569                *J = org_kind & need_j ? T(-sign(sp) * sign(Yv) * (Yv_scale != 0 ? sign(Yv_scale) : 1) * policies::raise_overflow_error<T>(function, nullptr, pol)) : T(0);
0570             else
0571                *J = cp * Jv - (sp == 0 ? T(0) : T((sp * Yv) / Yv_scale));     // reflection formula
0572             if((cp != 0) && (tools::max_value<T>() * fabs(Yv_scale) < fabs(cp * Yv)))
0573                *Y = org_kind & need_y ? T(-sign(cp) * sign(Yv) * (Yv_scale != 0 ? sign(Yv_scale) : 1) * policies::raise_overflow_error<T>(function, nullptr, pol)) : T(0);
0574             else
0575                *Y = (sp != 0 ? sp * Jv : T(0)) + (cp == 0 ? T(0) : T((cp * Yv) / Yv_scale));
0576          }
0577          else
0578          {
0579             *J = Jv;
0580             if(tools::max_value<T>() * fabs(Yv_scale) < fabs(Yv))
0581                *Y = org_kind & need_y ? T(sign(Yv) * sign(Yv_scale) * policies::raise_overflow_error<T>(function, nullptr, pol)) : T(0);
0582             else
0583                *Y = Yv / Yv_scale;
0584          }
0585 
0586          return 0;
0587       }
0588 
0589    } // namespace detail
0590 
0591 }} // namespaces
0592 
0593 #endif // BOOST_MATH_BESSEL_JY_HPP