xroms.utilities#

Functions that act on DataArrays or Datasets.

Functions

argsel2d(lons, lats, lon0, lat0)

Find the indices of coordinate pair closest to another point.

ddeta(var, xgrid[, z, hcoord, scoord, ...])

Calculate d/deta for a variable.

ddxi(var, xgrid[, z, hcoord, scoord, ...])

Calculate d/dxi for a variable.

ddz(var, xgrid[, hcoord, scoord, hboundary, ...])

Calculate d/dz for a variable.

grid_interp(xgrid, da, dim[, ...])

Interpolate da in dim

gridmean(var, xgrid, dim)

Calculate mean accounting for variable spatial grid.

gridsum(var, xgrid, dim)

Calculate sum accounting for variable spatial grid.

hgrad(q, xgrid[, which, z, hcoord, scoord, ...])

Return gradients of property q accounting for s coordinates.

order(var)

Reorder var to typical dimensional ordering.

sel2d(var, lons, lats, lon0, lat0)

Find the value of the var at closest location to lon0,lat0.

subset(ds[, X, Y])

Subset model output horizontally using isel, properly accounting for horizontal grids.

to_grid(var, xgrid[, hcoord, scoord, ...])

Implement grid changes.

to_psi(var, xgrid[, hboundary, hfill_value])

Change var to psi horizontal grid.

to_rho(var, xgrid[, hboundary, hfill_value])

Change var to rho horizontal grid.

to_s_rho(var, xgrid[, sboundary, sfill_value])

Change var to rho vertical grid.

to_s_w(var, xgrid[, sboundary, sfill_value])

Change var to w vertical grid.

to_u(var, xgrid[, hboundary, hfill_value])

Change var to u horizontal grid.

to_v(var, xgrid[, hboundary, hfill_value])

Change var to v horizontal grid.

xisoslice(iso_array, iso_value, ...)

Calculate an isosurface.

xroms.utilities.argsel2d(lons, lats, lon0, lat0)[source]#

Find the indices of coordinate pair closest to another point.

Parameters:
  • lons (DataArray, ndarray, list) – Longitudes of points to search through for closest point.

  • lats (DataArray, ndarray, list) – Latitudes of points to search through for closest point.

  • lon0 (float, int) – Longitude of comparison point.

  • lat0 (float, int) – Latitude of comparison point.

Returns:

  • Index or indices of location in coordinate pairs made up of lons, lats

  • that is closest to location lon0, lat0. Number of dimensions of

  • returned indices will correspond to the shape of input lons.

Notes

This function uses Great Circle distance to calculate distances assuming longitudes and latitudes as point coordinates. Uses cartopy function Geodesic: https://scitools.org.uk/cartopy/docs/latest/cartopy/geodesic.html

If searching for the closest grid node to a lon/lat location, be sure to use the correct horizontal grid (rho, u, v, or psi). This is accounted for if this function is used through the accessor.

Examples

>>> xroms.argsel2d(ds.lon_rho, ds.lat_rho, -96, 27)
xroms.utilities.ddeta(var, xgrid, z=None, hcoord=None, scoord=None, hboundary='extend', hfill_value=nan, sboundary='extend', sfill_value=nan, attrs=None)[source]#

Calculate d/deta for a variable.

Note that you need the 3D metrics for horizontal derivatives for ROMS, so include_3D_metrics=True in xroms.roms_dataset().

