summaryrefslogtreecommitdiff
path: root/numpy/lib
diff options
context:
space:
mode:
Diffstat (limited to 'numpy/lib')
-rw-r--r--numpy/lib/_iotools.py2
-rw-r--r--numpy/lib/arraysetops.py5
-rw-r--r--numpy/lib/function_base.py9
-rw-r--r--numpy/lib/index_tricks.py7
-rw-r--r--numpy/lib/nanfunctions.py6
-rw-r--r--numpy/lib/npyio.py7
-rw-r--r--numpy/lib/setup.py1
-rw-r--r--numpy/lib/tests/test_format.py32
-rw-r--r--numpy/lib/tests/test_function_base.py7
-rw-r--r--numpy/lib/tests/test_io.py26
-rw-r--r--numpy/lib/tests/test_twodim_base.py34
-rw-r--r--numpy/lib/twodim_base.py9
12 files changed, 102 insertions, 43 deletions
diff --git a/numpy/lib/_iotools.py b/numpy/lib/_iotools.py
index 1b1180893..9108b2e4c 100644
--- a/numpy/lib/_iotools.py
+++ b/numpy/lib/_iotools.py
@@ -687,7 +687,7 @@ class StringConverter(object):
def upgrade(self, value):
"""
- Rind the best converter for a given string, and return the result.
+ Find the best converter for a given string, and return the result.
The supplied string `value` is converted by testing different
converters in order. First the `func` method of the
diff --git a/numpy/lib/arraysetops.py b/numpy/lib/arraysetops.py
index 2d98c35d2..5ace7146a 100644
--- a/numpy/lib/arraysetops.py
+++ b/numpy/lib/arraysetops.py
@@ -204,8 +204,9 @@ def unique(ar, return_index=False, return_inverse=False, return_counts=False):
ret += (perm[flag],)
if return_inverse:
iflag = np.cumsum(flag) - 1
- iperm = perm.argsort()
- ret += (np.take(iflag, iperm),)
+ inv_idx = np.empty_like(ar, dtype=np.intp)
+ inv_idx[perm] = iflag
+ ret += (inv_idx,)
if return_counts:
idx = np.concatenate(np.nonzero(flag) + ([ar.size],))
ret += (np.diff(idx),)
diff --git a/numpy/lib/function_base.py b/numpy/lib/function_base.py
index ff69b07dc..3074a2f70 100644
--- a/numpy/lib/function_base.py
+++ b/numpy/lib/function_base.py
@@ -337,6 +337,11 @@ def histogramdd(sample, bins=10, range=None, normed=False, weights=None):
smin[i] = smin[i] - .5
smax[i] = smax[i] + .5
+ # avoid rounding issues for comparisons when dealing with inexact types
+ if np.issubdtype(sample.dtype, np.inexact):
+ edge_dt = sample.dtype
+ else:
+ edge_dt = float
# Create edge arrays
for i in arange(D):
if isscalar(bins[i]):
@@ -345,9 +350,9 @@ def histogramdd(sample, bins=10, range=None, normed=False, weights=None):
"Element at index %s in `bins` should be a positive "
"integer." % i)
nbin[i] = bins[i] + 2 # +2 for outlier bins
- edges[i] = linspace(smin[i], smax[i], nbin[i]-1)
+ edges[i] = linspace(smin[i], smax[i], nbin[i]-1, dtype=edge_dt)
else:
- edges[i] = asarray(bins[i], float)
+ edges[i] = asarray(bins[i], edge_dt)
nbin[i] = len(edges[i]) + 1 # +1 for outlier bins
dedges[i] = diff(edges[i])
if np.any(np.asarray(dedges[i]) <= 0):
diff --git a/numpy/lib/index_tricks.py b/numpy/lib/index_tricks.py
index 98c6b291b..f83024961 100644
--- a/numpy/lib/index_tricks.py
+++ b/numpy/lib/index_tricks.py
@@ -727,6 +727,7 @@ def fill_diagonal(a, val, wrap=False):
# tall matrices no wrap
>>> a = np.zeros((5, 3),int)
>>> fill_diagonal(a, 4)
+ >>> a
array([[4, 0, 0],
[0, 4, 0],
[0, 0, 4],
@@ -735,7 +736,8 @@ def fill_diagonal(a, val, wrap=False):
# tall matrices wrap
>>> a = np.zeros((5, 3),int)
- >>> fill_diagonal(a, 4)
+ >>> fill_diagonal(a, 4, wrap=True)
+ >>> a
array([[4, 0, 0],
[0, 4, 0],
[0, 0, 4],
@@ -744,7 +746,8 @@ def fill_diagonal(a, val, wrap=False):
# wide matrices
>>> a = np.zeros((3, 5),int)
- >>> fill_diagonal(a, 4)
+ >>> fill_diagonal(a, 4, wrap=True)
+ >>> a
array([[4, 0, 0, 0, 0],
[0, 4, 0, 0, 0],
[0, 0, 4, 0, 0]])
diff --git a/numpy/lib/nanfunctions.py b/numpy/lib/nanfunctions.py
index f5ac35e54..7260a35b8 100644
--- a/numpy/lib/nanfunctions.py
+++ b/numpy/lib/nanfunctions.py
@@ -33,6 +33,10 @@ def _replace_nan(a, val):
marking the locations where NaNs were present. If `a` is not of
inexact type, do nothing and return `a` together with a mask of None.
+ Note that scalars will end up as array scalars, which is important
+ for using the result as the value of the out argument in some
+ operations.
+
Parameters
----------
a : array-like
@@ -1037,7 +1041,7 @@ def nanvar(a, axis=None, dtype=None, out=None, ddof=0, keepdims=False):
avg = _divide_by_count(avg, cnt)
# Compute squared deviation from mean.
- arr -= avg
+ np.subtract(arr, avg, out=arr, casting='unsafe')
arr = _copyto(arr, 0, mask)
if issubclass(arr.dtype.type, np.complexfloating):
sqr = np.multiply(arr, arr.conj(), out=arr).real
diff --git a/numpy/lib/npyio.py b/numpy/lib/npyio.py
index e6b4d426f..641203f34 100644
--- a/numpy/lib/npyio.py
+++ b/numpy/lib/npyio.py
@@ -288,8 +288,7 @@ def load(file, mmap_mode=None):
Parameters
----------
file : file-like object or string
- The file to read. Compressed files with the filename extension
- ``.gz`` are acceptable. File-like objects must support the
+ The file to read. File-like objects must support the
``seek()`` and ``read()`` methods. Pickled files require that the
file-like object support the ``readline()`` method as well.
mmap_mode : {None, 'r+', 'r', 'w+', 'c'}, optional
@@ -1519,7 +1518,9 @@ def genfromtxt(fname, dtype=float, comments='#', delimiter=None,
# Process the filling_values ...............................
# Rename the input for convenience
- user_filling_values = filling_values or []
+ user_filling_values = filling_values
+ if user_filling_values is None:
+ user_filling_values = []
# Define the default
filling_values = [None] * nbcols
# We have a dictionary : update each entry individually
diff --git a/numpy/lib/setup.py b/numpy/lib/setup.py
index 68d99c33a..62d1dfbb8 100644
--- a/numpy/lib/setup.py
+++ b/numpy/lib/setup.py
@@ -13,7 +13,6 @@ def configuration(parent_package='',top_path=None):
sources=[join('src', '_compiled_base.c')]
)
- config.add_data_dir('benchmarks')
config.add_data_dir('tests')
return config
diff --git a/numpy/lib/tests/test_format.py b/numpy/lib/tests/test_format.py
index eea6f1e5c..c09386789 100644
--- a/numpy/lib/tests/test_format.py
+++ b/numpy/lib/tests/test_format.py
@@ -678,28 +678,28 @@ def test_bad_header():
def test_large_file_support():
from nose import SkipTest
+ if (sys.platform == 'win32' or sys.platform == 'cygwin'):
+ raise SkipTest("Unknown if Windows has sparse filesystems")
# try creating a large sparse file
- with tempfile.NamedTemporaryFile() as tf:
- try:
- # seek past end would work too, but linux truncate somewhat
- # increases the chances that we have a sparse filesystem and can
- # avoid actually writing 5GB
- import subprocess as sp
- sp.check_call(["truncate", "-s", "5368709120", tf.name])
- except:
- raise SkipTest("Could not create 5GB large file")
- # write a small array to the end
- f = open(tf.name, "wb")
+ tf_name = os.path.join(tempdir, 'sparse_file')
+ try:
+ # seek past end would work too, but linux truncate somewhat
+ # increases the chances that we have a sparse filesystem and can
+ # avoid actually writing 5GB
+ import subprocess as sp
+ sp.check_call(["truncate", "-s", "5368709120", tf_name])
+ except:
+ raise SkipTest("Could not create 5GB large file")
+ # write a small array to the end
+ with open(tf_name, "wb") as f:
f.seek(5368709120)
d = np.arange(5)
np.save(f, d)
- f.close()
- # read it back
- f = open(tf.name, "rb")
+ # read it back
+ with open(tf_name, "rb") as f:
f.seek(5368709120)
r = np.load(f)
- f.close()
- assert_array_equal(r, d)
+ assert_array_equal(r, d)
if __name__ == "__main__":
diff --git a/numpy/lib/tests/test_function_base.py b/numpy/lib/tests/test_function_base.py
index 6bc622f23..624b5f3eb 100644
--- a/numpy/lib/tests/test_function_base.py
+++ b/numpy/lib/tests/test_function_base.py
@@ -1080,6 +1080,13 @@ class TestHistogram(TestCase):
h, b = histogram(a, weights=np.ones(10, float))
assert_(issubdtype(h.dtype, float))
+ def test_f32_rounding(self):
+ # gh-4799, check that the rounding of the edges works with float32
+ x = np.array([276.318359 , -69.593948 , 21.329449], dtype=np.float32)
+ y = np.array([5005.689453, 4481.327637, 6010.369629], dtype=np.float32)
+ counts_hist, xedges, yedges = np.histogram2d(x, y, bins=100)
+ assert_equal(counts_hist.sum(), 3.)
+
def test_weights(self):
v = rand(100)
w = np.ones(100) * 5
diff --git a/numpy/lib/tests/test_io.py b/numpy/lib/tests/test_io.py
index 3d02c35ca..4038d6a7f 100644
--- a/numpy/lib/tests/test_io.py
+++ b/numpy/lib/tests/test_io.py
@@ -4,9 +4,7 @@ import sys
import gzip
import os
import threading
-import shutil
-import contextlib
-from tempfile import mkstemp, mkdtemp, NamedTemporaryFile
+from tempfile import mkstemp, NamedTemporaryFile
import time
import warnings
import gc
@@ -24,13 +22,7 @@ from numpy.ma.testutils import (
assert_raises, assert_raises_regex, run_module_suite
)
from numpy.testing import assert_warns, assert_, build_err_msg
-
-
-@contextlib.contextmanager
-def tempdir(change_dir=False):
- tmpdir = mkdtemp()
- yield tmpdir
- shutil.rmtree(tmpdir)
+from numpy.testing.utils import tempdir
class TextIO(BytesIO):
@@ -202,7 +194,7 @@ class TestSavezLoad(RoundtripTest, TestCase):
def test_big_arrays(self):
L = (1 << 31) + 100000
a = np.empty(L, dtype=np.uint8)
- with tempdir() as tmpdir:
+ with tempdir(prefix="numpy_test_big_arrays_") as tmpdir:
tmp = os.path.join(tmpdir, "file.npz")
np.savez(tmp, a=a)
del a
@@ -311,7 +303,7 @@ class TestSavezLoad(RoundtripTest, TestCase):
# Check that zipfile owns file and can close it.
# This needs to pass a file name to load for the
# test.
- with tempdir() as tmpdir:
+ with tempdir(prefix="numpy_test_closing_zipfile_after_load_") as tmpdir:
fd, tmp = mkstemp(suffix='.npz', dir=tmpdir)
os.close(fd)
np.savez(tmp, lab='place holder')
@@ -1323,6 +1315,16 @@ M 33 21.99
ctrl = np.array([(0, 3), (4, -999)], dtype=[(_, int) for _ in "ac"])
assert_equal(test, ctrl)
+ data2 = "1,2,*,4\n5,*,7,8\n"
+ test = np.genfromtxt(TextIO(data2), delimiter=',', dtype=int,
+ missing_values="*", filling_values=0)
+ ctrl = np.array([[1, 2, 0, 4], [5, 0, 7, 8]])
+ assert_equal(test, ctrl)
+ test = np.genfromtxt(TextIO(data2), delimiter=',', dtype=int,
+ missing_values="*", filling_values=-1)
+ ctrl = np.array([[1, 2, -1, 4], [5, -1, 7, 8]])
+ assert_equal(test, ctrl)
+
def test_withmissing_float(self):
data = TextIO('A,B\n0,1.5\n2,-999.00')
test = np.mafromtxt(data, dtype=None, delimiter=',',
diff --git a/numpy/lib/tests/test_twodim_base.py b/numpy/lib/tests/test_twodim_base.py
index e9dbef70f..739061a5d 100644
--- a/numpy/lib/tests/test_twodim_base.py
+++ b/numpy/lib/tests/test_twodim_base.py
@@ -311,6 +311,40 @@ def test_tril_triu_ndim3():
yield assert_equal, a_triu_observed.dtype, a.dtype
yield assert_equal, a_tril_observed.dtype, a.dtype
+def test_tril_triu_with_inf():
+ # Issue 4859
+ arr = np.array([[1, 1, np.inf],
+ [1, 1, 1],
+ [np.inf, 1, 1]])
+ out_tril = np.array([[1, 0, 0],
+ [1, 1, 0],
+ [np.inf, 1, 1]])
+ out_triu = out_tril.T
+ assert_array_equal(np.triu(arr), out_triu)
+ assert_array_equal(np.tril(arr), out_tril)
+
+
+def test_tril_triu_dtype():
+ # Issue 4916
+ # tril and triu should return the same dtype as input
+ for c in np.typecodes['All']:
+ if c == 'V':
+ continue
+ arr = np.zeros((3, 3), dtype=c)
+ assert_equal(np.triu(arr).dtype, arr.dtype)
+ assert_equal(np.tril(arr).dtype, arr.dtype)
+
+ # check special cases
+ arr = np.array([['2001-01-01T12:00', '2002-02-03T13:56'],
+ ['2004-01-01T12:00', '2003-01-03T13:45']],
+ dtype='datetime64')
+ assert_equal(np.triu(arr).dtype, arr.dtype)
+ assert_equal(np.tril(arr).dtype, arr.dtype)
+
+ arr = np.zeros((3,3), dtype='f4,f4')
+ assert_equal(np.triu(arr).dtype, arr.dtype)
+ assert_equal(np.tril(arr).dtype, arr.dtype)
+
def test_mask_indices():
# simple test without offset
diff --git a/numpy/lib/twodim_base.py b/numpy/lib/twodim_base.py
index 2861e1c4a..40a140b6b 100644
--- a/numpy/lib/twodim_base.py
+++ b/numpy/lib/twodim_base.py
@@ -387,7 +387,6 @@ def tri(N, M=None, k=0, dtype=float):
dtype : dtype, optional
Data type of the returned array. The default is float.
-
Returns
-------
tri : ndarray of shape (N, M)
@@ -452,7 +451,9 @@ def tril(m, k=0):
"""
m = asanyarray(m)
- return multiply(tri(*m.shape[-2:], k=k, dtype=bool), m, dtype=m.dtype)
+ mask = tri(*m.shape[-2:], k=k, dtype=bool)
+
+ return where(mask, m, zeros(1, m.dtype))
def triu(m, k=0):
@@ -478,7 +479,9 @@ def triu(m, k=0):
"""
m = asanyarray(m)
- return multiply(~tri(*m.shape[-2:], k=k-1, dtype=bool), m, dtype=m.dtype)
+ mask = tri(*m.shape[-2:], k=k-1, dtype=bool)
+
+ return where(mask, zeros(1, m.dtype), m)
# Originally borrowed from John Hunter and matplotlib