<p><b>dwj07@fsu.edu</b> 2012-07-23 12:06:01 -0600 (Mon, 23 Jul 2012)</p><p><br>
        -- BRANCH COMMIT --<br>
<br>
        Fixing branch to work with 0 blocks on processors.<br>
<br>
        Only tested with ocean core, but should work with other cores as well (assuming all block loops are written properly).<br>
        Should be able to test easily with config_number_of_blocks = 2 and mpirun -n 4.<br>
</p><hr noshade><pre><font color="gray">Modified: branches/omp_blocks/multiple_blocks/src/framework/mpas_block_creator.F
===================================================================
--- branches/omp_blocks/multiple_blocks/src/framework/mpas_block_creator.F        2012-07-23 17:39:48 UTC (rev 2041)
+++ branches/omp_blocks/multiple_blocks/src/framework/mpas_block_creator.F        2012-07-23 18:06:01 UTC (rev 2042)
@@ -60,61 +60,56 @@
      nBlocks = size(blockID)
      nHalos = config_num_halos
 
-     if(nBlocks &gt; 0) then
-       ! Setup first block
-       allocate(domain % blocklist)
-!      block =&gt; domain % blocklist
-       nullify(domain % blocklist % prev)
-       nullify(domain % blocklist % next)
-   
-       ! Setup first block field
-       allocate(indexToCellID)
-       nullify(indexToCellID % next)
+     ! Setup first block
+     allocate(domain % blocklist)
+     nullify(domain % blocklist % prev)
+     nullify(domain % blocklist % next)
   
-       ! Loop over blocks
-       blockCursor =&gt; domain % blocklist
-       fieldCursor =&gt; indexToCellID
-       do i = 1, nBlocks
-         ! Initialize block information
-         blockCursor % blockID = blockID(i)
-         blockCursor % localBlockID = i - 1
-         blockCursor % domain =&gt; domain
-   
-         ! Link to block, and setup array size
-         fieldCursor % block =&gt; blockCursor
-         fieldCursor % dimSizes(1) = blockCount(i)
-         nullify(fieldCursor % ioinfo)
+     ! Setup first block field
+     allocate(indexToCellID)
+     nullify(indexToCellID % next)

+     ! Loop over blocks
+     blockCursor =&gt; domain % blocklist
+     fieldCursor =&gt; indexToCellID
+     write(6,*) 'nBlocks = ', nBlocks
+     do i = 1, nBlocks
+       ! Initialize block information
+       blockCursor % blockID = blockID(i)
+       blockCursor % localBlockID = i - 1
+       blockCursor % domain =&gt; domain
   
-         ! Initialize exchange lists
-         call mpas_dmpar_init_mulithalo_exchange_list(fieldCursor % sendList, nHalos)
-         call mpas_dmpar_init_mulithalo_exchange_list(fieldCursor % recvList, nHalos)
-         call mpas_dmpar_init_mulithalo_exchange_list(fieldCursor % copyList, nHalos)
+       ! Link to block, and setup array size
+       fieldCursor % block =&gt; blockCursor
+       fieldCursor % dimSizes(1) = blockCount(i)
+       nullify(fieldCursor % ioinfo)

+       ! Initialize exchange lists
+       write(6,*) 'SETTING UP EXCH LISTS ON BLOCK', i, blockID(i)
+       call mpas_dmpar_init_mulithalo_exchange_list(fieldCursor % sendList, nHalos)
+       call mpas_dmpar_init_mulithalo_exchange_list(fieldCursor % recvList, nHalos)
+       call mpas_dmpar_init_mulithalo_exchange_list(fieldCursor % copyList, nHalos)

+       ! Allocate array, and copy indices into array
+       allocate(fieldCursor % array(fieldCursor % dimSizes(1)))
+       fieldCursor % array(:) = cellList(blockStart(i)+1:blockStart(i)+blockCount(i))
+       call mpas_quicksort(fieldCursor % dimSizes(1), fieldCursor % array)
   
