<p><b>mpetersen@lanl.gov</b> 2012-10-15 14:30:41 -0600 (Mon, 15 Oct 2012)</p><p>leith branch: Added Leith namelist and module.  Tests on 120km global match bit-for-bit with previous branch for both del2 and Leith.  With this implementation, Leith and del2 (as well as del4) can be run simultaniously.<br>
</p><hr noshade><pre><font color="gray">Modified: branches/ocean_projects/leith_mrp/namelist.input.ocean
===================================================================
--- branches/ocean_projects/leith_mrp/namelist.input.ocean        2012-10-11 22:58:41 UTC (rev 2208)
+++ branches/ocean_projects/leith_mrp/namelist.input.ocean        2012-10-15 20:30:41 UTC (rev 2209)
@@ -55,6 +55,12 @@
    config_h_tracer_eddy_diff2 = 1.0e5
    config_h_tracer_eddy_diff4 = 0.0
 /
+&amp;hmix_leith
+   config_use_leith_del2 = .false.
+   config_leith_parameter = 1.0
+   config_leith_dx = 15000.0
+   config_leith_visc2_max = 2.5e3
+/
 &amp;vmix
    config_vert_visc_type  = 'const'
    config_vert_diff_type  = 'const'

Modified: branches/ocean_projects/leith_mrp/src/core_ocean/Makefile
===================================================================
--- branches/ocean_projects/leith_mrp/src/core_ocean/Makefile        2012-10-11 22:58:41 UTC (rev 2208)
+++ branches/ocean_projects/leith_mrp/src/core_ocean/Makefile        2012-10-15 20:30:41 UTC (rev 2209)
@@ -10,6 +10,7 @@
            mpas_ocn_vel_vadv.o \
            mpas_ocn_vel_hmix.o \
            mpas_ocn_vel_hmix_del2.o \
+           mpas_ocn_vel_hmix_leith.o \
            mpas_ocn_vel_hmix_del4.o \
            mpas_ocn_vel_forcing.o \
            mpas_ocn_vel_forcing_windstress.o \
@@ -87,10 +88,12 @@
 
 mpas_ocn_vel_vadv.o:
 
-mpas_ocn_vel_hmix.o: mpas_ocn_vel_hmix_del2.o mpas_ocn_vel_hmix_del4.o
+mpas_ocn_vel_hmix.o: mpas_ocn_vel_hmix_del2.o mpas_ocn_vel_hmix_leith.o mpas_ocn_vel_hmix_del4.o
 
 mpas_ocn_vel_hmix_del2.o:
 
+mpas_ocn_vel_hmix_leith.o:
+
 mpas_ocn_vel_hmix_del4.o:
 
 mpas_ocn_vel_forcing.o: mpas_ocn_vel_forcing_windstress.o mpas_ocn_vel_forcing_bottomdrag.o mpas_ocn_vel_forcing_rayleigh.o
@@ -179,6 +182,7 @@
                                           mpas_ocn_vel_vadv.o \
                                           mpas_ocn_vel_hmix.o \
                                           mpas_ocn_vel_hmix_del2.o \
+                                          mpas_ocn_vel_hmix_leith.o \
                                           mpas_ocn_vel_hmix_del4.o \
                                           mpas_ocn_vel_forcing.o \
                                           mpas_ocn_vel_forcing_windstress.o \

Modified: branches/ocean_projects/leith_mrp/src/core_ocean/Registry
===================================================================
--- branches/ocean_projects/leith_mrp/src/core_ocean/Registry        2012-10-11 22:58:41 UTC (rev 2208)
+++ branches/ocean_projects/leith_mrp/src/core_ocean/Registry        2012-10-15 20:30:41 UTC (rev 2209)
@@ -58,10 +58,10 @@
 namelist logical   hmix     config_rayleigh_friction    false
 namelist real      hmix     config_rayleigh_damping_coeff 0.0
 namelist real      hmix     config_apvm_scale_factor      0.0
