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¶
Perform a fully implicit vertical diffusion of momentum. |
|
Calculate the fraction of momentum originally in a layer that remains after a time-step of viscosity, and 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 and CSa_v) 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. |
|
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 maxvel to truncvel. 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
[integer] :: Diagnostic identifiers.%
id_dv_dt_visc
[integer] :: Diagnostic identifiers.%
id_au_vv
[integer] :: Diagnostic identifiers.%
id_av_vv
[integer] :: Diagnostic identifiers.%
id_h_u
[integer] :: Diagnostic identifiers.%
id_h_v
[integer] :: Diagnostic identifiers.%
id_hml_u
[integer] :: Diagnostic identifiers.%
id_hml_v
[integer] :: Diagnostic identifiers.%
id_taux_bot
[integer] :: Diagnostic identifiers.%
id_tauy_bot
[integer] :: Diagnostic identifiers.%
id_kv_slow
[integer] :: Diagnostic identifiers.%
id_kv_u
[integer] :: Diagnostic identifiers.%
id_kv_v
[integer] :: Diagnostic identifiers.%
id_hf_du_dt_visc_2d
[integer] :: Diagnostic identifiers.%
id_hf_dv_dt_visc_2d
[integer] :: Diagnostic identifiers.%
hmix
[real] :: The mixed layer thickness in thickness units [H ~> m or kg m-2].%
hmix_stress
[real] :: The mixed layer thickness over which the wind stress is applied with direct_stress [H ~> m or kg m-2].%
kvml
[real] :: The mixed layer vertical viscosity [Z2 T-1 ~> m2 s-1].%
kv
[real] :: The interior vertical viscosity [Z2 T-1 ~> m2 s-1].%
hbbl
[real] :: The static bottom boundary layer thickness [H ~> m or kg m-2].%
kvbbl
[real] :: The vertical viscosity in the bottom boundary layer [Z2 T-1 ~> m2 s-1].%
maxvel
[real] :: Velocity components greater than maxvel are truncated [L T-1 ~> m s-1].%
vel_underflow
[real] :: Velocity components smaller than vel_underflow are set to 0 [L T-1 ~> m s-1].%
cfl_based_trunc
[logical] :: If true, base truncations on CFL numbers, not absolute velocities.%
cfl_trunc
[real] :: Velocity components will be truncated when they are large enough that the corresponding CFL number exceeds this value, nondim.%
cfl_report
[real] :: The value of the CFL number that will cause the accelerations to be reported, nondim. CFL_report will often equal CFL_trunc.%
truncramptime
[real] :: The time-scale over which to ramp up the value of CFL_trunc from CFL_truncS to CFL_truncE.%
cfl_truncs
[real] :: The start value of CFL_trunc.%
cfl_trunce
[real] :: The end/target value of CFL_trunc.%
cflrampingisactivated
[logical] :: True if the ramping has been initialized.%
rampstarttime
[type(time_type)] :: The time at which the ramping of CFL_trunc starts.%
real
(* h_v [*, *) :: The u-drag coefficient across an interface [Z T-1 ~> m s-1].%
real
:: The effective layer thickness at u-points [H ~> m or kg m-2].%
real
:: The v-drag coefficient across an interface [Z T-1 ~> m s-1].%
real
:: The effective layer thickness at v-points [H ~> m or kg m-2].%
a1_shelf_u
[real(:,:),pointer] :: The u-momentum coupling coefficient under ice shelves [Z T-1 ~> m s-1]. Retained to determine stress under shelves.%
a1_shelf_v
[real(:,:),pointer] :: The v-momentum coupling coefficient under ice shelves [Z T-1 ~> m s-1]. Retained to determine stress under shelves.%
split
[logical] :: If true, use the split time stepping scheme.%
bottomdraglaw
[logical] :: If 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.%
channel_drag
[logical] :: If true, the drag is exerted directly on each layer according to what fraction of the bottom they overlie.%
harmonic_visc
[logical] :: If 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
[real] :: A 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.%
direct_stress
[logical] :: If true, the wind stress is distributed over the topmost Hmix_stress of fluid and KVML may be very small.%
dynamic_viscous_ml
[logical] :: If true, use the results from a dynamic calculation, perhaps based on a bulk Richardson number criterion, to determine the mixed layer thickness for viscosity.%
answers_2018
[logical] :: If true, use the order of arithmetic and expressions that recover the answers from the end of 2018. Otherwise, use expressions that do not use an arbitary and hard-coded maximum viscous coupling coefficient between layers.%
debug
[logical] :: If true, write verbose checksums for debugging purposes.%
nkml
[integer] :: The number of layers in the mixed layer.%
ntrunc
[integer,pointer] :: The 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
[logical] :: If 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),pointer] :: A structure that is used to regulate the timing of diagnostic output.%
pointaccel_csp
[type(pointaccel_cs),pointer] :: A pointer to the control structure for recording accelerations leading to velocity truncations.
Function/Subroutine Documentation¶
-
subroutine
mom_vert_friction/
vertvisc
(u, v, h, forces, visc, dt, OBC, ADp, CDp, G, GV, US, CS, taux_bot, tauy_bot, 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 Rayleight 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 :: [inout] Zonal velocity [L T-1 ~> m s-1]
v :: [inout] Meridional velocity [L T-1 ~> m s-1]
h :: [in] Layer thickness [H ~> m or kg m-2]
forces :: [in] A structure with the driving mechanical forces
visc :: [inout] Viscosities and bottom drag
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 :: [out] Zonal bottom stress from ocean to
tauy_bot :: [out] Meridional bottom stress from ocean to
waves :: Container for wave/Stokes information
- Call to
- Called from
mom_dynamics_unsplit::step_mom_dyn_unsplit
mom_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 after a time-step of viscosity, and 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 :: [in] Viscosities and bottom drag
visc_rem_u :: [inout] Fraction of a time-step’s worth of a
visc_rem_v :: [inout] Fraction of a time-step’s worth of a
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, forces, visc, dt, G, GV, US, CS, OBC)¶ Calculate the coupling coefficients (CSa_u and CSa_v) 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 :: [in] Zonal velocity [L T-1 ~> m s-1]
v :: [in] Meridional velocity [L T-1 ~> m s-1]
h :: [in] Layer thickness [H ~> m or kg m-2]
forces :: [in] A structure with the driving mechanical forces
visc :: [in] Viscosities and bottom drag
dt :: [in] Time increment [T ~> s]
cs :: Vertical viscosity control structure
obc :: Open boundary condition structure
- Call to
find_coupling_coef
mom_error_handler::mom_error
mom_open_boundary::obc_direction_n
mom_open_boundary::obc_direction_s
mom_open_boundary::obc_direction_w
- Called from
mom_dynamics_unsplit::step_mom_dyn_unsplit
mom_dynamics_unsplit_rk2::step_mom_dyn_unsplit_rk2
-
subroutine
mom_vert_friction/
find_coupling_coef
(a_cpl, hvel, do_i, h_harm, bbl_thick, kv_bbl, z_i, h_ml, dt, j, G, GV, US, CS, visc, forces, 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 :: [out] Coupling coefficient across interfaces [Z T-1 ~> m s-1].
hvel :: [in] Thickness at velocity points [H ~> m or kg m-2]
do_i :: [in] If true, determine coupling coefficient for a column
h_harm :: [in] Harmonic mean of thicknesses around a velocity
bbl_thick :: [in] Bottom boundary layer thickness [H ~> m or kg m-2]
kv_bbl :: [in] Bottom boundary layer viscosity [Z2 T-1 ~> m2 s-1].
z_i :: [in] Estimate of interface heights above the bottom,
h_ml :: [out] Mixed layer depth [H ~> m or kg m-2]
j :: [in] j-index to find coupling coefficient for
dt :: [in] Time increment [T ~> s]
cs :: Vertical viscosity control structure
visc :: [in] Structure containing viscosities and bottom drag
forces :: [in] A structure with the driving mechanical forces
work_on_u :: [in] If true, u-points are being calculated, otherwise they are v-points
obc :: Open boundary condition structure
shelf :: [in] If present and true, use a surface boundary condition appropriate for an ice shelf.
- Call to
mom_open_boundary::obc_direction_n
mom_open_boundary::obc_direction_s
mom_open_boundary::obc_direction_w
- 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. 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 :: [inout] Zonal velocity [L T-1 ~> m s-1]
v :: [inout] Meridional velocity [L T-1 ~> m s-1]
h :: [in] Layer thickness [H ~> m or kg m-2]
adp :: [in] Acceleration diagnostic pointers
cdp :: [in] Continuity diagnostic pointers
forces :: [in] A structure with the driving mechanical forces
visc :: [in] Viscosities and bottom drag
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)¶ 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 :: [in] File to parse for parameters
diag :: [inout] Diagnostic control structure
adp :: [inout] Acceleration diagnostic pointers
dirs :: [in] Relevant directory paths
ntrunc :: [inout] Number of velocity truncations
cs :: Vertical viscosity control structure
- Call to
-
subroutine
mom_vert_friction/
updatecfltruncationvalue
(Time, CS, 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
activate :: [in] Specifiy whether to record the value of Time as the beginning of the ramp period
- Call to
- Called from
-
subroutine
mom_vert_friction/
vertvisc_end
(CS)¶ Clean up and deallocate the vertical friction module.
- Parameters
cs :: Vertical viscosity control structure that will be deallocated in this subroutine.