xroms.derived#

Variables derived from ROMS output are here.

Functions

EKE(ug, vg, xgrid[, hboundary, hfill_value])

Calculate EKE [m^2/s^2]

KE(rho0, speed)

Calculate kinetic energy [kg/(m*s^2)]

convergence(u, v, xgrid[, hboundary, ...])

Calculate 2D convergence from u and v [1/s].

dudz(u, xgrid[, sboundary, sfill_value])

Calculate the xi component of vertical shear [1/s]

dvdz(v, xgrid[, sboundary, sfill_value])

Calculate the eta component of vertical shear [1/s]

ertel(phi, u, v, f, xgrid[, hcoord, scoord, ...])

Calculate Ertel potential vorticity of phi.

omega(u, v)

Calculate s-grid vertical velocity from u and v [m/s]

relative_vorticity(u, v, xgrid[, hboundary, ...])

Calculate the vertical component of the relative vorticity [1/s]

speed(u, v, xgrid[, hboundary, hfill_value])

Calculate horizontal speed [m/s] from u and v components

uv_geostrophic(zeta, f, xgrid[, hboundary, ...])

Calculate geostrophic velocities from zeta [m/s]

vertical_shear(dudz, dvdz, xgrid[, ...])

Calculate the vertical shear [1/s]

w(u, v)

Calculate vertical velocity from u and v [m/s]

xroms.derived.EKE(ug, vg, xgrid, hboundary='extend', hfill_value=None)[source]#

Calculate EKE [m^2/s^2]

Parameters:
  • ug (DataArray) – Geostrophic or other xi component velocity [m/s]

  • vg (DataArray) – Geostrophic or other eta component velocity [m/s]

  • xgrid (xgcm.grid) – Grid object associated with ug, vg

  • hboundary (string, optional) –

    Passed to grid method calls; horizontal boundary selection for moving to rho grid. 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 for moving to rho grid. From xgcm documentation: The value to use in the boundary condition with boundary=’fill’.

Returns:

  • DataArray of eddy kinetic energy on rho grid.

  • Output is [T,Y,X].

Notes

EKE = 0.5*(ug^2 + vg^2)

Examples

>>> ug, vg = xroms.uv_geostrophic(ds.zeta, ds.f, xgrid)
>>> xroms.EKE(ug, vg, xgrid)
xroms.derived.KE(rho0, speed)[source]#

Calculate kinetic energy [kg/(m*s^2)]

Parameters:
  • rho0 (float) – background density of the water [kg/m^3]

  • speed (DataArray) – magnitude of horizontal velocity vector [m/s]

Returns:

  • DataArray of kinetic energy on rho/rho grids.

  • Output is [T,Z,Y,X].

Notes

KE = 0.5*rho*(u^2 + v^2)

Examples

>>> speed = xroms.speed(ds.u, ds.v, xgrid)
>>> xroms.KE(ds.rho0, speed)
xroms.derived.convergence(u, v, xgrid, hboundary='extend', hfill_value=None, sboundary='extend', sfill_value=None)[source]#

Calculate 2D convergence from u and v [1/s].

Parameters:
  • u (DataArray) – xi component of velocity [m/s]

  • v (DataArray) – eta component of velocity [m/s]

  • xgrid (xgcm.grid) – Grid object associated with u, v

  • hboundary (string, optional) –

    Passed to grid method calls; horizontal boundary selection for calculating horizontal derivatives of u and v. 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 u and v. 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 of 2D convergence of horizontal currents on rho/rho grids.

  • Output is [T,Z,Y,X].

Notes

2D convergence = u_x + v_y Resource for more information: https://uw.pressbooks.pub/ocean285/chapter/the-divergence/

Examples

>>> ds, xgrid = xroms.roms_dataset(ds)
>>> xroms.convergence(u, v, xgrid)
xroms.derived.dudz(u, xgrid, sboundary='extend', sfill_value=None)[source]#

Calculate the xi component of vertical shear [1/s]

Parameters:
  • u (DataArray) – xi component of velocity [m/s]

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

  • sboundary (string, optional) –

    Passed to grid method calls; vertical boundary selection for calculating z derivative. 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’.

Returns:

  • DataArray of xi component of vertical shear on u/w grids.

  • Output is [T,Z,Y,X].

Notes

u_z = ddz(u) Wrapper of ddz

Examples

>>> xroms.dudz(u, xgrid)
xroms.derived.dvdz(v, xgrid, sboundary='extend', sfill_value=None)[source]#

Calculate the eta component of vertical shear [1/s]

Parameters:
  • v (DataArray) – eta component of velocity [m/s]

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

  • sboundary (string, optional) –

    Passed to grid method calls; vertical boundary selection for calculating z derivative. 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’.

Returns:

  • DataArray of eta component of vertical shear on v/w grids.

  • Output is [T,Z,Y,X].

Notes

v_z = ddz(v) Wrapper of ddz

Examples

>>> xroms.dvdz(v, xgrid)
xroms.derived.ertel(phi, u, v, f, xgrid, hcoord='rho', scoord='s_rho', hboundary='extend', hfill_value=None, sboundary='extend', sfill_value=None)[source]#

Calculate Ertel potential vorticity of phi.