-namelist logical   hmix     config_use_leith_del2         false
-namelist real      hmix     config_leith_parameter        0.0
-namelist real      hmix     config_leith_dx               0.0
-namelist real      hmix     config_leith_visc2_max      1000000.0
+namelist logical   hmix_leith  config_use_leith_del2         false
+namelist real      hmix_leith  config_leith_parameter        0.0
+namelist real      hmix_leith  config_leith_dx               0.0
+namelist real      hmix_leith  config_leith_visc2_max      1000000.0
 namelist character vmix     config_vert_visc_type       const
 namelist character vmix     config_vert_diff_type       const
 namelist logical   vmix     config_implicit_vertical_mix  .true.

Modified: branches/ocean_projects/leith_mrp/src/core_ocean/mpas_ocn_vel_hmix.F
===================================================================
--- branches/ocean_projects/leith_mrp/src/core_ocean/mpas_ocn_vel_hmix.F        2012-10-11 22:58:41 UTC (rev 2208)
+++ branches/ocean_projects/leith_mrp/src/core_ocean/mpas_ocn_vel_hmix.F        2012-10-15 20:30:41 UTC (rev 2209)
@@ -20,6 +20,7 @@
    use mpas_configure
    use mpas_timer
    use ocn_vel_hmix_del2
+   use ocn_vel_hmix_leith
    use ocn_vel_hmix_del4
 
    implicit none
@@ -47,7 +48,7 @@
    !
    !--------------------------------------------------------------------
 
-   type (timer_node), pointer :: del2Timer, del4Timer
+   type (timer_node), pointer :: del2Timer, leithTimer, del4Timer
 
 
 !***********************************************************************
@@ -101,7 +102,6 @@
       real (kind=RKIND), dimension(:,:), intent(inout) :: &amp;
          viscosity     !&lt; Input: viscosity
 
-
       !-----------------------------------------------------------------
       !
       ! output variables
@@ -116,7 +116,7 @@
       !
       !-----------------------------------------------------------------
 
-      integer :: err1, err2
+      integer :: err1, err2, err3
 
       !-----------------------------------------------------------------
       !
@@ -127,13 +127,18 @@
       !-----------------------------------------------------------------
 
       call mpas_timer_start(&quot;del2&quot;, .false., del2Timer)
-      call ocn_vel_hmix_del2_tend(grid, divergence, vorticity, viscosity, tend, err1)
+      call ocn_vel_hmix_del2_tend(grid, divergence, vorticity, tend, err1)
       call mpas_timer_stop(&quot;del2&quot;, del2Timer)
+
+      call mpas_timer_start(&quot;leith&quot;, .false., leithTimer)
+      call ocn_vel_hmix_leith_tend(grid, divergence, vorticity, viscosity, tend, err2)
+      call mpas_timer_stop(&quot;leith&quot;, leithTimer)
+
       call mpas_timer_start(&quot;del4&quot;, .false., del4Timer)
-      call ocn_vel_hmix_del4_tend(grid, divergence, vorticity, tend, err2)
+      call ocn_vel_hmix_del4_tend(grid, divergence, vorticity, tend, err3)
       call mpas_timer_stop(&quot;del4&quot;, del4Timer)
 
-      err = ior(err1, err2)
+      err = ior(ior(err1, err2),err3)
 
    !--------------------------------------------------------------------
 
@@ -167,12 +172,13 @@
 
       integer, intent(out) :: err !&lt; Output: error flag
 
-      integer :: err1, err2
+      integer :: err1, err2, err3
 
       call ocn_vel_hmix_del2_init(err1)
-      call ocn_vel_hmix_del4_init(err2)
+      call ocn_vel_hmix_leith_init(err2)
+      call ocn_vel_hmix_del4_init(err3)
 
-      err = ior(err1, err2)
+      err = ior(ior(err1, err2),err3)
 
    !--------------------------------------------------------------------
 

Modified: branches/ocean_projects/leith_mrp/src/core_ocean/mpas_ocn_vel_hmix_del2.F
===================================================================
--- branches/ocean_projects/leith_mrp/src/core_ocean/mpas_ocn_vel_hmix_del2.F        2012-10-11 22:58:41 UTC (rev 2208)
+++ branches/ocean_projects/leith_mrp/src/core_ocean/mpas_ocn_vel_hmix_del2.F        2012-10-15 20:30:41 UTC (rev 2209)
@@ -70,7 +70,7 @@
 !
 !-----------------------------------------------------------------------
 
