mom_hor_visc module reference

Calculates horizontal viscosity and viscous stresses.

More…

Data Types

hor_visc_cs

Control structure for horizontal viscosity.

Functions/Subroutines

horizontal_viscosity()

Calculates the acceleration due to the horizontal viscosity.

hor_visc_init()

Allocates space for and calculates static variables used by horizontal_viscosity().

align_aniso_tensor_to_grid()

Calculates factors in the anisotropic orientation tensor to be align with the grid.

smooth_gme()

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

smooth_x9_h()

Apply a 9-point smoothing filter twice to a field staggered at a thickness point to reduce horizontal two-grid-point noise.

smooth_x9_uv()

Apply a 9-point smoothing filter twice to a pair of velocity components to reduce horizontal two-grid-point noise.

hor_visc_end()

Deallocates any variables allocated in hor_visc_init.

Detailed Description

This module contains the subroutine horizontal_viscosity() that calculates the effects of horizontal viscosity, including parameterizations of the value of the viscosity itself. 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.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

\[\begin{split}{\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}\end{split}\]

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

\[\begin{split}\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}\end{split}\]

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

\[\begin{split}{\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} .\end{split}\]

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

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

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{split}\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*}\end{split}\]

The form of the Laplacian viscosity in general coordinates is:

\[\begin{split}\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*}\end{split}\]

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{split}\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}\end{split}\]

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{split}\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*}\end{split}\]

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{split}\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*}\end{split}\]

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

Type Documentation

type mom_hor_visc/hor_visc_cs

Control structure for horizontal viscosity.

