MOM6
mom_hor_visc Module Reference

Detailed Description

Calculates horizontal viscosity and viscous stresses.

This module contains the subroutine horizontal_viscosity() that calculates the effects of horizontal viscosity, including parameterizations of the value of the viscosity itself. horizontal_viscosity() calculates the acceleration due to some combination of a biharmonic viscosity and a Laplacian viscosity. Either or both may use a coefficient that depends on the shear and strain of the flow. All metric terms are retained. The Laplacian is calculated as the divergence of a stress tensor, using the form suggested by Smagorinsky (1993). The biharmonic is calculated by twice applying the divergence of the stress tensor that is used to calculate the Laplacian, but without the dependence on thickness in the first pass. This form permits a variable viscosity, and indicates no acceleration for either resting fluid or solid body rotation.

The form of the viscous accelerations is discussed extensively in Griffies and Hallberg (2000), and the implementation here follows that discussion closely. We use the notation of Smith and McWilliams (2003) with the exception that the isotropic viscosity is \(\kappa_h\).

Horizontal viscosity in MOM

In general, the horizontal stress tensor can be written as

\[ {\bf \sigma} = \begin{pmatrix} \frac{1}{2} \left( \sigma_D + \sigma_T \right) & \frac{1}{2} \sigma_S \\\\ \frac{1}{2} \sigma_S & \frac{1}{2} \left( \sigma_D - \sigma_T \right) \end{pmatrix} \]

where \(\sigma_D\), \(\sigma_T\) and \(\sigma_S\) are stresses associated with invariant factors in the strain-rate tensor. For a Newtonian fluid, the stress tensor is usually linearly related to the strain-rate tensor. The horizontal strain-rate tensor is

\[ \dot{\bf e} = \begin{pmatrix} \frac{1}{2} \left( \dot{e}_D + \dot{e}_T \right) & \frac{1}{2} \dot{e}_S \\\\ \frac{1}{2} \dot{e}_S & \frac{1}{2} \left( \dot{e}_D - \dot{e}_T \right) \end{pmatrix} \]

where \(\dot{e}_D = \partial_x u + \partial_y v\) is the horizontal divergence, \(\dot{e}_T = \partial_x u - \partial_y v\) is the horizontal tension, and \(\dot{e}_S = \partial_y u + \partial_x v\) is the horizontal shear strain.

The trace of the stress tensor, \(tr(\bf \sigma) = \sigma_D\), is usually absorbed into the pressure and only the deviatoric stress tensor considered. From here on, we drop \(\sigma_D\). The trace of the strain tensor, \(tr(\bf e) = \dot{e}_D\) is non-zero for horizontally divergent flow but only enters the stress tensor through \(\sigma_D\) and so we will drop \(\sigma_D\) from calculations of the strain tensor in the code. Therefore the horizontal stress tensor can be considered to be

\[ {\bf \sigma} = \begin{pmatrix} \frac{1}{2} \sigma_T & \frac{1}{2} \sigma_S \\\\ \frac{1}{2} \sigma_S & - \frac{1}{2} \sigma_T \end{pmatrix} .\]

The stresses above are linearly related to the strain through a viscosity coefficient, \(\kappa_h\):

\begin{eqnarray*} \sigma_T & = & 2 \kappa_h \dot{e}_T \\\\ \sigma_S & = & 2 \kappa_h \dot{e}_S . \end{eqnarray*}

The viscosity \(\kappa_h\) may either be a constant or variable. For example, \(\kappa_h\) may vary with the shear, as proposed by Smagorinsky (1993).

The accelerations resulting form the divergence of the stress tensor are

\begin{eqnarray*} \hat{\bf x} \cdot \left( \nabla \cdot {\bf \sigma} \right) & = & \partial_x \left( \frac{1}{2} \sigma_T \right) + \partial_y \left( \frac{1}{2} \sigma_S \right) \\\\ & = & \partial_x \left( \kappa_h \dot{e}_T \right) + \partial_y \left( \kappa_h \dot{e}_S \right) \\\\ \hat{\bf y} \cdot \left( \nabla \cdot {\bf \sigma} \right) & = & \partial_x \left( \frac{1}{2} \sigma_S \right) + \partial_y \left( \frac{1}{2} \sigma_T \right) \\\\ & = & \partial_x \left( \kappa_h \dot{e}_S \right) + \partial_y \left( - \kappa_h \dot{e}_T \right) . \end{eqnarray*}

The form of the Laplacian viscosity in general coordinates is:

\begin{eqnarray*} \hat{\bf x} \cdot \left( \nabla \cdot \sigma \right) & = & \frac{1}{h} \left[ \partial_x \left( \kappa_h h \dot{e}_T \right) + \partial_y \left( \kappa_h h \dot{e}_S \right) \right] \\\\ \hat{\bf y} \cdot \left( \nabla \cdot \sigma \right) & = & \frac{1}{h} \left[ \partial_x \left( \kappa_h h \dot{e}_S \right) - \partial_y \left( \kappa_h h \dot{e}_T \right) \right] . \end{eqnarray*}

Laplacian viscosity coefficient

The horizontal viscosity coefficient, \(\kappa_h\), can have multiple components. The isotropic components are:

  • A uniform background component, \(\kappa_{bg}\).
  • A constant but spatially variable 2D map, \(\kappa_{2d}(x,y)\).
  • A ''MICOM'' viscosity, \(U_\nu \Delta(x,y)\), which uses a constant velocity scale, \(U_\nu\) and a measure of the grid-spacing \(\Delta(x,y)^2 = \frac{2 \Delta x^2 \Delta y^2}{\Delta x^2 + \Delta y^2}\).
  • A function of latitude, \(\kappa_{\phi}(x,y) = \kappa_{\pi/2} |\sin(\phi)|^n\).
  • A dynamic Smagorinsky viscosity, \(\kappa_{Sm}(x,y,t) = C_{Sm} \Delta^2 \sqrt{\dot{e}_T^2 + \dot{e}_S^2}\).
  • A dynamic Leith viscosity, \(\kappa_{Lth}(x,y,t) = C_{Lth} \Delta^3 \sqrt{|\nabla \zeta|^2 + |\nabla \dot{e}_D|^2}\).

A maximum stable viscosity, \(\kappa_{max}(x,y)\) is calculated based on the grid-spacing and time-step and used to clip calculated viscosities.

The static components of \(\kappa_h\) are first combined as follows:

\[ \kappa_{static} = \min \left[ \max\left( \kappa_{bg}, U_\nu \Delta(x,y), \kappa_{2d}(x,y), \kappa_\phi(x,y) \right) , \kappa_{max}(x,y) \right] \]

and stored in the module control structure as variables Kh_bg_xx and Kh_bg_xy for the tension (h-points) and shear (q-points) components respectively.

The full viscosity includes the dynamic components as follows:

\[ \kappa_h(x,y,t) = r(\Delta,L_d) \max \left( \kappa_{static}, \kappa_{Sm}, \kappa_{Lth} \right) \]

where \(r(\Delta,L_d)\) is a resolution function.

The dynamic Smagorinsky and Leith viscosity schemes are exclusive with each other.

Viscous boundary conditions

Free slip boundary conditions have been coded, although no slip boundary conditions can be used with the Laplacian viscosity based on the 2D land-sea mask. For a western boundary, for example, the boundary conditions with the biharmonic operator would be written as:

\[ \partial_x v = 0 ; \partial_x^3 v = 0 ; u = 0 ; \partial_x^2 u = 0 , \]

while for a Laplacian operator, they are simply

\[ \partial_x v = 0 ; u = 0 . \]

These boundary conditions are largely dictated by the use of an Arakawa C-grid and by the varying layer thickness.

Anisotropic viscosity

Large et al., 2001, proposed enhancing viscosity in a particular direction and the approach was generalized in Smith and McWilliams, 2003. We use the second form of their two coefficient anisotropic viscosity (section 4.3). We also replace their \(A^\prime\) and $D$ such that \(2A^\prime = 2 \kappa_h + D\) and \(\kappa_a = D\) so that \(\kappa_h\) can be considered the isotropic viscosity and \(\kappa_a=D\) can be consider the anisotropic viscosity. The direction of anisotropy is defined by a unit vector \(\hat{\bf n}=(n_1,n_2)\).

The contributions to the stress tensor are

\[ \begin{pmatrix} \sigma_T \\\\ \sigma_S \end{pmatrix} = \left[ \begin{pmatrix} 2 \kappa_h + \kappa_a & 0 \\\\ 0 & 2 \kappa_h \end{pmatrix} + 2 \kappa_a n_1 n_2 \begin{pmatrix} - 2 n_1 n_2 & n_1^2 - n_2^2 \\\\ n_1^2 - n_2^2 & 2 n_1 n_2 \end{pmatrix} \right] \begin{pmatrix} \dot{e}_T \\\\ \dot{e}_S \end{pmatrix} \]

Dissipation of kinetic energy requires \(\kappa_h \geq 0\) and \(2 \kappa_h + \kappa_a \geq 0\). Note that when anisotropy is aligned with the x-direction, \(n_1 = \pm 1\), then \(n_2 = 0\) and the cross terms vanish. The accelerations in this aligned limit with constant coefficients become

\begin{eqnarray*} \hat{\bf x} \cdot \nabla \cdot {\bf \sigma} & = & \partial_x \left( \left( \kappa_h + \frac{1}{2} \kappa_a \right) \dot{e}_T \right) + \partial_y \left( \kappa_h \dot{e}_S \right) \\\\ & = & \left( \kappa_h + \kappa_a \right) \partial_{xx} u + \kappa_h \partial_{yy} u - \frac{1}{2} \kappa_a \partial_x \left( \partial_x u + \partial_y v \right) \\\\ \hat{\bf y} \cdot \nabla \cdot {\bf \sigma} & = & \partial_x \left( \kappa_h \dot{e}_S \right) - \partial_y \left( \left( \kappa_h + \frac{1}{2} \kappa_a \right) \dot{e}_T \right) \\\\ & = & \kappa_h \partial_{xx} v + \left( \kappa_h + \kappa_a \right) \partial_{yy} v - \frac{1}{2} \kappa_a \partial_y \left( \partial_x u + \partial_y v \right) \end{eqnarray*}

which has contributions akin to a negative divergence damping (a divergence enhancement?) but which is weaker than the enhanced tension terms by half.

Discretization

The horizontal tension, \(\dot{e}_T\), is stored in variable sh_xx and discretized as

\[ \dot{e}_T = \frac{\Delta y}{\Delta x} \delta_i \left( \frac{1}{\Delta y} u \right) - \frac{\Delta x}{\Delta y} \delta_j \left( \frac{1}{\Delta x} v \right) . \]

The horizontal divergent strain, \(\dot{e}_D\), is stored in variable div_xx and discretized as

\[ \dot{e}_D = \frac{1}{h A} \left( \delta_i \left( \overline{h}^i \Delta y \, u \right) + \delta_j \left( \overline{h}^j\Delta x \, v \right) \right) . \]

Note that for expediency this is the exact discretization used in the continuity equation.

The horizontal shear strain, \(\dot{e}_S\), is stored in variable sh_xy and discretized as

\[ \dot{e}_S = v_x + u_y \]

where

\begin{align*} v_x &= \frac{\Delta y}{\Delta x} \delta_i \left( \frac{1}{\Delta y} v \right) \\\\ u_y &= \frac{\Delta x}{\Delta y} \delta_j \left( \frac{1}{\Delta x} u \right) \end{align*}

which are calculated separately so that no-slip or free-slip boundary conditions can be applied to \(v_x\) and \(u_y\) where appropriate.

The tendency for the x-component of the divergence of stress is stored in variable diffu and discretized as

\[ \hat{\bf x} \cdot \left( \nabla \cdot {\bf \sigma} \right) = \frac{1}{A \overline{h}^i} \left( \frac{1}{\Delta y} \delta_i \left( h \Delta y^2 \kappa_h \dot{e}_T \right) + \frac{1}{\Delta x} \delta_j \left( \tilde{h}^{ij} \Delta x^2 \kappa_h \dot{e}_S \right) \right) . \]

The tendency for the y-component of the divergence of stress is stored in variable diffv and discretized as

\[ \hat{\bf y} \cdot \left( \nabla \cdot {\bf \sigma} \right) = \frac{1}{A \overline{h}^j} \left( \frac{1}{\Delta y} \delta_i \left( \tilde{h}^{ij} \Delta y^2 A_M \dot{e}_S \right) - \frac{1}{\Delta x} \delta_j \left( h \Delta x^2 A_M \dot{e}_T \right) \right) . \]

References

Griffies, S.M., and Hallberg, R.W., 2000: Biharmonic friction with a Smagorinsky-like viscosity for use in large-scale eddy-permitting ocean models. Monthly Weather Review, 128(8), 2935-2946. https://doi.org/10.1175/1520-0493(2000)128%3C2935:BFWASL%3E2.0.CO;2

Large, W.G., Danabasoglu, G., McWilliams, J.C., Gent, P.R. and Bryan, F.O., 2001: Equatorial circulation of a global ocean climate model with anisotropic horizontal viscosity. Journal of Physical Oceanography, 31(2), pp.518-536. https://doi.org/10.1175/1520-0485(2001)031%3C0518:ECOAGO%3E2.0.CO;2

Smagorinsky, J., 1993: Some historical remarks on the use of nonlinear viscosities. Large eddy simulation of complex engineering and geophysical flows, 1, 69-106.

Smith, R.D., and McWilliams, J.C., 2003: Anisotropic horizontal viscosity for ocean models. Ocean Modelling, 5(2), 129-156. https://doi.org/10.1016/S1463-5003(02)00016-1

Data Types

type  hor_visc_cs
 Control structure for horizontal viscosity. More...
 

Functions/Subroutines

subroutine, public horizontal_viscosity (u, v, h, diffu, diffv, MEKE, VarMix, G, GV, US, CS, OBC, BT, TD, ADp)
 Calculates the acceleration due to the horizontal viscosity. More...
 
subroutine, public hor_visc_init (Time, G, US, param_file, diag, CS, MEKE, ADp)
 Allocates space for and calculates static variables used by horizontal_viscosity(). hor_visc_init calculates and stores the values of a number of metric functions that are used in horizontal_viscosity(). More...
 
subroutine align_aniso_tensor_to_grid (CS, n1, n2)
 Calculates factors in the anisotropic orientation tensor to be align with the grid. With n1=1 and n2=0, this recovers the approach of Large et al, 2001. More...
 
subroutine smooth_gme (CS, G, GME_flux_h, GME_flux_q)
 Apply a 1-1-4-1-1 Laplacian filter one time on GME diffusive flux to reduce any horizontal two-grid-point noise. More...
 
subroutine, public hor_visc_end (CS)
 Deallocates any variables allocated in hor_visc_init. More...
 

Function/Subroutine Documentation

◆ align_aniso_tensor_to_grid()

subroutine mom_hor_visc::align_aniso_tensor_to_grid ( type(hor_visc_cs), pointer  CS,
real, intent(in)  n1,
real, intent(in)  n2 
)
private

Calculates factors in the anisotropic orientation tensor to be align with the grid. With n1=1 and n2=0, this recovers the approach of Large et al, 2001.

Parameters
csControl structure for horizontal viscosity
[in]n1i-component of direction vector [nondim]
[in]n2j-component of direction vector [nondim]

Definition at line 2175 of file MOM_hor_visc.F90.

2176  type(hor_visc_CS), pointer :: CS !< Control structure for horizontal viscosity
2177  real, intent(in) :: n1 !< i-component of direction vector [nondim]
2178  real, intent(in) :: n2 !< j-component of direction vector [nondim]
2179  ! Local variables
2180  real :: recip_n2_norm
2181  ! For normalizing n=(n1,n2) in case arguments are not a unit vector
2182  recip_n2_norm = n1**2 + n2**2
2183  if (recip_n2_norm > 0.) recip_n2_norm = 1./recip_n2_norm
2184  cs%n1n2_h(:,:) = 2. * ( n1 * n2 ) * recip_n2_norm
2185  cs%n1n2_q(:,:) = 2. * ( n1 * n2 ) * recip_n2_norm
2186  cs%n1n1_m_n2n2_h(:,:) = ( n1 * n1 - n2 * n2 ) * recip_n2_norm
2187  cs%n1n1_m_n2n2_q(:,:) = ( n1 * n1 - n2 * n2 ) * recip_n2_norm

◆ hor_visc_end()

subroutine, public mom_hor_visc::hor_visc_end ( type(hor_visc_cs), pointer  CS)

Deallocates any variables allocated in hor_visc_init.

Parameters
csThe control structure returned by a previous call to hor_visc_init.

Definition at line 2252 of file MOM_hor_visc.F90.

2253  type(hor_visc_CS), pointer :: CS !< The control structure returned by a
2254  !! previous call to hor_visc_init.
2255  if (cs%Laplacian .or. cs%biharmonic) then
2256  dealloc_(cs%dx2h) ; dealloc_(cs%dx2q) ; dealloc_(cs%dy2h) ; dealloc_(cs%dy2q)
2257  dealloc_(cs%dx_dyT) ; dealloc_(cs%dy_dxT) ; dealloc_(cs%dx_dyBu) ; dealloc_(cs%dy_dxBu)
2258  dealloc_(cs%reduction_xx) ; dealloc_(cs%reduction_xy)
2259  endif
2260  if (cs%Laplacian) then
2261  dealloc_(cs%Kh_bg_xx) ; dealloc_(cs%Kh_bg_xy)
2262  dealloc_(cs%grid_sp_h2)
2263  if (cs%bound_Kh) then
2264  dealloc_(cs%Kh_Max_xx) ; dealloc_(cs%Kh_Max_xy)
2265  endif
2266  if (cs%Smagorinsky_Kh) then
2267  dealloc_(cs%Laplac2_const_xx) ; dealloc_(cs%Laplac2_const_xy)
2268  endif
2269  if (cs%Leith_Kh) then
2270  dealloc_(cs%Laplac3_const_xx) ; dealloc_(cs%Laplac3_const_xy)
2271  endif
2272  endif
2273  if (cs%biharmonic) then
2274  dealloc_(cs%grid_sp_h3)
2275  dealloc_(cs%Idx2dyCu) ; dealloc_(cs%Idx2dyCv)
2276  dealloc_(cs%Idxdy2u) ; dealloc_(cs%Idxdy2v)
2277  dealloc_(cs%Ah_bg_xx) ; dealloc_(cs%Ah_bg_xy)
2278  if (cs%bound_Ah) then
2279  dealloc_(cs%Ah_Max_xx) ; dealloc_(cs%Ah_Max_xy)
2280  endif
2281  if (cs%Smagorinsky_Ah) then
2282  dealloc_(cs%Biharm6_const_xx) ; dealloc_(cs%Biharm6_const_xy)
2283  endif
2284  if (cs%Leith_Ah) then
2285  dealloc_(cs%Biharm_const_xx) ; dealloc_(cs%Biharm_const_xy)
2286  endif
2287  if (cs%Re_Ah > 0.0) then
2288  dealloc_(cs%Re_Ah_const_xx) ; dealloc_(cs%Re_Ah_const_xy)
2289  endif
2290  endif
2291  if (cs%anisotropic) then
2292  dealloc_(cs%n1n2_h)
2293  dealloc_(cs%n1n2_q)
2294  dealloc_(cs%n1n1_m_n2n2_h)
2295  dealloc_(cs%n1n1_m_n2n2_q)
2296  endif
2297  deallocate(cs)

◆ hor_visc_init()

subroutine, public mom_hor_visc::hor_visc_init ( type(time_type), intent(in)  Time,
type(ocean_grid_type), intent(inout)  G,
type(unit_scale_type), intent(in)  US,
type(param_file_type), intent(in)  param_file,
type(diag_ctrl), intent(inout), target  diag,
type(hor_visc_cs), pointer  CS,
type(meke_type), pointer  MEKE,
type(accel_diag_ptrs), optional, pointer  ADp 
)

Allocates space for and calculates static variables used by horizontal_viscosity(). hor_visc_init calculates and stores the values of a number of metric functions that are used in horizontal_viscosity().

Parameters
[in]timeCurrent model time.
[in,out]gThe ocean's grid structure.
[in]usA dimensional unit scaling type
[in]param_fileA structure to parse for run-time parameters.
[in,out]diagStructure to regulate diagnostic output.
csPointer to the control structure for this module
mekeMEKE data
adpAcceleration diagnostic pointers

Definition at line 1417 of file MOM_hor_visc.F90.

1418  type(time_type), intent(in) :: Time !< Current model time.
1419  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure.
1420  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
1421  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time
1422  !! parameters.
1423  type(diag_ctrl), target, intent(inout) :: diag !< Structure to regulate diagnostic output.
1424  type(hor_visc_CS), pointer :: CS !< Pointer to the control structure for this module
1425  type(MEKE_type), pointer :: MEKE !< MEKE data
1426  type(accel_diag_ptrs), optional, pointer :: ADp !< Acceleration diagnostic pointers
1427  ! Local variables
1428  real, dimension(SZIB_(G),SZJ_(G)) :: u0u, u0v
1429  real, dimension(SZI_(G),SZJB_(G)) :: v0u, v0v
1430  ! u0v is the Laplacian sensitivities to the v velocities
1431  ! at u points [L-2 ~> m-2], with u0u, v0u, and v0v defined similarly.
1432  real :: grid_sp_h2 ! Harmonic mean of the squares of the grid [L2 ~> m2]
1433  real :: grid_sp_h3 ! Harmonic mean of the squares of the grid^(3/2) [L3 ~> m3]
1434  real :: grid_sp_q2 ! spacings at h and q points [L2 ~> m2]
1435  real :: grid_sp_q3 ! spacings at h and q points^(3/2) [L3 ~> m3]
1436  real :: min_grid_sp_h2 ! Minimum value of grid_sp_h2 [L2 ~> m2]
1437  real :: min_grid_sp_h4 ! Minimum value of grid_sp_h2**2 [L4 ~> m4]
1438  real :: Kh_Limit ! A coefficient [T-1 ~> s-1] used, along with the
1439  ! grid spacing, to limit Laplacian viscosity.
1440  real :: fmax ! maximum absolute value of f at the four
1441  ! vorticity points around a thickness point [T-1 ~> s-1]
1442  real :: BoundCorConst ! A constant used when using viscosity to bound the Coriolis accelerations
1443  ! [T2 L-2 ~> s2 m-2]
1444  real :: Ah_Limit ! coefficient [T-1 ~> s-1] used, along with the
1445  ! grid spacing, to limit biharmonic viscosity
1446  real :: Kh ! Lapacian horizontal viscosity [L2 T-1 ~> m2 s-1]
1447  real :: Ah ! biharmonic horizontal viscosity [L4 T-1 ~> m4 s-1]
1448  real :: Kh_vel_scale ! this speed [L T-1 ~> m s-1] times grid spacing gives Lap visc
1449  real :: Ah_vel_scale ! this speed [L T-1 ~> m s-1] times grid spacing cubed gives bih visc
1450  real :: Ah_time_scale ! damping time-scale for biharmonic visc [T ~> s]
1451  real :: Smag_Lap_const ! nondimensional Laplacian Smagorinsky constant
1452  real :: Smag_bi_const ! nondimensional biharmonic Smagorinsky constant
1453  real :: Leith_Lap_const ! nondimensional Laplacian Leith constant
1454  real :: Leith_bi_const ! nondimensional biharmonic Leith constant
1455  real :: dt ! The dynamics time step [T ~> s]
1456  real :: Idt ! The inverse of dt [T-1 ~> s-1]
1457  real :: denom ! work variable; the denominator of a fraction
1458  real :: maxvel ! largest permitted velocity components [m s-1]
1459  real :: bound_Cor_vel ! grid-scale velocity variations at which value
1460  ! the quadratically varying biharmonic viscosity
1461  ! balances Coriolis acceleration [L T-1 ~> m s-1]
1462  real :: Kh_sin_lat ! Amplitude of latitudinally dependent viscosity [L2 T-1 ~> m2 s-1]
1463  real :: Kh_pwr_of_sine ! Power used to raise sin(lat) when using Kh_sin_lat
1464  logical :: bound_Cor_def ! parameter setting of BOUND_CORIOLIS
1465  logical :: get_all ! If true, read and log all parameters, regardless of
1466  ! whether they are used, to enable spell-checking of
1467  ! valid parameters.
1468  logical :: split ! If true, use the split time stepping scheme.
1469  ! If false and USE_GME = True, issue a FATAL error.
1470  logical :: default_2018_answers
1471  character(len=64) :: inputdir, filename
1472  real :: deg2rad ! Converts degrees to radians
1473  real :: slat_fn ! sin(lat)**Kh_pwr_of_sine
1474  real :: aniso_grid_dir(2) ! Vector (n1,n2) for anisotropic direction
1475  integer :: aniso_mode ! Selects the mode for setting the anisotropic direction
1476  integer :: is, ie, js, je, Isq, Ieq, Jsq, Jeq, nz
1477  integer :: isd, ied, jsd, jed, IsdB, IedB, JsdB, JedB
1478  integer :: i, j
1479 ! This include declares and sets the variable "version".
1480 #include "version_variable.h"
1481  character(len=40) :: mdl = "MOM_hor_visc" ! module name
1482  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
1483  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
1484  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
1485  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
1486  if (associated(cs)) then
1487  call mom_error(warning, "hor_visc_init called with an associated "// &
1488  "control structure.")
1489  return
1490  endif
1491  allocate(cs)
1492  cs%diag => diag
1493  ! Read parameters and write them to the model log.
1494  call log_version(param_file, mdl, version, "")
1495  ! It is not clear whether these initialization lines are needed for the
1496  ! cases where the corresponding parameters are not read.
1497  cs%bound_Kh = .false. ; cs%better_bound_Kh = .false. ; cs%Smagorinsky_Kh = .false. ; cs%Leith_Kh = .false.
1498  cs%bound_Ah = .false. ; cs%better_bound_Ah = .false. ; cs%Smagorinsky_Ah = .false. ; cs%Leith_Ah = .false.
1499  cs%use_QG_Leith_visc = .false.
1500  cs%bound_Coriolis = .false.
1501  cs%Modified_Leith = .false.
1502  cs%anisotropic = .false.
1503  cs%dynamic_aniso = .false.
1504  kh = 0.0 ; ah = 0.0
1505  ! If GET_ALL_PARAMS is true, all parameters are read in all cases to enable
1506  ! parameter spelling checks.
1507  call get_param(param_file, mdl, "GET_ALL_PARAMS", get_all, default=.false.)
1508  call get_param(param_file, mdl, "DEFAULT_2018_ANSWERS", default_2018_answers, &
1509  "This sets the default value for the various _2018_ANSWERS parameters.", &
1510  default=.false.)
1511  call get_param(param_file, mdl, "HOR_VISC_2018_ANSWERS", cs%answers_2018, &
1512  "If true, use the order of arithmetic and expressions that recover the "//&
1513  "answers from the end of 2018. Otherwise, use updated and more robust "//&
1514  "forms of the same expressions.", default=default_2018_answers)
1515  call get_param(param_file, mdl, "DEBUG", cs%debug, default=.false.)
1516  call get_param(param_file, mdl, "LAPLACIAN", cs%Laplacian, &
1517  "If true, use a Laplacian horizontal viscosity.", &
1518  default=.false.)
1519  if (cs%Laplacian .or. get_all) then
1520  call get_param(param_file, mdl, "KH", kh, &
1521  "The background Laplacian horizontal viscosity.", &
1522  units = "m2 s-1", default=0.0, scale=us%m_to_L**2*us%T_to_s)
1523  call get_param(param_file, mdl, "KH_BG_MIN", cs%Kh_bg_min, &
1524  "The minimum value allowed for Laplacian horizontal viscosity, KH.", &
1525  units = "m2 s-1", default=0.0, scale=us%m_to_L**2*us%T_to_s)
1526  call get_param(param_file, mdl, "KH_VEL_SCALE", kh_vel_scale, &
1527  "The velocity scale which is multiplied by the grid "//&
1528  "spacing to calculate the Laplacian viscosity. "//&
1529  "The final viscosity is the largest of this scaled "//&
1530  "viscosity, the Smagorinsky and Leith viscosities, and KH.", &
1531  units="m s-1", default=0.0, scale=us%m_s_to_L_T)
1532  call get_param(param_file, mdl, "KH_SIN_LAT", kh_sin_lat, &
1533  "The amplitude of a latitudinally-dependent background "//&
1534  "viscosity of the form KH_SIN_LAT*(SIN(LAT)**KH_PWR_OF_SINE).", &
1535  units = "m2 s-1", default=0.0, scale=us%m_to_L**2*us%T_to_s)
1536  if (kh_sin_lat>0. .or. get_all) &
1537  call get_param(param_file, mdl, "KH_PWR_OF_SINE", kh_pwr_of_sine, &
1538  "The power used to raise SIN(LAT) when using a latitudinally "//&
1539  "dependent background viscosity.", &
1540  units = "nondim", default=4.0)
1541  call get_param(param_file, mdl, "SMAGORINSKY_KH", cs%Smagorinsky_Kh, &
1542  "If true, use a Smagorinsky nonlinear eddy viscosity.", &
1543  default=.false.)
1544  if (cs%Smagorinsky_Kh .or. get_all) &
1545  call get_param(param_file, mdl, "SMAG_LAP_CONST", smag_lap_const, &
1546  "The nondimensional Laplacian Smagorinsky constant, "//&
1547  "often 0.15.", units="nondim", default=0.0, &
1548  fail_if_missing = cs%Smagorinsky_Kh)
1549  call get_param(param_file, mdl, "LEITH_KH", cs%Leith_Kh, &
1550  "If true, use a Leith nonlinear eddy viscosity.", &
1551  default=.false.)
1552  call get_param(param_file, mdl, "MODIFIED_LEITH", cs%Modified_Leith, &
1553  "If true, add a term to Leith viscosity which is "//&
1554  "proportional to the gradient of divergence.", &
1555  default=.false.)
1556  call get_param(param_file, mdl, "RES_SCALE_MEKE_VISC", cs%res_scale_MEKE, &
1557  "If true, the viscosity contribution from MEKE is scaled by "//&
1558  "the resolution function.", default=.false.)
1559  if (cs%Leith_Kh .or. get_all) then
1560  call get_param(param_file, mdl, "LEITH_LAP_CONST", leith_lap_const, &
1561  "The nondimensional Laplacian Leith constant, "//&
1562  "often set to 1.0", units="nondim", default=0.0, &
1563  fail_if_missing = cs%Leith_Kh)
1564  call get_param(param_file, mdl, "USE_QG_LEITH_VISC", cs%use_QG_Leith_visc, &
1565  "If true, use QG Leith nonlinear eddy viscosity.", &
1566  default=.false.)
1567  if (cs%use_QG_Leith_visc .and. .not. cs%Leith_Kh) call mom_error(fatal, &
1568  "MOM_lateral_mixing_coeffs.F90, VarMix_init:"//&
1569  "LEITH_KH must be True when USE_QG_LEITH_VISC=True.")
1570  endif
1571  if (cs%Leith_Kh .or. cs%Leith_Ah .or. get_all) then
1572  call get_param(param_file, mdl, "USE_BETA_IN_LEITH", cs%use_beta_in_Leith, &
1573  "If true, include the beta term in the Leith nonlinear eddy viscosity.", &
1574  default=cs%Leith_Kh)
1575  call get_param(param_file, mdl, "MODIFIED_LEITH", cs%modified_Leith, &
1576  "If true, add a term to Leith viscosity which is \n"//&
1577  "proportional to the gradient of divergence.", &
1578  default=.false.)
1579  endif
1580  call get_param(param_file, mdl, "BOUND_KH", cs%bound_Kh, &
1581  "If true, the Laplacian coefficient is locally limited "//&
1582  "to be stable.", default=.true.)
1583  call get_param(param_file, mdl, "BETTER_BOUND_KH", cs%better_bound_Kh, &
1584  "If true, the Laplacian coefficient is locally limited "//&
1585  "to be stable with a better bounding than just BOUND_KH.", &
1586  default=cs%bound_Kh)
1587  call get_param(param_file, mdl, "ANISOTROPIC_VISCOSITY", cs%anisotropic, &
1588  "If true, allow anistropic viscosity in the Laplacian "//&
1589  "horizontal viscosity.", default=.false.)
1590  call get_param(param_file, mdl, "ADD_LES_VISCOSITY", cs%add_LES_viscosity, &
1591  "If true, adds the viscosity from Smagorinsky and Leith to the "//&
1592  "background viscosity instead of taking the maximum.", default=.false.)
1593  endif
1594  if (cs%anisotropic .or. get_all) then
1595  call get_param(param_file, mdl, "KH_ANISO", cs%Kh_aniso, &
1596  "The background Laplacian anisotropic horizontal viscosity.", &
1597  units = "m2 s-1", default=0.0, scale=us%m_to_L**2*us%T_to_s)
1598  call get_param(param_file, mdl, "ANISOTROPIC_MODE", aniso_mode, &
1599  "Selects the mode for setting the direction of anistropy.\n"//&
1600  "\t 0 - Points along the grid i-direction.\n"//&
1601  "\t 1 - Points towards East.\n"//&
1602  "\t 2 - Points along the flow direction, U/|U|.", &
1603  default=0)
1604  select case (aniso_mode)
1605  case (0)
1606  call get_param(param_file, mdl, "ANISO_GRID_DIR", aniso_grid_dir, &
1607  "The vector pointing in the direction of anistropy for "//&
1608  "horizont viscosity. n1,n2 are the i,j components relative "//&
1609  "to the grid.", units = "nondim", fail_if_missing=.true.)
1610  case (1)
1611  call get_param(param_file, mdl, "ANISO_GRID_DIR", aniso_grid_dir, &
1612  "The vector pointing in the direction of anistropy for "//&
1613  "horizont viscosity. n1,n2 are the i,j components relative "//&
1614  "to the spherical coordinates.", units = "nondim", fail_if_missing=.true.)
1615  end select
1616  endif
1617  call get_param(param_file, mdl, "BIHARMONIC", cs%biharmonic, &
1618  "If true, use a biharmonic horizontal viscosity. "//&
1619  "BIHARMONIC may be used with LAPLACIAN.", &
1620  default=.true.)
1621  if (cs%biharmonic .or. get_all) then
1622  call get_param(param_file, mdl, "AH", ah, &
1623  "The background biharmonic horizontal viscosity.", &
1624  units = "m4 s-1", default=0.0, scale=us%m_to_L**4*us%T_to_s)
1625  call get_param(param_file, mdl, "AH_VEL_SCALE", ah_vel_scale, &
1626  "The velocity scale which is multiplied by the cube of "//&
1627  "the grid spacing to calculate the biharmonic viscosity. "//&
1628  "The final viscosity is the largest of this scaled "//&
1629  "viscosity, the Smagorinsky and Leith viscosities, and AH.", &
1630  units="m s-1", default=0.0, scale=us%m_s_to_L_T)
1631  call get_param(param_file, mdl, "AH_TIME_SCALE", ah_time_scale, &
1632  "A time scale whose inverse is multiplied by the fourth "//&
1633  "power of the grid spacing to calculate biharmonic viscosity. "//&
1634  "The final viscosity is the largest of all viscosity "//&
1635  "formulations in use. 0.0 means that it's not used.", &
1636  units="s", default=0.0, scale=us%s_to_T)
1637  call get_param(param_file, mdl, "SMAGORINSKY_AH", cs%Smagorinsky_Ah, &
1638  "If true, use a biharmonic Smagorinsky nonlinear eddy "//&
1639  "viscosity.", default=.false.)
1640  call get_param(param_file, mdl, "LEITH_AH", cs%Leith_Ah, &
1641  "If true, use a biharmonic Leith nonlinear eddy "//&
1642  "viscosity.", default=.false.)
1643  call get_param(param_file, mdl, "BOUND_AH", cs%bound_Ah, &
1644  "If true, the biharmonic coefficient is locally limited "//&
1645  "to be stable.", default=.true.)
1646  call get_param(param_file, mdl, "BETTER_BOUND_AH", cs%better_bound_Ah, &
1647  "If true, the biharmonic coefficient is locally limited "//&
1648  "to be stable with a better bounding than just BOUND_AH.", &
1649  default=cs%bound_Ah)
1650  call get_param(param_file, mdl, "RE_AH", cs%Re_Ah, &
1651  "If nonzero, the biharmonic coefficient is scaled "//&
1652  "so that the biharmonic Reynolds number is equal to this.", &
1653  units="nondim", default=0.0)
1654 
1655  if (cs%Smagorinsky_Ah .or. get_all) then
1656  call get_param(param_file, mdl, "SMAG_BI_CONST",smag_bi_const, &
1657  "The nondimensional biharmonic Smagorinsky constant, "//&
1658  "typically 0.015 - 0.06.", units="nondim", default=0.0, &
1659  fail_if_missing = cs%Smagorinsky_Ah)
1660  call get_param(param_file, mdl, "BOUND_CORIOLIS", bound_cor_def, default=.false.)
1661  call get_param(param_file, mdl, "BOUND_CORIOLIS_BIHARM", cs%bound_Coriolis, &
1662  "If true use a viscosity that increases with the square "//&
1663  "of the velocity shears, so that the resulting viscous "//&
1664  "drag is of comparable magnitude to the Coriolis terms "//&
1665  "when the velocity differences between adjacent grid "//&
1666  "points is 0.5*BOUND_CORIOLIS_VEL. The default is the "//&
1667  "value of BOUND_CORIOLIS (or false).", default=bound_cor_def)
1668  if (cs%bound_Coriolis .or. get_all) then
1669  call get_param(param_file, mdl, "MAXVEL", maxvel, default=3.0e8)
1670  bound_cor_vel = maxvel
1671  call get_param(param_file, mdl, "BOUND_CORIOLIS_VEL", bound_cor_vel, &
1672  "The velocity scale at which BOUND_CORIOLIS_BIHARM causes "//&
1673  "the biharmonic drag to have comparable magnitude to the "//&
1674  "Coriolis acceleration. The default is set by MAXVEL.", &
1675  units="m s-1", default=maxvel, scale=us%m_s_to_L_T)
1676  endif
1677  endif
1678  if (cs%Leith_Ah .or. get_all) &
1679  call get_param(param_file, mdl, "LEITH_BI_CONST", leith_bi_const, &
1680  "The nondimensional biharmonic Leith constant, "//&
1681  "typical values are thus far undetermined.", units="nondim", default=0.0, &
1682  fail_if_missing = cs%Leith_Ah)
1683  endif
1684  call get_param(param_file, mdl, "USE_LAND_MASK_FOR_HVISC", cs%use_land_mask, &
1685  "If true, use Use the land mask for the computation of thicknesses "//&
1686  "at velocity locations. This eliminates the dependence on arbitrary "//&
1687  "values over land or outside of the domain.", default=.true.)
1688  if (cs%better_bound_Ah .or. cs%better_bound_Kh .or. get_all) &
1689  call get_param(param_file, mdl, "HORVISC_BOUND_COEF", cs%bound_coef, &
1690  "The nondimensional coefficient of the ratio of the "//&
1691  "viscosity bounds to the theoretical maximum for "//&
1692  "stability without considering other terms.", units="nondim", &
1693  default=0.8)
1694  call get_param(param_file, mdl, "NOSLIP", cs%no_slip, &
1695  "If true, no slip boundary conditions are used; otherwise "//&
1696  "free slip boundary conditions are assumed. The "//&
1697  "implementation of the free slip BCs on a C-grid is much "//&
1698  "cleaner than the no slip BCs. The use of free slip BCs "//&
1699  "is strongly encouraged, and no slip BCs are not used with "//&
1700  "the biharmonic viscosity.", default=.false.)
1701  call get_param(param_file, mdl, "USE_KH_BG_2D", cs%use_Kh_bg_2d, &
1702  "If true, read a file containing 2-d background harmonic "//&
1703  "viscosities. The final viscosity is the maximum of the other "//&
1704  "terms and this background value.", default=.false.)
1705 
1706  call get_param(param_file, mdl, "USE_GME", cs%use_GME, &
1707  "If true, use the GM+E backscatter scheme in association \n"//&
1708  "with the Gent and McWilliams parameterization.", default=.false.)
1709  if (cs%use_GME) then
1710  call get_param(param_file, mdl, "SPLIT", split, &
1711  "Use the split time stepping if true.", default=.true., &
1712  do_not_log=.true.)
1713  if (.not. split) call mom_error(fatal,"ERROR: Currently, USE_GME = True "// &
1714  "cannot be used with SPLIT=False.")
1715  call get_param(param_file, mdl, "GME_H0", cs%GME_h0, &
1716  "The strength of GME tapers quadratically to zero when the bathymetric "//&
1717  "depth is shallower than GME_H0.", units="m", scale=us%m_to_Z, &
1718  default=1000.0)
1719  call get_param(param_file, mdl, "GME_EFFICIENCY", cs%GME_efficiency, &
1720  "The nondimensional prefactor multiplying the GME coefficient.", &
1721  units="nondim", default=1.0)
1722  call get_param(param_file, mdl, "GME_LIMITER", cs%GME_limiter, &
1723  "The absolute maximum value the GME coefficient is allowed to take.", &
1724  units="m2 s-1", scale=us%m_to_L**2*us%T_to_s, default=1.0e7)
1725  endif
1726  if (cs%Laplacian .or. cs%biharmonic) then
1727  call get_param(param_file, mdl, "DT", dt, &
1728  "The (baroclinic) dynamics time step.", units="s", scale=us%s_to_T, &
1729  fail_if_missing=.true.)
1730  idt = 1.0 / dt
1731  endif
1732  if (cs%no_slip .and. cs%biharmonic) &
1733  call mom_error(fatal,"ERROR: NOSLIP and BIHARMONIC cannot be defined "// &
1734  "at the same time in MOM.")
1735  if (.not.(cs%Laplacian .or. cs%biharmonic)) then
1736  ! Only issue inviscid warning if not in single column mode (usually 2x2 domain)
1737  if ( max(g%domain%niglobal, g%domain%njglobal)>2 ) call mom_error(warning, &
1738  "hor_visc_init: It is usually a very bad idea not to use either "//&
1739  "LAPLACIAN or BIHARMONIC viscosity.")
1740  return ! We are not using either Laplacian or Bi-harmonic lateral viscosity
1741  endif
1742  deg2rad = atan(1.0) / 45.
1743  alloc_(cs%dx2h(isd:ied,jsd:jed)) ; cs%dx2h(:,:) = 0.0
1744  alloc_(cs%dy2h(isd:ied,jsd:jed)) ; cs%dy2h(:,:) = 0.0
1745  alloc_(cs%dx2q(isdb:iedb,jsdb:jedb)) ; cs%dx2q(:,:) = 0.0
1746  alloc_(cs%dy2q(isdb:iedb,jsdb:jedb)) ; cs%dy2q(:,:) = 0.0
1747  alloc_(cs%dx_dyT(isd:ied,jsd:jed)) ; cs%dx_dyT(:,:) = 0.0
1748  alloc_(cs%dy_dxT(isd:ied,jsd:jed)) ; cs%dy_dxT(:,:) = 0.0
1749  alloc_(cs%dx_dyBu(isdb:iedb,jsdb:jedb)) ; cs%dx_dyBu(:,:) = 0.0
1750  alloc_(cs%dy_dxBu(isdb:iedb,jsdb:jedb)) ; cs%dy_dxBu(:,:) = 0.0
1751  if (cs%Laplacian) then
1752  alloc_(cs%grid_sp_h2(isd:ied,jsd:jed)) ; cs%grid_sp_h2(:,:) = 0.0
1753  alloc_(cs%Kh_bg_xx(isd:ied,jsd:jed)) ; cs%Kh_bg_xx(:,:) = 0.0
1754  alloc_(cs%Kh_bg_xy(isdb:iedb,jsdb:jedb)) ; cs%Kh_bg_xy(:,:) = 0.0
1755  if (cs%bound_Kh .or. cs%better_bound_Kh) then
1756  alloc_(cs%Kh_Max_xx(isd:ied,jsd:jed)) ; cs%Kh_Max_xx(:,:) = 0.0
1757  alloc_(cs%Kh_Max_xy(isdb:iedb,jsdb:jedb)) ; cs%Kh_Max_xy(:,:) = 0.0
1758  endif
1759  if (cs%Smagorinsky_Kh) then
1760  alloc_(cs%Laplac2_const_xx(isd:ied,jsd:jed)) ; cs%Laplac2_const_xx(:,:) = 0.0
1761  alloc_(cs%Laplac2_const_xy(isdb:iedb,jsdb:jedb)) ; cs%Laplac2_const_xy(:,:) = 0.0
1762  endif
1763  if (cs%Leith_Kh) then
1764  alloc_(cs%Laplac3_const_xx(isd:ied,jsd:jed)) ; cs%Laplac3_const_xx(:,:) = 0.0
1765  alloc_(cs%Laplac3_const_xy(isdb:iedb,jsdb:jedb)) ; cs%Laplac3_const_xy(:,:) = 0.0
1766  endif
1767  endif
1768  alloc_(cs%reduction_xx(isd:ied,jsd:jed)) ; cs%reduction_xx(:,:) = 0.0
1769  alloc_(cs%reduction_xy(isdb:iedb,jsdb:jedb)) ; cs%reduction_xy(:,:) = 0.0
1770  if (cs%anisotropic) then
1771  alloc_(cs%n1n2_h(isd:ied,jsd:jed)) ; cs%n1n2_h(:,:) = 0.0
1772  alloc_(cs%n1n1_m_n2n2_h(isd:ied,jsd:jed)) ; cs%n1n1_m_n2n2_h(:,:) = 0.0
1773  alloc_(cs%n1n2_q(isdb:iedb,jsdb:jedb)) ; cs%n1n2_q(:,:) = 0.0
1774  alloc_(cs%n1n1_m_n2n2_q(isdb:iedb,jsdb:jedb)) ; cs%n1n1_m_n2n2_q(:,:) = 0.0
1775  select case (aniso_mode)
1776  case (0)
1777  call align_aniso_tensor_to_grid(cs, aniso_grid_dir(1), aniso_grid_dir(2))
1778  case (1)
1779  ! call align_aniso_tensor_to_grid(CS, aniso_grid_dir(1), aniso_grid_dir(2))
1780  case (2)
1781  cs%dynamic_aniso = .true.
1782  case default
1783  call mom_error(fatal, "MOM_hor_visc: "//&
1784  "Runtime parameter ANISOTROPIC_MODE is out of range.")
1785  end select
1786  endif
1787  if (cs%use_Kh_bg_2d) then
1788  alloc_(cs%Kh_bg_2d(isd:ied,jsd:jed)) ; cs%Kh_bg_2d(:,:) = 0.0
1789  call get_param(param_file, mdl, "KH_BG_2D_FILENAME", filename, &
1790  'The filename containing a 2d map of "Kh".', &
1791  default='KH_background_2d.nc')
1792  call get_param(param_file, mdl, "INPUTDIR", inputdir, default=".")
1793  inputdir = slasher(inputdir)
1794  call mom_read_data(trim(inputdir)//trim(filename), 'Kh', cs%Kh_bg_2d, &
1795  g%domain, timelevel=1, scale=us%m_to_L**2*us%T_to_s)
1796  call pass_var(cs%Kh_bg_2d, g%domain)
1797  endif
1798  if (cs%biharmonic) then
1799  alloc_(cs%Idx2dyCu(isdb:iedb,jsd:jed)) ; cs%Idx2dyCu(:,:) = 0.0
1800  alloc_(cs%Idx2dyCv(isd:ied,jsdb:jedb)) ; cs%Idx2dyCv(:,:) = 0.0
1801  alloc_(cs%Idxdy2u(isdb:iedb,jsd:jed)) ; cs%Idxdy2u(:,:) = 0.0
1802  alloc_(cs%Idxdy2v(isd:ied,jsdb:jedb)) ; cs%Idxdy2v(:,:) = 0.0
1803  alloc_(cs%Ah_bg_xx(isd:ied,jsd:jed)) ; cs%Ah_bg_xx(:,:) = 0.0
1804  alloc_(cs%Ah_bg_xy(isdb:iedb,jsdb:jedb)) ; cs%Ah_bg_xy(:,:) = 0.0
1805  alloc_(cs%grid_sp_h3(isd:ied,jsd:jed)) ; cs%grid_sp_h3(:,:) = 0.0
1806  if (cs%bound_Ah .or. cs%better_bound_Ah) then
1807  alloc_(cs%Ah_Max_xx(isd:ied,jsd:jed)) ; cs%Ah_Max_xx(:,:) = 0.0
1808  alloc_(cs%Ah_Max_xy(isdb:iedb,jsdb:jedb)) ; cs%Ah_Max_xy(:,:) = 0.0
1809  endif
1810  if (cs%Smagorinsky_Ah) then
1811  alloc_(cs%Biharm_const_xx(isd:ied,jsd:jed)) ; cs%Biharm_const_xx(:,:) = 0.0
1812  alloc_(cs%Biharm_const_xy(isdb:iedb,jsdb:jedb)) ; cs%Biharm_const_xy(:,:) = 0.0
1813  if (cs%bound_Coriolis) then
1814  alloc_(cs%Biharm_const2_xx(isd:ied,jsd:jed)) ; cs%Biharm_const2_xx(:,:) = 0.0
1815  alloc_(cs%Biharm_const2_xy(isdb:iedb,jsdb:jedb)) ; cs%Biharm_const2_xy(:,:) = 0.0
1816  endif
1817  endif
1818  if (cs%Leith_Ah) then
1819  alloc_(cs%biharm6_const_xx(isd:ied,jsd:jed)) ; cs%biharm6_const_xx(:,:) = 0.0
1820  alloc_(cs%biharm6_const_xy(isdb:iedb,jsdb:jedb)) ; cs%biharm6_const_xy(:,:) = 0.0
1821  endif
1822  if (cs%Re_Ah > 0.0) then
1823  alloc_(cs%Re_Ah_const_xx(isd:ied,jsd:jed)); cs%Re_Ah_const_xx(:,:) = 0.0
1824  alloc_(cs%Re_Ah_const_xy(isdb:iedb,jsdb:jedb)); cs%Re_Ah_const_xy(:,:) = 0.0
1825  endif
1826  endif
1827  do j=js-2,jeq+1 ; do i=is-2,ieq+1
1828  cs%dx2q(i,j) = g%dxBu(i,j)*g%dxBu(i,j) ; cs%dy2q(i,j) = g%dyBu(i,j)*g%dyBu(i,j)
1829  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)
1830  enddo ; enddo
1831  do j=jsq-1,jeq+2 ; do i=isq-1,ieq+2
1832  cs%dx2h(i,j) = g%dxT(i,j)*g%dxT(i,j) ; cs%dy2h(i,j) = g%dyT(i,j)*g%dyT(i,j)
1833  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)
1834  enddo ; enddo
1835  do j=jsq,jeq+1 ; do i=isq,ieq+1
1836  cs%reduction_xx(i,j) = 1.0
1837  if ((g%dy_Cu(i,j) > 0.0) .and. (g%dy_Cu(i,j) < g%dyCu(i,j)) .and. &
1838  (g%dy_Cu(i,j) < g%dyCu(i,j) * cs%reduction_xx(i,j))) &
1839  cs%reduction_xx(i,j) = g%dy_Cu(i,j) / (g%dyCu(i,j))
1840  if ((g%dy_Cu(i-1,j) > 0.0) .and. (g%dy_Cu(i-1,j) < g%dyCu(i-1,j)) .and. &
1841  (g%dy_Cu(i-1,j) < g%dyCu(i-1,j) * cs%reduction_xx(i,j))) &
1842  cs%reduction_xx(i,j) = g%dy_Cu(i-1,j) / (g%dyCu(i-1,j))
1843  if ((g%dx_Cv(i,j) > 0.0) .and. (g%dx_Cv(i,j) < g%dxCv(i,j)) .and. &
1844  (g%dx_Cv(i,j) < g%dxCv(i,j) * cs%reduction_xx(i,j))) &
1845  cs%reduction_xx(i,j) = g%dx_Cv(i,j) / (g%dxCv(i,j))
1846  if ((g%dx_Cv(i,j-1) > 0.0) .and. (g%dx_Cv(i,j-1) < g%dxCv(i,j-1)) .and. &
1847  (g%dx_Cv(i,j-1) < g%dxCv(i,j-1) * cs%reduction_xx(i,j))) &
1848  cs%reduction_xx(i,j) = g%dx_Cv(i,j-1) / (g%dxCv(i,j-1))
1849  enddo ; enddo
1850  do j=js-1,jeq ; do i=is-1,ieq
1851  cs%reduction_xy(i,j) = 1.0
1852  if ((g%dy_Cu(i,j) > 0.0) .and. (g%dy_Cu(i,j) < g%dyCu(i,j)) .and. &
1853  (g%dy_Cu(i,j) < g%dyCu(i,j) * cs%reduction_xy(i,j))) &
1854  cs%reduction_xy(i,j) = g%dy_Cu(i,j) / (g%dyCu(i,j))
1855  if ((g%dy_Cu(i,j+1) > 0.0) .and. (g%dy_Cu(i,j+1) < g%dyCu(i,j+1)) .and. &
1856  (g%dy_Cu(i,j+1) < g%dyCu(i,j+1) * cs%reduction_xy(i,j))) &
1857  cs%reduction_xy(i,j) = g%dy_Cu(i,j+1) / (g%dyCu(i,j+1))
1858  if ((g%dx_Cv(i,j) > 0.0) .and. (g%dx_Cv(i,j) < g%dxCv(i,j)) .and. &
1859  (g%dx_Cv(i,j) < g%dxCv(i,j) * cs%reduction_xy(i,j))) &
1860  cs%reduction_xy(i,j) = g%dx_Cv(i,j) / (g%dxCv(i,j))
1861  if ((g%dx_Cv(i+1,j) > 0.0) .and. (g%dx_Cv(i+1,j) < g%dxCv(i+1,j)) .and. &
1862  (g%dx_Cv(i+1,j) < g%dxCv(i+1,j) * cs%reduction_xy(i,j))) &
1863  cs%reduction_xy(i,j) = g%dx_Cv(i+1,j) / (g%dxCv(i+1,j))
1864  enddo ; enddo
1865  if (cs%Laplacian) then
1866  ! The 0.3 below was 0.4 in MOM1.10. The change in hq requires
1867  ! this to be less than 1/3, rather than 1/2 as before.
1868  if (cs%bound_Kh .or. cs%bound_Ah) kh_limit = 0.3 / (dt*4.0)
1869  ! Calculate and store the background viscosity at h-points
1870 
1871  min_grid_sp_h2 = huge(1.)
1872  do j=jsq,jeq+1 ; do i=isq,ieq+1
1873  ! Static factors in the Smagorinsky and Leith schemes
1874  grid_sp_h2 = (2.0*cs%dx2h(i,j)*cs%dy2h(i,j)) / (cs%dx2h(i,j) + cs%dy2h(i,j))
1875  cs%grid_sp_h2(i,j) = grid_sp_h2
1876  grid_sp_h3 = grid_sp_h2*sqrt(grid_sp_h2)
1877  if (cs%Smagorinsky_Kh) cs%Laplac2_const_xx(i,j) = smag_lap_const * grid_sp_h2
1878  if (cs%Leith_Kh) cs%Laplac3_const_xx(i,j) = leith_lap_const * grid_sp_h3
1879  ! Maximum of constant background and MICOM viscosity
1880  cs%Kh_bg_xx(i,j) = max(kh, kh_vel_scale * sqrt(grid_sp_h2))
1881  ! Use the larger of the above and values read from a file
1882  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))
1883  ! Use the larger of the above and a function of sin(latitude)
1884  if (kh_sin_lat>0.) then
1885  slat_fn = abs( sin( deg2rad * g%geoLatT(i,j) ) ) ** kh_pwr_of_sine
1886  cs%Kh_bg_xx(i,j) = max(kh_sin_lat * slat_fn, cs%Kh_bg_xx(i,j))
1887  endif
1888  if (cs%bound_Kh .and. .not.cs%better_bound_Kh) then
1889  ! Limit the background viscosity to be numerically stable
1890  cs%Kh_Max_xx(i,j) = kh_limit * grid_sp_h2
1891  cs%Kh_bg_xx(i,j) = min(cs%Kh_bg_xx(i,j), cs%Kh_Max_xx(i,j))
1892  endif
1893  min_grid_sp_h2 = min(grid_sp_h2, min_grid_sp_h2)
1894  enddo ; enddo
1895  call min_across_pes(min_grid_sp_h2)
1896 
1897  ! Calculate and store the background viscosity at q-points
1898  do j=js-1,jeq ; do i=is-1,ieq
1899  ! Static factors in the Smagorinsky and Leith schemes
1900  grid_sp_q2 = (2.0*cs%dx2q(i,j)*cs%dy2q(i,j)) / (cs%dx2q(i,j) + cs%dy2q(i,j))
1901  grid_sp_q3 = grid_sp_q2*sqrt(grid_sp_q2)
1902  if (cs%Smagorinsky_Kh) cs%Laplac2_const_xy(i,j) = smag_lap_const * grid_sp_q2
1903  if (cs%Leith_Kh) cs%Laplac3_const_xy(i,j) = leith_lap_const * grid_sp_q3
1904  ! Maximum of constant background and MICOM viscosity
1905  cs%Kh_bg_xy(i,j) = max(kh, kh_vel_scale * sqrt(grid_sp_q2))
1906  ! Use the larger of the above and values read from a file
1907  if (cs%use_Kh_bg_2d) then
1908  cs%Kh_bg_xy(i,j) = max(cs%Kh_bg_xy(i,j), &
1909  0.25*((cs%Kh_bg_2d(i,j) + cs%Kh_bg_2d(i+1,j+1)) + &
1910  (cs%Kh_bg_2d(i+1,j) + cs%Kh_bg_2d(i,j+1))) )
1911  endif
1912 
1913  ! Use the larger of the above and a function of sin(latitude)
1914  if (kh_sin_lat>0.) then
1915  slat_fn = abs( sin( deg2rad * g%geoLatBu(i,j) ) ) ** kh_pwr_of_sine
1916  cs%Kh_bg_xy(i,j) = max(kh_sin_lat * slat_fn, cs%Kh_bg_xy(i,j))
1917  endif
1918  if (cs%bound_Kh .and. .not.cs%better_bound_Kh) then
1919  ! Limit the background viscosity to be numerically stable
1920  cs%Kh_Max_xy(i,j) = kh_limit * grid_sp_q2
1921  cs%Kh_bg_xy(i,j) = min(cs%Kh_bg_xy(i,j), cs%Kh_Max_xy(i,j))
1922  endif
1923  enddo ; enddo
1924  endif
1925  if (cs%biharmonic) then
1926  do j=js-1,jeq+1 ; do i=isq-1,ieq+1
1927  cs%Idx2dyCu(i,j) = (g%IdxCu(i,j)*g%IdxCu(i,j)) * g%IdyCu(i,j)
1928  cs%Idxdy2u(i,j) = g%IdxCu(i,j) * (g%IdyCu(i,j)*g%IdyCu(i,j))
1929  enddo ; enddo
1930  do j=jsq-1,jeq+1 ; do i=is-1,ieq+1
1931  cs%Idx2dyCv(i,j) = (g%IdxCv(i,j)*g%IdxCv(i,j)) * g%IdyCv(i,j)
1932  cs%Idxdy2v(i,j) = g%IdxCv(i,j) * (g%IdyCv(i,j)*g%IdyCv(i,j))
1933  enddo ; enddo
1934  cs%Ah_bg_xy(:,:) = 0.0
1935  ! The 0.3 below was 0.4 in MOM1.10. The change in hq requires
1936  ! this to be less than 1/3, rather than 1/2 as before.
1937  if (cs%better_bound_Ah .or. cs%bound_Ah) ah_limit = 0.3 / (dt*64.0)
1938  if (cs%Smagorinsky_Ah .and. cs%bound_Coriolis) &
1939  boundcorconst = 1.0 / (5.0*(bound_cor_vel*bound_cor_vel))
1940 
1941  min_grid_sp_h4 = huge(1.)
1942  do j=jsq,jeq+1 ; do i=isq,ieq+1
1943  grid_sp_h2 = (2.0*cs%dx2h(i,j)*cs%dy2h(i,j)) / (cs%dx2h(i,j)+cs%dy2h(i,j))
1944  grid_sp_h3 = grid_sp_h2*sqrt(grid_sp_h2)
1945  cs%grid_sp_h3(i,j) = grid_sp_h3
1946  if (cs%Smagorinsky_Ah) then
1947  cs%Biharm_const_xx(i,j) = smag_bi_const * (grid_sp_h2 * grid_sp_h2)
1948  if (cs%bound_Coriolis) then
1949  fmax = max(abs(g%CoriolisBu(i-1,j-1)), abs(g%CoriolisBu(i,j-1)), &
1950  abs(g%CoriolisBu(i-1,j)), abs(g%CoriolisBu(i,j)))
1951  cs%Biharm_const2_xx(i,j) = (grid_sp_h2 * grid_sp_h2 * grid_sp_h2) * &
1952  (fmax * boundcorconst)
1953  endif
1954  endif
1955  if (cs%Leith_Ah) then
1956  cs%biharm6_const_xx(i,j) = leith_bi_const * (grid_sp_h3 * grid_sp_h3)
1957  endif
1958  cs%Ah_bg_xx(i,j) = max(ah, ah_vel_scale * grid_sp_h2 * sqrt(grid_sp_h2))
1959  if (cs%Re_Ah > 0.0) cs%Re_Ah_const_xx(i,j) = grid_sp_h3 / cs%Re_Ah
1960  if (ah_time_scale > 0.) cs%Ah_bg_xx(i,j) = &
1961  max(cs%Ah_bg_xx(i,j), (grid_sp_h2 * grid_sp_h2) / ah_time_scale)
1962  if (cs%bound_Ah .and. .not.cs%better_bound_Ah) then
1963  cs%Ah_Max_xx(i,j) = ah_limit * (grid_sp_h2 * grid_sp_h2)
1964  cs%Ah_bg_xx(i,j) = min(cs%Ah_bg_xx(i,j), cs%Ah_Max_xx(i,j))
1965  endif
1966  min_grid_sp_h4 = min(grid_sp_h2**2, min_grid_sp_h4)
1967  enddo ; enddo
1968  call min_across_pes(min_grid_sp_h4)
1969 
1970  do j=js-1,jeq ; do i=is-1,ieq
1971  grid_sp_q2 = (2.0*cs%dx2q(i,j)*cs%dy2q(i,j)) / (cs%dx2q(i,j)+cs%dy2q(i,j))
1972  grid_sp_q3 = grid_sp_q2*sqrt(grid_sp_q2)
1973  if (cs%Smagorinsky_Ah) then
1974  cs%Biharm_const_xy(i,j) = smag_bi_const * (grid_sp_q2 * grid_sp_q2)
1975  if (cs%bound_Coriolis) then
1976  cs%Biharm_const2_xy(i,j) = (grid_sp_q2 * grid_sp_q2 * grid_sp_q2) * &
1977  (abs(g%CoriolisBu(i,j)) * boundcorconst)
1978  endif
1979  endif
1980  if (cs%Leith_Ah) then
1981  cs%biharm6_const_xy(i,j) = leith_bi_const * (grid_sp_q3 * grid_sp_q3)
1982  endif
1983  cs%Ah_bg_xy(i,j) = max(ah, ah_vel_scale * grid_sp_q2 * sqrt(grid_sp_q2))
1984  if (cs%Re_Ah > 0.0) cs%Re_Ah_const_xy(i,j) = grid_sp_q3 / cs%Re_Ah
1985  if (ah_time_scale > 0.) cs%Ah_bg_xy(i,j) = &
1986  max(cs%Ah_bg_xy(i,j), (grid_sp_q2 * grid_sp_q2) / ah_time_scale)
1987  if (cs%bound_Ah .and. .not.cs%better_bound_Ah) then
1988  cs%Ah_Max_xy(i,j) = ah_limit * (grid_sp_q2 * grid_sp_q2)
1989  cs%Ah_bg_xy(i,j) = min(cs%Ah_bg_xy(i,j), cs%Ah_Max_xy(i,j))
1990  endif
1991  enddo ; enddo
1992  endif
1993  ! The Laplacian bounds should avoid overshoots when CS%bound_coef < 1.
1994  if (cs%Laplacian .and. cs%better_bound_Kh) then
1995  do j=jsq,jeq+1 ; do i=isq,ieq+1
1996  denom = max( &
1997  (cs%dy2h(i,j) * cs%DY_dxT(i,j) * (g%IdyCu(i,j) + g%IdyCu(i-1,j)) * &
1998  max(g%IdyCu(i,j)*g%IareaCu(i,j), g%IdyCu(i-1,j)*g%IareaCu(i-1,j)) ), &
1999  (cs%dx2h(i,j) * cs%DX_dyT(i,j) * (g%IdxCv(i,j) + g%IdxCv(i,j-1)) * &
2000  max(g%IdxCv(i,j)*g%IareaCv(i,j), g%IdxCv(i,j-1)*g%IareaCv(i,j-1)) ) )
2001  cs%Kh_Max_xx(i,j) = 0.0
2002  if (denom > 0.0) &
2003  cs%Kh_Max_xx(i,j) = cs%bound_coef * 0.25 * idt / denom
2004  enddo ; enddo
2005  do j=js-1,jeq ; do i=is-1,ieq
2006  denom = max( &
2007  (cs%dx2q(i,j) * cs%DX_dyBu(i,j) * (g%IdxCu(i,j+1) + g%IdxCu(i,j)) * &
2008  max(g%IdxCu(i,j)*g%IareaCu(i,j), g%IdxCu(i,j+1)*g%IareaCu(i,j+1)) ), &
2009  (cs%dy2q(i,j) * cs%DY_dxBu(i,j) * (g%IdyCv(i+1,j) + g%IdyCv(i,j)) * &
2010  max(g%IdyCv(i,j)*g%IareaCv(i,j), g%IdyCv(i+1,j)*g%IareaCv(i+1,j)) ) )
2011  cs%Kh_Max_xy(i,j) = 0.0
2012  if (denom > 0.0) &
2013  cs%Kh_Max_xy(i,j) = cs%bound_coef * 0.25 * idt / denom
2014  enddo ; enddo
2015  if (cs%debug) then
2016  call hchksum(cs%Kh_Max_xx, "Kh_Max_xx", g%HI, haloshift=0, scale=us%L_to_m**2*us%s_to_T)
2017  call bchksum(cs%Kh_Max_xy, "Kh_Max_xy", g%HI, haloshift=0, scale=us%L_to_m**2*us%s_to_T)
2018  endif
2019  endif
2020  ! The biharmonic bounds should avoid overshoots when CS%bound_coef < 0.5, but
2021  ! empirically work for CS%bound_coef <~ 1.0
2022  if (cs%biharmonic .and. cs%better_bound_Ah) then
2023  do j=js-1,jeq+1 ; do i=isq-1,ieq+1
2024  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)) + &
2025  cs%dy2h(i,j) * cs%DY_dxT(i,j) * (g%IdyCu(i,j) + g%IdyCu(i-1,j)) ) + &
2026  cs%Idx2dyCu(i,j)*(cs%dx2q(i,j) * cs%DX_dyBu(i,j) * (g%IdxCu(i,j+1) + g%IdxCu(i,j)) + &
2027  cs%dx2q(i,j-1)*cs%DX_dyBu(i,j-1)*(g%IdxCu(i,j) + g%IdxCu(i,j-1)) ) )
2028  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)) + &
2029  cs%dy2h(i,j) * cs%DX_dyT(i,j) * (g%IdxCv(i,j) + g%IdxCv(i,j-1)) ) + &
2030  cs%Idx2dyCu(i,j)*(cs%dx2q(i,j) * cs%DY_dxBu(i,j) * (g%IdyCv(i+1,j) + g%IdyCv(i,j)) + &
2031  cs%dx2q(i,j-1)*cs%DY_dxBu(i,j-1)*(g%IdyCv(i+1,j-1) + g%IdyCv(i,j-1)) ) )
2032  enddo ; enddo
2033  do j=jsq-1,jeq+1 ; do i=is-1,ieq+1
2034  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)) + &
2035  cs%dy2q(i-1,j)*cs%DX_dyBu(i-1,j)*(g%IdxCu(i-1,j+1) + g%IdxCu(i-1,j)) ) + &
2036  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)) + &
2037  cs%dx2h(i,j) * cs%DY_dxT(i,j) * (g%IdyCu(i,j) + g%IdyCu(i-1,j)) ) )
2038  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)) + &
2039  cs%dy2q(i-1,j)*cs%DY_dxBu(i-1,j)*(g%IdyCv(i,j) + g%IdyCv(i-1,j)) ) + &
2040  cs%Idx2dyCv(i,j)*(cs%dx2h(i,j+1)*cs%DX_dyT(i,j+1)*(g%IdxCv(i,j+1) + g%IdxCv(i,j)) + &
2041  cs%dx2h(i,j) * cs%DX_dyT(i,j) * (g%IdxCv(i,j) + g%IdxCv(i,j-1)) ) )
2042  enddo ; enddo
2043  do j=jsq,jeq+1 ; do i=isq,ieq+1
2044  denom = max( &
2045  (cs%dy2h(i,j) * &
2046  (cs%DY_dxT(i,j)*(g%IdyCu(i,j)*u0u(i,j) + g%IdyCu(i-1,j)*u0u(i-1,j)) + &
2047  cs%DX_dyT(i,j)*(g%IdxCv(i,j)*v0u(i,j) + g%IdxCv(i,j-1)*v0u(i,j-1))) * &
2048  max(g%IdyCu(i,j)*g%IareaCu(i,j), g%IdyCu(i-1,j)*g%IareaCu(i-1,j)) ), &
2049  (cs%dx2h(i,j) * &
2050  (cs%DY_dxT(i,j)*(g%IdyCu(i,j)*u0v(i,j) + g%IdyCu(i-1,j)*u0v(i-1,j)) + &
2051  cs%DX_dyT(i,j)*(g%IdxCv(i,j)*v0v(i,j) + g%IdxCv(i,j-1)*v0v(i,j-1))) * &
2052  max(g%IdxCv(i,j)*g%IareaCv(i,j), g%IdxCv(i,j-1)*g%IareaCv(i,j-1)) ) )
2053  cs%Ah_Max_xx(i,j) = 0.0
2054  if (denom > 0.0) &
2055  cs%Ah_Max_xx(i,j) = cs%bound_coef * 0.5 * idt / denom
2056  enddo ; enddo
2057  do j=js-1,jeq ; do i=is-1,ieq
2058  denom = max( &
2059  (cs%dx2q(i,j) * &
2060  (cs%DX_dyBu(i,j)*(u0u(i,j+1)*g%IdxCu(i,j+1) + u0u(i,j)*g%IdxCu(i,j)) + &
2061  cs%DY_dxBu(i,j)*(v0u(i+1,j)*g%IdyCv(i+1,j) + v0u(i,j)*g%IdyCv(i,j))) * &
2062  max(g%IdxCu(i,j)*g%IareaCu(i,j), g%IdxCu(i,j+1)*g%IareaCu(i,j+1)) ), &
2063  (cs%dy2q(i,j) * &
2064  (cs%DX_dyBu(i,j)*(u0v(i,j+1)*g%IdxCu(i,j+1) + u0v(i,j)*g%IdxCu(i,j)) + &
2065  cs%DY_dxBu(i,j)*(v0v(i+1,j)*g%IdyCv(i+1,j) + v0v(i,j)*g%IdyCv(i,j))) * &
2066  max(g%IdyCv(i,j)*g%IareaCv(i,j), g%IdyCv(i+1,j)*g%IareaCv(i+1,j)) ) )
2067  cs%Ah_Max_xy(i,j) = 0.0
2068  if (denom > 0.0) &
2069  cs%Ah_Max_xy(i,j) = cs%bound_coef * 0.5 * idt / denom
2070  enddo ; enddo
2071  if (cs%debug) then
2072  call hchksum(cs%Ah_Max_xx, "Ah_Max_xx", g%HI, haloshift=0, scale=us%L_to_m**4*us%s_to_T)
2073  call bchksum(cs%Ah_Max_xy, "Ah_Max_xy", g%HI, haloshift=0, scale=us%L_to_m**4*us%s_to_T)
2074  endif
2075  endif
2076  ! Register fields for output from this module.
2077  cs%id_diffu = register_diag_field('ocean_model', 'diffu', diag%axesCuL, time, &
2078  'Zonal Acceleration from Horizontal Viscosity', 'm s-2', conversion=us%L_T2_to_m_s2)
2079  cs%id_diffv = register_diag_field('ocean_model', 'diffv', diag%axesCvL, time, &
2080  'Meridional Acceleration from Horizontal Viscosity', 'm s-2', conversion=us%L_T2_to_m_s2)
2081 
2082  !CS%id_hf_diffu = register_diag_field('ocean_model', 'hf_diffu', diag%axesCuL, Time, &
2083  ! 'Fractional Thickness-weighted Zonal Acceleration from Horizontal Viscosity', 'm s-2', &
2084  ! v_extensive=.true., conversion=US%L_T2_to_m_s2)
2085  !if ((CS%id_hf_diffu > 0) .and. (present(ADp))) then
2086  ! call safe_alloc_ptr(CS%hf_diffu,G%IsdB,G%IedB,G%jsd,G%jed,G%ke)
2087  ! call safe_alloc_ptr(ADp%diag_hfrac_u,G%IsdB,G%IedB,G%jsd,G%jed,G%ke)
2088  !endif
2089 
2090  !CS%id_hf_diffv = register_diag_field('ocean_model', 'hf_diffv', diag%axesCvL, Time, &
2091  ! 'Fractional Thickness-weighted Meridional Acceleration from Horizontal Viscosity', 'm s-2', &
2092  ! v_extensive=.true., conversion=US%L_T2_to_m_s2)
2093  !if ((CS%id_hf_diffv > 0) .and. (present(ADp))) then
2094  ! call safe_alloc_ptr(CS%hf_diffv,G%isd,G%ied,G%JsdB,G%JedB,G%ke)
2095  ! call safe_alloc_ptr(ADp%diag_hfrac_v,G%isd,G%ied,G%JsdB,G%JedB,G%ke)
2096  !endif
2097 
2098  cs%id_hf_diffu_2d = register_diag_field('ocean_model', 'hf_diffu_2d', diag%axesCu1, time, &
2099  'Depth-sum Fractional Thickness-weighted Zonal Acceleration from Horizontal Viscosity', 'm s-2', &
2100  conversion=us%L_T2_to_m_s2)
2101  if ((cs%id_hf_diffu_2d > 0) .and. (present(adp))) then
2102  call safe_alloc_ptr(adp%diag_hfrac_u,g%IsdB,g%IedB,g%jsd,g%jed,g%ke)
2103  endif
2104 
2105  cs%id_hf_diffv_2d = register_diag_field('ocean_model', 'hf_diffv_2d', diag%axesCv1, time, &
2106  'Depth-sum Fractional Thickness-weighted Meridional Acceleration from Horizontal Viscosity', 'm s-2', &
2107  conversion=us%L_T2_to_m_s2)
2108  if ((cs%id_hf_diffv_2d > 0) .and. (present(adp))) then
2109  call safe_alloc_ptr(adp%diag_hfrac_v,g%isd,g%ied,g%JsdB,g%JedB,g%ke)
2110  endif
2111 
2112  if (cs%biharmonic) then
2113  cs%id_Ah_h = register_diag_field('ocean_model', 'Ahh', diag%axesTL, time, &
2114  'Biharmonic Horizontal Viscosity at h Points', 'm4 s-1', conversion=us%L_to_m**4*us%s_to_T, &
2115  cmor_field_name='difmxybo', &
2116  cmor_long_name='Ocean lateral biharmonic viscosity', &
2117  cmor_standard_name='ocean_momentum_xy_biharmonic_diffusivity')
2118  cs%id_Ah_q = register_diag_field('ocean_model', 'Ahq', diag%axesBL, time, &
2119  'Biharmonic Horizontal Viscosity at q Points', 'm4 s-1', conversion=us%L_to_m**4*us%s_to_T)
2120  cs%id_grid_Re_Ah = register_diag_field('ocean_model', 'grid_Re_Ah', diag%axesTL, time, &
2121  'Grid Reynolds number for the Biharmonic horizontal viscosity at h points', 'nondim')
2122 
2123  if (cs%id_grid_Re_Ah > 0) &
2124  ! Compute the smallest biharmonic viscosity capable of modifying the
2125  ! velocity at floating point precision.
2126  cs%min_grid_Ah = spacing(1.) * min_grid_sp_h4 * idt
2127  endif
2128  if (cs%Laplacian) then
2129  cs%id_Kh_h = register_diag_field('ocean_model', 'Khh', diag%axesTL, time, &
2130  'Laplacian Horizontal Viscosity at h Points', 'm2 s-1', conversion=us%L_to_m**2*us%s_to_T, &
2131  cmor_field_name='difmxylo', &
2132  cmor_long_name='Ocean lateral Laplacian viscosity', &
2133  cmor_standard_name='ocean_momentum_xy_laplacian_diffusivity')
2134  cs%id_Kh_q = register_diag_field('ocean_model', 'Khq', diag%axesBL, time, &
2135  'Laplacian Horizontal Viscosity at q Points', 'm2 s-1', conversion=us%L_to_m**2*us%s_to_T)
2136  cs%id_grid_Re_Kh = register_diag_field('ocean_model', 'grid_Re_Kh', diag%axesTL, time, &
2137  'Grid Reynolds number for the Laplacian horizontal viscosity at h points', 'nondim')
2138  cs%id_vort_xy_q = register_diag_field('ocean_model', 'vort_xy_q', diag%axesBL, time, &
2139  'Vertical vorticity at q Points', 's-1', conversion=us%s_to_T)
2140  cs%id_div_xx_h = register_diag_field('ocean_model', 'div_xx_h', diag%axesTL, time, &
2141  'Horizontal divergence at h Points', 's-1', conversion=us%s_to_T)
2142  cs%id_sh_xy_q = register_diag_field('ocean_model', 'sh_xy_q', diag%axesBL, time, &
2143  'Shearing strain at q Points', 's-1', conversion=us%s_to_T)
2144  cs%id_sh_xx_h = register_diag_field('ocean_model', 'sh_xx_h', diag%axesTL, time, &
2145  'Horizontal tension at h Points', 's-1', conversion=us%s_to_T)
2146 
2147  if (cs%id_grid_Re_Kh > 0) &
2148  ! Compute a smallest Laplacian viscosity capable of modifying the
2149  ! velocity at floating point precision.
2150  cs%min_grid_Kh = spacing(1.) * min_grid_sp_h2 * idt
2151  endif
2152  if (cs%use_GME) then
2153  cs%id_GME_coeff_h = register_diag_field('ocean_model', 'GME_coeff_h', diag%axesTL, time, &
2154  'GME coefficient at h Points', 'm2 s-1', conversion=us%L_to_m**2*us%s_to_T)
2155  cs%id_GME_coeff_q = register_diag_field('ocean_model', 'GME_coeff_q', diag%axesBL, time, &
2156  'GME coefficient at q Points', 'm2 s-1', conversion=us%L_to_m**2*us%s_to_T)
2157  cs%id_FrictWork_GME = register_diag_field('ocean_model','FrictWork_GME',diag%axesTL,time,&
2158  'Integral work done by lateral friction terms in GME (excluding diffusion of energy)', &
2159  'W m-2', conversion=us%RZ3_T3_to_W_m2*us%L_to_Z**2)
2160  endif
2161  cs%id_FrictWork = register_diag_field('ocean_model','FrictWork',diag%axesTL,time,&
2162  'Integral work done by lateral friction terms', &
2163  'W m-2', conversion=us%RZ3_T3_to_W_m2*us%L_to_Z**2)
2164  cs%id_FrictWorkIntz = register_diag_field('ocean_model','FrictWorkIntz',diag%axesT1,time, &
2165  'Depth integrated work done by lateral friction', &
2166  'W m-2', conversion=us%RZ3_T3_to_W_m2*us%L_to_Z**2, &
2167  cmor_field_name='dispkexyfo', &
2168  cmor_long_name='Depth integrated ocean kinetic energy dissipation due to lateral friction',&
2169  cmor_standard_name='ocean_kinetic_energy_dissipation_per_unit_area_due_to_xy_friction')
2170  if (cs%Laplacian .or. get_all) then
2171  endif

