<p><b>dwj07@fsu.edu</b> 2011-08-30 13:15:29 -0600 (Tue, 30 Aug 2011)</p><p><br>
        Modularizing and masking arrays in Vertical Mixing and horizontal pressure gradient.<br>
</p><hr noshade><pre><font color="gray">Modified: branches/ocean_projects/performance/mpas/src/core_ocean/Makefile
===================================================================
--- branches/ocean_projects/performance/mpas/src/core_ocean/Makefile        2011-08-30 18:03:51 UTC (rev 963)
+++ branches/ocean_projects/performance/mpas/src/core_ocean/Makefile        2011-08-30 19:15:29 UTC (rev 964)
@@ -16,8 +16,10 @@
        module_OcnVmixVelConst.o \
        module_OcnVmixVelTanh.o \
        module_OcnCoriolis.o \
+           module_OcnPressureGrad.o \
            module_OcnThickness.o \
            module_OcnVelocityForcing.o \
+           module_OcnVelVertAdv.o \
        module_time_integration.o \
        module_global_diagnostics.o
 
@@ -34,8 +36,12 @@
 
 module_advection.o:
 
+module_OcnPressureGrad.o:
+
 module_OcnCoriolis.o:
 
+module_OcnVelVertAdv.o:
+
 module_OcnVmixTracerConst.o:
 
 module_OcnVmixTracerTanh.o:
@@ -64,7 +70,7 @@
 
 module_global_diagnostics.o: 
 
-module_mpas_core.o: module_advection.o module_global_diagnostics.o module_test_cases.o module_OcnVelocityForcing.o module_OcnThickness.o module_OcnCoriolis.o module_OcnVmixVel.o module_OcnVmixVelConst.o module_OcnVmixVelTanh.o module_OcnVmixTracer.o module_OcnVmixTracerConst.o module_OcnVmixTracerTanh.o module_OcnHmixTracer.o module_OcnHmixTracerDel2.o module_OcnHmixTracerDel4.o module_OcnHmixVel.o module_OcnHmixVelDel2.o module_OcnHmixVelDel4.o module_time_integration.o
+module_mpas_core.o: module_advection.o module_global_diagnostics.o module_test_cases.o module_OcnVelVertAdv.o module_OcnPressureGrad.o module_OcnVelocityForcing.o module_OcnThickness.o module_OcnCoriolis.o module_OcnVmixVel.o module_OcnVmixVelConst.o module_OcnVmixVelTanh.o module_OcnVmixTracer.o module_OcnVmixTracerConst.o module_OcnVmixTracerTanh.o module_OcnHmixTracer.o module_OcnHmixTracerDel2.o module_OcnHmixTracerDel4.o module_OcnHmixVel.o module_OcnHmixVelDel2.o module_OcnHmixVelDel4.o module_time_integration.o
 
 clean:
         $(RM) *.o *.mod *.f90 libdycore.a