-   subroutine ocn_vel_hmix_del2_tend(grid, divergence, vorticity, viscosity, tend, err)!{{{
+   subroutine ocn_vel_hmix_del2_tend(grid, divergence, vorticity, tend, err)!{{{
 
       !-----------------------------------------------------------------
       !
@@ -96,8 +96,6 @@
       real (kind=RKIND), dimension(:,:), intent(inout) :: &amp;
          tend             !&lt; Input/Output: velocity tendency
 
-      real (kind=RKIND), dimension(:,:), intent(inout) :: &amp;
-         viscosity       !&lt; Input: viscosity
 
       !-----------------------------------------------------------------
       !
@@ -118,8 +116,8 @@
       integer, dimension(:), pointer :: maxLevelEdgeTop
       integer, dimension(:,:), pointer :: cellsOnEdge, verticesOnEdge, edgeMask
 
-      real (kind=RKIND) :: u_diffusion, invLength1, invLength2, visc2
-      real (kind=RKIND), dimension(:), pointer :: meshScalingDel2, meshScaling, &amp;
+      real (kind=RKIND) :: u_diffusion, invLength1, invLength2
+      real (kind=RKIND), dimension(:), pointer :: meshScalingDel2, &amp;
               dcEdge, dvEdge
 
       !-----------------------------------------------------------------
@@ -137,13 +135,10 @@
       cellsOnEdge =&gt; grid % cellsOnEdge % array
       verticesOnEdge =&gt; grid % verticesOnEdge % array
       meshScalingDel2 =&gt; grid % meshScalingDel2 % array
-      meshScaling =&gt; grid % meshScaling % array
       edgeMask =&gt; grid % edgeMask % array
       dcEdge =&gt; grid % dcEdge % array
       dvEdge =&gt; grid % dvEdge % array
 
-      viscosity = 1.0
-
       do iEdge=1,nEdgesSolve
          cell1 = cellsOnEdge(1,iEdge)
          cell2 = cellsOnEdge(2,iEdge)
@@ -163,19 +158,8 @@
                           -viscVortCoef &amp;
                           *( vorticity(k,vertex2) - vorticity(k,vertex1) ) * invLength2
 
-            if(config_use_leith_del2) then
-               visc2 = ( config_leith_parameter * config_leith_dx * meshScaling(iEdge) / 3.14)**3 &amp;
-                        * abs( vorticity(k,vertex2) - vorticity(k,vertex1) ) * invLength1 * sqrt(3.0)
-               visc2 = min(visc2, config_leith_visc2_max)
-               u_diffusion = visc2 * u_diffusion
-            else
-               visc2 = eddyVisc2
-               u_diffusion = meshScalingDel2(iEdge) * visc2 * u_diffusion
-            endif
+            u_diffusion = meshScalingDel2(iEdge) * eddyVisc2 * u_diffusion
 
-            viscosity(k,cell1) = visc2
-            viscosity(k,cell2) = visc2
-
             tend(k,iEdge) = tend(k,iEdge) + edgeMask(k, iEdge) * u_diffusion
 
          end do
@@ -214,7 +198,7 @@
 
    hmixDel2On = .false.
 
-   if ( config_h_mom_eddy_visc2 &gt; 0.0 .or. config_use_leith_del2 ) then
+   if ( config_h_mom_eddy_visc2 &gt; 0.0 ) then
       hmixDel2On = .true.
       eddyVisc2 = config_h_mom_eddy_visc2
 

Added: branches/ocean_projects/leith_mrp/src/core_ocean/mpas_ocn_vel_hmix_leith.F
===================================================================
--- branches/ocean_projects/leith_mrp/src/core_ocean/mpas_ocn_vel_hmix_leith.F                                (rev 0)
+++ branches/ocean_projects/leith_mrp/src/core_ocean/mpas_ocn_vel_hmix_leith.F        2012-10-15 20:30:41 UTC (rev 2209)
@@ -0,0 +1,231 @@
+!|||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||
+!
+!  ocn_vel_hmix_leith
+!
+!&gt; \brief Ocean horizontal mixing - Laplacian parameterization 
+!&gt; \author Phil Jones, Doug Jacobsen
+!&gt; \date   15 September 2011
+!&gt; \version SVN:$Id:$
+!&gt; \details
+!&gt;  This module contains routines for computing horizontal mixing 
+!&gt;  tendencies using a Laplacian formulation.
+!
+!-----------------------------------------------------------------------
+
+module ocn_vel_hmix_leith
+
+   use mpas_grid_types
+   use mpas_configure
+
+   implicit none
+   private
+   save
+
+   !--------------------------------------------------------------------
+   !
+   ! Public parameters
+   !
+   !--------------------------------------------------------------------
+
+   !--------------------------------------------------------------------
+   !
+   ! Public member functions
+   !
+   !--------------------------------------------------------------------
+
+   public :: ocn_vel_hmix_leith_tend, &amp;
+             ocn_vel_hmix_leith_init
+
+   !-------------------------------------------------------------------
+   !
+   ! Private module variables
+   !
+   !--------------------------------------------------------------------
+
+   logical ::  hmixLeithOn  !&lt; integer flag to determine whether leith chosen
+
+   real (kind=RKIND) :: &amp;
+      eddyVisc2,        &amp;!&lt; base eddy diffusivity for Laplacian
+      viscVortCoef
+
+!***********************************************************************
+
+contains
+
+!***********************************************************************
+!
+!  routine ocn_vel_hmix_leith_tend
+!
+!&gt; \brief   Computes tendency term for Laplacian horizontal momentum mixing
+!&gt; \author  Phil Jones, Doug Jacobsen
+!&gt; \date    22 August 2011
+!&gt; \version SVN:$Id$
+!&gt; \details 
+!&gt;  This routine computes the horizontal mixing tendency for momentum
+!&gt;  based on a Laplacian form for the mixing, \f$</font>
<font color="black">u_2 </font>
<font color="blue">abla^2 u\f$
+!&gt;  This tendency takes the
+!&gt;  form \f$</font>
<font color="black">u( </font>
<font color="black">abla divergence + k \times </font>
<font color="blue">abla vorticity )\f$,
+!&gt;  where \f$</font>
<font color="blue">u\f$ is a viscosity and \f$k\f$ is the vertical unit vector.
+!&gt;  This form is strictly only valid for constant \f$</font>
<font color="blue">u\f$ .
+!
+!-----------------------------------------------------------------------
+
+   subroutine ocn_vel_hmix_leith_tend(grid, divergence, vorticity, viscosity, tend, err)!{{{
+
+      !-----------------------------------------------------------------
+      !
+      ! input variables
+      !
+      !-----------------------------------------------------------------
+
+      real (kind=RKIND), dimension(:,:), intent(in) :: &amp;
+         divergence      !&lt; Input: velocity divergence
+
+      real (kind=RKIND), dimension(:,:), intent(in) :: &amp;
+         vorticity       !&lt; Input: vorticity
+
+      type (mesh_type), intent(in) :: &amp;
+         grid            !&lt; Input: grid information
+
+      !-----------------------------------------------------------------
+      !
+      ! input/output variables
+      !
+      !-----------------------------------------------------------------
+
+      real (kind=RKIND), dimension(:,:), intent(inout) :: &amp;
+         tend             !&lt; Input/Output: velocity tendency
+
+      real (kind=RKIND), dimension(:,:), intent(inout) :: &amp;
+         viscosity       !&lt; Input: viscosity
+
+      !-----------------------------------------------------------------
+      !
+      ! output variables
+      !
+      !-----------------------------------------------------------------
+
+      integer, intent(out) :: err !&lt; Output: error flag
+
+      !-----------------------------------------------------------------
+      !
+      ! local variables
+      !
+      !-----------------------------------------------------------------
+
+      integer :: iEdge, nEdgesSolve, cell1, cell2, vertex1, vertex2
+      integer :: k
+      integer, dimension(:), pointer :: maxLevelEdgeTop
+      integer, dimension(:,:), pointer :: cellsOnEdge, verticesOnEdge, edgeMask
+
+      real (kind=RKIND) :: u_diffusion, invLength1, invLength2, visc2
+      real (kind=RKIND), dimension(:), pointer :: meshScaling, &amp;
+              dcEdge, dvEdge
+
+      !-----------------------------------------------------------------
+      !
+      ! exit if this mixing is not selected
+      !
+      !-----------------------------------------------------------------
+
+      err = 0
+
+      if(.not.hmixLeithOn) return
+
+      nEdgesSolve = grid % nEdgesSolve
+      maxLevelEdgeTop =&gt; grid % maxLevelEdgeTop % array
+      cellsOnEdge =&gt; grid % cellsOnEdge % array
+      verticesOnEdge =&gt; grid % verticesOnEdge % array
+      meshScaling =&gt; grid % meshScaling % array
+      edgeMask =&gt; grid % edgeMask % array
+      dcEdge =&gt; grid % dcEdge % array
+      dvEdge =&gt; grid % dvEdge % array
+
+      viscosity = 1.0
+
+      do iEdge=1,nEdgesSolve
+         cell1 = cellsOnEdge(1,iEdge)
+         cell2 = cellsOnEdge(2,iEdge)
+         vertex1 = verticesOnEdge(1,iEdge)
+         vertex2 = verticesOnEdge(2,iEdge)
+
+         invLength1 = 1.0 / dcEdge(iEdge)
+         invLength2 = 1.0 / dvEdge(iEdge)
+
+         do k=1,maxLevelEdgeTop(iEdge)
+
+            ! Here -( vorticity(k,vertex2) - vorticity(k,vertex1) ) / dvEdge(iEdge)
+            ! is - </font>
<font color="blue">abla vorticity pointing from vertex 2 to vertex 1, or equivalently 
+            !    + k \times </font>
<font color="blue">abla vorticity pointing from cell1 to cell2.
+
+            u_diffusion = ( divergence(k,cell2)  - divergence(k,cell1) ) * invLength1 &amp;
+                          -viscVortCoef &amp;
+                          *( vorticity(k,vertex2) - vorticity(k,vertex1) ) * invLength2
+
+            visc2 = ( config_leith_parameter * config_leith_dx * meshScaling(iEdge) / 3.14)**3 &amp;
+                     * abs( vorticity(k,vertex2) - vorticity(k,vertex1) ) * invLength1 * sqrt(3.0)
+            visc2 = min(visc2, config_leith_visc2_max)
+            u_diffusion = visc2 * u_diffusion
+
+            viscosity(k,cell1) = visc2
+            viscosity(k,cell2) = visc2
+
+            tend(k,iEdge) = tend(k,iEdge) + edgeMask(k, iEdge) * u_diffusion
+
+         end do
+      end do
+
+   !--------------------------------------------------------------------
+
+   end subroutine ocn_vel_hmix_leith_tend!}}}
+
+!***********************************************************************
+!
+!  routine ocn_vel_hmix_leith_init
+!
+!&gt; \brief   Initializes ocean momentum Laplacian horizontal mixing
+!&gt; \author  Phil Jones, Doug Jacobsen
+!&gt; \date    15 September 2011
+!&gt; \version SVN:$Id$
+!&gt; \details 
+!&gt;  This routine initializes a variety of quantities related to 
+!&gt;  Laplacian horizontal momentum mixing in the ocean.  
+!
+!-----------------------------------------------------------------------
+
+   subroutine ocn_vel_hmix_leith_init(err)!{{{
+
+
+   integer, intent(out) :: err !&lt; Output: error flag
+
+   !--------------------------------------------------------------------
+   !
+   ! set some local module variables based on input config choices
+   !
+   !--------------------------------------------------------------------
+
+   err = 0
+
+   hmixLeithOn = .false.
+
+   if (config_use_leith_del2) then
+      hmixLeithOn = .true.
+
+      if (config_visc_vorticity_term) then
+         viscVortCoef = config_visc_vorticity_visc2_scale
+      else
+         viscVortCoef = 0.0
+      endif
+
+   endif
+
+   !--------------------------------------------------------------------
+
+   end subroutine ocn_vel_hmix_leith_init!}}}
+
+!***********************************************************************
+
+end module ocn_vel_hmix_leith
+
+!|||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||
+! vim: foldmethod=marker

</font>
</pre>