<p><b>dwj07@fsu.edu</b> 2011-10-27 09:02:51 -0600 (Thu, 27 Oct 2011)</p><p><br>
        -- BRANCH COMMIT --<br>
<br>
        Adding a routine to compute edgeMask, vertexMask, and cellMask.<br>
        These are based on boundaryEdge, boundaryCell, and boundaryVertex.<br>
<br>
        Will be used to remove if statements.<br>
</p><hr noshade><pre><font color="gray">Modified: branches/ocean_projects/performance/src/core_ocean/Registry
===================================================================
--- branches/ocean_projects/performance/src/core_ocean/Registry        2011-10-26 22:13:47 UTC (rev 1148)
+++ branches/ocean_projects/performance/src/core_ocean/Registry        2011-10-27 15:02:51 UTC (rev 1149)
@@ -180,6 +180,9 @@
 var persistent integer boundaryEdge ( nVertLevels nEdges ) 0 iro boundaryEdge mesh - -
 var persistent integer boundaryVertex ( nVertLevels nVertices ) 0 iro boundaryVertex mesh - -
 var persistent integer boundaryCell ( nVertLevels nCells ) 0 iro boundaryCell mesh - -
+var persistent integer edgeMask ( nVertLevels nEdges ) 0 o edgeMask mesh - -
+var persistent integer vertexMask ( nVertLevels nVertices ) 0 o vertexMask mesh - -
+var persistent integer cellMask ( nVertLevels nCells ) 0 o cellMask mesh - -
 var persistent real    u_src ( nVertLevels nEdges ) 0 ir u_src mesh - -
 var persistent real    temperatureRestore ( nCells ) 0 ir temperatureRestore mesh - -
 var persistent real    salinityRestore ( nCells ) 0 ir salinityRestore mesh - -

Modified: branches/ocean_projects/performance/src/core_ocean/mpas_ocn_mpas_core.F
===================================================================
--- branches/ocean_projects/performance/src/core_ocean/mpas_ocn_mpas_core.F        2011-10-26 22:13:47 UTC (rev 1148)
+++ branches/ocean_projects/performance/src/core_ocean/mpas_ocn_mpas_core.F        2011-10-27 15:02:51 UTC (rev 1149)
@@ -90,13 +90,13 @@
 
       if (.not. config_do_restart) call setup_sw_test_case(domain)
 
-      call compute_maxLevel(domain)
+      call ocn_compute_max_level(domain)
 
       if (config_vert_grid_type.eq.'isopycnal') then
          print *, ' Using isopycnal coordinates'
       elseif (config_vert_grid_type.eq.'zlevel') then
          print *, ' Using z-level coordinates'
-         call init_ZLevel(domain)
+         call ocn_init_z_level(domain)
       else 
          print *, ' Incorrect choice of config_vert_grid_type:',&amp;
            config_vert_grid_type
@@ -115,7 +115,7 @@
       !
       dt = config_dt
 
-      call simulation_clock_init(domain, dt, startTimeStamp)
+      call ocn_simulation_clock_init(domain, dt, startTimeStamp)
 
       block =&gt; domain % blocklist
       do while (associated(block))
@@ -133,7 +133,7 @@
           call ocn_compute_global_diagnostics(domain, 1 , 0, dt)
           call mpas_timer_stop(&quot;global diagnostics&quot;, globalDiagTimer)
 !         call mpas_output_state_init(output_obj, domain, &quot;OUTPUT&quot;)
-!         call write_output_frame(output_obj, output_frame, domain)
+!         call ocn_write_output_frame(output_obj, output_frame, domain)
       endif
 
       restart_frame = 1
@@ -141,7 +141,7 @@
 
    end subroutine mpas_core_init!}}}
 
-   subroutine simulation_clock_init(domain, dt, startTimeStamp)!{{{
+   subroutine ocn_simulation_clock_init(domain, dt, startTimeStamp)!{{{
 
       implicit none
 
@@ -197,7 +197,7 @@
 
       call mpas_get_time(curr_time=startTime, dateTimeString=startTimeStamp, ierr=ierr)
 
-   end subroutine simulation_clock_init!}}}
+   end subroutine ocn_simulation_clock_init!}}}
 
    subroutine mpas_init_block(block, mesh, dt)!{{{
    
@@ -217,7 +217,8 @@
       call ocn_diagnostic_solve(dt, block % state % time_levs(1) % state, mesh)
       call mpas_timer_stop(&quot;diagnostic solve&quot;, initDiagSolveTimer)
 
-      call compute_mesh_scaling(mesh)
+      call ocn_compute_mesh_scaling(mesh)
+      call build_boundary_masks(mesh) !dwj: need to fix this routine
  
       call mpas_rbf_interp_initialize(mesh)
       call mpas_init_reconstruct(mesh)
@@ -305,7 +306,7 @@
       call mpas_get_time(curr_time=currTime, dateTimeString=timeStamp, ierr=ierr)
       write(0,*) 'Initial time ', timeStamp
 
-      call write_output_frame(output_obj, output_frame, domain)
+      call ocn_write_output_frame(output_obj, output_frame, domain)
 
       ! During integration, time level 1 stores the model state at the beginning of the
       !   time step, and time level 2 stores the state advanced dt in time by timestep(...)
@@ -329,7 +330,7 @@
          if (mpas_is_alarm_ringing(clock, outputAlarmID, ierr=ierr)) then
             call mpas_reset_clock_alarm(clock, outputAlarmID, ierr=ierr)
             if(output_frame == 1) call mpas_output_state_init(output_obj, domain, &quot;OUTPUT&quot;, trim(timeStamp)) ! output_frame will always be &gt; 1 here unless it is reset after the output file is finalized
-            call write_output_frame(output_obj, output_frame, domain)
+            call ocn_write_output_frame(output_obj, output_frame, domain)
          end if
 
          if (mpas_is_alarm_ringing(clock, restartAlarmID, ierr=ierr)) then
@@ -343,7 +344,7 @@
 
    end subroutine mpas_core_run!}}}
    