Added: branches/ocean_projects/performance/mpas/src/core_ocean/module_OcnPressureGrad.F
===================================================================
--- branches/ocean_projects/performance/mpas/src/core_ocean/module_OcnPressureGrad.F                                (rev 0)
+++ branches/ocean_projects/performance/mpas/src/core_ocean/module_OcnPressureGrad.F        2011-08-30 19:15:29 UTC (rev 964)
@@ -0,0 +1,152 @@
+!|||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||
+!
+!  OcnPressureGrad
+!
+!&gt; \brief MPAS ocean horizontal pressure gradient tendency
+!&gt; \author Doug Jacobsen
+!&gt; \date   30 August 2011
+!&gt; \version SVN:$Id:$
+!&gt; \details
+!&gt;  This module contains the routine for computing 
+!&gt;  tendencies from horizontal pressure gradient
+!
+!-----------------------------------------------------------------------
+
+module OcnPressureGrad
+
+   use grid_types
+   use configure
+
+   implicit none
+   private
+   save
+
+   !--------------------------------------------------------------------
+   !
+   ! Public parameters
+   !
+   !--------------------------------------------------------------------
+
+   !--------------------------------------------------------------------
+   !
+   ! Public member functions
+   !
+   !--------------------------------------------------------------------
+
+   public :: OcnPressureGradTend
+
+   !--------------------------------------------------------------------
+   !
+   ! Private module variables
+   !
+   !--------------------------------------------------------------------
+
+
+!***********************************************************************
+
+contains
+
+!***********************************************************************
+!
+!  routine OcnPressureGradTend
+!
+!&gt; \brief   Computes tendency term from pressure gradient
+!&gt; \author  Doug Jacobsen
+!&gt; \date    30 August 2011
+!&gt; \version SVN:$Id$
+!&gt; \details 
+!&gt;  This routine computes the tendency from the pressure gradient for the
+!&gt;  momentum equation.
+!
+!-----------------------------------------------------------------------
+
+   subroutine OcnPressureGradTend(tend, MontPot, pressure, grid)
+
+      !-----------------------------------------------------------------
+      !
+      ! input variables
+      !
+      !-----------------------------------------------------------------
+
+      real (kind=RKIND), dimension(:,:), intent(in) :: &amp;
+         MontPot, pressure
+
+      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
+
+      !-----------------------------------------------------------------
+      !
+      ! output variables
+      !
+      !-----------------------------------------------------------------
+
+      !-----------------------------------------------------------------
+      !
+      ! local variables
+      !
+      !-----------------------------------------------------------------
+
+      integer :: iEdge, k, cell1, cell2
+      real (kind=RKIND) :: rho0Inv, invDCEdge
+
+      integer, dimension(:), pointer :: maxLeveLEdgeTop
+
+      integer, dimension(:,:), pointer :: cellsOnEdge
+
+      real (kind=RKIND), dimension(:), pointer :: dcEdge
+
+      cellsOnEdge =&gt; grid % cellsOnEdge % array
+      maxLevelEdgeTop =&gt; grid % maxLevelEdgeTop % array
+      dcEdge =&gt; grid % dcEdge % array
+
+      !
+      ! velocity tendency: pressure gradient
+      !
+      rho0Inv = 1.0/config_rho0
+      if (config_vert_grid_type.eq.'isopycnal') then
+        do iEdge=1,grid % nEdgesSolve
+          cell1 = cellsOnEdge(1,iEdge)
+          cell2 = cellsOnEdge(2,iEdge)
+
+          invDCEdge = 1.0 / dcEdge(iEdge)
+
+          do k=1,maxLevelEdgeTop(iEdge)
+             tend(k,iEdge) = tend(k,iEdge)     &amp;
+               - (MontPot(k,cell2) - MontPot(k,cell1))*invDCEdge
+           end do
+        enddo
+      elseif (config_vert_grid_type.eq.'zlevel') then
+        do iEdge=1,grid % nEdgesSolve
+          cell1 = cellsOnEdge(1,iEdge)
+          cell2 = cellsOnEdge(2,iEdge)
+
+          invDCEdge = 1.0 / dcEdge(iEdge)
+
+          do k=1,maxLevelEdgeTop(iEdge)
+             tend(k,iEdge) = tend(k,iEdge)     &amp;
+               - rho0Inv*(  pressure(k,cell2) &amp;
+                          - pressure(k,cell1) )*invDCEdge
+          end do
+        enddo
+      endif
+
+
+
+   !--------------------------------------------------------------------
+
+   end subroutine OcnPressureGradTend
+
+!***********************************************************************
+
+end module OcnPressureGrad
+
+!|||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||

