xroms.interp#
Interpolation functions.
Functions
|
Interpolate var to lons/lats positions. |
|
Interpolate var to iso_values. |
- xroms.interp.interpll(var, lons, lats, which='pairs', **kwargs)[source]#
Interpolate var to lons/lats positions.
Wraps xESMF to perform proper horizontal interpolation on non-flat Earth.
- Parameters:
var (DataArray) – Variable to operate on.
lons (list, ndarray) – Longitudes to interpolate to. Will be flattened upon input.
lats (list, ndarray) – Latitudes to interpolate to. Will be flattened upon input.
which (str, optional) –
Which type of interpolation to do:
”pairs”: lons/lats as unstructured coordinate pairs (in xESMF language, LocStream).
”grid”: 2D array of points with 1 dimension the lons and the other dimension the lats.
**kwargs – passed on to xESMF Regridder class
- Returns:
DataArray of var interpolated to lons/lats. Dimensionality will be the
same as var except the Y and X dimensions will be 1 dimension called
”locations” that lons.size if which==’pairs’, or 2 dimensions called
”lat” and “lon” if which==’grid’ that are of lats.size and lons.size,
respectively.
Notes
var cannot have chunks in the Y or X dimensions.
cf-xarray should still be usable after calling this function.
Examples
To return 1D pairs of points, in this case 3 points:
>>> xroms.interpll(var, [-96, -97, -96.5], [26.5, 27, 26.5], which='pairs')
To return 2D pairs of points, in this case a 3x3 array of points:
>>> xroms.interpll(var, [-96, -97, -96.5], [26.5, 27, 26.5], which='grid')
- xroms.interp.isoslice(var, iso_values, xgrid, iso_array=None, axis='Z')[source]#
Interpolate var to iso_values.
This wraps xgcm transform function for slice interpolation, though transform has additional functionality.
- Parameters:
var (DataArray) – Variable to operate on.
iso_values (list, ndarray) – Values to interpolate to. If calculating var at fixed depths, iso_values are the fixed depths, which should be negative if below mean sea level. If input as array, should be 1D.
xgrid (xgcm.grid, optional) – Grid object associated with var.
iso_array (DataArray, optional) – Array that var is interpolated onto (e.g., z coordinates or density). If calculating var on fixed depth slices, iso_array contains the depths [m] associated with var. In that case and if None, will use z coordinate attached to var. Also use this option if you want to interpolate with z depths constant in time and input the appropriate z coordinate.
dim (str, optional) – Dimension over which to calculate isoslice. If calculating var onto fixed depths, dim=’Z’. Options are ‘Z’, ‘Y’, and ‘X’.
- Returns:
DataArray of var interpolated to iso_values. Dimensionality will be the
same as var except with dim dimension of size of iso_values.
Notes
var cannot have chunks in the dimension dim.
cf-xarray should still be usable after calling this function.
Examples
To calculate temperature onto fixed depths:
>>> xroms.isoslice(ds.temp, np.linspace(0, -30, 50))
To calculate temperature onto salinity:
>>> xroms.isoslice(ds.temp, np.arange(0, 36), iso_array=ds.salt, axis='Z')
Calculate lat-z slice of salinity along a constant longitude value (-91.5):
>>> xroms.isoslice(ds.salt, -91.5, iso_array=ds.lon_rho, axis='X')
Calculate slice of salt at 28 deg latitude
>>> xroms.isoslice(ds.salt, 28, iso_array=ds.lat_rho, axis='Y')
Interpolate temp to salinity values between 0 and 36 in the X direction
>>> xroms.isoslice(ds.temp, np.linspace(0, 36, 50), iso_array=ds.salt, axis='X')
Interpolate temp to salinity values between 0 and 36 in the Z direction
>>> xroms.isoslice(ds.temp, np.linspace(0, 36, 50), iso_array=ds.salt, axis='Z')
Calculate the depth of a specific isohaline (33):
>>> xroms.isoslice(ds.salt, 33, iso_array=ds.z_rho, axis='Z')
Calculate dye 10 meters above seabed. Either do this on the vertical rho grid, or first change to the w grid and then use isoslice. You may prefer to do the latter if there is a possibility that the distance above the seabed you are interpolating to (10 m) could be below the deepest rho grid depth.
on rho grid directly:
>>> height_from_seabed = ds.z_rho + ds.h >>> height_from_seabed.name = 'z_rho' >>> xroms.isoslice(ds.dye_01, 10, iso_array=height_from_seabed, axis='Z')
on w grid:
>>> var_w = ds.dye_01.xroms.to_grid(scoord='w').chunk({'s_w': -1}) >>> ds['dye_01_w'] = var_w # currently this is the easiest way to reattached coords xgcm variables >>> height_from_seabed = ds.z_w + ds.h >>> height_from_seabed.name = 'z_w' >>> xroms.isoslice(ds['dye_01_w'], 10, iso_array=height_from_seabed, axis='Z')