-   subroutine write_output_frame(output_obj, output_frame, domain)!{{{
+   subroutine ocn_write_output_frame(output_obj, output_frame, domain)!{{{
    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
    ! Compute diagnostic fields for a domain and write model state to output file
    !
@@ -366,7 +367,7 @@
    
       block_ptr =&gt; domain % blocklist
       do while (associated(block_ptr))
-         call compute_output_diagnostics(block_ptr % state % time_levs(1) % state, block_ptr % mesh)
+         call ocn_compute_output_diagnostics(block_ptr % state % time_levs(1) % state, block_ptr % mesh)
          block_ptr =&gt; block_ptr % next
       end do
    
@@ -383,9 +384,9 @@
          end if
       end if
    
-   end subroutine write_output_frame!}}}
+   end subroutine ocn_write_output_frame!}}}
    
-   subroutine compute_output_diagnostics(state, grid)!{{{
+   subroutine ocn_compute_output_diagnostics(state, grid)!{{{
    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
    ! Compute diagnostic fields for a domain
    !
@@ -405,7 +406,7 @@
       integer :: i, eoe
       integer :: iEdge, k
    
-   end subroutine compute_output_diagnostics!}}}
+   end subroutine ocn_compute_output_diagnostics!}}}
    
    subroutine mpas_timestep(domain, itimestep, dt, timeStamp)!{{{
    
@@ -450,7 +451,7 @@
 
    end subroutine mpas_timestep!}}}
 
-subroutine init_ZLevel(domain)!{{{
+subroutine ocn_init_z_level(domain)!{{{
 ! Initialize maxLevel and bouncary grid variables.
 
    use mpas_grid_types
@@ -612,9 +613,9 @@
 
    end do
 
-end subroutine init_ZLevel!}}}
+end subroutine ocn_init_z_level!}}}
 
-subroutine compute_maxLevel(domain)!{{{
+subroutine ocn_compute_max_level(domain)!{{{
 ! Initialize maxLevel and bouncary grid variables.
 
    use mpas_grid_types
@@ -735,7 +736,7 @@
    ! Note: We do not update halos on maxLevel* variables.  I want the
    ! outside edge of a halo to be zero on each processor.
 
-end subroutine compute_maxLevel!}}}
+end subroutine ocn_compute_max_level!}}}
    
    subroutine mpas_core_finalize(domain)!{{{
    
@@ -753,7 +754,7 @@
 
    end subroutine mpas_core_finalize!}}}
 
-   subroutine compute_mesh_scaling(mesh)!{{{
+   subroutine ocn_compute_mesh_scaling(mesh)!{{{
 
       use mpas_grid_types
       use mpas_configure
@@ -783,8 +784,52 @@
          end do
       end if
 
-   end subroutine compute_mesh_scaling!}}}
+   end subroutine ocn_compute_mesh_scaling!}}}
 
+   subroutine build_boundary_masks(grid)!{{{
+
+      use mpas_grid_types
+
+      type (mesh_type), intent(inout) :: grid
+
+      integer :: i, k
+      integer :: nCells, nEdges, nVertices, nVertLevels
+      integer, dimension(:,:), pointer :: edgeMask, cellMask, vertexMask
+      integer, dimension(:,:), pointer :: boundaryEdge, boundaryCell, boundaryVertex
+
+      nCells = grid % nCells
+      nEdges = grid % nEdges
+      nVertices = grid % nVertices
+      nVertLevels = grid % nVertLevels
+
+      boundaryCell =&gt; grid % boundaryCell % array
+      boundaryEdge =&gt; grid % boundaryEdge % array
+      boundaryVertex =&gt; grid % boundaryVertex % array
+
+      edgeMask =&gt; grid % edgeMask % array
+      cellMask =&gt; grid % cellMask % array
+      vertexMask =&gt; grid % vertexMask % array
+
+      do i = 1, nCells
+        do k = 1, nVertLevels 
+           cellMask(k, i) = transfer(abs(.not.boundaryCell(k, i)), cellMask(k, i))
+        end do
+      end do
+
+      do i = 1, nEdges
+        do k = 1, nVertLevels
+           edgeMask(k, i) = transfer(abs(.not.boundaryEdge(k, i)), edgeMask(k, i))
+        end do
+      end do
+
+      do i = 1, nVertices
+        do k = 1, nVertLevels
+           vertexMask(k, i) = transfer(abs(.not.boundaryVertex(k, i)), vertexMask(k, i))
+        end do
+      end do
+
+   end subroutine build_boundary_masks!}}}
+
 end module mpas_core
 
 ! vim: foldmethod=marker

</font>
</pre>