/* Put in s an approximation of digamma(x). Assumes x >= 2. Assumes s does not overlap with x. Returns an integer e such that the error is bounded by 2^e ulps of the result s. */ static mpfr_exp_t mpfr_digamma_approx (mpfr_ptr s, mpfr_srcptr x) { mpfr_prec_t p = MPFR_PREC (s); mpfr_t t, u, invxx; mpfr_exp_t e, exps, f, expu; unsigned long n; MPFR_ASSERTN(MPFR_IS_POS(x) && (MPFR_EXP(x) >= 2)); mpfr_init2 (t, p); mpfr_init2 (u, p); mpfr_init2 (invxx, p); mpfr_log (s, x, MPFR_RNDN); /* error <= 1/2 ulp */ mpfr_ui_div (t, 1, x, MPFR_RNDN); /* error <= 1/2 ulp */ mpfr_div_2exp (t, t, 1, MPFR_RNDN); /* exact */ mpfr_sub (s, s, t, MPFR_RNDN); /* error <= 1/2 + 1/2*2^(EXP(olds)-EXP(s)) + 1/2*2^(EXP(t)-EXP(s)). For x >= 2, log(x) >= 2*(1/(2x)), thus olds >= 2t, and olds - t >= olds/2, thus 0 <= EXP(olds)-EXP(s) <= 1, and EXP(t)-EXP(s) <= 0, thus error <= 1/2 + 1/2*2 + 1/2 <= 2 ulps. */ e = 2; /* initial error */ mpfr_mul (invxx, x, x, MPFR_RNDZ); /* invxx = x^2 * (1 + theta) for |theta| <= 2^(-p) */ mpfr_ui_div (invxx, 1, invxx, MPFR_RNDU); /* invxx = 1/x^2 * (1 + theta)^2 */ /* in the following we note err=xxx when the ratio between the approximation and the exact result can be written (1 + theta)^xxx for |theta| <= 2^(-p), following Higham's method */ mpfr_set_ui (t, 1, MPFR_RNDN); /* err = 0 */ for (n = 1;; n++) { /* The main term is Bernoulli[2n]/(2n)/x^(2n) = B[n]/(2n+1)!(2n)/x^(2n) = B[n]*t[n]/(2n) where t[n]/t[n-1] = 1/(2n)/(2n+1)/x^2. */ mpfr_mul (t, t, invxx, MPFR_RNDU); /* err = err + 3 */ mpfr_div_ui (t, t, 2 * n, MPFR_RNDU); /* err = err + 1 */ mpfr_div_ui (t, t, 2 * n + 1, MPFR_RNDU); /* err = err + 1 */ /* we thus have err = 5n here */ mpfr_div_ui (u, t, 2 * n, MPFR_RNDU); /* err = 5n+1 */ mpfr_mul_z (u, u, mpfr_bernoulli_cache(n), MPFR_RNDU);/* err = 5n+2, and the absolute error is bounded by 10n+4 ulp(u) [Rule 11] */ /* if the terms 'u' are decreasing by a factor two at least, then the error coming from those is bounded by sum((10n+4)/2^n, n=1..infinity) = 24 */ exps = mpfr_get_exp (s); expu = mpfr_get_exp (u); if (expu < exps - (mpfr_exp_t) p) break; mpfr_sub (s, s, u, MPFR_RNDN); /* error <= 24 + n/2 */ if (mpfr_get_exp (s) < exps) e <<= exps - mpfr_get_exp (s); e ++; /* error in mpfr_sub */ f = 10 * n + 4; while (expu < exps) { f = (1 + f) / 2; expu ++; } e += f; /* total rounding error coming from 'u' term */ } mpfr_clear (t); mpfr_clear (u); mpfr_clear (invxx); f = 0; while (e > 1) { f++; e = (e + 1) / 2; /* Invariant: 2^f * e does not decrease */ } return f; }
/* Implements asymptotic expansion for jn or yn (formulae 9.2.5 and 9.2.6 from Abramowitz & Stegun). Assumes z > p log(2)/2, where p is the target precision. Return 0 if the expansion does not converge enough (the value 0 as inexact flag should not happen for normal input). */ static int FUNCTION (mpfr_ptr res, long n, mpfr_srcptr z, mp_rnd_t r) { mpfr_t s, c, P, Q, t, iz, err_t, err_s, err_u; mp_prec_t w; long k; int inex, stop, diverge = 0; mp_exp_t err2, err; MPFR_ZIV_DECL (loop); mpfr_init (c); w = MPFR_PREC(res) + MPFR_INT_CEIL_LOG2(MPFR_PREC(res)) + 4; MPFR_ZIV_INIT (loop, w); for (;;) { mpfr_set_prec (c, w); mpfr_init2 (s, w); mpfr_init2 (P, w); mpfr_init2 (Q, w); mpfr_init2 (t, w); mpfr_init2 (iz, w); mpfr_init2 (err_t, 31); mpfr_init2 (err_s, 31); mpfr_init2 (err_u, 31); /* Approximate sin(z) and cos(z). In the following, err <= k means that the approximate value y and the true value x are related by y = x * (1 + u)^k with |u| <= 2^(-w), following Higham's method. */ mpfr_sin_cos (s, c, z, GMP_RNDN); if (MPFR_IS_NEG(z)) mpfr_neg (s, s, GMP_RNDN); /* compute jn/yn(|z|), fix sign later */ /* The absolute error on s/c is bounded by 1/2 ulp(1/2) <= 2^(-w-1). */ mpfr_add (t, s, c, GMP_RNDN); mpfr_sub (c, s, c, GMP_RNDN); mpfr_swap (s, t); /* now s approximates sin(z)+cos(z), and c approximates sin(z)-cos(z), with total absolute error bounded by 2^(1-w). */ /* precompute 1/(8|z|) */ mpfr_si_div (iz, MPFR_IS_POS(z) ? 1 : -1, z, GMP_RNDN); /* err <= 1 */ mpfr_div_2ui (iz, iz, 3, GMP_RNDN); /* compute P and Q */ mpfr_set_ui (P, 1, GMP_RNDN); mpfr_set_ui (Q, 0, GMP_RNDN); mpfr_set_ui (t, 1, GMP_RNDN); /* current term */ mpfr_set_ui (err_t, 0, GMP_RNDN); /* error on t */ mpfr_set_ui (err_s, 0, GMP_RNDN); /* error on P and Q (sum of errors) */ for (k = 1, stop = 0; stop < 4; k++) { /* compute next term: t(k)/t(k-1) = (2n+2k-1)(2n-2k+1)/(8kz) */ mpfr_mul_si (t, t, 2 * (n + k) - 1, GMP_RNDN); /* err <= err_k + 1 */ mpfr_mul_si (t, t, 2 * (n - k) + 1, GMP_RNDN); /* err <= err_k + 2 */ mpfr_div_ui (t, t, k, GMP_RNDN); /* err <= err_k + 3 */ mpfr_mul (t, t, iz, GMP_RNDN); /* err <= err_k + 5 */ /* the relative error on t is bounded by (1+u)^(5k)-1, which is bounded by 6ku for 6ku <= 0.02: first |5 log(1+u)| <= |5.5u| for |u| <= 0.15, then |exp(5.5u)-1| <= 6u for |u| <= 0.02. */ mpfr_mul_ui (err_t, t, 6 * k, MPFR_IS_POS(t) ? GMP_RNDU : GMP_RNDD); mpfr_abs (err_t, err_t, GMP_RNDN); /* exact */ /* the absolute error on t is bounded by err_t * 2^(-w) */ mpfr_abs (err_u, t, GMP_RNDU); mpfr_mul_2ui (err_u, err_u, w, GMP_RNDU); /* t * 2^w */ mpfr_add (err_u, err_u, err_t, GMP_RNDU); /* max|t| * 2^w */ if (stop >= 2) { /* take into account the neglected terms: t * 2^w */ mpfr_div_2ui (err_s, err_s, w, GMP_RNDU); if (MPFR_IS_POS(t)) mpfr_add (err_s, err_s, t, GMP_RNDU); else mpfr_sub (err_s, err_s, t, GMP_RNDU); mpfr_mul_2ui (err_s, err_s, w, GMP_RNDU); stop ++; } /* if k is odd, add to Q, otherwise to P */ else if (k & 1) { /* if k = 1 mod 4, add, otherwise subtract */ if ((k & 2) == 0) mpfr_add (Q, Q, t, GMP_RNDN); else mpfr_sub (Q, Q, t, GMP_RNDN); /* check if the next term is smaller than ulp(Q): if EXP(err_u) <= EXP(Q), since the current term is bounded by err_u * 2^(-w), it is bounded by ulp(Q) */ if (MPFR_EXP(err_u) <= MPFR_EXP(Q)) stop ++; else stop = 0; } else { /* if k = 0 mod 4, add, otherwise subtract */ if ((k & 2) == 0) mpfr_add (P, P, t, GMP_RNDN); else mpfr_sub (P, P, t, GMP_RNDN); /* check if the next term is smaller than ulp(P) */ if (MPFR_EXP(err_u) <= MPFR_EXP(P)) stop ++; else stop = 0; } mpfr_add (err_s, err_s, err_t, GMP_RNDU); /* the sum of the rounding errors on P and Q is bounded by err_s * 2^(-w) */ /* stop when start to diverge */ if (stop < 2 && ((MPFR_IS_POS(z) && mpfr_cmp_ui (z, (k + 1) / 2) < 0) || (MPFR_IS_NEG(z) && mpfr_cmp_si (z, - ((k + 1) / 2)) > 0))) { /* if we have to stop the series because it diverges, then increasing the precision will most probably fail, since we will stop to the same point, and thus compute a very similar approximation */ diverge = 1; stop = 2; /* force stop */ } } /* the sum of the total errors on P and Q is bounded by err_s * 2^(-w) */ /* Now combine: the sum of the rounding errors on P and Q is bounded by err_s * 2^(-w), and the absolute error on s/c is bounded by 2^(1-w) */ if ((n & 1) == 0) /* n even: P * (sin + cos) + Q (cos - sin) for jn Q * (sin + cos) + P (sin - cos) for yn */ { #ifdef MPFR_JN mpfr_mul (c, c, Q, GMP_RNDN); /* Q * (sin - cos) */ mpfr_mul (s, s, P, GMP_RNDN); /* P * (sin + cos) */ #else mpfr_mul (c, c, P, GMP_RNDN); /* P * (sin - cos) */ mpfr_mul (s, s, Q, GMP_RNDN); /* Q * (sin + cos) */ #endif err = MPFR_EXP(c); if (MPFR_EXP(s) > err) err = MPFR_EXP(s); #ifdef MPFR_JN mpfr_sub (s, s, c, GMP_RNDN); #else mpfr_add (s, s, c, GMP_RNDN); #endif } else /* n odd: P * (sin - cos) + Q (cos + sin) for jn, Q * (sin - cos) - P (cos + sin) for yn */ { #ifdef MPFR_JN mpfr_mul (c, c, P, GMP_RNDN); /* P * (sin - cos) */ mpfr_mul (s, s, Q, GMP_RNDN); /* Q * (sin + cos) */ #else mpfr_mul (c, c, Q, GMP_RNDN); /* Q * (sin - cos) */ mpfr_mul (s, s, P, GMP_RNDN); /* P * (sin + cos) */ #endif err = MPFR_EXP(c); if (MPFR_EXP(s) > err) err = MPFR_EXP(s); #ifdef MPFR_JN mpfr_add (s, s, c, GMP_RNDN); #else mpfr_sub (s, c, s, GMP_RNDN); #endif } if ((n & 2) != 0) mpfr_neg (s, s, GMP_RNDN); if (MPFR_EXP(s) > err) err = MPFR_EXP(s); /* the absolute error on s is bounded by P*err(s/c) + Q*err(s/c) + err(P)*(s/c) + err(Q)*(s/c) + 3 * 2^(err - w - 1) <= (|P|+|Q|) * 2^(1-w) + err_s * 2^(1-w) + 2^err * 2^(1-w), since |c|, |old_s| <= 2. */ err2 = (MPFR_EXP(P) >= MPFR_EXP(Q)) ? MPFR_EXP(P) + 2 : MPFR_EXP(Q) + 2; /* (|P| + |Q|) * 2^(1 - w) <= 2^(err2 - w) */ err = MPFR_EXP(err_s) >= err ? MPFR_EXP(err_s) + 2 : err + 2; /* err_s * 2^(1-w) + 2^old_err * 2^(1-w) <= 2^err * 2^(-w) */ err2 = (err >= err2) ? err + 1 : err2 + 1; /* now the absolute error on s is bounded by 2^(err2 - w) */ /* multiply by sqrt(1/(Pi*z)) */ mpfr_const_pi (c, GMP_RNDN); /* Pi, err <= 1 */ mpfr_mul (c, c, z, GMP_RNDN); /* err <= 2 */ mpfr_si_div (c, MPFR_IS_POS(z) ? 1 : -1, c, GMP_RNDN); /* err <= 3 */ mpfr_sqrt (c, c, GMP_RNDN); /* err<=5/2, thus the absolute error is bounded by 3*u*|c| for |u| <= 0.25 */ mpfr_mul (err_t, c, s, MPFR_SIGN(c)==MPFR_SIGN(s) ? GMP_RNDU : GMP_RNDD); mpfr_abs (err_t, err_t, GMP_RNDU); mpfr_mul_ui (err_t, err_t, 3, GMP_RNDU); /* 3*2^(-w)*|old_c|*|s| [see below] is bounded by err_t * 2^(-w) */ err2 += MPFR_EXP(c); /* |old_c| * 2^(err2 - w) [see below] is bounded by 2^(err2-w) */ mpfr_mul (c, c, s, GMP_RNDN); /* the absolute error on c is bounded by 1/2 ulp(c) + 3*2^(-w)*|old_c|*|s| + |old_c| * 2^(err2 - w) */ /* compute err_t * 2^(-w) + 1/2 ulp(c) = (err_t + 2^EXP(c)) * 2^(-w) */ err = (MPFR_EXP(err_t) > MPFR_EXP(c)) ? MPFR_EXP(err_t) + 1 : MPFR_EXP(c) + 1; /* err_t * 2^(-w) + 1/2 ulp(c) <= 2^(err - w) */ /* now err_t * 2^(-w) bounds 1/2 ulp(c) + 3*2^(-w)*|old_c|*|s| */ err = (err >= err2) ? err + 1 : err2 + 1; /* the absolute error on c is bounded by 2^(err - w) */ mpfr_clear (s); mpfr_clear (P); mpfr_clear (Q); mpfr_clear (t); mpfr_clear (iz); mpfr_clear (err_t); mpfr_clear (err_s); mpfr_clear (err_u); err -= MPFR_EXP(c); if (MPFR_LIKELY (MPFR_CAN_ROUND (c, w - err, MPFR_PREC(res), r))) break; if (diverge != 0) { mpfr_set (c, z, r); /* will force inex=0 below, which means the asymptotic expansion failed */ break; } MPFR_ZIV_NEXT (loop, w); } MPFR_ZIV_FREE (loop); inex = MPFR_IS_POS(z) ? mpfr_set (res, c, r) : mpfr_neg (res, c, r); mpfr_clear (c); return inex; }
/* Don't need to save/restore exponent range: the cache does it. Catalan's constant is G = sum((-1)^k/(2*k+1)^2, k=0..infinity). We compute it using formula (31) of Victor Adamchik's page "33 representations for Catalan's constant" http://www-2.cs.cmu.edu/~adamchik/articles/catalan/catalan.htm G = Pi/8*log(2+sqrt(3)) + 3/8*sum(k!^2/(2k)!/(2k+1)^2,k=0..infinity) */ int mpfr_const_catalan_internal (mpfr_ptr g, mpfr_rnd_t rnd_mode) { mpfr_t x, y, z; mpz_t T, P, Q; mpfr_prec_t pg, p; int inex; MPFR_ZIV_DECL (loop); MPFR_GROUP_DECL (group); MPFR_LOG_FUNC (("rnd_mode=%d", rnd_mode), ("g[%Pu]=%.*Rg inex=%d", mpfr_get_prec (g), mpfr_log_prec, g, inex)); /* Here are the WC (max prec = 100.000.000) Once we have found a chain of 11, we only look for bigger chain. Found 3 '1' at 0 Found 5 '1' at 9 Found 6 '0' at 34 Found 9 '1' at 176 Found 11 '1' at 705 Found 12 '0' at 913 Found 14 '1' at 12762 Found 15 '1' at 152561 Found 16 '0' at 171725 Found 18 '0' at 525355 Found 20 '0' at 529245 Found 21 '1' at 6390133 Found 22 '0' at 7806417 Found 25 '1' at 11936239 Found 27 '1' at 51752950 */ pg = MPFR_PREC (g); p = pg + MPFR_INT_CEIL_LOG2 (pg) + 7; MPFR_GROUP_INIT_3 (group, p, x, y, z); mpz_init (T); mpz_init (P); mpz_init (Q); MPFR_ZIV_INIT (loop, p); for (;;) { mpfr_sqrt_ui (x, 3, MPFR_RNDU); mpfr_add_ui (x, x, 2, MPFR_RNDU); mpfr_log (x, x, MPFR_RNDU); mpfr_const_pi (y, MPFR_RNDU); mpfr_mul (x, x, y, MPFR_RNDN); S (T, P, Q, 0, (p - 1) / 2); mpz_mul_ui (T, T, 3); mpfr_set_z (y, T, MPFR_RNDU); mpfr_set_z (z, Q, MPFR_RNDD); mpfr_div (y, y, z, MPFR_RNDN); mpfr_add (x, x, y, MPFR_RNDN); mpfr_div_2ui (x, x, 3, MPFR_RNDN); if (MPFR_LIKELY (MPFR_CAN_ROUND (x, p - 5, pg, rnd_mode))) break; MPFR_ZIV_NEXT (loop, p); MPFR_GROUP_REPREC_3 (group, p, x, y, z); } MPFR_ZIV_FREE (loop); inex = mpfr_set (g, x, rnd_mode); MPFR_GROUP_CLEAR (group); mpz_clear (T); mpz_clear (P); mpz_clear (Q); return inex; }
/* This is the code of 'generic.c' slighty optimized for mpfr_atan. If x = p/2^r, put in y an approximation of atan(x)/x using 2^m terms for the series expansion, with an error of at most 1 ulp. Assumes |x| < 1. If X=x^2, we want 1 - X/3 + X^2/5 - ... + (-1)^k*X^k/(2k+1) + ... */ static void mpfr_atan_aux (mpfr_ptr y, mpz_ptr p, long r, int m, mpz_t *tab) { mpz_t *S, *Q, *ptoj; unsigned long n, i, k, j, l; mp_exp_t diff, expo; int im, done; mp_prec_t mult, *accu, *log2_nb_terms; mp_prec_t precy = MPFR_PREC(y); if (mpz_cmp_ui (p, 0) == 0) { mpfr_set_ui (y, 1, GMP_RNDN); /* limit(atan(x)/x, x=0) */ return; } accu = (mp_prec_t*) (*__gmp_allocate_func) ((2 * m + 2) * sizeof (mp_prec_t)); log2_nb_terms = accu + m + 1; /* Set Tables */ S = tab; /* S */ ptoj = S + 1*(m+1); /* p^2^j Precomputed table */ Q = S + 2*(m+1); /* Product of Odd integer table */ /* From p to p^2, and r to 2r */ mpz_mul (p, p, p); MPFR_ASSERTD (2 * r > r); r = 2 * r; /* Normalize p */ n = mpz_scan1 (p, 0); mpz_tdiv_q_2exp (p, p, n); /* exact */ MPFR_ASSERTD (r > n); r -= n; /* since |p/2^r| < 1, and p is a non-zero integer, necessarily r > 0 */ MPFR_ASSERTD (mpz_sgn (p) > 0); MPFR_ASSERTD (m > 0); /* check if p=1 (special case) */ l = 0; /* We compute by binary splitting, with X = x^2 = p/2^r: P(a,b) = p if a+1=b, P(a,c)*P(c,b) otherwise Q(a,b) = (2a+1)*2^r if a+1=b [except Q(0,1)=1], Q(a,c)*Q(c,b) otherwise S(a,b) = p*(2a+1) if a+1=b, Q(c,b)*S(a,c)+Q(a,c)*P(a,c)*S(c,b) otherwise Then atan(x)/x ~ S(0,i)/Q(0,i) for i so that (p/2^r)^i/i is small enough. The factor 2^(r*(b-a)) in Q(a,b) is implicit, thus we have to take it into account when we compute with Q. */ accu[0] = 0; /* accu[k] = Mult[0] + ... + Mult[k], where Mult[j] is the number of bits of the corresponding term S[j]/Q[j] */ if (mpz_cmp_ui (p, 1) != 0) { /* p <> 1: precompute ptoj table */ mpz_set (ptoj[0], p); for (im = 1 ; im <= m ; im ++) mpz_mul (ptoj[im], ptoj[im - 1], ptoj[im - 1]); /* main loop */ n = 1UL << m; /* the ith term being X^i/(2i+1) with X=p/2^r, we can stop when p^i/2^(r*i) < 2^(-precy), i.e. r*i > precy + log2(p^i) */ for (i = k = done = 0; (i < n) && (done == 0); i += 2, k ++) { /* initialize both S[k],Q[k] and S[k+1],Q[k+1] */ mpz_set_ui (Q[k+1], 2 * i + 3); /* Q(i+1,i+2) */ mpz_mul_ui (S[k+1], p, 2 * i + 1); /* S(i+1,i+2) */ mpz_mul_2exp (S[k], Q[k+1], r); mpz_sub (S[k], S[k], S[k+1]); /* S(i,i+2) */ mpz_mul_ui (Q[k], Q[k+1], 2 * i + 1); /* Q(i,i+2) */ log2_nb_terms[k] = 1; /* S[k]/Q[k] corresponds to 2 terms */ for (j = (i + 2) >> 1, l = 1; (j & 1) == 0; l ++, j >>= 1, k --) { /* invariant: S[k-1]/Q[k-1] and S[k]/Q[k] correspond to 2^l terms each. We combine them into S[k-1]/Q[k-1] */ MPFR_ASSERTD (k > 0); mpz_mul (S[k], S[k], Q[k-1]); mpz_mul (S[k], S[k], ptoj[l]); mpz_mul (S[k-1], S[k-1], Q[k]); mpz_mul_2exp (S[k-1], S[k-1], r << l); mpz_add (S[k-1], S[k-1], S[k]); mpz_mul (Q[k-1], Q[k-1], Q[k]); log2_nb_terms[k-1] = l + 1; /* now S[k-1]/Q[k-1] corresponds to 2^(l+1) terms */ MPFR_MPZ_SIZEINBASE2(mult, ptoj[l+1]); /* FIXME: precompute bits(ptoj[l+1]) outside the loop? */ mult = (r << (l + 1)) - mult - 1; accu[k-1] = (k == 1) ? mult : accu[k-2] + mult; if (accu[k-1] > precy) done = 1; } } }
mpfr_prec_t mpfr_get_prec (mpfr_srcptr x) { return MPFR_PREC(x); }
/* (y, z) <- (sin(x), cos(x)), return value is 0 iff both results are exact ie, iff x = 0 */ int mpfr_sin_cos (mpfr_ptr y, mpfr_ptr z, mpfr_srcptr x, mpfr_rnd_t rnd_mode) { mpfr_prec_t prec, m; int neg, reduce; mpfr_t c, xr; mpfr_srcptr xx; mpfr_exp_t err, expx; int inexy, inexz; MPFR_ZIV_DECL (loop); MPFR_SAVE_EXPO_DECL (expo); MPFR_ASSERTN (y != z); if (MPFR_UNLIKELY (MPFR_IS_SINGULAR (x))) { if (MPFR_IS_NAN(x) || MPFR_IS_INF(x)) { MPFR_SET_NAN (y); MPFR_SET_NAN (z); MPFR_RET_NAN; } else /* x is zero */ { MPFR_ASSERTD (MPFR_IS_ZERO (x)); MPFR_SET_ZERO (y); MPFR_SET_SAME_SIGN (y, x); /* y = 0, thus exact, but z is inexact in case of underflow or overflow */ inexy = 0; /* y is exact */ inexz = mpfr_set_ui (z, 1, rnd_mode); return INEX(inexy,inexz); } } MPFR_LOG_FUNC (("x[%#R]=%R rnd=%d", x, x, rnd_mode), ("sin[%#R]=%R cos[%#R]=%R", y, y, z, z)); MPFR_SAVE_EXPO_MARK (expo); prec = MAX (MPFR_PREC (y), MPFR_PREC (z)); m = prec + MPFR_INT_CEIL_LOG2 (prec) + 13; expx = MPFR_GET_EXP (x); /* When x is close to 0, say 2^(-k), then there is a cancellation of about 2k bits in 1-cos(x)^2. FIXME: in that case, it would be more efficient to compute sin(x) directly. VL: This is partly done by using MPFR_FAST_COMPUTE_IF_SMALL_INPUT from the mpfr_sin and mpfr_cos functions. Moreover, any overflow on m is avoided. */ if (expx < 0) { /* Warning: in case y = x, and the first call to MPFR_FAST_COMPUTE_IF_SMALL_INPUT succeeds but the second fails, we will have clobbered the original value of x. The workaround is to first compute z = cos(x) in that case, since y and z are different. */ if (y != x) /* y and x differ, thus we can safely try to compute y first */ { MPFR_FAST_COMPUTE_IF_SMALL_INPUT ( y, x, -2 * expx, 2, 0, rnd_mode, { inexy = _inexact; goto small_input; });
/* evaluates erf(x) using the expansion at x=0: erf(x) = 2/sqrt(Pi) * sum((-1)^k*x^(2k+1)/k!/(2k+1), k=0..infinity) Assumes x is neither NaN nor infinite nor zero. Assumes also that e*x^2 <= n (target precision). */ static int mpfr_erf_0 (mpfr_ptr res, mpfr_srcptr x, double xf2, mpfr_rnd_t rnd_mode) { mpfr_prec_t n, m; mpfr_exp_t nuk, sigmak; double tauk; mpfr_t y, s, t, u; unsigned int k; int log2tauk; int inex; MPFR_ZIV_DECL (loop); n = MPFR_PREC (res); /* target precision */ /* initial working precision */ m = n + (mpfr_prec_t) (xf2 / LOG2) + 8 + MPFR_INT_CEIL_LOG2 (n); mpfr_init2 (y, m); mpfr_init2 (s, m); mpfr_init2 (t, m); mpfr_init2 (u, m); MPFR_ZIV_INIT (loop, m); for (;;) { mpfr_mul (y, x, x, MPFR_RNDU); /* err <= 1 ulp */ mpfr_set_ui (s, 1, MPFR_RNDN); mpfr_set_ui (t, 1, MPFR_RNDN); tauk = 0.0; for (k = 1; ; k++) { mpfr_mul (t, y, t, MPFR_RNDU); mpfr_div_ui (t, t, k, MPFR_RNDU); mpfr_div_ui (u, t, 2 * k + 1, MPFR_RNDU); sigmak = MPFR_GET_EXP (s); if (k % 2) mpfr_sub (s, s, u, MPFR_RNDN); else mpfr_add (s, s, u, MPFR_RNDN); sigmak -= MPFR_GET_EXP(s); nuk = MPFR_GET_EXP(u) - MPFR_GET_EXP(s); if ((nuk < - (mpfr_exp_t) m) && ((double) k >= xf2)) break; /* tauk <- 1/2 + tauk * 2^sigmak + (1+8k)*2^nuk */ tauk = 0.5 + mul_2exp (tauk, sigmak) + mul_2exp (1.0 + 8.0 * (double) k, nuk); } mpfr_mul (s, x, s, MPFR_RNDU); MPFR_SET_EXP (s, MPFR_GET_EXP (s) + 1); mpfr_const_pi (t, MPFR_RNDZ); mpfr_sqrt (t, t, MPFR_RNDZ); mpfr_div (s, s, t, MPFR_RNDN); tauk = 4.0 * tauk + 11.0; /* final ulp-error on s */ log2tauk = __gmpfr_ceil_log2 (tauk); if (MPFR_LIKELY (MPFR_CAN_ROUND (s, m - log2tauk, n, rnd_mode))) break; /* Actualisation of the precision */ MPFR_ZIV_NEXT (loop, m); mpfr_set_prec (y, m); mpfr_set_prec (s, m); mpfr_set_prec (t, m); mpfr_set_prec (u, m); } MPFR_ZIV_FREE (loop); inex = mpfr_set (res, s, rnd_mode); mpfr_clear (y); mpfr_clear (t); mpfr_clear (u); mpfr_clear (s); return inex; }
/* Since MPFR-3.0, return the usual inexact value. The erange flag is set if an error occurred in the conversion (y is NaN, +Inf, or -Inf that have no equivalent in mpf) */ int mpfr_get_f (mpf_ptr x, mpfr_srcptr y, mpfr_rnd_t rnd_mode) { int inex; mp_size_t sx, sy; mpfr_prec_t precx, precy; mp_limb_t *xp; int sh; if (MPFR_UNLIKELY(MPFR_IS_SINGULAR(y))) { if (MPFR_IS_ZERO(y)) { mpf_set_ui (x, 0); return 0; } else if (MPFR_IS_NAN (y)) { MPFR_SET_ERANGEFLAG (); return 0; } else /* y is plus infinity (resp. minus infinity), set x to the maximum value (resp. the minimum value) in precision PREC(x) */ { int i; mp_limb_t *xp; MPFR_SET_ERANGEFLAG (); /* To this day, [mp_exp_t] and mp_size_t are #defined as the same type */ EXP (x) = MP_SIZE_T_MAX; sx = PREC (x); SIZ (x) = sx; xp = PTR (x); for (i = 0; i < sx; i++) xp[i] = MP_LIMB_T_MAX; if (MPFR_IS_POS (y)) return -1; else { mpf_neg (x, x); return +1; } } } sx = PREC(x); /* number of limbs of the mantissa of x */ precy = MPFR_PREC(y); precx = (mpfr_prec_t) sx * GMP_NUMB_BITS; sy = MPFR_LIMB_SIZE (y); xp = PTR (x); /* since mpf numbers are represented in base 2^GMP_NUMB_BITS, we loose -EXP(y) % GMP_NUMB_BITS bits in the most significant limb */ sh = MPFR_GET_EXP(y) % GMP_NUMB_BITS; sh = sh <= 0 ? - sh : GMP_NUMB_BITS - sh; MPFR_ASSERTD (sh >= 0); if (precy + sh <= precx) /* we can copy directly */ { mp_size_t ds; MPFR_ASSERTN (sx >= sy); ds = sx - sy; if (sh != 0) { mp_limb_t out; out = mpn_rshift (xp + ds, MPFR_MANT(y), sy, sh); MPFR_ASSERTN (ds > 0 || out == 0); if (ds > 0) xp[--ds] = out; } else MPN_COPY (xp + ds, MPFR_MANT (y), sy); if (ds > 0) MPN_ZERO (xp, ds); EXP(x) = (MPFR_GET_EXP(y) + sh) / GMP_NUMB_BITS; inex = 0; } else /* we have to round to precx - sh bits */ { mpfr_t z; mp_size_t sz; /* Recall that precx = (mpfr_prec_t) sx * GMP_NUMB_BITS, thus removing sh bits (sh < GMP_NUMB_BITSS) won't reduce the number of limbs. */ mpfr_init2 (z, precx - sh); sz = MPFR_LIMB_SIZE (z); MPFR_ASSERTN (sx == sz); inex = mpfr_set (z, y, rnd_mode); /* warning, sh may change due to rounding, but then z is a power of two, thus we can safely ignore its last bit which is 0 */ sh = MPFR_GET_EXP(z) % GMP_NUMB_BITS; sh = sh <= 0 ? - sh : GMP_NUMB_BITS - sh; MPFR_ASSERTD (sh >= 0); if (sh != 0) { mp_limb_t out; out = mpn_rshift (xp, MPFR_MANT(z), sz, sh); /* If sh hasn't changed, it is the number of the non-significant bits in the lowest limb of z. Therefore out == 0. */ MPFR_ASSERTD (out == 0); (void) out; /* avoid a warning */ } else MPN_COPY (xp, MPFR_MANT(z), sz); EXP(x) = (MPFR_GET_EXP(z) + sh) / GMP_NUMB_BITS; mpfr_clear (z); } /* set size and sign */ SIZ(x) = (MPFR_FROM_SIGN_TO_INT(MPFR_SIGN(y)) < 0) ? -sx : sx; return inex; }
static void test_generic (mpfr_prec_t p0, mpfr_prec_t p1, unsigned int N) { mpfr_prec_t prec, xprec, yprec; mpfr_t x, y, z, t; #ifdef TWO_ARGS mpfr_t u; #elif defined(DOUBLE_ARG1) || defined(DOUBLE_ARG2) mpfr_t u; double d; #endif mpfr_rnd_t rnd; int inexact, compare, compare2; unsigned int n; unsigned long ctrt = 0, ctrn = 0; mpfr_exp_t old_emin, old_emax; old_emin = mpfr_get_emin (); old_emax = mpfr_get_emax (); mpfr_init (x); mpfr_init (y); mpfr_init (z); mpfr_init (t); #if defined(TWO_ARGS) || defined(DOUBLE_ARG1) || defined(DOUBLE_ARG2) mpfr_init (u); #endif /* generic test */ for (prec = p0; prec <= p1; prec++) { mpfr_set_prec (z, prec); mpfr_set_prec (t, prec); yprec = prec + 10; mpfr_set_prec (y, yprec); /* Note: in precision p1, we test 4 special cases. */ for (n = 0; n < (prec == p1 ? N + 4 : N); n++) { xprec = prec; if (randlimb () & 1) { xprec *= (double) randlimb () / MP_LIMB_T_MAX; if (xprec < MPFR_PREC_MIN) xprec = MPFR_PREC_MIN; } mpfr_set_prec (x, xprec); #ifdef TWO_ARGS mpfr_set_prec (u, xprec); #elif defined(DOUBLE_ARG1) || defined(DOUBLE_ARG2) mpfr_set_prec (u, IEEE_DBL_MANT_DIG); #endif if (n > 3 || prec < p1) { #if defined(RAND_FUNCTION) RAND_FUNCTION (x); #if defined(TWO_ARGS) || defined(DOUBLE_ARG1) || defined(DOUBLE_ARG2) RAND_FUNCTION (u); #endif #else tests_default_random (x, TEST_RANDOM_POS, TEST_RANDOM_EMIN, TEST_RANDOM_EMAX); #if defined(TWO_ARGS) || defined(DOUBLE_ARG1) || defined(DOUBLE_ARG2) tests_default_random (u, TEST_RANDOM_POS2, TEST_RANDOM_EMIN, TEST_RANDOM_EMAX); #endif #endif } else { /* Special cases tested in precision p1 if n <= 3. They are useful really in the extended exponent range. */ set_emin (MPFR_EMIN_MIN); set_emax (MPFR_EMAX_MAX); if (n <= 1) { mpfr_set_si (x, n == 0 ? 1 : -1, MPFR_RNDN); mpfr_set_exp (x, mpfr_get_emin ()); #if defined(TWO_ARGS) || defined(DOUBLE_ARG1) || defined(DOUBLE_ARG2) mpfr_set_si (u, randlimb () % 2 == 0 ? 1 : -1, MPFR_RNDN); mpfr_set_exp (u, mpfr_get_emin ()); #endif } else /* 2 <= n <= 3 */ { if (getenv ("MPFR_CHECK_MAX") == NULL) goto next_n; mpfr_set_si (x, n == 0 ? 1 : -1, MPFR_RNDN); mpfr_setmax (x, REDUCE_EMAX); #if defined(TWO_ARGS) || defined(DOUBLE_ARG1) || defined(DOUBLE_ARG2) mpfr_set_si (u, randlimb () % 2 == 0 ? 1 : -1, MPFR_RNDN); mpfr_setmax (u, mpfr_get_emax ()); #endif } } rnd = RND_RAND (); mpfr_clear_flags (); #ifdef DEBUG_TGENERIC TGENERIC_INFO (TEST_FUNCTION, MPFR_PREC (y)); #endif #if defined(TWO_ARGS) compare = TEST_FUNCTION (y, x, u, rnd); #elif defined(DOUBLE_ARG1) d = mpfr_get_d (u, rnd); compare = TEST_FUNCTION (y, d, x, rnd); #elif defined(DOUBLE_ARG2) d = mpfr_get_d (u, rnd); compare = TEST_FUNCTION (y, x, d, rnd); #else compare = TEST_FUNCTION (y, x, rnd); #endif TGENERIC_CHECK ("Bad inexact flag", (compare != 0) ^ (mpfr_inexflag_p () == 0)); ctrt++; if (MPFR_IS_SINGULAR (y)) { if (MPFR_IS_NAN (y) || mpfr_nanflag_p ()) TGENERIC_CHECK ("Bad NaN flag", MPFR_IS_NAN (y) && mpfr_nanflag_p ()); else if (MPFR_IS_INF (y)) TGENERIC_CHECK ("Bad overflow flag", (compare != 0) ^ (mpfr_overflow_p () == 0)); else if (MPFR_IS_ZERO (y)) TGENERIC_CHECK ("Bad underflow flag", (compare != 0) ^ (mpfr_underflow_p () == 0)); } else if (mpfr_overflow_p ()) { TGENERIC_CHECK ("Bad compare value (overflow)", compare != 0); mpfr_nexttoinf (y); TGENERIC_CHECK ("Should have been max MPFR number", MPFR_IS_INF (y)); } else if (mpfr_underflow_p ()) { TGENERIC_CHECK ("Bad compare value (underflow)", compare != 0); mpfr_nexttozero (y); TGENERIC_CHECK ("Should have been min MPFR number", MPFR_IS_ZERO (y)); } else if (mpfr_can_round (y, yprec, rnd, rnd, prec)) { ctrn++; mpfr_set (t, y, rnd); /* Risk of failures are known when some flags are already set before the function call. Do not set the erange flag, as it will remain set after the function call and no checks are performed in such a case (see the mpfr_erangeflag_p test below). */ if (randlimb () & 1) __gmpfr_flags = MPFR_FLAGS_ALL ^ MPFR_FLAGS_ERANGE; #ifdef DEBUG_TGENERIC TGENERIC_INFO (TEST_FUNCTION, MPFR_PREC (z)); #endif /* Let's increase the precision of the inputs in a random way. In most cases, this doesn't make any difference, but for the mpfr_fmod bug fixed in r6230, this triggers the bug. */ mpfr_prec_round (x, mpfr_get_prec (x) + (randlimb () & 15), MPFR_RNDN); #if defined(TWO_ARGS) mpfr_prec_round (u, mpfr_get_prec (u) + (randlimb () & 15), MPFR_RNDN); inexact = TEST_FUNCTION (z, x, u, rnd); #elif defined(DOUBLE_ARG1) inexact = TEST_FUNCTION (z, d, x, rnd); #elif defined(DOUBLE_ARG2) inexact = TEST_FUNCTION (z, x, d, rnd); #else inexact = TEST_FUNCTION (z, x, rnd); #endif if (mpfr_erangeflag_p ()) goto next_n; if (mpfr_nan_p (z) || mpfr_cmp (t, z) != 0) { printf ("results differ for x="); mpfr_out_str (stdout, 2, xprec, x, MPFR_RNDN); #ifdef TWO_ARGS printf ("\nu="); mpfr_out_str (stdout, 2, xprec, u, MPFR_RNDN); #elif defined(DOUBLE_ARG1) || defined(DOUBLE_ARG2) printf ("\nu="); mpfr_out_str (stdout, 2, IEEE_DBL_MANT_DIG, u, MPFR_RNDN); #endif printf (" prec=%u rnd_mode=%s\n", (unsigned) prec, mpfr_print_rnd_mode (rnd)); printf ("got "); mpfr_out_str (stdout, 2, prec, z, MPFR_RNDN); puts (""); printf ("expected "); mpfr_out_str (stdout, 2, prec, t, MPFR_RNDN); puts (""); printf ("approx "); mpfr_print_binary (y); puts (""); exit (1); } compare2 = mpfr_cmp (t, y); /* if rounding to nearest, cannot know the sign of t - f(x) because of composed rounding: y = o(f(x)) and t = o(y) */ if (compare * compare2 >= 0) compare = compare + compare2; else compare = inexact; /* cannot determine sign(t-f(x)) */ if (((inexact == 0) && (compare != 0)) || ((inexact > 0) && (compare <= 0)) || ((inexact < 0) && (compare >= 0))) { printf ("Wrong inexact flag for rnd=%s: expected %d, got %d" "\n", mpfr_print_rnd_mode (rnd), compare, inexact); printf ("x="); mpfr_print_binary (x); puts (""); #if defined(TWO_ARGS) || defined(DOUBLE_ARG1) || defined(DOUBLE_ARG2) printf ("u="); mpfr_print_binary (u); puts (""); #endif printf ("y="); mpfr_print_binary (y); puts (""); printf ("t="); mpfr_print_binary (t); puts (""); exit (1); } } else if (getenv ("MPFR_SUSPICIOUS_OVERFLOW") != NULL) { /* For developers only! */ MPFR_ASSERTN (MPFR_IS_PURE_FP (y)); mpfr_nexttoinf (y); if (MPFR_IS_INF (y) && MPFR_IS_LIKE_RNDZ (rnd, MPFR_IS_NEG (y)) && !mpfr_overflow_p ()) { printf ("Possible bug! |y| is the maximum finite number " "and has been obtained when\nrounding toward zero" " (%s). Thus there is a very probable overflow,\n" "but the overflow flag is not set!\n", mpfr_print_rnd_mode (rnd)); printf ("x="); mpfr_print_binary (x); puts (""); #if defined(TWO_ARGS) || defined(DOUBLE_ARG1) || defined(DOUBLE_ARG2) printf ("u="); mpfr_print_binary (u); puts (""); #endif exit (1); } } next_n: /* In case the exponent range has been changed by tests_default_random() or for special values... */ mpfr_set_emin (old_emin); mpfr_set_emax (old_emax); } } #ifndef TGENERIC_NOWARNING if (3 * ctrn < 2 * ctrt) printf ("Warning! Too few normal cases in generic tests (%lu / %lu)\n", ctrn, ctrt); #endif mpfr_clear (x); mpfr_clear (y); mpfr_clear (z); mpfr_clear (t); #if defined(TWO_ARGS) || defined(DOUBLE_ARG1) || defined(DOUBLE_ARG2) mpfr_clear (u); #endif }
int mpfr_sqr (mpfr_ptr a, mpfr_srcptr b, mpfr_rnd_t rnd_mode) { int cc, inexact; mpfr_exp_t ax; mp_limb_t *tmp; mp_limb_t b1; mpfr_prec_t bq; mp_size_t bn, tn; MPFR_TMP_DECL(marker); MPFR_LOG_FUNC (("x[%#R]=%R rnd=%d", b, b, rnd_mode), ("y[%#R]=%R inexact=%d", a, a, inexact)); /* deal with special cases */ if (MPFR_UNLIKELY(MPFR_IS_SINGULAR(b))) { if (MPFR_IS_NAN(b)) { MPFR_SET_NAN(a); MPFR_RET_NAN; } MPFR_SET_POS (a); if (MPFR_IS_INF(b)) MPFR_SET_INF(a); else ( MPFR_ASSERTD(MPFR_IS_ZERO(b)), MPFR_SET_ZERO(a) ); MPFR_RET(0); } ax = 2 * MPFR_GET_EXP (b); bq = MPFR_PREC(b); MPFR_ASSERTD (2 * bq > bq); /* PREC_MAX is /2 so no integer overflow */ bn = MPFR_LIMB_SIZE(b); /* number of limbs of b */ tn = 1 + (2 * bq - 1) / GMP_NUMB_BITS; /* number of limbs of square, 2*bn or 2*bn-1 */ MPFR_TMP_MARK(marker); tmp = (mp_limb_t *) MPFR_TMP_ALLOC((size_t) 2 * bn * BYTES_PER_MP_LIMB); /* Multiplies the mantissa in temporary allocated space */ mpn_sqr_n (tmp, MPFR_MANT(b), bn); b1 = tmp[2 * bn - 1]; /* now tmp[0]..tmp[2*bn-1] contains the product of both mantissa, with tmp[2*bn-1]>=2^(GMP_NUMB_BITS-2) */ b1 >>= GMP_NUMB_BITS - 1; /* msb from the product */ /* if the mantissas of b and c are uniformly distributed in ]1/2, 1], then their product is in ]1/4, 1/2] with probability 2*ln(2)-1 ~ 0.386 and in [1/2, 1] with probability 2-2*ln(2) ~ 0.614 */ tmp += 2 * bn - tn; /* +0 or +1 */ if (MPFR_UNLIKELY(b1 == 0)) mpn_lshift (tmp, tmp, tn, 1); /* tn <= k, so no stack corruption */ cc = mpfr_round_raw (MPFR_MANT (a), tmp, 2 * bq, 0, MPFR_PREC (a), rnd_mode, &inexact); /* cc = 1 ==> result is a power of two */ if (MPFR_UNLIKELY(cc)) MPFR_MANT(a)[MPFR_LIMB_SIZE(a)-1] = MPFR_LIMB_HIGHBIT; MPFR_TMP_FREE(marker); { mpfr_exp_t ax2 = ax + (mpfr_exp_t) (b1 - 1 + cc); if (MPFR_UNLIKELY( ax2 > __gmpfr_emax)) return mpfr_overflow (a, rnd_mode, MPFR_SIGN_POS); if (MPFR_UNLIKELY( ax2 < __gmpfr_emin)) { /* In the rounding to the nearest mode, if the exponent of the exact result (i.e. before rounding, i.e. without taking cc into account) is < __gmpfr_emin - 1 or the exact result is a power of 2 (i.e. if both arguments are powers of 2), then round to zero. */ if (rnd_mode == MPFR_RNDN && (ax + (mpfr_exp_t) b1 < __gmpfr_emin || mpfr_powerof2_raw (b))) rnd_mode = MPFR_RNDZ; return mpfr_underflow (a, rnd_mode, MPFR_SIGN_POS); } MPFR_SET_EXP (a, ax2); MPFR_SET_POS (a); } MPFR_RET (inexact); }
/* Put in y an approximation of erfc(x) for large x, using formulae 7.1.23 and 7.1.24 from Abramowitz and Stegun. Returns e such that the error is bounded by 2^e ulp(y), or returns 0 in case of underflow. */ static mpfr_exp_t mpfr_erfc_asympt (mpfr_ptr y, mpfr_srcptr x) { mpfr_t t, xx, err; unsigned long k; mpfr_prec_t prec = MPFR_PREC(y); mpfr_exp_t exp_err; mpfr_init2 (t, prec); mpfr_init2 (xx, prec); mpfr_init2 (err, 31); /* let u = 2^(1-p), and let us represent the error as (1+u)^err with a bound for err */ mpfr_mul (xx, x, x, MPFR_RNDD); /* err <= 1 */ mpfr_ui_div (xx, 1, xx, MPFR_RNDU); /* upper bound for 1/(2x^2), err <= 2 */ mpfr_div_2ui (xx, xx, 1, MPFR_RNDU); /* exact */ mpfr_set_ui (t, 1, MPFR_RNDN); /* current term, exact */ mpfr_set (y, t, MPFR_RNDN); /* current sum */ mpfr_set_ui (err, 0, MPFR_RNDN); for (k = 1; ; k++) { mpfr_mul_ui (t, t, 2 * k - 1, MPFR_RNDU); /* err <= 4k-3 */ mpfr_mul (t, t, xx, MPFR_RNDU); /* err <= 4k */ /* for -1 < x < 1, and |nx| < 1, we have |(1+x)^n| <= 1+7/4|nx|. Indeed, for x>=0: log((1+x)^n) = n*log(1+x) <= n*x. Let y=n*x < 1, then exp(y) <= 1+7/4*y. For x<=0, let x=-x, we can prove by induction that (1-x)^n >= 1-n*x.*/ mpfr_mul_2si (err, err, MPFR_GET_EXP (y) - MPFR_GET_EXP (t), MPFR_RNDU); mpfr_add_ui (err, err, 14 * k, MPFR_RNDU); /* 2^(1-p) * t <= 2 ulp(t) */ mpfr_div_2si (err, err, MPFR_GET_EXP (y) - MPFR_GET_EXP (t), MPFR_RNDU); if (MPFR_GET_EXP (t) + (mpfr_exp_t) prec <= MPFR_GET_EXP (y)) { /* the truncation error is bounded by |t| < ulp(y) */ mpfr_add_ui (err, err, 1, MPFR_RNDU); break; } if (k & 1) mpfr_sub (y, y, t, MPFR_RNDN); else mpfr_add (y, y, t, MPFR_RNDN); } /* the error on y is bounded by err*ulp(y) */ mpfr_mul (t, x, x, MPFR_RNDU); /* rel. err <= 2^(1-p) */ mpfr_div_2ui (err, err, 3, MPFR_RNDU); /* err/8 */ mpfr_add (err, err, t, MPFR_RNDU); /* err/8 + xx */ mpfr_mul_2ui (err, err, 3, MPFR_RNDU); /* err + 8*xx */ mpfr_exp (t, t, MPFR_RNDU); /* err <= 1/2*ulp(t) + err(x*x)*t <= 1/2*ulp(t)+2*|x*x|*ulp(t) <= (2*|x*x|+1/2)*ulp(t) */ mpfr_mul (t, t, x, MPFR_RNDN); /* err <= 1/2*ulp(t) + (4*|x*x|+1)*ulp(t) <= (4*|x*x|+3/2)*ulp(t) */ mpfr_const_pi (xx, MPFR_RNDZ); /* err <= ulp(Pi) */ mpfr_sqrt (xx, xx, MPFR_RNDN); /* err <= 1/2*ulp(xx) + ulp(Pi)/2/sqrt(Pi) <= 3/2*ulp(xx) */ mpfr_mul (t, t, xx, MPFR_RNDN); /* err <= (8 |xx| + 13/2) * ulp(t) */ mpfr_div (y, y, t, MPFR_RNDN); /* the relative error on input y is bounded by (1+u)^err with u = 2^(1-p), that on t is bounded by (1+u)^(8 |xx| + 13/2), thus that on output y is bounded by 8 |xx| + 7 + err. */ if (MPFR_IS_ZERO(y)) { /* If y is zero, most probably we have underflow. We check it directly using the fact that erfc(x) <= exp(-x^2)/sqrt(Pi)/x for x >= 0. We compute an upper approximation of exp(-x^2)/sqrt(Pi)/x. */ mpfr_mul (t, x, x, MPFR_RNDD); /* t <= x^2 */ mpfr_neg (t, t, MPFR_RNDU); /* -x^2 <= t */ mpfr_exp (t, t, MPFR_RNDU); /* exp(-x^2) <= t */ mpfr_const_pi (xx, MPFR_RNDD); /* xx <= sqrt(Pi), cached */ mpfr_mul (xx, xx, x, MPFR_RNDD); /* xx <= sqrt(Pi)*x */ mpfr_div (y, t, xx, MPFR_RNDN); /* if y is zero, this means that the upper approximation of exp(-x^2)/sqrt(Pi)/x is nearer from 0 than from 2^(-emin-1), thus we have underflow. */ exp_err = 0; } else { mpfr_add_ui (err, err, 7, MPFR_RNDU); exp_err = MPFR_GET_EXP (err); } mpfr_clear (t); mpfr_clear (xx); mpfr_clear (err); return exp_err; }
int mpfr_fma (mpfr_ptr s, mpfr_srcptr x, mpfr_srcptr y, mpfr_srcptr z, mp_rnd_t rnd_mode) { int inexact; mpfr_t u; /* particular cases */ if (MPFR_IS_NAN(x) || MPFR_IS_NAN(y) || MPFR_IS_NAN(z)) { MPFR_SET_NAN(s); MPFR_RET_NAN; } if (MPFR_IS_INF(x) || MPFR_IS_INF(y)) { /* cases Inf*0+z, 0*Inf+z, Inf-Inf */ if ((MPFR_IS_FP(y) && MPFR_IS_ZERO(y)) || (MPFR_IS_FP(x) && MPFR_IS_ZERO(x)) || (MPFR_IS_INF(z) && ((MPFR_SIGN(x) * MPFR_SIGN(y)) != MPFR_SIGN(z)))) { MPFR_SET_NAN(s); MPFR_RET_NAN; } MPFR_CLEAR_NAN(s); if (MPFR_IS_INF(z)) /* case Inf-Inf already checked above */ { MPFR_SET_INF(s); MPFR_SET_SAME_SIGN(s, z); MPFR_RET(0); } else /* z is finite */ { MPFR_SET_INF(s); if (MPFR_SIGN(s) != (MPFR_SIGN(x) * MPFR_SIGN(y))) MPFR_CHANGE_SIGN(s); MPFR_RET(0); } } MPFR_CLEAR_NAN(s); /* now x and y are finite */ if (MPFR_IS_INF(z)) { MPFR_SET_INF(s); MPFR_SET_SAME_SIGN(s, z); MPFR_RET(0); } MPFR_CLEAR_INF(s); if (MPFR_IS_ZERO(x) || MPFR_IS_ZERO(y)) { if (MPFR_IS_ZERO(z)) { int sign_p, sign_z; sign_p = MPFR_SIGN(x) * MPFR_SIGN(y); sign_z = MPFR_SIGN(z); if (MPFR_SIGN(s) != (rnd_mode != GMP_RNDD ? ((sign_p < 0 && sign_z < 0) ? -1 : 1) : ((sign_p > 0 && sign_z > 0) ? 1 : -1))) { MPFR_CHANGE_SIGN(s); } MPFR_SET_ZERO(s); MPFR_RET(0); } else return mpfr_set (s, z, rnd_mode); } if (MPFR_IS_ZERO(z)) return mpfr_mul (s, x, y, rnd_mode); /* if we take prec(u) >= prec(x) + prec(y), the product u <- x*y is always exact */ mpfr_init2 (u, MPFR_PREC(x) + MPFR_PREC(y)); mpfr_mul (u, x, y, GMP_RNDN); /* always exact */ inexact = mpfr_add (s, z, u, rnd_mode); mpfr_clear(u); return inexact; }
int mpfr_zeta (mpfr_t z, mpfr_srcptr s, mpfr_rnd_t rnd_mode) { mpfr_t z_pre, s1, y, p; long add; mpfr_prec_t precz, prec1, precs, precs1; int inex; MPFR_GROUP_DECL (group); MPFR_ZIV_DECL (loop); MPFR_SAVE_EXPO_DECL (expo); MPFR_LOG_FUNC ( ("s[%Pu]=%.*Rg rnd=%d", mpfr_get_prec (s), mpfr_log_prec, s, rnd_mode), ("z[%Pu]=%.*Rg inexact=%d", mpfr_get_prec (z), mpfr_log_prec, z, inex)); /* Zero, Nan or Inf ? */ if (MPFR_UNLIKELY (MPFR_IS_SINGULAR (s))) { if (MPFR_IS_NAN (s)) { MPFR_SET_NAN (z); MPFR_RET_NAN; } else if (MPFR_IS_INF (s)) { if (MPFR_IS_POS (s)) return mpfr_set_ui (z, 1, MPFR_RNDN); /* Zeta(+Inf) = 1 */ MPFR_SET_NAN (z); /* Zeta(-Inf) = NaN */ MPFR_RET_NAN; } else /* s iz zero */ { MPFR_ASSERTD (MPFR_IS_ZERO (s)); return mpfr_set_si_2exp (z, -1, -1, rnd_mode); } } /* s is neither Nan, nor Inf, nor Zero */ /* check tiny s: we have zeta(s) = -1/2 - 1/2 log(2 Pi) s + ... around s=0, and for |s| <= 2^(-4), we have |zeta(s) + 1/2| <= |s|. EXP(s) + 1 < -PREC(z) is a sufficient condition to be able to round correctly, for any PREC(z) >= 1 (see algorithms.tex for details). */ if (MPFR_GET_EXP (s) + 1 < - (mpfr_exp_t) MPFR_PREC(z)) { int signs = MPFR_SIGN(s); MPFR_SAVE_EXPO_MARK (expo); mpfr_set_si_2exp (z, -1, -1, rnd_mode); /* -1/2 */ if (rnd_mode == MPFR_RNDA) rnd_mode = MPFR_RNDD; /* the result is around -1/2, thus negative */ if ((rnd_mode == MPFR_RNDU || rnd_mode == MPFR_RNDZ) && signs < 0) { mpfr_nextabove (z); /* z = -1/2 + epsilon */ inex = 1; } else if (rnd_mode == MPFR_RNDD && signs > 0) { mpfr_nextbelow (z); /* z = -1/2 - epsilon */ inex = -1; } else { if (rnd_mode == MPFR_RNDU) /* s > 0: z = -1/2 */ inex = 1; else if (rnd_mode == MPFR_RNDD) inex = -1; /* s < 0: z = -1/2 */ else /* (MPFR_RNDZ and s > 0) or MPFR_RNDN: z = -1/2 */ inex = (signs > 0) ? 1 : -1; } MPFR_SAVE_EXPO_FREE (expo); return mpfr_check_range (z, inex, rnd_mode); } /* Check for case s= -2n */ if (MPFR_IS_NEG (s)) { mpfr_t tmp; tmp[0] = *s; MPFR_EXP (tmp) = MPFR_GET_EXP (s) - 1; if (mpfr_integer_p (tmp)) { MPFR_SET_ZERO (z); MPFR_SET_POS (z); MPFR_RET (0); } } /* Check for case s=1 before changing the exponent range */ if (mpfr_cmp (s, __gmpfr_one) == 0) { MPFR_SET_INF (z); MPFR_SET_POS (z); MPFR_SET_DIVBY0 (); MPFR_RET (0); } MPFR_SAVE_EXPO_MARK (expo); /* Compute Zeta */ if (MPFR_IS_POS (s) && MPFR_GET_EXP (s) >= 0) /* Case s >= 1/2 */ inex = mpfr_zeta_pos (z, s, rnd_mode); else /* use reflection formula zeta(s) = 2^s*Pi^(s-1)*sin(Pi*s/2)*gamma(1-s)*zeta(1-s) */ { int overflow = 0; precz = MPFR_PREC (z); precs = MPFR_PREC (s); /* Precision precs1 needed to represent 1 - s, and s + 2, without any truncation */ precs1 = precs + 2 + MAX (0, - MPFR_GET_EXP (s)); /* Precision prec1 is the precision on elementary computations; it ensures a final precision prec1 - add for zeta(s) */ add = compute_add (s, precz); prec1 = precz + add; /* FIXME: To avoid that the working precision (prec1) depends on the input precision, one would need to take into account the error made when s1 is not exactly 1-s when computing zeta(s1) and gamma(s1) below, and also in the case y=Inf (i.e. when gamma(s1) overflows). Make sure that underflows do not occur in intermediate computations. Due to the limited precision, they are probably not possible in practice; add some MPFR_ASSERTN's to be sure that problems do not remain undetected? */ prec1 = MAX (prec1, precs1) + 10; MPFR_GROUP_INIT_4 (group, prec1, z_pre, s1, y, p); MPFR_ZIV_INIT (loop, prec1); for (;;) { mpfr_exp_t ey; mpfr_t z_up; mpfr_const_pi (p, MPFR_RNDD); /* p is Pi */ mpfr_sub (s1, __gmpfr_one, s, MPFR_RNDN); /* s1 = 1-s */ mpfr_gamma (y, s1, MPFR_RNDN); /* gamma(1-s) */ if (MPFR_IS_INF (y)) /* zeta(s) < 0 for -4k-2 < s < -4k, zeta(s) > 0 for -4k < s < -4k+2 */ { /* FIXME: An overflow in gamma(s1) does not imply that zeta(s) will overflow. A solution: 1. Compute log(|zeta(s)|/2) = (s-1)*log(2*pi) + lngamma(1-s) + log(abs(sin(Pi*s/2)) * zeta(1-s)) (possibly sharing computations with the normal case) with a rather good accuracy (see (2)). Memorize the sign of sin(...) for the final sign. 2. Take the exponential, ~= |zeta(s)|/2. If there is an overflow, then this means an overflow on the final result (due to the multiplication by 2, which has not been done yet). 3. Ziv test. 4. Correct the sign from the sign of sin(...). 5. Round then multiply by 2. Here, an overflow in either operation means a real overflow. */ mpfr_reflection_overflow (z_pre, s1, s, y, p, MPFR_RNDD); /* z_pre is a lower bound of |zeta(s)|/2, thus if it overflows, or has exponent emax, then |zeta(s)| overflows too. */ if (MPFR_IS_INF (z_pre) || MPFR_GET_EXP(z_pre) == __gmpfr_emax) { /* determine the sign of overflow */ mpfr_div_2ui (s1, s, 2, MPFR_RNDN); /* s/4, exact */ mpfr_frac (s1, s1, MPFR_RNDN); /* exact, -1 < s1 < 0 */ overflow = (mpfr_cmp_si_2exp (s1, -1, -1) > 0) ? -1 : 1; break; } else /* EXP(z_pre) < __gmpfr_emax */ { int ok = 0; mpfr_t z_down; mpfr_init2 (z_up, mpfr_get_prec (z_pre)); mpfr_reflection_overflow (z_up, s1, s, y, p, MPFR_RNDU); /* if the lower approximation z_pre does not overflow, but z_up does, we need more precision */ if (MPFR_IS_INF (z_up) || MPFR_GET_EXP(z_up) == __gmpfr_emax) goto next_loop; /* check if z_pre and z_up round to the same number */ mpfr_init2 (z_down, precz); mpfr_set (z_down, z_pre, rnd_mode); /* Note: it might be that EXP(z_down) = emax here, in that case we will have overflow below when we multiply by 2 */ mpfr_prec_round (z_up, precz, rnd_mode); ok = mpfr_cmp (z_down, z_up) == 0; mpfr_clear (z_up); mpfr_clear (z_down); if (ok) { /* get correct sign and multiply by 2 */ mpfr_div_2ui (s1, s, 2, MPFR_RNDN); /* s/4, exact */ mpfr_frac (s1, s1, MPFR_RNDN); /* exact, -1 < s1 < 0 */ if (mpfr_cmp_si_2exp (s1, -1, -1) > 0) mpfr_neg (z_pre, z_pre, rnd_mode); mpfr_mul_2ui (z_pre, z_pre, 1, rnd_mode); break; } else goto next_loop; } } mpfr_zeta_pos (z_pre, s1, MPFR_RNDN); /* zeta(1-s) */ mpfr_mul (z_pre, z_pre, y, MPFR_RNDN); /* gamma(1-s)*zeta(1-s) */ /* multiply z_pre by 2^s*Pi^(s-1) where p=Pi, s1=1-s */ mpfr_mul_2ui (y, p, 1, MPFR_RNDN); /* 2*Pi */ mpfr_neg (s1, s1, MPFR_RNDN); /* s-1 */ mpfr_pow (y, y, s1, MPFR_RNDN); /* (2*Pi)^(s-1) */ mpfr_mul (z_pre, z_pre, y, MPFR_RNDN); mpfr_mul_2ui (z_pre, z_pre, 1, MPFR_RNDN); /* multiply z_pre by sin(Pi*s/2) */ mpfr_mul (y, s, p, MPFR_RNDN); mpfr_div_2ui (p, y, 1, MPFR_RNDN); /* p = s*Pi/2 */ /* FIXME: sinpi will be available, we should replace the mpfr_sin call below by mpfr_sinpi(s/2), where s/2 will be exact. Can mpfr_sin underflow? Moreover, the code below should be improved so that the "if" condition becomes unlikely, e.g. by taking a slightly larger working precision. */ mpfr_sin (y, p, MPFR_RNDN); /* y = sin(Pi*s/2) */ ey = MPFR_GET_EXP (y); if (ey < 0) /* take account of cancellation in sin(p) */ { mpfr_t t; MPFR_ASSERTN (- ey < MPFR_PREC_MAX - prec1); mpfr_init2 (t, prec1 - ey); mpfr_const_pi (t, MPFR_RNDD); mpfr_mul (t, s, t, MPFR_RNDN); mpfr_div_2ui (t, t, 1, MPFR_RNDN); mpfr_sin (y, t, MPFR_RNDN); mpfr_clear (t); } mpfr_mul (z_pre, z_pre, y, MPFR_RNDN); if (MPFR_LIKELY (MPFR_CAN_ROUND (z_pre, prec1 - add, precz, rnd_mode))) break; next_loop: MPFR_ZIV_NEXT (loop, prec1); MPFR_GROUP_REPREC_4 (group, prec1, z_pre, s1, y, p); } MPFR_ZIV_FREE (loop); if (overflow != 0) { inex = mpfr_overflow (z, rnd_mode, overflow); MPFR_SAVE_EXPO_UPDATE_FLAGS (expo, MPFR_FLAGS_OVERFLOW); } else inex = mpfr_set (z, z_pre, rnd_mode); MPFR_GROUP_CLEAR (group); } MPFR_SAVE_EXPO_FREE (expo); return mpfr_check_range (z, inex, rnd_mode); }
/* Input: s - a floating-point number >= 1/2. rnd_mode - a rounding mode. Assumes s is neither NaN nor Infinite. Output: z - Zeta(s) rounded to the precision of z with direction rnd_mode */ static int mpfr_zeta_pos (mpfr_t z, mpfr_srcptr s, mpfr_rnd_t rnd_mode) { mpfr_t b, c, z_pre, f, s1; double beta, sd, dnep; mpfr_t *tc1; mpfr_prec_t precz, precs, d, dint; int p, n, l, add; int inex; MPFR_GROUP_DECL (group); MPFR_ZIV_DECL (loop); MPFR_ASSERTD (MPFR_IS_POS (s) && MPFR_GET_EXP (s) >= 0); precz = MPFR_PREC (z); precs = MPFR_PREC (s); /* Zeta(x) = 1+1/2^x+1/3^x+1/4^x+1/5^x+O(1/6^x) so with 2^(EXP(x)-1) <= x < 2^EXP(x) So for x > 2^3, k^x > k^8, so 2/k^x < 2/k^8 Zeta(x) = 1 + 1/2^x*(1+(2/3)^x+(2/4)^x+...) = 1 + 1/2^x*(1+sum((2/k)^x,k=3..infinity)) <= 1 + 1/2^x*(1+sum((2/k)^8,k=3..infinity)) And sum((2/k)^8,k=3..infinity) = -257+128*Pi^8/4725 ~= 0.0438035 So Zeta(x) <= 1 + 1/2^x*2 for x >= 8 The error is < 2^(-x+1) <= 2^(-2^(EXP(x)-1)+1) */ if (MPFR_GET_EXP (s) > 3) { mpfr_exp_t err; err = MPFR_GET_EXP (s) - 1; if (err > (mpfr_exp_t) (sizeof (mpfr_exp_t)*CHAR_BIT-2)) err = MPFR_EMAX_MAX; else err = ((mpfr_exp_t)1) << err; err = 1 - (-err+1); /* GET_EXP(one) - (-err+1) = err :) */ MPFR_FAST_COMPUTE_IF_SMALL_INPUT (z, __gmpfr_one, err, 0, 1, rnd_mode, {}); } d = precz + MPFR_INT_CEIL_LOG2(precz) + 10; /* we want that s1 = s-1 is exact, i.e. we should have PREC(s1) >= EXP(s) */ dint = (mpfr_uexp_t) MPFR_GET_EXP (s); mpfr_init2 (s1, MAX (precs, dint)); inex = mpfr_sub (s1, s, __gmpfr_one, MPFR_RNDN); MPFR_ASSERTD (inex == 0); /* case s=1 should have already been handled */ MPFR_ASSERTD (!MPFR_IS_ZERO (s1)); MPFR_GROUP_INIT_4 (group, MPFR_PREC_MIN, b, c, z_pre, f); MPFR_ZIV_INIT (loop, d); for (;;) { /* Principal loop: we compute, in z_pre, an approximation of Zeta(s), that we send to can_round */ if (MPFR_GET_EXP (s1) <= -(mpfr_exp_t) ((mpfr_prec_t) (d-3)/2)) /* Branch 1: when s-1 is very small, one uses the approximation Zeta(s)=1/(s-1)+gamma, where gamma is Euler's constant */ { dint = MAX (d + 3, precs); /* branch 1, with internal precision dint */ MPFR_GROUP_REPREC_4 (group, dint, b, c, z_pre, f); mpfr_div (z_pre, __gmpfr_one, s1, MPFR_RNDN); mpfr_const_euler (f, MPFR_RNDN); mpfr_add (z_pre, z_pre, f, MPFR_RNDN); } else /* Branch 2 */ { size_t size; /* branch 2 */ /* Computation of parameters n, p and working precision */ dnep = (double) d * LOG2; sd = mpfr_get_d (s, MPFR_RNDN); /* beta = dnep + 0.61 + sd * log (6.2832 / sd); but a larger value is OK */ #define LOG6dot2832 1.83787940484160805532 beta = dnep + 0.61 + sd * (LOG6dot2832 - LOG2 * __gmpfr_floor_log2 (sd)); if (beta <= 0.0) { p = 0; /* n = 1 + (int) (exp ((dnep - LOG2) / sd)); */ n = 1 + (int) __gmpfr_ceil_exp2 ((d - 1.0) / sd); } else { p = 1 + (int) beta / 2; n = 1 + (int) ((sd + 2.0 * (double) p - 1.0) / 6.2832); } /* add = 4 + floor(1.5 * log(d) / log (2)). We should have add >= 10, which is always fulfilled since d = precz + 11 >= 12, thus ceil(log2(d)) >= 4 */ add = 4 + (3 * MPFR_INT_CEIL_LOG2 (d)) / 2; MPFR_ASSERTD(add >= 10); dint = d + add; if (dint < precs) dint = precs; /* internal precision is dint */ size = (p + 1) * sizeof(mpfr_t); tc1 = (mpfr_t*) mpfr_allocate_func (size); for (l=1; l<=p; l++) mpfr_init2 (tc1[l], dint); MPFR_GROUP_REPREC_4 (group, dint, b, c, z_pre, f); /* precision of z is precz */ /* Computation of the coefficients c_k */ mpfr_zeta_c (p, tc1); /* Computation of the 3 parts of the function Zeta. */ mpfr_zeta_part_a (z_pre, s, n); mpfr_zeta_part_b (b, s, n, p, tc1); /* s1 = s-1 is already computed above */ mpfr_div (c, __gmpfr_one, s1, MPFR_RNDN); mpfr_ui_pow (f, n, s1, MPFR_RNDN); mpfr_div (c, c, f, MPFR_RNDN); mpfr_add (z_pre, z_pre, c, MPFR_RNDN); mpfr_add (z_pre, z_pre, b, MPFR_RNDN); for (l=1; l<=p; l++) mpfr_clear (tc1[l]); mpfr_free_func (tc1, size); /* End branch 2 */ } if (MPFR_LIKELY (MPFR_CAN_ROUND (z_pre, d-3, precz, rnd_mode))) break; MPFR_ZIV_NEXT (loop, d); } MPFR_ZIV_FREE (loop); inex = mpfr_set (z, z_pre, rnd_mode); MPFR_GROUP_CLEAR (group); mpfr_clear (s1); return inex; }
/* The test_one argument is seen a boolean. If it is true and rnd is a rounding mode toward infinity, then the function is tested in only one rounding mode (the one provided in rnd) and the variable rndnext is not used (due to the break). If it is true and rnd is a rounding mode toward or away from zero, then the function is tested twice, first with the provided rounding mode and second with the rounding mode toward the corresponding infinity (determined by the sign of the result). If it is false, then the function is tested in the 5 rounding modes, and rnd must initially be MPFR_RNDZ; thus rndnext will be initialized in the first iteration. If the test_one argument is 2, then this means that y is exact, and the ternary value is checked. As examples of use, see the calls to test5rm from the data_check and bad_cases functions. */ static void test5rm (int (*fct) (FLIST), mpfr_srcptr x, mpfr_ptr y, mpfr_ptr z, mpfr_rnd_t rnd, int test_one, const char *name) { mpfr_prec_t yprec = MPFR_PREC (y); mpfr_rnd_t rndnext = MPFR_RND_MAX; /* means uninitialized */ MPFR_ASSERTN (test_one || rnd == MPFR_RNDZ); mpfr_set_prec (z, yprec); while (1) { int inex; MPFR_ASSERTN (rnd != MPFR_RND_MAX); inex = fct (z, x, rnd); if (! (mpfr_equal_p (y, z) || (mpfr_nan_p (y) && mpfr_nan_p (z)))) { printf ("Error for %s with xprec=%lu, yprec=%lu, rnd=%s\nx = ", name, (unsigned long) MPFR_PREC (x), (unsigned long) yprec, mpfr_print_rnd_mode (rnd)); mpfr_out_str (stdout, 16, 0, x, MPFR_RNDN); printf ("\nexpected "); mpfr_out_str (stdout, 16, 0, y, MPFR_RNDN); printf ("\ngot "); mpfr_out_str (stdout, 16, 0, z, MPFR_RNDN); printf ("\n"); exit (1); } if (test_one == 2 && inex != 0) { printf ("Error for %s with xprec=%lu, yprec=%lu, rnd=%s\nx = ", name, (unsigned long) MPFR_PREC (x), (unsigned long) yprec, mpfr_print_rnd_mode (rnd)); mpfr_out_str (stdout, 16, 0, x, MPFR_RNDN); printf ("\nexact case, but non-zero ternary value (%d)\n", inex); exit (1); } if (rnd == MPFR_RNDN) break; if (test_one) { if (rnd == MPFR_RNDU || rnd == MPFR_RNDD) break; if (MPFR_IS_NEG (y)) rnd = (rnd == MPFR_RNDA) ? MPFR_RNDD : MPFR_RNDU; else rnd = (rnd == MPFR_RNDA) ? MPFR_RNDU : MPFR_RNDD; } else if (rnd == MPFR_RNDZ) { rnd = MPFR_IS_NEG (y) ? MPFR_RNDU : MPFR_RNDD; rndnext = MPFR_RNDA; } else { rnd = rndnext; if (rnd == MPFR_RNDA) { mpfr_nexttoinf (y); rndnext = (MPFR_IS_NEG (y)) ? MPFR_RNDD : MPFR_RNDU; } else if (rndnext != MPFR_RNDN) rndnext = MPFR_RNDN; else { if (yprec == MPFR_PREC_MIN) break; mpfr_prec_round (y, --yprec, MPFR_RNDZ); mpfr_set_prec (z, yprec); } } } }
/* computes tan(x) = sign(x)*sqrt(1/cos(x)^2-1) */ int mpfr_tan (mpfr_ptr y, mpfr_srcptr x, mpfr_rnd_t rnd_mode) { mpfr_prec_t precy, m; int inexact; mpfr_t s, c; MPFR_ZIV_DECL (loop); MPFR_SAVE_EXPO_DECL (expo); MPFR_GROUP_DECL (group); MPFR_LOG_FUNC (("x[%#R]=%R rnd=%d", x, x, rnd_mode), ("y[%#R]=%R inexact=%d", y, y, inexact)); if (MPFR_UNLIKELY(MPFR_IS_SINGULAR(x))) { if (MPFR_IS_NAN(x) || MPFR_IS_INF(x)) { MPFR_SET_NAN(y); MPFR_RET_NAN; } else /* x is zero */ { MPFR_ASSERTD(MPFR_IS_ZERO(x)); MPFR_SET_ZERO(y); MPFR_SET_SAME_SIGN(y, x); MPFR_RET(0); } } /* tan(x) = x + x^3/3 + ... so the error is < 2^(3*EXP(x)-1) */ MPFR_FAST_COMPUTE_IF_SMALL_INPUT (y, x, -2 * MPFR_GET_EXP (x), 1, 1, rnd_mode, {}); MPFR_SAVE_EXPO_MARK (expo); /* Compute initial precision */ precy = MPFR_PREC (y); m = precy + MPFR_INT_CEIL_LOG2 (precy) + 13; MPFR_ASSERTD (m >= 2); /* needed for the error analysis in algorithms.tex */ MPFR_GROUP_INIT_2 (group, m, s, c); MPFR_ZIV_INIT (loop, m); for (;;) { /* The only way to get an overflow is to get ~ Pi/2 But the result will be ~ 2^Prec(y). */ mpfr_sin_cos (s, c, x, MPFR_RNDN); /* err <= 1/2 ulp on s and c */ mpfr_div (c, s, c, MPFR_RNDN); /* err <= 4 ulps */ MPFR_ASSERTD (!MPFR_IS_SINGULAR (c)); if (MPFR_LIKELY (MPFR_CAN_ROUND (c, m - 2, precy, rnd_mode))) break; MPFR_ZIV_NEXT (loop, m); MPFR_GROUP_REPREC_2 (group, m, s, c); } MPFR_ZIV_FREE (loop); inexact = mpfr_set (y, c, rnd_mode); MPFR_GROUP_CLEAR (group); MPFR_SAVE_EXPO_FREE (expo); return mpfr_check_range (y, inexact, rnd_mode); }
int mpfr_cmp3 (mpfr_srcptr b, mpfr_srcptr c, int s) { mp_exp_t be, ce; mp_size_t bn, cn; mp_limb_t *bp, *cp; MPFR_ASSERTN(!MPFR_IS_NAN(b)); MPFR_ASSERTN(!MPFR_IS_NAN(c)); s *= MPFR_SIGN(c); if (MPFR_IS_INF(b)) return (MPFR_IS_INF(c) && s * MPFR_SIGN(b) > 0) ? 0 : MPFR_SIGN(b); if (MPFR_IS_INF(c)) return -s; /* b and c are real numbers */ if (MPFR_IS_ZERO(b)) return MPFR_IS_ZERO(c) ? 0 : -s; if (MPFR_IS_ZERO(c)) return MPFR_SIGN(b); if (s * MPFR_SIGN(b) < 0) return MPFR_SIGN(b); /* now signs are equal */ be = MPFR_EXP(b); ce = MPFR_EXP(c); if (be > ce) return s; if (be < ce) return -s; /* both signs and exponents are equal */ bn = (MPFR_PREC(b)-1)/BITS_PER_MP_LIMB; cn = (MPFR_PREC(c)-1)/BITS_PER_MP_LIMB; bp = MPFR_MANT(b); cp = MPFR_MANT(c); for ( ; bn >= 0 && cn >= 0; bn--, cn--) { if (bp[bn] > cp[cn]) return s; if (bp[bn] < cp[cn]) return -s; } for ( ; bn >= 0; bn--) if (bp[bn]) return s; for ( ; cn >= 0; cn--) if (cp[cn]) return -s; return 0; }
int mpfr_ui_pow_ui (mpfr_ptr x, unsigned long int y, unsigned long int n, mpfr_rnd_t rnd) { mpfr_exp_t err; unsigned long m; mpfr_t res; mpfr_prec_t prec; int size_n; int inexact; MPFR_ZIV_DECL (loop); MPFR_SAVE_EXPO_DECL (expo); if (MPFR_UNLIKELY (n <= 1)) { if (n == 1) return mpfr_set_ui (x, y, rnd); /* y^1 = y */ else return mpfr_set_ui (x, 1, rnd); /* y^0 = 1 for any y */ } else if (MPFR_UNLIKELY (y <= 1)) { if (y == 1) return mpfr_set_ui (x, 1, rnd); /* 1^n = 1 for any n > 0 */ else return mpfr_set_ui (x, 0, rnd); /* 0^n = 0 for any n > 0 */ } for (size_n = 0, m = n; m; size_n++, m >>= 1); MPFR_SAVE_EXPO_MARK (expo); prec = MPFR_PREC (x) + 3 + size_n; mpfr_init2 (res, prec); MPFR_ZIV_INIT (loop, prec); for (;;) { int i = size_n; inexact = mpfr_set_ui (res, y, MPFR_RNDU); err = 1; /* now 2^(i-1) <= n < 2^i: i=1+floor(log2(n)) */ for (i -= 2; i >= 0; i--) { inexact |= mpfr_mul (res, res, res, MPFR_RNDU); err++; if (n & (1UL << i)) inexact |= mpfr_mul_ui (res, res, y, MPFR_RNDU); } /* since the loop is executed floor(log2(n)) times, we have err = 1+floor(log2(n)). Since prec >= MPFR_PREC(x) + 4 + floor(log2(n)), prec > err */ err = prec - err; if (MPFR_LIKELY (inexact == 0 || MPFR_CAN_ROUND (res, err, MPFR_PREC (x), rnd))) break; /* Actualisation of the precision */ MPFR_ZIV_NEXT (loop, prec); mpfr_set_prec (res, prec); } MPFR_ZIV_FREE (loop); inexact = mpfr_set (x, res, rnd); mpfr_clear (res); MPFR_SAVE_EXPO_FREE (expo); return mpfr_check_range (x, inexact, rnd); }
int mpfr_grandom (mpfr_ptr rop1, mpfr_ptr rop2, gmp_randstate_t rstate, mpfr_rnd_t rnd) { int inex1, inex2, s1, s2; mpz_t x, y, xp, yp, t, a, b, s; mpfr_t sfr, l, r1, r2; mpfr_prec_t tprec, tprec0; inex2 = inex1 = 0; if (rop2 == NULL) /* only one output requested. */ { tprec0 = MPFR_PREC (rop1); } else { tprec0 = MAX (MPFR_PREC (rop1), MPFR_PREC (rop2)); } tprec0 += 11; /* We use "Marsaglia polar method" here (cf. George Marsaglia, Normal (Gaussian) random variables for supercomputers The Journal of Supercomputing, Volume 5, Number 1, 49–55 DOI: 10.1007/BF00155857). First we draw uniform x and y in [0,1] using mpz_urandomb (in fixed precision), and scale them to [-1, 1]. */ mpz_init (xp); mpz_init (yp); mpz_init (x); mpz_init (y); mpz_init (t); mpz_init (s); mpz_init (a); mpz_init (b); mpfr_init2 (sfr, MPFR_PREC_MIN); mpfr_init2 (l, MPFR_PREC_MIN); mpfr_init2 (r1, MPFR_PREC_MIN); if (rop2 != NULL) mpfr_init2 (r2, MPFR_PREC_MIN); mpz_set_ui (xp, 0); mpz_set_ui (yp, 0); for (;;) { tprec = tprec0; do { mpz_urandomb (xp, rstate, tprec); mpz_urandomb (yp, rstate, tprec); mpz_mul (a, xp, xp); mpz_mul (b, yp, yp); mpz_add (s, a, b); } while (mpz_sizeinbase (s, 2) > tprec * 2); /* x^2 + y^2 <= 2^{2tprec} */ for (;;) { /* FIXME: compute s as s += 2x + 2y + 2 */ mpz_add_ui (a, xp, 1); mpz_add_ui (b, yp, 1); mpz_mul (a, a, a); mpz_mul (b, b, b); mpz_add (s, a, b); if ((mpz_sizeinbase (s, 2) <= 2 * tprec) || ((mpz_sizeinbase (s, 2) == 2 * tprec + 1) && (mpz_scan1 (s, 0) == 2 * tprec))) goto yeepee; /* Extend by 32 bits */ mpz_mul_2exp (xp, xp, 32); mpz_mul_2exp (yp, yp, 32); mpz_urandomb (x, rstate, 32); mpz_urandomb (y, rstate, 32); mpz_add (xp, xp, x); mpz_add (yp, yp, y); tprec += 32; mpz_mul (a, xp, xp); mpz_mul (b, yp, yp); mpz_add (s, a, b); if (mpz_sizeinbase (s, 2) > tprec * 2) break; } } yeepee: /* FIXME: compute s with s -= 2x + 2y + 2 */ mpz_mul (a, xp, xp); mpz_mul (b, yp, yp); mpz_add (s, a, b); /* Compute the signs of the output */ mpz_urandomb (x, rstate, 2); s1 = mpz_tstbit (x, 0); s2 = mpz_tstbit (x, 1); for (;;) { /* s = xp^2 + yp^2 (loop invariant) */ mpfr_set_prec (sfr, 2 * tprec); mpfr_set_prec (l, tprec); mpfr_set_z (sfr, s, MPFR_RNDN); /* exact */ mpfr_mul_2si (sfr, sfr, -2 * tprec, MPFR_RNDN); /* exact */ mpfr_log (l, sfr, MPFR_RNDN); mpfr_neg (l, l, MPFR_RNDN); mpfr_mul_2si (l, l, 1, MPFR_RNDN); mpfr_div (l, l, sfr, MPFR_RNDN); mpfr_sqrt (l, l, MPFR_RNDN); mpfr_set_prec (r1, tprec); mpfr_mul_z (r1, l, xp, MPFR_RNDN); mpfr_div_2ui (r1, r1, tprec, MPFR_RNDN); /* exact */ if (s1) mpfr_neg (r1, r1, MPFR_RNDN); if (MPFR_CAN_ROUND (r1, tprec - 2, MPFR_PREC (rop1), rnd)) { if (rop2 != NULL) { mpfr_set_prec (r2, tprec); mpfr_mul_z (r2, l, yp, MPFR_RNDN); mpfr_div_2ui (r2, r2, tprec, MPFR_RNDN); /* exact */ if (s2) mpfr_neg (r2, r2, MPFR_RNDN); if (MPFR_CAN_ROUND (r2, tprec - 2, MPFR_PREC (rop2), rnd)) break; } else break; } /* Extend by 32 bits */ mpz_mul_2exp (xp, xp, 32); mpz_mul_2exp (yp, yp, 32); mpz_urandomb (x, rstate, 32); mpz_urandomb (y, rstate, 32); mpz_add (xp, xp, x); mpz_add (yp, yp, y); tprec += 32; mpz_mul (a, xp, xp); mpz_mul (b, yp, yp); mpz_add (s, a, b); } inex1 = mpfr_set (rop1, r1, rnd); if (rop2 != NULL) { inex2 = mpfr_set (rop2, r2, rnd); inex2 = mpfr_check_range (rop2, inex2, rnd); } inex1 = mpfr_check_range (rop1, inex1, rnd); if (rop2 != NULL) mpfr_clear (r2); mpfr_clear (r1); mpfr_clear (l); mpfr_clear (sfr); mpz_clear (b); mpz_clear (a); mpz_clear (s); mpz_clear (t); mpz_clear (y); mpz_clear (x); mpz_clear (yp); mpz_clear (xp); return INEX (inex1, inex2); }
int mpfr_sub1sp (mpfr_ptr a, mpfr_srcptr b, mpfr_srcptr c, mpfr_rnd_t rnd_mode) { mpfr_exp_t bx,cx; mpfr_uexp_t d; mpfr_prec_t p, sh, cnt; mp_size_t n; mp_limb_t *ap, *bp, *cp; mp_limb_t limb; int inexact; mp_limb_t bcp,bcp1; /* Cp and C'p+1 */ mp_limb_t bbcp = (mp_limb_t) -1, bbcp1 = (mp_limb_t) -1; /* Cp+1 and C'p+2, gcc claims that they might be used uninitialized. We fill them with invalid values, which should produce a failure if so. See README.dev file. */ MPFR_TMP_DECL(marker); MPFR_TMP_MARK(marker); MPFR_ASSERTD(MPFR_PREC(a) == MPFR_PREC(b) && MPFR_PREC(b) == MPFR_PREC(c)); MPFR_ASSERTD(MPFR_IS_PURE_FP(b)); MPFR_ASSERTD(MPFR_IS_PURE_FP(c)); /* Read prec and num of limbs */ p = MPFR_PREC(b); n = (p-1)/GMP_NUMB_BITS+1; /* Fast cmp of |b| and |c|*/ bx = MPFR_GET_EXP (b); cx = MPFR_GET_EXP (c); if (MPFR_UNLIKELY(bx == cx)) { mp_size_t k = n - 1; /* Check mantissa since exponent are equals */ bp = MPFR_MANT(b); cp = MPFR_MANT(c); while (k>=0 && MPFR_UNLIKELY(bp[k] == cp[k])) k--; if (MPFR_UNLIKELY(k < 0)) /* b == c ! */ { /* Return exact number 0 */ if (rnd_mode == MPFR_RNDD) MPFR_SET_NEG(a); else MPFR_SET_POS(a); MPFR_SET_ZERO(a); MPFR_RET(0); } else if (bp[k] > cp[k]) goto BGreater; else { MPFR_ASSERTD(bp[k]<cp[k]); goto CGreater; } } else if (MPFR_UNLIKELY(bx < cx)) { /* Swap b and c and set sign */ mpfr_srcptr t; mpfr_exp_t tx; CGreater: MPFR_SET_OPPOSITE_SIGN(a,b); t = b; b = c; c = t; tx = bx; bx = cx; cx = tx; } else { /* b > c */ BGreater: MPFR_SET_SAME_SIGN(a,b); } /* Now b > c */ MPFR_ASSERTD(bx >= cx); d = (mpfr_uexp_t) bx - cx; DEBUG (printf ("New with diff=%lu\n", (unsigned long) d)); if (MPFR_UNLIKELY(d <= 1)) { if (MPFR_LIKELY(d < 1)) { /* <-- b --> <-- c --> : exact sub */ ap = MPFR_MANT(a); mpn_sub_n (ap, MPFR_MANT(b), MPFR_MANT(c), n); /* Normalize */ ExactNormalize: limb = ap[n-1]; if (MPFR_LIKELY(limb)) { /* First limb is not zero. */ count_leading_zeros(cnt, limb); /* cnt could be == 0 <= SubD1Lose */ if (MPFR_LIKELY(cnt)) { mpn_lshift(ap, ap, n, cnt); /* Normalize number */ bx -= cnt; /* Update final expo */ } /* Last limb should be ok */ MPFR_ASSERTD(!(ap[0] & MPFR_LIMB_MASK((unsigned int) (-p) % GMP_NUMB_BITS))); } else { /* First limb is zero */ mp_size_t k = n-1, len; /* Find the first limb not equal to zero. FIXME:It is assume it exists (since |b| > |c| and same prec)*/ do { MPFR_ASSERTD( k > 0 ); limb = ap[--k]; } while (limb == 0); MPFR_ASSERTD(limb != 0); count_leading_zeros(cnt, limb); k++; len = n - k; /* Number of last limb */ MPFR_ASSERTD(k >= 0); if (MPFR_LIKELY(cnt)) mpn_lshift(ap+len, ap, k, cnt); /* Normalize the High Limb*/ else { /* Must use DECR since src and dest may overlap & dest>=src*/ MPN_COPY_DECR(ap+len, ap, k); } MPN_ZERO(ap, len); /* Zeroing the last limbs */ bx -= cnt + len*GMP_NUMB_BITS; /* Update Expo */ /* Last limb should be ok */ MPFR_ASSERTD(!(ap[len]&MPFR_LIMB_MASK((unsigned int) (-p) % GMP_NUMB_BITS))); } /* Check expo underflow */ if (MPFR_UNLIKELY(bx < __gmpfr_emin)) { MPFR_TMP_FREE(marker); /* inexact=0 */ DEBUG( printf("(D==0 Underflow)\n") ); if (rnd_mode == MPFR_RNDN && (bx < __gmpfr_emin - 1 || (/*inexact >= 0 &&*/ mpfr_powerof2_raw (a)))) rnd_mode = MPFR_RNDZ; return mpfr_underflow (a, rnd_mode, MPFR_SIGN(a)); } MPFR_SET_EXP (a, bx); /* No rounding is necessary since the result is exact */ MPFR_ASSERTD(ap[n-1] > ~ap[n-1]); MPFR_TMP_FREE(marker); return 0; } else /* if (d == 1) */ { /* | <-- b --> | <-- c --> */ mp_limb_t c0, mask; mp_size_t k; MPFR_UNSIGNED_MINUS_MODULO(sh, p); /* If we lose at least one bit, compute 2*b-c (Exact) * else compute b-c/2 */ bp = MPFR_MANT(b); cp = MPFR_MANT(c); k = n-1; limb = bp[k] - cp[k]/2; if (limb > MPFR_LIMB_HIGHBIT) { /* We can't lose precision: compute b-c/2 */ /* Shift c in the allocated temporary block */ SubD1NoLose: c0 = cp[0] & (MPFR_LIMB_ONE<<sh); cp = MPFR_TMP_LIMBS_ALLOC (n); mpn_rshift(cp, MPFR_MANT(c), n, 1); if (MPFR_LIKELY(c0 == 0)) { /* Result is exact: no need of rounding! */ ap = MPFR_MANT(a); mpn_sub_n (ap, bp, cp, n); MPFR_SET_EXP(a, bx); /* No expo overflow! */ /* No truncate or normalize is needed */ MPFR_ASSERTD(ap[n-1] > ~ap[n-1]); /* No rounding is necessary since the result is exact */ MPFR_TMP_FREE(marker); return 0; } ap = MPFR_MANT(a); mask = ~MPFR_LIMB_MASK(sh); cp[0] &= mask; /* Delete last bit of c */ mpn_sub_n (ap, bp, cp, n); MPFR_SET_EXP(a, bx); /* No expo overflow! */ MPFR_ASSERTD( !(ap[0] & ~mask) ); /* Check last bits */ /* No normalize is needed */ MPFR_ASSERTD(ap[n-1] > ~ap[n-1]); /* Rounding is necessary since c0 = 1*/ /* Cp =-1 and C'p+1=0 */ bcp = 1; bcp1 = 0; if (MPFR_LIKELY(rnd_mode == MPFR_RNDN)) { /* Even Rule apply: Check Ap-1 */ if (MPFR_LIKELY( (ap[0] & (MPFR_LIMB_ONE<<sh)) == 0) ) goto truncate; else goto sub_one_ulp; } MPFR_UPDATE_RND_MODE(rnd_mode, MPFR_IS_NEG(a)); if (rnd_mode == MPFR_RNDZ) goto sub_one_ulp; else goto truncate; } else if (MPFR_LIKELY(limb < MPFR_LIMB_HIGHBIT)) { /* We lose at least one bit of prec */ /* Calcul of 2*b-c (Exact) */ /* Shift b in the allocated temporary block */ SubD1Lose: bp = MPFR_TMP_LIMBS_ALLOC (n); mpn_lshift (bp, MPFR_MANT(b), n, 1); ap = MPFR_MANT(a); mpn_sub_n (ap, bp, cp, n); bx--; goto ExactNormalize; } else { /* Case: limb = 100000000000 */ /* Check while b[k] == c'[k] (C' is C shifted by 1) */ /* If b[k]<c'[k] => We lose at least one bit*/ /* If b[k]>c'[k] => We don't lose any bit */ /* If k==-1 => We don't lose any bit AND the result is 100000000000 0000000000 00000000000 */ mp_limb_t carry; do { carry = cp[k]&MPFR_LIMB_ONE; k--; } while (k>=0 && bp[k]==(carry=cp[k]/2+(carry<<(GMP_NUMB_BITS-1)))); if (MPFR_UNLIKELY(k<0)) { /*If carry then (sh==0 and Virtual c'[-1] > Virtual b[-1]) */ if (MPFR_UNLIKELY(carry)) /* carry = cp[0]&MPFR_LIMB_ONE */ { /* FIXME: Can be faster? */ MPFR_ASSERTD(sh == 0); goto SubD1Lose; } /* Result is a power of 2 */ ap = MPFR_MANT (a); MPN_ZERO (ap, n); ap[n-1] = MPFR_LIMB_HIGHBIT; MPFR_SET_EXP (a, bx); /* No expo overflow! */ /* No Normalize is needed*/ /* No Rounding is needed */ MPFR_TMP_FREE (marker); return 0; } /* carry = cp[k]/2+(cp[k-1]&1)<<(GMP_NUMB_BITS-1) = c'[k]*/ else if (bp[k] > carry) goto SubD1NoLose; else { MPFR_ASSERTD(bp[k]<carry); goto SubD1Lose; } } } } else if (MPFR_UNLIKELY(d >= p)) { ap = MPFR_MANT(a); MPFR_UNSIGNED_MINUS_MODULO(sh, p); /* We can't set A before since we use cp for rounding... */ /* Perform rounding: check if a=b or a=b-ulp(b) */ if (MPFR_UNLIKELY(d == p)) { /* cp == -1 and c'p+1 = ? */ bcp = 1; /* We need Cp+1 later for a very improbable case. */ bbcp = (MPFR_MANT(c)[n-1] & (MPFR_LIMB_ONE<<(GMP_NUMB_BITS-2))); /* We need also C'p+1 for an even more unprobable case... */ if (MPFR_LIKELY( bbcp )) bcp1 = 1; else { cp = MPFR_MANT(c); if (MPFR_UNLIKELY(cp[n-1] == MPFR_LIMB_HIGHBIT)) { mp_size_t k = n-1; do { k--; } while (k>=0 && cp[k]==0); bcp1 = (k>=0); } else bcp1 = 1; } DEBUG( printf("(D=P) Cp=-1 Cp+1=%d C'p+1=%d \n", bbcp!=0, bcp1!=0) ); bp = MPFR_MANT (b); /* Even if src and dest overlap, it is ok using MPN_COPY */ if (MPFR_LIKELY(rnd_mode == MPFR_RNDN)) { if (MPFR_UNLIKELY( bcp && bcp1==0 )) /* Cp=-1 and C'p+1=0: Even rule Apply! */ /* Check Ap-1 = Bp-1 */ if ((bp[0] & (MPFR_LIMB_ONE<<sh)) == 0) { MPN_COPY(ap, bp, n); goto truncate; } MPN_COPY(ap, bp, n); goto sub_one_ulp; } MPFR_UPDATE_RND_MODE(rnd_mode, MPFR_IS_NEG(a)); if (rnd_mode == MPFR_RNDZ) { MPN_COPY(ap, bp, n); goto sub_one_ulp; } else { MPN_COPY(ap, bp, n); goto truncate; } } else { /* Cp=0, Cp+1=-1 if d==p+1, C'p+1=-1 */ bcp = 0; bbcp = (d==p+1); bcp1 = 1; DEBUG( printf("(D>P) Cp=%d Cp+1=%d C'p+1=%d\n", bcp!=0,bbcp!=0,bcp1!=0) ); /* Need to compute C'p+2 if d==p+1 and if rnd_mode=NEAREST (Because of a very improbable case) */ if (MPFR_UNLIKELY(d==p+1 && rnd_mode==MPFR_RNDN)) { cp = MPFR_MANT(c); if (MPFR_UNLIKELY(cp[n-1] == MPFR_LIMB_HIGHBIT)) { mp_size_t k = n-1; do { k--; } while (k>=0 && cp[k]==0); bbcp1 = (k>=0); } else bbcp1 = 1; DEBUG( printf("(D>P) C'p+2=%d\n", bbcp1!=0) ); } /* Copy mantissa B in A */ MPN_COPY(ap, MPFR_MANT(b), n); /* Round */ if (MPFR_LIKELY(rnd_mode == MPFR_RNDN)) goto truncate; MPFR_UPDATE_RND_MODE(rnd_mode, MPFR_IS_NEG(a)); if (rnd_mode == MPFR_RNDZ) goto sub_one_ulp; else /* rnd_mode = AWAY */ goto truncate; } } else { mpfr_uexp_t dm; mp_size_t m; mp_limb_t mask; /* General case: 2 <= d < p */ MPFR_UNSIGNED_MINUS_MODULO(sh, p); cp = MPFR_TMP_LIMBS_ALLOC (n); /* Shift c in temporary allocated place */ dm = d % GMP_NUMB_BITS; m = d / GMP_NUMB_BITS; if (MPFR_UNLIKELY(dm == 0)) { /* dm = 0 and m > 0: Just copy */ MPFR_ASSERTD(m!=0); MPN_COPY(cp, MPFR_MANT(c)+m, n-m); MPN_ZERO(cp+n-m, m); } else if (MPFR_LIKELY(m == 0)) { /* dm >=2 and m == 0: just shift */ MPFR_ASSERTD(dm >= 2); mpn_rshift(cp, MPFR_MANT(c), n, dm); } else { /* dm > 0 and m > 0: shift and zero */ mpn_rshift(cp, MPFR_MANT(c)+m, n-m, dm); MPN_ZERO(cp+n-m, m); } DEBUG( mpfr_print_mant_binary("Before", MPFR_MANT(c), p) ); DEBUG( mpfr_print_mant_binary("B= ", MPFR_MANT(b), p) ); DEBUG( mpfr_print_mant_binary("After ", cp, p) ); /* Compute bcp=Cp and bcp1=C'p+1 */ if (MPFR_LIKELY(sh)) { /* Try to compute them from C' rather than C (FIXME: Faster?) */ bcp = (cp[0] & (MPFR_LIMB_ONE<<(sh-1))) ; if (MPFR_LIKELY( cp[0] & MPFR_LIMB_MASK(sh-1) )) bcp1 = 1; else { /* We can't compute C'p+1 from C'. Compute it from C */ /* Start from bit x=p-d+sh in mantissa C (+sh since we have already looked sh bits in C'!) */ mpfr_prec_t x = p-d+sh-1; if (MPFR_LIKELY(x>p)) /* We are already looked at all the bits of c, so C'p+1 = 0*/ bcp1 = 0; else { mp_limb_t *tp = MPFR_MANT(c); mp_size_t kx = n-1 - (x / GMP_NUMB_BITS); mpfr_prec_t sx = GMP_NUMB_BITS-1-(x%GMP_NUMB_BITS); DEBUG (printf ("(First) x=%lu Kx=%ld Sx=%lu\n", (unsigned long) x, (long) kx, (unsigned long) sx)); /* Looks at the last bits of limb kx (if sx=0 does nothing)*/ if (tp[kx] & MPFR_LIMB_MASK(sx)) bcp1 = 1; else { /*kx += (sx==0);*/ /*If sx==0, tp[kx] hasn't been checked*/ do { kx--; } while (kx>=0 && tp[kx]==0); bcp1 = (kx >= 0); } } } } else { /* Compute Cp and C'p+1 from C with sh=0 */ mp_limb_t *tp = MPFR_MANT(c); /* Start from bit x=p-d in mantissa C */ mpfr_prec_t x = p-d; mp_size_t kx = n-1 - (x / GMP_NUMB_BITS); mpfr_prec_t sx = GMP_NUMB_BITS-1-(x%GMP_NUMB_BITS); MPFR_ASSERTD(p >= d); bcp = (tp[kx] & (MPFR_LIMB_ONE<<sx)); /* Looks at the last bits of limb kx (If sx=0, does nothing)*/ if (tp[kx] & MPFR_LIMB_MASK(sx)) bcp1 = 1; else { /*kx += (sx==0);*/ /*If sx==0, tp[kx] hasn't been checked*/ do { kx--; } while (kx>=0 && tp[kx]==0); bcp1 = (kx>=0); } } DEBUG( printf("sh=%lu Cp=%d C'p+1=%d\n", sh, bcp!=0, bcp1!=0) ); /* Check if we can lose a bit, and if so compute Cp+1 and C'p+2 */ bp = MPFR_MANT(b); if (MPFR_UNLIKELY((bp[n-1]-cp[n-1]) <= MPFR_LIMB_HIGHBIT)) { /* We can lose a bit so we precompute Cp+1 and C'p+2 */ /* Test for trivial case: since C'p+1=0, Cp+1=0 and C'p+2 =0 */ if (MPFR_LIKELY(bcp1 == 0)) { bbcp = 0; bbcp1 = 0; } else /* bcp1 != 0 */ { /* We can lose a bit: compute Cp+1 and C'p+2 from mantissa C */ mp_limb_t *tp = MPFR_MANT(c); /* Start from bit x=(p+1)-d in mantissa C */ mpfr_prec_t x = p+1-d; mp_size_t kx = n-1 - (x/GMP_NUMB_BITS); mpfr_prec_t sx = GMP_NUMB_BITS-1-(x%GMP_NUMB_BITS); MPFR_ASSERTD(p > d); DEBUG (printf ("(pre) x=%lu Kx=%ld Sx=%lu\n", (unsigned long) x, (long) kx, (unsigned long) sx)); bbcp = (tp[kx] & (MPFR_LIMB_ONE<<sx)) ; /* Looks at the last bits of limb kx (If sx=0, does nothing)*/ /* If Cp+1=0, since C'p+1!=0, C'p+2=1 ! */ if (MPFR_LIKELY(bbcp==0 || (tp[kx]&MPFR_LIMB_MASK(sx)))) bbcp1 = 1; else { /*kx += (sx==0);*/ /*If sx==0, tp[kx] hasn't been checked*/ do { kx--; } while (kx>=0 && tp[kx]==0); bbcp1 = (kx>=0); DEBUG (printf ("(Pre) Scan done for %ld\n", (long) kx)); } } /*End of Bcp1 != 0*/ DEBUG( printf("(Pre) Cp+1=%d C'p+2=%d\n", bbcp!=0, bbcp1!=0) ); } /* End of "can lose a bit" */ /* Clean shifted C' */ mask = ~MPFR_LIMB_MASK (sh); cp[0] &= mask; /* Substract the mantissa c from b in a */ ap = MPFR_MANT(a); mpn_sub_n (ap, bp, cp, n); DEBUG( mpfr_print_mant_binary("Sub= ", ap, p) ); /* Normalize: we lose at max one bit*/ if (MPFR_UNLIKELY(MPFR_LIMB_MSB(ap[n-1]) == 0)) { /* High bit is not set and we have to fix it! */ /* Ap >= 010000xxx001 */ mpn_lshift(ap, ap, n, 1); /* Ap >= 100000xxx010 */ if (MPFR_UNLIKELY(bcp!=0)) /* Check if Cp = -1 */ /* Since Cp == -1, we have to substract one more */ { mpn_sub_1(ap, ap, n, MPFR_LIMB_ONE<<sh); MPFR_ASSERTD(MPFR_LIMB_MSB(ap[n-1]) != 0); } /* Ap >= 10000xxx001 */ /* Final exponent -1 since we have shifted the mantissa */ bx--; /* Update bcp and bcp1 */ MPFR_ASSERTN(bbcp != (mp_limb_t) -1); MPFR_ASSERTN(bbcp1 != (mp_limb_t) -1); bcp = bbcp; bcp1 = bbcp1; /* We dont't have anymore a valid Cp+1! But since Ap >= 100000xxx001, the final sub can't unnormalize!*/ } MPFR_ASSERTD( !(ap[0] & ~mask) ); /* Rounding */ if (MPFR_LIKELY(rnd_mode == MPFR_RNDN)) { if (MPFR_LIKELY(bcp==0)) goto truncate; else if ((bcp1) || ((ap[0] & (MPFR_LIMB_ONE<<sh)) != 0)) goto sub_one_ulp; else goto truncate; } /* Update rounding mode */ MPFR_UPDATE_RND_MODE(rnd_mode, MPFR_IS_NEG(a)); if (rnd_mode == MPFR_RNDZ && (MPFR_LIKELY(bcp || bcp1))) goto sub_one_ulp; goto truncate; } MPFR_RET_NEVER_GO_HERE (); /* Sub one ulp to the result */ sub_one_ulp: mpn_sub_1 (ap, ap, n, MPFR_LIMB_ONE << sh); /* Result should be smaller than exact value: inexact=-1 */ inexact = -1; /* Check normalisation */ if (MPFR_UNLIKELY(MPFR_LIMB_MSB(ap[n-1]) == 0)) { /* ap was a power of 2, and we lose a bit */ /* Now it is 0111111111111111111[00000 */ mpn_lshift(ap, ap, n, 1); bx--; /* And the lost bit x depends on Cp+1, and Cp */ /* Compute Cp+1 if it isn't already compute (ie d==1) */ /* FIXME: Is this case possible? */ if (MPFR_UNLIKELY(d == 1)) bbcp = 0; DEBUG( printf("(SubOneUlp)Cp=%d, Cp+1=%d C'p+1=%d\n", bcp!=0,bbcp!=0,bcp1!=0)); /* Compute the last bit (Since we have shifted the mantissa) we need one more bit!*/ MPFR_ASSERTN(bbcp != (mp_limb_t) -1); if ( (rnd_mode == MPFR_RNDZ && bcp==0) || (rnd_mode==MPFR_RNDN && bbcp==0) || (bcp && bcp1==0) ) /*Exact result*/ { ap[0] |= MPFR_LIMB_ONE<<sh; if (rnd_mode == MPFR_RNDN) inexact = 1; DEBUG( printf("(SubOneUlp) Last bit set\n") ); } /* Result could be exact if C'p+1 = 0 and rnd == Zero since we have had one more bit to the result */ /* Fixme: rnd_mode == MPFR_RNDZ needed ? */ if (bcp1==0 && rnd_mode==MPFR_RNDZ) { DEBUG( printf("(SubOneUlp) Exact result\n") ); inexact = 0; } } goto end_of_sub; truncate: /* Check if the result is an exact power of 2: 100000000000 in which cases, we could have to do sub_one_ulp due to some nasty reasons: If Result is a Power of 2: + If rnd = AWAY, | If Cp=-1 and C'p+1 = 0, SubOneUlp and the result is EXACT. If Cp=-1 and C'p+1 =-1, SubOneUlp and the result is above. Otherwise truncate + If rnd = NEAREST, If Cp= 0 and Cp+1 =-1 and C'p+2=-1, SubOneUlp and the result is above If cp=-1 and C'p+1 = 0, SubOneUlp and the result is exact. Otherwise truncate. X bit should always be set if SubOneUlp*/ if (MPFR_UNLIKELY(ap[n-1] == MPFR_LIMB_HIGHBIT)) { mp_size_t k = n-1; do { k--; } while (k>=0 && ap[k]==0); if (MPFR_UNLIKELY(k<0)) { /* It is a power of 2! */ /* Compute Cp+1 if it isn't already compute (ie d==1) */ /* FIXME: Is this case possible? */ if (d == 1) bbcp=0; DEBUG( printf("(Truncate) Cp=%d, Cp+1=%d C'p+1=%d C'p+2=%d\n", \ bcp!=0, bbcp!=0, bcp1!=0, bbcp1!=0) ); MPFR_ASSERTN(bbcp != (mp_limb_t) -1); MPFR_ASSERTN((rnd_mode != MPFR_RNDN) || (bcp != 0) || (bbcp == 0) || (bbcp1 != (mp_limb_t) -1)); if (((rnd_mode != MPFR_RNDZ) && bcp) || ((rnd_mode == MPFR_RNDN) && (bcp == 0) && (bbcp) && (bbcp1))) { DEBUG( printf("(Truncate) Do sub\n") ); mpn_sub_1 (ap, ap, n, MPFR_LIMB_ONE << sh); mpn_lshift(ap, ap, n, 1); ap[0] |= MPFR_LIMB_ONE<<sh; bx--; /* FIXME: Explain why it works (or why not)... */ inexact = (bcp1 == 0) ? 0 : (rnd_mode==MPFR_RNDN) ? -1 : 1; goto end_of_sub; } } } /* Calcul of Inexact flag.*/ inexact = MPFR_LIKELY(bcp || bcp1) ? 1 : 0; end_of_sub: /* Update Expo */ /* FIXME: Is this test really useful? If d==0 : Exact case. This is never called. if 1 < d < p : bx=MPFR_EXP(b) or MPFR_EXP(b)-1 > MPFR_EXP(c) > emin if d == 1 : bx=MPFR_EXP(b). If we could lose any bits, the exact normalisation is called. if d >= p : bx=MPFR_EXP(b) >= MPFR_EXP(c) + p > emin After SubOneUlp, we could have one bit less. if 1 < d < p : bx >= MPFR_EXP(b)-2 >= MPFR_EXP(c) > emin if d == 1 : bx >= MPFR_EXP(b)-1 = MPFR_EXP(c) > emin. if d >= p : bx >= MPFR_EXP(b)-1 > emin since p>=2. */ MPFR_ASSERTD( bx >= __gmpfr_emin); /* if (MPFR_UNLIKELY(bx < __gmpfr_emin)) { DEBUG( printf("(Final Underflow)\n") ); if (rnd_mode == MPFR_RNDN && (bx < __gmpfr_emin - 1 || (inexact >= 0 && mpfr_powerof2_raw (a)))) rnd_mode = MPFR_RNDZ; MPFR_TMP_FREE(marker); return mpfr_underflow (a, rnd_mode, MPFR_SIGN(a)); } */ MPFR_SET_EXP (a, bx); MPFR_TMP_FREE(marker); MPFR_RET (inexact * MPFR_INT_SIGN (a)); }
int mpfr_erf (mpfr_ptr y, mpfr_srcptr x, mpfr_rnd_t rnd_mode) { mpfr_t xf; int inex, large; MPFR_SAVE_EXPO_DECL (expo); MPFR_LOG_FUNC (("x[%Pu]=%.*Rg rnd=%d", mpfr_get_prec (x), mpfr_log_prec, x, rnd_mode), ("y[%Pu]=%.*Rg inexact=%d", mpfr_get_prec (y), mpfr_log_prec, y, inex)); if (MPFR_UNLIKELY (MPFR_IS_SINGULAR (x))) { if (MPFR_IS_NAN (x)) { MPFR_SET_NAN (y); MPFR_RET_NAN; } else if (MPFR_IS_INF (x)) /* erf(+inf) = +1, erf(-inf) = -1 */ return mpfr_set_si (y, MPFR_INT_SIGN (x), MPFR_RNDN); else /* erf(+0) = +0, erf(-0) = -0 */ { MPFR_ASSERTD (MPFR_IS_ZERO (x)); return mpfr_set (y, x, MPFR_RNDN); /* should keep the sign of x */ } } /* now x is neither NaN, Inf nor 0 */ /* first try expansion at x=0 when x is small, or asymptotic expansion where x is large */ MPFR_SAVE_EXPO_MARK (expo); /* around x=0, we have erf(x) = 2x/sqrt(Pi) (1 - x^2/3 + ...), with 1 - x^2/3 <= sqrt(Pi)*erf(x)/2/x <= 1 for x >= 0. This means that if x^2/3 < 2^(-PREC(y)-1) we can decide of the correct rounding, unless we have a worst-case for 2x/sqrt(Pi). */ if (MPFR_EXP(x) < - (mpfr_exp_t) (MPFR_PREC(y) / 2)) { /* we use 2x/sqrt(Pi) (1 - x^2/3) <= erf(x) <= 2x/sqrt(Pi) for x > 0 and 2x/sqrt(Pi) <= erf(x) <= 2x/sqrt(Pi) (1 - x^2/3) for x < 0. In both cases |2x/sqrt(Pi) (1 - x^2/3)| <= |erf(x)| <= |2x/sqrt(Pi)|. We will compute l and h such that l <= |2x/sqrt(Pi) (1 - x^2/3)| and |2x/sqrt(Pi)| <= h. If l and h round to the same value to precision PREC(y) and rounding rnd_mode, then we are done. */ mpfr_t l, h; /* lower and upper bounds for erf(x) */ int ok, inex2; mpfr_init2 (l, MPFR_PREC(y) + 17); mpfr_init2 (h, MPFR_PREC(y) + 17); /* first compute l */ mpfr_mul (l, x, x, MPFR_RNDU); mpfr_div_ui (l, l, 3, MPFR_RNDU); /* upper bound on x^2/3 */ mpfr_ui_sub (l, 1, l, MPFR_RNDZ); /* lower bound on 1 - x^2/3 */ mpfr_const_pi (h, MPFR_RNDU); /* upper bound of Pi */ mpfr_sqrt (h, h, MPFR_RNDU); /* upper bound on sqrt(Pi) */ mpfr_div (l, l, h, MPFR_RNDZ); /* lower bound on 1/sqrt(Pi) (1 - x^2/3) */ mpfr_mul_2ui (l, l, 1, MPFR_RNDZ); /* 2/sqrt(Pi) (1 - x^2/3) */ mpfr_mul (l, l, x, MPFR_RNDZ); /* |l| is a lower bound on |2x/sqrt(Pi) (1 - x^2/3)| */ /* now compute h */ mpfr_const_pi (h, MPFR_RNDD); /* lower bound on Pi */ mpfr_sqrt (h, h, MPFR_RNDD); /* lower bound on sqrt(Pi) */ mpfr_div_2ui (h, h, 1, MPFR_RNDD); /* lower bound on sqrt(Pi)/2 */ /* since sqrt(Pi)/2 < 1, the following should not underflow */ mpfr_div (h, x, h, MPFR_IS_POS(x) ? MPFR_RNDU : MPFR_RNDD); /* round l and h to precision PREC(y) */ inex = mpfr_prec_round (l, MPFR_PREC(y), rnd_mode); inex2 = mpfr_prec_round (h, MPFR_PREC(y), rnd_mode); /* Caution: we also need inex=inex2 (inex might be 0). */ ok = SAME_SIGN (inex, inex2) && mpfr_cmp (l, h) == 0; if (ok) mpfr_set (y, h, rnd_mode); mpfr_clear (l); mpfr_clear (h); if (ok) goto end; /* this test can still fail for small precision, for example for x=-0.100E-2 with a target precision of 3 bits, since the error term x^2/3 is not that small. */ } mpfr_init2 (xf, 53); mpfr_const_log2 (xf, MPFR_RNDU); mpfr_div (xf, x, xf, MPFR_RNDZ); /* round to zero ensures we get a lower bound of |x/log(2)| */ mpfr_mul (xf, xf, x, MPFR_RNDZ); large = mpfr_cmp_ui (xf, MPFR_PREC (y) + 1) > 0; mpfr_clear (xf); /* when x goes to infinity, we have erf(x) = 1 - 1/sqrt(Pi)/exp(x^2)/x + ... and |erf(x) - 1| <= exp(-x^2) is true for any x >= 0, thus if exp(-x^2) < 2^(-PREC(y)-1) the result is 1 or 1-epsilon. This rewrites as x^2/log(2) > p+1. */ if (MPFR_UNLIKELY (large)) /* |erf x| = 1 or 1- */ { mpfr_rnd_t rnd2 = MPFR_IS_POS (x) ? rnd_mode : MPFR_INVERT_RND(rnd_mode); if (rnd2 == MPFR_RNDN || rnd2 == MPFR_RNDU || rnd2 == MPFR_RNDA) { inex = MPFR_INT_SIGN (x); mpfr_set_si (y, inex, rnd2); } else /* round to zero */ { inex = -MPFR_INT_SIGN (x); mpfr_setmax (y, 0); /* warning: setmax keeps the old sign of y */ MPFR_SET_SAME_SIGN (y, x); } } else /* use Taylor */ { double xf2; /* FIXME: get rid of doubles/mpfr_get_d here */ xf2 = mpfr_get_d (x, MPFR_RNDN); xf2 = xf2 * xf2; /* xf2 ~ x^2 */ inex = mpfr_erf_0 (y, x, xf2, rnd_mode); } end: MPFR_SAVE_EXPO_FREE (expo); return mpfr_check_range (y, inex, rnd_mode); }
int mpfr_frac (mpfr_ptr r, mpfr_srcptr u, mp_rnd_t rnd_mode) { mp_exp_t re, ue; mp_prec_t uq, fq; mp_size_t un, tn, t0; mp_limb_t *up, *tp, k; int sh; mpfr_t tmp; mpfr_ptr t; /* Special cases */ if (MPFR_UNLIKELY(MPFR_IS_NAN(u))) { MPFR_SET_NAN(r); MPFR_RET_NAN; } else if (MPFR_UNLIKELY(MPFR_IS_INF(u) || mpfr_integer_p (u))) { MPFR_CLEAR_FLAGS(r); MPFR_SET_SAME_SIGN(r, u); MPFR_SET_ZERO(r); MPFR_RET(0); /* zero is exact */ } ue = MPFR_GET_EXP (u); if (ue <= 0) /* |u| < 1 */ return mpfr_set (r, u, rnd_mode); uq = MPFR_PREC(u); un = (uq - 1) / BITS_PER_MP_LIMB; /* index of most significant limb */ un -= (mp_size_t) (ue / BITS_PER_MP_LIMB); /* now the index of the MSL containing bits of the fractional part */ up = MPFR_MANT(u); sh = ue % BITS_PER_MP_LIMB; k = up[un] << sh; /* the first bit of the fractional part is the MSB of k */ if (k != 0) { int cnt; count_leading_zeros(cnt, k); /* first bit 1 of the fractional part -> MSB of the number */ re = -cnt; sh += cnt; MPFR_ASSERTN (sh < BITS_PER_MP_LIMB); k <<= cnt; } else { re = sh - BITS_PER_MP_LIMB; /* searching for the first bit 1 (exists since u isn't an integer) */ while (up[--un] == 0) re -= BITS_PER_MP_LIMB; MPFR_ASSERTN(un >= 0); k = up[un]; count_leading_zeros(sh, k); re -= sh; k <<= sh; } /* The exponent of r will be re */ /* un: index of the limb of u that contains the first bit 1 of the FP */ ue -= re; /* number of bits of u to discard */ fq = uq - ue; /* number of bits of the fractional part of u */ /* Temporary fix */ t = /* fq > MPFR_PREC(r) */ (mp_size_t) (MPFR_PREC(r) - 1) / BITS_PER_MP_LIMB < un ? (mpfr_init2 (tmp, (un + 1) * BITS_PER_MP_LIMB), tmp) : r; /* t has enough precision to contain the fractional part of u */ /* If we use a temporary variable, we take the non-significant bits of u into account, because of the mpn_lshift below. */ MPFR_CLEAR_FLAGS(t); MPFR_SET_SAME_SIGN(t, u); MPFR_SET_EXP (t, re); /* Put the fractional part of u into t */ tn = (MPFR_PREC(t) - 1) / BITS_PER_MP_LIMB; MPFR_ASSERTN(tn >= un); t0 = tn - un; tp = MPFR_MANT(t); if (sh == 0) MPN_COPY_DECR(tp + t0, up, un + 1); else /* warning: un may be 0 here */ tp[tn] = k | ((un) ? mpn_lshift (tp + t0, up, un, sh) : (mp_limb_t) 0); if (t0 > 0) MPN_ZERO(tp, t0); if (t != r) { /* t is tmp */ int inex; inex = mpfr_set (r, t, rnd_mode); mpfr_clear (t); return inex; } else MPFR_RET(0); }
int mpfr_rint (mpfr_ptr r, mpfr_srcptr u, mpfr_rnd_t rnd_mode) { int sign; int rnd_away; mpfr_exp_t exp; if (MPFR_UNLIKELY( MPFR_IS_SINGULAR(u) )) { if (MPFR_IS_NAN(u)) { MPFR_SET_NAN(r); MPFR_RET_NAN; } MPFR_SET_SAME_SIGN(r, u); if (MPFR_IS_INF(u)) { MPFR_SET_INF(r); MPFR_RET(0); /* infinity is exact */ } else /* now u is zero */ { MPFR_ASSERTD(MPFR_IS_ZERO(u)); MPFR_SET_ZERO(r); MPFR_RET(0); /* zero is exact */ } } MPFR_SET_SAME_SIGN (r, u); /* Does nothing if r==u */ sign = MPFR_INT_SIGN (u); exp = MPFR_GET_EXP (u); rnd_away = rnd_mode == MPFR_RNDD ? sign < 0 : rnd_mode == MPFR_RNDU ? sign > 0 : rnd_mode == MPFR_RNDZ ? 0 : rnd_mode == MPFR_RNDA ? 1 : -1; /* round to nearest-even (RNDN) or nearest-away (RNDNA) */ /* rnd_away: 1 if round away from zero, 0 if round to zero, -1 if not decided yet. */ if (MPFR_UNLIKELY (exp <= 0)) /* 0 < |u| < 1 ==> round |u| to 0 or 1 */ { /* Note: in the MPFR_RNDN mode, 0.5 must be rounded to 0. */ if (rnd_away != 0 && (rnd_away > 0 || (exp == 0 && (rnd_mode == MPFR_RNDNA || !mpfr_powerof2_raw (u))))) { mp_limb_t *rp; mp_size_t rm; rp = MPFR_MANT(r); rm = (MPFR_PREC(r) - 1) / GMP_NUMB_BITS; rp[rm] = MPFR_LIMB_HIGHBIT; MPN_ZERO(rp, rm); MPFR_SET_EXP (r, 1); /* |r| = 1 */ MPFR_RET(sign > 0 ? 2 : -2); } else { MPFR_SET_ZERO(r); /* r = 0 */ MPFR_RET(sign > 0 ? -2 : 2); } } else /* exp > 0, |u| >= 1 */ { mp_limb_t *up, *rp; mp_size_t un, rn, ui; int sh, idiff; int uflags; /* * uflags will contain: * _ 0 if u is an integer representable in r, * _ 1 if u is an integer not representable in r, * _ 2 if u is not an integer. */ up = MPFR_MANT(u); rp = MPFR_MANT(r); un = MPFR_LIMB_SIZE(u); rn = MPFR_LIMB_SIZE(r); MPFR_UNSIGNED_MINUS_MODULO (sh, MPFR_PREC (r)); MPFR_SET_EXP (r, exp); /* Does nothing if r==u */ if ((exp - 1) / GMP_NUMB_BITS >= un) { ui = un; idiff = 0; uflags = 0; /* u is an integer, representable or not in r */ } else { mp_size_t uj; ui = (exp - 1) / GMP_NUMB_BITS + 1; /* #limbs of the int part */ MPFR_ASSERTD (un >= ui); uj = un - ui; /* lowest limb of the integer part */ idiff = exp % GMP_NUMB_BITS; /* #int-part bits in up[uj] or 0 */ uflags = idiff == 0 || (up[uj] << idiff) == 0 ? 0 : 2; if (uflags == 0) while (uj > 0) if (up[--uj] != 0) { uflags = 2; break; } } if (ui > rn) { /* More limbs in the integer part of u than in r. Just round u with the precision of r. */ MPFR_ASSERTD (rp != up && un > rn); MPN_COPY (rp, up + (un - rn), rn); /* r != u */ if (rnd_away < 0) { /* This is a rounding to nearest mode (MPFR_RNDN or MPFR_RNDNA). Decide the rounding direction here. */ if (rnd_mode == MPFR_RNDN && (rp[0] & (MPFR_LIMB_ONE << sh)) == 0) { /* halfway cases rounded toward zero */ mp_limb_t a, b; /* a: rounding bit and some of the following bits */ /* b: boundary for a (weight of the rounding bit in a) */ if (sh != 0) { a = rp[0] & ((MPFR_LIMB_ONE << sh) - 1); b = MPFR_LIMB_ONE << (sh - 1); } else { a = up[un - rn - 1]; b = MPFR_LIMB_HIGHBIT; } rnd_away = a > b; if (a == b) { mp_size_t i; for (i = un - rn - 1 - (sh == 0); i >= 0; i--) if (up[i] != 0) { rnd_away = 1; break; } } } else /* halfway cases rounded away from zero */ rnd_away = /* rounding bit */ ((sh != 0 && (rp[0] & (MPFR_LIMB_ONE << (sh - 1))) != 0) || (sh == 0 && (up[un - rn - 1] & MPFR_LIMB_HIGHBIT) != 0)); } if (uflags == 0) { /* u is an integer; determine if it is representable in r */ if (sh != 0 && rp[0] << (GMP_NUMB_BITS - sh) != 0) uflags = 1; /* u is not representable in r */ else { mp_size_t i; for (i = un - rn - 1; i >= 0; i--) if (up[i] != 0) { uflags = 1; /* u is not representable in r */ break; } } } } else /* ui <= rn */ { mp_size_t uj, rj; int ush; uj = un - ui; /* lowest limb of the integer part in u */ rj = rn - ui; /* lowest limb of the integer part in r */ if (MPFR_LIKELY (rp != up)) MPN_COPY(rp + rj, up + uj, ui); /* Ignore the lowest rj limbs, all equal to zero. */ rp += rj; rn = ui; /* number of fractional bits in whole rp[0] */ ush = idiff == 0 ? 0 : GMP_NUMB_BITS - idiff; if (rj == 0 && ush < sh) { /* If u is an integer (uflags == 0), we need to determine if it is representable in r, i.e. if its sh - ush bits in the non-significant part of r are all 0. */ if (uflags == 0 && (rp[0] & ((MPFR_LIMB_ONE << sh) - (MPFR_LIMB_ONE << ush))) != 0) uflags = 1; /* u is an integer not representable in r */ } else /* The integer part of u fits in r, we'll round to it. */ sh = ush; if (rnd_away < 0) { /* This is a rounding to nearest mode. Decide the rounding direction here. */ if (uj == 0 && sh == 0) rnd_away = 0; /* rounding bit = 0 (not represented in u) */ else if (rnd_mode == MPFR_RNDN && (rp[0] & (MPFR_LIMB_ONE << sh)) == 0) { /* halfway cases rounded toward zero */ mp_limb_t a, b; /* a: rounding bit and some of the following bits */ /* b: boundary for a (weight of the rounding bit in a) */ if (sh != 0) { a = rp[0] & ((MPFR_LIMB_ONE << sh) - 1); b = MPFR_LIMB_ONE << (sh - 1); } else { MPFR_ASSERTD (uj >= 1); /* see above */ a = up[uj - 1]; b = MPFR_LIMB_HIGHBIT; } rnd_away = a > b; if (a == b) { mp_size_t i; for (i = uj - 1 - (sh == 0); i >= 0; i--) if (up[i] != 0) { rnd_away = 1; break; } } } else /* halfway cases rounded away from zero */ rnd_away = /* rounding bit */ ((sh != 0 && (rp[0] & (MPFR_LIMB_ONE << (sh - 1))) != 0) || (sh == 0 && (MPFR_ASSERTD (uj >= 1), up[uj - 1] & MPFR_LIMB_HIGHBIT) != 0)); } /* Now we can make the low rj limbs to 0 */ MPN_ZERO (rp-rj, rj); } if (sh != 0) rp[0] &= MP_LIMB_T_MAX << sh; /* If u is a representable integer, there is no rounding. */ if (uflags == 0) MPFR_RET(0); MPFR_ASSERTD (rnd_away >= 0); /* rounding direction is defined */ if (rnd_away && mpn_add_1(rp, rp, rn, MPFR_LIMB_ONE << sh)) { if (exp == __gmpfr_emax) return mpfr_overflow(r, rnd_mode, MPFR_SIGN(r)) >= 0 ? uflags : -uflags; else { MPFR_SET_EXP(r, exp + 1); rp[rn-1] = MPFR_LIMB_HIGHBIT; } } MPFR_RET (rnd_away ^ (sign < 0) ? uflags : -uflags); } /* exp > 0, |u| >= 1 */ }
int mpfr_cos (mpfr_ptr y, mpfr_srcptr x, mp_rnd_t rnd_mode) { int K0, K, precy, m, k, l; int inexact; mpfr_t r, s; mp_limb_t *rp, *sp; mp_size_t sm; mp_exp_t exps, cancel = 0; TMP_DECL (marker); if (MPFR_UNLIKELY(MPFR_IS_SINGULAR(x))) { if (MPFR_IS_NAN(x) || MPFR_IS_INF(x)) { MPFR_SET_NAN(y); MPFR_RET_NAN; } else { MPFR_ASSERTD(MPFR_IS_ZERO(x)); return mpfr_set_ui (y, 1, GMP_RNDN); } } mpfr_save_emin_emax (); precy = MPFR_PREC(y); K0 = __gmpfr_isqrt(precy / 2); /* Need K + log2(precy/K) extra bits */ m = precy + 3 * (K0 + 2 * MAX(MPFR_GET_EXP (x), 0)) + 3; TMP_MARK(marker); sm = (m + BITS_PER_MP_LIMB - 1) / BITS_PER_MP_LIMB; MPFR_TMP_INIT(rp, r, m, sm); MPFR_TMP_INIT(sp, s, m, sm); for (;;) { mpfr_mul (r, x, x, GMP_RNDU); /* err <= 1 ulp */ /* we need that |r| < 1 for mpfr_cos2_aux, i.e. up(x^2)/2^(2K) < 1 */ K = K0 + MAX (MPFR_GET_EXP (r), 0); mpfr_div_2ui (r, r, 2 * K, GMP_RNDN); /* r = (x/2^K)^2, err <= 1 ulp */ /* s <- 1 - r/2! + ... + (-1)^l r^l/(2l)! */ l = mpfr_cos2_aux (s, r); MPFR_SET_ONE (r); for (k = 0; k < K; k++) { mpfr_mul (s, s, s, GMP_RNDU); /* err <= 2*olderr */ mpfr_mul_2ui (s, s, 1, GMP_RNDU); /* err <= 4*olderr */ mpfr_sub (s, s, r, GMP_RNDN); } /* absolute error on s is bounded by (2l+1/3)*2^(2K-m) */ for (k = 2 * K, l = 2 * l + 1; l > 1; l = (l + 1) >> 1) k++; /* now the error is bounded by 2^(k-m) = 2^(EXP(s)-err) */ exps = MPFR_GET_EXP(s); if (MPFR_LIKELY(mpfr_can_round (s, exps + m - k, GMP_RNDN, GMP_RNDZ, precy + (rnd_mode == GMP_RNDN)))) break; m += BITS_PER_MP_LIMB; if (exps < cancel) { m += cancel - exps; cancel = exps; } sm = (m + BITS_PER_MP_LIMB - 1) / BITS_PER_MP_LIMB; MPFR_TMP_INIT(rp, r, m, sm); MPFR_TMP_INIT(sp, s, m, sm); } mpfr_restore_emin_emax (); inexact = mpfr_set (y, s, rnd_mode); /* FIXME: Dont' need check range? */ TMP_FREE(marker); return inexact; }
int mpfr_sinh_cosh (mpfr_ptr sh, mpfr_ptr ch, mpfr_srcptr xt, mpfr_rnd_t rnd_mode) { mpfr_t x; int inexact_sh, inexact_ch; MPFR_ASSERTN (sh != ch); MPFR_LOG_FUNC (("x[%Pu]=%.*Rg rnd=%d", mpfr_get_prec (xt), mpfr_log_prec, xt, rnd_mode), ("sh[%Pu]=%.*Rg ch[%Pu]=%.*Rg", mpfr_get_prec (sh), mpfr_log_prec, sh, mpfr_get_prec (ch), mpfr_log_prec, ch)); if (MPFR_UNLIKELY (MPFR_IS_SINGULAR (xt))) { if (MPFR_IS_NAN (xt)) { MPFR_SET_NAN (ch); MPFR_SET_NAN (sh); MPFR_RET_NAN; } else if (MPFR_IS_INF (xt)) { MPFR_SET_INF (sh); MPFR_SET_SAME_SIGN (sh, xt); MPFR_SET_INF (ch); MPFR_SET_POS (ch); MPFR_RET (0); } else /* xt is zero */ { MPFR_ASSERTD (MPFR_IS_ZERO (xt)); MPFR_SET_ZERO (sh); /* sinh(0) = 0 */ MPFR_SET_SAME_SIGN (sh, xt); inexact_sh = 0; inexact_ch = mpfr_set_ui (ch, 1, rnd_mode); /* cosh(0) = 1 */ return INEX(inexact_sh,inexact_ch); } } /* Warning: if we use MPFR_FAST_COMPUTE_IF_SMALL_INPUT here, make sure that the code also works in case of overlap (see sin_cos.c) */ MPFR_TMP_INIT_ABS (x, xt); { mpfr_t s, c, ti; mpfr_exp_t d; mpfr_prec_t N; /* Precision of the intermediary variables */ long int err; /* Precision of error */ MPFR_ZIV_DECL (loop); MPFR_SAVE_EXPO_DECL (expo); MPFR_GROUP_DECL (group); MPFR_SAVE_EXPO_MARK (expo); /* compute the precision of intermediary variable */ N = MPFR_PREC (ch); N = MAX (N, MPFR_PREC (sh)); /* the optimal number of bits : see algorithms.ps */ N = N + MPFR_INT_CEIL_LOG2 (N) + 4; /* initialise of intermediary variables */ MPFR_GROUP_INIT_3 (group, N, s, c, ti); /* First computation of sinh_cosh */ MPFR_ZIV_INIT (loop, N); for (;;) { MPFR_BLOCK_DECL (flags); /* compute sinh_cosh */ MPFR_BLOCK (flags, mpfr_exp (s, x, MPFR_RNDD)); if (MPFR_OVERFLOW (flags)) /* exp(x) does overflow */ { /* since cosh(x) >= exp(x), cosh(x) overflows too */ inexact_ch = mpfr_overflow (ch, rnd_mode, MPFR_SIGN_POS); /* sinh(x) may be representable */ inexact_sh = mpfr_sinh (sh, xt, rnd_mode); MPFR_SAVE_EXPO_UPDATE_FLAGS (expo, MPFR_FLAGS_OVERFLOW); break; } d = MPFR_GET_EXP (s); mpfr_ui_div (ti, 1, s, MPFR_RNDU); /* 1/exp(x) */ mpfr_add (c, s, ti, MPFR_RNDU); /* exp(x) + 1/exp(x) */ mpfr_sub (s, s, ti, MPFR_RNDN); /* exp(x) - 1/exp(x) */ mpfr_div_2ui (c, c, 1, MPFR_RNDN); /* 1/2(exp(x) + 1/exp(x)) */ mpfr_div_2ui (s, s, 1, MPFR_RNDN); /* 1/2(exp(x) - 1/exp(x)) */ /* it may be that s is zero (in fact, it can only occur when exp(x)=1, and thus ti=1 too) */ if (MPFR_IS_ZERO (s)) err = N; /* double the precision */ else { /* calculation of the error */ d = d - MPFR_GET_EXP (s) + 2; /* error estimate: err = N-(__gmpfr_ceil_log2(1+pow(2,d)));*/ err = N - (MAX (d, 0) + 1); if (MPFR_LIKELY (MPFR_CAN_ROUND (s, err, MPFR_PREC (sh), rnd_mode) && \ MPFR_CAN_ROUND (c, err, MPFR_PREC (ch), rnd_mode))) { inexact_sh = mpfr_set4 (sh, s, rnd_mode, MPFR_SIGN (xt)); inexact_ch = mpfr_set (ch, c, rnd_mode); break; } } /* actualisation of the precision */ N += err; MPFR_ZIV_NEXT (loop, N); MPFR_GROUP_REPREC_3 (group, N, s, c, ti); } MPFR_ZIV_FREE (loop); MPFR_GROUP_CLEAR (group); MPFR_SAVE_EXPO_FREE (expo); } /* now, let's raise the flags if needed */ inexact_sh = mpfr_check_range (sh, inexact_sh, rnd_mode); inexact_ch = mpfr_check_range (ch, inexact_ch, rnd_mode); return INEX(inexact_sh,inexact_ch); }
int mpfr_yn (mpfr_ptr res, long n, mpfr_srcptr z, mp_rnd_t r) { int inex; unsigned long absn; mp_prec_t prec; mp_exp_t err1, err2, err3; mpfr_t y, s1, s2, s3; MPFR_ZIV_DECL (loop); MPFR_LOG_FUNC (("x[%#R]=%R n=%d rnd=%d", z, z, n, r), ("y[%#R]=%R", res, res)); absn = SAFE_ABS (unsigned long, n); if (MPFR_UNLIKELY (MPFR_IS_SINGULAR (z))) { if (MPFR_IS_NAN (z)) { MPFR_SET_NAN (res); /* y(n,NaN) = NaN */ MPFR_RET_NAN; } /* y(n,z) tends to zero when z goes to +Inf, oscillating around 0. We choose to return +0 in that case. */ else if (MPFR_IS_INF (z)) { if (MPFR_SIGN(z) > 0) return mpfr_set_ui (res, 0, r); else /* y(n,-Inf) = NaN */ { MPFR_SET_NAN (res); MPFR_RET_NAN; } } else /* y(n,z) tends to -Inf for n >= 0 or n even, to +Inf otherwise, when z goes to zero */ { MPFR_SET_INF(res); if (n >= 0 || (n & 1) == 0) MPFR_SET_NEG(res); else MPFR_SET_POS(res); MPFR_RET(0); } } /* for z < 0, y(n,z) is imaginary except when j(n,|z|) = 0, which we assume does not happen for a rational z. */ if (MPFR_SIGN(z) < 0) { MPFR_SET_NAN (res); MPFR_RET_NAN; } /* now z is not singular, and z > 0 */ /* Deal with tiny arguments. We have: y0(z) = 2 log(z)/Pi + 2 (euler - log(2))/Pi + O(log(z)*z^2), more precisely for 0 <= z <= 1/2, with g(z) = 2/Pi + 2(euler-log(2))/Pi/log(z), g(z) - 0.41*z^2 < y0(z)/log(z) < g(z) thus since log(z) is negative: g(z)*log(z) < y0(z) < (g(z) - z^2/2)*log(z) and since |g(z)| >= 0.63 for 0 <= z <= 1/2, the relative error on y0(z)/log(z) is bounded by 0.41*z^2/0.63 <= 0.66*z^2. Note: we use both the main term in log(z) and the constant term, because otherwise the relative error would be only in 1/log(|log(z)|). */ if (n == 0 && MPFR_EXP(z) < - (mp_exp_t) (MPFR_PREC(res) / 2)) { mpfr_t l, h, t, logz; int ok, inex2; prec = MPFR_PREC(res) + 10; mpfr_init2 (l, prec); mpfr_init2 (h, prec); mpfr_init2 (t, prec); mpfr_init2 (logz, prec); /* first enclose log(z) + euler - log(2) = log(z/2) + euler */ mpfr_log (logz, z, GMP_RNDD); /* lower bound of log(z) */ mpfr_set (h, logz, GMP_RNDU); /* exact */ mpfr_nextabove (h); /* upper bound of log(z) */ mpfr_const_euler (t, GMP_RNDD); /* lower bound of euler */ mpfr_add (l, logz, t, GMP_RNDD); /* lower bound of log(z) + euler */ mpfr_nextabove (t); /* upper bound of euler */ mpfr_add (h, h, t, GMP_RNDU); /* upper bound of log(z) + euler */ mpfr_const_log2 (t, GMP_RNDU); /* upper bound of log(2) */ mpfr_sub (l, l, t, GMP_RNDD); /* lower bound of log(z/2) + euler */ mpfr_nextbelow (t); /* lower bound of log(2) */ mpfr_sub (h, h, t, GMP_RNDU); /* upper bound of log(z/2) + euler */ mpfr_const_pi (t, GMP_RNDU); /* upper bound of Pi */ mpfr_div (l, l, t, GMP_RNDD); /* lower bound of (log(z/2)+euler)/Pi */ mpfr_nextbelow (t); /* lower bound of Pi */ mpfr_div (h, h, t, GMP_RNDD); /* upper bound of (log(z/2)+euler)/Pi */ mpfr_mul_2ui (l, l, 1, GMP_RNDD); /* lower bound on g(z)*log(z) */ mpfr_mul_2ui (h, h, 1, GMP_RNDU); /* upper bound on g(z)*log(z) */ /* we now have l <= g(z)*log(z) <= h, and we need to add -z^2/2*log(z) to h */ mpfr_mul (t, z, z, GMP_RNDU); /* upper bound on z^2 */ /* since logz is negative, a lower bound corresponds to an upper bound for its absolute value */ mpfr_neg (t, t, GMP_RNDD); mpfr_div_2ui (t, t, 1, GMP_RNDD); mpfr_mul (t, t, logz, GMP_RNDU); /* upper bound on z^2/2*log(z) */ /* an underflow may happen in the above instructions, clear flag */ mpfr_clear_underflow (); mpfr_add (h, h, t, GMP_RNDU); inex = mpfr_prec_round (l, MPFR_PREC(res), r); inex2 = mpfr_prec_round (h, MPFR_PREC(res), r); /* we need h=l and inex=inex2 */ ok = (inex == inex2) && (mpfr_cmp (l, h) == 0); if (ok) mpfr_set (res, h, r); /* exact */ mpfr_clear (l); mpfr_clear (h); mpfr_clear (t); mpfr_clear (logz); if (ok) return inex; } /* small argument check for y1(z) = -2/Pi/z + O(log(z)): for 0 <= z <= 1, |y1(z) + 2/Pi/z| <= 0.25 */ if (n == 1 && MPFR_EXP(z) + 1 < - (mp_exp_t) MPFR_PREC(res)) { mpfr_t y; int ok; /* since 2/Pi > 0.5, and |y1(z)| >= |2/Pi/z|, if z <= 2^(-emax-1), then |y1(z)| > 2^emax */ prec = MPFR_PREC(res) + 10; mpfr_init2 (y, prec); mpfr_const_pi (y, GMP_RNDU); /* Pi*(1+u)^2, where here and below u represents a quantity <= 1/2^prec */ mpfr_mul (y, y, z, GMP_RNDU); /* Pi*z * (1+u)^4, upper bound */ mpfr_ui_div (y, 2, y, GMP_RNDZ); /* 2/Pi/z * (1+u)^6, lower bound */ mpfr_neg (y, y, GMP_RNDN); if (mpfr_overflow_p ()) { mpfr_clear (y); return mpfr_overflow (res, r, -1); } /* (1+u)^6 can be written 1+7u [for another value of u], thus the error on 2/Pi/z is less than 7ulp(y). The truncation error is less than 1/4, thus if ulp(y)>=1/4, the total error is less than 8ulp(y), otherwise it is less than 1/4+7/8 <= 2. */ if (MPFR_EXP(y) + 2 >= MPFR_PREC(y)) /* ulp(y) >= 1/4 */ err1 = 3; else /* ulp(y) <= 1/8 */ err1 = (mp_exp_t) MPFR_PREC(y) - MPFR_EXP(y) + 1; ok = MPFR_CAN_ROUND (y, prec - err1, MPFR_PREC(res), r); if (ok) inex = mpfr_set (res, y, r); mpfr_clear (y); if (ok) return inex; } /* we can use the asymptotic expansion as soon as z > p log(2)/2, but to get some margin we use it for z > p/2 */ if (mpfr_cmp_ui (z, MPFR_PREC(res) / 2 + 3) > 0) { inex = mpfr_yn_asympt (res, n, z, r); if (inex != 0) return inex; } mpfr_init (y); mpfr_init (s1); mpfr_init (s2); mpfr_init (s3); prec = MPFR_PREC(res) + 2 * MPFR_INT_CEIL_LOG2 (MPFR_PREC (res)) + 13; MPFR_ZIV_INIT (loop, prec); for (;;) { mpfr_set_prec (y, prec); mpfr_set_prec (s1, prec); mpfr_set_prec (s2, prec); mpfr_set_prec (s3, prec); mpfr_mul (y, z, z, GMP_RNDN); mpfr_div_2ui (y, y, 2, GMP_RNDN); /* z^2/4 */ /* store (z/2)^n temporarily in s2 */ mpfr_pow_ui (s2, z, absn, GMP_RNDN); mpfr_div_2si (s2, s2, absn, GMP_RNDN); /* compute S1 * (z/2)^(-n) */ if (n == 0) { mpfr_set_ui (s1, 0, GMP_RNDN); err1 = 0; } else err1 = mpfr_yn_s1 (s1, y, absn - 1); mpfr_div (s1, s1, s2, GMP_RNDN); /* (z/2)^(-n) * S1 */ /* See algorithms.tex: the relative error on s1 is bounded by (3n+3)*2^(e+1-prec). */ err1 = MPFR_INT_CEIL_LOG2 (3 * absn + 3) + err1 + 1; /* rel_err(s1) <= 2^(err1-prec), thus err(s1) <= 2^err1 ulps */ /* compute (z/2)^n * S3 */ mpfr_neg (y, y, GMP_RNDN); /* -z^2/4 */ err3 = mpfr_yn_s3 (s3, y, s2, absn); /* (z/2)^n * S3 */ /* the error on s3 is bounded by 2^err3 ulps */ /* add s1+s3 */ err1 += MPFR_EXP(s1); mpfr_add (s1, s1, s3, GMP_RNDN); /* the error is bounded by 1/2 + 2^err1*2^(- EXP(s1)) + 2^err3*2^(EXP(s3) - EXP(s1)) */ err3 += MPFR_EXP(s3); err1 = (err3 > err1) ? err3 + 1 : err1 + 1; err1 -= MPFR_EXP(s1); err1 = (err1 >= 0) ? err1 + 1 : 1; /* now the error on s1 is bounded by 2^err1*ulp(s1) */ /* compute S2 */ mpfr_div_2ui (s2, z, 1, GMP_RNDN); /* z/2 */ mpfr_log (s2, s2, GMP_RNDN); /* log(z/2) */ mpfr_const_euler (s3, GMP_RNDN); err2 = MPFR_EXP(s2) > MPFR_EXP(s3) ? MPFR_EXP(s2) : MPFR_EXP(s3); mpfr_add (s2, s2, s3, GMP_RNDN); /* log(z/2) + gamma */ err2 -= MPFR_EXP(s2); mpfr_mul_2ui (s2, s2, 1, GMP_RNDN); /* 2*(log(z/2) + gamma) */ mpfr_jn (s3, absn, z, GMP_RNDN); /* Jn(z) */ mpfr_mul (s2, s2, s3, GMP_RNDN); /* 2*(log(z/2) + gamma)*Jn(z) */ err2 += 4; /* the error on s2 is bounded by 2^err2 ulps, see algorithms.tex */ /* add all three sums */ err1 += MPFR_EXP(s1); /* the error on s1 is bounded by 2^err1 */ err2 += MPFR_EXP(s2); /* the error on s2 is bounded by 2^err2 */ mpfr_sub (s2, s2, s1, GMP_RNDN); /* s2 - (s1+s3) */ err2 = (err1 > err2) ? err1 + 1 : err2 + 1; err2 -= MPFR_EXP(s2); err2 = (err2 >= 0) ? err2 + 1 : 1; /* now the error on s2 is bounded by 2^err2*ulp(s2) */ mpfr_const_pi (y, GMP_RNDN); /* error bounded by 1 ulp */ mpfr_div (s2, s2, y, GMP_RNDN); /* error bounded by 2^(err2+1)*ulp(s2) */ err2 ++; if (MPFR_LIKELY (MPFR_CAN_ROUND (s2, prec - err2, MPFR_PREC(res), r))) break; MPFR_ZIV_NEXT (loop, prec); } MPFR_ZIV_FREE (loop); inex = (n >= 0 || (n & 1) == 0) ? mpfr_set (res, s2, r) : mpfr_neg (res, s2, r); mpfr_clear (y); mpfr_clear (s1); mpfr_clear (s2); mpfr_clear (s3); return inex; }
int mpfr_zeta_ui (mpfr_ptr z, unsigned long m, mpfr_rnd_t r) { MPFR_ZIV_DECL (loop); if (m == 0) { mpfr_set_ui (z, 1, r); mpfr_div_2ui (z, z, 1, r); MPFR_CHANGE_SIGN (z); MPFR_RET (0); } else if (m == 1) { MPFR_SET_INF (z); MPFR_SET_POS (z); mpfr_set_divby0 (); return 0; } else /* m >= 2 */ { mpfr_prec_t p = MPFR_PREC(z); unsigned long n, k, err, kbits; mpz_t d, t, s, q; mpfr_t y; int inex; if (r == MPFR_RNDA) r = MPFR_RNDU; /* since the result is always positive */ if (m >= p) /* 2^(-m) < ulp(1) = 2^(1-p). This means that 2^(-m) <= 1/2*ulp(1). We have 3^(-m)+4^(-m)+... < 2^(-m) i.e. zeta(m) < 1+2*2^(-m) for m >= 3 */ { if (m == 2) /* necessarily p=2 */ return mpfr_set_ui_2exp (z, 13, -3, r); else if (r == MPFR_RNDZ || r == MPFR_RNDD || (r == MPFR_RNDN && m > p)) { mpfr_set_ui (z, 1, r); return -1; } else { mpfr_set_ui (z, 1, r); mpfr_nextabove (z); return 1; } } /* now treat also the case where zeta(m) - (1+1/2^m) < 1/2*ulp(1), and the result is either 1+2^(-m) or 1+2^(-m)+2^(1-p). */ mpfr_init2 (y, 31); if (m >= p / 2) /* otherwise 4^(-m) > 2^(-p) */ { /* the following is a lower bound for log(3)/log(2) */ mpfr_set_str_binary (y, "1.100101011100000000011010001110"); mpfr_mul_ui (y, y, m, MPFR_RNDZ); /* lower bound for log2(3^m) */ if (mpfr_cmp_ui (y, p + 2) >= 0) { mpfr_clear (y); mpfr_set_ui (z, 1, MPFR_RNDZ); mpfr_div_2ui (z, z, m, MPFR_RNDZ); mpfr_add_ui (z, z, 1, MPFR_RNDZ); if (r != MPFR_RNDU) return -1; mpfr_nextabove (z); return 1; } } mpz_init (s); mpz_init (d); mpz_init (t); mpz_init (q); p += MPFR_INT_CEIL_LOG2(p); /* account of the n term in the error */ p += MPFR_INT_CEIL_LOG2(p) + 15; /* initial value */ MPFR_ZIV_INIT (loop, p); for(;;) { /* 0.39321985067869744 = log(2)/log(3+sqrt(8)) */ n = 1 + (unsigned long) (0.39321985067869744 * (double) p); err = n + 4; mpfr_set_prec (y, p); /* computation of the d[k] */ mpz_set_ui (s, 0); mpz_set_ui (t, 1); mpz_mul_2exp (t, t, 2 * n - 1); /* t[n] */ mpz_set (d, t); for (k = n; k > 0; k--) { count_leading_zeros (kbits, k); kbits = GMP_NUMB_BITS - kbits; /* if k^m is too large, use mpz_tdiv_q */ if (m * kbits > 2 * GMP_NUMB_BITS) { /* if we know in advance that k^m > d, then floor(d/k^m) will be zero below, so there is no need to compute k^m */ kbits = (kbits - 1) * m + 1; /* k^m has at least kbits bits */ if (kbits > mpz_sizeinbase (d, 2)) mpz_set_ui (q, 0); else { mpz_ui_pow_ui (q, k, m); mpz_tdiv_q (q, d, q); } } else /* use several mpz_tdiv_q_ui calls */ { unsigned long km = k, mm = m - 1; while (mm > 0 && km < ULONG_MAX / k) { km *= k; mm --; } mpz_tdiv_q_ui (q, d, km); while (mm > 0) { km = k; mm --; while (mm > 0 && km < ULONG_MAX / k) { km *= k; mm --; } mpz_tdiv_q_ui (q, q, km); } } if (k % 2) mpz_add (s, s, q); else mpz_sub (s, s, q); /* we have d[k] = sum(t[i], i=k+1..n) with t[i] = n*(n+i-1)!*4^i/(n-i)!/(2i)! t[k-1]/t[k] = k*(2k-1)/(n-k+1)/(n+k-1)/2 */ #if (GMP_NUMB_BITS == 32) #define KMAX 46341 /* max k such that k*(2k-1) < 2^32 */ #elif (GMP_NUMB_BITS == 64) #define KMAX 3037000500 #endif #ifdef KMAX if (k <= KMAX) mpz_mul_ui (t, t, k * (2 * k - 1)); else #endif { mpz_mul_ui (t, t, k); mpz_mul_ui (t, t, 2 * k - 1); } mpz_fdiv_q_2exp (t, t, 1); /* Warning: the test below assumes that an unsigned long has no padding bits. */ if (n < 1UL << ((sizeof(unsigned long) * CHAR_BIT) / 2)) /* (n - k + 1) * (n + k - 1) < n^2 */ mpz_divexact_ui (t, t, (n - k + 1) * (n + k - 1)); else { mpz_divexact_ui (t, t, n - k + 1); mpz_divexact_ui (t, t, n + k - 1); } mpz_add (d, d, t); } /* multiply by 1/(1-2^(1-m)) = 1 + 2^(1-m) + 2^(2-m) + ... */ mpz_fdiv_q_2exp (t, s, m - 1); do { err ++; mpz_add (s, s, t); mpz_fdiv_q_2exp (t, t, m - 1); } while (mpz_cmp_ui (t, 0) > 0); /* divide by d[n] */ mpz_mul_2exp (s, s, p); mpz_tdiv_q (s, s, d); mpfr_set_z (y, s, MPFR_RNDN); mpfr_div_2ui (y, y, p, MPFR_RNDN); err = MPFR_INT_CEIL_LOG2 (err); if (MPFR_LIKELY(MPFR_CAN_ROUND (y, p - err, MPFR_PREC(z), r))) break; MPFR_ZIV_NEXT (loop, p); } MPFR_ZIV_FREE (loop); mpz_clear (d); mpz_clear (t); mpz_clear (q); mpz_clear (s); inex = mpfr_set (z, y, r); mpfr_clear (y); return inex; } }
/* compute in s an approximation of S3 = c*sum((h(k)+h(n+k))*y^k/k!/(n+k)!,k=0..infinity) where h(k) = 1 + 1/2 + ... + 1/k k=0: h(n) k=1: 1+h(n+1) k=2: 3/2+h(n+2) Returns e such that the error is bounded by 2^e ulp(s). */ static mp_exp_t mpfr_yn_s3 (mpfr_ptr s, mpfr_srcptr y, mpfr_srcptr c, unsigned long n) { unsigned long k, zz; mpfr_t t, u; mpz_t p, q; /* p/q will store h(k)+h(n+k) */ mp_exp_t exps, expU; zz = mpfr_get_ui (y, GMP_RNDU); /* y = z^2/4 */ MPFR_ASSERTN (zz < ULONG_MAX - 2); zz += 2; /* z^2 <= 2^zz */ mpz_init_set_ui (p, 0); mpz_init_set_ui (q, 1); /* initialize p/q to h(n) */ for (k = 1; k <= n; k++) { /* p/q + 1/k = (k*p+q)/(q*k) */ mpz_mul_ui (p, p, k); mpz_add (p, p, q); mpz_mul_ui (q, q, k); } mpfr_init2 (t, MPFR_PREC(s)); mpfr_init2 (u, MPFR_PREC(s)); mpfr_fac_ui (t, n, GMP_RNDN); mpfr_div (t, c, t, GMP_RNDN); /* c/n! */ mpfr_mul_z (u, t, p, GMP_RNDN); mpfr_div_z (s, u, q, GMP_RNDN); exps = MPFR_EXP (s); expU = exps; for (k = 1; ;k ++) { /* update t */ mpfr_mul (t, t, y, GMP_RNDN); mpfr_div_ui (t, t, k, GMP_RNDN); mpfr_div_ui (t, t, n + k, GMP_RNDN); /* update p/q: p/q + 1/k + 1/(n+k) = [p*k*(n+k) + q*(n+k) + q*k]/(q*k*(n+k)) */ mpz_mul_ui (p, p, k); mpz_mul_ui (p, p, n + k); mpz_addmul_ui (p, q, n + 2 * k); mpz_mul_ui (q, q, k); mpz_mul_ui (q, q, n + k); mpfr_mul_z (u, t, p, GMP_RNDN); mpfr_div_z (u, u, q, GMP_RNDN); exps = MPFR_EXP (u); if (exps > expU) expU = exps; mpfr_add (s, s, u, GMP_RNDN); exps = MPFR_EXP (s); if (exps > expU) expU = exps; if (MPFR_EXP (u) + (mp_exp_t) MPFR_PREC (u) < MPFR_EXP (s) && zz / (2 * k) < k + n) break; } mpfr_clear (t); mpfr_clear (u); mpz_clear (p); mpz_clear (q); exps = expU - MPFR_EXP (s); /* the error is bounded by (6k^2+33/2k+11) 2^exps ulps <= 8*(k+2)^2 2^exps ulps */ return 3 + 2 * MPFR_INT_CEIL_LOG2(k + 2) + exps; }
/* Set iop to the integral part of op and fop to its fractional part */ int mpfr_modf (mpfr_ptr iop, mpfr_ptr fop, mpfr_srcptr op, mpfr_rnd_t rnd_mode) { mpfr_exp_t ope; mpfr_prec_t opq; int inexi, inexf; MPFR_LOG_FUNC (("op[%#R]=%R rnd=%d", op, op, rnd_mode), ("iop[%#R]=%R fop[%#R]=%R", iop, iop, fop, fop)); MPFR_ASSERTN (iop != fop); if ( MPFR_UNLIKELY (MPFR_IS_SINGULAR (op)) ) { if (MPFR_IS_NAN (op)) { MPFR_SET_NAN (iop); MPFR_SET_NAN (fop); MPFR_RET_NAN; } MPFR_SET_SAME_SIGN (iop, op); MPFR_SET_SAME_SIGN (fop, op); if (MPFR_IS_INF (op)) { MPFR_SET_INF (iop); MPFR_SET_ZERO (fop); MPFR_RET (0); } else /* op is zero */ { MPFR_ASSERTD (MPFR_IS_ZERO (op)); MPFR_SET_ZERO (iop); MPFR_SET_ZERO (fop); MPFR_RET (0); } } ope = MPFR_GET_EXP (op); opq = MPFR_PREC (op); if (ope <= 0) /* 0 < |op| < 1 */ { inexf = (fop != op) ? mpfr_set (fop, op, rnd_mode) : 0; MPFR_SET_SAME_SIGN (iop, op); MPFR_SET_ZERO (iop); MPFR_RET (INEX(0, inexf)); } else if (ope >= opq) /* op has no fractional part */ { inexi = (iop != op) ? mpfr_set (iop, op, rnd_mode) : 0; MPFR_SET_SAME_SIGN (fop, op); MPFR_SET_ZERO (fop); MPFR_RET (INEX(inexi, 0)); } else /* op has both integral and fractional parts */ { if (iop != op) { inexi = mpfr_rint_trunc (iop, op, rnd_mode); inexf = mpfr_frac (fop, op, rnd_mode); } else { MPFR_ASSERTN (fop != op); inexf = mpfr_frac (fop, op, rnd_mode); inexi = mpfr_rint_trunc (iop, op, rnd_mode); } MPFR_RET (INEX(inexi, inexf)); } }
int mpfr_digamma (mpfr_ptr y, mpfr_srcptr x, mpfr_rnd_t rnd_mode) { int inex; MPFR_SAVE_EXPO_DECL (expo); MPFR_LOG_FUNC (("x[%Pu]=%.*Rg rnd=%d", mpfr_get_prec(x), mpfr_log_prec, x, rnd_mode), ("y[%Pu]=%.*Rg inexact=%d", mpfr_get_prec(y), mpfr_log_prec, y, inex)); if (MPFR_UNLIKELY(MPFR_IS_SINGULAR(x))) { if (MPFR_IS_NAN(x)) { MPFR_SET_NAN(y); MPFR_RET_NAN; } else if (MPFR_IS_INF(x)) { if (MPFR_IS_POS(x)) /* Digamma(+Inf) = +Inf */ { MPFR_SET_SAME_SIGN(y, x); MPFR_SET_INF(y); MPFR_RET(0); } else /* Digamma(-Inf) = NaN */ { MPFR_SET_NAN(y); MPFR_RET_NAN; } } else /* Zero case */ { /* the following works also in case of overlap */ MPFR_SET_INF(y); MPFR_SET_OPPOSITE_SIGN(y, x); MPFR_SET_DIVBY0 (); MPFR_RET(0); } } /* Digamma is undefined for negative integers */ if (MPFR_IS_NEG(x) && mpfr_integer_p (x)) { MPFR_SET_NAN(y); MPFR_RET_NAN; } /* now x is a normal number */ MPFR_SAVE_EXPO_MARK (expo); /* for x very small, we have Digamma(x) = -1/x - gamma + O(x), more precisely -1 < Digamma(x) + 1/x < 0 for -0.2 < x < 0.2, thus: (i) either x is a power of two, then 1/x is exactly representable, and as long as 1/2*ulp(1/x) > 1, we can conclude; (ii) otherwise assume x has <= n bits, and y has <= n+1 bits, then |y + 1/x| >= 2^(-2n) ufp(y), where ufp means unit in first place. Since |Digamma(x) + 1/x| <= 1, if 2^(-2n) ufp(y) >= 2, then |y - Digamma(x)| >= 2^(-2n-1)ufp(y), and rounding -1/x gives the correct result. If x < 2^E, then y > 2^(-E), thus ufp(y) > 2^(-E-1). A sufficient condition is thus EXP(x) <= -2 MAX(PREC(x),PREC(Y)). */ if (MPFR_EXP(x) < -2) { if (MPFR_EXP(x) <= -2 * (mpfr_exp_t) MAX(MPFR_PREC(x), MPFR_PREC(y))) { int signx = MPFR_SIGN(x); inex = mpfr_si_div (y, -1, x, rnd_mode); if (inex == 0) /* x is a power of two */ { /* result always -1/x, except when rounding down */ if (rnd_mode == MPFR_RNDA) rnd_mode = (signx > 0) ? MPFR_RNDD : MPFR_RNDU; if (rnd_mode == MPFR_RNDZ) rnd_mode = (signx > 0) ? MPFR_RNDU : MPFR_RNDD; if (rnd_mode == MPFR_RNDU) inex = 1; else if (rnd_mode == MPFR_RNDD) { mpfr_nextbelow (y); inex = -1; } else /* nearest */ inex = 1; } MPFR_SAVE_EXPO_UPDATE_FLAGS (expo, __gmpfr_flags); goto end; } } if (MPFR_IS_NEG(x)) inex = mpfr_digamma_reflection (y, x, rnd_mode); /* if x < 1/2 we use the reflection formula */ else if (MPFR_EXP(x) < 0) inex = mpfr_digamma_reflection (y, x, rnd_mode); else inex = mpfr_digamma_positive (y, x, rnd_mode); end: MPFR_SAVE_EXPO_FREE (expo); return mpfr_check_range (y, inex, rnd_mode); }