mom_mixed_layer_restrat module reference

Parameterization of mixed layer restratification by unresolved mixed-layer eddies.

More…

Data Types

mixedlayer_restrat_cs

Control structure for mom_mixed_layer_restrat().

Functions/Subroutines

mixedlayer_restrat()

Driver for the mixed-layer restratification parameterization.

mixedlayer_restrat_om4()

Calculates a restratifying flow in the mixed layer, following the formulation used in OM4.

mu()

Stream function shape as a function of non-dimensional position within mixed-layer [nondim].

mixedlayer_restrat_bodner()

Calculates a restratifying flow in the mixed layer, following the formulation used in Bodner et al., 2023 (B22)

rmean2ts()

Two time-scale running mean [units of “signal” and “filtered”].

mixedlayer_restrat_bml()

Calculates a restratifying flow assuming a 2-layer bulk mixed layer.

growth_time()

Return the growth timescale for the submesoscale mixed layer eddies in [T ~> s].

mixedlayer_restrat_init()

Initialize the mixed layer restratification module.

mixedlayer_restrat_register_restarts()

Allocate and register fields in the mixed layer restratification structure for restarts.

mixedlayer_restrat_unit_tests()

test_answer()

Returns true if any cell of u and u_true are not identical.

Detailed Description

Mixed-layer eddy parameterization module

The subroutines in this module implement a parameterization of unresolved viscous mixed layer restratification of the mixed layer as described in Fox-Kemper et al., 2008, and whose impacts are described in Fox-Kemper et al., 2011. This is derived in part from the older parameterization that is described in Hallberg (Aha Hulikoa, 2003), which this new parameterization surpasses, which in turn is based on the sub-inertial mixed layer theory of Young (JPO, 1994). There is no net horizontal volume transport due to this parameterization, and no direct effect below the mixed layer. A revised of the parameterization by Bodner et al., 2023, is also available as an option.

This parameterization sets the restratification timescale to agree with high-resolution studies of mixed layer restratification.