Type fields:
  • % id_grid_re_ah [integer] :: Diagnostic id.

  • % id_grid_re_kh [integer] :: Diagnostic id.

  • % id_diffu [integer] :: Diagnostic id.

  • % id_diffv [integer] :: Diagnostic id.

  • % id_h_diffu [integer] :: Diagnostic id.

  • % id_h_diffv [integer] :: Diagnostic id.

  • % id_hf_diffu_2d [integer] :: Diagnostic id.

  • % id_hf_diffv_2d [integer] :: Diagnostic id.

  • % id_intz_diffu_2d [integer] :: Diagnostic id.

  • % id_intz_diffv_2d [integer] :: Diagnostic id.

  • % id_diffu_visc_rem [integer] :: Diagnostic id.

  • % id_diffv_visc_rem [integer] :: Diagnostic id.

  • % id_ah_h [integer] :: Diagnostic id.

  • % id_ah_q [integer] :: Diagnostic id.

  • % id_kh_h [integer] :: Diagnostic id.

  • % id_kh_q [integer] :: Diagnostic id.

  • % id_gme_coeff_h [integer] :: Diagnostic id.

  • % id_gme_coeff_q [integer] :: Diagnostic id.

  • % id_dudx_bt [integer] :: Diagnostic id.

  • % id_dvdy_bt [integer] :: Diagnostic id.

  • % id_dudy_bt [integer] :: Diagnostic id.

  • % id_dvdx_bt [integer] :: Diagnostic id.

  • % id_vort_xy_q [integer] :: Diagnostic id.

  • % id_div_xx_h [integer] :: Diagnostic id.

  • % id_sh_xy_q [integer] :: Diagnostic id.

  • % id_sh_xx_h [integer] :: Diagnostic id.

  • % id_frictwork [integer] :: Diagnostic id.

  • % id_frictworkintz [integer] :: Diagnostic id.

  • % id_frictwork_gme [integer] :: Diagnostic id.

  • % id_normstress [integer] :: Diagnostic id.

  • % id_shearstress [integer] :: Diagnostic id.

  • % initialized [logical] :: True if this control structure has been initialized.

  • % laplacian [logical] :: Use a Laplacian horizontal viscosity if true.

  • % biharmonic [logical] :: Use a biharmonic horizontal viscosity if true.

  • % debug [logical] :: If true, write verbose checksums for debugging purposes.

  • % no_slip [logical] :: If true, no slip boundary conditions are used. Otherwise free slip boundary conditions are assumed. The implementation of the free slip boundary conditions on a C-grid is much cleaner than the no slip boundary conditions. The use of free slip b.c.s is strongly encouraged. The no slip b.c.s are not implemented with the biharmonic viscosity.

  • % bound_kh [logical] :: If true, the Laplacian coefficient is locally limited to guarantee stability.

  • % better_bound_kh [logical] :: If true, use a more careful bounding of the Laplacian viscosity to guarantee stability.

  • % bound_ah [logical] :: If true, the biharmonic coefficient is locally limited to guarantee stability.

  • % better_bound_ah [logical] :: If true, use a more careful bounding of the biharmonic viscosity to guarantee stability.

  • % re_ah [real] :: so that the biharmonic Reynolds number is equal to this [nondim].

  • % bound_coef [real] :: The nondimensional coefficient of the ratio of the viscosity bounds to the theoretical maximum for stability without considering other terms [nondim]. The default is 0.8.

  • % smagorinsky_kh [logical] :: If true, use Smagorinsky nonlinear eddy viscosity. KH is the background value.

  • % smagorinsky_ah [logical] :: If true, use a biharmonic form of Smagorinsky nonlinear eddy viscosity. AH is the background.

  • % leith_kh [logical] :: If true, use 2D Leith nonlinear eddy viscosity. KH is the background value.

  • % modified_leith [logical] :: If true, use extra component of Leith viscosity to damp divergent flow. To use, still set Leith_Kh=.TRUE.

  • % use_beta_in_leith [logical] :: If true, includes the beta term in the Leith viscosity.

  • % leith_ah [logical] :: If true, use a biharmonic form of 2D Leith nonlinear eddy viscosity. AH is the background.

  • % use_leithy [logical] :: If true, use a biharmonic form of 2D Leith nonlinear eddy viscosity with harmonic backscatter. Ah is the background. Leithy = Leith+E.

  • % c_k [real] :: Fraction of energy dissipated by the biharmonic term that gets backscattered in the Leith+E scheme. [nondim].

  • % smooth_ah [logical] :: If true (default), then Ah and m_leithy are smoothed. This smoothing requires a lot of blocking communication.

  • % use_qg_leith_visc [logical] :: If true, use QG Leith nonlinear eddy viscosity. KH is the background value.

  • % bound_coriolis [logical] :: If true & SMAGORINSKY_AH is used, the biharmonic viscosity is modified to include a term that scales quadratically with the velocity shears.

  • % use_kh_bg_2d [logical] :: Read 2d background viscosity from a file.

  • % kh_bg_2d_bug [logical] :: If true, retain an answer-changing horizontal indexing bug in setting the corner-point viscosities when USE_KH_BG_2D=True.

  • % kh_bg_min [real] :: The minimum value allowed for Laplacian horizontal viscosity [L2 T-1 ~> m2 s-1]. The default is 0.0.

  • % use_land_mask [logical] :: Use the land mask for the computation of thicknesses at velocity locations. This eliminates the dependence on arbitrary values over land or outside of the domain. Default is False to maintain answers with legacy experiments but should be changed to True for new experiments.

  • % anisotropic [logical] :: If true, allow anisotropic component to the viscosity.

  • % add_les_viscosity [logical] :: If true, adds the viscosity from Smagorinsky and Leith to the background viscosity instead of taking the maximum.

  • % kh_aniso [real] :: The anisotropic viscosity [L2 T-1 ~> m2 s-1].

  • % dynamic_aniso [logical] :: If true, the anisotropic viscosity is recomputed as a function of state. This is set depending on ANISOTROPIC_MODE.

  • % res_scale_meke [logical] :: If true, the viscosity contribution from MEKE is scaled by the resolution function.

  • % use_gme [logical] :: If true, use GME backscatter scheme.

  • % answer_date [integer] :: The vintage of the order of arithmetic and expressions in the horizontal viscosity calculations. Values below 20190101 recover the answers from the end of 2018, while higher values use updated and more robust forms of the same expressions.

  • % gme_h0 [real] :: The strength of GME tapers quadratically to zero when the bathymetric total water column thickness is less than GME_H0 [H ~> m or kg m-2].

  • % gme_efficiency [real] :: The nondimensional prefactor multiplying the GME coefficient [nondim].

  • % gme_limiter [real] :: The absolute maximum value the GME coefficient is allowed to take [L2 T-1 ~> m2 s-1].

  • % min_grid_kh [real] :: Minimum horizontal Laplacian viscosity used to limit the grid Reynolds number [L2 T-1 ~> m2 s-1].

  • % min_grid_ah [real] :: Minimun horizontal biharmonic viscosity used to limit grid Reynolds number [L4 T-1 ~> m4 s-1].

  • % use_cont_thick [logical] :: If true, thickness at velocity points adopts h[uv] in BT_cont from continuity solver.

  • % zb2020 [type( zb2020_cs )] :: Zanna-Bolton 2020 control structure.

  • % use_zb2020 [logical] :: If true, use Zanna-Bolton 2020 parameterization.

  • % real (* re_ah_const_xy [*) :: The background Laplacian viscosity at h points [L2 T-1 ~> m2 s-1]. The actual viscosity may be the larger of this viscosity and the Smagorinsky and Leith viscosities.

  • % real :: The background Laplacian viscosity at h points [L2 T-1 ~> m2 s-1]. The actual viscosity may be the larger of this viscosity and the Smagorinsky and Leith viscosities.

  • % real :: The background biharmonic viscosity at h points [L4 T-1 ~> m4 s-1]. The actual viscosity may be the larger of this viscosity and the Smagorinsky and Leith viscosities.

  • % real :: The amount by which stresses through h points are reduced due to partial barriers [nondim].

  • % real :: The maximum permitted Laplacian viscosity [L2 T-1 ~> m2 s-1].

  • % real :: The maximum permitted biharmonic viscosity [L4 T-1 ~> m4 s-1].

  • % real :: Factor n1*n2 in the anisotropic direction tensor at h-points [nondim].

  • % real :: Factor n1**2-n2**2 in the anisotropic direction tensor at h-points [nondim].

  • % real :: Harmonic mean of the squares of the grid [L2 ~> m2].

  • % real :: Harmonic mean of the squares of the grid^(3/2) [L3 ~> m3].

  • % real :: The background Laplacian viscosity at q points [L2 T-1 ~> m2 s-1]. The actual viscosity may be the larger of this viscosity and the Smagorinsky and Leith viscosities.

  • % real :: The background biharmonic viscosity at q points [L4 T-1 ~> m4 s-1]. The actual viscosity may be the larger of this viscosity and the Smagorinsky and Leith viscosities.

  • % real :: The amount by which stresses through q points are reduced due to partial barriers [nondim].

  • % real :: The maximum permitted Laplacian viscosity [L2 T-1 ~> m2 s-1].

  • % real :: The maximum permitted biharmonic viscosity [L4 T-1 ~> m4 s-1].

  • % real :: Factor n1*n2 in the anisotropic direction tensor at q-points [nondim].

  • % real :: Factor n1**2-n2**2 in the anisotropic direction tensor at q-points [nondim].

  • % real :: Pre-calculated dx^2 at h points [L2 ~> m2].

  • % real :: Pre-calculated dy^2 at h points [L2 ~> m2].

  • % real :: Pre-calculated dx/dy at h points [nondim].

  • % real :: Pre-calculated dy/dx at h points [nondim].

  • % real :: Pre-calculated .5*sqrt(c_K)*max{dx,dy} [L ~> m].

  • % real :: Pre-calculated 4./max(dx,dy)^2 at h points [L-2 ~> m-2].

  • % real :: Pre-calculated dx^2 at q points [L2 ~> m2].

  • % real :: Pre-calculated dy^2 at q points [L2 ~> m2].

  • % real :: Pre-calculated dx/dy at q points [nondim].

  • % real :: Pre-calculated dy/dx at q points [nondim].

  • % real :: 1/(dx^2 dy) at u points [L-3 ~> m-3]

  • % real :: 1/(dx dy^2) at u points [L-3 ~> m-3]

  • % real :: 1/(dx^2 dy) at v points [L-3 ~> m-3]

  • % real :: 1/(dx dy^2) at v points [L-3 ~> m-3]

  • % real :: Laplacian metric-dependent constants [L2 ~> m2].

  • % real :: Biharmonic metric-dependent constants [L6 ~> m6].

  • % real :: Laplacian metric-dependent constants [L3 ~> m3].

  • % real :: Biharmonic metric-dependent constants [L4 ~> m4].

  • % real :: Biharmonic metric-dependent constants [T L4 ~> s m4].

  • % real :: Biharmonic metric-dependent constants [L3 ~> m3].

  • % real :: Laplacian metric-dependent constants [L2 ~> m2].

  • % real :: Biharmonic metric-dependent constants [L6 ~> m6].

  • % real :: Laplacian metric-dependent constants [L3 ~> m3].

  • % real :: Biharmonic metric-dependent constants [L4 ~> m4].

  • % real :: Biharmonic metric-dependent constants [T L4 ~> s m4].

  • % real :: Biharmonic metric-dependent constants [L3 ~> m3].

  • % diag [type( diag_ctrl ),pointer] :: structure to regulate diagnostics

  • % num_smooth_gme [integer] :: number of smoothing passes for the GME fluxes.

Function/Subroutine Documentation

subroutine mom_hor_visc/horizontal_viscosity(u, v, h, diffu, diffv, MEKE, VarMix, G, GV, US, CS, OBC, BT, TD, ADp, hu_cont, hv_cont)

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:
  • g :: [in] The ocean’s grid structure.

  • gv :: [in] The ocean’s vertical grid structure.

  • u :: [in] The zonal velocity [L T-1 ~> m s-1].

  • v :: [in] The meridional velocity [L T-1 ~> m s-1].

  • h :: [inout] Layer thicknesses [H ~> m or kg m-2].

  • diffu :: [out] Zonal acceleration due to convergence of

  • diffv :: [out] Meridional acceleration due to convergence

  • meke :: [inout] MEKE fields related to Mesoscale Eddy Kinetic Energy.

  • varmix :: [inout] Variable mixing control structure

  • us :: [in] A dimensional unit scaling type

  • cs :: [inout] Horizontal viscosity control structure

  • obc :: Pointer to an open boundary condition type

  • bt :: [in] Barotropic control structure

  • td :: [in] Thickness diffusion control structure

  • adp :: [in] Acceleration diagnostics

  • hu_cont :: [in] Layer thickness at u-points [H ~> m or kg m-2].

  • hv_cont :: [in] Layer thickness at v-points [H ~> m or kg m-2].

Call to:

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 smooth_gme smooth_x9_h smooth_x9_uv mom_thickness_diffuse::thickness_diffuse_get_kh mom_zanna_bolton::zb2020_copy_gradient_and_thickness mom_zanna_bolton::zb2020_lateral_stress

Called from:

mom_dynamics_unsplit::step_mom_dyn_unsplit mom_dynamics_unsplit_rk2::step_mom_dyn_unsplit_rk2

subroutine mom_hor_visc/hor_visc_init(Time, G, GV, US, param_file, diag, CS, ADp)

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

Parameters:
  • time :: [in] Current model time.

  • g :: [inout] The ocean’s grid structure.

  • gv :: [in] The ocean’s vertical grid structure

  • us :: [in] A dimensional unit scaling type

  • param_file :: [in] A structure to parse for run-time parameters.

  • diag :: [inout] Structure to regulate diagnostic output.

  • cs :: [inout] Horizontal viscosity control structure

  • adp :: [in] Acceleration diagnostics

Call to:

align_aniso_tensor_to_grid mom_error_handler::mom_error mom_zanna_bolton::zb2020_init

Called from:

mom_dynamics_split_rk2::initialize_dyn_split_rk2 mom_dynamics_split_rk2b::initialize_dyn_split_rk2b mom_dynamics_unsplit::initialize_dyn_unsplit mom_dynamics_unsplit_rk2::initialize_dyn_unsplit_rk2

subroutine mom_hor_visc/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.

Parameters:
  • cs :: [inout] Control structure for horizontal viscosity

  • n1 :: [in] i-component of direction vector [nondim]

  • n2 :: [in] j-component of direction vector [nondim]

Called from:

hor_visc_init

subroutine mom_hor_visc/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.

Parameters:
  • cs :: [in] Control structure

  • g :: [in] Ocean grid

  • gme_flux_h :: [inout] GME diffusive flux at h points [L2 T-2 ~> m2 s-2]

  • gme_flux_q :: [inout] GME diffusive flux at q points [L2 T-2 ~> m2 s-2]

Called from:

horizontal_viscosity

subroutine mom_hor_visc/smooth_x9_h(G, field_h, zero_land)

Apply a 9-point smoothing filter twice to a field staggered at a thickness point to reduce horizontal two-grid-point noise. Note that this subroutine does not conserve mass, so don’t use it in situations where you need conservation. Also note that it assumes that the input field has valid values in the first two halo points upon entry.

Parameters:
  • g :: [in] Ocean grid

  • field_h :: [inout] h-point field to be smoothed [arbitrary]

  • zero_land :: [in] If present and false, return the average of the surrounding ocean points when smoothing, otherwise use a value of 0 for land points and include them in the averages.

Called from:

horizontal_viscosity

subroutine mom_hor_visc/smooth_x9_uv(G, field_u, field_v, zero_land)

Apply a 9-point smoothing filter twice to a pair of velocity components to reduce horizontal two-grid-point noise. Note that this subroutine does not conserve angular momentum, so don’t use it in situations where you need conservation. Also note that it assumes that the input fields have valid values in the first two halo points upon entry.

Parameters:
  • g :: [in] Ocean grid

  • field_u :: [inout] u-point field to be smoothed [arbitrary]

  • field_v :: [inout] v-point field to be smoothed [arbitrary]

  • zero_land :: [in] If present and false, return the average of the surrounding ocean points when smoothing, otherwise use a value of 0 for land points and include them in the averages.

Called from:

horizontal_viscosity

subroutine mom_hor_visc/hor_visc_end(CS)

Deallocates any variables allocated in hor_visc_init.

Parameters:

cs :: [inout] Horizontal viscosity control structure

Call to:

mom_zanna_bolton::zb2020_end

Called from:

mom_dynamics_split_rk2::end_dyn_split_rk2 mom_dynamics_split_rk2b::end_dyn_split_rk2b