mom_barotropic module reference

Barotropic solver.

More…

Data Types

bt_obc_type

The barotropic stepping open boundary condition type.

barotropic_cs

The barotropic stepping control structure.

local_bt_cont_u_type

A description of the functional dependence of transport at a u-point.

local_bt_cont_v_type

A description of the functional dependence of transport at a v-point.

memory_size_type

A container for passing around active tracer point memory limits.

Functions/Subroutines

btstep()

This subroutine time steps the barotropic equations explicitly.

set_dtbt()

This subroutine automatically determines an optimal value for dtbt based on some state of the ocean.

apply_velocity_obcs()

The following 4 subroutines apply the open boundary conditions.

set_up_bt_obc()

This subroutine sets up the private structure used to apply the open boundary conditions, as developed by Mehmet Ilicak.

destroy_bt_obc()

Clean up the BT_OBC memory.

btcalc()

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.

find_uhbt()

The function find_uhbt determines the zonal transport for a given velocity, or with INTEGRAL_BT_CONT=True it determines the time-integrated zonal transport for a given time-integrated velocity.

find_duhbt_du()

The function find_duhbt_du determines the marginal zonal face area for a given velocity, or with INTEGRAL_BT_CONT=True for a given time-integrated velocity.

uhbt_to_ubt()

This function inverts the transport function to determine the barotopic velocity that is consistent with a given transport, or if INTEGRAL_BT_CONT=True this finds the time-integrated velocity that is consistent with a time-integrated transport.

find_vhbt()

The function find_vhbt determines the meridional transport for a given velocity, or with INTEGRAL_BT_CONT=True it determines the time-integrated meridional transport for a given time-integrated velocity.

find_dvhbt_dv()

The function find_dvhbt_dv determines the marginal meridional face area for a given velocity, or with INTEGRAL_BT_CONT=True for a given time-integrated velocity.

vhbt_to_vbt()

This function inverts the transport function to determine the barotopic velocity that is consistent with a given transport, or if INTEGRAL_BT_CONT=True this finds the time-integrated velocity that is consistent with a time-integrated transport.

set_local_bt_cont_types()

This subroutine sets up reordered versions of the BT_cont type in the local_BT_cont types, which have wide halos properly filled in.

adjust_local_bt_cont_types()

Adjust_local_BT_cont_types expands the range of velocities with a cubic curve translating velocities into transports to match the initial values of velocities and summed transports when the velocities are larger than the first guesses of the cubic transition velocities used to set up the local_BT_cont types.

bt_cont_to_face_areas()

This subroutine uses the BT_cont_type to find the maximum face areas, which can then be used for finding wave speeds, etc.

swap()

Swap the values of two real variables.

find_face_areas()

This subroutine determines the open face areas of cells for calculating the barotropic transport.

bt_mass_source()

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.

barotropic_init()

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.

barotropic_get_tav()

Copies ubtav and vbtav from private type into arrays.

barotropic_end()

Clean up the barotropic control structure.

register_barotropic_restarts()

This subroutine is used to register any fields from MOM_barotropic.F90 that should be written to or read from the restart file.

Detailed Description

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.

Type Documentation

type mom_barotropic/bt_obc_type

The barotropic stepping open boundary condition type.

Type fields:
  • % is_u_obc [integer,private] :: Index ranges for the open boundary conditions.

  • % ie_u_obc [integer,private] :: Index ranges for the open boundary conditions.

  • % js_u_obc [integer,private] :: Index ranges for the open boundary conditions.

  • % je_u_obc [integer,private] :: Index ranges for the open boundary conditions.

  • % is_v_obc [integer,private] :: Index ranges for the open boundary conditions.

  • % ie_v_obc [integer,private] :: Index ranges for the open boundary conditions.

  • % js_v_obc [integer,private] :: Index ranges for the open boundary conditions.

  • % je_v_obc [integer,private] :: Index ranges for the open boundary conditions.

  • % cg_u [real(:,:),allocatable, private] :: The external wave speed at u-points [L T-1 ~> m s-1].

  • % cg_v [real(:,:),allocatable, private] :: The external wave speed at u-points [L T-1 ~> m s-1].

  • % dz_u [real(:,:),allocatable, private] :: The total vertical column extent at the u-points [Z ~> m].

  • % dz_v [real(:,:),allocatable, private] :: The total vertical column extent at the v-points [Z ~> m].

  • % uhbt [real(:,:),allocatable, private] :: The zonal barotropic thickness fluxes specified for open boundary conditions (if any) [H L2 T-1 ~> m3 s-1 or kg s-1].

  • % vhbt [real(:,:),allocatable, private] :: The meridional barotropic thickness fluxes specified for open boundary conditions (if any) [H L2 T-1 ~> m3 s-1 or kg s-1].

  • % ubt_outer [real(:,:),allocatable, private] :: The zonal velocities just outside the domain, as set by the open boundary conditions [L T-1 ~> m s-1].

  • % vbt_outer [real(:,:),allocatable, private] :: The meridional velocities just outside the domain, as set by the open boundary conditions [L T-1 ~> m s-1].

  • % ssh_outer_u [real(:,:),allocatable, private] :: The surface height outside of the domain at a u-point with an open boundary condition [Z ~> m].

  • % ssh_outer_v [real(:,:),allocatable, private] :: The surface height outside of the domain at a v-point with an open boundary condition [Z ~> m].

  • % apply_u_obcs [logical,private] :: True if this PE has an open boundary at a u-point.

  • % apply_v_obcs [logical,private] :: True if this PE has an open boundary at a v-point.

  • % is_alloced [logical,private] :: True if BT_OBC is in use and has been allocated.

  • % pass_uv [type(group_pass_type),private] :: Structure for group halo pass.

  • % pass_uhvh [type(group_pass_type),private] :: Structure for group halo pass.

  • % pass_h [type(group_pass_type),private] :: Structure for group halo pass.

  • % pass_cg [type(group_pass_type),private] :: Structure for group halo pass.

  • % pass_eta_outer [type(group_pass_type),private] :: Structure for group halo pass.

type mom_barotropic/barotropic_cs

The barotropic stepping control structure.

