mom_vert_friction module reference
Implements vertical viscosity (vertvisc)
Data Types
The control structure with parameters and memory for the MOM_vert_friction module. |
Functions/Subroutines
Add nonlocal stress increments to ui^n and vi^n. |
|
Compute coupling coefficient associated with vertical viscosity parameterization as in Greatbatch and Lamb (1990), hereafter referred to as the GL90 vertical viscosity parameterization. |
|
Perform a fully implicit vertical diffusion of momentum. |
|
Calculate the fraction of momentum originally in a layer that remains in the water column after a time-step of viscosity, equivalently the fraction of a time-step's worth of barotropic acceleration that a layer experiences after viscosity is applied. |
|
Calculate the coupling coefficients (CSa_u, CSa_v, CSa_u_gl90, CSa_v_gl90) and effective layer thicknesses (CSh_u and CSh_v) for later use in the applying the implicit vertical viscosity via |
|
Calculate the 'coupling coefficient' (a_cpl) at the interfaces. |
|
Velocity components which exceed a threshold for physically reasonable values are truncated, and the running sum of the number of trunctionas within the non-symmetric memory computational domain is incremented. |
|
Initialize the vertical friction module. |
|
Update the CFL truncation value as a function of time. |
|
Clean up and deallocate the vertical friction module. |
Detailed Description
The vertical diffusion of momentum is fully implicit. This is necessary to allow for vanishingly small layers. The coupling is based on the distance between the centers of adjacent layers, except where a layer is close to the bottom compared with a bottom boundary layer thickness when a bottom drag law is used. A stress top b.c. and a no slip bottom b.c. are used. There is no limit on the time step for vertvisc.
Near the bottom, the horizontal thickness interpolation scheme changes to an upwind biased estimate to control the effect of spurious Montgomery potential gradients at the bottom where nearly massless layers layers ride over the topography. Within a few boundary layer depths of the bottom, the harmonic mean thickness (i.e. (2 h+ h-) / (h+ + h-) ) is used if the velocity is from the thinner side and the arithmetic mean thickness (i.e. (h+ + h-)/2) is used if the velocity is from the thicker side. Both of these thickness estimates are second order accurate. Above this the arithmetic mean thickness is used.
In addition, vertvisc truncates any velocity component that exceeds a maximum CFL number to a fraction of this value. This basically keeps instabilities spatially localized. The number of times the velocity is truncated is reported each time the energies are saved, and if exceeds CSMaxtrunc the model will stop itself and change the time to a large value. This has proven very useful in (1) diagnosing model failures and (2) letting the model settle down to a meaningful integration from a poorly specified initial condition.
The same code is used for the two velocity components, by indirectly referencing the velocities and defining a handful of direction-specific defined variables.
Macros written all in capital letters are defined in MOM_memory.h.
A small fragment of the grid is shown below:
j+1 x ^ x ^ x At x: q
j+1 > o > o > At ^: v, frhatv, tauy
j x ^ x ^ x At >: u, frhatu, taux
j > o > o > At o: h
j-1 x ^ x ^ x
i-1 i i+1 At x & ^:
i i+1 At > & o:
The boundaries always run through q grid points (x).
Type Documentation
- type mom_vert_friction/vertvisc_cs
The control structure with parameters and memory for the MOM_vert_friction module.
- Type fields:
% id_du_dt_visc ::
integerDiagnostic identifiers.% id_dv_dt_visc ::
integerDiagnostic identifiers.% id_du_dt_visc_gl90 ::
integerDiagnostic identifiers.% id_dv_dt_visc_gl90 ::
integerDiagnostic identifiers.% id_glwork ::
integerDiagnostic identifiers.% id_au_vv ::
integerDiagnostic identifiers.% id_av_vv ::
integerDiagnostic identifiers.% id_au_gl90_vv ::
integerDiagnostic identifiers.% id_av_gl90_vv ::
integerDiagnostic identifiers.% id_du_dt_str ::
integerDiagnostic identifiers.% id_dv_dt_str ::
integerDiagnostic identifiers.% id_h_u ::
integerDiagnostic identifiers.% id_h_v ::
integerDiagnostic identifiers.% id_hml_u ::
integerDiagnostic identifiers.% id_hml_v ::
integerDiagnostic identifiers.% id_omega_w2x ::
integerDiagnostic identifiers.% id_fptau2s ::
integerDiagnostic identifiers.% id_fptau2w ::
integerDiagnostic identifiers.% id_ue_h ::
integerDiagnostic identifiers.% id_ve_h ::
integerDiagnostic identifiers.% id_ustk ::
integerDiagnostic identifiers.% id_vstk ::
integerDiagnostic identifiers.% id_ustk0 ::
integerDiagnostic identifiers.% id_vstk0 ::
integerDiagnostic identifiers.% id_uinc_h ::
integerDiagnostic identifiers.% id_vinc_h ::
integerDiagnostic identifiers.% id_taux_bot ::
integerDiagnostic identifiers.% id_tauy_bot ::
integerDiagnostic identifiers.% id_kv_slow ::
integerDiagnostic identifiers.% id_kv_u ::
integerDiagnostic identifiers.% id_kv_v ::
integerDiagnostic identifiers.% id_kv_gl90_u ::
integerDiagnostic identifiers.% id_kv_gl90_v ::
integerDiagnostic identifiers.% id_h_du_dt_visc ::
integerDiagnostic identifiers.% id_h_dv_dt_visc ::
integerDiagnostic identifiers.% id_hf_du_dt_visc_2d ::
integerDiagnostic identifiers.% id_hf_dv_dt_visc_2d ::
integerDiagnostic identifiers.% id_h_du_dt_str ::
integerDiagnostic identifiers.% id_h_dv_dt_str ::
integerDiagnostic identifiers.% id_du_dt_str_visc_rem ::
integerDiagnostic identifiers.% id_dv_dt_str_visc_rem ::
integerDiagnostic identifiers.% initialized ::
logicalTrue if this control structure has been initialized.% hmix ::
realThe mixed layer thickness [Z ~> m].% hmix_stress ::
realThe mixed layer thickness over which the wind stress is applied with direct_stress [H ~> m or kg m-2].% kvml_invz2 ::
realThe extra vertical viscosity scale in [H Z T-1 ~> m2 s-1 or Pa s] in a surface mixed layer with a characteristic thickness given by Hmix, and scaling proportional to (Hmix/z)^2, where z is the distance from the surface; this can get very large with thin layers.% kv ::
realThe interior vertical viscosity [H Z T-1 ~> m2 s-1 or Pa s].% hbbl ::
realThe static bottom boundary layer thickness [Z ~> m].% hbbl_gl90 ::
realThe static bottom boundary layer thickness used for GL90 [Z ~> m].% kv_extra_bbl ::
realAn extra vertical viscosity in the bottom boundary layer of thickness Hbbl when there is not a bottom drag law in use [H Z T-1 ~> m2 s-1 or Pa s].% vonkar ::
realThe von Karman constant as used for mixed layer viscosity [nondim].% use_gl90_in_ssw ::
logicalIf true, use the GL90 parameterization in stacked shallow water mode (SSW). The calculation of the GL90 viscosity coefficient uses the fact that in SSW we simply have 1/N^2 = h/g^prime, where g^prime is the reduced gravity. This identity does not generalize to non-SSW setups.% use_gl90_n2 ::
logicalIf true, use GL90 vertical viscosity coefficient that is depth-independent; this corresponds to a kappa_GM that scales as N^2 with depth.% kappa_gl90 ::
realThe scalar diffusivity used in the GL90 vertical viscosity scheme [L2 H Z-1 T-1 ~> m2 s-1 or Pa s].% read_kappa_gl90 ::
logicalIf true, read a file containing the spatially varying kappa_gl90.% alpha_gl90 ::
realCoefficient used to compute a depth-independent GL90 vertical viscosity via Kv_gl90 = alpha_gl90 * f^2. Note that the implied Kv_gl90 corresponds to a kappa_gl90 that scales as N^2 with depth. [H Z T ~> m2 s or kg s m-1].% vel_underflow ::
realVelocity components smaller than vel_underflow are set to 0 [L T-1 ~> m s-1].% cfl_trunc ::
realVelocity components will be truncated when they are large enough that the corresponding CFL number exceeds this value [nondim].% cfl_report ::
realThe value of the CFL number that will cause the accelerations to be reported [nondim]. CFL_report will often equal CFL_trunc.% truncramptime ::
realThe time-scale over which to ramp up the value of CFL_trunc from CFL_truncS to CFL_truncE [T ~> s].% cfl_truncs ::
realThe start value of CFL_trunc [nondim].% cfl_trunce ::
realThe end/target value of CFL_trunc [nondim].% cflrampingisactivated ::
logicalTrue if the ramping has been initialized.% rampstarttime ::
type(time_type)The time at which the ramping of CFL_trunc starts.% a_u ::
real, dimension(:, :, :), allocatableThe u-drag coefficient across an interface [H T-1 ~> m s-1 or Pa s m-1].% a_u_gl90 ::
real, dimension(:, :, :), allocatableThe u-drag coefficient associated with GL90 across an interface [H T-1 ~> m s-1 or Pa s m-1].% h_u ::
real, dimension(:, :, :), allocatableThe effective layer thickness at u-points [H ~> m or kg m-2].% a_v ::
real, dimension(:, :, :), allocatableThe v-drag coefficient across an interface [H T-1 ~> m s-1 or Pa s m-1].% a_v_gl90 ::
real, dimension(:, :, :), allocatableThe v-drag coefficient associated with GL90 across an interface [H T-1 ~> m s-1 or Pa s m-1].% h_v ::
real, dimension(:, :, :), allocatableThe effective layer thickness at v-points [H ~> m or kg m-2].% a1_shelf_u ::
real, dimension(:,:), pointerThe u-momentum coupling coefficient under ice shelves [H T-1 ~> m s-1 or Pa s m-1]. Retained to determine stress under shelves.% a1_shelf_v ::
real, dimension(:,:), pointerThe v-momentum coupling coefficient under ice shelves [H T-1 ~> m s-1 or Pa s m-1]. Retained to determine stress under shelves.% split ::
logicalIf true, use the split time stepping scheme.% bottomdraglaw ::
logicalIf true, the bottom stress is calculated with a drag law c_drag*|u|*u. The velocity magnitude may be an assumed value or it may be based on the actual velocity in the bottommost HBBL, depending on whether linear_drag is true.% harmonic_visc ::
logicalIf true, the harmonic mean thicknesses are used to calculate the viscous coupling between layers except near the bottom. Otherwise the arithmetic mean thickness is used except near the bottom.% harm_bl_val ::
realA scale to determine when water is in the boundary layers based solely on harmonic mean thicknesses for the purpose of determining the extent to which the thicknesses used in the viscosities are upwinded [nondim].% direct_stress ::
logicalIf true, the wind stress is distributed over the topmost Hmix_stress of fluid, and an added mixed layer viscosity or a physically based boundary layer turbulence parameterization is not needed for stability.% dynamic_viscous_ml ::
logicalIf true, use the results from a dynamic calculation, perhaps based on a bulk Richardson number criterion, to determine the mixed layer thickness for viscosity.% fixed_lotw_ml ::
logicalIf true, use a Law-of-the-wall prescription for the mixed layer viscosity within a boundary layer that is the lesser of Hmix and the total depth of the ocean in a column.% apply_lotw_floor ::
logicalIf true, use a Law-of-the-wall prescription to set a lower bound on the viscous coupling between layers within the surface boundary layer, based the distance of interfaces from the surface. This only acts when there are large changes in the thicknesses of successive layers or when the viscosity is set externally and the wind stress has subsequently increased.% answer_date ::
integerThe vintage of the order of arithmetic and expressions in the viscous calculations. Values below 20190101 recover the answers from the end of 2018, while higher values use expressions that do not use an arbitrary and hard-coded maximum viscous coupling coefficient between layers. In non-Boussinesq cases, values below 20230601 recover a form of the viscosity within the mixed layer that breaks up the magnitude of the wind stress with BULKMIXEDLAYER, DYNAMIC_VISCOUS_ML or FIXED_DEPTH_LOTW_ML, but not LOTW_VISCOUS_ML_FLOOR.% debug ::
logicalIf true, write verbose checksums for debugging purposes.% nkml ::
integerThe number of layers in the mixed layer.% ntrunc ::
integer, pointerThe number of times the velocity has been truncated since the last call to write_energy.% u_trunc_file ::
character(len=200)The complete path to a file in which a column of u-accelerations are written if velocity truncations occur.% v_trunc_file ::
character(len=200)The complete path to a file in which a column of v-accelerations are written if velocity truncations occur.% stokesmixing ::
logicalIf true, do Stokes drift mixing via the Lagrangian current (Eulerian plus Stokes drift). False by default and set via STOKES_MIXING_COMBINED.% diag ::
type(diag_ctrl), pointerA structure that is used to regulate the timing of diagnostic output.% kappa_gl90_2d ::
real, dimension(:,:), allocatable2D kappa_gl90 at h-points [L2 H Z-1 T-1 ~> m2 s-1 or Pa s]% pointaccel_csp ::
type(pointaccel_cs), pointerA pointer to the control structure for recording accelerations leading to velocity truncations.% pass_ke_uv ::
type(group_pass_type)A handle used for group halo passes.
Function/Subroutine Documentation
- subroutine mom_vert_friction/vertfpmix(ui, vi, uold, vold, hbl_h, h, forces, dt, lpost, Cemp_NL, G, GV, US, CS, OBC, Waves)
Add nonlocal stress increments to ui^n and vi^n.
- Parameters:
g :: [in] Ocean grid structure
gv :: [in] Ocean vertical grid structure
ui ::
ui[inout] Zonal velocity after vertvisc [L T-1 ~> m s-1]vi ::
vi[inout] Meridional velocity after vertvisc [L T-1 ~> m s-1]uold ::
uold[inout] Old Zonal velocity [L T-1 ~> m s-1]vold ::
vold[inout] Old Meridional velocity [L T-1 ~> m s-1]hbl_h ::
hbl_h[inout] boundary layer depth [H ~> m]h ::
h[in] Layer thicknesses [H ~> m or kg m-2]forces ::
forces[in] A structure with the driving mechanical forcesdt ::
dt[in] Time increment [T ~> s]cemp_nl :: [in] empirical coefficient of non-local momentum mixing [nondim]
lpost ::
lpost[in] Compute and make available FPMix diagnosticsus :: [in] A dimensional unit scaling type
cs :: Vertical viscosity control structure
obc :: Open boundary condition structure
waves :: Container for wave/Stokes information
- Called from:
- subroutine mom_vert_friction/find_coupling_coef_gl90(a_cpl_gl90, hvel, i, j, z_i, G, GV, CS, VarMix, work_on_u)
Compute coupling coefficient associated with vertical viscosity parameterization as in Greatbatch and Lamb (1990), hereafter referred to as the GL90 vertical viscosity parameterization. This vertical viscosity scheme redistributes momentum in the vertical, and is the equivalent of the Gent & McWilliams (1990) parameterization, but in a TWA (thickness-weighted averaged) set of equations. The vertical viscosity coefficient nu is computed from kappa_GM via thermal wind balance, and the following relation: nu = kappa_GM * f^2 / N^2. In the following subroutine kappa_GM is assumed either (a) constant or (b) horizontally varying. In both cases, (a) and (b), one can additionally impose an EBT structure in the vertical for kappa_GM. A third possible formulation of nu is depth-independent: nu = f^2 * alpha The latter formulation would be equivalent to a kappa_GM that varies as N^2 with depth. The vertical viscosity del_z ( nu del_z u) is applied to the momentum equation with stress-free boundary conditions at the top and bottom.
In SSW mode, we have 1/N^2 = h/g’. The coupling coefficient is therefore equal to a_cpl_gl90 = nu / h = kappa_GM * f^2 / g’ or a_cpl_gl90 = nu / h = f^2 * alpha / h
- Parameters:
g :: [in] Grid structure.
gv :: [in] Vertical grid structure.
hvel ::
hvel[in] Distance between interfaces at velocity points [Z ~> m]i ::
i[in] Column i-indexj ::
j[in] Column j-indexz_i ::
z_i[in] Estimate of interface heights above the bottom, normalized by the GL90 bottom boundary layer thickness [nondim]a_cpl_gl90 ::
a_cpl_gl90[out] Coupling coefficient associated with GL90 across interfaces; is not included in a_cpl [H T-1 ~> m s-1 or Pa s m-1].cs :: [in] Vertical viscosity control structure
varmix :: [in] Variable mixing coefficients
work_on_u ::
work_on_u[in] If true, u-points are being calculated, otherwise they are v-points.
- Called from:
- subroutine mom_vert_friction/vertvisc(u, v, h, forces, visc, dt, OBC, ADp, CDp, G, GV, US, CS, taux_bot, tauy_bot, fpmix, Waves)
Perform a fully implicit vertical diffusion of momentum. Stress top and bottom boundary conditions are used.
This is solving the tridiagonal system
\[\left(h_k + a_{k + 1/2} + a_{k - 1/2} + r_k\right) u_k^{n+1} = h_k u_k^n + a_{k + 1/2} u_{k+1}^{n+1} + a_{k - 1/2} u_{k-1}^{n+1}\]where \(a_{k + 1/2} = \Delta t \nu_{k + 1/2} / h_{k + 1/2}\) is the interfacial coupling thickness per time step , encompassing background viscosity as well as contributions from enhanced mixed and bottom layer viscosities. $r_k$ is a Rayleigh drag term due to channel drag. There is an additional stress term on the right-hand side if DIRECT_STRESS is true, applied to the surface layer.
- Parameters:
g :: [in] Ocean grid structure
gv :: [in] Ocean vertical grid structure
us :: [in] A dimensional unit scaling type
u ::
u[inout] Zonal velocity [L T-1 ~> m s-1]v ::
v[inout] Meridional velocity [L T-1 ~> m s-1]h ::
h[in] Layer thickness [H ~> m or kg m-2]forces ::
forces[in] A structure with the driving mechanical forcesvisc ::
visc[inout] Viscosities and bottom dragdt ::
dt[in] Time increment [T ~> s]obc :: Open boundary condition structure
adp :: [inout] Accelerations in the momentum equations for diagnostics
cdp :: [inout] Continuity equation terms
cs :: Vertical viscosity control structure
taux_bot ::
taux_bot[out] Zonal bottom stress from ocean totauy_bot ::
tauy_bot[out] Meridional bottom stress from ocean tofpmix ::
fpmix[in] fpmix along Eulerian shearwaves :: Container for wave/Stokes information
- Call to:
mom_error_handler::mom_errormom_diag_mediator::query_averaging_enabledvertvisc_limit_vel- Called from:
mom_dynamics_unsplit::step_mom_dyn_unsplitmom_dynamics_unsplit_rk2::step_mom_dyn_unsplit_rk2
- subroutine mom_vert_friction/vertvisc_remnant(visc, visc_rem_u, visc_rem_v, dt, G, GV, US, CS)
Calculate the fraction of momentum originally in a layer that remains in the water column after a time-step of viscosity, equivalently the fraction of a time-step’s worth of barotropic acceleration that a layer experiences after viscosity is applied.
- Parameters:
g :: [in] Ocean grid structure
gv :: [in] Ocean vertical grid structure
visc ::
visc[in] Viscosities and bottom dragvisc_rem_u ::
visc_rem_u[inout] Fraction of a time-step’s worth of avisc_rem_v ::
visc_rem_v[inout] Fraction of a time-step’s worth of adt ::
dt[in] Time increment [T ~> s]us :: [in] A dimensional unit scaling type
cs :: Vertical viscosity control structure
- Call to:
- subroutine mom_vert_friction/vertvisc_coef(u, v, h, dz, forces, visc, tv, dt, G, GV, US, CS, OBC, VarMix)
Calculate the coupling coefficients (CSa_u, CSa_v, CSa_u_gl90, CSa_v_gl90) and effective layer thicknesses (CSh_u and CSh_v) for later use in the applying the implicit vertical viscosity via
vertvisc(). .- Parameters:
g :: [in] Ocean grid structure
gv :: [in] Ocean vertical grid structure
us :: [in] A dimensional unit scaling type
u ::
u[in] Zonal velocity [L T-1 ~> m s-1]v ::
v[in] Meridional velocity [L T-1 ~> m s-1]h ::
h[in] Layer thickness [H ~> m or kg m-2]dz ::
dz[in] Vertical distance across layers [Z ~> m]forces ::
forces[in] A structure with the driving mechanical forcesvisc ::
visc[in] Viscosities and bottom dragtv ::
tv[in] A structure containing pointers to any available thermodynamic fields.dt ::
dt[in] Time increment [T ~> s]cs :: [inout] Vertical viscosity control structure
obc :: Open boundary condition structure
varmix :: [in] Variable mixing coefficients
- Call to:
find_coupling_coeffind_coupling_coef_gl90mom_error_handler::mom_errormom_diag_mediator::query_averaging_enabled- Called from:
mom_dynamics_unsplit::step_mom_dyn_unsplitmom_dynamics_unsplit_rk2::step_mom_dyn_unsplit_rk2
- subroutine mom_vert_friction/find_coupling_coef(a_cpl, hvel, i, j, h_harm, bbl_thick, kv_bbl, z_i, h_ml, dt, G, GV, US, CS, visc, Ustar_2d, tv, work_on_u, OBC, shelf)
Calculate the ‘coupling coefficient’ (a_cpl) at the interfaces. If BOTTOMDRAGLAW is defined, the minimum of Hbbl and half the adjacent layer thicknesses are used to calculate a_cpl near the bottom.
- Parameters:
g :: [in] Ocean grid structure
gv :: [in] Ocean vertical grid structure
us :: [in] A dimensional unit scaling type
a_cpl ::
a_cpl[out] Coupling coefficient across interfaces [H T-1 ~> m s-1 or Pa s m-1]hvel ::
hvel[in] Distance between interfaces at velocity points [Z ~> m]i ::
i[in] Column i-indexj ::
j[in] Column j-indexh_harm ::
h_harm[in] Harmonic mean of thicknesses around a velocitybbl_thick ::
bbl_thick[in] Bottom boundary layer thickness [Z ~> m]kv_bbl ::
kv_bbl[in] Bottom boundary layer viscosity, exclusive of any depth-dependent contributions from viscKv_shear [H Z T-1 ~> m2 s-1 or Pa s]z_i ::
z_i[in] Estimate of interface heights above the bottom,h_ml ::
h_ml[out] Mixed layer depth [Z ~> m]dt ::
dt[in] Time increment [T ~> s]cs :: [in] Vertical viscosity control structure
visc ::
visc[in] Structure containing viscosities and bottom dragustar_2d :: [in] The wind friction velocity, calculated using
tv ::
tv[in] A structure containing pointers to any available thermodynamic fields.work_on_u ::
work_on_u[in] If true, u-points are being calculated, otherwise they are v-pointsobc :: Open boundary condition structure
shelf ::
shelf[in] If present and true, use a surface boundary condition appropriate for an ice shelf.
- Called from:
- subroutine mom_vert_friction/vertvisc_limit_vel(u, v, h, ADp, CDp, forces, visc, dt, G, GV, US, CS)
Velocity components which exceed a threshold for physically reasonable values are truncated, and the running sum of the number of trunctionas within the non-symmetric memory computational domain is incremented. Optionally, any column with excessive velocities may be sent to a diagnostic reporting subroutine.
- Parameters:
g :: [in] Ocean grid structure
gv :: [in] Ocean vertical grid structure
us :: [in] A dimensional unit scaling type
u ::
u[inout] Zonal velocity [L T-1 ~> m s-1]v ::
v[inout] Meridional velocity [L T-1 ~> m s-1]h ::
h[in] Layer thickness [H ~> m or kg m-2]adp :: [in] Acceleration diagnostic pointers
cdp :: [in] Continuity diagnostic pointers
forces ::
forces[in] A structure with the driving mechanical forcesvisc ::
visc[in] Viscosities and bottom dragdt ::
dt[in] Time increment [T ~> s]cs :: Vertical viscosity control structure
- Called from:
- subroutine mom_vert_friction/vertvisc_init(MIS, Time, G, GV, US, param_file, diag, ADp, dirs, ntrunc, CS, fpmix)
Initialize the vertical friction module.
- Parameters:
mis :: [in] The “MOM Internal State”, a set of pointers
time :: [in] Current model time
g :: [in] Ocean grid structure
gv :: [in] Ocean vertical grid structure
us :: [in] A dimensional unit scaling type
param_file ::
param_file[in] File to parse for parametersdiag ::
diag[inout] Diagnostic control structureadp :: [inout] Acceleration diagnostic pointers
dirs ::
dirs[in] Relevant directory pathsntrunc ::
ntrunc[inout] Number of velocity truncationscs :: Vertical viscosity control structure
fpmix ::
fpmix[in] Nonlocal momentum mixing
- Call to:
- Called from:
mom_dynamics_unsplit::initialize_dyn_unsplitmom_dynamics_unsplit_rk2::initialize_dyn_unsplit_rk2
- subroutine mom_vert_friction/updatecfltruncationvalue(Time, CS, US, activate)
Update the CFL truncation value as a function of time. If called with the optional argument activate=.true., record the value of Time as the beginning of the ramp period.
- Parameters:
time :: [in] Current model time
cs :: Vertical viscosity control structure
us :: [in] A dimensional unit scaling type
activate ::
activate[in] Specify whether to record the value of Time as the beginning of the ramp period
- Call to:
- Called from:
mom_dynamics_split_rk2::initialize_dyn_split_rk2mom_dynamics_split_rk2b::initialize_dyn_split_rk2bmom_dynamics_split_rk2::step_mom_dyn_split_rk2mom_dynamics_split_rk2b::step_mom_dyn_split_rk2b
- subroutine mom_vert_friction/vertvisc_end(CS)
Clean up and deallocate the vertical friction module.
- Parameters:
cs :: [inout] Vertical viscosity control structure that will be deallocated in this subroutine.