summaryrefslogtreecommitdiff
path: root/numpy/lib/twodim_base.py
diff options
context:
space:
mode:
Diffstat (limited to 'numpy/lib/twodim_base.py')
-rw-r--r--numpy/lib/twodim_base.py175
1 files changed, 103 insertions, 72 deletions
diff --git a/numpy/lib/twodim_base.py b/numpy/lib/twodim_base.py
index 12c0f9bb3..a8925592a 100644
--- a/numpy/lib/twodim_base.py
+++ b/numpy/lib/twodim_base.py
@@ -11,8 +11,23 @@ __all__ = ['diag', 'diagflat', 'eye', 'fliplr', 'flipud', 'rot90', 'tri',
from numpy.core.numeric import (
asanyarray, subtract, arange, zeros, greater_equal, multiply, ones,
- asarray, where,
+ asarray, where, less, int8, int16, int32, int64, empty, promote_types
)
+from numpy.core import iinfo
+
+
+i1 = iinfo(int8)
+i2 = iinfo(int16)
+i4 = iinfo(int32)
+def _min_int(low, high):
+ """ get small int that fits the range """
+ if high <= i1.max and low >= i1.min:
+ return int8
+ if high <= i2.max and low >= i2.min:
+ return int16
+ if high <= i4.max and low >= i4.min:
+ return int32
+ return int64
def fliplr(m):
@@ -25,7 +40,7 @@ def fliplr(m):
Parameters
----------
m : array_like
- Input array.
+ Input array, must be at least 2-D.
Returns
-------
@@ -40,8 +55,7 @@ def fliplr(m):
Notes
-----
- Equivalent to A[:,::-1]. Does not require the array to be
- two-dimensional.
+ Equivalent to A[:,::-1]. Requires the array to be at least 2-D.
Examples
--------
@@ -373,6 +387,7 @@ 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)
@@ -394,8 +409,14 @@ def tri(N, M=None, k=0, dtype=float):
"""
if M is None:
M = N
- m = greater_equal(subtract.outer(arange(N), arange(M)), -k)
- return m.astype(dtype)
+
+ m = greater_equal.outer(arange(N, dtype=_min_int(0, N)),
+ arange(-k, M-k, dtype=_min_int(-k, M - k)))
+
+ # Avoid making a copy if the requested type is already bool
+ m = m.astype(dtype, copy=False)
+
+ return m
def tril(m, k=0):
@@ -431,8 +452,7 @@ def tril(m, k=0):
"""
m = asanyarray(m)
- out = multiply(tri(m.shape[0], m.shape[1], k=k, dtype=m.dtype), m)
- return out
+ return multiply(tri(*m.shape[-2:], k=k, dtype=bool), m, dtype=m.dtype)
def triu(m, k=0):
@@ -458,22 +478,20 @@ def triu(m, k=0):
"""
m = asanyarray(m)
- out = multiply((1 - tri(m.shape[0], m.shape[1], k - 1, dtype=m.dtype)), m)
- return out
+ return multiply(~tri(*m.shape[-2:], k=k-1, dtype=bool), m, dtype=m.dtype)
# Originally borrowed from John Hunter and matplotlib
-def vander(x, N=None, order='decreasing'):
+def vander(x, N=None, increasing=False):
"""
Generate a Vandermonde matrix.
- The columns of the output matrix are powers of the input vector. The
- order of the powers is determined by the `order` argument, either
- "decreasing" (the default) or "increasing". Specifically, when
- `order` is "decreasing", the `i`-th output column is the input vector
- raised element-wise to the power of ``N - i - 1``. Such a matrix with
- a geometric progression in each row is named for Alexandre-Theophile
- Vandermonde.
+ The columns of the output matrix are powers of the input vector. The
+ order of the powers is determined by the `increasing` boolean argument.
+ Specifically, when `increasing` is False, the `i`-th output column is
+ the input vector raised element-wise to the power of ``N - i - 1``. Such
+ a matrix with a geometric progression in each row is named for Alexandre-
+ Theophile Vandermonde.
Parameters
----------
@@ -482,16 +500,18 @@ def vander(x, N=None, order='decreasing'):
N : int, optional
Number of columns in the output. If `N` is not specified, a square
array is returned (``N = len(x)``).
- order : str, optional
- Order of the powers of the columns. Must be either 'decreasing'
- (the default) or 'increasing'.
+ increasing : bool, optional
+ Order of the powers of the columns. If True, the powers increase
+ from left to right, if False (the default) they are reversed.
+
+ .. versionadded:: 1.9.0
Returns
-------
out : ndarray
- Vandermonde matrix. If `order` is "decreasing", the first column is
- ``x^(N-1)``, the second ``x^(N-2)`` and so forth. If `order` is
- "increasing", the columns are ``x^0, x^1, ..., x^(N-1)``.
+ Vandermonde matrix. If `increasing` is False, the first column is
+ ``x^(N-1)``, the second ``x^(N-2)`` and so forth. If `increasing` is
+ True, the columns are ``x^0, x^1, ..., x^(N-1)``.
See Also
--------
@@ -519,7 +539,7 @@ def vander(x, N=None, order='decreasing'):
[ 8, 4, 2, 1],
[ 27, 9, 3, 1],
[125, 25, 5, 1]])
- >>> np.vander(x, order='increasing')
+ >>> np.vander(x, increasing=True)
array([[ 1, 1, 1, 1],
[ 1, 2, 4, 8],
[ 1, 3, 9, 27],
@@ -534,22 +554,22 @@ def vander(x, N=None, order='decreasing'):
48
"""
- if order not in ['decreasing', 'increasing']:
- raise ValueError("Invalid order %r; order must be either "
- "'decreasing' or 'increasing'." % (order,))
x = asarray(x)
if x.ndim != 1:
raise ValueError("x must be a one-dimensional array or sequence.")
if N is None:
N = len(x)
- if order == "decreasing":
- powers = arange(N - 1, -1, -1)
- else:
- powers = arange(N)
- V = x.reshape(-1, 1) ** powers
+ v = empty((len(x), N), dtype=promote_types(x.dtype, int))
+ tmp = v[:, ::-1] if not increasing else v
+
+ if N > 0:
+ tmp[:, 0] = 1
+ if N > 1:
+ tmp[:, 1:] = x[:, None]
+ multiply.accumulate(tmp[:, 1:], out=tmp[:, 1:], axis=1)
- return V
+ return v
def histogram2d(x, y, bins=10, range=None, normed=False, weights=None):
@@ -559,9 +579,11 @@ def histogram2d(x, y, bins=10, range=None, normed=False, weights=None):
Parameters
----------
x : array_like, shape (N,)
- An array containing the x coordinates of the points to be histogrammed.
+ An array containing the x coordinates of the points to be
+ histogrammed.
y : array_like, shape (N,)
- An array containing the y coordinates of the points to be histogrammed.
+ An array containing the y coordinates of the points to be
+ histogrammed.
bins : int or [int, int] or array_like or [array, array], optional
The bin specification:
@@ -579,13 +601,13 @@ def histogram2d(x, y, bins=10, range=None, normed=False, weights=None):
``[[xmin, xmax], [ymin, ymax]]``. All values outside of this range
will be considered outliers and not tallied in the histogram.
normed : bool, optional
- If False, returns the number of samples in each bin. If True, returns
- the bin density, i.e. the bin count divided by the bin area.
+ If False, returns the number of samples in each bin. If True,
+ returns the bin density ``bin_count / sample_count / bin_area``.
weights : array_like, shape(N,), optional
- An array of values ``w_i`` weighing each sample ``(x_i, y_i)``. Weights
- are normalized to 1 if `normed` is True. If `normed` is False, the
- values of the returned histogram are equal to the sum of the weights
- belonging to the samples falling into each bin.
+ An array of values ``w_i`` weighing each sample ``(x_i, y_i)``.
+ Weights are normalized to 1 if `normed` is True. If `normed` is
+ False, the values of the returned histogram are equal to the sum of
+ the weights belonging to the samples falling into each bin.
Returns
-------
@@ -605,20 +627,15 @@ def histogram2d(x, y, bins=10, range=None, normed=False, weights=None):
Notes
-----
- When `normed` is True, then the returned histogram is the sample density,
- defined such that:
-
- .. math::
- \\sum_{i=0}^{nx-1} \\sum_{j=0}^{ny-1} H_{i,j} \\Delta x_i \\Delta y_j = 1
-
- where `H` is the histogram array and :math:`\\Delta x_i \\Delta y_i`
- the area of bin ``{i,j}``.
+ When `normed` is True, then the returned histogram is the sample
+ density, defined such that the sum over bins of the product
+ ``bin_value * bin_area`` is 1.
Please note that the histogram does not follow the Cartesian convention
- where `x` values are on the abcissa and `y` values on the ordinate axis.
- Rather, `x` is histogrammed along the first dimension of the array
- (vertical), and `y` along the second dimension of the array (horizontal).
- This ensures compatibility with `histogramdd`.
+ where `x` values are on the abscissa and `y` values on the ordinate
+ axis. Rather, `x` is histogrammed along the first dimension of the
+ array (vertical), and `y` along the second dimension of the array
+ (horizontal). This ensures compatibility with `histogramdd`.
Examples
--------
@@ -650,14 +667,14 @@ def histogram2d(x, y, bins=10, range=None, normed=False, weights=None):
>>> fig = plt.figure(figsize=(7, 3))
>>> ax = fig.add_subplot(131)
- >>> ax.set_title('imshow:\nequidistant')
+ >>> ax.set_title('imshow: equidistant')
>>> im = plt.imshow(H, interpolation='nearest', origin='low',
- extent=[xedges[0], xedges[-1], yedges[0], yedges[-1]])
+ extent=[xedges[0], xedges[-1], yedges[0], yedges[-1]])
pcolormesh can display exact bin edges:
>>> ax = fig.add_subplot(132)
- >>> ax.set_title('pcolormesh:\nexact bin edges')
+ >>> ax.set_title('pcolormesh: exact bin edges')
>>> X, Y = np.meshgrid(xedges, yedges)
>>> ax.pcolormesh(X, Y, H)
>>> ax.set_aspect('equal')
@@ -665,7 +682,7 @@ def histogram2d(x, y, bins=10, range=None, normed=False, weights=None):
NonUniformImage displays exact bin edges with interpolation:
>>> ax = fig.add_subplot(133)
- >>> ax.set_title('NonUniformImage:\ninterpolated')
+ >>> ax.set_title('NonUniformImage: interpolated')
>>> im = mpl.image.NonUniformImage(ax, interpolation='bilinear')
>>> xcenters = xedges[:-1] + 0.5 * (xedges[1:] - xedges[:-1])
>>> ycenters = yedges[:-1] + 0.5 * (yedges[1:] - yedges[:-1])
@@ -761,17 +778,24 @@ def mask_indices(n, mask_func, k=0):
return where(a != 0)
-def tril_indices(n, k=0):
+def tril_indices(n, k=0, m=None):
"""
- Return the indices for the lower-triangle of an (n, n) array.
+ Return the indices for the lower-triangle of an (n, m) array.
Parameters
----------
n : int
- The row dimension of the square arrays for which the returned
+ The row dimension of the arrays for which the returned
indices will be valid.
k : int, optional
Diagonal offset (see `tril` for details).
+ m : int, optional
+ .. versionadded:: 1.9.0
+
+ The column dimension of the arrays for which the returned
+ arrays will be valid.
+ By default `m` is taken equal to `n`.
+
Returns
-------
@@ -831,7 +855,7 @@ def tril_indices(n, k=0):
[-10, -10, -10, -10]])
"""
- return mask_indices(n, tril, k)
+ return where(tri(n, m, k=k, dtype=bool))
def tril_indices_from(arr, k=0):
@@ -857,14 +881,14 @@ def tril_indices_from(arr, k=0):
.. versionadded:: 1.4.0
"""
- if not (arr.ndim == 2 and arr.shape[0] == arr.shape[1]):
- raise ValueError("input array must be 2-d and square")
- return tril_indices(arr.shape[0], k)
+ if arr.ndim != 2:
+ raise ValueError("input array must be 2-d")
+ return tril_indices(arr.shape[-2], k=k, m=arr.shape[-1])
-def triu_indices(n, k=0):
+def triu_indices(n, k=0, m=None):
"""
- Return the indices for the upper-triangle of an (n, n) array.
+ Return the indices for the upper-triangle of an (n, m) array.
Parameters
----------
@@ -873,6 +897,13 @@ def triu_indices(n, k=0):
be valid.
k : int, optional
Diagonal offset (see `triu` for details).
+ m : int, optional
+ .. versionadded:: 1.9.0
+
+ The column dimension of the arrays for which the returned
+ arrays will be valid.
+ By default `m` is taken equal to `n`.
+
Returns
-------
@@ -934,12 +965,12 @@ def triu_indices(n, k=0):
[ 12, 13, 14, -1]])
"""
- return mask_indices(n, triu, k)
+ return where(~tri(n, m, k=k-1, dtype=bool))
def triu_indices_from(arr, k=0):
"""
- Return the indices for the upper-triangle of a (N, N) array.
+ Return the indices for the upper-triangle of arr.
See `triu_indices` for full details.
@@ -964,6 +995,6 @@ def triu_indices_from(arr, k=0):
.. versionadded:: 1.4.0
"""
- if not (arr.ndim == 2 and arr.shape[0] == arr.shape[1]):
- raise ValueError("input array must be 2-d and square")
- return triu_indices(arr.shape[0], k)
+ if arr.ndim != 2:
+ raise ValueError("input array must be 2-d")
+ return triu_indices(arr.shape[-2], k=k, m=arr.shape[-1])