namespace mom_barotropic¶
Overview¶
Baropotric solver. More…
namespace mom_barotropic { // global functions subroutine, public btstep( U_in U_in, V_in V_in, eta_in eta_in, dt dt, bc_accel_u bc_accel_u, bc_accel_v bc_accel_v, forces forces, pbce pbce, eta_PF_in eta_PF_in, U_Cor U_Cor, V_Cor V_Cor, accel_layer_u accel_layer_u, accel_layer_v accel_layer_v, eta_out eta_out, uhbtav uhbtav, vhbtav vhbtav, G G, GV GV, US US, CS CS, visc_rem_u visc_rem_u, visc_rem_v visc_rem_v, etaav etaav, OBC OBC, BT_cont BT_cont, eta_PF_start eta_PF_start, taux_bot taux_bot, tauy_bot tauy_bot, uh0 uh0, vh0 vh0, u_uh0 u_uh0, v_vh0 v_vh0 ); subroutine, public set_dtbt( G G, GV GV, US US, CS CS, eta eta, pbce pbce, BT_cont BT_cont, gtot_est gtot_est, SSH_add SSH_add ); subroutine, public btcalc( h h, G G, GV GV, CS CS, h_u h_u, h_v h_v, may_use_default may_use_default, OBC OBC ); subroutine, public bt_mass_source(h h, eta eta, set_cor set_cor, G G, GV GV, CS CS); subroutine, public barotropic_init( u u, v v, h h, eta eta, Time Time, G G, GV GV, US US, param_file param_file, diag diag, CS CS, restart_CS restart_CS, calc_dtbt calc_dtbt, BT_cont BT_cont, tides_CSp tides_CSp ); subroutine, public barotropic_get_tav(CS CS, ubtav ubtav, vbtav vbtav, G G, US US); subroutine, public barotropic_end(CS CS); subroutine, public register_barotropic_restarts( HI HI, GV GV, param_file param_file, CS CS, restart_CS restart_CS ); } // namespace mom_barotropic
Detailed Documentation¶
Baropotric solver.
By Robert Hallberg, April 1994 - January 2007
This program contains the subroutines that time steps the linearized barotropic equations. btstep is used to actually time step the barotropic equations, and contains most of the substance of this module.
btstep uses a forwards-backwards based scheme to time step the barotropic equations, returning the layers’ accelerations due to the barotropic changes in the ocean state, the final free surface height (or column mass), and the volume (or mass) fluxes summed through the layers and averaged over the baroclinic time step. As input, btstep takes the initial 3-D velocities, the inital free surface height, the 3-D accelerations of the layers, and the external forcing. Everything in btstep is cast in terms of anomalies, so if everything is in balance, there is explicitly no acceleration due to btstep.
The spatial discretization of the continuity equation is second order accurate. A flux conservative form is used to guarantee global conservation of volume. The spatial discretization of the momentum equation is second order accurate. The Coriolis force is written in a form which does not contribute to the energy tendency and which conserves linearized potential vorticity, f/D. These terms are exactly removed from the baroclinic momentum equations, so the linearization of vorticity advection will not degrade the overall solution.
btcalc calculates the fractional thickness of each layer at the velocity points, for later use in calculating the barotropic velocities and the averaged accelerations. Harmonic mean thicknesses (i.e. 2*h_L*h_R/(h_L + h_R)) are used to avoid overly strong weighting of overly thin layers. This may later be relaxed to use thicknesses determined from the continuity equations.
bt_mass_source determines the real mass sources for the barotropic solver, along with the corrective pseudo-fluxes that keep the barotropic and baroclinic estimates of the free surface height close to each other. Given the layer thicknesses and the free surface height that correspond to each other, it calculates a corrective mass source that is added to the barotropic continuity* equation, and optionally adjusts a slowly varying correction rate. Newer algorithmic changes have deemphasized the need for this, but it is still here to add net water sources to the barotropic solver.*
barotropic_init allocates and initializes any barotropic arrays that have not been read from a restart file, reads parameters from the inputfile, and sets up diagnostic fields.
barotropic_end deallocates anything allocated in barotropic_init or register_barotropic_restarts.
register_barotropic_restarts is used to indicate any fields that are private to the barotropic solver that need to be included in the restart files, and to ensure that they are read.
Global Functions¶
subroutine, public btstep( U_in U_in, V_in V_in, eta_in eta_in, dt dt, bc_accel_u bc_accel_u, bc_accel_v bc_accel_v, forces forces, pbce pbce, eta_PF_in eta_PF_in, U_Cor U_Cor, V_Cor V_Cor, accel_layer_u accel_layer_u, accel_layer_v accel_layer_v, eta_out eta_out, uhbtav uhbtav, vhbtav vhbtav, G G, GV GV, US US, CS CS, visc_rem_u visc_rem_u, visc_rem_v visc_rem_v, etaav etaav, OBC OBC, BT_cont BT_cont, eta_PF_start eta_PF_start, taux_bot taux_bot, tauy_bot tauy_bot, uh0 uh0, vh0 vh0, u_uh0 u_uh0, v_vh0 v_vh0 )
This subroutine time steps the barotropic equations explicitly. For gravity waves, anything between a forwards-backwards scheme and a simulated backwards Euler scheme is used, with bebt between 0.0 and 1.0 determining the scheme. In practice, bebt must be of order 0.2 or greater. A forwards-backwards treatment of the Coriolis terms is always used.
Parameters:
g |
The ocean’s grid structure. |
gv |
The ocean’s vertical grid structure. |
us |
A dimensional unit scaling type |
u_in |
The initial (3-D) zonal velocity [L T-1 ~> m s-1]. |
v_in |
The initial (3-D) meridional velocity [L T-1 ~> m s-1]. |
eta_in |
The initial barotropic free surface height anomaly or column mass anomaly [H ~> m or kg m-2]. |
dt |
The time increment to integrate over. |
bc_accel_u |
The zonal baroclinic accelerations [m s-2]. |
bc_accel_v |
The meridional baroclinic accelerations, [L T-2 ~> m s-2]. |
forces |
A structure with the driving mechanical forces |
pbce |
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]. |
eta_pf_in |
The 2-D eta field (either SSH anomaly or column mass anomaly) that was used to calculate the input pressure gradient accelerations (or its final value if eta_PF_start is provided [H ~> m or kg m-2]. Note: eta_in, pbce, and eta_PF_in must have up-to-date values in the first point of their halos. |
u_cor |
The (3-D) zonal velocities used to calculate the Coriolis terms in bc_accel_u [L T-1 ~> m s-1]. |
v_cor |
The (3-D) meridional velocities used to calculate the Coriolis terms in bc_accel_u [L T-1 ~> m s-1]. |
accel_layer_u |
The zonal acceleration of each layer due to the barotropic calculation [L T-2 ~> m s-2]. |
accel_layer_v |
The meridional acceleration of each layer due to the barotropic calculation [L T-2 ~> m s-2]. |
eta_out |
The final barotropic free surface height anomaly or column mass anomaly [H ~> m or kg m-2]. |
uhbtav |
the barotropic zonal volume or mass fluxes averaged through the barotropic steps [H L2 T-1 ~> m3 or kg s-1]. |
vhbtav |
the barotropic meridional volume or mass fluxes averaged through the barotropic steps [H L2 T-1 ~> m3 or kg s-1]. |
cs |
The control structure returned by a previous call to barotropic_init. |
visc_rem_u |
Both the fraction of the momentum originally in a layer that remains after a time-step of viscosity, and the fraction of a time-step’s worth of a barotropic acceleration that a layer experiences after viscosity is applied, in the zonal direction. Nondimensional between 0 (at the bottom) and 1 (far above the bottom). |
visc_rem_v |
Ditto for meridional direction [nondim]. |
etaav |
The free surface height or column mass averaged over the barotropic integration [H ~> m or kg m-2]. |
obc |
The open boundary condition structure. |
bt_cont |
A structure with elements that describe the effective open face areas as a function of barotropic flow. |
eta_pf_start |
The eta field consistent with the pressure gradient at the start of the barotropic stepping [H ~> m or kg m-2]. |
taux_bot |
The zonal bottom frictional stress from ocean to the seafloor [kg L Z T-2 m-3 ~> Pa]. |
tauy_bot |
The meridional bottom frictional stress from ocean to the seafloor [kg L Z T-2 m-3 ~> Pa]. |
uh0 |
The zonal layer transports at reference velocities [H L2 T-1 ~> m3 s-1 or kg s-1]. |
u_uh0 |
The velocities used to calculate uh0 [L T-1 ~> m s-1] |
vh0 |
The zonal layer transports at reference velocities [H L2 T-1 ~> m3 s-1 or kg s-1]. |
v_vh0 |
The velocities used to calculate vh0 [L T-1 ~> m s-1] |
subroutine, public set_dtbt( G G, GV GV, US US, CS CS, eta eta, pbce pbce, BT_cont BT_cont, gtot_est gtot_est, SSH_add SSH_add )
This subroutine automatically determines an optimal value for dtbt based on some state of the ocean.
Parameters:
g |
The ocean’s grid structure. |
gv |
The ocean’s vertical grid structure. |
us |
A dimensional unit scaling type |
cs |
Barotropic control structure. |
eta |
The barotropic free surface height anomaly or column mass anomaly [H ~> m or kg m-2]. |
pbce |
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]. |
bt_cont |
A structure with elements that describe the effective open face areas as a function of barotropic flow. |
gtot_est |
An estimate of the total gravitational acceleration [L2 Z-1 T-2 ~> m s-2]. |
ssh_add |
An additional contribution to SSH to provide a margin of error when calculating the external wave speed [Z ~> m]. |
subroutine, public btcalc( h h, G G, GV GV, CS CS, h_u h_u, h_v h_v, may_use_default may_use_default, OBC OBC )
btcalc calculates the barotropic velocities from the full velocity and thickness fields, determines the fraction of the total water column in each layer at velocity points, and determines a corrective fictitious mass source that will drive the barotropic estimate of the free surface height toward the baroclinic estimate.
Parameters:
g |
The ocean’s grid structure. |
gv |
The ocean’s vertical grid structure. |
h |
Layer thicknesses [H ~> m or kg m-2]. |
cs |
The control structure returned by a previous call to barotropic_init. |
h_u |
The specified thicknesses at u-points [H ~> m or kg m-2]. |
h_v |
The specified thicknesses at v-points [H ~> m or kg m-2]. |
may_use_default |
An optional logical argument to indicate that the default velocity point thicknesses may be used for this particular calculation, even though the setting of CShvel_scheme would usually require that h_u and h_v be passed in. |
obc |
Open boundary control structure. |
subroutine, public bt_mass_source( h h, eta eta, set_cor set_cor, G G, GV GV, CS CS )
bt_mass_source determines the appropriately limited mass source for the barotropic solver, along with a corrective fictitious mass source that will drive the barotropic estimate of the free surface height toward the baroclinic estimate.
Parameters:
g |
The ocean’s grid structure. |
gv |
The ocean’s vertical grid structure. |
h |
Layer thicknesses [H ~> m or kg m-2]. |
eta |
The free surface height that is to be corrected [H ~> m or kg m-2]. |
set_cor |
A flag to indicate whether to set the corrective fluxes (and update the slowly varying part of eta_cor) (.true.) or whether to incrementally update the corrective fluxes. |
cs |
The control structure returned by a previous call to barotropic_init. |
subroutine, public barotropic_init( u u, v v, h h, eta eta, Time Time, G G, GV GV, US US, param_file param_file, diag diag, CS CS, restart_CS restart_CS, calc_dtbt calc_dtbt, BT_cont BT_cont, tides_CSp tides_CSp )
barotropic_init initializes a number of time-invariant fields used in the barotropic calculation and initializes any barotropic fields that have not already been initialized.
Parameters:
g |
The ocean’s grid structure. |
gv |
The ocean’s vertical grid structure. |
us |
A dimensional unit scaling type |
u |
The zonal velocity [L T-1 ~> m s-1]. |
v |
The meridional velocity [L T-1 ~> m s-1]. |
h |
Layer thicknesses [H ~> m or kg m-2]. |
eta |
Free surface height or column mass anomaly |
time |
The current model time. |
param_file |
A structure to parse for run-time parameters. |
diag |
A structure that is used to regulate diagnostic output. |
cs |
A pointer to the control structure for this module that is set in register_barotropic_restarts. |
restart_cs |
A pointer to the restart control structure. |
calc_dtbt |
If true, the barotropic time step must be recalculated before stepping. |
bt_cont |
A structure with elements that describe the |
tides_csp |
A pointer to the control structure of the |
subroutine, public barotropic_get_tav( CS CS, ubtav ubtav, vbtav vbtav, G G, US US )
Copies ubtav and vbtav from private type into arrays.
Parameters:
cs |
Control structure for this module |
g |
Grid structure |
ubtav |
Zonal barotropic velocity averaged over a baroclinic timestep [L T-1 ~> m s-1] |
vbtav |
Meridional barotropic velocity averaged over a baroclinic timestep [L T-1 ~> m s-1] |
us |
A dimensional unit scaling type |
subroutine, public barotropic_end(CS CS)
Clean up the barotropic control structure.
Parameters:
cs |
Control structure to clear out. |
subroutine, public register_barotropic_restarts( HI HI, GV GV, param_file param_file, CS CS, restart_CS restart_CS )
This subroutine is used to register any fields from MOM_barotropic.F90
that should be written to or read from the restart file.
Parameters:
hi |
A horizontal index type structure. |
param_file |
A structure to parse for run-time parameters. |
cs |
A pointer that is set to point to the control structure for this module. |
gv |
The ocean’s vertical grid structure. |
restart_cs |
A pointer to the restart control structure. |