◆ horizontal_viscosity()

subroutine, public mom_hor_visc::horizontal_viscosity ( real, dimension( g %isdb: g %iedb, g %jsd: g %jed, g %ke), intent(in)  u,
real, dimension( g %isd: g %ied, g %jsdb: g %jedb, g %ke), intent(in)  v,
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(inout)  h,
real, dimension( g %isdb: g %iedb, g %jsd: g %jed, g %ke), intent(out)  diffu,
real, dimension( g %isd: g %ied, g %jsdb: g %jedb, g %ke), intent(out)  diffv,
type(meke_type), pointer  MEKE,
type(varmix_cs), pointer  VarMix,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in)  US,
type(hor_visc_cs), pointer  CS,
type(ocean_obc_type), optional, pointer  OBC,
type(barotropic_cs), optional, pointer  BT,
type(thickness_diffuse_cs), optional, pointer  TD,
type(accel_diag_ptrs), optional, pointer  ADp 
)

Calculates the acceleration due to the horizontal viscosity.

A combination of biharmonic and Laplacian forms can be used. The coefficient may either be a constant or a shear-dependent form. The biharmonic is determined by twice taking the divergence of an appropriately defined stress tensor. The Laplacian is determined by doing so once.

To work, the following fields must be set outside of the usual is:ie range before this subroutine is called: u[is-2:ie+2,js-2:je+2] v[is-2:ie+2,js-2:je+2] h[is-1:ie+1,js-1:je+1]

Parameters
[in]gThe ocean's grid structure.
[in]gvThe ocean's vertical grid structure.
[in]uThe zonal velocity [L T-1 ~> m s-1].
[in]vThe meridional velocity [L T-1 ~> m s-1].
[in,out]hLayer thicknesses [H ~> m or kg m-2].
[out]diffuZonal acceleration due to convergence of
[out]diffvMeridional acceleration due to convergence
mekePointer to a structure containing fields related to Mesoscale Eddy Kinetic Energy.
varmixPointer to a structure with fields that specify the spatially variable viscosities
[in]usA dimensional unit scaling type
csControl structure returned by a previous call to hor_visc_init.
obcPointer to an open boundary condition type
btPointer to a structure containing barotropic velocities.
tdPointer to a structure containing thickness diffusivities.
adpAcceleration diagnostic pointers

Definition at line 215 of file MOM_hor_visc.F90.

217  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
218  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
219  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
220  intent(in) :: u !< The zonal velocity [L T-1 ~> m s-1].
221  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
222  intent(in) :: v !< The meridional velocity [L T-1 ~> m s-1].
223  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
224  intent(inout) :: h !< Layer thicknesses [H ~> m or kg m-2].
225  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
226  intent(out) :: diffu !< Zonal acceleration due to convergence of
227  !! along-coordinate stress tensor [L T-2 ~> m s-2]
228  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
229  intent(out) :: diffv !< Meridional acceleration due to convergence
230  !! of along-coordinate stress tensor [L T-2 ~> m s-2].
231  type(MEKE_type), pointer :: MEKE !< Pointer to a structure containing fields
232  !! related to Mesoscale Eddy Kinetic Energy.
233  type(VarMix_CS), pointer :: VarMix !< Pointer to a structure with fields that
234  !! specify the spatially variable viscosities
235  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
236  type(hor_visc_CS), pointer :: CS !< Control structure returned by a previous
237  !! call to hor_visc_init.
238  type(ocean_OBC_type), optional, pointer :: OBC !< Pointer to an open boundary condition type
239  type(barotropic_CS), optional, pointer :: BT !< Pointer to a structure containing
240  !! barotropic velocities.
241  type(thickness_diffuse_CS), optional, pointer :: TD !< Pointer to a structure containing
242  !! thickness diffusivities.
243  type(accel_diag_ptrs), optional, pointer :: ADp !< Acceleration diagnostic pointers
244 
245  ! Local variables
246  real, dimension(SZIB_(G),SZJ_(G)) :: &
247  Del2u, & ! The u-compontent of the Laplacian of velocity [L-1 T-1 ~> m-1 s-1]
248  h_u, & ! Thickness interpolated to u points [H ~> m or kg m-2].
249  vort_xy_dy, & ! y-derivative of vertical vorticity (d/dy(dv/dx - du/dy)) [L-1 T-1 ~> m-1 s-1]
250  div_xx_dx, & ! x-derivative of horizontal divergence (d/dx(du/dx + dv/dy)) [L-1 T-1 ~> m-1 s-1]
251  ubtav ! zonal barotropic vel. ave. over baroclinic time-step [L T-1 ~> m s-1]
252  real, dimension(SZI_(G),SZJB_(G)) :: &
253  Del2v, & ! The v-compontent of the Laplacian of velocity [L-1 T-1 ~> m-1 s-1]
254  h_v, & ! Thickness interpolated to v points [H ~> m or kg m-2].
255  vort_xy_dx, & ! x-derivative of vertical vorticity (d/dx(dv/dx - du/dy)) [L-1 T-1 ~> m-1 s-1]
256  div_xx_dy, & ! y-derivative of horizontal divergence (d/dy(du/dx + dv/dy)) [L-1 T-1 ~> m-1 s-1]
257  vbtav ! meridional barotropic vel. ave. over baroclinic time-step [L T-1 ~> m s-1]
258  real, dimension(SZI_(G),SZJ_(G)) :: &
259  dudx_bt, dvdy_bt, & ! components in the barotropic horizontal tension [T-1 ~> s-1]
260  div_xx, & ! Estimate of horizontal divergence at h-points [T-1 ~> s-1]
261  sh_xx, & ! horizontal tension (du/dx - dv/dy) including metric terms [T-1 ~> s-1]
262  sh_xx_bt, & ! barotropic horizontal tension (du/dx - dv/dy) including metric terms [T-1 ~> s-1]
263  str_xx,& ! str_xx is the diagonal term in the stress tensor [H L2 T-2 ~> m3 s-2 or kg s-2]
264  str_xx_GME,& ! smoothed diagonal term in the stress tensor from GME [H L2 T-2 ~> m3 s-2 or kg s-2]
265  bhstr_xx, & ! A copy of str_xx that only contains the biharmonic contribution [H L2 T-2 ~> m3 s-2 or kg s-2]
266  FrictWorkIntz, & ! depth integrated energy dissipated by lateral friction [R L2 T-3 ~> W m-2]
267  ! Leith_Kh_h, & ! Leith Laplacian viscosity at h-points [L2 T-1 ~> m2 s-1]
268  ! Leith_Ah_h, & ! Leith bi-harmonic viscosity at h-points [L4 T-1 ~> m4 s-1]
269  grad_vort_mag_h, & ! Magnitude of vorticity gradient at h-points [L-1 T-1 ~> m-1 s-1]
270  grad_vort_mag_h_2d, & ! Magnitude of 2d vorticity gradient at h-points [L-1 T-1 ~> m-1 s-1]
271  del2vort_h, & ! Laplacian of vorticity at h-points [L-2 T-1 ~> m-2 s-1]
272  grad_div_mag_h, & ! Magnitude of divergence gradient at h-points [L-1 T-1 ~> m-1 s-1]
273  dudx, dvdy, & ! components in the horizontal tension [T-1 ~> s-1]
274  grad_vel_mag_h, & ! Magnitude of the velocity gradient tensor squared at h-points [T-2 ~> s-2]
275  grad_vel_mag_bt_h, & ! Magnitude of the barotropic velocity gradient tensor squared at h-points [T-2 ~> s-2]
276  grad_d2vel_mag_h, & ! Magnitude of the Laplacian of the velocity vector, squared [L-2 T-2 ~> m-2 s-2]
277  boundary_mask_h ! A mask that zeroes out cells with at least one land edge [nondim]
278 
279  real, allocatable, dimension(:,:) :: hf_diffu_2d ! Depth sum of hf_diffu [L T-2 ~> m s-2]
280  real, allocatable, dimension(:,:) :: hf_diffv_2d ! Depth sum of hf_diffv [L T-2 ~> m s-2]
281 
282  real, dimension(SZIB_(G),SZJB_(G)) :: &
283  dvdx, dudy, & ! components in the shearing strain [T-1 ~> s-1]
284  dDel2vdx, dDel2udy, & ! Components in the biharmonic equivalent of the shearing strain [L-2 T-1 ~> m-2 s-1]
285  dvdx_bt, dudy_bt, & ! components in the barotropic shearing strain [T-1 ~> s-1]
286  sh_xy, & ! horizontal shearing strain (du/dy + dv/dx) including metric terms [T-1 ~> s-1]
287  sh_xy_bt, & ! barotropic horizontal shearing strain (du/dy + dv/dx) inc. metric terms [T-1 ~> s-1]
288  str_xy, & ! str_xy is the cross term in the stress tensor [H L2 T-2 ~> m3 s-2 or kg s-2]
289  str_xy_GME, & ! smoothed cross term in the stress tensor from GME [H L2 T-2 ~> m3 s-2 or kg s-2]
290  bhstr_xy, & ! A copy of str_xy that only contains the biharmonic contribution [H L2 T-2 ~> m3 s-2 or kg s-2]
291  vort_xy, & ! Vertical vorticity (dv/dx - du/dy) including metric terms [T-1 ~> s-1]
292  Leith_Kh_q, & ! Leith Laplacian viscosity at q-points [L2 T-1 ~> m2 s-1]
293  Leith_Ah_q, & ! Leith bi-harmonic viscosity at q-points [L4 T-1 ~> m4 s-1]
294  grad_vort_mag_q, & ! Magnitude of vorticity gradient at q-points [L-1 T-1 ~> m-1 s-1]
295  grad_vort_mag_q_2d, & ! Magnitude of 2d vorticity gradient at q-points [L-1 T-1 ~> m-1 s-1]
296  Del2vort_q, & ! Laplacian of vorticity at q-points [L-2 T-1 ~> m-2 s-1]
297  grad_div_mag_q, & ! Magnitude of divergence gradient at q-points [L-1 T-1 ~> m-1 s-1]
298  grad_vel_mag_q, & ! Magnitude of the velocity gradient tensor squared at q-points [T-2 ~> s-2]
299  hq, & ! harmonic mean of the harmonic means of the u- & v point thicknesses [H ~> m or kg m-2]
300  ! This form guarantees that hq/hu < 4.
301  grad_vel_mag_bt_q, & ! Magnitude of the barotropic velocity gradient tensor squared at q-points [T-2 ~> s-2]
302  boundary_mask_q ! A mask that zeroes out cells with at least one land edge [nondim]
303 
304  real, dimension(SZIB_(G),SZJB_(G),SZK_(G)) :: &
305  Ah_q, & ! biharmonic viscosity at corner points [L4 T-1 ~> m4 s-1]
306  Kh_q, & ! Laplacian viscosity at corner points [L2 T-1 ~> m2 s-1]
307  vort_xy_q, & ! vertical vorticity at corner points [T-1 ~> s-1]
308  sh_xy_q, & ! horizontal shearing strain at corner points [T-1 ~> s-1]
309  GME_coeff_q, & !< GME coeff. at q-points [L2 T-1 ~> m2 s-1]
310  max_diss_rate_q ! maximum possible energy dissipated by lateral friction [L2 T-3 ~> m2 s-3]
311 
312  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)+1) :: &
313  KH_u_GME !< interface height diffusivities in u-columns [L2 T-1 ~> m2 s-1]
314  real, dimension(SZI_(G),SZJB_(G),SZK_(G)+1) :: &
315  KH_v_GME !< interface height diffusivities in v-columns [L2 T-1 ~> m2 s-1]
316  real, dimension(SZI_(G),SZJ_(G),SZK_(G)) :: &
317  Ah_h, & ! biharmonic viscosity at thickness points [L4 T-1 ~> m4 s-1]
318  Kh_h, & ! Laplacian viscosity at thickness points [L2 T-1 ~> m2 s-1]
319  max_diss_rate_h, & ! maximum possible energy dissipated by lateral friction [L2 T-3 ~> m2 s-3]
320  FrictWork, & ! work done by MKE dissipation mechanisms [R L2 T-3 ~> W m-2]
321  FrictWork_GME, & ! work done by GME [R L2 T-3 ~> W m-2]
322  div_xx_h, & ! horizontal divergence [T-1 ~> s-1]
323  sh_xx_h ! horizontal tension (du/dx - dv/dy) including metric terms [T-1 ~> s-1]
324  real, dimension(SZI_(G),SZJ_(G),SZK_(G)) :: &
325  grid_Re_Kh, & !< Grid Reynolds number for Laplacian horizontal viscosity at h points [nondim]
326  grid_Re_Ah, & !< Grid Reynolds number for Biharmonic horizontal viscosity at h points [nondim]
327  GME_coeff_h !< GME coeff. at h-points [L2 T-1 ~> m2 s-1]
328  real :: Ah ! biharmonic viscosity [L4 T-1 ~> m4 s-1]
329  real :: Kh ! Laplacian viscosity [L2 T-1 ~> m2 s-1]
330  real :: AhSm ! Smagorinsky biharmonic viscosity [L4 T-1 ~> m4 s-1]
331  real :: AhLth ! 2D Leith biharmonic viscosity [L4 T-1 ~> m4 s-1]
332  real :: mod_Leith ! nondimensional coefficient for divergence part of modified Leith
333  ! viscosity. Here set equal to nondimensional Laplacian Leith constant.
334  ! This is set equal to zero if modified Leith is not used.
335  real :: Shear_mag ! magnitude of the shear [T-1 ~> s-1]
336  real :: vert_vort_mag ! magnitude of the vertical vorticity gradient [L-1 T-1 ~> m-1 s-1]
337  real :: h2uq, h2vq ! temporary variables [H2 ~> m2 or kg2 m-4].
338  real :: hu, hv ! Thicknesses interpolated by arithmetic means to corner
339  ! points; these are first interpolated to u or v velocity
340  ! points where masks are applied [H ~> m or kg m-2].
341  real :: h_neglect ! thickness so small it can be lost in roundoff and so neglected [H ~> m or kg m-2]
342  real :: h_neglect3 ! h_neglect^3 [H3 ~> m3 or kg3 m-6]
343  real :: hrat_min ! minimum thicknesses at the 4 neighboring
344  ! velocity points divided by the thickness at the stress
345  ! point (h or q point) [nondim]
346  real :: visc_bound_rem ! fraction of overall viscous bounds that
347  ! remain to be applied [nondim]
348  real :: Kh_scale ! A factor between 0 and 1 by which the horizontal
349  ! Laplacian viscosity is rescaled [nondim]
350  real :: RoScl ! The scaling function for MEKE source term [nondim]
351  real :: FatH ! abs(f) at h-point for MEKE source term [T-1 ~> s-1]
352  real :: local_strain ! Local variable for interpolating computed strain rates [T-1 ~> s-1].
353  real :: meke_res_fn ! A copy of the resolution scaling factor if being applied to MEKE. Otherwise =1.
354  real :: GME_coeff ! The GME (negative) viscosity coefficient [L2 T-1 ~> m2 s-1]
355  real :: GME_coeff_limiter ! Maximum permitted value of the GME coefficient [L2 T-1 ~> m2 s-1]
356  real :: FWfrac ! Fraction of maximum theoretical energy transfer to use when scaling GME coefficient [nondim]
357  real :: DY_dxBu ! Ratio of meridional over zonal grid spacing at vertices [nondim]
358  real :: DX_dyBu ! Ratio of zonal over meridiononal grid spacing at vertices [nondim]
359  real :: DY_dxCv ! Ratio of meridional over zonal grid spacing at faces [nondim]
360  real :: DX_dyCu ! Ratio of zonal over meridional grid spacing at faces [nondim]
361  real :: Sh_F_pow ! The ratio of shear over the absolute value of f raised to some power and rescaled [nondim]
362  real :: backscat_subround ! The ratio of f over Shear_mag that is so small that the backscatter
363  ! calculation gives the same value as if f were 0 [nondim].
364  real :: H0_GME ! Depth used to scale down GME coefficient in shallow areas [Z ~> m]
365  real :: KE ! Local kinetic energy [L2 T-2 ~> m2 s-2]
366 
367  logical :: rescale_Kh, legacy_bound
368  logical :: find_FrictWork
369  logical :: apply_OBC = .false.
370  logical :: use_MEKE_Ku
371  logical :: use_MEKE_Au
372  integer :: is, ie, js, je, Isq, Ieq, Jsq, Jeq, nz
373  integer :: i, j, k, n
374  real :: inv_PI3, inv_PI2, inv_PI6
375  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
376  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
377 
378  h_neglect = gv%H_subroundoff
379  h_neglect3 = h_neglect**3
380  inv_pi3 = 1.0/((4.0*atan(1.0))**3)
381  inv_pi2 = 1.0/((4.0*atan(1.0))**2)
382  inv_pi6 = inv_pi3 * inv_pi3
383 
384  ah_h(:,:,:) = 0.0
385  kh_h(:,:,:) = 0.0
386 
387  if (present(obc)) then ; if (associated(obc)) then ; if (obc%OBC_pe) then
388  apply_obc = obc%Flather_u_BCs_exist_globally .or. obc%Flather_v_BCs_exist_globally
389  apply_obc = .true.
390  endif ; endif ; endif
391 
392  if (.not.associated(cs)) call mom_error(fatal, &
393  "MOM_hor_visc: Module must be initialized before it is used.")
394  if (.not.(cs%Laplacian .or. cs%biharmonic)) return
395 
396  find_frictwork = (cs%id_FrictWork > 0)
397  if (cs%id_FrictWorkIntz > 0) find_frictwork = .true.
398  if (associated(meke)) then
399  if (associated(meke%mom_src)) find_frictwork = .true.
400  backscat_subround = 0.0
401  if (find_frictwork .and. associated(meke%mom_src) .and. (meke%backscatter_Ro_c > 0.0) .and. &
402  (meke%backscatter_Ro_Pow /= 0.0)) &
403  backscat_subround = (1.0e-16/meke%backscatter_Ro_c)**(1.0/meke%backscatter_Ro_Pow)
404  endif
405 
406  rescale_kh = .false.
407  if (associated(varmix)) then
408  rescale_kh = varmix%Resoln_scaled_Kh
409  if (rescale_kh .and. &
410  (.not.associated(varmix%Res_fn_h) .or. .not.associated(varmix%Res_fn_q))) &
411  call mom_error(fatal, "MOM_hor_visc: VarMix%Res_fn_h and " //&
412  "VarMix%Res_fn_q both need to be associated with Resoln_scaled_Kh.")
413  endif
414  legacy_bound = (cs%Smagorinsky_Kh .or. cs%Leith_Kh) .and. &
415  (cs%bound_Kh .and. .not.cs%better_bound_Kh)
416 
417  ! Toggle whether to use a Laplacian viscosity derived from MEKE
418  use_meke_ku = associated(meke%Ku)
419  use_meke_au = associated(meke%Au)
420 
421  if (cs%use_GME) then
422  do j=jsq-1,jeq+2 ; do i=isq-1,ieq+2
423  boundary_mask_h(i,j) = (g%mask2dCu(i,j) * g%mask2dCv(i,j) * g%mask2dCu(i-1,j) * g%mask2dCv(i,j-1))
424  enddo ; enddo
425 
426  do j=js-2,jeq+1 ; do i=is-2,ieq+1
427  boundary_mask_q(i,j) = (g%mask2dCv(i,j) * g%mask2dCv(i+1,j) * g%mask2dCu(i,j) * g%mask2dCu(i,j-1))
428  enddo; enddo
429 
430  ! initialize diag. array with zeros
431  gme_coeff_h(:,:,:) = 0.0
432  gme_coeff_q(:,:,:) = 0.0
433  str_xx_gme(:,:) = 0.0
434  str_xy_gme(:,:) = 0.0
435 
436  ! Get barotropic velocities and their gradients
437  call barotropic_get_tav(bt, ubtav, vbtav, g, us)
438  call pass_vector(ubtav, vbtav, g%Domain)
439 
440  do j=js-1,je+2 ; do i=is-1,ie+2
441  dudx_bt(i,j) = cs%DY_dxT(i,j)*(g%IdyCu(i,j) * ubtav(i,j) - &
442  g%IdyCu(i-1,j) * ubtav(i-1,j))
443  dvdy_bt(i,j) = cs%DX_dyT(i,j)*(g%IdxCv(i,j) * vbtav(i,j) - &
444  g%IdxCv(i,j-1) * vbtav(i,j-1))
445  enddo; enddo
446 
447  do j=jsq,jeq+1 ; do i=isq,ieq+1
448  sh_xx_bt(i,j) = dudx_bt(i,j) - dvdy_bt(i,j)
449  enddo ; enddo
450 
451  ! Components for the barotropic shearing strain
452  do j=js-2,jeq+1 ; do i=is-2,ieq+1
453  dvdx_bt(i,j) = cs%DY_dxBu(i,j)*(vbtav(i+1,j)*g%IdyCv(i+1,j) &
454  - vbtav(i,j)*g%IdyCv(i,j))
455  dudy_bt(i,j) = cs%DX_dyBu(i,j)*(ubtav(i,j+1)*g%IdxCu(i,j+1) &
456  - ubtav(i,j)*g%IdxCu(i,j))
457  enddo ; enddo
458 
459  call pass_vector(dudx_bt, dvdy_bt, g%Domain, stagger=bgrid_ne)
460  call pass_vector(dvdx_bt, dudy_bt, g%Domain, stagger=agrid)
461 
462  if (cs%no_slip) then
463  do j=js-1,jeq ; do i=is-1,ieq
464  sh_xy_bt(i,j) = (2.0-g%mask2dBu(i,j)) * ( dvdx_bt(i,j) + dudy_bt(i,j) )
465  enddo ; enddo
466  else
467  do j=js-1,jeq ; do i=is-1,ieq
468  sh_xy_bt(i,j) = g%mask2dBu(i,j) * ( dvdx_bt(i,j) + dudy_bt(i,j) )
469  enddo ; enddo
470  endif
471 
472  do j=jsq-1,jeq+2 ; do i=isq-1,ieq+2
473  grad_vel_mag_bt_h(i,j) = boundary_mask_h(i,j) * (dudx_bt(i,j)**2 + dvdy_bt(i,j)**2 + &
474  (0.25*((dvdx_bt(i,j)+dvdx_bt(i-1,j-1))+(dvdx_bt(i,j-1)+dvdx_bt(i-1,j))))**2 + &
475  (0.25*((dudy_bt(i,j)+dudy_bt(i-1,j-1))+(dudy_bt(i,j-1)+dudy_bt(i-1,j))))**2)
476  enddo ; enddo
477 
478  do j=js-2,jeq+1 ; do i=is-2,ieq+1
479  grad_vel_mag_bt_q(i,j) = boundary_mask_q(i,j) * (dvdx_bt(i,j)**2 + dudy_bt(i,j)**2 + &
480  (0.25*((dudx_bt(i,j)+dudx_bt(i+1,j+1))+(dudx_bt(i,j+1)+dudx_bt(i+1,j))))**2 + &
481  (0.25*((dvdy_bt(i,j)+dvdy_bt(i+1,j+1))+(dvdy_bt(i,j+1)+dvdy_bt(i+1,j))))**2)
482  enddo ; enddo
483 
484  endif ! use_GME
485 
486  !$OMP parallel do default(none) &
487  !$OMP shared( &
488  !$OMP CS, G, GV, US, OBC, VarMix, MEKE, u, v, h, &
489  !$OMP is, ie, js, je, Isq, Ieq, Jsq, Jeq, nz, &
490  !$OMP apply_OBC, rescale_Kh, legacy_bound, find_FrictWork, &
491  !$OMP use_MEKE_Ku, use_MEKE_Au, boundary_mask_h, boundary_mask_q, &
492  !$OMP backscat_subround, GME_coeff_limiter, &
493  !$OMP h_neglect, h_neglect3, FWfrac, inv_PI3, inv_PI6, H0_GME, &
494  !$OMP diffu, diffv, max_diss_rate_h, max_diss_rate_q, &
495  !$OMP Kh_h, Kh_q, Ah_h, Ah_q, FrictWork, FrictWork_GME, &
496  !$OMP div_xx_h, sh_xx_h, vort_xy_q, sh_xy_q, GME_coeff_h, GME_coeff_q, &
497  !$OMP TD, KH_u_GME, KH_v_GME, grid_Re_Kh, grid_Re_Ah &
498  !$OMP ) &
499  !$OMP private( &
500  !$OMP i, j, k, n, &
501  !$OMP dudx, dudy, dvdx, dvdy, sh_xx, sh_xy, h_u, h_v, &
502  !$OMP Del2u, Del2v, DY_dxBu, DX_dyBu, sh_xx_bt, sh_xy_bt, &
503  !$OMP str_xx, str_xy, bhstr_xx, bhstr_xy, str_xx_GME, str_xy_GME, &
504  !$OMP vort_xy, vort_xy_dx, vort_xy_dy, div_xx, div_xx_dx, div_xx_dy, &
505  !$OMP grad_div_mag_h, grad_div_mag_q, grad_vort_mag_h, grad_vort_mag_q, &
506  !$OMP grad_vort_mag_h_2d, grad_vort_mag_q_2d, &
507  !$OMP grad_vel_mag_h, grad_vel_mag_q, &
508  !$OMP grad_vel_mag_bt_h, grad_vel_mag_bt_q, grad_d2vel_mag_h, &
509  !$OMP meke_res_fn, Shear_mag, vert_vort_mag, hrat_min, visc_bound_rem, &
510  !$OMP Kh, Ah, AhSm, AhLth, local_strain, Sh_F_pow, &
511  !$OMP dDel2vdx, dDel2udy, DY_dxCv, DX_dyCu, Del2vort_q, Del2vort_h, KE, &
512  !$OMP h2uq, h2vq, hu, hv, hq, FatH, RoScl, GME_coeff &
513  !$OMP )
514  do k=1,nz
515 
516  ! The following are the forms of the horizontal tension and horizontal
517  ! shearing strain advocated by Smagorinsky (1993) and discussed in
518  ! Griffies and Hallberg (2000).
519 
520  ! Calculate horizontal tension
521  do j=jsq-1,jeq+2 ; do i=isq-1,ieq+2
522  dudx(i,j) = cs%DY_dxT(i,j)*(g%IdyCu(i,j) * u(i,j,k) - &
523  g%IdyCu(i-1,j) * u(i-1,j,k))
524  dvdy(i,j) = cs%DX_dyT(i,j)*(g%IdxCv(i,j) * v(i,j,k) - &
525  g%IdxCv(i,j-1) * v(i,j-1,k))
526  sh_xx(i,j) = dudx(i,j) - dvdy(i,j)
527  enddo ; enddo
528 
529  ! Components for the shearing strain
530  do j=jsq-2,jeq+2 ; do i=isq-2,ieq+2
531  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))
532  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))
533  enddo ; enddo
534 
535  ! Interpolate the thicknesses to velocity points.
536  ! The extra wide halos are to accommodate the cross-corner-point projections
537  ! in OBCs, which are not ordinarily be necessary, and might not be necessary
538  ! even with OBCs if the accelerations are zeroed at OBC points, in which
539  ! case the j-loop for h_u could collapse to j=js=1,je+1. -RWH
540  if (cs%use_land_mask) then
541  do j=js-2,je+2 ; do i=isq-1,ieq+1
542  h_u(i,j) = 0.5 * (g%mask2dT(i,j)*h(i,j,k) + g%mask2dT(i+1,j)*h(i+1,j,k))
543  enddo ; enddo
544  do j=jsq-1,jeq+1 ; do i=is-2,ie+2
545  h_v(i,j) = 0.5 * (g%mask2dT(i,j)*h(i,j,k) + g%mask2dT(i,j+1)*h(i,j+1,k))
546  enddo ; enddo
547  else
548  do j=js-2,je+2 ; do i=isq-1,ieq+1
549  h_u(i,j) = 0.5 * (h(i,j,k) + h(i+1,j,k))
550  enddo ; enddo
551  do j=jsq-1,jeq+1 ; do i=is-2,ie+2
552  h_v(i,j) = 0.5 * (h(i,j,k) + h(i,j+1,k))
553  enddo ; enddo
554  endif
555 
556  ! Adjust contributions to shearing strain and interpolated values of
557  ! thicknesses on open boundaries.
558  if (apply_obc) then ; do n=1,obc%number_of_segments
559  j = obc%segment(n)%HI%JsdB ; i = obc%segment(n)%HI%IsdB
560  if (obc%zero_strain .or. obc%freeslip_strain .or. obc%computed_strain) then
561  if (obc%segment(n)%is_N_or_S .and. (j >= js-2) .and. (j <= jeq+1)) then
562  do i=obc%segment(n)%HI%IsdB,obc%segment(n)%HI%IedB
563  if (obc%zero_strain) then
564  dvdx(i,j) = 0. ; dudy(i,j) = 0.
565  elseif (obc%freeslip_strain) then
566  dudy(i,j) = 0.
567  elseif (obc%computed_strain) then
568  if (obc%segment(n)%direction == obc_direction_n) then
569  dudy(i,j) = 2.0*cs%DX_dyBu(i,j)* &
570  (obc%segment(n)%tangential_vel(i,j,k) - u(i,j,k))*g%IdxCu(i,j)
571  else
572  dudy(i,j) = 2.0*cs%DX_dyBu(i,j)* &
573  (u(i,j+1,k) - obc%segment(n)%tangential_vel(i,j,k))*g%IdxCu(i,j+1)
574  endif
575  elseif (obc%specified_strain) then
576  if (obc%segment(n)%direction == obc_direction_n) then
577  dudy(i,j) = cs%DX_dyBu(i,j)*obc%segment(n)%tangential_grad(i,j,k)*g%IdxCu(i,j)*g%dxBu(i,j)
578  else
579  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)
580  endif
581  endif
582  enddo
583  elseif (obc%segment(n)%is_E_or_W .and. (i >= is-2) .and. (i <= ieq+1)) then
584  do j=obc%segment(n)%HI%JsdB,obc%segment(n)%HI%JedB
585  if (obc%zero_strain) then
586  dvdx(i,j) = 0. ; dudy(i,j) = 0.
587  elseif (obc%freeslip_strain) then
588  dvdx(i,j) = 0.
589  elseif (obc%computed_strain) then
590  if (obc%segment(n)%direction == obc_direction_e) then
591  dvdx(i,j) = 2.0*cs%DY_dxBu(i,j)* &
592  (obc%segment(n)%tangential_vel(i,j,k) - v(i,j,k))*g%IdyCv(i,j)
593  else
594  dvdx(i,j) = 2.0*cs%DY_dxBu(i,j)* &
595  (v(i+1,j,k) - obc%segment(n)%tangential_vel(i,j,k))*g%IdyCv(i+1,j)
596  endif
597  elseif (obc%specified_strain) then
598  if (obc%segment(n)%direction == obc_direction_e) then
599  dvdx(i,j) = cs%DY_dxBu(i,j)*obc%segment(n)%tangential_grad(i,j,k)*g%IdyCv(i,j)*g%dxBu(i,j)
600  else
601  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)
602  endif
603  endif
604  enddo
605  endif
606  endif
607 
608 
609  if (obc%segment(n)%direction == obc_direction_n) then
610  ! There are extra wide halos here to accommodate the cross-corner-point
611  ! OBC projections, but they might not be necessary if the accelerations
612  ! are always zeroed out at OBC points, in which case the i-loop below
613  ! becomes do i=is-1,ie+1. -RWH
614  if ((j >= jsq-1) .and. (j <= jeq+1)) then
615  do i = max(is-2,obc%segment(n)%HI%isd), min(ie+2,obc%segment(n)%HI%ied)
616  h_v(i,j) = h(i,j,k)
617  enddo
618  endif
619  elseif (obc%segment(n)%direction == obc_direction_s) then
620  if ((j >= jsq-1) .and. (j <= jeq+1)) then
621  do i = max(is-2,obc%segment(n)%HI%isd), min(ie+2,obc%segment(n)%HI%ied)
622  h_v(i,j) = h(i,j+1,k)
623  enddo
624  endif
625  elseif (obc%segment(n)%direction == obc_direction_e) then
626  if ((i >= isq-1) .and. (i <= ieq+1)) then
627  do j = max(js-2,obc%segment(n)%HI%jsd), min(je+2,obc%segment(n)%HI%jed)
628  h_u(i,j) = h(i,j,k)
629  enddo
630  endif
631  elseif (obc%segment(n)%direction == obc_direction_w) then
632  if ((i >= isq-1) .and. (i <= ieq+1)) then
633  do j = max(js-2,obc%segment(n)%HI%jsd), min(je+2,obc%segment(n)%HI%jed)
634  h_u(i,j) = h(i+1,j,k)
635  enddo
636  endif
637  endif
638  enddo ; endif
639  ! Now project thicknesses across corner points on OBCs.
640  if (apply_obc) then ; do n=1,obc%number_of_segments
641  j = obc%segment(n)%HI%JsdB ; i = obc%segment(n)%HI%IsdB
642  if (obc%segment(n)%direction == obc_direction_n) then
643  if ((j >= js-2) .and. (j <= je)) then
644  do i = max(isq-1,obc%segment(n)%HI%IsdB), min(ieq+1,obc%segment(n)%HI%IedB)
645  h_u(i,j+1) = h_u(i,j)
646  enddo
647  endif
648  elseif (obc%segment(n)%direction == obc_direction_s) then
649  if ((j >= js-1) .and. (j <= je+1)) then
650  do i = max(isq-1,obc%segment(n)%HI%isd), min(ieq+1,obc%segment(n)%HI%ied)
651  h_u(i,j) = h_u(i,j+1)
652  enddo
653  endif
654  elseif (obc%segment(n)%direction == obc_direction_e) then
655  if ((i >= is-2) .and. (i <= ie)) then
656  do j = max(jsq-1,obc%segment(n)%HI%jsd), min(jeq+1,obc%segment(n)%HI%jed)
657  h_v(i+1,j) = h_v(i,j)
658  enddo
659  endif
660  elseif (obc%segment(n)%direction == obc_direction_w) then
661  if ((i >= is-1) .and. (i <= ie+1)) then
662  do j = max(jsq-1,obc%segment(n)%HI%jsd), min(jeq+1,obc%segment(n)%HI%jed)
663  h_v(i,j) = h_v(i+1,j)
664  enddo
665  endif
666  endif
667  enddo ; endif
668 
669  ! Shearing strain (including no-slip boundary conditions at the 2-D land-sea mask).
670  ! dudy and dvdx include modifications at OBCs from above.
671  if (cs%no_slip) then
672  do j=js-2,jeq+1 ; do i=is-2,ieq+1
673  sh_xy(i,j) = (2.0-g%mask2dBu(i,j)) * ( dvdx(i,j) + dudy(i,j) )
674  enddo ; enddo
675  else
676  do j=js-2,jeq+1 ; do i=is-2,ieq+1
677  sh_xy(i,j) = g%mask2dBu(i,j) * ( dvdx(i,j) + dudy(i,j) )
678  enddo ; enddo
679  endif
680 
681  ! Evaluate Del2u = x.Div(Grad u) and Del2v = y.Div( Grad u)
682  if (cs%biharmonic) then
683  do j=js-1,jeq+1 ; do i=isq-1,ieq+1
684  del2u(i,j) = cs%Idxdy2u(i,j)*(cs%dy2h(i+1,j)*sh_xx(i+1,j) - cs%dy2h(i,j)*sh_xx(i,j)) + &
685  cs%Idx2dyCu(i,j)*(cs%dx2q(i,j)*sh_xy(i,j) - cs%dx2q(i,j-1)*sh_xy(i,j-1))
686  enddo ; enddo
687  do j=jsq-1,jeq+1 ; do i=is-1,ieq+1
688  del2v(i,j) = cs%Idxdy2v(i,j)*(cs%dy2q(i,j)*sh_xy(i,j) - cs%dy2q(i-1,j)*sh_xy(i-1,j)) - &
689  cs%Idx2dyCv(i,j)*(cs%dx2h(i,j+1)*sh_xx(i,j+1) - cs%dx2h(i,j)*sh_xx(i,j))
690  enddo ; enddo
691  if (apply_obc) then; if (obc%zero_biharmonic) then
692  do n=1,obc%number_of_segments
693  i = obc%segment(n)%HI%IsdB ; j = obc%segment(n)%HI%JsdB
694  if (obc%segment(n)%is_N_or_S .and. (j >= jsq-1) .and. (j <= jeq+1)) then
695  do i=obc%segment(n)%HI%isd,obc%segment(n)%HI%ied
696  del2v(i,j) = 0.
697  enddo
698  elseif (obc%segment(n)%is_E_or_W .and. (i >= isq-1) .and. (i <= ieq+1)) then
699  do j=obc%segment(n)%HI%jsd,obc%segment(n)%HI%jed
700  del2u(i,j) = 0.
701  enddo
702  endif
703  enddo
704  endif; endif
705  endif
706 
707  ! Vorticity
708  if (cs%no_slip) then
709  do j=jsq-2,jeq+2 ; do i=isq-2,ieq+2
710  vort_xy(i,j) = (2.0-g%mask2dBu(i,j)) * ( dvdx(i,j) - dudy(i,j) )
711  enddo ; enddo
712  else
713  do j=jsq-2,jeq+2 ; do i=isq-2,ieq+2
714  vort_xy(i,j) = g%mask2dBu(i,j) * ( dvdx(i,j) - dudy(i,j) )
715  enddo ; enddo
716  endif
717 
718  ! Divergence
719  do j=jsq-1,jeq+2 ; do i=isq-1,ieq+2
720  div_xx(i,j) = dudx(i,j) + dvdy(i,j)
721  enddo ; enddo
722 
723  if ((cs%Leith_Kh) .or. (cs%Leith_Ah)) then
724 
725  ! Vorticity gradient
726  do j=jsq-1,jeq+1 ; do i=isq-1,ieq+2
727  dy_dxbu = g%dyBu(i,j) * g%IdxBu(i,j)
728  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))
729  enddo ; enddo
730 
731  do j=jsq-1,jeq+2 ; do i=isq-1,ieq+1
732  dx_dybu = g%dxBu(i,j) * g%IdyBu(i,j)
733  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))
734  enddo ; enddo
735 
736  ! Laplacian of vorticity
737  do j=jsq-1,jeq+1 ; do i=isq-1,ieq+1
738  dy_dxbu = g%dyBu(i,j) * g%IdxBu(i,j)
739  dx_dybu = g%dxBu(i,j) * g%IdyBu(i,j)
740 
741  del2vort_q(i,j) = dy_dxbu * (vort_xy_dx(i+1,j) * g%IdyCv(i+1,j) - vort_xy_dx(i,j) * g%IdyCv(i,j)) + &
742  dx_dybu * (vort_xy_dy(i,j+1) * g%IdyCu(i,j+1) - vort_xy_dy(i,j) * g%IdyCu(i,j))
743  enddo ; enddo
744  do j=jsq,jeq+1 ; do i=isq,ieq+1
745  del2vort_h(i,j) = 0.25*(del2vort_q(i,j) + del2vort_q(i-1,j) + del2vort_q(i,j-1) + del2vort_q(i-1,j-1))
746  enddo ; enddo
747 
748  if (cs%modified_Leith) then
749 
750  ! Divergence gradient
751  do j=jsq-1,jeq+2 ; do i=isq-1,ieq+1
752  div_xx_dx(i,j) = g%IdxCu(i,j)*(div_xx(i+1,j) - div_xx(i,j))
753  enddo ; enddo
754  do j=jsq-1,jeq+1 ; do i=isq-1,ieq+2
755  div_xx_dy(i,j) = g%IdyCv(i,j)*(div_xx(i,j+1) - div_xx(i,j))
756  enddo ; enddo
757 
758  ! Magnitude of divergence gradient
759  do j=jsq,jeq+1 ; do i=isq,ieq+1
760  grad_div_mag_h(i,j) =sqrt((0.5*(div_xx_dx(i,j) + div_xx_dx(i-1,j)))**2 + &
761  (0.5 * (div_xx_dy(i,j) + div_xx_dy(i,j-1)))**2)
762  enddo ; enddo
763  !do J=js-1,Jeq ; do I=is-1,Ieq
764  do j=jsq-1,jeq+1 ; do i=isq-1,ieq+1
765  grad_div_mag_q(i,j) =sqrt((0.5*(div_xx_dx(i,j) + div_xx_dx(i,j+1)))**2 + &
766  (0.5 * (div_xx_dy(i,j) + div_xx_dy(i+1,j)))**2)
767  enddo ; enddo
768 
769  else
770 
771  do j=jsq-1,jeq+2 ; do i=is-2,ieq+1
772  div_xx_dx(i,j) = 0.0
773  enddo ; enddo
774  do j=jsq-1,jeq+1 ; do i=isq-1,ieq+2
775  div_xx_dy(i,j) = 0.0
776  enddo ; enddo
777  do j=jsq,jeq+1 ; do i=isq,ieq+1
778  grad_div_mag_h(i,j) = 0.0
779  enddo ; enddo
780  do j=jsq-1,jeq+1 ; do i=isq-1,ieq+1
781  grad_div_mag_q(i,j) = 0.0
782  enddo ; enddo
783 
784  endif ! CS%modified_Leith
785 
786  ! Add in beta for the Leith viscosity
787  if (cs%use_beta_in_Leith) then
788  do j=js-2,jeq+1 ; do i=is-1,ieq+1
789  vort_xy_dx(i,j) = vort_xy_dx(i,j) + 0.5 * ( g%dF_dx(i,j) + g%dF_dx(i,j+1))
790  enddo ; enddo
791  do j=js-1,jeq+1 ; do i=is-2,ieq+1
792  vort_xy_dy(i,j) = vort_xy_dy(i,j) + 0.5 * ( g%dF_dy(i,j) + g%dF_dy(i+1,j))
793  enddo ; enddo
794  endif ! CS%use_beta_in_Leith
795 
796  if (cs%use_QG_Leith_visc) then
797 
798  do j=jsq,jeq+1 ; do i=isq,ieq+1
799  grad_vort_mag_h_2d(i,j) = sqrt((0.5*(vort_xy_dx(i,j) + vort_xy_dx(i,j-1)))**2 + &
800  (0.5*(vort_xy_dy(i,j) + vort_xy_dy(i-1,j)))**2 )
801  enddo ; enddo
802  do j=js-1,jeq ; do i=is-1,ieq
803  grad_vort_mag_q_2d(i,j) = sqrt((0.5*(vort_xy_dx(i,j) + vort_xy_dx(i+1,j)))**2 + &
804  (0.5*(vort_xy_dy(i,j) + vort_xy_dy(i,j+1)))**2 )
805  enddo ; enddo
806 
807  ! This accumulates terms, some of which are in VarMix, so rescaling can not be done here.
808  call calc_qg_leith_viscosity(varmix, g, gv, us, h, k, div_xx_dx, div_xx_dy, &
809  vort_xy_dx, vort_xy_dy)
810 
811  endif
812 
813  do j=jsq,jeq+1 ; do i=isq,ieq+1
814  grad_vort_mag_h(i,j) = sqrt((0.5*(vort_xy_dx(i,j) + vort_xy_dx(i,j-1)))**2 + &
815  (0.5*(vort_xy_dy(i,j) + vort_xy_dy(i-1,j)))**2 )
816  enddo ; enddo
817  do j=js-1,jeq ; do i=is-1,ieq
818  grad_vort_mag_q(i,j) = sqrt((0.5*(vort_xy_dx(i,j) + vort_xy_dx(i+1,j)))**2 + &
819  (0.5*(vort_xy_dy(i,j) + vort_xy_dy(i,j+1)))**2 )
820  enddo ; enddo
821 
822  endif ! CS%Leith_Kh
823 
824  meke_res_fn = 1.
825 
826  do j=jsq,jeq+1 ; do i=isq,ieq+1
827  if ((cs%Smagorinsky_Kh) .or. (cs%Smagorinsky_Ah)) then
828  shear_mag = sqrt(sh_xx(i,j)*sh_xx(i,j) + &
829  0.25*((sh_xy(i-1,j-1)*sh_xy(i-1,j-1) + sh_xy(i,j)*sh_xy(i,j)) + &
830  (sh_xy(i-1,j)*sh_xy(i-1,j) + sh_xy(i,j-1)*sh_xy(i,j-1))))
831  endif
832  if ((cs%Leith_Kh) .or. (cs%Leith_Ah)) then
833  if (cs%use_QG_Leith_visc) then
834  vert_vort_mag = min(grad_vort_mag_h(i,j) + grad_div_mag_h(i,j),3.*grad_vort_mag_h_2d(i,j))
835  else
836  vert_vort_mag = (grad_vort_mag_h(i,j) + grad_div_mag_h(i,j))
837  endif
838  endif
839  if (cs%better_bound_Ah .or. cs%better_bound_Kh) then
840  hrat_min = min(1.0, min(h_u(i,j), h_u(i-1,j), h_v(i,j), h_v(i,j-1)) / &
841  (h(i,j,k) + h_neglect) )
842  visc_bound_rem = 1.0
843  endif
844 
845  if (cs%Laplacian) then
846  ! Determine the Laplacian viscosity at h points, using the
847  ! largest value from several parameterizations.
848  kh = cs%Kh_bg_xx(i,j) ! Static (pre-computed) background viscosity
849  if (cs%add_LES_viscosity) then
850  if (cs%Smagorinsky_Kh) kh = kh + cs%Laplac2_const_xx(i,j) * shear_mag
851  if (cs%Leith_Kh) kh = kh + cs%Laplac3_const_xx(i,j) * vert_vort_mag*inv_pi3
852  else
853  if (cs%Smagorinsky_Kh) kh = max( kh, cs%Laplac2_const_xx(i,j) * shear_mag )
854  if (cs%Leith_Kh) kh = max( kh, cs%Laplac3_const_xx(i,j) * vert_vort_mag*inv_pi3)
855  endif
856  ! All viscosity contributions above are subject to resolution scaling
857  if (rescale_kh) kh = varmix%Res_fn_h(i,j) * kh
858  if (cs%res_scale_MEKE) meke_res_fn = varmix%Res_fn_h(i,j)
859  ! Older method of bounding for stability
860  if (legacy_bound) kh = min(kh, cs%Kh_Max_xx(i,j))
861  kh = max( kh, cs%Kh_bg_min ) ! Place a floor on the viscosity, if desired.
862  if (use_meke_ku) &
863  kh = kh + meke%Ku(i,j) * meke_res_fn ! *Add* the MEKE contribution (might be negative)
864  if (cs%anisotropic) kh = kh + cs%Kh_aniso * ( 1. - cs%n1n2_h(i,j)**2 ) ! *Add* the tension component
865  ! of anisotropic viscosity
866 
867  ! Newer method of bounding for stability
868  if (cs%better_bound_Kh) then
869  if (kh >= hrat_min*cs%Kh_Max_xx(i,j)) then
870  visc_bound_rem = 0.0
871  kh = hrat_min*cs%Kh_Max_xx(i,j)
872  else
873  visc_bound_rem = 1.0 - kh / (hrat_min*cs%Kh_Max_xx(i,j))
874  endif
875  endif
876 
877  if ((cs%id_Kh_h>0) .or. find_frictwork .or. cs%debug) kh_h(i,j,k) = kh
878 
879  if (cs%id_grid_Re_Kh>0) then
880  ke = 0.125*((u(i,j,k)+u(i-1,j,k))**2 + (v(i,j,k)+v(i,j-1,k))**2)
881  grid_re_kh(i,j,k) = (sqrt(ke) * sqrt(cs%grid_sp_h2(i,j))) &
882  / max(kh, cs%min_grid_Kh)
883  endif
884 
885  if (cs%id_div_xx_h>0) div_xx_h(i,j,k) = div_xx(i,j)
886  if (cs%id_sh_xx_h>0) sh_xx_h(i,j,k) = sh_xx(i,j)
887 
888  str_xx(i,j) = -kh * sh_xx(i,j)
889  else ! not Laplacian
890  str_xx(i,j) = 0.0
891  endif ! Laplacian
892 
893  if (cs%anisotropic) then
894  ! Shearing-strain averaged to h-points
895  local_strain = 0.25 * ( (sh_xy(i,j) + sh_xy(i-1,j-1)) + (sh_xy(i-1,j) + sh_xy(i,j-1)) )
896  ! *Add* the shear-strain contribution to the xx-component of stress
897  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
898  endif
899 
900  if (cs%biharmonic) then
901  ! Determine the biharmonic viscosity at h points, using the
902  ! largest value from several parameterizations.
903  ahsm = 0.0; ahlth = 0.0
904  if ((cs%Smagorinsky_Ah) .or. (cs%Leith_Ah)) then
905  if (cs%Smagorinsky_Ah) then
906  if (cs%bound_Coriolis) then
907  ahsm = shear_mag * (cs%Biharm_const_xx(i,j) + &
908  cs%Biharm_const2_xx(i,j)*shear_mag)
909  else
910  ahsm = cs%Biharm_const_xx(i,j) * shear_mag
911  endif
912  endif
913  if (cs%Leith_Ah) ahlth = cs%Biharm6_const_xx(i,j) * abs(del2vort_h(i,j)) * inv_pi6
914  ah = max(max(cs%Ah_bg_xx(i,j), ahsm), ahlth)
915  if (cs%bound_Ah .and. .not.cs%better_bound_Ah) &
916  ah = min(ah, cs%Ah_Max_xx(i,j))
917  else
918  ah = cs%Ah_bg_xx(i,j)
919  endif ! Smagorinsky_Ah or Leith_Ah
920 
921  if (use_meke_au) ah = ah + meke%Au(i,j) ! *Add* the MEKE contribution
922 
923  if (cs%Re_Ah > 0.0) then
924  ke = 0.125*((u(i,j,k)+u(i-1,j,k))**2 + (v(i,j,k)+v(i,j-1,k))**2)
925  ah = sqrt(ke) * cs%Re_Ah_const_xx(i,j)
926  endif
927 
928  if (cs%better_bound_Ah) then
929  ah = min(ah, visc_bound_rem*hrat_min*cs%Ah_Max_xx(i,j))
930  endif
931 
932  if ((cs%id_Ah_h>0) .or. find_frictwork .or. cs%debug) ah_h(i,j,k) = ah
933 
934  if (cs%id_grid_Re_Ah>0) then
935  ke = 0.125*((u(i,j,k)+u(i-1,j,k))**2 + (v(i,j,k)+v(i,j-1,k))**2)
936  grid_re_ah(i,j,k) = (sqrt(ke) * cs%grid_sp_h3(i,j)) &
937  / max(ah, cs%min_grid_Ah)
938  endif
939 
940  str_xx(i,j) = str_xx(i,j) + ah * &
941  (cs%DY_dxT(i,j) * (g%IdyCu(i,j)*del2u(i,j) - g%IdyCu(i-1,j)*del2u(i-1,j)) - &
942  cs%DX_dyT(i,j) * (g%IdxCv(i,j)*del2v(i,j) - g%IdxCv(i,j-1)*del2v(i,j-1)))
943 
944  ! Keep a copy of the biharmonic contribution for backscatter parameterization
945  bhstr_xx(i,j) = ah * &
946  (cs%DY_dxT(i,j) * (g%IdyCu(i,j)*del2u(i,j) - g%IdyCu(i-1,j)*del2u(i-1,j)) - &
947  cs%DX_dyT(i,j) * (g%IdxCv(i,j)*del2v(i,j) - g%IdxCv(i,j-1)*del2v(i,j-1)))
948  bhstr_xx(i,j) = bhstr_xx(i,j) * (h(i,j,k) * cs%reduction_xx(i,j))
949 
950  endif ! biharmonic
951 
952  enddo ; enddo
953 
954  if (cs%biharmonic) then
955  ! Gradient of Laplacian, for use in bi-harmonic term
956  do j=js-1,jeq ; do i=is-1,ieq
957  ddel2vdx(i,j) = cs%DY_dxBu(i,j)*(del2v(i+1,j)*g%IdyCv(i+1,j) - del2v(i,j)*g%IdyCv(i,j))
958  ddel2udy(i,j) = cs%DX_dyBu(i,j)*(del2u(i,j+1)*g%IdxCu(i,j+1) - del2u(i,j)*g%IdxCu(i,j))
959  enddo ; enddo
960  ! Adjust contributions to shearing strain on open boundaries.
961  if (apply_obc) then ; if (obc%zero_strain .or. obc%freeslip_strain) then
962  do n=1,obc%number_of_segments
963  j = obc%segment(n)%HI%JsdB ; i = obc%segment(n)%HI%IsdB
964  if (obc%segment(n)%is_N_or_S .and. (j >= js-1) .and. (j <= jeq)) then
965  do i=obc%segment(n)%HI%IsdB,obc%segment(n)%HI%IedB
966  if (obc%zero_strain) then
967  ddel2vdx(i,j) = 0. ; ddel2udy(i,j) = 0.
968  elseif (obc%freeslip_strain) then
969  ddel2udy(i,j) = 0.
970  endif
971  enddo
972  elseif (obc%segment(n)%is_E_or_W .and. (i >= is-1) .and. (i <= ieq)) then
973  do j=obc%segment(n)%HI%JsdB,obc%segment(n)%HI%JedB
974  if (obc%zero_strain) then
975  ddel2vdx(i,j) = 0. ; ddel2udy(i,j) = 0.
976  elseif (obc%freeslip_strain) then
977  ddel2vdx(i,j) = 0.
978  endif
979  enddo
980  endif
981  enddo
982  endif ; endif
983  endif
984 
985  meke_res_fn = 1.
986 
987  do j=js-1,jeq ; do i=is-1,ieq
988  if ((cs%Smagorinsky_Kh) .or. (cs%Smagorinsky_Ah)) then
989  shear_mag = sqrt(sh_xy(i,j)*sh_xy(i,j) + &
990  0.25*((sh_xx(i,j)*sh_xx(i,j) + sh_xx(i+1,j+1)*sh_xx(i+1,j+1)) + &
991  (sh_xx(i,j+1)*sh_xx(i,j+1) + sh_xx(i+1,j)*sh_xx(i+1,j))))
992  endif
993  if ((cs%Leith_Kh) .or. (cs%Leith_Ah)) then
994  if (cs%use_QG_Leith_visc) then
995  vert_vort_mag = min(grad_vort_mag_q(i,j) + grad_div_mag_q(i,j), 3.*grad_vort_mag_q_2d(i,j))
996  else
997  vert_vort_mag = (grad_vort_mag_q(i,j) + grad_div_mag_q(i,j))
998  endif
999  endif
1000  h2uq = 4.0 * h_u(i,j) * h_u(i,j+1)
1001  h2vq = 4.0 * h_v(i,j) * h_v(i+1,j)
1002  hq(i,j) = 2.0 * h2uq * h2vq / (h_neglect3 + (h2uq + h2vq) * &
1003  ((h_u(i,j) + h_u(i,j+1)) + (h_v(i,j) + h_v(i+1,j))))
1004 
1005  if (cs%better_bound_Ah .or. cs%better_bound_Kh) then
1006  hrat_min = min(1.0, min(h_u(i,j), h_u(i,j+1), h_v(i,j), h_v(i+1,j)) / &
1007  (hq(i,j) + h_neglect) )
1008  visc_bound_rem = 1.0
1009  endif
1010 
1011  if (cs%no_slip .and. (g%mask2dBu(i,j) < 0.5)) then
1012  if ((g%mask2dCu(i,j) + g%mask2dCu(i,j+1)) + &
1013  (g%mask2dCv(i,j) + g%mask2dCv(i+1,j)) > 0.0) then
1014  ! This is a coastal vorticity point, so modify hq and hrat_min.
1015 
1016  hu = g%mask2dCu(i,j) * h_u(i,j) + g%mask2dCu(i,j+1) * h_u(i,j+1)
1017  hv = g%mask2dCv(i,j) * h_v(i,j) + g%mask2dCv(i+1,j) * h_v(i+1,j)
1018  if ((g%mask2dCu(i,j) + g%mask2dCu(i,j+1)) * &
1019  (g%mask2dCv(i,j) + g%mask2dCv(i+1,j)) == 0.0) then
1020  ! Only one of hu and hv is nonzero, so just add them.
1021  hq(i,j) = hu + hv
1022  hrat_min = 1.0
1023  else
1024  ! Both hu and hv are nonzero, so take the harmonic mean.
1025  hq(i,j) = 2.0 * (hu * hv) / ((hu + hv) + h_neglect)
1026  hrat_min = min(1.0, min(hu, hv) / (hq(i,j) + h_neglect) )
1027  endif
1028  endif
1029  endif
1030 
1031  if (cs%Laplacian) then
1032  ! Determine the Laplacian viscosity at q points, using the
1033  ! largest value from several parameterizations.
1034  kh = cs%Kh_bg_xy(i,j) ! Static (pre-computed) background viscosity
1035  if (cs%add_LES_viscosity) then
1036  if (cs%Smagorinsky_Kh) kh = kh + cs%Laplac2_const_xx(i,j) * shear_mag
1037  if (cs%Leith_Kh) kh = kh + cs%Laplac3_const_xx(i,j) * vert_vort_mag*inv_pi3
1038  else
1039  if (cs%Smagorinsky_Kh) kh = max( kh, cs%Laplac2_const_xy(i,j) * shear_mag )
1040  if (cs%Leith_Kh) kh = max( kh, cs%Laplac3_const_xy(i,j) * vert_vort_mag*inv_pi3)
1041  endif
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 - kh / (hrat_min*cs%Kh_Max_xy(i,j))
1063  endif
1064  endif
1065 
1066  if (cs%id_Kh_q>0 .or. cs%debug) kh_q(i,j,k) = kh
1067  if (cs%id_vort_xy_q>0) vort_xy_q(i,j,k) = vort_xy(i,j)
1068  if (cs%id_sh_xy_q>0) sh_xy_q(i,j,k) = sh_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%Biharm6_const_xy(i,j) * abs(del2vort_q(i,j)) * inv_pi6
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%Re_Ah > 0.0) then
1109  ke = 0.125*((u(i,j,k)+u(i,j+1,k))**2 + (v(i,j,k)+v(i+1,j,k))**2)
1110  ah = sqrt(ke) * cs%Re_Ah_const_xy(i,j)
1111  endif
1112 
1113  if (cs%better_bound_Ah) then
1114  ah = min(ah, visc_bound_rem*hrat_min*cs%Ah_Max_xy(i,j))
1115  endif
1116 
1117  if (cs%id_Ah_q>0 .or. cs%debug) ah_q(i,j,k) = ah
1118 
1119  str_xy(i,j) = str_xy(i,j) + ah * ( ddel2vdx(i,j) + ddel2udy(i,j) )
1120 
1121  ! Keep a copy of the biharmonic contribution for backscatter parameterization
1122  bhstr_xy(i,j) = ah * ( ddel2vdx(i,j) + ddel2udy(i,j) ) * &
1123  (hq(i,j) * g%mask2dBu(i,j) * cs%reduction_xy(i,j))
1124 
1125  endif ! biharmonic
1126 
1127  enddo ; enddo
1128 
1129  if (cs%use_GME) then
1130  call thickness_diffuse_get_kh(td, kh_u_gme, kh_v_gme, g)
1131  call pass_vector(kh_u_gme, kh_v_gme, g%Domain)
1132 
1133  do j=jsq,jeq+1 ; do i=isq,ieq+1
1134  if (grad_vel_mag_bt_h(i,j)>0) then
1135  gme_coeff = cs%GME_efficiency * (min(g%bathyT(i,j)/cs%GME_h0,1.0)**2) * &
1136  (0.25*(kh_u_gme(i,j,k)+kh_u_gme(i-1,j,k)+kh_v_gme(i,j,k)+kh_v_gme(i,j-1,k)))
1137  else
1138  gme_coeff = 0.0
1139  endif
1140 
1141  ! apply mask
1142  gme_coeff = gme_coeff * boundary_mask_h(i,j)
1143  gme_coeff = min(gme_coeff, cs%GME_limiter)
1144 
1145  if ((cs%id_GME_coeff_h>0) .or. find_frictwork) gme_coeff_h(i,j,k) = gme_coeff
1146  str_xx_gme(i,j) = gme_coeff * sh_xx_bt(i,j)
1147 
1148  enddo ; enddo
1149 
1150  do j=js-1,jeq ; do i=is-1,ieq
1151  if (grad_vel_mag_bt_q(i,j)>0) then
1152  gme_coeff = cs%GME_efficiency * (min(g%bathyT(i,j)/cs%GME_h0,1.0)**2) * &
1153  (0.25*(kh_u_gme(i,j,k)+kh_u_gme(i,j+1,k)+kh_v_gme(i,j,k)+kh_v_gme(i+1,j,k)))
1154  else
1155  gme_coeff = 0.0
1156  endif
1157 
1158  ! apply mask
1159  gme_coeff = gme_coeff * boundary_mask_q(i,j)
1160  gme_coeff = min(gme_coeff, cs%GME_limiter)
1161 
1162  if (cs%id_GME_coeff_q>0) gme_coeff_q(i,j,k) = gme_coeff
1163  str_xy_gme(i,j) = gme_coeff * sh_xy_bt(i,j)
1164 
1165  enddo ; enddo
1166 
1167  ! Applying GME diagonal term. This is linear and the arguments can be rescaled.
1168  call smooth_gme(cs,g,gme_flux_h=str_xx_gme)
1169  call smooth_gme(cs,g,gme_flux_q=str_xy_gme)
1170 
1171  do j=jsq,jeq+1 ; do i=isq,ieq+1
1172  str_xx(i,j) = (str_xx(i,j) + str_xx_gme(i,j)) * (h(i,j,k) * cs%reduction_xx(i,j))
1173  enddo ; enddo
1174 
1175  do j=js-1,jeq ; do i=is-1,ieq
1176  ! GME is applied below
1177  if (cs%no_slip) then
1178  str_xy(i,j) = (str_xy(i,j) + str_xy_gme(i,j)) * (hq(i,j) * cs%reduction_xy(i,j))
1179  else
1180  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))
1181  endif
1182  enddo ; enddo
1183 
1184  if (associated(meke%GME_snk)) then
1185  do j=js,je ; do i=is,ie
1186  frictwork_gme(i,j,k) = gme_coeff_h(i,j,k) * h(i,j,k) * gv%H_to_kg_m2 * grad_vel_mag_bt_h(i,j)
1187  enddo ; enddo
1188  endif
1189 
1190  else ! use_GME
1191  do j=jsq,jeq+1 ; do i=isq,ieq+1
1192  str_xx(i,j) = str_xx(i,j) * (h(i,j,k) * cs%reduction_xx(i,j))
1193  enddo ; enddo
1194 
1195  do j=js-1,jeq ; do i=is-1,ieq
1196  if (cs%no_slip) then
1197  str_xy(i,j) = str_xy(i,j) * (hq(i,j) * cs%reduction_xy(i,j))
1198  else
1199  str_xy(i,j) = str_xy(i,j) * (hq(i,j) * g%mask2dBu(i,j) * cs%reduction_xy(i,j))
1200  endif
1201  enddo ; enddo
1202 
1203  endif ! use_GME
1204 
1205  ! Evaluate 1/h x.Div(h Grad u) or the biharmonic equivalent.
1206  do j=js,je ; do i=isq,ieq
1207  diffu(i,j,k) = ((g%IdyCu(i,j)*(cs%dy2h(i,j) *str_xx(i,j) - &
1208  cs%dy2h(i+1,j)*str_xx(i+1,j)) + &
1209  g%IdxCu(i,j)*(cs%dx2q(i,j-1)*str_xy(i,j-1) - &
1210  cs%dx2q(i,j) *str_xy(i,j))) * &
1211  g%IareaCu(i,j)) / (h_u(i,j) + h_neglect)
1212 
1213  enddo ; enddo
1214  if (apply_obc) then
1215  ! This is not the right boundary condition. If all the masking of tendencies are done
1216  ! correctly later then eliminating this block should not change answers.
1217  do n=1,obc%number_of_segments
1218  if (obc%segment(n)%is_E_or_W) then
1219  i = obc%segment(n)%HI%IsdB
1220  do j=obc%segment(n)%HI%jsd,obc%segment(n)%HI%jed
1221  diffu(i,j,k) = 0.
1222  enddo
1223  endif
1224  enddo
1225  endif
1226 
1227  ! Evaluate 1/h y.Div(h Grad u) or the biharmonic equivalent.
1228  do j=jsq,jeq ; do i=is,ie
1229  diffv(i,j,k) = ((g%IdyCv(i,j)*(cs%dy2q(i-1,j)*str_xy(i-1,j) - &
1230  cs%dy2q(i,j) *str_xy(i,j)) - &
1231  g%IdxCv(i,j)*(cs%dx2h(i,j) *str_xx(i,j) - &
1232  cs%dx2h(i,j+1)*str_xx(i,j+1))) * &
1233  g%IareaCv(i,j)) / (h_v(i,j) + h_neglect)
1234  enddo ; enddo
1235  if (apply_obc) then
1236  ! This is not the right boundary condition. If all the masking of tendencies are done
1237  ! correctly later then eliminating this block should not change answers.
1238  do n=1,obc%number_of_segments
1239  if (obc%segment(n)%is_N_or_S) then
1240  j = obc%segment(n)%HI%JsdB
1241  do i=obc%segment(n)%HI%isd,obc%segment(n)%HI%ied
1242  diffv(i,j,k) = 0.
1243  enddo
1244  endif
1245  enddo
1246  endif
1247 
1248  if (find_frictwork) then ; do j=js,je ; do i=is,ie
1249  ! Diagnose str_xx*d_x u - str_yy*d_y v + str_xy*(d_y u + d_x v)
1250  ! This is the old formulation that includes energy diffusion
1251  frictwork(i,j,k) = gv%H_to_RZ * ( &
1252  (str_xx(i,j)*(u(i,j,k)-u(i-1,j,k))*g%IdxT(i,j) &
1253  -str_xx(i,j)*(v(i,j,k)-v(i,j-1,k))*g%IdyT(i,j)) &
1254  +0.25*((str_xy(i,j)*( &
1255  (u(i,j+1,k)-u(i,j,k))*g%IdyBu(i,j) &
1256  +(v(i+1,j,k)-v(i,j,k))*g%IdxBu(i,j) ) &
1257  +str_xy(i-1,j-1)*( &
1258  (u(i-1,j,k)-u(i-1,j-1,k))*g%IdyBu(i-1,j-1) &
1259  +(v(i,j-1,k)-v(i-1,j-1,k))*g%IdxBu(i-1,j-1) )) &
1260  +(str_xy(i-1,j)*( &
1261  (u(i-1,j+1,k)-u(i-1,j,k))*g%IdyBu(i-1,j) &
1262  +(v(i,j,k)-v(i-1,j,k))*g%IdxBu(i-1,j) ) &
1263  +str_xy(i,j-1)*( &
1264  (u(i,j,k)-u(i,j-1,k))*g%IdyBu(i,j-1) &
1265  +(v(i+1,j-1,k)-v(i,j-1,k))*g%IdxBu(i,j-1) )) ) )
1266  enddo ; enddo ; endif
1267 
1268  ! Make a similar calculation as for FrictWork above but accumulating into
1269  ! the vertically integrated MEKE source term, and adjusting for any
1270  ! energy loss seen as a reduction in the (biharmonic) frictional source term.
1271  if (find_frictwork .and. associated(meke)) then ; if (associated(meke%mom_src)) then
1272  if (k==1) then
1273  do j=js,je ; do i=is,ie
1274  meke%mom_src(i,j) = 0.
1275  enddo ; enddo
1276  if (associated(meke%GME_snk)) then
1277  do j=js,je ; do i=is,ie
1278  meke%GME_snk(i,j) = 0.
1279  enddo ; enddo
1280  endif
1281  endif
1282  if (meke%backscatter_Ro_c /= 0.) then
1283  do j=js,je ; do i=is,ie
1284  fath = 0.25*( (abs(g%CoriolisBu(i-1,j-1)) + abs(g%CoriolisBu(i,j))) + &
1285  (abs(g%CoriolisBu(i-1,j)) + abs(g%CoriolisBu(i,j-1))) )
1286  shear_mag = sqrt(sh_xx(i,j)*sh_xx(i,j) + &
1287  0.25*((sh_xy(i-1,j-1)*sh_xy(i-1,j-1) + sh_xy(i,j)*sh_xy(i,j)) + &
1288  (sh_xy(i-1,j)*sh_xy(i-1,j) + sh_xy(i,j-1)*sh_xy(i,j-1))))
1289  if (cs%answers_2018) then
1290  fath = (us%s_to_T*fath)**meke%backscatter_Ro_pow ! f^n
1291  ! Note the hard-coded dimensional constant in the following line that can not
1292  ! be rescaled for dimensional consistency.
1293  shear_mag = ( ( (us%s_to_T*shear_mag)**meke%backscatter_Ro_pow ) + 1.e-30 ) &
1294  * meke%backscatter_Ro_c ! c * D^n
1295  ! The Rossby number function is g(Ro) = 1/(1+c.Ro^n)
1296  ! RoScl = 1 - g(Ro)
1297  roscl = shear_mag / ( fath + shear_mag ) ! = 1 - f^n/(f^n+c*D^n)
1298  else
1299  if (fath <= backscat_subround*shear_mag) then
1300  roscl = 1.0
1301  else
1302  sh_f_pow = meke%backscatter_Ro_c * (shear_mag / fath)**meke%backscatter_Ro_pow
1303  roscl = sh_f_pow / (1.0 + sh_f_pow) ! = 1 - f^n/(f^n+c*D^n)
1304  endif
1305  endif
1306 
1307 
1308  meke%mom_src(i,j) = meke%mom_src(i,j) + gv%H_to_RZ * ( &
1309  ((str_xx(i,j)-roscl*bhstr_xx(i,j))*(u(i,j,k)-u(i-1,j,k))*g%IdxT(i,j) &
1310  -(str_xx(i,j)-roscl*bhstr_xx(i,j))*(v(i,j,k)-v(i,j-1,k))*g%IdyT(i,j)) &
1311  +0.25*(((str_xy(i,j)-roscl*bhstr_xy(i,j))*( &
1312  (u(i,j+1,k)-u(i,j,k))*g%IdyBu(i,j) &
1313  +(v(i+1,j,k)-v(i,j,k))*g%IdxBu(i,j) ) &
1314  +(str_xy(i-1,j-1)-roscl*bhstr_xy(i-1,j-1))*( &
1315  (u(i-1,j,k)-u(i-1,j-1,k))*g%IdyBu(i-1,j-1) &
1316  +(v(i,j-1,k)-v(i-1,j-1,k))*g%IdxBu(i-1,j-1) )) &
1317  +((str_xy(i-1,j)-roscl*bhstr_xy(i-1,j))*( &
1318  (u(i-1,j+1,k)-u(i-1,j,k))*g%IdyBu(i-1,j) &
1319  +(v(i,j,k)-v(i-1,j,k))*g%IdxBu(i-1,j) ) &
1320  +(str_xy(i,j-1)-roscl*bhstr_xy(i,j-1))*( &
1321  (u(i,j,k)-u(i,j-1,k))*g%IdyBu(i,j-1) &
1322  +(v(i+1,j-1,k)-v(i,j-1,k))*g%IdxBu(i,j-1) )) ) )
1323  enddo ; enddo
1324  endif ! MEKE%backscatter_Ro_c
1325 
1326  do j=js,je ; do i=is,ie
1327  meke%mom_src(i,j) = meke%mom_src(i,j) + frictwork(i,j,k)
1328  enddo ; enddo
1329 
1330  if (cs%use_GME .and. associated(meke)) then ; if (associated(meke%GME_snk)) then
1331  do j=js,je ; do i=is,ie
1332  meke%GME_snk(i,j) = meke%GME_snk(i,j) + frictwork_gme(i,j,k)
1333  enddo ; enddo
1334  endif ; endif
1335 
1336  endif ; endif ! find_FrictWork and associated(mom_src)
1337 
1338  enddo ! end of k loop
1339 
1340  ! Offer fields for diagnostic averaging.
1341  if (cs%id_diffu>0) call post_data(cs%id_diffu, diffu, cs%diag)
1342  if (cs%id_diffv>0) call post_data(cs%id_diffv, diffv, cs%diag)
1343  if (cs%id_FrictWork>0) call post_data(cs%id_FrictWork, frictwork, cs%diag)
1344  if (cs%id_FrictWork_GME>0) call post_data(cs%id_FrictWork_GME, frictwork_gme, cs%diag)
1345  if (cs%id_Ah_h>0) call post_data(cs%id_Ah_h, ah_h, cs%diag)
1346  if (cs%id_grid_Re_Ah>0) call post_data(cs%id_grid_Re_Ah, grid_re_ah, cs%diag)
1347  if (cs%id_div_xx_h>0) call post_data(cs%id_div_xx_h, div_xx_h, cs%diag)
1348  if (cs%id_vort_xy_q>0) call post_data(cs%id_vort_xy_q, vort_xy_q, cs%diag)
1349  if (cs%id_sh_xx_h>0) call post_data(cs%id_sh_xx_h, sh_xx_h, cs%diag)
1350  if (cs%id_sh_xy_q>0) call post_data(cs%id_sh_xy_q, sh_xy_q, cs%diag)
1351  if (cs%id_Ah_q>0) call post_data(cs%id_Ah_q, ah_q, cs%diag)
1352  if (cs%id_Kh_h>0) call post_data(cs%id_Kh_h, kh_h, cs%diag)
1353  if (cs%id_grid_Re_Kh>0) call post_data(cs%id_grid_Re_Kh, grid_re_kh, cs%diag)
1354  if (cs%id_Kh_q>0) call post_data(cs%id_Kh_q, kh_q, cs%diag)
1355  if (cs%id_GME_coeff_h > 0) call post_data(cs%id_GME_coeff_h, gme_coeff_h, cs%diag)
1356  if (cs%id_GME_coeff_q > 0) call post_data(cs%id_GME_coeff_q, gme_coeff_q, cs%diag)
1357 
1358  if (cs%debug) then
1359  if (cs%Laplacian) then
1360  call hchksum(kh_h, "Kh_h", g%HI, haloshift=0, scale=us%L_to_m**2*us%s_to_T)
1361  call bchksum(kh_q, "Kh_q", g%HI, haloshift=0, scale=us%L_to_m**2*us%s_to_T)
1362  endif
1363  if (cs%biharmonic) call hchksum(ah_h, "Ah_h", g%HI, haloshift=0, scale=us%L_to_m**4*us%s_to_T)
1364  if (cs%biharmonic) call bchksum(ah_q, "Ah_q", g%HI, haloshift=0, scale=us%L_to_m**4*us%s_to_T)
1365  endif
1366 
1367  if (cs%id_FrictWorkIntz > 0) then
1368  do j=js,je
1369  do i=is,ie ; frictworkintz(i,j) = frictwork(i,j,1) ; enddo
1370  do k=2,nz ; do i=is,ie
1371  frictworkintz(i,j) = frictworkintz(i,j) + frictwork(i,j,k)
1372  enddo ; enddo
1373  enddo
1374  call post_data(cs%id_FrictWorkIntz, frictworkintz, cs%diag)
1375  endif
1376 
1377  ! Diagnostics for terms multiplied by fractional thicknesses
1378 
1379  ! 3D diagnostics hf_diffu(diffv) are commented because there is no clarity on proper remapping grid option.
1380  ! The code is retained for degugging purposes in the future.
1381  !if (present(ADp) .and. (CS%id_hf_diffu > 0)) then
1382  ! do k=1,nz ; do j=js,je ; do I=Isq,Ieq
1383  ! CS%hf_diffu(I,j,k) = diffu(I,j,k) * ADp%diag_hfrac_u(I,j,k)
1384  ! enddo ; enddo ; enddo
1385  ! call post_data(CS%id_hf_diffu, CS%hf_diffu, CS%diag)
1386  !endif
1387  !if (present(ADp) .and. (CS%id_hf_diffv > 0)) then
1388  ! do k=1,nz ; do J=Jsq,Jeq ; do i=is,ie
1389  ! CS%hf_diffv(i,J,k) = diffv(i,J,k) * ADp%diag_hfrac_v(i,J,k)
1390  ! enddo ; enddo ; enddo
1391  ! call post_data(CS%id_hf_diffv, CS%hf_diffv, CS%diag)
1392  !endif
1393  if (present(adp) .and. (cs%id_hf_diffu_2d > 0)) then
1394  allocate(hf_diffu_2d(g%IsdB:g%IedB,g%jsd:g%jed))
1395  hf_diffu_2d(:,:) = 0.0
1396  do k=1,nz ; do j=js,je ; do i=isq,ieq
1397  hf_diffu_2d(i,j) = hf_diffu_2d(i,j) + diffu(i,j,k) * adp%diag_hfrac_u(i,j,k)
1398  enddo ; enddo ; enddo
1399  call post_data(cs%id_hf_diffu_2d, hf_diffu_2d, cs%diag)
1400  deallocate(hf_diffu_2d)
1401  endif
1402  if (present(adp) .and. (cs%id_hf_diffv_2d > 0)) then
1403  allocate(hf_diffv_2d(g%isd:g%ied,g%JsdB:g%JedB))
1404  hf_diffv_2d(:,:) = 0.0
1405  do k=1,nz ; do j=jsq,jeq ; do i=is,ie
1406  hf_diffv_2d(i,j) = hf_diffv_2d(i,j) + diffv(i,j,k) * adp%diag_hfrac_v(i,j,k)
1407  enddo ; enddo ; enddo
1408  call post_data(cs%id_hf_diffv_2d, hf_diffv_2d, cs%diag)
1409  deallocate(hf_diffv_2d)
1410  endif
1411 

