mom_stochastics module reference¶
Top-level module for the MOM6 ocean model in coupled mode.
Data Types¶
This control structure holds parameters for the MOM_stochastics module. |
Functions/Subroutines¶
This subroutine initializes the stochastics physics control structure. |
|
Advances the stochastic patterns one time step. |
|
Adds a stochastic increment (backscatter) to the input velocity field. |
|
Apply a 9-point smoothing filter twice to a pair of velocity components to reduce horizontal two-grid-point noise. |
Detailed Description¶
This file contains subroutines that implement some stochastic parameterizations in MOM6. SPPT perturbations of the tendencies of S and T are turned on using DO_SPPT=True. Stochastic perturbations in ePBL are turned on using PERT_EPBL=True. Stochastic kinetic energy backscatter (SKEB) via the Stochastic GM+E scheme is turned on using DO_SKEB=True. For all three schemes the spatial and temporal correlation structure of the associated random fields is controlled from the nam_stochy namelist used by the external stochastic_physics package, which is called by subroutines in this module.
The SKEB backscatter can be set in a variety of ways. If SKEB_USE_GM=True then SKEB_GM_COEF times the GM work rate will be added to the backscatter rate. (The vertical structure for this component of backscatter is the so-called EBT struct.) If SKEB_USE_FRICT=True then SKEB_FRICT_COEF times the work rate from lateral viscosity will be added to the backscatter rate. The code uses the total contribution from Laplacian and biharmonic viscosities as computed within the horizontal viscosity module. If neither SKEB_USE_GM nor SKEB_USE_FRICT is true, then the code computes the dissipation rate as if it came from a lateral harmonic viscosity with coefficient 1 (MKS units). The only thoroughly tested SKEB option at this point is SKEB_USE_GM.
The contributions to the backscatter rate are smoothed before use. One smoothing pass uses a 3x3 moving average with weights proportional to the h-cell areas. The number of smoothing passes is controlled by SKEB_NPASS.
A taper is applied to the SKEB velocity increments (equivalently to the SKEB stochastic forcing). The taper zeros out the increments near land cells. The width of this taper can be controlled using SKEB_TAPER_WIDTH.
Type Documentation¶
-
type
mom_stochastics/stochastic_cs¶ This control structure holds parameters for the MOM_stochastics module.
- Type fields:
%do_sppt[logical] :: If true, stochastically perturb the diabatic.%do_skeb[logical] :: If true, stochastically perturb the horizontal velocity.%skeb_use_gm[logical] :: If true, adds GM work to the amplitude of SKEBS.%skeb_use_frict[logical] :: If true, adds viscous dissipation rate to the amplitude of SKEBS.%pert_epbl[logical] :: If true, then randomly perturb the KE dissipation and genration terms.%id_sppt_wts[integer] :: Diagnostic id for SPPT.%id_skeb_wts[integer] :: Diagnostic id for SKEB.%id_skebu[integer] :: Diagnostic id for SKEB.%id_skebv[integer] :: Diagnostic id for SKEB.%id_diss[integer] :: Diagnostic id for SKEB.%skeb_npass[integer] :: number of passes of the 9-point smoother for the dissipation estimate%id_psi[integer] :: Diagnostic id for SPPT.%id_epbl1_wts[integer] :: Diagnostic id for epbl generation perturbation.%id_epbl2_wts[integer] :: Diagnostic id for epbl dissipation perturbation.%id_skeb_taperu[integer] :: Diagnostic id for u taper of SKEB velocity increment.%id_skeb_taperv[integer] :: Diagnostic id for v taper of SKEB velocity increment.%skeb_gm_coef[real] :: If skeb_use_gm is true, then skeb_gm_coef * GM_work is added to the dissipation rate used to set the amplitude of SKEBS [nondim].%skeb_frict_coef[real] :: If skeb_use_frict is true, then skeb_gm_coef * GM_work is added to the dissipation rate used to set the amplitude of SKEBS [nondim].%skeb_diss[real, dimension (:,:,:), allocatable] :: Dissipation rate used to set amplitude of SKEBS [L2 T-3 ~> m2 s-2] Index into this at h points.%answer_date[integer] :: The vintage of the order of arithmetic in the stochastics calculations. Values below 20250701 recover the answers from early in 2025, while higher values use expressions that have been refactored for rotational symmetry, including with FMAs enabled.%sppt_wts[real, dimension (:,:), allocatable] :: Random pattern for ocean SPPT tendencies with a number between 0 and 2 [nondim].%skeb_wts[real, dimension (:,:), allocatable] :: Random pattern of lengthscales for ocean SKEB in mks units [m] Note that SKEB_wts is set via external code in mks units.%epbl1_wts[real, dimension (:,:), allocatable] :: Random pattern for K.E. generation [nondim].%epbl2_wts[real, dimension (:,:), allocatable] :: Random pattern for K.E. dissipation [nondim].%time[type(time_type),pointer] :: Pointer to model time (needed for sponges)%diag[type( diag_ctrl ),pointer] :: A structure that is used to regulate the.%tapercu[real, dimension (:,:), allocatable] :: Taper applied to u component of stochastic velocity increment range [0,1], [nondim].%tapercv[real, dimension (:,:), allocatable] :: Taper applied to v component of stochastic velocity increment range [0,1], [nondim].%allocable[real]%real(nimem dimension [njmem)%isum_area_wts[real] :: One over the 3x3 sum of area_wt [L-2 ~> m-2].%area_wt[real] :: Masked h cell areas. [L2 ~> m2].
Function/Subroutine Documentation¶
-
subroutine
mom_stochastics/stochastics_init(dt, grid, GV, US, CS, param_file, diag, Time)¶ This subroutine initializes the stochastics physics control structure.
- Parameters:
dt :: [in] time step [T ~> s]
grid :: [in] horizontal grid information
gv :: [in] vertical grid structure
us :: [in] A dimensional unit scaling type
cs :: [inout] stochastic control structure
param_file :: [in] A structure to parse for run-time parameters
diag :: [inout] structure to regulate diagnostic output
time :: model time
- Call to:
mom_error_handler::calltree_entermom_error_handler::calltree_leavemom_diag_mediator::register_static_fieldsmooth_x9_uv- Called from:
-
subroutine
mom_stochastics/update_stochastics(CS)¶ Advances the stochastic patterns one time step.
- Parameters:
cs :: [inout] diabatic control structure
- Call to:
mom_error_handler::calltree_entermom_error_handler::calltree_leave- Called from:
-
subroutine
mom_stochastics/apply_skeb(grid, GV, US, CS, uc, vc, thickness, tv, dt, Time_end)¶ Adds a stochastic increment (backscatter) to the input velocity field.
- Parameters:
grid :: [in] ocean grid structure
gv :: [in] ocean vertical grid
us :: [in] A dimensional unit scaling type
cs :: [inout] stochastic control structure
uc :: [inout] zonal velocity [L T-1 ~> m s-1]
vc :: [inout] meridional velocity [L T-1 ~> m s-1]
thickness :: [in] thickness [H ~> m or kg m-2]
tv :: [in] points to thermodynamic fields
dt :: [in] time increment [T ~> s]
time_end :: [in] Time at the end of the interval
- Call to:
mom_error_handler::calltree_entermom_error_handler::calltree_leavemom_diag_mediator::disable_averagingmom_diag_mediator::enable_averages- Called from:
-
subroutine
mom_stochastics/smooth_x9_uv(G, field_u, field_v, zero_land)¶ Apply a 9-point smoothing filter twice to a pair of velocity components to reduce horizontal two-grid-point noise. Note that this subroutine does not conserve angular momentum, so don’t use it in situations where you need conservation. Also note that it assumes that the input fields have valid values in the first two halo points upon entry.
- Parameters:
g :: [in] Ocean grid
field_u :: [inout] u-point field to be smoothed [arbitrary]
field_v :: [inout] v-point field to be smoothed [arbitrary]
zero_land :: [in] If present and false, return the average of the surrounding ocean points when smoothing, otherwise use a value of 0 for land points and include them in the averages.
- Called from: