MOM6
MOM_PressureForce.F90
1 !> A thin wrapper for Boussinesq/non-Boussinesq forms of the pressure force calculation.
3 
4 ! This file is part of MOM6. See LICENSE.md for the license.
5 
6 use mom_diag_mediator, only : diag_ctrl, time_type
7 use mom_error_handler, only : mom_error, mom_mesg, fatal, warning, is_root_pe
9 use mom_grid, only : ocean_grid_type
10 use mom_pressureforce_fv, only : pressureforce_fv_bouss, pressureforce_fv_nonbouss
11 use mom_pressureforce_fv, only : pressureforce_fv_init, pressureforce_fv_end
13 use mom_pressureforce_mont, only : pressureforce_mont_bouss, pressureforce_mont_nonbouss
14 use mom_pressureforce_mont, only : pressureforce_mont_init, pressureforce_mont_end
20 use mom_ale, only: ale_cs
21 implicit none ; private
22 
23 #include <MOM_memory.h>
24 
25 public pressureforce, pressureforce_init, pressureforce_end
26 
27 !> Pressure force control structure
28 type, public :: pressureforce_cs ; private
29  logical :: analytic_fv_pgf !< If true, use the analytic finite volume form
30  !! (Adcroft et al., Ocean Mod. 2008) of the PGF.
31  !> Control structure for the analytically integrated finite volume pressure force
32  type(pressureforce_fv_cs), pointer :: pressureforce_fv_csp => null()
33  !> Control structure for the Montgomery potential form of pressure force
34  type(pressureforce_mont_cs), pointer :: pressureforce_mont_csp => null()
35 end type pressureforce_cs
36 
37 contains
38 
39 !> A thin layer between the model and the Boussinesq and non-Boussinesq pressure force routines.
40 subroutine pressureforce(h, tv, PFu, PFv, G, GV, US, CS, ALE_CSp, p_atm, pbce, eta)
41  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
42  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
43  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
44  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
45  intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
46  type(thermo_var_ptrs), intent(in) :: tv !< A structure pointing to various thermodynamic variables
47  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
48  intent(out) :: PFu !< Zonal pressure force acceleration [L T-2 ~> m s-2]
49  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
50  intent(out) :: PFv !< Meridional pressure force acceleration [L T-2 ~> m s-2]
51  type(pressureforce_cs), pointer :: CS !< Pressure force control structure
52  type(ale_cs), pointer :: ALE_CSp !< ALE control structure
53  real, dimension(:,:), &
54  optional, pointer :: p_atm !< The pressure at the ice-ocean or
55  !! atmosphere-ocean interface [R L2 T-2 ~> Pa].
56  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
57  optional, intent(out) :: pbce !< The baroclinic pressure anomaly in each layer
58  !! due to eta anomalies [L2 T-2 H-1 ~> m s-2 or m4 s-2 kg-1].
59  real, dimension(SZI_(G),SZJ_(G)), &
60  optional, intent(out) :: eta !< The bottom mass used to calculate PFu and PFv,
61  !! [H ~> m or kg m-2], with any tidal contributions.
62 
63  if (cs%Analytic_FV_PGF) then
64  if (gv%Boussinesq) then
65  call pressureforce_fv_bouss(h, tv, pfu, pfv, g, gv, us, cs%PressureForce_FV_CSp, &
66  ale_csp, p_atm, pbce, eta)
67  else
68  call pressureforce_fv_nonbouss(h, tv, pfu, pfv, g, gv, us, cs%PressureForce_FV_CSp, &
69  ale_csp, p_atm, pbce, eta)
70  endif
71  else
72  if (gv%Boussinesq) then
73  call pressureforce_mont_bouss(h, tv, pfu, pfv, g, gv, us, cs%PressureForce_Mont_CSp, &
74  p_atm, pbce, eta)
75  else
76  call pressureforce_mont_nonbouss(h, tv, pfu, pfv, g, gv, us, cs%PressureForce_Mont_CSp, &
77  p_atm, pbce, eta)
78  endif
79  endif
80 
81 end subroutine pressureforce
82 
83 !> Initialize the pressure force control structure
84 subroutine pressureforce_init(Time, G, GV, US, param_file, diag, CS, tides_CSp)
85  type(time_type), target, intent(in) :: Time !< Current model time
86  type(ocean_grid_type), intent(in) :: G !< Ocean grid structure
87  type(verticalgrid_type), intent(in) :: GV !< Vertical grid structure
88  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
89  type(param_file_type), intent(in) :: param_file !< Parameter file handles
90  type(diag_ctrl), target, intent(inout) :: diag !< Diagnostics control structure
91  type(pressureforce_cs), pointer :: CS !< Pressure force control structure
92  type(tidal_forcing_cs), optional, pointer :: tides_CSp !< Tide control structure
93 #include "version_variable.h"
94  character(len=40) :: mdl = "MOM_PressureForce" ! This module's name.
95 
96  if (associated(cs)) then
97  call mom_error(warning, "PressureForce_init called with an associated "// &
98  "control structure.")
99  return
100  else ; allocate(cs) ; endif
101 
102  ! Read all relevant parameters and write them to the model log.
103  call log_version(param_file, mdl, version, "")
104  call get_param(param_file, mdl, "ANALYTIC_FV_PGF", cs%Analytic_FV_PGF, &
105  "If true the pressure gradient forces are calculated "//&
106  "with a finite volume form that analytically integrates "//&
107  "the equations of state in pressure to avoid any "//&
108  "possibility of numerical thermobaric instability, as "//&
109  "described in Adcroft et al., O. Mod. (2008).", default=.true.)
110 
111  if (cs%Analytic_FV_PGF) then
112  call pressureforce_fv_init(time, g, gv, us, param_file, diag, &
113  cs%PressureForce_FV_CSp, tides_csp)
114  else
115  call pressureforce_mont_init(time, g, gv, us, param_file, diag, &
116  cs%PressureForce_Mont_CSp, tides_csp)
117  endif
118 
119 end subroutine pressureforce_init
120 
121 !> Deallocate the pressure force control structure
122 subroutine pressureforce_end(CS)
123  type(pressureforce_cs), pointer :: CS !< Pressure force control structure
124 
125  if (cs%Analytic_FV_PGF) then
126  call pressureforce_fv_end(cs%PressureForce_FV_CSp)
127  else
128  call pressureforce_mont_end(cs%PressureForce_Mont_CSp)
129  endif
130 
131  if (associated(cs)) deallocate(cs)
132 end subroutine pressureforce_end
133 
134 !> \namespace mom_pressureforce
135 !!
136 !! This thin module provides a branch to two forms of the horizontal accelerations
137 !! due to pressure gradients. The two options currently available are a
138 !! Montgomery potential form (used in traditional isopycnal layer models), and the
139 !! analytic finite volume form.
140 
141 end module mom_pressureforce
This module contains the main regridding routines.
Definition: MOM_ALE.F90:9
Finite volume pressure gradient (integrated by quadrature or analytically)
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:26
The following data type a list of diagnostic fields an their variants, as well as variables that cont...
A structure that can be parsed to read and document run-time parameters.
Provides the ocean grid type.
Definition: MOM_grid.F90:2
A thin wrapper for Boussinesq/non-Boussinesq forms of the pressure force calculation.
Finite volume pressure gradient control structure.
The MOM6 facility to parse input files for runtime parameters.
Describes various unit conversion factors.
Provides the Montgomery potential form of pressure gradient.
Tidal contributions to geopotential.
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
The subroutines here provide convenient wrappers to the fms diag_manager interfaces with additional d...
Routines for error handling and I/O management.
An overloaded interface to log version information about modules.
Describes the vertical ocean grid, including unit conversion factors.
Pointers to an assortment of thermodynamic fields that may be available, including potential temperat...
The control structure for the MOM_tidal_forcing module.
Provides a transparent vertical ocean grid type and supporting routines.
ALE control structure.
Definition: MOM_ALE.F90:63
Provides transparent structures with groups of MOM6 variables and supporting routines.
Pressure force control structure.
Control structure for the Montgomery potential form of pressure gradient.
An overloaded interface to read and log the values of various types of parameters.