geocat.comp.rcm2points

geocat.comp.rcm2points(lat2d, lon2d, fi, lat1dPoints, lon1dPoints, opt=0, msg=None, meta=False)

Interpolates data on a curvilinear grid (i.e. RCM, WRF, NARR) to an unstructured grid.

Parameters
  • lat2d (numpy.ndarray) – A two-dimensional array that specifies the latitudes locations of fi. The latitude order must be south-to-north.

  • lon2d (numpy.ndarray) – A two-dimensional array that specifies the longitude locations of fi. The latitude order must be west-to-east.

  • fi (numpy.ndarray) – A multi-dimensional array to be interpolated. The rightmost two dimensions (latitude, longitude) are the dimensions to be interpolated.

  • lat1dPoints (numpy.ndarray) – A one-dimensional array that specifies the latitude coordinates of the output locations.

  • lon1dPoints (numpy.ndarray) – A one-dimensional array that specifies the longitude coordinates of the output locations.

  • opt (numpy.number) – opt=0 or 1 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 in fi. This argument allows a user to use a missing value scheme other than NaN or masked arrays, similar to what NCL allows.

  • meta (bool) –

  • set to True and the input array is an Xarray, the metadata (If) –

  • the input array will be copied to the output array; (from) –

  • is False. (default) –

  • Warning – this option is not currently supported.

Returns

The interpolated grid. A multi-dimensional array of the same size as fi except that the rightmost dimension sizes have been replaced by the number of coordinate pairs (lat1dPoints, lon1dPoints). Double if fi is double, otherwise float.

Return type

numpy.ndarray

Description:

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.

A inverse distance squared algorithm is used to perform the interpolation.

Missing values are allowed and no extrapolation is performed.

Examples

Example 1: Using rcm2points with xarray.DataArray input

import numpy as np
import xarray as xr
import geocat.comp

# Open a netCDF data file using xarray default engine and load the data stream
ds = xr.open_dataset("./ruc.nc")

# [INPUT] Grid & data info on the source curvilinear
ht_curv=ds.DIST_236_CBL[:]
lat2D_curv=ds.gridlat_236[:]
lon2D_curv=ds.gridlon_236[:]

# [OUTPUT] Grid on destination points grid (or read the 1D lat and lon from
#          an other .nc file.
newlat1D_points=np.linspace(lat2D_curv.min(), lat2D_curv.max(), 100)
newlon1D_points=np.linspace(lon2D_curv.min(), lon2D_curv.max(), 100)

ht_points = geocat.comp.rcm2points(lat2D_curv, lon2D_curv, ht_curv, newlat1D_points, newlon1D_points)