scipy.interpolate.interpn#

scipy.interpolate.interpn(points, values, xi, method='linear', bounds_error=True, fill_value=nan)[source]#

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
pointstuple 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.

valuesarray_like, shape (m1, …, mn, …)

The data on the regular grid in n dimensions. Complex data can be acceptable.

xindarray of shape (…, ndim)

The coordinates to sample the gridded data at

methodstr, optional

The method of interpolation to perform. Supported are “linear” and “nearest”, and “splinef2d”. “splinef2d” is only supported for 2-dimensional data.

bounds_errorbool, 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_valuenumber, 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_xndarray, shape xi.shape[:-1] + values.shape[ndim:]

Interpolated values at input coordinates.

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

New in version 0.14.

Examples

Evaluate a simple example function on the points of a regular 3-D grid:

>>> 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]