MOM6
MOM_PressureForce.F90
Go to the documentation of this file.
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
17 use mom_ale, only: ale_cs
18 implicit none ; private
19 
20 #include <MOM_memory.h>
21 
23 
24 ! Pressure force control structure
25 type, public :: pressureforce_cs ; private
26  logical :: analytic_fv_pgf !< If true, use the analytic finite volume form
27  !! (Adcroft et al., Ocean Mod. 2008) of the PGF.
28  !> Control structure for the analytically integrated finite volume pressure force
29  type(pressureforce_afv_cs), pointer :: pressureforce_afv_csp => null()
30  !> Control structure for the Montgomery potential form of pressure force
31  type(pressureforce_mont_cs), pointer :: pressureforce_mont_csp => null()
32 end type pressureforce_cs
33 
34 contains
35 
36 !> A thin layer between the model and the Boussinesq and non-Boussinesq pressure force routines.
37 subroutine pressureforce(h, tv, PFu, PFv, G, GV, CS, ALE_CSp, p_atm, pbce, eta)
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 
69 end subroutine pressureforce
70 
71 !> Initialize the pressure force control structure
72 subroutine pressureforce_init(Time, G, GV, param_file, diag, CS, tides_CSp)
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 
106 end subroutine pressureforce_init
107 
108 !> Deallocate the pressure force control structure
109 subroutine pressureforce_end(CS)
110  type(pressureforce_cs), pointer :: CS !< Pressure force control structure
111  if (associated(cs)) deallocate(cs)
112 end subroutine pressureforce_end
113 
114 !> \namespace mom_pressureforce
115 !!
116 !! This thin module provides a branch to two forms of the horizontal accelerations
117 !! due to pressure gradients. The two options currently available are a
118 !! Montgomery potential form (used in traditional isopycnal layer models), and the
119 !! analytic finite volume form.
120 
121 end module mom_pressureforce
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 ...
Definition: MOM_ALE.F90:7
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.
Definition: MOM_grid.F90:19
The following data type a list of diagnostic fields an their variants, as well as variables that cont...
Provides the ocean grid type.
Definition: MOM_grid.F90:2
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 ...
ALE control structure.
Definition: MOM_ALE.F90:61
Control structure for the Montgomery potential form of pressure gradient.
subroutine, public mom_error(level, message, all_print)