summaryrefslogtreecommitdiff
path: root/numpy/core
diff options
context:
space:
mode:
authorMarten van Kerkwijk <mhvk@astro.utoronto.ca>2017-02-21 15:58:08 -0500
committerGitHub <noreply@github.com>2017-02-21 15:58:08 -0500
commit2aabeafb97bea4e1bfa29d946fbf31e1104e7ae0 (patch)
tree2cd08a2211a3ec1f7403c17dd175aca73a93bccb /numpy/core
parent070b9660282288fa8bb376533667f31613373337 (diff)
parent8d6ec65c925ebef5e0567708de1d16df39077c9d (diff)
downloadnumpy-2aabeafb97bea4e1bfa29d946fbf31e1104e7ae0.tar.gz
Merge pull request #8584 from eric-wieser/resolve_axis
MAINT: Use the same exception for all bad axis requests
Diffstat (limited to 'numpy/core')
-rw-r--r--numpy/core/_internal.py3
-rw-r--r--numpy/core/numeric.py15
-rw-r--r--numpy/core/shape_base.py7
-rw-r--r--numpy/core/src/multiarray/common.h36
-rw-r--r--numpy/core/src/multiarray/conversion_utils.c20
-rw-r--r--numpy/core/src/multiarray/ctors.c8
-rw-r--r--numpy/core/src/multiarray/item_selection.c25
-rw-r--r--numpy/core/src/multiarray/multiarraymodule.c29
-rw-r--r--numpy/core/src/multiarray/shape.c7
-rw-r--r--numpy/core/src/umath/ufunc_object.c16
-rw-r--r--numpy/core/tests/test_multiarray.py16
-rw-r--r--numpy/core/tests/test_numeric.py16
-rw-r--r--numpy/core/tests/test_shape_base.py8
-rw-r--r--numpy/core/tests/test_ufunc.py8
14 files changed, 105 insertions, 109 deletions
diff --git a/numpy/core/_internal.py b/numpy/core/_internal.py
index 741c8bb5f..d73cdcc55 100644
--- a/numpy/core/_internal.py
+++ b/numpy/core/_internal.py
@@ -630,3 +630,6 @@ def _gcd(a, b):
# Exception used in shares_memory()
class TooHardError(RuntimeError):
pass
+
+class AxisError(ValueError, IndexError):
+ pass
diff --git a/numpy/core/numeric.py b/numpy/core/numeric.py
index 97d19f008..066697f3e 100644
--- a/numpy/core/numeric.py
+++ b/numpy/core/numeric.py
@@ -17,7 +17,7 @@ from .multiarray import (
inner, int_asbuffer, lexsort, matmul, may_share_memory,
min_scalar_type, ndarray, nditer, nested_iters, promote_types,
putmask, result_type, set_numeric_ops, shares_memory, vdot, where,
- zeros)
+ zeros, normalize_axis_index)
if sys.version_info[0] < 3:
from .multiarray import newbuffer, getbuffer
@@ -27,7 +27,7 @@ from .umath import (invert, sin, UFUNC_BUFSIZE_DEFAULT, ERR_IGNORE,
ERR_DEFAULT, PINF, NAN)
from . import numerictypes
from .numerictypes import longlong, intc, int_, float_, complex_, bool_
-from ._internal import TooHardError
+from ._internal import TooHardError, AxisError
bitwise_not = invert
ufunc = type(sin)
@@ -65,7 +65,7 @@ __all__ = [
'True_', 'bitwise_not', 'CLIP', 'RAISE', 'WRAP', 'MAXDIMS', 'BUFSIZE',
'ALLOW_THREADS', 'ComplexWarning', 'full', 'full_like', 'matmul',
'shares_memory', 'may_share_memory', 'MAY_SHARE_BOUNDS', 'MAY_SHARE_EXACT',
- 'TooHardError',
+ 'TooHardError', 'AxisError'
]
@@ -1527,15 +1527,12 @@ def rollaxis(a, axis, start=0):
"""
n = a.ndim
- if axis < 0:
- axis += n
+ axis = normalize_axis_index(axis, n)
if start < 0:
start += n
msg = "'%s' arg requires %d <= %s < %d, but %d was passed in"
- if not (0 <= axis < n):
- raise ValueError(msg % ('axis', -n, 'axis', n, axis))
if not (0 <= start < n + 1):
- raise ValueError(msg % ('start', -n, 'start', n + 1, start))
+ raise AxisError(msg % ('start', -n, 'start', n + 1, start))
if axis < start:
# it's been removed
start -= 1
@@ -1554,7 +1551,7 @@ def _validate_axis(axis, ndim, argname):
axis = list(axis)
axis = [a + ndim if a < 0 else a for a in axis]
if not builtins.all(0 <= a < ndim for a in axis):
- raise ValueError('invalid axis for this array in `%s` argument' %
+ raise AxisError('invalid axis for this array in `%s` argument' %
argname)
if len(set(axis)) != len(axis):
raise ValueError('repeated axis in `%s` argument' % argname)
diff --git a/numpy/core/shape_base.py b/numpy/core/shape_base.py
index 70afdb746..58b0dcaac 100644
--- a/numpy/core/shape_base.py
+++ b/numpy/core/shape_base.py
@@ -5,6 +5,7 @@ __all__ = ['atleast_1d', 'atleast_2d', 'atleast_3d', 'vstack', 'hstack',
from . import numeric as _nx
from .numeric import asanyarray, newaxis
+from .multiarray import normalize_axis_index
def atleast_1d(*arys):
"""
@@ -347,11 +348,7 @@ def stack(arrays, axis=0):
raise ValueError('all input arrays must have the same shape')
result_ndim = arrays[0].ndim + 1
- if not -result_ndim <= axis < result_ndim:
- msg = 'axis {0} out of bounds [-{1}, {1})'.format(axis, result_ndim)
- raise IndexError(msg)
- if axis < 0:
- axis += result_ndim
+ axis = normalize_axis_index(axis, result_ndim)
sl = (slice(None),) * axis + (_nx.newaxis,)
expanded_arrays = [arr[sl] for arr in arrays]
diff --git a/numpy/core/src/multiarray/common.h b/numpy/core/src/multiarray/common.h
index 5e14b80a7..625ca9d76 100644
--- a/numpy/core/src/multiarray/common.h
+++ b/numpy/core/src/multiarray/common.h
@@ -134,6 +134,42 @@ check_and_adjust_index(npy_intp *index, npy_intp max_item, int axis,
return 0;
}
+/*
+ * Returns -1 and sets an exception if *axis is an invalid axis for
+ * an array of dimension ndim, otherwise adjusts it in place to be
+ * 0 <= *axis < ndim, and returns 0.
+ */
+static NPY_INLINE int
+check_and_adjust_axis(int *axis, int ndim)
+{
+ /* Check that index is valid, taking into account negative indices */
+ if (NPY_UNLIKELY((*axis < -ndim) || (*axis >= ndim))) {
+ /*
+ * Load the exception type, if we don't already have it. Unfortunately
+ * we don't have access to npy_cache_import here
+ */
+ static PyObject *AxisError_cls = NULL;
+ if (AxisError_cls == NULL) {
+ PyObject *mod = PyImport_ImportModule("numpy.core._internal");
+
+ if (mod != NULL) {
+ AxisError_cls = PyObject_GetAttrString(mod, "AxisError");
+ Py_DECREF(mod);
+ }
+ }
+
+ PyErr_Format(AxisError_cls,
+ "axis %d is out of bounds for array of dimension %d",
+ *axis, ndim);
+ return -1;
+ }
+ /* adjust negative indices */
+ if (*axis < 0) {
+ *axis += ndim;
+ }
+ return 0;
+}
+
/*
* return true if pointer is aligned to 'alignment'
diff --git a/numpy/core/src/multiarray/conversion_utils.c b/numpy/core/src/multiarray/conversion_utils.c
index c016bb8d1..8ed08a366 100644
--- a/numpy/core/src/multiarray/conversion_utils.c
+++ b/numpy/core/src/multiarray/conversion_utils.c
@@ -259,17 +259,10 @@ PyArray_ConvertMultiAxis(PyObject *axis_in, int ndim, npy_bool *out_axis_flags)
PyObject *tmp = PyTuple_GET_ITEM(axis_in, i);
int axis = PyArray_PyIntAsInt_ErrMsg(tmp,
"integers are required for the axis tuple elements");
- int axis_orig = axis;
if (error_converting(axis)) {
return NPY_FAIL;
}
- if (axis < 0) {
- axis += ndim;
- }
- if (axis < 0 || axis >= ndim) {
- PyErr_Format(PyExc_ValueError,
- "'axis' entry %d is out of bounds [-%d, %d)",
- axis_orig, ndim, ndim);
+ if (check_and_adjust_axis(&axis, ndim) < 0) {
return NPY_FAIL;
}
if (out_axis_flags[axis]) {
@@ -284,20 +277,16 @@ PyArray_ConvertMultiAxis(PyObject *axis_in, int ndim, npy_bool *out_axis_flags)
}
/* Try to interpret axis as an integer */
else {
- int axis, axis_orig;
+ int axis;
memset(out_axis_flags, 0, ndim);
axis = PyArray_PyIntAsInt_ErrMsg(axis_in,
"an integer is required for the axis");
- axis_orig = axis;
if (error_converting(axis)) {
return NPY_FAIL;
}
- if (axis < 0) {
- axis += ndim;
- }
/*
* Special case letting axis={-1,0} slip through for scalars,
* for backwards compatibility reasons.
@@ -306,10 +295,7 @@ PyArray_ConvertMultiAxis(PyObject *axis_in, int ndim, npy_bool *out_axis_flags)
return NPY_SUCCEED;
}
- if (axis < 0 || axis >= ndim) {
- PyErr_Format(PyExc_ValueError,
- "'axis' entry %d is out of bounds [-%d, %d)",
- axis_orig, ndim, ndim);
+ if (check_and_adjust_axis(&axis, ndim) < 0) {
return NPY_FAIL;
}
diff --git a/numpy/core/src/multiarray/ctors.c b/numpy/core/src/multiarray/ctors.c
index 349b59c5f..ee6b66eef 100644
--- a/numpy/core/src/multiarray/ctors.c
+++ b/numpy/core/src/multiarray/ctors.c
@@ -2793,7 +2793,6 @@ PyArray_CheckAxis(PyArrayObject *arr, int *axis, int flags)
{
PyObject *temp1, *temp2;
int n = PyArray_NDIM(arr);
- int axis_orig = *axis;
if (*axis == NPY_MAXDIMS || n == 0) {
if (n != 1) {
@@ -2831,12 +2830,7 @@ PyArray_CheckAxis(PyArrayObject *arr, int *axis, int flags)
temp2 = (PyObject *)temp1;
}
n = PyArray_NDIM((PyArrayObject *)temp2);
- if (*axis < 0) {
- *axis += n;
- }
- if ((*axis < 0) || (*axis >= n)) {
- PyErr_Format(PyExc_ValueError,
- "axis(=%d) out of bounds", axis_orig);
+ if (check_and_adjust_axis(axis, n) < 0) {
Py_DECREF(temp2);
return NULL;
}
diff --git a/numpy/core/src/multiarray/item_selection.c b/numpy/core/src/multiarray/item_selection.c
index 08b9c5965..3c0f0782e 100644
--- a/numpy/core/src/multiarray/item_selection.c
+++ b/numpy/core/src/multiarray/item_selection.c
@@ -1101,16 +1101,12 @@ NPY_NO_EXPORT int
PyArray_Sort(PyArrayObject *op, int axis, NPY_SORTKIND which)
{
PyArray_SortFunc *sort;
- int axis_orig = axis;
- int n = PyArray_NDIM(op);
+ int n = PyArray_NDIM(op);
- if (axis < 0) {
- axis += n;
- }
- if (axis < 0 || axis >= n) {
- PyErr_Format(PyExc_ValueError, "axis(=%d) out of bounds", axis_orig);
+ if (check_and_adjust_axis(&axis, n) < 0) {
return -1;
}
+
if (PyArray_FailUnlessWriteable(op, "sort array") < 0) {
return -1;
}
@@ -1212,17 +1208,13 @@ PyArray_Partition(PyArrayObject *op, PyArrayObject * ktharray, int axis,
PyArrayObject *kthrvl;
PyArray_PartitionFunc *part;
PyArray_SortFunc *sort;
- int axis_orig = axis;
int n = PyArray_NDIM(op);
int ret;
- if (axis < 0) {
- axis += n;
- }
- if (axis < 0 || axis >= n) {
- PyErr_Format(PyExc_ValueError, "axis(=%d) out of bounds", axis_orig);
+ if (check_and_adjust_axis(&axis, n) < 0) {
return -1;
}
+
if (PyArray_FailUnlessWriteable(op, "partition array") < 0) {
return -1;
}
@@ -1455,12 +1447,7 @@ PyArray_LexSort(PyObject *sort_keys, int axis)
*((npy_intp *)(PyArray_DATA(ret))) = 0;
goto finish;
}
- if (axis < 0) {
- axis += nd;
- }
- if ((axis < 0) || (axis >= nd)) {
- PyErr_Format(PyExc_ValueError,
- "axis(=%d) out of bounds", axis);
+ if (check_and_adjust_axis(&axis, nd) < 0) {
goto fail;
}
diff --git a/numpy/core/src/multiarray/multiarraymodule.c b/numpy/core/src/multiarray/multiarraymodule.c
index f00de46c4..1c8d9b5e4 100644
--- a/numpy/core/src/multiarray/multiarraymodule.c
+++ b/numpy/core/src/multiarray/multiarraymodule.c
@@ -327,7 +327,6 @@ PyArray_ConcatenateArrays(int narrays, PyArrayObject **arrays, int axis)
PyArray_Descr *dtype = NULL;
PyArrayObject *ret = NULL;
PyArrayObject_fields *sliding_view = NULL;
- int orig_axis = axis;
if (narrays <= 0) {
PyErr_SetString(PyExc_ValueError,
@@ -345,13 +344,7 @@ PyArray_ConcatenateArrays(int narrays, PyArrayObject **arrays, int axis)
}
/* Handle standard Python negative indexing */
- if (axis < 0) {
- axis += ndim;
- }
-
- if (axis < 0 || axis >= ndim) {
- PyErr_Format(PyExc_IndexError,
- "axis %d out of bounds [0, %d)", orig_axis, ndim);
+ if (check_and_adjust_axis(&axis, ndim) < 0) {
return NULL;
}
@@ -4109,6 +4102,24 @@ array_may_share_memory(PyObject *NPY_UNUSED(ignored), PyObject *args, PyObject *
return array_shares_memory_impl(args, kwds, NPY_MAY_SHARE_BOUNDS, 0);
}
+static PyObject *
+normalize_axis_index(PyObject *NPY_UNUSED(self), PyObject *args, PyObject *kwds)
+{
+ static char *kwlist[] = {"axis", "ndim", NULL};
+ int axis;
+ int ndim;
+
+ if (!PyArg_ParseTupleAndKeywords(args, kwds, "ii", kwlist,
+ &axis, &ndim)) {
+ return NULL;
+ }
+
+ if(check_and_adjust_axis(&axis, ndim) < 0) {
+ return NULL;
+ }
+
+ return PyInt_FromLong(axis);
+}
static struct PyMethodDef array_module_methods[] = {
{"_get_ndarray_c_version",
@@ -4284,6 +4295,8 @@ static struct PyMethodDef array_module_methods[] = {
METH_VARARGS | METH_KEYWORDS, NULL},
{"unpackbits", (PyCFunction)io_unpack,
METH_VARARGS | METH_KEYWORDS, NULL},
+ {"normalize_axis_index", (PyCFunction)normalize_axis_index,
+ METH_VARARGS | METH_KEYWORDS, NULL},
{NULL, NULL, 0, NULL} /* sentinel */
};
diff --git a/numpy/core/src/multiarray/shape.c b/numpy/core/src/multiarray/shape.c
index 3bee562be..5207513bf 100644
--- a/numpy/core/src/multiarray/shape.c
+++ b/numpy/core/src/multiarray/shape.c
@@ -705,12 +705,7 @@ PyArray_Transpose(PyArrayObject *ap, PyArray_Dims *permute)
}
for (i = 0; i < n; i++) {
axis = axes[i];
- if (axis < 0) {
- axis = PyArray_NDIM(ap) + axis;
- }
- if (axis < 0 || axis >= PyArray_NDIM(ap)) {
- PyErr_SetString(PyExc_ValueError,
- "invalid axis for this array");
+ if (check_and_adjust_axis(&axis, PyArray_NDIM(ap)) < 0) {
return NULL;
}
if (reverse_permutation[axis] != -1) {
diff --git a/numpy/core/src/umath/ufunc_object.c b/numpy/core/src/umath/ufunc_object.c
index 0bae2d591..af4ce12db 100644
--- a/numpy/core/src/umath/ufunc_object.c
+++ b/numpy/core/src/umath/ufunc_object.c
@@ -4036,12 +4036,7 @@ PyUFunc_GenericReduction(PyUFuncObject *ufunc, PyObject *args,
Py_DECREF(mp);
return NULL;
}
- if (axis < 0) {
- axis += ndim;
- }
- if (axis < 0 || axis >= ndim) {
- PyErr_SetString(PyExc_ValueError,
- "'axis' entry is out of bounds");
+ if (check_and_adjust_axis(&axis, ndim) < 0) {
Py_XDECREF(otype);
Py_DECREF(mp);
return NULL;
@@ -4058,18 +4053,11 @@ PyUFunc_GenericReduction(PyUFuncObject *ufunc, PyObject *args,
Py_DECREF(mp);
return NULL;
}
- if (axis < 0) {
- axis += ndim;
- }
/* Special case letting axis={0 or -1} slip through for scalars */
if (ndim == 0 && (axis == 0 || axis == -1)) {
axis = 0;
}
- else if (axis < 0 || axis >= ndim) {
- PyErr_SetString(PyExc_ValueError,
- "'axis' entry is out of bounds");
- Py_XDECREF(otype);
- Py_DECREF(mp);
+ else if (check_and_adjust_axis(&axis, ndim) < 0) {
return NULL;
}
axes[0] = (int)axis;
diff --git a/numpy/core/tests/test_multiarray.py b/numpy/core/tests/test_multiarray.py
index 8229f1e1a..fa5051ba7 100644
--- a/numpy/core/tests/test_multiarray.py
+++ b/numpy/core/tests/test_multiarray.py
@@ -2013,13 +2013,13 @@ class TestMethods(TestCase):
d = np.array([2, 1])
d.partition(0, kind=k)
assert_raises(ValueError, d.partition, 2)
- assert_raises(ValueError, d.partition, 3, axis=1)
+ assert_raises(np.AxisError, d.partition, 3, axis=1)
assert_raises(ValueError, np.partition, d, 2)
- assert_raises(ValueError, np.partition, d, 2, axis=1)
+ assert_raises(np.AxisError, np.partition, d, 2, axis=1)
assert_raises(ValueError, d.argpartition, 2)
- assert_raises(ValueError, d.argpartition, 3, axis=1)
+ assert_raises(np.AxisError, d.argpartition, 3, axis=1)
assert_raises(ValueError, np.argpartition, d, 2)
- assert_raises(ValueError, np.argpartition, d, 2, axis=1)
+ assert_raises(np.AxisError, np.argpartition, d, 2, axis=1)
d = np.arange(10).reshape((2, 5))
d.partition(1, axis=0, kind=k)
d.partition(4, axis=1, kind=k)
@@ -3522,8 +3522,8 @@ class TestArgmin(TestCase):
class TestMinMax(TestCase):
def test_scalar(self):
- assert_raises(ValueError, np.amax, 1, 1)
- assert_raises(ValueError, np.amin, 1, 1)
+ assert_raises(np.AxisError, np.amax, 1, 1)
+ assert_raises(np.AxisError, np.amin, 1, 1)
assert_equal(np.amax(1, axis=0), 1)
assert_equal(np.amin(1, axis=0), 1)
@@ -3531,7 +3531,7 @@ class TestMinMax(TestCase):
assert_equal(np.amin(1, axis=None), 1)
def test_axis(self):
- assert_raises(ValueError, np.amax, [1, 2, 3], 1000)
+ assert_raises(np.AxisError, np.amax, [1, 2, 3], 1000)
assert_equal(np.amax([[1, 2, 3]], axis=1), 3)
def test_datetime(self):
@@ -3793,7 +3793,7 @@ class TestLexsort(TestCase):
def test_invalid_axis(self): # gh-7528
x = np.linspace(0., 1., 42*3).reshape(42, 3)
- assert_raises(ValueError, np.lexsort, x, axis=2)
+ assert_raises(np.AxisError, np.lexsort, x, axis=2)
class TestIO(object):
"""Test tofile, fromfile, tobytes, and fromstring"""
diff --git a/numpy/core/tests/test_numeric.py b/numpy/core/tests/test_numeric.py
index 4aa6bed33..906280e15 100644
--- a/numpy/core/tests/test_numeric.py
+++ b/numpy/core/tests/test_numeric.py
@@ -1010,7 +1010,7 @@ class TestNonzero(TestCase):
assert_raises(ValueError, np.count_nonzero, m, axis=(1, 1))
assert_raises(TypeError, np.count_nonzero, m, axis='foo')
- assert_raises(ValueError, np.count_nonzero, m, axis=3)
+ assert_raises(np.AxisError, np.count_nonzero, m, axis=3)
assert_raises(TypeError, np.count_nonzero,
m, axis=np.array([[1], [2]]))
@@ -2323,10 +2323,10 @@ class TestRollaxis(TestCase):
def test_exceptions(self):
a = np.arange(1*2*3*4).reshape(1, 2, 3, 4)
- assert_raises(ValueError, np.rollaxis, a, -5, 0)
- assert_raises(ValueError, np.rollaxis, a, 0, -5)
- assert_raises(ValueError, np.rollaxis, a, 4, 0)
- assert_raises(ValueError, np.rollaxis, a, 0, 5)
+ assert_raises(np.AxisError, np.rollaxis, a, -5, 0)
+ assert_raises(np.AxisError, np.rollaxis, a, 0, -5)
+ assert_raises(np.AxisError, np.rollaxis, a, 4, 0)
+ assert_raises(np.AxisError, np.rollaxis, a, 0, 5)
def test_results(self):
a = np.arange(1*2*3*4).reshape(1, 2, 3, 4).copy()
@@ -2413,11 +2413,11 @@ class TestMoveaxis(TestCase):
def test_errors(self):
x = np.random.randn(1, 2, 3)
- assert_raises_regex(ValueError, 'invalid axis .* `source`',
+ assert_raises_regex(np.AxisError, 'invalid axis .* `source`',
np.moveaxis, x, 3, 0)
- assert_raises_regex(ValueError, 'invalid axis .* `source`',
+ assert_raises_regex(np.AxisError, 'invalid axis .* `source`',
np.moveaxis, x, -4, 0)
- assert_raises_regex(ValueError, 'invalid axis .* `destination`',
+ assert_raises_regex(np.AxisError, 'invalid axis .* `destination`',
np.moveaxis, x, 0, 5)
assert_raises_regex(ValueError, 'repeated axis in `source`',
np.moveaxis, x, [0, 0], [0, 1])
diff --git a/numpy/core/tests/test_shape_base.py b/numpy/core/tests/test_shape_base.py
index ac8dc1eea..727608a17 100644
--- a/numpy/core/tests/test_shape_base.py
+++ b/numpy/core/tests/test_shape_base.py
@@ -184,8 +184,8 @@ class TestConcatenate(TestCase):
for ndim in [1, 2, 3]:
a = np.ones((1,)*ndim)
np.concatenate((a, a), axis=0) # OK
- assert_raises(IndexError, np.concatenate, (a, a), axis=ndim)
- assert_raises(IndexError, np.concatenate, (a, a), axis=-(ndim + 1))
+ assert_raises(np.AxisError, np.concatenate, (a, a), axis=ndim)
+ assert_raises(np.AxisError, np.concatenate, (a, a), axis=-(ndim + 1))
# Scalars cannot be concatenated
assert_raises(ValueError, concatenate, (0,))
@@ -294,8 +294,8 @@ def test_stack():
expected_shapes = [(10, 3), (3, 10), (3, 10), (10, 3)]
for axis, expected_shape in zip(axes, expected_shapes):
assert_equal(np.stack(arrays, axis).shape, expected_shape)
- assert_raises_regex(IndexError, 'out of bounds', stack, arrays, axis=2)
- assert_raises_regex(IndexError, 'out of bounds', stack, arrays, axis=-3)
+ assert_raises_regex(np.AxisError, 'out of bounds', stack, arrays, axis=2)
+ assert_raises_regex(np.AxisError, 'out of bounds', stack, arrays, axis=-3)
# all shapes for 2d input
arrays = [np.random.randn(3, 4) for _ in range(10)]
axes = [0, 1, 2, -1, -2, -3]
diff --git a/numpy/core/tests/test_ufunc.py b/numpy/core/tests/test_ufunc.py
index 3fea68700..f7b66f90c 100644
--- a/numpy/core/tests/test_ufunc.py
+++ b/numpy/core/tests/test_ufunc.py
@@ -703,14 +703,14 @@ class TestUfunc(TestCase):
def test_axis_out_of_bounds(self):
a = np.array([False, False])
- assert_raises(ValueError, a.all, axis=1)
+ assert_raises(np.AxisError, a.all, axis=1)
a = np.array([False, False])
- assert_raises(ValueError, a.all, axis=-2)
+ assert_raises(np.AxisError, a.all, axis=-2)
a = np.array([False, False])
- assert_raises(ValueError, a.any, axis=1)
+ assert_raises(np.AxisError, a.any, axis=1)
a = np.array([False, False])
- assert_raises(ValueError, a.any, axis=-2)
+ assert_raises(np.AxisError, a.any, axis=-2)
def test_scalar_reduction(self):
# The functions 'sum', 'prod', etc allow specifying axis=0