Type fields:
  • % id_pfu_bt [integer] :: Diagnostic IDs.

  • % id_pfv_bt [integer] :: Diagnostic IDs.

  • % id_coru_bt [integer] :: Diagnostic IDs.

  • % id_corv_bt [integer] :: Diagnostic IDs.

  • % id_ubtforce [integer] :: Diagnostic IDs.

  • % id_vbtforce [integer] :: Diagnostic IDs.

  • % id_uaccel [integer] :: Diagnostic IDs.

  • % id_vaccel [integer] :: Diagnostic IDs.

  • % id_visc_rem_u [integer] :: Diagnostic IDs.

  • % id_visc_rem_v [integer] :: Diagnostic IDs.

  • % id_eta_cor [integer] :: Diagnostic IDs.

  • % id_ubt [integer] :: Diagnostic IDs.

  • % id_vbt [integer] :: Diagnostic IDs.

  • % id_eta_bt [integer] :: Diagnostic IDs.

  • % id_ubtav [integer] :: Diagnostic IDs.

  • % id_vbtav [integer] :: Diagnostic IDs.

  • % id_ubt_st [integer] :: Diagnostic IDs.

  • % id_vbt_st [integer] :: Diagnostic IDs.

  • % id_eta_st [integer] :: Diagnostic IDs.

  • % id_ubtdt [integer] :: Diagnostic IDs.

  • % id_vbtdt [integer] :: Diagnostic IDs.

  • % id_ubt_hifreq [integer] :: Diagnostic IDs.

  • % id_vbt_hifreq [integer] :: Diagnostic IDs.

  • % id_eta_hifreq [integer] :: Diagnostic IDs.

  • % id_uhbt_hifreq [integer] :: Diagnostic IDs.

  • % id_vhbt_hifreq [integer] :: Diagnostic IDs.

  • % id_eta_pred_hifreq [integer] :: Diagnostic IDs.

  • % id_gtotn [integer] :: Diagnostic IDs.

  • % id_gtots [integer] :: Diagnostic IDs.

  • % id_gtote [integer] :: Diagnostic IDs.

  • % id_gtotw [integer] :: Diagnostic IDs.

  • % id_uhbt [integer] :: Diagnostic IDs.

  • % id_frhatu [integer] :: Diagnostic IDs.

  • % id_vhbt [integer] :: Diagnostic IDs.

  • % id_frhatv [integer] :: Diagnostic IDs.

  • % id_frhatu1 [integer] :: Diagnostic IDs.

  • % id_frhatv1 [integer] :: Diagnostic IDs.

  • % id_btc_fa_u_ee [integer] :: Diagnostic IDs.

  • % id_btc_fa_u_e0 [integer] :: Diagnostic IDs.

  • % id_btc_fa_u_w0 [integer] :: Diagnostic IDs.

  • % id_btc_fa_u_ww [integer] :: Diagnostic IDs.

  • % id_btc_ubt_ee [integer] :: Diagnostic IDs.

  • % id_btc_ubt_ww [integer] :: Diagnostic IDs.

  • % id_btc_fa_v_nn [integer] :: Diagnostic IDs.

  • % id_btc_fa_v_n0 [integer] :: Diagnostic IDs.

  • % id_btc_fa_v_s0 [integer] :: Diagnostic IDs.

  • % id_btc_fa_v_ss [integer] :: Diagnostic IDs.

  • % id_btc_vbt_nn [integer] :: Diagnostic IDs.

  • % id_btc_vbt_ss [integer] :: Diagnostic IDs.

  • % id_btc_fa_u_rat0 [integer] :: Diagnostic IDs.

  • % id_btc_fa_v_rat0 [integer] :: Diagnostic IDs.

  • % id_btc_fa_h_rat0 [integer] :: Diagnostic IDs.

  • % id_uhbt0 [integer] :: Diagnostic IDs.

  • % id_vhbt0 [integer] :: Diagnostic IDs.

  • % real (* q_d [*) :: The fraction of the total column thickness interpolated to u grid points in each layer [nondim].

  • % real :: The fraction of the total column thickness interpolated to v grid points in each layer [nondim].

  • % real :: Inverse of the total thickness at u grid points [H-1 ~> m-1 or m2 kg-1].

  • % real :: A spatially varying linear drag coefficient acting on the zonal barotropic flow [H T-1 ~> m s-1 or kg m-2 s-1].

  • % real :: The barotropic solvers estimate of the zonal velocity that will be the initial condition for the next call to btstep [L T-1 ~> m s-1].

  • % real :: The barotropic zonal velocity averaged over the baroclinic time step [L T-1 ~> m s-1].

  • % real :: Inverse of the basin depth at v grid points [Z-1 ~> m-1].

  • % real :: A spatially varying linear drag coefficient acting on the zonal barotropic flow [H T-1 ~> m s-1 or kg m-2 s-1].

  • % real :: The barotropic solvers estimate of the zonal velocity that will be the initial condition for the next call to btstep [L T-1 ~> m s-1].

  • % real :: The barotropic meridional velocity averaged over the baroclinic time step [L T-1 ~> m s-1].

  • % real :: The difference between the free surface height from the barotropic calculation and the sum of the layer thicknesses. This difference is imposed as a forcing term in the barotropic calculation over a baroclinic timestep [H ~> m or kg m-2].

  • % real :: A limit on the rate at which eta_cor can be applied while avoiding instability [H T-1 ~> m s-1 or kg m-2 s-1]. This is only used if CSbound_BT_corr is true.

  • % real :: Test vector components for checking grid polarity [nondim].

  • % real :: Test vector components for checking grid polarity [nondim].

  • % real :: A copy of bathyT (ocean bottom depth) with wide halos [Z ~> m].

  • % real :: This is a copy of GIareaT with wide halos, but will still utilize the macro IareaT when referenced, [L-2 ~> m-2].

  • % real :: A simply averaged depth at u points recast as a thickness [H ~> m or kg m-2].

  • % real :: A copy of Gdy_Cu with wide halos [L ~> m].

  • % real :: A copy of GIdxCu with wide halos [L-1 ~> m-1].

  • % real :: A simply averaged depth at v points recast as a thickness [H ~> m or kg m-2].

  • % real :: A copy of Gdx_Cv with wide halos [L ~> m].

  • % real :: A copy of GIdyCv with wide halos [L-1 ~> m-1].

  • % real :: f / D at PV points [Z-1 T-1 ~> m-1 s-1].

  • % frhatu1 [real(:,:,:),allocatable] :: Predictor step values of frhatu stored for diagnostics [nondim].

  • % frhatv1 [real(:,:,:),allocatable] :: Predictor step values of frhatv stored for diagnostics [nondim].

  • % bt_obc [type( bt_obc_type )] :: A structure with all of this modules fields for applying open boundary conditions.

  • % dtbt [real] :: The barotropic time step [T ~> s].

  • % dtbt_fraction [real] :: The fraction of the maximum time-step that should used [nondim]. The default is 0.98.

  • % dtbt_max [real] :: The maximum stable barotropic time step [T ~> s].

  • % dt_bt_filter [real] :: The time-scale over which the barotropic mode solutions are filtered [T ~> s] if positive, or as a fraction of DT if negative [nondim]. This can never be taken to be longer than 2*dt. Set this to 0 to apply no filtering.

  • % nstep_last [integer] :: The number of barotropic timesteps per baroclinic time step the last time btstep was called.

  • % bebt [real] :: A nondimensional number, from 0 to 1, that determines the gravity wave time stepping scheme [nondim]. 0.0 gives a forward-backward scheme, while 1.0 give backward Euler. In practice, bebt should be of order 0.2 or greater.

  • % rho_bt_lin [real] :: A density that is used to convert total water column thicknesses into mass in non-Boussinesq mode with linearized options in the barotropic solver or when estimating the stable barotropic timestep without access to the full baroclinic model state [R ~> kg m-3].

  • % split [logical] :: If true, use the split time stepping scheme.

  • % bound_bt_corr [logical] :: If true, the magnitude of the fake mass source in the barotropic equation that drives the two estimates of the free surface height toward each other is bounded to avoid driving corrective velocities that exceed MAXCFL_BT_CONT.

  • % gradual_bt_ics [logical] :: If true, adjust the initial conditions for the barotropic solver to the values from the layered solution over a whole timestep instead of instantly. This is a decent approximation to the inclusion of sum(u dh_dt) while also correcting for truncation errors.

  • % sadourny [logical] :: If true, the Coriolis terms are discretized with Sadourny’s energy conserving scheme, otherwise the Arakawa & Hsu scheme is used. If the deformation radius is not resolved Sadourny’s scheme should probably be used.

  • % integral_bt_cont [logical] :: If true, use the time-integrated velocity over the barotropic steps to determine the integrated transports used to update the continuity equation. Otherwise the transports are the sum of the transports based on a series of instantaneous velocities and the BT_CONT_TYPE for transports. This is only valid if a BT_CONT_TYPE is used.

  • % nonlinear_continuity [logical] :: If true, the barotropic continuity equation uses the full ocean thickness for transport.

  • % nonlin_cont_update_period [integer] :: The number of barotropic time steps between updates to the face area, or 0 only to update at the start of a call to btstep. The default is 1.

  • % bt_project_velocity [logical] :: If true, step the barotropic velocity first and project out the velocity tendency by 1+BEBT when calculating the transport. The default (false) is to use a predictor continuity step to find the pressure field, and then do a corrector continuity step using a weighted average of the old and new velocities, with weights of (1-BEBT) and BEBT.

  • % nonlin_stress [logical] :: If true, use the full depth of the ocean at the start of the barotropic step when calculating the surface stress contribution to the barotropic acclerations. Otherwise use the depth based on bathyT.

  • % bt_coriolis_scale [real] :: A factor by which the barotropic Coriolis acceleration anomaly terms are scaled [nondim].

  • % answer_date [integer] :: The vintage of the expressions in the barotropic solver. Values below 20190101 recover the answers from the end of 2018, while higher values use more efficient or general expressions.

  • % dynamic_psurf [logical] :: If true, add a dynamic pressure due to a viscous ice shelf, for instance.

  • % dmin_dyn_psurf [real] :: The minimum total thickness to use in limiting the size of the dynamic surface pressure for stability [H ~> m or kg m-2].

  • % ice_strength_length [real] :: The length scale at which the damping rate due to the ice strength should be the same as if a Laplacian were applied [L ~> m].

  • % const_dyn_psurf [real] :: The constant that scales the dynamic surface pressure [nondim]. Stable values are < ~1.0. The default is 0.9.

  • % calculate_sal [logical] :: If true, calculate self-attration and loading.

  • % tidal_sal_bug [logical] :: If true, the tidal self-attraction and loading anomaly in the barotropic solver has the wrong sign, replicating a long-standing bug.

  • % g_extra [real] :: A nondimensional factor by which gtot is enhanced [nondim].

  • % hvel_scheme [integer] :: An integer indicating how the thicknesses at velocity points are calculated. Valid values are given by the parameters defined below: HARMONIC, ARITHMETIC, HYBRID, and FROM_BT_CONT.

  • % strong_drag [logical] :: If true, use a stronger estimate of the retarding effects of strong bottom drag.

  • % linear_wave_drag [logical] :: If true, apply a linear drag to the barotropic velocities, using rates set by lin_drag_u & _v divided by the depth of the ocean.

  • % linearized_bt_pv [logical] :: If true, the PV and interface thicknesses used in the barotropic Coriolis calculation is time invariant and linearized.

  • % use_wide_halos [logical] :: If true, use wide halos and march in during the barotropic time stepping for efficiency.

  • % clip_velocity [logical] :: If true, limit any velocity components that are are large enough for a CFL number to exceed CFL_trunc. This should only be used as a desperate debugging measure.

  • % debug [logical] :: If true, write verbose checksums for debugging purposes.

  • % debug_bt [logical] :: If true, write verbose checksums for debugging purposes.

  • % vel_underflow [real] :: Velocity components smaller than vel_underflow are set to 0 [L T-1 ~> m s-1].

  • % maxvel [real] :: Velocity components greater than maxvel are truncated to maxvel [L T-1 ~> m s-1].

  • % cfl_trunc [real] :: If clip_velocity is true, velocity components will be truncated when they are large enough that the corresponding CFL number exceeds this value [nondim].

  • % maxcfl_bt_cont [real] :: The maximum permitted CFL number associated with the barotropic accelerations from the summed velocities times the time-derivatives of thicknesses [nondim]. The default is 0.1, and there will probably be real problems if this were set close to 1.

  • % bt_cont_bounds [logical] :: If true, use the BT_cont_type variables to set limits on the magnitude of the corrective mass fluxes.

  • % visc_rem_u_uh0 [logical] :: If true, use the viscous remnants when estimating the barotropic velocities that were used to calculate uh0 and vh0. False is probably the better choice.

  • % adjust_bt_cont [logical] :: If true, adjust the curve fit to the BT_cont type that is used by the barotropic solver to match the transport about which the flow is being linearized.

  • % use_old_coriolis_bracket_bug [logical] :: If True, use an order of operations that is not bitwise rotationally symmetric in the meridional Coriolis term of the barotropic solver.

  • % tidal_sal_flather [logical] :: Apply adjustment to external gravity wave speed consistent with tidal self-attraction and loading used within the barotropic solver.

  • % time [type(time_type),pointer] :: A pointer to the ocean models clock.

  • % diag [type( diag_ctrl ),pointer] :: A structure that is used to regulate the timing of diagnostic output.

  • % bt_domain [type(mom_domain_type),pointer] :: Barotropic MOM domain.

  • % debug_bt_hi [type( hor_index_type ),pointer] :: debugging copy of horizontal index_type

  • % sal_csp [type( sal_cs ),pointer] :: Control structure for SAL.

  • % module_is_initialized [logical] :: If true, module has been initialized.

  • % isdw [integer] :: The lower i-memory limit for the wide halo arrays.

  • % iedw [integer] :: The upper i-memory limit for the wide halo arrays.

  • % jsdw [integer] :: The lower j-memory limit for the wide halo arrays.

  • % jedw [integer] :: The upper j-memory limit for the wide halo arrays.

  • % pass_q_dcor [type(group_pass_type)] :: Handle for a group halo pass.

  • % pass_gtot [type(group_pass_type)] :: Handle for a group halo pass.

  • % pass_tmp_uv [type(group_pass_type)] :: Handle for a group halo pass.

  • % pass_eta_bt_rem [type(group_pass_type)] :: Handle for a group halo pass.

  • % pass_force_hbt0_cor_ref [type(group_pass_type)] :: Handle for a group halo pass.

  • % pass_dat_uv [type(group_pass_type)] :: Handle for a group halo pass.

  • % pass_eta_ubt [type(group_pass_type)] :: Handle for a group halo pass.

  • % pass_etaav [type(group_pass_type)] :: Handle for a group halo pass.

  • % pass_ubt_cor [type(group_pass_type)] :: Handle for a group halo pass.

  • % pass_ubta_uhbta [type(group_pass_type)] :: Handle for a group halo pass.

  • % pass_e_anom [type(group_pass_type)] :: Handle for a group halo pass.

  • % pass_spv_avg [type(group_pass_type)] :: Handle for a group halo pass.

type mom_barotropic/local_bt_cont_u_type

A description of the functional dependence of transport at a u-point.

Type fields:
  • % fa_u_ee [real,private] :: The effective open face area for zonal barotropic transport drawing from locations far to the east [H L ~> m2 or kg m-1].

  • % fa_u_e0 [real,private] :: The effective open face area for zonal barotropic transport drawing from nearby to the east [H L ~> m2 or kg m-1].

  • % fa_u_w0 [real,private] :: The effective open face area for zonal barotropic transport drawing from nearby to the west [H L ~> m2 or kg m-1].

  • % fa_u_ww [real,private] :: The effective open face area for zonal barotropic transport drawing from locations far to the west [H L ~> m2 or kg m-1].

  • % ubt_ww [real,private] :: uBT_WW is the barotropic velocity [L T-1 ~> m s-1], or with INTEGRAL_BT_CONTINUITY the time-integrated barotropic velocity [L ~> m], beyond which the marginal open face area is FA_u_WW. uBT_WW must be non-negative.

  • % ubt_ee [real,private] :: uBT_EE is a barotropic velocity [L T-1 ~> m s-1], or with INTEGRAL_BT_CONTINUITY the time-integrated barotropic velocity [L ~> m], beyond which the marginal open face area is FA_u_EE. uBT_EE must be non-positive.

  • % uh_crvw [real,private] :: The curvature of face area with velocity for flow from the west [H T2 L-1 ~> s2 or kg s2 m-3] or [H L-1 ~> nondim or kg m-3] with INTEGRAL_BT_CONTINUITY.

  • % uh_crve [real,private] :: The curvature of face area with velocity for flow from the east [H T2 L-1 ~> s2 or kg s2 m-3] or [H L-1 ~> nondim or kg m-3] with INTEGRAL_BT_CONTINUITY.

  • % uh_ww [real,private] :: The zonal transport when ubt=ubt_WW [H L2 T-1 ~> m3 s-1 or kg s-1], or the equivalent time-integrated transport with INTEGRAL_BT_CONTINUITY [H L2 ~> m3 or kg].

  • % uh_ee [real,private] :: The zonal transport when ubt=ubt_EE [H L2 T-1 ~> m3 s-1 or kg s-1], or the equivalent time-integrated transport with INTEGRAL_BT_CONTINUITY [H L2 ~> m3 or kg].

type mom_barotropic/local_bt_cont_v_type

A description of the functional dependence of transport at a v-point.

Type fields:
  • % fa_v_nn [real,private] :: The effective open face area for meridional barotropic transport drawing from locations far to the north [H L ~> m2 or kg m-1].

  • % fa_v_n0 [real,private] :: The effective open face area for meridional barotropic transport drawing from nearby to the north [H L ~> m2 or kg m-1].

  • % fa_v_s0 [real,private] :: The effective open face area for meridional barotropic transport drawing from nearby to the south [H L ~> m2 or kg m-1].

  • % fa_v_ss [real,private] :: The effective open face area for meridional barotropic transport drawing from locations far to the south [H L ~> m2 or kg m-1].

  • % vbt_ss [real,private] :: vBT_SS is the barotropic velocity [L T-1 ~> m s-1], or with INTEGRAL_BT_CONTINUITY the time-integrated barotropic velocity [L ~> m], beyond which the marginal open face area is FA_v_SS. vBT_SS must be non-negative.

  • % vbt_nn [real,private] :: vBT_NN is the barotropic velocity [L T-1 ~> m s-1], or with INTEGRAL_BT_CONTINUITY the time-integrated barotropic velocity [L ~> m], beyond which the marginal open face area is FA_v_NN. vBT_NN must be non-positive.

  • % vh_crvs [real,private] :: The curvature of face area with velocity for flow from the south [H T2 L-1 ~> s2 or kg s2 m-3] or [H L-1 ~> nondim or kg m-3] with INTEGRAL_BT_CONTINUITY.

  • % vh_crvn [real,private] :: The curvature of face area with velocity for flow from the north [H T2 L-1 ~> s2 or kg s2 m-3] or [H L-1 ~> nondim or kg m-3] with INTEGRAL_BT_CONTINUITY.

  • % vh_ss [real,private] :: The meridional transport when vbt=vbt_SS [H L2 T-1 ~> m3 s-1 or kg s-1], or the equivalent time-integrated transport with INTEGRAL_BT_CONTINUITY [H L2 ~> m3 or kg].

  • % vh_nn [real,private] :: The meridional transport when vbt=vbt_NN [H L2 T-1 ~> m3 s-1 or kg s-1], or the equivalent time-integrated transport with INTEGRAL_BT_CONTINUITY [H L2 ~> m3 or kg].

type mom_barotropic/memory_size_type

A container for passing around active tracer point memory limits.

Type fields:
  • % isdw [integer,private] :: Currently active memory limits.

  • % iedw [integer,private] :: Currently active memory limits.

  • % jsdw [integer,private] :: Currently active memory limits.

  • % jedw [integer,private] :: Currently active memory limits.

Function/Subroutine Documentation

subroutine mom_barotropic/btstep(U_in, V_in, eta_in, dt, bc_accel_u, bc_accel_v, forces, pbce, eta_PF_in, U_Cor, V_Cor, accel_layer_u, accel_layer_v, eta_out, uhbtav, vhbtav, G, GV, US, CS, visc_rem_u, visc_rem_v, SpV_avg, ADp, OBC, BT_cont, eta_PF_start, taux_bot, tauy_bot, uh0, vh0, u_uh0, v_vh0, etaav)

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 :: [inout] The ocean’s grid structure.

  • gv :: [in] The ocean’s vertical grid structure.

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

  • u_in :: [in] The initial (3-D) zonal velocity [L T-1 ~> m s-1].

  • v_in :: [in] The initial (3-D) meridional velocity [L T-1 ~> m s-1].

  • eta_in :: [in] The initial barotropic free surface height anomaly or column mass anomaly [H ~> m or kg m-2].

  • dt :: [in] The time increment to integrate over [T ~> s].

  • bc_accel_u :: [in] The zonal baroclinic accelerations, [L T-2 ~> m s-2].

  • bc_accel_v :: [in] The meridional baroclinic accelerations, [L T-2 ~> m s-2].

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

  • pbce :: [in] 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 :: [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 :: [in] The (3-D) zonal velocities used to calculate the Coriolis terms in bc_accel_u [L T-1 ~> m s-1].

  • v_cor :: [in] The (3-D) meridional velocities used to calculate the Coriolis terms in bc_accel_u [L T-1 ~> m s-1].

  • accel_layer_u :: [out] The zonal acceleration of each layer due to the barotropic calculation [L T-2 ~> m s-2].

  • accel_layer_v :: [out] The meridional acceleration of each layer due to the barotropic calculation [L T-2 ~> m s-2].

  • eta_out :: [out] The final barotropic free surface height anomaly or column mass anomaly [H ~> m or kg m-2].

  • uhbtav :: [out] the barotropic zonal volume or mass fluxes averaged through the barotropic steps [H L2 T-1 ~> m3 s-1 or kg s-1].

  • vhbtav :: [out] the barotropic meridional volume or mass fluxes averaged through the barotropic steps [H L2 T-1 ~> m3 s-1 or kg s-1].

  • cs :: [inout] Barotropic control structure

  • visc_rem_u :: [in] 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 [nondim]. Visc_rem_u is between 0 (at the bottom) and 1 (far above).

  • visc_rem_v :: [in] Ditto for meridional direction [nondim].

  • spv_avg :: [in] The column average specific volume, used in non-Boussinesq OBC calculations [R-1 ~> m3 kg-1]

  • adp :: Acceleration diagnostic pointers

  • 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 [R L Z T-2 ~> Pa].

  • tauy_bot :: The meridional bottom frictional stress from ocean to the seafloor [R L Z T-2 ~> 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]

  • etaav :: [out] The free surface height or column mass averaged over the barotropic integration [H ~> m or kg m-2].

Call to:

adjust_local_bt_cont_types apply_velocity_obcs bt_cont_to_face_areas mom_diag_mediator::enable_averages mom_diag_mediator::enable_averaging find_duhbt_du find_dvhbt_dv find_face_areas find_uhbt find_vhbt id_clock_calc id_clock_calc_post id_clock_calc_pre id_clock_pass_post id_clock_pass_pre id_clock_pass_step mom_error_handler::mom_error mom_error_handler::mom_mesg mom_open_boundary::obc_direction_n mom_open_boundary::obc_direction_s set_local_bt_cont_types set_up_bt_obc subroundoff swap

subroutine mom_barotropic/set_dtbt(G, GV, US, CS, eta, pbce, BT_cont, gtot_est, SSH_add)

This subroutine automatically determines an optimal value for dtbt based on some state of the ocean.

Parameters:
  • g :: [inout] The ocean’s grid structure.

  • gv :: [in] The ocean’s vertical grid structure.

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

  • cs :: [inout] Barotropic control structure

  • eta :: [in] The barotropic free surface height anomaly or column mass anomaly [H ~> m or kg m-2].

  • pbce :: [in] 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 :: [in] An estimate of the total gravitational acceleration [L2 H-1 T-2 ~> m s-2 or m4 kg-1 s-2].

  • ssh_add :: [in] An additional contribution to SSH to provide a margin of error when calculating the external wave speed [Z ~> m].

Call to:

bt_cont_to_face_areas mom_checksums::chksum0 find_face_areas id_clock_sync mom_error_handler::mom_error

Called from:

barotropic_init

subroutine mom_barotropic/apply_velocity_obcs(OBC, ubt, vbt, uhbt, vhbt, ubt_trans, vbt_trans, eta, SpV_avg, ubt_old, vbt_old, BT_OBC, G, MS, GV, US, CS, halo, dtbt, bebt, use_BT_cont, integral_BT_cont, dt_elapsed, Datu, Datv, BTCL_u, BTCL_v, uhbt0, vhbt0, ubt_int, vbt_int, uhbt_int, vhbt_int)

The following 4 subroutines apply the open boundary conditions. This subroutine applies the open boundary conditions on barotropic velocities and mass transports, as developed by Mehmet Ilicak.

Parameters:
  • obc :: An associated pointer to an OBC type.

  • g :: [in] The ocean’s grid structure.

  • ms :: [in] A type that describes the memory sizes of the argument arrays.

  • ubt :: [inout] the zonal barotropic velocity [L T-1 ~> m s-1].

  • uhbt :: [inout] the zonal barotropic transport [H L2 T-1 ~> m3 s-1 or kg s-1].

  • ubt_trans :: [inout] The zonal barotropic velocity used in transport [L T-1 ~> m s-1].

  • vbt :: [inout] The meridional barotropic velocity [L T-1 ~> m s-1].

  • vhbt :: [inout] the meridional barotropic transport [H L2 T-1 ~> m3 s-1 or kg s-1].

  • vbt_trans :: [inout] the meridional BT velocity used in transports [L T-1 ~> m s-1].

  • eta :: [in] The barotropic free surface height anomaly or column mass anomaly [H ~> m or kg m-2].

  • spv_avg :: [in] The column average specific volume [R-1 ~> m3 kg-1]

  • ubt_old :: [in] The starting value of ubt in a barotropic step [L T-1 ~> m s-1].

  • vbt_old :: [in] The starting value of vbt in a barotropic step [L T-1 ~> m s-1].

  • bt_obc :: [in] A structure with the private barotropic arrays related to the open boundary conditions, set by set_up_BT_OBC.

  • gv :: [in] The ocean’s vertical grid structure.

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

  • cs :: [in] Barotropic control structure

  • halo :: [in] The extra halo size to use here.

  • dtbt :: [in] The time step [T ~> s].

  • bebt :: [in] The fractional weighting of the future velocity in determining the transport [nondim]

  • use_bt_cont :: [in] If true, use the BT_cont_types to calculate transports.

  • integral_bt_cont :: [in] If true, update the barotropic continuity equation directly from the initial condition using the time-integrated barotropic velocity.

  • dt_elapsed :: [in] The amount of time in the barotropic stepping that will have elapsed [T ~> s].

  • datu :: [in] A fixed estimate of the face areas at u points [H L ~> m2 or kg m-1].

  • datv :: [in] A fixed estimate of the face areas at v points [H L ~> m2 or kg m-1].

  • btcl_u :: [in] Structure of information used for a dynamic estimate of the face areas at u-points.

  • btcl_v :: [in] Structure of information used for a dynamic estimate of the face areas at v-points.

  • uhbt0 :: [in] A correction to the zonal transport so that the barotropic functions agree with the sum of the layer transports [H L2 T-1 ~> m3 s-1 or kg s-1].

  • vhbt0 :: [in] A correction to the meridional transport so that the barotropic functions agree with the sum of the layer transports [H L2 T-1 ~> m3 s-1 or kg s-1].

  • ubt_int :: [in] The time-integrated zonal barotropic velocity before this update [L T-1 ~> m s-1].

  • uhbt_int :: [in] The time-integrated zonal barotropic transport [H L2 T-1 ~> m3 s-1 or kg s-1].

  • vbt_int :: [in] The time-integrated meridional barotropic velocity before this update [L T-1 ~> m s-1].

  • vhbt_int :: [in] The time-integrated meridional barotropic transport [H L2 T-1 ~> m3 s-1 or kg s-1].

Call to:

find_uhbt find_vhbt mom_open_boundary::obc_direction_n mom_open_boundary::obc_direction_s

Called from:

btstep

subroutine mom_barotropic/set_up_bt_obc(OBC, eta, SpV_avg, BT_OBC, BT_Domain, G, GV, US, CS, MS, halo, use_BT_cont, integral_BT_cont, dt_baroclinic, Datu, Datv, BTCL_u, BTCL_v, dgeo_de)

This subroutine sets up the private structure used to apply the open boundary conditions, as developed by Mehmet Ilicak.

Parameters:
  • obc :: [inout] An associated pointer to an OBC type.

  • ms :: [in] A type that describes the memory sizes of the argument arrays.

  • eta :: [in] The barotropic free surface height anomaly or column mass anomaly [H ~> m or kg m-2].

  • spv_avg :: [in] The column average specific volume [R-1 ~> m3 kg-1]

  • bt_obc :: [inout] A structure with the private barotropic arrays related to the open boundary conditions, set by set_up_BT_OBC.

  • bt_domain :: [inout] MOM_domain_type associated with wide arrays

  • g :: [inout] The ocean’s grid structure.

  • gv :: [in] The ocean’s vertical grid structure.

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

  • cs :: [in] Barotropic control structure

  • halo :: [in] The extra halo size to use here.

  • use_bt_cont :: [in] If true, use the BT_cont_types to calculate transports.

  • integral_bt_cont :: [in] If true, update the barotropic continuity equation directly from the initial condition using the time-integrated barotropic velocity.

  • dt_baroclinic :: [in] The baroclinic timestep for this cycle of updates to the barotropic solver [T ~> s]

  • datu :: [in] A fixed estimate of the face areas at u points [H L ~> m2 or kg m-1].

  • datv :: [in] A fixed estimate of the face areas at v points [H L ~> m2 or kg m-1].

  • btcl_u :: [in] Structure of information used for a dynamic estimate of the face areas at u-points.

  • btcl_v :: [in] Structure of information used for a dynamic estimate of the face areas at v-points.

  • dgeo_de :: [in] The constant of proportionality between geopotential and sea surface height [nondim].

Call to:

mom_error_handler::mom_error mom_open_boundary::obc_direction_n mom_open_boundary::obc_direction_s uhbt_to_ubt vhbt_to_vbt

Called from:

btstep

subroutine mom_barotropic/destroy_bt_obc(BT_OBC)

Clean up the BT_OBC memory.

Parameters:

bt_obc :: [inout] A structure with the private barotropic arrays related to the open boundary conditions, set by set_up_BT_OBC.

Called from:

barotropic_end

subroutine mom_barotropic/btcalc(h, G, GV, CS, h_u, h_v, may_use_default, 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 :: [inout] The ocean’s grid structure.

  • gv :: [in] The ocean’s vertical grid structure.

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

  • cs :: [inout] Barotropic control structure

  • h_u :: [in] The specified effective thicknesses at u-points,

  • h_v :: [in] The specified effective thicknesses at v-points,

  • may_use_default :: [in] 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.

Call to:

arithmetic harmonic hybrid mom_error_handler::mom_error mom_open_boundary::obc_direction_n mom_open_boundary::obc_direction_s

Called from:

barotropic_init

function mom_barotropic/find_uhbt(u, BTC) [real]

The function find_uhbt determines the zonal transport for a given velocity, or with INTEGRAL_BT_CONT=True it determines the time-integrated zonal transport for a given time-integrated velocity.

Parameters:
  • u :: [in] The local zonal velocity [L T-1 ~> m s-1] or time integrated velocity [L ~> m]

  • btc :: [in] A structure containing various fields that allow the barotropic transports to be calculated consistently with the layers’ continuity equations. The dimensions of some of the elements in this type vary depending on INTEGRAL_BT_CONT.

Return:

undefined :: The zonal barotropic transport [L2 H T-1 ~> m3 s-1] or time integrated transport [L2 H ~> m3]

Called from:

apply_velocity_obcs btstep

function mom_barotropic/find_duhbt_du(u, BTC) [real]

The function find_duhbt_du determines the marginal zonal face area for a given velocity, or with INTEGRAL_BT_CONT=True for a given time-integrated velocity.

Parameters:
  • u :: [in] The local zonal velocity [L T-1 ~> m s-1] or time integrated velocity [L ~> m]

  • btc :: [in] A structure containing various fields that allow the barotropic transports to be calculated consistently with the layers’ continuity equations. The dimensions of some of the elements in this type vary depending on INTEGRAL_BT_CONT.

Return:

undefined :: The zonal barotropic face area [L H ~> m2 or kg m-1]

Called from:

btstep

function mom_barotropic/uhbt_to_ubt(uhbt, BTC) [real]

This function inverts the transport function to determine the barotopic velocity that is consistent with a given transport, or if INTEGRAL_BT_CONT=True this finds the time-integrated velocity that is consistent with a time-integrated transport.

Parameters:
  • uhbt :: [in] The barotropic zonal transport that should be inverted for, [H L2 T-1 ~> m3 s-1 or kg s-1] or the time-integrated transport [H L2 ~> m3 or kg].

  • btc :: [in] A structure containing various fields that allow the barotropic transports to be calculated consistently with the layers’ continuity equations. The dimensions of some of the elements in this type vary depending on INTEGRAL_BT_CONT.

Return:

undefined :: The result - The velocity that gives uhbt transport [L T-1 ~> m s-1] or the time-integrated velocity [L ~> m].

Called from:

set_up_bt_obc

function mom_barotropic/find_vhbt(v, BTC) [real]

The function find_vhbt determines the meridional transport for a given velocity, or with INTEGRAL_BT_CONT=True it determines the time-integrated meridional transport for a given time-integrated velocity.

Parameters:
  • v :: [in] The local meridional velocity [L T-1 ~> m s-1] or time integrated velocity [L ~> m]

  • btc :: [in] A structure containing various fields that allow the barotropic transports to be calculated consistently with the layers’ continuity equations. The dimensions of some of the elements in this type vary depending on INTEGRAL_BT_CONT.

Return:

undefined :: The meridional barotropic transport [L2 H T-1 ~> m3 s-1] or time integrated transport [L2 H ~> m3]

Called from:

apply_velocity_obcs btstep

function mom_barotropic/find_dvhbt_dv(v, BTC) [real]

The function find_dvhbt_dv determines the marginal meridional face area for a given velocity, or with INTEGRAL_BT_CONT=True for a given time-integrated velocity.

Parameters:
  • v :: [in] The local meridional velocity [L T-1 ~> m s-1] or time integrated velocity [L ~> m]

  • btc :: [in] A structure containing various fields that allow the barotropic transports to be calculated consistently with the layers’ continuity equations. The dimensions of some of the elements in this type vary depending on INTEGRAL_BT_CONT.

Return:

undefined :: The meridional barotropic face area [L H ~> m2 or kg m-1]

Called from:

btstep

function mom_barotropic/vhbt_to_vbt(vhbt, BTC) [real]

This function inverts the transport function to determine the barotopic velocity that is consistent with a given transport, or if INTEGRAL_BT_CONT=True this finds the time-integrated velocity that is consistent with a time-integrated transport.

Parameters:
  • vhbt :: [in] The barotropic meridional transport that should be inverted for [H L2 T-1 ~> m3 s-1 or kg s-1] or the time-integrated transport [H L2 ~> m3 or kg].

  • btc :: [in] A structure containing various fields that allow the barotropic transports to be calculated consistently with the layers’ continuity equations. The dimensions of some of the elements in this type vary depending on INTEGRAL_BT_CONT.

Return:

undefined :: The result - The velocity that gives vhbt transport [L T-1 ~> m s-1] or the time-integrated velocity [L ~> m].

Called from:

set_up_bt_obc

subroutine mom_barotropic/set_local_bt_cont_types(BT_cont, BTCL_u, BTCL_v, G, US, MS, BT_Domain, halo, dt_baroclinic)

This subroutine sets up reordered versions of the BT_cont type in the local_BT_cont types, which have wide halos properly filled in.

Parameters:
  • bt_cont :: [inout] The BT_cont_type input to the barotropic solver

  • ms :: [in] A type that describes the memory sizes of the argument arrays

  • btcl_u :: [out] A structure with the u information from BT_cont

  • btcl_v :: [out] A structure with the v information from BT_cont

  • g :: [in] The ocean’s grid structure

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

  • bt_domain :: [inout] The domain to use for updating the halos of wide arrays

  • halo :: [in] The extra halo size to use here

  • dt_baroclinic :: [in] The baroclinic time step [T ~> s], which is provided if INTEGRAL_BT_CONTINUITY is true.

Call to:

id_clock_calc_pre id_clock_pass_pre swap

Called from:

btstep

subroutine mom_barotropic/adjust_local_bt_cont_types(ubt, uhbt, vbt, vhbt, BTCL_u, BTCL_v, G, US, MS, halo, dt_baroclinic)

Adjust_local_BT_cont_types expands the range of velocities with a cubic curve translating velocities into transports to match the initial values of velocities and summed transports when the velocities are larger than the first guesses of the cubic transition velocities used to set up the local_BT_cont types.

Parameters:
  • ms :: [in] A type that describes the memory sizes of the argument arrays.

  • ubt :: [in] The linearization zonal barotropic velocity [L T-1 ~> m s-1].

  • uhbt :: [in] The linearization zonal barotropic transport

  • vbt :: [in] The linearization meridional barotropic velocity [L T-1 ~> m s-1].

  • vhbt :: [in] The linearization meridional barotropic transport

  • btcl_u :: [out] A structure with the u information from BT_cont.

  • btcl_v :: [out] A structure with the v information from BT_cont.

  • g :: [in] The ocean’s grid structure.

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

  • halo :: [in] The extra halo size to use here.

  • dt_baroclinic :: [in] The baroclinic time step [T ~> s], which is provided if INTEGRAL_BT_CONTINUITY is true.

Called from:

btstep

subroutine mom_barotropic/bt_cont_to_face_areas(BT_cont, Datu, Datv, G, US, MS, halo)

This subroutine uses the BT_cont_type to find the maximum face areas, which can then be used for finding wave speeds, etc.

Parameters:
  • bt_cont :: [inout] The BT_cont_type input to the barotropic solver.

  • ms :: [in] A type that describes the memory sizes of the argument arrays.

  • datu :: [out] The effective zonal face area [H L ~> m2 or kg m-1].

  • datv :: [out] The effective meridional face area [H L ~> m2 or kg m-1].

  • g :: [in] The ocean’s grid structure.

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

  • halo :: [in] The extra halo size to use here.

Called from:

btstep set_dtbt

subroutine mom_barotropic/swap(a, b)

Swap the values of two real variables.

Parameters:
  • a :: [inout] The first variable to be swapped [arbitrary units]

  • b :: [inout] The second variable to be swapped [arbitrary units]

Called from:

btstep set_local_bt_cont_types

subroutine mom_barotropic/find_face_areas(Datu, Datv, G, GV, US, CS, MS, halo, eta, add_max)

This subroutine determines the open face areas of cells for calculating the barotropic transport.

Parameters:
  • ms :: [in] A type that describes the memory sizes of the argument arrays.

  • datu :: [out] The open zonal face area [H L ~> m2 or kg m-1].

  • datv :: [out] The open meridional face area [H L ~> m2 or kg m-1].

  • g :: [in] The ocean’s grid structure.

  • gv :: [in] The ocean’s vertical grid structure.

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

  • cs :: [in] Barotropic control structure

  • halo :: [in] The halo size to use, default = 1.

  • eta :: [in] The barotropic free surface height anomaly

  • add_max :: [in] A value to add to the maximum depth (used to overestimate the external wave speed) [Z ~> m].

Called from:

barotropic_init btstep set_dtbt

subroutine mom_barotropic/bt_mass_source(h, eta, set_cor, G, GV, 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 :: [in] The ocean’s grid structure.

  • gv :: [in] The ocean’s vertical grid structure.

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

  • eta :: [in] The free surface height that is to be corrected [H ~> m or kg m-2].

  • set_cor :: [in] 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 :: [inout] Barotropic control structure

Call to:

mom_error_handler::mom_error

subroutine mom_barotropic/barotropic_init(u, v, h, eta, Time, G, GV, US, param_file, diag, CS, restart_CS, calc_dtbt, BT_cont, SAL_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 :: [inout] The ocean’s grid structure.

  • gv :: [in] The ocean’s vertical grid structure.

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

  • u :: [in] The zonal velocity [L T-1 ~> m s-1].

  • v :: [in] The meridional velocity [L T-1 ~> m s-1].

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

  • eta :: [in] Free surface height or column mass anomaly

  • time :: [in] The current model time.

  • param_file :: [in] A structure to parse for run-time parameters.

  • diag :: [inout] A structure that is used to regulate diagnostic output.

  • cs :: [inout] Barotropic control structure

  • restart_cs :: [in] MOM restart control structure

  • calc_dtbt :: [out] If true, the barotropic time step must be recalculated before stepping.

  • bt_cont :: A structure with elements that describe the effective open face areas as a function of barotropic flow.

  • sal_csp :: A pointer to the control structure of the SAL module.

Call to:

arithmetic arithmetic_string bt_cont_string btcalc find_face_areas from_bt_cont harmonic harmonic_string hybrid hybrid_string id_clock_calc id_clock_calc_post id_clock_calc_pre id_clock_pass_post id_clock_pass_pre id_clock_pass_step id_clock_sync mom_error_handler::mom_error mom_error_handler::mom_mesg set_dtbt

subroutine mom_barotropic/barotropic_get_tav(CS, ubtav, vbtav, G, US)

Copies ubtav and vbtav from private type into arrays.

Parameters:
  • cs :: [in] Barotropic control structure

  • g :: [in] Grid structure

  • ubtav :: [inout] Zonal barotropic velocity averaged over a baroclinic timestep [L T-1 ~> m s-1]

  • vbtav :: [inout] Meridional barotropic velocity averaged over a baroclinic timestep [L T-1 ~> m s-1]

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

Called from:

mom_hor_visc::horizontal_viscosity

subroutine mom_barotropic/barotropic_end(CS)

Clean up the barotropic control structure.

Parameters:

cs :: [inout] Control structure to clear out.

Call to:

destroy_bt_obc

Called from:

mom_dynamics_split_rk2::end_dyn_split_rk2 mom_dynamics_split_rk2b::end_dyn_split_rk2b

subroutine mom_barotropic/register_barotropic_restarts(HI, GV, US, param_file, 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 :: [in] A horizontal index type structure.

  • gv :: [in] The ocean’s vertical grid structure.

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

  • param_file :: [in] A structure to parse for run-time parameters.

  • cs :: [inout] Barotropic control structure

  • restart_cs :: [inout] MOM restart control structure

Call to:

mom_io::var_desc