MOM6
MOM_thickness_diffuse.F90
1 !> Thickness diffusion (or Gent McWilliams)
3 
4 ! This file is part of MOM6. See LICENSE.md for the license.
5 
6 use mom_debugging, only : hchksum, uvchksum
7 use mom_diag_mediator, only : post_data, query_averaging_enabled, diag_ctrl
8 use mom_diag_mediator, only : register_diag_field, safe_alloc_ptr, time_type
9 use mom_diag_mediator, only : diag_update_remap_grids
10 use mom_domains, only : pass_var, corner, pass_vector
11 use mom_error_handler, only : mom_error, fatal, warning, is_root_pe
14 use mom_grid, only : ocean_grid_type
16 use mom_isopycnal_slopes, only : vert_fill_ts
18 use mom_meke_types, only : meke_type
22 
23 implicit none ; private
24 
25 #include <MOM_memory.h>
26 
27 public thickness_diffuse, thickness_diffuse_init, thickness_diffuse_end
28 ! public vert_fill_TS
29 public thickness_diffuse_get_kh
30 
31 ! A note on unit descriptions in comments: MOM6 uses units that can be rescaled for dimensional
32 ! consistency testing. These are noted in comments with units like Z, H, L, and T, along with
33 ! their mks counterparts with notation like "a velocity [Z T-1 ~> m s-1]". If the units
34 ! vary with the Boussinesq approximation, the Boussinesq variant is given first.
35 
36 !> Control structure for thickness diffusion
37 type, public :: thickness_diffuse_cs ; private
38  real :: khth !< Background interface depth diffusivity [m2 s-1]
39  real :: khth_slope_cff !< Slope dependence coefficient of Khth [m2 s-1]
40  real :: max_khth_cfl !< Maximum value of the diffusive CFL for thickness diffusion
41  real :: khth_min !< Minimum value of Khth [m2 s-1]
42  real :: khth_max !< Maximum value of Khth [m2 s-1], or 0 for no max
43  real :: slope_max !< Slopes steeper than slope_max are limited in some way [nondim].
44  real :: kappa_smooth !< Vertical diffusivity used to interpolate more
45  !! sensible values of T & S into thin layers [Z2 s-1 ~> m2 s-1].
46  logical :: thickness_diffuse !< If true, interfaces heights are diffused.
47  logical :: use_fgnv_streamfn !< If true, use the streamfunction formulation of
48  !! Ferrari et al., 2010, which effectively emphasizes
49  !! graver vertical modes by smoothing in the vertical.
50  real :: fgnv_scale !< A coefficient scaling the vertical smoothing term in the
51  !! Ferrari et al., 2010, streamfunction formulation [nondim].
52  real :: fgnv_c_min !< A minimum wave speed used in the Ferrari et al., 2010,
53  !! streamfunction formulation [m s-1].
54  real :: n2_floor !< A floor for Brunt-Vasaila frequency in the Ferrari et al., 2010,
55  !! streamfunction formulation [s-2].
56  logical :: detangle_interfaces !< If true, add 3-d structured interface height
57  !! diffusivities to horizontally smooth jagged layers.
58  real :: detangle_time !< If detangle_interfaces is true, this is the
59  !! timescale over which maximally jagged grid-scale
60  !! thickness variations are suppressed [s]. This must be
61  !! longer than DT, or 0 (the default) to use DT.
62  integer :: nkml !< number of layers within mixed layer
63  logical :: debug !< write verbose checksums for debugging purposes
64  logical :: use_gme_thickness_diffuse !< If true, passes GM coefficients to MOM_hor_visc for use
65  !! with GME closure.
66  logical :: meke_geometric !< If true, uses the GM coefficient formulation from the GEOMETRIC
67  !! framework (Marshall et al., 2012)
68  real :: meke_geometric_alpha!< The nondimensional coefficient governing the efficiency of
69  !! the GEOMETRIC thickness difussion [nondim]
70  real :: meke_geometric_epsilon !< Minimum Eady growth rate for the GEOMETRIC thickness
71  !! diffusivity [s-1].
72  logical :: use_kh_in_meke !< If true, uses the thickness diffusivity calculated here to diffuse MEKE.
73  logical :: gm_src_alt !< If true, use the GM energy conversion form S^2*N^2*kappa rather
74  !! than the streamfunction for the GM source term.
75  type(diag_ctrl), pointer :: diag => null() !< structure used to regulate timing of diagnostics
76  real, pointer :: gmwork(:,:) => null() !< Work by thickness diffusivity [W m-2]
77  real, pointer :: diagslopex(:,:,:) => null() !< Diagnostic: zonal neutral slope [nondim]
78  real, pointer :: diagslopey(:,:,:) => null() !< Diagnostic: zonal neutral slope [nondim]
79 
80  real, dimension(:,:,:), pointer :: &
81  kh_u_gme => null(), & !< interface height diffusivities in u-columns (m2 s-1)
82  kh_v_gme => null() !< interface height diffusivities in v-columns (m2 s-1)
83 
84  !>@{
85  !! Diagnostic identifier
86  integer :: id_uhgm = -1, id_vhgm = -1, id_gmwork = -1
87  integer :: id_kh_u = -1, id_kh_v = -1, id_kh_t = -1
88  integer :: id_kh_u1 = -1, id_kh_v1 = -1, id_kh_t1 = -1
89  integer :: id_slope_x = -1, id_slope_y = -1
90  integer :: id_sfn_unlim_x = -1, id_sfn_unlim_y = -1, id_sfn_x = -1, id_sfn_y = -1
91  !>@}
92 end type thickness_diffuse_cs
93 
94 contains
95 
96 !> Calculates thickness diffusion coefficients and applies thickness diffusion to layer
97 !! thicknesses, h. Diffusivities are limited to ensure stability.
98 !! Also returns along-layer mass fluxes used in the continuity equation.
99 subroutine thickness_diffuse(h, uhtr, vhtr, tv, dt, G, GV, US, MEKE, VarMix, CDp, CS)
100  type(ocean_grid_type), intent(in) :: g !< Ocean grid structure
101  type(verticalgrid_type), intent(in) :: gv !< Vertical grid structure
102  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
103  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(inout) :: h !< Layer thickness [H ~> m or kg m-2]
104  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), intent(inout) :: uhtr !< Accumulated zonal mass flux
105  !! [m2 H ~> m3 or kg]
106  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), intent(inout) :: vhtr !< Accumulated meridional mass flux
107  !! [m2 H ~> m3 or kg]
108  type(thermo_var_ptrs), intent(in) :: tv !< Thermodynamics structure
109  real, intent(in) :: dt !< Time increment [s]
110  type(meke_type), pointer :: meke !< MEKE control structure
111  type(varmix_cs), pointer :: varmix !< Variable mixing coefficients
112  type(cont_diag_ptrs), intent(inout) :: cdp !< Diagnostics for the continuity equation
113  type(thickness_diffuse_cs), pointer :: cs !< Control structure for thickness diffusion
114  ! Local variables
115  real :: e(szi_(g), szj_(g), szk_(g)+1) ! heights of interfaces, relative to mean
116  ! sea level [Z ~> m], positive up.
117  real :: uhd(szib_(g), szj_(g), szk_(g)) ! Diffusive u*h fluxes [m2 H s-1 ~> m3 s-1 or kg s-1]
118  real :: vhd(szi_(g), szjb_(g), szk_(g)) ! Diffusive v*h fluxes [m2 H s-1 ~> m3 s-1 or kg s-1]
119 
120  real, dimension(SZIB_(G), SZJ_(G), SZK_(G)+1) :: &
121  kh_u, & ! interface height diffusivities in u-columns [m2 s-1]
122  int_slope_u ! A nondimensional ratio from 0 to 1 that gives the relative
123  ! weighting of the interface slopes to that calculated also
124  ! using density gradients at u points. The physically correct
125  ! slopes occur at 0, while 1 is used for numerical closures.
126  real, dimension(SZI_(G), SZJB_(G), SZK_(G)+1) :: &
127  kh_v, & ! interface height diffusivities in v-columns [m2 s-1]
128  int_slope_v ! A nondimensional ratio from 0 to 1 that gives the relative
129  ! weighting of the interface slopes to that calculated also
130  ! using density gradients at v points. The physically correct
131  ! slopes occur at 0, while 1 is used for numerical closures.
132  real, dimension(SZI_(G), SZJ_(G), SZK_(G)) :: &
133  kh_t ! diagnosed diffusivity at tracer points [m2 s-1]
134 
135  real, dimension(SZIB_(G), SZJ_(G)) :: &
136  kh_u_cfl ! The maximum stable interface height diffusivity at u grid points [m2 s-1]
137  real, dimension(SZI_(G), SZJB_(G)) :: &
138  kh_v_cfl ! The maximum stable interface height diffusivity at v grid points [m2 s-1]
139  real :: khth_loc_u(szib_(g), szj_(g))
140  real :: khth_loc(szib_(g), szjb_(g)) ! locally calculated thickness diffusivity [m2 s-1]
141  real :: h_neglect ! A thickness that is so small it is usually lost
142  ! in roundoff and can be neglected [H ~> m or kg m-2].
143  real, dimension(:,:), pointer :: cg1 => null() !< Wave speed [m s-1]
144  logical :: use_varmix, resoln_scaled, use_stored_slopes, khth_use_ebt_struct, use_visbeck
145  logical :: use_qg_leith
146  integer :: i, j, k, is, ie, js, je, nz
147  real :: hu(szi_(g), szj_(g)) ! u-thickness [H ~> m or kg m-2]
148  real :: hv(szi_(g), szj_(g)) ! v-thickness [H ~> m or kg m-2]
149  real :: kh_u_lay(szi_(g), szj_(g)) ! layer ave thickness diffusivities [m2 s-1]
150  real :: kh_v_lay(szi_(g), szj_(g)) ! layer ave thickness diffusivities [m2 s-1]
151 
152  if (.not. associated(cs)) call mom_error(fatal, "MOM_thickness_diffuse:"// &
153  "Module must be initialized before it is used.")
154 
155  if ((.not.cs%thickness_diffuse) .or. &
156  .not.( cs%Khth > 0.0 .or. associated(varmix) .or. associated(meke) ) ) return
157 
158  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
159  h_neglect = gv%H_subroundoff
160 
161  if (associated(meke)) then
162  if (associated(meke%GM_src)) then
163  do j=js,je ; do i=is,ie ; meke%GM_src(i,j) = 0. ; enddo ; enddo
164  endif
165  endif
166 
167  use_varmix = .false. ; resoln_scaled = .false. ; use_stored_slopes = .false.
168  khth_use_ebt_struct = .false. ; use_visbeck = .false. ; use_qg_leith = .false.
169 
170  if (associated(varmix)) then
171  use_varmix = varmix%use_variable_mixing .and. (cs%KHTH_Slope_Cff > 0.)
172  resoln_scaled = varmix%Resoln_scaled_KhTh
173  use_stored_slopes = varmix%use_stored_slopes
174  khth_use_ebt_struct = varmix%khth_use_ebt_struct
175  use_visbeck = varmix%use_Visbeck
176  use_qg_leith = varmix%use_QG_Leith_GM
177  if (associated(varmix%cg1)) cg1 => varmix%cg1
178  else
179  cg1 => null()
180  endif
181 
182 
183 !$OMP parallel do default(none) shared(is,ie,js,je,KH_u_CFL,dt,G,CS)
184  do j=js,je ; do i=is-1,ie
185  kh_u_cfl(i,j) = (0.25*cs%max_Khth_CFL) / &
186  (dt*(g%IdxCu(i,j)*g%IdxCu(i,j) + g%IdyCu(i,j)*g%IdyCu(i,j)))
187  enddo ; enddo
188 !$OMP parallel do default(none) shared(is,ie,js,je,KH_v_CFL,dt,G,CS)
189  do j=js-1,je ; do i=is,ie
190  kh_v_cfl(i,j) = (0.25*cs%max_Khth_CFL) / &
191  (dt*(g%IdxCv(i,j)*g%IdxCv(i,j) + g%IdyCv(i,j)*g%IdyCv(i,j)))
192  enddo ; enddo
193 
194  ! Calculates interface heights, e, in [Z ~> m].
195  call find_eta(h, tv, g, gv, us, e, halo_size=1)
196 
197  ! Set the diffusivities.
198 !$OMP parallel default(none) shared(is,ie,js,je,Khth_Loc_u,CS,use_VarMix,VarMix, &
199 !$OMP MEKE,Resoln_scaled,KH_u, &
200 !$OMP KH_u_CFL,nz,Khth_Loc,KH_v,KH_v_CFL,int_slope_u, &
201 !$OMP int_slope_v,khth_use_ebt_struct)
202 !$OMP do
203  do j=js,je; do i=is-1,ie
204  khth_loc_u(i,j) = cs%Khth
205  enddo ; enddo
206 
207  if (use_varmix) then
208 !$OMP do
209  if (use_visbeck) then
210  do j=js,je ; do i=is-1,ie
211  khth_loc_u(i,j) = khth_loc_u(i,j) + cs%KHTH_Slope_Cff*varmix%L2u(i,j)*varmix%SN_u(i,j)
212  enddo ; enddo
213  endif
214  endif
215 
216  if (associated(meke)) then ; if (associated(meke%Kh)) then
217 !$OMP do
218  if (cs%MEKE_GEOMETRIC) then
219  do j=js,je ; do i=is-1,ie
220  khth_loc_u(i,j) = khth_loc_u(i,j) + &
221  g%mask2dCu(i,j) * cs%MEKE_GEOMETRIC_alpha * 0.5*(meke%MEKE(i,j)+meke%MEKE(i+1,j)) / &
222  (varmix%SN_u(i,j) + cs%MEKE_GEOMETRIC_epsilon)
223  enddo ; enddo
224  else
225  do j=js,je ; do i=is-1,ie
226  khth_loc_u(i,j) = khth_loc_u(i,j) + meke%KhTh_fac*sqrt(meke%Kh(i,j)*meke%Kh(i+1,j))
227  enddo ; enddo
228  endif
229  endif ; endif
230 
231  if (resoln_scaled) then
232 !$OMP do
233  do j=js,je; do i=is-1,ie
234  khth_loc_u(i,j) = khth_loc_u(i,j) * varmix%Res_fn_u(i,j)
235  enddo ; enddo
236  endif
237 
238  if (cs%Khth_Max > 0) then
239 !$OMP do
240  do j=js,je; do i=is-1,ie
241  khth_loc_u(i,j) = max(cs%Khth_min, min(khth_loc_u(i,j),cs%Khth_Max))
242  enddo ; enddo
243  else
244 !$OMP do
245  do j=js,je; do i=is-1,ie
246  khth_loc_u(i,j) = max(cs%Khth_min, khth_loc_u(i,j))
247  enddo ; enddo
248  endif
249 !$OMP do
250  do j=js,je; do i=is-1,ie
251  kh_u(i,j,1) = min(kh_u_cfl(i,j), khth_loc_u(i,j))
252  enddo ; enddo
253 
254  if (khth_use_ebt_struct) then
255 !$OMP do
256  do k=2,nz+1 ; do j=js,je ; do i=is-1,ie
257  kh_u(i,j,k) = kh_u(i,j,1) * 0.5 * ( varmix%ebt_struct(i,j,k-1) + varmix%ebt_struct(i+1,j,k-1) )
258  enddo ; enddo ; enddo
259  else
260 !$OMP do
261  do k=2,nz+1 ; do j=js,je ; do i=is-1,ie
262  kh_u(i,j,k) = kh_u(i,j,1)
263  enddo ; enddo ; enddo
264  endif
265 
266  if (use_varmix) then
267 !$OMP do
268  if (use_qg_leith) then
269  do k=1,nz ; do j=js,je ; do i=is-1,ie
270  kh_u(i,j,k) = varmix%KH_u_QG(i,j,k)
271  enddo ; enddo ; enddo
272  endif
273  endif
274 
275 !$OMP do
276  if (cs%use_GME_thickness_diffuse) then
277  do k=1,nz+1 ; do j=js,je ; do i=is-1,ie
278  cs%KH_u_GME(i,j,k) = kh_u(i,j,k)
279  enddo ; enddo ; enddo
280  endif
281 
282 !$OMP do
283  do j=js-1,je ; do i=is,ie
284  khth_loc(i,j) = cs%Khth
285  enddo ; enddo
286 
287  if (use_varmix) then
288 !$OMP do
289  if (use_visbeck) then
290  do j=js-1,je ; do i=is,ie
291  khth_loc(i,j) = khth_loc(i,j) + cs%KHTH_Slope_Cff*varmix%L2v(i,j)*varmix%SN_v(i,j)
292  enddo ; enddo
293  endif
294  endif
295  if (associated(meke)) then ; if (associated(meke%Kh)) then
296 !$OMP do
297  if (cs%MEKE_GEOMETRIC) then
298  do j=js-1,je ; do i=is,ie
299  khth_loc(i,j) = khth_loc(i,j) + &
300  g%mask2dCv(i,j) * cs%MEKE_GEOMETRIC_alpha * 0.5*(meke%MEKE(i,j)+meke%MEKE(i,j+1)) / &
301  (varmix%SN_v(i,j) + cs%MEKE_GEOMETRIC_epsilon)
302  enddo ; enddo
303  else
304  do j=js-1,je ; do i=is,ie
305  khth_loc(i,j) = khth_loc(i,j) + meke%KhTh_fac*sqrt(meke%Kh(i,j)*meke%Kh(i,j+1))
306  enddo ; enddo
307  endif
308  endif ; endif
309 
310  if (resoln_scaled) then
311 !$OMP do
312  do j=js-1,je ; do i=is,ie
313  khth_loc(i,j) = khth_loc(i,j) * varmix%Res_fn_v(i,j)
314  enddo ; enddo
315  endif
316 
317  if (cs%Khth_Max > 0) then
318 !$OMP do
319  do j=js-1,je ; do i=is,ie
320  khth_loc(i,j) = max(cs%Khth_min, min(khth_loc(i,j),cs%Khth_Max))
321  enddo ; enddo
322  else
323 !$OMP do
324  do j=js-1,je ; do i=is,ie
325  khth_loc(i,j) = max(cs%Khth_min, khth_loc(i,j))
326  enddo ; enddo
327  endif
328 
329  if (cs%max_Khth_CFL > 0.0) then
330 !$OMP do
331  do j=js-1,je ; do i=is,ie
332  kh_v(i,j,1) = min(kh_v_cfl(i,j), khth_loc(i,j))
333  enddo ; enddo
334  endif
335 
336  if (khth_use_ebt_struct) then
337 !$OMP do
338  do k=2,nz+1 ; do j=js-1,je ; do i=is,ie
339  kh_v(i,j,k) = kh_v(i,j,1) * 0.5 * ( varmix%ebt_struct(i,j,k-1) + varmix%ebt_struct(i,j+1,k-1) )
340  enddo ; enddo ; enddo
341  else
342 !$OMP do
343  do k=2,nz+1 ; do j=js-1,je ; do i=is,ie
344  kh_v(i,j,k) = kh_v(i,j,1)
345  enddo ; enddo ; enddo
346  endif
347 
348  if (use_varmix) then
349 !$OMP do
350  if (use_qg_leith) then
351  do k=1,nz ; do j=js-1,je ; do i=is,ie
352  kh_v(i,j,k) = varmix%KH_v_QG(i,j,k)
353  enddo ; enddo ; enddo
354  endif
355  endif
356 
357 !$OMP do
358  if (cs%use_GME_thickness_diffuse) then
359  do k=1,nz+1 ; do j=js-1,je ; do i=is,ie
360  cs%KH_v_GME(i,j,k) = kh_v(i,j,k)
361  enddo ; enddo ; enddo
362  endif
363 
364  if (associated(meke)) then ; if (associated(meke%Kh)) then
365 !$OMP do
366  if (cs%MEKE_GEOMETRIC) then
367  do j=js,je ; do i=is,ie
368  meke%Kh(i,j) = cs%MEKE_GEOMETRIC_alpha * meke%MEKE(i,j) / &
369  (0.25*(varmix%SN_u(i,j)+varmix%SN_u(i-1,j)+varmix%SN_v(i,j)+varmix%SN_v(i,j-1)) + &
370  cs%MEKE_GEOMETRIC_epsilon)
371  enddo ; enddo
372  endif
373  endif ; endif
374 
375 
376 !$OMP do
377  do k=1,nz+1 ; do j=js,je ; do i=is-1,ie ; int_slope_u(i,j,k) = 0.0 ; enddo ; enddo ; enddo
378 !$OMP do
379  do k=1,nz+1 ; do j=js-1,je ; do i=is,ie ; int_slope_v(i,j,k) = 0.0 ; enddo ; enddo ; enddo
380 !$OMP end parallel
381 
382  if (cs%detangle_interfaces) then
383  call add_detangling_kh(h, e, kh_u, kh_v, kh_u_cfl, kh_v_cfl, tv, dt, g, gv, us, &
384  cs, int_slope_u, int_slope_v)
385  endif
386 
387  if (cs%debug) then
388  call uvchksum("Kh_[uv]", kh_u, kh_v, g%HI,haloshift=0)
389  call uvchksum("int_slope_[uv]", int_slope_u, int_slope_v, g%HI, haloshift=0)
390  call hchksum(h, "thickness_diffuse_1 h", g%HI, haloshift=1, scale=gv%H_to_m)
391  call hchksum(e, "thickness_diffuse_1 e", g%HI, haloshift=1, scale=us%Z_to_m)
392  if (use_stored_slopes) then
393  call uvchksum("VarMix%slope_[xy]", varmix%slope_x, varmix%slope_y, &
394  g%HI, haloshift=0)
395  endif
396  if (associated(tv%eqn_of_state)) then
397  call hchksum(tv%T, "thickness_diffuse T", g%HI, haloshift=1)
398  call hchksum(tv%S, "thickness_diffuse S", g%HI, haloshift=1)
399  endif
400  endif
401 
402  ! Calculate uhD, vhD from h, e, KH_u, KH_v, tv%T/S
403  if (use_stored_slopes) then
404  call thickness_diffuse_full(h, e, kh_u, kh_v, tv, uhd, vhd, cg1, dt, g, gv, us, meke, cs, &
405  int_slope_u, int_slope_v, varmix%slope_x, varmix%slope_y)
406  else
407  call thickness_diffuse_full(h, e, kh_u, kh_v, tv, uhd, vhd, cg1, dt, g, gv, us, meke, cs, &
408  int_slope_u, int_slope_v)
409  endif
410 
411  if (associated(meke) .AND. associated(varmix)) then
412  if (associated(meke%Rd_dx_h) .and. associated(varmix%Rd_dx_h)) then
413 !$OMP parallel do default(none) shared(is,ie,js,je,MEKE,VarMix)
414  do j=js,je ; do i=is,ie
415  meke%Rd_dx_h(i,j) = varmix%Rd_dx_h(i,j)
416  enddo ; enddo
417  endif
418  endif
419 
420  ! offer diagnostic fields for averaging
421  if (query_averaging_enabled(cs%diag)) then
422  if (cs%id_uhGM > 0) call post_data(cs%id_uhGM, uhd, cs%diag)
423  if (cs%id_vhGM > 0) call post_data(cs%id_vhGM, vhd, cs%diag)
424  if (cs%id_GMwork > 0) call post_data(cs%id_GMwork, cs%GMwork, cs%diag)
425  if (cs%id_KH_u > 0) call post_data(cs%id_KH_u, kh_u, cs%diag)
426  if (cs%id_KH_v > 0) call post_data(cs%id_KH_v, kh_v, cs%diag)
427  if (cs%id_KH_u1 > 0) call post_data(cs%id_KH_u1, kh_u(:,:,1), cs%diag)
428  if (cs%id_KH_v1 > 0) call post_data(cs%id_KH_v1, kh_v(:,:,1), cs%diag)
429 
430  ! Diagnose diffusivity at T-cell point. Do simple average, rather than
431  ! thickness-weighted average, in order that KH_t is depth-independent
432  ! in the case where KH_u and KH_v are depth independent. Otherwise,
433  ! if use thickness weighted average, the variations of thickness with
434  ! depth will place a spurious depth dependence to the diagnosed KH_t.
435  if (cs%id_KH_t > 0 .or. cs%id_KH_t1 > 0 .or. cs%Use_KH_in_MEKE) then
436  do k=1,nz
437  ! thicknesses across u and v faces, converted to 0/1 mask
438  ! layer average of the interface diffusivities KH_u and KH_v
439  do j=js,je ; do i=is-1,ie
440  hu(i,j) = 2.0*h(i,j,k)*h(i+1,j,k)/(h(i,j,k)+h(i+1,j,k)+h_neglect)
441  if (hu(i,j) /= 0.0) hu(i,j) = 1.0
442  kh_u_lay(i,j) = 0.5*(kh_u(i,j,k)+kh_u(i,j,k+1))
443  enddo ; enddo
444  do j=js-1,je ; do i=is,ie
445  hv(i,j) = 2.0*h(i,j,k)*h(i,j+1,k)/(h(i,j,k)+h(i,j+1,k)+h_neglect)
446  if (hv(i,j) /= 0.0) hv(i,j) = 1.0
447  kh_v_lay(i,j) = 0.5*(kh_v(i,j,k)+kh_v(i,j,k+1))
448  enddo ; enddo
449  ! diagnose diffusivity at T-point
450  do j=js,je ; do i=is,ie
451  kh_t(i,j,k) = ((hu(i-1,j)*kh_u_lay(i-1,j)+hu(i,j)*kh_u_lay(i,j)) &
452  +(hv(i,j-1)*kh_v_lay(i,j-1)+hv(i,j)*kh_v_lay(i,j))) &
453  / (hu(i-1,j)+hu(i,j)+hv(i,j-1)+hv(i,j)+h_neglect)
454  enddo ; enddo
455  enddo
456 
457  if (cs%Use_KH_in_MEKE) then
458  meke%Kh_diff(:,:) = 0.0
459  do k=1,nz
460  do j=js,je ; do i=is,ie
461  meke%Kh_diff(i,j) = meke%Kh_diff(i,j) + kh_t(i,j,k) * h(i,j,k)
462  enddo; enddo
463  enddo
464 
465  do j=js,je ; do i=is,ie
466  meke%Kh_diff(i,j) = meke%Kh_diff(i,j) / max(1.0,g%bathyT(i,j))
467  enddo ; enddo
468  endif
469 
470  if (cs%id_KH_t > 0) call post_data(cs%id_KH_t, kh_t, cs%diag)
471  if (cs%id_KH_t1 > 0) call post_data(cs%id_KH_t1, kh_t(:,:,1), cs%diag)
472  endif
473 
474  endif
475 
476  !$OMP parallel do default(none) shared(is,ie,js,je,nz,uhtr,uhD,dt,vhtr,CDp,vhD,h,G,GV)
477  do k=1,nz
478  do j=js,je ; do i=is-1,ie
479  uhtr(i,j,k) = uhtr(i,j,k) + uhd(i,j,k)*dt
480  if (associated(cdp%uhGM)) cdp%uhGM(i,j,k) = uhd(i,j,k)
481  enddo ; enddo
482  do j=js-1,je ; do i=is,ie
483  vhtr(i,j,k) = vhtr(i,j,k) + vhd(i,j,k)*dt
484  if (associated(cdp%vhGM)) cdp%vhGM(i,j,k) = vhd(i,j,k)
485  enddo ; enddo
486  do j=js,je ; do i=is,ie
487  h(i,j,k) = h(i,j,k) - dt * g%IareaT(i,j) * &
488  ((uhd(i,j,k) - uhd(i-1,j,k)) + (vhd(i,j,k) - vhd(i,j-1,k)))
489  if (h(i,j,k) < gv%Angstrom_H) h(i,j,k) = gv%Angstrom_H
490  enddo ; enddo
491  enddo
492 
493  ! Whenever thickness changes let the diag manager know, target grids
494  ! for vertical remapping may need to be regenerated.
495  ! This needs to happen after the H update and before the next post_data.
496  call diag_update_remap_grids(cs%diag)
497 
498  if (cs%debug) then
499  call uvchksum("thickness_diffuse [uv]hD", uhd, vhd, &
500  g%HI, haloshift=0, scale=gv%H_to_m)
501  call uvchksum("thickness_diffuse [uv]htr", uhtr, vhtr, &
502  g%HI, haloshift=0, scale=gv%H_to_m)
503  call hchksum(h, "thickness_diffuse h", g%HI, haloshift=0, scale=gv%H_to_m)
504  endif
505 
506 end subroutine thickness_diffuse
507 
508 !> Calculates parameterized layer transports for use in the continuity equation.
509 !! Fluxes are limited to give positive definite thicknesses.
510 !! Called by thickness_diffuse().
511 subroutine thickness_diffuse_full(h, e, Kh_u, Kh_v, tv, uhD, vhD, cg1, dt, G, GV, US, MEKE, &
512  CS, int_slope_u, int_slope_v, slope_x, slope_y)
513  type(ocean_grid_type), intent(in) :: G !< Ocean grid structure
514  type(verticalgrid_type), intent(in) :: GV !< Vertical grid structure
515  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
516  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thickness [H ~> m or kg m-2]
517  real, dimension(SZI_(G),SZJ_(G),SZK_(G)+1), intent(in) :: e !< Interface positions [Z ~> m]
518  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)+1), intent(in) :: Kh_u !< Thickness diffusivity on interfaces
519  !! at u points [m2 s-1]
520  real, dimension(SZI_(G),SZJB_(G),SZK_(G)+1), intent(in) :: Kh_v !< Thickness diffusivity on interfaces
521  !! at v points [m2 s-1]
522  type(thermo_var_ptrs), intent(in) :: tv !< Thermodynamics structure
523  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), intent(out) :: uhD !< Zonal mass fluxes
524  !! [H m2 s-1 ~> m3 s-1 or kg s-1]
525  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), intent(out) :: vhD !< Meridional mass fluxes
526  !! [H m2 s-1 ~> m3 s-1 or kg s-1]
527  real, dimension(:,:), pointer :: cg1 !< Wave speed [m s-1]
528  real, intent(in) :: dt !< Time increment [s]
529  type(meke_type), pointer :: MEKE !< MEKE control structure
530  type(thickness_diffuse_cs), pointer :: CS !< Control structure for thickness diffusion
531  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)+1), optional, intent(in) :: int_slope_u !< Ratio that determine how much of
532  !! the isopycnal slopes are taken directly from
533  !! the interface slopes without consideration of
534  !! density gradients.
535  real, dimension(SZI_(G),SZJB_(G),SZK_(G)+1), optional, intent(in) :: int_slope_v !< Ratio that determine how much of
536  !! the isopycnal slopes are taken directly from
537  !! the interface slopes without consideration of
538  !! density gradients.
539  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)+1), optional, intent(in) :: slope_x !< Isopycnal slope at u-points
540  real, dimension(SZI_(G),SZJB_(G),SZK_(G)+1), optional, intent(in) :: slope_y !< Isopycnal slope at v-points
541  ! Local variables
542  real, dimension(SZI_(G), SZJ_(G), SZK_(G)) :: &
543  T, & ! The temperature (or density) [degC], with the values in
544  ! in massless layers filled vertically by diffusion.
545  s, & ! The filled salinity [ppt], with the values in
546  ! in massless layers filled vertically by diffusion.
547  rho, & ! Density itself [kg m-3], when a nonlinear equation of state is
548  ! not in use.
549  h_avail, & ! The mass available for diffusion out of each face, divided
550  ! by dt [H m2 s-1 ~> m3 s-1 or kg s-1].
551  h_frac ! The fraction of the mass in the column above the bottom
552  ! interface of a layer that is within a layer [nondim]. 0<h_frac<=1
553  real, dimension(SZI_(G), SZJB_(G), SZK_(G)+1) :: &
554  Slope_y_PE, & ! 3D array of neutral slopes at v-points, set equal to Slope (below, nondim)
555  hN2_y_PE ! thickness in m times Brunt-Vaisala freqeuncy at v-points [m s-2],
556  ! used for calculating PE release
557  real, dimension(SZIB_(G), SZJ_(G), SZK_(G)+1) :: &
558  Slope_x_PE, & ! 3D array of neutral slopes at u-points, set equal to Slope (below, nondim)
559  hN2_x_PE ! thickness in m times Brunt-Vaisala freqeuncy at u-points [m s-2]
560  ! used for calculating PE release
561  real, dimension(SZI_(G), SZJ_(G), SZK_(G)+1) :: &
562  pres, & ! The pressure at an interface [Pa].
563  h_avail_rsum ! The running sum of h_avail above an interface [H m2 s-1 ~> m3 s-1 or kg s-1].
564  real, dimension(SZIB_(G)) :: &
565  drho_dT_u, & ! The derivative of density with temperature at u points [kg m-3 degC-1]
566  drho_dS_u ! The derivative of density with salinity at u points [kg m-3 ppt-1].
567  real, dimension(SZI_(G)) :: &
568  drho_dT_v, & ! The derivative of density with temperature at v points [kg m-3 degC-1]
569  drho_dS_v ! The derivative of density with salinity at v points [kg m-3 ppt-1].
570  real :: uhtot(SZIB_(G), SZJ_(G)) ! The vertical sum of uhD [H m2 s-1 ~> m3 s-1 or kg s-1].
571  real :: vhtot(SZI_(G), SZJB_(G)) ! The vertical sum of vhD [H m2 s-1 ~> m3 s-1 or kg s-1].
572  real, dimension(SZIB_(G)) :: &
573  T_u, & ! Temperature on the interface at the u-point [degC].
574  S_u, & ! Salinity on the interface at the u-point [ppt].
575  pres_u ! Pressure on the interface at the u-point [Pa].
576  real, dimension(SZI_(G)) :: &
577  T_v, & ! Temperature on the interface at the v-point [degC].
578  S_v, & ! Salinity on the interface at the v-point [ppt].
579  pres_v ! Pressure on the interface at the v-point [Pa].
580  real :: Work_u(SZIB_(G), SZJ_(G)) ! The work being done by the thickness
581  real :: Work_v(SZI_(G), SZJB_(G)) ! diffusion integrated over a cell [W].
582  real :: Work_h ! The work averaged over an h-cell [W m-2].
583  real :: PE_release_h ! The amount of potential energy released by GM, averaged over an h-cell [m3 s-3].
584  ! The calculation is equal to h * S^2 * N^2 * kappa_GM.
585  real :: I4dt ! 1 / 4 dt [s-1].
586  real :: drdiA, drdiB ! Along layer zonal- and meridional- potential density
587  real :: drdjA, drdjB ! gradients in the layers above (A) and below(B) the
588  ! interface times the grid spacing [kg m-3].
589  real :: drdkL, drdkR ! Vertical density differences across an interface [kg m-3].
590  real :: drdi_u(SZIB_(G), SZK_(G)+1) ! Copy of drdi at u-points [kg m-3].
591  real :: drdj_v(SZI_(G), SZK_(G)+1) ! Copy of drdj at v-points [kg m-3].
592  real :: drdkDe_u(SZIB_(G),SZK_(G)+1) ! Lateral difference of product of drdk and e at u-points
593  ! [Z kg m-3 ~> kg m-2].
594  real :: drdkDe_v(SZI_(G),SZK_(G)+1) ! Lateral difference of product of drdk and e at v-points
595  ! [Z kg m-3 ~> kg m-2].
596  real :: hg2A, hg2B, hg2L, hg2R ! Squares of geometric mean thicknesses [H2 ~> m2 or kg2 m-4].
597  real :: haA, haB, haL, haR ! Arithmetic mean thicknesses [H ~> m or kg m-2].
598  real :: dzaL, dzaR ! Temporary thicknesses [Z ~> m].
599  real :: wtA, wtB, wtL, wtR ! Unscaled weights, with various units.
600  real :: drdx, drdy ! Zonal and meridional density gradients [kg m-4].
601  real :: drdz ! Vertical density gradient [kg m-3 Z-1 ~> kg m-4].
602  real :: h_harm ! Harmonic mean layer thickness [H ~> m or kg m-2].
603  real :: c2_h_u(SZIB_(G), SZK_(G)+1) ! Wave speed squared divided by h at u-points [m2 Z-1 s-2 ~> m s-2].
604  real :: c2_h_v(SZI_(G), SZK_(G)+1) ! Wave speed squared divided by h at v-points [m2 Z-1 s-2 ~> m s-2].
605  real :: hN2_u(SZIB_(G), SZK_(G)+1) ! Thickness in m times N2 at interfaces above u-points [m2 Z-1 s-2 ~> m s-2].
606  real :: hN2_v(SZI_(G), SZK_(G)+1) ! Thickness in m times N2 at interfaces above v-points [m2 Z-1 s-2 ~> m s-2].
607  real :: Sfn_est ! A preliminary estimate (before limiting) of the overturning
608  ! streamfunction [Z m2 s-1 ~> m3 s-1].
609  real :: Sfn_unlim_u(SZIB_(G), SZK_(G)+1) ! Streamfunction for u-points [Z m2 s-1 ~> m3 s-1].
610  real :: Sfn_unlim_v(SZI_(G), SZK_(G)+1) ! Streamfunction for v-points [Z m2 s-1 ~> m3 s-1].
611  real :: slope2_Ratio_u(SZIB_(G), SZK_(G)+1) ! The ratio of the slope squared to slope_max squared.
612  real :: slope2_Ratio_v(SZI_(G), SZK_(G)+1) ! The ratio of the slope squared to slope_max squared.
613  real :: Sfn_in_h ! The overturning streamfunction [H m2 s-1 ~> m3 s-1 or kg s-1] (note that
614  ! the units are different from other Sfn vars).
615  real :: Sfn_safe ! The streamfunction that goes linearly back to 0 at the surface. This is a
616  ! good thing to use when the slope is so large as to be meaningless [Z m2 s-1 ~> m3 s-1].
617  real :: Slope ! The slope of density surfaces, calculated in a way
618  ! that is always between -1 and 1, nondimensional.
619  real :: mag_grad2 ! The squared magnitude of the 3-d density gradient [kg2 m-8].
620  real :: I_slope_max2 ! The inverse of slope_max squared, nondimensional.
621  real :: h_neglect ! A thickness that is so small it is usually lost
622  ! in roundoff and can be neglected [H ~> m or kg m-2].
623  real :: h_neglect2 ! h_neglect^2 [H2 ~> m2 or kg2 m-4].
624  real :: dz_neglect ! A thickness [Z ~> m], that is so small it is usually lost
625  ! in roundoff and can be neglected [Z ~> m].
626  real :: G_scale ! The gravitational acceleration times some unit conversion
627  ! factors [m3 Z-1 H-1 s-2 ~> m s-2 or m4 kg-1 s-2].
628  logical :: use_EOS ! If true, density is calculated from T & S using an
629  ! equation of state.
630  logical :: find_work ! If true, find the change in energy due to the fluxes.
631  integer :: nk_linear ! The number of layers over which the streamfunction goes to 0.
632  real :: G_rho0 ! g/Rho0 [m5 Z-1 s-2 ~> m4 s-2].
633  real :: N2_floor ! A floor for N2 to avoid degeneracy in the elliptic solver
634  ! times unit conversion factors [s-2 m2 Z-2 ~> s-2]
635  real, dimension(SZIB_(G), SZJ_(G), SZK_(G)+1) :: diag_sfn_x, diag_sfn_unlim_x ! Diagnostics
636  real, dimension(SZI_(G), SZJB_(G), SZK_(G)+1) :: diag_sfn_y, diag_sfn_unlim_y ! Diagnostics
637  logical :: present_int_slope_u, present_int_slope_v
638  logical :: present_slope_x, present_slope_y, calc_derivatives
639  integer :: is, ie, js, je, nz, IsdB
640  integer :: i, j, k
641  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke ; isdb = g%IsdB
642 
643  i4dt = 0.25 / dt
644  i_slope_max2 = 1.0 / (cs%slope_max**2)
645  g_scale = gv%g_Earth*us%L_to_m**2*us%s_to_T**2 * gv%H_to_m
646  h_neglect = gv%H_subroundoff ; h_neglect2 = h_neglect**2
647  dz_neglect = gv%H_subroundoff*gv%H_to_Z
648  g_rho0 = gv%g_Earth*us%L_to_m**2*us%s_to_T**2 / gv%Rho0
649  n2_floor = cs%N2_floor*us%Z_to_m**2
650 
651  use_eos = associated(tv%eqn_of_state)
652  present_int_slope_u = PRESENT(int_slope_u)
653  present_int_slope_v = PRESENT(int_slope_v)
654  present_slope_x = PRESENT(slope_x)
655  present_slope_y = PRESENT(slope_y)
656 
657  nk_linear = max(gv%nkml, 1)
658 
659  slope_x_pe(:,:,:) = 0.0
660  slope_y_pe(:,:,:) = 0.0
661  hn2_x_pe(:,:,:) = 0.0
662  hn2_y_pe(:,:,:) = 0.0
663 
664  find_work = .false.
665  if (associated(meke)) find_work = associated(meke%GM_src)
666  find_work = (associated(cs%GMwork) .or. find_work)
667 
668  if (use_eos) then
669  call vert_fill_ts(h, tv%T, tv%S, cs%kappa_smooth*dt, t, s, g, gv, 1, larger_h_denom=.true.)
670  endif
671 
672  if (cs%use_FGNV_streamfn .and. .not. associated(cg1)) call mom_error(fatal, &
673  "cg1 must be associated when using FGNV streamfunction.")
674 
675 !$OMP parallel default(none) shared(is,ie,js,je,h_avail_rsum,pres,h_avail,I4dt, &
676 !$OMP G,GV,h,h_frac,nz,uhtot,Work_u,vhtot,Work_v, &
677 !$OMP diag_sfn_x, diag_sfn_y, diag_sfn_unlim_x, diag_sfn_unlim_y )
678  ! Find the maximum and minimum permitted streamfunction.
679 !$OMP do
680  do j=js-1,je+1 ; do i=is-1,ie+1
681  h_avail_rsum(i,j,1) = 0.0
682  pres(i,j,1) = 0.0 ! ### This should be atmospheric pressure.
683 
684  h_avail(i,j,1) = max(i4dt*g%areaT(i,j)*(h(i,j,1)-gv%Angstrom_H),0.0)
685  h_avail_rsum(i,j,2) = h_avail(i,j,1)
686  h_frac(i,j,1) = 1.0
687  pres(i,j,2) = pres(i,j,1) + gv%H_to_Pa*h(i,j,1)
688  enddo ; enddo
689 !$OMP do
690  do j=js-1,je+1
691  do k=2,nz ; do i=is-1,ie+1
692  h_avail(i,j,k) = max(i4dt*g%areaT(i,j)*(h(i,j,k)-gv%Angstrom_H),0.0)
693  h_avail_rsum(i,j,k+1) = h_avail_rsum(i,j,k) + h_avail(i,j,k)
694  h_frac(i,j,k) = 0.0 ; if (h_avail(i,j,k) > 0.0) &
695  h_frac(i,j,k) = h_avail(i,j,k) / h_avail_rsum(i,j,k+1)
696  pres(i,j,k+1) = pres(i,j,k) + gv%H_to_Pa*h(i,j,k)
697  enddo ; enddo
698  enddo
699 !$OMP do
700  do j=js,je ; do i=is-1,ie
701  uhtot(i,j) = 0.0 ; work_u(i,j) = 0.0
702  diag_sfn_x(i,j,1) = 0.0 ; diag_sfn_unlim_x(i,j,1) = 0.0
703  diag_sfn_x(i,j,nz+1) = 0.0 ; diag_sfn_unlim_x(i,j,nz+1) = 0.0
704  enddo ; enddo
705 !$OMP do
706  do j=js-1,je ; do i=is,ie
707  vhtot(i,j) = 0.0 ; work_v(i,j) = 0.0
708  diag_sfn_y(i,j,1) = 0.0 ; diag_sfn_unlim_y(i,j,1) = 0.0
709  diag_sfn_y(i,j,nz+1) = 0.0 ; diag_sfn_unlim_y(i,j,nz+1) = 0.0
710  enddo ; enddo
711 !$OMP end parallel
712 
713 !$OMP parallel do default(none) shared(nz,is,ie,js,je,find_work,use_EOS,G,GV,US,pres,T,S, &
714 !$OMP nk_linear,IsdB,tv,h,h_neglect,e,dz_neglect, &
715 !$OMP I_slope_max2,h_neglect2,present_int_slope_u, &
716 !$OMP int_slope_u,KH_u,uhtot,h_frac,h_avail_rsum, &
717 !$OMP uhD,h_avail,G_scale,work_u,CS,slope_x,cg1, &
718 !$OMP diag_sfn_x, diag_sfn_unlim_x,N2_floor, &
719 !$OMP present_slope_x,G_rho0) &
720 !$OMP private(drdiA,drdiB,drdkL,drdkR,pres_u,T_u,S_u, &
721 !$OMP drho_dT_u,drho_dS_u,hg2A,hg2B,hg2L,hg2R,haA, &
722 !$OMP haB,haL,haR,dzaL,dzaR,wtA,wtB,wtL,wtR,drdz, &
723 !$OMP drdx,mag_grad2,Slope,slope2_Ratio_u,hN2_u, &
724 !$OMP Sfn_unlim_u,drdi_u,drdkDe_u,h_harm,c2_h_u, &
725 !$OMP Sfn_safe,Sfn_est,Sfn_in_h,calc_derivatives)
726  do j=js,je
727  do i=is-1,ie ; hn2_u(i,1) = 0. ; hn2_u(i,nz+1) = 0. ; enddo
728  do k=nz,2,-1
729  if (find_work .and. .not.(use_eos)) then
730  drdia = 0.0 ; drdib = 0.0
731  drdkl = gv%Rlay(k)-gv%Rlay(k-1) ; drdkr = drdkl
732  endif
733 
734  calc_derivatives = use_eos .and. (k >= nk_linear) .and. &
735  (find_work .or. .not. present_slope_x .or. cs%use_FGNV_streamfn)
736 
737  ! Calculate the zonal fluxes and gradients.
738  if (calc_derivatives) then
739  do i=is-1,ie
740  pres_u(i) = 0.5*(pres(i,j,k) + pres(i+1,j,k))
741  t_u(i) = 0.25*((t(i,j,k) + t(i+1,j,k)) + (t(i,j,k-1) + t(i+1,j,k-1)))
742  s_u(i) = 0.25*((s(i,j,k) + s(i+1,j,k)) + (s(i,j,k-1) + s(i+1,j,k-1)))
743  enddo
744  call calculate_density_derivs(t_u, s_u, pres_u, drho_dt_u, &
745  drho_ds_u, (is-isdb+1)-1, ie-is+2, tv%eqn_of_state)
746  endif
747 
748  do i=is-1,ie
749  if (calc_derivatives) then
750  ! Estimate the horizontal density gradients along layers.
751  drdia = drho_dt_u(i) * (t(i+1,j,k-1)-t(i,j,k-1)) + &
752  drho_ds_u(i) * (s(i+1,j,k-1)-s(i,j,k-1))
753  drdib = drho_dt_u(i) * (t(i+1,j,k)-t(i,j,k)) + &
754  drho_ds_u(i) * (s(i+1,j,k)-s(i,j,k))
755 
756  ! Estimate the vertical density gradients times the grid spacing.
757  drdkl = (drho_dt_u(i) * (t(i,j,k)-t(i,j,k-1)) + &
758  drho_ds_u(i) * (s(i,j,k)-s(i,j,k-1)))
759  drdkr = (drho_dt_u(i) * (t(i+1,j,k)-t(i+1,j,k-1)) + &
760  drho_ds_u(i) * (s(i+1,j,k)-s(i+1,j,k-1)))
761  drdkde_u(i,k) = drdkr * e(i+1,j,k) - drdkl * e(i,j,k)
762  elseif (find_work) then ! This is used in pure stacked SW mode
763  drdkde_u(i,k) = drdkr * e(i+1,j,k) - drdkl * e(i,j,k)
764  endif
765 
766  if (find_work) drdi_u(i,k) = drdib
767 
768  if (k > nk_linear) then
769  if (use_eos) then
770  if (cs%use_FGNV_streamfn .or. find_work .or. .not.present_slope_x) then
771  hg2l = h(i,j,k-1)*h(i,j,k) + h_neglect2
772  hg2r = h(i+1,j,k-1)*h(i+1,j,k) + h_neglect2
773  hal = 0.5*(h(i,j,k-1) + h(i,j,k)) + h_neglect
774  har = 0.5*(h(i+1,j,k-1) + h(i+1,j,k)) + h_neglect
775  if (gv%Boussinesq) then
776  dzal = hal * gv%H_to_Z ; dzar = har * gv%H_to_Z
777  else
778  dzal = 0.5*(e(i,j,k-1) - e(i,j,k+1)) + dz_neglect
779  dzar = 0.5*(e(i+1,j,k-1) - e(i+1,j,k+1)) + dz_neglect
780  endif
781  ! Use the harmonic mean thicknesses to weight the horizontal gradients.
782  ! These unnormalized weights have been rearranged to minimize divisions.
783  wtl = hg2l*(har*dzar) ; wtr = hg2r*(hal*dzal)
784 
785  drdz = (wtl * drdkl + wtr * drdkr) / (dzal*wtl + dzar*wtr)
786  ! The expression for drdz above is mathematically equivalent to:
787  ! drdz = ((hg2L/haL) * drdkL/dzaL + (hg2R/haR) * drdkR/dzaR) / &
788  ! ((hg2L/haL) + (hg2R/haR))
789  hg2a = h(i,j,k-1)*h(i+1,j,k-1) + h_neglect2
790  hg2b = h(i,j,k)*h(i+1,j,k) + h_neglect2
791  haa = 0.5*(h(i,j,k-1) + h(i+1,j,k-1)) + h_neglect
792  hab = 0.5*(h(i,j,k) + h(i+1,j,k)) + h_neglect
793 
794  ! hN2_u is used with the FGNV streamfunction formulation
795  hn2_u(i,k) = (0.5 * gv%H_to_Z * ( hg2a / haa + hg2b / hab )) * &
796  max(drdz*g_rho0, n2_floor)
797  endif
798  if (present_slope_x) then
799  slope = slope_x(i,j,k)
800  slope2_ratio_u(i,k) = slope**2 * i_slope_max2
801  else
802  ! Use the harmonic mean thicknesses to weight the horizontal gradients.
803  ! These unnormalized weights have been rearranged to minimize divisions.
804  wta = hg2a*hab ; wtb = hg2b*haa
805  ! This is the gradient of density along geopotentials.
806  drdx = ((wta * drdia + wtb * drdib) / (wta + wtb) - &
807  drdz * (e(i,j,k)-e(i+1,j,k))) * g%IdxCu(i,j)
808 
809  ! This estimate of slope is accurate for small slopes, but bounded
810  ! to be between -1 and 1.
811  mag_grad2 = drdx**2 + (us%m_to_Z*drdz)**2
812  if (mag_grad2 > 0.0) then
813  slope = drdx / sqrt(mag_grad2)
814  slope2_ratio_u(i,k) = slope**2 * i_slope_max2
815  else ! Just in case mag_grad2 = 0 ever.
816  slope = 0.0
817  slope2_ratio_u(i,k) = 1.0e20 ! Force the use of the safe streamfunction.
818  endif
819  endif
820 
821  ! Adjust real slope by weights that bias towards slope of interfaces
822  ! that ignore density gradients along layers.
823  if (present_int_slope_u) then
824  slope = (1.0 - int_slope_u(i,j,k)) * slope + &
825  int_slope_u(i,j,k) * us%Z_to_m*((e(i+1,j,k)-e(i,j,k)) * g%IdxCu(i,j))
826  slope2_ratio_u(i,k) = (1.0 - int_slope_u(i,j,k)) * slope2_ratio_u(i,k)
827  endif
828 
829  slope_x_pe(i,j,k) = min(slope,cs%slope_max)
830  hn2_x_pe(i,j,k) = hn2_u(i,k) * us%m_to_Z
831  if (cs%id_slope_x > 0) cs%diagSlopeX(i,j,k) = slope
832 
833  ! Estimate the streamfunction at each interface [m3 s-1].
834  sfn_unlim_u(i,k) = -((kh_u(i,j,k)*g%dy_Cu(i,j))*us%m_to_Z*slope)
835 
836  ! Avoid moving dense water upslope from below the level of
837  ! the bottom on the receiving side.
838  if (sfn_unlim_u(i,k) > 0.0) then ! The flow below this interface is positive.
839  if (e(i,j,k) < e(i+1,j,nz+1)) then
840  sfn_unlim_u(i,k) = 0.0 ! This is not uhtot, because it may compensate for
841  ! deeper flow in very unusual cases.
842  elseif (e(i+1,j,nz+1) > e(i,j,k+1)) then
843  ! Scale the transport with the fraction of the donor layer above
844  ! the bottom on the receiving side.
845  sfn_unlim_u(i,k) = sfn_unlim_u(i,k) * ((e(i,j,k) - e(i+1,j,nz+1)) / &
846  ((e(i,j,k) - e(i,j,k+1)) + dz_neglect))
847  endif
848  else
849  if (e(i+1,j,k) < e(i,j,nz+1)) then ; sfn_unlim_u(i,k) = 0.0
850  elseif (e(i,j,nz+1) > e(i+1,j,k+1)) then
851  sfn_unlim_u(i,k) = sfn_unlim_u(i,k) * ((e(i+1,j,k) - e(i,j,nz+1)) / &
852  ((e(i+1,j,k) - e(i+1,j,k+1)) + dz_neglect))
853  endif
854  endif
855 
856  else ! .not. use_EOS
857  if (present_slope_x) then
858  slope = slope_x(i,j,k)
859  else
860  slope = us%Z_to_m*((e(i,j,k)-e(i+1,j,k))*g%IdxCu(i,j)) * g%mask2dCu(i,j)
861  endif
862  if (cs%id_slope_x > 0) cs%diagSlopeX(i,j,k) = slope
863  sfn_unlim_u(i,k) = ((kh_u(i,j,k)*g%dy_Cu(i,j))*us%m_to_Z*slope)
864  hn2_u(i,k) = us%L_to_m**2*us%s_to_T**2*gv%g_prime(k)
865  endif ! if (use_EOS)
866  else ! if (k > nk_linear)
867  hn2_u(i,k) = n2_floor * dz_neglect
868  sfn_unlim_u(i,k) = 0.
869  endif ! if (k > nk_linear)
870  if (cs%id_sfn_unlim_x>0) diag_sfn_unlim_x(i,j,k) = sfn_unlim_u(i,k)
871  enddo ! i-loop
872  enddo ! k-loop
873 
874  if (cs%use_FGNV_streamfn) then
875  do k=1,nz ; do i=is-1,ie ; if (g%mask2dCu(i,j)>0.) then
876  h_harm = max( h_neglect, &
877  2. * h(i,j,k) * h(i+1,j,k) / ( ( h(i,j,k) + h(i+1,j,k) ) + h_neglect ) )
878  c2_h_u(i,k) = cs%FGNV_scale * ( 0.5*( cg1(i,j) + cg1(i+1,j) ) )**2 / (gv%H_to_Z*h_harm)
879  endif ; enddo ; enddo
880 
881  ! Solve an elliptic equation for the streamfunction following Ferrari et al., 2010.
882  do i=is-1,ie
883  if (g%mask2dCu(i,j)>0.) then
884  sfn_unlim_u(i,:) = ( 1. + cs%FGNV_scale ) * sfn_unlim_u(i,:)
885  call streamfn_solver(nz, c2_h_u(i,:), hn2_u(i,:), sfn_unlim_u(i,:))
886  else
887  sfn_unlim_u(i,:) = 0.
888  endif
889  enddo
890  endif
891 
892  do k=nz,2,-1
893  do i=is-1,ie
894  if (k > nk_linear) then
895  if (use_eos) then
896 
897  if (uhtot(i,j) <= 0.0) then
898  ! The transport that must balance the transport below is positive.
899  sfn_safe = uhtot(i,j) * (1.0 - h_frac(i,j,k)) * gv%H_to_Z
900  else ! (uhtot(I,j) > 0.0)
901  sfn_safe = uhtot(i,j) * (1.0 - h_frac(i+1,j,k)) * gv%H_to_Z
902  endif
903 
904  ! The actual streamfunction at each interface.
905  sfn_est = (sfn_unlim_u(i,k) + slope2_ratio_u(i,k)*sfn_safe) / (1.0 + slope2_ratio_u(i,k))
906  else ! With .not.use_EOS, the layers are constant density.
907  sfn_est = sfn_unlim_u(i,k)
908  endif
909 
910  ! Make sure that there is enough mass above to allow the streamfunction
911  ! to satisfy the boundary condition of 0 at the surface.
912  sfn_in_h = min(max(sfn_est * gv%Z_to_H, -h_avail_rsum(i,j,k)), h_avail_rsum(i+1,j,k))
913 
914  ! The actual transport is limited by the mass available in the two
915  ! neighboring grid cells.
916  uhd(i,j,k) = max(min((sfn_in_h - uhtot(i,j)), h_avail(i,j,k)), &
917  -h_avail(i+1,j,k))
918 
919  if (cs%id_sfn_x>0) diag_sfn_x(i,j,k) = diag_sfn_x(i,j,k+1) + uhd(i,j,k)
920 ! sfn_x(I,j,K) = max(min(Sfn_in_h, uhtot(I,j)+h_avail(i,j,k)), &
921 ! uhtot(I,j)-h_avail(i+1,j,K))
922 ! sfn_slope_x(I,j,K) = max(uhtot(I,j)-h_avail(i+1,j,k), &
923 ! min(uhtot(I,j)+h_avail(i,j,k), &
924 ! min(h_avail_rsum(i+1,j,K), max(-h_avail_rsum(i,j,K), &
925 ! (KH_u(I,j,K)*G%dy_Cu(I,j)) * ((e(i,j,K)-e(i+1,j,K))*G%IdxCu(I,j)) )) ))
926  else ! k <= nk_linear
927  ! Balance the deeper flow with a return flow uniformly distributed
928  ! though the remaining near-surface layers. This is the same as
929  ! using Sfn_safe above. There is no need to apply the limiters in
930  ! this case.
931  if (uhtot(i,j) <= 0.0) then
932  uhd(i,j,k) = -uhtot(i,j) * h_frac(i,j,k)
933  else ! (uhtot(I,j) > 0.0)
934  uhd(i,j,k) = -uhtot(i,j) * h_frac(i+1,j,k)
935  endif
936 
937  if (cs%id_sfn_x>0) diag_sfn_x(i,j,k) = diag_sfn_x(i,j,k+1) + uhd(i,j,k)
938 ! if (sfn_slope_x(I,j,K+1) <= 0.0) then
939 ! sfn_slope_x(I,j,K) = sfn_slope_x(I,j,K+1) * (1.0 - h_frac(i,j,k))
940 ! else
941 ! sfn_slope_x(I,j,K) = sfn_slope_x(I,j,K+1) * (1.0 - h_frac(i+1,j,k))
942 ! endif
943  endif
944 
945  uhtot(i,j) = uhtot(i,j) + uhd(i,j,k)
946 
947  if (find_work) then
948  ! This is the energy tendency based on the original profiles, and does
949  ! not include any nonlinear terms due to a finite time step (which would
950  ! involve interactions between the fluxes through the different faces.
951  ! A second order centered estimate is used for the density transfered
952  ! between water columns.
953 
954  work_u(i,j) = work_u(i,j) + g_scale * &
955  ( uhtot(i,j) * drdkde_u(i,k) - &
956  (uhd(i,j,k) * drdi_u(i,k)) * 0.25 * &
957  ((e(i,j,k) + e(i,j,k+1)) + (e(i+1,j,k) + e(i+1,j,k+1))) )
958  endif
959 
960  enddo
961  enddo ! end of k-loop
962  enddo ! end of j-loop
963 
964  ! Calculate the meridional fluxes and gradients.
965 !$OMP parallel do default(none) shared(nz,is,ie,js,je,find_work,use_EOS,G,GV,US,pres,T,S, &
966 !$OMP nk_linear,IsdB,tv,h,h_neglect,e,dz_neglect, &
967 !$OMP I_slope_max2,h_neglect2,present_int_slope_v, &
968 !$OMP int_slope_v,KH_v,vhtot,h_frac,h_avail_rsum, &
969 !$OMP vhD,h_avail,G_scale,Work_v,CS,slope_y,cg1, &
970 !$OMP diag_sfn_y, diag_sfn_unlim_y,N2_floor, &
971 !$OMP present_slope_y,G_rho0) &
972 !$OMP private(drdjA,drdjB,drdkL,drdkR,pres_v,T_v,S_v, &
973 !$OMP drho_dT_v,drho_dS_v,hg2A,hg2B,hg2L,hg2R,haA, &
974 !$OMP haB,haL,haR,dzaL,dzaR,wtA,wtB,wtL,wtR,drdz, &
975 !$OMP drdy,mag_grad2,Slope,slope2_Ratio_v,hN2_v, &
976 !$OMP Sfn_unlim_v,drdj_v,drdkDe_v,h_harm,c2_h_v, &
977 !$OMP Sfn_safe,Sfn_est,Sfn_in_h,calc_derivatives)
978  do j=js-1,je
979  do k=nz,2,-1
980  if (find_work .and. .not.(use_eos)) then
981  drdja = 0.0 ; drdjb = 0.0
982  drdkl = gv%Rlay(k)-gv%Rlay(k-1) ; drdkr = drdkl
983  endif
984 
985  calc_derivatives = use_eos .and. (k >= nk_linear) .and. &
986  (find_work .or. .not. present_slope_y)
987 
988  if (calc_derivatives) then
989  do i=is,ie
990  pres_v(i) = 0.5*(pres(i,j,k) + pres(i,j+1,k))
991  t_v(i) = 0.25*((t(i,j,k) + t(i,j+1,k)) + (t(i,j,k-1) + t(i,j+1,k-1)))
992  s_v(i) = 0.25*((s(i,j,k) + s(i,j+1,k)) + (s(i,j,k-1) + s(i,j+1,k-1)))
993  enddo
994  call calculate_density_derivs(t_v, s_v, pres_v, drho_dt_v, &
995  drho_ds_v, is, ie-is+1, tv%eqn_of_state)
996  endif
997  do i=is,ie
998  if (calc_derivatives) then
999  ! Estimate the horizontal density gradients along layers.
1000  drdja = drho_dt_v(i) * (t(i,j+1,k-1)-t(i,j,k-1)) + &
1001  drho_ds_v(i) * (s(i,j+1,k-1)-s(i,j,k-1))
1002  drdjb = drho_dt_v(i) * (t(i,j+1,k)-t(i,j,k)) + &
1003  drho_ds_v(i) * (s(i,j+1,k)-s(i,j,k))
1004 
1005  ! Estimate the vertical density gradients times the grid spacing.
1006  drdkl = (drho_dt_v(i) * (t(i,j,k)-t(i,j,k-1)) + &
1007  drho_ds_v(i) * (s(i,j,k)-s(i,j,k-1)))
1008  drdkr = (drho_dt_v(i) * (t(i,j+1,k)-t(i,j+1,k-1)) + &
1009  drho_ds_v(i) * (s(i,j+1,k)-s(i,j+1,k-1)))
1010  drdkde_v(i,k) = drdkr * e(i,j+1,k) - drdkl * e(i,j,k)
1011  elseif (find_work) then ! This is used in pure stacked SW mode
1012  drdkde_v(i,k) = drdkr * e(i,j+1,k) - drdkl * e(i,j,k)
1013  endif
1014 
1015  if (find_work) drdj_v(i,k) = drdjb
1016 
1017  if (k > nk_linear) then
1018  if (use_eos) then
1019  if (cs%use_FGNV_streamfn .or. find_work .or. .not. present_slope_y) then
1020  hg2l = h(i,j,k-1)*h(i,j,k) + h_neglect2
1021  hg2r = h(i,j+1,k-1)*h(i,j+1,k) + h_neglect2
1022  hal = 0.5*(h(i,j,k-1) + h(i,j,k)) + h_neglect
1023  har = 0.5*(h(i,j+1,k-1) + h(i,j+1,k)) + h_neglect
1024  if (gv%Boussinesq) then
1025  dzal = hal * gv%H_to_Z ; dzar = har * gv%H_to_Z
1026  else
1027  dzal = 0.5*(e(i,j,k-1) - e(i,j,k+1)) + dz_neglect
1028  dzar = 0.5*(e(i,j+1,k-1) - e(i,j+1,k+1)) + dz_neglect
1029  endif
1030  ! Use the harmonic mean thicknesses to weight the horizontal gradients.
1031  ! These unnormalized weights have been rearranged to minimize divisions.
1032  wtl = hg2l*(har*dzar) ; wtr = hg2r*(hal*dzal)
1033 
1034  drdz = (wtl * drdkl + wtr * drdkr) / (dzal*wtl + dzar*wtr)
1035  ! The expression for drdz above is mathematically equivalent to:
1036  ! drdz = ((hg2L/haL) * drdkL/dzaL + (hg2R/haR) * drdkR/dzaR) / &
1037  ! ((hg2L/haL) + (hg2R/haR))
1038  hg2a = h(i,j,k-1)*h(i,j+1,k-1) + h_neglect2
1039  hg2b = h(i,j,k)*h(i,j+1,k) + h_neglect2
1040  haa = 0.5*(h(i,j,k-1) + h(i,j+1,k-1)) + h_neglect
1041  hab = 0.5*(h(i,j,k) + h(i,j+1,k)) + h_neglect
1042 
1043  ! hN2_v is used with the FGNV streamfunction formulation
1044  hn2_v(i,k) = (0.5 * gv%H_to_Z * ( hg2a / haa + hg2b / hab )) * &
1045  max(drdz*g_rho0, n2_floor)
1046  endif
1047  if (present_slope_y) then
1048  slope = slope_y(i,j,k)
1049  slope2_ratio_v(i,k) = slope**2 * i_slope_max2
1050  else
1051  ! Use the harmonic mean thicknesses to weight the horizontal gradients.
1052  ! These unnormalized weights have been rearranged to minimize divisions.
1053  wta = hg2a*hab ; wtb = hg2b*haa
1054  ! This is the gradient of density along geopotentials.
1055  drdy = ((wta * drdja + wtb * drdjb) / (wta + wtb) - &
1056  drdz * (e(i,j,k)-e(i,j+1,k))) * g%IdyCv(i,j)
1057 
1058  ! This estimate of slope is accurate for small slopes, but bounded
1059  ! to be between -1 and 1.
1060  mag_grad2 = drdy**2 + (us%m_to_Z*drdz)**2
1061  if (mag_grad2 > 0.0) then
1062  slope = drdy / sqrt(mag_grad2)
1063  slope2_ratio_v(i,k) = slope**2 * i_slope_max2
1064  else ! Just in case mag_grad2 = 0 ever.
1065  slope = 0.0
1066  slope2_ratio_v(i,k) = 1.0e20 ! Force the use of the safe streamfunction.
1067  endif
1068  endif
1069 
1070  ! Adjust real slope by weights that bias towards slope of interfaces
1071  ! that ignore density gradients along layers.
1072  if (present_int_slope_v) then
1073  slope = (1.0 - int_slope_v(i,j,k)) * slope + &
1074  int_slope_v(i,j,k) * us%Z_to_m*((e(i,j+1,k)-e(i,j,k)) * g%IdyCv(i,j))
1075  slope2_ratio_v(i,k) = (1.0 - int_slope_v(i,j,k)) * slope2_ratio_v(i,k)
1076  endif
1077 
1078  slope_y_pe(i,j,k) = min(slope,cs%slope_max)
1079  hn2_y_pe(i,j,k) = hn2_v(i,k) * us%m_to_Z
1080  if (cs%id_slope_y > 0) cs%diagSlopeY(i,j,k) = slope
1081 
1082  ! Estimate the streamfunction at each interface [m3 s-1].
1083  sfn_unlim_v(i,k) = -((kh_v(i,j,k)*g%dx_Cv(i,j))*us%m_to_Z*slope)
1084 
1085  ! Avoid moving dense water upslope from below the level of
1086  ! the bottom on the receiving side.
1087  if (sfn_unlim_v(i,k) > 0.0) then ! The flow below this interface is positive.
1088  if (e(i,j,k) < e(i,j+1,nz+1)) then
1089  sfn_unlim_v(i,k) = 0.0 ! This is not vhtot, because it may compensate for
1090  ! deeper flow in very unusual cases.
1091  elseif (e(i,j+1,nz+1) > e(i,j,k+1)) then
1092  ! Scale the transport with the fraction of the donor layer above
1093  ! the bottom on the receiving side.
1094  sfn_unlim_v(i,k) = sfn_unlim_v(i,k) * ((e(i,j,k) - e(i,j+1,nz+1)) / &
1095  ((e(i,j,k) - e(i,j,k+1)) + dz_neglect))
1096  endif
1097  else
1098  if (e(i,j+1,k) < e(i,j,nz+1)) then ; sfn_unlim_v(i,k) = 0.0
1099  elseif (e(i,j,nz+1) > e(i,j+1,k+1)) then
1100  sfn_unlim_v(i,k) = sfn_unlim_v(i,k) * ((e(i,j+1,k) - e(i,j,nz+1)) / &
1101  ((e(i,j+1,k) - e(i,j+1,k+1)) + dz_neglect))
1102  endif
1103  endif
1104 
1105  else ! .not. use_EOS
1106  if (present_slope_y) then
1107  slope = slope_y(i,j,k)
1108  else
1109  slope = us%Z_to_m*((e(i,j,k)-e(i,j+1,k))*g%IdyCv(i,j)) * g%mask2dCv(i,j)
1110  endif
1111  if (cs%id_slope_y > 0) cs%diagSlopeY(i,j,k) = slope
1112  sfn_unlim_v(i,k) = ((kh_v(i,j,k)*g%dx_Cv(i,j))*us%m_to_Z*slope)
1113  hn2_v(i,k) = us%L_to_m**2*us%s_to_T**2*gv%g_prime(k)
1114  endif ! if (use_EOS)
1115  else ! if (k > nk_linear)
1116  hn2_v(i,k) = n2_floor * dz_neglect
1117  sfn_unlim_v(i,k) = 0.
1118  endif ! if (k > nk_linear)
1119  if (cs%id_sfn_unlim_y>0) diag_sfn_unlim_y(i,j,k) = sfn_unlim_v(i,k)
1120  enddo ! i-loop
1121  enddo ! k-loop
1122 
1123  if (cs%use_FGNV_streamfn) then
1124  do k=1,nz ; do i=is,ie ; if (g%mask2dCv(i,j)>0.) then
1125  h_harm = max( h_neglect, &
1126  2. * h(i,j,k) * h(i,j+1,k) / ( ( h(i,j,k) + h(i,j+1,k) ) + h_neglect ) )
1127  c2_h_v(i,k) = cs%FGNV_scale * ( 0.5*( cg1(i,j) + cg1(i,j+1) ) )**2 / (gv%H_to_Z*h_harm)
1128  endif ; enddo ; enddo
1129 
1130  ! Solve an elliptic equation for the streamfunction following Ferrari et al., 2010.
1131  do i=is,ie
1132  if (g%mask2dCv(i,j)>0.) then
1133  sfn_unlim_v(i,:) = ( 1. + cs%FGNV_scale ) * sfn_unlim_v(i,:)
1134  call streamfn_solver(nz, c2_h_v(i,:), hn2_v(i,:), sfn_unlim_v(i,:))
1135  else
1136  sfn_unlim_v(i,:) = 0.
1137  endif
1138  enddo
1139  endif
1140 
1141  do k=nz,2,-1
1142  do i=is,ie
1143  if (k > nk_linear) then
1144  if (use_eos) then
1145 
1146  if (vhtot(i,j) <= 0.0) then
1147  ! The transport that must balance the transport below is positive.
1148  sfn_safe = vhtot(i,j) * (1.0 - h_frac(i,j,k)) * gv%H_to_Z
1149  else ! (vhtot(I,j) > 0.0)
1150  sfn_safe = vhtot(i,j) * (1.0 - h_frac(i,j+1,k)) * gv%H_to_Z
1151  endif
1152 
1153  ! The actual streamfunction at each interface.
1154  sfn_est = (sfn_unlim_v(i,k) + slope2_ratio_v(i,k)*sfn_safe) / (1.0 + slope2_ratio_v(i,k))
1155  else ! With .not.use_EOS, the layers are constant density.
1156  sfn_est = sfn_unlim_v(i,k)
1157  endif
1158 
1159  ! Make sure that there is enough mass above to allow the streamfunction
1160  ! to satisfy the boundary condition of 0 at the surface.
1161  sfn_in_h = min(max(sfn_est * gv%Z_to_H, -h_avail_rsum(i,j,k)), h_avail_rsum(i,j+1,k))
1162 
1163  ! The actual transport is limited by the mass available in the two
1164  ! neighboring grid cells.
1165  vhd(i,j,k) = max(min((sfn_in_h - vhtot(i,j)), h_avail(i,j,k)), &
1166  -h_avail(i,j+1,k))
1167 
1168  if (cs%id_sfn_y>0) diag_sfn_y(i,j,k) = diag_sfn_y(i,j,k+1) + vhd(i,j,k)
1169 ! sfn_y(i,J,K) = max(min(Sfn_in_h, vhtot(i,J)+h_avail(i,j,k)), &
1170 ! vhtot(i,J)-h_avail(i,j+1,k))
1171 ! sfn_slope_y(i,J,K) = max(vhtot(i,J)-h_avail(i,j+1,k), &
1172 ! min(vhtot(i,J)+h_avail(i,j,k), &
1173 ! min(h_avail_rsum(i,j+1,K), max(-h_avail_rsum(i,j,K), &
1174 ! (KH_v(i,J,K)*G%dx_Cv(i,J)) * ((e(i,j,K)-e(i,j+1,K))*G%IdyCv(i,J)) )) ))
1175  else ! k <= nk_linear
1176  ! Balance the deeper flow with a return flow uniformly distributed
1177  ! though the remaining near-surface layers. This is the same as
1178  ! using Sfn_safe above. There is no need to apply the limiters in
1179  ! this case.
1180  if (vhtot(i,j) <= 0.0) then
1181  vhd(i,j,k) = -vhtot(i,j) * h_frac(i,j,k)
1182  else ! (vhtot(i,J) > 0.0)
1183  vhd(i,j,k) = -vhtot(i,j) * h_frac(i,j+1,k)
1184  endif
1185 
1186  if (cs%id_sfn_y>0) diag_sfn_y(i,j,k) = diag_sfn_y(i,j,k+1) + vhd(i,j,k)
1187 ! if (sfn_slope_y(i,J,K+1) <= 0.0) then
1188 ! sfn_slope_y(i,J,K) = sfn_slope_y(i,J,K+1) * (1.0 - h_frac(i,j,k))
1189 ! else
1190 ! sfn_slope_y(i,J,K) = sfn_slope_y(i,J,K+1) * (1.0 - h_frac(i,j+1,k))
1191 ! endif
1192  endif
1193 
1194  vhtot(i,j) = vhtot(i,j) + vhd(i,j,k)
1195 
1196  if (find_work) then
1197  ! This is the energy tendency based on the original profiles, and does
1198  ! not include any nonlinear terms due to a finite time step (which would
1199  ! involve interactions between the fluxes through the different faces.
1200  ! A second order centered estimate is used for the density transfered
1201  ! between water columns.
1202 
1203  work_v(i,j) = work_v(i,j) + g_scale * &
1204  ( vhtot(i,j) * drdkde_v(i,k) - &
1205  (vhd(i,j,k) * drdj_v(i,k)) * 0.25 * &
1206  ((e(i,j,k) + e(i,j,k+1)) + (e(i,j+1,k) + e(i,j+1,k+1))) )
1207  endif
1208 
1209  enddo
1210  enddo ! end of k-loop
1211  enddo ! end of j-loop
1212 
1213  ! In layer 1, enforce the boundary conditions that Sfn(z=0) = 0.0
1214  if (.not.find_work .or. .not.(use_eos)) then
1215  do j=js,je ; do i=is-1,ie ; uhd(i,j,1) = -uhtot(i,j) ; enddo ; enddo
1216  do j=js-1,je ; do i=is,ie ; vhd(i,j,1) = -vhtot(i,j) ; enddo ; enddo
1217  else
1218  !$OMP parallel do default(shared) private(pres_u,T_u,S_u,drho_dT_u,drho_dS_u,drdiB)
1219  do j=js,je
1220  if (use_eos) then
1221  do i=is-1,ie
1222  pres_u(i) = 0.5*(pres(i,j,1) + pres(i+1,j,1))
1223  t_u(i) = 0.5*(t(i,j,1) + t(i+1,j,1))
1224  s_u(i) = 0.5*(s(i,j,1) + s(i+1,j,1))
1225  enddo
1226  call calculate_density_derivs(t_u, s_u, pres_u, drho_dt_u, &
1227  drho_ds_u, (is-isdb+1)-1, ie-is+2, tv%eqn_of_state)
1228  endif
1229  do i=is-1,ie
1230  uhd(i,j,1) = -uhtot(i,j)
1231 
1232  if (use_eos) then
1233  drdib = drho_dt_u(i) * (t(i+1,j,1)-t(i,j,1)) + &
1234  drho_ds_u(i) * (s(i+1,j,1)-s(i,j,1))
1235  endif
1236  work_u(i,j) = work_u(i,j) + g_scale * &
1237  ( (uhd(i,j,1) * drdib) * 0.25 * &
1238  ((e(i,j,1) + e(i,j,2)) + (e(i+1,j,1) + e(i+1,j,2))) )
1239 
1240  enddo
1241  enddo
1242 
1243  !$OMP parallel do default(shared) private(pres_v,T_v,S_v,drho_dT_v,drho_dS_v,drdjB)
1244  do j=js-1,je
1245  if (use_eos) then
1246  do i=is,ie
1247  pres_v(i) = 0.5*(pres(i,j,1) + pres(i,j+1,1))
1248  t_v(i) = 0.5*(t(i,j,1) + t(i,j+1,1))
1249  s_v(i) = 0.5*(s(i,j,1) + s(i,j+1,1))
1250  enddo
1251  call calculate_density_derivs(t_v, s_v, pres_v, drho_dt_v, &
1252  drho_ds_v, is, ie-is+1, tv%eqn_of_state)
1253  endif
1254  do i=is,ie
1255  vhd(i,j,1) = -vhtot(i,j)
1256 
1257  if (use_eos) then
1258  drdjb = drho_dt_v(i) * (t(i,j+1,1)-t(i,j,1)) + &
1259  drho_ds_v(i) * (s(i,j+1,1)-s(i,j,1))
1260  endif
1261  work_v(i,j) = work_v(i,j) - g_scale * &
1262  ( (vhd(i,j,1) * drdjb) * 0.25 * &
1263  ((e(i,j,1) + e(i,j,2)) + (e(i,j+1,1) + e(i,j+1,2))) )
1264  enddo
1265  enddo
1266  endif
1267 
1268 
1269  !if (find_work) then ; do j=js,je ; do i=is,ie ; do k=nz,1,-1
1270  if (find_work) then ; do j=js,je ; do i=is,ie
1271  ! Note that the units of Work_v and Work_u are W, while Work_h is W m-2.
1272  work_h = 0.5 * g%IareaT(i,j) * &
1273  ((work_u(i-1,j) + work_u(i,j)) + (work_v(i,j-1) + work_v(i,j)))
1274  pe_release_h = -0.25*(kh_u(i,j,k)*(slope_x_pe(i,j,k)**2) * hn2_x_pe(i,j,k) + &
1275  kh_u(i-1,j,k)*(slope_x_pe(i-1,j,k)**2) * hn2_x_pe(i-1,j,k) + &
1276  kh_v(i,j,k)*(slope_y_pe(i,j,k)**2) * hn2_y_pe(i,j,k) + &
1277  kh_v(i,j-1,k)*(slope_y_pe(i,j-1,k)**2) * hn2_y_pe(i,j-1,k))
1278  if (associated(cs%GMwork)) cs%GMwork(i,j) = work_h
1279  if (associated(meke)) then ; if (associated(meke%GM_src)) then
1280  if (cs%GM_src_alt) then
1281  meke%GM_src(i,j) = meke%GM_src(i,j) + pe_release_h
1282  else
1283  meke%GM_src(i,j) = meke%GM_src(i,j) + work_h
1284  endif
1285  endif ; endif
1286  !enddo ; enddo ; enddo ; endif
1287  enddo ; enddo ; endif
1288 
1289  if (cs%id_slope_x > 0) call post_data(cs%id_slope_x, cs%diagSlopeX, cs%diag)
1290  if (cs%id_slope_y > 0) call post_data(cs%id_slope_y, cs%diagSlopeY, cs%diag)
1291  if (cs%id_sfn_x > 0) call post_data(cs%id_sfn_x, diag_sfn_x, cs%diag)
1292  if (cs%id_sfn_y > 0) call post_data(cs%id_sfn_y, diag_sfn_y, cs%diag)
1293  if (cs%id_sfn_unlim_x > 0) call post_data(cs%id_sfn_unlim_x, diag_sfn_unlim_x, cs%diag)
1294  if (cs%id_sfn_unlim_y > 0) call post_data(cs%id_sfn_unlim_y, diag_sfn_unlim_y, cs%diag)
1295 
1296 end subroutine thickness_diffuse_full
1297 
1298 !> Tridiagonal solver for streamfunction at interfaces
1299 subroutine streamfn_solver(nk, c2_h, hN2, sfn)
1300  integer, intent(in) :: nk !< Number of layers
1301  real, dimension(nk), intent(in) :: c2_h !< Wave speed squared over thickness in layers [m s-2]
1302  real, dimension(nk+1), intent(in) :: hN2 !< Thickness times N2 at interfaces [m s-2]
1303  real, dimension(nk+1), intent(inout) :: sfn !< Streamfunction [Z m2 s-1 ~> m3 s-1] or arbitrary units
1304  !! On entry, equals diffusivity times slope.
1305  !! On exit, equals the streamfunction.
1306  ! Local variables
1307  integer :: k
1308 
1309  real :: b_denom, beta, d1, c1(nk)
1310 
1311  sfn(1) = 0.
1312  b_denom = hn2(2) + c2_h(1)
1313  beta = 1.0 / ( b_denom + c2_h(2) )
1314  d1 = beta * b_denom
1315  sfn(2) = ( beta * hn2(2) )*sfn(2)
1316  do k=3,nk
1317  c1(k-1) = beta * c2_h(k-1)
1318  b_denom = hn2(k) + d1*c2_h(k-1)
1319  beta = 1.0 / (b_denom + c2_h(k))
1320  d1 = beta * b_denom
1321  sfn(k) = beta * (hn2(k)*sfn(k) + c2_h(k-1)*sfn(k-1))
1322  enddo
1323  c1(nk) = beta * c2_h(nk)
1324  sfn(nk+1) = 0.
1325  do k=nk,2,-1
1326  sfn(k) = sfn(k) + c1(k)*sfn(k+1)
1327  enddo
1328 
1329 end subroutine streamfn_solver
1330 
1331 !> Modifies thickness diffusivities to untangle layer structures
1332 subroutine add_detangling_kh(h, e, Kh_u, Kh_v, KH_u_CFL, KH_v_CFL, tv, dt, G, GV, US, CS, &
1333  int_slope_u, int_slope_v)
1334  type(ocean_grid_type), intent(in) :: G !< Ocean grid structure
1335  type(verticalgrid_type), intent(in) :: GV !< Vertical grid structure
1336  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
1337  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thickness [H ~> m or kg m-2]
1338  real, dimension(SZI_(G),SZJ_(G),SZK_(G)+1), intent(in) :: e !< Interface positions [Z ~> m]
1339  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)+1), intent(inout) :: Kh_u !< Thickness diffusivity on interfaces
1340  !! at u points [m2 s-1]
1341  real, dimension(SZI_(G),SZJB_(G),SZK_(G)+1), intent(inout) :: Kh_v !< Thickness diffusivity on interfaces
1342  !! at v points [m2 s-1]
1343  real, dimension(SZIB_(G),SZJ_(G)), intent(in) :: Kh_u_CFL !< Maximum stable thickness diffusivity
1344  !! at u points [m2 s-1]
1345  real, dimension(SZI_(G),SZJB_(G)), intent(in) :: Kh_v_CFL !< Maximum stable thickness diffusivity
1346  !! at v points [m2 s-1]
1347  type(thermo_var_ptrs), intent(in) :: tv !< Thermodynamics structure
1348  real, intent(in) :: dt !< Time increment [s]
1349  type(thickness_diffuse_cs), pointer :: CS !< Control structure for thickness diffusion
1350  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)+1), intent(inout) :: int_slope_u !< Ratio that determine how much of
1351  !! the isopycnal slopes are taken directly from
1352  !! the interface slopes without consideration
1353  !! of density gradients.
1354  real, dimension(SZI_(G),SZJB_(G),SZK_(G)+1), intent(inout) :: int_slope_v !< Ratio that determine how much of
1355  !! the isopycnal slopes are taken directly from
1356  !! the interface slopes without consideration
1357  !! of density gradients.
1358  ! Local variables
1359  real, dimension(SZI_(G),SZJ_(G),SZK_(G)) :: &
1360  de_top ! The distances between the top of a layer and the top of the
1361  ! region where the detangling is applied [H ~> m or kg m-2].
1362  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)) :: &
1363  Kh_lay_u ! The tentative interface height diffusivity for each layer at
1364  ! u points [m2 s-1].
1365  real, dimension(SZI_(G),SZJB_(G),SZK_(G)) :: &
1366  Kh_lay_v ! The tentative interface height diffusivity for each layer at
1367  ! v points [m2 s-1].
1368  real, dimension(SZI_(G),SZJ_(G)) :: &
1369  de_bot ! The distances from the bottom of the region where the
1370  ! detangling is applied [H ~> m or kg m-2].
1371  real :: h1, h2 ! The thinner and thicker surrounding thicknesses [H ~> m or kg m-2],
1372  ! with the thinner modified near the boundaries to mask out
1373  ! thickness variations due to topography, etc.
1374  real :: jag_Rat ! The nondimensional jaggedness ratio for a layer, going
1375  ! from 0 (smooth) to 1 (jagged). This is the difference
1376  ! between the arithmetic and harmonic mean thicknesses
1377  ! normalized by the arithmetic mean thickness.
1378  real :: Kh_scale ! A ratio by which Kh_u_CFL is scaled for maximally jagged
1379  ! layers [nondim].
1380  real :: Kh_det ! The detangling diffusivity [m2 s-1].
1381  real :: h_neglect ! A thickness that is so small it is usually lost
1382  ! in roundoff and can be neglected [H ~> m or kg m-2].
1383 
1384  real :: I_sl ! The absolute value of the larger in magnitude of the slopes
1385  ! above and below.
1386  real :: Rsl ! The ratio of the smaller magnitude slope to the larger
1387  ! magnitude one [nondim]. 0 <= Rsl <1.
1388  real :: IRsl ! The (limited) inverse of Rsl [nondim]. 1 < IRsl <= 1e9.
1389  real :: dH ! The thickness gradient divided by the damping timescale
1390  ! and the ratio of the face length to the adjacent cell
1391  ! areas for comparability with the diffusivities [m2 s-1].
1392  real :: adH ! The absolute value of dH [m2 s-1].
1393  real :: sign ! 1 or -1, with the same sign as the layer thickness gradient.
1394  real :: sl_K ! The sign-corrected slope of the interface above [nondim].
1395  real :: sl_Kp1 ! The sign-corrected slope of the interface below [nondim].
1396  real :: I_sl_K ! The (limited) inverse of sl_K [nondim].
1397  real :: I_sl_Kp1 ! The (limited) inverse of sl_Kp1 [nondim].
1398  real :: I_4t ! A quarter of a unit conversion factor divided by
1399  ! the damping timescale [s-1].
1400  real :: Fn_R ! A function of Rsl, such that Rsl < Fn_R < 1.
1401  real :: denom, I_denom ! A denominator and its inverse, various units.
1402  real :: Kh_min ! A local floor on the diffusivity [m2 s-1].
1403  real :: Kh_max ! A local ceiling on the diffusivity [m2 s-1].
1404  real :: wt1, wt2 ! Nondimensional weights.
1405  ! Variables used only in testing code.
1406  ! real, dimension(SZK_(G)) :: uh_here
1407  ! real, dimension(SZK_(G)+1) :: Sfn
1408  real :: dKh ! An increment in the diffusivity [m2 s-1].
1409 
1410  real, dimension(SZIB_(G),SZK_(G)+1) :: &
1411  Kh_bg, & ! The background (floor) value of Kh [m2 s-1].
1412  Kh, & ! The tentative value of Kh [m2 s-1].
1413  Kh_detangle, & ! The detangling diffusivity that could be used [m2 s-1].
1414  Kh_min_max_p, & ! The smallest ceiling that can be placed on Kh(I,K)
1415  ! based on the value of Kh(I,K+1) [m2 s-1].
1416  kh_min_max_m, & ! The smallest ceiling that can be placed on Kh(I,K)
1417  ! based on the value of Kh(I,K-1) [m2 s-1].
1418  ! The following are variables that define the relationships between
1419  ! successive values of Kh.
1420  ! Search for Kh that satisfy...
1421  ! Kh(I,K) >= Kh_min_m(I,K)*Kh(I,K-1) + Kh0_min_m(I,K)
1422  ! Kh(I,K) >= Kh_min_p(I,K)*Kh(I,K+1) + Kh0_min_p(I,K)
1423  ! Kh(I,K) <= Kh_max_m(I,K)*Kh(I,K-1) + Kh0_max_m(I,K)
1424  ! Kh(I,K) <= Kh_max_p(I,K)*Kh(I,K+1) + Kh0_max_p(I,K)
1425  kh_min_m , & ! See above [nondim].
1426  kh0_min_m , & ! See above [m2 s-1].
1427  kh_max_m , & ! See above [nondim].
1428  kh0_max_m, & ! See above [m2 s-1].
1429  kh_min_p , & ! See above [nondim].
1430  kh0_min_p , & ! See above [m2 s-1].
1431  kh_max_p , & ! See above [nondim].
1432  kh0_max_p ! See above [m2 s-1].
1433  real, dimension(SZIB_(G)) :: &
1434  Kh_max_max ! The maximum diffusivity permitted in a column.
1435  logical, dimension(SZIB_(G)) :: &
1436  do_i ! If true, work on a column.
1437  integer :: i, j, k, n, ish, jsh, is, ie, js, je, nz, k_top
1438  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
1439 
1440  k_top = gv%nk_rho_varies + 1
1441  h_neglect = gv%H_subroundoff
1442  ! The 0.5 is because we are not using uniform weightings, but are
1443  ! distributing the diffusivities more effectively (with wt1 & wt2), but this
1444  ! means that the additions to a single interface can be up to twice as large.
1445  kh_scale = 0.5
1446  if (cs%detangle_time > dt) kh_scale = 0.5 * dt / cs%detangle_time
1447 
1448  do j=js-1,je+1 ; do i=is-1,ie+1
1449  de_top(i,j,k_top) = 0.0 ; de_bot(i,j) = 0.0
1450  enddo ; enddo
1451  do k=k_top+1,nz ; do j=js-1,je+1 ; do i=is-1,ie+1
1452  de_top(i,j,k) = de_top(i,j,k-1) + h(i,j,k-1)
1453  enddo ; enddo ; enddo
1454 
1455  do j=js,je ; do i=is-1,ie
1456  kh_lay_u(i,j,nz) = 0.0 ; kh_lay_u(i,j,k_top) = 0.0
1457  enddo ; enddo
1458  do j=js-1,je ; do i=is,ie
1459  kh_lay_v(i,j,nz) = 0.0 ; kh_lay_v(i,j,k_top) = 0.0
1460  enddo ; enddo
1461 
1462  do k=nz-1,k_top+1,-1
1463  ! Find the diffusivities associated with each layer.
1464  do j=js-1,je+1 ; do i=is-1,ie+1
1465  de_bot(i,j) = de_bot(i,j) + h(i,j,k+1)
1466  enddo ; enddo
1467 
1468  do j=js,je ; do i=is-1,ie ; if (g%mask2dCu(i,j) > 0.0) then
1469  if (h(i,j,k) > h(i+1,j,k)) then
1470  h2 = h(i,j,k)
1471  h1 = max( h(i+1,j,k), h2 - min(de_bot(i+1,j), de_top(i+1,j,k)) )
1472  else
1473  h2 = h(i+1,j,k)
1474  h1 = max( h(i,j,k), h2 - min(de_bot(i,j), de_top(i,j,k)) )
1475  endif
1476  jag_rat = (h2 - h1)**2 / (h2 + h1 + h_neglect)**2
1477  kh_lay_u(i,j,k) = (kh_scale * kh_u_cfl(i,j)) * jag_rat**2
1478  endif ; enddo ; enddo
1479 
1480  do j=js-1,je ; do i=is,ie ; if (g%mask2dCv(i,j) > 0.0) then
1481  if (h(i,j,k) > h(i,j+1,k)) then
1482  h2 = h(i,j,k)
1483  h1 = max( h(i,j+1,k), h2 - min(de_bot(i,j+1), de_top(i,j+1,k)) )
1484  else
1485  h2 = h(i,j+1,k)
1486  h1 = max( h(i,j,k), h2 - min(de_bot(i,j), de_top(i,j,k)) )
1487  endif
1488  jag_rat = (h2 - h1)**2 / (h2 + h1 + h_neglect)**2
1489  kh_lay_v(i,j,k) = (kh_scale * kh_v_cfl(i,j)) * jag_rat**2
1490  endif ; enddo ; enddo
1491  enddo
1492 
1493  ! Limit the diffusivities
1494 
1495  i_4t = us%Z_to_m*kh_scale / (4.0*dt)
1496 
1497  do n=1,2
1498  if (n==1) then ; jsh = js ; ish = is-1
1499  else ; jsh = js-1 ; ish = is ; endif
1500 
1501  do j=jsh,je
1502 
1503  ! First, populate the diffusivities
1504  if (n==1) then ! This is a u-column.
1505  do i=ish,ie
1506  do_i(i) = (g%mask2dCu(i,j) > 0.0)
1507  kh_max_max(i) = kh_u_cfl(i,j)
1508  enddo
1509  do k=1,nz+1 ; do i=ish,ie
1510  kh_bg(i,k) = kh_u(i,j,k) ; kh(i,k) = kh_bg(i,k)
1511  kh_min_max_p(i,k) = kh_bg(i,k) ; kh_min_max_m(i,k) = kh_bg(i,k)
1512  kh_detangle(i,k) = 0.0
1513  enddo ; enddo
1514  else ! This is a v-column.
1515  do i=ish,ie
1516  do_i(i) = (g%mask2dCv(i,j) > 0.0) ; kh_max_max(i) = kh_v_cfl(i,j)
1517  enddo
1518  do k=1,nz+1 ; do i=ish,ie
1519  kh_bg(i,k) = kh_v(i,j,k) ; kh(i,k) = kh_bg(i,k)
1520  kh_min_max_p(i,k) = kh_bg(i,k) ; kh_min_max_m(i,k) = kh_bg(i,k)
1521  kh_detangle(i,k) = 0.0
1522  enddo ; enddo
1523  endif
1524 
1525  ! Determine the limits on the diffusivities.
1526  do k=k_top,nz ; do i=ish,ie ; if (do_i(i)) then
1527  if (n==1) then ! This is a u-column.
1528  dh = 0.0
1529  denom = ((g%IareaT(i+1,j) + g%IareaT(i,j))*g%dy_Cu(i,j))
1530  ! This expression uses differences in e in place of h for better
1531  ! consistency with the slopes.
1532  if (denom > 0.0) &
1533  dh = i_4t * ((e(i+1,j,k) - e(i+1,j,k+1)) - &
1534  (e(i,j,k) - e(i,j,k+1))) / denom
1535  ! dH = I_4t * (h(i+1,j,k) - h(i,j,k)) / denom
1536 
1537  adh = abs(dh)
1538  sign = 1.0*us%Z_to_m ; if (dh < 0) sign = -1.0*us%Z_to_m
1539  sl_k = sign * (e(i+1,j,k)-e(i,j,k)) * g%IdxCu(i,j)
1540  sl_kp1 = sign * (e(i+1,j,k+1)-e(i,j,k+1)) * g%IdxCu(i,j)
1541 
1542  ! Add the incremental diffusivites to the surrounding interfaces.
1543  ! Adding more to the more steeply sloping layers (as below) makes
1544  ! the diffusivities more than twice as effective.
1545  denom = (sl_k**2 + sl_kp1**2)
1546  wt1 = 0.5 ; wt2 = 0.5
1547  if (denom > 0.0) then
1548  wt1 = sl_k**2 / denom ; wt2 = sl_kp1**2 / denom
1549  endif
1550  kh_detangle(i,k) = kh_detangle(i,k) + wt1*kh_lay_u(i,j,k)
1551  kh_detangle(i,k+1) = kh_detangle(i,k+1) + wt2*kh_lay_u(i,j,k)
1552  else ! This is a v-column.
1553  dh = 0.0
1554  denom = ((g%IareaT(i,j+1) + g%IareaT(i,j))*g%dx_Cv(i,j))
1555  if (denom > 0.0) &
1556  dh = i_4t * ((e(i,j+1,k) - e(i,j+1,k+1)) - &
1557  (e(i,j,k) - e(i,j,k+1))) / denom
1558  ! dH = I_4t * (h(i,j+1,k) - h(i,j,k)) / denom
1559 
1560  adh = abs(dh)
1561  sign = 1.0*us%Z_to_m ; if (dh < 0) sign = -1.0*us%Z_to_m
1562  sl_k = sign * (e(i,j+1,k)-e(i,j,k)) * g%IdyCv(i,j)
1563  sl_kp1 = sign * (e(i,j+1,k+1)-e(i,j,k+1)) * g%IdyCv(i,j)
1564 
1565  ! Add the incremental diffusviites to the surrounding interfaces.
1566  ! Adding more to the more steeply sloping layers (as below) makes
1567  ! the diffusivities more than twice as effective.
1568  denom = (sl_k**2 + sl_kp1**2)
1569  wt1 = 0.5 ; wt2 = 0.5
1570  if (denom > 0.0) then
1571  wt1 = sl_k**2 / denom ; wt2 = sl_kp1**2 / denom
1572  endif
1573  kh_detangle(i,k) = kh_detangle(i,k) + wt1*kh_lay_v(i,j,k)
1574  kh_detangle(i,k+1) = kh_detangle(i,k+1) + wt2*kh_lay_v(i,j,k)
1575  endif
1576 
1577  if (adh == 0.0) then
1578  kh_min_m(i,k+1) = 1.0 ; kh0_min_m(i,k+1) = 0.0
1579  kh_max_m(i,k+1) = 1.0 ; kh0_max_m(i,k+1) = 0.0
1580  kh_min_p(i,k) = 1.0 ; kh0_min_p(i,k) = 0.0
1581  kh_max_p(i,k) = 1.0 ; kh0_max_p(i,k) = 0.0
1582  elseif (adh > 0.0) then
1583  if (sl_k <= sl_kp1) then
1584  ! This case should only arise from nonlinearities in the equation of state.
1585  ! Treat it as though dedx(K) = dedx(K+1) & dH = 0.
1586  kh_min_m(i,k+1) = 1.0 ; kh0_min_m(i,k+1) = 0.0
1587  kh_max_m(i,k+1) = 1.0 ; kh0_max_m(i,k+1) = 0.0
1588  kh_min_p(i,k) = 1.0 ; kh0_min_p(i,k) = 0.0
1589  kh_max_p(i,k) = 1.0 ; kh0_max_p(i,k) = 0.0
1590  elseif (sl_k <= 0.0) then ! Both slopes are opposite to dH
1591  i_sl = -1.0 / sl_kp1
1592  rsl = -sl_k * i_sl ! 0 <= Rsl < 1
1593  irsl = 1e9 ; if (rsl > 1e-9) irsl = 1.0/rsl ! 1 < IRsl <= 1e9
1594 
1595  fn_r = rsl
1596  if (kh_max_max(i) > 0) &
1597  fn_r = min(sqrt(rsl), rsl + (adh * i_sl) / kh_max_max(i))
1598 
1599  kh_min_m(i,k+1) = fn_r ; kh0_min_m(i,k+1) = 0.0
1600  kh_max_m(i,k+1) = rsl ; kh0_max_m(i,k+1) = adh * i_sl
1601  kh_min_p(i,k) = irsl ; kh0_min_p(i,k) = -adh * (i_sl*irsl)
1602  kh_max_p(i,k) = 1.0/(fn_r + 1.0e-30) ; kh0_max_p(i,k) = 0.0
1603  elseif (sl_kp1 < 0.0) then ! Opposite (nonzero) signs of slopes.
1604  i_sl_k = 1e18 ; if (sl_k > 1e-18) i_sl_k = 1.0 / sl_k
1605  i_sl_kp1 = 1e18 ; if (-sl_kp1 > 1e-18) i_sl_kp1 = -1.0 / sl_kp1
1606 
1607  kh_min_m(i,k+1) = 0.0 ; kh0_min_m(i,k+1) = 0.0
1608  kh_max_m(i,k+1) = - sl_k*i_sl_kp1 ; kh0_max_m(i,k+1) = adh*i_sl_kp1
1609  kh_min_p(i,k) = 0.0 ; kh0_min_p(i,k) = 0.0
1610  kh_max_p(i,k) = sl_kp1*i_sl_k ; kh0_max_p(i,k) = adh*i_sl_k
1611 
1612  ! This limit does not use the slope weighting so that potentially
1613  ! sharp gradients in diffusivities are not forced to occur.
1614  kh_max = adh / (sl_k - sl_kp1)
1615  kh_min_max_p(i,k) = max(kh_min_max_p(i,k), kh_max)
1616  kh_min_max_m(i,k+1) = max(kh_min_max_m(i,k+1), kh_max)
1617  else ! Both slopes are of the same sign as dH.
1618  i_sl = 1.0 / sl_k
1619  rsl = sl_kp1 * i_sl ! 0 <= Rsl < 1
1620  irsl = 1e9 ; if (rsl > 1e-9) irsl = 1.0/rsl ! 1 < IRsl <= 1e9
1621 
1622  ! Rsl <= Fn_R <= 1
1623  fn_r = rsl
1624  if (kh_max_max(i) > 0) &
1625  fn_r = min(sqrt(rsl), rsl + (adh * i_sl) / kh_max_max(i))
1626 
1627  kh_min_m(i,k+1) = irsl ; kh0_min_m(i,k+1) = -adh * (i_sl*irsl)
1628  kh_max_m(i,k+1) = 1.0/(fn_r + 1.0e-30) ; kh0_max_m(i,k+1) = 0.0
1629  kh_min_p(i,k) = fn_r ; kh0_min_p(i,k) = 0.0
1630  kh_max_p(i,k) = rsl ; kh0_max_p(i,k) = adh * i_sl
1631  endif
1632  endif
1633  endif ; enddo ; enddo ! I-loop & k-loop
1634 
1635  do k=k_top,nz+1,nz+1-k_top ; do i=ish,ie ; if (do_i(i)) then
1636  ! The diffusivities at k_top and nz+1 are both fixed.
1637  kh_min_m(i,k) = 0.0 ; kh0_min_m(i,k) = 0.0
1638  kh_max_m(i,k) = 0.0 ; kh0_max_m(i,k) = 0.0
1639  kh_min_p(i,k) = 0.0 ; kh0_min_p(i,k) = 0.0
1640  kh_max_p(i,k) = 0.0 ; kh0_max_p(i,k) = 0.0
1641  kh_min_max_p(i,k) = kh_bg(i,k)
1642  kh_min_max_m(i,k) = kh_bg(i,k)
1643  endif ; enddo ; enddo ! I-loop and k_top/nz+1 loop
1644 
1645  ! Search for Kh that satisfy...
1646  ! Kh(I,K) >= Kh_min_m(I,K)*Kh(I,K-1) + Kh0_min_m(I,K)
1647  ! Kh(I,K) >= Kh_min_p(I,K)*Kh(I,K+1) + Kh0_min_p(I,K)
1648  ! Kh(I,K) <= Kh_max_m(I,K)*Kh(I,K-1) + Kh0_max_m(I,K)
1649  ! Kh(I,K) <= Kh_max_p(I,K)*Kh(I,K+1) + Kh0_max_p(I,K)
1650 
1651  ! Increase the diffusivies to satisfy the min constraints.
1652  ! All non-zero min constraints on one diffusivity are max constraints on another.
1653  do k=k_top+1,nz ; do i=ish,ie ; if (do_i(i)) then
1654  kh(i,k) = max(kh_bg(i,k), kh_detangle(i,k), &
1655  min(kh_min_m(i,k)*kh(i,k-1) + kh0_min_m(i,k), kh(i,k-1)))
1656 
1657  if (kh0_max_m(i,k) > kh_bg(i,k)) kh(i,k) = min(kh(i,k), kh0_max_m(i,k))
1658  if (kh0_max_p(i,k) > kh_bg(i,k)) kh(i,k) = min(kh(i,k), kh0_max_p(i,k))
1659  endif ; enddo ; enddo ! I-loop & k-loop
1660  ! This is still true... do i=ish,ie ; Kh(I,nz+1) = Kh_bg(I,nz+1) ; enddo
1661  do k=nz,k_top+1,-1 ; do i=ish,ie ; if (do_i(i)) then
1662  kh(i,k) = max(kh(i,k), min(kh_min_p(i,k)*kh(i,k+1) + kh0_min_p(i,k), kh(i,k+1)))
1663 
1664  kh_max = max(kh_min_max_p(i,k), kh_max_p(i,k)*kh(i,k+1) + kh0_max_p(i,k))
1665  kh(i,k) = min(kh(i,k), kh_max)
1666  endif ; enddo ; enddo ! I-loop & k-loop
1667  ! All non-zero min constraints on one diffusivity are max constraints on
1668  ! another layer, so the min constraints can now be discounted.
1669 
1670  ! Decrease the diffusivities to satisfy the max constraints.
1671  do k=k_top+1,nz ; do i=ish,ie ; if (do_i(i)) then
1672  kh_max = max(kh_min_max_m(i,k), kh_max_m(i,k)*kh(i,k-1) + kh0_max_m(i,k))
1673  if (kh(i,k) > kh_max) kh(i,k) = kh_max
1674  endif ; enddo ; enddo ! i- and K-loops
1675 
1676  ! This code tests the solutions...
1677 ! do i=ish,ie
1678 ! Sfn(:) = 0.0 ; uh_here(:) = 0.0
1679 ! do K=k_top,nz
1680 ! if ((Kh(i,K) > Kh_bg(i,K)) .or. (Kh(i,K+1) > Kh_bg(i,K+1))) then
1681 ! if (n==1) then ! u-point.
1682 ! if ((h(i+1,j,k) - h(i,j,k)) * &
1683 ! ((e(i+1,j,K)-e(i+1,j,K+1)) - (e(i,j,K)-e(i,j,K+1))) > 0.0) then
1684 ! Sfn(K) = -Kh(i,K) * (e(i+1,j,K)-e(i,j,K)) * G%IdxCu(I,j)
1685 ! Sfn(K+1) = -Kh(i,K+1) * (e(i+1,j,K+1)-e(i,j,K+1)) * G%IdxCu(I,j)
1686 ! uh_here(k) = (Sfn(K) - Sfn(K+1))*G%dy_Cu(I,j)
1687 ! if (abs(uh_here(k))*min(G%IareaT(i,j), G%IareaT(i+1,j)) > &
1688 ! (1e-10*GV%m_to_H)) then
1689 ! if (uh_here(k) * (h(i+1,j,k) - h(i,j,k)) > 0.0) then
1690 ! call MOM_error(WARNING, &
1691 ! "Corrective u-transport is up the thickness gradient.", .true.)
1692 ! endif
1693 ! if (((h(i,j,k) - 4.0*dt*G%IareaT(i,j)*uh_here(k)) - &
1694 ! (h(i+1,j,k) + 4.0*dt*G%IareaT(i+1,j)*uh_here(k))) * &
1695 ! (h(i,j,k) - h(i+1,j,k)) < 0.0) then
1696 ! call MOM_error(WARNING, &
1697 ! "Corrective u-transport is too large.", .true.)
1698 ! endif
1699 ! endif
1700 ! endif
1701 ! else ! v-point
1702 ! if ((h(i,j+1,k) - h(i,j,k)) * &
1703 ! ((e(i,j+1,K)-e(i,j+1,K+1)) - (e(i,j,K)-e(i,j,K+1))) > 0.0) then
1704 ! Sfn(K) = -Kh(i,K) * (e(i,j+1,K)-e(i,j,K)) * G%IdyCv(i,J)
1705 ! Sfn(K+1) = -Kh(i,K+1) * (e(i,j+1,K+1)-e(i,j,K+1)) * G%IdyCv(i,J)
1706 ! uh_here(k) = (Sfn(K) - Sfn(K+1))*G%dx_Cv(i,J)
1707 ! if (abs(uh_here(K))*min(G%IareaT(i,j), G%IareaT(i,j+1)) > &
1708 ! (1e-10*GV%m_to_H)) then
1709 ! if (uh_here(K) * (h(i,j+1,k) - h(i,j,k)) > 0.0) then
1710 ! call MOM_error(WARNING, &
1711 ! "Corrective v-transport is up the thickness gradient.", .true.)
1712 ! endif
1713 ! if (((h(i,j,k) - 4.0*dt*G%IareaT(i,j)*uh_here(K)) - &
1714 ! (h(i,j+1,k) + 4.0*dt*G%IareaT(i,j+1)*uh_here(K))) * &
1715 ! (h(i,j,k) - h(i,j+1,k)) < 0.0) then
1716 ! call MOM_error(WARNING, &
1717 ! "Corrective v-transport is too large.", .true.)
1718 ! endif
1719 ! endif
1720 ! endif
1721 ! endif ! u- or v- selection.
1722 ! ! de_dx(I,K) = (e(i+1,j,K)-e(i,j,K)) * G%IdxCu(I,j)
1723 ! endif
1724 ! enddo
1725 ! enddo
1726 
1727  if (n==1) then ! This is a u-column.
1728  do k=k_top+1,nz ; do i=ish,ie
1729  if (kh(i,k) > kh_u(i,j,k)) then
1730  dkh = (kh(i,k) - kh_u(i,j,k))
1731  int_slope_u(i,j,k) = dkh / kh(i,k)
1732  kh_u(i,j,k) = kh(i,k)
1733  endif
1734  enddo ; enddo
1735  else ! This is a v-column.
1736  do k=k_top+1,nz ; do i=ish,ie
1737  if (kh(i,k) > kh_v(i,j,k)) then
1738  dkh = kh(i,k) - kh_v(i,j,k)
1739  int_slope_v(i,j,k) = dkh / kh(i,k)
1740  kh_v(i,j,k) = kh(i,k)
1741  endif
1742  enddo ; enddo
1743  endif
1744 
1745  enddo ! j-loop
1746  enddo ! n-loop over u- and v- directions.
1747 
1748 end subroutine add_detangling_kh
1749 
1750 !> Initialize the thickness diffusion module/structure
1751 subroutine thickness_diffuse_init(Time, G, GV, US, param_file, diag, CDp, CS)
1752  type(time_type), intent(in) :: time !< Current model time
1753  type(ocean_grid_type), intent(in) :: g !< Ocean grid structure
1754  type(verticalgrid_type), intent(in) :: gv !< Vertical grid structure
1755  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
1756  type(param_file_type), intent(in) :: param_file !< Parameter file handles
1757  type(diag_ctrl), target, intent(inout) :: diag !< Diagnostics control structure
1758  type(cont_diag_ptrs), intent(inout) :: cdp !< Continuity equation diagnostics
1759  type(thickness_diffuse_cs), pointer :: cs !< Control structure for thickness diffusion
1760 
1761 ! This include declares and sets the variable "version".
1762 #include "version_variable.h"
1763  character(len=40) :: mdl = "MOM_thickness_diffuse" ! This module's name.
1764  real :: omega, strat_floor, flux_to_kg_per_s
1765 
1766  if (associated(cs)) then
1767  call mom_error(warning, &
1768  "Thickness_diffuse_init called with an associated control structure.")
1769  return
1770  else ; allocate(cs) ; endif
1771 
1772  cs%diag => diag
1773 
1774  ! Read all relevant parameters and write them to the model log.
1775  call log_version(param_file, mdl, version, "")
1776  call get_param(param_file, mdl, "THICKNESSDIFFUSE", cs%thickness_diffuse, &
1777  "If true, interface heights are diffused with a "//&
1778  "coefficient of KHTH.", default=.false.)
1779  call get_param(param_file, mdl, "KHTH", cs%Khth, &
1780  "The background horizontal thickness diffusivity.", &
1781  units = "m2 s-1", default=0.0)
1782  call get_param(param_file, mdl, "KHTH_SLOPE_CFF", cs%KHTH_Slope_Cff, &
1783  "The nondimensional coefficient in the Visbeck formula "//&
1784  "for the interface depth diffusivity", units="nondim", &
1785  default=0.0)
1786  call get_param(param_file, mdl, "KHTH_MIN", cs%KHTH_Min, &
1787  "The minimum horizontal thickness diffusivity.", &
1788  units = "m2 s-1", default=0.0)
1789  call get_param(param_file, mdl, "KHTH_MAX", cs%KHTH_Max, &
1790  "The maximum horizontal thickness diffusivity.", &
1791  units = "m2 s-1", default=0.0)
1792  call get_param(param_file, mdl, "KHTH_MAX_CFL", cs%max_Khth_CFL, &
1793  "The maximum value of the local diffusive CFL ratio that "//&
1794  "is permitted for the thickness diffusivity. 1.0 is the "//&
1795  "marginally unstable value in a pure layered model, but "//&
1796  "much smaller numbers (e.g. 0.1) seem to work better for "//&
1797  "ALE-based models.", units = "nondimensional", default=0.8)
1798 
1799 ! call get_param(param_file, mdl, "USE_QG_LEITH_GM", CS%QG_Leith_GM, &
1800 ! "If true, use the QG Leith viscosity as the GM coefficient.", &
1801 ! default=.false.)
1802 
1803  if (cs%max_Khth_CFL < 0.0) cs%max_Khth_CFL = 0.0
1804  call get_param(param_file, mdl, "DETANGLE_INTERFACES", cs%detangle_interfaces, &
1805  "If defined add 3-d structured enhanced interface height "//&
1806  "diffusivities to horizontally smooth jagged layers.", &
1807  default=.false.)
1808  cs%detangle_time = 0.0
1809  if (cs%detangle_interfaces) &
1810  call get_param(param_file, mdl, "DETANGLE_TIMESCALE", cs%detangle_time, &
1811  "A timescale over which maximally jagged grid-scale "//&
1812  "thickness variations are suppressed. This must be "//&
1813  "longer than DT, or 0 to use DT.", units = "s", default=0.0)
1814  call get_param(param_file, mdl, "KHTH_SLOPE_MAX", cs%slope_max, &
1815  "A slope beyond which the calculated isopycnal slope is "//&
1816  "not reliable and is scaled away.", units="nondim", default=0.01)
1817  call get_param(param_file, mdl, "KD_SMOOTH", cs%kappa_smooth, &
1818  "A diapycnal diffusivity that is used to interpolate "//&
1819  "more sensible values of T & S into thin layers.", &
1820  default=1.0e-6, scale=us%m_to_Z**2)
1821  call get_param(param_file, mdl, "KHTH_USE_FGNV_STREAMFUNCTION", cs%use_FGNV_streamfn, &
1822  "If true, use the streamfunction formulation of "//&
1823  "Ferrari et al., 2010, which effectively emphasizes "//&
1824  "graver vertical modes by smoothing in the vertical.", &
1825  default=.false.)
1826  call get_param(param_file, mdl, "FGNV_FILTER_SCALE", cs%FGNV_scale, &
1827  "A coefficient scaling the vertical smoothing term in the "//&
1828  "Ferrari et al., 2010, streamfunction formulation.", &
1829  default=1., do_not_log=.not.cs%use_FGNV_streamfn)
1830  call get_param(param_file, mdl, "FGNV_C_MIN", cs%FGNV_c_min, &
1831  "A minium wave speed used in the Ferrari et al., 2010, "//&
1832  "streamfunction formulation.", &
1833  default=0., units="m s-1", do_not_log=.not.cs%use_FGNV_streamfn)
1834  call get_param(param_file, mdl, "FGNV_STRAT_FLOOR", strat_floor, &
1835  "A floor for Brunt-Vasaila frequency in the Ferrari et al., 2010, "//&
1836  "streamfunction formulation, expressed as a fraction of planetary "//&
1837  "rotation, OMEGA. This should be tiny but non-zero to avoid degeneracy.", &
1838  default=1.e-15, units="nondim", do_not_log=.not.cs%use_FGNV_streamfn)
1839  call get_param(param_file, mdl, "OMEGA",omega, &
1840  "The rotation rate of the earth.", units="s-1", &
1841  default=7.2921e-5, do_not_log=.not.cs%use_FGNV_streamfn)
1842  if (cs%use_FGNV_streamfn) cs%N2_floor = (strat_floor*omega)**2
1843  call get_param(param_file, mdl, "DEBUG", cs%debug, &
1844  "If true, write out verbose debugging data.", &
1845  default=.false., debuggingparam=.true.)
1846 
1847  call get_param(param_file, mdl, "MEKE_GM_SRC_ALT", cs%GM_src_alt, &
1848  "If true, use the GM energy conversion form S^2*N^2*kappa rather \n"//&
1849  "than the streamfunction for the GM source term.", default=.false.)
1850  call get_param(param_file, mdl, "MEKE_GEOMETRIC", cs%MEKE_GEOMETRIC, &
1851  "If true, uses the GM coefficient formulation \n"//&
1852  "from the GEOMETRIC framework (Marshall et al., 2012).", default=.false.)
1853  if (cs%MEKE_GEOMETRIC) then
1854 
1855  call get_param(param_file, mdl, "MEKE_GEOMETRIC_EPSILON", cs%MEKE_GEOMETRIC_epsilon, &
1856  "Minimum Eady growth rate used in the calculation of \n"//&
1857  "GEOMETRIC thickness diffusivity.", units="s-1", default=1.0e-7)
1858 
1859  call get_param(param_file, mdl, "MEKE_GEOMETRIC_ALPHA", cs%MEKE_GEOMETRIC_alpha, &
1860  "The nondimensional coefficient governing the efficiency of the GEOMETRIC \n"//&
1861  "thickness diffusion.", units="nondim", default=0.05)
1862  endif
1863 
1864  call get_param(param_file, mdl, "USE_KH_IN_MEKE", cs%Use_KH_in_MEKE, &
1865  "If true, uses the thickness diffusivity calculated here to diffuse \n"//&
1866  "MEKE.", default=.false.)
1867 
1868  call get_param(param_file, mdl, "USE_GME", cs%use_GME_thickness_diffuse, &
1869  "If true, use the GM+E backscatter scheme in association \n"//&
1870  "with the Gent and McWilliams parameterization.", default=.false.)
1871 
1872  if (cs%use_GME_thickness_diffuse) then
1873  call safe_alloc_ptr(cs%KH_u_GME,g%IsdB,g%IedB,g%jsd,g%jed,g%ke+1)
1874  call safe_alloc_ptr(cs%KH_v_GME,g%isd,g%ied,g%JsdB,g%JedB,g%ke+1)
1875  endif
1876 
1877  if (gv%Boussinesq) then ; flux_to_kg_per_s = gv%Rho0
1878  else ; flux_to_kg_per_s = 1. ; endif
1879 
1880  cs%id_uhGM = register_diag_field('ocean_model', 'uhGM', diag%axesCuL, time, &
1881  'Time Mean Diffusive Zonal Thickness Flux', 'kg s-1', &
1882  y_cell_method='sum', v_extensive=.true., conversion=flux_to_kg_per_s)
1883  if (cs%id_uhGM > 0) call safe_alloc_ptr(cdp%uhGM,g%IsdB,g%IedB,g%jsd,g%jed,g%ke)
1884  cs%id_vhGM = register_diag_field('ocean_model', 'vhGM', diag%axesCvL, time, &
1885  'Time Mean Diffusive Meridional Thickness Flux', 'kg s-1', &
1886  x_cell_method='sum', v_extensive=.true., conversion=flux_to_kg_per_s)
1887  if (cs%id_vhGM > 0) call safe_alloc_ptr(cdp%vhGM,g%isd,g%ied,g%JsdB,g%JedB,g%ke)
1888 
1889  cs%id_GMwork = register_diag_field('ocean_model', 'GMwork', diag%axesT1, time, &
1890  'Integrated Tendency of Ocean Mesoscale Eddy KE from Parameterized Eddy Advection', &
1891  'W m-2', cmor_field_name='tnkebto', &
1892  cmor_long_name='Integrated Tendency of Ocean Mesoscale Eddy KE from Parameterized Eddy Advection',&
1893  cmor_standard_name='tendency_of_ocean_eddy_kinetic_energy_content_due_to_parameterized_eddy_advection')
1894  if (cs%id_GMwork > 0) call safe_alloc_ptr(cs%GMwork,g%isd,g%ied,g%jsd,g%jed)
1895 
1896  cs%id_KH_u = register_diag_field('ocean_model', 'KHTH_u', diag%axesCui, time, &
1897  'Parameterized mesoscale eddy advection diffusivity at U-point', 'm2 s-1')
1898  cs%id_KH_v = register_diag_field('ocean_model', 'KHTH_v', diag%axesCvi, time, &
1899  'Parameterized mesoscale eddy advection diffusivity at V-point', 'm2 s-1')
1900  cs%id_KH_t = register_diag_field('ocean_model', 'KHTH_t', diag%axesTL, time, &
1901  'Ocean Tracer Diffusivity due to Parameterized Mesoscale Advection', 'm2 s-1',&
1902  cmor_field_name='diftrblo', &
1903  cmor_long_name='Ocean Tracer Diffusivity due to Parameterized Mesoscale Advection', &
1904  cmor_units='m2 s-1', &
1905  cmor_standard_name='ocean_tracer_diffusivity_due_to_parameterized_mesoscale_advection')
1906 
1907  cs%id_KH_u1 = register_diag_field('ocean_model', 'KHTH_u1', diag%axesCu1, time, &
1908  'Parameterized mesoscale eddy advection diffusivity at U-points (2-D)', 'm2 s-1')
1909  cs%id_KH_v1 = register_diag_field('ocean_model', 'KHTH_v1', diag%axesCv1, time, &
1910  'Parameterized mesoscale eddy advection diffusivity at V-points (2-D)', 'm2 s-1')
1911  cs%id_KH_t1 = register_diag_field('ocean_model', 'KHTH_t1', diag%axesT1, time,&
1912  'Parameterized mesoscale eddy advection diffusivity at T-points (2-D)', 'm2 s-1')
1913 
1914  cs%id_slope_x = register_diag_field('ocean_model', 'neutral_slope_x', diag%axesCui, time, &
1915  'Zonal slope of neutral surface', 'nondim')
1916  if (cs%id_slope_x > 0) call safe_alloc_ptr(cs%diagSlopeX,g%IsdB,g%IedB,g%jsd,g%jed,g%ke+1)
1917  cs%id_slope_y = register_diag_field('ocean_model', 'neutral_slope_y', diag%axesCvi, time, &
1918  'Meridional slope of neutral surface', 'nondim')
1919  if (cs%id_slope_y > 0) call safe_alloc_ptr(cs%diagSlopeY,g%isd,g%ied,g%JsdB,g%JedB,g%ke+1)
1920  cs%id_sfn_x = register_diag_field('ocean_model', 'GM_sfn_x', diag%axesCui, time, &
1921  'Parameterized Zonal Overturning Streamfunction', 'm3 s-1')
1922  cs%id_sfn_y = register_diag_field('ocean_model', 'GM_sfn_y', diag%axesCvi, time, &
1923  'Parameterized Meridional Overturning Streamfunction', 'm3 s-1')
1924  cs%id_sfn_unlim_x = register_diag_field('ocean_model', 'GM_sfn_unlim_x', diag%axesCui, time, &
1925  'Parameterized Zonal Overturning Streamfunction before limiting/smoothing', &
1926  'm3 s-1', conversion=us%Z_to_m)
1927  cs%id_sfn_unlim_y = register_diag_field('ocean_model', 'GM_sfn_unlim_y', diag%axesCvi, time, &
1928  'Parameterized Meridional Overturning Streamfunction before limiting/smoothing', &
1929  'm3 s-1', conversion=us%Z_to_m)
1930 
1931 end subroutine thickness_diffuse_init
1932 
1933 !> Copies ubtav and vbtav from private type into arrays
1934 subroutine thickness_diffuse_get_kh(CS, KH_u_GME, KH_v_GME, G)
1935  type(thickness_diffuse_cs), pointer :: cs !< Control structure for
1936  !! this module
1937  type(ocean_grid_type), intent(in) :: g !< Grid structure
1938  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)+1), intent(inout) :: kh_u_gme!< interface height
1939  !! diffusivities in u-columns [m2 s-1]
1940  real, dimension(SZI_(G),SZJB_(G),SZK_(G)+1), intent(inout) :: kh_v_gme!< interface height
1941  !! diffusivities in v-columns [m2 s-1]
1942  ! Local variables
1943  integer :: i,j,k
1944 
1945  do k=1,g%ke+1 ; do j = g%jsc, g%jec ; do i = g%isc-1, g%iec
1946  kh_u_gme(i,j,k) = cs%KH_u_GME(i,j,k)
1947  enddo ; enddo ; enddo
1948 
1949  do k=1,g%ke+1 ; do j = g%jsc-1, g%jec ; do i = g%isc, g%iec
1950  kh_v_gme(i,j,k) = cs%KH_v_GME(i,j,k)
1951  enddo ; enddo ; enddo
1952 
1953 end subroutine thickness_diffuse_get_kh
1954 
1955 !> Deallocate the thickness diffusion control structure
1956 subroutine thickness_diffuse_end(CS)
1957  type(thickness_diffuse_cs), pointer :: cs !< Control structure for thickness diffusion
1958 
1959  if (associated(cs)) deallocate(cs)
1960 end subroutine thickness_diffuse_end
1961 
1962 !> \namespace mom_thickness_diffuse
1963 !!
1964 !! \section section_gm Thickness diffusion (aka Gent-McWilliams)
1965 !!
1966 !! Thickness diffusion is implemented via along-layer mass fluxes
1967 !! \f[
1968 !! h^\dagger \leftarrow h^n - \Delta t \nabla \cdot ( \vec{uh}^* )
1969 !! \f]
1970 !! where the mass fluxes are cast as the difference in vector streamfunction
1971 !!
1972 !! \f[
1973 !! \vec{uh}^* = \delta_k \vec{\psi} .
1974 !! \f]
1975 !!
1976 !! The GM implementation of thickness diffusion made the streamfunction proportional to the potential density slope
1977 !! \f[
1978 !! \vec{\psi} = - \kappa_h \frac{\nabla_z \rho}{\partial_z \rho}
1979 !! = \frac{g\kappa_h}{\rho_o} \frac{\nabla \rho}{N^2} = \kappa_h \frac{M^2}{N^2}
1980 !! \f]
1981 !! but for robustness the scheme is implemented as
1982 !! \f[
1983 !! \vec{\psi} = \kappa_h \frac{M^2}{\sqrt{N^4 + M^4}}
1984 !! \f]
1985 !! since the quantity \f$\frac{M^2}{\sqrt{N^2 + M^2}}\f$ is bounded between $-1$ and $1$ and does not change sign
1986 !! if \f$N^2<0\f$.
1987 !!
1988 !! Optionally, the method of Ferrari et al, 2010, can be used to obtain the streamfunction which solves the
1989 !! vertically elliptic equation:
1990 !! \f[
1991 !! \gamma_F \partial_z c^2 \partial_z \psi - N_*^2 \psi = ( 1 + \gamma_F ) \kappa_h N_*^2 \frac{M^2}{\sqrt{N^4+M^4}}
1992 !! \f]
1993 !! which recovers the previous streamfunction relation in the limit that \f$ c \rightarrow 0 \f$.
1994 !! Here, \f$c=\max(c_{min},c_g)\f$ is the maximum of either \f$c_{min}\f$ and either the first baroclinic mode
1995 !! wave-speed or the equivalent barotropic mode wave-speed.
1996 !! \f$N_*^2 = \max(N^2,0)\f$ is a non-negative form of the square of the Brunt-Vaisala frequency.
1997 !! The parameter \f$\gamma_F\f$ is used to reduce the vertical smoothing length scale.
1998 !! This elliptic form for \f$ \psi \f$ is turned on with the logical <code>KHTH_USE_FGNV_STREAMFUNCTION</code>.
1999 !!
2000 !! Thickness diffusivities are calculated independently at u- and v-points using the following expression
2001 !! \f[
2002 !! \kappa_h = \left( \kappa_o + \alpha_{s} L_{s}^2 < S N > + \alpha_{M} \kappa_{M} \right) r(\Delta x,L_d)
2003 !! \f]
2004 !! where \f$ S \f$ is the isoneutral slope magnitude, \f$ N \f$ is the square root of Brunt-Vaisala frequency,
2005 !! \f$\kappa_{M}\f$ is the diffusivity calculated by the MEKE parameterization (mom_meke module) and
2006 !! \f$ r(\Delta x,L_d) \f$ is a function of the local resolution (ratio of grid-spacing, \f$\Delta x\f$,
2007 !! to deformation radius, \f$L_d\f$). The length \f$L_s\f$ is provided by the mom_lateral_mixing_coeffs module
2008 !! (enabled with <code>USE_VARIABLE_MIXING=True</code> and the term \f$<SN>\f$ is the vertical average slope
2009 !! times the Brunt-Vaisala frequency prescribed by Visbeck et al., 1996.
2010 !!
2011 !! The result of the above expression is subsequently bounded by minimum and maximum values, including an upper
2012 !! diffusivity consistent with numerical stability (\f$ \kappa_{cfl} \f$ is calculated internally).
2013 !! \f[
2014 !! \kappa_h \leftarrow \min{\left( \kappa_{max}, \kappa_{cfl}, \max{\left( \kappa_{min}, \kappa_h \right)} \right)}
2015 !! f(c_g,z)
2016 !! \f]
2017 !!
2018 !! where \f$f(c_g,z)\f$ is a vertical structure function.
2019 !! \f$f(c_g,z)\f$ is calculated in module mom_lateral_mixing_coeffs.
2020 !! If <code>KHTH_USE_EBT_STRUCT=True</code> then \f$f(c_g,z)\f$ is set to look like the equivalent barotropic
2021 !! modal velocity structure. Otherwise \f$f(c_g,z)=1\f$ and the diffusivity is independent of depth.
2022 !!
2023 !! In order to calculate meaningful slopes in vanished layers, temporary copies of the thermodynamic variables
2024 !! are passed through a vertical smoother, function vert_fill_ts():
2025 !! \f{eqnarray*}{
2026 !! \left[ 1 + \Delta t \kappa_{smth} \frac{\partial^2}{\partial_z^2} \right] \theta & \leftarrow & \theta \\
2027 !! \left[ 1 + \Delta t \kappa_{smth} \frac{\partial^2}{\partial_z^2} \right] s & \leftarrow & s
2028 !! \f}
2029 !!
2030 !! \subsection section_khth_module_parameters Module mom_thickness_diffuse parameters
2031 !!
2032 !! | Symbol | Module parameter |
2033 !! | ------ | --------------- |
2034 !! | - | <code>THICKNESSDIFFUSE</code> |
2035 !! | \f$ \kappa_o \f$ | <code>KHTH</code> |
2036 !! | \f$ \alpha_{s} \f$ | <code>KHTH_SLOPE_CFF</code> |
2037 !! | \f$ \kappa_{min} \f$ | <code>KHTH_MIN</code> |
2038 !! | \f$ \kappa_{max} \f$ | <code>KHTH_MAX</code> |
2039 !! | - | <code>KHTH_MAX_CFL</code> |
2040 !! | \f$ \kappa_{smth} \f$ | <code>KD_SMOOTH</code> |
2041 !! | \f$ \alpha_{M} \f$ | <code>MEKE_KHTH_FAC</code> (from mom_meke module) |
2042 !! | - | <code>KHTH_USE_EBT_STRUCT</code> (from mom_lateral_mixing_coeffs module) |
2043 !! | - | <code>KHTH_USE_FGNV_STREAMFUNCTION</code> |
2044 !! | \f$ \gamma_F \f$ | <code>FGNV_FILTER_SCALE</code> |
2045 !! | \f$ c_{min} \f$ | <code>FGNV_C_MIN</code> |
2046 !!
2047 !! \subsection section_khth_module_reference References
2048 !!
2049 !! Ferrari, R., S.M. Griffies, A.J.G. Nurser and G.K. Vallis, 2010:
2050 !! A boundary-value problem for the parameterized mesoscale eddy transport.
2051 !! Ocean Modelling, 32, 143-156. http://doi.org/10.1016/j.ocemod.2010.01.004
2052 !!
2053 !! Viscbeck, M., J.C. Marshall, H. Jones, 1996:
2054 !! On he dynamics of convective "chimneys" in the ocean. J. Phys. Oceangr., 26, 1721-1734.
2055 !! http://dx.doi.org/10.1175/1520-0485(1996)026%3C1721:DOICRI%3E2.0.CO;2
2056 
2057 end module mom_thickness_diffuse
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_domains::pass_var
Do a halo update on an array.
Definition: MOM_domains.F90:49
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_isopycnal_slopes
Calculations of isoneutral slopes and stratification.
Definition: MOM_isopycnal_slopes.F90:2
mom_diag_mediator::post_data
Make a diagnostic available for averaging or output.
Definition: MOM_diag_mediator.F90:70
mom_domains::pass_vector
Do a halo update on a pair of arrays representing the two components of a vector.
Definition: MOM_domains.F90:54
mom_thickness_diffuse::thickness_diffuse_cs
Control structure for thickness diffusion.
Definition: MOM_thickness_diffuse.F90:37
mom_verticalgrid::verticalgrid_type
Describes the vertical ocean grid, including unit conversion factors.
Definition: MOM_verticalGrid.F90:24
mom_meke_types::meke_type
This type is used to exchange information related to the MEKE calculations.
Definition: MOM_MEKE_types.F90:8
mom_domains
Describes the decomposed MOM domain and has routines for communications across PEs.
Definition: MOM_domains.F90:2
mom_lateral_mixing_coeffs::varmix_cs
Variable mixing coefficients.
Definition: MOM_lateral_mixing_coeffs.F90:26
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_variables::cont_diag_ptrs
Pointers to arrays with transports, which can later be used for derived diagnostics,...
Definition: MOM_variables.F90:185
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_interface_heights::find_eta
Calculates the heights of sruface or all interfaces from layer thicknesses.
Definition: MOM_interface_heights.F90:21
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_debugging
Provides checksumming functions for debugging.
Definition: MOM_debugging.F90:7
mom_thickness_diffuse
Thickness diffusion (or Gent McWilliams)
Definition: MOM_thickness_diffuse.F90:2
mom_lateral_mixing_coeffs
Variable mixing coefficients.
Definition: MOM_lateral_mixing_coeffs.F90:2
mom_error_handler
Routines for error handling and I/O management.
Definition: MOM_error_handler.F90:2
mom_interface_heights
Functions for calculating interface heights, including free surface height.
Definition: MOM_interface_heights.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