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)
 Calculates the acceleration due to the horizontal viscosity. More...
 
subroutine, public hor_visc_init (Time, G, US, param_file, diag, CS, MEKE)
 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 2123 of file MOM_hor_visc.F90.

2123  type(hor_visc_CS), pointer :: CS !< Control structure for horizontal viscosity
2124  real, intent(in) :: n1 !< i-component of direction vector [nondim]
2125  real, intent(in) :: n2 !< j-component of direction vector [nondim]
2126  ! Local variables
2127  real :: recip_n2_norm
2128 
2129  ! For normalizing n=(n1,n2) in case arguments are not a unit vector
2130  recip_n2_norm = n1**2 + n2**2
2131  if (recip_n2_norm > 0.) recip_n2_norm = 1./recip_n2_norm
2132 
2133  cs%n1n2_h(:,:) = 2. * ( n1 * n2 ) * recip_n2_norm
2134  cs%n1n2_q(:,:) = 2. * ( n1 * n2 ) * recip_n2_norm
2135  cs%n1n1_m_n2n2_h(:,:) = ( n1 * n1 - n2 * n2 ) * recip_n2_norm
2136  cs%n1n1_m_n2n2_q(:,:) = ( n1 * n1 - n2 * n2 ) * recip_n2_norm
2137 

Referenced by hor_visc_init().

Here is the caller graph for this function:

◆ 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 2215 of file MOM_hor_visc.F90.

2215  type(hor_visc_CS), pointer :: CS !< The control structure returned by a
2216  !! previous call to hor_visc_init.
2217 
2218  if (cs%Laplacian .or. cs%biharmonic) then
2219  dealloc_(cs%dx2h) ; dealloc_(cs%dx2q) ; dealloc_(cs%dy2h) ; dealloc_(cs%dy2q)
2220  dealloc_(cs%dx_dyT) ; dealloc_(cs%dy_dxT) ; dealloc_(cs%dx_dyBu) ; dealloc_(cs%dy_dxBu)
2221  dealloc_(cs%reduction_xx) ; dealloc_(cs%reduction_xy)
2222  endif
2223 
2224  if (cs%Laplacian) then
2225  dealloc_(cs%Kh_bg_xx) ; dealloc_(cs%Kh_bg_xy)
2226  if (cs%bound_Kh) then
2227  dealloc_(cs%Kh_Max_xx) ; dealloc_(cs%Kh_Max_xy)
2228  endif
2229  if (cs%Smagorinsky_Kh) then
2230  dealloc_(cs%Laplac2_const_xx) ; dealloc_(cs%Laplac2_const_xy)
2231  endif
2232  if (cs%Leith_Kh) then
2233  dealloc_(cs%Laplac3_const_xx) ; dealloc_(cs%Laplac3_const_xy)
2234  endif
2235  endif
2236 
2237  if (cs%biharmonic) then
2238  dealloc_(cs%Idx2dyCu) ; dealloc_(cs%Idx2dyCv)
2239  dealloc_(cs%Idxdy2u) ; dealloc_(cs%Idxdy2v)
2240  dealloc_(cs%Ah_bg_xx) ; dealloc_(cs%Ah_bg_xy)
2241  if (cs%bound_Ah) then
2242  dealloc_(cs%Ah_Max_xx) ; dealloc_(cs%Ah_Max_xy)
2243  endif
2244  if (cs%Smagorinsky_Ah) then
2245  dealloc_(cs%Biharm5_const_xx) ; dealloc_(cs%Biharm5_const_xy)
2246  ! if (CS%bound_Coriolis) then
2247  ! DEALLOC_(CS%Biharm5_const2_xx) ; DEALLOC_(CS%Biharm5_const2_xy)
2248  ! endif
2249  endif
2250  if (cs%Leith_Ah) then
2251  dealloc_(cs%Biharm_const_xx) ; dealloc_(cs%Biharm_const_xy)
2252  endif
2253  endif
2254  if (cs%anisotropic) then
2255  dealloc_(cs%n1n2_h)
2256  dealloc_(cs%n1n2_q)
2257  dealloc_(cs%n1n1_m_n2n2_h)
2258  dealloc_(cs%n1n1_m_n2n2_q)
2259  endif
2260  deallocate(cs)
2261 

◆ 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 
)

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

Definition at line 1349 of file MOM_hor_visc.F90.

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