-         ! Allocate array, and copy indices into array
-         allocate(fieldCursor % array(fieldCursor % dimSizes(1)))
-         fieldCursor % array(:) = cellList(blockStart(i)+1:blockStart(i)+blockCount(i))
-         call mpas_quicksort(fieldCursor % dimSizes(1), fieldCursor % array)
-   
-         ! Advance cursors, and create new blocks as needed
-         if(i &lt; nBlocks) then
-           allocate(blockCursor % next)
-           allocate(fieldCursor % next)
+       ! Advance cursors, and create new blocks as needed
+       if(i &lt; nBlocks) then
+         allocate(blockCursor % next)
+         allocate(fieldCursor % next)

+         blockCursor % next % prev =&gt; blockCursor
   
-           blockCursor % next % prev =&gt; blockCursor
-   
-           blockCursor =&gt; blockCursor % next
-           fieldCursor =&gt; fieldCursor % next
-         end if
-  
-         ! Nullify next pointers
-         nullify(blockCursor % next)
-         nullify(fieldCursor % next)
-       end do
-     else
-       ! If no blocks, nullify block and field pointers
-       nullify(domain % blocklist)
-       nullify(indexToCellID)
-     end if
+         blockCursor =&gt; blockCursor % next
+         fieldCursor =&gt; fieldCursor % next
+       end if

+       ! Nullify next pointers
+       nullify(blockCursor % next)
+       nullify(fieldCursor % next)
+     end do
    end subroutine mpas_block_creator_setup_blocks_and_0halo_cells!}}}
 
 !***********************************************************************
@@ -167,26 +162,18 @@
      call mpas_dmpar_get_exch_list(1, indexToCellIDBlock, indexToCellID_0Halo)
 
      ! Setup header fields if at least 1 block exists
-     if(associated(indexToCellID_0Halo)) then
-       allocate(nEdgesOnCell_0Halo)
-       nullify(nEdgesOncell_0Halo % next)
+     allocate(nEdgesOnCell_0Halo)
+     nullify(nEdgesOncell_0Halo % next)
 
-       allocate(cellsOnCell_0Halo)
-       nullify(cellsOnCell_0Halo % next)
+     allocate(cellsOnCell_0Halo)
+     nullify(cellsOnCell_0Halo % next)
   
-       allocate(verticesOnCell_0Halo)
-       nullify(verticesOnCell_0Halo % next)
+     allocate(verticesOnCell_0Halo)
+     nullify(verticesOnCell_0Halo % next)
   
-       allocate(edgesOnCell_0Halo)
-       nullify(edgesOnCell_0Halo % next)
-     else
-       nullify(nEdgesOnCell_0Halo)
-       nullify(cellsOnCell_0Halo)
-       nullify(verticesOnCell_0Halo)
-       nullify(edgesOnCell_0Halo)
-     end if
+     allocate(edgesOnCell_0Halo)
+     nullify(edgesOnCell_0Halo % next)
 
-
      ! Loop over blocks
      indexCursor =&gt; indexToCellID_0Halo
      nEdgesCursor =&gt; nEdgesOnCell_0Halo
@@ -317,16 +304,11 @@
      nHalos = config_num_halos
 
      ! Setup initial block for each field
-     if(associated(indexToCellID_0Halo)) then
-       allocate(cellsOnEdge_0Halo)
-       allocate(indexToEdgeID_0Halo)
+     allocate(cellsOnEdge_0Halo)
+     allocate(indexToEdgeID_0Halo)
 
-       nullify(cellsOnEdge_0Halo % next)
-       nullify(indexToEdgeID_0Halo % next)
-     else
-       nullify(cellsOnEdge_0Halo)
-       nullify(indexToEdgeID_0Halo)
-     end if
+     nullify(cellsOnEdge_0Halo % next)
+     nullify(indexToEdgeID_0Halo % next)
 
      ! Loop over blocks
      indexToCellCursor =&gt; indexToCellID_0Halo
@@ -538,19 +520,13 @@
      allocate(sendingHaloLayers(1))
 
      ! Setup header fields
