<p><b>dwj07@fsu.edu</b> 2011-08-30 14:18:32 -0600 (Tue, 30 Aug 2011)</p><p><br>
        Various cleanup of modules.<br>
</p><hr noshade><pre><font color="gray">Modified: branches/ocean_projects/performance/mpas/src/core_ocean/module_OcnCoriolis.F
===================================================================
--- branches/ocean_projects/performance/mpas/src/core_ocean/module_OcnCoriolis.F        2011-08-30 19:23:25 UTC (rev 965)
+++ branches/ocean_projects/performance/mpas/src/core_ocean/module_OcnCoriolis.F        2011-08-30 20:18:32 UTC (rev 966)
@@ -95,7 +95,7 @@
       !-----------------------------------------------------------------
 
       integer :: iEdge, j, k, cell1, cell2, eoe, nEdgesSolve
-      real (kind=RKIND) :: q, workpv
+      real (kind=RKIND) :: q, workpv, invDCEdge
       integer, dimension(:,:), pointer :: cellsOnEdge, edgesOnEdge
       integer, dimension(:), pointer :: maxLeveLEdgeTop, nEdgesOnEdge
 
@@ -120,6 +120,8 @@
          cell1 = cellsOnEdge(1,iEdge)
          cell2 = cellsOnEdge(2,iEdge)
 
+         invDCEdge = 1.d0 / dcEdge(iEdge)
+
          do k=1,maxLevelEdgeTop(iEdge)
 
             q = 0.0
@@ -130,7 +132,7 @@
             end do
             tend(k,iEdge) = tend(k,iEdge)     &amp;
                    + q     &amp;
-                   - (   ke(k,cell2) - ke(k,cell1) ) / dcEdge(iEdge)
+                   - (   ke(k,cell2) - ke(k,cell1) ) * invDCEdge
 
          end do
       end do

Modified: branches/ocean_projects/performance/mpas/src/core_ocean/module_OcnHmixTracerDel2.F
===================================================================
--- branches/ocean_projects/performance/mpas/src/core_ocean/module_OcnHmixTracerDel2.F        2011-08-30 19:23:25 UTC (rev 965)
+++ branches/ocean_projects/performance/mpas/src/core_ocean/module_OcnHmixTracerDel2.F        2011-08-30 20:18:32 UTC (rev 966)
@@ -108,7 +108,7 @@
 
       real (kind=RKIND) :: &amp;
          invAreaCell1, invAreaCell2, &amp;! reciprocal of cell areas
-         flux                         ! flux across edge
+         invDCEdge, flux              ! reciprocal of dcedge and flux across edge
 
       real (kind=RKIND), dimension(:), pointer :: &amp;
          dvEdge, dcEdge, areaCell
@@ -145,8 +145,6 @@
       maxLevelEdgeTop   =&gt; grid % maxLevelEdgeTop % array
 
       nEdges      = grid % nEdges
-      nCells      = grid % nCells
-      nVertLevels = size(tracers,dim=2)
       nTracers    = size(tracers,dim=1)
 
       !-----------------------------------------------------------------
@@ -164,6 +162,7 @@
          cell2 = grid % cellsOnEdge % array(2,iEdge)
          invAreaCell1 = 1.d0/areaCell(cell1)
          invAreaCell2 = 1.d0/areaCell(cell2)
+         invDCEdge = 1.d0/dcEdge(iEdge)
 
          !*** compute \kappa_2 </font>
<font color="gray">abla \phi on edge
 
@@ -171,11 +170,7 @@
          !*** compute a boundary mask at each vertical level
          !*** on this edge
 
-            if (boundaryEdge(k,iEdge) == 1) then
-               boundaryMask = 0.d0
-            else
-               boundaryMask = 1.d0
-            endif
+            boundaryMask = boundaryEdge(k,iEdge) * 1.d0
 
             do iTracer=1,nTracers
 
@@ -183,13 +178,13 @@
 
                 flux = dvEdge(iEdge) * hEdge(k,iEdge) * eddyDiff2 * &amp;
                        ( tracers(iTracer,k,cell2)                   &amp;
-                       - tracers(iTracer,k,cell1))/dcEdge(iEdge)*   &amp;
+                       - tracers(iTracer,k,cell1))*invDCEdge*   &amp;
                          boundaryMask
 
                 tend(iTracer,k,cell1) = &amp; 
-                tend(iTracer,k,cell1) + flux * invAreaCell1
+                  tend(iTracer,k,cell1) + flux * invAreaCell1
                 tend(iTracer,k,cell2) = &amp;
-                tend(iTracer,k,cell2) - flux * invAreaCell2
+                  tend(iTracer,k,cell2) - flux * invAreaCell2
 
              end do
          end do
@@ -224,11 +219,8 @@
    hmixDel2On = .false.
 
    if (config_OcnHmixTracerDel2Diff &gt; 0.0) then
-
       hmixDel2On = .true.
-
       eddyDiff2 = config_OcnHmixTracerDel2Diff
-
    endif
 
    !--------------------------------------------------------------------

Modified: branches/ocean_projects/performance/mpas/src/core_ocean/module_OcnHmixTracerDel4.F
===================================================================
--- branches/ocean_projects/performance/mpas/src/core_ocean/module_OcnHmixTracerDel4.F        2011-08-30 19:23:25 UTC (rev 965)
+++ branches/ocean_projects/performance/mpas/src/core_ocean/module_OcnHmixTracerDel4.F        2011-08-30 20:18:32 UTC (rev 966)
@@ -104,23 +104,19 @@
       !-----------------------------------------------------------------
 
       integer :: &amp;
-         i, k, iCell, iEdge, iTracer, &amp;! loop counters
-         cell1, cell2,                &amp;! cell indices
-         nEdges, nCells, nVertLevels, nTracers ! array sizes
+         i, k, iEdge, iTracer, &amp;! loop counters
+         cell1, cell2,         &amp;! cell indices
+         nTracers, nVertLevels
 
       real (kind=RKIND) :: &amp;
          invAreaCell1, invAreaCell2, &amp;! reciprocal of cell areas
-         flux                         ! flux across edge
+         invDCEdge, flux              ! reciprocal of dcedge and flux across edge
 
       real (kind=RKIND), dimension(:), pointer :: &amp;
          dvEdge, dcEdge, areaCell
 
-      integer, dimension(:,:), pointer :: boundaryEdge
+      integer, dimension(:), pointer :: maxLevelEdgeTop
 
-      integer, dimension(:), pointer :: maxLevelCell, maxLevelEdgeTop
-
-      integer, dimension(:,:), pointer :: cellsOnEdge
-
       real (kind=RKIND) :: boundaryMask
 
       real (kind=RKIND), dimension(:,:,:), allocatable:: del2tracer
@@ -140,17 +136,12 @@
       !-----------------------------------------------------------------
 
       areaCell          =&gt; grid % areaCell % array
-      cellsOnEdge       =&gt; grid % cellsOnEdge % array
       dvEdge            =&gt; grid % dvEdge % array
       dcEdge            =&gt; grid % dcEdge % array
-      boundaryEdge      =&gt; grid % boundaryEdge % array
-      maxLevelCell      =&gt; grid % maxLevelCell % array
       maxLevelEdgeTop   =&gt; grid % maxLevelEdgeTop % array
 
-      nEdges      = grid % nEdges
-      nCells      = grid % nCells
-      nVertLevels = size(tracers,dim=2)
-      nTracers    = size(tracers,dim=1)
+      nTracers = size(tracers, dim=1)
+      nVertLEvels = size(tracers, dim=2)
 
       !-----------------------------------------------------------------
       !
@@ -160,86 +151,72 @@
       !
       !-----------------------------------------------------------------
 
-      allocate(del2tracer(nTracers,nVertLevels,nCells+1))
+      allocate(del2tracer(nTracers,nVertLevels,grid % nCells+1))
 
       del2tracer(:,:,:) = 0.d0
 
       !*** first application of Laplacian
       !*** del2: div(h </font>
<font color="red">abla \phi) at cell center
 
-      do iEdge=1,nEdges
+      do iEdge=1,grid % nEdges
 
          !*** get cell pair that share edge
 
          cell1 = grid % cellsOnEdge % array(1,iEdge)
          cell2 = grid % cellsOnEdge % array(2,iEdge)
 
+         invAreaCell1 = 1.d0/areaCell(cell1)
+         invAreaCell2 = 1.d0/areaCell(cell2)
 
+         invDCEdge = 1.d0/dcEdge(iEdge)
 
          do k=1,maxLevelEdgeTop(iEdge)
          !*** compute a boundary mask at each vertical level
          !*** on this edge
 
-            if (boundaryEdge(k,iEdge) == 1) then
-               boundaryMask = 0.d0
-            else
-               boundaryMask = 1.d0
-            endif
+            boundaryMask = grid % boundaryEdge % array(k,iEdge) * 1.d0
+
             do iTracer=1,nTracers
 
                 flux = dvEdge(iEdge)*hEdge(k,iEdge) *            &amp;
                        (tracers(iTracer,k,cell2) -               &amp;
-                        tracers(iTracer,k,cell1))/dcEdge(iEdge)* &amp;
+                        tracers(iTracer,k,cell1))*invDCEdge* &amp;
                        boundaryMask
 
                 del2tracer(iTracer,k,cell1) = &amp;
-                del2tracer(iTracer,k,cell1) + flux
+                del2tracer(iTracer,k,cell1) + flux * invAreaCell1
     
                 del2tracer(iTracer,k,cell2) = &amp;
-                del2tracer(iTracer,k,cell2) - flux
+                del2tracer(iTracer,k,cell2) - flux * invAreaCell2
 
              end do
          end do
       end do
 
-      do iCell = 1,nCells
-         invAreaCell1 = 1.d0/areaCell(iCell)
-         do k=1,maxLevelCell(iCell)
-         do iTracer=1,nTracers
-            del2tracer(iTracer,k,iCell) = &amp;
-            del2tracer(iTracer,k,iCell) * invAreaCell1
-         end do
-         end do
-      end do
-
       !*** second application of Laplacian
       !*** del2: div(h </font>
<font color="gray">abla [delsq_tracer]) at cell center
 
-      do iEdge=1,nEdges
+      do iEdge=1,grid % nEdges
 
          !*** get cell pair that share edge
 
          cell1 = grid % cellsOnEdge % array(1,iEdge)
          cell2 = grid % cellsOnEdge % array(2,iEdge)
-         invAreaCell1 = 1.0 / areaCell(cell1)
-         invAreaCell2 = 1.0 / areaCell(cell2)
+         invAreaCell1 = 1.d0 / areaCell(cell1)
+         invAreaCell2 = 1.d0 / areaCell(cell2)
+         invDCEdge = 1.d0/dcEdge(iEdge)
 
-
-
          do k=1,maxLevelEdgeTop(iEdge)
          !*** compute a boundary mask at each vertical level
          !*** on this edge
 
-            if (boundaryEdge(k,iEdge) == 1) then
-               boundaryMask = 0.d0
-            else
-               boundaryMask = 1.d0
-            endif
+            boundaryMask = grid % boundaryEdge % array(k,iEdge) * 1.d0
+
             do iTracer=1,nTracers
 
                flux = dvEdge(iEdge)* eddyDiff4 *                    &amp;
                       (del2tracer(iTracer,k,cell2) -                &amp;
-                       del2tracer(iTracer,k,cell1))/dcEdge(iEdge) * &amp;
+                       del2tracer(iTracer,k,cell1))*invDCEdge * &amp;
                        boundaryMask
 
                tend(iTracer,k,cell1) = &amp;

Modified: branches/ocean_projects/performance/mpas/src/core_ocean/module_OcnHmixVelDel2.F
===================================================================
--- branches/ocean_projects/performance/mpas/src/core_ocean/module_OcnHmixVelDel2.F        2011-08-30 19:23:25 UTC (rev 965)
+++ branches/ocean_projects/performance/mpas/src/core_ocean/module_OcnHmixVelDel2.F        2011-08-30 20:18:32 UTC (rev 966)
@@ -105,9 +105,10 @@
       integer ::               &amp;
          iEdge, k,             &amp;! loop counters
          cell1, cell2,         &amp;! cell indices
-         vertex1, vertex2,     &amp;! vertex indices
-         nEdgesSolve            ! array sizes
+         vertex1, vertex2       ! vertex indices
 
+      real (kind=RKIND) :: invDCEdge, invDVEdge
+
       real (kind=RKIND), dimension(:), pointer :: &amp;
          dvEdge, dcEdge
 
@@ -137,8 +138,6 @@
       dcEdge            =&gt; grid % dcEdge % array
       maxLevelEdgeTop   =&gt; grid % maxLevelEdgeTop % array
 
-      nEdgesSolve       =  grid % nEdgesSolve
-
       !-----------------------------------------------------------------
       !
       ! del2 horizontal momentum diffusion
@@ -150,20 +149,23 @@
       !
       !-----------------------------------------------------------------
 
-      do iEdge=1,nEdgesSolve
+      do iEdge=1,grid % nEdgesSolve
 
          cell1   = cellsOnEdge(1,iEdge)
          cell2   = cellsOnEdge(2,iEdge)
          vertex1 = verticesOnEdge(1,iEdge)
          vertex2 = verticesOnEdge(2,iEdge)
 
+         invDCEdge = 1.d0 / dcEdge(iEdge)
+         invDVEdge = 1.d0 / dvEdge(iEdge)
+
          do k=1,maxLevelEdgeTop(iEdge)
 
             tend(k,iEdge) = tend(k,iEdge) + eddyVisc2 *(               &amp;
                             ( divergence(k,cell2  ) -                  &amp;
-                              divergence(k,cell1  ) ) / dcEdge(iEdge)  &amp;
+                              divergence(k,cell1  ) ) * invDCEdge      &amp;
                            -( vorticity (k,vertex2) -                  &amp;
-                              vorticity (k,vertex1) ) / dvEdge(iEdge))
+                              vorticity (k,vertex1) ) * invDVEdge)
 
          end do
       end do

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 19:23:25 UTC (rev 965)
+++ branches/ocean_projects/performance/mpas/src/core_ocean/module_time_integration.F        2011-08-30 20:18:32 UTC (rev 966)
@@ -278,16 +278,13 @@
       type (mesh_type), intent(in) :: grid
 
       real (kind=RKIND), dimension(:,:), pointer :: &amp;
-        h_edge, h, u, v, pressure, tend_h, tend_u, circulation,      &amp;
-        vorticity, ke, ke_edge, pv_edge, MontPot, wTop, divergence
+        h_edge, u, pressure, tend_h, tend_u, vorticity,  &amp;
+        ke, ke_edge, pv_edge, MontPot, wTop, divergence
       type (dm_info) :: dminfo
 
-      h           =&gt; s % h % array
       u           =&gt; s % u % array
-      v           =&gt; s % v % array
       wTop        =&gt; s % wTop % array
       h_edge      =&gt; s % h_edge % array
-      circulation =&gt; s % circulation % array
       vorticity   =&gt; s % vorticity % array
       divergence  =&gt; s % divergence % array
       ke          =&gt; s % ke % array
@@ -307,7 +304,6 @@
 
       call OcnThicknessTend(tend_h, h_edge, u, wTop, grid)
 
-
       call timer_stop(&quot;height&quot;)
       call timer_start(&quot;velocity&quot;)
       !

</font>
</pre>