References align_aniso_tensor_to_grid(), and mom_error_handler::mom_error().

Referenced by mom_dynamics_split_rk2::initialize_dyn_split_rk2(), mom_dynamics_unsplit::initialize_dyn_unsplit(), and mom_dynamics_unsplit_rk2::initialize_dyn_unsplit_rk2().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ 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 
)

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.

Definition at line 207 of file MOM_hor_visc.F90.

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

References mom_barotropic::barotropic_get_tav(), mom_lateral_mixing_coeffs::calc_qg_leith_viscosity(), mom_error_handler::mom_error(), mom_open_boundary::obc_direction_n, mom_open_boundary::obc_direction_s, and smooth_gme().

Referenced by mom_dynamics_split_rk2::initialize_dyn_split_rk2(), mom_dynamics_split_rk2::step_mom_dyn_split_rk2(), mom_dynamics_unsplit::step_mom_dyn_unsplit(), and mom_dynamics_unsplit_rk2::step_mom_dyn_unsplit_rk2().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ 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 2143 of file MOM_hor_visc.F90.

2143  ! Arguments
2144  type(hor_visc_CS), pointer :: CS !< Control structure
2145  type(ocean_grid_type), intent(in) :: G !< Ocean grid
2146  real, dimension(SZI_(G),SZJ_(G)), optional, intent(inout) :: GME_flux_h!< GME diffusive flux
2147  !! at h points
2148  real, dimension(SZIB_(G),SZJB_(G)), optional, intent(inout) :: GME_flux_q!< GME diffusive flux
2149  !! at q points
2150 
2151  ! local variables
2152  real, dimension(SZI_(G),SZJ_(G)) :: GME_flux_h_original
2153  real, dimension(SZIB_(G),SZJB_(G)) :: GME_flux_q_original
2154  real :: wc, ww, we, wn, ws ! averaging weights for smoothing
2155  integer :: i, j, k, s
2156 
2157  do s=1,1
2158 
2159  ! Update halos
2160  if (present(gme_flux_h)) then
2161  call pass_var(gme_flux_h, g%Domain)
2162  gme_flux_h_original = gme_flux_h
2163  ! apply smoothing on GME
2164  do j = g%jsc, g%jec
2165  do i = g%isc, g%iec
2166  ! skip land points
2167  if (g%mask2dT(i,j)==0.) cycle
2168 
2169  ! compute weights
2170  ww = 0.125 * g%mask2dT(i-1,j)
2171  we = 0.125 * g%mask2dT(i+1,j)
2172  ws = 0.125 * g%mask2dT(i,j-1)
2173  wn = 0.125 * g%mask2dT(i,j+1)
2174  wc = 1.0 - (ww+we+wn+ws)
2175 
2176  gme_flux_h(i,j) = wc * gme_flux_h_original(i,j) &
2177  + ww * gme_flux_h_original(i-1,j) &
2178  + we * gme_flux_h_original(i+1,j) &
2179  + ws * gme_flux_h_original(i,j-1) &
2180  + wn * gme_flux_h_original(i,j+1)
2181  enddo; enddo
2182  endif
2183 
2184  ! Update halos
2185  if (present(gme_flux_q)) then
2186  call pass_var(gme_flux_q, g%Domain, position=corner, complete=.true.)
2187  gme_flux_q_original = gme_flux_q
2188  ! apply smoothing on GME
2189  do j = g%JscB, g%JecB
2190  do i = g%IscB, g%IecB
2191  ! skip land points
2192  if (g%mask2dBu(i,j)==0.) cycle
2193 
2194  ! compute weights
2195  ww = 0.125 * g%mask2dBu(i-1,j)
2196  we = 0.125 * g%mask2dBu(i+1,j)
2197  ws = 0.125 * g%mask2dBu(i,j-1)
2198  wn = 0.125 * g%mask2dBu(i,j+1)
2199  wc = 1.0 - (ww+we+wn+ws)
2200 
2201  gme_flux_q(i,j) = wc * gme_flux_q_original(i,j) &
2202  + ww * gme_flux_q_original(i-1,j) &
2203  + we * gme_flux_q_original(i+1,j) &
2204  + ws * gme_flux_q_original(i,j-1) &
2205  + wn * gme_flux_q_original(i,j+1)
2206  enddo; enddo
2207  endif
2208 
2209  enddo ! s-loop
2210 

Referenced by horizontal_viscosity().

Here is the caller graph for this function: