summaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authorvlefevre <vlefevre@280ebfd0-de03-0410-8827-d642c229c3f4>2018-02-01 13:32:53 +0000
committervlefevre <vlefevre@280ebfd0-de03-0410-8827-d642c229c3f4>2018-02-01 13:32:53 +0000
commit50158da09ed01cfe1e3d26d6d2449312297865b1 (patch)
tree96b0494e13691e638744d49516bef4f1027cc8fe
parent42bc83502a1552f6f4af7d92faa7f2f36dcd9aee (diff)
downloadmpfr-50158da09ed01cfe1e3d26d6d2449312297865b1.tar.gz
[src/div_ui.c] Various bug fixes.
(merged changesets r10704,11014,11141-11155,12139-12163 from the trunk) [tests/tdiv_ui.c] Added tests, in particular triggering the above bugs. (merged changesets r11141-11144,12139-12159,12167-12172 from the trunk) git-svn-id: svn://scm.gforge.inria.fr/svn/mpfr/branches/3.1@12173 280ebfd0-de03-0410-8827-d642c229c3f4
-rw-r--r--src/div_ui.c192
-rw-r--r--tests/tdiv_ui.c277
2 files changed, 386 insertions, 83 deletions
diff --git a/src/div_ui.c b/src/div_ui.c
index 3976d9b26..1cf20c852 100644
--- a/src/div_ui.c
+++ b/src/div_ui.c
@@ -25,14 +25,16 @@ http://www.gnu.org/licenses/ or write to the Free Software Foundation, Inc.,
/* returns 0 if result exact, non-zero otherwise */
int
-mpfr_div_ui (mpfr_ptr y, mpfr_srcptr x, unsigned long int u, mpfr_rnd_t rnd_mode)
+mpfr_div_ui (mpfr_ptr y, mpfr_srcptr x, unsigned long int u,
+ mpfr_rnd_t rnd_mode)
{
- long i;
int sh;
- mp_size_t xn, yn, dif;
+ mp_size_t i, xn, yn, dif;
mp_limb_t *xp, *yp, *tmp, c, d;
mpfr_exp_t exp;
- int inexact, middle = 1, nexttoinf;
+ int inexact;
+ mp_limb_t rb; /* round bit */
+ mp_limb_t sb; /* sticky bit */
MPFR_TMP_DECL(marker);
MPFR_LOG_FUNC
@@ -56,10 +58,10 @@ mpfr_div_ui (mpfr_ptr y, mpfr_srcptr x, unsigned long int u, mpfr_rnd_t rnd_mode
}
else
{
- MPFR_ASSERTD (MPFR_IS_ZERO(x));
+ MPFR_ASSERTD (MPFR_IS_ZERO (x));
if (u == 0) /* 0/0 is NaN */
{
- MPFR_SET_NAN(y);
+ MPFR_SET_NAN (y);
MPFR_RET_NAN;
}
else
@@ -74,7 +76,7 @@ mpfr_div_ui (mpfr_ptr y, mpfr_srcptr x, unsigned long int u, mpfr_rnd_t rnd_mode
{
if (u < 1)
{
- /* x/0 is Inf since x != 0*/
+ /* x/0 is Inf since x != 0 */
MPFR_SET_INF (y);
MPFR_SET_SAME_SIGN (y, x);
mpfr_set_divby0 ();
@@ -89,6 +91,7 @@ mpfr_div_ui (mpfr_ptr y, mpfr_srcptr x, unsigned long int u, mpfr_rnd_t rnd_mode
MPFR_SET_SAME_SIGN (y, x);
MPFR_TMP_MARK (marker);
+
xn = MPFR_LIMB_SIZE (x);
yn = MPFR_LIMB_SIZE (y);
@@ -98,105 +101,136 @@ mpfr_div_ui (mpfr_ptr y, mpfr_srcptr x, unsigned long int u, mpfr_rnd_t rnd_mode
dif = yn + 1 - xn;
- /* we need to store yn+1 = xn + dif limbs of the quotient */
- /* don't use tmp=yp since the mpn_lshift call below requires yp >= tmp+1 */
+ /* we need to store yn + 1 = xn + dif limbs of the quotient */
tmp = MPFR_TMP_LIMBS_ALLOC (yn + 1);
- c = (mp_limb_t) u;
- MPFR_ASSERTN (u == c);
- if (dif >= 0)
- c = mpn_divrem_1 (tmp, dif, xp, xn, c); /* used all the dividend */
- else /* dif < 0 i.e. xn > yn, don't use the (-dif) low limbs from x */
- c = mpn_divrem_1 (tmp, 0, xp - dif, yn + 1, c);
-
- inexact = (c != 0);
+ /* Notation: {p, n} denotes the integer formed by the n limbs
+ from p[0] to p[n-1]. Let B = 2^GMP_NUMB_BITS.
+ One has: 0 <= {p, n} < B^n. */
- /* First pass in estimating next bit of the quotient, in case of RNDN *
- * In case we just have the right number of bits (postpone this ?), *
- * we need to check whether the remainder is more or less than half *
- * the divisor. The test must be performed with a subtraction, so as *
- * to prevent carries. */
-
- if (MPFR_LIKELY (rnd_mode == MPFR_RNDN))
+ MPFR_ASSERTN (u == (mp_limb_t) u);
+ if (dif >= 0)
{
- if (c < (mp_limb_t) u - c) /* We have u > c */
- middle = -1;
- else if (c > (mp_limb_t) u - c)
- middle = 1;
- else
- middle = 0; /* exactly in the middle */
+ c = mpn_divrem_1 (tmp, dif, xp, xn, u); /* used all the dividend */
+ /* {xp, xn} = ({tmp, xn+dif} * u + c) * B^(-dif)
+ = ({tmp, yn+1} * u + c) * B^(-dif) */
+ }
+ else /* dif < 0, i.e. xn > yn+1; ignore the (-dif) low limbs from x */
+ {
+ c = mpn_divrem_1 (tmp, 0, xp - dif, yn + 1, u);
+ /* {xp-dif, yn+1} = {tmp, yn+1} * u + c
+ thus
+ {xp, xn} = {xp, -dif} + {xp-dif, yn+1} * B^(-dif)
+ = {xp, -dif} + ({tmp, yn+1} * u + c) * B^(-dif) */
}
- /* If we believe that we are right in the middle or exact, we should check
- that we did not neglect any word of x (division large / 1 -> small). */
+ /* Let r = {xp, -dif} / B^(-dif) if dif < 0, r = 0 otherwise; 0 <= r < 1.
+ Then {xp, xn} = ({tmp, yn+1} * u + c + r) * B^(-dif).
+ x / u = ({xp, xn} / u) * B^(-xn) * 2^exp
+ = ({tmp, yn+1} + (c + r) / u) * B^(-(yn+1)) * 2^exp
+ where 0 <= (c + r) / u < 1. */
- for (i=0; ((inexact == 0) || (middle == 0)) && (i < -dif); i++)
+ for (sb = 0, i = 0; sb == 0 && i < -dif; i++)
if (xp[i])
- inexact = middle = 1; /* larger than middle */
+ sb = 1;
+ /* sb != 0 iff r != 0 */
/*
- If the high limb of the result is 0 (xp[xn-1] < u), remove it.
+ If the highest limb of the result is 0 (xp[xn-1] < u), remove it.
Otherwise, compute the left shift to be performed to normalize.
In the latter case, we discard some low bits computed. They
contain information useful for the rounding, hence the updating
of middle and inexact.
*/
+ MPFR_UNSIGNED_MINUS_MODULO (sh, MPFR_PREC (y));
+ /* sh: number of the trailing bits of y */
+
if (tmp[yn] == 0)
{
MPN_COPY(yp, tmp, yn);
exp -= GMP_NUMB_BITS;
+ if (sh == 0) /* round bit is 1 iff (c + r) / u >= 1/2 */
+ {
+ /* In this case tmp[yn]=0 and sh=0, the round bit is not in
+ {tmp,yn+1}. It is 1 iff 2*(c+r) - u >= 0. This means that in
+ some cases, we should look at the most significant bit of r. */
+ if (c >= u - c) /* i.e. 2c >= u: round bit is always 1 */
+ {
+ rb = 1;
+ /* The sticky bit is 1 unless 2c-u = 0 and r = 0. */
+ sb |= 2 * c - u;
+ }
+ else /* 2*c < u */
+ {
+ /* The round bit is 1 iff r >= 1/2 and 2*(c+1/2) = u. */
+ rb = (c == u/2) && (dif < 0) && (xp[-dif-1] & MPFR_LIMB_HIGHBIT);
+ /* If rb is set, we need to recompute sb, since it might have
+ taken into account the msb of xp[-dif-1]. */
+ if (rb)
+ {
+ sb = xp[-dif-1] << 1; /* discard the most significant bit */
+ for (i = 0; sb == 0 && i < -dif-1; i++)
+ if (xp[i])
+ sb = 1;
+ }
+ else
+ sb |= c;
+ }
+ }
+ else
+ {
+ /* round bit is in tmp[0] */
+ rb = tmp[0] & (MPFR_LIMB_ONE << (sh - 1));
+ sb |= (tmp[0] & MPFR_LIMB_MASK(sh - 1)) | c;
+ }
}
- else
+ else /* tmp[yn] != 0 */
{
int shlz;
+ mp_limb_t w;
+ MPFR_ASSERTD (tmp[yn] != 0);
count_leading_zeros (shlz, tmp[yn]);
- /* shift left to normalize */
- if (MPFR_LIKELY (shlz != 0))
- {
- mp_limb_t w = tmp[0] << shlz;
-
- mpn_lshift (yp, tmp + 1, yn, shlz);
- yp[0] += tmp[0] >> (GMP_NUMB_BITS - shlz);
+ MPFR_ASSERTD (u >= 2); /* see special cases at the beginning */
+ MPFR_ASSERTD (shlz > 0); /* since u >= 2 */
- if (w > (MPFR_LIMB_ONE << (GMP_NUMB_BITS - 1)))
- { middle = 1; }
- else if (w < (MPFR_LIMB_ONE << (GMP_NUMB_BITS - 1)))
- { middle = -1; }
- else
- { middle = (c != 0); }
+ /* shift left to normalize */
+ w = tmp[0] << shlz;
+ mpn_lshift (yp, tmp + 1, yn, shlz);
+ yp[0] |= tmp[0] >> (GMP_NUMB_BITS - shlz);
+ /* now {yp, yn} is the approximate quotient, w is the next limb */
- inexact = inexact || (w != 0);
- exp -= shlz;
+ if (sh == 0) /* round bit is upper bit from w */
+ {
+ rb = w & MPFR_LIMB_HIGHBIT;
+ sb |= (w - rb) | c;
}
else
- { /* this happens only if u == 1 and xp[xn-1] >=
- 1<<(GMP_NUMB_BITS-1). It might be better to handle the
- u == 1 case separately?
- */
-
- MPN_COPY (yp, tmp + 1, yn);
+ {
+ rb = yp[0] & (MPFR_LIMB_ONE << (sh - 1));
+ sb |= (yp[0] & MPFR_LIMB_MASK(sh - 1)) | w | c;
}
- }
- MPFR_UNSIGNED_MINUS_MODULO (sh, MPFR_PREC (y));
- /* it remains sh bits in less significant limb of y */
+ exp -= shlz;
+ }
- d = *yp & MPFR_LIMB_MASK (sh);
- *yp ^= d; /* set to zero lowest sh bits */
+ d = yp[0] & MPFR_LIMB_MASK (sh);
+ yp[0] ^= d; /* clear the lowest sh bits */
MPFR_TMP_FREE (marker);
- if (exp < __gmpfr_emin - 1)
+ if (MPFR_UNLIKELY (exp < __gmpfr_emin - 1))
return mpfr_underflow (y, rnd_mode == MPFR_RNDN ? MPFR_RNDZ : rnd_mode,
MPFR_SIGN (y));
- if (MPFR_UNLIKELY (d == 0 && inexact == 0))
- nexttoinf = 0; /* result is exact */
+ if (MPFR_UNLIKELY (rb == 0 && sb == 0))
+ inexact = 0; /* result is exact */
else
{
+ int nexttoinf;
+
MPFR_UPDATE2_RND_MODE(rnd_mode, MPFR_SIGN (y));
switch (rnd_mode)
{
@@ -213,26 +247,19 @@ mpfr_div_ui (mpfr_ptr y, mpfr_srcptr x, unsigned long int u, mpfr_rnd_t rnd_mode
default: /* should be MPFR_RNDN */
MPFR_ASSERTD (rnd_mode == MPFR_RNDN);
/* We have one more significant bit in yn. */
- if (sh && d < (MPFR_LIMB_ONE << (sh - 1)))
+ if (rb == 0)
{
inexact = - MPFR_INT_SIGN (y);
nexttoinf = 0;
}
- else if (sh && d > (MPFR_LIMB_ONE << (sh - 1)))
+ else if (sb != 0) /* necessarily rb != 0 */
{
inexact = MPFR_INT_SIGN (y);
nexttoinf = 1;
}
- else /* sh = 0 or d = 1 << (sh-1) */
+ else /* middle case */
{
- /* The first case is "false" even rounding (significant bits
- indicate even rounding, but the result is inexact, so up) ;
- The second case is the case where middle should be used to
- decide the direction of rounding (no further bit computed) ;
- The third is the true even rounding.
- */
- if ((sh && inexact) || (!sh && middle > 0) ||
- (!inexact && *yp & (MPFR_LIMB_ONE << sh)))
+ if (yp[0] & (MPFR_LIMB_ONE << sh))
{
inexact = MPFR_INT_SIGN (y);
nexttoinf = 1;
@@ -244,13 +271,12 @@ mpfr_div_ui (mpfr_ptr y, mpfr_srcptr x, unsigned long int u, mpfr_rnd_t rnd_mode
}
}
}
- }
-
- if (nexttoinf &&
- MPFR_UNLIKELY (mpn_add_1 (yp, yp, yn, MPFR_LIMB_ONE << sh)))
- {
- exp++;
- yp[yn-1] = MPFR_LIMB_HIGHBIT;
+ if (nexttoinf &&
+ MPFR_UNLIKELY (mpn_add_1 (yp, yp, yn, MPFR_LIMB_ONE << sh)))
+ {
+ exp++;
+ yp[yn-1] = MPFR_LIMB_HIGHBIT;
+ }
}
/* Set the exponent. Warning! One may still have an underflow. */
diff --git a/tests/tdiv_ui.c b/tests/tdiv_ui.c
index 9bfa10bcb..34792d260 100644
--- a/tests/tdiv_ui.c
+++ b/tests/tdiv_ui.c
@@ -196,6 +196,275 @@ check_inexact (void)
mpfr_clear (z);
}
+#if GMP_NUMB_BITS == 64
+/* With r11140, on a 64-bit machine with GMP_CHECK_RANDOMIZE=1484406128:
+ Consistency error for i = 2577
+*/
+static void
+test_20170105 (void)
+{
+ mpfr_t x,z, t;
+
+ if (sizeof (unsigned long) * CHAR_BIT == 64)
+ {
+ mpfr_init2 (x, 138);
+ mpfr_init2 (z, 128);
+ mpfr_init2 (t, 128);
+ mpfr_set_str_binary (x, "100110111001001000101111010010011101111110111111110001110100000001110111010100111010100011101010110000010100000011100100110101101011000000");
+ /* up to exponents, x/y is exactly
+ 367625553447399614694201910705139062483, which has 129 bits,
+ thus we are in the round-to-nearest-even case, and since the
+ penultimate bit of x/y is 1, we should round upwards */
+ mpfr_set_str_binary (t, "10001010010010010000110110010110111111111100011011101010000000000110101000010001011110011011010000111010000000001100101101101010E-53");
+ mpfr_div_ui (z, x, 36UL << 58, MPFR_RNDN);
+ MPFR_ASSERTN(mpfr_equal_p (z, t));
+
+ mpfr_set_prec (x, 189);
+ mpfr_set_prec (z, 185);
+ mpfr_set_prec (t, 185);
+ mpfr_set_str_binary (x, "100001010000111100011110111010000011110000000110100010001010101011110001110000110111101000100100001101010011000111110100011111110110011011101000000000001010010010111011001100111111111101001");
+ mpfr_set_str_binary (t, "10011000000100010100011111100100110101101110001011100101010101011010011010010110010000100111001010000101111011111111001011011010101111101100000000000000101111000100001110101001001001000E-60");
+ mpfr_div_ui (z, x, 7UL << 61, MPFR_RNDN);
+ MPFR_ASSERTN(mpfr_equal_p (z, t));
+
+ mpfr_clears (x, z, t, (mpfr_ptr) 0);
+ }
+}
+#endif
+
+static void
+bug20180126 (void)
+{
+ mpfr_t w, x, y, z, t;
+ unsigned long u, v;
+ int i, k, m, n, p, pmax, q, r;
+ int inex;
+
+ /* Let m = n * q + r, with 0 <= r < v.
+ (2^m-1) / (2^n-1) = 2^r * (2^(n*q)-1) / (2^n-1) + (2^r-1) / (2^n-1)
+ = sum(i=0,q-1,2^(r+n*i)) + sum(i=1,inf,(2^r-1)*2^(-n*i))
+ */
+ n = 1;
+ for (u = 1; u != ULONG_MAX; u = (u << 1) + 1)
+ n++;
+ pmax = 6 * n;
+ mpfr_init2 (t, n);
+ for (m = n; m < 4 * n; m++)
+ {
+ q = m / n;
+ r = m % n;
+ mpfr_init2 (w, pmax + n + 1);
+ mpfr_set_zero (w, 1);
+ for (i = 0; i < q; i++)
+ {
+ inex = mpfr_set_ui_2exp (t, 1, r + n * i, MPFR_RNDN);
+ MPFR_ASSERTN (inex == 0);
+ inex = mpfr_add (w, w, t, MPFR_RNDN);
+ MPFR_ASSERTN (inex == 0);
+ }
+ v = (1UL << r) - 1;
+ for (i = 1; n * (q - 1 + i) <= MPFR_PREC (w); i++)
+ {
+ inex = mpfr_set_ui_2exp (t, v, - n * i, MPFR_RNDN);
+ MPFR_ASSERTN (inex == 0);
+ mpfr_add (w, w, t, MPFR_RNDN);
+ }
+ for (p = pmax; p >= MPFR_PREC_MIN; p--)
+ {
+ mpfr_inits2 (p, y, z, (mpfr_ptr) 0);
+ mpfr_set (z, w, MPFR_RNDN); /* the sticky bit is not 0 */
+ mpfr_init2 (x, m);
+ inex = mpfr_set_ui_2exp (x, 1, m, MPFR_RNDN);
+ MPFR_ASSERTN (inex == 0);
+ inex = mpfr_sub_ui (x, x, 1, MPFR_RNDN); /* x = 2^m-1 */
+ MPFR_ASSERTN (inex == 0);
+ for (k = 0; k < 2; k++)
+ {
+ if (k)
+ {
+ inex = mpfr_prec_round (x, 6 * n, MPFR_RNDN);
+ MPFR_ASSERTN (inex == 0);
+ }
+ inex = mpfr_div_ui (y, x, u, MPFR_RNDN);
+ if (! mpfr_equal_p (y, z))
+ {
+ printf ("Error in bug20180126 for (2^%d-1)/(2^%d-1)"
+ " with px=%d py=%d\n", m, n,
+ (int) MPFR_PREC (x), p);
+ printf ("Expected ");
+ mpfr_dump (z);
+ printf ("Got ");
+ mpfr_dump (y);
+ exit (1);
+ }
+ }
+ mpfr_clears (x, y, z, (mpfr_ptr) 0);
+ }
+ mpfr_clear (w);
+ }
+ mpfr_clear (t);
+
+ /* This test expects that a limb fits in an unsigned long.
+ One failing case from function bug20180126() in tdiv.c,
+ for GMP_NUMB_BITS == 64. */
+ if (GMP_NUMB_BITS == 64 && MP_LIMB_T_MAX <= ULONG_MAX)
+ {
+ mpfr_init2 (x, 133);
+ mpfr_init2 (y, 64);
+ mpfr_set_ui (x, 1, MPFR_RNDN);
+ mpfr_nextbelow (x); /* 1 - 2^(-133) = (2^133-1)/2^133 */
+ u = MP_LIMB_T_MAX; /* 2^64 - 1 */
+ inex = mpfr_div_ui (y, x, u, MPFR_RNDN);
+ /* 2^133*x/u = (2^133-1)/(2^64-1) = (2^64+1)*2^5 + 31/(2^64-1)
+ and should be rounded to 2^69+2^6, thus x/u should be rounded
+ to 2^(-133)*(2^69+2^6). */
+ MPFR_ASSERTN (inex > 0);
+ mpfr_nextbelow (y);
+ MPFR_ASSERTN (mpfr_cmp_ui_2exp (y, 1, -64) == 0);
+
+ mpfr_set_prec (x, 49);
+ mpfr_set_str_binary (x, "0.1000000000000000000111111111111111111111100000000E0");
+ /* x = 281476050452224/2^49 */
+ /* let X = 2^256*x = q*u+r, then q has 192 bits, and
+ r = 8222597979955926678 > u/2 thus we should round to (q+1)/2^256 */
+ mpfr_set_prec (y, 192);
+ /* The cast below avoid spurious warnings from GCC with a 32-bit ABI. */
+ u = (mp_limb_t) 10865468317030705979U;
+ inex = mpfr_div_ui (y, x, u, MPFR_RNDN);
+ mpfr_init2 (z, 192);
+ mpfr_set_str_binary (z, "0.110110010100111111000100101011011110010101010010001101100110101111001010100011010111010011100001101000110100011101001010000001010000001001011100000100000110101111110100100101011000000110011111E-64");
+ MPFR_ASSERTN (inex > 0);
+ MPFR_ASSERTN (mpfr_equal_p (y, z));
+ mpfr_clear (x);
+ mpfr_clear (y);
+ mpfr_clear (z);
+ }
+}
+
+/* check corner cases where the round bit is located in the upper bit of r */
+static void
+corner_cases (int n)
+{
+ mpfr_t x, y, t;
+ unsigned long u, v;
+ int i, xn;
+
+ if (MP_LIMB_T_MAX <= ULONG_MAX)
+ {
+ /* We need xn > yn + 1, thus we take xn = 3 and yn = 1.
+ Also take xn = 4 to 6 to cover more code. */
+ for (xn = 3; xn < 6; xn++)
+ {
+ mpfr_init2 (x, xn * GMP_NUMB_BITS);
+ mpfr_init2 (y, GMP_NUMB_BITS);
+ mpfr_init2 (t, 2 * GMP_NUMB_BITS);
+ for (i = 0; i < n; i++)
+ {
+ u = randlimb ();
+ do
+ v = randlimb ();
+ while (v <= MPFR_LIMB_HIGHBIT);
+ mpfr_set_ui (t, v, MPFR_RNDN);
+ mpfr_sub_d (t, t, 0.5, MPFR_RNDN);
+ /* t = v-1/2 */
+ mpfr_mul_ui (x, t, u, MPFR_RNDN);
+
+ /* when x = (v-1/2)*u, x/u should give v-1/2, which should give
+ either v (if v is even) or v-1 (if v is odd) */
+ mpfr_div_ui (y, x, u, MPFR_RNDN);
+ MPFR_ASSERTN(mpfr_cmp_ui (y, v - (v & 1)) == 0);
+
+ /* when x = (v-1/2)*u - epsilon, x/u should give v-1 */
+ mpfr_nextbelow (x);
+ mpfr_div_ui (y, x, u, MPFR_RNDN);
+ MPFR_ASSERTN(mpfr_cmp_ui (y, v - 1) == 0);
+
+ /* when x = (v-1/2)*u + epsilon, x/u should give v */
+ mpfr_nextabove (x);
+ mpfr_nextabove (x);
+ mpfr_div_ui (y, x, u, MPFR_RNDN);
+ MPFR_ASSERTN(mpfr_cmp_ui (y, v) == 0);
+ }
+ mpfr_clear (x);
+ mpfr_clear (y);
+ mpfr_clear (t);
+ }
+ }
+}
+
+static void
+midpoint_exact (void)
+{
+ mpfr_t x, y1, y2;
+ unsigned long j;
+ int i, kx, ky, px, pxmin, py, pymin, r;
+ int inex1, inex2;
+
+ pymin = 1;
+ for (i = 3; i < 32; i += 2)
+ {
+ if ((i & (i-2)) == 1)
+ pymin++;
+ for (j = 1; j != 0; j++)
+ {
+ if (j == 31)
+ j = ULONG_MAX;
+ /* Test of (i*j) / j with various precisions. The target precisions
+ include: large, length(i), and length(i)-1; the latter case
+ corresponds to a midpoint. */
+ mpfr_init2 (x, 5 + sizeof(long) * CHAR_BIT);
+ inex1 = mpfr_set_ui (x, j, MPFR_RNDN);
+ MPFR_ASSERTN (inex1 == 0);
+ inex1 = mpfr_mul_ui (x, x, i, MPFR_RNDN);
+ MPFR_ASSERTN (inex1 == 0);
+ /* x = (i*j) */
+ pxmin = mpfr_min_prec (x);
+ if (pxmin < MPFR_PREC_MIN)
+ pxmin = MPFR_PREC_MIN;
+ for (kx = 0; kx < 8; kx++)
+ {
+ px = pxmin;
+ if (kx != 0)
+ px += randlimb () % (4 * GMP_NUMB_BITS);
+ inex1 = mpfr_prec_round (x, px, MPFR_RNDN);
+ MPFR_ASSERTN (inex1 == 0);
+ for (ky = 0; ky < 8; ky++)
+ {
+ py = pymin;
+ if (ky == 0)
+ py--;
+ else if (ky > 1)
+ py += randlimb () % (4 * GMP_NUMB_BITS);
+ if (py < MPFR_PREC_MIN)
+ break;
+ mpfr_inits2 (py, y1, y2, (mpfr_ptr) 0);
+ RND_LOOP (r)
+ {
+ inex1 = mpfr_set_ui (y1, i, (mpfr_rnd_t) r);
+ inex2 = mpfr_div_ui (y2, x, j, (mpfr_rnd_t) r);
+ if (! mpfr_equal_p (y1, y2) ||
+ ! SAME_SIGN (inex1, inex2))
+ {
+ printf ("Error in midpoint_exact for "
+ "i=%d j=%lu px=%d py=%d %s\n", i, j, px, py,
+ mpfr_print_rnd_mode ((mpfr_rnd_t) r));
+ printf ("Expected ");
+ mpfr_dump (y1);
+ printf ("with inex = %d\n", inex1);
+ printf ("Got ");
+ mpfr_dump (y2);
+ printf ("with inex = %d\n", inex2);
+ exit (1);
+ }
+ }
+ mpfr_clears (y1, y2, (mpfr_ptr) 0);
+ }
+ }
+ mpfr_clear (x);
+ }
+ }
+}
+
#define TEST_FUNCTION mpfr_div_ui
#define ULONG_ARG2
#define RAND_FUNCTION(x) mpfr_random2(x, MPFR_LIMB_SIZE (x), 1, RANDS)
@@ -208,6 +477,10 @@ main (int argc, char **argv)
tests_start_mpfr ();
+ corner_cases (100);
+
+ bug20180126 ();
+
special ();
check_inexact ();
@@ -218,6 +491,9 @@ main (int argc, char **argv)
check("1.0", 3, MPFR_RNDD, "3.3333333333333331483e-1");
check("1.0", 2116118, MPFR_RNDN, "4.7256343927890600483e-7");
check("1.098612288668109782", 5, MPFR_RNDN, "0.21972245773362195087");
+#if GMP_NUMB_BITS == 64
+ test_20170105 ();
+#endif
mpfr_init2 (x, 53);
mpfr_set_ui (x, 3, MPFR_RNDD);
@@ -231,6 +507,7 @@ main (int argc, char **argv)
mpfr_clear (x);
test_generic (2, 200, 100);
+ midpoint_exact ();
tests_end_mpfr ();
return 0;