<p><b>dwj07@fsu.edu</b> 2012-09-19 15:47:35 -0600 (Wed, 19 Sep 2012)</p><p><br>
        -- BRANCH COMMIT --<br>
<br>
        Bit reproducibile version of ocean core.<br>
        This is for implementing OpenMP on element loops.<br>
</p><hr noshade><pre><font color="gray">Modified: branches/omp_blocks/openmp_test/src/core_ocean_elements/Registry
===================================================================
--- branches/omp_blocks/openmp_test/src/core_ocean_elements/Registry        2012-09-18 20:02:55 UTC (rev 2164)
+++ branches/omp_blocks/openmp_test/src/core_ocean_elements/Registry        2012-09-19 21:47:35 UTC (rev 2165)
@@ -209,9 +209,9 @@
 var persistent real    salinityRestore ( nCells ) 0 iro salinityRestore mesh - -
 
 % mrp trying to figure out why these do not appear
-var persistent real    windStressMonthly ( nMonths nEdges ) 0 iro windStressMonthly mesh - -
-var persistent real    temperatureRestoreMonthly ( nMonths nCells ) 0 iro temperatureRestoreMonthly mesh - -
-var persistent real    salinityRestoreMonthly ( nMonths nCells ) 0 iro salinityRestoreMonthly mesh - -
+var persistent real    windStressMonthly ( nMonths nEdges ) 0 - windStressMonthly mesh - -
+var persistent real    temperatureRestoreMonthly ( nMonths nCells ) 0 - temperatureRestoreMonthly mesh - -
+var persistent real    salinityRestoreMonthly ( nMonths nCells ) 0 - salinityRestoreMonthly mesh - -
 
 % Prognostic variables: read from input, saved in restart, and written to output
 var persistent real    u ( nVertLevels nEdges Time ) 2 ir u state - -
@@ -307,3 +307,8 @@
 var persistent real    acc_uReconstructMeridionalVar ( nVertLevels nCells Time ) 2 o acc_uReconstructMeridionalVar state - -
 var persistent real    acc_u ( nVertLevels nEdges Time ) 2 o acc_u state - -
 var persistent real    acc_uVar ( nVertLevels nEdges Time ) 2 o acc_uVar state - -
+
+% Sign fields, for openmp and bit reproducibility without branching statements.
+var persistent integer edgeSignOnCell ( maxEdges nCells ) 0 - edgeSignOnCell mesh - -
+var persistent integer edgeSignOnVertex ( maxEdges nVertices ) 0 - edgeSignOnVertex mesh - -
+var persistent integer kiteIndexOnCell ( maxEdges nCells ) 0 - kiteIndexOnCell mesh - -

Modified: branches/omp_blocks/openmp_test/src/core_ocean_elements/mpas_ocn_mpas_core.F
===================================================================
--- branches/omp_blocks/openmp_test/src/core_ocean_elements/mpas_ocn_mpas_core.F        2012-09-18 20:02:55 UTC (rev 2164)
+++ branches/omp_blocks/openmp_test/src/core_ocean_elements/mpas_ocn_mpas_core.F        2012-09-19 21:47:35 UTC (rev 2165)
@@ -249,6 +249,7 @@
       integer :: i, iEdge, iCell, k
       integer :: err1
    
+      call ocn_setup_sign_and_index_fields(mesh)
       call ocn_initialize_advection_rk(mesh, err)
       call mpas_ocn_tracer_advection_coefficients(mesh, err1)
       err = ior(err, err1)
@@ -960,6 +961,72 @@
 
    end subroutine ocn_compute_mesh_scaling!}}}
 
