xroms.accessor#

This is an accessor to xarray. It is basically a convenient way to use some of the xroms functions, which has bookkeeping in the background where possible. No functions are available only here; this connects to functions in other files.

Classes

xromsDataArrayAccessor(da)

Accessor for DataArrays.

xromsDatasetAccessor(ds)

Accessor for Datasets.

class xroms.accessor.xromsDataArrayAccessor(da)[source]#

Bases: object

Accessor for DataArrays.

Methods

argsel2d(lon0, lat0)

Find the indices of coordinate pair closest to another point.

ddeta(xgrid[, hcoord, scoord, hboundary, ...])

Calculate d/deta for a variable.

ddxi(xgrid[, hcoord, scoord, hboundary, ...])

Calculate d/dxi for variable.

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

Calculate d/dz for a variable.

gridmean(xgrid, dim)

Calculate mean accounting for variable spatial grid.

gridsum(xgrid, dim)

Calculate sum accounting for variable spatial grid.

interpll(lons, lats[, which])

Interpolate var to lons/lats positions.

order()

Reorder self to typical dimensional ordering.

sel2d(lon0, lat0)

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

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

Implement grid changes.

zslice(xgrid, depths[, z])

Interpolate var to depths.

argsel2d(lon0, lat0)[source]#

Find the indices of coordinate pair closest to another point.

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

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

Return type:

Indices in eta, xi of closest location to lon0, lat0.

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

Examples

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

Calculate d/deta for a variable.

Parameters:
  • xgrid – xgcm grid

  • 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 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 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

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

Calculate d/dxi for variable.

Parameters:
  • xgrid – xgcm grid

  • 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

>>> ds.salt.xroms.ddxi(xgrid)
ddz(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:
  • xgrid – xgcm grid

  • 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 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

>>> ds.salt.xroms.ddz(xgrid)
gridmean(xgrid, dim)[source]#

Calculate mean accounting for variable spatial grid.

Parameters:
  • xgrid – xgcm grid

  • 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 = ds.u.xroms.gridmean(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)

gridsum(xgrid, dim)[source]#

Calculate sum accounting for variable spatial grid.

Parameters:
  • xgrid – xgcm grid

  • 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 = ds.u.xroms.gridsum(xgrid, (‘Z’,’X’)) >>> app2 = (ds.u*ds.dz_u * ds.dx_u).sum((‘s_rho’,’xi_u’)) >>> np.allclose(app1, app2)

interpll(lons, lats, which='pairs', **kwargs)[source]#

Interpolate var to lons/lats positions.

Wraps xESMF to perform proper horizontal interpolation on non-flat Earth.

Parameters:
  • 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’)

order()[source]#

Reorder self to typical dimensional ordering.

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

>>> ds.temp.xroms.order()
sel2d(lon0, lat0)[source]#

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

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

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

Return type:

DataArray value(s) of closest location to lon0/lat0.

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

This wraps argsel2d.

Examples

>>> ds.temp.xroms.sel2d(-96, 27)
to_grid(xgrid, hcoord=None, scoord=None, hboundary='extend', hfill_value=None, sboundary='extend', sfill_value=None)[source]#

Implement grid changes.

Parameters:
  • xgrid – xgcm grid

  • 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 interpolated onto hcoord horizontal and scoord

  • vertical grids.

Notes

If var is already on selected grid, nothing happens.

Examples

>>> ds.salt.xroms.to_grid(xgrid, hcoord='rho', scoord='w')
zslice(xgrid, depths, z=None)[source]#

Interpolate var to depths.

This wraps xgcm transform function for slice interpolation, though transform has additional functionality. See xroms.isoslice for full docs.

Parameters:
  • xgrid – xgcm grid

  • depths (list, ndarray) – Values to interpolate to (called iso_values in other functions). Should be negative if below mean sea level. If input as array, should be 1D.

  • z (DataArray, optional) – Array that var is interpolated onto (e.g., z coordinates or density). The “vertical” coordinate is selected by default. Use this option if you want to interpolate with z depths constant in time and input the appropriate z coordinate (e.g. z_rho0).

