mom_set_diffusivity module reference¶
Calculate vertical diffusivity from all mixing processes.
Data Types¶
This structure has memory for used in calculating diagnostics of diffusivity. |
|
This control structure contains parameters for MOM_set_diffusivity. |
Functions/Subroutines¶
Sets the interior vertical diffusion of scalars due to the following processes: |
|
Convert turbulent kinetic energy to diffusivity. |
|
Calculate Brunt-Vaisala frequency, N^2. |
|
This subroutine sets the additional diffusivities of temperature and salinity due to double diffusion, using the same functional form as is used in MOM4.1, and taken from an NCAR technical note (REF?) that updates what was in Large et al. |
|
This routine adds diffusion sustained by flow energy extracted by bottom drag. |
|
Calculates a BBL diffusivity use a Prandtl number 1 diffusivity with a law of the wall turbulent viscosity, up to a BBL height where the energy used for mixing has consumed the mechanical TKE input. |
|
This routine adds effects of mixed layer radiation to the layer diffusivities. |
|
This subroutine calculates several properties related to bottom boundary layer turbulence. |
|
Clear pointers and dealocate memory. |
Detailed Description¶
Calculate vertical diffusivity from all mixing processes.
Type Documentation¶
-
type
mom_set_diffusivity/
diffusivity_diags
¶ This structure has memory for used in calculating diagnostics of diffusivity.
- Type fields
%
n2_3d
[real(:,:,:),pointer, private] :: squared buoyancy frequency at interfaces [T-2 ~> s-2]%
kd_user
[real(:,:,:),pointer, private] :: user-added diffusivity at interfaces [Z2 T-1 ~> m2 s-1]%
kd_bbl
[real(:,:,:),pointer, private] :: BBL diffusivity at interfaces [Z2 T-1 ~> m2 s-1].%
kd_work
[real(:,:,:),pointer, private] :: layer integrated work by diapycnal mixing [R Z3 T-3 ~> W m-2]%
maxtke
[real(:,:,:),pointer, private] :: energy required to entrain to h_max [Z3 T-3 ~> m3 s-3]%
kd_bkgnd
[real(:,:,:),pointer, private] :: Background diffusivity at interfaces [Z2 T-1 ~> m2 s-1].%
kv_bkgnd
[real(:,:,:),pointer, private] :: Viscosity from ackground diffusivity at interfaces [Z2 T-1 ~> m2 s-1].%
kt_extra
[real(:,:,:),pointer, private] :: double diffusion diffusivity for temp [Z2 T-1 ~> m2 s-1].%
ks_extra
[real(:,:,:),pointer, private] :: double diffusion diffusivity for saln [Z2 T-1 ~> m2 s-1].%
drho_rat
[real(:,:,:),pointer, private] :: The density difference ratio used in double diffusion [nondim].%
tke_to_kd
[real(:,:,:),pointer, private] :: conversion rate (~1.0 / (G_Earth + dRho_lay)) between TKE dissipated within a layer and Kd in that layer [Z2 T-1 / Z3 T-3 = T2 Z-1 ~> s2 m-1]
-
type
mom_set_diffusivity/
set_diffusivity_cs
¶ This control structure contains parameters for MOM_set_diffusivity.
- Type fields
%
id_maxtke
[integer] :: Diagnostic IDs.%
id_tke_to_kd
[integer] :: Diagnostic IDs.%
id_kd_user
[integer] :: Diagnostic IDs.%
id_kd_layer
[integer] :: Diagnostic IDs.%
id_kd_bbl
[integer] :: Diagnostic IDs.%
id_n2
[integer] :: Diagnostic IDs.%
id_kd_work
[integer] :: Diagnostic IDs.%
id_kt_extra
[integer] :: Diagnostic IDs.%
id_ks_extra
[integer] :: Diagnostic IDs.%
id_r_rho
[integer] :: Diagnostic IDs.%
id_kd_bkgnd
[integer] :: Diagnostic IDs.%
id_kv_bkgnd
[integer] :: Diagnostic IDs.%
debug
[logical] :: If true, write verbose checksums for debugging.%
bulkmixedlayer
[logical] :: If true, a refined bulk mixed layer is used with GVnk_rho_varies variable density mixed & buffer layers.%
fluxri_max
[real] :: The flux Richardson number where the stratification is large enough that N2 > omega2 [nondim]. The full expression for the Flux Richardson number is usually FLUX_RI_MAX*N2/(N2+OMEGA2). The default is 0.2.%
bottomdraglaw
[logical] :: If true, the bottom stress is calculated with a drag law c_drag*|u|*u.%
bbl_mixing_as_max
[logical] :: If true, take the maximum of the diffusivity from the BBL mixing and the other diffusivities. Otherwise, diffusivities from the BBL_mixing is added.%
use_lotw_bbl_diffusivity
[logical] :: If true, use simpler/less precise, BBL diffusivity.%
lotw_bbl_use_omega
[logical] :: If true, use simpler/less precise, BBL diffusivity.%
bbl_effic
[real] :: efficiency with which the energy extracted by bottom drag drives BBL diffusion [nondim]%
cdrag
[real] :: quadratic drag coefficient [nondim]%
imax_decay
[real] :: inverse of a maximum decay scale for bottom-drag driven turbulence [Z-1 ~> m-1].%
kv
[real] :: The interior vertical viscosity [Z2 T-1 ~> m2 s-1].%
kd
[real] :: interior diapycnal diffusivity [Z2 T-1 ~> m2 s-1].%
kd_min
[real] :: minimum diapycnal diffusivity [Z2 T-1 ~> m2 s-1].%
kd_max
[real] :: maximum increment for diapycnal diffusivity [Z2 T-1 ~> m2 s-1]. Set to a negative value to have no limit.%
kd_add
[real] :: uniform diffusivity added everywhere without filtering or scaling [Z2 T-1 ~> m2 s-1].%
kd_smooth
[real] :: Vertical diffusivity used to interpolate more sensible values of T & S into thin layers [Z2 T-1 ~> m2 s-1].%
diag
[type(diag_ctrl),pointer] :: structure to regulate diagnostic output timing%
limit_dissipation
[logical] :: If enabled, dissipation is limited to be larger than the following:%
dissip_min
[real] :: Minimum dissipation [R Z2 T-3 ~> W m-3].%
dissip_n0
[real] :: Coefficient a in minimum dissipation = a+b*N [R Z2 T-3 ~> W m-3].%
dissip_n1
[real] :: Coefficient b in minimum dissipation = a+b*N [R Z2 T-2 ~> J m-3].%
dissip_n2
[real] :: Coefficient c in minimum dissipation = c*N2 [R Z2 T-1 ~> J s m-3].%
dissip_kd_min
[real] :: Minimum Kd [Z2 T-1 ~> m2 s-1], with dissipation Rho0*Kd_min*N^2.%
omega
[real] :: Earth’s rotation frequency [T-1 ~> s-1].%
ml_radiation
[logical] :: allow a fraction of TKE available from wind work to penetrate below mixed layer base with a vertical decay scale determined by the minimum of (1) The depth of the mixed layer, or (2) An Ekman length scale. Energy available to drive mixing below the mixed layer is given by E = ML_RAD_COEFF*MSTAR*USTAR**3. Optionally, if ML_rad_TKE_decay is true, this is further reduced by a factor of exp(-h_ML*Idecay_len_TkE), where Idecay_len_TKE is calculated the same way as in the mixed layer code. The diapycnal diffusivity is KD(k) = E/(N2(k)+OMEGA2), where N2 is the squared buoyancy frequency [T-2 ~> s-2] and OMEGA2 is the rotation rate of the earth squared.%
ml_rad_kd_max
[real] :: Maximum diapycnal diffusivity due to turbulence radiated from the base of the mixed layer [Z2 T-1 ~> m2 s-1].%
ml_rad_efold_coeff
[real] :: non-dim coefficient to scale penetration depth%
ml_rad_coeff
[real] :: coefficient, which scales MSTAR*USTAR^3 to obtain energy available for mixing below mixed layer base [nondim]%
ml_rad_bug
[logical] :: If true use code with a bug that reduces the energy available in the transition layer by a factor of the inverse of the energy deposition lenthscale (in m).%
ml_rad_tke_decay
[logical] :: If true, apply same exponential decay to ML_rad as applied to the other surface sources of TKE in the mixed layer code.%
ustar_min
[real] :: A minimum value of ustar to avoid numerical problems [Z T-1 ~> m s-1]. If the value is small enough, this parameter should not affect the solution.%
tke_decay
[real] :: ratio of natural Ekman depth to TKE decay scale [nondim]%
mstar
[real] :: ratio of friction velocity cubed to TKE input to the mixed layer [nondim]%
ml_use_omega
[logical] :: If true, use absolute rotation rate instead of the vertical component of rotation when setting the decay scale for mixed layer turbulence.%
ml_omega_frac
[real] :: When setting the decay scale for turbulence, use this fraction of the absolute rotation rate blended with the local value of f, as f^2 ~= (1-of)*f^2 + of*4*omega^2.%
user_change_diff
[logical] :: If true, call user-defined code to change diffusivity.%
usekappashear
[logical] :: If true, use the kappa_shear module to find the shear-driven diapycnal diffusivity.%
vertex_shear
[logical] :: If true, do the calculations of the shear-driven mixing at the cell vertices (i.e., the vorticity points).%
use_cvmix_shear
[logical] :: If true, use one of the CVMix modules to find shear-driven diapycnal diffusivity.%
double_diffusion
[logical] :: If true, enable double-diffusive mixing using an old method.%
use_cvmix_ddiff
[logical] :: If true, enable double-diffusive mixing via CVMix.%
use_tidal_mixing
[logical] :: If true, activate tidal mixing diffusivity.%
simple_tke_to_kd
[logical] :: If true, uses a simple estimate of Kd/TKE that does not rely on a layer-formulation.%
max_rrho_salt_fingers
[real] :: max density ratio for salt fingering%
max_salt_diff_salt_fingers
[real] :: max salt diffusivity for salt fingers [Z2 T-1 ~> m2 s-1]%
kv_molecular
[real] :: molecular visc for double diff convect [Z2 T-1 ~> m2 s-1]%
answers_2018
[logical] :: If true, use the order of arithmetic and expressions that recover the answers from the end of 2018. Otherwise, use updated and more robust forms of the same expressions.%
inputdir
[character (len=200)] :: The directory in which input files are found.%
user_change_diff_csp
[type(user_change_diff_cs),pointer] :: Control structure for a child module.%
kappashear_csp
[type(kappa_shear_cs),pointer] :: Control structure for a child module.%
cvmix_shear_csp
[type(cvmix_shear_cs),pointer] :: Control structure for a child module.%
cvmix_ddiff_csp
[type(cvmix_ddiff_cs),pointer] :: Control structure for a child module.%
bkgnd_mixing_csp
[type(bkgnd_mixing_cs),pointer] :: Control structure for a child module.%
int_tide_csp
[type(int_tide_cs),pointer] :: Control structure for a child module.%
tidal_mixing_csp
[type(tidal_mixing_cs),pointer] :: Control structure for a child module.
Function/Subroutine Documentation¶
-
subroutine
mom_set_diffusivity/
set_diffusivity
(u, v, h, u_h, v_h, tv, fluxes, optics, visc, dt, G, GV, US, CS, Kd_lay, Kd_int)¶ Sets the interior vertical diffusion of scalars due to the following processes:
Shear-driven mixing: two options, Jackson et at. and KPP interior;
Background mixing via CVMix (Bryan-Lewis profile) or the scheme described by Harrison & Hallberg, JPO 2008;
Double-diffusion, old method and new method via CVMix;
Tidal mixing: many options available, see
MOM_tidal_mixing.F90
; In addition, this subroutine has the option to set the interior vertical viscosity associated with processes 1,2 and 4 listed above, which is stored in viscKv_slow. Vertical viscosity due to shear-driven mixing is passed via viscKv_shear
- Parameters
g :: [in] The ocean’s grid structure.
gv :: [in] The ocean’s vertical grid structure.
us :: [in] A dimensional unit scaling type
u :: [in] The zonal velocity [L T-1 ~> m s-1].
v :: [in] The meridional velocity [L T-1 ~> m s-1].
h :: [in] Layer thicknesses [H ~> m or kg m-2].
u_h :: [in] Zonal velocity interpolated to h points [L T-1 ~> m s-1].
v_h :: [in] Meridional velocity interpolated to h points [L T-1 ~> m s-1].
tv :: [inout] Structure with pointers to thermodynamic fields. Out is for tvTempxPmE.
fluxes :: [in] A structure of thermodynamic surface fluxes
optics :: A structure describing the optical properties of the ocean.
visc :: [inout] Structure containing vertical viscosities, bottom boundary layer properies, and related fields.
dt :: [in] Time increment [T ~> s].
cs :: Module control structure.
kd_lay :: [out] Diapycnal diffusivity of each layer [Z2 T-1 ~> m2 s-1].
kd_int :: [out] Diapycnal diffusivity at each interface [Z2 T-1 ~> m2 s-1].
- Call to
add_drag_diffusivity
add_lotw_bbl_diffusivity
add_mlrad_diffusivity
mom_kappa_shear::calc_kappa_shear_vertex
mom_error_handler::calltree_enter
mom_error_handler::calltree_leave
mom_error_handler::calltree_waypoint
mom_cvmix_ddiff::compute_ddiff_coeffs
double_diffusion
find_n2
find_tke_to_kd
mom_full_convection::full_convection
id_clock_cvmix_ddiff
id_clock_kappashear
mom_isopycnal_slopes::vert_fill_ts
-
subroutine
mom_set_diffusivity/
find_tke_to_kd
(h, tv, dRho_int, N2_lay, j, dt, G, GV, US, CS, TKE_to_Kd, maxTKE, kb)¶ Convert turbulent kinetic energy to diffusivity.
- Parameters
g :: [in] The ocean’s grid structure
gv :: [in] The ocean’s vertical grid structure
us :: [in] A dimensional unit scaling type
h :: [in] Layer thicknesses [H ~> m or kg m-2]
tv :: [in] Structure containing pointers to any available thermodynamic fields.
drho_int :: [in] Change in locally referenced potential density across each interface [R ~> kg m-3].
n2_lay :: [in] The squared buoyancy frequency of the layers [T-2 ~> s-2].
j :: [in] j-index of row to work on
dt :: [in] Time increment [T ~> s].
cs :: Diffusivity control structure
tke_to_kd :: [out] The conversion rate between the TKE dissipated within a layer and the diapycnal diffusivity witin that layer, usually (~Rho_0 / (G_Earth * dRho_lay)) [Z2 T-1 / Z3 T-3 = T2 Z-1 ~> s2 m-1]
maxtke :: [out] The energy required to for a layer to entrain to its maximum realizable thickness [Z3 T-3 ~> m3 s-3]
kb :: [out] Index of lightest layer denser than the buffer layer, or -1 without a bulk mixed layer.
- Call to
- Called from
-
subroutine
mom_set_diffusivity/
find_n2
(h, tv, T_f, S_f, fluxes, j, G, GV, US, CS, dRho_int, N2_lay, N2_int, N2_bot)¶ Calculate Brunt-Vaisala frequency, N^2.
- Parameters
g :: [in] The ocean’s grid structure
gv :: [in] The ocean’s vertical grid structure
us :: [in] A dimensional unit scaling type
h :: [in] Layer thicknesses [H ~> m or kg m-2]
tv :: [in] Structure containing pointers to any available thermodynamic fields.
t_f :: [in] layer temperature with the values in massless layers
s_f :: [in] Layer salinities with values in massless
fluxes :: [in] A structure of thermodynamic surface fluxes
j :: [in] j-index of row to work on
cs :: Diffusivity control structure
drho_int :: [out] Change in locally referenced potential density
n2_int :: [out] The squared buoyancy frequency at the interfaces [T-2 ~> s-2].
n2_lay :: [out] The squared buoyancy frequency of the layers [T-2 ~> s-2].
n2_bot :: [out] The near-bottom squared buoyancy frequency [T-2 ~> s-2].
- Call to
- Called from
-
subroutine
mom_set_diffusivity/
double_diffusion
(tv, h, T_f, S_f, j, G, GV, US, CS, Kd_T_dd, Kd_S_dd)¶ This subroutine sets the additional diffusivities of temperature and salinity due to double diffusion, using the same functional form as is used in MOM4.1, and taken from an NCAR technical note (REF?) that updates what was in Large et al. (1994). All the coefficients here should probably be made run-time variables rather than hard-coded constants.
Todo
Find reference for NCAR tech note above.
- Parameters
g :: [in] The ocean’s grid structure.
gv :: [in] The ocean’s vertical grid structure.
us :: [in] A dimensional unit scaling type
tv :: [in] Structure containing pointers to any available thermodynamic fields; absent fields have NULL ptrs.
h :: [in] Layer thicknesses [H ~> m or kg m-2].
t_f :: [in] layer temperatures with the values in massless layers
s_f :: [in] Layer salinities with values in massless
j :: [in] Meridional index upon which to work.
cs :: Module control structure.
kd_t_dd :: [out] Interface double diffusion diapycnal
kd_s_dd :: [out] Interface double diffusion diapycnal
- Call to
- Called from
-
subroutine
mom_set_diffusivity/
add_drag_diffusivity
(h, u, v, tv, fluxes, visc, j, TKE_to_Kd, maxTKE, kb, G, GV, US, CS, Kd_lay, Kd_int, Kd_BBL)¶ This routine adds diffusion sustained by flow energy extracted by bottom drag.
- Parameters
g :: [in] The ocean’s grid structure
gv :: [in] The ocean’s vertical grid structure
us :: [in] A dimensional unit scaling type
u :: [in] The zonal velocity [L T-1 ~> m s-1]
v :: [in] The meridional velocity [L T-1 ~> m s-1]
h :: [in] Layer thicknesses [H ~> m or kg m-2]
tv :: [in] Structure containing pointers to any available thermodynamic fields.
fluxes :: [in] A structure of thermodynamic surface fluxes
visc :: [in] Structure containing vertical viscosities, bottom boundary layer properies, and related fields
j :: [in] j-index of row to work on
tke_to_kd :: [in] The conversion rate between the TKE TKE dissipated within a layer and the diapycnal diffusivity witin that layer, usually (~Rho_0 / (G_Earth * dRho_lay)) [Z2 T-1 / Z3 T-3 = T2 Z-1 ~> s2 m-1]
maxtke :: [in] The energy required to for a layer to entrain to its maximum-realizable thickness [Z3 T-3 ~> m3 s-3]
kb :: [in] Index of lightest layer denser than the buffer layer, or -1 without a bulk mixed layer
cs :: Diffusivity control structure
kd_lay :: [inout] The diapycnal diffusivity in layers, [Z2 T-1 ~> m2 s-1].
kd_int :: [inout] The diapycnal diffusivity at interfaces,
kd_bbl :: Interface BBL diffusivity [Z2 T-1 ~> m2 s-1].
- Called from
-
subroutine
mom_set_diffusivity/
add_lotw_bbl_diffusivity
(h, u, v, tv, fluxes, visc, j, N2_int, G, GV, US, CS, Kd_BBL, Kd_lay, Kd_int)¶ Calculates a BBL diffusivity use a Prandtl number 1 diffusivity with a law of the wall turbulent viscosity, up to a BBL height where the energy used for mixing has consumed the mechanical TKE input.
- Parameters
g :: [in] Grid structure
gv :: [in] Vertical grid structure
us :: [in] A dimensional unit scaling type
u :: [in] u component of flow [L T-1 ~> m s-1]
v :: [in] v component of flow [L T-1 ~> m s-1]
h :: [in] Layer thickness [H ~> m or kg m-2]
tv :: [in] Structure containing pointers to any available thermodynamic fields.
fluxes :: [in] Surface fluxes structure
visc :: [in] Structure containing vertical viscosities, bottom boundary layer properies, and related fields.
j :: [in] j-index of row to work on
n2_int :: [in] Square of Brunt-Vaisala at interfaces [T-2 ~> s-2]
cs :: Diffusivity control structure
kd_bbl :: Interface BBL diffusivity [Z2 T-1 ~> m2 s-1]
kd_lay :: [inout] Layer net diffusivity [Z2 T-1 ~> m2 s-1]
kd_int :: [inout] Interface net diffusivity [Z2 T-1 ~> m2 s-1]
- Called from
-
subroutine
mom_set_diffusivity/
add_mlrad_diffusivity
(h, fluxes, j, G, GV, US, CS, TKE_to_Kd, Kd_lay, Kd_int)¶ This routine adds effects of mixed layer radiation to the layer diffusivities.
- Parameters
g :: [in] The ocean’s grid structure
gv :: [in] The ocean’s vertical grid structure
us :: [in] A dimensional unit scaling type
h :: [in] Layer thicknesses [H ~> m or kg m-2]
fluxes :: [in] Surface fluxes structure
j :: [in] The j-index to work on
cs :: Diffusivity control structure
tke_to_kd :: [in] The conversion rate between the TKE TKE dissipated within a layer and the diapycnal diffusivity witin that layer, usually (~Rho_0 / (G_Earth * dRho_lay)) [Z2 T-1 / Z3 T-3 = T2 Z-1 ~> s2 m-1]
kd_lay :: [inout] The diapycnal diffusivity in layers [Z2 T-1 ~> m2 s-1].
kd_int :: [inout] The diapycnal diffusivity at interfaces
- Called from
-
subroutine
mom_set_diffusivity/
set_bbl_tke
(u, v, h, fluxes, visc, G, GV, US, CS, OBC)¶ This subroutine calculates several properties related to bottom boundary layer turbulence.
- Parameters
g :: [in] The ocean’s grid structure
gv :: [in] The ocean’s vertical grid structure
us :: [in] A dimensional unit scaling type
u :: [in] The zonal velocity [L T-1 ~> m s-1]
v :: [in] The meridional velocity [L T-1 ~> m s-1]
h :: [in] Layer thicknesses [H ~> m or kg m-2]
fluxes :: [in] A structure of thermodynamic surface fluxes
visc :: [in] Structure containing vertical viscosities, bottom boundary layer properies, and related fields.
cs :: Diffusivity control structure
obc :: Open boundaries control structure.
- Call to
mom_open_boundary::obc_direction_e
mom_open_boundary::obc_direction_n
-
subroutine
mom_set_diffusivity/
set_density_ratios
(h, tv, kb, G, GV, US, CS, j, ds_dsp1, rho_0)¶ - Parameters
g :: [in] The ocean’s grid structure.
gv :: [in] The ocean’s vertical grid structure.
h :: [in] Layer thicknesses [H ~> m or kg m-2].
tv :: [in] Structure containing pointers to any available thermodynamic fields; absent fields have NULL ptrs.
kb :: [in] Index of lightest layer denser than the buffer layer, or -1 without a bulk mixed layer.
us :: [in] A dimensional unit scaling type
cs :: Control structure returned by previous call to diabatic_entrain_init.
j :: [in] Meridional index upon which to work.
ds_dsp1 :: [out] Coordinate variable (sigma-2) difference across an interface divided by the difference across the interface below it [nondim]
rho_0 :: [in] Layer potential densities relative to
- Call to
- Called from
-
subroutine
mom_set_diffusivity/
set_diffusivity_init
(Time, G, GV, US, param_file, diag, CS, int_tide_CSp, halo_TS)¶ - Parameters
time :: [in] The current model time
g :: [inout] The ocean’s grid structure.
gv :: [in] The ocean’s vertical grid structure.
us :: [in] A dimensional unit scaling type
param_file :: [in] A structure to parse for run-time parameters.
diag :: [inout] A structure used to regulate diagnostic output.
cs :: pointer set to point to the module control structure.
int_tide_csp :: A pointer to the internal tides control structure
halo_ts :: [out] The halo size of tracer points that must be valid for the calculations in set_diffusivity.
- Call to
id_clock_cvmix_ddiff
id_clock_kappashear
mom_kappa_shear::kappa_shear_at_vertex
mom_diag_mediator::register_diag_field
mom_tidal_mixing::tidal_mixing_init
-
subroutine
mom_set_diffusivity/
set_diffusivity_end
(CS)¶ Clear pointers and dealocate memory.
- Parameters
cs :: Control structure for this module
- Call to
mom_bkgnd_mixing::bkgnd_mixing_end
mom_cvmix_shear::cvmix_shear_end
mom_tidal_mixing::tidal_mixing_end
user_change_diffusivity::user_change_diff_end