mom_lateral_mixing_coeffs module reference

Variable mixing coefficients.

More…

Data Types

varmix_cs

Variable mixing coefficients.

Functions/Subroutines

calc_depth_function()

Calculates the non-dimensional depth functions.

calc_resoln_function()

Calculates and stores the non-dimensional resolution functions.

calc_slope_functions()

Calculates and stores functions of isopycnal slopes, e.g.

calc_visbeck_coeffs()

Calculates factors used when setting diffusivity coefficients similar to Visbeck et al.

calc_slope_functions_using_just_e()

The original calc_slope_function() that calculated slopes using interface positions only, not accounting for density variations.

calc_qg_leith_viscosity()

Calculates the Leith Laplacian and bi-harmonic viscosity coefficients.

varmix_init()

Initializes the variables mixing coefficients container.

Detailed Description

This module provides a container for various factors used in prescribing diffusivities, that are a function of the state (in particular the stratification and isoneutral slopes).

The resolution function

The resolution function is expressed in terms of the ratio of grid-spacing to deformation radius. The square of the resolution parameter is

\[R^2 = \frac{L_d^2}{\Delta^2} = \frac{ c_g^2 }{ f^2 \Delta^2 + c_g \beta \Delta^2 }\]

where the grid spacing is calculated as

\[\Delta^2 = \Delta x^2 + \Delta y^2 .\]

Todo

Check this reference to Bob on/off paper. The resolution function used in scaling diffusivities (Hallberg, 2010) is

\[r(\Delta,L_d) = \frac{1}{1+(\alpha R)^p}\]

The resolution function can be applied independently to thickness diffusion (module mom_thickness_diffuse()), tracer diffusion (mom_tracer_hordiff) lateral viscosity ( ), tracer diffusion (mom_tracer_hordiff) lateral viscosity (mom_hor_visc()).).

Robert Hallberg, 2013: Using a resolution function to regulate parameterizations of oceanic mesoscale eddy effects. Ocean Modelling, 71, pp 92-103. http://dx.doi.org/10.1016/j.ocemod.2013.08.007

Symbol

Module parameter

USE_VARIABLE_MIXING

RESOLN_SCALED_KH

RESOLN_SCALED_KHTH

RESOLN_SCALED_KHTR

\(\alpha\)

KH_RES_SCALE_COEF (for thickness and tracer diffusivity)

\(p\)

KH_RES_FN_POWER (for thickness and tracer diffusivity)

\(\alpha\)

VISC_RES_SCALE_COEF (for lateral viscosity)

\(p\)

VISC_RES_FN_POWER (for lateral viscosity)

GILL_EQUATORIAL_LD

Visbeck diffusivity

This module also calculates factors used in setting the thickness diffusivity similar to a Visbeck et al., 1997, scheme. The factors are combined in mom_thickness_diffuse::thickness_diffuse() but calculated in this module. but calculated in this module.

\[\kappa_h = \alpha_s L_s^2 S N\]

where \(S\) is the magnitude of the isoneutral slope and \(N\) is the Brunt-Vaisala frequency.

Visbeck, Marshall, Haine and Spall, 1997: Specification of Eddy Transfer Coefficients in Coarse-Resolution Ocean Circulation Models. J. Phys. Oceanogr. http://dx.doi.org/10.1175/1520-0485(1997)027%3C0381:SOETCI%3E2.0.CO;2

Symbol

Module parameter

USE_VARIABLE_MIXING

\(\alpha_s\)

KHTH_SLOPE_CFF (for mom_thickness_diffuse() module) module)

\(\alpha_s\)

KHTR_SLOPE_CFF (for mom_tracer_hordiff module)

\(L_{s}\)

VISBECK_L_SCALE

\(S_{max}\)

VISBECK_MAX_SLOPE

Vertical structure function for KhTh

The thickness diffusivity can be prescribed a vertical distribution with the shape of the equivalent barotropic velocity mode. The structure function is stored in the control structure for thie module (varmix_cs()) but is calculated using subroutines in ) but is calculated using subroutines in mom_wave_speed()..

Symbol

Module parameter

KHTH_USE_EBT_STRUCT

