MOM6
MOM_PressureForce_Montgomery.F90
1 !> Provides the Montgomery potential form of pressure gradient
3 
4 ! This file is part of MOM6. See LICENSE.md for the license.
5 
6 use mom_diag_mediator, only : post_data, register_diag_field
7 use mom_diag_mediator, only : safe_alloc_ptr, diag_ctrl, time_type
8 use mom_error_handler, only : mom_error, mom_mesg, fatal, warning, is_root_pe
10 use mom_grid, only : ocean_grid_type
11 use mom_tidal_forcing, only : calc_tidal_forcing, tidal_forcing_cs
16 use mom_eos, only : int_specific_vol_dp, query_compressible
17 
18 implicit none ; private
19 
20 #include <MOM_memory.h>
21 
22 public pressureforce_mont_bouss, pressureforce_mont_nonbouss, set_pbce_bouss
23 public set_pbce_nonbouss, pressureforce_mont_init, pressureforce_mont_end
24 
25 ! A note on unit descriptions in comments: MOM6 uses units that can be rescaled for dimensional
26 ! consistency testing. These are noted in comments with units like Z, H, L, and T, along with
27 ! their mks counterparts with notation like "a velocity [Z T-1 ~> m s-1]". If the units
28 ! vary with the Boussinesq approximation, the Boussinesq variant is given first.
29 
30 !> Control structure for the Montgomery potential form of pressure gradient
31 type, public :: pressureforce_mont_cs ; private
32  logical :: tides !< If true, apply tidal momentum forcing.
33  real :: rho0 !< The density used in the Boussinesq
34  !! approximation [kg m-3].
35  real :: rho_atm !< The assumed atmospheric density [kg m-3].
36  !! By default, Rho_atm is 0.
37  real :: gfs_scale !< Ratio between gravity applied to top interface and the
38  !! gravitational acceleration of the planet [nondim].
39  !! Usually this ratio is 1.
40  type(time_type), pointer :: time => null() !< A pointer to the ocean model's clock.
41  type(diag_ctrl), pointer :: diag => null() !< A structure that is used to regulate
42  !! the timing of diagnostic output.
43  real, pointer :: pfu_bc(:,:,:) => null() !< Accelerations due to pressure
44  real, pointer :: pfv_bc(:,:,:) => null() !< gradients deriving from density
45  !! gradients within layers [m s-2].
46  !>@{ Diagnostic IDs
47  integer :: id_pfu_bc = -1, id_pfv_bc = -1, id_e_tidal = -1
48  !!@}
49  type(tidal_forcing_cs), pointer :: tides_csp => null() !< The tidal forcing control structure
50 end type pressureforce_mont_cs
51 
52 contains
53 
54 !> \brief Non-Boussinesq Montgomery-potential form of pressure gradient
55 !!
56 !! Determines the acceleration due to pressure forces in a
57 !! non-Boussinesq fluid using the compressibility compensated (if appropriate)
58 !! Montgomery-potential form described in Hallberg (Ocean Mod., 2005).
59 !!
60 !! To work, the following fields must be set outside of the usual (is:ie,js:je)
61 !! range before this subroutine is called:
62 !! h(isB:ie+1,jsB:je+1), T(isB:ie+1,jsB:je+1), and S(isB:ie+1,jsB:je+1).
63 subroutine pressureforce_mont_nonbouss(h, tv, PFu, PFv, G, GV, US, CS, p_atm, pbce, eta)
64  type(ocean_grid_type), intent(in) :: g !< Ocean grid structure.
65  type(verticalgrid_type), intent(in) :: gv !< Vertical grid structure.
66  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
67  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thickness, [H ~> kg m-2].
68  type(thermo_var_ptrs), intent(in) :: tv !< Thermodynamic variables.
69  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), intent(out) :: pfu !< Zonal acceleration due to pressure gradients
70  !! (equal to -dM/dx) [m s-2].
71  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), intent(out) :: pfv !< Meridional acceleration due to pressure gradients
72  !! (equal to -dM/dy) [m s-2].
73  type(pressureforce_mont_cs), pointer :: cs !< Control structure for Montgomery potential PGF
74  real, dimension(:,:), optional, pointer :: p_atm !< The pressure at the ice-ocean or
75  !! atmosphere-ocean [Pa].
76  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
77  optional, intent(out) :: pbce !< The baroclinic pressure anomaly in
78  !! each layer due to free surface height anomalies,
79  !! [m2 s-2 H-1 ~> m s-2 or m4 kg-1 s-2].
80  real, dimension(SZI_(G),SZJ_(G)), optional, intent(out) :: eta !< Free surface height [H ~> kg m-1].
81 
82  ! Local variables
83  real, dimension(SZI_(G),SZJ_(G),SZK_(G)) :: &
84  m, & ! The Montgomery potential, M = (p/rho + gz) [m2 s-2].
85  alpha_star, & ! Compression adjusted specific volume [m3 kg-1].
86  dz_geo ! The change in geopotential across a layer [m2 s-2].
87  real, dimension(SZI_(G),SZJ_(G),SZK_(G)+1) :: p ! Interface pressure [Pa].
88  ! p may be adjusted (with a nonlinear equation of state) so that
89  ! its derivative compensates for the adiabatic compressibility
90  ! in seawater, but p will still be close to the pressure.
91  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), target :: &
92  t_tmp, & ! Temporary array of temperatures where layers that are lighter
93  ! than the mixed layer have the mixed layer's properties [degC].
94  s_tmp ! Temporary array of salinities where layers that are lighter
95  ! than the mixed layer have the mixed layer's properties [ppt].
96 
97  real, dimension(SZI_(G)) :: rho_cv_bl ! The coordinate potential density in the
98  ! deepest variable density near-surface layer [kg m-3].
99 
100  real, dimension(SZI_(G),SZJ_(G)) :: &
101  dm, & ! A barotropic correction to the Montgomery potentials to
102  ! enable the use of a reduced gravity form of the equations
103  ! [m2 s-2].
104  dp_star, & ! Layer thickness after compensation for compressibility [Pa].
105  ssh, & ! The sea surface height anomaly, in depth units [Z ~> m].
106  e_tidal, & ! Bottom geopotential anomaly due to tidal forces from
107  ! astronomical sources and self-attraction and loading [Z ~> m].
108  geopot_bot ! Bottom geopotential relative to time-mean sea level,
109  ! including any tidal contributions [m2 s-2].
110  real :: p_ref(szi_(g)) ! The pressure used to calculate the coordinate
111  ! density [Pa] (usually 2e7 Pa = 2000 dbar).
112  real :: rho_in_situ(szi_(g)) !In-situ density of a layer [kg m-3].
113  real :: pfu_bc, pfv_bc ! The pressure gradient force due to along-layer
114  ! compensated density gradients [m s-2]
115  real :: dp_neglect ! A thickness that is so small it is usually lost
116  ! in roundoff and can be neglected [Pa].
117  logical :: use_p_atm ! If true, use the atmospheric pressure.
118  logical :: use_eos ! If true, density is calculated from T & S using
119  ! an equation of state.
120  logical :: is_split ! A flag indicating whether the pressure
121  ! gradient terms are to be split into
122  ! barotropic and baroclinic pieces.
123  type(thermo_var_ptrs) :: tv_tmp! A structure of temporary T & S.
124 
125  real :: i_gearth ! The inverse of g_Earth [s2 Z m-2 ~> s2 m-1]
126 ! real :: dalpha
127  real :: pa_to_h ! A factor to convert from Pa to the thicknesss units (H).
128  real :: alpha_lay(szk_(g)) ! The specific volume of each layer [kg m-3].
129  real :: dalpha_int(szk_(g)+1) ! The change in specific volume across each
130  ! interface [kg m-3].
131  integer :: is, ie, js, je, isq, ieq, jsq, jeq, nz, nkmb
132  integer :: i, j, k
133  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
134  nkmb=gv%nk_rho_varies
135  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
136 
137  use_p_atm = .false.
138  if (present(p_atm)) then ; if (associated(p_atm)) use_p_atm = .true. ; endif
139  is_split = .false. ; if (present(pbce)) is_split = .true.
140  use_eos = associated(tv%eqn_of_state)
141 
142  if (.not.associated(cs)) call mom_error(fatal, &
143  "MOM_PressureForce_Mont: Module must be initialized before it is used.")
144  if (use_eos) then
145  if (query_compressible(tv%eqn_of_state)) call mom_error(fatal, &
146  "PressureForce_Mont_nonBouss: The Montgomery form of the pressure force "//&
147  "can no longer be used with a compressible EOS. Use #define ANALYTIC_FV_PGF.")
148  endif
149 
150  i_gearth = 1.0 / (us%L_T_to_m_s**2 * gv%g_Earth)
151  dp_neglect = gv%H_to_Pa * gv%H_subroundoff
152  do k=1,nz ; alpha_lay(k) = 1.0 / gv%Rlay(k) ; enddo
153  do k=2,nz ; dalpha_int(k) = alpha_lay(k-1) - alpha_lay(k) ; enddo
154 
155  if (use_p_atm) then
156  !$OMP parallel do default(shared)
157  do j=jsq,jeq+1 ; do i=isq,ieq+1 ; p(i,j,1) = p_atm(i,j) ; enddo ; enddo
158  else
159  !$OMP parallel do default(shared)
160  do j=jsq,jeq+1 ; do i=isq,ieq+1 ; p(i,j,1) = 0.0 ; enddo ; enddo
161  endif
162  !$OMP parallel do default(shared)
163  do j=jsq,jeq+1 ; do k=1,nz ; do i=isq,ieq+1
164  p(i,j,k+1) = p(i,j,k) + gv%H_to_Pa * h(i,j,k)
165  enddo ; enddo ; enddo
166 
167  if (present(eta)) then
168  pa_to_h = 1.0 / gv%H_to_Pa
169  if (use_p_atm) then
170  !$OMP parallel do default(shared)
171  do j=jsq,jeq+1 ; do i=isq,ieq+1
172  eta(i,j) = (p(i,j,nz+1) - p_atm(i,j))*pa_to_h ! eta has the same units as h.
173  enddo ; enddo
174  else
175  !$OMP parallel do default(shared)
176  do j=jsq,jeq+1 ; do i=isq,ieq+1
177  eta(i,j) = p(i,j,nz+1)*pa_to_h ! eta has the same units as h.
178  enddo ; enddo
179  endif
180  endif
181 
182  if (cs%tides) then
183  ! Determine the sea surface height anomalies, to enable the calculation
184  ! of self-attraction and loading.
185  !$OMP parallel do default(shared)
186  do j=jsq,jeq+1 ; do i=isq,ieq+1
187  ssh(i,j) = -g%bathyT(i,j)
188  enddo ; enddo
189  if (use_eos) then
190  !$OMP parallel do default(shared)
191  do k=1,nz
192  call int_specific_vol_dp(tv%T(:,:,k), tv%S(:,:,k), p(:,:,k), p(:,:,k+1), &
193  0.0, g%HI, tv%eqn_of_state, dz_geo(:,:,k), halo_size=1)
194  enddo
195  !$OMP parallel do default(shared)
196  do j=jsq,jeq+1 ; do k=1,nz; do i=isq,ieq+1
197  ssh(i,j) = ssh(i,j) + i_gearth * dz_geo(i,j,k)
198  enddo ; enddo ; enddo
199  else
200  !$OMP parallel do default(shared)
201  do j=jsq,jeq+1 ; do k=1,nz ; do i=isq,ieq+1
202  ssh(i,j) = ssh(i,j) + (us%m_to_Z*gv%H_to_kg_m2)*h(i,j,k)*alpha_lay(k)
203  enddo ; enddo ; enddo
204  endif
205 
206  call calc_tidal_forcing(cs%Time, ssh, e_tidal, g, cs%tides_CSp, m_to_z=us%m_to_Z)
207  !$OMP parallel do default(shared)
208  do j=jsq,jeq+1 ; do i=isq,ieq+1
209  geopot_bot(i,j) = -us%L_T_to_m_s**2 * gv%g_Earth*(e_tidal(i,j) + g%bathyT(i,j))
210  enddo ; enddo
211  else
212  !$OMP parallel do default(shared)
213  do j=jsq,jeq+1 ; do i=isq,ieq+1
214  geopot_bot(i,j) = -us%L_T_to_m_s**2 * gv%g_Earth*g%bathyT(i,j)
215  enddo ; enddo
216  endif
217 
218  if (use_eos) then
219  ! Calculate in-situ specific volumes (alpha_star).
220 
221  ! With a bulk mixed layer, replace the T & S of any layers that are
222  ! lighter than the the buffer layer with the properties of the buffer
223  ! layer. These layers will be massless anyway, and it avoids any
224  ! formal calculations with hydrostatically unstable profiles.
225  if (nkmb>0) then
226  tv_tmp%T => t_tmp ; tv_tmp%S => s_tmp
227  tv_tmp%eqn_of_state => tv%eqn_of_state
228  do i=isq,ieq+1 ; p_ref(i) = tv%P_Ref ; enddo
229  !$OMP parallel do default(shared) private(Rho_cv_BL)
230  do j=jsq,jeq+1
231  do k=1,nkmb ; do i=isq,ieq+1
232  tv_tmp%T(i,j,k) = tv%T(i,j,k) ; tv_tmp%S(i,j,k) = tv%S(i,j,k)
233  enddo ; enddo
234  call calculate_density(tv%T(:,j,nkmb), tv%S(:,j,nkmb), p_ref, &
235  rho_cv_bl(:), isq, ieq-isq+2, tv%eqn_of_state)
236  do k=nkmb+1,nz ; do i=isq,ieq+1
237  if (gv%Rlay(k) < rho_cv_bl(i)) then
238  tv_tmp%T(i,j,k) = tv%T(i,j,nkmb) ; tv_tmp%S(i,j,k) = tv%S(i,j,nkmb)
239  else
240  tv_tmp%T(i,j,k) = tv%T(i,j,k) ; tv_tmp%S(i,j,k) = tv%S(i,j,k)
241  endif
242  enddo ; enddo
243  enddo
244  else
245  tv_tmp%T => tv%T ; tv_tmp%S => tv%S
246  tv_tmp%eqn_of_state => tv%eqn_of_state
247  do i=isq,ieq+1 ; p_ref(i) = 0 ; enddo
248  endif
249  !$OMP parallel do default(shared) private(rho_in_situ)
250  do k=1,nz ; do j=jsq,jeq+1
251  call calculate_density(tv_tmp%T(:,j,k),tv_tmp%S(:,j,k),p_ref, &
252  rho_in_situ,isq,ieq-isq+2,tv%eqn_of_state)
253  do i=isq,ieq+1 ; alpha_star(i,j,k) = 1.0 / rho_in_situ(i) ; enddo
254  enddo ; enddo
255  endif ! use_EOS
256 
257  if (use_eos) then
258  !$OMP parallel do default(shared)
259  do j=jsq,jeq+1
260  do i=isq,ieq+1
261  m(i,j,nz) = geopot_bot(i,j) + p(i,j,nz+1) * alpha_star(i,j,nz)
262  enddo
263  do k=nz-1,1,-1 ; do i=isq,ieq+1
264  m(i,j,k) = m(i,j,k+1) + p(i,j,k+1) * (alpha_star(i,j,k) - alpha_star(i,j,k+1))
265  enddo ; enddo
266  enddo
267  else ! not use_EOS
268  !$OMP parallel do default(shared)
269  do j=jsq,jeq+1
270  do i=isq,ieq+1
271  m(i,j,nz) = geopot_bot(i,j) + p(i,j,nz+1) * alpha_lay(nz)
272  enddo
273  do k=nz-1,1,-1 ; do i=isq,ieq+1
274  m(i,j,k) = m(i,j,k+1) + p(i,j,k+1) * dalpha_int(k+1)
275  enddo ; enddo
276  enddo
277  endif ! use_EOS
278 
279  if (cs%GFS_scale < 1.0) then
280  ! Adjust the Montgomery potential to make this a reduced gravity model.
281  !$OMP parallel do default(shared)
282  do j=jsq,jeq+1 ; do i=isq,ieq+1
283  dm(i,j) = (cs%GFS_scale - 1.0) * m(i,j,1)
284  enddo ; enddo
285  !$OMP parallel do default(shared)
286  do k=1,nz ; do j=jsq,jeq+1 ; do i=isq,ieq+1
287  m(i,j,k) = m(i,j,k) + dm(i,j)
288  enddo ; enddo ; enddo
289 
290  ! Could instead do the following, to avoid taking small differences
291  ! of large numbers...
292 ! do j=Jsq,Jeq+1 ; do i=Isq,Ieq+1
293 ! M(i,j,1) = CS%GFS_scale * M(i,j,1)
294 ! enddo ; enddo
295 ! if (use_EOS) then
296 ! do k=2,nz ; do j=Jsq,Jeq+1 ; do i=Isq,Ieq+1
297 ! M(i,j,k) = M(i,j,k-1) - p(i,j,K) * (alpha_star(i,j,k-1) - alpha_star(i,j,k))
298 ! enddo ; enddo ; enddo
299 ! else ! not use_EOS
300 ! do k=2,nz ; do j=Jsq,Jeq+1 ; do i=Isq,Ieq+1
301 ! M(i,j,k) = M(i,j,k-1) - p(i,j,K) * dalpha_int(K)
302 ! enddo ; enddo ; enddo
303 ! endif ! use_EOS
304 
305  endif
306 
307  ! Note that ddM/dPb = alpha_star(i,j,1)
308  if (present(pbce)) then
309  call set_pbce_nonbouss(p, tv_tmp, g, gv, us, cs%GFS_scale, pbce, alpha_star)
310  endif
311 
312 ! Calculate the pressure force. On a Cartesian grid,
313 ! PFu = - dM/dx and PFv = - dM/dy.
314  if (use_eos) then
315  !$OMP parallel do default(shared) private(dp_star,PFu_bc,PFv_bc)
316  do k=1,nz
317  do j=jsq,jeq+1 ; do i=isq,ieq+1
318  dp_star(i,j) = (p(i,j,k+1) - p(i,j,k)) + dp_neglect
319  enddo ; enddo
320  do j=js,je ; do i=isq,ieq
321  ! PFu_bc = p* grad alpha*
322  pfu_bc = (alpha_star(i+1,j,k) - alpha_star(i,j,k)) * (g%IdxCu(i,j) * &
323  ((dp_star(i,j) * dp_star(i+1,j) + (p(i,j,k) * dp_star(i+1,j) + &
324  p(i+1,j,k) * dp_star(i,j))) / (dp_star(i,j) + dp_star(i+1,j))))
325  pfu(i,j,k) = -(m(i+1,j,k) - m(i,j,k)) * g%IdxCu(i,j) + pfu_bc
326  if (associated(cs%PFu_bc)) cs%PFu_bc(i,j,k) = pfu_bc
327  enddo ; enddo
328  do j=jsq,jeq ; do i=is,ie
329  pfv_bc = (alpha_star(i,j+1,k) - alpha_star(i,j,k)) * (g%IdyCv(i,j) * &
330  ((dp_star(i,j) * dp_star(i,j+1) + (p(i,j,k) * dp_star(i,j+1) + &
331  p(i,j+1,k) * dp_star(i,j))) / (dp_star(i,j) + dp_star(i,j+1))))
332  pfv(i,j,k) = -(m(i,j+1,k) - m(i,j,k)) * g%IdyCv(i,j) + pfv_bc
333  if (associated(cs%PFv_bc)) cs%PFv_bc(i,j,k) = pfv_bc
334  enddo ; enddo
335  enddo ! k-loop
336  else ! .not. use_EOS
337  !$OMP parallel do default(shared)
338  do k=1,nz
339  do j=js,je ; do i=isq,ieq
340  pfu(i,j,k) = -(m(i+1,j,k) - m(i,j,k)) * g%IdxCu(i,j)
341  enddo ; enddo
342  do j=jsq,jeq ; do i=is,ie
343  pfv(i,j,k) = -(m(i,j+1,k) - m(i,j,k)) * g%IdyCv(i,j)
344  enddo ; enddo
345  enddo
346  endif ! use_EOS
347 
348  if (cs%id_PFu_bc>0) call post_data(cs%id_PFu_bc, cs%PFu_bc, cs%diag)
349  if (cs%id_PFv_bc>0) call post_data(cs%id_PFv_bc, cs%PFv_bc, cs%diag)
350  if (cs%id_e_tidal>0) call post_data(cs%id_e_tidal, e_tidal, cs%diag)
351 
352 end subroutine pressureforce_mont_nonbouss
353 
354 !> \brief Boussinesq Montgomery-potential form of pressure gradient
355 !!
356 !! Determines the acceleration due to pressure forces.
357 !!
358 !! To work, the following fields must be set outside of the usual (is:ie,js:je)
359 !! range before this subroutine is called:
360 !! h(isB:ie+1,jsB:je+1), T(isB:ie+1,jsB:je+1), and S(isB:ie+1,jsB:je+1).
361 subroutine pressureforce_mont_bouss(h, tv, PFu, PFv, G, GV, US, CS, p_atm, pbce, eta)
362  type(ocean_grid_type), intent(in) :: g !< Ocean grid structure.
363  type(verticalgrid_type), intent(in) :: gv !< Vertical grid structure.
364  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
365  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thickness [H ~> m].
366  type(thermo_var_ptrs), intent(in) :: tv !< Thermodynamic variables.
367  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), intent(out) :: pfu !< Zonal acceleration due to pressure gradients
368  !! (equal to -dM/dx) [m s-2].
369  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), intent(out) :: pfv !< Meridional acceleration due to pressure gradients
370  !! (equal to -dM/dy) [m s2].
371  type(pressureforce_mont_cs), pointer :: cs !< Control structure for Montgomery potential PGF
372  real, dimension(:,:), optional, pointer :: p_atm !< The pressure at the ice-ocean or
373  !! atmosphere-ocean [Pa].
374  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), optional, intent(out) :: pbce !< The baroclinic pressure anomaly in
375  !! each layer due to free surface height anomalies
376  !! [m2 s-2 H-1 ~> m s-2].
377  real, dimension(SZI_(G),SZJ_(G)), optional, intent(out) :: eta !< Free surface height [H ~> m].
378  ! Local variables
379  real, dimension(SZI_(G),SZJ_(G),SZK_(G)) :: &
380  m, & ! The Montgomery potential, M = (p/rho + gz) [m2 s-2].
381  rho_star ! In-situ density divided by the derivative with depth of the
382  ! corrected e times (G_Earth/Rho0) [m2 Z-1 s-2 ~> m s-2].
383  real, dimension(SZI_(G),SZJ_(G),SZK_(G)+1) :: e ! Interface height in m.
384  ! e may be adjusted (with a nonlinear equation of state) so that
385  ! its derivative compensates for the adiabatic compressibility
386  ! in seawater, but e will still be close to the interface depth.
387  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), target :: &
388  t_tmp, & ! Temporary array of temperatures where layers that are lighter
389  ! than the mixed layer have the mixed layer's properties [degC].
390  s_tmp ! Temporary array of salinities where layers that are lighter
391  ! than the mixed layer have the mixed layer's properties [ppt].
392 
393  real :: rho_cv_bl(szi_(g)) ! The coordinate potential density in
394  ! the deepest variable density near-surface layer [kg m-3].
395  real :: h_star(szi_(g),szj_(g)) ! Layer thickness after compensation
396  ! for compressibility [Z ~> m].
397  real :: e_tidal(szi_(g),szj_(g)) ! Bottom geopotential anomaly due to tidal
398  ! forces from astronomical sources and self-
399  ! attraction and loading, in depth units [Z ~> m].
400  real :: p_ref(szi_(g)) ! The pressure used to calculate the coordinate
401  ! density [Pa] (usually 2e7 Pa = 2000 dbar).
402  real :: i_rho0 ! 1/Rho0 [m3 kg-1].
403  real :: g_rho0 ! G_Earth / Rho0 [m5 Z-1 s-2 kg-1 ~> m4 s-2 kg-1].
404  real :: pfu_bc, pfv_bc ! The pressure gradient force due to along-layer
405  ! compensated density gradients [m s-2]
406 ! real :: dr ! Temporary variables.
407  real :: h_neglect ! A thickness that is so small it is usually lost
408  ! in roundoff and can be neglected [Z ~> m].
409  logical :: use_p_atm ! If true, use the atmospheric pressure.
410  logical :: use_eos ! If true, density is calculated from T & S using
411  ! an equation of state.
412  logical :: is_split ! A flag indicating whether the pressure
413  ! gradient terms are to be split into
414  ! barotropic and baroclinic pieces.
415  type(thermo_var_ptrs) :: tv_tmp! A structure of temporary T & S.
416  integer :: is, ie, js, je, isq, ieq, jsq, jeq, nz, nkmb
417  integer :: i, j, k
418 
419  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
420  nkmb=gv%nk_rho_varies
421  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
422 
423  use_p_atm = .false.
424  if (present(p_atm)) then ; if (associated(p_atm)) use_p_atm = .true. ; endif
425  is_split = .false. ; if (present(pbce)) is_split = .true.
426  use_eos = associated(tv%eqn_of_state)
427 
428  if (.not.associated(cs)) call mom_error(fatal, &
429  "MOM_PressureForce_Mont: Module must be initialized before it is used.")
430  if (use_eos) then
431  if (query_compressible(tv%eqn_of_state)) call mom_error(fatal, &
432  "PressureForce_Mont_Bouss: The Montgomery form of the pressure force "//&
433  "can no longer be used with a compressible EOS. Use #define ANALYTIC_FV_PGF.")
434  endif
435 
436  h_neglect = gv%H_subroundoff * gv%H_to_Z
437  i_rho0 = 1.0/cs%Rho0
438  g_rho0 = us%L_T_to_m_s**2 * gv%g_Earth/gv%Rho0
439 
440  if (cs%tides) then
441  ! Determine the surface height anomaly for calculating self attraction
442  ! and loading. This should really be based on bottom pressure anomalies,
443  ! but that is not yet implemented, and the current form is correct for
444  ! barotropic tides.
445  !$OMP parallel do default(shared)
446  do j=jsq,jeq+1
447  do i=isq,ieq+1 ; e(i,j,1) = -g%bathyT(i,j) ; enddo
448  do k=1,nz ; do i=isq,ieq+1
449  e(i,j,1) = e(i,j,1) + h(i,j,k)*gv%H_to_Z
450  enddo ; enddo
451  enddo
452  call calc_tidal_forcing(cs%Time, e(:,:,1), e_tidal, g, cs%tides_CSp, m_to_z=us%m_to_Z)
453  endif
454 
455 ! Here layer interface heights, e, are calculated.
456  if (cs%tides) then
457  !$OMP parallel do default(shared)
458  do j=jsq,jeq+1 ; do i=isq,ieq+1
459  e(i,j,nz+1) = -(g%bathyT(i,j) + e_tidal(i,j))
460  enddo ; enddo
461  else
462  !$OMP parallel do default(shared)
463  do j=jsq,jeq+1 ; do i=isq,ieq+1
464  e(i,j,nz+1) = -g%bathyT(i,j)
465  enddo ; enddo
466  endif
467  !$OMP parallel do default(shared)
468  do j=jsq,jeq+1 ; do k=nz,1,-1 ; do i=isq,ieq+1
469  e(i,j,k) = e(i,j,k+1) + h(i,j,k)*gv%H_to_Z
470  enddo ; enddo ; enddo
471 
472  if (use_eos) then
473 ! Calculate in-situ densities (rho_star).
474 
475 ! With a bulk mixed layer, replace the T & S of any layers that are
476 ! lighter than the the buffer layer with the properties of the buffer
477 ! layer. These layers will be massless anyway, and it avoids any
478 ! formal calculations with hydrostatically unstable profiles.
479 
480  if (nkmb>0) then
481  tv_tmp%T => t_tmp ; tv_tmp%S => s_tmp
482  tv_tmp%eqn_of_state => tv%eqn_of_state
483 
484  do i=isq,ieq+1 ; p_ref(i) = tv%P_Ref ; enddo
485  !$OMP parallel do default(shared) private(Rho_cv_BL)
486  do j=jsq,jeq+1
487  do k=1,nkmb ; do i=isq,ieq+1
488  tv_tmp%T(i,j,k) = tv%T(i,j,k) ; tv_tmp%S(i,j,k) = tv%S(i,j,k)
489  enddo ; enddo
490  call calculate_density(tv%T(:,j,nkmb), tv%S(:,j,nkmb), p_ref, &
491  rho_cv_bl(:), isq, ieq-isq+2, tv%eqn_of_state)
492 
493  do k=nkmb+1,nz ; do i=isq,ieq+1
494  if (gv%Rlay(k) < rho_cv_bl(i)) then
495  tv_tmp%T(i,j,k) = tv%T(i,j,nkmb) ; tv_tmp%S(i,j,k) = tv%S(i,j,nkmb)
496  else
497  tv_tmp%T(i,j,k) = tv%T(i,j,k) ; tv_tmp%S(i,j,k) = tv%S(i,j,k)
498  endif
499  enddo ; enddo
500  enddo
501  else
502  tv_tmp%T => tv%T ; tv_tmp%S => tv%S
503  tv_tmp%eqn_of_state => tv%eqn_of_state
504  do i=isq,ieq+1 ; p_ref(i) = 0.0 ; enddo
505  endif
506 
507  ! This no longer includes any pressure dependency, since this routine
508  ! will come down with a fatal error if there is any compressibility.
509  !$OMP parallel do default(shared)
510  do k=1,nz+1 ; do j=jsq,jeq+1
511  call calculate_density(tv_tmp%T(:,j,k), tv_tmp%S(:,j,k), p_ref, rho_star(:,j,k), &
512  isq,ieq-isq+2,tv%eqn_of_state)
513  do i=isq,ieq+1 ; rho_star(i,j,k) = g_rho0*rho_star(i,j,k) ; enddo
514  enddo ; enddo
515  endif ! use_EOS
516 
517 ! Here the layer Montgomery potentials, M, are calculated.
518  if (use_eos) then
519  !$OMP parallel do default(shared)
520  do j=jsq,jeq+1
521  do i=isq,ieq+1
522  m(i,j,1) = cs%GFS_scale * (rho_star(i,j,1) * e(i,j,1))
523  if (use_p_atm) m(i,j,1) = m(i,j,1) + p_atm(i,j) * i_rho0
524  enddo
525  do k=2,nz ; do i=isq,ieq+1
526  m(i,j,k) = m(i,j,k-1) + (rho_star(i,j,k) - rho_star(i,j,k-1)) * e(i,j,k)
527  enddo ; enddo
528  enddo
529  else ! not use_EOS
530  !$OMP parallel do default(shared)
531  do j=jsq,jeq+1
532  do i=isq,ieq+1
533  m(i,j,1) = us%L_to_m**2*us%s_to_T**2*gv%g_prime(1) * e(i,j,1)
534  if (use_p_atm) m(i,j,1) = m(i,j,1) + p_atm(i,j) * i_rho0
535  enddo
536  do k=2,nz ; do i=isq,ieq+1
537  m(i,j,k) = m(i,j,k-1) + us%L_to_m**2*us%s_to_T**2*gv%g_prime(k) * e(i,j,k)
538  enddo ; enddo
539  enddo
540  endif ! use_EOS
541 
542  if (present(pbce)) then
543  call set_pbce_bouss(e, tv_tmp, g, gv, us, cs%Rho0, cs%GFS_scale, pbce, rho_star)
544  endif
545 
546 ! Calculate the pressure force. On a Cartesian grid,
547 ! PFu = - dM/dx and PFv = - dM/dy.
548  if (use_eos) then
549  !$OMP parallel do default(shared) private(h_star,PFu_bc,PFv_bc)
550  do k=1,nz
551  do j=jsq,jeq+1 ; do i=isq,ieq+1
552  h_star(i,j) = (e(i,j,k) - e(i,j,k+1)) + h_neglect
553  enddo ; enddo
554  do j=js,je ; do i=isq,ieq
555  pfu_bc = -1.0*(rho_star(i+1,j,k) - rho_star(i,j,k)) * (g%IdxCu(i,j) * &
556  ((h_star(i,j) * h_star(i+1,j) - (e(i,j,k) * h_star(i+1,j) + &
557  e(i+1,j,k) * h_star(i,j))) / (h_star(i,j) + h_star(i+1,j))))
558  pfu(i,j,k) = -(m(i+1,j,k) - m(i,j,k)) * g%IdxCu(i,j) + pfu_bc
559  if (associated(cs%PFu_bc)) cs%PFu_bc(i,j,k) = pfu_bc
560  enddo ; enddo
561  do j=jsq,jeq ; do i=is,ie
562  pfv_bc = -1.0*(rho_star(i,j+1,k) - rho_star(i,j,k)) * (g%IdyCv(i,j) * &
563  ((h_star(i,j) * h_star(i,j+1) - (e(i,j,k) * h_star(i,j+1) + &
564  e(i,j+1,k) * h_star(i,j))) / (h_star(i,j) + h_star(i,j+1))))
565  pfv(i,j,k) = -(m(i,j+1,k) - m(i,j,k)) * g%IdyCv(i,j) + pfv_bc
566  if (associated(cs%PFv_bc)) cs%PFv_bc(i,j,k) = pfv_bc
567  enddo ; enddo
568  enddo ! k-loop
569  else ! .not. use_EOS
570  !$OMP parallel do default(shared)
571  do k=1,nz
572  do j=js,je ; do i=isq,ieq
573  pfu(i,j,k) = -(m(i+1,j,k) - m(i,j,k)) * g%IdxCu(i,j)
574  enddo ; enddo
575  do j=jsq,jeq ; do i=is,ie
576  pfv(i,j,k) = -(m(i,j+1,k) - m(i,j,k)) * g%IdyCv(i,j)
577  enddo ; enddo
578  enddo
579  endif ! use_EOS
580 
581  if (present(eta)) then
582  if (cs%tides) then
583  ! eta is the sea surface height relative to a time-invariant geoid, for
584  ! comparison with what is used for eta in btstep. See how e was calculated
585  ! about 200 lines above.
586  !$OMP parallel do default(shared)
587  do j=jsq,jeq+1 ; do i=isq,ieq+1
588  eta(i,j) = e(i,j,1)*gv%Z_to_H + e_tidal(i,j)*gv%Z_to_H
589  enddo ; enddo
590  else
591  !$OMP parallel do default(shared)
592  do j=jsq,jeq+1 ; do i=isq,ieq+1
593  eta(i,j) = e(i,j,1)*gv%Z_to_H
594  enddo ; enddo
595  endif
596  endif
597 
598  if (cs%id_PFu_bc>0) call post_data(cs%id_PFu_bc, cs%PFu_bc, cs%diag)
599  if (cs%id_PFv_bc>0) call post_data(cs%id_PFv_bc, cs%PFv_bc, cs%diag)
600  if (cs%id_e_tidal>0) call post_data(cs%id_e_tidal, e_tidal, cs%diag)
601 
602 end subroutine pressureforce_mont_bouss
603 
604 !> Determines the partial derivative of the acceleration due
605 !! to pressure forces with the free surface height.
606 subroutine set_pbce_bouss(e, tv, G, GV, US, Rho0, GFS_scale, pbce, rho_star)
607  type(ocean_grid_type), intent(in) :: g !< Ocean grid structure
608  type(verticalgrid_type), intent(in) :: gv !< Vertical grid structure
609  real, dimension(SZI_(G),SZJ_(G),SZK_(G)+1), intent(in) :: e !< Interface height [Z ~> m].
610  type(thermo_var_ptrs), intent(in) :: tv !< Thermodynamic variables
611  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
612  real, intent(in) :: rho0 !< The "Boussinesq" ocean density [kg m-3].
613  real, intent(in) :: gfs_scale !< Ratio between gravity applied to top
614  !! interface and the gravitational acceleration of
615  !! the planet [nondim]. Usually this ratio is 1.
616  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
617  intent(out) :: pbce !< The baroclinic pressure anomaly in each layer due
618  !! to free surface height anomalies
619  !! [m2 H-1 s-2 ~> m4 kg-2 s-2].
620  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
621  optional, intent(in) :: rho_star !< The layer densities (maybe compressibility
622  !! compensated), times g/rho_0 [m2 Z-1 s-2 ~> m s-2].
623 
624  ! Local variables
625  real :: ihtot(szi_(g)) ! The inverse of the sum of the layer thicknesses [H-1 ~> m-1 or m2 kg-1].
626  real :: press(szi_(g)) ! Interface pressure [Pa].
627  real :: t_int(szi_(g)) ! Interface temperature [degC].
628  real :: s_int(szi_(g)) ! Interface salinity [ppt].
629  real :: dr_dt(szi_(g)) ! Partial derivative of density with temperature [kg m-3 degC-1].
630  real :: dr_ds(szi_(g)) ! Partial derivative of density with salinity [kg m-3 ppt-1].
631  real :: rho_in_situ(szi_(g)) !In-situ density at the top of a layer [kg m-3].
632  real :: g_rho0 ! A scaled version of g_Earth / Rho0 [L2 m3 Z-1 T-2 kg-1 ~> m4 s-2 kg-1]
633  real :: rho0xg ! g_Earth * Rho0 [kg s-2 m-1 Z-1 ~> kg s-2 m-2]
634  logical :: use_eos ! If true, density is calculated from T & S using
635  ! an equation of state.
636  real :: z_neglect ! A thickness that is so small it is usually lost
637  ! in roundoff and can be neglected [Z ~> m].
638  integer :: isq, ieq, jsq, jeq, nz, i, j, k
639 
640  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB ; nz = g%ke
641 
642  rho0xg = rho0*us%L_T_to_m_s**2 * gv%g_Earth
643  g_rho0 = gv%g_Earth / gv%Rho0
644  use_eos = associated(tv%eqn_of_state)
645  z_neglect = gv%H_subroundoff*gv%H_to_Z
646 
647  if (use_eos) then
648  if (present(rho_star)) then
649  !$OMP parallel do default(shared) private(Ihtot)
650  do j=jsq,jeq+1
651  do i=isq,ieq+1
652  ihtot(i) = gv%H_to_Z / ((e(i,j,1)-e(i,j,nz+1)) + z_neglect)
653  pbce(i,j,1) = gfs_scale * us%m_s_to_L_T**2*rho_star(i,j,1) * gv%H_to_Z
654  enddo
655  do k=2,nz ; do i=isq,ieq+1
656  pbce(i,j,k) = pbce(i,j,k-1) + us%m_s_to_L_T**2*(rho_star(i,j,k)-rho_star(i,j,k-1)) * &
657  ((e(i,j,k) - e(i,j,nz+1)) * ihtot(i))
658  enddo ; enddo
659  enddo ! end of j loop
660  else
661  !$OMP parallel do default(shared) private(Ihtot,press,rho_in_situ,T_int,S_int,dR_dT,dR_dS)
662  do j=jsq,jeq+1
663  do i=isq,ieq+1
664  ihtot(i) = gv%H_to_Z / ((e(i,j,1)-e(i,j,nz+1)) + z_neglect)
665  press(i) = -rho0xg*e(i,j,1)
666  enddo
667  call calculate_density(tv%T(:,j,1), tv%S(:,j,1), press, rho_in_situ, &
668  isq, ieq-isq+2, tv%eqn_of_state)
669  do i=isq,ieq+1
670  pbce(i,j,1) = g_rho0*(gfs_scale * rho_in_situ(i)) * gv%H_to_Z
671  enddo
672  do k=2,nz
673  do i=isq,ieq+1
674  press(i) = -rho0xg*e(i,j,k)
675  t_int(i) = 0.5*(tv%T(i,j,k-1)+tv%T(i,j,k))
676  s_int(i) = 0.5*(tv%S(i,j,k-1)+tv%S(i,j,k))
677  enddo
678  call calculate_density_derivs(t_int, s_int, press, dr_dt, dr_ds, &
679  isq, ieq-isq+2, tv%eqn_of_state)
680  do i=isq,ieq+1
681  pbce(i,j,k) = pbce(i,j,k-1) + g_rho0 * &
682  ((e(i,j,k) - e(i,j,nz+1)) * ihtot(i)) * &
683  (dr_dt(i)*(tv%T(i,j,k)-tv%T(i,j,k-1)) + &
684  dr_ds(i)*(tv%S(i,j,k)-tv%S(i,j,k-1)))
685  enddo
686  enddo
687  enddo ! end of j loop
688  endif
689  else ! not use_EOS
690  !$OMP parallel do default(shared) private(Ihtot)
691  do j=jsq,jeq+1
692  do i=isq,ieq+1
693  ihtot(i) = 1.0 / ((e(i,j,1)-e(i,j,nz+1)) + z_neglect)
694  pbce(i,j,1) = gv%g_prime(1) * gv%H_to_Z
695  enddo
696  do k=2,nz ; do i=isq,ieq+1
697  pbce(i,j,k) = pbce(i,j,k-1) + &
698  (gv%g_prime(k)*gv%H_to_Z) * ((e(i,j,k) - e(i,j,nz+1)) * ihtot(i))
699  enddo ; enddo
700  enddo ! end of j loop
701  endif ! use_EOS
702 
703 end subroutine set_pbce_bouss
704 
705 !> Determines the partial derivative of the acceleration due
706 !! to pressure forces with the column mass.
707 subroutine set_pbce_nonbouss(p, tv, G, GV, US, GFS_scale, pbce, alpha_star)
708  type(ocean_grid_type), intent(in) :: g !< Ocean grid structure
709  type(verticalgrid_type), intent(in) :: gv !< Vertical grid structure
710  real, dimension(SZI_(G),SZJ_(G),SZK_(G)+1), intent(in) :: p !< Interface pressures [Pa].
711  type(thermo_var_ptrs), intent(in) :: tv !< Thermodynamic variables
712  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
713  real, intent(in) :: gfs_scale !< Ratio between gravity applied to top
714  !! interface and the gravitational acceleration of
715  !! the planet [nondim]. Usually this ratio is 1.
716  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(out) :: pbce !< The baroclinic pressure anomaly in each layer due
717  !! to free surface height anomalies
718  !! [L2 H-1 T-2 ~> m4 kg-1 s-2].
719  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), optional, intent(in) :: alpha_star !< The layer specific volumes
720  !! (maybe compressibility compensated) [m3 kg-1].
721  ! Local variables
722  real, dimension(SZI_(G),SZJ_(G)) :: &
723  dpbce, & ! A barotropic correction to the pbce to enable the use of
724  ! a reduced gravity form of the equations [L2 H-1 T-2 ~> m4 kg-1 s-2].
725  c_htot ! dP_dH divided by the total ocean pressure [Z2 s2 m-2 T-2 H-1 ~> m2 kg-1].
726  real :: t_int(szi_(g)) ! Interface temperature [degC].
727  real :: s_int(szi_(g)) ! Interface salinity [ppt].
728  real :: dr_dt(szi_(g)) ! Partial derivative of density with temperature [kg m-3 degC-1].
729  real :: dr_ds(szi_(g)) ! Partial derivative of density with salinity [kg m-3 ppt-1].
730  real :: rho_in_situ(szi_(g)) ! In-situ density at an interface [kg m-3].
731  real :: alpha_lay(szk_(g)) ! The specific volume of each layer [kg m-3].
732  real :: dalpha_int(szk_(g)+1) ! The change in specific volume across each interface [kg m-3].
733  real :: dp_dh ! A factor that converts from thickness to pressure times other dimensional
734  ! conversion factors [Z2 s2 Pa m-2 T-2 H-1 ~> Pa m2 kg-1].
735  real :: dp_neglect ! A thickness that is so small it is usually lost
736  ! in roundoff and can be neglected [Pa].
737  logical :: use_eos ! If true, density is calculated from T & S using
738  ! an equation of state.
739  integer :: isq, ieq, jsq, jeq, nz, i, j, k
740 
741  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB ; nz = g%ke
742 
743  use_eos = associated(tv%eqn_of_state)
744 
745  dp_dh = us%m_s_to_L_T**2*gv%H_to_Pa
746  dp_neglect = gv%H_to_Pa * gv%H_subroundoff
747 
748  do k=1,nz ; alpha_lay(k) = 1.0 / gv%Rlay(k) ; enddo
749  do k=2,nz ; dalpha_int(k) = alpha_lay(k-1) - alpha_lay(k) ; enddo
750 
751  if (use_eos) then
752  if (present(alpha_star)) then
753  !$OMP parallel do default(shared)
754  do j=jsq,jeq+1
755  do i=isq,ieq+1
756  c_htot(i,j) = dp_dh / ((p(i,j,nz+1)-p(i,j,1)) + dp_neglect)
757  pbce(i,j,nz) = dp_dh * alpha_star(i,j,nz)
758  enddo
759  do k=nz-1,1,-1 ; do i=isq,ieq+1
760  pbce(i,j,k) = pbce(i,j,k+1) + ((p(i,j,k+1)-p(i,j,1)) * c_htot(i,j)) * &
761  (alpha_star(i,j,k) - alpha_star(i,j,k+1))
762  enddo ; enddo
763  enddo
764  else
765  !$OMP parallel do default(shared) private(T_int,S_int,dR_dT,dR_dS,rho_in_situ)
766  do j=jsq,jeq+1
767  call calculate_density(tv%T(:,j,nz), tv%S(:,j,nz), p(:,j,nz+1), &
768  rho_in_situ, isq, ieq-isq+2, tv%eqn_of_state)
769  do i=isq,ieq+1
770  c_htot(i,j) = dp_dh / ((p(i,j,nz+1)-p(i,j,1)) + dp_neglect)
771  pbce(i,j,nz) = dp_dh / rho_in_situ(i)
772  enddo
773  do k=nz-1,1,-1
774  do i=isq,ieq+1
775  t_int(i) = 0.5*(tv%T(i,j,k)+tv%T(i,j,k+1))
776  s_int(i) = 0.5*(tv%S(i,j,k)+tv%S(i,j,k+1))
777  enddo
778  call calculate_density(t_int, s_int, p(:,j,k+1), rho_in_situ, &
779  isq, ieq-isq+2, tv%eqn_of_state)
780  call calculate_density_derivs(t_int, s_int, p(:,j,k+1), dr_dt, dr_ds, &
781  isq, ieq-isq+2, tv%eqn_of_state)
782  do i=isq,ieq+1
783  pbce(i,j,k) = pbce(i,j,k+1) + ((p(i,j,k+1)-p(i,j,1))*c_htot(i,j)) * &
784  ((dr_dt(i)*(tv%T(i,j,k+1)-tv%T(i,j,k)) + &
785  dr_ds(i)*(tv%S(i,j,k+1)-tv%S(i,j,k))) / rho_in_situ(i)**2)
786  enddo
787  enddo
788  enddo
789  endif
790  else ! not use_EOS
791  !$OMP parallel do default(shared)
792  do j=jsq,jeq+1
793  do i=isq,ieq+1
794  c_htot(i,j) = dp_dh / ((p(i,j,nz+1)-p(i,j,1)) + dp_neglect)
795  pbce(i,j,nz) = dp_dh * alpha_lay(nz)
796  enddo
797  do k=nz-1,1,-1 ; do i=isq,ieq+1
798  pbce(i,j,k) = pbce(i,j,k+1) + ((p(i,j,k+1)-p(i,j,1))*c_htot(i,j)) * &
799  dalpha_int(k+1)
800  enddo ; enddo
801  enddo
802  endif ! use_EOS
803 
804  if (gfs_scale < 1.0) then
805  ! Adjust the Montgomery potential to make this a reduced gravity model.
806  !$OMP parallel do default(shared)
807  do j=jsq,jeq+1 ; do i=isq,ieq+1
808  dpbce(i,j) = (gfs_scale - 1.0) * pbce(i,j,1)
809  enddo ; enddo
810  !$OMP parallel do default(shared)
811  do k=1,nz ; do j=jsq,jeq+1 ; do i=isq,ieq+1
812  pbce(i,j,k) = pbce(i,j,k) + dpbce(i,j)
813  enddo ; enddo ; enddo
814  endif
815 
816 end subroutine set_pbce_nonbouss
817 
818 !> Initialize the Montgomery-potential form of PGF control structure
819 subroutine pressureforce_mont_init(Time, G, GV, US, param_file, diag, CS, tides_CSp)
820  type(time_type), target, intent(in) :: time !< Current model time
821  type(ocean_grid_type), intent(in) :: g !< ocean grid structure
822  type(verticalgrid_type), intent(in) :: gv !< Vertical grid structure
823  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
824  type(param_file_type), intent(in) :: param_file !< Parameter file handles
825  type(diag_ctrl), target, intent(inout) :: diag !< Diagnostics control structure
826  type(pressureforce_mont_cs), pointer :: cs !< Montgomery PGF control structure
827  type(tidal_forcing_cs), optional, pointer :: tides_csp !< Tides control structure
828  ! Local variables
829  logical :: use_temperature, use_eos
830 ! This include declares and sets the variable "version".
831 #include "version_variable.h"
832  character(len=40) :: mdl ! This module's name.
833 
834  if (associated(cs)) then
835  call mom_error(warning, "PressureForce_init called with an associated "// &
836  "control structure.")
837  return
838  else ; allocate(cs) ; endif
839 
840  cs%diag => diag ; cs%Time => time
841  if (present(tides_csp)) then
842  if (associated(tides_csp)) cs%tides_CSp => tides_csp
843  endif
844 
845  mdl = "MOM_PressureForce_Mont"
846  call log_version(param_file, mdl, version, "")
847  call get_param(param_file, mdl, "RHO_0", cs%Rho0, &
848  "The mean ocean density used with BOUSSINESQ true to "//&
849  "calculate accelerations and the mass for conservation "//&
850  "properties, or with BOUSSINSEQ false to convert some "//&
851  "parameters from vertical units of m to kg m-2.", &
852  units="kg m-3", default=1035.0)
853  call get_param(param_file, mdl, "TIDES", cs%tides, &
854  "If true, apply tidal momentum forcing.", default=.false.)
855  call get_param(param_file, mdl, "USE_EOS", use_eos, default=.true., &
856  do_not_log=.true.) ! Input for diagnostic use only.
857 
858  if (use_eos) then
859  cs%id_PFu_bc = register_diag_field('ocean_model', 'PFu_bc', diag%axesCuL, time, &
860  'Density Gradient Zonal Pressure Force Accel.', "meter second-2")
861  cs%id_PFv_bc = register_diag_field('ocean_model', 'PFv_bc', diag%axesCvL, time, &
862  'Density Gradient Meridional Pressure Force Accel.', "meter second-2")
863  if (cs%id_PFu_bc > 0) then
864  call safe_alloc_ptr(cs%PFu_bc,g%IsdB,g%IedB,g%jsd,g%jed,g%ke)
865  cs%PFu_bc(:,:,:) = 0.0
866  endif
867  if (cs%id_PFv_bc > 0) then
868  call safe_alloc_ptr(cs%PFv_bc,g%isd,g%ied,g%JsdB,g%JedB,g%ke)
869  cs%PFv_bc(:,:,:) = 0.0
870  endif
871  endif
872 
873  if (cs%tides) then
874  cs%id_e_tidal = register_diag_field('ocean_model', 'e_tidal', diag%axesT1, &
875  time, 'Tidal Forcing Astronomical and SAL Height Anomaly', 'meter', conversion=us%Z_to_m)
876  endif
877 
878  cs%GFS_scale = 1.0
879  if (gv%g_prime(1) /= gv%g_Earth) cs%GFS_scale = gv%g_prime(1) / gv%g_Earth
880 
881  call log_param(param_file, mdl, "GFS / G_EARTH", cs%GFS_scale)
882 
883 end subroutine pressureforce_mont_init
884 
885 !> Deallocates the Montgomery-potential form of PGF control structure
886 subroutine pressureforce_mont_end(CS)
887  type(pressureforce_mont_cs), pointer :: cs !< Control structure for Montgomery potential PGF
888  if (associated(cs)) deallocate(cs)
889 end subroutine pressureforce_mont_end
890 
891 !>\namespace mom_pressureforce_mont
892 !!
893 !! Provides the Boussunesq and non-Boussinesq forms of the horizontal
894 !! accelerations due to pressure gradients using the Montgomery potential. A
895 !! second-order accurate, centered scheme is used. If a split time stepping
896 !! scheme is used, the vertical decomposition into barotropic and baroclinic
897 !! contributions described by Hallberg (J Comp Phys 1997) is used. With a
898 !! nonlinear equation of state, compressibility is added along the lines proposed
899 !! by Sun et al. (JPO 1999), but with compressibility coefficients based on a fit
900 !! to a user-provided reference profile.
901 
902 end module mom_pressureforce_mont
mom_verticalgrid
Provides a transparent vertical ocean grid type and supporting routines.
Definition: MOM_verticalGrid.F90:2
mom_file_parser::log_version
An overloaded interface to log version information about modules.
Definition: MOM_file_parser.F90:109
mom_diag_mediator
The subroutines here provide convenient wrappers to the fms diag_manager interfaces with additional d...
Definition: MOM_diag_mediator.F90:3
mom_variables::thermo_var_ptrs
Pointers to an assortment of thermodynamic fields that may be available, including potential temperat...
Definition: MOM_variables.F90:82
mom_eos
Provides subroutines for quantities specific to the equation of state.
Definition: MOM_EOS.F90:2
mom_file_parser::param_file_type
A structure that can be parsed to read and document run-time parameters.
Definition: MOM_file_parser.F90:54
mom_file_parser::get_param
An overloaded interface to read and log the values of various types of parameters.
Definition: MOM_file_parser.F90:102
mom_unit_scaling::unit_scale_type
Describes various unit conversion factors.
Definition: MOM_unit_scaling.F90:14
mom_diag_mediator::post_data
Make a diagnostic available for averaging or output.
Definition: MOM_diag_mediator.F90:70
mom_tidal_forcing
Tidal contributions to geopotential.
Definition: MOM_tidal_forcing.F90:2
mom_tidal_forcing::tidal_forcing_cs
The control structure for the MOM_tidal_forcing module.
Definition: MOM_tidal_forcing.F90:26
mom_verticalgrid::verticalgrid_type
Describes the vertical ocean grid, including unit conversion factors.
Definition: MOM_verticalGrid.F90:24
mom_variables
Provides transparent structures with groups of MOM6 variables and supporting routines.
Definition: MOM_variables.F90:2
mom_file_parser
The MOM6 facility to parse input files for runtime parameters.
Definition: MOM_file_parser.F90:2
mom_pressureforce_mont
Provides the Montgomery potential form of pressure gradient.
Definition: MOM_PressureForce_Montgomery.F90:2
mom_eos::calculate_density_derivs
Calculate the derivatives of density with temperature and salinity from T, S, and P.
Definition: MOM_EOS.F90:70
mom_grid
Provides the ocean grid type.
Definition: MOM_grid.F90:2
mom_unit_scaling
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
Definition: MOM_unit_scaling.F90:2
mom_pressureforce_mont::pressureforce_mont_cs
Control structure for the Montgomery potential form of pressure gradient.
Definition: MOM_PressureForce_Montgomery.F90:31
mom_file_parser::log_param
An overloaded interface to log the values of various types of parameters.
Definition: MOM_file_parser.F90:96
mom_error_handler
Routines for error handling and I/O management.
Definition: MOM_error_handler.F90:2
mom_grid::ocean_grid_type
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:25
mom_diag_mediator::diag_ctrl
The following data type a list of diagnostic fields an their variants, as well as variables that cont...
Definition: MOM_diag_mediator.F90:239
mom_eos::calculate_density
Calculates density of sea water from T, S and P.
Definition: MOM_EOS.F90:60