double qf(double p, double df1, double df2, int lower_tail, int log_p) { #ifdef IEEE_754 if (ISNAN(p) || ISNAN(df1) || ISNAN(df2)) return p + df1 + df2; #endif if (df1 <= 0. || df2 <= 0.) ML_ERR_return_NAN; R_Q_P01_boundaries(p, 0, ML_POSINF); /* fudge the extreme DF cases -- qbeta doesn't do this well. But we still need to fudge the infinite ones. */ if (df1 <= df2 && df2 > 4e5) { if(!R_FINITE(df1)) /* df1 == df2 == Inf : */ return 1.; /* else */ return qchisq(p, df1, lower_tail, log_p) / df1; } if (df1 > 4e5) { /* and so df2 < df1 */ return df2 / qchisq(p, df2, !lower_tail, log_p); } p = (1. / qbeta(p, df2/2, df1/2, !lower_tail, log_p) - 1.) * (df2 / df1); return ML_VALID(p) ? p : ML_NAN; }
double qlnorm(double p, double meanlog, double sdlog, int lower_tail, int log_p) { #ifdef IEEE_754 if (ISNAN(p) || ISNAN(meanlog) || ISNAN(sdlog)) return p + meanlog + sdlog; #endif R_Q_P01_boundaries(p, 0, ML_POSINF); return exp(qnorm(p, meanlog, sdlog, lower_tail, log_p)); }
double qpois(double p, double lambda, int lower_tail, int log_p) { double mu, sigma, gamma, z, y; #ifdef IEEE_754 if (ISNAN(p) || ISNAN(lambda)) return p + lambda; #endif if(!R_FINITE(lambda)) ML_ERR_return_NAN; if(lambda < 0) ML_ERR_return_NAN; if(lambda == 0) return 0; R_Q_P01_boundaries(p, 0, ML_POSINF); mu = lambda; sigma = sqrt(lambda); /* gamma = sigma; PR#8058 should be kurtosis which is mu^-0.5 */ gamma = 1.0/sigma; /* Note : "same" code in qpois.c, qbinom.c, qnbinom.c -- * FIXME: This is far from optimal [cancellation for p ~= 1, etc]: */ if(!lower_tail || log_p) { p = R_DT_qIv(p); /* need check again (cancellation!): */ if (p == 0.) return 0; if (p == 1.) return ML_POSINF; } /* temporary hack --- FIXME --- */ if (p + 1.01*DBL_EPSILON >= 1.) return ML_POSINF; /* y := approx.value (Cornish-Fisher expansion) : */ z = qnorm(p, 0., 1., /*lower_tail*/TRUE, /*log_p*/FALSE); #ifdef HAVE_NEARBYINT y = nearbyint(mu + sigma * (z + gamma * (z*z - 1) / 6)); #else y = round(mu + sigma * (z + gamma * (z*z - 1) / 6)); #endif z = ppois(y, lambda, /*lower_tail*/TRUE, /*log_p*/FALSE); /* fuzz to ensure left continuity; 1 - 1e-7 may lose too much : */ p *= 1 - 64*DBL_EPSILON; /* If the mean is not too large a simple search is OK */ if(lambda < 1e5) return do_search(y, &z, p, lambda, 1); /* Otherwise be a bit cleverer in the search */ { double incr = floor(y * 0.001), oldincr; do { oldincr = incr; y = do_search(y, &z, p, lambda, incr); incr = fmax2(1, floor(incr/100)); } while(oldincr > 1 && incr > lambda*1e-15); return y; } }
double qnbinom(double p, double size, double prob, int lower_tail, int log_p) { double P, Q, mu, sigma, gamma, z, y; #ifdef IEEE_754 if (ISNAN(p) || ISNAN(size) || ISNAN(prob)) return p + size + prob; #endif if (prob <= 0 || prob > 1 || size <= 0) ML_ERR_return_NAN; /* FIXME: size = 0 is well defined ! */ if (prob == 1) return 0; R_Q_P01_boundaries(p, 0, ML_POSINF); Q = 1.0 / prob; P = (1.0 - prob) * Q; mu = size * P; sigma = sqrt(size * P * Q); gamma = (Q + P)/sigma; /* Note : "same" code in qpois.c, qbinom.c, qnbinom.c -- * FIXME: This is far from optimal [cancellation for p ~= 1, etc]: */ if(!lower_tail || log_p) { p = R_DT_qIv(p); /* need check again (cancellation!): */ if (p == R_DT_0) return 0; if (p == R_DT_1) return ML_POSINF; } /* temporary hack --- FIXME --- */ if (p + 1.01*DBL_EPSILON >= 1.) return ML_POSINF; /* y := approx.value (Cornish-Fisher expansion) : */ z = qnorm(p, 0., 1., /*lower_tail*/TRUE, /*log_p*/FALSE); y = floor(mu + sigma * (z + gamma * (z*z - 1) / 6) + 0.5); z = pnbinom(y, size, prob, /*lower_tail*/TRUE, /*log_p*/FALSE); /* fuzz to ensure left continuity: */ p *= 1 - 64*DBL_EPSILON; /* If the C-F value is not too large a simple search is OK */ if(y < 1e5) return do_search(y, &z, p, size, prob, 1); /* Otherwise be a bit cleverer in the search */ { double incr = floor(y * 0.001), oldincr; do { oldincr = incr; y = do_search(y, &z, p, size, prob, incr); incr = fmax2(1, floor(incr/100)); } while(oldincr > 1 && incr > y*1e-15); return y; } }
double qgeom(double p, double prob, int lower_tail, int log_p) { if (prob <= 0 || prob > 1) ML_ERR_return_NAN; R_Q_P01_boundaries(p, 0, ML_POSINF); #ifdef IEEE_754 if (ISNAN(p) || ISNAN(prob)) return p + prob; #endif if (prob == 1) return(0); /* add a fuzz to ensure left continuity, but value must be >= 0 */ return fmax2(0, ceil(R_DT_Clog(p) / log1p(- prob) - 1 - 1e-12)); }
double qnt(double p, double df, double ncp, int lower_tail, int log_p) { const static double accu = 1e-13; const static double Eps = 1e-11; /* must be > accu */ double ux, lx, nx, pp; #ifdef IEEE_754 if (ISNAN(p) || ISNAN(df) || ISNAN(ncp)) return p + df + ncp; #endif if (!R_FINITE(df)) ML_ERR_return_NAN; /* Was * df = floor(df + 0.5); * if (df < 1 || ncp < 0) ML_ERR_return_NAN; */ if (df <= 0.0) ML_ERR_return_NAN; if(ncp == 0.0 && df >= 1.0) return qt(p, df, lower_tail, log_p); R_Q_P01_boundaries(p, ML_NEGINF, ML_POSINF); p = R_DT_qIv(p); /* Invert pnt(.) : * 1. finding an upper and lower bound */ if(p > 1 - DBL_EPSILON) return ML_POSINF; pp = fmin2(1 - DBL_EPSILON, p * (1 + Eps)); for(ux = fmax2(1., ncp); ux < DBL_MAX && pnt(ux, df, ncp, TRUE, FALSE) < pp; ux *= 2); pp = p * (1 - Eps); for(lx = fmin2(-1., -ncp); lx > -DBL_MAX && pnt(lx, df, ncp, TRUE, FALSE) > pp; lx *= 2); /* 2. interval (lx,ux) halving : */ do { nx = 0.5 * (lx + ux); if (pnt(nx, df, ncp, TRUE, FALSE) > p) ux = nx; else lx = nx; } while ((ux - lx) / fabs(nx) > accu); return 0.5 * (lx + ux); }
double qnf(double p, double df1, double df2, double ncp, int lower_tail, int log_p) { double y; #ifdef IEEE_754 if (ISNAN(p) || ISNAN(df1) || ISNAN(df2) || ISNAN(ncp)) return p + df1 + df2 + ncp; #endif if (df1 <= 0. || df2 <= 0. || ncp < 0) ML_ERR_return_NAN; if (!R_FINITE(ncp)) ML_ERR_return_NAN; if (!R_FINITE(df1) && !R_FINITE(df2)) ML_ERR_return_NAN; R_Q_P01_boundaries(p, 0, ML_POSINF); if (df2 > 1e8) /* avoid problems with +Inf and loss of accuracy */ return qnchisq(p, df1, ncp, lower_tail, log_p)/df1; y = qnbeta(p, df1 / 2., df2 / 2., ncp, lower_tail, log_p); return y/(1-y) * (df2/df1); }
double qnbeta(double p, double a, double b, double ncp, int lower_tail, int log_p) { const static double accu = 1e-15; const static double Eps = 1e-14; /* must be > accu */ double ux, lx, nx, pp; #ifdef IEEE_754 if (ISNAN(p) || ISNAN(a) || ISNAN(b) || ISNAN(ncp)) return p + a + b + ncp; #endif if (!R_FINITE(a)) ML_ERR_return_NAN; if (ncp < 0. || a <= 0. || b <= 0.) ML_ERR_return_NAN; R_Q_P01_boundaries(p, 0, 1); p = R_DT_qIv(p); /* Invert pnbeta(.) : * 1. finding an upper and lower bound */ if(p > 1 - DBL_EPSILON) return 1.0; pp = fmin2(1 - DBL_EPSILON, p * (1 + Eps)); for(ux = 0.5; ux < 1 - DBL_EPSILON && pnbeta(ux, a, b, ncp, TRUE, FALSE) < pp; ux = 0.5*(1+ux)); pp = p * (1 - Eps); for(lx = 0.5; lx > DBL_MIN && pnbeta(lx, a, b, ncp, TRUE, FALSE) > pp; lx *= 0.5); /* 2. interval (lx,ux) halving : */ do { nx = 0.5 * (lx + ux); if (pnbeta(nx, a, b, ncp, TRUE, FALSE) > p) ux = nx; else lx = nx; } while ((ux - lx) / nx > accu); return 0.5 * (ux + lx); }
double qbinom(double p, double n, double pr, int lower_tail, int log_p) { double q, mu, sigma, gamma, z, y; #ifdef IEEE_754 if (ISNAN(p) || ISNAN(n) || ISNAN(pr)) return p + n + pr; #endif if(!R_FINITE(n) || !R_FINITE(pr)) ML_ERR_return_NAN; /* if log_p is true, p = -Inf is a legitimate value */ if(!R_FINITE(p) && !log_p) ML_ERR_return_NAN; if(n != floor(n + 0.5)) ML_ERR_return_NAN; if (pr < 0 || pr > 1 || n < 0) ML_ERR_return_NAN; R_Q_P01_boundaries(p, 0, n); if (pr == 0. || n == 0) return 0.; q = 1 - pr; if(q == 0.) return n; /* covers the full range of the distribution */ mu = n * pr; sigma = sqrt(n * pr * q); gamma = (q - pr) / sigma; #ifdef DEBUG_qbinom REprintf("qbinom(p=%7g, n=%g, pr=%7g, l.t.=%d, log=%d): sigm=%g, gam=%g\n", p,n,pr, lower_tail, log_p, sigma, gamma); #endif /* Note : "same" code in qpois.c, qbinom.c, qnbinom.c -- * FIXME: This is far from optimal [cancellation for p ~= 1, etc]: */ if(!lower_tail || log_p) { p = R_DT_qIv(p); /* need check again (cancellation!): */ if (p == 0.) return 0.; if (p == 1.) return n; } /* temporary hack --- FIXME --- */ if (p + 1.01*DBL_EPSILON >= 1.) return n; /* y := approx.value (Cornish-Fisher expansion) : */ z = qnorm(p, 0., 1., /*lower_tail*/TRUE, /*log_p*/FALSE); y = floor(mu + sigma * (z + gamma * (z*z - 1) / 6) + 0.5); if(y > n) /* way off */ y = n; #ifdef DEBUG_qbinom REprintf(" new (p,1-p)=(%7g,%7g), z=qnorm(..)=%7g, y=%5g\n", p, 1-p, z, y); #endif z = pbinom(y, n, pr, /*lower_tail*/TRUE, /*log_p*/FALSE); /* fuzz to ensure left continuity: */ p *= 1 - 64*DBL_EPSILON; if(n < 1e5) return do_search(y, &z, p, n, pr, 1); /* Otherwise be a bit cleverer in the search */ { double incr = floor(n * 0.001), oldincr; do { oldincr = incr; y = do_search(y, &z, p, n, pr, incr); incr = fmax2(1, floor(incr/100)); } while(oldincr > 1 && incr > n*1e-15); return y; } }
double qnorm5(double p, double mu, double sigma, int lower_tail, int log_p) { double p_, q, r, val; #ifdef IEEE_754 if (ISNAN(p) || ISNAN(mu) || ISNAN(sigma)) return p + mu + sigma; #endif R_Q_P01_boundaries(p, ML_NEGINF, ML_POSINF); if(sigma < 0) ML_ERR_return_NAN; if(sigma == 0) return mu; p_ = R_DT_qIv(p);/* real lower_tail prob. p */ q = p_ - 0.5; #ifdef DEBUG_qnorm REprintf("qnorm(p=%10.7g, m=%g, s=%g, l.t.= %d, log= %d): q = %g\n", p,mu,sigma, lower_tail, log_p, q); #endif /*-- use AS 241 --- */ /* double ppnd16_(double *p, long *ifault)*/ /* ALGORITHM AS241 APPL. STATIST. (1988) VOL. 37, NO. 3 Produces the normal deviate Z corresponding to a given lower tail area of P; Z is accurate to about 1 part in 10**16. (original fortran code used PARAMETER(..) for the coefficients and provided hash codes for checking them...) */ if (fabs(q) <= .425) {/* 0.075 <= p <= 0.925 */ r = .180625 - q * q; val = q * (((((((r * 2509.0809287301226727 + 33430.575583588128105) * r + 67265.770927008700853) * r + 45921.953931549871457) * r + 13731.693765509461125) * r + 1971.5909503065514427) * r + 133.14166789178437745) * r + 3.387132872796366608) / (((((((r * 5226.495278852854561 + 28729.085735721942674) * r + 39307.89580009271061) * r + 21213.794301586595867) * r + 5394.1960214247511077) * r + 687.1870074920579083) * r + 42.313330701600911252) * r + 1.); } else { /* closer than 0.075 from {0,1} boundary */ /* r = min(p, 1-p) < 0.075 */ if (q > 0) r = R_DT_CIv(p);/* 1-p */ else r = p_;/* = R_DT_Iv(p) ^= p */ r = sqrt(- ((log_p && ((lower_tail && q <= 0) || (!lower_tail && q > 0))) ? p : /* else */ log(r))); /* r = sqrt(-log(r)) <==> min(p, 1-p) = exp( - r^2 ) */ #ifdef DEBUG_qnorm REprintf("\t close to 0 or 1: r = %7g\n", r); #endif if (r <= 5.) { /* <==> min(p,1-p) >= exp(-25) ~= 1.3888e-11 */ r += -1.6; val = (((((((r * 7.7454501427834140764e-4 + .0227238449892691845833) * r + .24178072517745061177) * r + 1.27045825245236838258) * r + 3.64784832476320460504) * r + 5.7694972214606914055) * r + 4.6303378461565452959) * r + 1.42343711074968357734) / (((((((r * 1.05075007164441684324e-9 + 5.475938084995344946e-4) * r + .0151986665636164571966) * r + .14810397642748007459) * r + .68976733498510000455) * r + 1.6763848301838038494) * r + 2.05319162663775882187) * r + 1.); } else { /* very close to 0 or 1 */ r += -5.; val = (((((((r * 2.01033439929228813265e-7 + 2.71155556874348757815e-5) * r + .0012426609473880784386) * r + .026532189526576123093) * r + .29656057182850489123) * r + 1.7848265399172913358) * r + 5.4637849111641143699) * r + 6.6579046435011037772) / (((((((r * 2.04426310338993978564e-15 + 1.4215117583164458887e-7)* r + 1.8463183175100546818e-5) * r + 7.868691311456132591e-4) * r + .0148753612908506148525) * r + .13692988092273580531) * r + .59983220655588793769) * r + 1.); } if(q < 0.0) val = -val; /* return (q >= 0.)? r : -r ;*/ } return mu + sigma * val; }
double qnchisq(double p, double df, double ncp, int lower_tail, int log_p) { const static double accu = 1e-13; const static double racc = 4*DBL_EPSILON; /* these two are for the "search" loops, can have less accuracy: */ const static double Eps = 1e-11; /* must be > accu */ const static double rEps= 1e-10; /* relative tolerance ... */ double ux, lx, ux0, nx, pp; #ifdef IEEE_754 if (ISNAN(p) || ISNAN(df) || ISNAN(ncp)) return p + df + ncp; #endif if (!R_FINITE(df)) ML_ERR_return_NAN; /* Was * df = floor(df + 0.5); * if (df < 1 || ncp < 0) ML_ERR_return_NAN; */ if (df < 0 || ncp < 0) ML_ERR_return_NAN; R_Q_P01_boundaries(p, 0, ML_POSINF); /* Invert pnchisq(.) : * 1. finding an upper and lower bound */ { /* This is Pearson's (1959) approximation, which is usually good to 4 figs or so. */ double b, c, ff; b = (ncp*ncp)/(df + 3*ncp); c = (df + 3*ncp)/(df + 2*ncp); ff = (df + 2 * ncp)/(c*c); ux = b + c * qchisq(p, ff, lower_tail, log_p); if(ux < 0) ux = 1; ux0 = ux; } p = R_D_qIv(p); if(!lower_tail && ncp >= 80) { /* pnchisq is only for lower.tail = TRUE */ if(p < 1e-10) ML_ERROR(ME_PRECISION, "qnchisq"); p = 1. - p; lower_tail = TRUE; } if(lower_tail) { if(p > 1 - DBL_EPSILON) return ML_POSINF; pp = fmin2(1 - DBL_EPSILON, p * (1 + Eps)); for(; ux < DBL_MAX && pnchisq_raw(ux, df, ncp, Eps, rEps, 10000, TRUE) < pp; ux *= 2); pp = p * (1 - Eps); for(lx = fmin2(ux0, DBL_MAX); lx > DBL_MIN && pnchisq_raw(lx, df, ncp, Eps, rEps, 10000, TRUE) > pp; lx *= 0.5); } else { if(p > 1 - DBL_EPSILON) return 0.0; pp = fmin2(1 - DBL_EPSILON, p * (1 + Eps)); for(; ux < DBL_MAX && pnchisq_raw(ux, df, ncp, Eps, rEps, 10000, FALSE) > pp; ux *= 2); pp = p * (1 - Eps); for(lx = fmin2(ux0, DBL_MAX); lx > DBL_MIN && pnchisq_raw(lx, df, ncp, Eps, rEps, 10000, FALSE) < pp; lx *= 0.5); } /* 2. interval (lx,ux) halving : */ if(lower_tail) { do { nx = 0.5 * (lx + ux); if (pnchisq_raw(nx, df, ncp, accu, racc, 100000, TRUE) > p) ux = nx; else lx = nx; } while ((ux - lx) / nx > accu); } else { do { nx = 0.5 * (lx + ux); if (pnchisq_raw(nx, df, ncp, accu, racc, 100000, FALSE) < p) ux = nx; else lx = nx; } while ((ux - lx) / nx > accu); } return 0.5 * (ux + lx); }
double qbeta(double alpha, double p, double q, int lower_tail, int log_p) { int swap_tail, i_pb, i_inn; double a, adj, logbeta, g, h, pp, p_, prev, qq, r, s, t, tx, w, y, yprev; double acu; volatile double xinbta; /* test for admissibility of parameters */ if (isnan(p) || isnan(q) || isnan(alpha)){ return p + q + alpha; } if(p < 0. || q < 0.){ report_error("shape parameters for qbeta must be > 0."); } R_Q_P01_boundaries(alpha, 0, 1); p_ = R_DT_qIv(alpha);/* lower_tail prob (in any case) */ if(log_p && (p_ == 0. || p_ == 1.)) return p_; /* better than NaN or infinite loop; FIXME: suboptimal, since -Inf < alpha ! */ /* initialize */ logbeta = lbeta(p, q); /* change tail if necessary; afterwards 0 < a <= 1/2 */ if (p_ <= 0.5) { a = p_; pp = p; qq = q; swap_tail = 0; } else { /* change tail, swap p <-> q :*/ a = (!lower_tail && !log_p)? alpha : 1 - p_; pp = q; qq = p; swap_tail = 1; } /* calculate the initial approximation */ /* y := {fast approximation of} qnorm(1 - a) :*/ r = sqrt(-2 * log(a)); y = r - (const1 + const2 * r) / (1. + (const3 + const4 * r) * r); if (pp > 1 && qq > 1) { r = (y * y - 3.) / 6.; s = 1. / (pp + pp - 1.); t = 1. / (qq + qq - 1.); h = 2. / (s + t); w = y * sqrt(h + r) / h - (t - s) * (r + 5. / 6. - 2. / (3. * h)); xinbta = pp / (pp + qq * exp(w + w)); } else { r = qq + qq; t = 1. / (9. * qq); t = r * pow(1. - t + y * sqrt(t), 3.0); if (t <= 0.) xinbta = 1. - exp((log1p(-a)+ log(qq) + logbeta) / qq); else { t = (4. * pp + r - 2.) / t; if (t <= 1.) xinbta = exp((log(a * pp) + logbeta) / pp); else xinbta = 1. - 2. / (t + 1.); } } /* solve for x by a modified newton-raphson method, */ /* using the function pbeta_raw */ r = 1 - pp; t = 1 - qq; yprev = 0.; adj = 1; /* Sometimes the approximation is negative! */ if (xinbta < lower) xinbta = 0.5; else if (xinbta > upper) xinbta = 0.5; /* Desired accuracy should depend on (a,p) * This is from Remark .. on AS 109, adapted. * However, it's not clear if this is "optimal" for IEEE double prec. * acu = std::max<double>(acu_min, pow(10., -25. - 5./(pp * pp) - 1./(a * a))); * NEW: 'acu' accuracy NOT for squared adjustment, but simple; * ---- i.e., "new acu" = sqrt(old acu) */ acu = std::max<double>(acu_min, pow(10., -13 - 2.5/(pp * pp) - 0.5/(a * a))); tx = prev = 0.; /* keep -Wall happy */ for (i_pb=0; i_pb < 1000; i_pb++) { y = pbeta_raw(xinbta, pp, qq, /*lower_tail = */ true, false); if(!std::isfinite(y)){ report_error("algorithm blew up ni qbeta"); } y = (y - a) * exp(logbeta + r * log(xinbta) + t * log1p(-xinbta)); if (y * yprev <= 0.) prev = std::max<double>(fabs(adj),fpu); g = 1; for (i_inn=0; i_inn < 1000;i_inn++) { adj = g * y; if (fabs(adj) < prev) { tx = xinbta - adj; /* trial new x */ if (tx >= 0. && tx <= 1) { if (prev <= acu) goto L_converged; if (fabs(y) <= acu) goto L_converged; if (tx != 0. && tx != 1) break; } } g /= 3; } if (fabs(tx - xinbta) < 1e-15*xinbta) goto L_converged; xinbta = tx; yprev = y; } /*-- NOT converged: Iteration count --*/ report_error("algorithm did not converge in qbeta"); L_converged: return swap_tail ? 1 - xinbta : xinbta; }
double qtukey(double p, double rr, double cc, double df, int lower_tail, int log_p) { const static double eps = 0.0001; const int maxiter = 50; double ans = 0.0, valx0, valx1, x0, x1, xabs; int iter; #ifdef IEEE_754 if (ISNAN(p) || ISNAN(rr) || ISNAN(cc) || ISNAN(df)) { ML_ERROR(ME_DOMAIN, "qtukey"); return p + rr + cc + df; } #endif /* df must be > 1 ; there must be at least two values */ if (df < 2 || rr < 1 || cc < 2) ML_ERR_return_NAN; R_Q_P01_boundaries(p, 0, ML_POSINF); p = R_DT_qIv(p); /* lower_tail,non-log "p" */ /* Initial value */ x0 = qinv(p, cc, df); /* Find prob(value < x0) */ valx0 = ptukey(x0, rr, cc, df, /*LOWER*/TRUE, /*LOG_P*/FALSE) - p; /* Find the second iterate and prob(value < x1). */ /* If the first iterate has probability value */ /* exceeding p then second iterate is 1 less than */ /* first iterate; otherwise it is 1 greater. */ if (valx0 > 0.0) x1 = fmax2(0.0, x0 - 1.0); else x1 = x0 + 1.0; valx1 = ptukey(x1, rr, cc, df, /*LOWER*/TRUE, /*LOG_P*/FALSE) - p; /* Find new iterate */ for(iter=1 ; iter < maxiter ; iter++) { ans = x1 - ((valx1 * (x1 - x0)) / (valx1 - valx0)); valx0 = valx1; /* New iterate must be >= 0 */ x0 = x1; if (ans < 0.0) { ans = 0.0; valx1 = -p; } /* Find prob(value < new iterate) */ valx1 = ptukey(ans, rr, cc, df, /*LOWER*/TRUE, /*LOG_P*/FALSE) - p; x1 = ans; /* If the difference between two successive */ /* iterates is less than eps, stop */ xabs = fabs(x1 - x0); if (xabs < eps) return ans; } /* The process did not converge in 'maxiter' iterations */ ML_ERROR(ME_NOCONV, "qtukey"); return ans; }
double qt(double p, double ndf, int lower_tail, int log_p) { const static double eps = 1.e-12; double P, q; Rboolean neg; #ifdef IEEE_754 if (ISNAN(p) || ISNAN(ndf)) return p + ndf; #endif R_Q_P01_boundaries(p, ML_NEGINF, ML_POSINF); if (ndf <= 0) ML_ERR_return_NAN; if (ndf < 1) { /* based on qnt */ const static double accu = 1e-13; const static double Eps = 1e-11; /* must be > accu */ double ux, lx, nx, pp; int iter = 0; p = R_DT_qIv(p); /* Invert pt(.) : * 1. finding an upper and lower bound */ if(p > 1 - DBL_EPSILON) return ML_POSINF; pp = fmin2(1 - DBL_EPSILON, p * (1 + Eps)); for(ux = 1.; ux < DBL_MAX && pt(ux, ndf, TRUE, FALSE) < pp; ux *= 2); pp = p * (1 - Eps); for(lx =-1.; lx > -DBL_MAX && pt(lx, ndf, TRUE, FALSE) > pp; lx *= 2); /* 2. interval (lx,ux) halving regula falsi failed on qt(0.1, 0.1) */ do { nx = 0.5 * (lx + ux); if (pt(nx, ndf, TRUE, FALSE) > p) ux = nx; else lx = nx; } while ((ux - lx) / fabs(nx) > accu && ++iter < 1000); if(iter >= 1000) ML_ERROR(ME_PRECISION, "qt"); return 0.5 * (lx + ux); } /* Old comment: * FIXME: "This test should depend on ndf AND p !! * ----- and in fact should be replaced by * something like Abramowitz & Stegun 26.7.5 (p.949)" * * That would say that if the qnorm value is x then * the result is about x + (x^3+x)/4df + (5x^5+16x^3+3x)/96df^2 * The differences are tiny even if x ~ 1e5, and qnorm is not * that accurate in the extreme tails. */ if (ndf > 1e20) return qnorm(p, 0., 1., lower_tail, log_p); P = R_D_qIv(p); /* if exp(p) underflows, we fix below */ neg = (!lower_tail || P < 0.5) && (lower_tail || P > 0.5); if(neg) P = 2 * (log_p ? (lower_tail ? P : -expm1(p)) : R_D_Lval(p)); else P = 2 * (log_p ? (lower_tail ? -expm1(p) : P) : R_D_Cval(p)); /* 0 <= P <= 1 ; P = 2*min(P', 1 - P') in all cases */ /* Use this if(log_p) only : */ #define P_is_exp_2p (lower_tail == neg) /* both TRUE or FALSE == !xor */ if (fabs(ndf - 2) < eps) { /* df ~= 2 */ if(P > DBL_MIN) { if(3* P < DBL_EPSILON) /* P ~= 0 */ q = 1 / sqrt(P); else if (P > 0.9) /* P ~= 1 */ q = (1 - P) * sqrt(2 /(P * (2 - P))); else /* eps/3 <= P <= 0.9 */ q = sqrt(2 / (P * (2 - P)) - 2); } else { /* P << 1, q = 1/sqrt(P) = ... */ if(log_p) q = P_is_exp_2p ? exp(- p/2) / M_SQRT2 : 1/sqrt(-expm1(p)); else q = ML_POSINF; } } else if (ndf < 1 + eps) { /* df ~= 1 (df < 1 excluded above): Cauchy */ if(P > 0) q = 1/tan(P * M_PI_2);/* == - tan((P+1) * M_PI_2) -- suffers for P ~= 0 */ else { /* P = 0, but maybe = 2*exp(p) ! */ if(log_p) /* 1/tan(e) ~ 1/e */ q = P_is_exp_2p ? M_1_PI * exp(-p) : -1./(M_PI * expm1(p)); else q = ML_POSINF; } } else { /*-- usual case; including, e.g., df = 1.1 */ double x = 0., y, log_P2 = 0./* -Wall */, a = 1 / (ndf - 0.5), b = 48 / (a * a), c = ((20700 * a / b - 98) * a - 16) * a + 96.36, d = ((94.5 / (b + c) - 3) / b + 1) * sqrt(a * M_PI_2) * ndf; Rboolean P_ok1 = P > DBL_MIN || !log_p, P_ok = P_ok1; if(P_ok1) { y = pow(d * P, 2 / ndf); P_ok = (y >= DBL_EPSILON); } if(!P_ok) { /* log_p && P very small */ log_P2 = P_is_exp_2p ? p : R_Log1_Exp(p); /* == log(P / 2) */ x = (log(d) + M_LN2 + log_P2) / ndf; y = exp(2 * x); } if ((ndf < 2.1 && P > 0.5) || y > 0.05 + a) { /* P > P0(df) */ /* Asymptotic inverse expansion about normal */ if(P_ok) x = qnorm(0.5 * P, 0., 1., /*lower_tail*/TRUE, /*log_p*/FALSE); else /* log_p && P underflowed */ x = qnorm(log_P2, 0., 1., lower_tail, /*log_p*/ TRUE); y = x * x; if (ndf < 5) c += 0.3 * (ndf - 4.5) * (x + 0.6); c = (((0.05 * d * x - 5) * x - 7) * x - 2) * x + b + c; y = (((((0.4 * y + 6.3) * y + 36) * y + 94.5) / c - y - 3) / b + 1) * x; y = expm1(a * y * y); q = sqrt(ndf * y); } else { /* re-use 'y' from above */ if(!P_ok && x < - M_LN2 * DBL_MANT_DIG) {/* 0.5* log(DBL_EPSILON) */ /* y above might have underflown */ q = sqrt(ndf) * exp(-x); } else { y = ((1 / (((ndf + 6) / (ndf * y) - 0.089 * d - 0.822) * (ndf + 2) * 3) + 0.5 / (ndf + 4)) * y - 1) * (ndf + 1) / (ndf + 2) + 1 / y; q = sqrt(ndf * y); } } /* Now apply 2-term Taylor expansion improvement (1-term = Newton): * as by Hill (1981) [ref.above] */ /* FIXME: This can be far from optimal when log_p = TRUE * but is still needed, e.g. for qt(-2, df=1.01, log=TRUE). * Probably also improvable when lower_tail = FALSE */ if(P_ok1) { int it=0; while(it++ < 10 && (y = dt(q, ndf, FALSE)) > 0 && R_FINITE(x = (pt(q, ndf, FALSE, FALSE) - P/2) / y) && fabs(x) > 1e-14*fabs(q)) /* Newton (=Taylor 1 term): * q += x; * Taylor 2-term : */ q += x * (1. + x * q * (ndf + 1) / (2 * (q * q + ndf))); } } if(neg) q = -q; return q; }