Parameters:
  • phi (DataArray) – Conservative tracer. Usually this would be the buoyancy but could be another approximately conservative tracer. The buoyancy can be calculated as: >>> xroms.buoyancy(temp, salt, 0) and then input as phi.

  • u (DataArray) – xi component of velocity [m/s]

  • v (DataArray) – eta component of velocity [m/s]

  • f (DataArray) – Coriolis parameter [1/s]

  • xgrid (xgcm.grid) – Grid object associated with u, v

  • 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 phi and for calculating relative vorticity. 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 and vertical derivatives of phi, and for calculating relative vorticity. 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’.

Returns:

  • DataArray of the Ertel potential vorticity for the input tracer.

  • Output is [T,Z,Y,X].

Notes

epv = -v_z * phi_x + u_z * phi_y + (f + v_x - u_y) * phi_z

This is not set up to accept different boundary choices for different variables.

Example usage: >>> xroms.ertel(ds.dye_01, ds.u, ds.v, ds.f, xgrid, scoord=’s_w’);

xroms.derived.omega(u, v)[source]#

Calculate s-grid vertical velocity from u and v [m/s]

TO BE INPUT BY VRX.

Parameters:
  • u (DataArray) – xi component of velocity [m/s]

  • v (DataArray) – eta component of velocity [m/s]

Returns:

  • DataArray of vertical component of velocity with respect to the s grid

  • on [horizontal]/[vertical] grids.

  • Output is [T,Z,Y,X].

Notes

[Give calculation]

Examples

>>> xroms.omega(u, v)
xroms.derived.relative_vorticity(u, v, xgrid, hboundary='extend', hfill_value=None, sboundary='extend', sfill_value=None)[source]#

Calculate the vertical component of the relative vorticity [1/s]

Parameters:
  • u (DataArray) – xi component of velocity [m/s]

  • v (DataArray) – eta component of velocity [m/s]

  • xgrid (xgcm.grid) – Grid object associated with u, v

  • hboundary (string, optional) –

    Passed to grid method calls; horizontal boundary selection for calculating horizontal derivatives of u and v. 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 u and v. 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 of vertical component of relative vorticity psi/w grids.

  • Output is [T,Z,Y,X].

Notes

relative_vorticity = v_x - u_y

Examples

>>> xroms.relative_vorticity(u, v, xgrid)
xroms.derived.speed(u, v, xgrid, hboundary='extend', hfill_value=None)[source]#

Calculate horizontal speed [m/s] from u and v components

Parameters:
  • u (DataArray) – xi component of velocity [m/s]

  • v (DataArray) – eta component of velocity [m/s]

  • xgrid (xgcm.grid) – Grid object associated with u, v

  • hboundary (string, optional) –

    Passed to grid method calls; horizontal boundary selection for moving to rho grid. 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 fill value selection for moving to rho grid. From xgcm documentation: The value to use in the boundary condition with boundary=’fill’.

Returns:

  • DataArray of speed calculated on rho/rho grids.

  • Output is [T,Z,Y,X].

Notes

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

xroms.derived.uv_geostrophic(zeta, f, xgrid, hboundary='extend', hfill_value=None, which='both')[source]#

Calculate geostrophic velocities from zeta [m/s]

Parameters:
  • zeta (DataArray) – sea surface height [m]

  • f (DataArray or ndarray) – Coriolis parameter [1/s]

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

  • hboundary (string, optional) –

    Passed to grid method calls; horizontal boundary selection for moving f to rho grid. 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 for moving f to rho grid. From xgcm documentation: The value to use in the boundary condition with boundary=’fill’.

  • which (string, optional) –

    Which components of geostrophic velocity to return.

    • ’both’: return both components of hgrad

    • ’xi’: return only xi-direction.

    • ’eta’: return only eta-direction.

Returns:

  • DataArrays of components of geostrophic velocity

  • calculated on their respective grids.

  • Output is [T,Y,X].

Notes

ug = -g * zeta_eta / (d eta * f) # on u grid

vg = g * zeta_xi / (d xi * f) # on v grid

Translation to Python of Matlab copy of surf_geostr_vel of IRD Roms_Tools.

Good resourcefor more information: https://uw.pressbooks.pub/ocean285/chapter/geostrophic-balance/

Examples

>>> xroms.uv_geostrophic(ds.zeta, ds.f, xgrid)
xroms.derived.vertical_shear(dudz, dvdz, xgrid, hboundary='extend', hfill_value=None)[source]#

Calculate the vertical shear [1/s]

Parameters:
  • dudz (DataArray) – xi component of vertical shear [1/s]

  • dvdz (DataArray) – eta compoenent of vertical shear [1/s]

  • xgrid (xgcm.grid) – Grid object associated with dudz, dvdz

  • hboundary (string, optional) –

    Passed to grid method calls; horizontal boundary selection for moving dudz and dvdz to rho grid. 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 for moving to rho grid. From xgcm documentation: The value to use in the boundary condition with boundary=’fill’.

Returns:

  • DataArray of vertical shear on rho/w grids.

  • Output is [T,Z,Y,X].

Notes

vertical_shear = np.sqrt(u_z^2 + v_z^2)

Examples

>>> xroms.vertical_shear(dudz, dvdz, xgrid)
xroms.derived.w(u, v)[source]#

Calculate vertical velocity from u and v [m/s]

TO BE INPUT BY VRX.

Parameters:
  • u (DataArray) – xi component of velocity [m/s]

  • v (DataArray) – eta component of velocity [m/s]

Returns:

  • DataArray of vertical component of velocity on [horizontal]/[vertical] grids.

  • Output is [T,Z,Y,X].

Notes

[Give calculation]

Examples

>>> xroms.w(u, v)