+   subroutine ocn_setup_sign_and_index_fields(mesh)!{{{
+
+       type (mesh_type), intent(inout) :: mesh
+
+       integer, dimension(:), pointer :: nEdgesOnCell
+       integer, dimension(:,:), pointer :: edgesOnCell, edgesOnVertex, cellsOnVertex, cellsOnEdge, verticesOnCell, verticesOnEdge
+       integer, dimension(:,:), pointer :: edgeSignOnCell, edgeSignOnVertex, kiteIndexOnCell
+
+       integer :: nCells, nEdges, nVertices, vertexDegree
+       integer :: iCell, iEdge, iVertex, i, j, k
+
+       nCells = mesh % nCells
+       nEdges = mesh % nEdges
+       nVertices = mesh % nVertices
+       vertexDegree = mesh % vertexDegree
+
+       nEdgesOnCell =&gt; mesh % nEdgesOnCell % array
+       edgesOnCell =&gt; mesh % edgeSOnCell % array
+       edgesOnVertex =&gt; mesh % edgesOnVertex % array
+       cellsOnVertex =&gt; mesh % cellsOnVertex % array
+       cellsOnEdge =&gt; mesh % cellsOnEdge % array
+       verticesOnCell =&gt; mesh % verticesOnCell % array
+       verticesOnEdge =&gt; mesh % verticesOnEdge % array
+       edgeSignOnCell =&gt; mesh % edgeSignOnCell % array
+       edgeSignOnVertex =&gt; mesh % edgeSignOnVertex % array
+       kiteIndexOnCell =&gt; mesh % kiteIndexOnCell % array
+
+       edgeSignOnCell = 0.0_RKIND
+       edgeSignOnVertex = 0.0_RKIND
+       kiteIndexOnCell = 0.0_RKIND
+
+       do iCell = 1, nCells
+         do i = 1, nEdgesOnCell(iCell) 
+           iEdge = edgesOnCell(i, iCell)
+           iVertex = verticesOnCell(i, iCell)
+
+           ! Vector points from cell 1 to cell 2
+           if(iCell == cellsOnEdge(1, iEdge)) then
+             edgeSignOnCell(i, iCell) = -1
+           else
+             edgeSignOnCell(i, iCell) =  1
+           end if
+
+           do j = 1, vertexDegree
+             if(cellsOnVertex(j, iVertex) == iCell) then
+               kiteIndexOnCell(i, iCell) = j
+             end if
+           end do
+         end do
+       end do
+
+       do iVertex = 1, nVertices
+         do i = 1, vertexDegree
+           iEdge = edgesOnVertex(i, iVertex)
+
+           ! Vector points from vertex 1 to vertex 2
+           if(iVertex == verticesOnEdge(1, iVertex)) then
+             edgeSignOnVertex(i, iVertex) = -1
+           else
+             edgeSignOnVertex(i, iVertex) =  1
+           end if
+         end do
+       end do
+
+   end subroutine ocn_setup_sign_and_index_fields!}}}
+
 end module mpas_core
 
 ! vim: foldmethod=marker

Modified: branches/omp_blocks/openmp_test/src/core_ocean_elements/mpas_ocn_time_integration_split.F
===================================================================
--- branches/omp_blocks/openmp_test/src/core_ocean_elements/mpas_ocn_time_integration_split.F        2012-09-18 20:02:55 UTC (rev 2164)
+++ branches/omp_blocks/openmp_test/src/core_ocean_elements/mpas_ocn_time_integration_split.F        2012-09-19 21:47:35 UTC (rev 2165)
@@ -452,6 +452,8 @@
       
                 block =&gt; domain % blocklist
                 do while (associated(block))
+                   allocate(uTemp(block % mesh % nEdges+1))
+                   uTemp(:) = block % state % time_levs(newBtrSubcycleTime) % state % uBtrSubcycle % array(:)
                    do iEdge=1,block % mesh % nEdges 
                      cell1 = block % mesh % cellsOnEdge % array(1,iEdge)
                      cell2 = block % mesh % cellsOnEdge % array(2,iEdge)
@@ -460,8 +462,8 @@
                      CoriolisTerm = 0.0
                      do i = 1,block % mesh % nEdgesOnEdge % array(iEdge)
                          eoe = block % mesh % edgesOnEdge % array(i,iEdge)
-                       CoriolisTerm = CoriolisTerm + block % mesh % weightsOnEdge % array(i,iEdge) &amp;
-                             * block % state % time_levs(uPerpTime) % state % uBtrSubcycle % array(eoe) &amp;
+                         CoriolisTerm = CoriolisTerm + block % mesh % weightsOnEdge % array(i,iEdge) &amp;
+                             * uTemp(eoe) &amp;
                              * block % mesh % fEdge  % array(eoe) 
                      end do
       
@@ -478,6 +480,7 @@
                          + dt/config_n_btr_subcycles *(CoriolisTerm - gravity *(sshCell2 - sshCell1) /block % mesh % dcEdge % array(iEdge) &amp;
                          + block % state % time_levs(1) % state % GBtrForcing % array(iEdge))) * block % mesh % edgeMask % array(1,iEdge)
                    end do
+                   deallocate(uTemp)
       
                    block =&gt; block % next
                 end do  ! block

</font>
</pre>