-     if(associated(indexToCellID)) then
-       allocate(nCellsSolve)
-       allocate(cellLimitField)
-       allocate(offSetField)
+     allocate(nCellsSolve)
+     allocate(cellLimitField)
+     allocate(offSetField)
 
-       nullify(nCellsSolve % next)
-       nullify(cellLimitField % next)
-       nullify(offSetField % next)
-     else
-       nullify(nCellsSolve)
-       nullify(cellLimitField)
-       nullify(offSetField)
-     end if
+     nullify(nCellsSolve % next)
+     nullify(cellLimitField % next)
+     nullify(offSetField % next)
 
      ! Loop over blocks
      offSetCursor =&gt; offsetField
@@ -813,20 +789,15 @@
 
      ! Allocate some needed arrays and fields
      allocate(sendingHaloLayers(1))
-     if(associated(indexToEdgeID)) then
-       allocate(haloIndices)
-       allocate(offSetField)
-       allocate(edgeLimitField)
 
-       nullify(haloIndices % next)
-       nullify(offSetField % next)
-       nullify(edgeLimitField % next)
-     else
-       nullify(haloIndices)
-       nullify(offSetField)
-       nullify(edgeLimitField)
-     end if
+     allocate(haloIndices)
+     allocate(offSetField)
+     allocate(edgeLimitField)
 
+     nullify(haloIndices % next)
+     nullify(offSetField % next)
+     nullify(edgeLimitField % next)
+
      ! Determine number of blocks, and setup field lists
      ! Loop over blocks
      nBlocks = 0

Modified: branches/omp_blocks/multiple_blocks/src/framework/mpas_block_decomp.F
===================================================================
--- branches/omp_blocks/multiple_blocks/src/framework/mpas_block_decomp.F        2012-07-23 17:39:48 UTC (rev 2041)
+++ branches/omp_blocks/multiple_blocks/src/framework/mpas_block_decomp.F        2012-07-23 18:06:01 UTC (rev 2042)
@@ -5,6 +5,7 @@
    use mpas_hash
    use mpas_sort
    use mpas_grid_types
+   use mpas_configure
 
    type graph
       integer :: nVerticesTotal
@@ -47,6 +48,10 @@
       integer, dimension(:), pointer :: local_nvertices
       character (len=StrKIND) :: filename
 
+      logical :: no_blocks
+
+      no_blocks = .false.
+
       if(config_number_of_blocks == 0) then
         total_blocks = dminfo % nProcs
       else
@@ -142,71 +147,83 @@
                                    global_start, local_nvertices, global_block_list, local_block_list)
         end if
 
-        allocate(sorted_local_cell_list(2, local_nvertices(dminfo % my_proc_id + 1)))
-        allocate(block_id(blocks_per_proc))
-        allocate(block_start(blocks_per_proc))
-        allocate(block_count(blocks_per_proc))
+        if(blocks_per_proc == 0) then
+           no_blocks = .true.
+           blocks_per_proc = 1
+        end if
 
-        do i = 1, blocks_per_proc
-          block_start = 0
-          block_count = 0
-        end do
+        if(no_blocks) then
+           allocate(block_id(blocks_per_proc))
+           allocate(block_start(blocks_per_proc))
+           allocate(block_count(blocks_per_proc))
 