Type Documentation

type mom_lateral_mixing_coeffs/varmix_cs

Variable mixing coefficients.

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

  • % id_sn_v [integer] :: Diagnostic identifier.

  • % id_l2u [integer] :: Diagnostic identifier.

  • % id_l2v [integer] :: Diagnostic identifier.

  • % id_res_fn [integer] :: Diagnostic identifier.

  • % id_n2_u [integer] :: Diagnostic identifier.

  • % id_n2_v [integer] :: Diagnostic identifier.

  • % id_s2_u [integer] :: Diagnostic identifier.

  • % id_s2_v [integer] :: Diagnostic identifier.

  • % id_rd_dx [integer] :: Diagnostic identifier.

  • % id_kh_u_qg [integer] :: Diagnostic identifier.

  • % id_kh_v_qg [integer] :: Diagnostic identifier.

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

  • % use_variable_mixing [logical] :: If true, use the variable mixing.

  • % resoln_scaled_kh [logical] :: If true, scale away the Laplacian viscosity when the deformation radius is well resolved.

  • % resoln_scaled_khth [logical] :: If true, scale away the thickness diffusivity when the deformation radius is well resolved.

  • % depth_scaled_khth [logical] :: If true, KHTH is scaled away when the depth is shallower than a reference depth.

  • % resoln_scaled_khtr [logical] :: If true, scale away the tracer diffusivity when the deformation radius is well resolved.

  • % interpolate_res_fn [logical] :: If true, interpolate the resolution function to the velocity points from the thickness points; otherwise interpolate the wave speed and calculate the resolution function independently at each point.

  • % use_stored_slopes [logical] :: If true, stores isopycnal slopes in this structure.

  • % resoln_use_ebt [logical] :: If true, uses the equivalent barotropic wave speed instead of first baroclinic wave for calculating the resolution fn.

  • % khth_use_ebt_struct [logical] :: If true, uses the equivalent barotropic structure as the vertical structure of thickness diffusivity.

  • % calculate_cg1 [logical] :: If true, calls

  • % calculate_rd_dx [logical] :: If true, calculates Rd/dx and populate CSRd_dx_h. This parameter is set depending on other parameters.

  • % calculate_res_fns [logical] :: If true, calculate all the resolution factors. This parameter is set depending on other parameters.

  • % calculate_depth_fns [logical] :: If true, calculate all the depth factors. This parameter is set depending on other parameters.

  • % calculate_eady_growth_rate [logical] :: If true, calculate all the Eady growth rate. This parameter is set depending on other parameters.

  • % sn_u [real(:,:),pointer] :: S*N at u-points [T-1 ~> s-1].

  • % sn_v [real(:,:),pointer] :: S*N at v-points [T-1 ~> s-1].

  • % l2u [real(:,:),pointer] :: Length scale^2 at u-points [L2 ~> m2].

  • % l2v [real(:,:),pointer] :: Length scale^2 at v-points [L2 ~> m2].

  • % cg1 [real(:,:),pointer] :: The first baroclinic gravity wave speed [L T-1 ~> m s-1].

  • % res_fn_h [real(:,:),pointer] :: Non-dimensional function of the ratio the first baroclinic.

  • % res_fn_q [real(:,:),pointer] :: Non-dimensional function of the ratio the first baroclinic.

  • % res_fn_u [real(:,:),pointer] :: Non-dimensional function of the ratio the first baroclinic.

  • % res_fn_v [real(:,:),pointer] :: Non-dimensional function of the ratio the first baroclinic.

  • % depth_fn_u [real(:,:),pointer] :: Non-dimensional function of the ratio of the depth to.

  • % depth_fn_v [real(:,:),pointer] :: Non-dimensional function of the ratio of the depth to.

  • % beta_dx2_h [real(:,:),pointer] :: The magnitude of the gradient of the Coriolis parameter.

  • % beta_dx2_q [real(:,:),pointer] :: The magnitude of the gradient of the Coriolis parameter.

  • % beta_dx2_u [real(:,:),pointer] :: The magnitude of the gradient of the Coriolis parameter.

  • % beta_dx2_v [real(:,:),pointer] :: The magnitude of the gradient of the Coriolis parameter.

  • % f2_dx2_h [real(:,:),pointer] :: The Coriolis parameter squared times the grid.

  • % f2_dx2_q [real(:,:),pointer] :: The Coriolis parameter squared times the grid.

  • % f2_dx2_u [real(:,:),pointer] :: The Coriolis parameter squared times the grid.

  • % f2_dx2_v [real(:,:),pointer] :: The Coriolis parameter squared times the grid.

  • % rd_dx_h [real(:,:),pointer] :: Deformation radius over grid spacing [nondim].

  • % slope_x [real(:,:,:),pointer] :: Zonal isopycnal slope [nondim].

  • % slope_y [real(:,:,:),pointer] :: Meridional isopycnal slope [nondim].

  • % ebt_struct [real(:,:,:),pointer] :: Vertical structure function to scale diffusivities with [nondim].

  • % real (* kh_v_qg [*, *) :: Laplacian metric-dependent constants [L3 ~> m3].

  • % real :: Laplacian metric-dependent constants [L3 ~> m3].

  • % real :: QG Leith GM coefficient at u-points [L2 T-1 ~> m2 s-1].

  • % real :: QG Leith GM coefficient at v-points [L2 T-1 ~> m2 s-1].

  • % use_visbeck [logical] :: Use Visbeck formulation for thickness diffusivity.

  • % varmix_ktop [integer] :: Top layer to start downward integrals.

  • % visbeck_l_scale [real] :: Fixed length scale in Visbeck formula.

  • % res_coef_khth [real] :: A non-dimensional number that determines the function of resolution, used for thickness and tracer mixing, as: F = 1 / (1 + (Res_coef_khth*Ld/dx)^Res_fn_power)

  • % res_coef_visc [real] :: A non-dimensional number that determines the function of resolution, used for lateral viscosity, as: F = 1 / (1 + (Res_coef_visc*Ld/dx)^Res_fn_power)

  • % depth_scaled_khth_h0 [real] :: The depth above which KHTH is linearly scaled away [Z ~> m].

  • % depth_scaled_khth_exp [real] :: The exponent used in the depth dependent scaling function for KHTH [nondim].

  • % kappa_smooth [real] :: A diffusivity for smoothing T/S in vanished layers [Z2 T-1 ~> m2 s-1].

  • % res_fn_power_khth [integer] :: The power of dx/Ld in the KhTh resolution function. Any positive integer power may be used, but even powers and especially 2 are coded to be more efficient.

  • % res_fn_power_visc [integer] :: The power of dx/Ld in the Kh resolution function. Any positive integer power may be used, but even powers and especially 2 are coded to be more efficient.

  • % visbeck_s_max [real] :: Upper bound on slope used in Eady growth rate [nondim].

  • % use_qg_leith_gm [logical] :: If true, uses the QG Leith viscosity as the GM coefficient.

  • % use_beta_in_qg_leith [logical] :: If true, includes the beta term in the QG Leith GM coefficient.

  • % wave_speed_csp [type(wave_speed_cs),pointer] :: Wave speed control structure.

  • % pass_cg1 [type(group_pass_type)] :: For group halo pass.

  • % debug [logical] :: If true, write out checksums of data for debugging.

Function/Subroutine Documentation

subroutine mom_lateral_mixing_coeffs/calc_depth_function(G, CS)

Calculates the non-dimensional depth functions.

Parameters
  • g :: [in] Ocean grid structure

  • cs :: Variable mixing coefficients

Call to

mom_error_handler::mom_error

Called from

mom::step_mom

subroutine mom_lateral_mixing_coeffs/calc_resoln_function(h, tv, G, GV, US, CS)

Calculates and stores the non-dimensional resolution functions.

Parameters
  • g :: [inout] Ocean grid structure

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

  • tv :: [in] Thermodynamic variables

  • gv :: [in] Vertical grid structure

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

  • cs :: Variable mixing coefficients

Call to

mom_error_handler::mom_error mom_diag_mediator::query_averaging_enabled mom_wave_speed::wave_speed

Called from

mom::step_mom

subroutine mom_lateral_mixing_coeffs/calc_slope_functions(h, tv, dt, G, GV, US, CS, OBC)

Calculates and stores functions of isopycnal slopes, e.g. Sx, Sy, S*N, mostly used in the Visbeck et al. style scaling of diffusivity.

Parameters
  • g :: [inout] Ocean grid structure

  • gv :: [in] Vertical grid structure

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

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

  • tv :: [in] Thermodynamic variables

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

  • cs :: Variable mixing coefficients

  • obc :: Open boundaries control structure.

Call to

mom_isopycnal_slopes::calc_isoneutral_slopes calc_slope_functions_using_just_e calc_visbeck_coeffs mom_error_handler::mom_error mom_diag_mediator::query_averaging_enabled

subroutine mom_lateral_mixing_coeffs/calc_visbeck_coeffs(h, slope_x, slope_y, N2_u, N2_v, G, GV, US, CS, OBC)

Calculates factors used when setting diffusivity coefficients similar to Visbeck et al.

Parameters
  • g :: [inout] Ocean grid structure

  • gv :: [in] Vertical grid structure

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

  • slope_x :: [in] Zonal isoneutral slope

  • n2_u :: [in] Buoyancy (Brunt-Vaisala) frequency at u-points [T-2 ~> s-2]

  • slope_y :: [in] Meridional isoneutral slope

  • n2_v :: [in] Buoyancy (Brunt-Vaisala) frequency at v-points [T-2 ~> s-2]

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

  • cs :: Variable mixing coefficients

  • obc :: Open boundaries control structure.

Call to

mom_error_handler::mom_error mom_open_boundary::obc_none mom_diag_mediator::query_averaging_enabled

Called from

calc_slope_functions

subroutine mom_lateral_mixing_coeffs/calc_slope_functions_using_just_e(h, G, GV, US, CS, e, calculate_slopes, OBC)

The original calc_slope_function() that calculated slopes using interface positions only, not accounting for density variations.

Parameters
  • g :: [inout] Ocean grid structure

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

  • gv :: [in] Vertical grid structure

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

  • cs :: Variable mixing coefficients

  • e :: [in] Interface position [Z ~> m]

  • calculate_slopes :: [in] If true, calculate slopes internally otherwise use slopes stored in CS

  • obc :: Open boundaries control structure.

Call to

mom_error_handler::mom_error mom_open_boundary::obc_none

Called from

calc_slope_functions

subroutine mom_lateral_mixing_coeffs/calc_qg_leith_viscosity(CS, G, GV, US, h, k, div_xx_dx, div_xx_dy, vort_xy_dx, vort_xy_dy)

Calculates the Leith Laplacian and bi-harmonic viscosity coefficients.

Parameters
  • cs :: Variable mixing coefficients

  • g :: [in] Ocean grid structure

  • gv :: [in] The ocean’s vertical grid structure.

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

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

  • k :: [in] Layer for which to calculate vorticity magnitude

  • div_xx_dx :: [in] x-derivative of horizontal divergence (d/dx(du/dx + dv/dy)) [L-1 T-1 ~> m-1 s-1]

  • div_xx_dy :: [in] y-derivative of horizontal divergence (d/dy(du/dx + dv/dy)) [L-1 T-1 ~> m-1 s-1]

  • vort_xy_dx :: [inout] x-derivative of vertical vorticity (d/dx(dv/dx - du/dy)) [L-1 T-1 ~> m-1 s-1]

  • vort_xy_dy :: [inout] y-derivative of vertical vorticity (d/dy(dv/dx - du/dy)) [L-1 T-1 ~> m-1 s-1]

subroutine mom_lateral_mixing_coeffs/varmix_init(Time, G, GV, US, param_file, diag, CS)

Initializes the variables mixing coefficients container.

Parameters
  • time :: [in] Current model time

  • g :: [in] Ocean grid structure

  • gv :: [in] The ocean’s vertical grid structure

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

  • param_file :: [in] Parameter file handles

  • diag :: [inout] Diagnostics control structure

  • cs :: Variable mixing coefficients

Call to

mom_error_handler::mom_error