The run-time parameter FOX_KEMPER_ML_RESTRAT_COEF is a non-dimensional number of order a few tens, proportional to the ratio of the deformation radius or the grid scale (whichever is smaller to the dominant horizontal length-scale of the sub-meso-scale mixed layer instabilities.

“Sub-meso” in a nutshell

The parameterization is colloquially referred to as “sub-meso”.

The original Fox-Kemper et al., (2008b) paper proposed a quasi-Stokes advection described by the stream function (eq. 5 of Fox-Kemper et al., 2011):

\[{\bf \Psi}_o = C_e \frac{ H^2 \nabla \bar{b} \times \hat{\bf z} }{ |f| } \mu(z)\]

where the vertical profile function is

\[\mu(z) = \max \left\{ 0, \left[ 1 - \left(\frac{2z}{H}+1\right)^2 \right] \left[ 1 + \frac{5}{21} \left(\frac{2z}{H}+1\right)^2 \right] \right\}\]

and \(H\) is the mixed-layer depth, \(f\) is the local Coriolis parameter, \(C_e \sim 0.06-0.08\) and \(\nabla \bar{b}\) is a depth mean buoyancy gradient averaged over the mixed layer.

For use in coarse-resolution models, an upscaling of the buoyancy gradients and adaption for the equator leads to the following parameterization (eq. 6 of Fox-Kemper et al., 2011):

\[{\bf \Psi} = C_e \Gamma_\Delta \frac{\Delta s}{l_f} \frac{ H^2 \nabla \bar{b} \times \hat{\bf z} } { \sqrt{ f^2 + \tau^{-2}} } \mu(z)\]

where \(\Delta s\) is the minimum of grid-scale and deformation radius, \(l_f\) is the width of the mixed-layer fronts, and \(\Gamma_\Delta=1\). \(\tau\) is a time-scale for mixing momentum across the mixed layer. \(l_f\) is thought to be of order hundreds of meters.

The upscaling factor \(\frac{\Delta s}{l_f}\) can be a global constant, model parameter FOX_KEMPER_ML_RESTRAT, so that in practice the parameterization is:

\[{\bf \Psi} = C_e \Gamma_\Delta \frac{ H^2 \nabla \bar{b} \times \hat{\bf z} }{ \sqrt{ f^2 + \tau^{-2}} } \mu(z)\]

with non-unity \(\Gamma_\Delta\).

\(C_e\) is hard-coded as 0.0625. \(\tau\) is calculated from the surface friction velocity \(u^*\). .. admonition:: Todo

Explain expression for momentum mixing time-scale.

Symbol

Module parameter

\(\Gamma_\Delta\)

FOX_KEMPER_ML_RESTRAT

\(l_f\)

MLE_FRONT_LENGTH

\(\Delta \rho\)

MLE_DENSITY_DIFF

Time-filtering of mixed-layer depth

Using the instantaneous mixed-layer depth is inconsistent with the finite life-time of mixed-layer instabilities. We provide a one-sided running-mean filter of mixed-layer depth, \(H\), of the form:

\[\bar{H} \leftarrow \max \left( H, \frac{ \Delta t H + \tau_h \bar{H} }{ \Delta t + \tau_h } \right)\]

which allows the effective mixed-layer depth seen by the parameterization, \(\bar{H}\), to instantaneously deepen but to decay with time-scale \(\tau_h\). \(\bar{H}\) is substituted for \(H\) in the above equations.

Symbol

Module parameter

\(\tau_h\)

MLE_MLD_DECAY_TIME

Defining the mixed-layer-depth

If the parameter MLE_USE_PBL_MLD=True then the mixed-layer depth is defined/diagnosed by the boundary-layer parameterization (e.g. ePBL, KPP, etc.).

If the parameter MLE_USE_PBL_MLD=False then the mixed-layer depth is diagnosed in this module as the depth of a given density difference, \(\Delta \rho\), with the surface where the density difference is the parameter MLE_DENSITY_DIFF.

Bodner (2023) modification

To use this variant of the parameterization, set MLE%USE_BODNER23=True which then changes the available parameters. MLE_USE_PBL_MLD must be True to use the B23 modification.

Bodner et al., 2023, (B23) use an expression for the frontal width which changes the scaling from \(H^2\) to \(h H^2\):

\[{\bf \Psi} = C_r \frac{\Delta s |f| \bar{h} \bar{H}^2 \nabla \bar{b} \times \hat{\bf z} } { \left( m_*u_*^3 + n_* w_*^3 \right)^{2/3} } \mu(z)\]

(see eq. 27 of B23). Here, the \(h\) is the activate boundary layer depth, and \(H\) is the mixed layer depth. The denominator is an approximation of the vertical turbulent momentum flux \(\overline{w'u'}\) (see eq. 18 of B23) calculated from the surface friction velocity \(u_*\), and from the surface buoyancy flux, \(B\), using the relation \(w_*^3 \sim -B h\). An advantage of this form of “sub-meso” is the denominator is well behaved at the equator but we apply a lower bound of \(w_{min}^2\) to avoid division by zero under zero forcing. As for the original Fox-Kemper parameterization, \(\nabla \bar{b}\) is the buoyancy gradient averaged over the mixed-layer.

The instantaneous boundary layer depth, \(h\), is time filtered primarily to remove the diurnal cycle:

\[\bar{h} \leftarrow \max \left( \min \left( h, \frac{ \Delta t h + \tau_{h+} \bar{h} }{ \Delta t + \tau_{h+} } \right), \frac{ \Delta t h + \tau_{h-} \bar{h} }{ \Delta t + \tau_{h-} } \right)\]

Setting \(\tau_{h+}=0\) means that when \(h>\bar{h}\) then \(\bar{h}\leftarrow h\), i.e. the effective (filtered) depth, \(\bar{h}\), is instantly deepened. When \(h<\bar{h}\) then the effective depth shoals with time-scale \(\tau_{h-}\).

A second filter is applied to \(\bar{h}\) to yield and effective “mixed layer depth”, \(\bar{H}\), defined as the deepest the boundary layer over some time-scale \(\tau_{H-}\):

\[\bar{H} \leftarrow \max \left( \min \left( \bar{h}, \frac{ \Delta t \bar{h} + \tau_{H+} \bar{H} }{ \Delta t + \tau_{H+} } \right), \frac{ \Delta t \bar{h} + \tau_{h-} \bar{H} }{ \Delta t + \tau_{H-} } \right)\]

Again, setting \(\tau_{H+}=0\) allows the effective mixed layer to instantly deepend to \(\bar{h}\).

Symbol

Module parameter

\(C_r\)

MLE%CR

\(n_*\)

MLE%BODNER_NSTAR

\(m_*\)

MLE%BODNER_MSTAR

\(w_*\)

MLE%BODNER_MSTAR

\(w_{min}^2\)

MLE%MIN_WSTAR2

\(\tau_{h+}\)

MLE%BLD_GROWING_TFILTER

\(\tau_{h-}\)

MLE%BLD_DECAYING_TFILTER

\(\tau_{H+}\)

MLE%MLD_GROWING_TFILTER

\(\tau_{H-}\)

MLE%BLD_DECAYING_TFILTER

References

Fox-Kemper, B., Ferrari, R. and Hallberg, R., 2008: Parameterization of Mixed Layer Eddies. Part I: Theory and Diagnosis J. Phys. Oceangraphy, 38 (6), p1145-1165. https://doi.org/10.1175/2007JPO3792.1

Fox-Kemper, B. and Ferrari, R. 2008: Parameterization of Mixed Layer Eddies. Part II: Prognosis and Impact J. Phys. Oceangraphy, 38 (6), p1166-1179. https://doi.org/10.1175/2007JPO3788.1

B. Fox-Kemper, G. Danabasoglu, R. Ferrari, S.M. Griffies, R.W. Hallberg, M.M. Holland, M.E. Maltrud, S. Peacock, and B.L. Samuels, 2011: Parameterization of mixed layer eddies. III: Implementation and impact in global ocean climate simulations. Ocean Modell., 39(1), p61-78. https://doi.org/10.1016/j.ocemod.2010.09.002

A.S. Bodner, B. Fox-Kemper, L. Johnson, L. P. Van Roekel, J. C. McWilliams, P. P. Sullivan, P. S. Hall, and J. Dong, 2023: Modifying the Mixed Layer Eddy Parameterization to Include Frontogenesis Arrest by Boundary Layer Turbulence. J. Phys. Oceanogr., 53(1), p323-339. https://doi.org/10.1175/JPO-D-21-0297.1

Type Documentation

type mom_mixed_layer_restrat/mixedlayer_restrat_cs

Control structure for mom_mixed_layer_restrat(). .

Type fields:
  • % id_urestrat_time [integer] :: Diagnostic identifier.

  • % id_vrestrat_time [integer] :: Diagnostic identifier.

  • % id_uhml [integer] :: Diagnostic identifier.

  • % id_vhml [integer] :: Diagnostic identifier.

  • % id_mld [integer] :: Diagnostic identifier.

  • % id_bld [integer] :: Diagnostic identifier.

  • % id_rml [integer] :: Diagnostic identifier.

  • % id_udml [integer] :: Diagnostic identifier.

  • % id_vdml [integer] :: Diagnostic identifier.

  • % id_uml [integer] :: Diagnostic identifier.

  • % id_vml [integer] :: Diagnostic identifier.

  • % id_wpup [integer] :: Diagnostic identifier.

  • % id_ustar [integer] :: Diagnostic identifier.

  • % id_bflux [integer] :: Diagnostic identifier.

  • % initialized [logical] :: True if this control structure has been initialized.

  • % ml_restrat_coef [real] :: A non-dimensional factor by which the instability is enhanced over what would be predicted based on the resolved gradients [nondim]. This increases with grid spacing^2, up to something of order 500.

  • % ml_restrat_coef2 [real] :: As for ml_restrat_coef but using the slow filtered MLD [nondim].

  • % front_length [real] :: If non-zero, is the frontal-length scale [L ~> m] used to calculate the upscaling of buoyancy gradients that is otherwise represented by the parameter FOX_KEMPER_ML_RESTRAT_COEF. If MLE_FRONT_LENGTH is non-zero, it is recommended to set FOX_KEMPER_ML_RESTRAT_COEF=1.0.

  • % mle_use_pbl_mld [logical] :: If true, use the MLD provided by the PBL parameterization. if false, MLE will calculate a MLD based on a density difference based on the parameter MLE_DENSITY_DIFF.

  • % vonkar [real] :: The von Karman constant as used for mixed layer viscosity [nondim].

  • % mle_mld_decay_time [real] :: Time-scale to use in a running-mean when MLD is retreating [T ~> s].

  • % mle_mld_decay_time2 [real] :: Time-scale to use in a running-mean when filtered MLD is retreating [T ~> s].

  • % mle_density_diff [real] :: Density difference used in detecting mixed-layer depth [R ~> kg m-3].

  • % mle_tail_dh [real] :: Fraction by which to extend the mixed-layer restratification depth used for a smoother stream function at the base of the mixed-layer [nondim].

  • % mle_mld_stretch [real] :: A scaling coefficient for stretching/shrinking the MLD used in the MLE scheme [nondim]. This simply multiplies MLD wherever used.

  • % use_bodner [logical] :: If true, use the Bodner et al., 2023, parameterization.

  • % cr [real] :: Efficiency coefficient from Bodner et al., 2023 [nondim].

  • % mstar [real] :: The m* value used to estimate the turbulent vertical momentum flux [nondim].

  • % nstar [real] :: The n* value used to estimate the turbulent vertical momentum flux [nondim].

  • % min_wstar2 [real] :: The minimum lower bound to apply to the vertical momentum flux, w’u’, in the Bodner et al., restratification parameterization [m2 s-2]. This avoids a division-by-zero in the limit when u* and the buoyancy flux are zero.

  • % bld_growing_tfilt [real] :: The time-scale for a running-mean filter applied to the boundary layer depth (BLD) when the BLD is deeper than the running mean [T ~> s]. A value of 0 instantaneously sets the running mean to the current value of BLD.

  • % bld_decaying_tfilt [real] :: The time-scale for a running-mean filter applied to the boundary layer depth (BLD) when the BLD is shallower than the running mean [T ~> s]. A value of 0 instantaneously sets the running mean to the current value of BLD.

  • % mld_decaying_tfilt [real] :: The time-scale for a running-mean filter applied to the time-filtered MLD, when the latter is shallower than the running mean [T ~> s]. A value of 0 instantaneously sets the running mean to the current value of MLD.

  • % mld_growing_tfilt [real] :: The time-scale for a running-mean filter applied to the time-filtered MLD, when the latter is deeper than the running mean [T ~> s]. A value of 0 instantaneously sets the running mean to the current value of MLD.

  • % debug [logical] :: If true, calculate checksums of fields for debugging.

  • % diag [type( diag_ctrl ),pointer] :: A structure that is used to regulate the timing of diagnostic output.

  • % use_stanley_ml [logical] :: If true, use the Stanley parameterization of SGS T variance.

  • % ustar_min [real] :: A minimum value of ustar to avoid numerical problems [Z T-1 ~> m s-1].

  • % kv_restrat [real] :: A viscosity that sets a floor on the momentum mixing rate during restratification [Z2 T-1 ~> m2 s-1].

  • % mld_filtered [real(:,:),allocatable] :: Time-filtered MLD [H ~> m or kg m-2].

  • % mld_filtered_slow [real(:,:),allocatable] :: Slower time-filtered MLD [H ~> m or kg m-2].

  • % wpup_filtered [real(:,:),allocatable] :: Time-filtered vertical momentum flux [Z2 T-2 ~> m2 s-2].

Function/Subroutine Documentation

subroutine mom_mixed_layer_restrat/mixedlayer_restrat(h, uhtr, vhtr, tv, forces, dt, MLD, bflux, VarMix, G, GV, US, CS)

Driver for the mixed-layer restratification parameterization. The code branches between two different implementations depending on whether the bulk-mixed layer or a general coordinate are in use.

Parameters:
  • g :: [inout] Ocean grid structure

  • gv :: [in] Ocean vertical grid structure

  • us :: [in] A dimensional unit scaling type

  • h :: [inout] Layer thickness [H ~> m or kg m-2]

  • uhtr :: [inout] Accumulated zonal mass flux [H L2 ~> m3 or kg]

  • vhtr :: [inout] Accumulated meridional mass flux [H L2 ~> m3 or kg]

  • tv :: [in] Thermodynamic variables structure

  • forces :: [in] A structure with the driving mechanical forces

  • dt :: [in] Time increment [T ~> s]

  • mld :: Mixed layer depth provided by the planetary boundary layer scheme [Z ~> m]

  • bflux :: Surface buoyancy flux provided by the PBL scheme [Z2 T-3 ~> m2 s-3]

  • varmix :: [in] Variable mixing control structure

  • cs :: [inout] Module control structure

Call to:

mixedlayer_restrat_bml mixedlayer_restrat_bodner mixedlayer_restrat_om4 mom_error_handler::mom_error

subroutine mom_mixed_layer_restrat/mixedlayer_restrat_om4(h, uhtr, vhtr, tv, forces, dt, MLD_in, VarMix, G, GV, US, CS)

Calculates a restratifying flow in the mixed layer, following the formulation used in OM4.

Parameters:
  • g :: [inout] Ocean grid structure

  • gv :: [in] Ocean vertical grid structure

  • us :: [in] A dimensional unit scaling type

  • h :: [inout] Layer thickness [H ~> m or kg m-2]

  • uhtr :: [inout] Accumulated zonal mass flux [H L2 ~> m3 or kg]

  • vhtr :: [inout] Accumulated meridional mass flux [H L2 ~> m3 or kg]

  • tv :: [in] Thermodynamic variables structure

  • forces :: [in] A structure with the driving mechanical forces

  • dt :: [in] Time increment [T ~> s]

  • mld_in :: Mixed layer depth provided by the PBL scheme [Z ~> m] (not H)

  • varmix :: [in] Variable mixing control structure

  • cs :: [inout] Module control structure

Call to:

mom_diag_mediator::diag_update_remap_grids mom_eos::eos_domain mom_error_handler::mom_error mu

Called from:

mixedlayer_restrat

function mom_mixed_layer_restrat/mu(sigma, dh) [real]

Stream function shape as a function of non-dimensional position within mixed-layer [nondim].

Parameters:
  • sigma :: [in] Fractional position within mixed layer [nondim] z=0 is surface, z=-1 is the bottom of the mixed layer

  • dh :: [in] Non-dimensional distance over which to extend stream function to smooth transport at base [nondim]

Called from:

mixedlayer_restrat_bodner mixedlayer_restrat_om4 mixedlayer_restrat_unit_tests

subroutine mom_mixed_layer_restrat/mixedlayer_restrat_bodner(CS, G, GV, US, h, uhtr, vhtr, tv, forces, dt, BLD, bflux)

Calculates a restratifying flow in the mixed layer, following the formulation used in Bodner et al., 2023 (B22)

Parameters:
  • cs :: [inout] Module control structure

  • g :: [inout] Ocean grid structure

  • gv :: [in] Ocean vertical grid structure

  • us :: [in] A dimensional unit scaling type

  • h :: [inout] Layer thickness [H ~> m or kg m-2]

  • uhtr :: [inout] Accumulated zonal mass flux [H L2 ~> m3 or kg]

  • vhtr :: [inout] Accumulated meridional mass flux [H L2 ~> m3 or kg]

  • tv :: [in] Thermodynamic variables structure

  • forces :: [in] A structure with the driving mechanical forces

  • dt :: [in] Time increment [T ~> s]

  • bld :: Active boundary layer depth provided by the PBL scheme [Z ~> m] (not H)

  • bflux :: Surface buoyancy flux provided by the PBL scheme [Z2 T-3 ~> m2 s-3]

Call to:

mom_diag_mediator::diag_update_remap_grids mom_eos::eos_domain mom_error_handler::mom_error mu rmean2ts

Called from:

mixedlayer_restrat

function mom_mixed_layer_restrat/rmean2ts(signal, filtered, tau_growing, tau_decaying, dt) [elemental]

Two time-scale running mean [units of “signal” and “filtered”].

If signal > filtered, returns running-mean with time scale “tau_growing”. If signal <= filtered, returns running-mean with time scale “tau_decaying”.

The running mean of \(s\) with time scale “of \(\tau\) is:

\[\bar{s} <- ( \Delta t * s + \tau * \bar{s} ) / ( \Delta t + \tau )\]

Note that if \(tau=0\), then the running mean equals the signal. Thus, rmean2ts with tau_growing=0 recovers the “resetting running mean” used in OM4.

Called from:

mixedlayer_restrat_bodner mixedlayer_restrat_unit_tests

subroutine mom_mixed_layer_restrat/mixedlayer_restrat_bml(h, uhtr, vhtr, tv, forces, dt, G, GV, US, CS)

Calculates a restratifying flow assuming a 2-layer bulk mixed layer.

Parameters:
  • g :: [in] Ocean grid structure

  • gv :: [in] Ocean vertical grid structure

  • us :: [in] A dimensional unit scaling type

  • h :: [inout] Layer thickness [H ~> m or kg m-2]

  • uhtr :: [inout] Accumulated zonal mass flux [H L2 ~> m3 or kg]

  • vhtr :: [inout] Accumulated meridional mass flux [H L2 ~> m3 or kg]

  • tv :: [in] Thermodynamic variables structure

  • forces :: [in] A structure with the driving mechanical forces

  • dt :: [in] Time increment [T ~> s]

  • cs :: [inout] Module control structure

Called from:

mixedlayer_restrat

function mom_mixed_layer_restrat/growth_time(u_star, hBL, absf, h_neg, vonKar, Kv_rest, restrat_coef) [real]

Return the growth timescale for the submesoscale mixed layer eddies in [T ~> s].

Parameters:
  • u_star :: [in] Surface friction velocity [Z T-1 ~> m s-1]

  • hbl :: [in] Boundary layer thickness including at least a neglible value to keep it positive definite [Z ~> m]

  • absf :: [in] Absolute value of the Coriolis parameter [T-1 ~> s-1]

  • h_neg :: [in] A tiny thickness that is usually lost in roundoff so can be neglected [Z ~> m]

  • kv_rest :: [in] The background laminar vertical viscosity used for restratification [Z2 T-1 ~> m2 s-1]

  • vonkar :: [in] The von Karman constant, used to scale the turbulent limits on the restratification timescales [nondim]

  • restrat_coef :: [in] An overall scaling factor for the restratification timescale [nondim]

function mom_mixed_layer_restrat/mixedlayer_restrat_init(Time, G, GV, US, param_file, diag, CS, restart_CS) [logical]

Initialize the mixed layer restratification module.

Parameters:
  • time :: [in] Current model time

  • g :: [inout] Ocean grid structure

  • gv :: [in] Ocean vertical grid structure

  • us :: [in] A dimensional unit scaling type

  • param_file :: [in] Parameter file to parse

  • diag :: [inout] Regulate diagnostics

  • cs :: [inout] Module control structure

  • restart_cs :: [in] MOM restart control structure

Call to:

mdl

Called from:

mixedlayer_restrat_register_restarts

subroutine mom_mixed_layer_restrat/mixedlayer_restrat_register_restarts(HI, GV, US, param_file, CS, restart_CS)

Allocate and register fields in the mixed layer restratification structure for restarts.

Parameters:
  • hi :: [in] Horizontal index structure

  • gv :: [in] Ocean vertical grid structure

  • us :: [in] A dimensional unit scaling type

  • param_file :: [in] Parameter file to parse

  • cs :: [inout] Module control structure

  • restart_cs :: [inout] MOM restart control structure

Call to:

mdl mixedlayer_restrat_init

Called from:

mom::initialize_mom

function mom_mixed_layer_restrat/mixedlayer_restrat_unit_tests(verbose) [logical]
Parameters:

verbose :: [in] If true, write results to stdout

Call to:

mu rmean2ts test_answer

Called from:

mom_unit_tests::unit_tests

function mom_mixed_layer_restrat/test_answer(verbose, u, u_true, label, tol) [logical]

Returns true if any cell of u and u_true are not identical. Returns false otherwise.

Parameters:
  • verbose :: [in] If true, write results to stdout

  • u :: [in] Values to test

  • u_true :: [in] Values to test against (correct answer)

  • label :: [in] Message

  • tol :: [in] The tolerance for differences between u and u_true

Called from:

mixedlayer_restrat_unit_tests