geocat.comp.rcm2rgrid#

geocat.comp.rcm2rgrid(lat2d, lon2d, fi, lat1d, lon1d, msg=None, meta=False)#

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

Parameters
  • lat2d (xarray.DataArray, numpy.ndarray:) – A two-dimensional array that specifies the latitudes locations of the input (fi). Because this array is two-dimensional, it is not an associated coordinate variable of fi; therefore, it always needs to be explicitly provided. The latitude order must be south-to-north.

  • lon2d (xarray.DataArray, numpy.ndarray:) – A two-dimensional array that specifies the longitudes locations of the input (fi). Because this array is two-dimensional, it is not an associated coordinate variable of fi; therefore, it always needs to be explicitly provided. The latitude order must be west-to-east.

  • 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 regular grid. Must be monotonically increasing.

  • lon1d (xarray.DataArray, numpy.ndarray:) – A one-dimensional array that specifies the longitude coordinates of the regular grid. Must be monotonically increasing.

  • 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:) – 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 as fi except that the rightmost dimension sizes have been replaced by the sizes of lat1d and lon1d, respectively.

  • Description

  • ----------- – Interpolates RCM (Regional Climate Model), WRF (Weather Research and Forecasting) and NARR (North American Regional Reanalysis) grids to a rectilinear grid. Actually, this function will interpolate most grids that use curvilinear latitude/longitude grids. No extrapolation is performed beyond the range of the input coordinates. Missing values are allowed but ignored.

    The weighting method used is simple inverse distance squared. Missing values are allowed but ignored.

    The code searches the input curvilinear grid latitudes and longitudes for the four grid points that surround a specified output grid coordinate. Because one or more of these input points could contain missing values, fewer than four points could be used in the interpolation.

    Curvilinear grids which have two-dimensional latitude and longitude coordinate axes present some issues because the coordinates are not necessarily monotonically increasing. The simple search algorithm used by rcm2rgrid is not capable of handling all cases. The result is that, sometimes, there are small gaps in the interpolated grids. Any interior points not interpolated in the initial interpolation pass will be filled using linear interpolation.

    In some cases, edge points may not be filled.

Examples

Example 1: Using rcm2rgrid 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 rectilinear grid (or read the 1D lat and lon from
#          an other .nc file.
newlat1D_rect=np.linspace(lat2D_curv.min(), lat2D_curv.max(), 100)
newlon1D_rect=np.linspace(lon2D_curv.min(), lon2D_curv.max(), 100)

ht_rect = geocat.comp.rcm2rgrid(lat2D_curv, lon2D_curv, ht_curv, newlat1D_rect, newlon1D_rect)