mom_pressureforce_mont module reference¶
Provides the Montgomery potential form of pressure gradient.
Data Types¶
Control structure for the Montgomery potential form of pressure gradient. |
Functions/Subroutines¶
Non-Boussinesq Montgomery-potential form of pressure gradient. |
|
Boussinesq Montgomery-potential form of pressure gradient. |
|
Determines the partial derivative of the acceleration due to pressure forces with the free surface height. |
|
Determines the partial derivative of the acceleration due to pressure forces with the column mass. |
|
Initialize the Montgomery-potential form of PGF control structure. |
Detailed Description¶
Provides the Boussunesq and non-Boussinesq forms of the horizontal accelerations due to pressure gradients using the Montgomery potential. A second-order accurate, centered scheme is used. If a split time stepping scheme is used, the vertical decomposition into barotropic and baroclinic contributions described by Hallberg (J Comp Phys 1997) is used. With a nonlinear equation of state, compressibility is added along the lines proposed by Sun et al. (JPO 1999), but with compressibility coefficients based on a fit to a user-provided reference profile.
Type Documentation¶
-
type
mom_pressureforce_mont/
pressureforce_mont_cs
¶ Control structure for the Montgomery potential form of pressure gradient.
- Type fields:
%
id_pfu_bc
[integer] :: Diagnostic IDs.%
id_pfv_bc
[integer] :: Diagnostic IDs.%
id_e_sal
[integer] :: Diagnostic IDs.%
id_e_tide
[integer] :: Diagnostic IDs.%
id_e_tide_eq
[integer] :: Diagnostic IDs.%
id_e_tide_sal
[integer] :: Diagnostic IDs.%
initialized
[logical] :: True if this control structure has been initialized.%
calculate_sal
[logical] :: If true, calculate self-attraction and loading.%
tides
[logical] :: If true, apply tidal momentum forcing.%
rho0
[real] :: The density used in the Boussinesq approximation [R ~> kg m-3].%
gfs_scale
[real] :: Ratio between gravity applied to top interface and the gravitational acceleration of the planet [nondim]. Usually this ratio is 1.%
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.%
pfu_bc
[real(:,:,:),allocatable] :: Zonal accelerations due to pressure gradients deriving from density gradients within layers [L T-2 ~> m s-2].%
pfv_bc
[real(:,:,:),allocatable] :: Meridional accelerations due to pressure gradients deriving from density gradients within layers [L T-2 ~> m s-2].%
sal_csp
[type( sal_cs ),pointer] :: SAL control structure.%
tides_csp
[type( tidal_forcing_cs ),pointer] :: The tidal forcing control structure.
Function/Subroutine Documentation¶
-
subroutine
mom_pressureforce_mont/
pressureforce_mont_nonbouss
(h, tv, PFu, PFv, G, GV, US, CS, p_atm, pbce, eta)¶ Non-Boussinesq Montgomery-potential form of pressure gradient.
Determines the acceleration due to pressure forces in a non-Boussinesq fluid using the compressibility compensated (if appropriate) Montgomery-potential form described in Hallberg (Ocean Mod., 2005).
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 due to pressure gradients (equal to -dM/dx) [L T-2 ~> m s-2].
pfv :: [out] Meridional acceleration due to pressure gradients (equal to -dM/dy) [L T-2 ~> m s-2].
cs :: [inout] Control structure for Montgomery potential PGF
p_atm :: The pressure at the ice-ocean or atmosphere-ocean [R L2 T-2 ~> Pa].
pbce :: [out] The baroclinic pressure anomaly in
eta :: [out] The total column mass used to calculate PFu and PFv [H ~> kg m-2].
- Call to:
mom_self_attr_load::calc_sal
mom_tidal_forcing::calc_tidal_forcing
mom_density_integrals::int_specific_vol_dp
mom_error_handler::mom_error
mom_eos::query_compressible
set_pbce_nonbouss
-
subroutine
mom_pressureforce_mont/
pressureforce_mont_bouss
(h, tv, PFu, PFv, G, GV, US, CS, p_atm, pbce, eta)¶ Boussinesq Montgomery-potential form of pressure gradient.
Determines the acceleration due to pressure forces.
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 due to pressure gradients (equal to -dM/dx) [L T-2 ~> m s-2].
pfv :: [out] Meridional acceleration due to pressure gradients (equal to -dM/dy) [L T-2 ~> m s-2].
cs :: [inout] Control structure for Montgomery potential PGF
p_atm :: The pressure at the ice-ocean or atmosphere-ocean [R L2 T-2 ~> Pa].
pbce :: [out] The baroclinic pressure anomaly in each layer due to free surface height anomalies [L2 T-2 H-1 ~> m s-2].
eta :: [out] Free surface height [H ~> m].
- Call to:
mom_self_attr_load::calc_sal
mom_tidal_forcing::calc_tidal_forcing
mom_error_handler::mom_error
mom_eos::query_compressible
set_pbce_bouss
-
subroutine
mom_pressureforce_mont/
set_pbce_bouss
(e, tv, G, GV, US, Rho0, GFS_scale, pbce, rho_star)¶ Determines the partial derivative of the acceleration due to pressure forces with the free surface height.
- Parameters:
g :: [in] Ocean grid structure
gv :: [in] Vertical grid structure
e :: [in] Interface height [Z ~> m].
tv :: [in] Thermodynamic variables
us :: [in] A dimensional unit scaling type
rho0 :: [in] The “Boussinesq” ocean density [R ~> kg m-3].
gfs_scale :: [in] Ratio between gravity applied to top interface and the gravitational acceleration of the planet [nondim]. Usually this ratio is 1.
pbce :: [out] The baroclinic pressure anomaly in each layer due
rho_star :: [in] The layer densities (maybe compressibility
- Called from:
mom_pressureforce_fv::pressureforce_fv_bouss
pressureforce_mont_bouss
-
subroutine
mom_pressureforce_mont/
set_pbce_nonbouss
(p, tv, G, GV, US, GFS_scale, pbce, alpha_star)¶ Determines the partial derivative of the acceleration due to pressure forces with the column mass.
- Parameters:
g :: [in] Ocean grid structure
gv :: [in] Vertical grid structure
p :: [in] Interface pressures [R L2 T-2 ~> Pa].
tv :: [in] Thermodynamic variables
us :: [in] A dimensional unit scaling type
gfs_scale :: [in] Ratio between gravity applied to top interface and the gravitational acceleration of the planet [nondim]. Usually this ratio is 1.
pbce :: [out] The baroclinic pressure anomaly in each layer due to free surface height anomalies [L2 H-1 T-2 ~> m4 kg-1 s-2].
alpha_star :: [in] The layer specific volumes (maybe compressibility compensated) [R-1 ~> m3 kg-1].
- Called from:
mom_pressureforce_fv::pressureforce_fv_nonbouss
pressureforce_mont_nonbouss
-
subroutine
mom_pressureforce_mont/
pressureforce_mont_init
(Time, G, GV, US, param_file, diag, CS, SAL_CSp, tides_CSp)¶ Initialize the Montgomery-potential form of PGF 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] Montgomery PGF control structure
sal_csp :: [in] SAL control structure
tides_csp :: [in] Tides control structure