geocat.comp.linint2#
- geocat.comp.linint2(fi, xo, yo, xi=None, yi=None, icycx=0, msg_py=None)#
Interpolates a regular grid to a rectilinear one using bi-linear interpolation. The input grid may be cyclic in the x-direction. The interpolation is first performed in the x-direction, and then in the y-direction.
- Parameters:
fi (
xarray.DataArray
,numpy.ndarray
) – An array of two or more dimensions. Ifxi
is passed in as an argument, then the size of the rightmost dimension offi
must match the rightmost dimension ofxi
. Similarly, ifyi
is passed in as an argument, then the size of the second-rightmost dimension offi
must match the rightmost dimension ofyi
.If missing values are present, then linint2 will perform the bilinear interpolation at all points possible, but will return missing values at coordinates which could not be used.
Note: This variable must be supplied as a
xarray.DataArray
in order to copy the dimension names to the output. Otherwise, default names will be used.xo (
xarray.DataArray
,numpy.ndarray
) – A one-dimensional array that specifies the X-coordinates of the return array. It must be strictly monotonically increasing, but may be unequally spaced. For geo-referenced data,xo
is generally the longitude array.If the output coordinates
xo
are outside those of the input coordinatesxi
, then thefo
values at those coordinates will be set to missing (i.e. no extrapolation is performed).yo (
xarray.DataArray
,numpy.ndarray
) – A one-dimensional array that specifies the Y coordinates of the return array. It must be strictly monotonically increasing, but may be unequally spaced. For geo-referenced data,yo
is typically the latitude array.If the output coordinates
yo
are outside those of the input coordinatesyi
, then thefo
values at those coordinates will be set to missing (i.e. no extrapolation is performed).xi (
xarray.DataArray
,numpy.ndarray
) – An array that specifies the X-coordinates of thefi
array. Most frequently, this is a 1D strictly monotonically increasing array that may be unequally spaced. In some cases,xi
can be a multi-dimensional array (see next paragraph). The rightmost dimension (call itnxi
) must have at least two elements, and is the last (fastest varying) dimension offi
.If
xi
is a multi-dimensional array, then eachnxi
subsection ofxi
must be strictly monotonically increasing, but may be unequally spaced. All but its rightmost dimension must be the same size as all butfi
’s rightmost two dimensions. For geo-referenced data,xi
is generally the longitude array.Note: If
fi
is of typexarray.DataArray
andxi
is left unspecified, then the rightmost coordinate dimension offi
will be used. Iffi
is not of typexarray.DataArray
, thenxi
becomes a mandatory parameter. This parameter must be specified as a keyword argument.yi (
xarray.DataArray
,numpy.ndarray
) – An array that specifies the Y-coordinates of thefi
array. Most frequently, this is a 1D strictly monotonically increasing array that may be unequally spaced. In some cases,yi
can be a multi-dimensional array (see next paragraph). The rightmost dimension (call it nyi) must have at least two elements, and is the second-to-last dimension offi
.If
yi
is a multi-dimensional array, then each nyi subsection ofyi
must be strictly monotonically increasing, but may be unequally spaced. All but its rightmost dimension must be the same size as all butfi
’s rightmost two dimensions. For geo-referenced data,yi
is generally the latitude array.Note: If
fi
is of typexarray.DataArray
andxi
is left unspecified, then the second-to-rightmost coordinate dimension offi
will be used. Iffi
is not of typexarray.DataArray
, thenxi
becomes a mandatory parameter. This parameter must be specified as a keyword argument.icycx (
bool
) – An option to indicate whether the rightmost dimension offi
is cyclic. This should be set to True only if you have global data, but your longitude values don’t quite wrap all the way around the globe. For example, if your longitude values go from, say, -179.75 to 179.75, or 0.5 to 359.5, then you would set this to True.msg_py (
numpy.number
) – A numpy scalar value that represent a missing value infi
. This argument allows a user to use a missing value scheme other than NaN or masked arrays, similar to what NCL allows.
- Returns:
fo (
xarray.DataArray
,numpy.ndarray
) – The interpolated grid. If the meta parameter is True, then the result will include named dimensions matching the input array. The returned value will have the same dimensions asfi
, except for the rightmost two dimensions which will have the same dimension sizes as the lengths ofyo
andxo
. The return type will be double iffi
is double, and float otherwise.
Examples
Example 1: Using
linint2
withxarray.DataArray
inputimport numpy as np import xarray as xr import geocat.comp fi_np = np.random.rand(30, 80) # random 30x80 array # xi and yi do not have to be equally spaced, but they are # in this example xi = np.arange(80) yi = np.arange(30) # create target coordinate arrays, in this case use the same # min/max values as xi and yi, but with different spacing xo = np.linspace(xi.min(), xi.max(), 100) yo = np.linspace(yi.min(), yi.max(), 50) # create :class:`xarray.DataArray` and chunk it using the # full shape of the original array. # note that xi and yi are attached as coordinate arrays fi = xr.DataArray(fi_np, dims=['lat', 'lon'], coords={'lat': yi, 'lon': xi} ).chunk(fi_np.shape) fo = geocat.comp.linint2(fi, xo, yo, icycx=0)