MOM6
mom_pressureforce Module Reference

Detailed Description

A thin wrapper for Boussinesq/non-Boussinesq forms of the pressure force calculation.

This thin module provides a branch to two forms of the horizontal accelerations due to pressure gradients. The two options currently available are a Montgomery potential form (used in traditional isopycnal layer models), and the analytic finite volume form.

Data Types

type  pressureforce_cs
 

Functions/Subroutines

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. More...
 
subroutine, public pressureforce_init (Time, G, GV, param_file, diag, CS, tides_CSp)
 Initialize the pressure force control structure. More...
 
subroutine, public pressureforce_end (CS)
 Deallocate the pressure force control structure. More...
 

Function/Subroutine Documentation

◆ pressureforce()

subroutine, public mom_pressureforce::pressureforce ( real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in)  h,
type(thermo_var_ptrs), intent(in)  tv,
real, dimension( g %isdb: g %iedb, g %jsd: g %jed, g %ke), intent(out)  PFu,
real, dimension( g %isd: g %ied, g %jsdb: g %jedb, g %ke), intent(out)  PFv,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(pressureforce_cs), pointer  CS,
type(ale_cs), pointer  ALE_CSp,
real, dimension(:,:), optional, pointer  p_atm,
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(out), optional  pbce,
real, dimension( g %isd: g %ied, g %jsd: g %jed), intent(out), optional  eta 
)

A thin layer between the model and the Boussinesq and non-Boussinesq pressure force routines.

Parameters
[in]gThe ocean's grid structure
[in]gvThe ocean's vertical grid structure
[in]hLayer thicknesses, in H (usually m or kg m-2)
[in]tvA structure pointing to various thermodynamic variables

Definition at line 38 of file MOM_PressureForce.F90.

Referenced by mom_dynamics_legacy_split::step_mom_dyn_legacy_split(), mom_dynamics_split_rk2::step_mom_dyn_split_rk2(), mom_dynamics_unsplit::step_mom_dyn_unsplit(), and mom_dynamics_unsplit_rk2::step_mom_dyn_unsplit_rk2().

38  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
39  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure
40  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thicknesses, in H (usually m or kg m-2)
41  type(thermo_var_ptrs), intent(in) :: tv !< A structure pointing to various thermodynamic variables
42  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), intent(out) :: pfu
43  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), intent(out) :: pfv
44  type(pressureforce_cs), pointer :: cs
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
49 
50 
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)
55  else
56  call pressureforce_afv_nonbouss(h, tv, pfu, pfv, g, gv, cs%PressureForce_AFV_CSp, &
57  p_atm, pbce, eta)
58  endif
59  else
60  if (gv%Boussinesq) then
61  call pressureforce_mont_bouss(h, tv, pfu, pfv, g, gv, cs%PressureForce_Mont_CSp, &
62  p_atm, pbce, eta)
63  else
64  call pressureforce_mont_nonbouss(h, tv, pfu, pfv, g, gv, cs%PressureForce_Mont_CSp, &
65  p_atm, pbce, eta)
66  endif
67  endif
68 
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:19
Here is the caller graph for this function:

◆ pressureforce_end()

subroutine, public mom_pressureforce::pressureforce_end ( type(pressureforce_cs), pointer  CS)

Deallocate the pressure force control structure.

Parameters
csPressure force control structure

Definition at line 110 of file MOM_PressureForce.F90.

110  type(pressureforce_cs), pointer :: cs !< Pressure force control structure
111  if (associated(cs)) deallocate(cs)

◆ pressureforce_init()

subroutine, public mom_pressureforce::pressureforce_init ( type(time_type), intent(in), target  Time,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(param_file_type), intent(in)  param_file,
type(diag_ctrl), intent(inout), target  diag,
type(pressureforce_cs), pointer  CS,
type(tidal_forcing_cs), optional, pointer  tides_CSp 
)

Initialize the pressure force control structure.

Parameters
[in]timeCurrent model time
[in]gOcean grid structure
[in]gvVertical grid structure
[in]param_fileParameter file handles
[in,out]diagDiagnostics control structure
csPressure force control structure
tides_cspTide control structure

Definition at line 73 of file MOM_PressureForce.F90.

References mom_error_handler::mom_error(), mom_pressureforce_afv::pressureforce_afv_init(), and mom_pressureforce_mont::pressureforce_mont_init().

73  type(time_type), target, intent(in) :: time !< Current model time
74  type(ocean_grid_type), intent(in) :: g !< Ocean grid structure
75  type(verticalgrid_type), intent(in) :: gv !< Vertical grid structure
76  type(param_file_type), intent(in) :: param_file !< Parameter file handles
77  type(diag_ctrl), target, intent(inout) :: diag !< Diagnostics control structure
78  type(pressureforce_cs), pointer :: cs !< Pressure force control structure
79  type(tidal_forcing_cs), optional, pointer :: tides_csp !< Tide control structure
80 #include "version_variable.h"
81  character(len=40) :: mdl = "MOM_PressureForce" ! This module's name.
82 
83  if (associated(cs)) then
84  call mom_error(warning, "PressureForce_init called with an associated "// &
85  "control structure.")
86  return
87  else ; allocate(cs) ; endif
88 
89  ! Read all relevant parameters and write them to the model log.
90  call log_version(param_file, mdl, version, "")
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.)
97 
98  if (cs%Analytic_FV_PGF) then
99  call pressureforce_afv_init(time, g, gv, param_file, diag, &
100  cs%PressureForce_AFV_CSp, tides_csp)
101  else
102  call pressureforce_mont_init(time, g, gv, param_file, diag, &
103  cs%PressureForce_Mont_CSp, tides_csp)
104  endif
105 
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:19
Here is the call graph for this function: