geocat.comp.grid_to_triple#
- geocat.comp.grid_to_triple(data, x_in=None, y_in=None, msg_py=None)#
Converts a two-dimensional grid with one-dimensional coordinate variables to an array where each grid value is associated with its coordinates.
- Parameters:
data (
xarray.DataArray,numpy.ndarray) – Two-dimensional input array of size ny x mx containing the data values. Missing values may be present indata, but they are ignored.x_in (
xarray.DataArray,numpy.ndarray) – A one-dimensional array that specifies the the right dimension coordinates of the input (data).Note: It should only be explicitly provided when the input (
fi) isnumpy.ndarray; otherwise, it should come fromfi.coords.y_in (
xarray.DataArray,numpy.ndarray) – A one-dimensional array that specifies the the left dimension coordinates of the input (data).Note: It should only be explicitly provided when the input (
fi) isnumpy.ndarray; otherwise, it should come fromfi.coords.msg_py (
numpy.number) – A numpy scalar value that represent a missing value indata. This argument allows a user to use a missing value scheme other than NaN or masked arrays, similar to what NCL allows.
- Returns:
out (
xarray.DataArray,numpy.ndarray) – The maximum size of the returned array will be3 x ld, whereld <= ny x mx. If no missing values are encountered indata, thenld = ny x mx. If missing values are encountered indata, they are not returned and henceldwill be equal tony x mxminus the number of missing values found indata.
Examples
Example 1: Using grid_to_triple with
xarray.DataArrayinputimport 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("./NETCDF_FILE.nc") # [INPUT] Grid & data info on the source curvilinear data = ds.DIST_236_CBL[:] x_in = ds.gridlat_236[:] y_in = ds.gridlon_236[:] output = geocat.comp.grid_to_triple(data, x_in, y_in)