Calculates horizontal viscosity and viscous stresses
hor_visc_CS |
Control structure for horizontal viscosity |
horizontal_viscosity |
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] |
hor_visc_init |
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(). |
align_aniso_tensor_to_grid |
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. |
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 |
hor_visc_end |
Deallocates any variables allocated in hor_visc_init. |
namespace mom_hor_visc
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 f$kappa_hf$.
section section_horizontal_viscosity Horizontal viscosity in MOM
In general, the horizontal stress tensor can be written as f[ {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} f] where f$sigma_Df$, f$sigma_Tf$ and f$sigma_Sf$ 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 f[ 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} f] where f$dot{e}_D = partial_x u + partial_y vf$ is the horizontal divergence, f$dot{e}_T = partial_x u - partial_y vf$ is the horizontal tension, and f$dot{e}_S = partial_y u + partial_x vf$ is the horizontal shear strain.
The trace of the stress tensor, f$tr(bf sigma) = sigma_Df$, is usually absorbed into the pressure and only the deviatoric stress tensor considered. From here on, we drop f$sigma_Df$. The trace of the strain tensor, f$tr(bf e) = dot{e}_Df$ is non-zero for horizontally divergent flow but only enters the stress tensor through f$sigma_Df$ and so we will drop f$sigma_Df$ from calculations of the strain tensor in the code. Therefore the horizontal stress tensor can be considered to be f[ {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} .f]
The stresses above are linearly related to the strain through a viscosity coefficient, f$kappa_hf$: f{eqnarray*}{ sigma_T & = & 2 kappa_h dot{e}_T \\ sigma_S & = & 2 kappa_h dot{e}_S . f}
The viscosity f$kappa_hf$ may either be a constant or variable. For example, f$kappa_hf$ may vary with the shear, as proposed by Smagorinsky (1993).
The accelerations resulting form the divergence of the stress tensor are f{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) . f}
The form of the Laplacian viscosity in general coordinates is: f{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] . f}
subsection section_laplacian_viscosity_coefficient Laplacian viscosity coefficient
The horizontal viscosity coefficient, f$kappa_hf$, can have multiple components. The isotropic components are:
- A uniform background component, f$kappa_{bg}f$.
- A constant but spatially variable 2D map, f$kappa_{2d}(x,y)f$.
- A ''MICOM'' viscosity, f$U_nu Delta(x,y)f$, which uses a constant
velocity scale, f$U_nuf$ and a measure of the grid-spacing f$Delta(x,y)^2 = frac{2 Delta x^2 Delta y^2}{Delta x^2 + Delta y^2}f$.
- A function of
- latitude, f$kappa_{phi}(x,y) = kappa_{pi/2} |sin(phi)|^nf$.
A maximum stable viscosity, f$kappa_{max}(x,y)f$ is calculated based on the grid-spacing and time-step and used to clip calculated viscosities.
The static components of f$kappa_hf$ are first combined as follows: f[ kappa_{static} = min left[ maxleft( kappa_{bg}, U_nu Delta(x,y), kappa_{2d}(x,y), kappa_phi(x,y) right) , kappa_{max}(x,y) right] f] and stored in the module control structure as variables <code>Kh_bg_xx</code> and <code>Kh_bg_xy</code> for the tension (h-points) and shear (q-points) components respectively.
The full viscosity includes the dynamic components as follows: f[ kappa_h(x,y,t) = r(Delta,L_d) max left( kappa_{static}, kappa_{Sm}, kappa_{Lth} right) f] where f$r(Delta,L_d)f$ is a resolution function.
The dynamic Smagorinsky and Leith viscosity schemes are exclusive with each other.
subsection section_viscous_boundary_conditions 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: f[
partial_x v = 0 ; partial_x^3 v = 0 ; u = 0 ; partial_x^2 u = 0 ,
f] while for a Laplacian operator, they are simply f[
partial_x v = 0 ; u = 0 .
f] These boundary conditions are largely dictated by the use of an Arakawa C-grid and by the varying layer thickness.
subsection section_anisotropic_viscosity 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
f$A^primef$ and
The contributions to the stress tensor are f[ 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} f] Dissipation of kinetic energy requires f$kappa_h geq 0f$ and f$2 kappa_h + kappa_a geq 0f$. Note that when anisotropy is aligned with the x-direction, f$n_1 = pm 1f$, then f$n_2 = 0f$ and the cross terms vanish. The accelerations in this aligned limit with constant coefficients become f{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) f} which has contributions akin to a negative divergence damping (a divergence enhancement?) but which is weaker than the enhanced tension terms by half.
subsection section_viscous_discretization Discretization
The horizontal tension, f$dot{e}_Tf$, is stored in variable <code>sh_xx</code> and discretized as f[ 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) . f] The horizontal divergent strain, f$dot{e}_Df$, is stored in variable <code>div_xx</code> and discretized as f[ dot{e}_D = frac{1}{h A} left( delta_i left( overline{h}^i Delta y , u right) + delta_j left( overline{h}^jDelta x , v right) right) . f] Note that for expediency this is the exact discretization used in the continuity equation.
The horizontal shear strain, f$dot{e}_Sf$, is stored in variable <code>sh_xy</code> and discretized as f[ dot{e}_S = v_x + u_y f] where f{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) f} which are calculated separately so that no-slip or free-slip boundary conditions can be applied to f$v_xf$ and f$u_yf$ where appropriate.
The tendency for the x-component of the divergence of stress is stored in variable <code>diffu</code> and discretized as f[ 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) . f]
The tendency for the y-component of the divergence of stress is stored in variable <code>diffv</code> and discretized as f[ 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) . f]
subsection section_viscous_refs 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 , public :: hor_visc_CS
Control structure for horizontal viscosity
- Parameters
- Laplacian :: Use a Laplacian horizontal viscosity if true.
- biharmonic :: Use a biharmonic horizontal viscosity if true.
- debug :: If true, write verbose checksums for debugging purposes.
- no_slip :: 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 :: If true, the Laplacian coefficient is locally limited to guarantee stability.
- better_bound_Kh :: If true, use a more careful bounding of the Laplacian viscosity to guarantee stability.
- bound_Ah :: If true, the biharmonic coefficient is locally limited to guarantee stability.
- better_bound_Ah :: If true, use a more careful bounding of the biharmonic viscosity to guarantee stability.
- bound_coef :: 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 :: If true, use Smagorinsky nonlinear eddy viscosity. KH is the background value.
- Smagorinsky_Ah :: If true, use a biharmonic form of Smagorinsky nonlinear eddy viscosity. AH is the background.
- Leith_Kh :: If true, use 2D Leith nonlinear eddy viscosity. KH is the background value.
- Modified_Leith :: If true, use extra component of Leith viscosity to damp divergent flow. To use, still set Leith_Kh=.TRUE.
- use_beta_in_Leith :: If true, includes the beta term in the Leith viscosity
- Leith_Ah :: If true, use a biharmonic form of 2D Leith nonlinear eddy viscosity. AH is the background.
- use_QG_Leith_visc :: If true, use QG Leith nonlinear eddy viscosity. KH is the background value.
- bound_Coriolis :: 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 :: Read 2d background viscosity from a file.
- Kh_bg_2d_bug :: If true, retain an answer-changing horizontal indexing bug in setting the corner-point viscosities when USE_KH_BG_2D=True.
- Kh_bg_min :: The minimum value allowed for Laplacian horizontal viscosity [L2 T-1 ~> m2 s-1]. The default is 0.0.
- use_land_mask :: 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 :: If true, allow anisotropic component to the viscosity.
- add_LES_viscosity :: If true, adds the viscosity from Smagorinsky and Leith to the background viscosity instead of taking the maximum.
- Kh_aniso :: The anisotropic viscosity [L2 T-1 ~> m2 s-1].
- dynamic_aniso :: If true, the anisotropic viscosity is recomputed as a function of state. This is set depending on ANISOTROPIC_MODE.
- res_scale_MEKE :: If true, the viscosity contribution from MEKE is scaled by the resolution function.
- use_GME :: If true, use GME backscatter scheme.
- answers_2018 :: If true, use the order of arithmetic and expressions that recover the answers from the end of 2018. Otherwise, use updated and more robust forms of the same expressions.
- GME_h0 :: The strength of GME tapers quadratically to zero when the bathymetric depth is shallower than GME_H0 [Z ~> m]
- GME_efficiency :: The nondimensional prefactor multiplying the GME coefficient [nondim]
- GME_limiter :: The absolute maximum value the GME coefficient is allowed to take [L2 T-1 ~> m2 s-1].
- Kh_bg_xx :: 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.
- Kh_bg_2d :: 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.
- Ah_bg_xx :: 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.
- reduction_xx :: The amount by which stresses through h points are reduced due to partial barriers [nondim].
- Kh_Max_xx :: The maximum permitted Laplacian viscosity [L2 T-1 ~> m2 s-1].
- Ah_Max_xx :: The maximum permitted biharmonic viscosity [L4 T-1 ~> m4 s-1].
- n1n2_h :: Factor n1*n2 in the anisotropic direction tensor at h-points
- n1n1_m_n2n2_h :: Factor n1**2-n2**2 in the anisotropic direction tensor at h-points
- Kh_bg_xy :: 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.
- Ah_bg_xy :: 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.
- reduction_xy :: The amount by which stresses through q points are reduced due to partial barriers [nondim].
- Kh_Max_xy :: The maximum permitted Laplacian viscosity [L2 T-1 ~> m2 s-1].
- Ah_Max_xy :: The maximum permitted biharmonic viscosity [L4 T-1 ~> m4 s-1].
- n1n2_q :: Factor n1*n2 in the anisotropic direction tensor at q-points
- n1n1_m_n2n2_q :: Factor n1**2-n2**2 in the anisotropic direction tensor at q-points
- dx2h :: Pre-calculated dx^2 at h points [L2 ~> m2]
- dy2h :: Pre-calculated dy^2 at h points [L2 ~> m2]
- dx_dyT :: Pre-calculated dx/dy at h points [nondim]
- dy_dxT :: Pre-calculated dy/dx at h points [nondim]
- dx2q :: Pre-calculated dx^2 at q points [L2 ~> m2]
- dy2q :: Pre-calculated dy^2 at q points [L2 ~> m2]
- dx_dyBu :: Pre-calculated dx/dy at q points [nondim]
- dy_dxBu :: Pre-calculated dy/dx at q points [nondim]
- Idx2dyCu :: 1/(dx^2 dy) at u points [L-3 ~> m-3]
- Idxdy2u :: 1/(dx dy^2) at u points [L-3 ~> m-3]
- Idx2dyCv :: 1/(dx^2 dy) at v points [L-3 ~> m-3]
- Idxdy2v :: 1/(dx dy^2) at v points [L-3 ~> m-3]
- Laplac2_const_xx :: Laplacian metric-dependent constants [L2 ~> m2]
- Biharm5_const_xx :: Biharmonic metric-dependent constants [L5 ~> m5]
- Laplac3_const_xx :: Laplacian metric-dependent constants [L3 ~> m3]
- Biharm_const_xx :: Biharmonic metric-dependent constants [L4 ~> m4]
- Biharm_const2_xx :: Biharmonic metric-dependent constants [T L4 ~> s m4]
- Laplac2_const_xy :: Laplacian metric-dependent constants [L2 ~> m2]
- Biharm5_const_xy :: Biharmonic metric-dependent constants [L5 ~> m5]
- Laplac3_const_xy :: Laplacian metric-dependent constants [L3 ~> m3]
- Biharm_const_xy :: Biharmonic metric-dependent constants [L4 ~> m4]
- Biharm_const2_xy :: Biharmonic metric-dependent constants [T L4 ~> s m4]
- diag :: structure to regulate diagnostics !>@{ Diagnostic id
- subroutine horizontal_viscosity ( u , v , h , diffu , diffv , MEKE , VarMix , G , GV , US , CS , OBC , BT , TD )
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 along-coordinate stress tensor [L T-2 ~> m s-2]
- diffv :: [out] Meridional acceleration due to convergence of along-coordinate stress tensor [L T-2 ~> m s-2].
- MEKE :: Pointer to a structure containing fields related to Mesoscale Eddy Kinetic Energy.
- VarMix :: Pointer to a structure with fields that specify the spatially variable viscosities
- US :: [in] A dimensional unit scaling type
- CS :: Control structure returned by a previous call to hor_visc_init.
- OBC :: Pointer to an open boundary condition type
- BT :: Pointer to a structure containing barotropic velocities.
- TD :: Pointer to a structure containing thickness diffusivities.
- Ah_q :: GME coeff. at q-points [L2 T-1 ~> m2 s-1]
- KH_u_GME :: interface height diffusivities in u-columns [L2 T-1 ~> m2 s-1]
- KH_v_GME :: interface height diffusivities in v-columns [L2 T-1 ~> m2 s-1]
- GME_coeff_h :: GME coeff. at h-points [L2 T-1 ~> m2 s-1]
- Calls into
- |post_data|_ |calc_QG_Leith_viscosity|_ |thickness_diffuse_get_KH|_
smooth_GME
|barotropic_get_tav|_ |Bchksum|_ |hchksum|_ |MOM_error|_ |pass_vector|_ - Called by
- |MOM_dynamics_split_RK2::initialize_dyn_split_RK2|_ |MOM_dynamics_unsplit_RK2::step_MOM_dyn_unsplit_RK2|_ |MOM_dynamics_unsplit::step_MOM_dyn_unsplit|_ |MOM_dynamics_split_RK2::step_MOM_dyn_split_RK2|_
- subroutine 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().
- Parameters
- Time :: [in] Current model time.
- G :: [inout] The ocean's 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 :: Pointer to the control structure for this module
- MEKE :: MEKE data
- Calls into
- |register_diag_field|_ |slasher|_ |log_version|_ |MOM_read_data|_ |get_param|_
align_aniso_tensor_to_grid
|hchksum|_ |Bchksum|_ |MOM_error|_ |pass_var|_ |case|_ - Called by
- |MOM_dynamics_split_RK2::initialize_dyn_split_RK2|_ |MOM_dynamics_unsplit_RK2::initialize_dyn_unsplit_RK2|_ |MOM_dynamics_unsplit::initialize_dyn_unsplit|_
- 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.
- Parameters
- CS :: Control structure for horizontal viscosity
- n1 :: [in] i-component of direction vector [nondim]
- n2 :: [in] j-component of direction vector [nondim]
- Called by
hor_visc_init
- 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
- Parameters
- CS :: Control structure
- G :: [in] Ocean grid
- GME_flux_h :: [inout] GME diffusive flux at h points
- GME_flux_q :: [inout] GME diffusive flux at q points
- Calls into
- |pass_var|_
- Called by
horizontal_viscosity
- subroutine hor_visc_end ( CS )
Deallocates any variables allocated in hor_visc_init.
- Parameters
- CS :: The control structure returned by a previous call to hor_visc_init.