/* the full definition of C_L(eta) for any valid L and eta * [Abramowitz and Stegun 14.1.7] * This depends on the complex gamma function. For large * arguments the phase of the complex gamma function is not * very accurately determined. However the modulus is, and that * is all that we need to calculate C_L. * * This is not valid for L <= -3/2 or L = -1. */ static int CLeta(double L, double eta, gsl_sf_result * result) { gsl_sf_result ln1; /* log of numerator Gamma function */ gsl_sf_result ln2; /* log of denominator Gamma function */ double sgn = 1.0; double arg_val, arg_err; if(fabs(eta/(L+1.0)) < GSL_DBL_EPSILON) { gsl_sf_lngamma_e(L+1.0, &ln1); } else { gsl_sf_result p1; /* phase of numerator Gamma -- not used */ gsl_sf_lngamma_complex_e(L+1.0, eta, &ln1, &p1); /* should be ok */ } gsl_sf_lngamma_e(2.0*(L+1.0), &ln2); if(L < -1.0) sgn = -sgn; arg_val = L*M_LN2 - 0.5*eta*M_PI + ln1.val - ln2.val; arg_err = ln1.err + ln2.err; arg_err += GSL_DBL_EPSILON * (fabs(L*M_LN2) + fabs(0.5*eta*M_PI)); return gsl_sf_exp_err_e(arg_val, arg_err, result); }
static int pow_omx(const double x, const double p, gsl_sf_result * result) { double ln_omx; double ln_result; if(fabs(x) < GSL_ROOT5_DBL_EPSILON) { ln_omx = -x*(1.0 + x*(1.0/2.0 + x*(1.0/3.0 + x/4.0 + x*x/5.0))); } else { ln_omx = log(1.0-x); } ln_result = p * ln_omx; return gsl_sf_exp_err_e(ln_result, GSL_DBL_EPSILON * fabs(ln_result), result); }
/* asymptotic expansion * j + 2.0 > 0.0 */ static int fd_asymp(const double j, const double x, gsl_sf_result * result) { const int j_integer = ( fabs(j - floor(j+0.5)) < 100.0*GSL_DBL_EPSILON ); const int itmax = 200; gsl_sf_result lg; int stat_lg = gsl_sf_lngamma_e(j + 2.0, &lg); double seqn_val = 0.5; double seqn_err = 0.0; double xm2 = (1.0/x)/x; double xgam = 1.0; double add = GSL_DBL_MAX; double cos_term; double ln_x; double ex_term_1; double ex_term_2; gsl_sf_result fneg; gsl_sf_result ex_arg; gsl_sf_result ex; int stat_fneg; int stat_e; int n; for(n=1; n<=itmax; n++) { double add_previous = add; gsl_sf_result eta; gsl_sf_eta_int_e(2*n, &eta); xgam = xgam * xm2 * (j + 1.0 - (2*n-2)) * (j + 1.0 - (2*n-1)); add = eta.val * xgam; if(!j_integer && fabs(add) > fabs(add_previous)) break; if(fabs(add/seqn_val) < GSL_DBL_EPSILON) break; seqn_val += add; seqn_err += 2.0 * GSL_DBL_EPSILON * fabs(add); } seqn_err += fabs(add); stat_fneg = fd_neg(j, -x, &fneg); ln_x = log(x); ex_term_1 = (j+1.0)*ln_x; ex_term_2 = lg.val; ex_arg.val = ex_term_1 - ex_term_2; /*(j+1.0)*ln_x - lg.val; */ ex_arg.err = GSL_DBL_EPSILON*(fabs(ex_term_1) + fabs(ex_term_2)) + lg.err; stat_e = gsl_sf_exp_err_e(ex_arg.val, ex_arg.err, &ex); cos_term = cos(j*M_PI); result->val = cos_term * fneg.val + 2.0 * seqn_val * ex.val; result->err = fabs(2.0 * ex.err * seqn_val); result->err += fabs(2.0 * ex.val * seqn_err); result->err += fabs(cos_term) * fneg.err; result->err += 4.0 * GSL_DBL_EPSILON * fabs(result->val); return GSL_ERROR_SELECT_3(stat_e, stat_fneg, stat_lg); }
/* normalization for hydrogenic wave functions */ static int R_norm(const int n, const int l, const double Z, gsl_sf_result * result) { double A = 2.0*Z/n; double pre = sqrt(A*A*A /(2.0*n)); gsl_sf_result ln_a, ln_b; gsl_sf_result ex; int stat_a = gsl_sf_lnfact_e(n+l, &ln_a); int stat_b = gsl_sf_lnfact_e(n-l-1, &ln_b); double diff_val = 0.5*(ln_b.val - ln_a.val); double diff_err = 0.5*(ln_b.err + ln_a.err) + GSL_DBL_EPSILON * fabs(diff_val); int stat_e = gsl_sf_exp_err_e(diff_val, diff_err, &ex); result->val = pre * ex.val; result->err = pre * ex.err; result->err += 2.0 * GSL_DBL_EPSILON * fabs(result->val); return GSL_ERROR_SELECT_3(stat_e, stat_a, stat_b); }
int gsl_sf_poch_e(const double a, const double x, gsl_sf_result * result) { /* CHECK_POINTER(result) */ if(x == 0.0) { result->val = 1.0; result->err = 0.0; return GSL_SUCCESS; } else { gsl_sf_result lnpoch; double sgn; int stat_lnpoch = gsl_sf_lnpoch_sgn_e(a, x, &lnpoch, &sgn); int stat_exp = gsl_sf_exp_err_e(lnpoch.val, lnpoch.err, result); result->val *= sgn; result->err += 2.0 * GSL_DBL_EPSILON * fabs(result->val); return GSL_ERROR_SELECT_2(stat_exp, stat_lnpoch); } }
int gsl_sf_beta_e(const double x, const double y, gsl_sf_result * result) { if((x > 0 && y > 0) && x < 50.0 && y < 50.0) { /* Handle the easy case */ gsl_sf_result gx, gy, gxy; gsl_sf_gamma_e(x, &gx); gsl_sf_gamma_e(y, &gy); gsl_sf_gamma_e(x+y, &gxy); result->val = (gx.val*gy.val)/gxy.val; result->err = gx.err * fabs(gy.val/gxy.val); result->err += gy.err * fabs(gx.val/gxy.val); result->err += fabs((gx.val*gy.val)/(gxy.val*gxy.val)) * gxy.err; result->err += 2.0 * GSL_DBL_EPSILON * fabs(result->val); return GSL_SUCCESS; } else if (isnegint(x) || isnegint(y)) { DOMAIN_ERROR(result); } else if (isnegint(x+y)) { /* infinity in the denominator */ result->val = 0.0; result->err = 0.0; return GSL_SUCCESS; } else { gsl_sf_result lb; double sgn; int stat_lb = gsl_sf_lnbeta_sgn_e(x, y, &lb, &sgn); if(stat_lb == GSL_SUCCESS) { int status = gsl_sf_exp_err_e(lb.val, lb.err, result); result->val *= sgn; return status; } else { result->val = 0.0; result->err = 0.0; return stat_lb; } } }
int gsl_sf_exprel_n_e(const int N, const double x, gsl_sf_result * result) { if(N < 0) { DOMAIN_ERROR(result); } else if(x == 0.0) { result->val = 1.0; result->err = 0.0; return GSL_SUCCESS; } else if(fabs(x) < GSL_ROOT3_DBL_EPSILON * N) { result->val = 1.0 + x/(N+1) * (1.0 + x/(N+2)); result->err = 2.0 * GSL_DBL_EPSILON; return GSL_SUCCESS; } else if(N == 0) { return gsl_sf_exp_e(x, result); } else if(N == 1) { return gsl_sf_exprel_e(x, result); } else if(N == 2) { return gsl_sf_exprel_2_e(x, result); } else { if(x > N && (-x + N*(1.0 + log(x/N)) < GSL_LOG_DBL_EPSILON)) { /* x is much larger than n. * Ignore polynomial part, so * exprel_N(x) ~= e^x N!/x^N */ gsl_sf_result lnf_N; double lnr_val; double lnr_err; double lnterm; gsl_sf_lnfact_e(N, &lnf_N); lnterm = N*log(x); lnr_val = x + lnf_N.val - lnterm; lnr_err = GSL_DBL_EPSILON * (fabs(x) + fabs(lnf_N.val) + fabs(lnterm)); lnr_err += lnf_N.err; return gsl_sf_exp_err_e(lnr_val, lnr_err, result); } else if(x > N) { /* Write the identity * exprel_n(x) = e^x n! / x^n (1 - Gamma[n,x]/Gamma[n]) * then use the asymptotic expansion * Gamma[n,x] ~ x^(n-1) e^(-x) (1 + (n-1)/x + (n-1)(n-2)/x^2 + ...) */ double ln_x = log(x); gsl_sf_result lnf_N; double lg_N; double lnpre_val; double lnpre_err; gsl_sf_lnfact_e(N, &lnf_N); /* log(N!) */ lg_N = lnf_N.val - log(N); /* log(Gamma(N)) */ lnpre_val = x + lnf_N.val - N*ln_x; lnpre_err = GSL_DBL_EPSILON * (fabs(x) + fabs(lnf_N.val) + fabs(N*ln_x)); lnpre_err += lnf_N.err; if(lnpre_val < GSL_LOG_DBL_MAX - 5.0) { int stat_eG; gsl_sf_result bigG_ratio; gsl_sf_result pre; int stat_ex = gsl_sf_exp_err_e(lnpre_val, lnpre_err, &pre); double ln_bigG_ratio_pre = -x + (N-1)*ln_x - lg_N; double bigGsum = 1.0; double term = 1.0; int k; for(k=1; k<N; k++) { term *= (N-k)/x; bigGsum += term; } stat_eG = gsl_sf_exp_mult_e(ln_bigG_ratio_pre, bigGsum, &bigG_ratio); if(stat_eG == GSL_SUCCESS) { result->val = pre.val * (1.0 - bigG_ratio.val); result->err = pre.val * (2.0*GSL_DBL_EPSILON + bigG_ratio.err); result->err += pre.err * fabs(1.0 - bigG_ratio.val); result->err += 2.0 * GSL_DBL_EPSILON * fabs(result->val); return stat_ex; } else { result->val = 0.0; result->err = 0.0; return stat_eG; } } else { OVERFLOW_ERROR(result); } } else if(x > -10.0*N) { return exprel_n_CF(N, x, result); } else { /* x -> -Inf asymptotic: * exprel_n(x) ~ e^x n!/x^n - n/x (1 + (n-1)/x + (n-1)(n-2)/x + ...) * ~ - n/x (1 + (n-1)/x + (n-1)(n-2)/x + ...) */ double sum = 1.0; double term = 1.0; int k; for(k=1; k<N; k++) { term *= (N-k)/x; sum += term; } result->val = -N/x * sum; result->err = 2.0 * GSL_DBL_EPSILON * fabs(result->val); return GSL_SUCCESS; } } }
int gsl_sf_bessel_IJ_taylor_e(const double nu, const double x, const int sign, const int kmax, const double threshold, gsl_sf_result * result ) { /* CHECK_POINTER(result) */ if(nu < 0.0 || x < 0.0) { DOMAIN_ERROR(result); } else if(x == 0.0) { if(nu == 0.0) { result->val = 1.0; result->err = 0.0; } else { result->val = 0.0; result->err = 0.0; } return GSL_SUCCESS; } else { gsl_sf_result prefactor; /* (x/2)^nu / Gamma(nu+1) */ gsl_sf_result sum; int stat_pre; int stat_sum; int stat_mul; if(nu == 0.0) { prefactor.val = 1.0; prefactor.err = 0.0; stat_pre = GSL_SUCCESS; } else if(nu < INT_MAX-1) { /* Separate the integer part and use * y^nu / Gamma(nu+1) = y^N /N! y^f / (N+1)_f, * to control the error. */ const int N = (int)floor(nu + 0.5); const double f = nu - N; gsl_sf_result poch_factor; gsl_sf_result tc_factor; const int stat_poch = gsl_sf_poch_e(N+1.0, f, &poch_factor); const int stat_tc = gsl_sf_taylorcoeff_e(N, 0.5*x, &tc_factor); const double p = pow(0.5*x,f); prefactor.val = tc_factor.val * p / poch_factor.val; prefactor.err = tc_factor.err * p / poch_factor.val; prefactor.err += fabs(prefactor.val) / poch_factor.val * poch_factor.err; prefactor.err += 2.0 * GSL_DBL_EPSILON * fabs(prefactor.val); stat_pre = GSL_ERROR_SELECT_2(stat_tc, stat_poch); } else { gsl_sf_result lg; const int stat_lg = gsl_sf_lngamma_e(nu+1.0, &lg); const double term1 = nu*log(0.5*x); const double term2 = lg.val; const double ln_pre = term1 - term2; const double ln_pre_err = GSL_DBL_EPSILON * (fabs(term1)+fabs(term2)) + lg.err; const int stat_ex = gsl_sf_exp_err_e(ln_pre, ln_pre_err, &prefactor); stat_pre = GSL_ERROR_SELECT_2(stat_ex, stat_lg); } /* Evaluate the sum. * [Abramowitz+Stegun, 9.1.10] * [Abramowitz+Stegun, 9.6.7] */ { const double y = sign * 0.25 * x*x; double sumk = 1.0; double term = 1.0; int k; for(k=1; k<=kmax; k++) { term *= y/((nu+k)*k); sumk += term; if(fabs(term/sumk) < threshold) break; } sum.val = sumk; sum.err = threshold * fabs(sumk); stat_sum = ( k >= kmax ? GSL_EMAXITER : GSL_SUCCESS ); } stat_mul = gsl_sf_multiply_err_e(prefactor.val, prefactor.err, sum.val, sum.err, result); return GSL_ERROR_SELECT_3(stat_mul, stat_pre, stat_sum); } }
int gsl_sf_legendre_sphPlm_e(const int l, int m, const double x, gsl_sf_result * result) { /* CHECK_POINTER(result) */ if(m < 0 || l < m || x < -1.0 || x > 1.0) { DOMAIN_ERROR(result); } else if(m == 0) { gsl_sf_result P; int stat_P = gsl_sf_legendre_Pl_e(l, x, &P); double pre = sqrt((2.0*l + 1.0)/(4.0*M_PI)); result->val = pre * P.val; result->err = pre * P.err; result->err += 2.0 * GSL_DBL_EPSILON * fabs(result->val); return stat_P; } else if(x == 1.0 || x == -1.0) { /* m > 0 here */ result->val = 0.0; result->err = 0.0; return GSL_SUCCESS; } else { /* m > 0 and |x| < 1 here */ /* Starting value for recursion. * Y_m^m(x) = sqrt( (2m+1)/(4pi m) gamma(m+1/2)/gamma(m) ) (-1)^m (1-x^2)^(m/2) / pi^(1/4) */ gsl_sf_result lncirc; gsl_sf_result lnpoch; double lnpre_val; double lnpre_err; gsl_sf_result ex_pre; double sr; const double sgn = ( GSL_IS_ODD(m) ? -1.0 : 1.0); const double y_mmp1_factor = x * sqrt(2.0*m + 3.0); double y_mm, y_mm_err; double y_mmp1; gsl_sf_log_1plusx_e(-x*x, &lncirc); gsl_sf_lnpoch_e(m, 0.5, &lnpoch); /* Gamma(m+1/2)/Gamma(m) */ lnpre_val = -0.25*M_LNPI + 0.5 * (lnpoch.val + m*lncirc.val); lnpre_err = 0.25*M_LNPI*GSL_DBL_EPSILON + 0.5 * (lnpoch.err + fabs(m)*lncirc.err); gsl_sf_exp_err_e(lnpre_val, lnpre_err, &ex_pre); sr = sqrt((2.0+1.0/m)/(4.0*M_PI)); y_mm = sgn * sr * ex_pre.val; y_mmp1 = y_mmp1_factor * y_mm; y_mm_err = 2.0 * GSL_DBL_EPSILON * fabs(y_mm) + sr * ex_pre.err; y_mm_err *= 1.0 + 1.0/(GSL_DBL_EPSILON + fabs(1.0-x)); if(l == m){ result->val = y_mm; result->err = y_mm_err; result->err += 2.0 * GSL_DBL_EPSILON * fabs(y_mm); return GSL_SUCCESS; } else if(l == m + 1) { result->val = y_mmp1; result->err = fabs(y_mmp1_factor) * y_mm_err; result->err += 2.0 * GSL_DBL_EPSILON * fabs(y_mmp1); return GSL_SUCCESS; } else{ double y_ell = 0.0; int ell; /* Compute Y_l^m, l > m+1, upward recursion on l. */ for(ell=m+2; ell <= l; ell++){ const double rat1 = (double)(ell-m)/(double)(ell+m); const double rat2 = (ell-m-1.0)/(ell+m-1.0); const double factor1 = sqrt(rat1*(2*ell+1)*(2*ell-1)); const double factor2 = sqrt(rat1*rat2*(2*ell+1)/(2*ell-3)); y_ell = (x*y_mmp1*factor1 - (ell+m-1)*y_mm*factor2) / (ell-m); y_mm = y_mmp1; y_mmp1 = y_ell; } result->val = y_ell; result->err = (0.5*(l-m) + 1.0) * GSL_DBL_EPSILON * fabs(y_ell); result->err += fabs(y_mm_err/y_mm) * fabs(y_ell); return GSL_SUCCESS; } } }
int gsl_sf_beta_inc_e( const double a, const double b, const double x, gsl_sf_result * result ) { if(a <= 0.0 || b <= 0.0 || x < 0.0 || x > 1.0) { DOMAIN_ERROR(result); } else if(x == 0.0) { result->val = 0.0; result->err = 0.0; return GSL_SUCCESS; } else if(x == 1.0) { result->val = 1.0; result->err = 0.0; return GSL_SUCCESS; } else { gsl_sf_result ln_beta; gsl_sf_result ln_x; gsl_sf_result ln_1mx; gsl_sf_result prefactor; const int stat_ln_beta = gsl_sf_lnbeta_e(a, b, &ln_beta); const int stat_ln_1mx = gsl_sf_log_1plusx_e(-x, &ln_1mx); const int stat_ln_x = gsl_sf_log_e(x, &ln_x); const int stat_ln = GSL_ERROR_SELECT_3(stat_ln_beta, stat_ln_1mx, stat_ln_x); const double ln_pre_val = -ln_beta.val + a * ln_x.val + b * ln_1mx.val; const double ln_pre_err = ln_beta.err + fabs(a*ln_x.err) + fabs(b*ln_1mx.err); const int stat_exp = gsl_sf_exp_err_e(ln_pre_val, ln_pre_err, &prefactor); if(stat_ln != GSL_SUCCESS) { result->val = 0.0; result->err = 0.0; GSL_ERROR ("error", GSL_ESANITY); } if(x < (a + 1.0)/(a+b+2.0)) { /* Apply continued fraction directly. */ gsl_sf_result cf; const int stat_cf = beta_cont_frac(a, b, x, &cf); int stat; result->val = prefactor.val * cf.val / a; result->err = (fabs(prefactor.err * cf.val) + fabs(prefactor.val * cf.err))/a; stat = GSL_ERROR_SELECT_2(stat_exp, stat_cf); if(stat == GSL_SUCCESS) { CHECK_UNDERFLOW(result); } return stat; } else { /* Apply continued fraction after hypergeometric transformation. */ gsl_sf_result cf; const int stat_cf = beta_cont_frac(b, a, 1.0-x, &cf); int stat; const double term = prefactor.val * cf.val / b; result->val = 1.0 - term; result->err = fabs(prefactor.err * cf.val)/b; result->err += fabs(prefactor.val * cf.err)/b; result->err += 2.0 * GSL_DBL_EPSILON * (1.0 + fabs(term)); stat = GSL_ERROR_SELECT_2(stat_exp, stat_cf); if(stat == GSL_SUCCESS) { CHECK_UNDERFLOW(result); } return stat; } } }
int gsl_sf_hyperg_U_large_b_e(const double a, const double b, const double x, gsl_sf_result * result, double * ln_multiplier ) { double N = floor(b); /* b = N + eps */ double eps = b - N; if(fabs(eps) < GSL_SQRT_DBL_EPSILON) { double lnpre_val; double lnpre_err; gsl_sf_result M; if(b > 1.0) { double tmp = (1.0-b)*log(x); gsl_sf_result lg_bm1; gsl_sf_result lg_a; gsl_sf_lngamma_e(b-1.0, &lg_bm1); gsl_sf_lngamma_e(a, &lg_a); lnpre_val = tmp + x + lg_bm1.val - lg_a.val; lnpre_err = lg_bm1.err + lg_a.err + GSL_DBL_EPSILON * (fabs(x) + fabs(tmp)); gsl_sf_hyperg_1F1_large_b_e(1.0-a, 2.0-b, -x, &M); } else { gsl_sf_result lg_1mb; gsl_sf_result lg_1pamb; gsl_sf_lngamma_e(1.0-b, &lg_1mb); gsl_sf_lngamma_e(1.0+a-b, &lg_1pamb); lnpre_val = lg_1mb.val - lg_1pamb.val; lnpre_err = lg_1mb.err + lg_1pamb.err; gsl_sf_hyperg_1F1_large_b_e(a, b, x, &M); } if(lnpre_val > GSL_LOG_DBL_MAX-10.0) { result->val = M.val; result->err = M.err; *ln_multiplier = lnpre_val; GSL_ERROR ("overflow", GSL_EOVRFLW); } else { gsl_sf_result epre; int stat_e = gsl_sf_exp_err_e(lnpre_val, lnpre_err, &epre); result->val = epre.val * M.val; result->err = epre.val * M.err + epre.err * fabs(M.val); result->err += 2.0 * GSL_DBL_EPSILON * fabs(result->val); *ln_multiplier = 0.0; return stat_e; } } else { double omb_lnx = (1.0-b)*log(x); gsl_sf_result lg_1mb; double sgn_1mb; gsl_sf_result lg_1pamb; double sgn_1pamb; gsl_sf_result lg_bm1; double sgn_bm1; gsl_sf_result lg_a; double sgn_a; gsl_sf_result M1, M2; double lnpre1_val, lnpre2_val; double lnpre1_err, lnpre2_err; double sgpre1, sgpre2; gsl_sf_hyperg_1F1_large_b_e( a, b, x, &M1); gsl_sf_hyperg_1F1_large_b_e(1.0-a, 2.0-b, x, &M2); gsl_sf_lngamma_sgn_e(1.0-b, &lg_1mb, &sgn_1mb); gsl_sf_lngamma_sgn_e(1.0+a-b, &lg_1pamb, &sgn_1pamb); gsl_sf_lngamma_sgn_e(b-1.0, &lg_bm1, &sgn_bm1); gsl_sf_lngamma_sgn_e(a, &lg_a, &sgn_a); lnpre1_val = lg_1mb.val - lg_1pamb.val; lnpre1_err = lg_1mb.err + lg_1pamb.err; lnpre2_val = lg_bm1.val - lg_a.val - omb_lnx - x; lnpre2_err = lg_bm1.err + lg_a.err + GSL_DBL_EPSILON * (fabs(omb_lnx)+fabs(x)); sgpre1 = sgn_1mb * sgn_1pamb; sgpre2 = sgn_bm1 * sgn_a; if(lnpre1_val > GSL_LOG_DBL_MAX-10.0 || lnpre2_val > GSL_LOG_DBL_MAX-10.0) { double max_lnpre_val = GSL_MAX(lnpre1_val,lnpre2_val); double max_lnpre_err = GSL_MAX(lnpre1_err,lnpre2_err); double lp1 = lnpre1_val - max_lnpre_val; double lp2 = lnpre2_val - max_lnpre_val; double t1 = sgpre1*exp(lp1); double t2 = sgpre2*exp(lp2); result->val = t1*M1.val + t2*M2.val; result->err = fabs(t1)*M1.err + fabs(t2)*M2.err; result->err += GSL_DBL_EPSILON * exp(max_lnpre_err) * (fabs(t1*M1.val) + fabs(t2*M2.val)); result->err += 2.0 * GSL_DBL_EPSILON * fabs(result->val); *ln_multiplier = max_lnpre_val; GSL_ERROR ("overflow", GSL_EOVRFLW); } else { double t1 = sgpre1*exp(lnpre1_val); double t2 = sgpre2*exp(lnpre2_val); result->val = t1*M1.val + t2*M2.val; result->err = fabs(t1) * M1.err + fabs(t2)*M2.err; result->err += GSL_DBL_EPSILON * (exp(lnpre1_err)*fabs(t1*M1.val) + exp(lnpre2_err)*fabs(t2*M2.val)); result->err += 2.0 * GSL_DBL_EPSILON * fabs(result->val); *ln_multiplier = 0.0; return GSL_SUCCESS; } } }
/** * C++ version of gsl_sf_exp_err_e(). * Exponentiate a quantity with an associated error. * @param x A real number * @param dx A real number * @param result The result as a @c gsl::sf::result object * @return Error code on failure */ inline int exp_err_e( double const x, double const dx, result& result ){ return gsl_sf_exp_err_e( x, dx, &result ); }
int gsl_sf_hyperg_2F1_e(double a, double b, const double c, const double x, gsl_sf_result * result) { const double d = c - a - b; const double rinta = floor(a + 0.5); const double rintb = floor(b + 0.5); const double rintc = floor(c + 0.5); const int a_neg_integer = ( a < 0.0 && fabs(a - rinta) < locEPS ); const int b_neg_integer = ( b < 0.0 && fabs(b - rintb) < locEPS ); const int c_neg_integer = ( c < 0.0 && fabs(c - rintc) < locEPS ); result->val = 0.0; result->err = 0.0; /* Handle x == 1.0 RJM */ if (fabs (x - 1.0) < locEPS && (c - a - b) > 0 && c != 0 && !c_neg_integer) { gsl_sf_result lngamc, lngamcab, lngamca, lngamcb; double lngamc_sgn, lngamca_sgn, lngamcb_sgn; int status; int stat1 = gsl_sf_lngamma_sgn_e (c, &lngamc, &lngamc_sgn); int stat2 = gsl_sf_lngamma_e (c - a - b, &lngamcab); int stat3 = gsl_sf_lngamma_sgn_e (c - a, &lngamca, &lngamca_sgn); int stat4 = gsl_sf_lngamma_sgn_e (c - b, &lngamcb, &lngamcb_sgn); if (stat1 != GSL_SUCCESS || stat2 != GSL_SUCCESS || stat3 != GSL_SUCCESS || stat4 != GSL_SUCCESS) { DOMAIN_ERROR (result); } status = gsl_sf_exp_err_e (lngamc.val + lngamcab.val - lngamca.val - lngamcb.val, lngamc.err + lngamcab.err + lngamca.err + lngamcb.err, result); result->val *= lngamc_sgn / (lngamca_sgn * lngamcb_sgn); return status; } if(x < -1.0 || 1.0 <= x) { DOMAIN_ERROR(result); } if(c_neg_integer) { /* If c is a negative integer, then either a or b must be a negative integer of smaller magnitude than c to ensure cancellation of the series. */ if(! (a_neg_integer && a > c + 0.1) && ! (b_neg_integer && b > c + 0.1)) { DOMAIN_ERROR(result); } } if(fabs(c-b) < locEPS || fabs(c-a) < locEPS) { return pow_omx(x, d, result); /* (1-x)^(c-a-b) */ } if(a >= 0.0 && b >= 0.0 && c >=0.0 && x >= 0.0 && x < 0.995) { /* Series has all positive definite * terms and x is not close to 1. */ return hyperg_2F1_series(a, b, c, x, result); } if(fabs(a) < 10.0 && fabs(b) < 10.0) { /* a and b are not too large, so we attempt * variations on the series summation. */ if(a_neg_integer) { return hyperg_2F1_series(rinta, b, c, x, result); } if(b_neg_integer) { return hyperg_2F1_series(a, rintb, c, x, result); } if(x < -0.25) { return hyperg_2F1_luke(a, b, c, x, result); } else if(x < 0.5) { return hyperg_2F1_series(a, b, c, x, result); } else { if(fabs(c) > 10.0) { return hyperg_2F1_series(a, b, c, x, result); } else { return hyperg_2F1_reflect(a, b, c, x, result); } } } else { /* Either a or b or both large. * Introduce some new variables ap,bp so that bp is * the larger in magnitude. */ double ap, bp; if(fabs(a) > fabs(b)) { bp = a; ap = b; } else { bp = b; ap = a; } if(x < 0.0) { /* What the hell, maybe Luke will converge. */ return hyperg_2F1_luke(a, b, c, x, result); } if(GSL_MAX_DBL(fabs(a),1.0)*fabs(bp)*fabs(x) < 2.0*fabs(c)) { /* If c is large enough or x is small enough, * we can attempt the series anyway. */ return hyperg_2F1_series(a, b, c, x, result); } if(fabs(bp*bp*x*x) < 0.001*fabs(bp) && fabs(a) < 10.0) { /* The famous but nearly worthless "large b" asymptotic. */ int stat = gsl_sf_hyperg_1F1_e(a, c, bp*x, result); result->err = 0.001 * fabs(result->val); return stat; } /* We give up. */ result->val = 0.0; result->err = 0.0; GSL_ERROR ("error", GSL_EUNIMPL); } }
/* Do the reflection described in [Moshier, p. 334]. * Assumes a,b,c != neg integer. */ static int hyperg_2F1_reflect(const double a, const double b, const double c, const double x, gsl_sf_result * result) { const double d = c - a - b; const int intd = floor(d+0.5); const int d_integer = ( fabs(d - intd) < locEPS ); if(d_integer) { const double ln_omx = log(1.0 - x); const double ad = fabs(d); int stat_F2 = GSL_SUCCESS; double sgn_2; gsl_sf_result F1; gsl_sf_result F2; double d1, d2; gsl_sf_result lng_c; gsl_sf_result lng_ad2; gsl_sf_result lng_bd2; int stat_c; int stat_ad2; int stat_bd2; if(d >= 0.0) { d1 = d; d2 = 0.0; } else { d1 = 0.0; d2 = d; } stat_ad2 = gsl_sf_lngamma_e(a+d2, &lng_ad2); stat_bd2 = gsl_sf_lngamma_e(b+d2, &lng_bd2); stat_c = gsl_sf_lngamma_e(c, &lng_c); /* Evaluate F1. */ if(ad < GSL_DBL_EPSILON) { /* d = 0 */ F1.val = 0.0; F1.err = 0.0; } else { gsl_sf_result lng_ad; gsl_sf_result lng_ad1; gsl_sf_result lng_bd1; int stat_ad = gsl_sf_lngamma_e(ad, &lng_ad); int stat_ad1 = gsl_sf_lngamma_e(a+d1, &lng_ad1); int stat_bd1 = gsl_sf_lngamma_e(b+d1, &lng_bd1); if(stat_ad1 == GSL_SUCCESS && stat_bd1 == GSL_SUCCESS && stat_ad == GSL_SUCCESS) { /* Gamma functions in the denominator are ok. * Proceed with evaluation. */ int i; double sum1 = 1.0; double term = 1.0; double ln_pre1_val = lng_ad.val + lng_c.val + d2*ln_omx - lng_ad1.val - lng_bd1.val; double ln_pre1_err = lng_ad.err + lng_c.err + lng_ad1.err + lng_bd1.err + GSL_DBL_EPSILON * fabs(ln_pre1_val); int stat_e; /* Do F1 sum. */ for(i=1; i<ad; i++) { int j = i-1; term *= (a + d2 + j) * (b + d2 + j) / (1.0 + d2 + j) / i * (1.0-x); sum1 += term; } stat_e = gsl_sf_exp_mult_err_e(ln_pre1_val, ln_pre1_err, sum1, GSL_DBL_EPSILON*fabs(sum1), &F1); if(stat_e == GSL_EOVRFLW) { OVERFLOW_ERROR(result); } } else { /* Gamma functions in the denominator were not ok. * So the F1 term is zero. */ F1.val = 0.0; F1.err = 0.0; } } /* end F1 evaluation */ /* Evaluate F2. */ if(stat_ad2 == GSL_SUCCESS && stat_bd2 == GSL_SUCCESS) { /* Gamma functions in the denominator are ok. * Proceed with evaluation. */ const int maxiter = 2000; double psi_1 = -M_EULER; gsl_sf_result psi_1pd; gsl_sf_result psi_apd1; gsl_sf_result psi_bpd1; int stat_1pd = gsl_sf_psi_e(1.0 + ad, &psi_1pd); int stat_apd1 = gsl_sf_psi_e(a + d1, &psi_apd1); int stat_bpd1 = gsl_sf_psi_e(b + d1, &psi_bpd1); int stat_dall = GSL_ERROR_SELECT_3(stat_1pd, stat_apd1, stat_bpd1); double psi_val = psi_1 + psi_1pd.val - psi_apd1.val - psi_bpd1.val - ln_omx; double psi_err = psi_1pd.err + psi_apd1.err + psi_bpd1.err + GSL_DBL_EPSILON*fabs(psi_val); double fact = 1.0; double sum2_val = psi_val; double sum2_err = psi_err; double ln_pre2_val = lng_c.val + d1*ln_omx - lng_ad2.val - lng_bd2.val; double ln_pre2_err = lng_c.err + lng_ad2.err + lng_bd2.err + GSL_DBL_EPSILON * fabs(ln_pre2_val); int stat_e; int j; /* Do F2 sum. */ for(j=1; j<maxiter; j++) { /* values for psi functions use recurrence; Abramowitz+Stegun 6.3.5 */ double term1 = 1.0/(double)j + 1.0/(ad+j); double term2 = 1.0/(a+d1+j-1.0) + 1.0/(b+d1+j-1.0); double delta = 0.0; psi_val += term1 - term2; psi_err += GSL_DBL_EPSILON * (fabs(term1) + fabs(term2)); fact *= (a+d1+j-1.0)*(b+d1+j-1.0)/((ad+j)*j) * (1.0-x); delta = fact * psi_val; sum2_val += delta; sum2_err += fabs(fact * psi_err) + GSL_DBL_EPSILON*fabs(delta); if(fabs(delta) < GSL_DBL_EPSILON * fabs(sum2_val)) break; } if(j == maxiter) stat_F2 = GSL_EMAXITER; if(sum2_val == 0.0) { F2.val = 0.0; F2.err = 0.0; } else { stat_e = gsl_sf_exp_mult_err_e(ln_pre2_val, ln_pre2_err, sum2_val, sum2_err, &F2); if(stat_e == GSL_EOVRFLW) { result->val = 0.0; result->err = 0.0; GSL_ERROR ("error", GSL_EOVRFLW); } } stat_F2 = GSL_ERROR_SELECT_2(stat_F2, stat_dall); } else { /* Gamma functions in the denominator not ok. * So the F2 term is zero. */ F2.val = 0.0; F2.err = 0.0; } /* end F2 evaluation */ sgn_2 = ( GSL_IS_ODD(intd) ? -1.0 : 1.0 ); result->val = F1.val + sgn_2 * F2.val; result->err = F1.err + F2. err; result->err += 2.0 * GSL_DBL_EPSILON * (fabs(F1.val) + fabs(F2.val)); result->err += 2.0 * GSL_DBL_EPSILON * fabs(result->val); return stat_F2; } else { /* d not an integer */ gsl_sf_result pre1, pre2; double sgn1, sgn2; gsl_sf_result F1, F2; int status_F1, status_F2; /* These gamma functions appear in the denominator, so we * catch their harmless domain errors and set the terms to zero. */ gsl_sf_result ln_g1ca, ln_g1cb, ln_g2a, ln_g2b; double sgn_g1ca, sgn_g1cb, sgn_g2a, sgn_g2b; int stat_1ca = gsl_sf_lngamma_sgn_e(c-a, &ln_g1ca, &sgn_g1ca); int stat_1cb = gsl_sf_lngamma_sgn_e(c-b, &ln_g1cb, &sgn_g1cb); int stat_2a = gsl_sf_lngamma_sgn_e(a, &ln_g2a, &sgn_g2a); int stat_2b = gsl_sf_lngamma_sgn_e(b, &ln_g2b, &sgn_g2b); int ok1 = (stat_1ca == GSL_SUCCESS && stat_1cb == GSL_SUCCESS); int ok2 = (stat_2a == GSL_SUCCESS && stat_2b == GSL_SUCCESS); gsl_sf_result ln_gc, ln_gd, ln_gmd; double sgn_gc, sgn_gd, sgn_gmd; gsl_sf_lngamma_sgn_e( c, &ln_gc, &sgn_gc); gsl_sf_lngamma_sgn_e( d, &ln_gd, &sgn_gd); gsl_sf_lngamma_sgn_e(-d, &ln_gmd, &sgn_gmd); sgn1 = sgn_gc * sgn_gd * sgn_g1ca * sgn_g1cb; sgn2 = sgn_gc * sgn_gmd * sgn_g2a * sgn_g2b; if(ok1 && ok2) { double ln_pre1_val = ln_gc.val + ln_gd.val - ln_g1ca.val - ln_g1cb.val; double ln_pre2_val = ln_gc.val + ln_gmd.val - ln_g2a.val - ln_g2b.val + d*log(1.0-x); double ln_pre1_err = ln_gc.err + ln_gd.err + ln_g1ca.err + ln_g1cb.err; double ln_pre2_err = ln_gc.err + ln_gmd.err + ln_g2a.err + ln_g2b.err; if(ln_pre1_val < GSL_LOG_DBL_MAX && ln_pre2_val < GSL_LOG_DBL_MAX) { gsl_sf_exp_err_e(ln_pre1_val, ln_pre1_err, &pre1); gsl_sf_exp_err_e(ln_pre2_val, ln_pre2_err, &pre2); pre1.val *= sgn1; pre2.val *= sgn2; } else { OVERFLOW_ERROR(result); } } else if(ok1 && !ok2) { double ln_pre1_val = ln_gc.val + ln_gd.val - ln_g1ca.val - ln_g1cb.val; double ln_pre1_err = ln_gc.err + ln_gd.err + ln_g1ca.err + ln_g1cb.err; if(ln_pre1_val < GSL_LOG_DBL_MAX) { gsl_sf_exp_err_e(ln_pre1_val, ln_pre1_err, &pre1); pre1.val *= sgn1; pre2.val = 0.0; pre2.err = 0.0; } else { OVERFLOW_ERROR(result); } } else if(!ok1 && ok2) { double ln_pre2_val = ln_gc.val + ln_gmd.val - ln_g2a.val - ln_g2b.val + d*log(1.0-x); double ln_pre2_err = ln_gc.err + ln_gmd.err + ln_g2a.err + ln_g2b.err; if(ln_pre2_val < GSL_LOG_DBL_MAX) { pre1.val = 0.0; pre1.err = 0.0; gsl_sf_exp_err_e(ln_pre2_val, ln_pre2_err, &pre2); pre2.val *= sgn2; } else { OVERFLOW_ERROR(result); } } else { pre1.val = 0.0; pre2.val = 0.0; UNDERFLOW_ERROR(result); } status_F1 = hyperg_2F1_series( a, b, 1.0-d, 1.0-x, &F1); status_F2 = hyperg_2F1_series(c-a, c-b, 1.0+d, 1.0-x, &F2); result->val = pre1.val*F1.val + pre2.val*F2.val; result->err = fabs(pre1.val*F1.err) + fabs(pre2.val*F2.err); result->err += fabs(pre1.err*F1.val) + fabs(pre2.err*F2.val); result->err += 2.0 * GSL_DBL_EPSILON * (fabs(pre1.val*F1.val) + fabs(pre2.val*F2.val)); result->err += 2.0 * GSL_DBL_EPSILON * fabs(result->val); return GSL_SUCCESS; } }