mom_vert_friction module reference

Implements vertical viscosity (vertvisc)

More…

Data Types

vertvisc_cs

The control structure with parameters and memory for the MOM_vert_friction module.

Functions/Subroutines

vertfpmix()

Add nonlocal stress increments to u^n (uold) and v^n (vold) using ui and vi.

g_sig()

Returns the empirical shape-function given sigma.

find_coupling_coef_gl90()

Compute coupling coefficient associated with vertical viscosity parameterization as in Greatbatch and Lamb (1990), hereafter referred to as the GL90 vertical viscosity parameterization.

vertvisc()

Perform a fully implicit vertical diffusion of momentum.

vertvisc_remnant()

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.

vertvisc_coef()

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().

find_coupling_coef()

Calculate the ‘coupling coefficient’ (a_cpl) at the interfaces.

vertvisc_limit_vel()

Velocity components which exceed a threshold for physically reasonable values are truncated.

vertvisc_init()

Initialize the vertical friction module.

updatecfltruncationvalue()

Update the CFL truncation value as a function of time.

vertvisc_end()

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_du_dt_visc_gl90 [integer] :: Diagnostic identifiers.

  • % id_dv_dt_visc_gl90 [integer] :: Diagnostic identifiers.

  • % id_glwork [integer] :: Diagnostic identifiers.

  • % id_au_vv [integer] :: Diagnostic identifiers.

  • % id_av_vv [integer] :: Diagnostic identifiers.

  • % id_au_gl90_vv [integer] :: Diagnostic identifiers.

  • % id_av_gl90_vv [integer] :: Diagnostic identifiers.

  • % id_du_dt_str [integer] :: Diagnostic identifiers.

  • % id_dv_dt_str [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_fpw2x [integer] :: Diagnostic identifiers.

  • % id_taufp_u [integer] :: Diagnostic identifiers.

  • % id_taufp_v [integer] :: Diagnostic identifiers.

  • % id_fptau2s_u [integer] :: Diagnostic identifiers.

  • % id_fptau2s_v [integer] :: Diagnostic identifiers.

  • % id_fptau2w_u [integer] :: Diagnostic identifiers.

  • % id_fptau2w_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_kv_gl90_u [integer] :: Diagnostic identifiers.

  • % id_kv_gl90_v [integer] :: Diagnostic identifiers.

  • % id_h_du_dt_visc [integer] :: Diagnostic identifiers.

  • % id_h_dv_dt_visc [integer] :: Diagnostic identifiers.

  • % id_hf_du_dt_visc_2d [integer] :: Diagnostic identifiers.

  • % id_hf_dv_dt_visc_2d [integer] :: Diagnostic identifiers.

  • % id_h_du_dt_str [integer] :: Diagnostic identifiers.

  • % id_h_dv_dt_str [integer] :: Diagnostic identifiers.

  • % id_du_dt_str_visc_rem [integer] :: Diagnostic identifiers.

  • % id_dv_dt_str_visc_rem [integer] :: Diagnostic identifiers.

  • % initialized [logical] :: True if this control structure has been initialized.

  • % hmix [real] :: The mixed layer thickness [Z ~> m].

  • % hmix_stress [real] :: The mixed layer thickness over which the wind stress is applied with direct_stress [H ~> m or kg m-2].

  • % kvml_invz2 [real] :: The 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 [real] :: The interior vertical viscosity [H Z T-1 ~> m2 s-1 or Pa s].

  • % hbbl [real] :: The static bottom boundary layer thickness [Z ~> m].

  • % hbbl_gl90 [real] :: The static bottom boundary layer thickness used for GL90 [Z ~> m].

  • % kv_extra_bbl [real] :: An 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 [real] :: The von Karman constant as used for mixed layer viscosity [nondim].

  • % use_gl90_in_ssw [logical] :: If 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 [logical] :: If 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 [real] :: The scalar diffusivity used in the GL90 vertical viscosity scheme [L2 H Z-1 T-1 ~> m2 s-1 or Pa s].

  • % read_kappa_gl90 [logical] :: If true, read a file containing the spatially varying kappa_gl90.

  • % alpha_gl90 [real] :: Coefficient 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].

  • % 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 [T ~> s].

  • % cfl_truncs [real] :: The start value of CFL_trunc [nondim].

  • % cfl_trunce [real] :: The end/target value of CFL_trunc [nondim].

  • % 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 [H T-1 ~> m s-1 or Pa s m-1].

  • % real :: The u-drag coefficient associated with GL90 across an interface [H T-1 ~> m s-1 or Pa s m-1].

  • % real :: The effective layer thickness at u-points [H ~> m or kg m-2].

  • % real :: The v-drag coefficient across an interface [H T-1 ~> m s-1 or Pa s m-1].

  • % real :: The v-drag coefficient associated with GL90 across an interface [H T-1 ~> m s-1 or Pa s m-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 [H T-1 ~> m s-1 or Pa s m-1]. Retained to determine stress under shelves.

  • % a1_shelf_v [real(:,:),pointer] :: The 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 [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.

  • % 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 [nondim].

  • % direct_stress [logical] :: If 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 [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.

  • % fixed_lotw_ml [logical] :: If 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 [logical] :: If 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 [integer] :: The 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 [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.

  • % kappa_gl90_2d [real(:,:),allocatable] :: 2D kappa_gl90 at h-points [L2 H Z-1 T-1 ~> m2 s-1 or Pa s]

  • % pointaccel_csp [type( pointaccel_cs ),pointer] :: A 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, G, GV, US, CS, OBC)

Add nonlocal stress increments to u^n (uold) and v^n (vold) using ui and vi.

Parameters:
  • g :: [in] Ocean grid structure

  • gv :: [in] Ocean vertical grid structure

  • ui :: [inout] Zonal velocity after vertvisc [L T-1 ~> m s-1]

  • vi :: [inout] Meridional velocity after vertvisc [L T-1 ~> m s-1]

  • uold :: [inout] Old Zonal velocity [L T-1 ~> m s-1]

  • vold :: [inout] Old Meridional velocity [L T-1 ~> m s-1]

  • hbl_h :: [inout] boundary layer depth [H ~> m]

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

  • forces :: [in] A structure with the driving mechanical forces

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

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

  • cs :: Vertical viscosity control structure

  • obc :: Open boundary condition structure

Call to:

g_sig mom_set_visc::set_u_at_v mom_set_visc::set_v_at_u

Called from:

mom_dynamics_split_rk2::step_mom_dyn_split_rk2

function mom_vert_friction/g_sig(sigma) [real]

Returns the empirical shape-function given sigma.

Parameters:

sigma :: [in] non-dimensional normalized boundary layer depth [m]

Called from:

vertfpmix

subroutine mom_vert_friction/find_coupling_coef_gl90(a_cpl_gl90, hvel, do_i, z_i, j, 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 :: [in] Distance between interfaces at velocity points [Z ~> m]

  • do_i :: [in] If true, determine coupling coefficient for a column

  • z_i :: [in] Estimate of interface heights above the bottom, normalized by the GL90 bottom boundary layer thickness [nondim]

  • a_cpl_gl90 :: [inout] Coupling coefficient associated with GL90 across interfaces; is not included in a_cpl [H T-1 ~> m s-1 or Pa s m-1].

  • j :: [in] j-index to find coupling coefficient for

  • cs :: Vertical viscosity control structure

  • varmix :: [in] Variable mixing coefficients

  • work_on_u :: [in] If true, u-points are being calculated, otherwise they are v-points.

Called from:

vertvisc_coef

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 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 :: [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:

mom_error_handler::mom_error mom_diag_mediator::query_averaging_enabled vertvisc_limit_vel

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 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 :: [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:

mom_error_handler::mom_error

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 :: [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]

  • dz :: [in] Vertical distance across layers [Z ~> m]

  • forces :: [in] A structure with the driving mechanical forces

  • visc :: [in] Viscosities and bottom drag

  • tv :: [in] A structure containing pointers to any available thermodynamic fields.

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

  • cs :: Vertical viscosity control structure

  • obc :: Open boundary condition structure

  • varmix :: [in] Variable mixing coefficients

Call to:

find_coupling_coef find_coupling_coef_gl90 mom_error_handler::mom_error mom_open_boundary::obc_direction_n mom_open_boundary::obc_direction_s mom_open_boundary::obc_direction_w mom_diag_mediator::query_averaging_enabled

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, 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 :: [out] Coupling coefficient across interfaces [H T-1 ~> m s-1 or Pa s m-1]

  • hvel :: [in] Distance between interfaces at velocity points [Z ~> m]

  • 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 [Z ~> m]

  • 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 :: [in] Estimate of interface heights above the bottom,

  • h_ml :: [out] Mixed layer depth [Z ~> m]

  • 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

  • ustar_2d :: [in] The wind friction velocity, calculated using

  • tv :: [in] A structure containing pointers to any available thermodynamic fields.

  • 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:

vertvisc_coef

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:

vertvisc

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:

mom_error_handler::mom_error

Called from:

mom_dynamics_unsplit::initialize_dyn_unsplit mom_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 :: [in] Specify whether to record the value of Time as the beginning of the ramp period

Call to:

mom_error_handler::mom_error

Called from:

mom_dynamics_split_rk2::initialize_dyn_split_rk2 mom_dynamics_split_rk2b::initialize_dyn_split_rk2b mom_dynamics_split_rk2::step_mom_dyn_split_rk2 mom_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.