Returns:

  • DataArray of var interpolated to depths. Dimensionality will be the

  • same as var except with dim dimension of size of depths.

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:

>>> ds.temp.xroms.zslice(depths)

To calculate temperature onto fixed depths without considering time for z coord:

>>> ds.temp.xroms.zslice(depths, z=ds.temp.z_rho0)
class xroms.accessor.xromsDatasetAccessor(ds)[source]#

Bases: object

Accessor for Datasets.

Attributes:
EKE

Calculate EKE [m^2/s^2], on rho grid.

KE

Calculate kinetic energy [kg/(m*s^2)], on rho/rho grids.

M2

Calculate the horizontal buoyancy gradient on rho/w grids.

N2

Calculate buoyancy frequency squared on rho/w grids.

buoyancy

Calculate buoyancy on rho/rho grids.

convergence

Calculate convergence, rho/rho grid.

convergence_norm

Calculate normalized surface convergence, rho/rho grid.

dudz

Calculate dudz [1/s] on u/w grids.

dvdz

Calculate dvdz [1/s] on v/w grids.

east

Rotate grid-aligned u velocity to be eastward.

ertel

Calculate Ertel potential vorticity of buoyancy on rho/rho grids.

north

Rotate grid-aligned v velocity to be northward.

omega

Calculate s-grid vertical velocity on [horizontal]/[vertical] grids.

rho

Return existing rho or calculate, on rho/rho grids.

sig0

Calculate potential density referenced to z=0, on rho/rho grids.

speed

Calculate horizontal speed [m/s] from u and v components, on rho/rho grids.

u

Rotate eastward velocity to be grid-aligned u velocity

ug

Calculate geostrophic u velocity from zeta, on u grid.

v

Rotate northward velocity to be grid-aligned v velocity

vertical_shear

Calculate vertical shear [1/s], rho/w grids.

vg

Calculate geostrophic v velocity from zeta, on v grid.

vort

Calculate vertical relative vorticity, psi/w grids.

w

Calculate vertical velocity on [horizontal]/[vertical] grids.

xgrid

Methods

ddeta(varname[, hcoord, scoord, hboundary, ...])

Calculate d/deta for a variable.

ddxi(varname[, hcoord, scoord, hboundary, ...])

Calculate d/dxi for a variable.

ddz(varname[, hcoord, scoord, hboundary, ...])

Calculate d/dz for a variable.

east_rotated(angle[, name])

Rotate eastward velocity by angle.

mld([thresh])

Calculate mixed layer depth [m] on rho grid.

north_rotated(angle[, name])

Rotate northward velocity by angle.

set_grid(xgrid)

If you already have a xgrid object and don't want to rerun

subset([X, Y])

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

to_grid(varname[, hcoord, scoord, ...])

Implement grid changes.

zslice(varname, depths[, z])

Interpolate var to depths.

find_horizontal_velocities

property EKE#

Calculate EKE [m^2/s^2], on rho grid.

Notes

EKE = 0.5*(ug^2 + vg^2) Puts geostrophic speed on rho grid.

See xroms.EKE for full docstring.

Examples

>>> ds.xroms.EKE
property KE#

Calculate kinetic energy [kg/(m*s^2)], on rho/rho grids.

Notes

Uses speed that has been extended out to the rho grid and rho0.

See xroms.KE for full docstring.

Examples

>>> ds.xroms.KE
property M2#

Calculate the horizontal buoyancy gradient on rho/w grids.

Notes

See xroms.M2 for full docstring.

hboundary set to ‘extend’ and sboundary=’fill’ with sfill_value=np.nan.

Examples

>>> ds.xroms.M2
property N2#

Calculate buoyancy frequency squared on rho/w grids.

Notes

See xroms.N2 for full docstring.

sboundary set to ‘fill’ with sfill_value=np.nan.

Examples

>>> ds.xroms.N2
_eastnorth2uv()[source]#

