geocat.comp.interpolation.interp_multidim#

geocat.comp.interpolation.interp_multidim(data_in, lat_out, lon_out, lat_in=None, lon_in=None, cyclic=False, missing_val=None, method='linear', fill_value=nan)#

Multidimensional interpolation of variables. Uses xarray.interp to perform linear interpolation. Will not perform extrapolation, returns missing values if any surrounding points contain missing values.

Parameters
  • data_in (xarray.DataArray, numpy.ndarray) – Data array with data to be interpolated and associated coords. If it is a np array, then lat_in and lon_in must be provided. Length must be coordinated with given coordinates.

  • lat_out (numpy.ndarray) – List of latitude coordinates to be interpolated to.

  • lon_out (numpy.ndarray) – List of longitude coordinates to be interpolated to.

  • lat_in (numpy.ndarray) – List of latitude coordinates corresponding to data_in. Must be given if data_in is not an xarray.

  • lon_in (numpy.ndarray) – List of longitude coordinates corresponding to data_in. Must be given if data_in is not an xarray.

  • cyclic (bool, optional) – Set as true if lon values are cyclical but do not fully wrap around the globe (0, 1.5, 3, …, 354, 355.5) Default is false

  • missing_val (np.number, optional) – Provide a number to represent missing data. Alternative to using np.nan

  • method (str, optional) – Provide specific method of interpolation. Default is “linear” “linear” or “nearest” for multidimensional array

  • fill_value (str, optional) – Set as ‘extrapolate’ to allow extrapoltion of data. Default is no extrapolation.

Returns

data_out (numpy.ndarray, xarray.DataArray) – Returns same data type as input data_in. Shape will be the same as input array except for last two dimensions which will be equal to coordinates given in data_out

Examples

>>> import xarray as xr
>>> import numpy as np
>>> import geocat.comp
>>> data = np.asarray([[1, 2, 3, 4, 5, 99], [2, 4, 6, 8, 10, 12]])
>>> lat_in = [0, 1]
>>> lon_in = [0, 50, 100, 250, 300, 350]
>>> data_in = xr.DataArray(data, dims=['lat', 'lon'], coords={'lat':
lat_in, 'lon': lon_in})
>>> data_out = xr.DataArray(dims=['lat', 'lon'], coords={'lat': [0, 1],
'lon': [0, 50, 360]})
>>> do = interp_multidim(data_in, [0, 1], [0, 50, 360], cyclic=True, missing_val=99)
>>> print(do)
<xarray.DataArray (lat: 2, lon: 3)>
array([[ 1.,  2., 99.],
   [ 2.,  4., 99.]])
Coordinates:
  * lat      (lat) int64 0 1
  * lon      (lon) int64 0 50 360