summaryrefslogtreecommitdiff
path: root/numpy/lib/histograms.py
diff options
context:
space:
mode:
authorEric Wieser <wieser.eric@gmail.com>2018-03-26 00:36:56 -0700
committerEric Wieser <wieser.eric@gmail.com>2018-04-06 16:25:27 -0700
commit992163e51d48d26df05430690052e42dbdd7fe9b (patch)
treea2956938f13c2b1861d25ea08cee5fab27a5e29e /numpy/lib/histograms.py
parentfd92d02152d50735b5ffafcce3a7a722e366b192 (diff)
downloadnumpy-992163e51d48d26df05430690052e42dbdd7fe9b.tar.gz
MAINT/ENH: Reuse range computation code from np.histogram in np.histogramdd
This also adds support for inferring the range along a subset of the axes, rather than an all or nothing approach.
Diffstat (limited to 'numpy/lib/histograms.py')
-rw-r--r--numpy/lib/histograms.py44
1 files changed, 15 insertions, 29 deletions
diff --git a/numpy/lib/histograms.py b/numpy/lib/histograms.py
index a8cd9acb8..c0cfa07c7 100644
--- a/numpy/lib/histograms.py
+++ b/numpy/lib/histograms.py
@@ -791,9 +791,12 @@ def histogramdd(sample, bins=10, range=None, normed=False, weights=None):
* The number of bins for all dimensions (nx=ny=...=bins).
range : sequence, optional
- A sequence of lower and upper bin edges to be used if the edges are
- not given explicitly in `bins`. Defaults to the minimum and maximum
- values along each dimension.
+ A sequence of length D, each an optional (lower, upper) tuple giving
+ the outer bin edges to be used if the edges are not given explicitly in
+ `bins`.
+ An entry of None in the sequence results in the minimum and maximum
+ values being used for the corresponding dimension.
+ The default, None, is equivalent to passing a tuple of D None values.
normed : bool, optional
If False, returns the number of samples in each bin. If True,
returns the bin density ``bin_count / sample_count / bin_volume``.
@@ -849,43 +852,26 @@ def histogramdd(sample, bins=10, range=None, normed=False, weights=None):
# bins is an integer
bins = D*[bins]
- # Select range for each dimension
- # Used only if number of bins is given.
- if range is None:
- # Handle empty input. Range can't be determined in that case, use 0-1.
- if N == 0:
- smin = np.zeros(D)
- smax = np.ones(D)
- else:
- smin = np.atleast_1d(np.array(sample.min(0), float))
- smax = np.atleast_1d(np.array(sample.max(0), float))
- else:
- if not np.all(np.isfinite(range)):
- raise ValueError(
- 'range parameter must be finite.')
- smin = np.zeros(D)
- smax = np.zeros(D)
- for i in np.arange(D):
- smin[i], smax[i] = range[i]
-
- # Make sure the bins have a finite width.
- for i in np.arange(len(smin)):
- if smin[i] == smax[i]:
- 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
+
+ # normalize the range argument
+ if range is None:
+ range = (None,) * D
+ elif len(range) != D:
+ raise ValueError('range argument must have one entry per dimension')
+
# Create edge arrays
for i in np.arange(D):
if np.ndim(bins[i]) == 0:
if bins[i] < 1:
raise ValueError(
'`bins[{}]` must be positive, when an integer'.format(i))
- edges[i] = np.linspace(smin[i], smax[i], bins[i] + 1, dtype=edge_dt)
+ smin, smax = _get_outer_edges(sample[:,i], range[i])
+ edges[i] = np.linspace(smin, smax, bins[i] + 1, dtype=edge_dt)
elif np.ndim(bins[i]) == 1:
edges[i] = np.asarray(bins[i], edge_dt)
# not just monotonic, due to the use of mindiff below