geocat.comp.linint2pts#
- geocat.comp.linint2pts(fi, xo, yo, icycx=False, msg_py=None, xi=None, yi=None)#
Interpolates from a rectilinear grid to an unstructured grid or locations using bilinear interpolation.
The
linint2ptsfunction uses bilinear interpolation to interpolate from a rectilinear grid to an unstructured grid.If missing values are present, then
linint2ptswill perform the piecewise linear interpolation at all points possible, but will return missing values at coordinates which could not be used.If one or more of the four closest grid points to a particular (
xo,yo) coordinate pair are missing, then the return value for this coordinate pair will be missing.If the user inadvertently specifies output coordinates (
xo,yo) that are outside those of the input coordinates (xi,yi), the output value at this coordinate pair will be set to missing as no extrapolation is performed.linint2ptsis different fromlinint2in thatxoandyoare coordinate pairs, and need not be monotonically increasing. It is also different in the dimensioning of the return array. This function could be used if the user wanted to interpolate gridded data to, say, the location of rawinsonde sites or buoy/xbt locations.Warning: if
xicontains longitudes, then thexovalues must be in the same range. In addition, if thexi`values span 0 to 360, then thexovalues must also be specified in this range (i.e. -180 to 180 will not work).- Parameters:
fi (
xarray.DataArray,numpy.ndarray) – An array of two or more dimensions. The two rightmost dimensions (nyixnxi) are the dimensions to be used in the interpolation. If user-defined missing values are present (other than NaNs), the value ofmsg_pymust be set appropriately.xo (
xarray.DataArray,numpy.ndarray) – A one-dimensional array that specifies the X-coordinates of the unstructured grid.yo (
xarray.DataArray,numpy.ndarray) – A one-dimensional array that specifies the Y-coordinates of the unstructured grid. It must be the same length asxo.icycx (
bool) – An option to indicate whether the rightmost dimension offiis cyclic. Default valus is 0. 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.xi (
xarray.DataArray,numpy.ndarray) – A strictly monotonically increasing array that specifies the X-coordinates of thefiarray.ximight be defined as the coordinates offiwhenfiis of typexarray.DataArray; in this caseximay not be explicitly given as a function argument.yi (
xarray.DataArray,numpy.ndarray) – A strictly monotonically increasing array that specifies the Y [latitude] coordinates of thefiarray.yimight be defined as the coordinates offiwhenfiis of typexarray.DataArray; in this caseyimay not be explicitly given as a function argument.
- Returns:
fo (
xarray.DataArray,numpy.ndarray) – The returned value will have the same dimensions asfi, except for the rightmost dimension which will have the same dimension size as the length ofyoandxo.
Examples
Example 1: Using
linint2ptswithxarray.DataArrayinputimport 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.linint2pts(fi, xo, yo, 0)