pyinterp.backends.xarray.RegularGridInterpolator.__call__

pyinterp.backends.xarray.RegularGridInterpolator.__call__#

RegularGridInterpolator.__call__(coords, method='bilinear', **kwargs)[source]#

Interpolate at coordinates.

Perform interpolation at the specified coordinates using the chosen method and parameters.

Parameters:
  • coords (dict) – Mapping from dimension names to the new coordinates. Coordinates can be scalars or array-like. For temporal axes, provide datetime64 arrays.

  • method (Literal['bilinear', 'idw', 'nearest', 'akima', 'akima_periodic', 'bicubic', 'c_spline', 'c_spline_not_a_knot', 'c_spline_periodic', 'linear', 'polynomial', 'steffen']) –

    The method of interpolation to perform. Supported methods depend on the grid type. Common methods include:

    • Geometric methods: nearest, bilinear, idw

    • Windowed methods: akima, akima_periodic, bicubic, bilinear, c_spline, c_spline_not_a_knot, c_spline_periodic, linear, polynomial, steffen.

  • **kwargs (Any) –

    Additional keyword arguments passed to the interpolation function. Common options include:

    • bounds_error (bool): Raise error if coordinates are out of bounds. Default: False (returns NaN).

    • num_threads (int): Number of threads for parallel computation. 0 uses all CPUs. Default: 0.

    For windowed methods (bicubic, c_spline, etc.), additional options include:

    • half_window_size_x (int): Half window size in X direction

    • half_window_size_y (int): Half window size in Y direction

    • boundary_mode (str): Boundary handling mode ("shrink", "undef")

    For 3D/4D grids:

    • third_axis (str): Method for 3rd axis ("linear", "nearest")

    • fourth_axis (str): Method for 4th axis ("linear", "nearest")

Returns:

Interpolated values as numpy array with same shape as input coordinate arrays.

Raises:
  • ValueError – If bounds_error=True and coordinates are out of bounds.

  • IndexError – If coordinate dimensions don’t match grid dimensions.

Return type:

ndarray

Examples

Simple bilinear interpolation

>>> result = interp(
...     {"lon": [10.5], "lat": [45.2]}, method="bilinear"
... )

Bicubic with custom window size

>>> result = interp(
...     {"lon": [10.5], "lat": [45.2]},
...     method="bicubic",
...     half_window_size_x=10,
...     half_window_size_y=10,
... )

With bounds checking

>>> result = interp(
...     {"lon": [10.5], "lat": [45.2]},
...     method="bilinear",
...     bounds_error=True,
... )

Multi-threaded

>>> result = interp(
...     {"lon": lon_array, "lat": lat_array},
...     method="bilinear",
...     num_threads=4,
... )

3D with temporal axis

>>> result = interp(
...     {
...         "lon": [10.5],
...         "lat": [45.2],
...         "time": np.array(["2020-01-01"], dtype="datetime64"),
...     },
...     method="bilinear",
... )