summaryrefslogtreecommitdiff
path: root/numpy/core
diff options
context:
space:
mode:
Diffstat (limited to 'numpy/core')
-rw-r--r--numpy/core/numeric.py89
-rw-r--r--numpy/core/src/multiarray/multiarraymodule.c8
-rw-r--r--numpy/core/tests/test_numeric.py112
3 files changed, 198 insertions, 11 deletions
diff --git a/numpy/core/numeric.py b/numpy/core/numeric.py
index b3eed9714..8db4e1302 100644
--- a/numpy/core/numeric.py
+++ b/numpy/core/numeric.py
@@ -6,6 +6,7 @@ import operator
import sys
import warnings
+import numpy as np
from . import multiarray
from .multiarray import (
_fastCopyAndTranspose as fastCopyAndTranspose, ALLOW_THREADS,
@@ -376,6 +377,89 @@ def extend_all(module):
__all__.append(a)
+def count_nonzero(a, axis=None):
+ """
+ Counts the number of non-zero values in the array ``a``.
+
+ The word "non-zero" is in reference to the Python 2.x
+ built-in method ``__nonzero__()`` (renamed ``__bool__()``
+ in Python 3.x) of Python objects that tests an object's
+ "truthfulness". For example, any number is considered
+ truthful if it is nonzero, whereas any string is considered
+ truthful if it is not the empty string. Thus, this function
+ (recursively) counts how many elements in ``a`` (and in
+ sub-arrays thereof) have their ``__nonzero__()`` or ``__bool__()``
+ method evaluated to ``True``.
+
+ Parameters
+ ----------
+ a : array_like
+ The array for which to count non-zeros.
+ axis : int or tuple, optional
+ Axis or tuple of axes along which to count non-zeros.
+ Default is None, meaning that non-zeros will be counted
+ along a flattened version of ``a``.
+
+ .. versionadded:: 1.12.0
+
+ Returns
+ -------
+ count : int or array of int
+ Number of non-zero values in the array along a given axis.
+ Otherwise, the total number of non-zero values in the array
+ is returned.
+
+ See Also
+ --------
+ nonzero : Return the coordinates of all the non-zero values.
+
+ Examples
+ --------
+ >>> np.count_nonzero(np.eye(4))
+ 4
+ >>> np.count_nonzero([[0,1,7,0,0],[3,0,0,2,19]])
+ 5
+ >>> np.count_nonzero([[0,1,7,0,0],[3,0,0,2,19]], axis=0)
+ array([1, 1, 1, 1, 1])
+ >>> np.count_nonzero([[0,1,7,0,0],[3,0,0,2,19]], axis=1)
+ array([2, 3])
+
+ """
+ if axis is None or axis == ():
+ return multiarray.count_nonzero(a)
+
+ a = asanyarray(a)
+
+ if a.dtype == bool:
+ return a.sum(axis=axis, dtype=np.intp)
+
+ if issubdtype(a.dtype, np.number):
+ return (a != 0).sum(axis=axis, dtype=np.intp)
+
+ if (issubdtype(a.dtype, np.string_) or
+ issubdtype(a.dtype, np.unicode_)):
+ nullstr = a.dtype.type('')
+ return (a != nullstr).sum(axis=axis, dtype=np.intp)
+
+ axis = asarray(_validate_axis(axis, a.ndim, 'axis'))
+ counts = np.apply_along_axis(multiarray.count_nonzero, axis[0], a)
+
+ if axis.size == 1:
+ return counts
+ else:
+ # for subsequent axis numbers, that number decreases
+ # by one in this new 'counts' array if it was larger
+ # than the first axis upon which 'count_nonzero' was
+ # applied but remains unchanged if that number was
+ # smaller than that first axis
+ #
+ # this trick enables us to perform counts on object-like
+ # elements across multiple axes very quickly because integer
+ # addition is very well optimized
+ return counts.sum(axis=tuple(axis[1:] - (
+ axis[1:] > axis[0])), dtype=np.intp)
+
+
def asarray(a, dtype=None, order=None):
"""Convert the input to an array.
@@ -891,7 +975,7 @@ def correlate(a, v, mode='valid'):
return multiarray.correlate2(a, v, mode)
-def convolve(a,v,mode='full'):
+def convolve(a, v, mode='full'):
"""
Returns the discrete, linear convolution of two one-dimensional sequences.
@@ -1752,7 +1836,7 @@ def cross(a, b, axisa=-1, axisb=-1, axisc=-1, axis=None):
return rollaxis(cp, -1, axisc)
-#Use numarray's printing function
+# Use numarray's printing function
from .arrayprint import array2string, get_printoptions, set_printoptions
@@ -2283,6 +2367,7 @@ def load(file):
# These are all essentially abbreviations
# These might wind up in a special abbreviations module
+
def _maketup(descr, val):
dt = dtype(descr)
# Place val in all scalar tuples:
diff --git a/numpy/core/src/multiarray/multiarraymodule.c b/numpy/core/src/multiarray/multiarraymodule.c
index 62b562856..d90debc6f 100644
--- a/numpy/core/src/multiarray/multiarraymodule.c
+++ b/numpy/core/src/multiarray/multiarraymodule.c
@@ -1980,16 +1980,10 @@ fail:
static PyObject *
array_count_nonzero(PyObject *NPY_UNUSED(self), PyObject *args, PyObject *kwds)
{
- PyObject *array_in;
PyArrayObject *array;
npy_intp count;
- if (!PyArg_ParseTuple(args, "O", &array_in)) {
- return NULL;
- }
-
- array = (PyArrayObject *)PyArray_FromAny(array_in, NULL, 0, 0, 0, NULL);
- if (array == NULL) {
+ if (!PyArg_ParseTuple(args, "O&", PyArray_Converter, &array)) {
return NULL;
}
diff --git a/numpy/core/tests/test_numeric.py b/numpy/core/tests/test_numeric.py
index 7fdbbc930..7830c5c5a 100644
--- a/numpy/core/tests/test_numeric.py
+++ b/numpy/core/tests/test_numeric.py
@@ -67,6 +67,13 @@ class TestNonarrayArgs(TestCase):
out = np.compress([0, 1], arr, axis=0)
assert_equal(out, tgt)
+ def test_count_nonzero(self):
+ arr = [[0, 1, 7, 0, 0],
+ [3, 0, 0, 2, 19]]
+ tgt = np.array([2, 3])
+ out = np.count_nonzero(arr, axis=1)
+ assert_equal(out, tgt)
+
def test_cumproduct(self):
A = [[1, 2, 3], [4, 5, 6]]
assert_(np.all(np.cumproduct(A) == np.array([1, 2, 6, 24, 120, 720])))
@@ -991,9 +998,110 @@ class TestNonzero(TestCase):
assert_(type(nzx_i) is np.ndarray)
assert_(nzx_i.flags.writeable)
- # Tests that the array method
- # call works
+ def test_count_nonzero_axis(self):
+ # Basic check of functionality
+ m = np.array([[0, 1, 7, 0, 0], [3, 0, 0, 2, 19]])
+
+ expected = np.array([1, 1, 1, 1, 1])
+ assert_equal(np.count_nonzero(m, axis=0), expected)
+
+ expected = np.array([2, 3])
+ assert_equal(np.count_nonzero(m, axis=1), expected)
+
+ 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(TypeError, np.count_nonzero,
+ m, axis=np.array([[1], [2]]))
+
+ def test_count_nonzero_axis_all_dtypes(self):
+ # More thorough test that the axis argument is respected
+ # for all dtypes and responds correctly when presented with
+ # either integer or tuple arguments for axis
+ msg = "Mismatch for dtype: %s"
+
+ for dt in np.typecodes['All']:
+ err_msg = msg % (np.dtype(dt).name,)
+
+ if dt != 'V':
+ if dt != 'M':
+ m = np.zeros((3, 3), dtype=dt)
+ n = np.ones(1, dtype=dt)
+
+ m[0, 0] = n[0]
+ m[1, 0] = n[0]
+
+ else: # np.zeros doesn't work for np.datetime64
+ m = np.array(['1970-01-01'] * 9)
+ m = m.reshape((3, 3))
+
+ m[0, 0] = '1970-01-12'
+ m[1, 0] = '1970-01-12'
+ m = m.astype(dt)
+
+ expected = np.array([2, 0, 0])
+ assert_equal(np.count_nonzero(m, axis=0),
+ expected, err_msg=err_msg)
+
+ expected = np.array([1, 1, 0])
+ assert_equal(np.count_nonzero(m, axis=1),
+ expected, err_msg=err_msg)
+
+ expected = np.array(2)
+ assert_equal(np.count_nonzero(m, axis=(0, 1)),
+ expected, err_msg=err_msg)
+ assert_equal(np.count_nonzero(m, axis=None),
+ expected, err_msg=err_msg)
+ assert_equal(np.count_nonzero(m),
+ expected, err_msg=err_msg)
+
+ if dt == 'V':
+ # There are no 'nonzero' objects for np.void, so the testing
+ # setup is slightly different for this dtype
+ m = np.array([np.void(1)] * 6).reshape((2, 3))
+
+ expected = np.array([0, 0, 0])
+ assert_equal(np.count_nonzero(m, axis=0),
+ expected, err_msg=err_msg)
+
+ expected = np.array([0, 0])
+ assert_equal(np.count_nonzero(m, axis=1),
+ expected, err_msg=err_msg)
+
+ expected = np.array(0)
+ assert_equal(np.count_nonzero(m, axis=(0, 1)),
+ expected, err_msg=err_msg)
+ assert_equal(np.count_nonzero(m, axis=None),
+ expected, err_msg=err_msg)
+ assert_equal(np.count_nonzero(m),
+ expected, err_msg=err_msg)
+
+ def test_count_nonzero_axis_consistent(self):
+ # Check that the axis behaviour for valid axes in
+ # non-special cases is consistent (and therefore
+ # correct) by checking it against an integer array
+ # that is then casted to the generic object dtype
+ from itertools import combinations, permutations
+
+ axis = (0, 1, 2, 3)
+ size = (5, 5, 5, 5)
+ msg = "Mismatch for axis: %s"
+
+ rng = np.random.RandomState(1234)
+ m = rng.randint(-100, 100, size=size)
+ n = m.astype(np.object)
+
+ for length in range(len(axis)):
+ for combo in combinations(axis, length):
+ for perm in permutations(combo):
+ assert_equal(
+ np.count_nonzero(m, axis=perm),
+ np.count_nonzero(n, axis=perm),
+ err_msg=msg % (perm,))
+
def test_array_method(self):
+ # Tests that the array method
+ # call to nonzero works
m = np.array([[1, 0, 0], [4, 0, 6]])
tgt = [[0, 1, 1], [0, 0, 2]]