767 lines
31 KiB
Python
767 lines
31 KiB
Python
|
__all__ = ['RegularGridInterpolator', 'interpn']
|
||
|
|
||
|
import itertools
|
||
|
import warnings
|
||
|
|
||
|
import numpy as np
|
||
|
|
||
|
import scipy.sparse.linalg as ssl
|
||
|
|
||
|
from .interpnd import _ndim_coords_from_arrays
|
||
|
from ._cubic import PchipInterpolator
|
||
|
from ._rgi_cython import evaluate_linear_2d, find_indices
|
||
|
from ._bsplines import make_interp_spline
|
||
|
from ._fitpack2 import RectBivariateSpline
|
||
|
from ._ndbspline import make_ndbspl
|
||
|
|
||
|
|
||
|
def _check_points(points):
|
||
|
descending_dimensions = []
|
||
|
grid = []
|
||
|
for i, p in enumerate(points):
|
||
|
# early make points float
|
||
|
# see https://github.com/scipy/scipy/pull/17230
|
||
|
p = np.asarray(p, dtype=float)
|
||
|
if not np.all(p[1:] > p[:-1]):
|
||
|
if np.all(p[1:] < p[:-1]):
|
||
|
# input is descending, so make it ascending
|
||
|
descending_dimensions.append(i)
|
||
|
p = np.flip(p)
|
||
|
else:
|
||
|
raise ValueError(
|
||
|
"The points in dimension %d must be strictly "
|
||
|
"ascending or descending" % i)
|
||
|
# see https://github.com/scipy/scipy/issues/17716
|
||
|
p = np.ascontiguousarray(p)
|
||
|
grid.append(p)
|
||
|
return tuple(grid), tuple(descending_dimensions)
|
||
|
|
||
|
|
||
|
def _check_dimensionality(points, values):
|
||
|
if len(points) > values.ndim:
|
||
|
raise ValueError("There are %d point arrays, but values has %d "
|
||
|
"dimensions" % (len(points), values.ndim))
|
||
|
for i, p in enumerate(points):
|
||
|
if not np.asarray(p).ndim == 1:
|
||
|
raise ValueError("The points in dimension %d must be "
|
||
|
"1-dimensional" % i)
|
||
|
if not values.shape[i] == len(p):
|
||
|
raise ValueError("There are %d points and %d values in "
|
||
|
"dimension %d" % (len(p), values.shape[i], i))
|
||
|
|
||
|
|
||
|
class RegularGridInterpolator:
|
||
|
"""
|
||
|
Interpolator on a regular or rectilinear grid in arbitrary dimensions.
|
||
|
|
||
|
The data must be defined on a rectilinear grid; that is, a rectangular
|
||
|
grid with even or uneven spacing. Linear, nearest-neighbor, spline
|
||
|
interpolations are supported. After setting up the interpolator object,
|
||
|
the interpolation method may be chosen at each evaluation.
|
||
|
|
||
|
Parameters
|
||
|
----------
|
||
|
points : tuple of ndarray of float, with shapes (m1, ), ..., (mn, )
|
||
|
The points defining the regular grid in n dimensions. The points in
|
||
|
each dimension (i.e. every elements of the points tuple) must be
|
||
|
strictly ascending or descending.
|
||
|
|
||
|
values : array_like, shape (m1, ..., mn, ...)
|
||
|
The data on the regular grid in n dimensions. Complex data is
|
||
|
accepted.
|
||
|
|
||
|
.. deprecated:: 1.13.0
|
||
|
Complex data is deprecated with ``method="pchip"`` and will raise an
|
||
|
error in SciPy 1.15.0. This is because ``PchipInterpolator`` only
|
||
|
works with real values. If you are trying to use the real components of
|
||
|
the passed array, use ``np.real`` on ``values``.
|
||
|
|
||
|
method : str, optional
|
||
|
The method of interpolation to perform. Supported are "linear",
|
||
|
"nearest", "slinear", "cubic", "quintic" and "pchip". This
|
||
|
parameter will become the default for the object's ``__call__``
|
||
|
method. Default is "linear".
|
||
|
|
||
|
bounds_error : bool, optional
|
||
|
If True, when interpolated values are requested outside of the
|
||
|
domain of the input data, a ValueError is raised.
|
||
|
If False, then `fill_value` is used.
|
||
|
Default is True.
|
||
|
|
||
|
fill_value : float or None, optional
|
||
|
The value to use for points outside of the interpolation domain.
|
||
|
If None, values outside the domain are extrapolated.
|
||
|
Default is ``np.nan``.
|
||
|
|
||
|
solver : callable, optional
|
||
|
Only used for methods "slinear", "cubic" and "quintic".
|
||
|
Sparse linear algebra solver for construction of the NdBSpline instance.
|
||
|
Default is the iterative solver `scipy.sparse.linalg.gcrotmk`.
|
||
|
|
||
|
.. versionadded:: 1.13
|
||
|
|
||
|
solver_args: dict, optional
|
||
|
Additional arguments to pass to `solver`, if any.
|
||
|
|
||
|
.. versionadded:: 1.13
|
||
|
|
||
|
Methods
|
||
|
-------
|
||
|
__call__
|
||
|
|
||
|
Attributes
|
||
|
----------
|
||
|
grid : tuple of ndarrays
|
||
|
The points defining the regular grid in n dimensions.
|
||
|
This tuple defines the full grid via
|
||
|
``np.meshgrid(*grid, indexing='ij')``
|
||
|
values : ndarray
|
||
|
Data values at the grid.
|
||
|
method : str
|
||
|
Interpolation method.
|
||
|
fill_value : float or ``None``
|
||
|
Use this value for out-of-bounds arguments to `__call__`.
|
||
|
bounds_error : bool
|
||
|
If ``True``, out-of-bounds argument raise a ``ValueError``.
|
||
|
|
||
|
Notes
|
||
|
-----
|
||
|
Contrary to `LinearNDInterpolator` and `NearestNDInterpolator`, this class
|
||
|
avoids expensive triangulation of the input data by taking advantage of the
|
||
|
regular grid structure.
|
||
|
|
||
|
In other words, this class assumes that the data is defined on a
|
||
|
*rectilinear* grid.
|
||
|
|
||
|
.. versionadded:: 0.14
|
||
|
|
||
|
The 'slinear'(k=1), 'cubic'(k=3), and 'quintic'(k=5) methods are
|
||
|
tensor-product spline interpolators, where `k` is the spline degree,
|
||
|
If any dimension has fewer points than `k` + 1, an error will be raised.
|
||
|
|
||
|
.. versionadded:: 1.9
|
||
|
|
||
|
If the input data is such that dimensions have incommensurate
|
||
|
units and differ by many orders of magnitude, the interpolant may have
|
||
|
numerical artifacts. Consider rescaling the data before interpolating.
|
||
|
|
||
|
**Choosing a solver for spline methods**
|
||
|
|
||
|
Spline methods, "slinear", "cubic" and "quintic" involve solving a
|
||
|
large sparse linear system at instantiation time. Depending on data,
|
||
|
the default solver may or may not be adequate. When it is not, you may
|
||
|
need to experiment with an optional `solver` argument, where you may
|
||
|
choose between the direct solver (`scipy.sparse.linalg.spsolve`) or
|
||
|
iterative solvers from `scipy.sparse.linalg`. You may need to supply
|
||
|
additional parameters via the optional `solver_args` parameter (for instance,
|
||
|
you may supply the starting value or target tolerance). See the
|
||
|
`scipy.sparse.linalg` documentation for the full list of available options.
|
||
|
|
||
|
Alternatively, you may instead use the legacy methods, "slinear_legacy",
|
||
|
"cubic_legacy" and "quintic_legacy". These methods allow faster construction
|
||
|
but evaluations will be much slower.
|
||
|
|
||
|
Examples
|
||
|
--------
|
||
|
**Evaluate a function on the points of a 3-D grid**
|
||
|
|
||
|
As a first example, we evaluate a simple example function on the points of
|
||
|
a 3-D grid:
|
||
|
|
||
|
>>> from scipy.interpolate import RegularGridInterpolator
|
||
|
>>> import numpy as np
|
||
|
>>> def f(x, y, z):
|
||
|
... return 2 * x**3 + 3 * y**2 - z
|
||
|
>>> x = np.linspace(1, 4, 11)
|
||
|
>>> y = np.linspace(4, 7, 22)
|
||
|
>>> z = np.linspace(7, 9, 33)
|
||
|
>>> xg, yg ,zg = np.meshgrid(x, y, z, indexing='ij', sparse=True)
|
||
|
>>> data = f(xg, yg, zg)
|
||
|
|
||
|
``data`` is now a 3-D array with ``data[i, j, k] = f(x[i], y[j], z[k])``.
|
||
|
Next, define an interpolating function from this data:
|
||
|
|
||
|
>>> interp = RegularGridInterpolator((x, y, z), data)
|
||
|
|
||
|
Evaluate the interpolating function at the two points
|
||
|
``(x,y,z) = (2.1, 6.2, 8.3)`` and ``(3.3, 5.2, 7.1)``:
|
||
|
|
||
|
>>> pts = np.array([[2.1, 6.2, 8.3],
|
||
|
... [3.3, 5.2, 7.1]])
|
||
|
>>> interp(pts)
|
||
|
array([ 125.80469388, 146.30069388])
|
||
|
|
||
|
which is indeed a close approximation to
|
||
|
|
||
|
>>> f(2.1, 6.2, 8.3), f(3.3, 5.2, 7.1)
|
||
|
(125.54200000000002, 145.894)
|
||
|
|
||
|
**Interpolate and extrapolate a 2D dataset**
|
||
|
|
||
|
As a second example, we interpolate and extrapolate a 2D data set:
|
||
|
|
||
|
>>> x, y = np.array([-2, 0, 4]), np.array([-2, 0, 2, 5])
|
||
|
>>> def ff(x, y):
|
||
|
... return x**2 + y**2
|
||
|
|
||
|
>>> xg, yg = np.meshgrid(x, y, indexing='ij')
|
||
|
>>> data = ff(xg, yg)
|
||
|
>>> interp = RegularGridInterpolator((x, y), data,
|
||
|
... bounds_error=False, fill_value=None)
|
||
|
|
||
|
>>> import matplotlib.pyplot as plt
|
||
|
>>> fig = plt.figure()
|
||
|
>>> ax = fig.add_subplot(projection='3d')
|
||
|
>>> ax.scatter(xg.ravel(), yg.ravel(), data.ravel(),
|
||
|
... s=60, c='k', label='data')
|
||
|
|
||
|
Evaluate and plot the interpolator on a finer grid
|
||
|
|
||
|
>>> xx = np.linspace(-4, 9, 31)
|
||
|
>>> yy = np.linspace(-4, 9, 31)
|
||
|
>>> X, Y = np.meshgrid(xx, yy, indexing='ij')
|
||
|
|
||
|
>>> # interpolator
|
||
|
>>> ax.plot_wireframe(X, Y, interp((X, Y)), rstride=3, cstride=3,
|
||
|
... alpha=0.4, color='m', label='linear interp')
|
||
|
|
||
|
>>> # ground truth
|
||
|
>>> ax.plot_wireframe(X, Y, ff(X, Y), rstride=3, cstride=3,
|
||
|
... alpha=0.4, label='ground truth')
|
||
|
>>> plt.legend()
|
||
|
>>> plt.show()
|
||
|
|
||
|
Other examples are given
|
||
|
:ref:`in the tutorial <tutorial-interpolate_regular_grid_interpolator>`.
|
||
|
|
||
|
See Also
|
||
|
--------
|
||
|
NearestNDInterpolator : Nearest neighbor interpolator on *unstructured*
|
||
|
data in N dimensions
|
||
|
|
||
|
LinearNDInterpolator : Piecewise linear interpolator on *unstructured* data
|
||
|
in N dimensions
|
||
|
|
||
|
interpn : a convenience function which wraps `RegularGridInterpolator`
|
||
|
|
||
|
scipy.ndimage.map_coordinates : interpolation on grids with equal spacing
|
||
|
(suitable for e.g., N-D image resampling)
|
||
|
|
||
|
References
|
||
|
----------
|
||
|
.. [1] Python package *regulargrid* by Johannes Buchner, see
|
||
|
https://pypi.python.org/pypi/regulargrid/
|
||
|
.. [2] Wikipedia, "Trilinear interpolation",
|
||
|
https://en.wikipedia.org/wiki/Trilinear_interpolation
|
||
|
.. [3] Weiser, Alan, and Sergio E. Zarantonello. "A note on piecewise linear
|
||
|
and multilinear table interpolation in many dimensions." MATH.
|
||
|
COMPUT. 50.181 (1988): 189-196.
|
||
|
https://www.ams.org/journals/mcom/1988-50-181/S0025-5718-1988-0917826-0/S0025-5718-1988-0917826-0.pdf
|
||
|
:doi:`10.1090/S0025-5718-1988-0917826-0`
|
||
|
|
||
|
"""
|
||
|
# this class is based on code originally programmed by Johannes Buchner,
|
||
|
# see https://github.com/JohannesBuchner/regulargrid
|
||
|
|
||
|
_SPLINE_DEGREE_MAP = {"slinear": 1, "cubic": 3, "quintic": 5, 'pchip': 3,
|
||
|
"slinear_legacy": 1, "cubic_legacy": 3, "quintic_legacy": 5,}
|
||
|
_SPLINE_METHODS_recursive = {"slinear_legacy", "cubic_legacy",
|
||
|
"quintic_legacy", "pchip"}
|
||
|
_SPLINE_METHODS_ndbspl = {"slinear", "cubic", "quintic"}
|
||
|
_SPLINE_METHODS = list(_SPLINE_DEGREE_MAP.keys())
|
||
|
_ALL_METHODS = ["linear", "nearest"] + _SPLINE_METHODS
|
||
|
|
||
|
def __init__(self, points, values, method="linear", bounds_error=True,
|
||
|
fill_value=np.nan, *, solver=None, solver_args=None):
|
||
|
if method not in self._ALL_METHODS:
|
||
|
raise ValueError("Method '%s' is not defined" % method)
|
||
|
elif method in self._SPLINE_METHODS:
|
||
|
self._validate_grid_dimensions(points, method)
|
||
|
self.method = method
|
||
|
self.bounds_error = bounds_error
|
||
|
self.grid, self._descending_dimensions = _check_points(points)
|
||
|
self.values = self._check_values(values)
|
||
|
self._check_dimensionality(self.grid, self.values)
|
||
|
self.fill_value = self._check_fill_value(self.values, fill_value)
|
||
|
if self._descending_dimensions:
|
||
|
self.values = np.flip(values, axis=self._descending_dimensions)
|
||
|
if self.method == "pchip" and np.iscomplexobj(self.values):
|
||
|
msg = ("`PchipInterpolator` only works with real values. Passing "
|
||
|
"complex-dtyped `values` with `method='pchip'` is deprecated "
|
||
|
"and will raise an error in SciPy 1.15.0. If you are trying to "
|
||
|
"use the real components of the passed array, use `np.real` on "
|
||
|
"the array before passing to `RegularGridInterpolator`.")
|
||
|
warnings.warn(msg, DeprecationWarning, stacklevel=2)
|
||
|
if method in self._SPLINE_METHODS_ndbspl:
|
||
|
if solver_args is None:
|
||
|
solver_args = {}
|
||
|
self._spline = self._construct_spline(method, solver, **solver_args)
|
||
|
else:
|
||
|
if solver is not None or solver_args:
|
||
|
raise ValueError(
|
||
|
f"{method =} does not accept the 'solver' argument. Got "
|
||
|
f" {solver = } and with arguments {solver_args}."
|
||
|
)
|
||
|
|
||
|
def _construct_spline(self, method, solver=None, **solver_args):
|
||
|
if solver is None:
|
||
|
solver = ssl.gcrotmk
|
||
|
spl = make_ndbspl(
|
||
|
self.grid, self.values, self._SPLINE_DEGREE_MAP[method],
|
||
|
solver=solver, **solver_args
|
||
|
)
|
||
|
return spl
|
||
|
|
||
|
def _check_dimensionality(self, grid, values):
|
||
|
_check_dimensionality(grid, values)
|
||
|
|
||
|
def _check_points(self, points):
|
||
|
return _check_points(points)
|
||
|
|
||
|
def _check_values(self, values):
|
||
|
if not hasattr(values, 'ndim'):
|
||
|
# allow reasonable duck-typed values
|
||
|
values = np.asarray(values)
|
||
|
|
||
|
if hasattr(values, 'dtype') and hasattr(values, 'astype'):
|
||
|
if not np.issubdtype(values.dtype, np.inexact):
|
||
|
values = values.astype(float)
|
||
|
|
||
|
return values
|
||
|
|
||
|
def _check_fill_value(self, values, fill_value):
|
||
|
if fill_value is not None:
|
||
|
fill_value_dtype = np.asarray(fill_value).dtype
|
||
|
if (hasattr(values, 'dtype') and not
|
||
|
np.can_cast(fill_value_dtype, values.dtype,
|
||
|
casting='same_kind')):
|
||
|
raise ValueError("fill_value must be either 'None' or "
|
||
|
"of a type compatible with values")
|
||
|
return fill_value
|
||
|
|
||
|
def __call__(self, xi, method=None, *, nu=None):
|
||
|
"""
|
||
|
Interpolation at coordinates.
|
||
|
|
||
|
Parameters
|
||
|
----------
|
||
|
xi : ndarray of shape (..., ndim)
|
||
|
The coordinates to evaluate the interpolator at.
|
||
|
|
||
|
method : str, optional
|
||
|
The method of interpolation to perform. Supported are "linear",
|
||
|
"nearest", "slinear", "cubic", "quintic" and "pchip". Default is
|
||
|
the method chosen when the interpolator was created.
|
||
|
|
||
|
nu : sequence of ints, length ndim, optional
|
||
|
If not None, the orders of the derivatives to evaluate.
|
||
|
Each entry must be non-negative.
|
||
|
Only allowed for methods "slinear", "cubic" and "quintic".
|
||
|
|
||
|
.. versionadded:: 1.13
|
||
|
|
||
|
Returns
|
||
|
-------
|
||
|
values_x : ndarray, shape xi.shape[:-1] + values.shape[ndim:]
|
||
|
Interpolated values at `xi`. See notes for behaviour when
|
||
|
``xi.ndim == 1``.
|
||
|
|
||
|
Notes
|
||
|
-----
|
||
|
In the case that ``xi.ndim == 1`` a new axis is inserted into
|
||
|
the 0 position of the returned array, values_x, so its shape is
|
||
|
instead ``(1,) + values.shape[ndim:]``.
|
||
|
|
||
|
Examples
|
||
|
--------
|
||
|
Here we define a nearest-neighbor interpolator of a simple function
|
||
|
|
||
|
>>> import numpy as np
|
||
|
>>> x, y = np.array([0, 1, 2]), np.array([1, 3, 7])
|
||
|
>>> def f(x, y):
|
||
|
... return x**2 + y**2
|
||
|
>>> data = f(*np.meshgrid(x, y, indexing='ij', sparse=True))
|
||
|
>>> from scipy.interpolate import RegularGridInterpolator
|
||
|
>>> interp = RegularGridInterpolator((x, y), data, method='nearest')
|
||
|
|
||
|
By construction, the interpolator uses the nearest-neighbor
|
||
|
interpolation
|
||
|
|
||
|
>>> interp([[1.5, 1.3], [0.3, 4.5]])
|
||
|
array([2., 9.])
|
||
|
|
||
|
We can however evaluate the linear interpolant by overriding the
|
||
|
`method` parameter
|
||
|
|
||
|
>>> interp([[1.5, 1.3], [0.3, 4.5]], method='linear')
|
||
|
array([ 4.7, 24.3])
|
||
|
"""
|
||
|
method = self.method if method is None else method
|
||
|
is_method_changed = self.method != method
|
||
|
if method not in self._ALL_METHODS:
|
||
|
raise ValueError("Method '%s' is not defined" % method)
|
||
|
if is_method_changed and method in self._SPLINE_METHODS_ndbspl:
|
||
|
self._spline = self._construct_spline(method)
|
||
|
|
||
|
if nu is not None and method not in self._SPLINE_METHODS_ndbspl:
|
||
|
raise ValueError(
|
||
|
f"Can only compute derivatives for methods "
|
||
|
f"{self._SPLINE_METHODS_ndbspl}, got {method =}."
|
||
|
)
|
||
|
|
||
|
xi, xi_shape, ndim, nans, out_of_bounds = self._prepare_xi(xi)
|
||
|
|
||
|
if method == "linear":
|
||
|
indices, norm_distances = self._find_indices(xi.T)
|
||
|
if (ndim == 2 and hasattr(self.values, 'dtype') and
|
||
|
self.values.ndim == 2 and self.values.flags.writeable and
|
||
|
self.values.dtype in (np.float64, np.complex128) and
|
||
|
self.values.dtype.byteorder == '='):
|
||
|
# until cython supports const fused types, the fast path
|
||
|
# cannot support non-writeable values
|
||
|
# a fast path
|
||
|
out = np.empty(indices.shape[1], dtype=self.values.dtype)
|
||
|
result = evaluate_linear_2d(self.values,
|
||
|
indices,
|
||
|
norm_distances,
|
||
|
self.grid,
|
||
|
out)
|
||
|
else:
|
||
|
result = self._evaluate_linear(indices, norm_distances)
|
||
|
elif method == "nearest":
|
||
|
indices, norm_distances = self._find_indices(xi.T)
|
||
|
result = self._evaluate_nearest(indices, norm_distances)
|
||
|
elif method in self._SPLINE_METHODS:
|
||
|
if is_method_changed:
|
||
|
self._validate_grid_dimensions(self.grid, method)
|
||
|
if method in self._SPLINE_METHODS_recursive:
|
||
|
result = self._evaluate_spline(xi, method)
|
||
|
else:
|
||
|
result = self._spline(xi, nu=nu)
|
||
|
|
||
|
if not self.bounds_error and self.fill_value is not None:
|
||
|
result[out_of_bounds] = self.fill_value
|
||
|
|
||
|
# f(nan) = nan, if any
|
||
|
if np.any(nans):
|
||
|
result[nans] = np.nan
|
||
|
return result.reshape(xi_shape[:-1] + self.values.shape[ndim:])
|
||
|
|
||
|
def _prepare_xi(self, xi):
|
||
|
ndim = len(self.grid)
|
||
|
xi = _ndim_coords_from_arrays(xi, ndim=ndim)
|
||
|
if xi.shape[-1] != len(self.grid):
|
||
|
raise ValueError("The requested sample points xi have dimension "
|
||
|
f"{xi.shape[-1]} but this "
|
||
|
f"RegularGridInterpolator has dimension {ndim}")
|
||
|
|
||
|
xi_shape = xi.shape
|
||
|
xi = xi.reshape(-1, xi_shape[-1])
|
||
|
xi = np.asarray(xi, dtype=float)
|
||
|
|
||
|
# find nans in input
|
||
|
nans = np.any(np.isnan(xi), axis=-1)
|
||
|
|
||
|
if self.bounds_error:
|
||
|
for i, p in enumerate(xi.T):
|
||
|
if not np.logical_and(np.all(self.grid[i][0] <= p),
|
||
|
np.all(p <= self.grid[i][-1])):
|
||
|
raise ValueError("One of the requested xi is out of bounds "
|
||
|
"in dimension %d" % i)
|
||
|
out_of_bounds = None
|
||
|
else:
|
||
|
out_of_bounds = self._find_out_of_bounds(xi.T)
|
||
|
|
||
|
return xi, xi_shape, ndim, nans, out_of_bounds
|
||
|
|
||
|
def _evaluate_linear(self, indices, norm_distances):
|
||
|
# slice for broadcasting over trailing dimensions in self.values
|
||
|
vslice = (slice(None),) + (None,)*(self.values.ndim - len(indices))
|
||
|
|
||
|
# Compute shifting up front before zipping everything together
|
||
|
shift_norm_distances = [1 - yi for yi in norm_distances]
|
||
|
shift_indices = [i + 1 for i in indices]
|
||
|
|
||
|
# The formula for linear interpolation in 2d takes the form:
|
||
|
# values = self.values[(i0, i1)] * (1 - y0) * (1 - y1) + \
|
||
|
# self.values[(i0, i1 + 1)] * (1 - y0) * y1 + \
|
||
|
# self.values[(i0 + 1, i1)] * y0 * (1 - y1) + \
|
||
|
# self.values[(i0 + 1, i1 + 1)] * y0 * y1
|
||
|
# We pair i with 1 - yi (zipped1) and i + 1 with yi (zipped2)
|
||
|
zipped1 = zip(indices, shift_norm_distances)
|
||
|
zipped2 = zip(shift_indices, norm_distances)
|
||
|
|
||
|
# Take all products of zipped1 and zipped2 and iterate over them
|
||
|
# to get the terms in the above formula. This corresponds to iterating
|
||
|
# over the vertices of a hypercube.
|
||
|
hypercube = itertools.product(*zip(zipped1, zipped2))
|
||
|
value = np.array([0.])
|
||
|
for h in hypercube:
|
||
|
edge_indices, weights = zip(*h)
|
||
|
weight = np.array([1.])
|
||
|
for w in weights:
|
||
|
weight = weight * w
|
||
|
term = np.asarray(self.values[edge_indices]) * weight[vslice]
|
||
|
value = value + term # cannot use += because broadcasting
|
||
|
return value
|
||
|
|
||
|
def _evaluate_nearest(self, indices, norm_distances):
|
||
|
idx_res = [np.where(yi <= .5, i, i + 1)
|
||
|
for i, yi in zip(indices, norm_distances)]
|
||
|
return self.values[tuple(idx_res)]
|
||
|
|
||
|
def _validate_grid_dimensions(self, points, method):
|
||
|
k = self._SPLINE_DEGREE_MAP[method]
|
||
|
for i, point in enumerate(points):
|
||
|
ndim = len(np.atleast_1d(point))
|
||
|
if ndim <= k:
|
||
|
raise ValueError(f"There are {ndim} points in dimension {i},"
|
||
|
f" but method {method} requires at least "
|
||
|
f" {k+1} points per dimension.")
|
||
|
|
||
|
def _evaluate_spline(self, xi, method):
|
||
|
# ensure xi is 2D list of points to evaluate (`m` is the number of
|
||
|
# points and `n` is the number of interpolation dimensions,
|
||
|
# ``n == len(self.grid)``.)
|
||
|
if xi.ndim == 1:
|
||
|
xi = xi.reshape((1, xi.size))
|
||
|
m, n = xi.shape
|
||
|
|
||
|
# Reorder the axes: n-dimensional process iterates over the
|
||
|
# interpolation axes from the last axis downwards: E.g. for a 4D grid
|
||
|
# the order of axes is 3, 2, 1, 0. Each 1D interpolation works along
|
||
|
# the 0th axis of its argument array (for 1D routine it's its ``y``
|
||
|
# array). Thus permute the interpolation axes of `values` *and keep
|
||
|
# trailing dimensions trailing*.
|
||
|
axes = tuple(range(self.values.ndim))
|
||
|
axx = axes[:n][::-1] + axes[n:]
|
||
|
values = self.values.transpose(axx)
|
||
|
|
||
|
if method == 'pchip':
|
||
|
_eval_func = self._do_pchip
|
||
|
else:
|
||
|
_eval_func = self._do_spline_fit
|
||
|
k = self._SPLINE_DEGREE_MAP[method]
|
||
|
|
||
|
# Non-stationary procedure: difficult to vectorize this part entirely
|
||
|
# into numpy-level operations. Unfortunately this requires explicit
|
||
|
# looping over each point in xi.
|
||
|
|
||
|
# can at least vectorize the first pass across all points in the
|
||
|
# last variable of xi.
|
||
|
last_dim = n - 1
|
||
|
first_values = _eval_func(self.grid[last_dim],
|
||
|
values,
|
||
|
xi[:, last_dim],
|
||
|
k)
|
||
|
|
||
|
# the rest of the dimensions have to be on a per point-in-xi basis
|
||
|
shape = (m, *self.values.shape[n:])
|
||
|
result = np.empty(shape, dtype=self.values.dtype)
|
||
|
for j in range(m):
|
||
|
# Main process: Apply 1D interpolate in each dimension
|
||
|
# sequentially, starting with the last dimension.
|
||
|
# These are then "folded" into the next dimension in-place.
|
||
|
folded_values = first_values[j, ...]
|
||
|
for i in range(last_dim-1, -1, -1):
|
||
|
# Interpolate for each 1D from the last dimensions.
|
||
|
# This collapses each 1D sequence into a scalar.
|
||
|
folded_values = _eval_func(self.grid[i],
|
||
|
folded_values,
|
||
|
xi[j, i],
|
||
|
k)
|
||
|
result[j, ...] = folded_values
|
||
|
|
||
|
return result
|
||
|
|
||
|
@staticmethod
|
||
|
def _do_spline_fit(x, y, pt, k):
|
||
|
local_interp = make_interp_spline(x, y, k=k, axis=0)
|
||
|
values = local_interp(pt)
|
||
|
return values
|
||
|
|
||
|
@staticmethod
|
||
|
def _do_pchip(x, y, pt, k):
|
||
|
local_interp = PchipInterpolator(x, y, axis=0)
|
||
|
values = local_interp(pt)
|
||
|
return values
|
||
|
|
||
|
def _find_indices(self, xi):
|
||
|
return find_indices(self.grid, xi)
|
||
|
|
||
|
def _find_out_of_bounds(self, xi):
|
||
|
# check for out of bounds xi
|
||
|
out_of_bounds = np.zeros((xi.shape[1]), dtype=bool)
|
||
|
# iterate through dimensions
|
||
|
for x, grid in zip(xi, self.grid):
|
||
|
out_of_bounds += x < grid[0]
|
||
|
out_of_bounds += x > grid[-1]
|
||
|
return out_of_bounds
|
||
|
|
||
|
|
||
|
def interpn(points, values, xi, method="linear", bounds_error=True,
|
||
|
fill_value=np.nan):
|
||
|
"""
|
||
|
Multidimensional interpolation on regular or rectilinear grids.
|
||
|
|
||
|
Strictly speaking, not all regular grids are supported - this function
|
||
|
works on *rectilinear* grids, that is, a rectangular grid with even or
|
||
|
uneven spacing.
|
||
|
|
||
|
Parameters
|
||
|
----------
|
||
|
points : tuple of ndarray of float, with shapes (m1, ), ..., (mn, )
|
||
|
The points defining the regular grid in n dimensions. The points in
|
||
|
each dimension (i.e. every elements of the points tuple) must be
|
||
|
strictly ascending or descending.
|
||
|
|
||
|
values : array_like, shape (m1, ..., mn, ...)
|
||
|
The data on the regular grid in n dimensions. Complex data is
|
||
|
accepted.
|
||
|
|
||
|
.. deprecated:: 1.13.0
|
||
|
Complex data is deprecated with ``method="pchip"`` and will raise an
|
||
|
error in SciPy 1.15.0. This is because ``PchipInterpolator`` only
|
||
|
works with real values. If you are trying to use the real components of
|
||
|
the passed array, use ``np.real`` on ``values``.
|
||
|
|
||
|
xi : ndarray of shape (..., ndim)
|
||
|
The coordinates to sample the gridded data at
|
||
|
|
||
|
method : str, optional
|
||
|
The method of interpolation to perform. Supported are "linear",
|
||
|
"nearest", "slinear", "cubic", "quintic", "pchip", and "splinef2d".
|
||
|
"splinef2d" is only supported for 2-dimensional data.
|
||
|
|
||
|
bounds_error : bool, optional
|
||
|
If True, when interpolated values are requested outside of the
|
||
|
domain of the input data, a ValueError is raised.
|
||
|
If False, then `fill_value` is used.
|
||
|
|
||
|
fill_value : number, optional
|
||
|
If provided, the value to use for points outside of the
|
||
|
interpolation domain. If None, values outside
|
||
|
the domain are extrapolated. Extrapolation is not supported by method
|
||
|
"splinef2d".
|
||
|
|
||
|
Returns
|
||
|
-------
|
||
|
values_x : ndarray, shape xi.shape[:-1] + values.shape[ndim:]
|
||
|
Interpolated values at `xi`. See notes for behaviour when
|
||
|
``xi.ndim == 1``.
|
||
|
|
||
|
See Also
|
||
|
--------
|
||
|
NearestNDInterpolator : Nearest neighbor interpolation on unstructured
|
||
|
data in N dimensions
|
||
|
LinearNDInterpolator : Piecewise linear interpolant on unstructured data
|
||
|
in N dimensions
|
||
|
RegularGridInterpolator : interpolation on a regular or rectilinear grid
|
||
|
in arbitrary dimensions (`interpn` wraps this
|
||
|
class).
|
||
|
RectBivariateSpline : Bivariate spline approximation over a rectangular mesh
|
||
|
scipy.ndimage.map_coordinates : interpolation on grids with equal spacing
|
||
|
(suitable for e.g., N-D image resampling)
|
||
|
|
||
|
Notes
|
||
|
-----
|
||
|
|
||
|
.. versionadded:: 0.14
|
||
|
|
||
|
In the case that ``xi.ndim == 1`` a new axis is inserted into
|
||
|
the 0 position of the returned array, values_x, so its shape is
|
||
|
instead ``(1,) + values.shape[ndim:]``.
|
||
|
|
||
|
If the input data is such that input dimensions have incommensurate
|
||
|
units and differ by many orders of magnitude, the interpolant may have
|
||
|
numerical artifacts. Consider rescaling the data before interpolation.
|
||
|
|
||
|
Examples
|
||
|
--------
|
||
|
Evaluate a simple example function on the points of a regular 3-D grid:
|
||
|
|
||
|
>>> import numpy as np
|
||
|
>>> from scipy.interpolate import interpn
|
||
|
>>> def value_func_3d(x, y, z):
|
||
|
... return 2 * x + 3 * y - z
|
||
|
>>> x = np.linspace(0, 4, 5)
|
||
|
>>> y = np.linspace(0, 5, 6)
|
||
|
>>> z = np.linspace(0, 6, 7)
|
||
|
>>> points = (x, y, z)
|
||
|
>>> values = value_func_3d(*np.meshgrid(*points, indexing='ij'))
|
||
|
|
||
|
Evaluate the interpolating function at a point
|
||
|
|
||
|
>>> point = np.array([2.21, 3.12, 1.15])
|
||
|
>>> print(interpn(points, values, point))
|
||
|
[12.63]
|
||
|
|
||
|
"""
|
||
|
# sanity check 'method' kwarg
|
||
|
if method not in ["linear", "nearest", "cubic", "quintic", "pchip",
|
||
|
"splinef2d", "slinear",
|
||
|
"slinear_legacy", "cubic_legacy", "quintic_legacy"]:
|
||
|
raise ValueError("interpn only understands the methods 'linear', "
|
||
|
"'nearest', 'slinear', 'cubic', 'quintic', 'pchip', "
|
||
|
f"and 'splinef2d'. You provided {method}.")
|
||
|
|
||
|
if not hasattr(values, 'ndim'):
|
||
|
values = np.asarray(values)
|
||
|
|
||
|
ndim = values.ndim
|
||
|
if ndim > 2 and method == "splinef2d":
|
||
|
raise ValueError("The method splinef2d can only be used for "
|
||
|
"2-dimensional input data")
|
||
|
if not bounds_error and fill_value is None and method == "splinef2d":
|
||
|
raise ValueError("The method splinef2d does not support extrapolation.")
|
||
|
|
||
|
# sanity check consistency of input dimensions
|
||
|
if len(points) > ndim:
|
||
|
raise ValueError("There are %d point arrays, but values has %d "
|
||
|
"dimensions" % (len(points), ndim))
|
||
|
if len(points) != ndim and method == 'splinef2d':
|
||
|
raise ValueError("The method splinef2d can only be used for "
|
||
|
"scalar data with one point per coordinate")
|
||
|
|
||
|
grid, descending_dimensions = _check_points(points)
|
||
|
_check_dimensionality(grid, values)
|
||
|
|
||
|
# sanity check requested xi
|
||
|
xi = _ndim_coords_from_arrays(xi, ndim=len(grid))
|
||
|
if xi.shape[-1] != len(grid):
|
||
|
raise ValueError("The requested sample points xi have dimension "
|
||
|
"%d, but this RegularGridInterpolator has "
|
||
|
"dimension %d" % (xi.shape[-1], len(grid)))
|
||
|
|
||
|
if bounds_error:
|
||
|
for i, p in enumerate(xi.T):
|
||
|
if not np.logical_and(np.all(grid[i][0] <= p),
|
||
|
np.all(p <= grid[i][-1])):
|
||
|
raise ValueError("One of the requested xi is out of bounds "
|
||
|
"in dimension %d" % i)
|
||
|
|
||
|
# perform interpolation
|
||
|
if method in RegularGridInterpolator._ALL_METHODS:
|
||
|
interp = RegularGridInterpolator(points, values, method=method,
|
||
|
bounds_error=bounds_error,
|
||
|
fill_value=fill_value)
|
||
|
return interp(xi)
|
||
|
elif method == "splinef2d":
|
||
|
xi_shape = xi.shape
|
||
|
xi = xi.reshape(-1, xi.shape[-1])
|
||
|
|
||
|
# RectBivariateSpline doesn't support fill_value; we need to wrap here
|
||
|
idx_valid = np.all((grid[0][0] <= xi[:, 0], xi[:, 0] <= grid[0][-1],
|
||
|
grid[1][0] <= xi[:, 1], xi[:, 1] <= grid[1][-1]),
|
||
|
axis=0)
|
||
|
result = np.empty_like(xi[:, 0])
|
||
|
|
||
|
# make a copy of values for RectBivariateSpline
|
||
|
interp = RectBivariateSpline(points[0], points[1], values[:])
|
||
|
result[idx_valid] = interp.ev(xi[idx_valid, 0], xi[idx_valid, 1])
|
||
|
result[np.logical_not(idx_valid)] = fill_value
|
||
|
|
||
|
return result.reshape(xi_shape[:-1])
|
||
|
else:
|
||
|
raise ValueError(f"unknown {method = }")
|