namespace mom_pressureforce_mont¶
Overview¶
Provides the Montgomery potential form of pressure gradient. More…
namespace mom_pressureforce_mont { // global functions subroutine, public pressureforce_mont_nonbouss( h h, tv tv, PFu PFu, PFv PFv, G G, GV GV, US US, CS CS, p_atm p_atm, pbce pbce, eta eta ); subroutine, public pressureforce_mont_bouss( h h, tv tv, PFu PFu, PFv PFv, G G, GV GV, US US, CS CS, p_atm p_atm, pbce pbce, eta eta ); subroutine, public set_pbce_bouss( e e, tv tv, G G, GV GV, US US, Rho0 Rho0, GFS_scale GFS_scale, pbce pbce, rho_star rho_star ); subroutine, public set_pbce_nonbouss( p p, tv tv, G G, GV GV, US US, GFS_scale GFS_scale, pbce pbce, alpha_star alpha_star ); subroutine, public pressureforce_mont_init( Time Time, G G, GV GV, US US, param_file param_file, diag diag, CS CS, tides_CSp tides_CSp ); subroutine, public pressureforce_mont_end(CS CS); } // namespace mom_pressureforce_mont
Detailed Documentation¶
Provides the Montgomery potential form of pressure gradient.
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.
Global Functions¶
subroutine, public pressureforce_mont_nonbouss( h h, tv tv, PFu PFu, PFv PFv, G G, GV GV, US US, CS CS, p_atm p_atm, pbce pbce, eta 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 |
Ocean grid structure. |
gv |
Vertical grid structure. |
us |
A dimensional unit scaling type |
h |
Layer thickness, [H ~> kg m-2]. |
tv |
Thermodynamic variables. |
pfu |
Zonal acceleration due to pressure gradients (equal to -dM/dx) [L T-2 ~> m s-2]. |
pfv |
Meridional acceleration due to pressure gradients (equal to -dM/dy) [L T-2 ~> m s-2]. |
cs |
Control structure for Montgomery potential PGF |
p_atm |
The pressure at the ice-ocean or atmosphere-ocean [Pa]. |
pbce |
The baroclinic pressure anomaly in |
eta |
Free surface height [H ~> kg m-1]. |
subroutine, public pressureforce_mont_bouss( h h, tv tv, PFu PFu, PFv PFv, G G, GV GV, US US, CS CS, p_atm p_atm, pbce pbce, eta 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 |
Ocean grid structure. |
gv |
Vertical grid structure. |
us |
A dimensional unit scaling type |
h |
Layer thickness [H ~> m]. |
tv |
Thermodynamic variables. |
pfu |
Zonal acceleration due to pressure gradients (equal to -dM/dx) [L T-2 ~> m s-2]. |
pfv |
Meridional acceleration due to pressure gradients (equal to -dM/dy) [L T-2 ~> m s2]. |
cs |
Control structure for Montgomery potential PGF |
p_atm |
The pressure at the ice-ocean or atmosphere-ocean [Pa]. |
pbce |
The baroclinic pressure anomaly in each layer due to free surface height anomalies [m2 s-2 H-1 ~> m s-2]. |
eta |
Free surface height [H ~> m]. |
subroutine, public set_pbce_bouss( e e, tv tv, G G, GV GV, US US, Rho0 Rho0, GFS_scale GFS_scale, pbce pbce, rho_star rho_star )
Determines the partial derivative of the acceleration due to pressure forces with the free surface height.
Parameters:
g |
Ocean grid structure |
gv |
Vertical grid structure |
e |
Interface height [Z ~> m]. |
tv |
Thermodynamic variables |
us |
A dimensional unit scaling type |
rho0 |
The “Boussinesq” ocean density [kg m-3]. |
gfs_scale |
Ratio between gravity applied to top interface and the gravitational acceleration of the planet [nondim]. Usually this ratio is 1. |
pbce |
The baroclinic pressure anomaly in each layer due |
rho_star |
The layer densities (maybe compressibility |
subroutine, public set_pbce_nonbouss( p p, tv tv, G G, GV GV, US US, GFS_scale GFS_scale, pbce pbce, alpha_star alpha_star )
Determines the partial derivative of the acceleration due to pressure forces with the column mass.
Parameters:
g |
Ocean grid structure |
gv |
Vertical grid structure |
p |
Interface pressures [Pa]. |
tv |
Thermodynamic variables |
us |
A dimensional unit scaling type |
gfs_scale |
Ratio between gravity applied to top interface and the gravitational acceleration of the planet [nondim]. Usually this ratio is 1. |
pbce |
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 |
The layer specific volumes (maybe compressibility compensated) [m3 kg-1]. |
subroutine, public pressureforce_mont_init( Time Time, G G, GV GV, US US, param_file param_file, diag diag, CS CS, tides_CSp tides_CSp )
Initialize the Montgomery-potential form of PGF control structure.
Parameters:
time |
Current model time |
g |
ocean grid structure |
gv |
Vertical grid structure |
us |
A dimensional unit scaling type |
param_file |
Parameter file handles |
diag |
Diagnostics control structure |
cs |
Montgomery PGF control structure |
tides_csp |
Tides control structure |
subroutine, public pressureforce_mont_end(CS CS)
Deallocates the Montgomery-potential form of PGF control structure.
Parameters:
cs |
Control structure for Montgomery potential PGF |