xroms.roms_seawater#
Functions related to density of seawater.
Functions
|
Calculate the horizontal buoyancy gradient. |
|
Calculate buoyancy frequency squared (vertical buoyancy gradient). |
|
Calculate buoyancy [m/s^2] based on potential density. |
|
Calculate the density [kg/m^3] as calculated in ROMS. |
|
Calculate the mixed layer depth [m], return positive and as depth if no value calculated. |
|
Calculate potential density [kg/m^3] with constant depth reference. |
- xroms.roms_seawater.M2(rho, xgrid, rho0=1025.0, hboundary='extend', hfill_value=None, sboundary='fill', sfill_value=nan, z=None)[source]#
Calculate the horizontal buoyancy gradient.
- Parameters:
rho (DataArray) – Density [kg/m^3]
xgrid (xgcm.grid) – Grid object associated with rho
rho0 (int, float, optional) – Reference density [kg/m^3].
hboundary (string, optional) –
Passed to grid method calls; horizontal boundary selection for calculating horizontal derivatives of rho. 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 rho. 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’.
z (DataArray, optional) – Depths [m] associated with rho. If None, use z coordinate attached to temperature.
- Returns:
DataArray of the horizontal buoyancy gradient on rho/w grids.
Output is [T,Z,Y,X].
Notes
M2 = g/rho0 * sqrt(d(rho)/dxi^2 + d(rho)deta^2)
g=9.81 [m/s^2]
Examples
>>> xroms.M2(rho, xgrid)
- xroms.roms_seawater.N2(rho, xgrid, rho0=1025.0, sboundary='fill', sfill_value=nan)[source]#
Calculate buoyancy frequency squared (vertical buoyancy gradient).
- Parameters:
rho (DataArray) – Density [kg/m^3]
xgrid (xgcm.grid) – Grid object associated with rho
rho0 (int, float) – Reference density [kg/m^3].
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 buoyancy frequency squared on rho/w grids.
Output is [T,Z,Y,X].
Notes
N2 = -g d(rho)/dz / rho0
Examples
>>> xroms.N2(rho, xgrid)
- xroms.roms_seawater.buoyancy(sig0, rho0=1025.0)[source]#
Calculate buoyancy [m/s^2] based on potential density.
- Parameters:
sig0 (DataArray, ndarray) – Potential density [kg/m^3]
rho0 (int, float, optional) – Reference density [kg/m^3].
- Returns:
DataArray or ndarray of calculated buoyancy on rho/rho grids.
Output is [T,Z,Y,X].
Notes
buoyancy = -g * rho / rho0
Uses equation of state based on ROMS Nonlinear/rho_eos.F
g=9.81 [m/s^2]
Examples
>>> xroms.potential_density(ds.temp, ds.salt)
- xroms.roms_seawater.density(temp, salt, z=None)[source]#
Calculate the density [kg/m^3] as calculated in ROMS.
- Parameters:
temp (DataArray, ndarray) – Temperature [Celsius]
salt (DataArray, ndarray) – Salinity
z (DataArray, ndarray, int, float, optional) – Depth [m]. To specify a reference depth, use a constant. If None, use z coordinate attached to temperature.
- Returns:
DataArray or ndarray of calculated density on rho/rho grids.
Output is [T,Z,Y,X].
Notes
Equation of state based on ROMS Nonlinear/rho_eos.F
Examples
>>> xroms.density(ds.temp, ds.salt)
- xroms.roms_seawater.mld(sig0, xgrid, h, mask, z=None, thresh=0.03)[source]#
Calculate the mixed layer depth [m], return positive and as depth if no value calculated.
- Parameters:
sig0 (DataArray) – Potential density [kg/m^3]
xgrid – xgcm grid
h (DataArray, ndarray) – Depths [m].
mask (DataArray, ndarray) – mask to match sig0
z (DataArray, ndarray, optional) – The vertical depths associated with sig0. Should be on ‘rho’ grid horizontally and vertically. Use z coords associated with DataArray sig0 if not input.
thresh (float, optional) – For detection of mixed layer [kg/m^3]
- Returns:
DataArray of mixed layer depth on rho horizontal grid.
Output is [T,Y,X].
Notes
Mixed layer depth is based on the fixed Potential Density (PD) threshold.
Converted to xroms by K. Thyng Aug 2020 from:
Update history: v1.0 DL 2020Jun07
References: ncl mixed_layer_depth function at https://github.com/NCAR/ncl/blob/ed6016bf579f8c8e8f77341503daef3c532f1069/ni/src/lib/nfpfort/ocean.f de Boyer Montégut, C., Madec, G., Fischer, A. S., Lazar, A., & Iudicone, D. (2004). Mixed layer depth over the global ocean: An examination of profile data and a profile‐based climatology. Journal of Geophysical Research: Oceans, 109(C12).
Useful resources:
Climate Data Toolbox documentation: https://www.chadagreene.com/CDT/mld_documentation.html
MLD calculation from MDTF: https://github.com/NOAA-GFDL/MDTF-diagnostics/blob/437d30590c45e8b7dd0cd01a3dc67066a2137115/diagnostics/mixed_layer_depth/mixed_layer_depth.py#L147
Examples
>>> xroms.mld(sig0, h, mask)
- xroms.roms_seawater.potential_density(temp, salt, z=0)[source]#
Calculate potential density [kg/m^3] with constant depth reference.
- Parameters:
temp (DataArray, ndarray) – Temperature [Celsius]
salt (DataArray, ndarray) – Salinity
z (int, float, optional) – Reference depth [m].
- Returns:
DataArray or ndarray of calculated potential density on rho/rho grids.
Output is [T,Z,Y,X].
Notes
Uses equation of state based on ROMS Nonlinear/rho_eos.F
Examples
>>> xroms.potential_density(ds.temp, ds.salt)