pycontrails.core.interpolation

Interpolation utilities.

Functions

interp(longitude, latitude, level, time, da, ...)

Interpolate over a grid with localize option.

Classes

EmissionsProfileInterpolator(xp, fp[, ...])

Support for interpolating a profile on a linear or logarithmic scale.

PycontrailsRegularGridInterpolator(points, ...)

Support for performant interpolation over a regular grid.

RGIArtifacts(xi_indices, norm_distances, ...)

An interface to intermediate RGI interpolation artifacts.

class pycontrails.core.interpolation.EmissionsProfileInterpolator(xp, fp, drop_duplicates=True)

Bases: object

Support for interpolating a profile on a linear or logarithmic scale.

This class simply wraps numpy.interp() with fixed values for the xp and fp arguments. Unlike xarray.DataArray interpolation, the numpy.interp automatically clips values outside the range of the xp array.

Parameters:
  • xp (npt.NDarray[np.float64]) – Array of x-values. These must be strictly increasing and free from any nan values. Passed to numpy.interp().

  • fp (npt.NDarray[np.float64]) – Array of y-values. Passed to numpy.interp().

  • drop_duplicates (bool, optional) – Whether to drop duplicate values in xp. Defaults to True.

Examples

>>> xp = np.array([3, 7, 10, 30], dtype=float)
>>> fp = np.array([0.1, 0.2, 0.3, 0.4], dtype=float)
>>> epi = EmissionsProfileInterpolator(xp, fp)
>>> # Interpolate a single value
>>> epi.interp(5)
np.float64(0.150000...)
>>> # Interpolate a single value on a logarithmic scale
>>> epi.log_interp(5)
np.float64(1.105171...)
>>> # Demonstrate speed up compared with xarray.DataArray interpolation
>>> import time, xarray as xr
>>> da = xr.DataArray(fp, dims=["x"], coords={"x": xp})
>>> inputs = [np.random.uniform(0, 31, size=200) for _ in range(1000)]
>>> t0 = time.perf_counter()
>>> xr_out = [da.interp(x=x.clip(3, 30)).values for x in inputs]
>>> t1 = time.perf_counter()
>>> np_out = [epi.interp(x) for x in inputs]
>>> t2 = time.perf_counter()
>>> assert np.allclose(xr_out, np_out)
>>> # We see a 100 fold speed up (more like 500x faster, but we don't
>>> # want the test to fail!)
>>> assert t2 - t1 < (t1 - t0) / 100
interp(x)

Interpolate x against xp and fp.

Parameters:

x (npt.NDArray[np.float64]) – Array of x-values to interpolate.

Returns:

npt.NDArray[np.float64] – Array of interpolated y-values arising from the x-values. The dtype of the output array is the same as the dtype of x.

log_interp(x)

Interpolate x against xp and fp on a logarithmic scale.

This method composes the following three functions.
  1. numpy.log()

  2. interp()

  3. numpy.exp()

Parameters:

x (npt.NDArray[np.float64]) – Array of x-values to interpolate.

Returns:

npt.NDArray[np.float64] – Array of interpolated y-values arising from the x-values.

class pycontrails.core.interpolation.PycontrailsRegularGridInterpolator(points, values, method, bounds_error, fill_value)

Bases: RegularGridInterpolator

Support for performant interpolation over a regular grid.

This class is a thin wrapper around the scipy.interpolate.RegularGridInterpolator in order to make typical pycontrails use-cases more efficient.

#. Avoid RegularGridInterpolator constructor validation. In interp(), parameters are carefully crafted to fit into the intended form, thereby making validation unnecessary. #. Override the _evaluate_linear() method with a faster implementation. See the _evaluate_linear() docstring for more information.

This class should not be used directly. Instead, use the interp() function.

Changed in version 0.40.0: The _evaluate_linear() method now uses a Cython implementation. The dtype of the output is now consistent with the dtype of the underlying values

Parameters:
class pycontrails.core.interpolation.RGIArtifacts(xi_indices, norm_distances, out_of_bounds)

Bases: object

An interface to intermediate RGI interpolation artifacts.

norm_distances
out_of_bounds
xi_indices
pycontrails.core.interpolation.interp(longitude, latitude, level, time, da, method, bounds_error, fill_value, localize, *, indices=None, return_indices=False)

Interpolate over a grid with localize option.

Changed in version 0.25.6: Utilize scipy 1.9 upgrades to remove singleton dimensions.

Changed in version 0.26.0: Include indices and return_indices experimental parameters. Currently, nan values in longitude, latitude, level, or time are always propagated through to the output, regardless of bounds_error. In other words, a ValueError for an out of bounds coordinate is only raised if a non-nan value is out of bounds.

Changed in version 0.40.0: When return_indices is True, an instance of RGIArtifacts is used to store the indices artifacts.

Parameters:
  • longitude, latitude, level, time (numpy.ndarray) – Coordinates of points to be interpolated. These parameters have the same meaning as x in analogy with numpy.interp(). All four of these arrays must be 1 dimensional of the same size.

  • da (xarray.DataArray) – Gridded data interpolated over. Must adhere to MetDataArray conventions. In particular, the dimensions of da must be longitude, latitude, level, and time. The three spatial dimensions must be monotonically increasing with float64 dtype. The time dimension must be monotonically increasing with datetime64 dtype. Assumed to be cheap to load into memory (xr.DataArray.values is used without hesitation).

  • method (str) – Passed into scipy.interpolate.RegularGridInterpolator.

  • bounds_error (bool) – Passed into scipy.interpolate.RegularGridInterpolator.

  • fill_value (float | np.float64 | None) – Passed into scipy.interpolate.RegularGridInterpolator.

  • localize (bool) – If True, clip da to the smallest bounding box that contains all of coords.

  • indices (tuple | None, optional) – Experimental. Provide intermediate artifacts computed by :meth:scipy.interpolate.RegularGridInterpolator._find_indices` to avoid redundant computation. If known and provided, this can speed up interpolation by avoiding an unnecessary call to ``_find_indices. By default, None. Must be used precisely.

  • return_indices (bool, optional) – If True, return output of scipy.interpolate.RegularGridInterpolator._find_indices() in addition to interpolated values.

Returns:

npt.NDArray[np.float64] | tuple[npt.NDArray[np.float64], RGIArtifacts] – Interpolated values with same size as longitude. If return_indices is True, return intermediate indices artifacts as well.

See also

-

meth:MetDataArray.interpolate

-

func:scipy.interpolate.interpn

-

class:scipy.interpolate.RegularGridInterpolator