mom_dynamics_split_rk2b module reference
Time step the adiabatic dynamic core of MOM using RK2 method with greater use of the time-filtered velocities and less inheritance of tedencies from the previous step in the predictor step than in the original MOM_dyanmics_split_RK2.
Data Types
MOM_dynamics_split_RK2b module control structure. |
Functions/Subroutines
RK2 splitting for time stepping MOM adiabatic dynamics. |
|
This subroutine sets up any auxiliary restart variables that are specific to the split-explicit time stepping scheme. |
|
This subroutine does remapping for the auxiliary restart variables that are used with the split RK2 time stepping scheme. |
|
This subroutine initializes all of the variables that are used by this dynamic core, including diagnostics and the cpu clocks. |
|
Close the dyn_split_RK2b module. |
Detailed Description
This file time steps the adiabatic dynamic core by splitting between baroclinic and barotropic modes. It uses a pseudo-second order Runge-Kutta time stepping scheme for the baroclinic momentum equation and a forward-backward coupling between the baroclinic momentum and continuity equations. This split time-stepping scheme is described in detail in Hallberg (JCP, 1997). Additional issues related to exact tracer conservation and how to ensure consistency between the barotropic and layered estimates of the free surface height are described in Hallberg and Adcroft (Ocean Modelling, 2009). This was the time stepping code that is used for most GOLD applications, including GFDL’s ESM2G Earth system model, and all of the examples provided with the MOM code (although several of these solutions are routinely verified by comparison with the slower unsplit schemes).
The subroutine step_MOM_dyn_split_RK2b actually does the time stepping, while register_restarts_dyn_split_RK2b sets the fields that are found in a full restart file with this scheme, and initialize_dyn_split_RK2b initializes the cpu clocks that are used in this module. For largely historical reasons, this module does not have its own control structure, but shares the same control structure with MOM.F90 and the other MOM_dynamics_… modules.
Type Documentation
- type mom_dynamics_split_rk2b/mom_dyn_split_rk2b_cs
MOM_dynamics_split_RK2b module control structure.
- Type fields:
% id_uh ::
integerDiagnostic IDs.% id_vh ::
integerDiagnostic IDs.% id_umo ::
integerDiagnostic IDs.% id_vmo ::
integerDiagnostic IDs.% id_umo_2d ::
integerDiagnostic IDs.% id_vmo_2d ::
integerDiagnostic IDs.% id_pfu ::
integerDiagnostic IDs.% id_pfv ::
integerDiagnostic IDs.% id_cau ::
integerDiagnostic IDs.% id_cav ::
integerDiagnostic IDs.% id_ueffa ::
integerDiagnostic IDs.% id_veffa ::
integerDiagnostic IDs.% id_h_pfu ::
integerDiagnostic IDs.% id_h_pfv ::
integerDiagnostic IDs.% id_hf_pfu_2d ::
integerDiagnostic IDs.% id_hf_pfv_2d ::
integerDiagnostic IDs.% id_intz_pfu_2d ::
integerDiagnostic IDs.% id_intz_pfv_2d ::
integerDiagnostic IDs.% id_pfu_visc_rem ::
integerDiagnostic IDs.% id_pfv_visc_rem ::
integerDiagnostic IDs.% id_h_cau ::
integerDiagnostic IDs.% id_h_cav ::
integerDiagnostic IDs.% id_hf_cau_2d ::
integerDiagnostic IDs.% id_hf_cav_2d ::
integerDiagnostic IDs.% id_intz_cau_2d ::
integerDiagnostic IDs.% id_intz_cav_2d ::
integerDiagnostic IDs.% id_cau_visc_rem ::
integerDiagnostic IDs.% id_cav_visc_rem ::
integerDiagnostic IDs.% id_deta_dt ::
integerDiagnostic IDs.% id_uav ::
integerDiagnostic IDs.% id_vav ::
integerDiagnostic IDs.% id_u_bt_accel ::
integerDiagnostic IDs.% id_v_bt_accel ::
integerDiagnostic IDs.% id_h_u_bt_accel ::
integerDiagnostic IDs.% id_h_v_bt_accel ::
integerDiagnostic IDs.% id_hf_u_bt_accel_2d ::
integerDiagnostic IDs.% id_hf_v_bt_accel_2d ::
integerDiagnostic IDs.% id_intz_u_bt_accel_2d ::
integerDiagnostic IDs.% id_intz_v_bt_accel_2d ::
integerDiagnostic IDs.% id_u_bt_accel_visc_rem ::
integerDiagnostic IDs.% id_v_bt_accel_visc_rem ::
integerDiagnostic IDs.% cau ::
real, dimension(:, :, :), allocatableCAu = f*v - u.grad(u) [L T-2 ~> m s-2].% cau_pred ::
real, dimension(:, :, :), allocatableThe predictor step value of CAu = f*v - u.grad(u) [L T-2 ~> m s-2].% pfu ::
real, dimension(:, :, :), allocatablePFu = -dM/dx [L T-2 ~> m s-2].% pfu_stokes ::
real, dimension(:, :, :), allocatablePFu_Stokes = -d/dx int_r (u_L*duS/dr) [L T-2 ~> m s-2].% diffu ::
real, dimension(:, :, :), allocatableZonal acceleration due to convergence of the along-isopycnal stress tensor [L T-2 ~> m s-2].% cav ::
real, dimension(:, :, :), allocatableCAv = -f*u - u.grad(v) [L T-2 ~> m s-2].% cav_pred ::
real, dimension(:, :, :), allocatableThe predictor step value of CAv = -f*u - u.grad(v) [L T-2 ~> m s-2].% pfv ::
real, dimension(:, :, :), allocatablePFv = -dM/dy [L T-2 ~> m s-2].% pfv_stokes ::
real, dimension(:, :, :), allocatablePFv_Stokes = -d/dy int_r (v_L*dvS/dr) [L T-2 ~> m s-2].% diffv ::
real, dimension(:, :, :), allocatableMeridional acceleration due to convergence of the along-isopycnal stress tensor [L T-2 ~> m s-2].% visc_rem_u ::
real, dimension(:, :, :), allocatableBoth the fraction of the zonal momentum originally in a layer that remains after a time-step of viscosity, and the fraction of a time-step worth of a barotropic acceleration that a layer experiences after viscosity is applied [nondim]. Nondimensional between 0 (at the bottom) and 1 (far above).% u_accel_bt ::
real, dimension(:, :, :), allocatableThe zonal layer accelerations due to the difference between the barotropic accelerations and the baroclinic accelerations that were fed into the barotopic calculation [L T-2 ~> m s-2].% visc_rem_v ::
real, dimension(:, :, :), allocatableBoth the fraction of the meridional momentum originally in a layer that remains after a time-step of viscosity, and the fraction of a time-step worth of a barotropic acceleration that a layer experiences after viscosity is applied [nondim]. Nondimensional between 0 (at the bottom) and 1 (far above).% v_accel_bt ::
real, dimension(:, :, :), allocatableThe meridional layer accelerations due to the difference between the barotropic accelerations and the baroclinic accelerations that were fed into the barotopic calculation [L T-2 ~> m s-2].% eta ::
real, dimension(:, :), allocatableInstantaneous free surface height (in Boussinesq mode) or column mass anomaly (in non-Boussinesq mode) [H ~> m or kg m-2].% u_av ::
real, dimension(:, :, :), allocatablelayer x-velocity with vertical mean replaced by time-mean barotropic velocity over a baroclinic timestep [L T-1 ~> m s-1]% v_av ::
real, dimension(:, :, :), allocatablelayer y-velocity with vertical mean replaced by time-mean barotropic velocity over a baroclinic timestep [L T-1 ~> m s-1]% h_av ::
real, dimension(:, :, :), allocatablearithmetic mean of two successive layer thicknesses [H ~> m or kg m-2]% eta_pf ::
real, dimension(:, :), allocatableinstantaneous SSH used in calculating PFu and PFv [H ~> m or kg m-2]% uhbt ::
real, dimension(:, :), allocatableaverage x-volume or mass flux determined by the barotropic solver [H L2 T-1 ~> m3 s-1 or kg s-1]. uhbt is roughly equal to the vertical sum of uh.% vhbt ::
real, dimension(:, :), allocatableaverage y-volume or mass flux determined by the barotropic solver [H L2 T-1 ~> m3 s-1 or kg s-1]. vhbt is roughly equal to vertical sum of vh.% pbce ::
real, dimension(:, :, :), allocatablepbce times eta gives the baroclinic pressure anomaly in each layer due to free surface height anomalies [L2 H-1 T-2 ~> m s-2 or m4 kg-1 s-2].% du_av_inst ::
real, dimension(:, :), allocatableThe barotropic zonal velocity increment between filtered and instantaneous velocities [L T-1 ~> m s-1].% dv_av_inst ::
real, dimension(:, :), allocatableThe barotropic meridional velocity increment between filtered and instantaneous velocities [L T-1 ~> m s-1].% taux_bot ::
real, dimension(:,:), pointerfrictional x-bottom stress from the ocean to the seafloor [R L Z T-2 ~> Pa]% tauy_bot ::
real, dimension(:,:), pointerfrictional y-bottom stress from the ocean to the seafloor [R L Z T-2 ~> Pa]% bt_cont ::
type(bt_cont_type), pointerA structure with elements that describe the effective summed open face areas as a function of barotropic flow.% bt_adj_corr_mass_src ::
logicalIf true, recalculates the barotropic mass source after predictor step. This should make little difference in the deep ocean but appears to help for vanished layers.% split_bottom_stress ::
logicalIf true, provide the bottom stress calculated by the vertical viscosity to the barotropic solver.% dtbt_use_bt_cont ::
logicalIf true, use BT_cont to calculate DTBT.% calculate_sal ::
logicalIf true, calculate self-attraction and loading.% use_tides ::
logicalIf true, tidal forcing is enabled.% use_ha ::
logicalIf true, perform inline harmonic analysis.% remap_aux ::
logicalIf true, apply ALE remapping to all of the auxiliary 3-D variables that are needed to reproduce across restarts, similarly to what is done with the primary state variables.% be ::
realA nondimensional number from 0.5 to 1 that controls the backward weighting of the time stepping scheme [nondim].% begw ::
realA nondimensional number from 0 to 1 that controls the extent to which the treatment of gravity waves is forward-backward (0) or simulated backward Euler (1) [nondim]. 0 is often used.% debug ::
logicalIf true, write verbose checksums for debugging purposes.% debug_obc ::
logicalIf true, do additional calls resetting values to help verify the correctness of the open boundary condition code.% fpmix ::
logicalIf true, applies profiles of momentum flux magnitude and direction.% module_is_initialized ::
logicalRecord whether this module has been initialized.% visc_rem_dt_bug ::
logicalIf true, recover a bug that uses dt_pred rather than dt for vertvisc_rem at the end of predictor.% diag ::
type(diag_ctrl), pointerA structure that is used to regulate the timing of diagnostic output.% adp ::
type(accel_diag_ptrs), pointerA structure pointing to the various accelerations in the momentum equations, which can later be used to calculate derived diagnostics like energy budgets.% ad_pred ::
type(accel_diag_ptrs), pointerA structure pointing to the various predictor step accelerations in the momentum equations, which can be used to debug truncations.% cdp ::
type(cont_diag_ptrs), pointerA structure with pointers to various terms in the continuity equations, which can later be used to calculate derived diagnostics like energy budgets.% hor_visc ::
type(hor_visc_cs)A pointer to the horizontal viscosity control structure.% continuity_csp ::
type(continuity_cs)A pointer to the continuity control structure.% coriolisadv ::
type(coriolisadv_cs)The CoriolisAdv control structure.% pressureforce_csp ::
type(pressureforce_cs)A pointer to the PressureForce control structure.% vertvisc_csp ::
type(vertvisc_cs), pointerA pointer to a structure containing interface height diffusivities.% set_visc_csp ::
type(set_visc_cs), pointerA pointer to the set_visc control structure.% barotropic_csp ::
type(barotropic_cs)A pointer to the barotropic stepping control structure.% sal_csp ::
type(sal_cs)A pointer to the SAL control structure.% tides_csp ::
type(tidal_forcing_cs)A pointer to the tidal forcing control structure.% ha_csp ::
type(harmonic_analysis_cs)A pointer to the harmonic analysis control structure.% ale_csp ::
type(ale_cs), pointerA pointer to the ALE control structure.% obc ::
type(ocean_obc_type), pointerA pointer to an open boundary condition type that specifies whether, where, and what open boundary conditions are used. If no open BCs are used, this pointer stays nullified. Flather OBCs use open boundary_CS as well.% update_obc_csp ::
type(update_obc_cs), pointerA pointer to the update_OBC control structure.% pass_eta ::
type(group_pass_type)Structure for group halo pass.% pass_visc_rem ::
type(group_pass_type)Structure for group halo pass.% pass_uvp ::
type(group_pass_type)Structure for group halo pass.% pass_uv_inst ::
type(group_pass_type)Structure for group halo pass.% pass_hp_uv ::
type(group_pass_type)Structure for group halo pass.% pass_hp_uhvh ::
type(group_pass_type)Structure for group halo pass.% pass_h_uv ::
type(group_pass_type)Structure for group halo pass.
Function/Subroutine Documentation
- subroutine mom_dynamics_split_rk2b/step_mom_dyn_split_rk2b(u_av, v_av, h, tv, visc, Time_local, dt, forces, p_surf_begin, p_surf_end, uh, vh, uhtr, vhtr, eta_av, G, GV, US, CS, calc_dtbt, VarMix, MEKE, thickness_diffuse_CSp, pbv, Waves)
RK2 splitting for time stepping MOM adiabatic dynamics.
- Parameters:
g :: [inout] Ocean grid structure
gv :: [in] Ocean vertical grid structure
us :: [in] A dimensional unit scaling type
u_av ::
u_av[inout] Zonal velocity [L T-1 ~> m s-1]v_av ::
v_av[inout] Meridional velocity [L T-1 ~> m s-1]h ::
h[inout] Layer thickness [H ~> m or kg m-2]tv ::
tv[in] Thermodynamic typevisc ::
visc[inout] Vertical visc, bottom drag, and relatedtime_local :: [in] Model time at end of time step
dt ::
dt[in] Baroclinic dynamics time step [T ~> s]forces ::
forces[in] A structure with the driving mechanical forcesp_surf_begin ::
p_surf_beginSurface pressure at the start of this dynamic time step [R L2 T-2 ~> Pa]p_surf_end ::
p_surf_endSurface pressure at the end of this dynamic time step [R L2 T-2 ~> Pa]uh ::
uh[inout] Zonal volume or mass transportvh ::
vh[inout] Meridional volume or mass transportuhtr ::
uhtr[inout] Accumulated zonal volume or mass transportvhtr ::
vhtr[inout] Accumulated meridional volume or mass transporteta_av ::
eta_av[out] Free surface height or column mass averaged over time step [H ~> m or kg m-2]cs :: Module control structure
calc_dtbt ::
calc_dtbt[in] If true, recalculate the barotropic time stepvarmix :: [inout] Variable mixing control structure
meke :: [inout] MEKE fields
thickness_diffuse_csp :: [inout] Pointer to a structure containing interface height diffusivities
pbv ::
pbv[in] porous barrier fractional cell metricswaves :: A pointer to a structure containing fields related to the surface wave conditions
- Call to:
mom_error_handler::calltree_entermom_error_handler::calltree_leavemom_error_handler::calltree_waypointmom_coriolisadv::coriolisadv_stencilmom_diag_mediator::diag_update_remap_gridsmom_interface_heights::find_col_avg_spvid_clock_btcalcid_clock_btforceid_clock_btstepid_clock_continuityid_clock_corid_clock_horviscid_clock_mom_updateid_clock_passid_clock_presid_clock_vertviscmom_checksum_packages::mom_accel_chksummom_set_visc::set_viscous_mlmom_wave_interface::stokes_pgfmom_boundary_update::update_obc_datamom_open_boundary::update_segment_thickness_reservoirsmom_vert_friction::updatecfltruncationvalue
- subroutine mom_dynamics_split_rk2b/register_restarts_dyn_split_rk2b(HI, GV, US, param_file, CS, restart_CS, uh, vh)
This subroutine sets up any auxiliary restart variables that are specific to the split-explicit time stepping scheme. All variables registered here should have the ability to be recreated if they are not present in a restart file.
- Parameters:
hi :: [in] Horizontal index structure
gv :: [in] ocean vertical grid structure
us :: [in] A dimensional unit scaling type
param_file ::
param_file[in] parameter filecs :: module control structure
restart_cs :: [inout] MOM restart control structure
uh ::
uh[inout] zonal volume or mass transport [H L2 T-1 ~> m3 s-1 or kg s-1]vh ::
vh[inout] merid volume or mass transport [H L2 T-1 ~> m3 s-1 or kg s-1]
- Call to:
- subroutine mom_dynamics_split_rk2b/remap_dyn_split_rk2b_aux_vars(G, GV, CS, h_old_u, h_old_v, h_new_u, h_new_v, ALE_CSp)
This subroutine does remapping for the auxiliary restart variables that are used with the split RK2 time stepping scheme.
- Parameters:
g :: [inout] ocean grid structure
gv :: [in] ocean vertical grid structure
cs :: module control structure
h_old_u ::
h_old_u[in] Source grid thickness at zonalh_old_v ::
h_old_v[in] Source grid thickness at meridionalh_new_u ::
h_new_u[in] Destination grid thickness at zonalh_new_v ::
h_new_v[in] Destination grid thickness at meridionalale_csp :: ALE control structure to use when remapping
- Called from:
- subroutine mom_dynamics_split_rk2b/initialize_dyn_split_rk2b(u, v, h, tv, uh, vh, eta, Time, G, GV, US, param_file, diag, CS, HA_CSp, restart_CS, dt, Accel_diag, Cont_diag, MIS, VarMix, MEKE, thickness_diffuse_CSp, OBC, update_OBC_CSp, ALE_CSp, set_visc, visc, dirs, ntrunc, pbv, calc_dtbt, cont_stencil, dyn_h_stencil)
This subroutine initializes all of the variables that are used by this dynamic core, including diagnostics and the cpu clocks.
- Parameters:
g :: [inout] 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] merid velocity [L T-1 ~> m s-1]h ::
h[inout] layer thickness [H ~> m or kg m-2]tv ::
tv[in] Thermodynamic typeuh ::
uh[inout] zonal volume/mass transport [H L2 T-1 ~> m3 s-1 or kg s-1]vh ::
vh[inout] merid volume/mass transport [H L2 T-1 ~> m3 s-1 or kg s-1]eta ::
eta[inout] free surface height or column mass [H ~> m or kg m-2]time :: [in] current model time
param_file ::
param_file[in] parameter file for parsingdiag ::
diag[inout] to control diagnosticscs :: module control structure
ha_csp :: A pointer to the control structure of the harmonic analysis module
restart_cs :: [inout] MOM restart control structure
dt ::
dt[in] time step [T ~> s]accel_diag :: [inout] points to momentum equation terms for budget analysis
cont_diag :: [inout] points to terms in continuity equation
mis :: [inout] “MOM6 internal state” used to pass diagnostic pointers
varmix :: [inout] points to spatially variable viscosities
meke :: [inout] MEKE fields
thickness_diffuse_csp :: [inout] Pointer to the control structure used for the isopycnal height diffusive transport.
obc :: points to OBC related fields
update_obc_csp :: points to OBC update related fields
ale_csp :: points to ALE control structure
set_visc ::
set_visc[in] set_visc control structurevisc ::
visc[inout] vertical viscosities, bottom drag, and relateddirs ::
dirs[in] contains directory pathsntrunc ::
ntrunc[inout] A target for the variable that records the number of times the velocity is truncated (this should be 0).calc_dtbt ::
calc_dtbt[out] If true, recalculate the barotropic time steppbv ::
pbv[in] porous barrier fractional cell metricscont_stencil ::
cont_stencil[out] The stencil for thickness from the continuity solver.dyn_h_stencil ::
dyn_h_stencil[out] The stencil for thickness for the dynamics based on the continuity solver and Coriolis scheme.
- Call to:
mom_coriolisadv::coriolisadv_initmom_coriolisadv::coriolisadv_stencilmom_verticalgrid::get_flux_unitsmom_harmonic_analysis::ha_initmom_hor_visc::hor_visc_initid_clock_btcalcid_clock_btforceid_clock_btstepid_clock_continuityid_clock_corid_clock_horviscid_clock_mom_updateid_clock_passid_clock_presid_clock_vertviscmom_restart::is_new_runmom_pressureforce::pressureforce_initmom_self_attr_load::sal_initmom_tidal_forcing::tidal_forcing_initmom_vert_friction::updatecfltruncationvalue
- subroutine mom_dynamics_split_rk2b/end_dyn_split_rk2b(CS)
Close the dyn_split_RK2b module.
- Parameters:
cs :: module control structure
- Call to:
mom_barotropic::barotropic_endmom_coriolisadv::coriolisadv_endmom_hor_visc::hor_visc_endmom_self_attr_load::sal_endmom_tidal_forcing::tidal_forcing_end