Actual source code: ex1f.F90

  1: !    Description: A star forest is a simple tree with one root and zero or more leaves.
  2: !    Many common communication patterns can be expressed as updates of rootdata using leafdata and vice-versa.
  3: !     This example creates a star forest, communicates values using the graph  views the graph, then destroys it.
  4: !
  5: !     This is a copy of ex1.c but currently only tests the broadcast operation
  6: #include <petsc/finclude/petscvec.h>
  7: program main
  8:   use petscvec
  9:   implicit none

 11:   PetscErrorCode ierr
 12:   PetscInt i, nroots, nrootsalloc, nleaves, nleavesalloc, mine(6)
 13:   PetscInt, parameter :: stride = 2
 14:   PetscSFNode remote(6)
 15:   PetscMPIInt rank, size
 16:   PetscSF sf
 17:   PetscInt rootdata(6), leafdata(6)

 19: ! used with PetscSFGetGraph()
 20:   PetscSFNode, pointer ::       gremote(:)
 21:   PetscInt, pointer ::          gmine(:)
 22:   PetscInt gnroots, gnleaves

 24:   PetscMPIInt niranks, nranks
 25:   PetscMPIInt, pointer ::       iranks(:), ranks(:)
 26:   PetscInt, pointer ::          ioffset(:), irootloc(:), roffset(:), rmine(:), rremote(:)

 28:   PetscCallA(PetscInitialize(ierr))
 29:   PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD, rank, ierr))
 30:   PetscCallMPIA(MPI_Comm_size(PETSC_COMM_WORLD, size, ierr))

 32:   if (rank == 0) then
 33:     nroots = 3
 34:   else
 35:     nroots = 2
 36:   end if
 37:   nrootsalloc = nroots*stride
 38:   if (rank > 0) then
 39:     nleaves = 3
 40:   else
 41:     nleaves = 2
 42:   end if
 43:   nleavesalloc = nleaves*stride
 44:   if (stride > 1) then
 45:     do i = 1, nleaves
 46:       mine(i) = stride*(i - 1)
 47:     end do
 48:   end if

 50: ! Left periodic neighbor
 51:   remote(1)%rank = modulo(rank + size - 1, size)
 52:   remote(1)%index = 1*stride
 53: ! Right periodic neighbor
 54:   remote(2)%rank = modulo(rank + 1, size)
 55:   remote(2)%index = 0*stride
 56:   if (rank > 0) then !               All processes reference rank 0, index
 57:     remote(3)%rank = 0
 58:     remote(3)%index = 2*stride
 59:   end if

 61: ! Create a star forest for communication
 62:   PetscCallA(PetscSFCreate(PETSC_COMM_WORLD, sf, ierr))
 63:   PetscCallA(PetscSFSetFromOptions(sf, ierr))
 64:   PetscCallA(PetscSFSetGraph(sf, nrootsalloc, nleaves, mine, PETSC_COPY_VALUES, remote, PETSC_COPY_VALUES, ierr))
 65:   PetscCallA(PetscSFSetUp(sf, ierr))

 67: ! View graph, mostly useful for debugging purposes.
 68:   PetscCallA(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_INFO_DETAIL, ierr))
 69:   PetscCallA(PetscSFView(sf, PETSC_VIEWER_STDOUT_WORLD, ierr))
 70:   PetscCallA(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD, ierr))

 72: ! Allocate space for send and receive buffers. This example communicates PetscInt, but other types, including
 73: !   * user-defined structures, could also be used.
 74: ! Set rootdata buffer to be broadcast
 75:   rootdata(1:nrootsalloc) = -1
 76:   do i = 1, nroots
 77:     rootdata(1 + (i - 1)*stride) = 100*(rank + 1) + i - 1
 78:   end do

 80: ! Initialize local buffer, these values are never used.
 81:   leafdata(1:nleavesalloc) = -1

 83: ! Broadcast entries from rootdata to leafdata. Computation or other communication can be performed between the begin and end calls.
 84:   PetscCallA(PetscSFBcastBegin(sf, MPIU_INTEGER, rootdata, leafdata, MPI_REPLACE, ierr))
 85:   PetscCallA(PetscSFBcastEnd(sf, MPIU_INTEGER, rootdata, leafdata, MPI_REPLACE, ierr))
 86:   PetscCallA(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, '## Bcast Rootdata\n', ierr))
 87:   PetscCallA(PetscIntView(nrootsalloc, rootdata, PETSC_VIEWER_STDOUT_WORLD, ierr))
 88:   PetscCallA(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, '## Bcast Leafdata\n', ierr))
 89:   PetscCallA(PetscIntView(nleavesalloc, leafdata, PETSC_VIEWER_STDOUT_WORLD, ierr))

 91: ! Reduce entries from leafdata to rootdata. Computation or other communication can be performed between the begin and end calls.
 92:   PetscCallA(PetscSFReduceBegin(sf, MPIU_INTEGER, leafdata, rootdata, MPI_SUM, ierr))
 93:   PetscCallA(PetscSFReduceEnd(sf, MPIU_INTEGER, leafdata, rootdata, MPI_SUM, ierr))
 94:   PetscCallA(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, '## Reduce Leafdata\n', ierr))
 95:   PetscCallA(PetscIntView(nleavesalloc, leafdata, PETSC_VIEWER_STDOUT_WORLD, ierr))
 96:   PetscCallA(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, '## Reduce Rootdata\n', ierr))
 97:   PetscCallA(PetscIntView(nrootsalloc, rootdata, PETSC_VIEWER_STDOUT_WORLD, ierr))

 99:   PetscCallA(PetscSFGetGraph(sf, gnroots, gnleaves, gmine, gremote, ierr))