Added: branches/ocean_projects/performance/mpas/src/core_ocean/module_OcnVelVertAdv.F
===================================================================
--- branches/ocean_projects/performance/mpas/src/core_ocean/module_OcnVelVertAdv.F                                (rev 0)
+++ branches/ocean_projects/performance/mpas/src/core_ocean/module_OcnVelVertAdv.F        2011-08-30 19:15:29 UTC (rev 964)
@@ -0,0 +1,186 @@
+!|||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||
+!
+!  OcnVelVertAdv
+!
+!&gt; \brief MPAS Ocean vertical advection for momentum
+!&gt; \author Doug Jacobsen
+!&gt; \date   30 August 2011
+!&gt; \version SVN:$Id:$
+!&gt; \details
+!&gt;  This module contains routines for computing vertical advection 
+!&gt;  tendencies 
+!
+!-----------------------------------------------------------------------
+
+module OcnVelVertAdv
+
+   use grid_types
+   use configure
+
+   implicit none
+   private
+   save
+
+   !--------------------------------------------------------------------
+   !
+   ! Public parameters
+   !
+   !--------------------------------------------------------------------
+
+   !--------------------------------------------------------------------
+   !
+   ! Public member functions
+   !
+   !--------------------------------------------------------------------
+
+   public :: OcnVelVertAdvTend, &amp;
+             OcnVelVertAdvInit
+
+   !-------------------------------------------------------------------
+   !
+   ! Private module variables
+   !
+   !--------------------------------------------------------------------
+
+   logical :: &amp;
+      vertAdvOn         !&lt; local flag to determine whether Const chosen
+
+!***********************************************************************
+
+contains
+
+!***********************************************************************
+!
+!  routine OcnVelVertAdvTend
+!
+!&gt; \brief   Computes tendency term from vertical advection
+!&gt; \author  Doug Jacobsen
+!&gt; \date    30 August 2011
+!&gt; \version SVN:$Id$
+!&gt; \details 
+!&gt;  This routine computes the vertical advection tendency for momentum
+!
+!-----------------------------------------------------------------------
+
+   subroutine OcnVelVertAdvTend(tend, u, wTop, grid)
+
+      !-----------------------------------------------------------------
+      !
+      ! input variables
+      !
+      !-----------------------------------------------------------------
+
+      real (kind=RKIND), dimension(:,:), intent(in) :: &amp;
+                                 u,  &amp; !&lt; Input: current velocity 
+                                 wTop  !&lt; Input: top layer vertical velocity
+
+      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
+
+      !-----------------------------------------------------------------
+      !
+      ! local variables
+      !
+      !-----------------------------------------------------------------
+
+      integer ::               &amp;
+         iEdge, k,             &amp;! loop counters
+         cell1, cell2
+
+      real (kind=RKIND) :: wTopEdge
+
+      integer, dimension(:), pointer :: maxLevelEdgeTop
+      integer, dimension(:,:), pointer :: cellsOnEdge
+
+      real (kind=RKIND), dimension(:), pointer :: zMidZLevel, w_dudzTopEdge
+
+      !-----------------------------------------------------------------
+      !
+      ! exit if this mixing is not selected
+      !
+      !-----------------------------------------------------------------
+
+      if (.not. vertAdvOn) return
+
+      maxLevelEdgeTop =&gt; grid % maxLevelEdgeTop % array
+      zMidZLevel =&gt; grid % zMidZLevel % array
+      cellsOnEdge =&gt; grid % cellsOnEdge % array
+
+      !
+      ! velocity tendency: vertical advection term -w du/dz
+      !
+
+      allocate(w_dudzTopEdge(grid % nVertLevels+1))
+      w_dudzTopEdge(1) = 0.0
+      do iEdge=1,grid % nEdgesSolve
+        cell1 = cellsOnEdge(1,iEdge)
+        cell2 = cellsOnEdge(2,iEdge)
+
+        do k=2,maxLevelEdgeTop(iEdge)
+          ! Average w from cell center to edge
+          wTopEdge = 0.5*(wTop(k,cell1)+wTop(k,cell2))
+
+          ! compute dudz at vertical interface with first order derivative.
+          w_dudzTopEdge(k) = wTopEdge * (u(k-1,iEdge)-u(k,iEdge)) &amp;
+                       / (zMidZLevel(k-1) - zMidZLevel(k))
+        end do
+        w_dudzTopEdge(maxLevelEdgeTop(iEdge)+1) = 0.0
+
+        ! Average w*du/dz from vertical interface to vertical middle of cell
+        do k=1,maxLevelEdgeTop(iEdge)
+           tend(k,iEdge) = - 0.5 * (w_dudzTopEdge(k) + w_dudzTopEdge(k+1))
+        enddo
+      enddo
+      deallocate(w_dudzTopEdge)
+
+
+   !--------------------------------------------------------------------
+
+   end subroutine OcnVelVertAdvTend
+
+!***********************************************************************
+!
+!  routine OcnVelVertAdvInit
+!
+!&gt; \brief   Initializes ocean momentum constant vertical mixing
+!&gt; \author  Doug Jacobsen
+!&gt; \date    30 August 2011
+!&gt; \version SVN:$Id$
+!&gt; \details 
+!&gt;  This routine initializes a variety of quantities related to 
+!&gt;  constant vertical momentum mixing in the ocean.  
+!
+!-----------------------------------------------------------------------
+
+   subroutine OcnVelVertAdvInit
+
+   !--------------------------------------------------------------------
+   !
+   ! set some local module variables based on input config choices
+   !
+   !--------------------------------------------------------------------
+
+   vertAdvOn = .false.
+
+   if (config_vert_grid_type .eq. 'zlevel') then 
+      vertAdvOn = .true.
+   endif
+
+   !--------------------------------------------------------------------
+
+   end subroutine OcnVelVertAdvInit
+
+!***********************************************************************
+
+end module OcnVelVertAdv
+
+!|||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||

