MOM6
|
Implements vertical viscosity (vertvisc)
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).
Data Types | |
type | vertvisc_cs |
Functions/Subroutines | |
subroutine, public | vertvisc (u, v, h, fluxes, visc, dt, OBC, ADp, CDp, G, GV, CS, taux_bot, tauy_bot) |
Perform a fully implicit vertical diffusion of momentum. Stress top and bottom boundary conditions are used. More... | |
subroutine, public | vertvisc_remnant (visc, visc_rem_u, visc_rem_v, dt, G, GV, 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. More... | |
subroutine, public | vertvisc_coef (u, v, h, fluxes, visc, dt, G, GV, 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(). More... | |
subroutine | find_coupling_coef (a, hvel, do_i, h_harm, bbl_thick, kv_bbl, z_i, h_ml, dt, j, G, GV, CS, visc, fluxes, work_on_u, OBC, shelf) |
Calculate the 'coupling coefficient' (a[k]) at the interfaces. If BOTTOMDRAGLAW is defined, the minimum of Hbbl and half the adjacent layer thicknesses are used to calculate a[k] near the bottom. More... | |
subroutine, public | vertvisc_limit_vel (u, v, h, ADp, CDp, fluxes, visc, dt, G, GV, 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. More... | |
subroutine, public | vertvisc_init (MIS, Time, G, GV, param_file, diag, ADp, dirs, ntrunc, CS) |
Initialise the vertical friction module. More... | |
subroutine, public | 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. More... | |
subroutine, public | vertvisc_end (CS) |
Clean up and deallocate the vertical friction module. More... | |
|
private |
Calculate the 'coupling coefficient' (a[k]) at the interfaces. If BOTTOMDRAGLAW is defined, the minimum of Hbbl and half the adjacent layer thicknesses are used to calculate a[k] near the bottom.
[in] | g | Ocean grid structure |
[in] | gv | Ocean vertical grid structure |
[out] | a | Coupling coefficient across interfaces, in m s-1 |
[in] | hvel | Thickness at velocity points, in H |
[in] | do_i | If true, determine coupling coefficient for a column |
[in] | h_harm | Harmonic mean of thicknesses around a velocity grid point, in H |
[in] | bbl_thick | Bottom boundary layer thickness, in H |
[in] | kv_bbl | Bottom boundary layer viscosity, in m2 s-1 |
[in] | z_i | Estimate of interface heights above the bottom, normalised by the bottom boundary layer thickness |
[out] | h_ml | Mixed layer depth, in H |
[in] | j | j-index to find coupling coefficient for |
[in] | dt | Time increment, in s |
cs | Vertical viscosity control structure | |
[in] | visc | Structure containing viscosities and bottom drag |
[in] | fluxes | Structure containing forcing fields |
[in] | work_on_u | If true, u-points are being calculated, otherwise v-points |
obc | Open boundary condition structure | |
[in] | shelf | If present and true, use a surface boundary condition appropriate for an ice shelf. |
Definition at line 981 of file MOM_vert_friction.F90.
References mom_open_boundary::obc_direction_n, mom_open_boundary::obc_direction_s, and mom_open_boundary::obc_direction_w.
Referenced by vertvisc_coef().
subroutine, public mom_vert_friction::updatecfltruncationvalue | ( | type(time_type), intent(in), target | Time, |
type(vertvisc_cs), pointer | CS, | ||
logical, intent(in), optional | 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.
[in] | time | Current model time |
cs | Vertical viscosity control structure | |
[in] | activate | Whether to record the value of Time as the beginning of the ramp period |
Definition at line 1637 of file MOM_vert_friction.F90.
References mom_error_handler::mom_error().
Referenced by mom_dynamics_split_rk2::step_mom_dyn_split_rk2().
subroutine, public mom_vert_friction::vertvisc | ( | real, dimension(szib_(g),szj_(g),szk_(gv)), intent(inout) | u, |
real, dimension(szi_(g),szjb_(g),szk_(gv)), intent(inout) | v, | ||
real, dimension(szi_(g),szj_(g),szk_(gv)), intent(in) | h, | ||
type(forcing), intent(in) | fluxes, | ||
type(vertvisc_type), intent(inout) | visc, | ||
real, intent(in) | dt, | ||
type(ocean_obc_type), pointer | OBC, | ||
type(accel_diag_ptrs), intent(inout) | ADp, | ||
type(cont_diag_ptrs), intent(inout) | CDp, | ||
type(ocean_grid_type), intent(in) | G, | ||
type(verticalgrid_type), intent(in) | GV, | ||
type(vertvisc_cs), pointer | CS, | ||
real, dimension(szib_(g),szj_(g)), intent(out), optional | taux_bot, | ||
real, dimension(szi_(g),szjb_(g)), intent(out), optional | tauy_bot | ||
) |
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.
[in] | g | Ocean grid structure |
[in] | gv | Ocean vertical grid structure |
[in,out] | u | Zonal velocity in m s-1 |
[in,out] | v | Meridional velocity in m s-1 |
[in] | h | Layer thickness in H |
[in] | fluxes | Structure containing forcing fields |
[in,out] | visc | Viscosities and bottom drag |
[in] | dt | Time increment in s |
obc | Open boundary condition structure | |
[in,out] | adp | Accelerations in the momentum equations for diagnostics |
[in,out] | cdp | Continuity equation terms |
cs | Vertical viscosity control structure | |
[out] | taux_bot | Zonal bottom stress from ocean to rock in Pa |
[out] | tauy_bot | Meridional bottom stress from ocean to rock in Pa |
Definition at line 140 of file MOM_vert_friction.F90.
References mom_error_handler::mom_error(), and vertvisc_limit_vel().
subroutine, public mom_vert_friction::vertvisc_coef | ( | real, dimension( g %isdb: g %iedb, g %jsd: g %jed, gv %ke), intent(in) | u, |
real, dimension( g %isd: g %ied, g %jsdb: g %jedb, gv %ke), intent(in) | v, | ||
real, dimension( g %isd: g %ied, g %jsd: g %jed, gv %ke), intent(in) | h, | ||
type(forcing), intent(in) | fluxes, | ||
type(vertvisc_type), intent(in) | visc, | ||
real, intent(in) | dt, | ||
type(ocean_grid_type), intent(in) | G, | ||
type(verticalgrid_type), intent(in) | GV, | ||
type(vertvisc_cs), pointer | CS, | ||
type(ocean_obc_type), pointer | 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().
[in] | g | Ocean grid structure |
[in] | gv | Ocean vertical grid structure |
[in] | u | Zonal velocity in m s-1 |
[in] | v | Meridional velocity in m s-1 |
[in] | h | Layer thickness in H |
[in] | fluxes | Structure containing forcing fields |
[in] | visc | Viscosities and bottom drag |
[in] | dt | Time increment in s |
cs | Vertical viscosity control structure | |
obc | Open boundary condition structure |
Definition at line 534 of file MOM_vert_friction.F90.
References find_coupling_coef(), mom_error_handler::mom_error(), mom_open_boundary::obc_direction_n, mom_open_boundary::obc_direction_s, and mom_open_boundary::obc_direction_w.
subroutine, public mom_vert_friction::vertvisc_end | ( | type(vertvisc_cs), pointer | CS | ) |
Clean up and deallocate the vertical friction module.
Definition at line 1675 of file MOM_vert_friction.F90.
subroutine, public mom_vert_friction::vertvisc_init | ( | type(ocean_internal_state), intent(in), target | MIS, |
type(time_type), intent(in), target | Time, | ||
type(ocean_grid_type), intent(in) | G, | ||
type(verticalgrid_type), intent(in) | GV, | ||
type(param_file_type), intent(in) | param_file, | ||
type(diag_ctrl), intent(inout), target | diag, | ||
type(accel_diag_ptrs), intent(inout) | ADp, | ||
type(directories), intent(in) | dirs, | ||
integer, intent(inout), target | ntrunc, | ||
type(vertvisc_cs), pointer | CS | ||
) |
Initialise the vertical friction module.
[in] | mis | "MOM Internal State", a set of pointers to the fields and accelerations that make up the ocean's physical state |
[in] | time | Current model time |
[in] | g | Ocean grid structure |
[in] | gv | Ocean vertical grid structure |
[in] | param_file | File to parse for parameters |
[in,out] | diag | Diagnostic control structure |
[in,out] | adp | Acceleration diagnostic pointers |
[in] | dirs | Relevant directory paths |
[in,out] | ntrunc | Number of velocity truncations |
cs | Vertical viscosity control structure |
Definition at line 1458 of file MOM_vert_friction.F90.
References mom_error_handler::mom_error().
subroutine, public mom_vert_friction::vertvisc_limit_vel | ( | real, dimension( g %isdb: g %iedb, g %jsd: g %jed, gv %ke), intent(inout) | u, |
real, dimension( g %isd: g %ied, g %jsdb: g %jedb, gv %ke), intent(inout) | v, | ||
real, dimension( g %isd: g %ied, g %jsd: g %jed, gv %ke), intent(in) | h, | ||
type(accel_diag_ptrs), intent(in) | ADp, | ||
type(cont_diag_ptrs), intent(in) | CDp, | ||
type(forcing), intent(in) | fluxes, | ||
type(vertvisc_type), intent(in) | visc, | ||
real, intent(in) | dt, | ||
type(ocean_grid_type), intent(in) | G, | ||
type(verticalgrid_type), intent(in) | GV, | ||
type(vertvisc_cs), pointer | 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.
[in] | g | Ocean grid structure |
[in] | gv | Ocean vertical grid structure |
[in,out] | u | Zonal velocity in m s-1 |
[in,out] | v | Meridional velocity in m s-1 |
[in] | h | Layer thickness in H |
[in] | adp | Acceleration diagnostic pointers |
[in] | cdp | Continuity diagnostic pointers |
[in] | fluxes | Forcing fields |
[in] | visc | Viscosities and bottom drag |
[in] | dt | Time increment in s |
cs | Vertical viscosity control structure |
Definition at line 1258 of file MOM_vert_friction.F90.
Referenced by vertvisc().
subroutine, public mom_vert_friction::vertvisc_remnant | ( | type(vertvisc_type), intent(in) | visc, |
real, dimension(szib_(g),szj_(g),szk_(gv)), intent(inout) | visc_rem_u, | ||
real, dimension( g %isd: g %ied, g %jsdb: g %jedb, gv %ke), intent(inout) | visc_rem_v, | ||
real, intent(in) | dt, | ||
type(ocean_grid_type), intent(in) | G, | ||
type(verticalgrid_type), intent(in) | GV, | ||
type(vertvisc_cs), pointer | 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.
[in] | g | Ocean grid structure |
[in] | gv | Ocean vertical grid structure |
[in] | visc | Viscosities and bottom drag |
[in,out] | visc_rem_u | Fraction of a time-step's worth of a barotopic acceleration that a layer experiences after viscosity is applied in the zonal direction |
[in,out] | visc_rem_v | Fraction of a time-step's worth of a barotopic acceleration that a layer experiences after viscosity is applied in the meridional direction |
[in] | dt | Time increment in s |
cs | Vertical viscosity control structure |
Definition at line 427 of file MOM_vert_friction.F90.
References mom_error_handler::mom_error().