diff options
author | Charles Harris <charlesr.harris@gmail.com> | 2016-02-15 18:40:14 -0700 |
---|---|---|
committer | Charles Harris <charlesr.harris@gmail.com> | 2016-02-19 09:25:49 -0700 |
commit | e8302fd038147f4689ae021d3b5c1a02f2a8cd76 (patch) | |
tree | d6300cd8569ad88971de2094e75065e13a8524e6 | |
parent | 735174b88768999a2247f960669185978587408a (diff) | |
download | numpy-e8302fd038147f4689ae021d3b5c1a02f2a8cd76.tar.gz |
TST: Add tests for '//' and '%' (floor_divide, remainder).
Add tests for some corner cases involving inf, zero, and nan.
Check that small integers are handled exactly.
-rw-r--r-- | numpy/core/tests/test_multiarray.py | 36 | ||||
-rw-r--r-- | numpy/core/tests/test_scalarmath.py | 81 | ||||
-rw-r--r-- | numpy/core/tests/test_umath.py | 97 |
3 files changed, 170 insertions, 44 deletions
diff --git a/numpy/core/tests/test_multiarray.py b/numpy/core/tests/test_multiarray.py index d57e7c106..16df605f4 100644 --- a/numpy/core/tests/test_multiarray.py +++ b/numpy/core/tests/test_multiarray.py @@ -2448,42 +2448,6 @@ class TestMethods(TestCase): assert_raises(AttributeError, lambda: a.conj()) assert_raises(AttributeError, lambda: a.conjugate()) - def test_divmod_basic(self): - dt = np.typecodes['AllInteger'] + np.typecodes['Float'] - for dt1, dt2 in itertools.product(dt, dt): - for sg1, sg2 in itertools.product((+1, -1), (+1, -1)): - if sg1 == -1 and dt1 in np.typecodes['UnsignedInteger']: - continue - if sg2 == -1 and dt2 in np.typecodes['UnsignedInteger']: - continue - fmt = 'dt1: %s, dt2: %s, sg1: %s, sg2: %s' - msg = fmt % (dt1, dt2, sg1, sg2) - a = np.array(sg1*71, dtype=dt1) - b = np.array(sg2*19, dtype=dt2) - div, rem = divmod(a, b) - assert_allclose(div*b + rem, a, err_msg=msg) - if sg2 == -1: - assert_(b < rem <= 0, msg) - else: - assert_(b > rem >= 0, msg) - - def test_divmod_roundoff(self): - # gh-6127 - dt = 'fdg' - for dt1, dt2 in itertools.product(dt, dt): - for sg1, sg2 in itertools.product((+1, -1), (+1, -1)): - fmt = 'dt1: %s, dt2: %s, sg1: %s, sg2: %s' - msg = fmt % (dt1, dt2, sg1, sg2) - a = np.array(sg1*78*6e-8, dtype=dt1) - b = np.array(sg2*6e-8, dtype=dt2) - div, rem = divmod(a, b) - assert_allclose(div*b + rem, a, err_msg=msg) - if sg2 == -1: - assert_(b < rem <= 0, msg) - else: - assert_(b > rem >= 0, msg) - - class TestBinop(object): def test_inplace(self): # test refcount 1 inplace conversion diff --git a/numpy/core/tests/test_scalarmath.py b/numpy/core/tests/test_scalarmath.py index 17f70f6c9..12b1a0fe3 100644 --- a/numpy/core/tests/test_scalarmath.py +++ b/numpy/core/tests/test_scalarmath.py @@ -2,6 +2,8 @@ from __future__ import division, absolute_import, print_function import sys import itertools +import warnings +import operator import numpy as np from numpy.testing.utils import _gen_alignment_data @@ -137,8 +139,12 @@ class TestPower(TestCase): assert_almost_equal(result, 9, err_msg=msg) -class TestDivmod(TestCase): - def test_divmod_basic(self): +class TestModulus(TestCase): + + floordiv = operator.floordiv + mod = operator.mod + + def test_modulus_basic(self): dt = np.typecodes['AllInteger'] + np.typecodes['Float'] for dt1, dt2 in itertools.product(dt, dt): for sg1, sg2 in itertools.product((+1, -1), (+1, -1)): @@ -150,29 +156,88 @@ class TestDivmod(TestCase): msg = fmt % (dt1, dt2, sg1, sg2) a = np.array(sg1*71, dtype=dt1)[()] b = np.array(sg2*19, dtype=dt2)[()] - div, rem = divmod(a, b) - assert_allclose(div*b + rem, a, err_msg=msg) + div = self.floordiv(a, b) + rem = self.mod(a, b) + assert_equal(div*b + rem, a, err_msg=msg) if sg2 == -1: assert_(b < rem <= 0, msg) else: assert_(b > rem >= 0, msg) - def test_divmod_roundoff(self): + def test_float_modulus_exact(self): + # test that float results are exact for small integers. This also + # holds for the same integers scaled by powers of two. + nlst = list(range(-127, 0)) + plst = list(range(1, 128)) + dividend = nlst + [0] + plst + divisor = nlst + plst + arg = list(itertools.product(dividend, divisor)) + tgt = list(divmod(*t) for t in arg) + + a, b = np.array(arg, dtype=int).T + # convert exact integer results from Python to float so that + # signed zero can be used, it is checked. + tgtdiv, tgtrem = np.array(tgt, dtype=float).T + tgtdiv = np.where((tgtdiv == 0.0) & ((b < 0) ^ (a < 0)), -0.0, tgtdiv) + tgtrem = np.where((tgtrem == 0.0) & (b < 0), -0.0, tgtrem) + + for dt in np.typecodes['Float']: + msg = 'dtype: %s' % (dt,) + fa = a.astype(dt) + fb = b.astype(dt) + # use list comprehension so a_ and b_ are scalars + div = [self.floordiv(a_, b_) for a_, b_ in zip(fa, fb)] + rem = [self.mod(a_, b_) for a_, b_ in zip(fa, fb)] + assert_equal(div, tgtdiv, err_msg=msg) + assert_equal(rem, tgtrem, err_msg=msg) + + def test_float_modulus_roundoff(self): # gh-6127 - dt = 'fdg' + dt = np.typecodes['Float'] for dt1, dt2 in itertools.product(dt, dt): for sg1, sg2 in itertools.product((+1, -1), (+1, -1)): fmt = 'dt1: %s, dt2: %s, sg1: %s, sg2: %s' msg = fmt % (dt1, dt2, sg1, sg2) a = np.array(sg1*78*6e-8, dtype=dt1)[()] b = np.array(sg2*6e-8, dtype=dt2)[()] - div, rem = divmod(a, b) - assert_allclose(div*b + rem, a, err_msg=msg) + div = self.floordiv(a, b) + rem = self.mod(a, b) + # Equal assertion should hold when fmod is used + assert_equal(div*b + rem, a, err_msg=msg) if sg2 == -1: assert_(b < rem <= 0, msg) else: assert_(b > rem >= 0, msg) + def test_float_modulus_corner_cases(self): + # Check remainder magnitude. + for dt in np.typecodes['Float']: + b = np.array(1.0, dtype=dt) + a = np.nextafter(np.array(0.0, dtype=dt), -b) + rem = self.mod(a, b) + assert_(rem <= b, 'dt: %s' % dt) + rem = self.mod(-a, -b) + assert_(rem >= -b, 'dt: %s' % dt) + + # Check nans, inf + with warnings.catch_warnings(): + warnings.simplefilter('always') + warnings.simplefilter('ignore', RuntimeWarning) + for dt in np.typecodes['Float']: + fone = np.array(1.0, dtype=dt) + fzer = np.array(0.0, dtype=dt) + finf = np.array(np.inf, dtype=dt) + fnan = np.array(np.nan, dtype=dt) + rem = self.mod(fone, fzer) + assert_(np.isnan(rem), 'dt: %s' % dt) + # MSVC 2008 returns NaN here, so disable the check. + #rem = self.mod(fone, finf) + #assert_(rem == fone, 'dt: %s' % dt) + rem = self.mod(fone, fnan) + assert_(np.isnan(rem), 'dt: %s' % dt) + rem = self.mod(finf, fone) + assert_(np.isnan(rem), 'dt: %s' % dt) + class TestComplexDivision(TestCase): def test_zero_division(self): diff --git a/numpy/core/tests/test_umath.py b/numpy/core/tests/test_umath.py index 917e05e6a..da52e0dde 100644 --- a/numpy/core/tests/test_umath.py +++ b/numpy/core/tests/test_umath.py @@ -3,6 +3,7 @@ from __future__ import division, absolute_import, print_function import sys import platform import warnings +import itertools from numpy.testing.utils import _gen_alignment_data import numpy.core.umath as ncu @@ -222,6 +223,102 @@ class TestDivision(TestCase): assert_equal(y, [1.e+110, 0], err_msg=msg) +class TestRemainder(TestCase): + + def test_remainder_basic(self): + dt = np.typecodes['AllInteger'] + np.typecodes['Float'] + for dt1, dt2 in itertools.product(dt, dt): + for sg1, sg2 in itertools.product((+1, -1), (+1, -1)): + if sg1 == -1 and dt1 in np.typecodes['UnsignedInteger']: + continue + if sg2 == -1 and dt2 in np.typecodes['UnsignedInteger']: + continue + fmt = 'dt1: %s, dt2: %s, sg1: %s, sg2: %s' + msg = fmt % (dt1, dt2, sg1, sg2) + a = np.array(sg1*71, dtype=dt1) + b = np.array(sg2*19, dtype=dt2) + div = np.floor_divide(a, b) + rem = np.remainder(a, b) + assert_equal(div*b + rem, a, err_msg=msg) + if sg2 == -1: + assert_(b < rem <= 0, msg) + else: + assert_(b > rem >= 0, msg) + + def test_float_remainder_exact(self): + # test that float results are exact for small integers. This also + # holds for the same integers scaled by powers of two. + nlst = list(range(-127, 0)) + plst = list(range(1, 128)) + dividend = nlst + [0] + plst + divisor = nlst + plst + arg = list(itertools.product(dividend, divisor)) + tgt = list(divmod(*t) for t in arg) + + a, b = np.array(arg, dtype=int).T + # convert exact integer results from Python to float so that + # signed zero can be used, it is checked. + tgtdiv, tgtrem = np.array(tgt, dtype=float).T + tgtdiv = np.where((tgtdiv == 0.0) & ((b < 0) ^ (a < 0)), -0.0, tgtdiv) + tgtrem = np.where((tgtrem == 0.0) & (b < 0), -0.0, tgtrem) + + for dt in np.typecodes['Float']: + msg = 'dtype: %s' % (dt,) + fa = a.astype(dt) + fb = b.astype(dt) + div = np.floor_divide(fa, fb) + rem = np.remainder(fa, fb) + assert_equal(div, tgtdiv, err_msg=msg) + assert_equal(rem, tgtrem, err_msg=msg) + + def test_float_remainder_roundoff(self): + # gh-6127 + dt = np.typecodes['Float'] + for dt1, dt2 in itertools.product(dt, dt): + for sg1, sg2 in itertools.product((+1, -1), (+1, -1)): + fmt = 'dt1: %s, dt2: %s, sg1: %s, sg2: %s' + msg = fmt % (dt1, dt2, sg1, sg2) + a = np.array(sg1*78*6e-8, dtype=dt1) + b = np.array(sg2*6e-8, dtype=dt2) + div = np.floor_divide(a, b) + rem = np.remainder(a, b) + # Equal assertion should hold when fmod is used + assert_equal(div*b + rem, a, err_msg=msg) + if sg2 == -1: + assert_(b < rem <= 0, msg) + else: + assert_(b > rem >= 0, msg) + + def test_float_remainder_corner_cases(self): + # Check remainder magnitude. + for dt in np.typecodes['Float']: + b = np.array(1.0, dtype=dt) + a = np.nextafter(np.array(0.0, dtype=dt), -b) + rem = np.remainder(a, b) + assert_(rem <= b, 'dt: %s' % dt) + rem = np.remainder(-a, -b) + assert_(rem >= -b, 'dt: %s' % dt) + + # Check nans, inf + with warnings.catch_warnings(): + warnings.simplefilter('always') + warnings.simplefilter('ignore', RuntimeWarning) + for dt in np.typecodes['Float']: + fone = np.array(1.0, dtype=dt) + fzer = np.array(0.0, dtype=dt) + finf = np.array(np.inf, dtype=dt) + fnan = np.array(np.nan, dtype=dt) + rem = np.remainder(fone, fzer) + assert_(np.isnan(rem), 'dt: %s, rem: %s' % (dt, rem)) + # MSVC 2008 returns NaN here, so disable the check. + #rem = np.remainder(fone, finf) + #assert_(rem == fone, 'dt: %s, rem: %s' % (dt, rem)) + rem = np.remainder(fone, fnan) + assert_(np.isnan(rem), 'dt: %s, rem: %s' % (dt, rem)) + rem = np.remainder(finf, fone) + assert_(np.isnan(rem), 'dt: %s, rem: %s' % (dt, rem)) + + class TestCbrt(TestCase): def test_cbrt_scalar(self): assert_almost_equal((np.cbrt(np.float32(-2.5)**3)), -2.5) |