Call the velocity rotation for accessor.

_eastnorth_rotated(angle, include_vars_adcp=False, **kwargs)[source]#

Call the velocity rotation for accessor.

include_vars_adcpbool

If True, include all variables that might be compared with ADCP data and ways to convert between: east_rotated, north_rotated, angle, east, north, grid_angle.

_uv2eastnorth()[source]#

Call the velocity rotation for accessor.

property buoyancy#

Calculate buoyancy on rho/rho grids.

Notes

See xroms.buoyancy for full docstring.

Examples

>>> ds.xroms.buoyancy
property convergence#

Calculate convergence, rho/rho grid.

Notes

See xroms.convergence for full docstring.

hboundary and sboundary both set to ‘extend’.

Examples

>>> ds.xroms.convergence
property convergence_norm#

Calculate normalized surface convergence, rho/rho grid.

The surface currents are selected for this calculation, so return is [T,Y,X]. The convergence is normalized by $f$. It is dimensionless.

Notes

See xroms.convergence for full docstring.

hboundary and sboundary both set to ‘extend’.

Examples

>>> ds.xroms.convergence_norm
ddeta(varname, hcoord=None, scoord=None, hboundary='extend', hfill_value=None, sboundary='extend', sfill_value=None, attrs=None)[source]#

Calculate d/deta for a variable.

Parameters:
  • varname (str) – Name of variable in Dataset to operate on.

  • 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 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 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

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

Calculate d/dxi for a variable.

Parameters:
  • varname (str) – Name of variable in Dataset to operate on.

  • 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

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

Calculate d/dz for a variable.

Parameters:
  • varname (str) – Name of variable in Dataset to operate on.

  • 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 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

>>> ds.xroms.ddz('salt')
property dudz#

Calculate dudz [1/s] on u/w grids.

Notes

See xroms.dudz for full docstring.

sboundary is set to ‘extend’.

Examples

>>> ds.xroms.dudz
property dvdz#

Calculate dvdz [1/s] on v/w grids.

Notes

See xroms.dvdz for full docstring.

sboundary is set to ‘extend’.

Examples

>>> ds.xroms.dvdz
property east#

Rotate grid-aligned u velocity to be eastward.

Notes

See xroms.rotate_vectors for full docstring.

Examples

>>> ds.xroms.east
east_rotated(angle, name=None, **kwargs)[source]#

Rotate eastward velocity by angle.

Parameters:
  • angle (float,xr.DataArray) – Angle to rotate eastward, northward velocities by to get x component of rotated velocities.

  • name (str, optional) – If input, will be used for output array name.

  • kwargs (optional) – will be input to xroms.rotate_vectors().

Notes

See xroms.rotate_vectors() for full docstring.

Examples

>>> ds.xroms.east_rotated(angle, reference="compass", isradians=False, name="along_channel")
property ertel#

Calculate Ertel potential vorticity of buoyancy on rho/rho grids.

Notes

See xroms.ertel for full docstring.

hboundary and sboundary both set to ‘extend’.

Examples

>>> ds.xroms.ertel
find_horizontal_velocities()[source]#
mld(thresh=0.03)[source]#

Calculate mixed layer depth [m] on rho grid.

property north#

Rotate grid-aligned v velocity to be northward.

Notes

See xroms.rotate_vectors for full docstring.

Examples

>>> ds.xroms.north
north_rotated(angle, name=None, **kwargs)[source]#

Rotate northward velocity by angle.

Parameters:
  • angle (float,xr.DataArray) – Angle to rotate eastward, northward velocities by to get y component of rotated velocities.

  • name (str, optional) – If input, will be used for output array name.

  • kwargs (optional) – will be input to xroms.rotate_vectors().

Notes

See xroms.rotate_vectors() for full docstring.

Examples

>>> ds.xroms.north_rotated(angle, reference="compass", isradians=False, name="across_channel")
property omega#

Calculate s-grid vertical velocity on [horizontal]/[vertical] grids.

Notes

