geocat.f2py
¶
Subpackages¶
Submodules¶
Package Contents¶
Functions¶
|
Interpolates from one series to another using piecewise linear |
|
Interpolates a regular grid to a rectilinear one using bi-linear |
|
|
|
Interpolates from a rectilinear grid to an unstructured grid or |
|
|
|
|
|
Facilitates calculating the meridional overturning circulation for the |
|
Interpolates data on a curvilinear grid (i.e. RCM, WRF, NARR) to an |
|
Interpolates data on a curvilinear grid (i.e. RCM, WRF, NARR) to a |
|
Interpolates data on a rectilinear lat/lon grid to a curvilinear grid |
|
|
|
Converts a two-dimensional grid with one-dimensional coordinate |
|
|
|
Places unstructured (randomly-spaced) data onto the nearest locations of |
- exception geocat.f2py.AttributeError¶
Bases:
Error
Exception raised when the arguments of GeoCAT-comp functions argument has a mismatch of attributes with other arguments.
- exception geocat.f2py.ChunkError¶
Bases:
Error
Exception raised when a Dask array is chunked in a way that is incompatible with an _ncomp function.
- exception geocat.f2py.CoordinateError¶
Bases:
Error
Exception raised when a GeoCAT-comp function is passed a NumPy array as an argument without a required coordinate array being passed separately.
- exception geocat.f2py.DimensionError¶
Bases:
Error
Exception raised when the arguments of GeoCAT-comp functions argument has a mismatch of the necessary dimensionality.
- exception geocat.f2py.MetaError¶
Bases:
Error
Exception raised when the support for the retention of metadata is not supported.
- geocat.f2py.linint1(fi: supported_types, xo: supported_types, xi: supported_types = None, icycx: numpy.number = 0, msg_py: numpy.number = None) supported_types ¶
Interpolates from one series to another using piecewise linear interpolation across the rightmost dimension. The series may be cyclic in the X direction.
If missing values are present, then linint1 will perform the piecewise linear interpolation at all points possible, but will return missing values at coordinates which could not be used.
If any of the output coordinates
xo
are outside those of the input coordinatesxi
, thefo
values at those coordinates will be set to missing (i.e. no extrapolation is performed).- Parameters
fi (
xarray.DataArray
,numpy.ndarray
) –An array of one or more dimensions. If
xi
is passed in as an argument, then the size of the rightmost dimension offi
must match the rightmost dimension ofxi
.If missing values are present, then
linint1
will perform the piecewise linear interpolation at all points possible, but will return missing values at coordinates which could not be used.Note: This variable must be supplied as a
xarray.DataArray
in order to copy the dimension names to the output. Otherwise, default names will be used.xo (
xarray.DataArray
,numpy.ndarray
) –A one-dimensional array that specifies the X coordinates of the return array. It must be strictly monotonically increasing or decreasing, but may be unequally spaced.
If the output coordinates
xo
are outside those of the input coordinatesxi
, then the fo values at those coordinates will be set to missing (i.e. no extrapolation is performed).xi (
xarray.DataArray
,numpy.ndarray
) –An array that specifies the X coordinates of the
fi
array. Most frequently, this array is one-dimensional. It must be strictly monotonically increasing or decreasing, but can be unequally spaced. Ifxi
is multi-dimensional, then itsimensions must be the same asfi
’s dimensions. If it is one-dimensional, its length must be the same as the rightmost (fastest varying) dimension offi
.Note: If
fi
is of typexarray.DataArray
andxi
is left unspecified, then the rightmost coordinate dimension offi
will be used. Iffi
is not of typexarray.DataArray
, thenxi
becomes a mandatory parameter. This parameter must be specified as a keyword argument.icycx (
bool
) – An option to indicate whether the rightmost dimension offi
is cyclic. 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.
- Returns
fo – The interpolated series. The returned value will have the same dimensions as
fi
, except for the rightmost dimension which will have the same dimension size as the length of xo. The return type will be double iffi
is double, and float otherwise.- Return type
Examples
Example 1: Using
linint1
withxarray.DataArray
inputimport numpy as np import xarray as xr import geocat.comp fi_np = np.random.rand(80) # random 80-element array # xi does not have to be equally spaced, but it is # in this example xi = np.arange(80) # create target coordinate array, in this case use the same # min/max values as xi, but with different spacing xo = np.linspace(xi.min(), xi.max(), 100) # create :class:`xarray.DataArray` and chunk it using the # full shape of the original array. # note that xi is attached as a coordinate array fi = xr.DataArray(fi_np, dims=['x'], coords={'x': xi} ).chunk(fi_np.shape) fo = geocat.comp.linint1(fi, xo, icycx=0)
- geocat.f2py.linint2(fi: supported_types, xo: supported_types, yo: supported_types, xi: supported_types = None, yi: supported_types = None, icycx: bool = 0, msg_py: numpy.number = None) supported_types ¶
Interpolates a regular grid to a rectilinear one using bi-linear interpolation. The input grid may be cyclic in the x-direction. The interpolation is first performed in the x-direction, and then in the y-direction.
- Parameters
fi (
xarray.DataArray
,numpy.ndarray
) –An array of two or more dimensions. If
xi
is passed in as an argument, then the size of the rightmost dimension offi
must match the rightmost dimension ofxi
. Similarly, ifyi
is passed in as an argument, then the size of the second-rightmost dimension offi
must match the rightmost dimension ofyi
.If missing values are present, then linint2 will perform the bilinear interpolation at all points possible, but will return missing values at coordinates which could not be used.
Note: This variable must be supplied as a
xarray.DataArray
in order to copy the dimension names to the output. Otherwise, default names will be used.xo (
xarray.DataArray
,numpy.ndarray
) –A one-dimensional array that specifies the X-coordinates of the return array. It must be strictly monotonically increasing, but may be unequally spaced. For geo-referenced data,
xo
is generally the longitude array.If the output coordinates
xo
are outside those of the input coordinatesxi
, then thefo
values at those coordinates will be set to missing (i.e. no extrapolation is performed).yo (
xarray.DataArray
,numpy.ndarray
) –A one-dimensional array that specifies the Y coordinates of the return array. It must be strictly monotonically increasing, but may be unequally spaced. For geo-referenced data,
yo
is typically the latitude array.If the output coordinates
yo
are outside those of the input coordinatesyi
, then thefo
values at those coordinates will be set to missing (i.e. no extrapolation is performed).xi (
xarray.DataArray
,numpy.ndarray
) –An array that specifies the X-coordinates of the
fi
array. Most frequently, this is a 1D strictly monotonically increasing array that may be unequally spaced. In some cases,xi
can be a multi-dimensional array (see next paragraph). The rightmost dimension (call itnxi
) must have at least two elements, and is the last (fastest varying) dimension offi
.If
xi
is a multi-dimensional array, then eachnxi
subsection ofxi
must be strictly monotonically increasing, but may be unequally spaced. All but its rightmost dimension must be the same size as all butfi
’s rightmost two dimensions. For geo-referenced data,xi
is generally the longitude array.Note: If
fi
is of typexarray.DataArray
andxi
is left unspecified, then the rightmost coordinate dimension offi
will be used. Iffi
is not of typexarray.DataArray
, thenxi
becomes a mandatory parameter. This parameter must be specified as a keyword argument.yi (
xarray.DataArray
,numpy.ndarray
) –An array that specifies the Y-coordinates of the
fi
array. Most frequently, this is a 1D strictly monotonically increasing array that may be unequally spaced. In some cases,yi
can be a multi-dimensional array (see next paragraph). The rightmost dimension (call it nyi) must have at least two elements, and is the second-to-last dimension offi
.If
yi
is a multi-dimensional array, then each nyi subsection ofyi
must be strictly monotonically increasing, but may be unequally spaced. All but its rightmost dimension must be the same size as all butfi
’s rightmost two dimensions. For geo-referenced data,yi
is generally the latitude array.Note: If
fi
is of typexarray.DataArray
andxi
is left unspecified, then the second-to-rightmost coordinate dimension offi
will be used. Iffi
is not of typexarray.DataArray
, thenxi
becomes a mandatory parameter. This parameter must be specified as a keyword argument.icycx (
bool
) – An option to indicate whether the rightmost dimension offi
is cyclic. 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.
- Returns
fo – The interpolated grid. If the meta parameter is True, then the result will include named dimensions matching the input array. The returned value will have the same dimensions as
fi
, except for the rightmost two dimensions which will have the same dimension sizes as the lengths ofyo
andxo
. The return type will be double iffi
is double, and float otherwise.- Return type
Examples
Example 1: Using
linint2
withxarray.DataArray
inputimport 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.linint2(fi, xo, yo, icycx=0)
- geocat.f2py.linint2_points(fi, xo, yo, icycx, msg=None, meta=False, xi=None, yi=None)¶
- geocat.f2py.linint2pts(fi: supported_types, xo: supported_types, yo: supported_types, icycx: bool = False, msg_py: numpy.number = None, xi: supported_types = None, yi: supported_types = None) supported_types ¶
Interpolates from a rectilinear grid to an unstructured grid or locations using bilinear interpolation.
The
linint2pts
function uses bilinear interpolation to interpolate from a rectilinear grid to an unstructured grid.If missing values are present, then
linint2pts
will 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.linint2pts
is different fromlinint2
in thatxo
andyo
are 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
xi
contains longitudes, then thexo
values must be in the same range. In addition, if thexi`
values span 0 to 360, then thexo
values 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 (nyi
xnxi
) are the dimensions to be used in the interpolation. If user-defined missing values are present (other than NaNs), the value ofmsg_py
must 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 offi
is 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 thefi
array.xi
might be defined as the coordinates offi
whenfi
is of typexarray.DataArray
; in this casexi
may 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 thefi
array.yi
might be defined as the coordinates offi
whenfi
is of typexarray.DataArray
; in this caseyi
may not be explicitly given as a function argument.
- Returns
fo – The returned value will have the same dimensions as
fi
, except for the rightmost dimension which will have the same dimension size as the length ofyo
andxo
.- Return type
Examples
Example 1: Using
linint2pts
withxarray.DataArray
inputimport 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)
- geocat.f2py.fort2py_msg(ndarray, msg_fort=None, msg_py=None)¶
- geocat.f2py.py2fort_msg(ndarray, msg_py=None, msg_fort=None)¶
- geocat.f2py.moc_globe_atl(lat_aux_grid: supported_types, a_wvel: supported_types, a_bolus: supported_types, a_submeso: supported_types, t_lat: supported_types, rmlak: supported_types, msg: numpy.number = None, meta: bool = False) supported_types ¶
Facilitates calculating the meridional overturning circulation for the globe and Atlantic.
- Parameters
lat_aux_grid (
xarray.DataArray
,numpy.ndarray
) – Latitude grid for transport diagnostics.a_wvel (
xarray.DataArray
,numpy.ndarray
) – Area weighted Eulerian-mean vertical velocity [TAREA x WVEL
].a_bolus (
xarray.DataArray
,numpy.ndarray
) – Area weighted Eddy-induced (bolus) vertical velocity [TAREA x WISOP
].a_submeso (
xarray.DataArray
,numpy.ndarray
) – Area weighted submeso vertical velocity [TAREA x WSUBM
].tlat (
xarray.DataArray
,numpy.ndarray
) – Array of t-grid latitudes.rmlak (
xarray.DataArray
,numpy.ndarray
) – Basin index number: [0]=Globe, [1]=Atlanticmsg (
numpy.number
) – A numpy scalar value that represent a missing value. 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 – A multi-dimensional array of size [
moc_comp
] x [n_transport_reg
] x [kdepth
] x [nyaux
] where:moc_comp
refers to the three components returnedn_transport_reg
refers to the Globe and Atlantickdepth
is the the number of vertical levels of the work arraysnyaux
is the size of thelat_aux_grid
- Return type
Examples
# Input data can be read from a data set as follows: import xarray as xr ds = xr.open_dataset("input.nc") lat_aux_grid = ds.lat_aux_grid a_wvel = ds.a_wvel a_bolus = ds.a_bolus a_submeso = ds.a_submeso tlat = ds.tlat rmlak = ds.rmlak # (1) Calling with xArray inputs and default arguments (Missing value = np.nan, NO meta information) out_arr = moc_globe_atl(lat_aux_grid, a_wvel, a_bolus, a_submeso, tlat, rmlak) # (2) Calling with Numpy inputs and default arguments (Missing value = np.nan, NO meta information) out_arr = moc_globe_atl(lat_aux_grid.values, a_wvel.values, a_bolus.values, a_submeso.values, tlat.values, rmlak.values) # (3) Calling with xArray inputs and user-defined arguments (Missing value = np.nan, NO meta information) out_arr = moc_globe_atl(lat_aux_grid, a_wvel, a_bolus, a_submeso, tlat, rmlak, msg=-99.0, meta=True) # (4) Calling with Numpy inputs and user-defined arguments (Missing value = np.nan, NO meta information) out_arr = moc_globe_atl(lat_aux_grid.values, a_wvel.values, a_bolus.values, a_submeso.values, tlat.values, rmlak.values, msg=-99.0, meta=True)
- geocat.f2py.rcm2points(lat2d: supported_types, lon2d: supported_types, fi: supported_types, lat1d: supported_types, lon1d: supported_types, opt: numpy.number = 0, msg: numpy.number = None, meta: bool = False) supported_types ¶
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 – 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 (lat1d
,lon1d
).- Return type
- geocat.f2py.rcm2rgrid(lat2d: supported_types, lon2d: supported_types, fi: supported_types, lat1d: supported_types, lon1d: supported_types, msg: numpy.number = None, meta: bool = False) supported_types ¶
Interpolates data on a curvilinear grid (i.e. RCM, WRF, NARR) to a rectilinear grid.
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.
- 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 offi
; 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 offi
; 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 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 – 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 oflat1d
andlon1d
, respectively.- Return type
Examples
Example 1: Using rcm2rgrid with
xarray.DataArray
inputimport 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)
- geocat.f2py.rgrid2rcm(lat1d: supported_types, lon1d: supported_types, fi: supported_types, lat2d: supported_types, lon2d: supported_types, msg: numpy.number = None, meta: bool = False) supported_types ¶
Interpolates data on a rectilinear lat/lon grid to a curvilinear grid like those used by the RCM, WRF and NARR models/datasets.
Interpolates data on a rectilinear lat/lon grid to 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. No extrapolation is performed beyond the range of the input coordinates. The method used is simple inverse distance weighting. Missing values are allowed but ignored.
- Parameters
lat1d (
xarray.DataArray
,numpy.ndarray
) –A one-dimensional array that specifies the latitude coordinates of the regular grid. Must be monotonically increasing.
Note: It should only be explicitly provided when the input (
fi
) isnumpy.ndarray
; otherwise, it should come fromfi.coords
.lon1d (
xarray.DataArray
,numpy.ndarray
) –A one-dimensional array that specifies the longitude coordinates of the regular grid. Must be monotonically increasing.
Note: It should only be explicitly provided when the input (
fi
) isnumpy.ndarray
; otherwise, it should come fromfi.coords
.fi (
xarray.DataArray
,numpy.ndarray
) – A multi-dimensional array to be interpolated. The rightmost two dimensions (latitude, longitude) are the dimensions to be interpolated.lat2d (
xarray.DataArray
,numpy.ndarray
) – A two-dimensional array that specifies the latitude locations of the input (fi
). Because this array is two-dimensional it is not an associated coordinate variable offi
.lon2d (
xarray.DataArray
,numpy.ndarray
) – A two-dimensional array that specifies the longitude locations of the input (fi
). Because this array is two-dimensional it is not an associated coordinate variable offi
.msg (
numpy.number
) – A numpy scalar value that represents 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 – 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 oflat2d
(orlon2d
).- Return type
Examples
Example 1: Using rgrid2rcm with
xarray.DataArray
inputimport numpy as np import xarray as xr import geocat.comp # Open a netCDF data file using xarray default engine and load the data stream # input grid and data ds_rect = xr.open_dataset("./DATAFILE_RECT.nc") # [INPUT] Grid & data info on the source rectilinear ht_rect =ds_rect.SOME_FIELD[:] lat1D_rect=ds_rect.gridlat_[:] lon1D_rect=ds_rect.gridlon_[:] # Open a netCDF data file using xarray default engine and load the data stream # for output grid ds_curv = xr.open_dataset("./DATAFILE_CURV.nc") # [OUTPUT] Grid on destination curvilinear grid (or read the 2D lat and lon from # an other .nc file newlat2D_rect=ds_curv.gridlat2D_[:] newlon2D_rect=ds_curv.gridlat2D_[:] ht_curv = geocat.comp.rgrid2rcm(lat1D_rect, lon1D_rect, ht_rect, newlat2D_curv, newlon2D_curv)
- geocat.f2py.grid2triple(x_in, y_in, data, msg_py)¶
- geocat.f2py.grid_to_triple(data: supported_types, x_in: supported_types = None, y_in: supported_types = None, msg_py: supported_types = None) supported_types ¶
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 – The maximum size of the returned array will be
3 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 henceld
will be equal tony x mx
minus the number of missing values found indata
.- Return type
Examples
Example 1: Using grid_to_triple with
xarray.DataArray
inputimport 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)
- geocat.f2py.triple2grid(x_in, y_in, data, x_out, y_out, **kwargs)¶
- geocat.f2py.triple_to_grid(data: supported_types, x_in: supported_types, y_in: supported_types, x_out: supported_types, y_out: supported_types, method: int = 1, domain: float = 1.0, distmx: float = None, missing_value: numpy.number = None, meta: bool = False) supported_types ¶
Places unstructured (randomly-spaced) data onto the nearest locations of a rectilinear grid.
This function puts unstructured data (randomly-spaced) onto the nearest locations of a rectilinear grid. A default value of
domain
option is now set to 1.0 instead of 0.0.This function does not perform interpolation; rather, each individual data point is assigned to the nearest grid point. It is possible that upon return, grid will contain grid points set to missing value if no
x_in(n)
,y_in(n)
are nearby.- Parameters
data (
xarray.DataArray
,numpy.ndarray
) – A multi-dimensional array, whose rightmost dimension is the same length asx_in
andy_in
, containing the values associated with the “x” and “y” coordinates. Missing values may be present but will be ignored.x_in (
xarray.DataArray
,numpy.ndarray
) – A one-dimensional array that specifies the x-coordinate associated with the input (data
).y_in (
xarray.DataArray
,numpy.ndarray
) – A one-dimensional array that specifies the y-coordinate associated with the input (data
).x_out (
xarray.DataArray
,numpy.ndarray
) – A one-dimensional array of lengthM
containing the x-coordinates associated with the returned two-dimensional grid. The coordinate values must be monotonically increasing.y_out (
xarray.DataArray
ornumpy.ndarray
) – A one-dimensional array of lengthN
containing the y-coordinates associated with the returned two-dimensional grid. The coordinate values must be monotonically increasing.method (
int
, optional) – An integer value that can be 0 or 1. The default value is 1. A value of 1 means to use the great circle distance formula for distance calculations. Warning:method = 0
, together withdomain = 1.0
, could result in many of the target grid points to be set to the missing value if the number of grid points is large (ie: a high resolution grid) and the number of observations relatively small.domain (
float
, optional) – A float value that should be set to a value>= 0
. The default value is 1.0. If present, the larger this factor, the wider the spatial domain allowed to influence grid boundary points. Typically,domain
is 1.0 or 2.0. Ifdomain <= 0.0
, then values located outside the grid domain specified byx_out
andy_out
arguments will not be used.distmx (
float
, optional) – Settingdistmx
allows the user to specify a search radius (km) beyond which observations are not considered for nearest neighbor. Only applicable whenmethod = 1
. The defaultdistmx=1e20 (km)
means that every grid point will have a nearest neighbor. It is suggested that users specify a reasonable value fordistmx
.missing_value (
numpy.number
, optional) – A numpy scalar value that represent a missing value indata
. The default value isnp.nan
. If specified explicitly, this argument allows the user to use a missing value scheme other than NaN or masked arrays.meta (
bool
, optional) – 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 yet supported for this function.
- Returns
grid – The returned array will be
K
xN
xM
, whereK
represents the leftmost dimensions ofdata
, N represent the size ofy_out
, and M represent the size ofx_out
coordinate vectors.- Return type
Examples
Example 1: Using triple_to_grid with
xarray.DataArray
inputimport 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 data = ds.DIST_236_CBL[:] x_in = ds.gridlat_236[:] y_in = ds.gridlon_236[:] x_out = ds.gridlat_236[:] y_out = ds.gridlon_236[:] # [OUTPUT] Grid on destination points grid (or read the 1D lat and lon from # another .nc file. newlat1D_points=np.linspace(lat2D_curv.min(), lat2D_curv.max(), 100) newlon1D_points=np.linspace(lon2D_curv.min(), lon2D_curv.max(), 100) output = geocat.comp.triple_to_grid(data, x_in, y_in, x_out, y_out)