geocat.comp.rcm2points#
- geocat.comp.rcm2points(lat2d, lon2d, fi, lat1d, lon1d, opt=0, msg=None, meta=False)#
Interpolates data on a curvilinear grid (i.e. RCM, WRF, NARR) to an unstructured grid.
Interpolates data on a curvilinear grid, such as those used by the RCM (Regional Climate Model), WRF (Weather Research and Forecasting) and NARR (North American Regional Reanalysis) models/datasets to an unstructured grid. All of these have latitudes that are oriented south-to-north. An inverse distance squared algorithm is used to perform the interpolation. Missing values are allowed and no extrapolation is performed.
- Parameters:
lat2d (
xarray.DataArray
,numpy.ndarray
) – A two-dimensional array that specifies the latitudes locations offi
. The latitude order must be south-to-north. Because this array is two-dimensional it is not an associated coordinate variable offi
, so it should always be explicitly provided.lon2d (
xarray.DataArray
,numpy.ndarray
) – A two-dimensional array that specifies the longitude locations offi
. The longitude order must be west-to-east. Because this array is two-dimensional it is not an associated coordinate variable offi
, so it should always be explicitly provided.fi (
xarray.DataArray
,numpy.ndarray
) – A multi-dimensional array to be interpolated. The rightmost two dimensions (latitude, longitude) are the dimensions to be interpolated.lat1d (
xarray.DataArray
,numpy.ndarray
) – A one-dimensional array that specifies the latitude coordinates of the output locations.lon1d (
xarray.DataArray
,numpy.ndarray
) – A one-dimensional array that specifies the longitude coordinates of the output locations.opt (
numpy.number
) –opt=0
or1
means use an inverse distance weight interpolation.opt=2
means use a bilinear interpolation.msg (
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.meta (
bool
) – If set to True and the input array is an Xarray, the metadata from the input array will be copied to the output array; default is False. Warning: This option is not currently supported.
- Returns:
fo (
xarray.DataArray
,numpy.ndarray
) – The interpolated grid. A multi-dimensional array of the same size asfi
except that the rightmost dimension sizes have been replaced by the number of coordinate pairs (lat1d
,lon1d
).