🍾 Xarray is now 10 years old! 🎉

xarray.DataArray.interp

xarray.DataArray.interp#

DataArray.interp(coords=None, method='linear', assume_sorted=False, kwargs=None, **coords_kwargs)[source]#

Interpolate a DataArray onto new coordinates

Performs univariate or multivariate interpolation of a DataArray onto new coordinates using scipy’s interpolation routines. If interpolating along an existing dimension, scipy.interpolate.interp1d is called. When interpolating along multiple existing dimensions, an attempt is made to decompose the interpolation into multiple 1-dimensional interpolations. If this is possible, scipy.interpolate.interp1d is called. Otherwise, scipy.interpolate.interpn() is called.

Parameters:
  • coords (dict, optional) – Mapping from dimension names to the new coordinates. New coordinate can be a scalar, array-like or DataArray. If DataArrays are passed as new coordinates, their dimensions are used for the broadcasting. Missing values are skipped.

  • method ({"linear", "nearest", "zero", "slinear", "quadratic", "cubic", "polynomial"}, default: "linear") – The method used to interpolate. The method should be supported by the scipy interpolator:

    • interp1d: {“linear”, “nearest”, “zero”, “slinear”, “quadratic”, “cubic”, “polynomial”}

    • interpn: {“linear”, “nearest”}

    If "polynomial" is passed, the order keyword argument must also be provided.

  • assume_sorted (bool, default: False) – If False, values of x can be in any order and they are sorted first. If True, x has to be an array of monotonically increasing values.

  • kwargs (dict-like or None, default: None) – Additional keyword arguments passed to scipy’s interpolator. Valid options and their behavior depend whether interp1d or interpn is used.

  • **coords_kwargs ({dim: coordinate, ...}, optional) – The keyword arguments form of coords. One of coords or coords_kwargs must be provided.

Returns:

interpolated (DataArray) – New dataarray on the new coordinates.

Notes

scipy is required.

See also

scipy.interpolate.interp1d scipy.interpolate.interpn

Manipulating Dimensions (Data Resolution)

Tutorial material on manipulating data resolution using interp()

Examples

>>> da = xr.DataArray(
...     data=[[1, 4, 2, 9], [2, 7, 6, np.nan], [6, np.nan, 5, 8]],
...     dims=("x", "y"),
...     coords={"x": [0, 1, 2], "y": [10, 12, 14, 16]},
... )
>>> da
<xarray.DataArray (x: 3, y: 4)> Size: 96B
array([[ 1.,  4.,  2.,  9.],
       [ 2.,  7.,  6., nan],
       [ 6., nan,  5.,  8.]])
Coordinates:
  * x        (x) int64 24B 0 1 2
  * y        (y) int64 32B 10 12 14 16

1D linear interpolation (the default):

>>> da.interp(x=[0, 0.75, 1.25, 1.75])
<xarray.DataArray (x: 4, y: 4)> Size: 128B
array([[1.  , 4.  , 2.  ,  nan],
       [1.75, 6.25, 5.  ,  nan],
       [3.  ,  nan, 5.75,  nan],
       [5.  ,  nan, 5.25,  nan]])
Coordinates:
  * y        (y) int64 32B 10 12 14 16
  * x        (x) float64 32B 0.0 0.75 1.25 1.75

1D nearest interpolation:

>>> da.interp(x=[0, 0.75, 1.25, 1.75], method="nearest")
<xarray.DataArray (x: 4, y: 4)> Size: 128B
array([[ 1.,  4.,  2.,  9.],
       [ 2.,  7.,  6., nan],
       [ 2.,  7.,  6., nan],
       [ 6., nan,  5.,  8.]])
Coordinates:
  * y        (y) int64 32B 10 12 14 16
  * x        (x) float64 32B 0.0 0.75 1.25 1.75

1D linear extrapolation:

>>> da.interp(
...     x=[1, 1.5, 2.5, 3.5],
...     method="linear",
...     kwargs={"fill_value": "extrapolate"},
... )
<xarray.DataArray (x: 4, y: 4)> Size: 128B
array([[ 2. ,  7. ,  6. ,  nan],
       [ 4. ,  nan,  5.5,  nan],
       [ 8. ,  nan,  4.5,  nan],
       [12. ,  nan,  3.5,  nan]])
Coordinates:
  * y        (y) int64 32B 10 12 14 16
  * x        (x) float64 32B 1.0 1.5 2.5 3.5

2D linear interpolation:

>>> da.interp(x=[0, 0.75, 1.25, 1.75], y=[11, 13, 15], method="linear")
<xarray.DataArray (x: 4, y: 3)> Size: 96B
array([[2.5  , 3.   ,   nan],
       [4.   , 5.625,   nan],
       [  nan,   nan,   nan],
       [  nan,   nan,   nan]])
Coordinates:
  * x        (x) float64 32B 0.0 0.75 1.25 1.75
  * y        (y) int64 24B 11 13 15