18 implicit none ;
private 20 #include <MOM_memory.h> 26 logical :: analytic_fv_pgf
37 subroutine pressureforce(h, tv, PFu, PFv, G, GV, CS, ALE_CSp, p_atm, pbce, eta)
40 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: h
42 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)),
intent(out) :: PFu
43 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)),
intent(out) :: PFv
45 type(
ale_cs),
pointer :: ALE_CSp
46 real,
dimension(:,:),
optional,
pointer :: p_atm
47 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
optional,
intent(out) :: pbce
48 real,
dimension(SZI_(G),SZJ_(G)),
optional,
intent(out) :: eta
51 if (cs%Analytic_FV_PGF)
then 52 if (gv%Boussinesq)
then 53 call pressureforce_afv_bouss(h, tv, pfu, pfv, g, gv, cs%PressureForce_AFV_CSp, &
54 ale_csp, p_atm, pbce, eta)
56 call pressureforce_afv_nonbouss(h, tv, pfu, pfv, g, gv, cs%PressureForce_AFV_CSp, &
60 if (gv%Boussinesq)
then 61 call pressureforce_mont_bouss(h, tv, pfu, pfv, g, gv, cs%PressureForce_Mont_CSp, &
64 call pressureforce_mont_nonbouss(h, tv, pfu, pfv, g, gv, cs%PressureForce_Mont_CSp, &
73 type(time_type),
target,
intent(in) :: Time
77 type(
diag_ctrl),
target,
intent(inout) :: diag
80 #include "version_variable.h" 81 character(len=40) :: mdl =
"MOM_PressureForce" 83 if (
associated(cs))
then 84 call mom_error(warning,
"PressureForce_init called with an associated "// &
87 else ;
allocate(cs) ;
endif 91 call get_param(param_file, mdl,
"ANALYTIC_FV_PGF", cs%Analytic_FV_PGF, &
92 "If true the pressure gradient forces are calculated \n"//&
93 "with a finite volume form that analytically integrates \n"//&
94 "the equations of state in pressure to avoid any \n"//&
95 "possibility of numerical thermobaric instability, as \n"//&
96 "described in Adcroft et al., O. Mod. (2008).", default=.true.)
98 if (cs%Analytic_FV_PGF)
then 100 cs%PressureForce_AFV_CSp, tides_csp)
103 cs%PressureForce_Mont_CSp, tides_csp)
111 if (
associated(cs))
deallocate(cs)
subroutine, public pressureforce_mont_bouss(h, tv, PFu, PFv, G, GV, CS, p_atm, pbce, eta)
Boussinesq Montgomery-potential form of pressure gradient.
subroutine, public calc_tidal_forcing(Time, eta, eta_tidal, G, CS, deta_tidal_deta)
This subroutine calculates the geopotential anomalies that drive the tides, including self-attraction...
This module contains the main regridding routines. Regridding comprises two steps: (1) Interpolation ...
subroutine, public pressureforce_mont_init(Time, G, GV, param_file, diag, CS, tides_CSp)
Initialize the Montgomery-potential form of PGF control structure.
Ocean grid type. See mom_grid for details.
Provides the ocean grid type.
subroutine, public pressureforce_init(Time, G, GV, param_file, diag, CS, tides_CSp)
Initialize the pressure force control structure.
A thin wrapper for Boussinesq/non-Boussinesq forms of the pressure force calculation.
subroutine, public pressureforce_afv_nonbouss(h, tv, PFu, PFv, G, GV, CS, p_atm, pbce, eta)
Non-Boussinesq analytically-integrated finite volume form of pressure gradient.
Analytically integrated finite volume pressure gradient.
Provides the Montgomery potential form of pressure gradient.
subroutine, public pressureforce_afv_init(Time, G, GV, param_file, diag, CS, tides_CSp)
Initializes the finite volume pressure gradient control structure.
Finite volume pressure gradient control structure.
logical function, public is_root_pe()
subroutine, public pressureforce(h, tv, PFu, PFv, G, GV, CS, ALE_CSp, p_atm, pbce, eta)
A thin layer between the model and the Boussinesq and non-Boussinesq pressure force routines...
subroutine, public pressureforce_mont_nonbouss(h, tv, PFu, PFv, G, GV, CS, p_atm, pbce, eta)
Non-Boussinesq Montgomery-potential form of pressure gradient.
subroutine, public mom_mesg(message, verb, all_print)
subroutine, public pressureforce_end(CS)
Deallocate the pressure force control structure.
subroutine, public pressureforce_afv_bouss(h, tv, PFu, PFv, G, GV, CS, ALE_CSp, p_atm, pbce, eta)
Boussinesq analytically-integrated finite volume form of pressure gradient.
The thermo_var_ptrs structure contains pointers to an assortment of thermodynamic fields that may be ...
Control structure for the Montgomery potential form of pressure gradient.
subroutine, public mom_error(level, message, all_print)