Index: lib/msun/Makefile =================================================================== --- lib/msun/Makefile +++ lib/msun/Makefile @@ -59,7 +59,7 @@ SHLIB_MAJOR= 5 WARNS?= 1 IGNORE_PRAGMA= -COMMON_SRCS= b_exp.c b_log.c b_tgamma.c \ +COMMON_SRCS= b_tgamma.c \ e_acos.c e_acosf.c e_acosh.c e_acoshf.c e_asin.c e_asinf.c \ e_atan2.c e_atan2f.c e_atanh.c e_atanhf.c e_cosh.c e_coshf.c e_exp.c \ e_expf.c e_fmod.c e_fmodf.c e_gamma.c e_gamma_r.c e_gammaf.c \ @@ -69,7 +69,6 @@ e_pow.c e_powf.c e_rem_pio2.c \ e_rem_pio2f.c e_remainder.c e_remainderf.c e_scalb.c e_scalbf.c \ e_sinh.c e_sinhf.c e_sqrt.c e_sqrtf.c fenv.c \ - imprecise.c \ k_cos.c k_cosf.c k_exp.c k_expf.c k_rem_pio2.c k_sin.c k_sinf.c \ k_tan.c k_tanf.c \ s_asinh.c s_asinhf.c s_atan.c s_atanf.c s_carg.c s_cargf.c s_cargl.c \ @@ -116,7 +115,7 @@ COMMON_SRCS+= s_copysignl.c s_fabsl.c s_llrintl.c s_lrintl.c s_modfl.c .if ${LDBL_PREC} != 53 # If long double != double use these; otherwise, we alias the double versions. -COMMON_SRCS+= catrigl.c \ +COMMON_SRCS+= b_tgammal.c catrigl.c \ e_acoshl.c e_acosl.c e_asinl.c e_atan2l.c e_atanhl.c \ e_coshl.c e_fmodl.c e_hypotl.c \ e_lgammal.c e_lgammal_r.c e_powl.c \ Index: lib/msun/bsdsrc/b_exp.c =================================================================== --- lib/msun/bsdsrc/b_exp.c +++ lib/msun/bsdsrc/b_exp.c @@ -33,7 +33,6 @@ #include __FBSDID("$FreeBSD$"); - /* EXP(X) * RETURN THE EXPONENTIAL OF X * DOUBLE PRECISION (IEEE 53 bits, VAX D FORMAT 56 BITS) @@ -41,14 +40,14 @@ * REVISED BY K.C. NG on 2/6/85, 2/15/85, 3/7/85, 3/24/85, 4/16/85, 6/14/86. * * Required system supported functions: - * scalb(x,n) + * ldexp(x,n) * copysign(x,y) - * finite(x) + * isfinite(x) * * Method: * 1. Argument Reduction: given the input x, find r and integer k such * that - * x = k*ln2 + r, |r| <= 0.5*ln2 . + * x = k*ln2 + r, |r| <= 0.5*ln2. * r will be represented as r := z+c for better accuracy. * * 2. Compute exp(r) by @@ -69,105 +68,59 @@ * with 1,156,000 random arguments on a VAX, the maximum observed * error was 0.869 ulps (units in the last place). */ - -#include "mathimpl.h" - -static const double p1 = 0x1.555555555553ep-3; -static const double p2 = -0x1.6c16c16bebd93p-9; -static const double p3 = 0x1.1566aaf25de2cp-14; -static const double p4 = -0x1.bbd41c5d26bf1p-20; -static const double p5 = 0x1.6376972bea4d0p-25; -static const double ln2hi = 0x1.62e42fee00000p-1; -static const double ln2lo = 0x1.a39ef35793c76p-33; -static const double lnhuge = 0x1.6602b15b7ecf2p9; -static const double lntiny = -0x1.77af8ebeae354p9; -static const double invln2 = 0x1.71547652b82fep0; - -#if 0 -double exp(x) -double x; -{ - double z,hi,lo,c; - int k; - -#if !defined(vax)&&!defined(tahoe) - if(x!=x) return(x); /* x is NaN */ -#endif /* !defined(vax)&&!defined(tahoe) */ - if( x <= lnhuge ) { - if( x >= lntiny ) { - - /* argument reduction : x --> x - k*ln2 */ - - k=invln2*x+copysign(0.5,x); /* k=NINT(x/ln2) */ - - /* express x-k*ln2 as hi-lo and let x=hi-lo rounded */ - - hi=x-k*ln2hi; - x=hi-(lo=k*ln2lo); - - /* return 2^k*[1+x+x*c/(2+c)] */ - z=x*x; - c= x - z*(p1+z*(p2+z*(p3+z*(p4+z*p5)))); - return scalb(1.0+(hi-(lo-(x*c)/(2.0-c))),k); - - } - /* end of x > lntiny */ - - else - /* exp(-big#) underflows to zero */ - if(finite(x)) return(scalb(1.0,-5000)); - - /* exp(-INF) is zero */ - else return(0.0); - } - /* end of x < lnhuge */ - - else - /* exp(INF) is INF, exp(+big#) overflows to INF */ - return( finite(x) ? scalb(1.0,5000) : x); -} -#endif +static const double + p1 = 1.6666666666666660e-01, /* 0x3fc55555, 0x55555553 */ + p2 = -2.7777777777564776e-03, /* 0xbf66c16c, 0x16c0ac3c */ + p3 = 6.6137564717940088e-05, /* 0x3f11566a, 0xb5c2ba0d */ + p4 = -1.6534060280704225e-06, /* 0xbebbbd53, 0x273e8fb7 */ + p5 = 4.1437773411069054e-08; /* 0x3e663f2a, 0x09c94b6c */ + +static const double + ln2hi = 0x1.62e42fee00000p-1, /* High 32 bits round-down. */ + ln2lo = 0x1.a39ef35793c76p-33; /* Next 53 bits round-to-nearst. */ + +static const double + lnhuge = 0x1.6602b15b7ecf2p9, /* (DBL_MAX_EXP + 9) * log(2.) */ + lntiny = -0x1.77af8ebeae354p9, /* (DBL_MIN_EXP - 53 - 10) * log(2.) */ + invln2 = 0x1.71547652b82fep0; /* 1 / log(2.) */ /* returns exp(r = x + c) for |c| < |x| with no overlap. */ -double __exp__D(x, c) -double x, c; +static double +__exp__D(double x, double c) { - double z,hi,lo; + double hi, lo, z; int k; - if (x != x) /* x is NaN */ + if (x != x) /* x is NaN. */ return(x); - if ( x <= lnhuge ) { - if ( x >= lntiny ) { - /* argument reduction : x --> x - k*ln2 */ - z = invln2*x; - k = z + copysign(.5, x); - - /* express (x+c)-k*ln2 as hi-lo and let x=hi-lo rounded */ - - hi=(x-k*ln2hi); /* Exact. */ - x= hi - (lo = k*ln2lo-c); - /* return 2^k*[1+x+x*c/(2+c)] */ - z=x*x; - c= x - z*(p1+z*(p2+z*(p3+z*(p4+z*p5)))); - c = (x*c)/(2.0-c); - - return scalb(1.+(hi-(lo - c)), k); + if (x <= lnhuge) { + if (x >= lntiny) { + /* argument reduction: x --> x - k*ln2 */ + z = invln2 * x; + k = z + copysign(0.5, x); + + /* + * Express (x + c) - k * ln2 as hi - lo. + * Let x = hi - lo rounded. + */ + hi = x - k * ln2hi; /* Exact. */ + lo = k * ln2lo - c; + x = hi - lo; + + /* Return 2^k*[1+x+x*c/(2+c)] */ + z = x * x; + c = x - z * (p1 + z * (p2 + z * (p3 + z * (p4 + + z * p5)))); + c = (x * c) / (2 - c); + + return (ldexp(1 + (hi - (lo - c)), k)); + } else { + /* exp(-INF) is 0. exp(-big) underflows to 0. */ + return (isfinite(x) ? ldexp(1., -5000) : 0); } - /* end of x > lntiny */ - - else - /* exp(-big#) underflows to zero */ - if(finite(x)) return(scalb(1.0,-5000)); - - /* exp(-INF) is zero */ - else return(0.0); - } - /* end of x < lnhuge */ - - else + } else /* exp(INF) is INF, exp(+big#) overflows to INF */ - return( finite(x) ? scalb(1.0,5000) : x); + return (isfinite(x) ? ldexp(1., 5000) : x); } Index: lib/msun/bsdsrc/b_log.c =================================================================== --- lib/msun/bsdsrc/b_log.c +++ lib/msun/bsdsrc/b_log.c @@ -33,10 +33,6 @@ #include __FBSDID("$FreeBSD$"); -#include - -#include "mathimpl.h" - /* Table-driven natural logarithm. * * This code was derived, with minor modifications, from: @@ -44,25 +40,27 @@ * Logarithm in IEEE Floating-Point arithmetic." ACM Trans. * Math Software, vol 16. no 4, pp 378-400, Dec 1990). * - * Calculates log(2^m*F*(1+f/F)), |f/j| <= 1/256, + * Calculates log(2^m*F*(1+f/F)), |f/F| <= 1/256, * where F = j/128 for j an integer in [0, 128]. * * log(2^m) = log2_hi*m + log2_tail*m - * since m is an integer, the dominant term is exact. + * The leading term is exact, because m is an integer, * m has at most 10 digits (for subnormal numbers), * and log2_hi has 11 trailing zero bits. * - * log(F) = logF_hi[j] + logF_lo[j] is in tabular form in log_table.h + * log(F) = logF_hi[j] + logF_lo[j] is in table below. * logF_hi[] + 512 is exact. * * log(1+f/F) = 2*f/(2*F + f) + 1/12 * (2*f/(2*F + f))**3 + ... - * the leading term is calculated to extra precision in two + * + * The leading term is calculated to extra precision in two * parts, the larger of which adds exactly to the dominant * m and F terms. + * * There are two cases: - * 1. when m, j are non-zero (m | j), use absolute + * 1. When m and j are non-zero (m | j), use absolute * precision for the leading term. - * 2. when m = j = 0, |1-x| < 1/256, and log(x) ~= (x-1). + * 2. When m = j = 0, |1-x| < 1/256, and log(x) ~= (x-1). * In this case, use a relative precision of 24 bits. * (This is done differently in the original paper) * @@ -70,11 +68,21 @@ * 0 return signalling -Inf * neg return signalling NaN * +Inf return +Inf -*/ + */ #define N 128 -/* Table of log(Fj) = logF_head[j] + logF_tail[j], for Fj = 1+j/128. +/* + * Coefficients in the polynomial approximation of log(1+f/F). + * Domain of x is [0,1./256] with 2**(-64.187) precision. + */ +static const double + A1 = 8.3333333333333329e-02, /* 0x3fb55555, 0x55555555 */ + A2 = 1.2499999999943598e-02, /* 0x3f899999, 0x99991a98 */ + A3 = 2.2321527525957776e-03; /* 0x3f624929, 0xe24e70be */ + +/* + * Table of log(Fj) = logF_head[j] + logF_tail[j], for Fj = 1+j/128. * Used for generation of extend precision logarithms. * The constant 35184372088832 is 2^45, so the divide is exact. * It ensures correct reading of logF_head, even for inaccurate @@ -82,12 +90,7 @@ * right answer for integers less than 2^53.) * Values for log(F) were generated using error < 10^-57 absolute * with the bc -l package. -*/ -static double A1 = .08333333333333178827; -static double A2 = .01250000000377174923; -static double A3 = .002232139987919447809; -static double A4 = .0004348877777076145742; - + */ static double logF_head[N+1] = { 0., .007782140442060381246, @@ -351,118 +354,51 @@ .00000000000025144230728376072, -.00000000000017239444525614834 }; - -#if 0 -double -#ifdef _ANSI_SOURCE -log(double x) -#else -log(x) double x; -#endif -{ - int m, j; - double F, f, g, q, u, u2, v, zero = 0.0, one = 1.0; - volatile double u1; - - /* Catch special cases */ - if (x <= 0) - if (x == zero) /* log(0) = -Inf */ - return (-one/zero); - else /* log(neg) = NaN */ - return (zero/zero); - else if (!finite(x)) - return (x+x); /* x = NaN, Inf */ - - /* Argument reduction: 1 <= g < 2; x/2^m = g; */ - /* y = F*(1 + f/F) for |f| <= 2^-8 */ - - m = logb(x); - g = ldexp(x, -m); - if (m == -1022) { - j = logb(g), m += j; - g = ldexp(g, -j); - } - j = N*(g-1) + .5; - F = (1.0/N) * j + 1; /* F*128 is an integer in [128, 512] */ - f = g - F; - - /* Approximate expansion for log(1+f/F) ~= u + q */ - g = 1/(2*F+f); - u = 2*f*g; - v = u*u; - q = u*v*(A1 + v*(A2 + v*(A3 + v*A4))); - - /* case 1: u1 = u rounded to 2^-43 absolute. Since u < 2^-8, - * u1 has at most 35 bits, and F*u1 is exact, as F has < 8 bits. - * It also adds exactly to |m*log2_hi + log_F_head[j] | < 750 - */ - if (m | j) - u1 = u + 513, u1 -= 513; - - /* case 2: |1-x| < 1/256. The m- and j- dependent terms are zero; - * u1 = u to 24 bits. - */ - else - u1 = u, TRUNC(u1); - u2 = (2.0*(f - F*u1) - u1*f) * g; - /* u1 + u2 = 2f/(2F+f) to extra precision. */ - - /* log(x) = log(2^m*F*(1+f/F)) = */ - /* (m*log2_hi+logF_head[j]+u1) + (m*log2_lo+logF_tail[j]+q); */ - /* (exact) + (tiny) */ - - u1 += m*logF_head[N] + logF_head[j]; /* exact */ - u2 = (u2 + logF_tail[j]) + q; /* tiny */ - u2 += logF_tail[N]*m; - return (u1 + u2); -} -#endif - /* * Extra precision variant, returning struct {double a, b;}; - * log(x) = a+b to 63 bits, with a rounded to 26 bits. + * log(x) = a+b to 63 bits, with 'a' rounded to 24 bits. */ -struct Double -#ifdef _ANSI_SOURCE +static struct Double __log__D(double x) -#else -__log__D(x) double x; -#endif { int m, j; - double F, f, g, q, u, v, u2; - volatile double u1; + double F, f, g, q, u, v, u1, u2; struct Double r; - /* Argument reduction: 1 <= g < 2; x/2^m = g; */ - /* y = F*(1 + f/F) for |f| <= 2^-8 */ - - m = logb(x); - g = ldexp(x, -m); + /* + * Argument reduction: 1 <= g < 2; x/2^m = g; + * y = F*(1 + f/F) for |f| <= 2^-8 + */ + g = frexp(x, &m); + g *= 2; + m--; if (m == -1022) { - j = logb(g), m += j; + j = ilogb(g); + m += j; g = ldexp(g, -j); } - j = N*(g-1) + .5; - F = (1.0/N) * j + 1; + j = N * (g - 1) + 0.5; + F = (1. / N) * j + 1; f = g - F; - g = 1/(2*F+f); - u = 2*f*g; - v = u*u; - q = u*v*(A1 + v*(A2 + v*(A3 + v*A4))); - if (m | j) - u1 = u + 513, u1 -= 513; - else - u1 = u, TRUNC(u1); - u2 = (2.0*(f - F*u1) - u1*f) * g; + g = 1 / (2 * F + f); + u = 2 * f * g; + v = u * u; + q = u * v * (A1 + v * (A2 + v * A3)); + if (m | j) { + u1 = u + 513; + u1 -= 513; + } else { + u1 = (float)u; + } + u2 = (2 * (f - F * u1) - u1 * f) * g; - u1 += m*logF_head[N] + logF_head[j]; + u1 += m * logF_head[N] + logF_head[j]; - u2 += logF_tail[j]; u2 += q; - u2 += logF_tail[N]*m; - r.a = u1 + u2; /* Only difference is here */ - TRUNC(r.a); + u2 += logF_tail[j]; + u2 += q; + u2 += logF_tail[N] * m; + r.a = (float)(u1 + u2); /* Only difference is here. */ r.b = (u1 - r.a) + u2; return (r); } Index: lib/msun/bsdsrc/b_tgamma.c =================================================================== --- lib/msun/bsdsrc/b_tgamma.c +++ lib/msun/bsdsrc/b_tgamma.c @@ -29,37 +29,46 @@ * SUCH DAMAGE. */ +/* + * The original code, FreeBSD's old svn r93211, contained the following + * attribution: + * + * This code by P. McIlroy, Oct 1992; + * + * The financial support of UUNET Communications Services is greatfully + * acknowledged. + * + * The algorithm remains, but the code has been re-arranged to facilitate + * porting to other precisions. + */ + /* @(#)gamma.c 8.1 (Berkeley) 6/4/93 */ #include __FBSDID("$FreeBSD$"); -/* - * This code by P. McIlroy, Oct 1992; - * - * The financial support of UUNET Communications Services is greatfully - * acknowledged. - */ +#include -#include -#include "mathimpl.h" +#include "math.h" +#include "math_private.h" -/* METHOD: - * x < 0: Use reflection formula, G(x) = pi/(sin(pi*x)*x*G(x)) - * At negative integers, return NaN and raise invalid. - * - * x < 6.5: - * Use argument reduction G(x+1) = xG(x) to reach the - * range [1.066124,2.066124]. Use a rational - * approximation centered at the minimum (x0+1) to - * ensure monotonicity. - * - * x >= 6.5: Use the asymptotic approximation (Stirling's formula) - * adjusted for equal-ripples: - * - * log(G(x)) ~= (x-.5)*(log(x)-1) + .5(log(2*pi)-1) + 1/x*P(1/(x*x)) +/* Used in b_log.c and below. */ +struct Double { + double a; + double b; +}; + +#include "b_log.c" +#include "b_exp.c" + +/* + * The range is broken into several subranges. Each is handled by its + * helper functions. * - * Keep extra precision in multiplying (x-.5)(log(x)-1), to - * avoid premature round-off. + * x >= 6.0: large_gam(x) + * 6.0 > x >= xleft: small_gam(x) where xleft = 1 + left + x0. + * xleft > x > iota: smaller_gam(x) where iota = 1e-17. + * iota > x > -itoa: Handle x near 0. + * -iota > x : neg_gam * * Special values: * -Inf: return NaN and raise invalid; @@ -77,201 +86,224 @@ * Maximum observed error < 4ulp in 1,000,000 trials. */ -static double neg_gam(double); -static double small_gam(double); -static double smaller_gam(double); -static struct Double large_gam(double); -static struct Double ratfun_gam(double, double); - -/* - * Rational approximation, A0 + x*x*P(x)/Q(x), on the interval - * [1.066.., 2.066..] accurate to 4.25e-19. - */ -#define LEFT -.3955078125 /* left boundary for rat. approx */ -#define x0 .461632144968362356785 /* xmin - 1 */ - -#define a0_hi 0.88560319441088874992 -#define a0_lo -.00000000000000004996427036469019695 -#define P0 6.21389571821820863029017800727e-01 -#define P1 2.65757198651533466104979197553e-01 -#define P2 5.53859446429917461063308081748e-03 -#define P3 1.38456698304096573887145282811e-03 -#define P4 2.40659950032711365819348969808e-03 -#define Q0 1.45019531250000000000000000000e+00 -#define Q1 1.06258521948016171343454061571e+00 -#define Q2 -2.07474561943859936441469926649e-01 -#define Q3 -1.46734131782005422506287573015e-01 -#define Q4 3.07878176156175520361557573779e-02 -#define Q5 5.12449347980666221336054633184e-03 -#define Q6 -1.76012741431666995019222898833e-03 -#define Q7 9.35021023573788935372153030556e-05 -#define Q8 6.13275507472443958924745652239e-06 /* * Constants for large x approximation (x in [6, Inf]) * (Accurate to 2.8*10^-19 absolute) */ -#define lns2pi_hi 0.418945312500000 -#define lns2pi_lo -.000006779295327258219670263595 -#define Pa0 8.33333333333333148296162562474e-02 -#define Pa1 -2.77777777774548123579378966497e-03 -#define Pa2 7.93650778754435631476282786423e-04 -#define Pa3 -5.95235082566672847950717262222e-04 -#define Pa4 8.41428560346653702135821806252e-04 -#define Pa5 -1.89773526463879200348872089421e-03 -#define Pa6 5.69394463439411649408050664078e-03 -#define Pa7 -1.44705562421428915453880392761e-02 - -static const double zero = 0., one = 1.0, tiny = 1e-300; -double -tgamma(x) - double x; +static const double zero = 0.; +static const volatile double tiny = 1e-300; +/* + * x >= 6 + * + * Use the asymptotic approximation (Stirling's formula) adjusted fof + * equal-ripples: + * + * log(G(x)) ~= (x-0.5)*(log(x)-1) + 0.5(log(2*pi)-1) + 1/x*P(1/(x*x)) + * + * Keep extra precision in multiplying (x-.5)(log(x)-1), to avoid + * premature round-off. + * + * Accurate to max(ulp(1/128) absolute, 2^-66 relative) error. + */ +static const double + ln2pi_hi = 0.41894531250000000, + ln2pi_lo = -6.7792953272582197e-6, + Pa0 = 8.3333333333333329e-02, /* 0x3fb55555, 0x55555555 */ + Pa1 = -2.7777777777735404e-03, /* 0xbf66c16c, 0x16c145ec */ + Pa2 = 7.9365079044114095e-04, /* 0x3f4a01a0, 0x183de82d */ + Pa3 = -5.9523715464225254e-04, /* 0xbf438136, 0x0e681f62 */ + Pa4 = 8.4161391899445698e-04, /* 0x3f4b93f8, 0x21042a13 */ + Pa5 = -1.9065246069191080e-03, /* 0xbf5f3c8b, 0x357cb64e */ + Pa6 = 5.9047708485785158e-03, /* 0x3f782f99, 0xdaf5d65f */ + Pa7 = -1.6484018705183290e-02; /* 0xbf90e12f, 0xc4fb4df0 */ + +static struct Double +large_gam(double x) { + double p, z, thi, tlo, xhi, xlo; struct Double u; - if (x >= 6) { - if(x > 171.63) - return (x / zero); - u = large_gam(x); - return(__exp__D(u.a, u.b)); - } else if (x >= 1.0 + LEFT + x0) - return (small_gam(x)); - else if (x > 1.e-17) - return (smaller_gam(x)); - else if (x > -1.e-17) { - if (x != 0.0) - u.a = one - tiny; /* raise inexact */ - return (one/x); - } else if (!finite(x)) - return (x - x); /* x is NaN or -Inf */ - else - return (neg_gam(x)); + z = 1 / (x * x); + p = Pa0 + z * (Pa1 + z * (Pa2 + z * (Pa3 + z * (Pa4 + z * (Pa5 + + z * (Pa6 + z * Pa7)))))); + p = p / x; + + u = __log__D(x); + u.a -= 1; + + /* Split (x - 0.5) in high and low parts. */ + x -= 0.5; + xhi = (float)x; + xlo = x - xhi; + + /* Compute t = (x-.5)*(log(x)-1) in extra precision. */ + thi = xhi * u.a; + tlo = xlo * u.a + x * u.b; + + /* Compute thi + tlo + ln2pi_hi + ln2pi_lo + p. */ + tlo += ln2pi_lo; + tlo += p; + u.a = ln2pi_hi + tlo; + u.a += thi; + u.b = thi - u.a; + u.b += ln2pi_hi; + u.b += tlo; + return (u); } /* - * Accurate to max(ulp(1/128) absolute, 2^-66 relative) error. + * Rational approximation, A0 + x * x * P(x) / Q(x), on the interval + * [1.066.., 2.066..] accurate to 4.25e-19. + * + * Returns r.a + r.b = a0 + (z + c)^2 * p / q, with r.a truncated. */ +static const double +#if 0 + a0_hi = 8.8560319441088875e-1, + a0_lo = -4.9964270364690197e-17, +#else + a0_hi = 8.8560319441088875e-01, /* 0x3fec56dc, 0x82a74aef */ + a0_lo = -4.9642368725563397e-17, /* 0xbc8c9deb, 0xaa64afc3 */ +#endif + P0 = 6.2138957182182086e-1, + P1 = 2.6575719865153347e-1, + P2 = 5.5385944642991746e-3, + P3 = 1.3845669830409657e-3, + P4 = 2.4065995003271137e-3, + Q0 = 1.4501953125000000e+0, + Q1 = 1.0625852194801617e+0, + Q2 = -2.0747456194385994e-1, + Q3 = -1.4673413178200542e-1, + Q4 = 3.0787817615617552e-2, + Q5 = 5.1244934798066622e-3, + Q6 = -1.7601274143166700e-3, + Q7 = 9.3502102357378894e-5, + Q8 = 6.1327550747244396e-6; + static struct Double -large_gam(x) - double x; +ratfun_gam(double z, double c) { - double z, p; - struct Double t, u, v; + double p, q, thi, tlo; + struct Double r; - z = one/(x*x); - p = Pa0+z*(Pa1+z*(Pa2+z*(Pa3+z*(Pa4+z*(Pa5+z*(Pa6+z*Pa7)))))); - p = p/x; + q = Q0 + z * (Q1 + z * (Q2 + z * (Q3 + z * (Q4 + z * (Q5 + + z * (Q6 + z * (Q7 + z * Q8))))))); + p = P0 + z * (P1 + z * (P2 + z * (P3 + z * P4))); + p = p / q; - u = __log__D(x); - u.a -= one; - v.a = (x -= .5); - TRUNC(v.a); - v.b = x - v.a; - t.a = v.a*u.a; /* t = (x-.5)*(log(x)-1) */ - t.b = v.b*u.a + x*u.b; - /* return t.a + t.b + lns2pi_hi + lns2pi_lo + p */ - t.b += lns2pi_lo; t.b += p; - u.a = lns2pi_hi + t.b; u.a += t.a; - u.b = t.a - u.a; - u.b += lns2pi_hi; u.b += t.b; - return (u); + /* Split z into high and low parts. */ + thi = (float)z; + tlo = (z - thi) + c; + tlo *= (thi + z); + + /* Split (z+c)^2 into high and low parts. */ + thi *= thi; + q = thi; + thi = (float)thi; + tlo += (q - thi); + + /* Split p/q into high and low parts. */ + r.a = (float)p; + r.b = p - r.a; + + tlo = tlo * p + thi * r.b + a0_lo; + thi *= r.a; /* t = (z+c)^2*(P/Q) */ + r.a = (float)(thi + a0_hi); + r.b = ((a0_hi - r.a) + thi) + tlo; + return (r); /* r = a0 + t */ } /* + * x < 6 + * + * Use argument reduction G(x+1) = xG(x) to reach the range [1.066124, + * 2.066124]. Use a rational approximation centered at the minimum + * (x0+1) to ensure monotonicity. + * * Good to < 1 ulp. (provably .90 ulp; .87 ulp on 1,000,000 runs.) * It also has correct monotonicity. */ +static const double + left = -0.3955078125, /* left boundary for rat. approx */ + x0 = 4.6163214496836236e-1; /* xmin - 1 */ + static double -small_gam(x) - double x; +small_gam(double x) { - double y, ym1, t; + double t, y, ym1; struct Double yy, r; - y = x - one; - ym1 = y - one; - if (y <= 1.0 + (LEFT + x0)) { + + y = x - 1; + if (y <= 1 + (left + x0)) { yy = ratfun_gam(y - x0, 0); return (yy.a + yy.b); } - r.a = y; - TRUNC(r.a); - yy.a = r.a - one; - y = ym1; - yy.b = r.b = y - yy.a; + + r.a = (float)y; + yy.a = r.a - 1; + y = y - 1 ; + r.b = yy.b = y - yy.a; + /* Argument reduction: G(x+1) = x*G(x) */ - for (ym1 = y-one; ym1 > LEFT + x0; y = ym1--, yy.a--) { - t = r.a*yy.a; - r.b = r.a*yy.b + y*r.b; - r.a = t; - TRUNC(r.a); + for (ym1 = y - 1; ym1 > left + x0; y = ym1--, yy.a--) { + t = r.a * yy.a; + r.b = r.a * yy.b + y * r.b; + r.a = (float)t; r.b += (t - r.a); } + /* Return r*tgamma(y). */ yy = ratfun_gam(y - x0, 0); - y = r.b*(yy.a + yy.b) + r.a*yy.b; - y += yy.a*r.a; + y = r.b * (yy.a + yy.b) + r.a * yy.b; + y += yy.a * r.a; return (y); } /* - * Good on (0, 1+x0+LEFT]. Accurate to 1ulp. + * Good on (0, 1+x0+left]. Accurate to 1 ulp. */ static double -smaller_gam(x) - double x; +smaller_gam(double x) { - double t, d; - struct Double r, xx; - if (x < x0 + LEFT) { - t = x, TRUNC(t); - d = (t+x)*(x-t); + double d, rhi, rlo, t, xhi, xlo; + struct Double r; + + if (x < x0 + left) { + t = (float)x; + d = (t + x) * (x - t); t *= t; - xx.a = (t + x), TRUNC(xx.a); - xx.b = x - xx.a; xx.b += t; xx.b += d; - t = (one-x0); t += x; - d = (one-x0); d -= t; d += x; - x = xx.a + xx.b; + xhi = (float)(t + x); + xlo = x - xhi; + xlo += t; + xlo += d; + t = 1 - x0; + t += x; + d = 1 - x0; + d -= t; + d += x; + x = xhi + xlo; } else { - xx.a = x, TRUNC(xx.a); - xx.b = x - xx.a; + xhi = (float)x; + xlo = x - xhi; t = x - x0; - d = (-x0 -t); d += x; + d = - x0 - t; + d += x; } + r = ratfun_gam(t, d); - d = r.a/x, TRUNC(d); - r.a -= d*xx.a; r.a -= d*xx.b; r.a += r.b; - return (d + r.a/x); + d = (float)(r.a / x); + r.a -= d * xhi; + r.a -= d * xlo; + r.a += r.b; + + return (d + r.a / x); } /* - * returns (z+c)^2 * P(z)/Q(z) + a0 + * x < 0 + * + * Use reflection formula, G(x) = pi/(sin(pi*x)*x*G(x)). + * At negative integers, return NaN and raise invalid. */ -static struct Double -ratfun_gam(z, c) - double z, c; -{ - double p, q; - struct Double r, t; - - q = Q0 +z*(Q1+z*(Q2+z*(Q3+z*(Q4+z*(Q5+z*(Q6+z*(Q7+z*Q8))))))); - p = P0 + z*(P1 + z*(P2 + z*(P3 + z*P4))); - - /* return r.a + r.b = a0 + (z+c)^2*p/q, with r.a truncated to 26 bits. */ - p = p/q; - t.a = z, TRUNC(t.a); /* t ~= z + c */ - t.b = (z - t.a) + c; - t.b *= (t.a + z); - q = (t.a *= t.a); /* t = (z+c)^2 */ - TRUNC(t.a); - t.b += (q - t.a); - r.a = p, TRUNC(r.a); /* r = P/Q */ - r.b = p - r.a; - t.b = t.b*p + t.a*r.b + a0_lo; - t.a *= r.a; /* t = (z+c)^2*(P/Q) */ - r.a = t.a + a0_hi, TRUNC(r.a); - r.b = ((a0_hi-r.a) + t.a) + t.b; - return (r); /* r = a0 + t */ -} - static double -neg_gam(x) - double x; +neg_gam(double x) { int sgn = 1; struct Double lg, lsine; @@ -280,23 +312,29 @@ y = ceil(x); if (y == x) /* Negative integer. */ return ((x - x) / zero); + z = y - x; if (z > 0.5) - z = one - z; - y = 0.5 * y; + z = 1 - z; + + y = y / 2; if (y == ceil(y)) sgn = -1; - if (z < .25) - z = sin(M_PI*z); + + if (z < 0.25) + z = sinpi(z); else - z = cos(M_PI*(0.5-z)); + z = cospi(0.5 - z); + /* Special case: G(1-x) = Inf; G(x) may be nonzero. */ if (x < -170) { + if (x < -190) - return ((double)sgn*tiny*tiny); - y = one - x; /* exact: 128 < |x| < 255 */ + return (sgn * tiny * tiny); + + y = 1 - x; /* exact: 128 < |x| < 255 */ lg = large_gam(y); - lsine = __log__D(M_PI/z); /* = TRUNC(log(u)) + small */ + lsine = __log__D(M_PI / z); /* = TRUNC(log(u)) + small */ lg.a -= lsine.a; /* exact (opposite signs) */ lg.b -= lsine.b; y = -(lg.a + lg.b); @@ -305,11 +343,58 @@ if (sgn < 0) y = -y; return (y); } - y = one-x; - if (one-y == x) + + y = 1 - x; + if (1 - y == x) y = tgamma(y); else /* 1-x is inexact */ - y = -x*tgamma(-x); + y = - x * tgamma(-x); + if (sgn < 0) y = -y; - return (M_PI / (y*z)); + return (M_PI / (y * z)); +} +/* + * xmax comes from lgamma(xmax) - emax * log(2) = 0. + * static const float xmax = 35.040095f + * static const double xmax = 171.624376956302725; + * ld80: LD80C(0xdb718c066b352e20, 10, 1.75554834290446291689e+03L), + * ld128: 1.75554834290446291700388921607020320e+03L, + * + * iota is a sloppy threshold to isolate x = 0. + */ +static const double xmax = 171.624376956302725; +static const double iota = 0x1p-56; + +double +tgamma(double x) +{ + struct Double u; + + if (x >= 6) { + if (x > xmax) + return (x / zero); + u = large_gam(x); + return (__exp__D(u.a, u.b)); + } + + if (x >= 1 + left + x0) + return (small_gam(x)); + + if (x > iota) + return (smaller_gam(x)); + + if (x > -iota) { + if (x != 0.) + u.a = 1 - tiny; /* raise inexact */ + return (1 / x); + } + + if (!isfinite(x)) + return (x - x); /* x is NaN or -Inf */ + + return (neg_gam(x)); } + +#if (LDBL_MANT_DIG == 53) +__weak_reference(tgamma, tgammal); +#endif Index: lib/msun/ld80/b_expl.c =================================================================== --- /dev/null +++ lib/msun/ld80/b_expl.c @@ -0,0 +1,113 @@ +/*- + * SPDX-License-Identifier: BSD-3-Clause + * + * Copyright (c) 1985, 1993 + * The Regents of the University of California. All rights reserved. + * + * Redistribution and use in source and binary forms, with or without + * modification, are permitted provided that the following conditions + * are met: + * 1. Redistributions of source code must retain the above copyright + * notice, this list of conditions and the following disclaimer. + * 2. Redistributions in binary form must reproduce the above copyright + * notice, this list of conditions and the following disclaimer in the + * documentation and/or other materials provided with the distribution. + * 3. Neither the name of the University nor the names of its contributors + * may be used to endorse or promote products derived from this software + * without specific prior written permission. + * + * THIS SOFTWARE IS PROVIDED BY THE REGENTS AND CONTRIBUTORS ``AS IS'' AND + * ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE + * IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE + * ARE DISCLAIMED. IN NO EVENT SHALL THE REGENTS OR CONTRIBUTORS BE LIABLE + * FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL + * DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS + * OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) + * HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT + * LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY + * OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF + * SUCH DAMAGE. + */ + +/* + * See bsdsrc/b_exp.c for implementation details. + * + * bsdrc/b_exp.c converted to long double by Steven G. Kargl. + */ + +#include "fpmath.h" +#include "math_private.h" + +static const union IEEEl2bits + p0u = LD80C(0xaaaaaaaaaaaaaaab, -3, 1.66666666666666666671e-01L), + p1u = LD80C(0xb60b60b60b60b59a, -9, -2.77777777777777775377e-03L), + p2u = LD80C(0x8ab355e008a3cfce, -14, 6.61375661375629297465e-05L), + p3u = LD80C(0xddebbc994b0c1376, -20, -1.65343915327882529784e-06L), + p4u = LD80C(0xb354784cb4ef4c41, -25, 4.17535101591534118469e-08L), + p5u = LD80C(0x913e8a718382ce75, -30, -1.05679137034774806475e-09L), + p6u = LD80C(0xe8f0042aa134502e, -36, 2.64819349895429516863e-11L); +#define p1 (p0u.e) +#define p2 (p1u.e) +#define p3 (p2u.e) +#define p4 (p3u.e) +#define p5 (p4u.e) +#define p6 (p5u.e) +#define p7 (p6u.e) + +/* + * lnhuge = (LDBL_MAX_EXP + 9) * log(2.) + * lntiny = (LDBL_MIN_EXP - 64 - 10) * log(2.) + * invln2 = 1 / log(2.) + */ +static const union IEEEl2bits +ln2hiu = LD80C(0xb17217f700000000, -1, 6.93147180369123816490e-01L), +ln2lou = LD80C(0xd1cf79abc9e3b398, -33, 1.90821492927058781614e-10L), +lnhugeu = LD80C(0xb18b0c0330a8fad9, 13, 1.13627617309191834574e+04L), +lntinyu = LD80C(0xb236f28a68bc3bd7, 13, -1.14057368561139000667e+04L), +invln2u = LD80C(0xb8aa3b295c17f0bc, 0, 1.44269504088896340739e+00L); +#define ln2hi (ln2hiu.e) +#define ln2lo (ln2lou.e) +#define lnhuge (lnhugeu.e) +#define lntiny (lntinyu.e) +#define invln2 (invln2u.e) + +/* returns exp(r = x + c) for |c| < |x| with no overlap. */ + +static long double +__exp__D(long double x, long double c) +{ + long double hi, lo, z; + int k; + + if (x != x) /* x is NaN. */ + return(x); + + if (x <= lnhuge) { + if (x >= lntiny) { + /* argument reduction: x --> x - k*ln2 */ + z = invln2 * x; + k = z + copysignl(0.5L, x); + + /* + * Express (x + c) - k * ln2 as hi - lo. + * Let x = hi - lo rounded. + */ + hi = x - k * ln2hi; /* Exact. */ + lo = k * ln2lo - c; + x = hi - lo; + + /* Return 2^k*[1+x+x*c/(2+c)] */ + z = x * x; + c = x - z * (p1 + z * (p2 + z * (p3 + z * (p4 + + z * (p5 + z * (p6 + z * p7)))))); + c = (x * c) / (2 - c); + + return (ldexpl(1 + (hi - (lo - c)), k)); + } else { + /* exp(-INF) is 0. exp(-big) underflows to 0. */ + return (isfinite(x) ? ldexpl(1., -5000) : 0); + } + } else + /* exp(INF) is INF, exp(+big#) overflows to INF */ + return (isfinite(x) ? ldexpl(1., 5000) : x); +} Index: lib/msun/ld80/b_logl.c =================================================================== --- lib/msun/ld80/b_logl.c +++ lib/msun/ld80/b_logl.c @@ -29,52 +29,30 @@ * SUCH DAMAGE. */ -/* @(#)log.c 8.2 (Berkeley) 11/30/93 */ -#include -__FBSDID("$FreeBSD$"); - -#include - -#include "mathimpl.h" - -/* Table-driven natural logarithm. - * - * This code was derived, with minor modifications, from: - * Peter Tang, "Table-Driven Implementation of the - * Logarithm in IEEE Floating-Point arithmetic." ACM Trans. - * Math Software, vol 16. no 4, pp 378-400, Dec 1990). - * - * Calculates log(2^m*F*(1+f/F)), |f/j| <= 1/256, - * where F = j/128 for j an integer in [0, 128]. - * - * log(2^m) = log2_hi*m + log2_tail*m - * since m is an integer, the dominant term is exact. - * m has at most 10 digits (for subnormal numbers), - * and log2_hi has 11 trailing zero bits. - * - * log(F) = logF_hi[j] + logF_lo[j] is in tabular form in log_table.h - * logF_hi[] + 512 is exact. - * - * log(1+f/F) = 2*f/(2*F + f) + 1/12 * (2*f/(2*F + f))**3 + ... - * the leading term is calculated to extra precision in two - * parts, the larger of which adds exactly to the dominant - * m and F terms. - * There are two cases: - * 1. when m, j are non-zero (m | j), use absolute - * precision for the leading term. - * 2. when m = j = 0, |1-x| < 1/256, and log(x) ~= (x-1). - * In this case, use a relative precision of 24 bits. - * (This is done differently in the original paper) +/* + * See bsdsrc/b_log.c for implementation details. * - * Special cases: - * 0 return signalling -Inf - * neg return signalling NaN - * +Inf return +Inf -*/ + * bsdrc/b_log.c converted to long double by Steven G. Kargl. + */ #define N 128 -/* Table of log(Fj) = logF_head[j] + logF_tail[j], for Fj = 1+j/128. +/* + * Coefficients in the polynomial approximation of log(1+f/F). + * Domain of x is [0,1./256] with 2**(-84.48) precision. + */ +static const union IEEEl2bits + a1u = LD80C(0xaaaaaaaaaaaaaaab, -4, 8.33333333333333333356e-02L), + a2u = LD80C(0xcccccccccccccd29, -7, 1.25000000000000000781e-02L), + a3u = LD80C(0x9249249241ed3764, -9, 2.23214285711721994134e-03L), + a4u = LD80C(0xe38e959e1e7e01cf, -12, 4.34030476540000360640e-04L); +#define A1 (a1u.e) +#define A2 (a2u.e) +#define A3 (a3u.e) +#define A4 (a4u.e) + +/* + * Table of log(Fj) = logF_head[j] + logF_tail[j], for Fj = 1+j/128. * Used for generation of extend precision logarithms. * The constant 35184372088832 is 2^45, so the divide is exact. * It ensures correct reading of logF_head, even for inaccurate @@ -82,11 +60,7 @@ * right answer for integers less than 2^53.) * Values for log(F) were generated using error < 10^-57 absolute * with the bc -l package. -*/ -static double A1 = .08333333333333178827; -static double A2 = .01250000000377174923; -static double A3 = .002232139987919447809; -static double A4 = .0004348877777076145742; + */ static double logF_head[N+1] = { 0., @@ -351,118 +325,51 @@ .00000000000025144230728376072, -.00000000000017239444525614834 }; - -#if 0 -double -#ifdef _ANSI_SOURCE -log(double x) -#else -log(x) double x; -#endif -{ - int m, j; - double F, f, g, q, u, u2, v, zero = 0.0, one = 1.0; - volatile double u1; - - /* Catch special cases */ - if (x <= 0) - if (x == zero) /* log(0) = -Inf */ - return (-one/zero); - else /* log(neg) = NaN */ - return (zero/zero); - else if (!finite(x)) - return (x+x); /* x = NaN, Inf */ - - /* Argument reduction: 1 <= g < 2; x/2^m = g; */ - /* y = F*(1 + f/F) for |f| <= 2^-8 */ - - m = logb(x); - g = ldexp(x, -m); - if (m == -1022) { - j = logb(g), m += j; - g = ldexp(g, -j); - } - j = N*(g-1) + .5; - F = (1.0/N) * j + 1; /* F*128 is an integer in [128, 512] */ - f = g - F; - - /* Approximate expansion for log(1+f/F) ~= u + q */ - g = 1/(2*F+f); - u = 2*f*g; - v = u*u; - q = u*v*(A1 + v*(A2 + v*(A3 + v*A4))); - - /* case 1: u1 = u rounded to 2^-43 absolute. Since u < 2^-8, - * u1 has at most 35 bits, and F*u1 is exact, as F has < 8 bits. - * It also adds exactly to |m*log2_hi + log_F_head[j] | < 750 - */ - if (m | j) - u1 = u + 513, u1 -= 513; - - /* case 2: |1-x| < 1/256. The m- and j- dependent terms are zero; - * u1 = u to 24 bits. - */ - else - u1 = u, TRUNC(u1); - u2 = (2.0*(f - F*u1) - u1*f) * g; - /* u1 + u2 = 2f/(2F+f) to extra precision. */ - - /* log(x) = log(2^m*F*(1+f/F)) = */ - /* (m*log2_hi+logF_head[j]+u1) + (m*log2_lo+logF_tail[j]+q); */ - /* (exact) + (tiny) */ - - u1 += m*logF_head[N] + logF_head[j]; /* exact */ - u2 = (u2 + logF_tail[j]) + q; /* tiny */ - u2 += logF_tail[N]*m; - return (u1 + u2); -} -#endif - /* * Extra precision variant, returning struct {double a, b;}; - * log(x) = a+b to 63 bits, with a rounded to 26 bits. + * log(x) = a + b to 63 bits, with 'a' rounded to 24 bits. */ -struct Double -#ifdef _ANSI_SOURCE -__log__D(double x) -#else -__log__D(x) double x; -#endif +static struct Double +__log__D(long double x) { int m, j; - double F, f, g, q, u, v, u2; - volatile double u1; + long double F, f, g, q, u, v, u1, u2; struct Double r; - /* Argument reduction: 1 <= g < 2; x/2^m = g; */ - /* y = F*(1 + f/F) for |f| <= 2^-8 */ - - m = logb(x); - g = ldexp(x, -m); - if (m == -1022) { - j = logb(g), m += j; - g = ldexp(g, -j); + /* + * Argument reduction: 1 <= g < 2; x/2^m = g; + * y = F*(1 + f/F) for |f| <= 2^-8 + */ + g = frexpl(x, &m); + g *= 2; + m--; + if (m == DBL_MIN_EXP - 1) { + j = ilogbl(g); + m += j; + g = ldexpl(g, -j); } - j = N*(g-1) + .5; - F = (1.0/N) * j + 1; + j = N * (g - 1) + 0.5L; + F = (1.L / N) * j + 1; f = g - F; - g = 1/(2*F+f); - u = 2*f*g; - v = u*u; - q = u*v*(A1 + v*(A2 + v*(A3 + v*A4))); - if (m | j) - u1 = u + 513, u1 -= 513; - else - u1 = u, TRUNC(u1); - u2 = (2.0*(f - F*u1) - u1*f) * g; + g = 1 / (2 * F + f); + u = 2 * f * g; + v = u * u; + q = u * v * (A1 + v * (A2 + v * (A3 + v * A4))); + if (m | j) { + u1 = u + 513; + u1 -= 513; + } else { + u1 = (float)u; + } + u2 = (2 * (f - F * u1) - u1 * f) * g; - u1 += m*logF_head[N] + logF_head[j]; + u1 += m * (long double)logF_head[N] + logF_head[j]; - u2 += logF_tail[j]; u2 += q; - u2 += logF_tail[N]*m; - r.a = u1 + u2; /* Only difference is here */ - TRUNC(r.a); + u2 += logF_tail[j]; + u2 += q; + u2 += logF_tail[N] * m; + r.a = (float)(u1 + u2); /* Only difference is here. */ r.b = (u1 - r.a) + u2; return (r); } Index: lib/msun/ld80/b_tgammal.c =================================================================== --- /dev/null +++ lib/msun/ld80/b_tgammal.c @@ -0,0 +1,419 @@ +/*- + * SPDX-License-Identifier: BSD-3-Clause + * + * Copyright (c) 1992, 1993 + * The Regents of the University of California. All rights reserved. + * + * Redistribution and use in source and binary forms, with or without + * modification, are permitted provided that the following conditions + * are met: + * 1. Redistributions of source code must retain the above copyright + * notice, this list of conditions and the following disclaimer. + * 2. Redistributions in binary form must reproduce the above copyright + * notice, this list of conditions and the following disclaimer in the + * documentation and/or other materials provided with the distribution. + * 3. Neither the name of the University nor the names of its contributors + * may be used to endorse or promote products derived from this software + * without specific prior written permission. + * + * THIS SOFTWARE IS PROVIDED BY THE REGENTS AND CONTRIBUTORS ``AS IS'' AND + * ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE + * IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE + * ARE DISCLAIMED. IN NO EVENT SHALL THE REGENTS OR CONTRIBUTORS BE LIABLE + * FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL + * DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS + * OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) + * HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT + * LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY + * OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF + * SUCH DAMAGE. + */ + +/* + * The original code, FreeBSD's old svn r93211, contain the following + * attribution: + * + * This code by P. McIlroy, Oct 1992; + * + * The financial support of UUNET Communications Services is greatfully + * acknowledged. + * + * bsdrc/b_tgamma.c converted to long double by Steven G. Kargl. + */ + +/* + * See bsdsrc/t_tgamma.c for implementation details. + */ + +#include + +#if LDBL_MAX_EXP != 0x4000 +#error "Unsupported long double format" +#endif + +#ifdef __i386__ +#include +#endif + +#include "fpmath.h" +#include "math.h" +#include "math_private.h" + +/* Used in b_log.c and below. */ +struct Double { + long double a; + long double b; +}; + +#include "b_logl.c" +#include "b_expl.c" + +static const double zero = 0.; +static const volatile double tiny = 1e-300; +/* + * x >= 6 + * + * Use the asymptotic approximation (Stirling's formula) adjusted for + * equal-ripples: + * + * log(G(x)) ~= (x-0.5)*(log(x)-1) + 0.5(log(2*pi)-1) + 1/x*P(1/(x*x)) + * + * Keep extra precision in multiplying (x-.5)(log(x)-1), to avoid + * premature round-off. + * + * Accurate to max(ulp(1/128) absolute, 2^-66 relative) error. + */ + +/* + * The following is a decomposition of 0.5 * (log(2*pi) - 1) into the + * first 12 bits in ln2pi_hi and the trailing 64 bits in ln2pi_lo. The + * variables are clearly misnamed. + */ +static const union IEEEl2bits +ln2pi_hiu = LD80C(0xd680000000000000, -2, 4.18945312500000000000e-01L), +ln2pi_lou = LD80C(0xe379b414b596d687, -18, -6.77929532725821967032e-06L); +#define ln2pi_hi (ln2pi_hiu.e) +#define ln2pi_lo (ln2pi_lou.e) + +static const union IEEEl2bits + Pa0u = LD80C(0xaaaaaaaaaaaaaaaa, -4, 8.33333333333333333288e-02L), + Pa1u = LD80C(0xb60b60b60b5fcd59, -9, -2.77777777777776516326e-03L), + Pa2u = LD80C(0xd00d00cffbb47014, -11, 7.93650793635429639018e-04L), + Pa3u = LD80C(0x9c09c07c0805343e, -11, -5.95238087960599252215e-04L), + Pa4u = LD80C(0xdca8d31f8e6e5e8f, -11, 8.41749082509607342883e-04L), + Pa5u = LD80C(0xfb4d4289632f1638, -10, -1.91728055205541624556e-03L), + Pa6u = LD80C(0xd15a4ba04078d3f8, -8, 6.38893788027752396194e-03L), + Pa7u = LD80C(0xe877283110bcad95, -6, -2.83771309846297590312e-02L), + Pa8u = LD80C(0x8da97eed13717af8, -3, 1.38341887683837576925e-01L), + Pa9u = LD80C(0xf093b1c1584e30ce, -2, -4.69876818515470146031e-01L); +#define Pa0 (Pa0u.e) +#define Pa1 (Pa1u.e) +#define Pa2 (Pa2u.e) +#define Pa3 (Pa3u.e) +#define Pa4 (Pa4u.e) +#define Pa5 (Pa5u.e) +#define Pa6 (Pa6u.e) +#define Pa7 (Pa7u.e) +#define Pa8 (Pa8u.e) +#define Pa9 (Pa9u.e) + +static struct Double +large_gam(long double x) +{ + long double p, z, thi, tlo, xhi, xlo; + long double logx; + struct Double u; + + z = 1 / (x * x); + p = Pa0 + z * (Pa1 + z * (Pa2 + z * (Pa3 + z * (Pa4 + z * (Pa5 + + z * (Pa6 + z * (Pa7 + z * (Pa8 + z * Pa9)))))))); + p = p / x; + + u = __log__D(x); + u.a -= 1; + + /* Split (x - 0.5) in high and low parts. */ + x -= 0.5L; + xhi = (float)x; + xlo = x - xhi; + + /* Compute t = (x-.5)*(log(x)-1) in extra precision. */ + thi = xhi * u.a; + tlo = xlo * u.a + x * u.b; + + /* Compute thi + tlo + ln2pi_hi + ln2pi_lo + p. */ + tlo += ln2pi_lo; + tlo += p; + u.a = ln2pi_hi + tlo; + u.a += thi; + u.b = thi - u.a; + u.b += ln2pi_hi; + u.b += tlo; + return (u); +} +/* + * Rational approximation, A0 + x * x * P(x) / Q(x), on the interval + * [1.066.., 2.066..] accurate to 4.25e-19. + * + * Returns r.a + r.b = a0 + (z + c)^2 * p / q, with r.a truncated. + */ +static const union IEEEl2bits + a0_hiu = LD80C(0xe2b6e4153a57746c, -1, 8.85603194410888700265e-01L), + a0_lou = LD80C(0x851566d40f32c76d, -66, 1.40907742727049706207e-20L); +#define a0_hi (a0_hiu.e) +#define a0_lo (a0_lou.e) + +static const union IEEEl2bits +P0u = LD80C(0xdb629fb9bbdc1c1d, -2, 4.28486815855585429733e-01L), +P1u = LD80C(0xe6f4f9f5641aa6be, -3, 2.25543885805587730552e-01L), +P2u = LD80C(0xead1bd99fdaf7cc1, -6, 2.86644652514293482381e-02L), +P3u = LD80C(0x9ccc8b25838ab1e0, -8, 4.78512567772456362048e-03L), +P4u = LD80C(0x8f0c4383ef9ce72a, -9, 2.18273781132301146458e-03L), +P5u = LD80C(0xe732ab2c0a2778da, -13, 2.20487522485636008928e-04L), +P6u = LD80C(0xce70b27ca822b297, -16, 2.46095923774929264284e-05L), +P7u = LD80C(0xa309e2e16fb63663, -19, 2.42946473022376182921e-06L), +P8u = LD80C(0xaf9c110efb2c633d, -23, 1.63549217667765869987e-07L), +Q1u = LD80C(0xd4d7422719f48f15, -1, 8.31409582658993993626e-01L), +Q2u = LD80C(0xe13138ea404f1268, -5, -5.49785826915643198508e-02L), +Q3u = LD80C(0xd1c6cc91989352c0, -4, -1.02429960435139887683e-01L), +Q4u = LD80C(0xa7e9435a84445579, -7, 1.02484853505908820524e-02L), +Q5u = LD80C(0x83c7c34db89b7bda, -8, 4.02161632832052872697e-03L), +Q6u = LD80C(0xbed06bf6e1c14e5b, -11, -7.27898206351223022157e-04L), +Q7u = LD80C(0xef05bf841d4504c0, -18, 7.12342421869453515194e-06L), +Q8u = LD80C(0xf348d08a1ff53cb1, -19, 3.62522053809474067060e-06L); +#define P0 (P0u.e) +#define P1 (P1u.e) +#define P2 (P2u.e) +#define P3 (P3u.e) +#define P4 (P4u.e) +#define P5 (P5u.e) +#define P6 (P6u.e) +#define P7 (P7u.e) +#define P8 (P8u.e) +#define Q1 (Q1u.e) +#define Q2 (Q2u.e) +#define Q3 (Q3u.e) +#define Q4 (Q4u.e) +#define Q5 (Q5u.e) +#define Q6 (Q6u.e) +#define Q7 (Q7u.e) +#define Q8 (Q8u.e) + +static struct Double +ratfun_gam(long double z, long double c) +{ + long double p, q, thi, tlo; + struct Double r; + + q = 1 + z * (Q1 + z * (Q2 + z * (Q3 + z * (Q4 + z * (Q5 + + z * (Q6 + z * (Q7 + z * Q8))))))); + p = P0 + z * (P1 + z * (P2 + z * (P3 + z * (P4 + z * (P5 + + z * (P6 + z * (P7 + z * P8))))))); + p = p / q; + + /* Split z into high and low parts. */ + thi = (float)z; + tlo = (z - thi) + c; + tlo *= (thi + z); + + /* Split (z+c)^2 into high and low parts. */ + thi *= thi; + q = thi; + thi = (float)thi; + tlo += (q - thi); + + /* Split p/q into high and low parts. */ + r.a = (float)p; + r.b = p - r.a; + + tlo = tlo * p + thi * r.b + a0_lo; + thi *= r.a; /* t = (z+c)^2*(P/Q) */ + r.a = (float)(thi + a0_hi); + r.b = ((a0_hi - r.a) + thi) + tlo; + return (r); /* r = a0 + t */ +} +/* + * x < 6 + * + * Use argument reduction G(x+1) = xG(x) to reach the range [1.066124, + * 2.066124]. Use a rational approximation centered at the minimum + * (x0+1) to ensure monotonicity. + * + * Good to < 1 ulp. (provably .90 ulp; .87 ulp on 1,000,000 runs.) + * It also has correct monotonicity. + */ +static const union IEEEl2bits + xm1u = LD80C(0xec5b0c6ad7c7edc3, -2, 4.61632144968362341254e-01L); +#define x0 (xm1u.e) + +static const double + left = -0.3955078125; /* left boundary for rat. approx */ + +static long double +small_gam(long double x) +{ + long double t, y, ym1; + struct Double yy, r; + + y = x - 1; + + if (y <= 1 + (left + x0)) { + yy = ratfun_gam(y - x0, 0); + return (yy.a + yy.b); + } + + r.a = (float)y; + yy.a = r.a - 1; + y = y - 1 ; + r.b = yy.b = y - yy.a; + + /* Argument reduction: G(x+1) = x*G(x) */ + for (ym1 = y - 1; ym1 > left + x0; y = ym1--, yy.a--) { + t = r.a * yy.a; + r.b = r.a * yy.b + y * r.b; + r.a = (float)t; + r.b += (t - r.a); + } + + /* Return r*tgamma(y). */ + yy = ratfun_gam(y - x0, 0); + y = r.b * (yy.a + yy.b) + r.a * yy.b; + y += yy.a * r.a; + return (y); +} +/* + * Good on (0, 1+x0+left]. Accurate to 1 ulp. + */ +static long double +smaller_gam(long double x) +{ + long double d, rhi, rlo, t, xhi, xlo; + struct Double r; + + if (x < x0 + left) { + t = (float)x; + d = (t + x) * (x - t); + t *= t; + xhi = (float)(t + x); + xlo = x - xhi; + xlo += t; + xlo += d; + t = 1 - x0; + t += x; + d = 1 - x0; + d -= t; + d += x; + x = xhi + xlo; + } else { + xhi = (float)x; + xlo = x - xhi; + t = x - x0; + d = - x0 - t; + d += x; + } + + r = ratfun_gam(t, d); + d = (float)(r.a / x); + r.a -= d * xhi; + r.a -= d * xlo; + r.a += r.b; + + return (d + r.a / x); +} +/* + * x < 0 + * + * Use reflection formula, G(x) = pi/(sin(pi*x)*x*G(x)). + * At negative integers, return NaN and raise invalid. + */ +static const union IEEEl2bits +piu = LD80C(0xc90fdaa22168c235, 1, 3.14159265358979323851e+00L); +#define pi (piu.e) + +static long double +neg_gam(long double x) +{ + int sgn = 1; + struct Double lg, lsine; + long double y, z; + + y = ceill(x); + if (y == x) /* Negative integer. */ + return ((x - x) / zero); + + z = y - x; + if (z > 0.5) + z = 1 - z; + + y = y / 2; + if (y == ceill(y)) + sgn = -1; + + if (z < 0.25) + z = sinpil(z); + else + z = cospil(0.5 - z); + + /* Special case: G(1-x) = Inf; G(x) may be nonzero. */ + if (x < -1753) { + + if (x < -1760) + return (sgn * tiny * tiny); + y = expl(lgammal(x) / 2); + y *= y; + return (sgn < 0 ? -y : y); + } + + + y = 1 - x; + if (1 - y == x) + y = tgammal(y); + else /* 1-x is inexact */ + y = - x * tgammal(-x); + + if (sgn < 0) y = -y; + return (pi / (y * z)); +} +/* + * xmax comes from lgamma(xmax) - emax * log(2) = 0. + * static const float xmax = 35.040095f + * static const double xmax = 171.624376956302725; + * ld80: LD80C(0xdb718c066b352e20, 10, 1.75554834290446291689e+03L), + * ld128: 1.75554834290446291700388921607020320e+03L, + * + * iota is a sloppy threshold to isolate x = 0. + */ +static const double xmax = 1755.54834290446291689; +static const double iota = 0x1p-116; + +long double +tgammal(long double x) +{ + struct Double u; + + ENTERI(); + + if (x >= 6) { + if (x > xmax) + RETURNI(x / zero); + u = large_gam(x); + RETURNI(__exp__D(u.a, u.b)); + } + + if (x >= 1 + left + x0) + RETURNI(small_gam(x)); + + if (x > iota) + RETURNI(smaller_gam(x)); + + if (x > -iota) { + if (x != 0) + u.a = 1 - tiny; /* raise inexact */ + RETURNI(1 / x); + } + + if (!isfinite(x)) + RETURNI(x - x); /* x is NaN or -Inf */ + + RETURNI(neg_gam(x)); +} Index: lib/msun/src/imprecise.c =================================================================== --- /dev/null +++ lib/msun/src/imprecise.c @@ -1,57 +0,0 @@ -/*- - * SPDX-License-Identifier: BSD-2-Clause-FreeBSD - * - * Copyright (c) 2013 David Chisnall - * All rights reserved. - * - * Redistribution and use in source and binary forms, with or without - * modification, are permitted provided that the following conditions - * are met: - * 1. Redistributions of source code must retain the above copyright - * notice, this list of conditions and the following disclaimer. - * 2. Redistributions in binary form must reproduce the above copyright - * notice, this list of conditions and the following disclaimer in the - * documentation and/or other materials provided with the distribution. - * - * THIS SOFTWARE IS PROVIDED BY THE AUTHOR AND CONTRIBUTORS ``AS IS'' AND - * ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE - * IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE - * ARE DISCLAIMED. IN NO EVENT SHALL THE AUTHOR OR CONTRIBUTORS BE LIABLE - * FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL - * DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS - * OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) - * HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT - * LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY - * OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF - * SUCH DAMAGE. - * - * $FreeBSD$ - */ - -#include -#include - -/* - * If long double is not the same size as double, then these will lose - * precision and we should emit a warning whenever something links against - * them. - */ -#if (LDBL_MANT_DIG > 53) -#define WARN_IMPRECISE(x) \ - __warn_references(x, # x " has lower than advertised precision"); -#else -#define WARN_IMPRECISE(x) -#endif -/* - * Declare the functions as weak variants so that other libraries providing - * real versions can override them. - */ -#define DECLARE_WEAK(x)\ - __weak_reference(imprecise_## x, x);\ - WARN_IMPRECISE(x) - -#define DECLARE_IMPRECISE(f) \ - long double imprecise_ ## f ## l(long double v) { return f(v); }\ - DECLARE_WEAK(f ## l) - -DECLARE_IMPRECISE(tgamma);