Modified: branches/ocean_projects/performance/mpas/src/core_ocean/module_OcnVmixVelConst.F
===================================================================
--- branches/ocean_projects/performance/mpas/src/core_ocean/module_OcnVmixVelConst.F        2011-08-30 18:03:51 UTC (rev 963)
+++ branches/ocean_projects/performance/mpas/src/core_ocean/module_OcnVmixVelConst.F        2011-08-30 19:15:29 UTC (rev 964)
@@ -4,7 +4,7 @@
 !
 !&gt; \brief Ocean vertical mixing - constant parameterization 
 !&gt; \author Doug Jacobsen
-!&gt; \date   22 August 2011
+!&gt; \date   30 August 2011
 !&gt; \version SVN:$Id:$
 !&gt; \details
 !&gt;  This module contains routines for computing vertical mixing 
@@ -58,7 +58,7 @@
 !
 !&gt; \brief   Computes tendency term for constant vertical momentum mixing
 !&gt; \author  Doug Jacobsen
-!&gt; \date    22 August 2011
+!&gt; \date    30 August 2011
 !&gt; \version SVN:$Id$
 !&gt; \details 
 !&gt;  This routine computes the vertical mixing tendency for momentum

Modified: branches/ocean_projects/performance/mpas/src/core_ocean/module_OcnVmixVelTanh.F
===================================================================
--- branches/ocean_projects/performance/mpas/src/core_ocean/module_OcnVmixVelTanh.F        2011-08-30 18:03:51 UTC (rev 963)
+++ branches/ocean_projects/performance/mpas/src/core_ocean/module_OcnVmixVelTanh.F        2011-08-30 19:15:29 UTC (rev 964)
@@ -4,7 +4,7 @@
 !
 !&gt; \brief Ocean vertical mixing - Tanhant parameterization 
 !&gt; \author Doug Jacobsen
-!&gt; \date   22 August 2011
+!&gt; \date   30 August 2011
 !&gt; \version SVN:$Id:$
 !&gt; \details
 !&gt;  This module contains routines for computing vertical mixing 
@@ -58,7 +58,7 @@
 !
 !&gt; \brief   Computes tendency term for Tanhant vertical momentum mixing
 !&gt; \author  Doug Jacobsen
-!&gt; \date    22 August 2011
+!&gt; \date    30 August 2011
 !&gt; \version SVN:$Id$
 !&gt; \details 
 !&gt;  This routine computes the vertical mixing tendency for momentum

Modified: branches/ocean_projects/performance/mpas/src/core_ocean/module_mpas_core.F
===================================================================
--- branches/ocean_projects/performance/mpas/src/core_ocean/module_mpas_core.F        2011-08-30 18:03:51 UTC (rev 963)
+++ branches/ocean_projects/performance/mpas/src/core_ocean/module_mpas_core.F        2011-08-30 19:15:29 UTC (rev 964)
@@ -7,6 +7,7 @@
    use OcnHmixTracer
    use OcnVmixVel
    use OcnVmixTracer
+   use OcnVelVertAdv
 
    type (io_output_object) :: restart_obj
    integer :: restart_frame
@@ -85,6 +86,7 @@
       endif
 
       call OcnVmixVelInit
+      call OcnVelVertAdvInit
 
    ! mrp 100316 In order for this to work, we need to pass domain % dminfo as an 
    ! input arguement into mpas_init.  Ask about that later.  For now, there will be

Modified: branches/ocean_projects/performance/mpas/src/core_ocean/module_time_integration.F
===================================================================
--- branches/ocean_projects/performance/mpas/src/core_ocean/module_time_integration.F        2011-08-30 18:03:51 UTC (rev 963)
+++ branches/ocean_projects/performance/mpas/src/core_ocean/module_time_integration.F        2011-08-30 19:15:29 UTC (rev 964)
@@ -14,8 +14,10 @@
    use OcnVmixTracer
    use OcnVmixVel
    use OcnCoriolis
+   use OcnPressureGrad
    use OcnThickness
    use OcnVelocityForcing
+   use OcnVelVertAdv
 
    contains
 
@@ -371,70 +373,17 @@
       !
       ! velocity tendency: start accumulating tendency terms
       !
-      tend_u(:,:) = 0.0
+      tend_u = 0.0
 
       call timer_start(&quot;vel vert advect&quot;)
-      !
-      ! velocity tendency: vertical advection term -w du/dz
-      !
-      if (config_vert_grid_type.eq.'zlevel') then
-        allocate(w_dudzTopEdge(nVertLevels+1))
-        w_dudzTopEdge(1) = 0.0
-        do iEdge=1,grid % nEdgesSolve
-          cell1 = cellsOnEdge(1,iEdge)
-          cell2 = cellsOnEdge(2,iEdge)
 
-          do k=2,maxLevelEdgeTop(iEdge)
-            ! Average w from cell center to edge
-            wTopEdge = 0.5*(wTop(k,cell1)+wTop(k,cell2))
+      call OcnVelVertAdvTend(tend_u, u, wTop, grid)
 
-            ! compute dudz at vertical interface with first order derivative.
-            w_dudzTopEdge(k) = wTopEdge * (u(k-1,iEdge)-u(k,iEdge)) &amp;
-                         / (zMidZLevel(k-1) - zMidZLevel(k))
-          end do
-          w_dudzTopEdge(maxLevelEdgeTop(iEdge)+1) = 0.0
-
-          ! Average w*du/dz from vertical interface to vertical middle of cell
-          do k=1,maxLevelEdgeTop(iEdge)
-             tend_u(k,iEdge) = - 0.5 * (w_dudzTopEdge(k) + w_dudzTopEdge(k+1))
-          enddo
-        enddo
-        deallocate(w_dudzTopEdge)
-      endif
-
       call timer_stop(&quot;vel vert advect&quot;)
       call timer_start(&quot;press grad&quot;)
-      !
-      ! velocity tendency: pressure gradient
-      !
-      rho0Inv = 1.0/config_rho0
-      if (config_vert_grid_type.eq.'isopycnal') then
-        do iEdge=1,grid % nEdgesSolve
-          cell1 = cellsOnEdge(1,iEdge)
-          cell2 = cellsOnEdge(2,iEdge)
 
-          invDCEdge = 1.0 / dcEdge(iEdge)
+      call OcnPressureGradTend(tend_u, MontPot, pressure, grid)
 
-          do k=1,maxLevelEdgeTop(iEdge)
-             tend_u(k,iEdge) = tend_u(k,iEdge)     &amp;
-               - (MontPot(k,cell2) - MontPot(k,cell1))*invDCEdge
-           end do
-        enddo
-      elseif (config_vert_grid_type.eq.'zlevel') then
-        do iEdge=1,grid % nEdgesSolve
-          cell1 = cellsOnEdge(1,iEdge)
-          cell2 = cellsOnEdge(2,iEdge)
-
-          invDCEdge = 1.0 / dcEdge(iEdge)
-
-          do k=1,maxLevelEdgeTop(iEdge)
-             tend_u(k,iEdge) = tend_u(k,iEdge)     &amp;
-               - rho0Inv*(  pressure(k,cell2) &amp;
-                          - pressure(k,cell1) )*invDCEdge
-          end do
-        enddo
-      endif
-
       call timer_stop(&quot;press grad&quot;)
       call timer_start(&quot;vel dissipation&quot;)
 
@@ -733,6 +682,8 @@
             ! namelist, if desired.
             flux3Coef = 1.0
             do iCell=1,nCellsSolve 
+               ! DWJ : Efficency note - Could fuse the k=2 and k = maxLevelCell
+               ! loops to reduce total number of loops.
                k=2
                do iTracer=1,num_tracers
                  tracerTop(iTracer,k,iCell) = &amp;
@@ -764,6 +715,8 @@
             ! Compute tracerTop using 4rd order stencil [-1 7 7 -1]
 
             do iCell=1,nCellsSolve 
+               ! DWJ : Efficency note - Could fuse the k=2 and k = maxLevelCell
+               ! loops to reduce total number of loops.
                k=2
                   do iTracer=1,num_tracers
                     tracerTop(iTracer,k,iCell) = &amp;

</font>
</pre>