mom_diagnose_mld module reference
Provides functions for some diabatic processes such as fraxil, brine rejection, tendency due to surface flux divergence.
Functions/Subroutines
Diagnose a mixed layer depth (MLD) determined by a given density difference with the surface. |
|
Diagnose a mixed layer depth (MLD) determined by the depth a given energy value would mix. |
Detailed Description
This module contains subroutines that apply various diabatic processes. Usually these subroutines are called from the MOM_diabatic module. All of these routines use appropriate limiters or logic to work properly with arbitrary layer thicknesses (including massless layers) and an arbitrarily large timestep.
The subroutine diagnoseMLDbyDensityDifference diagnoses a mixed layer depth based on a density difference criterion, and may also estimate the stratification of the water below this diagnosed mixed layer.
The subroutine diagnoseMLDbyEnergy diagnoses a mixed layer depth based on a mixing-energy criterion, as described by Reichl et al., 2022, JGR: Oceans, doi:10.1029/2021JC018140.
Function/Subroutine Documentation
- subroutine mom_diagnose_mld/diagnosemldbydensitydifference(id_MLD, h, tv, densityDiff, G, GV, US, diagPtr, ref_h_mld, id_ref_z, id_ref_rho, id_N2subML, id_MLDsq, dz_subML, MLD_out)
Diagnose a mixed layer depth (MLD) determined by a given density difference with the surface. This routine is appropriate in MOM_diabatic_aux due to its position within the time stepping.
- Parameters:
g :: [in] Grid type
gv :: [in] ocean vertical grid structure
us :: [in] A dimensional unit scaling type
id_mld :: [in] Handle (ID) of MLD diagnostic
h ::
h[in] Layer thickness [H ~> m or kg m-2]tv ::
tv[in] Structure containing pointers to any available thermodynamic fields.densitydiff :: [in] Density difference to determine MLD [R ~> kg m-3]
diagptr :: Diagnostics structure
ref_h_mld ::
ref_h_mld[in] Depth of the calculated “surface” densisty [Z ~> m]id_ref_z ::
id_ref_z[in] Handle (ID) of reference depth diagnosticid_ref_rho ::
id_ref_rho[in] Handle (ID) of reference density diagnosticid_n2subml :: [in] Optional handle (ID) of subML stratification
id_mldsq :: [in] Optional handle (ID) of squared MLD
dz_subml :: [in] The distance over which to calculate N2subML or 50 m if missing [Z ~> m]
mld_out :: [out] Send MLD to other routines [Z ~> m]
- Call to:
- Called from:
- subroutine mom_diagnose_mld/diagnosemldbyenergy(id_MLD, h, tv, G, GV, US, Mixing_Energy, k_bounds, diagPtr, OM4_iteration, MLD_out)
Diagnose a mixed layer depth (MLD) determined by the depth a given energy value would mix. This routine is appropriate in MOM_diabatic_aux due to its position within the time stepping.
- Parameters:
id_mld :: [in] Energy output diagnostic IDs
g :: [in] Grid type
gv :: [in] ocean vertical grid structure
us :: [in] A dimensional unit scaling type
mixing_energy :: [in] Energy values for up to 3 MLDs [R Z3 T-2 ~> J m-2]
h ::
h[in] Layer thickness [H ~> m or kg m-2]tv ::
tv[in] Structure containing pointers to any available thermodynamic fields.diagptr :: Diagnostics structure
k_bounds ::
k_bounds[in] vertical interface bounds to apply calculationsom4_iteration :: [in] Uses a legacy version of the MLD iteration it is kept to reproduce OM4 output
mld_out :: [out] Send MLD to other routines [Z ~> m]
- Called from: