<p><b>mmwolf@sandia.gov</b> 2010-01-05 17:06:40 -0700 (Tue, 05 Jan 2010)</p><p>Added correct Zoltan query functions for hsfc ordering<br>
</p><hr noshade><pre><font color="gray">Modified: trunk/swmodel/src/module_io_input.F
===================================================================
--- trunk/swmodel/src/module_io_input.F        2010-01-05 21:48:08 UTC (rev 102)
+++ trunk/swmodel/src/module_io_input.F        2010-01-06 00:06:40 UTC (rev 103)
@@ -63,6 +63,13 @@
type (field2dInteger) :: verticesOnCellField
type (field2dInteger) :: cellsOnEdgeField
type (field2dInteger) :: cellsOnVertexField
+
+#ifdef HAVE_ZOLTAN
+#ifdef _MPI
+ type (field1dReal) :: xCellField, yCellField, zCellField
+#endif
+#endif
+
type (field1dReal) :: xtime
integer, dimension(:), pointer :: indexToCellID_0Halo
@@ -77,6 +84,12 @@
integer, dimension(:,:), pointer :: cellIDSorted
integer, dimension(:,:), pointer :: edgeIDSorted
integer, dimension(:,:), pointer :: vertexIDSorted
+
+#ifdef HAVE_ZOLTAN
+#ifdef _MPI
+ real (kind=RKIND), dimension(:), pointer :: xCell, yCell, zCell
+#endif
+#endif
integer, dimension(:), pointer :: local_cell_list, local_edge_list, local_vertex_list
integer, dimension(:), pointer :: local_vertlevel_list, needed_vertlevel_list
@@ -138,6 +151,35 @@
allocate(indexToCellIDField % array(nReadCells))
call io_input_field(input_obj, indexToCellIDField)
+#ifdef HAVE_ZOLTAN
+#ifdef _MPI
+ ! Cell x-coordinates (in 3d Cartesian space)
+ allocate(xCellField % ioinfo)
+ xCellField % ioinfo % fieldName = 'xCell'
+ xCellField % ioinfo % start(1) = readCellStart
+ xCellField % ioinfo % count(1) = nReadCells
+ allocate(xCellField % array(nReadCells))
+ call io_input_field(input_obj, xCellField)
+
+ ! Cell y-coordinates (in 3d Cartesian space)
+ allocate(yCellField % ioinfo)
+ yCellField % ioinfo % fieldName = 'yCell'
+ yCellField % ioinfo % start(1) = readCellStart
+ yCellField % ioinfo % count(1) = nReadCells
+ allocate(yCellField % array(nReadCells))
+ call io_input_field(input_obj, yCellField)
+
+ ! Cell z-coordinates (in 3d Cartesian space)
+ allocate(zCellField % ioinfo)
+ zCellField % ioinfo % fieldName = 'zCell'
+ zCellField % ioinfo % start(1) = readCellStart
+ zCellField % ioinfo % count(1) = nReadCells
+ allocate(zCellField % array(nReadCells))
+ call io_input_field(input_obj, zCellField)
+#endif
+#endif
+
+
! Global edge indices
allocate(indexToEdgeIDField % ioinfo)
indexToEdgeIDField % ioinfo % fieldName = 'indexToEdgeID'
@@ -254,6 +296,13 @@
allocate(nEdgesOnCell_0Halo(size(local_cell_list)))
allocate(cellsOnCell_0Halo(maxEdges, size(local_cell_list)))
+#ifdef HAVE_ZOLTAN
+#ifdef _MPI
+ allocate(xCell(size(local_cell_list)))
+ allocate(yCell(size(local_cell_list)))
+ allocate(zCell(size(local_cell_list)))
+#endif
+#endif
!
! Now that each process has a list of cells that it owns, exchange cell connectivity
@@ -276,6 +325,23 @@
size(cellsOnCellField % array, 1), size(indexToCellIDField % array), size(local_cell_list), &
sendCellList, recvCellList)
+#ifdef HAVE_ZOLTAN
+#ifdef _MPI
+ call dmpar_alltoall_field(domain % dminfo, xCellField % array, xCell, &
+ size(xCellField % array), size(local_cell_list), &
+ sendCellList, recvCellList)
+
+ call dmpar_alltoall_field(domain % dminfo, yCellField % array, yCell, &
+ size(yCellField % array), size(local_cell_list), &
+ sendCellList, recvCellList)
+
+ call dmpar_alltoall_field(domain % dminfo, zCellField % array, zCell, &
+ size(zCellField % array), size(local_cell_list), &
+ sendCellList, recvCellList)
+#endif
+#endif
+
+
deallocate(sendCellList % list)
deallocate(sendCellList)
deallocate(recvCellList % list)
@@ -342,7 +408,7 @@
call zoltanStart()
!! Zoltan hook for cells
- call zoltanOrderLocHSFC_Cells(block_graph_2Halo%nVertices,block_graph_2Halo%VertexID)
+ call zoltanOrderLocHSFC_Cells(block_graph_2Halo%nVertices,block_graph_2Halo%VertexID,3,xCell,yCell,zCell)
#endif
#endif
@@ -718,6 +784,16 @@
!
deallocate(indexToCellIDField % ioinfo)
deallocate(indexToCellIDField % array)
+#ifdef HAVE_ZOLTAN
+#ifdef _MPI
+ deallocate(xCellField % ioinfo)
+ deallocate(xCellField % array)
+ deallocate(yCellField % ioinfo)
+ deallocate(yCellField % array)
+ deallocate(zCellField % ioinfo)
+ deallocate(zCellField % array)
+#endif
+#endif
deallocate(indexToEdgeIDField % ioinfo)
deallocate(indexToEdgeIDField % array)
deallocate(indexToVertexIDField % ioinfo)
@@ -742,7 +818,13 @@
deallocate(block_graph_0Halo % vertexID)
deallocate(block_graph_0Halo % nAdjacent)
deallocate(block_graph_0Halo % adjacencyList)
-
+#ifdef HAVE_ZOLTAN
+#ifdef _MPI
+ deallocate(xCell)
+ deallocate(yCell)
+ deallocate(zCell)
+#endif
+#endif
end subroutine input_state_for_domain
Modified: trunk/swmodel/src/module_zoltan_interface.F
===================================================================
--- trunk/swmodel/src/module_zoltan_interface.F        2010-01-05 21:48:08 UTC (rev 102)
+++ trunk/swmodel/src/module_zoltan_interface.F        2010-01-06 00:06:40 UTC (rev 103)
@@ -6,6 +6,8 @@
integer :: numCells
integer, dimension(:), pointer :: cellIDs
+ integer :: geomDim
+ real (kind=RKIND), dimension(:), pointer :: cellCoordX, cellCoordY, cellCoordZ
contains
@@ -29,12 +31,15 @@
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
- subroutine zoltanOrderLocHSFC_Cells(in_numcells,in_cellIDs)
+ subroutine zoltanOrderLocHSFC_Cells(in_numcells,in_cellIDs,in_geomDim,in_cellX, &
+ in_cellY, in_cellZ)
implicit none
integer :: in_numcells
integer, dimension(:), pointer :: in_cellIDs
+ integer :: in_geomDim
+ real (kind=RKIND), dimension(:), pointer :: in_cellX, in_cellY, in_cellZ
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!! local variables
@@ -43,7 +48,7 @@
type(Zoltan_Struct), pointer :: zz_obj
integer(ZOLTAN_INT) :: ierr
- integer :: numGidEntries, myrank, locNumObj, i
+ integer :: numGidEntries, i
integer(ZOLTAN_INT), allocatable :: global_ids(:)
integer(ZOLTAN_INT), allocatable :: permGIDs(:)
integer(ZOLTAN_INT), allocatable :: permIndices(:)
@@ -54,6 +59,10 @@
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
numCells = in_numcells
cellIDs => in_cellIDs
+ geomDim = in_geomDim
+ cellCoordX => in_cellX
+ cellCoordY => in_cellY
+ cellCoordZ => in_cellZ
nullify(zz_obj)
zz_obj => Zoltan_Create(MPI_COMM_SELF)
@@ -71,31 +80,28 @@
ierr = Zoltan_Set_Fn(zz_obj, ZOLTAN_NUM_GEOM_FN_TYPE,zqfGeomDim)
ierr = Zoltan_Set_Fn(zz_obj, ZOLTAN_GEOM_FN_TYPE, zqfGetCellGeom)
- locNumObj=16
numGidEntries=1
- allocate(global_ids(locNumObj))
- allocate(permGIDs(locNumObj))
- allocate(permIndices(locNumObj))
+ allocate(global_ids(numCells))
+ allocate(permGIDs(numCells))
+ allocate(permIndices(numCells))
- call MPI_Comm_rank(MPI_COMM_WORLD, myrank, ierr)
-
- do i=1,locNumObj
- global_ids(i) = i
+ do i=1,numCells
+ global_ids(i) = cellIDs(i)
end do
- ierr = Zoltan_Order(zz_obj, numGidEntries, locNumObj, global_ids, permIndices);
+ ierr = Zoltan_Order(zz_obj, numGidEntries, numCells, global_ids, permIndices);
!!!!!!!!!!!!!!!!!!!!!!!!!!
!! This is necessary for now until we fix a small bug in Zoltan_Order
!!!!!!!!!!!!!!!!!!!!!!!!!!
- do i=1,locNumObj
+ do i=1,numCells
permGIDs(i) = global_ids(permIndices(i)+1)
end do
-!! do i=1,locNumObj
-!! write(*,*) global_ids(i), permGIDs(i)
-!! end do
+ !!do i=1,numCells
+ !! write(*,*) global_ids(i), permGIDs(i)
+ !!end do
deallocate(global_ids)
deallocate(permGIDs)
@@ -119,7 +125,7 @@
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!! MMW: Should be set to numCells
- zqfNumCells = 16
+ zqfNumCells = numCells
ierr = ZOLTAN_OK
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
end function zqfNumCells
@@ -149,8 +155,8 @@
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
!! MMW: cellIDs should have global cell IDs
- do i= 1, 16
- global_ids(i) = i
+ do i= 1, numCells
+ global_ids(i) = cellIDs(i)
local_ids(i) = i
end do
@@ -170,7 +176,8 @@
integer(ZOLTAN_INT) :: ierr
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
- zqfGeomDim = 2
+ !! MMW should set to geomDim
+ zqfGeomDim = geomDim
ierr = ZOLTAN_OK
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
@@ -194,9 +201,15 @@
integer(ZOLTAN_INT), intent(out) :: ierr
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
- geom_vec(1) = 4*( (global_id-1)/32 ) + MODULO(global_id-1,4) + 1;
- geom_vec(2) = MODULO(global_id-1,32) /4 + 1;
+ !!geom_vec(1) = 4*( (global_id-1)/32 ) + MODULO(global_id-1,4) + 1;
+ !!geom_vec(2) = MODULO(global_id-1,32) /4 + 1;
+ !! MMW write(*,*) cellCoordX(global_id), cellCoordY(global_id), cellCoordZ(global_id)
+ !! Assuming geom_dim is 3
+ geom_vec(1) = cellCoordX(local_id)
+ geom_vec(2) = cellCoordY(local_id)
+ geom_vec(3) = cellCoordZ(local_id)
+
ierr = ZOLTAN_OK
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
</font>
</pre>