MOM6
MOM_hor_visc.F90
1 !> Calculates horizontal viscosity and viscous stresses
3 
4 ! This file is part of MOM6. See LICENSE.md for the license.
5 
6 use mom_diag_mediator, only : post_data, register_diag_field, safe_alloc_ptr
7 use mom_diag_mediator, only : diag_ctrl, time_type
8 use mom_domains, only : pass_var, corner, pass_vector
9 use mom_error_handler, only : mom_error, fatal, warning, is_root_pe
11 use mom_grid, only : ocean_grid_type
12 use mom_lateral_mixing_coeffs, only : varmix_cs, calc_qg_leith_viscosity
13 use mom_barotropic, only : barotropic_cs, barotropic_get_tav
14 use mom_thickness_diffuse, only : thickness_diffuse_cs, thickness_diffuse_get_kh
15 use mom_io, only : mom_read_data, slasher
16 use mom_meke_types, only : meke_type
17 use mom_open_boundary, only : ocean_obc_type, obc_direction_e, obc_direction_w
18 use mom_open_boundary, only : obc_direction_n, obc_direction_s, obc_none
21 
22 implicit none ; private
23 
24 #include <MOM_memory.h>
25 
26 public horizontal_viscosity, hor_visc_init, hor_visc_end
27 
28 !> Control structure for horizontal viscosity
29 type, public :: hor_visc_cs ; private
30  logical :: laplacian !< Use a Laplacian horizontal viscosity if true.
31  logical :: biharmonic !< Use a biharmonic horizontal viscosity if true.
32  logical :: no_slip !< If true, no slip boundary conditions are used.
33  !! Otherwise free slip boundary conditions are assumed.
34  !! The implementation of the free slip boundary
35  !! conditions on a C-grid is much cleaner than the
36  !! no slip boundary conditions. The use of free slip
37  !! b.c.s is strongly encouraged. The no slip b.c.s
38  !! are not implemented with the biharmonic viscosity.
39  logical :: bound_kh !< If true, the Laplacian coefficient is locally
40  !! limited to guarantee stability.
41  logical :: better_bound_kh !< If true, use a more careful bounding of the
42  !! Laplacian viscosity to guarantee stability.
43  logical :: bound_ah !< If true, the biharmonic coefficient is locally
44  !! limited to guarantee stability.
45  logical :: better_bound_ah !< If true, use a more careful bounding of the
46  !! biharmonic viscosity to guarantee stability.
47  real :: bound_coef !< The nondimensional coefficient of the ratio of
48  !! the viscosity bounds to the theoretical maximum
49  !! for stability without considering other terms [nondim].
50  !! The default is 0.8.
51  logical :: smagorinsky_kh !< If true, use Smagorinsky nonlinear eddy
52  !! viscosity. KH is the background value.
53  logical :: smagorinsky_ah !< If true, use a biharmonic form of Smagorinsky
54  !! nonlinear eddy viscosity. AH is the background.
55  logical :: leith_kh !< If true, use 2D Leith nonlinear eddy
56  !! viscosity. KH is the background value.
57  logical :: modified_leith !< If true, use extra component of Leith viscosity
58  !! to damp divergent flow. To use, still set Leith_Kh=.TRUE.
59  logical :: use_beta_in_leith !< If true, includes the beta term in the Leith viscosity
60  logical :: leith_ah !< If true, use a biharmonic form of 2D Leith
61  !! nonlinear eddy viscosity. AH is the background.
62  logical :: use_qg_leith_visc !< If true, use QG Leith nonlinear eddy viscosity.
63  !! KH is the background value.
64  logical :: bound_coriolis !< If true & SMAGORINSKY_AH is used, the biharmonic
65  !! viscosity is modified to include a term that
66  !! scales quadratically with the velocity shears.
67  logical :: use_kh_bg_2d !< Read 2d background viscosity from a file.
68  real :: kh_bg_min !< The minimum value allowed for Laplacian horizontal
69  !! viscosity [m2 T-1 ~> m2 s-1]. The default is 0.0
70  logical :: use_land_mask !< Use the land mask for the computation of thicknesses
71  !! at velocity locations. This eliminates the dependence on
72  !! arbitrary values over land or outside of the domain.
73  !! Default is False to maintain answers with legacy experiments
74  !! but should be changed to True for new experiments.
75  logical :: anisotropic !< If true, allow anisotropic component to the viscosity.
76  real :: kh_aniso !< The anisotropic viscosity [m2 T-1 ~> m2 s-1].
77  logical :: dynamic_aniso !< If true, the anisotropic viscosity is recomputed as a function
78  !! of state. This is set depending on ANISOTROPIC_MODE.
79  logical :: res_scale_meke !< If true, the viscosity contribution from MEKE is scaled by
80  !! the resolution function.
81  logical :: use_gme !< If true, use GME backscatter scheme.
82  logical :: answers_2018 !< If true, use the order of arithmetic and expressions that recover the
83  !! answers from the end of 2018. Otherwise, use updated and more robust
84  !! forms of the same expressions.
85 
86  real allocable_, dimension(NIMEM_,NJMEM_) :: kh_bg_xx
87  !< The background Laplacian viscosity at h points [m2 T-1 ~> m2 s-1].
88  !! The actual viscosity may be the larger of this
89  !! viscosity and the Smagorinsky and Leith viscosities.
90  real allocable_, dimension(NIMEM_,NJMEM_) :: kh_bg_2d
91  !< The background Laplacian viscosity at h points [m2 T-1 ~> m2 s-1].
92  !! The actual viscosity may be the larger of this
93  !! viscosity and the Smagorinsky and Leith viscosities.
94  real allocable_, dimension(NIMEM_,NJMEM_) :: ah_bg_xx
95  !< The background biharmonic viscosity at h points [m4 T-1 ~> m4 s-1].
96  !! The actual viscosity may be the larger of this
97  !! viscosity and the Smagorinsky and Leith viscosities.
98  real allocable_, dimension(NIMEM_,NJMEM_) :: biharm5_const2_xx
99  !< A constant relating the biharmonic viscosity to the
100  !! square of the velocity shear [m4 s]. This value is
101  !! set to be the magnitude of the Coriolis terms once the
102  !! velocity differences reach a value of order 1/2 MAXVEL.
103  real allocable_, dimension(NIMEM_,NJMEM_) :: reduction_xx
104  !< The amount by which stresses through h points are reduced
105  !! due to partial barriers. Nondimensional.
106  real allocable_, dimension(NIMEM_,NJMEM_) :: &
107  kh_max_xx, & !< The maximum permitted Laplacian viscosity [m2 T-1 ~> m2 s-1].
108  ah_max_xx, & !< The maximum permitted biharmonic viscosity [m4 T-1 ~> m4 s-1].
109  n1n2_h, & !< Factor n1*n2 in the anisotropic direction tensor at h-points
110  n1n1_m_n2n2_h !< Factor n1**2-n2**2 in the anisotropic direction tensor at h-points
111 
112  real allocable_, dimension(NIMEMB_PTR_,NJMEMB_PTR_) :: kh_bg_xy
113  !< The background Laplacian viscosity at q points [m2 T-1 ~> m2 s-1].
114  !! The actual viscosity may be the larger of this
115  !! viscosity and the Smagorinsky and Leith viscosities.
116  real allocable_, dimension(NIMEMB_PTR_,NJMEMB_PTR_) :: ah_bg_xy
117  !< The background biharmonic viscosity at q points [m4 T-1 ~> m4 s-1].
118  !! The actual viscosity may be the larger of this
119  !! viscosity and the Smagorinsky and Leith viscosities.
120  real allocable_, dimension(NIMEMB_PTR_,NJMEMB_PTR_) :: biharm5_const2_xy
121  !< A constant relating the biharmonic viscosity to the
122  !! square of the velocity shear [m4 s]. This value is
123  !! set to be the magnitude of the Coriolis terms once the
124  !! velocity differences reach a value of order 1/2 MAXVEL.
125  real allocable_, dimension(NIMEMB_PTR_,NJMEMB_PTR_) :: reduction_xy
126  !< The amount by which stresses through q points are reduced
127  !! due to partial barriers [nondim].
128  real allocable_, dimension(NIMEMB_PTR_,NJMEMB_PTR_) :: &
129  kh_max_xy, & !< The maximum permitted Laplacian viscosity [m2 T-1 ~> m2 s-1].
130  ah_max_xy, & !< The maximum permitted biharmonic viscosity [m4 T-1 ~> m4 s-1].
131  n1n2_q, & !< Factor n1*n2 in the anisotropic direction tensor at q-points
132  n1n1_m_n2n2_q !< Factor n1**2-n2**2 in the anisotropic direction tensor at q-points
133 
134  real allocable_, dimension(NIMEM_,NJMEM_) :: &
135  dx2h, & !< Pre-calculated dx^2 at h points [m2]
136  dy2h, & !< Pre-calculated dy^2 at h points [m2]
137  dx_dyt, & !< Pre-calculated dx/dy at h points [nondim]
138  dy_dxt !< Pre-calculated dy/dx at h points [nondim]
139  real allocable_, dimension(NIMEMB_PTR_,NJMEMB_PTR_) :: &
140  dx2q, & !< Pre-calculated dx^2 at q points [m2]
141  dy2q, & !< Pre-calculated dy^2 at q points [m2]
142  dx_dybu, & !< Pre-calculated dx/dy at q points [nondim]
143  dy_dxbu !< Pre-calculated dy/dx at q points [nondim]
144  real allocable_, dimension(NIMEMB_PTR_,NJMEM_) :: &
145  idx2dycu, & !< 1/(dx^2 dy) at u points [m-3]
146  idxdy2u !< 1/(dx dy^2) at u points [m-3]
147  real allocable_, dimension(NIMEM_,NJMEMB_PTR_) :: &
148  idx2dycv, & !< 1/(dx^2 dy) at v points [m-3]
149  idxdy2v !< 1/(dx dy^2) at v points [m-3]
150 
151  ! The following variables are precalculated time-invariant combinations of
152  ! parameters and metric terms.
153  real allocable_, dimension(NIMEM_,NJMEM_) :: &
154  laplac2_const_xx, & !< Laplacian metric-dependent constants [m2]
155  biharm5_const_xx, & !< Biharmonic metric-dependent constants [m5]
156  laplac3_const_xx, & !< Laplacian metric-dependent constants [m3]
157  biharm_const_xx, & !< Biharmonic metric-dependent constants [m4]
158  biharm_const2_xx !< Biharmonic metric-dependent constants [T m4 ~> s m4]
159 
160  real allocable_, dimension(NIMEMB_PTR_,NJMEMB_PTR_) :: &
161  laplac2_const_xy, & !< Laplacian metric-dependent constants [m2]
162  biharm5_const_xy, & !< Biharmonic metric-dependent constants [m5]
163  laplac3_const_xy, & !< Laplacian metric-dependent constants [m3]
164  biharm_const_xy, & !< Biharmonic metric-dependent constants [m4]
165  biharm_const2_xy !< Biharmonic metric-dependent constants [T m4 ~> s m4]
166 
167  type(diag_ctrl), pointer :: diag => null() !< structure to regulate diagnostics
168 
169  !>@{
170  !! Diagnostic id
171  integer :: id_diffu = -1, id_diffv = -1
172  integer :: id_ah_h = -1, id_ah_q = -1
173  integer :: id_kh_h = -1, id_kh_q = -1
174  integer :: id_gme_coeff_h = -1, id_gme_coeff_q = -1
175  integer :: id_vort_xy_q = -1, id_div_xx_h = -1
176  integer :: id_frictwork = -1, id_frictworkintz = -1
177  integer :: id_frictworkmax = -1
178  integer :: id_frictwork_diss = -1, id_frictwork_gme = -1
179  !!@}
180 
181 
182 end type hor_visc_cs
183 
184 contains
185 
186 !> Calculates the acceleration due to the horizontal viscosity.
187 !!
188 !! A combination of biharmonic and Laplacian forms can be used. The coefficient
189 !! may either be a constant or a shear-dependent form. The biharmonic is
190 !! determined by twice taking the divergence of an appropriately defined stress
191 !! tensor. The Laplacian is determined by doing so once.
192 !!
193 !! To work, the following fields must be set outside of the usual
194 !! is:ie range before this subroutine is called:
195 !! u[is-2:ie+2,js-2:je+2]
196 !! v[is-2:ie+2,js-2:je+2]
197 !! h[is-1:ie+1,js-1:je+1]
198 subroutine horizontal_viscosity(u, v, h, diffu, diffv, MEKE, VarMix, G, GV, US, &
199  CS, OBC, BT)
200  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure.
201  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure.
202  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
203  intent(in) :: u !< The zonal velocity [m s-1].
204  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
205  intent(in) :: v !< The meridional velocity [m s-1].
206  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
207  intent(inout) :: h !< Layer thicknesses [H ~> m or kg m-2].
208  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
209  intent(out) :: diffu !< Zonal acceleration due to convergence of
210  !! along-coordinate stress tensor [m s-1 T-1 ~> m s-2]
211  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
212  intent(out) :: diffv !< Meridional acceleration due to convergence
213  !! of along-coordinate stress tensor [m s-1 T-1 ~> m s-2].
214  type(meke_type), pointer :: meke !< Pointer to a structure containing fields
215  !! related to Mesoscale Eddy Kinetic Energy.
216  type(varmix_cs), pointer :: varmix !< Pointer to a structure with fields that
217  !! specify the spatially variable viscosities
218  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
219  type(hor_visc_cs), pointer :: cs !< Control structure returned by a previous
220  !! call to hor_visc_init.
221  type(ocean_obc_type), optional, pointer :: obc !< Pointer to an open boundary condition type
222  type(barotropic_cs), optional, pointer :: bt !< Pointer to a structure containing
223  !! barotropic velocities.
224 
225  ! Local variables
226  real, dimension(SZIB_(G),SZJ_(G)) :: &
227  u0, & ! Laplacian of u [m-1 s-1]
228  h_u, & ! Thickness interpolated to u points [H ~> m or kg m-2].
229  vort_xy_dy, & ! y-derivative of vertical vorticity (d/dy(dv/dx - du/dy)) [m-1 s-1]
230  div_xx_dx, & ! x-derivative of horizontal divergence (d/dx(du/dx + dv/dy)) [m-1 s-1]
231  ubtav ! zonal barotropic vel. ave. over baroclinic time-step [m s-1]
232  real, dimension(SZI_(G),SZJB_(G)) :: &
233  v0, & ! Laplacian of v [m-1 s-1]
234  h_v, & ! Thickness interpolated to v points [H ~> m or kg m-2].
235  vort_xy_dx, & ! x-derivative of vertical vorticity (d/dx(dv/dx - du/dy)) [m-1 s-1]
236  div_xx_dy, & ! y-derivative of horizontal divergence (d/dy(du/dx + dv/dy)) [m-1 s-1]
237  vbtav ! meridional barotropic vel. ave. over baroclinic time-step [m s-1]
238  real, dimension(SZI_(G),SZJ_(G)) :: &
239  dudx_bt, dvdy_bt, & ! components in the barotropic horizontal tension [s-1]
240  div_xx, & ! Estimate of horizontal divergence at h-points [s-1]
241  sh_xx, & ! horizontal tension (du/dx - dv/dy) including metric terms [s-1]
242  sh_xx_bt, & ! barotropic horizontal tension (du/dx - dv/dy) including metric terms [s-1]
243  str_xx,& ! str_xx is the diagonal term in the stress tensor [H m2 s-1 T-1 ~> m3 s-2 or kg s-2]
244  str_xx_gme,& ! smoothed diagonal term in the stress tensor from GME [H m2 s-1 T-1 ~> m3 s-2 or kg s-2]
245  bhstr_xx,& ! A copy of str_xx that only contains the biharmonic contribution
246  ! [H m2 T-1 s-1 ~> m3 s-2 or kg s-2]
247  frictworkintz, & ! depth integrated energy dissipated by lateral friction [W m-2]
248  leith_kh_h, & ! Leith Laplacian viscosity at h-points [m2 s-1]
249  leith_ah_h, & ! Leith bi-harmonic viscosity at h-points [m4 s-1]
250  beta_h, & ! Gradient of planetary vorticity at h-points [m-1 s-1]
251  grad_vort_mag_h, & ! Magnitude of vorticity gradient at h-points [m-1 s-1]
252  grad_vort_mag_h_2d, & ! Magnitude of 2d vorticity gradient at h-points [m-1 s-1]
253  grad_div_mag_h, & ! Magnitude of divergence gradient at h-points [m-1 s-1]
254  dudx, dvdy, & ! components in the horizontal tension [s-1]
255  grad_vel_mag_h, & ! Magnitude of the velocity gradient tensor squared at h-points [s-2]
256  grad_vel_mag_bt_h, & ! Magnitude of the barotropic velocity gradient tensor squared at h-points [s-2]
257  grad_d2vel_mag_h, & ! Magnitude of the Laplacian of the velocity vector, squared [m-2 s-2]
258  max_diss_rate_bt, & ! maximum possible energy dissipated by barotropic lateral friction [m2 s-3]
259  boundary_mask ! A mask that zeroes out cells with at least one land edge
260 
261  real, dimension(SZIB_(G),SZJB_(G)) :: &
262  dvdx, dudy, & ! components in the shearing strain [s-1]
263  dvdx_bt, dudy_bt, & ! components in the barotropic shearing strain [s-1]
264  sh_xy, & ! horizontal shearing strain (du/dy + dv/dx) including metric terms [s-1]
265  sh_xy_bt, & ! barotropic horizontal shearing strain (du/dy + dv/dx) inc. metric terms [s-1]
266  str_xy, & ! str_xy is the cross term in the stress tensor [H m2 s-2 ~> m3 s-2 or kg s-2]
267  str_xy_gme, & ! smoothed cross term in the stress tensor from GME [H m2 s-2]
268  bhstr_xy, & ! A copy of str_xy that only contains the biharmonic contribution
269  ! [H m2 s-2 ~> m3 s-2 or kg s-2]
270  vort_xy, & ! Vertical vorticity (dv/dx - du/dy) including metric terms [s-1]
271  leith_kh_q, & ! Leith Laplacian viscosity at q-points [m2 s-1]
272  leith_ah_q, & ! Leith bi-harmonic viscosity at q-points [m4 s-1]
273  beta_q, & ! Gradient of planetary vorticity at q-points [m-1 s-1]
274  grad_vort_mag_q, & ! Magnitude of vorticity gradient at q-points [m-1 s-1]
275  grad_vort_mag_q_2d, & ! Magnitude of 2d vorticity gradient at q-points [m-1 s-1]
276  grad_div_mag_q, & ! Magnitude of divergence gradient at q-points [m-1 s-1]
277  grad_vel_mag_q, & ! Magnitude of the velocity gradient tensor squared at q-points [s-2]
278  hq, & ! harmonic mean of the harmonic means of the u- & v point thicknesses [H ~> m or kg m-2]
279  ! This form guarantees that hq/hu < 4.
280  grad_vel_mag_bt_q ! Magnitude of the barotropic velocity gradient tensor squared at q-points [s-2]
281 
282  real, dimension(SZIB_(G),SZJB_(G),SZK_(G)) :: &
283  ah_q, & ! biharmonic viscosity at corner points [m4 T-1 ~> m4 s-1]
284  kh_q, & ! Laplacian viscosity at corner points [m2 s-1]
285  vort_xy_q, & ! vertical vorticity at corner points [s-1]
286  gme_coeff_q !< GME coeff. at q-points [m2 T-1 ~> m2 s-1]
287 
288  ! These 3-d arrays are unused.
289  ! real, dimension(SZIB_(G),SZJ_(G),SZK_(G)+1) :: &
290  ! KH_u_GME !< interface height diffusivities in u-columns [m2 s-1]
291  ! real, dimension(SZI_(G),SZJB_(G),SZK_(G)+1) :: &
292  ! KH_v_GME !< interface height diffusivities in v-columns [m2 s-1]
293  real, dimension(SZI_(G),SZJ_(G),SZK_(G)) :: &
294  ah_h, & ! biharmonic viscosity at thickness points [m4 T-1 ~> m4 s-1]
295  kh_h, & ! Laplacian viscosity at thickness points [m2 T-1 ~> m2 s-1]
296  diss_rate, & ! MKE dissipated by parameterized shear production [m2 s-3]
297  max_diss_rate, & ! maximum possible energy dissipated by lateral friction [m2 s-3]
298  target_diss_rate_gme, & ! the maximum theoretical dissipation plus the amount spuriously dissipated
299  ! by friction [m2 s-3]
300  frictwork, & ! work done by MKE dissipation mechanisms [W m-2]
301  frictwork_diss, & ! negative definite work done by MKE dissipation mechanisms [W m-2]
302  frictworkmax, & ! maximum possible work done by MKE dissipation mechanisms [W m-2]
303  frictwork_gme, & ! work done by GME [W m-2]
304  div_xx_h ! horizontal divergence [s-1]
305  ! real, dimension(SZI_(G),SZJ_(G),SZK_(G)+1) :: &
306  real, dimension(SZI_(G),SZJ_(G),SZK_(G)) :: &
307  ! KH_t_GME, & !< interface height diffusivities in t-columns [m2 s-1]
308  gme_coeff_h !< GME coeff. at h-points [m2 T-1 ~> m2 s-1]
309  real :: ah ! biharmonic viscosity [m4 T-1 ~> m4 s-1]
310  real :: kh ! Laplacian viscosity [m2 T-1 ~> m2 s-1]
311  real :: ahsm ! Smagorinsky biharmonic viscosity [m4 T-1 ~> m4 s-1]
312 ! real :: KhSm ! Smagorinsky Laplacian viscosity [m2 T-1 ~> m2 s-1]
313  real :: ahlth ! 2D Leith biharmonic viscosity [m4 T-1 ~> m4 s-1]
314 ! real :: KhLth ! 2D Leith Laplacian viscosity [m2 s-1]
315  real :: mod_leith ! nondimensional coefficient for divergence part of modified Leith
316  ! viscosity. Here set equal to nondimensional Laplacian Leith constant.
317  ! This is set equal to zero if modified Leith is not used.
318  real :: shear_mag ! magnitude of the shear [T-1 ~> s-1]
319  real :: vert_vort_mag ! magnitude of the vertical vorticity gradient [m-1 T-1 ~> m-1 s-1]
320  real :: h2uq, h2vq ! temporary variables [H2 ~> m2 or kg2 m-4].
321  real :: hu, hv ! Thicknesses interpolated by arithmetic means to corner
322  ! points; these are first interpolated to u or v velocity
323  ! points where masks are applied [H ~> m or kg m-2].
324 ! real :: hq ! harmonic mean of the harmonic means of the u- & v-
325 ! ! point thicknesses, in H; This form guarantees that hq/hu < 4.
326  real :: h_neglect ! thickness so small it can be lost in roundoff and so neglected [H ~> m or kg m-2]
327  real :: h_neglect3 ! h_neglect^3 [H3 ~> m3 or kg3 m-6]
328  real :: hrat_min ! minimum thicknesses at the 4 neighboring
329  ! velocity points divided by the thickness at the stress
330  ! point (h or q point) [nondim]
331  real :: visc_bound_rem ! fraction of overall viscous bounds that
332  ! remain to be applied [nondim]
333  real :: kh_scale ! A factor between 0 and 1 by which the horizontal
334  ! Laplacian viscosity is rescaled [nondim]
335  real :: roscl ! The scaling function for MEKE source term [nondim]
336  real :: fath ! abs(f) at h-point for MEKE source term [T-1 ~> s-1]
337  real :: local_strain ! Local variable for interpolating computed strain rates [s-1].
338  real :: meke_res_fn ! A copy of the resolution scaling factor if being applied to MEKE. Otherwise =1.
339  real :: gme_coeff ! The GME (negative) viscosity coefficient [m2 T-1 ~> m2 s-1]
340  real :: gme_coeff_limiter ! Maximum permitted value of the GME coefficient [m2 T-1 ~> m2 s-1]
341  real :: fwfrac ! Fraction of maximum theoretical energy transfer to use when scaling GME coefficient [nondim]
342  real :: dy_dxbu ! Ratio of meridional over zonal grid spacing at vertices [nondim]
343  real :: dx_dybu ! Ratio of zonal over meridiononal grid spacing at vertices [nondim]
344  real :: sh_f_pow ! The ratio of shear over the absolute value of f raised to some power and rescaled [nondim]
345  real :: backscat_subround ! The ratio of f over Shear_mag that is so small that the backscatter
346  ! calculation gives the same value as if f were 0 [nondim].
347  real :: h0_gme ! Depth used to scale down GME coefficient in shallow areas [Z ~> m]
348  logical :: rescale_kh, legacy_bound
349  logical :: find_frictwork
350  logical :: apply_obc = .false.
351  logical :: use_meke_ku
352  logical :: use_meke_au
353  integer :: is, ie, js, je, isq, ieq, jsq, jeq, nz
354  integer :: i, j, k, n
355  real :: inv_pi3, inv_pi2, inv_pi5
356  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
357  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
358 
359  h_neglect = gv%H_subroundoff
360  h_neglect3 = h_neglect**3
361  inv_pi3 = 1.0/((4.0*atan(1.0))**3)
362  inv_pi2 = 1.0/((4.0*atan(1.0))**2)
363  inv_pi5 = inv_pi3 * inv_pi2
364 
365  ah_h(:,:,:) = 0.0
366  kh_h(:,:,:) = 0.0
367 
368  if (present(obc)) then ; if (associated(obc)) then ; if (obc%OBC_pe) then
369  apply_obc = obc%Flather_u_BCs_exist_globally .or. obc%Flather_v_BCs_exist_globally
370  apply_obc = .true.
371  endif ; endif ; endif
372 
373  if (.not.associated(cs)) call mom_error(fatal, &
374  "MOM_hor_visc: Module must be initialized before it is used.")
375  if (.not.(cs%Laplacian .or. cs%biharmonic)) return
376 
377  find_frictwork = (cs%id_FrictWork > 0)
378  if (cs%id_FrictWorkIntz > 0) find_frictwork = .true.
379  if (associated(meke)) then
380  if (associated(meke%mom_src)) find_frictwork = .true.
381  backscat_subround = 0.0
382  if (find_frictwork .and. associated(meke%mom_src) .and. (meke%backscatter_Ro_c > 0.0) .and. &
383  (meke%backscatter_Ro_Pow /= 0.0)) &
384  backscat_subround = (1.0e-16/meke%backscatter_Ro_c)**(1.0/meke%backscatter_Ro_Pow)
385  endif
386 
387  rescale_kh = .false.
388  if (associated(varmix)) then
389  rescale_kh = varmix%Resoln_scaled_Kh
390  if (rescale_kh .and. &
391  (.not.associated(varmix%Res_fn_h) .or. .not.associated(varmix%Res_fn_q))) &
392  call mom_error(fatal, "MOM_hor_visc: VarMix%Res_fn_h and " //&
393  "VarMix%Res_fn_q both need to be associated with Resoln_scaled_Kh.")
394  endif
395  legacy_bound = (cs%Smagorinsky_Kh .or. cs%Leith_Kh) .and. &
396  (cs%bound_Kh .and. .not.cs%better_bound_Kh)
397 
398  ! Toggle whether to use a Laplacian viscosity derived from MEKE
399  use_meke_ku = associated(meke%Ku)
400  use_meke_au = associated(meke%Au)
401 
402  do j=js,je ; do i=is,ie
403  boundary_mask(i,j) = (g%mask2dCu(i,j) * g%mask2dCv(i,j) * g%mask2dCu(i-1,j) * g%mask2dCv(i,j-1))
404  enddo ; enddo
405 
406  if (cs%use_GME) then
407  ! GME tapers off above this depth
408  h0_gme = 1000.0*us%m_to_Z
409  fwfrac = 1.0
410  gme_coeff_limiter = 1e7*us%T_to_s
411 
412  ! initialize diag. array with zeros
413  gme_coeff_h(:,:,:) = 0.0
414  gme_coeff_q(:,:,:) = 0.0
415  str_xx_gme(:,:) = 0.0
416  str_xy_gme(:,:) = 0.0
417 
418 ! call pass_var(boundary_mask, G%Domain, complete=.true.)
419 
420  ! Get barotropic velocities and their gradients
421  call barotropic_get_tav(bt, ubtav, vbtav, g, us)
422  call pass_vector(ubtav, vbtav, g%Domain)
423 
424  !#GME# The following loop range should be: do j=js-1,je+1 ; do i=is-1,ie+1
425  do j=js,je ; do i=is,ie
426  dudx_bt(i,j) = cs%DY_dxT(i,j)*(g%IdyCu(i,j) * ubtav(i,j) - &
427  g%IdyCu(i-1,j) * ubtav(i-1,j))
428  dvdy_bt(i,j) = cs%DX_dyT(i,j)*(g%IdxCv(i,j) * vbtav(i,j) - &
429  g%IdxCv(i,j-1) * vbtav(i,j-1))
430  enddo; enddo
431 
432  !#GME# These should be combined into a vactor pass
433  call pass_var(dudx_bt, g%Domain, complete=.true.)
434  call pass_var(dvdy_bt, g%Domain, complete=.true.)
435 
436  !#GME# These loop bounds should be:
437  !#GME# do j=Jsq,Jeq+1 ; do i=Isq,Ieq+1
438  do j=jsq-1,jeq+2 ; do i=isq-1,ieq+2
439  sh_xx_bt(i,j) = dudx_bt(i,j) - dvdy_bt(i,j)
440  enddo ; enddo
441 
442  ! Components for the barotropic shearing strain
443  do j=js-2,jeq+1 ; do i=is-2,ieq+1
444  dvdx_bt(i,j) = cs%DY_dxBu(i,j)*(vbtav(i+1,j)*g%IdyCv(i+1,j) &
445  - vbtav(i,j)*g%IdyCv(i,j))
446  dudy_bt(i,j) = cs%DX_dyBu(i,j)*(ubtav(i,j+1)*g%IdxCu(i,j+1) &
447  - ubtav(i,j)*g%IdxCu(i,j))
448  enddo ; enddo
449 
450  !#GME# These should be combined into a vactor pass
451  call pass_var(dvdx_bt, g%Domain, position=corner, complete=.true.)
452  call pass_var(dudy_bt, g%Domain, position=corner, complete=.true.)
453 
454  if (cs%no_slip) then
455  !#GME# These loop bounds should be
456  !#GME# do J=js-1,Jeq ; do I=is-1,Ieq
457  do j=js-2,jeq+1 ; do i=is-2,ieq+1
458  sh_xy_bt(i,j) = (2.0-g%mask2dBu(i,j)) * ( dvdx_bt(i,j) + dudy_bt(i,j) )
459  enddo ; enddo
460  else
461  !#GME# These loop bounds should be
462  !#GME# do J=js-1,Jeq ; do I=is-1,Ieq
463  do j=js-2,jeq+1 ; do i=is-2,ieq+1
464  sh_xy_bt(i,j) = g%mask2dBu(i,j) * ( dvdx_bt(i,j) + dudy_bt(i,j) )
465  enddo ; enddo
466  endif
467 
468  ! Get thickness diffusivity for use in GME
469 ! call thickness_diffuse_get_KH(thickness_diffuse, KH_u_GME, KH_v_GME, G)
470 
471  !#GME# These loops bounds should probably be: do j=js-1,je+1 ; do i=is-1,is+1
472  !#GME# Group the 4-point sums so they are rotationally invariant.`
473  do j=jsq-1,jeq+2 ; do i=isq-1,ieq+2
474  grad_vel_mag_bt_h(i,j) = boundary_mask(i,j) * (dudx_bt(i,j)**2 + dvdy_bt(i,j)**2 + &
475  (0.25*(dvdx_bt(i,j)+dvdx_bt(i-1,j)+dvdx_bt(i,j-1)+dvdx_bt(i-1,j-1)) )**2 + &
476  (0.25*(dudy_bt(i,j)+dudy_bt(i-1,j)+dudy_bt(i,j-1)+dudy_bt(i-1,j-1)) )**2)
477  enddo ; enddo
478 
479  !#GME# max_diss_rate_bt is not used.
480  if (associated(meke)) then ; if (associated(meke%mom_src)) then
481  !#GME# These loops bounds should be: do j=js-1,je+1 ; do i=is-1,is+1
482  do j=jsq-1,jeq+2 ; do i=isq-1,ieq+2
483  max_diss_rate_bt(i,j) = 2.0*meke%MEKE(i,j) * grad_vel_mag_bt_h(i,j)
484  enddo ; enddo
485  endif ; endif
486 
487  !#GME# boundary_mask is defined at h points, not q points as used here.
488  !#GME# boundary_mask has only been defined over the range is:ie, js:je.
489  !#GME# Group the 4-point sums so they are rotationally invariant.`
490  !#GME# The following loop range should be: do J=js-1,Jeq ; do I=is-1,Ieq
491  do j=js-2,jeq+1 ; do i=is-2,ieq+1
492  grad_vel_mag_bt_q(i,j) = boundary_mask(i,j) * (dvdx_bt(i,j)**2 + dudy_bt(i,j)**2 + &
493  (0.25*(dudx_bt(i,j)+dudx_bt(i+1,j)+dudx_bt(i,j+1)+dudx_bt(i+1,j+1)))**2 + &
494  (0.25*(dvdy_bt(i,j)+dvdy_bt(i+1,j)+dvdy_bt(i,j+1)+dvdy_bt(i+1,j+1)) )**2)
495  enddo ; enddo
496 
497  endif ! use_GME
498 
499  !$OMP parallel do default(none) shared(Isq,Ieq,Jsq,Jeq,nz,CS,G,GV,u,v,is,js,ie,je,h, &
500  !$OMP rescale_Kh,VarMix,h_neglect,h_neglect3, &
501  !$OMP Kh_h,Ah_h,Kh_q,Ah_q,diffu,diffv,apply_OBC,OBC, &
502  !$OMP find_FrictWork,FrictWork,use_MEKE_Ku, &
503  !$OMP use_MEKE_Au, MEKE, hq, &
504  !$OMP mod_Leith, legacy_bound, div_xx_h, vort_xy_q) &
505  !$OMP private(u0, v0, sh_xx, str_xx, visc_bound_rem, &
506  !$OMP sh_xy, str_xy, Ah, Kh, AhSm, dvdx, dudy, &
507  !$OMP sh_xx_bt, sh_xy_bt, dvdx_bt, dudy_bt, &
508  !$OMP bhstr_xx, bhstr_xy,FatH,RoScl, hu, hv, h_u, h_v, &
509  !$OMP vort_xy,vort_xy_dx,vort_xy_dy,Vort_mag,AhLth, &
510  !$OMP div_xx, div_xx_dx, div_xx_dy, local_strain, &
511  !$OMP meke_res_fn,Sh_F_pow, &
512  !$OMP Shear_mag, h2uq, h2vq, hq, Kh_scale, hrat_min)
513  do k=1,nz
514 
515  ! The following are the forms of the horizontal tension and horizontal
516  ! shearing strain advocated by Smagorinsky (1993) and discussed in
517  ! Griffies and Hallberg (2000).
518 
519  ! Calculate horizontal tension
520  do j=jsq-1,jeq+2 ; do i=isq-1,ieq+2
521  dudx(i,j) = cs%DY_dxT(i,j)*(g%IdyCu(i,j) * u(i,j,k) - &
522  g%IdyCu(i-1,j) * u(i-1,j,k))
523  dvdy(i,j) = cs%DX_dyT(i,j)*(g%IdxCv(i,j) * v(i,j,k) - &
524  g%IdxCv(i,j-1) * v(i,j-1,k))
525  sh_xx(i,j) = dudx(i,j) - dvdy(i,j)
526  enddo ; enddo
527 
528  ! Components for the shearing strain
529  do j=js-2,jeq+1 ; do i=is-2,ieq+1
530  dvdx(i,j) = cs%DY_dxBu(i,j)*(v(i+1,j,k)*g%IdyCv(i+1,j) - v(i,j,k)*g%IdyCv(i,j))
531  dudy(i,j) = cs%DX_dyBu(i,j)*(u(i,j+1,k)*g%IdxCu(i,j+1) - u(i,j,k)*g%IdxCu(i,j))
532  enddo ; enddo
533 
534  ! Interpolate the thicknesses to velocity points.
535  ! The extra wide halos are to accommodate the cross-corner-point projections
536  ! in OBCs, which are not ordinarily be necessary, and might not be necessary
537  ! even with OBCs if the accelerations are zeroed at OBC points, in which
538  ! case the j-loop for h_u could collapse to j=js=1,je+1. -RWH
539  if (cs%use_land_mask) then
540  do j=js-2,je+2 ; do i=isq-1,ieq+1
541  h_u(i,j) = 0.5 * (g%mask2dT(i,j)*h(i,j,k) + g%mask2dT(i+1,j)*h(i+1,j,k))
542  enddo ; enddo
543  do j=jsq-1,jeq+1 ; do i=is-2,ie+2
544  h_v(i,j) = 0.5 * (g%mask2dT(i,j)*h(i,j,k) + g%mask2dT(i,j+1)*h(i,j+1,k))
545  enddo ; enddo
546  else
547  do j=js-2,je+2 ; do i=isq-1,ieq+1
548  h_u(i,j) = 0.5 * (h(i,j,k) + h(i+1,j,k))
549  enddo ; enddo
550  do j=jsq-1,jeq+1 ; do i=is-2,ie+2
551  h_v(i,j) = 0.5 * (h(i,j,k) + h(i,j+1,k))
552  enddo ; enddo
553  endif
554 
555  ! Adjust contributions to shearing strain and interpolated values of
556  ! thicknesses on open boundaries.
557  if (apply_obc) then ; do n=1,obc%number_of_segments
558  j = obc%segment(n)%HI%JsdB ; i = obc%segment(n)%HI%IsdB
559  if (obc%zero_strain .or. obc%freeslip_strain .or. obc%computed_strain) then
560  if (obc%segment(n)%is_N_or_S .and. (j >= js-2) .and. (j <= jeq+1)) then
561  do i=obc%segment(n)%HI%IsdB,obc%segment(n)%HI%IedB
562  if (obc%zero_strain) then
563  dvdx(i,j) = 0. ; dudy(i,j) = 0.
564  elseif (obc%freeslip_strain) then
565  dudy(i,j) = 0.
566  elseif (obc%computed_strain) then
567  if (obc%segment(n)%direction == obc_direction_n) then
568  dudy(i,j) = 2.0*cs%DX_dyBu(i,j)*(obc%segment(n)%tangential_vel(i,j,k) - u(i,j,k))*g%IdxCu(i,j)
569  else
570  dudy(i,j) = 2.0*cs%DX_dyBu(i,j)*(u(i,j+1,k) - obc%segment(n)%tangential_vel(i,j,k))*g%IdxCu(i,j+1)
571  endif
572  elseif (obc%specified_strain) then
573  if (obc%segment(n)%direction == obc_direction_n) then
574  dudy(i,j) = cs%DX_dyBu(i,j)*obc%segment(n)%tangential_grad(i,j,k)*g%IdxCu(i,j)*g%dxBu(i,j)
575  else
576  dudy(i,j) = cs%DX_dyBu(i,j)*obc%segment(n)%tangential_grad(i,j,k)*g%IdxCu(i,j+1)*g%dxBu(i,j)
577  endif
578  endif
579  enddo
580  elseif (obc%segment(n)%is_E_or_W .and. (i >= is-2) .and. (i <= ieq+1)) then
581  do j=obc%segment(n)%HI%JsdB,obc%segment(n)%HI%JedB
582  if (obc%zero_strain) then
583  dvdx(i,j) = 0. ; dudy(i,j) = 0.
584  elseif (obc%freeslip_strain) then
585  dvdx(i,j) = 0.
586  elseif (obc%computed_strain) then
587  if (obc%segment(n)%direction == obc_direction_e) then
588  dvdx(i,j) = 2.0*cs%DY_dxBu(i,j)*(obc%segment(n)%tangential_vel(i,j,k) - v(i,j,k))*g%IdyCv(i,j)
589  else
590  dvdx(i,j) = 2.0*cs%DY_dxBu(i,j)*(v(i+1,j,k) - obc%segment(n)%tangential_vel(i,j,k))*g%IdyCv(i+1,j)
591  endif
592  elseif (obc%specified_strain) then
593  if (obc%segment(n)%direction == obc_direction_e) then
594  dvdx(i,j) = cs%DY_dxBu(i,j)*obc%segment(n)%tangential_grad(i,j,k)*g%IdyCv(i,j)*g%dxBu(i,j)
595  else
596  dvdx(i,j) = cs%DY_dxBu(i,j)*obc%segment(n)%tangential_grad(i,j,k)*g%IdyCv(i+1,j)*g%dxBu(i,j)
597  endif
598  endif
599  enddo
600  endif
601  endif
602  if (obc%segment(n)%direction == obc_direction_n) then
603  ! There are extra wide halos here to accommodate the cross-corner-point
604  ! OBC projections, but they might not be necessary if the accelerations
605  ! are always zeroed out at OBC points, in which case the i-loop below
606  ! becomes do i=is-1,ie+1. -RWH
607  if ((j >= jsq-1) .and. (j <= jeq+1)) then
608  do i = max(is-2,obc%segment(n)%HI%isd), min(ie+2,obc%segment(n)%HI%ied)
609  h_v(i,j) = h(i,j,k)
610  enddo
611  endif
612  elseif (obc%segment(n)%direction == obc_direction_s) then
613  if ((j >= jsq-1) .and. (j <= jeq+1)) then
614  do i = max(is-2,obc%segment(n)%HI%isd), min(ie+2,obc%segment(n)%HI%ied)
615  h_v(i,j) = h(i,j+1,k)
616  enddo
617  endif
618  elseif (obc%segment(n)%direction == obc_direction_e) then
619  if ((i >= isq-1) .and. (i <= ieq+1)) then
620  do j = max(js-2,obc%segment(n)%HI%jsd), min(je+2,obc%segment(n)%HI%jed)
621  h_u(i,j) = h(i,j,k)
622  enddo
623  endif
624  elseif (obc%segment(n)%direction == obc_direction_w) then
625  if ((i >= isq-1) .and. (i <= ieq+1)) then
626  do j = max(js-2,obc%segment(n)%HI%jsd), min(je+2,obc%segment(n)%HI%jed)
627  h_u(i,j) = h(i+1,j,k)
628  enddo
629  endif
630  endif
631  enddo ; endif
632  ! Now project thicknesses across corner points on OBCs.
633  if (apply_obc) then ; do n=1,obc%number_of_segments
634  j = obc%segment(n)%HI%JsdB ; i = obc%segment(n)%HI%IsdB
635  if (obc%segment(n)%direction == obc_direction_n) then
636  if ((j >= js-2) .and. (j <= je)) then
637  do i = max(isq-1,obc%segment(n)%HI%IsdB), min(ieq+1,obc%segment(n)%HI%IedB)
638  h_u(i,j+1) = h_u(i,j)
639  enddo
640  endif
641  elseif (obc%segment(n)%direction == obc_direction_s) then
642  if ((j >= js-1) .and. (j <= je+1)) then
643  do i = max(isq-1,obc%segment(n)%HI%isd), min(ieq+1,obc%segment(n)%HI%ied)
644  h_u(i,j) = h_u(i,j+1)
645  enddo
646  endif
647  elseif (obc%segment(n)%direction == obc_direction_e) then
648  if ((i >= is-2) .and. (i <= ie)) then
649  do j = max(jsq-1,obc%segment(n)%HI%jsd), min(jeq+1,obc%segment(n)%HI%jed)
650  h_v(i+1,j) = h_v(i,j)
651  enddo
652  endif
653  elseif (obc%segment(n)%direction == obc_direction_w) then
654  if ((i >= is-1) .and. (i <= ie+1)) then
655  do j = max(jsq-1,obc%segment(n)%HI%jsd), min(jeq+1,obc%segment(n)%HI%jed)
656  h_v(i,j) = h_v(i+1,j)
657  enddo
658  endif
659  endif
660  enddo ; endif
661 
662  ! Shearing strain (including no-slip boundary conditions at the 2-D land-sea mask).
663  ! dudy and dvdx include modifications at OBCs from above.
664  if (cs%no_slip) then
665  do j=js-2,jeq+1 ; do i=is-2,ieq+1
666  sh_xy(i,j) = (2.0-g%mask2dBu(i,j)) * ( dvdx(i,j) + dudy(i,j) )
667  enddo ; enddo
668  else
669  do j=js-2,jeq+1 ; do i=is-2,ieq+1
670  sh_xy(i,j) = g%mask2dBu(i,j) * ( dvdx(i,j) + dudy(i,j) )
671  enddo ; enddo
672  endif
673 
674  ! Evaluate u0 = x.Div(Grad u) and v0 = y.Div( Grad u)
675  if (cs%biharmonic) then
676  do j=js-1,jeq+1 ; do i=isq-1,ieq+1
677  u0(i,j) = cs%IDXDY2u(i,j)*(cs%DY2h(i+1,j)*sh_xx(i+1,j) - cs%DY2h(i,j)*sh_xx(i,j)) + &
678  cs%IDX2dyCu(i,j)*(cs%DX2q(i,j)*sh_xy(i,j) - cs%DX2q(i,j-1)*sh_xy(i,j-1))
679  enddo ; enddo
680  do j=jsq-1,jeq+1 ; do i=is-1,ieq+1
681  v0(i,j) = cs%IDXDY2v(i,j)*(cs%DY2q(i,j)*sh_xy(i,j) - cs%DY2q(i-1,j)*sh_xy(i-1,j)) - &
682  cs%IDX2dyCv(i,j)*(cs%DX2h(i,j+1)*sh_xx(i,j+1) - cs%DX2h(i,j)*sh_xx(i,j))
683  enddo ; enddo
684  if (apply_obc) then; if (obc%zero_biharmonic) then
685  do n=1,obc%number_of_segments
686  i = obc%segment(n)%HI%IsdB ; j = obc%segment(n)%HI%JsdB
687  if (obc%segment(n)%is_N_or_S .and. (j >= jsq-1) .and. (j <= jeq+1)) then
688  do i=obc%segment(n)%HI%isd,obc%segment(n)%HI%ied
689  v0(i,j) = 0.
690  enddo
691  elseif (obc%segment(n)%is_E_or_W .and. (i >= isq-1) .and. (i <= ieq+1)) then
692  do j=obc%segment(n)%HI%jsd,obc%segment(n)%HI%jed
693  u0(i,j) = 0.
694  enddo
695  endif
696  enddo
697  endif; endif
698  endif
699 
700  if ((cs%Leith_Kh) .or. (cs%Leith_Ah)) then
701 
702  ! Components for the vertical vorticity
703  ! Note this a simple re-calculation of shearing components using the same discretization.
704  ! We will consider using a circulation based calculation of vorticity later.
705  ! Also note this will need OBC boundary conditions re-applied...
706  do j=js-2,jeq+1 ; do i=is-2,ieq+1
707  dy_dxbu = g%dyBu(i,j) * g%IdxBu(i,j)
708  dvdx(i,j) = dy_dxbu * (v(i+1,j,k) * g%IdyCv(i+1,j) - v(i,j,k) * g%IdyCv(i,j))
709  dx_dybu = g%dxBu(i,j) * g%IdyBu(i,j)
710  dudy(i,j) = dx_dybu * (u(i,j+1,k) * g%IdxCu(i,j+1) - u(i,j,k) * g%IdxCu(i,j))
711  enddo ; enddo
712 
713  ! Vorticity
714  if (cs%no_slip) then
715  do j=js-2,jeq+1 ; do i=is-2,ieq+1
716  vort_xy(i,j) = (2.0-g%mask2dBu(i,j)) * ( dvdx(i,j) - dudy(i,j) )
717  dudy(i,j) = (2.0-g%mask2dBu(i,j)) * dudy(i,j)
718  dvdx(i,j) = (2.0-g%mask2dBu(i,j)) * dvdx(i,j)
719  enddo ; enddo
720  else
721  do j=js-2,jeq+1 ; do i=is-2,ieq+1
722  vort_xy(i,j) = g%mask2dBu(i,j) * ( dvdx(i,j) - dudy(i,j) )
723  dudy(i,j) = g%mask2dBu(i,j) * dudy(i,j)
724  dvdx(i,j) = g%mask2dBu(i,j) * dvdx(i,j)
725  enddo ; enddo
726  endif
727 
728  call pass_var(vort_xy, g%Domain, position=corner, complete=.true.)
729 
730  ! Vorticity gradient
731  do j=js-2,jeq+1 ; do i=is-1,ieq+1
732  dy_dxbu = g%dyBu(i,j) * g%IdxBu(i,j)
733  vort_xy_dx(i,j) = dy_dxbu * (vort_xy(i,j) * g%IdyCu(i,j) - vort_xy(i-1,j) * g%IdyCu(i-1,j))
734  enddo ; enddo
735 
736  do j=js-1,jeq+1 ; do i=is-2,ieq+1
737  dx_dybu = g%dxBu(i,j) * g%IdyBu(i,j)
738  vort_xy_dy(i,j) = dx_dybu * (vort_xy(i,j) * g%IdxCv(i,j) - vort_xy(i,j-1) * g%IdxCv(i,j-1))
739  enddo ; enddo
740 
741  call pass_vector(vort_xy_dy, vort_xy_dx, g%Domain)
742 
743  if (cs%modified_Leith) then
744  ! Divergence
745  do j=jsq-1,jeq+2 ; do i=isq-1,ieq+2
746  div_xx(i,j) = 0.5*((g%dyCu(i,j) * u(i,j,k) * (h(i+1,j,k)+h(i,j,k)) - &
747  g%dyCu(i-1,j) * u(i-1,j,k) * (h(i-1,j,k)+h(i,j,k)) ) + &
748  (g%dxCv(i,j) * v(i,j,k) * (h(i,j,k)+h(i,j+1,k)) - &
749  g%dxCv(i,j-1)*v(i,j-1,k)*(h(i,j,k)+h(i,j-1,k))))*g%IareaT(i,j)/ &
750  (h(i,j,k) + gv%H_subroundoff)
751  enddo ; enddo
752 
753  !#GME# Adding so many halo updates will make this code very slow!
754  !#GME# With the correct index range, this halo update is unnecessary.
755  call pass_var(div_xx, g%Domain, complete=.true.)
756 
757  ! Divergence gradient
758  !#GME# This index range should be: do j=Jsq,Jeq+1 ; do I=Isq-1,Ieq+1
759  do j=jsq-1,jeq+2 ; do i=is-2,ieq+1
760  div_xx_dx(i,j) = g%IdxCu(i,j)*(div_xx(i+1,j) - div_xx(i,j))
761  enddo ; enddo
762  !#GME# This index range should be: do j=Jsq-1,Jeq+1 ; do i=Isq,Ieq+1
763  do j=js-2,jeq+1 ; do i=isq-1,ieq+2
764  div_xx_dy(i,j) = g%IdyCv(i,j)*(div_xx(i,j+1) - div_xx(i,j))
765  enddo ; enddo
766 
767  !#GME# With the correct index ranges, this halo update is unnecessary.
768  call pass_vector(div_xx_dx, div_xx_dy, g%Domain)
769 
770  ! Magnitude of divergence gradient
771  ! Why use the magnitude of the average instead of the average magnitude?
772  !#GME# This index range should be: do j=Jsq,Jeq+1 ; do i=Isq,Ieq+1
773  do j=jsq-1,jeq+2 ; do i=isq-1,ieq+2
774  grad_div_mag_h(i,j) = sqrt((0.5*(div_xx_dx(i,j) + div_xx_dx(i-1,j)))**2 + &
775  (0.5*(div_xx_dy(i,j) + div_xx_dy(i,j-1)))**2)
776  enddo ; enddo
777  !#GME# This index range should be: do J=js-1,Jeq ; do I=is-1,Ieq
778  do j=js-2,jeq+1 ; do i=is-2,ieq+1
779  grad_div_mag_q(i,j) = sqrt((0.5*(div_xx_dx(i,j) + div_xx_dx(i,j+1)))**2 + &
780  (0.5*(div_xx_dy(i,j) + div_xx_dy(i+1,j)))**2)
781  enddo ; enddo
782 
783  else
784 
785  do j=jsq-1,jeq+2 ; do i=is-2,ieq+1
786  div_xx_dx(i,j) = 0.0
787  enddo ; enddo
788  do j=js-2,jeq+1 ; do i=isq-1,ieq+2
789  div_xx_dy(i,j) = 0.0
790  enddo ; enddo
791  !#GME# This index range should be: do j=Jsq,Jeq+1 ; do i=Isq,Ieq+1
792  do j=jsq-1,jeq+2 ; do i=isq-1,ieq+2
793  grad_div_mag_h(i,j) = 0.0
794  enddo ; enddo
795  !#GME# This index range should be: do J=js-1,Jeq ; do I=is-1,Ieq
796  do j=js-2,jeq+1 ; do i=is-2,ieq+1
797  grad_div_mag_q(i,j) = 0.0
798  enddo ; enddo
799 
800  endif ! CS%modified_Leith
801 
802  ! Add in beta for the Leith viscosity
803  if (cs%use_beta_in_Leith) then
804  !#GME# beta_h and beta_q are never used.
805  do j=jsq-1,jeq+2 ; do i=isq-1,ieq+2
806  beta_h(i,j) = sqrt( g%dF_dx(i,j)**2 + g%dF_dy(i,j)**2 )
807  enddo; enddo
808  do j=js-2,jeq+1 ; do i=is-2,ieq+1
809  beta_q(i,j) = sqrt( (0.25*(g%dF_dx(i,j)+g%dF_dx(i+1,j)+g%dF_dx(i,j+1)+g%dF_dx(i+1,j+1))**2) + &
810  (0.25*(g%dF_dy(i,j)+g%dF_dy(i+1,j)+g%dF_dy(i,j+1)+g%dF_dy(i+1,j+1))**2) )
811  enddo ; enddo
812 
813  do j=js-2,jeq+1 ; do i=is-1,ieq+1
814  vort_xy_dx(i,j) = vort_xy_dx(i,j) + 0.5 * ( g%dF_dx(i,j) + g%dF_dx(i,j+1))
815  enddo ; enddo
816  do j=js-1,jeq+1 ; do i=is-2,ieq+1
817  vort_xy_dy(i,j) = vort_xy_dy(i,j) + 0.5 * ( g%dF_dy(i,j) + g%dF_dy(i+1,j))
818  enddo ; enddo
819  endif ! CS%use_beta_in_Leith
820 
821  if (cs%use_QG_Leith_visc) then
822 
823  !#GME# This should be do j=Jsq,Jeq+1 ; do i=Isq,Ieq+1
824  do j=jsq-1,jeq+2 ; do i=isq-1,ieq+2
825  grad_vort_mag_h_2d(i,j) = sqrt((0.5*(vort_xy_dx(i,j) + vort_xy_dx(i,j-1)))**2 + &
826  (0.5*(vort_xy_dy(i,j) + vort_xy_dy(i-1,j)))**2 )
827  enddo ; enddo
828  !#GME# This index range should be: do J=js-1,Jeq ; do I=is-1,Ieq
829  do j=js-2,jeq+1 ; do i=is-2,ieq+1
830  grad_vort_mag_q_2d(i,j) = sqrt((0.5*(vort_xy_dx(i,j) + vort_xy_dx(i+1,j)))**2 + &
831  (0.5*(vort_xy_dy(i,j) + vort_xy_dy(i,j+1)))**2 )
832  enddo ; enddo
833 
834  call calc_qg_leith_viscosity(varmix, g, gv, h, k, div_xx_dx, div_xx_dy, &
835  vort_xy_dx, vort_xy_dy)
836 
837  endif
838 
839  !#GME# This should be do j=Jsq,Jeq+1 ; do i=Isq,Ieq+1
840  do j=jsq-1,jeq+2 ; do i=isq-1,ieq+2
841  grad_vort_mag_h(i,j) = sqrt((0.5*(vort_xy_dx(i,j) + vort_xy_dx(i,j-1)))**2 + &
842  (0.5*(vort_xy_dy(i,j) + vort_xy_dy(i-1,j)))**2 )
843  enddo ; enddo
844  !#GME# This index range should be: do J=js-1,Jeq ; do I=is-1,Ieq
845  do j=js-2,jeq+1 ; do i=is-2,ieq+1
846  grad_vort_mag_q(i,j) = sqrt((0.5*(vort_xy_dx(i,j) + vort_xy_dx(i+1,j)))**2 + &
847  (0.5*(vort_xy_dy(i,j) + vort_xy_dy(i,j+1)))**2 )
848  enddo ; enddo
849 
850  endif ! CS%Leith_Kh
851 
852  meke_res_fn = 1.
853 
854  do j=jsq,jeq+1 ; do i=isq,ieq+1
855  if ((cs%Smagorinsky_Kh) .or. (cs%Smagorinsky_Ah)) then
856  shear_mag = us%T_to_s * sqrt(sh_xx(i,j)*sh_xx(i,j) + &
857  0.25*((sh_xy(i-1,j-1)*sh_xy(i-1,j-1) + sh_xy(i,j)*sh_xy(i,j)) + &
858  (sh_xy(i-1,j)*sh_xy(i-1,j) + sh_xy(i,j-1)*sh_xy(i,j-1))))
859  endif
860  if ((cs%Leith_Kh) .or. (cs%Leith_Ah)) then
861  if (cs%use_QG_Leith_visc) then
862  vert_vort_mag = us%T_to_s*min(grad_vort_mag_h(i,j) + grad_div_mag_h(i,j),3*grad_vort_mag_h_2d(i,j))
863  else
864  vert_vort_mag = us%T_to_s*(grad_vort_mag_h(i,j) + grad_div_mag_h(i,j))
865  endif
866  endif
867  if (cs%better_bound_Ah .or. cs%better_bound_Kh) then
868  hrat_min = min(1.0, min(h_u(i,j), h_u(i-1,j), h_v(i,j), h_v(i,j-1)) / &
869  (h(i,j,k) + h_neglect) )
870  visc_bound_rem = 1.0
871  endif
872 
873  if (cs%Laplacian) then
874  ! Determine the Laplacian viscosity at h points, using the
875  ! largest value from several parameterizations.
876  kh = cs%Kh_bg_xx(i,j) ! Static (pre-computed) background viscosity
877  if (cs%Smagorinsky_Kh) kh = max( kh, cs%Laplac2_const_xx(i,j) * shear_mag )
878  if (cs%Leith_Kh) kh = max( kh, cs%Laplac3_const_xx(i,j) * vert_vort_mag*inv_pi3)
879  ! All viscosity contributions above are subject to resolution scaling
880  if (rescale_kh) kh = varmix%Res_fn_h(i,j) * kh
881  if (cs%res_scale_MEKE) meke_res_fn = varmix%Res_fn_h(i,j)
882  ! Older method of bounding for stability
883  if (legacy_bound) kh = min(kh, cs%Kh_Max_xx(i,j))
884  kh = max( kh, cs%Kh_bg_min ) ! Place a floor on the viscosity, if desired.
885  if (use_meke_ku) &
886  kh = kh + meke%Ku(i,j) * meke_res_fn ! *Add* the MEKE contribution (might be negative)
887  if (cs%anisotropic) kh = kh + cs%Kh_aniso * ( 1. - cs%n1n2_h(i,j)**2 ) ! *Add* the tension component
888  ! of anisotropic viscosity
889 
890  ! Newer method of bounding for stability
891  if (cs%better_bound_Kh) then
892  if (kh >= hrat_min*cs%Kh_Max_xx(i,j)) then
893  visc_bound_rem = 0.0
894  kh = hrat_min*cs%Kh_Max_xx(i,j)
895  else
896  !visc_bound_rem = 1.0 - abs(Kh) / (hrat_min*CS%Kh_Max_xx(i,j))
897  visc_bound_rem = 1.0 - kh / (hrat_min*cs%Kh_Max_xx(i,j))
898  endif
899  endif
900 
901  if ((cs%id_Kh_h>0) .or. find_frictwork) kh_h(i,j,k) = kh
902  if (cs%id_div_xx_h>0) div_xx_h(i,j,k) = div_xx(i,j)
903 
904  str_xx(i,j) = -kh * sh_xx(i,j)
905  else ! not Laplacian
906  str_xx(i,j) = 0.0
907  endif ! Laplacian
908 
909  if (cs%anisotropic) then
910  ! Shearing-strain averaged to h-points
911  local_strain = 0.25 * ( (sh_xy(i,j) + sh_xy(i-1,j-1)) + (sh_xy(i-1,j) + sh_xy(i,j-1)) )
912  ! *Add* the shear-strain contribution to the xx-component of stress
913  str_xx(i,j) = str_xx(i,j) - cs%Kh_aniso * cs%n1n2_h(i,j) * cs%n1n1_m_n2n2_h(i,j) * local_strain
914  endif
915 
916  if (cs%biharmonic) then
917  ! Determine the biharmonic viscosity at h points, using the
918  ! largest value from several parameterizations.
919  ahsm = 0.0; ahlth = 0.0
920  if ((cs%Smagorinsky_Ah) .or. (cs%Leith_Ah)) then
921  if (cs%Smagorinsky_Ah) then
922  if (cs%bound_Coriolis) then
923  ahsm = shear_mag * (cs%Biharm_const_xx(i,j) + &
924  cs%Biharm_const2_xx(i,j)*shear_mag)
925  else
926  ahsm = cs%Biharm_const_xx(i,j) * shear_mag
927  endif
928  endif
929  if (cs%Leith_Ah) ahlth = cs%biharm5_const_xx(i,j) * vert_vort_mag * inv_pi5
930  ah = max(max(cs%Ah_bg_xx(i,j), ahsm), ahlth)
931  if (cs%bound_Ah .and. .not.cs%better_bound_Ah) &
932  ah = min(ah, cs%Ah_Max_xx(i,j))
933  else
934  ah = cs%Ah_bg_xx(i,j)
935  endif ! Smagorinsky_Ah or Leith_Ah
936 
937  if (use_meke_au) ah = ah + meke%Au(i,j) ! *Add* the MEKE contribution
938 
939  if (cs%better_bound_Ah) then
940  ah = min(ah, visc_bound_rem*hrat_min*cs%Ah_Max_xx(i,j))
941  endif
942 
943  if ((cs%id_Ah_h>0) .or. find_frictwork) ah_h(i,j,k) = ah
944 
945  str_xx(i,j) = str_xx(i,j) + ah * &
946  (cs%DY_dxT(i,j)*(g%IdyCu(i,j)*u0(i,j) - g%IdyCu(i-1,j)*u0(i-1,j)) - &
947  cs%DX_dyT(i,j) *(g%IdxCv(i,j)*v0(i,j) - g%IdxCv(i,j-1)*v0(i,j-1)))
948 
949  ! Keep a copy of the biharmonic contribution for backscatter parameterization
950  bhstr_xx(i,j) = ah * &
951  (cs%DY_dxT(i,j)*(g%IdyCu(i,j)*u0(i,j) - g%IdyCu(i-1,j)*u0(i-1,j)) - &
952  cs%DX_dyT(i,j) *(g%IdxCv(i,j)*v0(i,j) - g%IdxCv(i,j-1)*v0(i,j-1)))
953  bhstr_xx(i,j) = bhstr_xx(i,j) * (h(i,j,k) * cs%reduction_xx(i,j))
954 
955  endif ! biharmonic
956 
957  enddo ; enddo
958 
959  if (cs%biharmonic) then
960  ! Gradient of Laplacian, for use in bi-harmonic term
961  do j=js-1,jeq ; do i=is-1,ieq
962  dvdx(i,j) = cs%DY_dxBu(i,j)*(v0(i+1,j)*g%IdyCv(i+1,j) - v0(i,j)*g%IdyCv(i,j))
963  dudy(i,j) = cs%DX_dyBu(i,j)*(u0(i,j+1)*g%IdxCu(i,j+1) - u0(i,j)*g%IdxCu(i,j))
964  enddo ; enddo
965  ! Adjust contributions to shearing strain on open boundaries.
966  if (apply_obc) then ; if (obc%zero_strain .or. obc%freeslip_strain) then
967  do n=1,obc%number_of_segments
968  j = obc%segment(n)%HI%JsdB ; i = obc%segment(n)%HI%IsdB
969  if (obc%segment(n)%is_N_or_S .and. (j >= js-1) .and. (j <= jeq)) then
970  do i=obc%segment(n)%HI%IsdB,obc%segment(n)%HI%IedB
971  if (obc%zero_strain) then
972  dvdx(i,j) = 0. ; dudy(i,j) = 0.
973  elseif (obc%freeslip_strain) then
974  dudy(i,j) = 0.
975  endif
976  enddo
977  elseif (obc%segment(n)%is_E_or_W .and. (i >= is-1) .and. (i <= ieq)) then
978  do j=obc%segment(n)%HI%JsdB,obc%segment(n)%HI%JedB
979  if (obc%zero_strain) then
980  dvdx(i,j) = 0. ; dudy(i,j) = 0.
981  elseif (obc%freeslip_strain) then
982  dvdx(i,j) = 0.
983  endif
984  enddo
985  endif
986  enddo
987  endif ; endif
988  endif
989 
990  meke_res_fn = 1.
991 
992  do j=js-1,jeq ; do i=is-1,ieq
993  if ((cs%Smagorinsky_Kh) .or. (cs%Smagorinsky_Ah)) then
994  shear_mag = us%T_to_s * sqrt(sh_xy(i,j)*sh_xy(i,j) + &
995  0.25*((sh_xx(i,j)*sh_xx(i,j) + sh_xx(i+1,j+1)*sh_xx(i+1,j+1)) + &
996  (sh_xx(i,j+1)*sh_xx(i,j+1) + sh_xx(i+1,j)*sh_xx(i+1,j))))
997  endif
998  if ((cs%Leith_Kh) .or. (cs%Leith_Ah)) then
999  if (cs%use_QG_Leith_visc) then
1000  vert_vort_mag = us%T_to_s*min(grad_vort_mag_q(i,j) + grad_div_mag_q(i,j), 3*grad_vort_mag_q_2d(i,j))
1001  else
1002  vert_vort_mag = us%T_to_s*(grad_vort_mag_q(i,j) + grad_div_mag_q(i,j))
1003  endif
1004  endif
1005  h2uq = 4.0 * h_u(i,j) * h_u(i,j+1)
1006  h2vq = 4.0 * h_v(i,j) * h_v(i+1,j)
1007  hq(i,j) = 2.0 * h2uq * h2vq / (h_neglect3 + (h2uq + h2vq) * &
1008  ((h_u(i,j) + h_u(i,j+1)) + (h_v(i,j) + h_v(i+1,j))))
1009 
1010  if (cs%better_bound_Ah .or. cs%better_bound_Kh) then
1011  hrat_min = min(1.0, min(h_u(i,j), h_u(i,j+1), h_v(i,j), h_v(i+1,j)) / &
1012  (hq(i,j) + h_neglect) )
1013  visc_bound_rem = 1.0
1014  endif
1015 
1016  if (cs%no_slip .and. (g%mask2dBu(i,j) < 0.5)) then
1017  if ((g%mask2dCu(i,j) + g%mask2dCu(i,j+1)) + &
1018  (g%mask2dCv(i,j) + g%mask2dCv(i+1,j)) > 0.0) then
1019  ! This is a coastal vorticity point, so modify hq and hrat_min.
1020 
1021  hu = g%mask2dCu(i,j) * h_u(i,j) + g%mask2dCu(i,j+1) * h_u(i,j+1)
1022  hv = g%mask2dCv(i,j) * h_v(i,j) + g%mask2dCv(i+1,j) * h_v(i+1,j)
1023  if ((g%mask2dCu(i,j) + g%mask2dCu(i,j+1)) * &
1024  (g%mask2dCv(i,j) + g%mask2dCv(i+1,j)) == 0.0) then
1025  ! Only one of hu and hv is nonzero, so just add them.
1026  hq(i,j) = hu + hv
1027  hrat_min = 1.0
1028  else
1029  ! Both hu and hv are nonzero, so take the harmonic mean.
1030  hq(i,j) = 2.0 * (hu * hv) / ((hu + hv) + h_neglect)
1031  hrat_min = min(1.0, min(hu, hv) / (hq(i,j) + h_neglect) )
1032  endif
1033  endif
1034  endif
1035 
1036  if (cs%Laplacian) then
1037  ! Determine the Laplacian viscosity at q points, using the
1038  ! largest value from several parameterizations.
1039  kh = cs%Kh_bg_xy(i,j) ! Static (pre-computed) background viscosity
1040  if (cs%Smagorinsky_Kh) kh = max( kh, cs%Laplac2_const_xy(i,j) * shear_mag )
1041  if (cs%Leith_Kh) kh = max( kh, cs%Laplac3_const_xy(i,j) * vert_vort_mag*inv_pi3)
1042  ! All viscosity contributions above are subject to resolution scaling
1043  if (rescale_kh) kh = varmix%Res_fn_q(i,j) * kh
1044  if (cs%res_scale_MEKE) meke_res_fn = varmix%Res_fn_q(i,j)
1045  ! Older method of bounding for stability
1046  if (legacy_bound) kh = min(kh, cs%Kh_Max_xy(i,j))
1047  kh = max( kh, cs%Kh_bg_min ) ! Place a floor on the viscosity, if desired.
1048  if (use_meke_ku) then ! *Add* the MEKE contribution (might be negative)
1049  kh = kh + 0.25*( (meke%Ku(i,j) + meke%Ku(i+1,j+1)) + &
1050  (meke%Ku(i+1,j) + meke%Ku(i,j+1)) ) * meke_res_fn
1051  endif
1052  ! Older method of bounding for stability
1053  if (cs%anisotropic) kh = kh + cs%Kh_aniso * cs%n1n2_q(i,j)**2 ! *Add* the shear component
1054  ! of anisotropic viscosity
1055 
1056  ! Newer method of bounding for stability
1057  if (cs%better_bound_Kh) then
1058  if (kh >= hrat_min*cs%Kh_Max_xy(i,j)) then
1059  visc_bound_rem = 0.0
1060  kh = hrat_min*cs%Kh_Max_xy(i,j)
1061  elseif (cs%Kh_Max_xy(i,j)>0.) then
1062  !visc_bound_rem = 1.0 - abs(Kh) / (hrat_min*CS%Kh_Max_xy(I,J))
1063  visc_bound_rem = 1.0 - kh / (hrat_min*cs%Kh_Max_xy(i,j))
1064  endif
1065  endif
1066 
1067  if (cs%id_Kh_q>0) kh_q(i,j,k) = kh
1068  if (cs%id_vort_xy_q>0) vort_xy_q(i,j,k) = vort_xy(i,j)
1069 
1070  str_xy(i,j) = -kh * sh_xy(i,j)
1071  else ! not Laplacian
1072  str_xy(i,j) = 0.0
1073  endif ! Laplacian
1074 
1075  if (cs%anisotropic) then
1076  ! Horizontal-tension averaged to q-points
1077  local_strain = 0.25 * ( (sh_xx(i,j) + sh_xx(i+1,j+1)) + (sh_xx(i+1,j) + sh_xx(i,j+1)) )
1078  ! *Add* the tension contribution to the xy-component of stress
1079  str_xy(i,j) = str_xy(i,j) - cs%Kh_aniso * cs%n1n2_q(i,j) * cs%n1n1_m_n2n2_q(i,j) * local_strain
1080  endif
1081 
1082  if (cs%biharmonic) then
1083  ! Determine the biharmonic viscosity at q points, using the
1084  ! largest value from several parameterizations.
1085  ahsm = 0.0 ; ahlth = 0.0
1086  if (cs%Smagorinsky_Ah .or. cs%Leith_Ah) then
1087  if (cs%Smagorinsky_Ah) then
1088  if (cs%bound_Coriolis) then
1089  ahsm = shear_mag * (cs%Biharm_const_xy(i,j) + &
1090  cs%Biharm_const2_xy(i,j)*shear_mag)
1091  else
1092  ahsm = cs%Biharm_const_xy(i,j) * shear_mag
1093  endif
1094  endif
1095  if (cs%Leith_Ah) ahlth = cs%Biharm5_const_xy(i,j) * vert_vort_mag * inv_pi5
1096  ah = max(max(cs%Ah_bg_xy(i,j), ahsm), ahlth)
1097  if (cs%bound_Ah .and. .not.cs%better_bound_Ah) &
1098  ah = min(ah, cs%Ah_Max_xy(i,j))
1099  else
1100  ah = cs%Ah_bg_xy(i,j)
1101  endif ! Smagorinsky_Ah or Leith_Ah
1102 
1103  if (use_meke_au) then ! *Add* the MEKE contribution
1104  ah = ah + 0.25*( (meke%Au(i,j) + meke%Au(i+1,j+1)) + &
1105  (meke%Au(i+1,j) + meke%Au(i,j+1)) )
1106  endif
1107 
1108  if (cs%better_bound_Ah) then
1109  ah = min(ah, visc_bound_rem*hrat_min*cs%Ah_Max_xy(i,j))
1110  endif
1111 
1112  if (cs%id_Ah_q>0) ah_q(i,j,k) = ah
1113 
1114  str_xy(i,j) = str_xy(i,j) + ah * ( dvdx(i,j) + dudy(i,j) )
1115 
1116  ! Keep a copy of the biharmonic contribution for backscatter parameterization
1117  bhstr_xy(i,j) = ah * ( dvdx(i,j) + dudy(i,j) ) * &
1118  (hq(i,j) * g%mask2dBu(i,j) * cs%reduction_xy(i,j))
1119 
1120  endif ! biharmonic
1121 
1122  enddo ; enddo
1123 
1124 
1125  if (find_frictwork) then
1126  if (cs%Laplacian) then
1127  if (cs%answers_2018) then
1128  do j=js,je ; do i=is,ie
1129  grad_vel_mag_h(i,j) = boundary_mask(i,j) * (dudx(i,j)**2 + dvdy(i,j)**2 + &
1130  (0.25*(dvdx(i,j)+dvdx(i-1,j)+dvdx(i,j-1)+dvdx(i-1,j-1)) )**2 + &
1131  (0.25*(dudy(i,j)+dudy(i-1,j)+dudy(i,j-1)+dudy(i-1,j-1)) )**2)
1132  enddo ; enddo
1133  else
1134  do j=js,je ; do i=is,ie
1135  grad_vel_mag_h(i,j) = boundary_mask(i,j) * ((dudx(i,j)**2 + dvdy(i,j)**2) + &
1136  ((0.25*((dvdx(i,j) + dvdx(i-1,j-1)) + (dvdx(i-1,j) + dvdx(i,j-1))) )**2 + &
1137  (0.25*((dudy(i,j) + dudy(i-1,j-1)) + (dudy(i-1,j) + dudy(i,j-1))) )**2))
1138  enddo ; enddo
1139  endif
1140  else
1141  do j=js,je ; do i=is,ie
1142  grad_vel_mag_h(i,j) = 0.0
1143  enddo ; enddo
1144  endif
1145 
1146  if (cs%biharmonic) then
1147  do j=js,je ; do i=is,ie
1148  grad_d2vel_mag_h(i,j) = boundary_mask(i,j) * ((0.5*(u0(i,j) + u0(i-1,j)))**2 + &
1149  (0.5*(v0(i,j) + v0(i,j-1)))**2)
1150  enddo ; enddo
1151  else
1152  do j=js,je ; do i=is,ie
1153  grad_d2vel_mag_h(i,j) = 0.0
1154  enddo ; enddo
1155  endif
1156 
1157  do j=js,je ; do i=is,ie
1158  ! Diagnose -Kh * |del u|^2 - Ah * |del^2 u|^2
1159  diss_rate(i,j,k) = -us%s_to_T*kh_h(i,j,k) * grad_vel_mag_h(i,j) - &
1160  us%s_to_T*ah_h(i,j,k) * grad_d2vel_mag_h(i,j)
1161 
1162  if (associated(meke)) then ; if (associated(meke%mom_src)) then
1163  ! This is the maximum possible amount of energy that can be converted
1164  ! per unit time, according to theory (multiplied by h)
1165  max_diss_rate(i,j,k) = 2.0*meke%MEKE(i,j) * sqrt(grad_vel_mag_h(i,j))
1166  frictwork_diss(i,j,k) = diss_rate(i,j,k) * h(i,j,k) * gv%H_to_kg_m2
1167  frictworkmax(i,j,k) = -max_diss_rate(i,j,k) * h(i,j,k) * gv%H_to_kg_m2
1168 
1169  ! Determine how much work GME needs to do to reach the "target" ratio between
1170  ! the amount of work actually done and the maximum allowed by theory. Note that
1171  ! we need to add the FrictWork done by the dissipation operators, since this work
1172  ! is done only for numerical stability and is therefore spurious
1173  if (cs%use_GME) then
1174  target_diss_rate_gme(i,j,k) = fwfrac * max_diss_rate(i,j,k) - diss_rate(i,j,k)
1175  endif
1176 
1177  else
1178  frictwork_diss(i,j,k) = diss_rate(i,j,k) * h(i,j,k) * gv%H_to_kg_m2
1179  endif ; endif
1180 
1181  enddo ; enddo
1182  endif
1183 
1184  if (cs%use_GME) then
1185  if (.not. (associated(meke))) call mom_error(fatal, "MEKE must be enabled for GME to be used.")
1186 
1187  do j=jsq,jeq+1 ; do i=isq,ieq+1
1188  gme_coeff = 0.0
1189  if ((max_diss_rate(i,j,k) > 0) .and. (grad_vel_mag_bt_h(i,j)>0) ) then
1190  gme_coeff = fwfrac*us%T_to_s*max_diss_rate(i,j,k) / grad_vel_mag_bt_h(i,j)
1191 ! GME_coeff = FWfrac*US%T_to_s*target_diss_rate_GME(i,j,k) / grad_vel_mag_bt_h(i,j)
1192 
1193  if ((g%bathyT(i,j) < h0_gme) .and. (h0_gme > 0.0)) &
1194  gme_coeff = (g%bathyT(i,j) / h0_gme)**2 * gme_coeff
1195 
1196  ! apply mask and limiter
1197  gme_coeff = min(gme_coeff * boundary_mask(i,j), gme_coeff_limiter)
1198  endif
1199 
1200  if ((cs%id_GME_coeff_h>0) .or. find_frictwork) gme_coeff_h(i,j,k) = gme_coeff
1201 
1202  str_xx_gme(i,j) = gme_coeff * sh_xx_bt(i,j)
1203 
1204  enddo ; enddo
1205 
1206  do j=js-1,jeq ; do i=is-1,ieq
1207  gme_coeff = 0.0
1208  if ((max_diss_rate(i,j,k) > 0) .and. (grad_vel_mag_bt_q(i,j)>0) ) then
1209  !#GME# target_diss_rate_GME and max_diss_rate are defined at h points, not q points as used here.
1210  gme_coeff = fwfrac*us%T_to_s*max_diss_rate(i,j,k) / grad_vel_mag_bt_q(i,j)
1211 ! GME_coeff = FWfrac*US%T_to_s*target_diss_rate_GME(i,j,k) / grad_vel_mag_bt_q(I,J)
1212  if ((g%bathyT(i,j) < h0_gme) .and. (h0_gme > 0.0)) &
1213  gme_coeff = (g%bathyT(i,j) / h0_gme)**2 * gme_coeff
1214 
1215  !#GME# boundary_mask is defined at h points, not q points as used here.
1216  ! apply mask and limiter
1217  gme_coeff = min(gme_coeff * boundary_mask(i,j), gme_coeff_limiter)
1218  endif
1219 
1220  if (cs%id_GME_coeff_q>0) gme_coeff_q(i,j,k) = gme_coeff
1221  str_xy_gme(i,j) = gme_coeff * sh_xy_bt(i,j)
1222 
1223  enddo ; enddo
1224 
1225  ! applying GME diagonal term
1226  call smooth_gme(cs,g,gme_flux_h=str_xx_gme)
1227  call smooth_gme(cs,g,gme_flux_q=str_xy_gme)
1228 
1229  do j=jsq,jeq+1 ; do i=isq,ieq+1
1230  str_xx(i,j) = (str_xx(i,j) + str_xx_gme(i,j)) * (h(i,j,k) * cs%reduction_xx(i,j))
1231  enddo ; enddo
1232 
1233  do j=js-1,jeq ; do i=is-1,ieq
1234  ! GME is applied below
1235  if (cs%no_slip) then
1236  str_xy(i,j) = (str_xy(i,j) + str_xy_gme(i,j)) * (hq(i,j) * cs%reduction_xy(i,j))
1237  else
1238  str_xy(i,j) = (str_xy(i,j) + str_xy_gme(i,j)) * (hq(i,j) * g%mask2dBu(i,j) * cs%reduction_xy(i,j))
1239  endif
1240  enddo ; enddo
1241 
1242  if (associated(meke%GME_snk)) then
1243  do j=js,je ; do i=is,ie
1244  frictwork_gme(i,j,k) = us%s_to_T*gme_coeff_h(i,j,k) * h(i,j,k) * gv%H_to_kg_m2 * grad_vel_mag_bt_h(i,j)
1245  enddo ; enddo
1246  endif
1247 
1248  else ! use_GME
1249  do j=jsq,jeq+1 ; do i=isq,ieq+1
1250  str_xx(i,j) = str_xx(i,j) * (h(i,j,k) * cs%reduction_xx(i,j))
1251  enddo ; enddo
1252 
1253  do j=js-1,jeq ; do i=is-1,ieq
1254  if (cs%no_slip) then
1255  str_xy(i,j) = str_xy(i,j) * (hq(i,j) * cs%reduction_xy(i,j))
1256  else
1257  str_xy(i,j) = str_xy(i,j) * (hq(i,j) * g%mask2dBu(i,j) * cs%reduction_xy(i,j))
1258  endif
1259  enddo ; enddo
1260 
1261  endif ! use_GME
1262 
1263  ! Evaluate 1/h x.Div(h Grad u) or the biharmonic equivalent.
1264  do j=js,je ; do i=isq,ieq
1265  diffu(i,j,k) = ((g%IdyCu(i,j)*(cs%DY2h(i,j) *str_xx(i,j) - &
1266  cs%DY2h(i+1,j)*str_xx(i+1,j)) + &
1267  g%IdxCu(i,j)*(cs%DX2q(i,j-1)*str_xy(i,j-1) - &
1268  cs%DX2q(i,j) *str_xy(i,j))) * &
1269  g%IareaCu(i,j)) / (h_u(i,j) + h_neglect)
1270 
1271  enddo ; enddo
1272  if (apply_obc) then
1273  ! This is not the right boundary condition. If all the masking of tendencies are done
1274  ! correctly later then eliminating this block should not change answers.
1275  do n=1,obc%number_of_segments
1276  if (obc%segment(n)%is_E_or_W) then
1277  i = obc%segment(n)%HI%IsdB
1278  do j=obc%segment(n)%HI%jsd,obc%segment(n)%HI%jed
1279  diffu(i,j,k) = 0.
1280  enddo
1281  endif
1282  enddo
1283  endif
1284 
1285  ! Evaluate 1/h y.Div(h Grad u) or the biharmonic equivalent.
1286  do j=jsq,jeq ; do i=is,ie
1287  diffv(i,j,k) = ((g%IdyCv(i,j)*(cs%DY2q(i-1,j)*str_xy(i-1,j) - &
1288  cs%DY2q(i,j) *str_xy(i,j)) - &
1289  g%IdxCv(i,j)*(cs%DX2h(i,j) *str_xx(i,j) - &
1290  cs%DX2h(i,j+1)*str_xx(i,j+1))) * &
1291  g%IareaCv(i,j)) / (h_v(i,j) + h_neglect)
1292  enddo ; enddo
1293  if (apply_obc) then
1294  ! This is not the right boundary condition. If all the masking of tendencies are done
1295  ! correctly later then eliminating this block should not change answers.
1296  do n=1,obc%number_of_segments
1297  if (obc%segment(n)%is_N_or_S) then
1298  j = obc%segment(n)%HI%JsdB
1299  do i=obc%segment(n)%HI%isd,obc%segment(n)%HI%ied
1300  diffv(i,j,k) = 0.
1301  enddo
1302  endif
1303  enddo
1304  endif
1305 
1306  if (find_frictwork) then ; do j=js,je ; do i=is,ie
1307  ! Diagnose str_xx*d_x u - str_yy*d_y v + str_xy*(d_y u + d_x v)
1308  ! This is the old formulation that includes energy diffusion
1309  frictwork(i,j,k) = us%s_to_T*gv%H_to_kg_m2 * ( &
1310  (str_xx(i,j)*(u(i,j,k)-u(i-1,j,k))*g%IdxT(i,j) &
1311  -str_xx(i,j)*(v(i,j,k)-v(i,j-1,k))*g%IdyT(i,j)) &
1312  +0.25*((str_xy(i,j)*( &
1313  (u(i,j+1,k)-u(i,j,k))*g%IdyBu(i,j) &
1314  +(v(i+1,j,k)-v(i,j,k))*g%IdxBu(i,j) ) &
1315  +str_xy(i-1,j-1)*( &
1316  (u(i-1,j,k)-u(i-1,j-1,k))*g%IdyBu(i-1,j-1) &
1317  +(v(i,j-1,k)-v(i-1,j-1,k))*g%IdxBu(i-1,j-1) )) &
1318  +(str_xy(i-1,j)*( &
1319  (u(i-1,j+1,k)-u(i-1,j,k))*g%IdyBu(i-1,j) &
1320  +(v(i,j,k)-v(i-1,j,k))*g%IdxBu(i-1,j) ) &
1321  +str_xy(i,j-1)*( &
1322  (u(i,j,k)-u(i,j-1,k))*g%IdyBu(i,j-1) &
1323  +(v(i+1,j-1,k)-v(i,j-1,k))*g%IdxBu(i,j-1) )) ) )
1324  enddo ; enddo ; endif
1325 
1326  ! Make a similar calculation as for FrictWork above but accumulating into
1327  ! the vertically integrated MEKE source term, and adjusting for any
1328  ! energy loss seen as a reduction in the [biharmonic] frictional source term.
1329  if (find_frictwork .and. associated(meke)) then ; if (associated(meke%mom_src)) then
1330  if (k==1) then
1331  do j=js,je ; do i=is,ie
1332  meke%mom_src(i,j) = 0.
1333  meke%GME_snk(i,j) = 0.
1334  enddo ; enddo
1335  endif
1336  if (meke%backscatter_Ro_c /= 0.) then
1337  do j=js,je ; do i=is,ie
1338  fath = 0.25*( (abs(g%CoriolisBu(i-1,j-1)) + abs(g%CoriolisBu(i,j))) + &
1339  (abs(g%CoriolisBu(i-1,j)) + abs(g%CoriolisBu(i,j-1))) )
1340  shear_mag = us%T_to_s * sqrt(sh_xx(i,j)*sh_xx(i,j) + &
1341  0.25*((sh_xy(i-1,j-1)*sh_xy(i-1,j-1) + sh_xy(i,j)*sh_xy(i,j)) + &
1342  (sh_xy(i-1,j)*sh_xy(i-1,j) + sh_xy(i,j-1)*sh_xy(i,j-1))))
1343  if (cs%answers_2018) then
1344  fath = (us%s_to_T*fath)**meke%backscatter_Ro_pow ! f^n
1345  ! Note the hard-coded dimensional constant in the following line that can not
1346  ! be rescaled for dimensional consistency.
1347  shear_mag = ( ( (us%s_to_T*shear_mag)**meke%backscatter_Ro_pow ) + 1.e-30 ) &
1348  * meke%backscatter_Ro_c ! c * D^n
1349  ! The Rossby number function is g(Ro) = 1/(1+c.Ro^n)
1350  ! RoScl = 1 - g(Ro)
1351  roscl = shear_mag / ( fath + shear_mag ) ! = 1 - f^n/(f^n+c*D^n)
1352  else
1353  if (fath <= backscat_subround*shear_mag) then
1354  roscl = 1.0
1355  else
1356  sh_f_pow = meke%backscatter_Ro_c * (shear_mag / fath)**meke%backscatter_Ro_pow
1357  roscl = sh_f_pow / (1.0 + sh_f_pow) ! = 1 - f^n/(f^n+c*D^n)
1358  endif
1359  endif
1360  meke%mom_src(i,j) = meke%mom_src(i,j) + us%s_to_T*gv%H_to_kg_m2 * ( &
1361  ((str_xx(i,j)-roscl*bhstr_xx(i,j))*(u(i,j,k)-u(i-1,j,k))*g%IdxT(i,j) &
1362  -(str_xx(i,j)-roscl*bhstr_xx(i,j))*(v(i,j,k)-v(i,j-1,k))*g%IdyT(i,j)) &
1363  +0.25*(((str_xy(i,j)-roscl*bhstr_xy(i,j))*( &
1364  (u(i,j+1,k)-u(i,j,k))*g%IdyBu(i,j) &
1365  +(v(i+1,j,k)-v(i,j,k))*g%IdxBu(i,j) ) &
1366  +(str_xy(i-1,j-1)-roscl*bhstr_xy(i-1,j-1))*( &
1367  (u(i-1,j,k)-u(i-1,j-1,k))*g%IdyBu(i-1,j-1) &
1368  +(v(i,j-1,k)-v(i-1,j-1,k))*g%IdxBu(i-1,j-1) )) &
1369  +((str_xy(i-1,j)-roscl*bhstr_xy(i-1,j))*( &
1370  (u(i-1,j+1,k)-u(i-1,j,k))*g%IdyBu(i-1,j) &
1371  +(v(i,j,k)-v(i-1,j,k))*g%IdxBu(i-1,j) ) &
1372  +(str_xy(i,j-1)-roscl*bhstr_xy(i,j-1))*( &
1373  (u(i,j,k)-u(i,j-1,k))*g%IdyBu(i,j-1) &
1374  +(v(i+1,j-1,k)-v(i,j-1,k))*g%IdxBu(i,j-1) )) ) )
1375  enddo ; enddo
1376  else
1377  do j=js,je ; do i=is,ie
1378  ! MEKE%mom_src now is sign definite because it only uses the dissipation
1379  meke%mom_src(i,j) = meke%mom_src(i,j) + max(frictwork_diss(i,j,k), frictworkmax(i,j,k))
1380  enddo ; enddo
1381  endif ! MEKE%backscatter
1382 
1383  if (cs%use_GME .and. associated(meke)) then
1384  if (associated(meke%GME_snk)) then
1385  do j=js,je ; do i=is,ie
1386  meke%GME_snk(i,j) = meke%GME_snk(i,j) + frictwork_gme(i,j,k)
1387  enddo ; enddo
1388  endif
1389  endif
1390 
1391  endif ; endif ! find_FrictWork and associated(mom_src)
1392 
1393  enddo ! end of k loop
1394 
1395  ! Offer fields for diagnostic averaging.
1396  if (cs%id_diffu>0) call post_data(cs%id_diffu, diffu, cs%diag)
1397  if (cs%id_diffv>0) call post_data(cs%id_diffv, diffv, cs%diag)
1398  if (cs%id_FrictWork>0) call post_data(cs%id_FrictWork, frictwork, cs%diag)
1399  if (cs%id_FrictWorkMax>0) call post_data(cs%id_FrictWorkMax, frictworkmax, cs%diag)
1400  if (cs%id_FrictWork_diss>0) call post_data(cs%id_FrictWork_diss, frictwork_diss, cs%diag)
1401  if (cs%id_FrictWork_GME>0) call post_data(cs%id_FrictWork_GME, frictwork_gme, cs%diag)
1402  if (cs%id_Ah_h>0) call post_data(cs%id_Ah_h, ah_h, cs%diag)
1403  if (cs%id_div_xx_h>0) call post_data(cs%id_div_xx_h, div_xx_h, cs%diag)
1404  if (cs%id_vort_xy_q>0) call post_data(cs%id_vort_xy_q, vort_xy_q, cs%diag)
1405  if (cs%id_Ah_q>0) call post_data(cs%id_Ah_q, ah_q, cs%diag)
1406  if (cs%id_Kh_h>0) call post_data(cs%id_Kh_h, kh_h, cs%diag)
1407  if (cs%id_Kh_q>0) call post_data(cs%id_Kh_q, kh_q, cs%diag)
1408  if (cs%id_GME_coeff_h > 0) call post_data(cs%id_GME_coeff_h, gme_coeff_h, cs%diag)
1409  if (cs%id_GME_coeff_q > 0) call post_data(cs%id_GME_coeff_q, gme_coeff_q, cs%diag)
1410 
1411  if (cs%id_FrictWorkIntz > 0) then
1412  do j=js,je
1413  do i=is,ie ; frictworkintz(i,j) = frictwork(i,j,1) ; enddo
1414  do k=2,nz ; do i=is,ie
1415  frictworkintz(i,j) = frictworkintz(i,j) + frictwork(i,j,k)
1416  enddo ; enddo
1417  enddo
1418  call post_data(cs%id_FrictWorkIntz, frictworkintz, cs%diag)
1419  endif
1420 
1421 end subroutine horizontal_viscosity
1422 
1423 !> Allocates space for and calculates static variables used by horizontal_viscosity().
1424 !! hor_visc_init calculates and stores the values of a number of metric functions that
1425 !! are used in horizontal_viscosity().
1426 subroutine hor_visc_init(Time, G, US, param_file, diag, CS, MEKE)
1427  type(time_type), intent(in) :: time !< Current model time.
1428  type(ocean_grid_type), intent(inout) :: g !< The ocean's grid structure.
1429  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
1430  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time
1431  !! parameters.
1432  type(diag_ctrl), target, intent(inout) :: diag !< Structure to regulate diagnostic output.
1433  type(hor_visc_cs), pointer :: cs !< Pointer to the control structure for this module
1434  type(meke_type), pointer :: meke !< MEKE data
1435  ! Local variables
1436  real, dimension(SZIB_(G),SZJ_(G)) :: u0u, u0v
1437  real, dimension(SZI_(G),SZJB_(G)) :: v0u, v0v
1438  ! u0v is the Laplacian sensitivities to the v velocities
1439  ! at u points [m-2], with u0u, v0u, and v0v defined similarly.
1440  real :: grid_sp_h2 ! Harmonic mean of the squares of the grid [m2]
1441  real :: grid_sp_h3 ! Harmonic mean of the squares of the grid^(3/2) [m3]
1442  real :: grid_sp_q2 ! spacings at h and q points [m2]
1443  real :: grid_sp_q3 ! spacings at h and q points^(3/2) [m3]
1444  real :: kh_limit ! A coefficient [T-1 ~> s-1] used, along with the
1445  ! grid spacing, to limit Laplacian viscosity.
1446  real :: fmax ! maximum absolute value of f at the four
1447  ! vorticity points around a thickness point [T-1 ~> s-1]
1448  real :: boundcorconst ! A constant used when using viscosity to bound the Coriolis accelerations
1449  ! [T2 L-2 ~> s2 m-2]
1450  real :: ah_limit ! coefficient [T-1 ~> s-1] used, along with the
1451  ! grid spacing, to limit biharmonic viscosity
1452  real :: kh ! Lapacian horizontal viscosity [m2 s-1]
1453  real :: ah ! biharmonic horizontal viscosity [m4 s-1]
1454  real :: kh_vel_scale ! this speed [m T-1 ~> m s-1] times grid spacing gives Lap visc
1455  real :: ah_vel_scale ! this speed [m T-1 ~> m s-1] times grid spacing cubed gives bih visc
1456  real :: ah_time_scale ! damping time-scale for biharmonic visc [T ~> s]
1457  real :: smag_lap_const ! nondimensional Laplacian Smagorinsky constant
1458  real :: smag_bi_const ! nondimensional biharmonic Smagorinsky constant
1459  real :: leith_lap_const ! nondimensional Laplacian Leith constant
1460  real :: leith_bi_const ! nondimensional biharmonic Leith constant
1461  real :: dt ! The dynamics time step [T ~> s]
1462  real :: idt ! The inverse of dt [T-1 ~> s-1]
1463  real :: denom ! work variable; the denominator of a fraction
1464  real :: maxvel ! largest permitted velocity components [m s-1]
1465  real :: bound_cor_vel ! grid-scale velocity variations at which value
1466  ! the quadratically varying biharmonic viscosity
1467  ! balances Coriolis acceleration [L T-1 ~> m s-1]
1468  real :: kh_sin_lat ! Amplitude of latitudinally dependent viscosity [m2 T-1 ~> m2 s-1]
1469  real :: kh_pwr_of_sine ! Power used to raise sin(lat) when using Kh_sin_lat
1470  logical :: bound_cor_def ! parameter setting of BOUND_CORIOLIS
1471  logical :: get_all ! If true, read and log all parameters, regardless of
1472  ! whether they are used, to enable spell-checking of
1473  ! valid parameters.
1474  logical :: split ! If true, use the split time stepping scheme.
1475  ! If false and USE_GME = True, issue a FATAL error.
1476  logical :: default_2018_answers
1477  character(len=64) :: inputdir, filename
1478  real :: deg2rad ! Converts degrees to radians
1479  real :: slat_fn ! sin(lat)**Kh_pwr_of_sine
1480  real :: aniso_grid_dir(2) ! Vector (n1,n2) for anisotropic direction
1481  integer :: aniso_mode ! Selects the mode for setting the anisotropic direction
1482  integer :: is, ie, js, je, isq, ieq, jsq, jeq, nz
1483  integer :: isd, ied, jsd, jed, isdb, iedb, jsdb, jedb
1484  integer :: i, j
1485 
1486 ! This include declares and sets the variable "version".
1487 #include "version_variable.h"
1488  character(len=40) :: mdl = "MOM_hor_visc" ! module name
1489 
1490  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
1491  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
1492  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
1493  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
1494 
1495  if (associated(cs)) then
1496  call mom_error(warning, "hor_visc_init called with an associated "// &
1497  "control structure.")
1498  return
1499  endif
1500  allocate(cs)
1501 
1502  cs%diag => diag
1503 
1504  ! Read parameters and write them to the model log.
1505  call log_version(param_file, mdl, version, "")
1506 
1507  ! It is not clear whether these initialization lines are needed for the
1508  ! cases where the corresponding parameters are not read.
1509  cs%bound_Kh = .false. ; cs%better_bound_Kh = .false. ; cs%Smagorinsky_Kh = .false. ; cs%Leith_Kh = .false.
1510  cs%bound_Ah = .false. ; cs%better_bound_Ah = .false. ; cs%Smagorinsky_Ah = .false. ; cs%Leith_Ah = .false.
1511  cs%use_QG_Leith_visc = .false.
1512  cs%bound_Coriolis = .false.
1513  cs%Modified_Leith = .false.
1514  cs%anisotropic = .false.
1515  cs%dynamic_aniso = .false.
1516 
1517  kh = 0.0 ; ah = 0.0
1518 
1519  ! If GET_ALL_PARAMS is true, all parameters are read in all cases to enable
1520  ! parameter spelling checks.
1521  call get_param(param_file, mdl, "GET_ALL_PARAMS", get_all, default=.false.)
1522 
1523  call get_param(param_file, mdl, "DEFAULT_2018_ANSWERS", default_2018_answers, &
1524  "This sets the default value for the various _2018_ANSWERS parameters.", &
1525  default=.true.)
1526  call get_param(param_file, mdl, "HOR_VISC_2018_ANSWERS", cs%answers_2018, &
1527  "If true, use the order of arithmetic and expressions that recover the "//&
1528  "answers from the end of 2018. Otherwise, use updated and more robust "//&
1529  "forms of the same expressions.", default=default_2018_answers)
1530  call get_param(param_file, mdl, "LAPLACIAN", cs%Laplacian, &
1531  "If true, use a Laplacian horizontal viscosity.", &
1532  default=.false.)
1533  if (cs%Laplacian .or. get_all) then
1534  call get_param(param_file, mdl, "KH", kh, &
1535  "The background Laplacian horizontal viscosity.", &
1536  units = "m2 s-1", default=0.0, scale=us%T_to_s)
1537  call get_param(param_file, mdl, "KH_BG_MIN", cs%Kh_bg_min, &
1538  "The minimum value allowed for Laplacian horizontal viscosity, KH.", &
1539  units = "m2 s-1", default=0.0, scale=us%T_to_s)
1540  call get_param(param_file, mdl, "KH_VEL_SCALE", kh_vel_scale, &
1541  "The velocity scale which is multiplied by the grid "//&
1542  "spacing to calculate the Laplacian viscosity. "//&
1543  "The final viscosity is the largest of this scaled "//&
1544  "viscosity, the Smagorinsky and Leith viscosities, and KH.", &
1545  units="m s-1", default=0.0, scale=us%T_to_s)
1546  call get_param(param_file, mdl, "KH_SIN_LAT", kh_sin_lat, &
1547  "The amplitude of a latitudinally-dependent background "//&
1548  "viscosity of the form KH_SIN_LAT*(SIN(LAT)**KH_PWR_OF_SINE).", &
1549  units = "m2 s-1", default=0.0, scale=us%T_to_s)
1550  if (kh_sin_lat>0. .or. get_all) &
1551  call get_param(param_file, mdl, "KH_PWR_OF_SINE", kh_pwr_of_sine, &
1552  "The power used to raise SIN(LAT) when using a latitudinally "//&
1553  "dependent background viscosity.", &
1554  units = "nondim", default=4.0)
1555 
1556  call get_param(param_file, mdl, "SMAGORINSKY_KH", cs%Smagorinsky_Kh, &
1557  "If true, use a Smagorinsky nonlinear eddy viscosity.", &
1558  default=.false.)
1559  if (cs%Smagorinsky_Kh .or. get_all) &
1560  call get_param(param_file, mdl, "SMAG_LAP_CONST", smag_lap_const, &
1561  "The nondimensional Laplacian Smagorinsky constant, "//&
1562  "often 0.15.", units="nondim", default=0.0, &
1563  fail_if_missing = cs%Smagorinsky_Kh)
1564 
1565  call get_param(param_file, mdl, "LEITH_KH", cs%Leith_Kh, &
1566  "If true, use a Leith nonlinear eddy viscosity.", &
1567  default=.false.)
1568 
1569  call get_param(param_file, mdl, "MODIFIED_LEITH", cs%Modified_Leith, &
1570  "If true, add a term to Leith viscosity which is "//&
1571  "proportional to the gradient of divergence.", &
1572  default=.false.)
1573  call get_param(param_file, mdl, "RES_SCALE_MEKE_VISC", cs%res_scale_MEKE, &
1574  "If true, the viscosity contribution from MEKE is scaled by "//&
1575  "the resolution function.", default=.false.)
1576 
1577  if (cs%Leith_Kh .or. get_all) then
1578  call get_param(param_file, mdl, "LEITH_LAP_CONST", leith_lap_const, &
1579  "The nondimensional Laplacian Leith constant, "//&
1580  "often set to 1.0", units="nondim", default=0.0, &
1581  fail_if_missing = cs%Leith_Kh)
1582  call get_param(param_file, mdl, "USE_QG_LEITH_VISC", cs%use_QG_Leith_visc, &
1583  "If true, use QG Leith nonlinear eddy viscosity.", &
1584  default=.false.)
1585  if (cs%use_QG_Leith_visc .and. .not. cs%Leith_Kh) call mom_error(fatal, &
1586  "MOM_lateral_mixing_coeffs.F90, VarMix_init:"//&
1587  "LEITH_KH must be True when USE_QG_LEITH_VISC=True.")
1588  endif
1589  if (cs%Leith_Kh .or. cs%Leith_Ah .or. get_all) then
1590  call get_param(param_file, mdl, "USE_BETA_IN_LEITH", cs%use_beta_in_Leith, &
1591  "If true, include the beta term in the Leith nonlinear eddy viscosity.", &
1592  default=cs%Leith_Kh)
1593  call get_param(param_file, mdl, "MODIFIED_LEITH", cs%modified_Leith, &
1594  "If true, add a term to Leith viscosity which is \n"//&
1595  "proportional to the gradient of divergence.", &
1596  default=.false.)
1597  endif
1598  call get_param(param_file, mdl, "BOUND_KH", cs%bound_Kh, &
1599  "If true, the Laplacian coefficient is locally limited "//&
1600  "to be stable.", default=.true.)
1601  call get_param(param_file, mdl, "BETTER_BOUND_KH", cs%better_bound_Kh, &
1602  "If true, the Laplacian coefficient is locally limited "//&
1603  "to be stable with a better bounding than just BOUND_KH.", &
1604  default=cs%bound_Kh)
1605  call get_param(param_file, mdl, "ANISOTROPIC_VISCOSITY", cs%anisotropic, &
1606  "If true, allow anistropic viscosity in the Laplacian "//&
1607  "horizontal viscosity.", default=.false.)
1608  endif
1609  if (cs%anisotropic .or. get_all) then
1610  call get_param(param_file, mdl, "KH_ANISO", cs%Kh_aniso, &
1611  "The background Laplacian anisotropic horizontal viscosity.", &
1612  units = "m2 s-1", default=0.0, scale=us%T_to_s)
1613  call get_param(param_file, mdl, "ANISOTROPIC_MODE", aniso_mode, &
1614  "Selects the mode for setting the direction of anistropy.\n"//&
1615  "\t 0 - Points along the grid i-direction.\n"//&
1616  "\t 1 - Points towards East.\n"//&
1617  "\t 2 - Points along the flow direction, U/|U|.", &
1618  default=0)
1619  select case (aniso_mode)
1620  case (0)
1621  call get_param(param_file, mdl, "ANISO_GRID_DIR", aniso_grid_dir, &
1622  "The vector pointing in the direction of anistropy for "//&
1623  "horizont viscosity. n1,n2 are the i,j components relative "//&
1624  "to the grid.", units = "nondim", fail_if_missing=.true.)
1625  case (1)
1626  call get_param(param_file, mdl, "ANISO_GRID_DIR", aniso_grid_dir, &
1627  "The vector pointing in the direction of anistropy for "//&
1628  "horizont viscosity. n1,n2 are the i,j components relative "//&
1629  "to the spherical coordinates.", units = "nondim", fail_if_missing=.true.)
1630  end select
1631  endif
1632 
1633  call get_param(param_file, mdl, "BIHARMONIC", cs%biharmonic, &
1634  "If true, use a biharmonic horizontal viscosity. "//&
1635  "BIHARMONIC may be used with LAPLACIAN.", &
1636  default=.true.)
1637  if (cs%biharmonic .or. get_all) then
1638  call get_param(param_file, mdl, "AH", ah, &
1639  "The background biharmonic horizontal viscosity.", &
1640  units = "m4 s-1", default=0.0, scale=us%T_to_s)
1641  call get_param(param_file, mdl, "AH_VEL_SCALE", ah_vel_scale, &
1642  "The velocity scale which is multiplied by the cube of "//&
1643  "the grid spacing to calculate the biharmonic viscosity. "//&
1644  "The final viscosity is the largest of this scaled "//&
1645  "viscosity, the Smagorinsky and Leith viscosities, and AH.", &
1646  units="m s-1", default=0.0, scale=us%T_to_s)
1647  call get_param(param_file, mdl, "AH_TIME_SCALE", ah_time_scale, &
1648  "A time scale whose inverse is multiplied by the fourth "//&
1649  "power of the grid spacing to calculate biharmonic viscosity. "//&
1650  "The final viscosity is the largest of all viscosity "//&
1651  "formulations in use. 0.0 means that it's not used.", &
1652  units="s", default=0.0, scale=us%s_to_T)
1653  call get_param(param_file, mdl, "SMAGORINSKY_AH", cs%Smagorinsky_Ah, &
1654  "If true, use a biharmonic Smagorinsky nonlinear eddy "//&
1655  "viscosity.", default=.false.)
1656  call get_param(param_file, mdl, "LEITH_AH", cs%Leith_Ah, &
1657  "If true, use a biharmonic Leith nonlinear eddy "//&
1658  "viscosity.", default=.false.)
1659 
1660  call get_param(param_file, mdl, "BOUND_AH", cs%bound_Ah, &
1661  "If true, the biharmonic coefficient is locally limited "//&
1662  "to be stable.", default=.true.)
1663  call get_param(param_file, mdl, "BETTER_BOUND_AH", cs%better_bound_Ah, &
1664  "If true, the biharmonic coefficient is locally limited "//&
1665  "to be stable with a better bounding than just BOUND_AH.", &
1666  default=cs%bound_Ah)
1667 
1668  if (cs%Smagorinsky_Ah .or. get_all) then
1669  call get_param(param_file, mdl, "SMAG_BI_CONST",smag_bi_const, &
1670  "The nondimensional biharmonic Smagorinsky constant, "//&
1671  "typically 0.015 - 0.06.", units="nondim", default=0.0, &
1672  fail_if_missing = cs%Smagorinsky_Ah)
1673 
1674  call get_param(param_file, mdl, "BOUND_CORIOLIS", bound_cor_def, default=.false.)
1675  call get_param(param_file, mdl, "BOUND_CORIOLIS_BIHARM", cs%bound_Coriolis, &
1676  "If true use a viscosity that increases with the square "//&
1677  "of the velocity shears, so that the resulting viscous "//&
1678  "drag is of comparable magnitude to the Coriolis terms "//&
1679  "when the velocity differences between adjacent grid "//&
1680  "points is 0.5*BOUND_CORIOLIS_VEL. The default is the "//&
1681  "value of BOUND_CORIOLIS (or false).", default=bound_cor_def)
1682  if (cs%bound_Coriolis .or. get_all) then
1683  call get_param(param_file, mdl, "MAXVEL", maxvel, default=3.0e8)
1684  bound_cor_vel = maxvel
1685  call get_param(param_file, mdl, "BOUND_CORIOLIS_VEL", bound_cor_vel, &
1686  "The velocity scale at which BOUND_CORIOLIS_BIHARM causes "//&
1687  "the biharmonic drag to have comparable magnitude to the "//&
1688  "Coriolis acceleration. The default is set by MAXVEL.", &
1689  units="m s-1", default=maxvel, scale=us%m_s_to_L_T)
1690  endif
1691  endif
1692 
1693  if (cs%Leith_Ah .or. get_all) &
1694  call get_param(param_file, mdl, "LEITH_BI_CONST", leith_bi_const, &
1695  "The nondimensional biharmonic Leith constant, "//&
1696  "typical values are thus far undetermined.", units="nondim", default=0.0, &
1697  fail_if_missing = cs%Leith_Ah)
1698 
1699  endif
1700 
1701  call get_param(param_file, mdl, "USE_LAND_MASK_FOR_HVISC", cs%use_land_mask, &
1702  "If true, use Use the land mask for the computation of thicknesses "//&
1703  "at velocity locations. This eliminates the dependence on arbitrary "//&
1704  "values over land or outside of the domain. Default is False in order to "//&
1705  "maintain answers with legacy experiments but should be changed to True "//&
1706  "for new experiments.", default=.false.)
1707 
1708  if (cs%better_bound_Ah .or. cs%better_bound_Kh .or. get_all) &
1709  call get_param(param_file, mdl, "HORVISC_BOUND_COEF", cs%bound_coef, &
1710  "The nondimensional coefficient of the ratio of the "//&
1711  "viscosity bounds to the theoretical maximum for "//&
1712  "stability without considering other terms.", units="nondim", &
1713  default=0.8)
1714 
1715  call get_param(param_file, mdl, "NOSLIP", cs%no_slip, &
1716  "If true, no slip boundary conditions are used; otherwise "//&
1717  "free slip boundary conditions are assumed. The "//&
1718  "implementation of the free slip BCs on a C-grid is much "//&
1719  "cleaner than the no slip BCs. The use of free slip BCs "//&
1720  "is strongly encouraged, and no slip BCs are not used with "//&
1721  "the biharmonic viscosity.", default=.false.)
1722 
1723  call get_param(param_file, mdl, "USE_KH_BG_2D", cs%use_Kh_bg_2d, &
1724  "If true, read a file containing 2-d background harmonic "//&
1725  "viscosities. The final viscosity is the maximum of the other "//&
1726  "terms and this background value.", default=.false.)
1727 
1728  call get_param(param_file, mdl, "USE_GME", cs%use_GME, &
1729  "If true, use the GM+E backscatter scheme in association \n"//&
1730  "with the Gent and McWilliams parameterization.", default=.false.)
1731 
1732  if (cs%use_GME) then
1733  call get_param(param_file, mdl, "SPLIT", split, &
1734  "Use the split time stepping if true.", default=.true., &
1735  do_not_log=.true.)
1736  if (.not. split) call mom_error(fatal,"ERROR: Currently, USE_GME = True "// &
1737  "cannot be used with SPLIT=False.")
1738  endif
1739 
1740  if (cs%bound_Kh .or. cs%bound_Ah .or. cs%better_bound_Kh .or. cs%better_bound_Ah) &
1741  call get_param(param_file, mdl, "DT", dt, &
1742  "The (baroclinic) dynamics time step.", units="s", scale=us%s_to_T, &
1743  fail_if_missing=.true.)
1744 
1745  if (cs%no_slip .and. cs%biharmonic) &
1746  call mom_error(fatal,"ERROR: NOSLIP and BIHARMONIC cannot be defined "// &
1747  "at the same time in MOM.")
1748 
1749  if (.not.(cs%Laplacian .or. cs%biharmonic)) then
1750  ! Only issue inviscid warning if not in single column mode (usually 2x2 domain)
1751  if ( max(g%domain%niglobal, g%domain%njglobal)>2 ) call mom_error(warning, &
1752  "hor_visc_init: It is usually a very bad idea not to use either "//&
1753  "LAPLACIAN or BIHARMONIC viscosity.")
1754  return ! We are not using either Laplacian or Bi-harmonic lateral viscosity
1755  endif
1756 
1757  deg2rad = atan(1.0) / 45.
1758 
1759  alloc_(cs%dx2h(isd:ied,jsd:jed)) ; cs%dx2h(:,:) = 0.0
1760  alloc_(cs%dy2h(isd:ied,jsd:jed)) ; cs%dy2h(:,:) = 0.0
1761  alloc_(cs%dx2q(isdb:iedb,jsdb:jedb)) ; cs%dx2q(:,:) = 0.0
1762  alloc_(cs%dy2q(isdb:iedb,jsdb:jedb)) ; cs%dy2q(:,:) = 0.0
1763  alloc_(cs%dx_dyT(isd:ied,jsd:jed)) ; cs%dx_dyT(:,:) = 0.0
1764  alloc_(cs%dy_dxT(isd:ied,jsd:jed)) ; cs%dy_dxT(:,:) = 0.0
1765  alloc_(cs%dx_dyBu(isdb:iedb,jsdb:jedb)) ; cs%dx_dyBu(:,:) = 0.0
1766  alloc_(cs%dy_dxBu(isdb:iedb,jsdb:jedb)) ; cs%dy_dxBu(:,:) = 0.0
1767 
1768  if (cs%Laplacian) then
1769  alloc_(cs%Kh_bg_xx(isd:ied,jsd:jed)) ; cs%Kh_bg_xx(:,:) = 0.0
1770  alloc_(cs%Kh_bg_xy(isdb:iedb,jsdb:jedb)) ; cs%Kh_bg_xy(:,:) = 0.0
1771  if (cs%bound_Kh .or. cs%better_bound_Kh) then
1772  alloc_(cs%Kh_Max_xx(isdb:iedb,jsdb:jedb)) ; cs%Kh_Max_xx(:,:) = 0.0
1773  alloc_(cs%Kh_Max_xy(isdb:iedb,jsdb:jedb)) ; cs%Kh_Max_xy(:,:) = 0.0
1774  endif
1775  if (cs%Smagorinsky_Kh) then
1776  alloc_(cs%Laplac2_const_xx(isd:ied,jsd:jed)) ; cs%Laplac2_const_xx(:,:) = 0.0
1777  alloc_(cs%Laplac2_const_xy(isdb:iedb,jsdb:jedb)) ; cs%Laplac2_const_xy(:,:) = 0.0
1778  endif
1779  if (cs%Leith_Kh) then
1780  alloc_(cs%Laplac3_const_xx(isd:ied,jsd:jed)) ; cs%Laplac3_const_xx(:,:) = 0.0
1781  alloc_(cs%Laplac3_const_xy(isdb:iedb,jsdb:jedb)) ; cs%Laplac3_const_xy(:,:) = 0.0
1782  endif
1783  endif
1784  alloc_(cs%reduction_xx(isd:ied,jsd:jed)) ; cs%reduction_xx(:,:) = 0.0
1785  alloc_(cs%reduction_xy(isdb:iedb,jsdb:jedb)) ; cs%reduction_xy(:,:) = 0.0
1786 
1787  if (cs%anisotropic) then
1788  alloc_(cs%n1n2_h(isd:ied,jsd:jed)) ; cs%n1n2_h(:,:) = 0.0
1789  alloc_(cs%n1n1_m_n2n2_h(isd:ied,jsd:jed)) ; cs%n1n1_m_n2n2_h(:,:) = 0.0
1790  alloc_(cs%n1n2_q(isdb:iedb,jsdb:jedb)) ; cs%n1n2_q(:,:) = 0.0
1791  alloc_(cs%n1n1_m_n2n2_q(isdb:iedb,jsdb:jedb)) ; cs%n1n1_m_n2n2_q(:,:) = 0.0
1792  select case (aniso_mode)
1793  case (0)
1794  call align_aniso_tensor_to_grid(cs, aniso_grid_dir(1), aniso_grid_dir(2))
1795  case (1)
1796  ! call align_aniso_tensor_to_grid(CS, aniso_grid_dir(1), aniso_grid_dir(2))
1797  case (2)
1798  cs%dynamic_aniso = .true.
1799  case default
1800  call mom_error(fatal, "MOM_hor_visc: "//&
1801  "Runtime parameter ANISOTROPIC_MODE is out of range.")
1802  end select
1803  endif
1804 
1805  if (cs%use_Kh_bg_2d) then
1806  alloc_(cs%Kh_bg_2d(isd:ied,jsd:jed)) ; cs%Kh_bg_2d(:,:) = 0.0
1807  call get_param(param_file, mdl, "KH_BG_2D_FILENAME", filename, &
1808  'The filename containing a 2d map of "Kh".', &
1809  default='KH_background_2d.nc')
1810  call get_param(param_file, mdl, "INPUTDIR", inputdir, default=".")
1811  inputdir = slasher(inputdir)
1812  call mom_read_data(trim(inputdir)//trim(filename), 'Kh', cs%Kh_bg_2d, &
1813  g%domain, timelevel=1, scale=us%T_to_s)
1814  call pass_var(cs%Kh_bg_2d, g%domain)
1815  endif
1816 
1817  if (cs%biharmonic) then
1818  alloc_(cs%Idx2dyCu(isdb:iedb,jsd:jed)) ; cs%Idx2dyCu(:,:) = 0.0
1819  alloc_(cs%Idx2dyCv(isd:ied,jsdb:jedb)) ; cs%Idx2dyCv(:,:) = 0.0
1820  alloc_(cs%Idxdy2u(isdb:iedb,jsd:jed)) ; cs%Idxdy2u(:,:) = 0.0
1821  alloc_(cs%Idxdy2v(isd:ied,jsdb:jedb)) ; cs%Idxdy2v(:,:) = 0.0
1822 
1823  alloc_(cs%Ah_bg_xx(isd:ied,jsd:jed)) ; cs%Ah_bg_xx(:,:) = 0.0
1824  alloc_(cs%Ah_bg_xy(isdb:iedb,jsdb:jedb)) ; cs%Ah_bg_xy(:,:) = 0.0
1825  if (cs%bound_Ah .or. cs%better_bound_Ah) then
1826  alloc_(cs%Ah_Max_xx(isd:ied,jsd:jed)) ; cs%Ah_Max_xx(:,:) = 0.0
1827  alloc_(cs%Ah_Max_xy(isdb:iedb,jsdb:jedb)) ; cs%Ah_Max_xy(:,:) = 0.0
1828  endif
1829  if (cs%Smagorinsky_Ah) then
1830  alloc_(cs%Biharm_const_xx(isd:ied,jsd:jed)) ; cs%Biharm_const_xx(:,:) = 0.0
1831  alloc_(cs%Biharm_const_xy(isdb:iedb,jsdb:jedb)) ; cs%Biharm_const_xy(:,:) = 0.0
1832  if (cs%bound_Coriolis) then
1833  alloc_(cs%Biharm_const2_xx(isd:ied,jsd:jed)) ; cs%Biharm_const2_xx(:,:) = 0.0
1834  alloc_(cs%Biharm_const2_xy(isdb:iedb,jsdb:jedb)) ; cs%Biharm_const2_xy(:,:) = 0.0
1835  endif
1836  endif
1837  if (cs%Leith_Ah) then
1838  alloc_(cs%biharm5_const_xx(isd:ied,jsd:jed)) ; cs%biharm5_const_xx(:,:) = 0.0
1839  alloc_(cs%biharm5_const_xy(isdb:iedb,jsdb:jedb)) ; cs%biharm5_const_xy(:,:) = 0.0
1840  endif
1841  endif
1842 
1843  do j=js-2,jeq+1 ; do i=is-2,ieq+1
1844  cs%DX2q(i,j) = g%dxBu(i,j)*g%dxBu(i,j) ; cs%DY2q(i,j) = g%dyBu(i,j)*g%dyBu(i,j)
1845  cs%DX_dyBu(i,j) = g%dxBu(i,j)*g%IdyBu(i,j) ; cs%DY_dxBu(i,j) = g%dyBu(i,j)*g%IdxBu(i,j)
1846  enddo ; enddo
1847  do j=jsq-1,jeq+2 ; do i=isq-1,ieq+2
1848  cs%DX2h(i,j) = g%dxT(i,j)*g%dxT(i,j) ; cs%DY2h(i,j) = g%dyT(i,j)*g%dyT(i,j)
1849  cs%DX_dyT(i,j) = g%dxT(i,j)*g%IdyT(i,j) ; cs%DY_dxT(i,j) = g%dyT(i,j)*g%IdxT(i,j)
1850  enddo ; enddo
1851 
1852  do j=jsq,jeq+1 ; do i=isq,ieq+1
1853  cs%reduction_xx(i,j) = 1.0
1854  if ((g%dy_Cu(i,j) > 0.0) .and. (g%dy_Cu(i,j) < g%dyCu(i,j)) .and. &
1855  (g%dy_Cu(i,j) < g%dyCu(i,j) * cs%reduction_xx(i,j))) &
1856  cs%reduction_xx(i,j) = g%dy_Cu(i,j) / g%dyCu(i,j)
1857  if ((g%dy_Cu(i-1,j) > 0.0) .and. (g%dy_Cu(i-1,j) < g%dyCu(i-1,j)) .and. &
1858  (g%dy_Cu(i-1,j) < g%dyCu(i-1,j) * cs%reduction_xx(i,j))) &
1859  cs%reduction_xx(i,j) = g%dy_Cu(i-1,j) / g%dyCu(i-1,j)
1860  if ((g%dx_Cv(i,j) > 0.0) .and. (g%dx_Cv(i,j) < g%dxCv(i,j)) .and. &
1861  (g%dx_Cv(i,j) < g%dxCv(i,j) * cs%reduction_xx(i,j))) &
1862  cs%reduction_xx(i,j) = g%dx_Cv(i,j) / g%dxCv(i,j)
1863  if ((g%dx_Cv(i,j-1) > 0.0) .and. (g%dx_Cv(i,j-1) < g%dxCv(i,j-1)) .and. &
1864  (g%dx_Cv(i,j-1) < g%dxCv(i,j-1) * cs%reduction_xx(i,j))) &
1865  cs%reduction_xx(i,j) = g%dx_Cv(i,j-1) / g%dxCv(i,j-1)
1866  enddo ; enddo
1867 
1868  do j=js-1,jeq ; do i=is-1,ieq
1869  cs%reduction_xy(i,j) = 1.0
1870  if ((g%dy_Cu(i,j) > 0.0) .and. (g%dy_Cu(i,j) < g%dyCu(i,j)) .and. &
1871  (g%dy_Cu(i,j) < g%dyCu(i,j) * cs%reduction_xy(i,j))) &
1872  cs%reduction_xy(i,j) = g%dy_Cu(i,j) / g%dyCu(i,j)
1873  if ((g%dy_Cu(i,j+1) > 0.0) .and. (g%dy_Cu(i,j+1) < g%dyCu(i,j+1)) .and. &
1874  (g%dy_Cu(i,j+1) < g%dyCu(i,j+1) * cs%reduction_xy(i,j))) &
1875  cs%reduction_xy(i,j) = g%dy_Cu(i,j+1) / g%dyCu(i,j+1)
1876  if ((g%dx_Cv(i,j) > 0.0) .and. (g%dx_Cv(i,j) < g%dxCv(i,j)) .and. &
1877  (g%dx_Cv(i,j) < g%dxCv(i,j) * cs%reduction_xy(i,j))) &
1878  cs%reduction_xy(i,j) = g%dx_Cv(i,j) / g%dxCv(i,j)
1879  if ((g%dx_Cv(i+1,j) > 0.0) .and. (g%dx_Cv(i+1,j) < g%dxCv(i+1,j)) .and. &
1880  (g%dx_Cv(i+1,j) < g%dxCv(i+1,j) * cs%reduction_xy(i,j))) &
1881  cs%reduction_xy(i,j) = g%dx_Cv(i+1,j) / g%dxCv(i+1,j)
1882  enddo ; enddo
1883 
1884  if (cs%Laplacian) then
1885  ! The 0.3 below was 0.4 in MOM1.10. The change in hq requires
1886  ! this to be less than 1/3, rather than 1/2 as before.
1887  if (cs%bound_Kh .or. cs%bound_Ah) kh_limit = 0.3 / (dt*4.0)
1888 
1889  ! Calculate and store the background viscosity at h-points
1890  do j=jsq,jeq+1 ; do i=isq,ieq+1
1891  ! Static factors in the Smagorinsky and Leith schemes
1892  grid_sp_h2 = (2.0*cs%DX2h(i,j)*cs%DY2h(i,j)) / (cs%DX2h(i,j) + cs%DY2h(i,j))
1893  grid_sp_h3 = grid_sp_h2*sqrt(grid_sp_h2)
1894  if (cs%Smagorinsky_Kh) cs%Laplac2_const_xx(i,j) = smag_lap_const * grid_sp_h2
1895  if (cs%Leith_Kh) cs%Laplac3_const_xx(i,j) = leith_lap_const * grid_sp_h3
1896  ! Maximum of constant background and MICOM viscosity
1897  cs%Kh_bg_xx(i,j) = max(kh, kh_vel_scale * sqrt(grid_sp_h2))
1898 
1899  ! Use the larger of the above and values read from a file
1900  if (cs%use_Kh_bg_2d) cs%Kh_bg_xx(i,j) = max(cs%Kh_bg_2d(i,j), cs%Kh_bg_xx(i,j))
1901 
1902  ! Use the larger of the above and a function of sin(latitude)
1903  if (kh_sin_lat>0.) then
1904  slat_fn = abs( sin( deg2rad * g%geoLatT(i,j) ) ) ** kh_pwr_of_sine
1905  cs%Kh_bg_xx(i,j) = max(kh_sin_lat * slat_fn, cs%Kh_bg_xx(i,j))
1906  endif
1907 
1908  if (cs%bound_Kh .and. .not.cs%better_bound_Kh) then
1909  ! Limit the background viscosity to be numerically stable
1910  cs%Kh_Max_xx(i,j) = kh_limit * grid_sp_h2
1911  cs%Kh_bg_xx(i,j) = min(cs%Kh_bg_xx(i,j), cs%Kh_Max_xx(i,j))
1912  endif
1913  enddo ; enddo
1914 
1915  ! Calculate and store the background viscosity at q-points
1916  do j=js-1,jeq ; do i=is-1,ieq
1917  ! Static factors in the Smagorinsky and Leith schemes
1918  grid_sp_q2 = (2.0*cs%DX2q(i,j)*cs%DY2q(i,j)) / (cs%DX2q(i,j) + cs%DY2q(i,j))
1919  grid_sp_q3 = grid_sp_q2*sqrt(grid_sp_q2)
1920  if (cs%Smagorinsky_Kh) cs%Laplac2_const_xy(i,j) = smag_lap_const * grid_sp_q2
1921  if (cs%Leith_Kh) cs%Laplac3_const_xy(i,j) = leith_lap_const * grid_sp_q3
1922  ! Maximum of constant background and MICOM viscosity
1923  cs%Kh_bg_xy(i,j) = max(kh, kh_vel_scale * sqrt(grid_sp_q2))
1924 
1925  ! Use the larger of the above and values read from a file
1926  !### This expression uses inconsistent staggering
1927  if (cs%use_Kh_bg_2d) cs%Kh_bg_xy(i,j) = max(cs%Kh_bg_2d(i,j), cs%Kh_bg_xy(i,j))
1928 
1929  ! Use the larger of the above and a function of sin(latitude)
1930  if (kh_sin_lat>0.) then
1931  slat_fn = abs( sin( deg2rad * g%geoLatBu(i,j) ) ) ** kh_pwr_of_sine
1932  cs%Kh_bg_xy(i,j) = max(kh_sin_lat * slat_fn, cs%Kh_bg_xy(i,j))
1933  endif
1934 
1935  if (cs%bound_Kh .and. .not.cs%better_bound_Kh) then
1936  ! Limit the background viscosity to be numerically stable
1937  cs%Kh_Max_xy(i,j) = kh_limit * grid_sp_q2
1938  cs%Kh_bg_xy(i,j) = min(cs%Kh_bg_xy(i,j), cs%Kh_Max_xy(i,j))
1939  endif
1940  enddo ; enddo
1941  endif
1942 
1943  if (cs%biharmonic) then
1944 
1945  do j=js-1,jeq+1 ; do i=isq-1,ieq+1
1946  cs%IDX2dyCu(i,j) = (g%IdxCu(i,j)*g%IdxCu(i,j)) * g%IdyCu(i,j)
1947  cs%IDXDY2u(i,j) = g%IdxCu(i,j) * (g%IdyCu(i,j)*g%IdyCu(i,j))
1948  enddo ; enddo
1949  do j=jsq-1,jeq+1 ; do i=is-1,ieq+1
1950  cs%IDX2dyCv(i,j) = (g%IdxCv(i,j)*g%IdxCv(i,j)) * g%IdyCv(i,j)
1951  cs%IDXDY2v(i,j) = g%IdxCv(i,j) * (g%IdyCv(i,j)*g%IdyCv(i,j))
1952  enddo ; enddo
1953 
1954  cs%Ah_bg_xy(:,:) = 0.0
1955  ! The 0.3 below was 0.4 in MOM1.10. The change in hq requires
1956  ! this to be less than 1/3, rather than 1/2 as before.
1957  if (cs%better_bound_Ah .or. cs%bound_Ah) ah_limit = 0.3 / (dt*64.0)
1958  if (cs%Smagorinsky_Ah .and. cs%bound_Coriolis) &
1959  boundcorconst = 1.0 / (5.0*(bound_cor_vel*bound_cor_vel))
1960  do j=jsq,jeq+1 ; do i=isq,ieq+1
1961  grid_sp_h2 = (2.0*cs%DX2h(i,j)*cs%DY2h(i,j)) / (cs%DX2h(i,j)+cs%DY2h(i,j))
1962  grid_sp_h3 = grid_sp_h2*sqrt(grid_sp_h2)
1963 
1964  if (cs%Smagorinsky_Ah) then
1965  cs%Biharm_const_xx(i,j) = smag_bi_const * (grid_sp_h2 * grid_sp_h2)
1966  if (cs%bound_Coriolis) then
1967  fmax = max(abs(g%CoriolisBu(i-1,j-1)), abs(g%CoriolisBu(i,j-1)), &
1968  abs(g%CoriolisBu(i-1,j)), abs(g%CoriolisBu(i,j)))
1969  cs%Biharm_const2_xx(i,j) = us%m_to_L**2*(grid_sp_h2 * grid_sp_h2 * grid_sp_h2) * &
1970  (fmax * boundcorconst)
1971  endif
1972  endif
1973  if (cs%Leith_Ah) then
1974  cs%biharm5_const_xx(i,j) = leith_bi_const * (grid_sp_h3 * grid_sp_h2)
1975  endif
1976  cs%Ah_bg_xx(i,j) = max(ah, ah_vel_scale * grid_sp_h2 * sqrt(grid_sp_h2))
1977  if (ah_time_scale > 0.) cs%Ah_bg_xx(i,j) = &
1978  max(cs%Ah_bg_xx(i,j), (grid_sp_h2 * grid_sp_h2) / ah_time_scale)
1979  if (cs%bound_Ah .and. .not.cs%better_bound_Ah) then
1980  cs%Ah_Max_xx(i,j) = ah_limit * (grid_sp_h2 * grid_sp_h2)
1981  cs%Ah_bg_xx(i,j) = min(cs%Ah_bg_xx(i,j), cs%Ah_Max_xx(i,j))
1982  endif
1983  enddo ; enddo
1984  do j=js-1,jeq ; do i=is-1,ieq
1985  grid_sp_q2 = (2.0*cs%DX2q(i,j)*cs%DY2q(i,j)) / (cs%DX2q(i,j)+cs%DY2q(i,j))
1986  grid_sp_q3 = grid_sp_q2*sqrt(grid_sp_q2)
1987 
1988  if (cs%Smagorinsky_Ah) then
1989  cs%Biharm_const_xy(i,j) = smag_bi_const * (grid_sp_q2 * grid_sp_q2)
1990  if (cs%bound_Coriolis) then
1991  cs%Biharm_const2_xy(i,j) = us%m_to_L**2*(grid_sp_q2 * grid_sp_q2 * grid_sp_q2) * &
1992  (abs(g%CoriolisBu(i,j)) * boundcorconst)
1993  endif
1994  endif
1995  if (cs%Leith_Ah) then
1996  cs%biharm5_const_xy(i,j) = leith_bi_const * (grid_sp_q3 * grid_sp_q2)
1997  endif
1998 
1999  cs%Ah_bg_xy(i,j) = max(ah, ah_vel_scale * grid_sp_q2 * sqrt(grid_sp_q2))
2000  if (ah_time_scale > 0.) cs%Ah_bg_xy(i,j) = &
2001  max(cs%Ah_bg_xy(i,j), (grid_sp_q2 * grid_sp_q2) / ah_time_scale)
2002  if (cs%bound_Ah .and. .not.cs%better_bound_Ah) then
2003  cs%Ah_Max_xy(i,j) = ah_limit * (grid_sp_q2 * grid_sp_q2)
2004  cs%Ah_bg_xy(i,j) = min(cs%Ah_bg_xy(i,j), cs%Ah_Max_xy(i,j))
2005  endif
2006  enddo ; enddo
2007  endif
2008 
2009  ! The Laplacian bounds should avoid overshoots when CS%bound_coef < 1.
2010  if (cs%Laplacian .and. cs%better_bound_Kh) then
2011  idt = 1.0 / dt
2012  do j=jsq,jeq+1 ; do i=isq,ieq+1
2013  denom = max( &
2014  (cs%DY2h(i,j) * cs%DY_dxT(i,j) * (g%IdyCu(i,j) + g%IdyCu(i-1,j)) * &
2015  max(g%IdyCu(i,j)*g%IareaCu(i,j), g%IdyCu(i-1,j)*g%IareaCu(i-1,j)) ), &
2016  (cs%DX2h(i,j) * cs%DX_dyT(i,j) * (g%IdxCv(i,j) + g%IdxCv(i,j-1)) * &
2017  max(g%IdxCv(i,j)*g%IareaCv(i,j), g%IdxCv(i,j-1)*g%IareaCv(i,j-1)) ) )
2018  cs%Kh_Max_xx(i,j) = 0.0
2019  if (denom > 0.0) &
2020  cs%Kh_Max_xx(i,j) = cs%bound_coef * 0.25 * idt / denom
2021  enddo ; enddo
2022  do j=js-1,jeq ; do i=is-1,ieq
2023  denom = max( &
2024  (cs%DX2q(i,j) * cs%DX_dyBu(i,j) * (g%IdxCu(i,j+1) + g%IdxCu(i,j)) * &
2025  max(g%IdxCu(i,j)*g%IareaCu(i,j), g%IdxCu(i,j+1)*g%IareaCu(i,j+1)) ), &
2026  (cs%DY2q(i,j) * cs%DY_dxBu(i,j) * (g%IdyCv(i+1,j) + g%IdyCv(i,j)) * &
2027  max(g%IdyCv(i,j)*g%IareaCv(i,j), g%IdyCv(i+1,j)*g%IareaCv(i+1,j)) ) )
2028  cs%Kh_Max_xy(i,j) = 0.0
2029  if (denom > 0.0) &
2030  cs%Kh_Max_xy(i,j) = cs%bound_coef * 0.25 * idt / denom
2031  enddo ; enddo
2032  endif
2033 
2034  ! The biharmonic bounds should avoid overshoots when CS%bound_coef < 0.5, but
2035  ! empirically work for CS%bound_coef <~ 1.0
2036  if (cs%biharmonic .and. cs%better_bound_Ah) then
2037  idt = 1.0 / dt
2038  do j=js-1,jeq+1 ; do i=isq-1,ieq+1
2039  u0u(i,j) = cs%IDXDY2u(i,j)*(cs%DY2h(i+1,j)*cs%DY_dxT(i+1,j)*(g%IdyCu(i+1,j) + g%IdyCu(i,j)) + &
2040  cs%DY2h(i,j) * cs%DY_dxT(i,j) * (g%IdyCu(i,j) + g%IdyCu(i-1,j)) ) + &
2041  cs%IDX2dyCu(i,j)*(cs%DX2q(i,j) * cs%DX_dyBu(i,j) * (g%IdxCu(i,j+1) + g%IdxCu(i,j)) + &
2042  cs%DX2q(i,j-1)*cs%DX_dyBu(i,j-1)*(g%IdxCu(i,j) + g%IdxCu(i,j-1)) )
2043 
2044  u0v(i,j) = cs%IDXDY2u(i,j)*(cs%DY2h(i+1,j)*cs%DX_dyT(i+1,j)*(g%IdxCv(i+1,j) + g%IdxCv(i+1,j-1)) + &
2045  cs%DY2h(i,j) * cs%DX_dyT(i,j) * (g%IdxCv(i,j) + g%IdxCv(i,j-1)) ) + &
2046  cs%IDX2dyCu(i,j)*(cs%DX2q(i,j) * cs%DY_dxBu(i,j) * (g%IdyCv(i+1,j) + g%IdyCv(i,j)) + &
2047  cs%DX2q(i,j-1)*cs%DY_dxBu(i,j-1)*(g%IdyCv(i+1,j-1) + g%IdyCv(i,j-1)) )
2048  enddo ; enddo
2049  do j=jsq-1,jeq+1 ; do i=is-1,ieq+1
2050  v0u(i,j) = cs%IDXDY2v(i,j)*(cs%DY2q(i,j) * cs%DX_dyBu(i,j) * (g%IdxCu(i,j+1) + g%IdxCu(i,j)) + &
2051  cs%DY2q(i-1,j)*cs%DX_dyBu(i-1,j)*(g%IdxCu(i-1,j+1) + g%IdxCu(i-1,j)) ) + &
2052  cs%IDX2dyCv(i,j)*(cs%DX2h(i,j+1)*cs%DY_dxT(i,j+1)*(g%IdyCu(i,j+1) + g%IdyCu(i-1,j+1)) + &
2053  cs%DX2h(i,j) * cs%DY_dxT(i,j) * (g%IdyCu(i,j) + g%IdyCu(i-1,j)) )
2054 
2055  v0v(i,j) = cs%IDXDY2v(i,j)*(cs%DY2q(i,j) * cs%DY_dxBu(i,j) * (g%IdyCv(i+1,j) + g%IdyCv(i,j)) + &
2056  cs%DY2q(i-1,j)*cs%DY_dxBu(i-1,j)*(g%IdyCv(i,j) + g%IdyCv(i-1,j)) ) + &
2057  cs%IDX2dyCv(i,j)*(cs%DX2h(i,j+1)*cs%DX_dyT(i,j+1)*(g%IdxCv(i,j+1) + g%IdxCv(i,j)) + &
2058  cs%DX2h(i,j) * cs%DX_dyT(i,j) * (g%IdxCv(i,j) + g%IdxCv(i,j-1)) )
2059  enddo ; enddo
2060 
2061  do j=jsq,jeq+1 ; do i=isq,ieq+1
2062  denom = max( &
2063  (cs%DY2h(i,j) * &
2064  (cs%DY_dxT(i,j)*(g%IdyCu(i,j)*u0u(i,j) + g%IdyCu(i-1,j)*u0u(i-1,j)) + &
2065  cs%DX_dyT(i,j)*(g%IdxCv(i,j)*v0u(i,j) + g%IdxCv(i,j-1)*v0u(i,j-1))) * &
2066  max(g%IdyCu(i,j)*g%IareaCu(i,j), g%IdyCu(i-1,j)*g%IareaCu(i-1,j)) ), &
2067  (cs%DX2h(i,j) * &
2068  (cs%DY_dxT(i,j)*(g%IdyCu(i,j)*u0v(i,j) + g%IdyCu(i-1,j)*u0v(i-1,j)) + &
2069  cs%DX_dyT(i,j)*(g%IdxCv(i,j)*v0v(i,j) + g%IdxCv(i,j-1)*v0v(i,j-1))) * &
2070  max(g%IdxCv(i,j)*g%IareaCv(i,j), g%IdxCv(i,j-1)*g%IareaCv(i,j-1)) ) )
2071  cs%Ah_Max_xx(i,j) = 0.0
2072  if (denom > 0.0) &
2073  cs%Ah_Max_xx(i,j) = cs%bound_coef * 0.5 * idt / denom
2074  enddo ; enddo
2075 
2076  do j=js-1,jeq ; do i=is-1,ieq
2077  denom = max( &
2078  (cs%DX2q(i,j) * &
2079  (cs%DX_dyBu(i,j)*(u0u(i,j+1)*g%IdxCu(i,j+1) + u0u(i,j)*g%IdxCu(i,j)) + &
2080  cs%DY_dxBu(i,j)*(v0u(i+1,j)*g%IdyCv(i+1,j) + v0u(i,j)*g%IdyCv(i,j))) * &
2081  max(g%IdxCu(i,j)*g%IareaCu(i,j), g%IdxCu(i,j+1)*g%IareaCu(i,j+1)) ), &
2082  (cs%DY2q(i,j) * &
2083  (cs%DX_dyBu(i,j)*(u0v(i,j+1)*g%IdxCu(i,j+1) + u0v(i,j)*g%IdxCu(i,j)) + &
2084  cs%DY_dxBu(i,j)*(v0v(i+1,j)*g%IdyCv(i+1,j) + v0v(i,j)*g%IdyCv(i,j))) * &
2085  max(g%IdyCv(i,j)*g%IareaCv(i,j), g%IdyCv(i+1,j)*g%IareaCv(i+1,j)) ) )
2086  cs%Ah_Max_xy(i,j) = 0.0
2087  if (denom > 0.0) &
2088  cs%Ah_Max_xy(i,j) = cs%bound_coef * 0.5 * idt / denom
2089  enddo ; enddo
2090  endif
2091 
2092  ! Register fields for output from this module.
2093 
2094  cs%id_diffu = register_diag_field('ocean_model', 'diffu', diag%axesCuL, time, &
2095  'Zonal Acceleration from Horizontal Viscosity', 'm s-2', conversion=us%s_to_T)
2096 
2097  cs%id_diffv = register_diag_field('ocean_model', 'diffv', diag%axesCvL, time, &
2098  'Meridional Acceleration from Horizontal Viscosity', 'm s-2', conversion=us%s_to_T)
2099 
2100  if (cs%biharmonic) then
2101  cs%id_Ah_h = register_diag_field('ocean_model', 'Ahh', diag%axesTL, time, &
2102  'Biharmonic Horizontal Viscosity at h Points', 'm4 s-1', conversion=us%s_to_T, &
2103  cmor_field_name='difmxybo', &
2104  cmor_long_name='Ocean lateral biharmonic viscosity', &
2105  cmor_standard_name='ocean_momentum_xy_biharmonic_diffusivity')
2106 
2107  cs%id_Ah_q = register_diag_field('ocean_model', 'Ahq', diag%axesBL, time, &
2108  'Biharmonic Horizontal Viscosity at q Points', 'm4 s-1', conversion=us%s_to_T)
2109  endif
2110 
2111  if (cs%Laplacian) then
2112  cs%id_Kh_h = register_diag_field('ocean_model', 'Khh', diag%axesTL, time, &
2113  'Laplacian Horizontal Viscosity at h Points', 'm2 s-1', conversion=us%s_to_T, &
2114  cmor_field_name='difmxylo', &
2115  cmor_long_name='Ocean lateral Laplacian viscosity', &
2116  cmor_standard_name='ocean_momentum_xy_laplacian_diffusivity')
2117 
2118  cs%id_Kh_q = register_diag_field('ocean_model', 'Khq', diag%axesBL, time, &
2119  'Laplacian Horizontal Viscosity at q Points', 'm2 s-1', conversion=us%s_to_T)
2120 
2121  if (cs%Leith_Kh) then
2122  cs%id_vort_xy_q = register_diag_field('ocean_model', 'vort_xy_q', diag%axesBL, time, &
2123  'Vertical vorticity at q Points', 's-1')
2124 
2125  cs%id_div_xx_h = register_diag_field('ocean_model', 'div_xx_h', diag%axesTL, time, &
2126  'Horizontal divergence at h Points', 's-1')
2127  endif
2128 
2129  endif
2130 
2131  if (cs%use_GME) then
2132  cs%id_GME_coeff_h = register_diag_field('ocean_model', 'GME_coeff_h', diag%axesTL, time, &
2133  'GME coefficient at h Points', 'm2 s-1', conversion=us%s_to_T)
2134 
2135  cs%id_GME_coeff_q = register_diag_field('ocean_model', 'GME_coeff_q', diag%axesBL, time, &
2136  'GME coefficient at q Points', 'm2 s-1', conversion=us%s_to_T)
2137 
2138  cs%id_FrictWork_GME = register_diag_field('ocean_model','FrictWork_GME',diag%axesTL,time,&
2139  'Integral work done by lateral friction terms in GME (excluding diffusion of energy)', 'W m-2')
2140  endif
2141 
2142  cs%id_FrictWork = register_diag_field('ocean_model','FrictWork',diag%axesTL,time,&
2143  'Integral work done by lateral friction terms', 'W m-2')
2144 
2145  cs%id_FrictWork_diss = register_diag_field('ocean_model','FrictWork_diss',diag%axesTL,time,&
2146  'Integral work done by lateral friction terms (excluding diffusion of energy)', 'W m-2')
2147 
2148  if (associated(meke)) then
2149  if (associated(meke%mom_src)) then
2150  cs%id_FrictWorkMax = register_diag_field('ocean_model', 'FrictWorkMax', diag%axesTL, time,&
2151  'Maximum possible integral work done by lateral friction terms', 'W m-2')
2152  endif
2153  endif
2154 
2155  cs%id_FrictWorkIntz = register_diag_field('ocean_model','FrictWorkIntz',diag%axesT1,time, &
2156  'Depth integrated work done by lateral friction', 'W m-2', &
2157  cmor_field_name='dispkexyfo', &
2158  cmor_long_name='Depth integrated ocean kinetic energy dissipation due to lateral friction',&
2159  cmor_standard_name='ocean_kinetic_energy_dissipation_per_unit_area_due_to_xy_friction')
2160 
2161  if (cs%Laplacian .or. get_all) then
2162  endif
2163 
2164 end subroutine hor_visc_init
2165 
2166 !> Calculates factors in the anisotropic orientation tensor to be align with the grid.
2167 !! With n1=1 and n2=0, this recovers the approach of Large et al, 2001.
2168 subroutine align_aniso_tensor_to_grid(CS, n1, n2)
2169  type(hor_visc_cs), pointer :: CS !< Control structure for horizontal viscosity
2170  real, intent(in) :: n1 !< i-component of direction vector [nondim]
2171  real, intent(in) :: n2 !< j-component of direction vector [nondim]
2172  ! Local variables
2173  real :: recip_n2_norm
2174 
2175  ! For normalizing n=(n1,n2) in case arguments are not a unit vector
2176  recip_n2_norm = n1**2 + n2**2
2177  if (recip_n2_norm > 0.) recip_n2_norm = 1./recip_n2_norm
2178 
2179  cs%n1n2_h(:,:) = 2. * ( n1 * n2 ) * recip_n2_norm
2180  cs%n1n2_q(:,:) = 2. * ( n1 * n2 ) * recip_n2_norm
2181  cs%n1n1_m_n2n2_h(:,:) = ( n1 * n1 - n2 * n2 ) * recip_n2_norm
2182  cs%n1n1_m_n2n2_q(:,:) = ( n1 * n1 - n2 * n2 ) * recip_n2_norm
2183 
2184 end subroutine align_aniso_tensor_to_grid
2185 
2186 !> Apply a 1-1-4-1-1 Laplacian filter one time on GME diffusive flux to reduce any
2187 !! horizontal two-grid-point noise
2188 subroutine smooth_gme(CS,G,GME_flux_h,GME_flux_q)
2189  ! Arguments
2190  type(hor_visc_cs), pointer :: CS !< Control structure
2191  type(ocean_grid_type), intent(in) :: G !< Ocean grid
2192  real, dimension(SZI_(G),SZJ_(G)), optional, intent(inout) :: GME_flux_h!< GME diffusive flux
2193  !! at h points
2194  real, dimension(SZIB_(G),SZJB_(G)), optional, intent(inout) :: GME_flux_q!< GME diffusive flux
2195  !! at q points
2196 
2197  ! local variables
2198  real, dimension(SZI_(G),SZJ_(G)) :: GME_flux_h_original
2199  real, dimension(SZIB_(G),SZJB_(G)) :: GME_flux_q_original
2200  real :: wc, ww, we, wn, ws ! averaging weights for smoothing
2201  integer :: i, j, k, s
2202 
2203  !do s=1,CS%n_smooth
2204  do s=1,1
2205 
2206  ! Update halos
2207  if (present(gme_flux_h)) then
2208  call pass_var(gme_flux_h, g%Domain)
2209  gme_flux_h_original = gme_flux_h
2210  ! apply smoothing on GME
2211  do j = g%jsc, g%jec
2212  do i = g%isc, g%iec
2213  ! skip land points
2214  if (g%mask2dT(i,j)==0.) cycle
2215 
2216  ! compute weights
2217  ww = 0.125 * g%mask2dT(i-1,j)
2218  we = 0.125 * g%mask2dT(i+1,j)
2219  ws = 0.125 * g%mask2dT(i,j-1)
2220  wn = 0.125 * g%mask2dT(i,j+1)
2221  wc = 1.0 - (ww+we+wn+ws)
2222 
2223  gme_flux_h(i,j) = wc * gme_flux_h_original(i,j) &
2224  + ww * gme_flux_h_original(i-1,j) &
2225  + we * gme_flux_h_original(i+1,j) &
2226  + ws * gme_flux_h_original(i,j-1) &
2227  + wn * gme_flux_h_original(i,j+1)
2228  enddo; enddo
2229  endif
2230 
2231  ! Update halos
2232  if (present(gme_flux_q)) then
2233  call pass_var(gme_flux_q, g%Domain, position=corner, complete=.true.)
2234  gme_flux_q_original = gme_flux_q
2235  ! apply smoothing on GME
2236  do j = g%JscB, g%JecB
2237  do i = g%IscB, g%IecB
2238  ! skip land points
2239  if (g%mask2dBu(i,j)==0.) cycle
2240 
2241  ! compute weights
2242  ww = 0.125 * g%mask2dBu(i-1,j)
2243  we = 0.125 * g%mask2dBu(i+1,j)
2244  ws = 0.125 * g%mask2dBu(i,j-1)
2245  wn = 0.125 * g%mask2dBu(i,j+1)
2246  wc = 1.0 - (ww+we+wn+ws)
2247 
2248  gme_flux_q(i,j) = wc * gme_flux_q_original(i,j) &
2249  + ww * gme_flux_q_original(i-1,j) &
2250  + we * gme_flux_q_original(i+1,j) &
2251  + ws * gme_flux_q_original(i,j-1) &
2252  + wn * gme_flux_q_original(i,j+1)
2253  enddo; enddo
2254  endif
2255 
2256  enddo ! s-loop
2257 
2258 end subroutine smooth_gme
2259 
2260 !> Deallocates any variables allocated in hor_visc_init.
2261 subroutine hor_visc_end(CS)
2262  type(hor_visc_cs), pointer :: cs !< The control structure returned by a
2263  !! previous call to hor_visc_init.
2264 
2265  if (cs%Laplacian .or. cs%biharmonic) then
2266  dealloc_(cs%dx2h) ; dealloc_(cs%dx2q) ; dealloc_(cs%dy2h) ; dealloc_(cs%dy2q)
2267  dealloc_(cs%dx_dyT) ; dealloc_(cs%dy_dxT) ; dealloc_(cs%dx_dyBu) ; dealloc_(cs%dy_dxBu)
2268  dealloc_(cs%reduction_xx) ; dealloc_(cs%reduction_xy)
2269  endif
2270 
2271  if (cs%Laplacian) then
2272  dealloc_(cs%Kh_bg_xx) ; dealloc_(cs%Kh_bg_xy)
2273  if (cs%bound_Kh) then
2274  dealloc_(cs%Kh_Max_xx) ; dealloc_(cs%Kh_Max_xy)
2275  endif
2276  if (cs%Smagorinsky_Kh) then
2277  dealloc_(cs%Laplac2_const_xx) ; dealloc_(cs%Laplac2_const_xy)
2278  endif
2279  if (cs%Leith_Kh) then
2280  dealloc_(cs%Laplac3_const_xx) ; dealloc_(cs%Laplac3_const_xy)
2281  endif
2282  endif
2283 
2284  if (cs%biharmonic) then
2285  dealloc_(cs%Idx2dyCu) ; dealloc_(cs%Idx2dyCv)
2286  dealloc_(cs%Idxdy2u) ; dealloc_(cs%Idxdy2v)
2287  dealloc_(cs%Ah_bg_xx) ; dealloc_(cs%Ah_bg_xy)
2288  if (cs%bound_Ah) then
2289  dealloc_(cs%Ah_Max_xx) ; dealloc_(cs%Ah_Max_xy)
2290  endif
2291  if (cs%Smagorinsky_Ah) then
2292  dealloc_(cs%Biharm5_const_xx) ; dealloc_(cs%Biharm5_const_xy)
2293  if (cs%bound_Coriolis) then
2294  dealloc_(cs%Biharm5_const2_xx) ; dealloc_(cs%Biharm5_const2_xy)
2295  endif
2296  endif
2297  if (cs%Leith_Ah) then
2298  dealloc_(cs%Biharm_const_xx) ; dealloc_(cs%Biharm_const_xy)
2299  endif
2300  endif
2301  if (cs%anisotropic) then
2302  dealloc_(cs%n1n2_h)
2303  dealloc_(cs%n1n2_q)
2304  dealloc_(cs%n1n1_m_n2n2_h)
2305  dealloc_(cs%n1n1_m_n2n2_q)
2306  endif
2307  deallocate(cs)
2308 
2309 end subroutine hor_visc_end
2310 
2311 
2312 !> \namespace mom_hor_visc
2313 !!
2314 !! This module contains the subroutine horizontal_viscosity() that calculates the
2315 !! effects of horizontal viscosity, including parameterizations of the value of
2316 !! the viscosity itself. horizontal_viscosity() calculates the acceleration due to
2317 !! some combination of a biharmonic viscosity and a Laplacian viscosity. Either or
2318 !! both may use a coefficient that depends on the shear and strain of the flow.
2319 !! All metric terms are retained. The Laplacian is calculated as the divergence of
2320 !! a stress tensor, using the form suggested by Smagorinsky (1993). The biharmonic
2321 !! is calculated by twice applying the divergence of the stress tensor that is
2322 !! used to calculate the Laplacian, but without the dependence on thickness in the
2323 !! first pass. This form permits a variable viscosity, and indicates no
2324 !! acceleration for either resting fluid or solid body rotation.
2325 !!
2326 !! The form of the viscous accelerations is discussed extensively in Griffies and
2327 !! Hallberg (2000), and the implementation here follows that discussion closely.
2328 !! We use the notation of Smith and McWilliams (2003) with the exception that the
2329 !! isotropic viscosity is \f$\kappa_h\f$.
2330 !!
2331 !! \section section_horizontal_viscosity Horizontal viscosity in MOM
2332 !!
2333 !! In general, the horizontal stress tensor can be written as
2334 !! \f[
2335 !! {\bf \sigma} =
2336 !! \begin{pmatrix}
2337 !! \frac{1}{2} \left( \sigma_D + \sigma_T \right) & \frac{1}{2} \sigma_S \\\\
2338 !! \frac{1}{2} \sigma_S & \frac{1}{2} \left( \sigma_D - \sigma_T \right)
2339 !! \end{pmatrix}
2340 !! \f]
2341 !! where \f$\sigma_D\f$, \f$\sigma_T\f$ and \f$\sigma_S\f$ are stresses associated with
2342 !! invariant factors in the strain-rate tensor. For a Newtonian fluid, the stress
2343 !! tensor is usually linearly related to the strain-rate tensor. The horizontal
2344 !! strain-rate tensor is
2345 !! \f[
2346 !! \dot{\bf e} =
2347 !! \begin{pmatrix}
2348 !! \frac{1}{2} \left( \dot{e}_D + \dot{e}_T \right) & \frac{1}{2} \dot{e}_S \\\\
2349 !! \frac{1}{2} \dot{e}_S & \frac{1}{2} \left( \dot{e}_D - \dot{e}_T \right)
2350 !! \end{pmatrix}
2351 !! \f]
2352 !! where \f$\dot{e}_D = \partial_x u + \partial_y v\f$ is the horizontal divergence,
2353 !! \f$\dot{e}_T = \partial_x u - \partial_y v\f$ is the horizontal tension, and
2354 !! \f$\dot{e}_S = \partial_y u + \partial_x v\f$ is the horizontal shear strain.
2355 !!
2356 !! The trace of the stress tensor, \f$tr(\bf \sigma) = \sigma_D\f$, is usually
2357 !! absorbed into the pressure and only the deviatoric stress tensor considered.
2358 !! From here on, we drop \f$\sigma_D\f$. The trace of the strain tensor, \f$tr(\bf e) =
2359 !! \dot{e}_D\f$ is non-zero for horizontally divergent flow but only enters the
2360 !! stress tensor through \f$\sigma_D\f$ and so we will drop \f$\sigma_D\f$ from
2361 !! calculations of the strain tensor in the code. Therefore the horizontal stress
2362 !! tensor can be considered to be
2363 !! \f[
2364 !! {\bf \sigma} =
2365 !! \begin{pmatrix}
2366 !! \frac{1}{2} \sigma_T & \frac{1}{2} \sigma_S \\\\
2367 !! \frac{1}{2} \sigma_S & - \frac{1}{2} \sigma_T
2368 !! \end{pmatrix}
2369 !! .\f]
2370 !!
2371 !! The stresses above are linearly related to the strain through a viscosity
2372 !! coefficient, \f$\kappa_h\f$:
2373 !! \f{eqnarray*}{
2374 !! \sigma_T & = & 2 \kappa_h \dot{e}_T \\\\
2375 !! \sigma_S & = & 2 \kappa_h \dot{e}_S
2376 !! .
2377 !! \f}
2378 !!
2379 !! The viscosity \f$\kappa_h\f$ may either be a constant or variable. For example,
2380 !! \f$\kappa_h\f$ may vary with the shear, as proposed by Smagorinsky (1993).
2381 !!
2382 !! The accelerations resulting form the divergence of the stress tensor are
2383 !! \f{eqnarray*}{
2384 !! \hat{\bf x} \cdot \left( \nabla \cdot {\bf \sigma} \right)
2385 !! & = &
2386 !! \partial_x \left( \frac{1}{2} \sigma_T \right)
2387 !! + \partial_y \left( \frac{1}{2} \sigma_S \right)
2388 !! \\\\
2389 !! & = &
2390 !! \partial_x \left( \kappa_h \dot{e}_T \right)
2391 !! + \partial_y \left( \kappa_h \dot{e}_S \right)
2392 !! \\\\
2393 !! \hat{\bf y} \cdot \left( \nabla \cdot {\bf \sigma} \right)
2394 !! & = &
2395 !! \partial_x \left( \frac{1}{2} \sigma_S \right)
2396 !! + \partial_y \left( \frac{1}{2} \sigma_T \right)
2397 !! \\\\
2398 !! & = &
2399 !! \partial_x \left( \kappa_h \dot{e}_S \right)
2400 !! + \partial_y \left( - \kappa_h \dot{e}_T \right)
2401 !! .
2402 !! \f}
2403 !!
2404 !! The form of the Laplacian viscosity in general coordinates is:
2405 !! \f{eqnarray*}{
2406 !! \hat{\bf x} \cdot \left( \nabla \cdot \sigma \right)
2407 !! & = &
2408 !! \frac{1}{h} \left[ \partial_x \left( \kappa_h h \dot{e}_T \right)
2409 !! + \partial_y \left( \kappa_h h \dot{e}_S \right) \right]
2410 !! \\\\
2411 !! \hat{\bf y} \cdot \left( \nabla \cdot \sigma \right)
2412 !! & = &
2413 !! \frac{1}{h} \left[ \partial_x \left( \kappa_h h \dot{e}_S \right)
2414 !! - \partial_y \left( \kappa_h h \dot{e}_T \right) \right]
2415 !! .
2416 !! \f}
2417 !!
2418 !! \subsection section_laplacian_viscosity_coefficient Laplacian viscosity coefficient
2419 !!
2420 !! The horizontal viscosity coefficient, \f$\kappa_h\f$, can have multiple components.
2421 !! The isotropic components are:
2422 !! - A uniform background component, \f$\kappa_{bg}\f$.
2423 !! - A constant but spatially variable 2D map, \f$\kappa_{2d}(x,y)\f$.
2424 !! - A ''MICOM'' viscosity, \f$U_\nu \Delta(x,y)\f$, which uses a constant
2425 !! velocity scale, \f$U_\nu\f$ and a measure of the grid-spacing \f$\Delta(x,y)^2 =
2426 !! \frac{2 \Delta x^2 \Delta y^2}{\Delta x^2 + \Delta y^2}\f$.
2427 !! - A function of
2428 !! latitude, \f$\kappa_{\phi}(x,y) = \kappa_{\pi/2} |\sin(\phi)|^n\f$.
2429 !! - A dynamic Smagorinsky viscosity, \f$\kappa_{Sm}(x,y,t) = C_{Sm} \Delta^2 \sqrt{\dot{e}_T^2 + \dot{e}_S^2}\f$.
2430 !! - A dynamic Leith viscosity, \f$\kappa_{Lth}(x,y,t) =
2431 !! C_{Lth} \Delta^3 \sqrt{|\nabla \zeta|^2 + |\nabla \dot{e}_D|^2}\f$.
2432 !!
2433 !! A maximum stable viscosity, \f$\kappa_{max}(x,y)\f$ is calculated based on the
2434 !! grid-spacing and time-step and used to clip calculated viscosities.
2435 !!
2436 !! The static components of \f$\kappa_h\f$ are first combined as follows:
2437 !! \f[
2438 !! \kappa_{static} = \min \left[ \max\left(
2439 !! \kappa_{bg},
2440 !! U_\nu \Delta(x,y),
2441 !! \kappa_{2d}(x,y),
2442 !! \kappa_\phi(x,y)
2443 !! \right)
2444 !! , \kappa_{max}(x,y) \right]
2445 !! \f]
2446 !! and stored in the module control structure as variables <code>Kh_bg_xx</code> and
2447 !! <code>Kh_bg_xy</code> for the tension (h-points) and shear (q-points) components
2448 !! respectively.
2449 !!
2450 !! The full viscosity includes the dynamic components as follows:
2451 !! \f[
2452 !! \kappa_h(x,y,t) = r(\Delta,L_d)
2453 !! \max \left( \kappa_{static}, \kappa_{Sm}, \kappa_{Lth} \right)
2454 !! \f]
2455 !! where \f$r(\Delta,L_d)\f$ is a resolution function.
2456 !!
2457 !! The dynamic Smagorinsky and Leith viscosity schemes are exclusive with each
2458 !! other.
2459 !!
2460 !! \subsection section_viscous_boundary_conditions Viscous boundary conditions
2461 !!
2462 !! Free slip boundary conditions have been coded, although no slip boundary
2463 !! conditions can be used with the Laplacian viscosity based on the 2D land-sea
2464 !! mask. For a western boundary, for example, the boundary conditions with the
2465 !! biharmonic operator would be written as:
2466 !! \f[
2467 !! \partial_x v = 0 ; \partial_x^3 v = 0 ; u = 0 ; \partial_x^2 u = 0 ,
2468 !! \f]
2469 !! while for a Laplacian operator, they are simply
2470 !! \f[
2471 !! \partial_x v = 0 ; u = 0 .
2472 !! \f]
2473 !! These boundary conditions are largely dictated by the use of an Arakawa
2474 !! C-grid and by the varying layer thickness.
2475 !!
2476 !! \subsection section_anisotropic_viscosity Anisotropic viscosity
2477 !!
2478 !! Large et al., 2001, proposed enhancing viscosity in a particular direction and the
2479 !! approach was generalized in Smith and McWilliams, 2003. We use the second form of their
2480 !! two coefficient anisotropic viscosity (section 4.3). We also replace their
2481 !! \f$A^\prime\f$ and $D$ such that \f$2A^\prime = 2 \kappa_h + D\f$ and
2482 !! \f$\kappa_a = D\f$ so that \f$\kappa_h\f$ can be considered the isotropic
2483 !! viscosity and \f$\kappa_a=D\f$ can be consider the anisotropic viscosity. The
2484 !! direction of anisotropy is defined by a unit vector \f$\hat{\bf
2485 !! n}=(n_1,n_2)\f$.
2486 !!
2487 !! The contributions to the stress tensor are
2488 !! \f[
2489 !! \begin{pmatrix}
2490 !! \sigma_T \\\\ \sigma_S
2491 !! \end{pmatrix}
2492 !! =
2493 !! \left[
2494 !! \begin{pmatrix}
2495 !! 2 \kappa_h + \kappa_a & 0 \\\\
2496 !! 0 & 2 \kappa_h
2497 !! \end{pmatrix}
2498 !! + 2 \kappa_a n_1 n_2
2499 !! \begin{pmatrix}
2500 !! - 2 n_1 n_2 & n_1^2 - n_2^2 \\\\
2501 !! n_1^2 - n_2^2 & 2 n_1 n_2
2502 !! \end{pmatrix}
2503 !! \right]
2504 !! \begin{pmatrix}
2505 !! \dot{e}_T \\\\ \dot{e}_S
2506 !! \end{pmatrix}
2507 !! \f]
2508 !! Dissipation of kinetic energy requires \f$\kappa_h \geq 0\f$ and \f$2 \kappa_h + \kappa_a \geq 0\f$.
2509 !! Note that when anisotropy is aligned with the x-direction, \f$n_1 = \pm 1\f$, then
2510 !! \f$n_2 = 0\f$ and the cross terms vanish. The accelerations in this aligned limit
2511 !! with constant coefficients become
2512 !! \f{eqnarray*}{
2513 !! \hat{\bf x} \cdot \nabla \cdot {\bf \sigma}
2514 !! & = &
2515 !! \partial_x \left( \left( \kappa_h + \frac{1}{2} \kappa_a \right) \dot{e}_T \right)
2516 !! + \partial_y \left( \kappa_h \dot{e}_S \right)
2517 !! \\\\
2518 !! & = &
2519 !! \left( \kappa_h + \kappa_a \right) \partial_{xx} u
2520 !! + \kappa_h \partial_{yy} u
2521 !! - \frac{1}{2} \kappa_a \partial_x \left( \partial_x u + \partial_y v \right)
2522 !! \\\\
2523 !! \hat{\bf y} \cdot \nabla \cdot {\bf \sigma}
2524 !! & = &
2525 !! \partial_x \left( \kappa_h \dot{e}_S \right)
2526 !! - \partial_y \left( \left( \kappa_h + \frac{1}{2} \kappa_a \right) \dot{e}_T \right)
2527 !! \\\\
2528 !! & = &
2529 !! \kappa_h \partial_{xx} v
2530 !! + \left( \kappa_h + \kappa_a \right) \partial_{yy} v
2531 !! - \frac{1}{2} \kappa_a \partial_y \left( \partial_x u + \partial_y v \right)
2532 !! \f}
2533 !! which has contributions akin to a negative divergence damping (a divergence
2534 !! enhancement?) but which is weaker than the enhanced tension terms by half.
2535 !!
2536 !! \subsection section_viscous_discretization Discretization
2537 !!
2538 !! The horizontal tension, \f$\dot{e}_T\f$, is stored in variable <code>sh_xx</code> and
2539 !! discretized as
2540 !! \f[
2541 !! \dot{e}_T
2542 !! = \frac{\Delta y}{\Delta x} \delta_i \left( \frac{1}{\Delta y} u \right)
2543 !! - \frac{\Delta x}{\Delta y} \delta_j \left( \frac{1}{\Delta x} v \right)
2544 !! .
2545 !! \f]
2546 !! The horizontal divergent strain, \f$\dot{e}_D\f$, is stored in variable
2547 !! <code>div_xx</code> and discretized as
2548 !! \f[
2549 !! \dot{e}_D
2550 !! = \frac{1}{h A} \left( \delta_i \left( \overline{h}^i \Delta y \, u \right)
2551 !! + \delta_j \left( \overline{h}^j\Delta x \, v \right) \right)
2552 !! .
2553 !! \f]
2554 !! Note that for expediency this is the exact discretization used in the
2555 !! continuity equation.
2556 !!
2557 !! The horizontal shear strain, \f$\dot{e}_S\f$, is stored in variable <code>sh_xy</code>
2558 !! and discretized as
2559 !! \f[
2560 !! \dot{e}_S = v_x + u_y
2561 !! \f]
2562 !! where
2563 !! \f{align*}{
2564 !! v_x &= \frac{\Delta y}{\Delta x} \delta_i \left( \frac{1}{\Delta y} v \right) \\\\
2565 !! u_y &= \frac{\Delta x}{\Delta y} \delta_j \left( \frac{1}{\Delta x} u \right)
2566 !! \f}
2567 !! which are calculated separately so that no-slip or free-slip boundary
2568 !! conditions can be applied to \f$v_x\f$ and \f$u_y\f$ where appropriate.
2569 !!
2570 !! The tendency for the x-component of the divergence of stress is stored in
2571 !! variable <code>diffu</code> and discretized as
2572 !! \f[
2573 !! \hat{\bf x} \cdot \left( \nabla \cdot {\bf \sigma} \right) =
2574 !! \frac{1}{A \overline{h}^i} \left(
2575 !! \frac{1}{\Delta y} \delta_i \left( h \Delta y^2 \kappa_h \dot{e}_T \right) +
2576 !! \frac{1}{\Delta x} \delta_j \left( \tilde{h}^{ij} \Delta x^2 \kappa_h \dot{e}_S \right)
2577 !! \right)
2578 !! .
2579 !! \f]
2580 !!
2581 !! The tendency for the y-component of the divergence of stress is stored in
2582 !! variable <code>diffv</code> and discretized as
2583 !! \f[
2584 !! \hat{\bf y} \cdot \left( \nabla \cdot {\bf \sigma} \right) =
2585 !! \frac{1}{A \overline{h}^j} \left(
2586 !! \frac{1}{\Delta y} \delta_i \left( \tilde{h}^{ij} \Delta y^2 A_M \dot{e}_S \right)
2587 !! - \frac{1}{\Delta x} \delta_j \left( h \Delta x^2 A_M \dot{e}_T \right)
2588 !! \right)
2589 !! .
2590 !! \f]
2591 !!
2592 !! \subsection section_viscous_refs References
2593 !!
2594 !! Griffies, S.M., and Hallberg, R.W., 2000: Biharmonic friction with a
2595 !! Smagorinsky-like viscosity for use in large-scale eddy-permitting ocean models.
2596 !! Monthly Weather Review, 128(8), 2935-2946.
2597 !! https://doi.org/10.1175/1520-0493(2000)128%3C2935:BFWASL%3E2.0.CO;2
2598 !!
2599 !! Large, W.G., Danabasoglu, G., McWilliams, J.C., Gent, P.R. and Bryan, F.O.,
2600 !! 2001: Equatorial circulation of a global ocean climate model with
2601 !! anisotropic horizontal viscosity.
2602 !! Journal of Physical Oceanography, 31(2), pp.518-536.
2603 !! https://doi.org/10.1175/1520-0485(2001)031%3C0518:ECOAGO%3E2.0.CO;2
2604 !!
2605 !! Smagorinsky, J., 1993: Some historical remarks on the use of nonlinear
2606 !! viscosities. Large eddy simulation of complex engineering and geophysical
2607 !! flows, 1, 69-106.
2608 !!
2609 !! Smith, R.D., and McWilliams, J.C., 2003: Anisotropic horizontal viscosity for
2610 !! ocean models. Ocean Modelling, 5(2), 129-156.
2611 !! https://doi.org/10.1016/S1463-5003(02)00016-1
2612 
2613 end module mom_hor_visc
mom_verticalgrid
Provides a transparent vertical ocean grid type and supporting routines.
Definition: MOM_verticalGrid.F90:2
mom_hor_visc::hor_visc_cs
Control structure for horizontal viscosity.
Definition: MOM_hor_visc.F90:29
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_barotropic
Baropotric solver.
Definition: MOM_barotropic.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_io
This module contains I/O framework code.
Definition: MOM_io.F90:2
mom_unit_scaling::unit_scale_type
Describes various unit conversion factors.
Definition: MOM_unit_scaling.F90:14
mom_diag_mediator::post_data
Make a diagnostic available for averaging or output.
Definition: MOM_diag_mediator.F90:70
mom_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_hor_visc
Calculates horizontal viscosity and viscous stresses.
Definition: MOM_hor_visc.F90:2
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_io::mom_read_data
Read a data field from a file.
Definition: MOM_io.F90:74
mom_open_boundary
Controls where open boundary conditions are applied.
Definition: MOM_open_boundary.F90:2
mom_file_parser
The MOM6 facility to parse input files for runtime parameters.
Definition: MOM_file_parser.F90:2
mom_grid
Provides the ocean grid type.
Definition: MOM_grid.F90:2
mom_open_boundary::ocean_obc_type
Open-boundary data.
Definition: MOM_open_boundary.F90:186
mom_unit_scaling
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
Definition: MOM_unit_scaling.F90:2
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_barotropic::barotropic_cs
The barotropic stepping control stucture.
Definition: MOM_barotropic.F90:100
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