mom_pressureforce_fv module reference¶
Finite volume pressure gradient (integrated by quadrature or analytically)
Data Types¶
Finite volume pressure gradient control structure. |
Functions/Subroutines¶
Non-Boussinesq analytically-integrated finite volume form of pressure gradient. |
|
Boussinesq analytically-integrated finite volume form of pressure gradient. |
|
Initializes the finite volume pressure gradient control structure. |
Detailed Description¶
Provides the Boussinesq and non-Boussinesq forms of horizontal accelerations due to pressure gradients using a vertically integrated finite volume form, as described by Adcroft et al., 2008. Integration in the vertical is made either by quadrature or analytically.
This form eliminates the thermobaric instabilities that had been a problem with previous forms of the pressure gradient force calculation, as described by Hallberg, 2005.
Adcroft, A., R. Hallberg, and M. Harrison, 2008: A finite volume discretization of the pressure gradient force using analytic integration. Ocean Modelling, 22, 106-113. http://doi.org/10.1016/j.ocemod.2008.02.001
Hallberg, 2005: A thermobaric instability of Lagrangian vertical coordinate ocean models. Ocean Modelling, 8, 279-300. http://dx.doi.org/10.1016/j.ocemod.2004.01.001
Type Documentation¶
-
type
mom_pressureforce_fv/pressureforce_fv_cs¶ Finite volume pressure gradient control structure.
- Type fields:
%initialized[logical] :: True if this control structure has been initialized.%calculate_sal[logical] :: If true, calculate self-attraction and loading.%sal_use_bpa[logical] :: If true, use bottom pressure anomaly instead of SSH to calculate SAL.%tides[logical] :: If true, apply tidal momentum forcing.%rho_ref[real] :: The reference density that is subtracted off when calculating pressure gradient forces [R ~> kg m-3].%rho_ref_bug[logical] :: If true, recover a bug that mixes GVRho0 and CSrho_ref in Boussinesq mode.%gfs_scale[real] :: A scaling of the surface pressure gradients to allow the use of a reduced gravity model [nondim].%time[type(time_type),pointer] :: A pointer to the ocean model’s clock.%diag[type( diag_ctrl ),pointer] :: A structure that is used to regulate the timing of diagnostic output.%masswghtinterp[integer] :: A flag indicating whether and how to use mass weighting in T/S interpolation.%correction_intxpa[logical] :: If true, apply a correction to the value of intxpa at a selected interface under ice, using matching at the end values along with a 5-point quadrature integral of the hydrostatic pressure or height changes along that interface. The selected interface is either at the ocean’s surface or in the interior, depending on reset_intxpa_integral.%reset_intxpa_integral[logical] :: If true and the surface displacement between adjacent cells exceeds the vertical grid spacing, reset intxpa at the interface below a trusted interior cell. (This often applies in ice shelf cavities.)%masswghtinterpvanonly[logical] :: If true, don’t do mass weighting of T/S interpolation unless vanished.%reset_intxpa_flattest[logical] :: If true, use flattest interface rather than top for reset integral in cases where no best nonvanished interface.%h_nonvanished[real] :: A minimal layer thickness that indicates that a layer is thick enough to usefully reestimate the pressure integral across the interface below it [H ~> m or kg m-2].%use_inaccurate_pgf_rho_anom[logical] :: If true, uses the older and less accurate method to calculate density anomalies, as used prior to March 2018.%boundary_extrap[logical] :: Indicate whether high-order boundary extrapolation should be used within boundary cells.%reconstruct[logical] :: If true, polynomial profiles of T & S will be reconstructed and used in the integrals for the finite volume pressure gradient calculation. The default depends on whether regridding is being used.%recon_scheme[integer] :: Order of the polynomial of the reconstruction of T & S for the finite volume pressure gradient calculation. By the default (1) is for a piecewise linear method.%debug[logical] :: If true, write verbose checksums for debugging purposes.%use_ssh_in_z0p[logical] :: If true, adjust the height at which the pressure used in the equation of state is 0 to account for the displacement of the sea surface including adjustments for atmospheric or sea-ice pressure.%use_stanley_pgf[logical] :: If true, turn on Stanley parameterization in the PGF.%bq_sal_tides[logical] :: If true, use an alternative method for SAL and tides in Boussinesq mode.%tides_answer_date[integer] :: Recover old answers with tides.%id_e_tide[integer] :: Diagnostic identifier.%id_e_tidal_eq[integer] :: Diagnostic identifier.%id_e_tidal_sal[integer] :: Diagnostic identifier.%id_e_sal[integer] :: Diagnostic identifier.%id_rho_pgf[integer] :: Diagnostic identifier.%id_rho_stanley_pgf[integer] :: Diagnostic identifier.%id_p_stanley[integer] :: Diagnostic identifier.%id_masswt_u[integer] :: Diagnostic identifier.%id_masswt_v[integer] :: Diagnostic identifier.%id_sal_u[integer] :: Diagnostic identifier.%id_sal_v[integer] :: Diagnostic identifier.%id_tides_u[integer] :: Diagnostic identifier.%id_tides_v[integer] :: Diagnostic identifier.%sal_csp[type( sal_cs ),pointer] :: SAL control structure.%tides_csp[type( tidal_forcing_cs ),pointer] :: Tides control structure.
Function/Subroutine Documentation¶
-
subroutine
mom_pressureforce_fv/pressureforce_fv_nonbouss(h, tv, PFu, PFv, G, GV, US, CS, ALE_CSp, ADp, p_atm, pbce, eta)¶ Non-Boussinesq analytically-integrated finite volume form of pressure gradient.
Determines the acceleration due to hydrostatic pressure forces, using the analytic finite volume form of the Pressure gradient, and does not make the Boussinesq approximation.
To work, the following fields must be set outside of the usual (is:ie,js:je) range before this subroutine is called: h(isB:ie+1,jsB:je+1), T(isB:ie+1,jsB:je+1), and S(isB:ie+1,jsB:je+1).
- Parameters:
g :: [in] Ocean grid structure
gv :: [in] Vertical grid structure
us :: [in] A dimensional unit scaling type
h :: [in] Layer thickness [H ~> kg m-2]
tv :: [in] Thermodynamic variables
pfu :: [out] Zonal acceleration [L T-2 ~> m s-2]
pfv :: [out] Meridional acceleration [L T-2 ~> m s-2]
cs :: [in] Finite volume PGF control structure
ale_csp :: ALE control structure
adp :: Acceleration diagnostic pointers
p_atm :: The pressure at the ice-ocean or atmosphere-ocean interface [R L2 T-2 ~> Pa].
pbce :: [out] The baroclinic pressure anomaly in each layer due to eta anomalies [L2 T-2 H-1 ~> m4 s-2 kg-1].
eta :: [out] The total column mass used to calculate PFu and PFv [H ~> kg m-2].
- Call to:
mom_self_attr_load::calc_salmom_tidal_forcing::calc_tidal_forcing_legacymom_density_integrals::diagnose_mass_weight_pmom_error_handler::mom_errormom_pressureforce_mont::set_pbce_nonboussmom_ale::ts_plm_edge_valuesmom_ale::ts_ppm_edge_values
-
subroutine
mom_pressureforce_fv/pressureforce_fv_bouss(h, tv, PFu, PFv, G, GV, US, CS, ALE_CSp, ADp, p_atm, pbce, eta)¶ Boussinesq analytically-integrated finite volume form of pressure gradient.
Determines the acceleration due to hydrostatic pressure forces, using the finite volume form of the terms and analytic integrals in depth.
To work, the following fields must be set outside of the usual (is:ie,js:je) range before this subroutine is called: h(isB:ie+1,jsB:je+1), T(isB:ie+1,jsB:je+1), and S(isB:ie+1,jsB:je+1).
- Parameters:
g :: [in] Ocean grid structure
gv :: [in] Vertical grid structure
us :: [in] A dimensional unit scaling type
h :: [in] Layer thickness [H ~> m]
tv :: [in] Thermodynamic variables
pfu :: [out] Zonal acceleration [L T-2 ~> m s-2]
pfv :: [out] Meridional acceleration [L T-2 ~> m s-2]
cs :: [in] Finite volume PGF control structure
ale_csp :: ALE control structure
adp :: Acceleration diagnostic pointers
p_atm :: The pressure at the ice-ocean or atmosphere-ocean interface [R L2 T-2 ~> Pa].
pbce :: [out] The baroclinic pressure anomaly in each layer due to eta anomalies [L2 T-2 H-1 ~> m s-2].
eta :: [out] The sea-surface height used to calculate PFu and PFv [H ~> m], with any tidal contributions.
- Call to:
mom_self_attr_load::calc_salmom_tidal_forcing::calc_tidal_forcing_legacymom_density_integrals::diagnose_mass_weight_zmom_eos::eos_domainmom_error_handler::mom_errormom_pressureforce_mont::set_pbce_boussmom_ale::ts_plm_edge_valuesmom_ale::ts_ppm_edge_values
-
subroutine
mom_pressureforce_fv/pressureforce_fv_init(Time, G, GV, US, param_file, diag, CS, ADp, SAL_CSp, tides_CSp)¶ Initializes the finite volume pressure gradient control structure.
- Parameters:
time :: [in] Current model time
g :: [in] Ocean grid structure
gv :: [in] Vertical grid structure
us :: [in] A dimensional unit scaling type
param_file :: [in] Parameter file handles
diag :: [inout] Diagnostics control structure
cs :: [inout] Finite volume PGF control structure
adp :: Acceleration diagnostic pointers
sal_csp :: [in] SAL control structure
tides_csp :: [in] Tides control structure
- Call to: