bundles / scipy 1.17.1 / scipy / interpolate / _rgi / RegularGridInterpolator
class
scipy.interpolate._rgi:RegularGridInterpolator
source: /scipy/interpolate/_rgi.py :59
Signature
class RegularGridInterpolator ( points , values , method = linear , bounds_error = True , fill_value = nan , * , solver = None , solver_args = None ) Members
Summary
Interpolator of specified order on a rectilinear grid in N ≥ 1 dimensions.
Extended Summary
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.
method: str, optionalThe 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, optionalIf True, when interpolated values are requested outside of the domain of the input data, a ValueError is raised. If False, then
fill_valueis used. Default is True.fill_value: float or None, optionalThe value to use for points outside of the interpolation domain. If None, values outside the domain are extrapolated. Default is
np.nan.solver: callable, optionalOnly 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.
solver_args: dict, optionalAdditional arguments to pass to
solver, if any.
Attributes
grid: tuple of ndarraysThe points defining the regular grid in n dimensions. This tuple defines the full grid via
np.meshgrid(*grid, indexing='ij')values: ndarrayData values at the grid.
method: strInterpolation method.
fill_value: float or ``None``Use this value for out-of-bounds arguments to __call__.
bounds_error: boolIf
True, out-of-bounds argument raise aValueError.
Methods
__call__
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.
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.
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.
Rounding rule at half points with `nearest` method
The rounding rule with the nearest method at half points is rounding down.
Array API Standard Support
RegularGridInterpolator has experimental support for Python Array API Standard compatible backends in addition to NumPy. Please consider testing these features by setting an environment variable SCIPY_ARRAY_API=1 and providing CuPy, PyTorch, JAX, or Dask arrays as array arguments. The following combinations of backend and device (or other capability) are supported.
==================== ==================== ==================== Library CPU GPU ==================== ==================== ==================== NumPy ✅ n/a CuPy n/a ⛔ PyTorch ✅ ⛔ JAX ⚠️ no JIT ⛔ Dask ⛔ n/a ==================== ==================== ====================
See
dev-arrayapifor more information.
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)✓
interp = RegularGridInterpolator((x, y, z), data)
✓pts = np.array([[2.1, 6.2, 8.3], [3.3, 5.2, 7.1]])✓
interp(pts)
✗f(2.1, 6.2, 8.3), f(3.3, 5.2, 7.1)
✓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')✗
xx = np.linspace(-4, 9, 31) yy = np.linspace(-4, 9, 31) X, Y = np.meshgrid(xx, yy, indexing='ij')✓
ax.plot_wireframe(X, Y, interp((X, Y)), rstride=3, cstride=3, alpha=0.4, color='m', label='linear interp')✗
ax.plot_wireframe(X, Y, ff(X, Y), rstride=3, cstride=3, alpha=0.4, label='ground truth') plt.legend()✗
plt.show()
✓
See also
- LinearNDInterpolator
Piecewise linear interpolator on unstructured data in N dimensions
- NearestNDInterpolator
Nearest neighbor 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)
Aliases
-
scipy.interpolate.RegularGridInterpolator