namespace mom_pressureforce_afv¶
Overview¶
Analytically integrated finite volume pressure gradient. More…
namespace mom_pressureforce_afv { // global functions subroutine, public pressureforce_afv( h h, tv tv, PFu PFu, PFv PFv, G G, GV GV, US US, CS CS, ALE_CSp ALE_CSp, p_atm p_atm, pbce pbce, eta eta ); subroutine, public pressureforce_afv_nonbouss( h h, tv tv, PFu PFu, PFv PFv, G G, GV GV, US US, CS CS, ALE_CSp ALE_CSp, p_atm p_atm, pbce pbce, eta eta ); subroutine, public pressureforce_afv_bouss( h h, tv tv, PFu PFu, PFv PFv, G G, GV GV, US US, CS CS, ALE_CSp ALE_CSp, p_atm p_atm, pbce pbce, eta eta ); subroutine, public pressureforce_afv_init( Time Time, G G, GV GV, US US, param_file param_file, diag diag, CS CS, tides_CSp tides_CSp ); subroutine, public pressureforce_afv_end(CS CS); } // namespace mom_pressureforce_afv
Detailed Documentation¶
Analytically integrated finite volume pressure gradient.
Provides the Boussinesq and non-Boussinesq forms of horizontal accelerations due to pressure gradients using a 2nd-order analytically vertically integrated finite volume form, as described by Adcroft et al., 2008.
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
Global Functions¶
subroutine, public pressureforce_afv( h h, tv tv, PFu PFu, PFv PFv, G G, GV GV, US US, CS CS, ALE_CSp ALE_CSp, p_atm p_atm, pbce pbce, eta eta )
Thin interface between the model and the Boussinesq and non-Boussinesq pressure force routines.
Parameters:
g |
Ocean grid structure |
gv |
Vertical grid structure |
us |
A dimensional unit scaling type |
h |
Layer thickness [H ~> m or kg m-2] |
tv |
Thermodynamic variables |
pfu |
Zonal acceleration [L T-2 ~> m s-2] |
pfv |
Meridional acceleration [L T-2 ~> m s-2] |
cs |
Finite volume PGF control structure |
ale_csp |
ALE control structure |
p_atm |
The pressure at the ice-ocean or atmosphere-ocean interface [Pa]. |
pbce |
The baroclinic pressure anomaly in each layer due to eta anomalies [m2 s-2 H-1 ~> m s-2 or m4 s-2 kg-1]. |
eta |
The bottom mass used to calculate PFu and PFv [H ~> m or kg m-2], with any tidal contributions or compressibility compensation. |
subroutine, public pressureforce_afv_nonbouss( h h, tv tv, PFu PFu, PFv PFv, G G, GV GV, US US, CS CS, ALE_CSp ALE_CSp, p_atm p_atm, pbce pbce, eta 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 |
Ocean grid structure |
gv |
Vertical grid structure |
us |
A dimensional unit scaling type |
h |
Layer thickness [H ~> kg/m2] |
tv |
Thermodynamic variables |
pfu |
Zonal acceleration [L T-2 ~> m s-2] |
pfv |
Meridional acceleration [L T-2 ~> m s-2] |
cs |
Finite volume PGF control structure |
ale_csp |
ALE control structure |
p_atm |
The pressure at the ice-ocean or atmosphere-ocean interface [Pa]. |
pbce |
The baroclinic pressure anomaly in each layer due to eta anomalies [m2 s-2 H-1 ~> m s-2 or m4 s-2 kg-1]. |
eta |
The bottom mass used to calculate PFu and PFv [H ~> m or kg m-2], with any tidal contributions or compressibility compensation. |
subroutine, public pressureforce_afv_bouss( h h, tv tv, PFu PFu, PFv PFv, G G, GV GV, US US, CS CS, ALE_CSp ALE_CSp, p_atm p_atm, pbce pbce, eta 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 |
Ocean grid structure |
gv |
Vertical grid structure |
us |
A dimensional unit scaling type |
h |
Layer thickness [H ~> m] |
tv |
Thermodynamic variables |
pfu |
Zonal acceleration [L T-2 ~> m s-2] |
pfv |
Meridional acceleration [L T-2 ~> m s-2] |
cs |
Finite volume PGF control structure |
ale_csp |
ALE control structure |
p_atm |
The pressure at the ice-ocean or atmosphere-ocean interface [Pa]. |
pbce |
The baroclinic pressure anomaly in each layer due to eta anomalies [m2 s-2 H-1 ~> m s-2 or m4 s-2 kg-1]. |
eta |
The bottom mass used to calculate PFu and PFv [H ~> m or kg m-2], with any tidal contributions or compressibility compensation. |
subroutine, public pressureforce_afv_init( Time Time, G G, GV GV, US US, param_file param_file, diag diag, CS CS, tides_CSp tides_CSp )
Initializes the finite volume pressure gradient 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 |
Finite volume PGF control structure |
tides_csp |
Tides control structure |
subroutine, public pressureforce_afv_end(CS CS)
Deallocates the finite volume pressure gradient control structure.
Parameters:
cs |
Finite volume pressure control structure that will be deallocated in this subroutine. |