100:   PetscCheckA(gnleaves == nleaves, PETSC_COMM_WORLD, PETSC_ERR_PLIB, 'nleaves returned from PetscSFGetGraph() does not match that set with PetscSFSetGraph()')
101:   do i = 1, nleaves
102:     PetscCheckA(gmine(i) == mine(i), PETSC_COMM_WORLD, PETSC_ERR_PLIB, 'Root from PetscSFGetGraph() does not match that set with PetscSFSetGraph()')
103:   end do
104:   do i = 1, nleaves
105:     PetscCheckA(gremote(i)%index == remote(i)%index, PETSC_COMM_WORLD, PETSC_ERR_PLIB, 'Leaf from PetscSFGetGraph() does not match that set with PetscSFSetGraph()')
106:   end do
107:   PetscCallA(PetscSFRestoreGraph(sf, gnroots, gnleaves, gmine, gremote, ierr))

109: ! Test PetscSFGet{Leaf,Root}Ranks
110:   PetscCallA(PetscSFGetLeafRanks(sf, niranks, iranks, ioffset, irootloc, ierr))
111:   PetscCallA(PetscSFGetRootRanks(sf, nranks, ranks, roffset, rmine, rremote, ierr))

113: !    Clean storage for star forest.
114:   PetscCallA(PetscSFDestroy(sf, ierr))

116: !  Create a star forest with continuous leaves and hence no buffer
117:   PetscCallA(PetscSFCreate(PETSC_COMM_WORLD, sf, ierr))
118:   PetscCallA(PetscSFSetFromOptions(sf, ierr))
119:   PetscCallA(PetscSFSetGraph(sf, nrootsalloc, nleaves, PETSC_NULL_INTEGER_ARRAY, PETSC_COPY_VALUES, remote, PETSC_COPY_VALUES, ierr))
120:   PetscCallA(PetscSFSetUp(sf, ierr))

122: !   View graph, mostly useful for debugging purposes.
123:   PetscCallA(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_INFO_DETAIL, ierr))
124:   PetscCallA(PetscSFView(sf, PETSC_VIEWER_STDOUT_WORLD, ierr))
125:   PetscCallA(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD, ierr))

127:   PetscCallA(PetscSFGetGraph(sf, gnroots, gnleaves, gmine, gremote, ierr))
128:   PetscCheckA(loc(gmine) == loc(PETSC_NULL_INTEGER), PETSC_COMM_WORLD, PETSC_ERR_PLIB, 'Leaves from PetscSFGetGraph() not null as expected')
129:   PetscCallA(PetscSFRestoreGraph(sf, gnroots, gnleaves, gmine, gremote, ierr))
130:   PetscCallA(PetscSFDestroy(sf, ierr))
131:   PetscCallA(PetscFinalize(ierr))
132: end

134: !/*TEST
135: !
136: !  test:
137: !    nsize: 3
138: !
139: !TEST*/