◆ smooth_gme()

subroutine mom_hor_visc::smooth_gme ( type(hor_visc_cs), pointer  CS,
type(ocean_grid_type), intent(in)  G,
real, dimension(szi_(g),szj_(g)), intent(inout), optional  GME_flux_h,
real, dimension(szib_(g),szjb_(g)), intent(inout), optional  GME_flux_q 
)
private

Apply a 1-1-4-1-1 Laplacian filter one time on GME diffusive flux to reduce any horizontal two-grid-point noise.

Parameters
csControl structure
[in]gOcean grid
[in,out]gme_flux_hGME diffusive flux at h points
[in,out]gme_flux_qGME diffusive flux at q points

Definition at line 2191 of file MOM_hor_visc.F90.

2192  ! Arguments
2193  type(hor_visc_CS), pointer :: CS !< Control structure
2194  type(ocean_grid_type), intent(in) :: G !< Ocean grid
2195  real, dimension(SZI_(G),SZJ_(G)), optional, intent(inout) :: GME_flux_h!< GME diffusive flux
2196  !! at h points
2197  real, dimension(SZIB_(G),SZJB_(G)), optional, intent(inout) :: GME_flux_q!< GME diffusive flux
2198  !! at q points
2199  ! local variables
2200  real, dimension(SZI_(G),SZJ_(G)) :: GME_flux_h_original
2201  real, dimension(SZIB_(G),SZJB_(G)) :: GME_flux_q_original
2202  real :: wc, ww, we, wn, ws ! averaging weights for smoothing
2203  integer :: i, j, k, s
2204  do s=1,1
2205  ! Update halos
2206  if (present(gme_flux_h)) then
2207  call pass_var(gme_flux_h, g%Domain)
2208  gme_flux_h_original = gme_flux_h
2209  ! apply smoothing on GME
2210  do j = g%jsc, g%jec
2211  do i = g%isc, g%iec
2212  ! skip land points
2213  if (g%mask2dT(i,j)==0.) cycle
2214  ! compute weights
2215  ww = 0.125 * g%mask2dT(i-1,j)
2216  we = 0.125 * g%mask2dT(i+1,j)
2217  ws = 0.125 * g%mask2dT(i,j-1)
2218  wn = 0.125 * g%mask2dT(i,j+1)
2219  wc = 1.0 - (ww+we+wn+ws)
2220  gme_flux_h(i,j) = wc * gme_flux_h_original(i,j) &
2221  + ww * gme_flux_h_original(i-1,j) &
2222  + we * gme_flux_h_original(i+1,j) &
2223  + ws * gme_flux_h_original(i,j-1) &
2224  + wn * gme_flux_h_original(i,j+1)
2225  enddo; enddo
2226  endif
2227  ! Update halos
2228  if (present(gme_flux_q)) then
2229  call pass_var(gme_flux_q, g%Domain, position=corner, complete=.true.)
2230  gme_flux_q_original = gme_flux_q
2231  ! apply smoothing on GME
2232  do j = g%JscB, g%JecB
2233  do i = g%IscB, g%IecB
2234  ! skip land points
2235  if (g%mask2dBu(i,j)==0.) cycle
2236  ! compute weights
2237  ww = 0.125 * g%mask2dBu(i-1,j)
2238  we = 0.125 * g%mask2dBu(i+1,j)
2239  ws = 0.125 * g%mask2dBu(i,j-1)
2240  wn = 0.125 * g%mask2dBu(i,j+1)
2241  wc = 1.0 - (ww+we+wn+ws)
2242  gme_flux_q(i,j) = wc * gme_flux_q_original(i,j) &
2243  + ww * gme_flux_q_original(i-1,j) &
2244  + we * gme_flux_q_original(i+1,j) &
2245  + ws * gme_flux_q_original(i,j-1) &
2246  + wn * gme_flux_q_original(i,j+1)
2247  enddo; enddo
2248  endif
2249  enddo ! s-loop