-        do i = 1,local_nvertices(dminfo % my_proc_id +1)
-          call mpas_get_local_block_id(dminfo, local_block_list(i), local_block_id)
-  
-          block_id(local_block_id+1) = local_block_list(i)
-  
-          sorted_local_cell_list(1, i) = local_block_list(i)
-          sorted_local_cell_list(2, i) = local_cell_list(i)
-  
-          block_count(local_block_id+1) = block_count(local_block_id+1) + 1
-        end do
+           block_id(1) = config_number_of_blocks + 1
+           block_start(1) = 0
+           block_count(1) = 0
+        else
+           allocate(sorted_local_cell_list(2, local_nvertices(dminfo % my_proc_id + 1)))
+           allocate(block_id(blocks_per_proc))
+           allocate(block_start(blocks_per_proc))
+           allocate(block_count(blocks_per_proc))
+   
+           do i = 1, blocks_per_proc
+             block_start = 0
+             block_count = 0
+           end do
+   
+           do i = 1,local_nvertices(dminfo % my_proc_id +1)
+             call mpas_get_local_block_id(dminfo, local_block_list(i), local_block_id)
+     
+             block_id(local_block_id+1) = local_block_list(i)
+     
+             sorted_local_cell_list(1, i) = local_block_list(i)
+             sorted_local_cell_list(2, i) = local_cell_list(i)
+     
+             block_count(local_block_id+1) = block_count(local_block_id+1) + 1
+           end do
+   
+           call mpas_quicksort(local_nvertices(dminfo % my_proc_id + 1), sorted_local_cell_list)
+   
+           do i = 1, local_nvertices(dminfo % my_proc_id+1)
+             local_cell_list(i) = sorted_local_cell_list(2, i)
+           end do
+   
+           do i = 2,blocks_per_proc
+             block_start(i) = block_start(i-1) + block_count(i-1)
+           end do
 
-        call mpas_quicksort(local_nvertices(dminfo % my_proc_id + 1), sorted_local_cell_list)
+           deallocate(sorted_local_cell_list)
+           deallocate(local_block_list)
+           deallocate(local_nvertices)
+           deallocate(global_start)
+           deallocate(global_cell_list)
+           deallocate(global_block_list)
+        end if
+      else
 
-        do i = 1, local_nvertices(dminfo % my_proc_id+1)
-          local_cell_list(i) = sorted_local_cell_list(2, i)
-        end do
-
-        do i = 2,blocks_per_proc
-          block_start(i) = block_start(i-1) + block_count(i-1)
-        end do
-
-        !dwj 01/31/12 debugging multiple blocks
-!       do i=1,local_nvertices(dminfo % my_proc_id +1)
-!         call mpas_get_local_block_id(dminfo, sorted_local_cell_list(1, i), local_block_id)
-!         write(*,*) sorted_local_cell_list(1, i), local_block_id, sorted_local_cell_list(2,i)
-!       end do
-  
-        deallocate(sorted_local_cell_list)
-        deallocate(local_block_list)
-        deallocate(local_nvertices)
-        deallocate(global_start)
-        deallocate(global_cell_list)
-        deallocate(global_block_list)
-      else
-        allocate(local_cell_list(partial_global_graph_info % nVerticesTotal))
-        allocate(block_id(1))
-        allocate(block_start(1))
-        allocate(block_count(1))
-        block_id(1) = 0
-        block_start(1) = 0
-        block_count(1) = size(local_cell_list)
-        do i=1,size(local_cell_list)
-          local_cell_list(i) = i
-        end do
+        if (dminfo % my_proc_id == IO_NODE) then
+           allocate(local_cell_list(partial_global_graph_info % nVerticesTotal))
+           allocate(block_id(1))
+           allocate(block_start(1))
+           allocate(block_count(1))
+           block_id(1) = 0
+           block_start(1) = 0
+           block_count(1) = size(local_cell_list)
+           do i=1,size(local_cell_list)
+             local_cell_list(i) = i
+           end do
+        else
+           allocate(local_cell_list(1))
+           allocate(block_id(1))
+           allocate(block_start(1))
+           allocate(block_count(1))
+           local_cell_list(1) = 0
+           block_id(1) = config_number_of_blocks + 1
+           block_start(1) = 0
+           block_count(1) = 0
+        end if
       end if
 
-      !dwj 01/31/12 debugging multiple blocks
-!     write(*,*) 'Blocks per proc = ', blocks_per_proc, 'total_blocks = ', total_blocks
-
-!     do i=1,blocks_per_proc
-!       write(*,*) block_id(i), block_start(i), block_count(i)
-!     end do
-
-!     call mpas_dmpar_abort(dminfo)
-
    end subroutine mpas_block_decomp_cells_for_proc!}}}
 
    subroutine mpas_block_decomp_partitioned_edge_list(nCells, cellIDList, maxCells, nEdges, cellsOnEdge, edgeIDList, ghostEdgeStart)!{{{

</font>
</pre>