* lib/msun/Makefile:

. Disconnect b_exp.c and b_log.c from the build.

* lib/msun/bsdsrc/b_exp.c:
  . Replace scalb() usage with C99's ldexp().
  . Replace finite(x) usage with C99's isfinite().
  . Whitespace changes towards style(9).
  . Remove include of "mathimpl.h".  It is no longer needed.
  . Remove #if 0 ... #endif code, which has been present since svn r93211
    (2002-03-26).
  . New minimax polynomial coefficients.
  . Add comments to explain origins of some constants.
  . Use ansi-C prototype.  Remove K&R prototype.  Add static to prototype.

* lib/msun/bsdsrc/b_log.c:
  . Remove include of "mathimpl.h".  It is no longer needed.
  . Fix comments to actually describe the code.
  . Reduce minimax polynomial from degree 4 to degree 3.
    This uses newly computed coefficients.
  . Use ansi-C prototype.  Remove K&R prototype.  Add static to prototype.
  . Remove volatile in declaration of u1.
  . Alphabetize decalaration list.
  . Whitespace changes towards style(9).
  . In argument reduction of x to g and m, replace use of logb() and
    ldexp() with a single call to frexp().  Add code to get 1 <= g < 2.
  . Remove #if 0 ... #endif code, which has been present since svn r93211
    (2002-03-26).
  . The special case m == -1022, replace logb() with ilogb().

* lib/msun/bsdsrc/b_tgamma.c:
  . Update comments.  Fix comments where needed.
  . Add float.h to get LDBL_MANT_DIG for weak reference of tgammal to tgamma.
  . Remove include of "mathimpl.h".  It is no longer needed.
  . Use "math.h" instead of <math.h>.
  . Add '#include math_private.h"
  . Add struct Double from mathimpl.h and include b_log.c and b_exp.c.
  . Remove forward declarations of neg_gam(), small_gam(), smaller_gam,
    large_gam() and ratfun_gam() by re-arranging the code to move these
    function above their first reference.
  . New minimax coefficients for polynomial in large_gam().
  . New splitting of a0 into a0hi nd a0lo, which include additional
    bits of precision.
  . Use ansi-C prototype.  Remove K&R prototype.
  . Replace the TRUNC() macro with a simple cast of a double entities
    to float before assignment (functional changes).
  . Replace sin(M_PI*z) with sinpi(z) and cos(M_PI*(0.5-z)) with cospi(0.5-z).

Submitted by:		Steve Kargl
Differential Revision:	https://reviews.freebsd.org/D33444
Reviewed by:		pfg
This commit is contained in:
Mark Murray 2021-12-14 09:02:45 +00:00
parent a46856c3f9
commit 455b2ccda3
4 changed files with 369 additions and 395 deletions

View file

@ -59,7 +59,7 @@ SHLIBDIR?= /lib
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 \

View file

@ -33,7 +33,6 @@
#include <sys/cdefs.h>
__FBSDID("$FreeBSD$");
/* EXP(X)
* RETURN THE EXPONENTIAL OF X
* DOUBLE PRECISION (IEEE 53 bits, VAX D FORMAT 56 BITS)
@ -41,14 +40,14 @@ __FBSDID("$FreeBSD$");
* 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 @@ __FBSDID("$FreeBSD$");
* with 1,156,000 random arguments on a VAX, the maximum observed
* error was 0.869 ulps (units in the last place).
*/
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 */
#include "mathimpl.h"
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 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
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);
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 and let x=hi-lo rounded */
/*
* 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;
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 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 scalb(1.+(hi-(lo - c)), k);
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);
}

View file

@ -33,10 +33,6 @@
#include <sys/cdefs.h>
__FBSDID("$FreeBSD$");
#include <math.h>
#include "mathimpl.h"
/* Table-driven natural logarithm.
*
* This code was derived, with minor modifications, from:
@ -44,25 +40,27 @@ __FBSDID("$FreeBSD$");
* 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 @@ __FBSDID("$FreeBSD$");
* 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 @@ __FBSDID("$FreeBSD$");
* 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 @@ static double logF_tail[N+1] = {
.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);
}

View file

@ -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 <sys/cdefs.h>
__FBSDID("$FreeBSD$");
#include <float.h>
#include "math.h"
#include "math_private.h"
/* Used in b_log.c and below. */
struct Double {
double a;
double b;
};
#include "b_log.c"
#include "b_exp.c"
/*
* This code by P. McIlroy, Oct 1992;
* The range is broken into several subranges. Each is handled by its
* helper functions.
*
* The financial support of UUNET Communications Services is greatfully
* acknowledged.
*/
#include <math.h>
#include "mathimpl.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))
*
* 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 @@ __FBSDID("$FreeBSD$");
* 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;
{
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));
}
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 struct Double
large_gam(x)
double x;
{
double z, p;
struct Double t, u, v;
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 */
z = one/(x*x);
p = Pa0+z*(Pa1+z*(Pa2+z*(Pa3+z*(Pa4+z*(Pa5+z*(Pa6+z*Pa7))))));
p = p/x;
static struct Double
large_gam(double x)
{
double p, z, thi, tlo, xhi, xlo;
struct Double u;
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 -= 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;
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);
}
/*
* 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
ratfun_gam(double z, double c)
{
double p, q, thi, tlo;
struct Double r;
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;
/* 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 @@ neg_gam(x)
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 @@ neg_gam(x)
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