/* * Abramowitz and Stegun 6.5.29 [right] */ static double pgamma_smallx (double x, double alph, int lower_tail, int log_p) { double sum = 0, c = alph, n = 0, term; #ifdef DEBUG_p REprintf (" pg_smallx(x=%.12g, alph=%.12g): ", x, alph); #endif /* * Relative to 6.5.29 all terms have been multiplied by alph * and the first, thus being 1, is omitted. */ do { n++; c *= -x / n; term = c / (alph + n); sum += term; } while (fabs (term) > DBL_EPSILON * fabs (sum)); #ifdef DEBUG_p REprintf ("%5.0f terms --> conv.sum=%g;", n, sum); #endif if (lower_tail) { double f1 = log_p ? log1p (sum) : 1 + sum; double f2; if (alph > 1) { f2 = dpois_raw (alph, x, log_p); f2 = log_p ? f2 + x : f2 * exp (x); } else if (log_p) f2 = alph * log (x) - lgamma1p (alph); else f2 = pow (x, alph) / exp (lgamma1p (alph)); #ifdef DEBUG_p REprintf (" (f1,f2)= (%g,%g)\n", f1,f2); #endif return log_p ? f1 + f2 : f1 * f2; } else { double lf2 = alph * log (x) - lgamma1p (alph); #ifdef DEBUG_p REprintf (" 1:%.14g 2:%.14g\n", alph * log (x), lgamma1p (alph)); REprintf (" sum=%.14g log(1+sum)=%.14g lf2=%.14g\n", sum, log1p (sum), lf2); #endif if (log_p) return R_Log1_Exp (log1p (sum) + lf2); else { double f1m1 = sum; double f2m1 = expm1 (lf2); return -(f1m1 + f2m1 + f1m1 * f2m1); } } } /* pgamma_smallx() */
double lgamma_r (double x, int *signp) { *signp = +1; if (gnm_isnan (x)) return gnm_nan; if (x > 0) { gnm_float f = 1; while (x < 10) { f *= x; x++; } return (M_LN_SQRT_2PI + (x - 0.5) * gnm_log(x) - x + lgammacor(x)) - gnm_log (f); } else { gnm_float axm2 = gnm_fmod (-x, 2.0); gnm_float y = gnm_sinpi (axm2) / M_PIgnum; *signp = axm2 > 1.0 ? +1 : -1; return y == 0 ? gnm_nan : - gnm_log (gnm_abs (y)) - lgamma1p (-x); } }
gnm_float stirlerr(gnm_float n) { #define S0 GNM_const(0.083333333333333333333) /* 1/12 */ #define S1 GNM_const(0.00277777777777777777778) /* 1/360 */ #define S2 GNM_const(0.00079365079365079365079365) /* 1/1260 */ #define S3 GNM_const(0.000595238095238095238095238) /* 1/1680 */ #define S4 GNM_const(0.0008417508417508417508417508)/* 1/1188 */ /* error for 0, 0.5, 1.0, 1.5, ..., 14.5, 15.0. */ static const gnm_float sferr_halves[31] = { 0.0, /* n=0 - wrong, place holder only */ GNM_const(0.1534264097200273452913848), /* 0.5 */ GNM_const(0.0810614667953272582196702), /* 1.0 */ GNM_const(0.0548141210519176538961390), /* 1.5 */ GNM_const(0.0413406959554092940938221), /* 2.0 */ GNM_const(0.03316287351993628748511048), /* 2.5 */ GNM_const(0.02767792568499833914878929), /* 3.0 */ GNM_const(0.02374616365629749597132920), /* 3.5 */ GNM_const(0.02079067210376509311152277), /* 4.0 */ GNM_const(0.01848845053267318523077934), /* 4.5 */ GNM_const(0.01664469118982119216319487), /* 5.0 */ GNM_const(0.01513497322191737887351255), /* 5.5 */ GNM_const(0.01387612882307074799874573), /* 6.0 */ GNM_const(0.01281046524292022692424986), /* 6.5 */ GNM_const(0.01189670994589177009505572), /* 7.0 */ GNM_const(0.01110455975820691732662991), /* 7.5 */ GNM_const(0.010411265261972096497478567), /* 8.0 */ GNM_const(0.009799416126158803298389475), /* 8.5 */ GNM_const(0.009255462182712732917728637), /* 9.0 */ GNM_const(0.008768700134139385462952823), /* 9.5 */ GNM_const(0.008330563433362871256469318), /* 10.0 */ GNM_const(0.007934114564314020547248100), /* 10.5 */ GNM_const(0.007573675487951840794972024), /* 11.0 */ GNM_const(0.007244554301320383179543912), /* 11.5 */ GNM_const(0.006942840107209529865664152), /* 12.0 */ GNM_const(0.006665247032707682442354394), /* 12.5 */ GNM_const(0.006408994188004207068439631), /* 13.0 */ GNM_const(0.006171712263039457647532867), /* 13.5 */ GNM_const(0.005951370112758847735624416), /* 14.0 */ GNM_const(0.005746216513010115682023589), /* 14.5 */ GNM_const(0.005554733551962801371038690) /* 15.0 */ }; gnm_float nn; if (n <= 15.0) { nn = n + n; if (nn == (int)nn) return(sferr_halves[(int)nn]); return(lgamma1p (n ) - (n + 0.5)*gnm_log(n) + n - M_LN_SQRT_2PI); } nn = n*n; if (n>500) return((S0-S1/nn)/n); if (n> 80) return((S0-(S1-S2/nn)/nn)/n); if (n> 35) return((S0-(S1-(S2-S3/nn)/nn)/nn)/n); /* 15 < n <= 35 : */ return((S0-(S1-(S2-(S3-S4/nn)/nn)/nn)/nn)/n); }
// MM_R attribute_hidden double qchisq_appr(double p, double nu, double g /* = log Gamma(nu/2) */, logical lower_tail, logical log_p, double tol /* EPS1 */) { #define C7 4.67 #define C8 6.66 #define C9 6.73 #define C10 13.32 double alpha, a, c, ch, p1; double p2, q, t, x; /* test arguments and initialise */ #ifdef IEEE_754 if (ISNAN(p) || ISNAN(nu)) return p + nu; #endif R_Q_P01_check(p); if (nu <= 0) ML_ERR_return_NAN; alpha = 0.5 * nu;/* = [pq]gamma() shape */ c = alpha-1; if(nu < (-1.24)*(p1 = R_DT_log(p))) { /* for small chi-squared */ /* log(alpha) + g = log(alpha) + log(gamma(alpha)) = * = log(alpha*gamma(alpha)) = lgamma(alpha+1) suffers from * catastrophic cancellation when alpha << 1 */ double lgam1pa = (alpha < 0.5) ? lgamma1p(alpha) : (log(alpha) + g); ch = exp((lgam1pa + p1)/alpha + M_LN2); #ifdef DEBUG_qgamma REprintf(" small chi-sq., ch0 = %g\n", ch); #endif } else if(nu > 0.32) { /* using Wilson and Hilferty estimate */ x = qnorm(p, 0, 1, lower_tail, log_p); p1 = 2./(9*nu); ch = nu*pow(x*sqrt(p1) + 1-p1, 3); #ifdef DEBUG_qgamma REprintf(" nu > .32: Wilson-Hilferty; x = %7g\n", x); #endif /* approximation for p tending to 1: */ if( ch > 2.2*nu + 6 ) ch = -2*(R_DT_Clog(p) - c*log(0.5*ch) + g); } else { /* "small nu" : 1.24*(-log(p)) <= nu <= 0.32 */ ch = 0.4; a = R_DT_Clog(p) + g + c*M_LN2; #ifdef DEBUG_qgamma REprintf(" nu <= .32: a = %7g\n", a); #endif do { q = ch; p1 = 1. / (1+ch*(C7+ch)); p2 = ch*(C9+ch*(C8+ch)); t = -0.5 +(C7+2*ch)*p1 - (C9+ch*(C10+3*ch))/p2; ch -= (1- exp(a+0.5*ch)*p2*p1)/t; } while(fabs(q - ch) > tol * fabs(ch)); } return ch; }