Parameters:
  • var (DataArray, ndarray) – Variable to operate on.

  • xgrid (xgcm.grid) – Grid object associated with var.

  • z (DataArray, ndarray, optional) – Depth [m]. If None, use z coordinate attached to var.

  • hcoord (string, optional) – Name of horizontal grid to interpolate output to. Options are ‘rho’, ‘psi’, ‘u’, ‘v’.

  • scoord (string, optional) – Name of vertical grid to interpolate output to. Options are ‘s_rho’, ‘s_w’, ‘rho’, ‘w’.

  • hboundary (string, optional) –

    Passed to grid method calls; horizontal boundary selection for calculating horizontal derivative of var. This same value will be used for grid changes too. From xgcm documentation: A flag indicating how to handle boundaries:

    • None: Do not apply any boundary conditions. Raise an error if boundary conditions are required for the operation.

    • ’fill’: Set values outside the array boundary to fill_value (i.e. a Neumann boundary condition.)

    • ’extend’: Set values outside the array to the nearest array value. (i.e. a limited form of Dirichlet boundary condition.

  • hfill_value (float, optional) – Passed to grid method calls; horizontal boundary selection fill value. From xgcm documentation: The value to use in the boundary condition with boundary=’fill’.

  • sboundary (string, optional) –

    Passed to grid method calls; vertical boundary selection for calculating horizontal derivative of var. This same value will be used for vertical grid changes too. From xgcm documentation: A flag indicating how to handle boundaries:

    • None: Do not apply any boundary conditions. Raise an error if boundary conditions are required for the operation.

    • ’fill’: Set values outside the array boundary to fill_value (i.e. a Neumann boundary condition.)

    • ’extend’: Set values outside the array to the nearest array value. (i.e. a limited form of Dirichlet boundary condition.

  • sfill_value (float, optional) – Passed to grid method calls; vertical boundary selection fill value. From xgcm documentation: The value to use in the boundary condition with boundary=’fill’.

  • attrs (dict, optional) – Dictionary of attributes to add to resultant arrays. Requires that q is DataArray. For example: attrs={‘name’: ‘varname’, ‘long_name’: ‘longvarname’, ‘units’: ‘units’}

Returns:

  • DataArray or ndarray of dqdeta, the gradient of q in the eta-direction with

  • attributes altered to reflect calculation.

Notes

dqdeta = dqdy*dzdz - dqdz*dzdy

Derivatives are taken in the ROMS curvilinear grid native eta-direction.

These derivatives properly account for the fact that ROMS vertical coordinates are s coordinates and therefore can vary in time and space.

This will alter the number of points in the eta and s dimensions.

Examples

>>> xroms.ddeta(ds.salt, xgrid)
xroms.utilities.ddxi(var, xgrid, z=None, hcoord=None, scoord=None, hboundary='extend', hfill_value=nan, sboundary='extend', sfill_value=nan, attrs=None)[source]#

Calculate d/dxi for a variable.

Note that you need the 3D metrics for horizontal derivatives for ROMS, so include_3D_metrics=True in xroms.roms_dataset().

Parameters:
  • var (DataArray) – Variable to operate on.

  • xgrid (xgcm.grid) – Grid object associated with var.

  • z (DataArray, ndarray, optional) – Depth [m]. If None, use z coordinate attached to var.

  • hcoord (string, optional) – Name of horizontal grid to interpolate output to. Options are ‘rho’, ‘psi’, ‘u’, ‘v’.

  • scoord (string, optional) – Name of vertical grid to interpolate output to. Options are ‘s_rho’, ‘s_w’, ‘rho’, ‘w’.

  • hboundary (string, optional) –

    Passed to grid method calls; horizontal boundary selection for calculating horizontal derivative of var. This same value will be used for all horizontal grid changes too. From xgcm documentation: A flag indicating how to handle boundaries:

    • None: Do not apply any boundary conditions. Raise an error if boundary conditions are required for the operation.

    • ’fill’: Set values outside the array boundary to fill_value (i.e. a Neumann boundary condition.)

    • ’extend’: Set values outside the array to the nearest array value. (i.e. a limited form of Dirichlet boundary condition.

  • hfill_value (float, optional) – Passed to grid method calls; horizontal boundary selection fill value. From xgcm documentation: The value to use in the boundary condition with boundary=’fill’.

  • sboundary (string, optional) –

    Passed to grid method calls; vertical boundary selection for calculating horizontal derivative of var. This same value will be used for all vertical grid changes too. From xgcm documentation: A flag indicating how to handle boundaries:

    • None: Do not apply any boundary conditions. Raise an error if boundary conditions are required for the operation.

    • ’fill’: Set values outside the array boundary to fill_value (i.e. a Neumann boundary condition.)

    • ’extend’: Set values outside the array to the nearest array value. (i.e. a limited form of Dirichlet boundary condition.

  • sfill_value (float, optional) – Passed to grid method calls; vertical boundary selection fill value. From xgcm documentation: The value to use in the boundary condition with boundary=’fill’.

  • attrs (dict, optional) – Dictionary of attributes to add to resultant arrays. Requires that q is DataArray. For example: attrs={‘name’: ‘varname’, ‘long_name’: ‘longvarname’, ‘units’: ‘units’}

Returns:

  • DataArray of dqdxi, the gradient of q in the xi-direction with

  • attributes altered to reflect calculation.

Notes

dqdxi = dqdx*dzdz - dqdz*dzdx

Derivatives are taken in the ROMS curvilinear grid native xi-direction.

These derivatives properly account for the fact that ROMS vertical coordinates are s coordinates and therefore can vary in time and space.

This will alter the number of points in the xi and s dimensions.

Examples

>>> xroms.ddxi(ds.salt, xgrid)
xroms.utilities.ddz(var, xgrid, hcoord=None, scoord=None, hboundary='extend', hfill_value=None, sboundary='extend', sfill_value=None, attrs=None)[source]#

Calculate d/dz for a variable.

Parameters:
  • var (DataArray) – Variable to operate on.

  • xgrid (xgcm.grid) – Grid object associated with var

  • hcoord (string, optional.) – Name of horizontal grid to interpolate output to. Options are ‘rho’, ‘psi’, ‘u’, ‘v’.

  • scoord (string, optional.) – Name of vertical grid to interpolate output to. Options are ‘s_rho’, ‘s_w’, ‘rho’, ‘w’.

  • hboundary (string, optional) –

    Passed to grid method calls; horizontal boundary selection for calculating horizontal derivative of var. This same value will be used for grid changes too. From xgcm documentation: A flag indicating how to handle boundaries:

    • None: Do not apply any boundary conditions. Raise an error if boundary conditions are required for the operation.

    • ’fill’: Set values outside the array boundary to fill_value (i.e. a Neumann boundary condition.)

    • ’extend’: Set values outside the array to the nearest array value. (i.e. a limited form of Dirichlet boundary condition.

  • hfill_value (float, optional) – Passed to grid method calls; horizontal boundary selection fill value. From xgcm documentation: The value to use in the boundary condition with boundary=’fill’.

  • sboundary (string, optional) –

    Passed to grid method calls; vertical boundary selection for calculating z derivative. This same value will be used for grid changes too. From xgcm documentation: A flag indicating how to handle boundaries:

    • None: Do not apply any boundary conditions. Raise an error if boundary conditions are required for the operation.

    • ’fill’: Set values outside the array boundary to fill_value (i.e. a Neumann boundary condition.)

    • ’extend’: Set values outside the array to the nearest array value. (i.e. a limited form of Dirichlet boundary condition.

  • sfill_value (float, optional) – Passed to grid method calls; vertical boundary fill value associated with sboundary input. From xgcm documentation: The value to use in the boundary condition with boundary=’fill’.

  • attrs (dict, optional) – Dictionary of attributes to add to resultant arrays. Requires that q is DataArray. For example: attrs={‘name’: ‘varname’, ‘long_name’: ‘longvarname’, ‘units’: ‘units’}

Returns:

  • DataArray of vertical derivative of variable with

  • attributes altered to reflect calculation.

Notes

This will alter the number of points in the s dimension.

Examples

>>> xroms.ddz(ds.salt, xgrid)
xroms.utilities.grid_interp(xgrid, da, dim, which_xgcm_function='interp', **kwargs)[source]#

Interpolate da in dim

This function is necessary because of weirdness with chunking with xgcm. More info: https://github.com/xgcm/xgcm/issues/522

Parameters:
  • xgrid (xgcm grid object) – _description_

  • da (DataArray) – interpolating from this dataarray

  • dim (str) – interpolating grids in this dimension

  • which_xgcm_function ("interp") – But could instead be “integrate”

Returns:

interpolated down one dimension in dim

Return type:

DataArray

xroms.utilities.gridmean(var, xgrid, dim)[source]#

Calculate mean accounting for variable spatial grid.

Parameters:
  • var (DataArray or ndarray) – Variable to operate on.

  • xgrid (xgcm.grid) – Grid object associated with var

  • dim (str, list, tuple) – Spatial dimension names to average over. In the xgcm convention, the allowable names are ‘Z’, ‘Y’, or ‘X’.

Returns:

  • DataArray or ndarray of average calculated over dim accounting

  • for variable spatial grid.

Notes

If result is DataArray, long name attribute is modified to describe calculation.

Examples

Note that the following two approaches are equivalent:

>>> app1 = xroms.gridmean(ds.u, xgrid, ('Y','X'))
>>> app2 = (ds.u*ds.dy_u*ds.dx_u).sum(('eta_rho','xi_u'))/(ds.dy_u*ds.dx_u).sum(('eta_rho','xi_u'))
>>> np.allclose(app1, app2)
xroms.utilities.gridsum(var, xgrid, dim)[source]#

Calculate sum accounting for variable spatial grid.

Parameters:
  • var (DataArray or ndarray) – Variable to operate on.

  • xgrid (xgcm.grid) – Grid object associated with var

  • dim (str, list, tuple) – Spatial dimension names to sum over. In the xgcm convention, the allowable names are ‘Z’, ‘Y’, or ‘X’.

Returns:

  • DataArray or ndarray of sum calculated over dim accounting

  • for variable spatial grid.

Notes

If result is DataArray, long name attribute is modified to describe calculation.

Examples

Note that the following two approaches are equivalent:

>>> app1 = xroms.gridsum(ds.u, xgrid, ('Z','X'))
>>> app2 = (ds.u*ds.dz_u * ds.dx_u).sum(('s_rho','xi_u'))
>>> np.allclose(app1, app2)
xroms.utilities.hgrad(q, xgrid, which='both', z=None, hcoord=None, scoord=None, hboundary='extend', hfill_value=None, sboundary='extend', sfill_value=None, attrs=None)[source]#

Return gradients of property q accounting for s coordinates.

Note that you need the 3D metrics for horizontal derivatives for ROMS, so include_3D_metrics=True in xroms.roms_dataset().

Parameters:
  • q (DataArray) – Property to take gradients of.

  • xgrid (xgcm.grid) – Grid object associated with q.

  • which (string, optional) –

    Which components of gradient to return.

    • ’both’: return both components of hgrad.

    • ’xi’: return only xi-direction.

    • ’eta’: return only eta-direction.

  • z (DataArray, ndarray, optional) – Depth [m]. If None, use z coordinate attached to q.

  • hcoord (string, optional) – Name of horizontal grid to interpolate output to. Options are ‘rho’, ‘psi’, ‘u’, ‘v’.

  • scoord (string, optional) – Name of vertical grid to interpolate output to. Options are ‘s_rho’, ‘s_w’, ‘rho’, ‘w’.

  • hboundary (string, optional) –

    Passed to grid method calls; horizontal boundary selection for calculating horizontal derivatives of q. This same value will be used for all horizontal grid changes too. From xgcm documentation: A flag indicating how to handle boundaries:

    • None: Do not apply any boundary conditions. Raise an error if boundary conditions are required for the operation.

    • ’fill’: Set values outside the array boundary to fill_value (i.e. a Neumann boundary condition.)

    • ’extend’: Set values outside the array to the nearest array value. (i.e. a limited form of Dirichlet boundary condition.

  • hfill_value (float, optional) – Passed to grid method calls; horizontal boundary selection fill value. From xgcm documentation: The value to use in the boundary condition with boundary=’fill’.

  • sboundary (string, optional) –

    Passed to grid method calls; vertical boundary selection for calculating horizontal derivatives of q. This same value will be used for all vertical grid changes too. From xgcm documentation: A flag indicating how to handle boundaries:

    • None: Do not apply any boundary conditions. Raise an error if boundary conditions are required for the operation.

    • ’fill’: Set values outside the array boundary to fill_value (i.e. a Neumann boundary condition.)

    • ’extend’: Set values outside the array to the nearest array value. (i.e. a limited form of Dirichlet boundary condition.

  • sfill_value (float, optional) – Passed to grid method calls; vertical boundary selection fill value. From xgcm documentation: The value to use in the boundary condition with boundary=’fill’.

  • attrs (dict, optional) – Dictionary of attributes to add to resultant arrays. Requires that q is DataArray.

Returns:

  • DataArray(s) of dqdxi and/or dqdeta, the gradients of q

  • in the xi- and eta-directions with attributes altered to reflect calculation.

Notes

dqdxi = dqdx*dzdz - dqdz*dzdx

dqdeta = dqdy*dzdz - dqdz*dzdy

Derivatives are taken in the ROMS curvilinear grid native xi- and eta- directions.

These derivatives properly account for the fact that ROMS vertical coordinates are s coordinates and therefore can vary in time and space.

The xi derivative will alter the number of points in the xi and s dimensions. The eta derivative will alter the number of points in the eta and s dimensions.

Examples

>>> dtempdxi, dtempdeta = xroms.hgrad(ds.temp, xgrid)
xroms.utilities.order(var)[source]#

Reorder var to typical dimensional ordering.

Parameters:

var (DataArray) – Variable to operate on.

Returns:

  • DataArray with dimensional order [‘T’, ‘Z’, ‘Y’, ‘X’], or whatever subset of

  • dimensions are present in var.

Notes

Do not consider previously-selected dimensions that are kept on as coordinates but cannot be transposed anymore. This is accomplished with .reset_coords(drop=True).

Examples

>>> xroms.order(var)
xroms.utilities.sel2d(var, lons, lats, lon0, lat0)[source]#

Find the value of the var at closest location to lon0,lat0.

Parameters:
  • var (DataArray, ndarray) – Variable to operate on.

  • lons (DataArray, ndarray, list) – Longitudes of points to search through for closest point.

  • lats (DataArray, ndarray, list) – Latitudes of points to search through for closest point.

  • lon0 (float, int) – Longitude of comparison point.

  • lat0 (float, int) – Latitude of comparison point.

Returns:

  • Value in var of location in coordinate pairs made up of lons, lats

  • that is closest to location lon0, lat0. If var has other

  • dimensions, they are brought along.

Notes

This function uses Great Circle distance to calculate distances assuming longitudes and latitudes as point coordinates. Uses cartopy function Geodesic: https://scitools.org.uk/cartopy/docs/latest/cartopy/geodesic.html

If searching for the closest grid node to a lon/lat location, be sure to use the correct horizontal grid (rho, u, v, or psi). This is accounted for if this function is used through the accessor.

This is meant to be used by the accessor to conveniently wrap argsel2d.

Examples

>>> xroms.sel2d(ds.temp, ds.lon_rho, ds.lat_rho, -96, 27)
xroms.utilities.subset(ds, X=None, Y=None)[source]#

Subset model output horizontally using isel, properly accounting for horizontal grids.

Parameters:
  • ds (xarray Dataset) – Dataset of ROMS model output. Assumes that full regular grid setup is available and has been read in using xroms so that dimension names have been updated.

  • X (slice, optional) –

    Slice in X dimension using form X=slice(start, stop, step). For example,

    >>> X=slice(20,40,2)
    

    Indices are used for rho grid, and psi grid is reduced accordingly.

  • Y (slice, optional) –

    Slice in Y dimension using form Y=slice(start, stop, step). For example,

    >>> Y=slice(20,40,2)
    

    Indices are used for rho grid, and psi grid is reduced accordingly.

Returns:

  • Dataset with form as if model had been run at the subsetted size. That is, the outermost

  • cells of the rho grid are like ghost cells and the psi grid is one inward from this size

  • in each direction.

Notes

X and Y must be slices, not single numbers.

Examples

Subset only in Y direction:

>>> xroms.subset(ds, Y=slice(50,100))

Subset in X and Y:

>>> xroms.subset(ds, X=slice(20,40), Y=slice(50,100))
xroms.utilities.to_grid(var, xgrid, hcoord=None, scoord=None, hboundary='extend', hfill_value=None, sboundary='extend', sfill_value=None, attrs=None)[source]#

Implement grid changes.

Parameters:
  • var (DataArray or ndarray) – Variable to operate on.

  • xgrid (xgcm.grid) – Grid object associated with var

  • hcoord (string, optional.) – Name of horizontal grid to interpolate output to. Options are ‘rho’, ‘psi’, ‘u’, ‘v’.

  • scoord (string, optional.) – Name of vertical grid to interpolate output to. Options are ‘s_rho’, ‘s_w’, ‘rho’, ‘w’.

  • hboundary (string, optional) –

    Passed to grid method calls; horizontal boundary selection for grid changes. From xgcm documentation:

    A flag indicating how to handle boundaries: * None: Do not apply any boundary conditions. Raise an error if

    boundary conditions are required for the operation.

    • ’fill’: Set values outside the array boundary to fill_value (i.e. a Neumann boundary condition.)

    • ’extend’: Set values outside the array to the nearest array value. (i.e. a limited form of Dirichlet boundary condition.

  • hfill_value (float, optional) – Passed to grid method calls; horizontal boundary selection fill value. From xgcm documentation: The value to use in the boundary condition with boundary=’fill’.

  • sboundary (string, optional) –

    Passed to grid method calls; vertical boundary selection for grid changes. From xgcm documentation: A flag indicating how to handle boundaries:

    • None: Do not apply any boundary conditions. Raise an error if boundary conditions are required for the operation.

    • ’fill’: Set values outside the array boundary to fill_value (i.e. a Neumann boundary condition.)

    • ’extend’: Set values outside the array to the nearest array value. (i.e. a limited form of Dirichlet boundary condition.

  • sfill_value (float, optional) – Passed to grid method calls; vertical boundary selection fill value. From xgcm documentation: The value to use in the boundary condition with boundary=’fill’.

Returns:

  • DataArray or ndarray interpolated onto hcoord horizontal and scoord

  • vertical grids.

Notes

If var is already on selected grid, nothing happens.

Examples

>>> xroms.to_grid(ds.salt, xgrid, hcoord='rho', scoord='w')
xroms.utilities.to_psi(var, xgrid, hboundary='extend', hfill_value=None)[source]#

Change var to psi horizontal grid.

Parameters:
  • var (DataArray or ndarray) – Variable to operate on.

  • xgrid (xgcm.grid) – Grid object associated with var

  • hboundary (string, optional) –

    Passed to grid method calls; horizontal boundary selection for grid changes. From xgcm documentation: A flag indicating how to handle boundaries:

    • None: Do not apply any boundary conditions. Raise an error if boundary conditions are required for the operation.

    • ’fill’: Set values outside the array boundary to fill_value (i.e. a Neumann boundary condition.)

    • ’extend’: Set values outside the array to the nearest array value. (i.e. a limited form of Dirichlet boundary condition.

  • hfill_value (float, optional) – Passed to grid method calls; horizontal boundary selection fill value. From xgcm documentation: The value to use in the boundary condition with boundary=’fill’.

Return type:

DataArray or ndarray interpolated onto psi horizontal grid.

Notes

If var is already on psi grid, nothing happens.

to_grid function wraps all of the to_* functions so one function can be used for all grid changes.

Examples

>>> xroms.to_psi('salt', xgrid)
xroms.utilities.to_rho(var, xgrid, hboundary='extend', hfill_value=None)[source]#

Change var to rho horizontal grid.

Parameters:
  • var (DataArray or ndarray) – Variable to operate on.

  • xgrid (xgcm.grid) – Grid object associated with var

  • hboundary (string, optional) –

    Passed to grid method calls; horizontal boundary selection for grid changes. From xgcm documentation: A flag indicating how to handle boundaries:

    • None: Do not apply any boundary conditions. Raise an error if boundary conditions are required for the operation.

    • ’fill’: Set values outside the array boundary to fill_value (i.e. a Neumann boundary condition.)

    • ’extend’: Set values outside the array to the nearest array value. (i.e. a limited form of Dirichlet boundary condition.

  • hfill_value (float, optional) – Passed to grid method calls; horizontal boundary selection fill value. From xgcm documentation: The value to use in the boundary condition with boundary=’fill’.

Return type:

DataArray or ndarray interpolated onto rho horizontal grid.

Notes

If var is already on rho grid, nothing happens.

to_grid function wraps all of the to_* functions so one function can be used for all grid changes.

Examples

>>> xroms.to_rho('salt', xgrid)
xroms.utilities.to_s_rho(var, xgrid, sboundary='extend', sfill_value=None)[source]#

Change var to rho vertical grid.

Parameters:
  • var (DataArray or ndarray) – Variable to operate on.

  • xgrid (xgcm.grid) – Grid object associated with var

  • sboundary (string, optional) –

    Passed to grid method calls; vertical boundary selection for grid changes. From xgcm documentation: A flag indicating how to handle boundaries:

    • None: Do not apply any boundary conditions. Raise an error if boundary conditions are required for the operation.

    • ’fill’: Set values outside the array boundary to fill_value (i.e. a Neumann boundary condition.)

    • ’extend’: Set values outside the array to the nearest array value. (i.e. a limited form of Dirichlet boundary condition.

  • sfill_value (float, optional) – Passed to grid method calls; vertical boundary selection fill value. From xgcm documentation: The value to use in the boundary condition with boundary=’fill’.

Return type:

DataArray or ndarray interpolated onto rho vertical grid.

Notes

If var is already on rho grid, nothing happens.

to_grid function wraps all of the to_* functions so one function can be used for all grid changes.

Examples

>>> xroms.to_s_rho('salt', xgrid)
xroms.utilities.to_s_w(var, xgrid, sboundary='extend', sfill_value=None)[source]#

Change var to w vertical grid.

Parameters:
  • var (DataArray or ndarray) – Variable to operate on.

  • xgrid (xgcm.grid) – Grid object associated with var

  • sboundary (string, optional) –

    Passed to grid method calls; vertical boundary selection for grid changes. From xgcm documentation: A flag indicating how to handle boundaries:

    • None: Do not apply any boundary conditions. Raise an error if boundary conditions are required for the operation.

    • ’fill’: Set values outside the array boundary to fill_value (i.e. a Neumann boundary condition.)

    • ’extend’: Set values outside the array to the nearest array value. (i.e. a limited form of Dirichlet boundary condition.

  • sfill_value (float, optional) – Passed to grid method calls; vertical boundary selection fill value. From xgcm documentation: The value to use in the boundary condition with boundary=’fill’.

Return type:

DataArray or ndarray interpolated onto rho vertical grid.

Notes

If var is already on w grid, nothing happens.

to_grid function wraps all of the to_* functions so one function can be used for all grid changes.

Examples

>>> xroms.to_s_w('salt', xgrid)
xroms.utilities.to_u(var, xgrid, hboundary='extend', hfill_value=None)[source]#

Change var to u horizontal grid.

Parameters:
  • var (DataArray or ndarray) – Variable to operate on.

  • xgrid (xgcm.grid) – Grid object associated with var

  • hboundary (string, optional) –

    Passed to grid method calls; horizontal boundary selection for grid changes. From xgcm documentation: A flag indicating how to handle boundaries:

    • None: Do not apply any boundary conditions. Raise an error if boundary conditions are required for the operation.

    • ’fill’: Set values outside the array boundary to fill_value (i.e. a Neumann boundary condition.)

    • ’extend’: Set values outside the array to the nearest array value. (i.e. a limited form of Dirichlet boundary condition.

  • hfill_value (float, optional) – Passed to grid method calls; horizontal boundary selection fill value. From xgcm documentation: The value to use in the boundary condition with boundary=’fill’.

Return type:

DataArray or ndarray interpolated onto u horizontal grid.

Notes

If var is already on u grid, nothing happens.

to_grid function wraps all of the to_* functions so one function can be used for all grid changes.

Examples

>>> xroms.to_u('salt', xgrid)
xroms.utilities.to_v(var, xgrid, hboundary='extend', hfill_value=None)[source]#

Change var to v horizontal grid.

Parameters:
  • var (DataArray or ndarray) – Variable to operate on.

  • xgrid (xgcm.grid) – Grid object associated with var

  • hboundary (string, optional) –

    Passed to grid method calls; horizontal boundary selection for grid changes. From xgcm documentation: A flag indicating how to handle boundaries:

    • None: Do not apply any boundary conditions. Raise an error if boundary conditions are required for the operation.

    • ’fill’: Set values outside the array boundary to fill_value (i.e. a Neumann boundary condition.)

    • ’extend’: Set values outside the array to the nearest array value. (i.e. a limited form of Dirichlet boundary condition.

  • hfill_value (float, optional) – Passed to grid method calls; horizontal boundary selection fill value. From xgcm documentation: The value to use in the boundary condition with boundary=’fill’.

Return type:

DataArray or ndarray interpolated onto v horizontal grid.

Notes

If var is already on v grid, nothing happens.

to_grid function wraps all of the to_* functions so one function can be used for all grid changes.

Examples

>>> xroms.to_v('salt', xgrid)
xroms.utilities.xisoslice(iso_array, iso_value, projected_array, coord)[source]#

Calculate an isosurface.

This function has been possibly superseded by isoslice that wraps xgcm.grid.transform for the following reasons, but more testing is needed:

  • The implementation of xgcm.grid.transform is more robust than xisoslice which has extra code for in case iso_value is exactly in iso_array.

  • For a 5-day model file, the run time for the same call for was approximately the same for xisolice and isoslice.

  • isoslice might be more computationally robust for not breaking mid-way, but this is still unclear.

This function calculates the value of projected_array on an isosurface in the array iso_array defined by iso_value.

Parameters:
  • iso_array (DataArray, ndarray) – Array in which the isosurface is defined

  • iso_value (float) – Value of the isosurface in iso_array

  • projected_array (DataArray, ndarray) – Array in which to project values on the isosurface. This can have multiple time outputs. Needs to be broadcastable from iso_array?

  • coord (string) – Name of coordinate associated with the dimension along which to project

Return type:

DataArray or ndarray of values of projected_array on the isosurface

Notes

Performs lazy evaluation.

xisoslice requires that iso_array be monotonic. If iso_value is not monotonic it will still run but values may be incorrect where not monotonic. If iso_value is exactly in iso_array or the value is passed twice in iso_array, a message will be printed. iso_value is changed a tiny amount in this case to account for it being in iso_array exactly. The latter case is not deal with.

Examples

Calculate lat-z slice of salinity along a constant longitude value (-91.5):

>>> sl = xroms.utilities.xisoslice(ds.lon_rho, -91.5, ds.salt, 'xi_rho')

Calculate a lon-lat slice at a constant z value (-10):

>>> sl = xroms.utilities.xisoslice(ds.z_rho, -10, ds.temp, 's_rho')

Calculate a lon-lat slice at a constant z value (-10) but without zeta changing in time:

(use ds.z_rho0 which is relative to mean sea level and does not vary in time) >>> sl = xroms.utilities.xisoslice(ds.z_rho0, -10, ds.temp, ‘s_rho’)

Calculate the depth of a specific isohaline (33):

>>> sl = xroms.utilities.xisoslice(ds.salt, 33, ds.z_rho, 's_rho')

Calculate the salt 10 meters above the seabed. Either do this on the vertical rho grid, or first change to the w grid and then use xisoslice. 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:

    >>> sl = xroms.xisoslice(ds.z_rho + ds.h, 10., ds.salt, 's_rho')
    
  • on w grid:

    >>> var_w = xroms.to_s_w(ds.salt, ds.xroms.xgrid)
    >>> sl = xroms.xisoslice(ds.z_w + ds.h, 10., var_w, 's_w')
    

In addition to calculating the slices themselves, you may need to calculate related coordinates for plotting. For example, to accompany the lat-z slice, you may want the following:

calculate z values (s_rho)

>>> slz = xroms.utilities.xisoslice(ds.lon_rho, -91.5, ds.z_rho, 'xi_rho')

calculate latitude values (eta_rho)

>>> sllat = xroms.utilities.xisoslice(ds.lon_rho, -91.5, ds.lat_rho, 'xi_rho')

assign these as coords to be used in plot

>>> sl = sl.assign_coords(z=slz, lat=sllat)

points that should be masked

>>> slmask = xroms.utilities.xisoslice(ds.lon_rho, -91.5, ds.mask_rho, 'xi_rho')

drop masked values

>>> sl = sl.where(slmask==1, drop=True)