summaryrefslogtreecommitdiff
path: root/numpy/lib
diff options
context:
space:
mode:
Diffstat (limited to 'numpy/lib')
-rw-r--r--numpy/lib/function_base.py13
-rw-r--r--numpy/lib/polynomial.py3
-rw-r--r--numpy/lib/shape_base.py10
-rw-r--r--numpy/lib/tests/test__datasource.py8
-rw-r--r--numpy/lib/tests/test_format.py3
-rw-r--r--numpy/lib/tests/test_function_base.py57
-rw-r--r--numpy/lib/tests/test_packbits.py3
7 files changed, 81 insertions, 16 deletions
diff --git a/numpy/lib/function_base.py b/numpy/lib/function_base.py
index fef69dff3..9261dba22 100644
--- a/numpy/lib/function_base.py
+++ b/numpy/lib/function_base.py
@@ -336,8 +336,12 @@ def histogram(a, bins=10, range=None, normed=False, weights=None,
if (range is not None):
mn, mx = range
if (mn > mx):
- raise AttributeError(
+ raise ValueError(
'max must be larger than min in range parameter.')
+ if not np.all(np.isfinite([mn, mx])):
+ raise ValueError(
+ 'range parameter must be finite.')
+
if isinstance(bins, basestring):
bins = _hist_optim_numbins_estimator(a, bins)
@@ -422,7 +426,7 @@ def histogram(a, bins=10, range=None, normed=False, weights=None,
else:
bins = asarray(bins)
if (np.diff(bins) < 0).any():
- raise AttributeError(
+ raise ValueError(
'bins must increase monotonically.')
# Initialize empty histogram
@@ -533,7 +537,7 @@ def histogramdd(sample, bins=10, range=None, normed=False, weights=None):
try:
M = len(bins)
if M != D:
- raise AttributeError(
+ raise ValueError(
'The dimension of bins must be equal to the dimension of the '
' sample x.')
except TypeError:
@@ -551,6 +555,9 @@ def histogramdd(sample, bins=10, range=None, normed=False, weights=None):
smin = atleast_1d(array(sample.min(0), float))
smax = atleast_1d(array(sample.max(0), float))
else:
+ if not np.all(np.isfinite(range)):
+ raise ValueError(
+ 'range parameter must be finite.')
smin = zeros(D)
smax = zeros(D)
for i in arange(D):
diff --git a/numpy/lib/polynomial.py b/numpy/lib/polynomial.py
index de9376300..2f677438b 100644
--- a/numpy/lib/polynomial.py
+++ b/numpy/lib/polynomial.py
@@ -427,7 +427,8 @@ def polyfit(x, y, deg, rcond=None, full=False, w=None, cov=False):
default) just the coefficients are returned, when True diagnostic
information from the singular value decomposition is also returned.
w : array_like, shape (M,), optional
- weights to apply to the y-coordinates of the sample points.
+ Weights to apply to the y-coordinates of the sample points. For
+ gaussian uncertainties, use 1/sigma (not 1/sigma**2).
cov : bool, optional
Return the estimate and the covariance matrix of the estimate
If full is True, then cov is not returned.
diff --git a/numpy/lib/shape_base.py b/numpy/lib/shape_base.py
index 615cf88f4..ffbe56721 100644
--- a/numpy/lib/shape_base.py
+++ b/numpy/lib/shape_base.py
@@ -799,6 +799,9 @@ def tile(A, reps):
Thus for an `A` of shape (2, 3, 4, 5), a `reps` of (2, 2) is treated as
(1, 1, 2, 2).
+ Note : Although tile may be used for broadcasting, it is strongly
+ recommended to use numpy's broadcasting operations and functions.
+
Parameters
----------
A : array_like
@@ -814,6 +817,7 @@ def tile(A, reps):
See Also
--------
repeat : Repeat elements of an array.
+ broadcast_to : Broadcast an array to a new shape
Examples
--------
@@ -837,6 +841,12 @@ def tile(A, reps):
[1, 2],
[3, 4]])
+ >>> c = np.array([1,2,3,4])
+ >>> np.tile(c,(4,1))
+ array([[1, 2, 3, 4],
+ [1, 2, 3, 4],
+ [1, 2, 3, 4],
+ [1, 2, 3, 4]])
"""
try:
tup = tuple(reps)
diff --git a/numpy/lib/tests/test__datasource.py b/numpy/lib/tests/test__datasource.py
index 090f71f67..f4bece352 100644
--- a/numpy/lib/tests/test__datasource.py
+++ b/numpy/lib/tests/test__datasource.py
@@ -7,7 +7,7 @@ from shutil import rmtree
from numpy.compat import asbytes
from numpy.testing import (
- run_module_suite, TestCase, assert_
+ run_module_suite, TestCase, assert_, SkipTest
)
import numpy.lib._datasource as datasource
@@ -137,8 +137,7 @@ class TestDataSourceOpen(TestCase):
import gzip
except ImportError:
# We don't have the gzip capabilities to test.
- import nose
- raise nose.SkipTest
+ raise SkipTest
# Test datasource's internal file_opener for Gzip files.
filepath = os.path.join(self.tmpdir, 'foobar.txt.gz')
fp = gzip.open(filepath, 'w')
@@ -154,8 +153,7 @@ class TestDataSourceOpen(TestCase):
import bz2
except ImportError:
# We don't have the bz2 capabilities to test.
- import nose
- raise nose.SkipTest
+ raise SkipTest
# Test datasource's internal file_opener for BZip2 files.
filepath = os.path.join(self.tmpdir, 'foobar.txt.bz2')
fp = bz2.BZ2File(filepath, 'w')
diff --git a/numpy/lib/tests/test_format.py b/numpy/lib/tests/test_format.py
index 4f8a65148..1bf65fa61 100644
--- a/numpy/lib/tests/test_format.py
+++ b/numpy/lib/tests/test_format.py
@@ -287,7 +287,7 @@ import numpy as np
from numpy.compat import asbytes, asbytes_nested, sixu
from numpy.testing import (
run_module_suite, assert_, assert_array_equal, assert_raises, raises,
- dec
+ dec, SkipTest
)
from numpy.lib import format
@@ -812,7 +812,6 @@ 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
diff --git a/numpy/lib/tests/test_function_base.py b/numpy/lib/tests/test_function_base.py
index cc53c2b8e..a5ac78e33 100644
--- a/numpy/lib/tests/test_function_base.py
+++ b/numpy/lib/tests/test_function_base.py
@@ -1267,6 +1267,13 @@ class TestHistogram(TestCase):
assert_array_equal(a, np.array([0]))
assert_array_equal(b, np.array([0, 1]))
+ def test_finite_range(self):
+ # Normal ranges should be fine
+ vals = np.linspace(0.0, 1.0, num=100)
+ histogram(vals, range=[0.25,0.75])
+ assert_raises(ValueError, histogram, vals, range=[np.nan,0.75])
+ assert_raises(ValueError, histogram, vals, range=[0.25,np.inf])
+
class TestHistogramOptimBinNums(TestCase):
"""
@@ -1489,6 +1496,16 @@ class TestHistogramdd(TestCase):
assert_(hist[0] == 0.0)
assert_(hist[1] == 0.0)
+ def test_finite_range(self):
+ vals = np.random.random((100,3))
+ histogramdd(vals, range=[[0.0,1.0],[0.25,0.75],[0.25,0.5]])
+ assert_raises(ValueError, histogramdd, vals,
+ range=[[0.0,1.0],[0.25,0.75],[0.25,np.inf]])
+ assert_raises(ValueError, histogramdd, vals,
+ range=[[0.0,1.0],[np.nan,0.75],[0.25,0.5]])
+
+
+
class TestUnique(TestCase):
@@ -1957,10 +1974,42 @@ class TestInterp(TestCase):
assert_almost_equal(np.interp(x0, x, y), x0)
def test_right_left_behavior(self):
- assert_equal(interp([-1, 0, 1], [0], [1]), [1, 1, 1])
- assert_equal(interp([-1, 0, 1], [0], [1], left=0), [0, 1, 1])
- assert_equal(interp([-1, 0, 1], [0], [1], right=0), [1, 1, 0])
- assert_equal(interp([-1, 0, 1], [0], [1], left=0, right=0), [0, 1, 0])
+ # Needs range of sizes to test different code paths.
+ # size ==1 is special cased, 1 < size < 5 is linear search, and
+ # size >= 5 goes through local search and possibly binary search.
+ for size in range(1, 10):
+ xp = np.arange(size, dtype=np.double)
+ yp = np.ones(size, dtype=np.double)
+ incpts = np.array([-1, 0, size - 1, size], dtype=np.double)
+ decpts = incpts[::-1]
+
+ incres = interp(incpts, xp, yp)
+ decres = interp(decpts, xp, yp)
+ inctgt = np.array([1, 1, 1, 1], dtype=np.float)
+ dectgt = inctgt[::-1]
+ assert_equal(incres, inctgt)
+ assert_equal(decres, dectgt)
+
+ incres = interp(incpts, xp, yp, left=0)
+ decres = interp(decpts, xp, yp, left=0)
+ inctgt = np.array([0, 1, 1, 1], dtype=np.float)
+ dectgt = inctgt[::-1]
+ assert_equal(incres, inctgt)
+ assert_equal(decres, dectgt)
+
+ incres = interp(incpts, xp, yp, right=2)
+ decres = interp(decpts, xp, yp, right=2)
+ inctgt = np.array([1, 1, 1, 2], dtype=np.float)
+ dectgt = inctgt[::-1]
+ assert_equal(incres, inctgt)
+ assert_equal(decres, dectgt)
+
+ incres = interp(incpts, xp, yp, left=0, right=2)
+ decres = interp(decpts, xp, yp, left=0, right=2)
+ inctgt = np.array([0, 1, 1, 2], dtype=np.float)
+ dectgt = inctgt[::-1]
+ assert_equal(incres, inctgt)
+ assert_equal(decres, dectgt)
def test_scalar_interpolation_point(self):
x = np.linspace(0, 1, 5)
diff --git a/numpy/lib/tests/test_packbits.py b/numpy/lib/tests/test_packbits.py
index 186e8960d..0de084ef9 100644
--- a/numpy/lib/tests/test_packbits.py
+++ b/numpy/lib/tests/test_packbits.py
@@ -1,5 +1,6 @@
-import numpy as np
+from __future__ import division, absolute_import, print_function
+import numpy as np
from numpy.testing import assert_array_equal, assert_equal, assert_raises