Back to home page

EIC code displayed by LXR

 
 

    


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

0001 //  Copyright John Maddock 2006.
0002 //  Copyright Paul A. Bristow 2007
0003 //  Use, modification and distribution are subject to the
0004 //  Boost Software License, Version 1.0. (See accompanying file
0005 //  LICENSE_1_0.txt or copy at http://www.boost.org/LICENSE_1_0.txt)
0006 
0007 #ifndef BOOST_MATH_SPECIAL_FUNCTIONS_IBETA_INVERSE_HPP
0008 #define BOOST_MATH_SPECIAL_FUNCTIONS_IBETA_INVERSE_HPP
0009 
0010 #ifdef _MSC_VER
0011 #pragma once
0012 #endif
0013 
0014 #include <boost/math/special_functions/beta.hpp>
0015 #include <boost/math/special_functions/erf.hpp>
0016 #include <boost/math/tools/roots.hpp>
0017 #include <boost/math/special_functions/detail/t_distribution_inv.hpp>
0018 
0019 namespace boost{ namespace math{ namespace detail{
0020 
0021 //
0022 // Helper object used by root finding
0023 // code to convert eta to x.
0024 //
0025 template <class T>
0026 struct temme_root_finder
0027 {
0028    temme_root_finder(const T t_, const T a_) : t(t_), a(a_) {
0029       const T x_extrema = 1 / (1 + a);
0030       BOOST_MATH_ASSERT(0 < x_extrema && x_extrema < 1);
0031    }
0032 
0033    boost::math::tuple<T, T> operator()(T x)
0034    {
0035       BOOST_MATH_STD_USING // ADL of std names
0036 
0037       T y = 1 - x;
0038       T f = log(x) + a * log(y) + t;
0039       T f1 = (1 / x) - (a / (y));
0040       return boost::math::make_tuple(f, f1);
0041    }
0042 private:
0043    T t, a;
0044 };
0045 //
0046 // See:
0047 // "Asymptotic Inversion of the Incomplete Beta Function"
0048 // N.M. Temme
0049 // Journal of Computation and Applied Mathematics 41 (1992) 145-157.
0050 // Section 2.
0051 //
0052 template <class T, class Policy>
0053 T temme_method_1_ibeta_inverse(T a, T b, T z, const Policy& pol)
0054 {
0055    BOOST_MATH_STD_USING // ADL of std names
0056 
0057    const T r2 = sqrt(T(2));
0058    //
0059    // get the first approximation for eta from the inverse
0060    // error function (Eq: 2.9 and 2.10).
0061    //
0062    T eta0 = boost::math::erfc_inv(2 * z, pol);
0063    eta0 /= -sqrt(a / 2);
0064 
0065    T terms[4] = { eta0 };
0066    T workspace[7];
0067    //
0068    // calculate powers:
0069    //
0070    T B = b - a;
0071    T B_2 = B * B;
0072    T B_3 = B_2 * B;
0073    //
0074    // Calculate correction terms:
0075    //
0076 
0077    // See eq following 2.15:
0078    workspace[0] = -B * r2 / 2;
0079    workspace[1] = (1 - 2 * B) / 8;
0080    workspace[2] = -(B * r2 / 48);
0081    workspace[3] = T(-1) / 192;
0082    workspace[4] = -B * r2 / 3840;
0083    terms[1] = tools::evaluate_polynomial(workspace, eta0, 5);
0084    // Eq Following 2.17:
0085    workspace[0] = B * r2 * (3 * B - 2) / 12;
0086    workspace[1] = (20 * B_2 - 12 * B + 1) / 128;
0087    workspace[2] = B * r2 * (20 * B - 1) / 960;
0088    workspace[3] = (16 * B_2 + 30 * B - 15) / 4608;
0089    workspace[4] = B * r2 * (21 * B + 32) / 53760;
0090    workspace[5] = (-32 * B_2 + 63) / 368640;
0091    workspace[6] = -B * r2 * (120 * B + 17) / 25804480;
0092    terms[2] = tools::evaluate_polynomial(workspace, eta0, 7);
0093    // Eq Following 2.17:
0094    workspace[0] = B * r2 * (-75 * B_2 + 80 * B - 16) / 480;
0095    workspace[1] = (-1080 * B_3 + 868 * B_2 - 90 * B - 45) / 9216;
0096    workspace[2] = B * r2 * (-1190 * B_2 + 84 * B + 373) / 53760;
0097    workspace[3] = (-2240 * B_3 - 2508 * B_2 + 2100 * B - 165) / 368640;
0098    terms[3] = tools::evaluate_polynomial(workspace, eta0, 4);
0099    //
0100    // Bring them together to get a final estimate for eta:
0101    //
0102    T eta = tools::evaluate_polynomial(terms, T(1/a), 4);
0103    //
0104    // now we need to convert eta to x, by solving the appropriate
0105    // quadratic equation:
0106    //
0107    T eta_2 = eta * eta;
0108    T c = -exp(-eta_2 / 2);
0109    T x;
0110    if(eta_2 == 0)
0111       x = static_cast<T>(0.5f);
0112    else
0113       x = (1 + eta * sqrt((1 + c) / eta_2)) / 2;
0114    //
0115    // These are post-conditions of the method, but the addition above
0116    // may result in us being out by 1ulp either side of the boundary,
0117    // so just check that we're in bounds and adjust as needed.
0118    // See https://github.com/boostorg/math/issues/961
0119    //
0120    if (x < 0)
0121       x = 0;
0122    else if (x > 1)
0123       x = 1;
0124    
0125    BOOST_MATH_ASSERT(eta * (x - 0.5) >= 0);
0126 #ifdef BOOST_INSTRUMENT
0127    std::cout << "Estimating x with Temme method 1: " << x << std::endl;
0128 #endif
0129    return x;
0130 }
0131 //
0132 // See:
0133 // "Asymptotic Inversion of the Incomplete Beta Function"
0134 // N.M. Temme
0135 // Journal of Computation and Applied Mathematics 41 (1992) 145-157.
0136 // Section 3.
0137 //
0138 template <class T, class Policy>
0139 T temme_method_2_ibeta_inverse(T /*a*/, T /*b*/, T z, T r, T theta, const Policy& pol)
0140 {
0141    BOOST_MATH_STD_USING // ADL of std names
0142 
0143    //
0144    // Get first estimate for eta, see Eq 3.9 and 3.10,
0145    // but note there is a typo in Eq 3.10:
0146    //
0147    T eta0 = boost::math::erfc_inv(2 * z, pol);
0148    eta0 /= -sqrt(r / 2);
0149 
0150    T s = sin(theta);
0151    T c = cos(theta);
0152    //
0153    // Now we need to perturb eta0 to get eta, which we do by
0154    // evaluating the polynomial in 1/r at the bottom of page 151,
0155    // to do this we first need the error terms e1, e2 e3
0156    // which we'll fill into the array "terms".  Since these
0157    // terms are themselves polynomials, we'll need another
0158    // array "workspace" to calculate those...
0159    //
0160    T terms[4] = { eta0 };
0161    T workspace[6];
0162    //
0163    // some powers of sin(theta)cos(theta) that we'll need later:
0164    //
0165    T sc = s * c;
0166    T sc_2 = sc * sc;
0167    T sc_3 = sc_2 * sc;
0168    T sc_4 = sc_2 * sc_2;
0169    T sc_5 = sc_2 * sc_3;
0170    T sc_6 = sc_3 * sc_3;
0171    T sc_7 = sc_4 * sc_3;
0172    //
0173    // Calculate e1 and put it in terms[1], see the middle of page 151:
0174    //
0175    workspace[0] = (2 * s * s - 1) / (3 * s * c);
0176    static const BOOST_MATH_INT_TABLE_TYPE(T, int) co1[] = { -1, -5, 5 };
0177    workspace[1] = -tools::evaluate_even_polynomial(co1, s, 3) / (36 * sc_2);
0178    static const BOOST_MATH_INT_TABLE_TYPE(T, int) co2[] = { 1, 21, -69, 46 };
0179    workspace[2] = tools::evaluate_even_polynomial(co2, s, 4) / (1620 * sc_3);
0180    static const BOOST_MATH_INT_TABLE_TYPE(T, int) co3[] = { 7, -2, 33, -62, 31 };
0181    workspace[3] = -tools::evaluate_even_polynomial(co3, s, 5) / (6480 * sc_4);
0182    static const BOOST_MATH_INT_TABLE_TYPE(T, int) co4[] = { 25, -52, -17, 88, -115, 46 };
0183    workspace[4] = tools::evaluate_even_polynomial(co4, s, 6) / (90720 * sc_5);
0184    terms[1] = tools::evaluate_polynomial(workspace, eta0, 5);
0185    //
0186    // Now evaluate e2 and put it in terms[2]:
0187    //
0188    static const BOOST_MATH_INT_TABLE_TYPE(T, int) co5[] = { 7, 12, -78, 52 };
0189    workspace[0] = -tools::evaluate_even_polynomial(co5, s, 4) / (405 * sc_3);
0190    static const BOOST_MATH_INT_TABLE_TYPE(T, int) co6[] = { -7, 2, 183, -370, 185 };
0191    workspace[1] = tools::evaluate_even_polynomial(co6, s, 5) / (2592 * sc_4);
0192    static const BOOST_MATH_INT_TABLE_TYPE(T, int) co7[] = { -533, 776, -1835, 10240, -13525, 5410 };
0193    workspace[2] = -tools::evaluate_even_polynomial(co7, s, 6) / (204120 * sc_5);
0194    static const BOOST_MATH_INT_TABLE_TYPE(T, int) co8[] = { -1579, 3747, -3372, -15821, 45588, -45213, 15071 };
0195    workspace[3] = -tools::evaluate_even_polynomial(co8, s, 7) / (2099520 * sc_6);
0196    terms[2] = tools::evaluate_polynomial(workspace, eta0, 4);
0197    //
0198    // And e3, and put it in terms[3]:
0199    //
0200    static const BOOST_MATH_INT_TABLE_TYPE(T, int) co9[] = {449, -1259, -769, 6686, -9260, 3704 };
0201    workspace[0] = tools::evaluate_even_polynomial(co9, s, 6) / (102060 * sc_5);
0202    static const BOOST_MATH_INT_TABLE_TYPE(T, int) co10[] = { 63149, -151557, 140052, -727469, 2239932, -2251437, 750479 };
0203    workspace[1] = -tools::evaluate_even_polynomial(co10, s, 7) / (20995200 * sc_6);
0204    static const BOOST_MATH_INT_TABLE_TYPE(T, int) co11[] = { 29233, -78755, 105222, 146879, -1602610, 3195183, -2554139, 729754 };
0205    workspace[2] = tools::evaluate_even_polynomial(co11, s, 8) / (36741600 * sc_7);
0206    terms[3] = tools::evaluate_polynomial(workspace, eta0, 3);
0207    //
0208    // Bring the correction terms together to evaluate eta,
0209    // this is the last equation on page 151:
0210    //
0211    T eta = tools::evaluate_polynomial(terms, T(1/r), 4);
0212    //
0213    // Now that we have eta we need to back solve for x,
0214    // we seek the value of x that gives eta in Eq 3.2.
0215    // The two methods used are described in section 5.
0216    //
0217    // Begin by defining a few variables we'll need later:
0218    //
0219    T x;
0220    T s_2 = s * s;
0221    T c_2 = c * c;
0222    T alpha = c / s;
0223    alpha *= alpha;
0224    T lu = (-(eta * eta) / (2 * s_2) + log(s_2) + c_2 * log(c_2) / s_2);
0225    //
0226    // Temme doesn't specify what value to switch on here,
0227    // but this seems to work pretty well:
0228    //
0229    if(fabs(eta) < 0.7)
0230    {
0231       //
0232       // Small eta use the expansion Temme gives in the second equation
0233       // of section 5, it's a polynomial in eta:
0234       //
0235       workspace[0] = s * s;
0236       workspace[1] = s * c;
0237       workspace[2] = (1 - 2 * workspace[0]) / 3;
0238       static const BOOST_MATH_INT_TABLE_TYPE(T, int) co12[] = { 1, -13, 13 };
0239       workspace[3] = tools::evaluate_polynomial(co12, workspace[0], 3) / (36 * s * c);
0240       static const BOOST_MATH_INT_TABLE_TYPE(T, int) co13[] = { 1, 21, -69, 46 };
0241       workspace[4] = tools::evaluate_polynomial(co13, workspace[0], 4) / (270 * workspace[0] * c * c);
0242       x = tools::evaluate_polynomial(workspace, eta, 5);
0243 #ifdef BOOST_INSTRUMENT
0244       std::cout << "Estimating x with Temme method 2 (small eta): " << x << std::endl;
0245 #endif
0246    }
0247    else
0248    {
0249       //
0250       // If eta is large we need to solve Eq 3.2 more directly,
0251       // begin by getting an initial approximation for x from
0252       // the last equation on page 155, this is a polynomial in u:
0253       //
0254       T u = exp(lu);
0255       workspace[0] = u;
0256       workspace[1] = alpha;
0257       workspace[2] = 0;
0258       workspace[3] = 3 * alpha * (3 * alpha + 1) / 6;
0259       workspace[4] = 4 * alpha * (4 * alpha + 1) * (4 * alpha + 2) / 24;
0260       workspace[5] = 5 * alpha * (5 * alpha + 1) * (5 * alpha + 2) * (5 * alpha + 3) / 120;
0261       x = tools::evaluate_polynomial(workspace, u, 6);
0262       //
0263       // At this point we may or may not have the right answer, Eq-3.2 has
0264       // two solutions for x for any given eta, however the mapping in 3.2
0265       // is 1:1 with the sign of eta and x-sin^2(theta) being the same.
0266       // So we can check if we have the right root of 3.2, and if not
0267       // switch x for 1-x.  This transformation is motivated by the fact
0268       // that the distribution is *almost* symmetric so 1-x will be in the right
0269       // ball park for the solution:
0270       //
0271       if((x - s_2) * eta < 0)
0272          x = 1 - x;
0273 #ifdef BOOST_INSTRUMENT
0274       std::cout << "Estimating x with Temme method 2 (large eta): " << x << std::endl;
0275 #endif
0276    }
0277    //
0278    // The final step is a few Newton-Raphson iterations to
0279    // clean up our approximation for x, this is pretty cheap
0280    // in general, and very cheap compared to an incomplete beta
0281    // evaluation.  The limits set on x come from the observation
0282    // that the sign of eta and x-sin^2(theta) are the same.
0283    //
0284    T lower, upper;
0285    if(eta < 0)
0286    {
0287       lower = 0;
0288       upper = s_2;
0289    }
0290    else
0291    {
0292       lower = s_2;
0293       upper = 1;
0294    }
0295    //
0296    // If our initial approximation is out of bounds then bisect:
0297    //
0298    if((x < lower) || (x > upper))
0299       x = (lower+upper) / 2;
0300    //
0301    // And iterate:
0302    //
0303    x = tools::newton_raphson_iterate(
0304       temme_root_finder<T>(-lu, alpha), x, lower, upper, policies::digits<T, Policy>() / 2);
0305 
0306    return x;
0307 }
0308 //
0309 // See:
0310 // "Asymptotic Inversion of the Incomplete Beta Function"
0311 // N.M. Temme
0312 // Journal of Computation and Applied Mathematics 41 (1992) 145-157.
0313 // Section 4.
0314 //
0315 template <class T, class Policy>
0316 T temme_method_3_ibeta_inverse(T a, T b, T p, T q, const Policy& pol)
0317 {
0318    BOOST_MATH_STD_USING // ADL of std names
0319 
0320    //
0321    // Begin by getting an initial approximation for the quantity
0322    // eta from the dominant part of the incomplete beta:
0323    //
0324    T eta0;
0325    if(p < q)
0326       eta0 = boost::math::gamma_q_inv(b, p, pol);
0327    else
0328       eta0 = boost::math::gamma_p_inv(b, q, pol);
0329    eta0 /= a;
0330    //
0331    // Define the variables and powers we'll need later on:
0332    //
0333    T mu = b / a;
0334    T w = sqrt(1 + mu);
0335    T w_2 = w * w;
0336    T w_3 = w_2 * w;
0337    T w_4 = w_2 * w_2;
0338    T w_5 = w_3 * w_2;
0339    T w_6 = w_3 * w_3;
0340    T w_7 = w_4 * w_3;
0341    T w_8 = w_4 * w_4;
0342    T w_9 = w_5 * w_4;
0343    T w_10 = w_5 * w_5;
0344    T d = eta0 - mu;
0345    T d_2 = d * d;
0346    T d_3 = d_2 * d;
0347    T d_4 = d_2 * d_2;
0348    T w1 = w + 1;
0349    T w1_2 = w1 * w1;
0350    T w1_3 = w1 * w1_2;
0351    T w1_4 = w1_2 * w1_2;
0352    //
0353    // Now we need to compute the perturbation error terms that
0354    // convert eta0 to eta, these are all polynomials of polynomials.
0355    // Probably these should be re-written to use tabulated data
0356    // (see examples above), but it's less of a win in this case as we
0357    // need to calculate the individual powers for the denominator terms
0358    // anyway, so we might as well use them for the numerator-polynomials
0359    // as well....
0360    //
0361    // Refer to p154-p155 for the details of these expansions:
0362    //
0363    T e1 = (w + 2) * (w - 1) / (3 * w);
0364    e1 += (w_3 + 9 * w_2 + 21 * w + 5) * d / (36 * w_2 * w1);
0365    e1 -= (w_4 - 13 * w_3 + 69 * w_2 + 167 * w + 46) * d_2 / (1620 * w1_2 * w_3);
0366    e1 -= (7 * w_5 + 21 * w_4 + 70 * w_3 + 26 * w_2 - 93 * w - 31) * d_3 / (6480 * w1_3 * w_4);
0367    e1 -= (75 * w_6 + 202 * w_5 + 188 * w_4 - 888 * w_3 - 1345 * w_2 + 118 * w + 138) * d_4 / (272160 * w1_4 * w_5);
0368 
0369    T e2 = (28 * w_4 + 131 * w_3 + 402 * w_2 + 581 * w + 208) * (w - 1) / (1620 * w1 * w_3);
0370    e2 -= (35 * w_6 - 154 * w_5 - 623 * w_4 - 1636 * w_3 - 3983 * w_2 - 3514 * w - 925) * d / (12960 * w1_2 * w_4);
0371    e2 -= (2132 * w_7 + 7915 * w_6 + 16821 * w_5 + 35066 * w_4 + 87490 * w_3 + 141183 * w_2 + 95993 * w + 21640) * d_2  / (816480 * w_5 * w1_3);
0372    e2 -= (11053 * w_8 + 53308 * w_7 + 117010 * w_6 + 163924 * w_5 + 116188 * w_4 - 258428 * w_3 - 677042 * w_2 - 481940 * w - 105497) * d_3 / (T(14696640) * w1_4 * w_6);
0373 
0374    T e3 = -((3592 * w_7 + 8375 * w_6 - 1323 * w_5 - 29198 * w_4 - 89578 * w_3 - 154413 * w_2 - 116063 * w - 29632) * (w - 1)) / (816480 * w_5 * w1_2);
0375    e3 -= (442043 * w_9 + T(2054169) * w_8 + T(3803094) * w_7 + T(3470754) * w_6 + T(2141568) * w_5 - T(2393568) * w_4 - T(19904934) * w_3 - T(34714674) * w_2 - T(23128299) * w - T(5253353)) * d / (T(146966400) * w_6 * w1_3);
0376    e3 -= (116932 * w_10 + 819281 * w_9 + T(2378172) * w_8 + T(4341330) * w_7 + T(6806004) * w_6 + T(10622748) * w_5 + T(18739500) * w_4 + T(30651894) * w_3 + T(30869976) * w_2 + T(15431867) * w + T(2919016)) * d_2 / (T(146966400) * w1_4 * w_7);
0377    //
0378    // Combine eta0 and the error terms to compute eta (Second equation p155):
0379    //
0380    T eta = eta0 + e1 / a + e2 / (a * a) + e3 / (a * a * a);
0381    //
0382    // Now we need to solve Eq 4.2 to obtain x.  For any given value of
0383    // eta there are two solutions to this equation, and since the distribution
0384    // may be very skewed, these are not related by x ~ 1-x we used when
0385    // implementing section 3 above.  However we know that:
0386    //
0387    //  cross < x <= 1       ; iff eta < mu
0388    //          x == cross   ; iff eta == mu
0389    //     0 <= x < cross    ; iff eta > mu
0390    //
0391    // Where cross == 1 / (1 + mu)
0392    // Many thanks to Prof Temme for clarifying this point.
0393    //
0394    // Therefore we'll just jump straight into Newton iterations
0395    // to solve Eq 4.2 using these bounds, and simple bisection
0396    // as the first guess, in practice this converges pretty quickly
0397    // and we only need a few digits correct anyway:
0398    //
0399    if(eta <= 0)
0400       eta = tools::min_value<T>();
0401    T u = eta - mu * log(eta) + (1 + mu) * log(1 + mu) - mu;
0402    T cross = 1 / (1 + mu);
0403    T lower = eta < mu ? cross : 0;
0404    T upper = eta < mu ? 1 : cross;
0405    T x = (lower + upper) / 2;
0406 
0407    // Early exit for cases with numerical precision issues.
0408    if (cross == 0 || cross == 1) { return cross; }
0409    
0410    x = tools::newton_raphson_iterate(
0411       temme_root_finder<T>(u, mu), x, lower, upper, policies::digits<T, Policy>() / 2);
0412 #ifdef BOOST_INSTRUMENT
0413    std::cout << "Estimating x with Temme method 3: " << x << std::endl;
0414 #endif
0415    return x;
0416 }
0417 
0418 template <class T, class Policy>
0419 struct ibeta_roots
0420 {
0421    ibeta_roots(T _a, T _b, T t, bool inv = false)
0422       : a(_a), b(_b), target(t), invert(inv) {}
0423 
0424    boost::math::tuple<T, T, T> operator()(T x)
0425    {
0426       BOOST_MATH_STD_USING // ADL of std names
0427 
0428       BOOST_FPU_EXCEPTION_GUARD
0429 
0430       T f1;
0431       T y = 1 - x;
0432       T f = ibeta_imp(a, b, x, Policy(), invert, true, &f1) - target;
0433       if(invert)
0434          f1 = -f1;
0435       if(y == 0)
0436          y = tools::min_value<T>() * 64;
0437       if(x == 0)
0438          x = tools::min_value<T>() * 64;
0439 
0440       T f2 = f1 * (-y * a + (b - 2) * x + 1);
0441       if(fabs(f2) < y * x * tools::max_value<T>())
0442          f2 /= (y * x);
0443       if(invert)
0444          f2 = -f2;
0445 
0446       // make sure we don't have a zero derivative:
0447       if(f1 == 0)
0448          f1 = (invert ? -1 : 1) * tools::min_value<T>() * 64;
0449 
0450       return boost::math::make_tuple(f, f1, f2);
0451    }
0452 private:
0453    T a, b, target;
0454    bool invert;
0455 };
0456 
0457 template <class T, class Policy>
0458 T ibeta_inv_imp(T a, T b, T p, T q, const Policy& pol, T* py)
0459 {
0460    BOOST_MATH_STD_USING  // For ADL of math functions.
0461 
0462    //
0463    // The flag invert is set to true if we swap a for b and p for q,
0464    // in which case the result has to be subtracted from 1:
0465    //
0466    bool invert = false;
0467    //
0468    // Handle trivial cases first:
0469    //
0470    if(q == 0)
0471    {
0472       if(py) *py = 0;
0473       return 1;
0474    }
0475    else if(p == 0)
0476    {
0477       if(py) *py = 1;
0478       return 0;
0479    }
0480    else if(a == 1)
0481    {
0482       if(b == 1)
0483       {
0484          if(py) *py = 1 - p;
0485          return p;
0486       }
0487       // Change things around so we can handle as b == 1 special case below:
0488       std::swap(a, b);
0489       std::swap(p, q);
0490       invert = true;
0491    }
0492    //
0493    // Depending upon which approximation method we use, we may end up
0494    // calculating either x or y initially (where y = 1-x):
0495    //
0496    T x = 0; // Set to a safe zero to avoid a
0497    // MSVC 2005 warning C4701: potentially uninitialized local variable 'x' used
0498    // But code inspection appears to ensure that x IS assigned whatever the code path.
0499    T y;
0500 
0501    // For some of the methods we can put tighter bounds
0502    // on the result than simply [0,1]:
0503    //
0504    T lower = 0;
0505    T upper = 1;
0506    //
0507    // Student's T with b = 0.5 gets handled as a special case, swap
0508    // around if the arguments are in the "wrong" order:
0509    //
0510    if(a == 0.5f)
0511    {
0512       if(b == 0.5f)
0513       {
0514          x = sin(p * constants::half_pi<T>());
0515          x *= x;
0516          if(py)
0517          {
0518             *py = sin(q * constants::half_pi<T>());
0519             *py *= *py;
0520          }
0521          return x;
0522       }
0523       else if(b > 0.5f)
0524       {
0525          std::swap(a, b);
0526          std::swap(p, q);
0527          invert = !invert;
0528       }
0529    }
0530    //
0531    // Select calculation method for the initial estimate:
0532    //
0533    if((b == 0.5f) && (a >= 0.5f) && (p != 1))
0534    {
0535       //
0536       // We have a Student's T distribution:
0537       x = find_ibeta_inv_from_t_dist(a, p, q, &y, pol);
0538    }
0539    else if(b == 1)
0540    {
0541       if(p < q)
0542       {
0543          if(a > 1)
0544          {
0545             x = pow(p, 1 / a);
0546             y = -boost::math::expm1(log(p) / a, pol);
0547          }
0548          else
0549          {
0550             x = pow(p, 1 / a);
0551             y = 1 - x;
0552          }
0553       }
0554       else
0555       {
0556          x = exp(boost::math::log1p(-q, pol) / a);
0557          y = -boost::math::expm1(boost::math::log1p(-q, pol) / a, pol);
0558       }
0559       if(invert)
0560          std::swap(x, y);
0561       if(py)
0562          *py = y;
0563       return x;
0564    }
0565    else if(a + b > 5)
0566    {
0567       //
0568       // When a+b is large then we can use one of Prof Temme's
0569       // asymptotic expansions, begin by swapping things around
0570       // so that p < 0.5, we do this to avoid cancellations errors
0571       // when p is large.
0572       //
0573       if(p > 0.5)
0574       {
0575          std::swap(a, b);
0576          std::swap(p, q);
0577          invert = !invert;
0578       }
0579       T minv = (std::min)(a, b);
0580       T maxv = (std::max)(a, b);
0581       if((sqrt(minv) > (maxv - minv)) && (minv > 5))
0582       {
0583          //
0584          // When a and b differ by a small amount
0585          // the curve is quite symmetrical and we can use an error
0586          // function to approximate the inverse. This is the cheapest
0587          // of the three Temme expansions, and the calculated value
0588          // for x will never be much larger than p, so we don't have
0589          // to worry about cancellation as long as p is small.
0590          //
0591          x = temme_method_1_ibeta_inverse(a, b, p, pol);
0592          y = 1 - x;
0593       }
0594       else
0595       {
0596          T r = a + b;
0597          T theta = asin(sqrt(a / r));
0598          T lambda = minv / r;
0599          if((lambda >= 0.2) && (lambda <= 0.8) && (r >= 10))
0600          {
0601             //
0602             // The second error function case is the next cheapest
0603             // to use, it brakes down when the result is likely to be
0604             // very small, if a+b is also small, but we can use a
0605             // cheaper expansion there in any case.  As before x won't
0606             // be much larger than p, so as long as p is small we should
0607             // be free of cancellation error.
0608             //
0609             T ppa = pow(p, 1/a);
0610             if((ppa < 0.0025) && (a + b < 200))
0611             {
0612                x = ppa * pow(a * boost::math::beta(a, b, pol), 1/a);
0613             }
0614             else
0615                x = temme_method_2_ibeta_inverse(a, b, p, r, theta, pol);
0616             y = 1 - x;
0617          }
0618          else
0619          {
0620             //
0621             // If we get here then a and b are very different in magnitude
0622             // and we need to use the third of Temme's methods which
0623             // involves inverting the incomplete gamma.  This is much more
0624             // expensive than the other methods.  We also can only use this
0625             // method when a > b, which can lead to cancellation errors
0626             // if we really want y (as we will when x is close to 1), so
0627             // a different expansion is used in that case.
0628             //
0629             if(a < b)
0630             {
0631                std::swap(a, b);
0632                std::swap(p, q);
0633                invert = !invert;
0634             }
0635             //
0636             // Try and compute the easy way first:
0637             //
0638             T bet = 0;
0639             if (b < 2)
0640             {
0641 #ifndef BOOST_NO_EXCEPTIONS
0642                try
0643 #endif
0644                {
0645                   bet = boost::math::beta(a, b, pol);
0646 
0647                   typedef typename Policy::overflow_error_type overflow_type;
0648 
0649                   BOOST_IF_CONSTEXPR(overflow_type::value != boost::math::policies::throw_on_error)
0650                      if(bet > tools::max_value<T>())
0651                         bet = tools::max_value<T>();
0652                }
0653 #ifndef BOOST_NO_EXCEPTIONS
0654                catch (const std::overflow_error&)
0655                {
0656                   bet = tools::max_value<T>();
0657                }
0658 #endif
0659             }
0660             if(bet != 0)
0661             {
0662                y = pow(b * q * bet, 1/b);
0663                x = 1 - y;
0664             }
0665             else
0666                y = 1;
0667             if(y > 1e-5)
0668             {
0669                x = temme_method_3_ibeta_inverse(a, b, p, q, pol);
0670                y = 1 - x;
0671             }
0672          }
0673       }
0674    }
0675    else if((a < 1) && (b < 1))
0676    {
0677       //
0678       // Both a and b less than 1,
0679       // there is a point of inflection at xs:
0680       //
0681       T xs = (1 - a) / (2 - a - b);
0682       //
0683       // Now we need to ensure that we start our iteration from the
0684       // right side of the inflection point:
0685       //
0686       T fs = boost::math::ibeta(a, b, xs, pol) - p;
0687       if(fabs(fs) / p < tools::epsilon<T>() * 3)
0688       {
0689          // The result is at the point of inflection, best just return it:
0690          *py = invert ? xs : 1 - xs;
0691          return invert ? 1-xs : xs;
0692       }
0693       if(fs < 0)
0694       {
0695          std::swap(a, b);
0696          std::swap(p, q);
0697          invert = !invert;
0698          xs = 1 - xs;
0699       }
0700       if ((a < tools::min_value<T>()) && (b > tools::min_value<T>()))
0701       {
0702          if (py)
0703          {
0704             *py = invert ? 0 : 1;
0705          }
0706          return invert ? 1 : 0; // nothing interesting going on here.
0707       }
0708       //
0709       // The call to beta may overflow, plus the alternative using lgamma may do the same
0710       // if T is a type where 1/T is infinite for small values (denorms for example).
0711       //
0712       T bet = 0;
0713       T xg;
0714       bool overflow = false;
0715 #ifndef BOOST_NO_EXCEPTIONS
0716       try {
0717 #endif
0718          bet = boost::math::beta(a, b, pol);
0719 #ifndef BOOST_NO_EXCEPTIONS
0720       }
0721       catch (const std::runtime_error&)
0722       {
0723          overflow = true;
0724       }
0725 #endif
0726       if (overflow || !(boost::math::isfinite)(bet))
0727       {
0728          xg = exp((boost::math::lgamma(a + 1, pol) + boost::math::lgamma(b, pol) - boost::math::lgamma(a + b, pol) + log(p)) / a);
0729          if (xg > 2 / tools::epsilon<T>())
0730             xg = 2 / tools::epsilon<T>();
0731       }
0732       else
0733          xg = pow(a * p * bet, 1/a);
0734       x = xg / (1 + xg);
0735       y = 1 / (1 + xg);
0736       //
0737       // And finally we know that our result is below the inflection
0738       // point, so set an upper limit on our search:
0739       //
0740       if(x > xs)
0741          x = xs;
0742       upper = xs;
0743    }
0744    else if((a > 1) && (b > 1))
0745    {
0746       //
0747       // Small a and b, both greater than 1,
0748       // there is a point of inflection at xs,
0749       // and it's complement is xs2, we must always
0750       // start our iteration from the right side of the
0751       // point of inflection.
0752       //
0753       T xs = (a - 1) / (a + b - 2);
0754       T xs2 = (b - 1) / (a + b - 2);
0755       T ps = boost::math::ibeta(a, b, xs, pol) - p;
0756 
0757       if(ps < 0)
0758       {
0759          std::swap(a, b);
0760          std::swap(p, q);
0761          std::swap(xs, xs2);
0762          invert = !invert;
0763       }
0764       //
0765       // Estimate x and y, using expm1 to get a good estimate
0766       // for y when it's very small:
0767       //
0768       T lx = log(p * a * boost::math::beta(a, b, pol)) / a;
0769       x = exp(lx);
0770       y = x < 0.9 ? T(1 - x) : (T)(-boost::math::expm1(lx, pol));
0771 
0772       if((b < a) && (x < 0.2))
0773       {
0774          //
0775          // Under a limited range of circumstances we can improve
0776          // our estimate for x, frankly it's clear if this has much effect!
0777          //
0778          T ap1 = a - 1;
0779          T bm1 = b - 1;
0780          T a_2 = a * a;
0781          T a_3 = a * a_2;
0782          T b_2 = b * b;
0783          T terms[5] = { 0, 1 };
0784          terms[2] = bm1 / ap1;
0785          ap1 *= ap1;
0786          terms[3] = bm1 * (3 * a * b + 5 * b + a_2 - a - 4) / (2 * (a + 2) * ap1);
0787          ap1 *= (a + 1);
0788          terms[4] = bm1 * (33 * a * b_2 + 31 * b_2 + 8 * a_2 * b_2 - 30 * a * b - 47 * b + 11 * a_2 * b + 6 * a_3 * b + 18 + 4 * a - a_3 + a_2 * a_2 - 10 * a_2)
0789                     / (3 * (a + 3) * (a + 2) * ap1);
0790          x = tools::evaluate_polynomial(terms, x, 5);
0791       }
0792       //
0793       // And finally we know that our result is below the inflection
0794       // point, so set an upper limit on our search:
0795       //
0796       if(x > xs)
0797          x = xs;
0798       upper = xs;
0799    }
0800    else /*if((a <= 1) != (b <= 1))*/
0801    {
0802       //
0803       // If all else fails we get here, only one of a and b
0804       // is above 1, and a+b is small.  Start by swapping
0805       // things around so that we have a concave curve with b > a
0806       // and no points of inflection in [0,1].  As long as we expect
0807       // x to be small then we can use the simple (and cheap) power
0808       // term to estimate x, but when we expect x to be large then
0809       // this greatly underestimates x and leaves us trying to
0810       // iterate "round the corner" which may take almost forever...
0811       //
0812       // We could use Temme's inverse gamma function case in that case,
0813       // this works really rather well (albeit expensively) even though
0814       // strictly speaking we're outside it's defined range.
0815       //
0816       // However it's expensive to compute, and an alternative approach
0817       // which models the curve as a distorted quarter circle is much
0818       // cheaper to compute, and still keeps the number of iterations
0819       // required down to a reasonable level.  With thanks to Prof Temme
0820       // for this suggestion.
0821       //
0822       if(b < a)
0823       {
0824          std::swap(a, b);
0825          std::swap(p, q);
0826          invert = !invert;
0827       }
0828       if (a < tools::min_value<T>())
0829       {
0830          // Avoid spurious overflows for denorms:
0831          if (p < 1)
0832          {
0833             x = 1;
0834             y = 0;
0835          }
0836          else
0837          {
0838             x = 0;
0839             y = 1;
0840          }
0841       }
0842       else if(pow(p, 1/a) < 0.5)
0843       {
0844 #ifndef BOOST_NO_EXCEPTIONS
0845          try 
0846          {
0847 #endif
0848             x = pow(p * a * boost::math::beta(a, b, pol), 1 / a);
0849             if ((x > 1) || !(boost::math::isfinite)(x))
0850                x = 1;
0851 #ifndef BOOST_NO_EXCEPTIONS
0852          }
0853          catch (const std::overflow_error&)
0854          {
0855             x = 1;
0856          }
0857 #endif
0858          if(x == 0)
0859             x = boost::math::tools::min_value<T>();
0860          y = 1 - x;
0861       }
0862       else /*if(pow(q, 1/b) < 0.1)*/
0863       {
0864          // model a distorted quarter circle:
0865 #ifndef BOOST_NO_EXCEPTIONS
0866          try 
0867          {
0868 #endif
0869             y = pow(1 - pow(p, b * boost::math::beta(a, b, pol)), 1/b);
0870             if ((y > 1) || !(boost::math::isfinite)(y))
0871                y = 1;
0872 #ifndef BOOST_NO_EXCEPTIONS
0873          }
0874          catch (const std::overflow_error&)
0875          {
0876             y = 1;
0877          }
0878 #endif
0879          if(y == 0)
0880             y = boost::math::tools::min_value<T>();
0881          x = 1 - y;
0882       }
0883    }
0884 
0885    //
0886    // Now we have a guess for x (and for y) we can set things up for
0887    // iteration.  If x > 0.5 it pays to swap things round:
0888    //
0889    if(x > 0.5)
0890    {
0891       std::swap(a, b);
0892       std::swap(p, q);
0893       std::swap(x, y);
0894       invert = !invert;
0895       T l = 1 - upper;
0896       T u = 1 - lower;
0897       lower = l;
0898       upper = u;
0899    }
0900    //
0901    // lower bound for our search:
0902    //
0903    // We're not interested in denormalised answers as these tend to
0904    // these tend to take up lots of iterations, given that we can't get
0905    // accurate derivatives in this area (they tend to be infinite).
0906    //
0907    if(lower == 0)
0908    {
0909       if(invert && (py == 0))
0910       {
0911          //
0912          // We're not interested in answers smaller than machine epsilon:
0913          //
0914          lower = boost::math::tools::epsilon<T>();
0915          if(x < lower)
0916             x = lower;
0917       }
0918       else
0919          lower = boost::math::tools::min_value<T>();
0920       if(x < lower)
0921          x = lower;
0922    }
0923    std::uintmax_t max_iter = policies::get_max_root_iterations<Policy>();
0924    std::uintmax_t max_iter_used = 0;
0925    //
0926    // Figure out how many digits to iterate towards:
0927    //
0928    int digits = boost::math::policies::digits<T, Policy>() / 2;
0929    if((x < 1e-50) && ((a < 1) || (b < 1)))
0930    {
0931       //
0932       // If we're in a region where the first derivative is very
0933       // large, then we have to take care that the root-finder
0934       // doesn't terminate prematurely.  We'll bump the precision
0935       // up to avoid this, but we have to take care not to set the
0936       // precision too high or the last few iterations will just
0937       // thrash around and convergence may be slow in this case.
0938       // Try 3/4 of machine epsilon:
0939       //
0940       digits *= 3;
0941       digits /= 2;
0942    }
0943    //
0944    // Now iterate, we can use either p or q as the target here
0945    // depending on which is smaller:
0946    //
0947    x = boost::math::tools::halley_iterate(
0948       boost::math::detail::ibeta_roots<T, Policy>(a, b, (p < q ? p : q), (p < q ? false : true)), x, lower, upper, digits, max_iter);
0949    policies::check_root_iterations<T>("boost::math::ibeta<%1%>(%1%, %1%, %1%)", max_iter + max_iter_used, pol);
0950    //
0951    // We don't really want these asserts here, but they are useful for sanity
0952    // checking that we have the limits right, uncomment if you suspect bugs *only*.
0953    //
0954    //BOOST_MATH_ASSERT(x != upper);
0955    //BOOST_MATH_ASSERT((x != lower) || (x == boost::math::tools::min_value<T>()) || (x == boost::math::tools::epsilon<T>()));
0956    //
0957    // Tidy up, if we "lower" was too high then zero is the best answer we have:
0958    //
0959    if(x == lower)
0960       x = 0;
0961    if(py)
0962       *py = invert ? x : 1 - x;
0963    return invert ? 1-x : x;
0964 }
0965 
0966 } // namespace detail
0967 
0968 template <class T1, class T2, class T3, class T4, class Policy>
0969 inline typename tools::promote_args<T1, T2, T3, T4>::type
0970    ibeta_inv(T1 a, T2 b, T3 p, T4* py, const Policy& pol)
0971 {
0972    static const char* function = "boost::math::ibeta_inv<%1%>(%1%,%1%,%1%)";
0973    BOOST_FPU_EXCEPTION_GUARD
0974    typedef typename tools::promote_args<T1, T2, T3, T4>::type result_type;
0975    typedef typename policies::evaluation<result_type, Policy>::type value_type;
0976    typedef typename policies::normalise<
0977       Policy,
0978       policies::promote_float<false>,
0979       policies::promote_double<false>,
0980       policies::discrete_quantile<>,
0981       policies::assert_undefined<> >::type forwarding_policy;
0982 
0983    if(a <= 0)
0984       return policies::raise_domain_error<result_type>(function, "The argument a to the incomplete beta function inverse must be greater than zero (got a=%1%).", a, pol);
0985    if(b <= 0)
0986       return policies::raise_domain_error<result_type>(function, "The argument b to the incomplete beta function inverse must be greater than zero (got b=%1%).", b, pol);
0987    if((p < 0) || (p > 1))
0988       return policies::raise_domain_error<result_type>(function, "Argument p outside the range [0,1] in the incomplete beta function inverse (got p=%1%).", p, pol);
0989 
0990    value_type rx, ry;
0991 
0992    rx = detail::ibeta_inv_imp(
0993          static_cast<value_type>(a),
0994          static_cast<value_type>(b),
0995          static_cast<value_type>(p),
0996          static_cast<value_type>(1 - p),
0997          forwarding_policy(), &ry);
0998 
0999    if(py) *py = policies::checked_narrowing_cast<T4, forwarding_policy>(ry, function);
1000    return policies::checked_narrowing_cast<result_type, forwarding_policy>(rx, function);
1001 }
1002 
1003 template <class T1, class T2, class T3, class T4>
1004 inline typename tools::promote_args<T1, T2, T3, T4>::type
1005    ibeta_inv(T1 a, T2 b, T3 p, T4* py)
1006 {
1007    return ibeta_inv(a, b, p, py, policies::policy<>());
1008 }
1009 
1010 template <class T1, class T2, class T3>
1011 inline typename tools::promote_args<T1, T2, T3>::type
1012    ibeta_inv(T1 a, T2 b, T3 p)
1013 {
1014    typedef typename tools::promote_args<T1, T2, T3>::type result_type;
1015    return ibeta_inv(a, b, p, static_cast<result_type*>(nullptr), policies::policy<>());
1016 }
1017 
1018 template <class T1, class T2, class T3, class Policy>
1019 inline typename tools::promote_args<T1, T2, T3>::type
1020    ibeta_inv(T1 a, T2 b, T3 p, const Policy& pol)
1021 {
1022    typedef typename tools::promote_args<T1, T2, T3>::type result_type;
1023    return ibeta_inv(a, b, p, static_cast<result_type*>(nullptr), pol);
1024 }
1025 
1026 template <class T1, class T2, class T3, class T4, class Policy>
1027 inline typename tools::promote_args<T1, T2, T3, T4>::type
1028    ibetac_inv(T1 a, T2 b, T3 q, T4* py, const Policy& pol)
1029 {
1030    static const char* function = "boost::math::ibetac_inv<%1%>(%1%,%1%,%1%)";
1031    BOOST_FPU_EXCEPTION_GUARD
1032    typedef typename tools::promote_args<T1, T2, T3, T4>::type result_type;
1033    typedef typename policies::evaluation<result_type, Policy>::type value_type;
1034    typedef typename policies::normalise<
1035       Policy,
1036       policies::promote_float<false>,
1037       policies::promote_double<false>,
1038       policies::discrete_quantile<>,
1039       policies::assert_undefined<> >::type forwarding_policy;
1040 
1041    if(a <= 0)
1042       return policies::raise_domain_error<result_type>(function, "The argument a to the incomplete beta function inverse must be greater than zero (got a=%1%).", a, pol);
1043    if(b <= 0)
1044       return policies::raise_domain_error<result_type>(function, "The argument b to the incomplete beta function inverse must be greater than zero (got b=%1%).", b, pol);
1045    if((q < 0) || (q > 1))
1046       return policies::raise_domain_error<result_type>(function, "Argument q outside the range [0,1] in the incomplete beta function inverse (got q=%1%).", q, pol);
1047 
1048    value_type rx, ry;
1049 
1050    rx = detail::ibeta_inv_imp(
1051          static_cast<value_type>(a),
1052          static_cast<value_type>(b),
1053          static_cast<value_type>(1 - q),
1054          static_cast<value_type>(q),
1055          forwarding_policy(), &ry);
1056 
1057    if(py) *py = policies::checked_narrowing_cast<T4, forwarding_policy>(ry, function);
1058    return policies::checked_narrowing_cast<result_type, forwarding_policy>(rx, function);
1059 }
1060 
1061 template <class T1, class T2, class T3, class T4>
1062 inline typename tools::promote_args<T1, T2, T3, T4>::type
1063    ibetac_inv(T1 a, T2 b, T3 q, T4* py)
1064 {
1065    return ibetac_inv(a, b, q, py, policies::policy<>());
1066 }
1067 
1068 template <class RT1, class RT2, class RT3>
1069 inline typename tools::promote_args<RT1, RT2, RT3>::type
1070    ibetac_inv(RT1 a, RT2 b, RT3 q)
1071 {
1072    typedef typename tools::promote_args<RT1, RT2, RT3>::type result_type;
1073    return ibetac_inv(a, b, q, static_cast<result_type*>(nullptr), policies::policy<>());
1074 }
1075 
1076 template <class RT1, class RT2, class RT3, class Policy>
1077 inline typename tools::promote_args<RT1, RT2, RT3>::type
1078    ibetac_inv(RT1 a, RT2 b, RT3 q, const Policy& pol)
1079 {
1080    typedef typename tools::promote_args<RT1, RT2, RT3>::type result_type;
1081    return ibetac_inv(a, b, q, static_cast<result_type*>(nullptr), pol);
1082 }
1083 
1084 } // namespace math
1085 } // namespace boost
1086 
1087 #endif // BOOST_MATH_SPECIAL_FUNCTIONS_IGAMMA_INVERSE_HPP
1088 
1089 
1090 
1091