See xroms.omega for full docstring.

Examples

>>> ds.xroms.omega
property rho#

Return existing rho or calculate, on rho/rho grids.

Notes

See xroms.density for full docstring.

Examples

>>> ds.xroms.rho
set_grid(xgrid)[source]#

If you already have a xgrid object and don’t want to rerun

Or, you want to have more options in the xgrid setup, input it to the xroms accessor this way.

Examples

>>> ds.xroms.set_grid(xgrid)
property sig0#

Calculate potential density referenced to z=0, on rho/rho grids.

Notes

See xroms.potential_density for full docstring.

Examples

>>> ds.xroms.sig0
property speed#

Calculate horizontal speed [m/s] from u and v components, on rho/rho grids.

Notes

speed = np.sqrt(u^2 + v^2)

Uses ‘extend’ for horizontal boundary.

See xroms.speed for full docstring.

Examples

>>> ds.xroms.speed
subset(X=None, Y=None)[source]#

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

Parameters:
  • 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: >>> ds.xroms.subset(Y=slice(50,100)) Subset in X and Y: >>> ds.xroms.subset(X=slice(20,40), Y=slice(50,100))

to_grid(varname, hcoord=None, scoord=None, hboundary='extend', hfill_value=None, sboundary='extend', sfill_value=None)[source]#

Implement grid changes.

Parameters:
  • varname (str) – Name of variable in Dataset to operate on.

  • 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 interpolated onto hcoord horizontal and scoord

  • vertical grids.

Notes

If var is already on selected grid, nothing happens.

Examples

>>> ds.xroms.to_grid('salt', hcoord='rho', scoord='w')
property u#

Rotate eastward velocity to be grid-aligned u velocity

Notes

See xroms.rotate_vectors for full docstring.

Examples

>>> ds.xroms.u
property ug#

Calculate geostrophic u velocity from zeta, on u grid.

Notes

ug = -g * zeta_xi / (d xi * f) # on u grid

See xroms.uv_geostrophic for full docstring.

Examples

>>> ds.xroms.ug
property v#

Rotate northward velocity to be grid-aligned v velocity

Notes

See xroms.rotate_vectors for full docstring.

Examples

>>> ds.xroms.v
property vertical_shear#

Calculate vertical shear [1/s], rho/w grids.

Notes

See xroms.vertical_shear for full docstring.

hboundary is set to ‘extend’.

Examples

>>> ds.xroms.vertical_shear
property vg#

Calculate geostrophic v velocity from zeta, on v grid.

Notes

vg = g * zeta_eta / (d eta * f) # on v grid

See xroms.uv_geostrophic for full docstring.

Examples

>>> ds.xroms.vg
property vort#

Calculate vertical relative vorticity, psi/w grids.

Notes

See xroms.relative_vorticity for full docstring.

hboundary and sboundary both set to ‘extend’.

Examples

>>> ds.xroms.vort
property w#

Calculate vertical velocity on [horizontal]/[vertical] grids.

Notes

See xroms.w for full docstring.

Examples

>>> ds.xroms.w
property xgrid#
zslice(varname, depths, z=None)[source]#

Interpolate var to depths.

This wraps xgcm transform function for slice interpolation, though transform has additional functionality. See xroms.isoslice for full docs.

Parameters:
  • depths (list, ndarray) – Values to interpolate to (called iso_values in other functions). Should be negative if below mean sea level. If input as array, should be 1D.

  • z (DataArray, optional) – Array that var is interpolated onto (e.g., z coordinates or density). The “vertical” coordinate is selected by default. Use this option if you want to interpolate with z depths constant in time and input the appropriate z coordinate (e.g. z_rho0).

Returns:

  • DataArray of var interpolated to depths. Dimensionality will be the

  • same as var except with dim dimension of size of depths.

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:

>>> ds.temp.xroms.zslice(depths)

To calculate temperature onto fixed depths without considering time for z coord:

>>> ds.temp.xroms.zslice(depths, z=ds.temp.z_rho0)