1c4762a1bSJed Brown 2c4762a1bSJed Brown! Description: A star forest is a simple tree with one root and zero or more leaves. 3c4762a1bSJed Brown! Many common communication patterns can be expressed as updates of rootdata using leafdata and vice-versa. 4c4762a1bSJed Brown! This example creates a star forest, communicates values using the graph views the graph, then destroys it. 5c4762a1bSJed Brown! 6c4762a1bSJed Brown! This is a copy of ex1.c but currently only tests the broadcast operation 7c4762a1bSJed Brown 8c4762a1bSJed Brown program main 9c4762a1bSJed Brown#include <petsc/finclude/petscvec.h> 108c8af28eSPedro Ricardo C. Souza use petscmpi ! or mpi or mpi_f08 11c4762a1bSJed Brown use petscvec 12c4762a1bSJed Brown implicit none 13c4762a1bSJed Brown 14c4762a1bSJed Brown PetscErrorCode ierr 15c4762a1bSJed Brown PetscInt i,nroots,nrootsalloc,nleaves,nleavesalloc,mine(6),stride 16c4762a1bSJed Brown type(PetscSFNode) remote(6) 17c4762a1bSJed Brown PetscMPIInt rank,size 18c4762a1bSJed Brown PetscSF sf 19c4762a1bSJed Brown PetscInt rootdata(6),leafdata(6) 20c4762a1bSJed Brown 21c4762a1bSJed Brown! used with PetscSFGetGraph() 22c4762a1bSJed Brown type(PetscSFNode), pointer :: gremote(:) 23c4762a1bSJed Brown PetscInt, pointer :: gmine(:) 24c4762a1bSJed Brown PetscInt gnroots,gnleaves; 25c4762a1bSJed Brown 2694a885e8SJunchao Zhang PetscInt niranks,nranks 2794a885e8SJunchao Zhang PetscMPIInt, pointer :: iranks(:), ranks(:) 2894a885e8SJunchao Zhang PetscInt, pointer :: ioffset(:),irootloc(:),roffset(:),rmine(:),rremote(:) 2994a885e8SJunchao Zhang 30d8606c27SBarry Smith PetscCallA(PetscInitialize(ierr)) 31c4762a1bSJed Brown stride = 2 32d8606c27SBarry Smith PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)) 33d8606c27SBarry Smith PetscCallMPIA(MPI_Comm_size(PETSC_COMM_WORLD,size,ierr)) 34c4762a1bSJed Brown 35c4762a1bSJed Brown if (rank == 0) then 36c4762a1bSJed Brown nroots = 3 37c4762a1bSJed Brown else 38c4762a1bSJed Brown nroots = 2 39c4762a1bSJed Brown endif 40c4762a1bSJed Brown nrootsalloc = nroots * stride; 41c4762a1bSJed Brown if (rank > 0) then 42c4762a1bSJed Brown nleaves = 3 43c4762a1bSJed Brown else 44c4762a1bSJed Brown nleaves = 2 45c4762a1bSJed Brown endif 46c4762a1bSJed Brown nleavesalloc = nleaves * stride 47c4762a1bSJed Brown if (stride > 1) then 48c4762a1bSJed Brown do i=1,nleaves 49c4762a1bSJed Brown mine(i) = stride * (i-1) 50c4762a1bSJed Brown enddo 51c4762a1bSJed Brown endif 52c4762a1bSJed Brown 53c4762a1bSJed Brown! Left periodic neighbor 54c4762a1bSJed Brown remote(1)%rank = modulo(rank+size-1,size) 55c4762a1bSJed Brown remote(1)%index = 1 * stride 56c4762a1bSJed Brown! Right periodic neighbor 57c4762a1bSJed Brown remote(2)%rank = modulo(rank+1,size) 58c4762a1bSJed Brown remote(2)%index = 0 * stride 59c4762a1bSJed Brown if (rank > 0) then ! All processes reference rank 0, index 60c4762a1bSJed Brown remote(3)%rank = 0 61c4762a1bSJed Brown remote(3)%index = 2 * stride 62c4762a1bSJed Brown endif 63c4762a1bSJed Brown 64c4762a1bSJed Brown! Create a star forest for communication 65d8606c27SBarry Smith PetscCallA(PetscSFCreate(PETSC_COMM_WORLD,sf,ierr)) 66d8606c27SBarry Smith PetscCallA(PetscSFSetFromOptions(sf,ierr)) 67d8606c27SBarry Smith PetscCallA(PetscSFSetGraph(sf,nrootsalloc,nleaves,mine,PETSC_COPY_VALUES,remote,PETSC_COPY_VALUES,ierr)) 68d8606c27SBarry Smith PetscCallA(PetscSFSetUp(sf,ierr)) 69c4762a1bSJed Brown 70c4762a1bSJed Brown! View graph, mostly useful for debugging purposes. 71d8606c27SBarry Smith PetscCallA(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL,ierr)) 72d8606c27SBarry Smith PetscCallA(PetscSFView(sf,PETSC_VIEWER_STDOUT_WORLD,ierr)) 73d8606c27SBarry Smith PetscCallA(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD,ierr)) 74c4762a1bSJed Brown 752d4ee042Sprj-! Allocate space for send and receive buffers. This example communicates PetscInt, but other types, including 76c4762a1bSJed Brown! * user-defined structures, could also be used. 77c4762a1bSJed Brown! Set rootdata buffer to be broadcast 78c4762a1bSJed Brown do i=1,nrootsalloc 79c4762a1bSJed Brown rootdata(i) = -1 80c4762a1bSJed Brown enddo 81c4762a1bSJed Brown do i=1,nroots 82c4762a1bSJed Brown rootdata(1 + (i-1)*stride) = 100*(rank+1) + i - 1 83c4762a1bSJed Brown enddo 84c4762a1bSJed Brown 85c4762a1bSJed Brown! Initialize local buffer, these values are never used. 86c4762a1bSJed Brown do i=1,nleavesalloc 87c4762a1bSJed Brown leafdata(i) = -1 88c4762a1bSJed Brown enddo 89c4762a1bSJed Brown 90c4762a1bSJed Brown! Broadcast entries from rootdata to leafdata. Computation or other communication can be performed between the begin and end calls. 91d8606c27SBarry Smith PetscCallA(PetscSFBcastBegin(sf,MPIU_INTEGER,rootdata,leafdata,MPI_REPLACE,ierr)) 92d8606c27SBarry Smith PetscCallA(PetscSFBcastEnd(sf,MPIU_INTEGER,rootdata,leafdata,MPI_REPLACE,ierr)) 93*dcb3e689SBarry Smith PetscCallA(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,'## Bcast Rootdata\n',ierr)) 94d8606c27SBarry Smith PetscCallA(PetscIntView(nrootsalloc,rootdata,PETSC_VIEWER_STDOUT_WORLD,ierr)) 95*dcb3e689SBarry Smith PetscCallA(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,'## Bcast Leafdata\n',ierr)) 96d8606c27SBarry Smith PetscCallA(PetscIntView(nleavesalloc,leafdata,PETSC_VIEWER_STDOUT_WORLD,ierr)) 97c4762a1bSJed Brown 989037d788SNicholas Arnold-Medabalimi! Reduce entries from leafdata to rootdata. Computation or other communication can be performed between the begin and end calls. 999037d788SNicholas Arnold-Medabalimi PetscCallA(PetscSFReduceBegin(sf,MPIU_INTEGER,leafdata,rootdata,MPI_SUM,ierr)) 1009037d788SNicholas Arnold-Medabalimi PetscCallA(PetscSFReduceEnd(sf,MPIU_INTEGER,leafdata,rootdata,MPI_SUM,ierr)) 101*dcb3e689SBarry Smith PetscCallA(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,'## Reduce Leafdata\n',ierr)) 1029037d788SNicholas Arnold-Medabalimi PetscCallA(PetscIntView(nleavesalloc,leafdata,PETSC_VIEWER_STDOUT_WORLD,ierr)) 103*dcb3e689SBarry Smith PetscCallA(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,'## Reduce Rootdata\n',ierr)) 1049037d788SNicholas Arnold-Medabalimi PetscCallA(PetscIntView(nrootsalloc,rootdata,PETSC_VIEWER_STDOUT_WORLD,ierr)) 1059037d788SNicholas Arnold-Medabalimi 106d8606c27SBarry Smith PetscCallA(PetscSFGetGraph(sf,gnroots,gnleaves,gmine,gremote,ierr)) 107*dcb3e689SBarry Smith PetscCheckA(gnleaves .eq. nleaves,PETSC_COMM_WORLD,PETSC_ERR_PLIB,'nleaves returned from PetscSFGetGraph() does not match that set with PetscSFSetGraph()') 108c4762a1bSJed Brown do i=1,nleaves 109*dcb3e689SBarry Smith PetscCheckA(gmine(i) .eq. mine(i),PETSC_COMM_WORLD,PETSC_ERR_PLIB,'Root from PetscSFGetGraph() does not match that set with PetscSFSetGraph()') 110c4762a1bSJed Brown enddo 111c4762a1bSJed Brown do i=1,nleaves 112*dcb3e689SBarry Smith PetscCheckA(gremote(i)%index .eq. remote(i)%index,PETSC_COMM_WORLD,PETSC_ERR_PLIB,'Leaf from PetscSFGetGraph() does not match that set with PetscSFSetGraph()') 113c4762a1bSJed Brown enddo 114c4762a1bSJed Brown 115c4762a1bSJed Brown deallocate(gremote) 11694a885e8SJunchao Zhang 11794a885e8SJunchao Zhang! Test PetscSFGet{Leaf,Root}Ranks 11894a885e8SJunchao Zhang PetscCallA(PetscSFGetLeafRanks(sf,niranks,iranks,ioffset,irootloc,ierr)) 11994a885e8SJunchao Zhang PetscCallA(PetscSFGetRootRanks(sf,nranks,ranks,roffset,rmine,rremote,ierr)) 12094a885e8SJunchao Zhang 121c4762a1bSJed Brown! Clean storage for star forest. 122d8606c27SBarry Smith PetscCallA(PetscSFDestroy(sf,ierr)) 1238dbb0df6SBarry Smith 124da81f932SPierre Jolivet! Create a star forest with continuous leaves and hence no buffer 125d8606c27SBarry Smith PetscCallA(PetscSFCreate(PETSC_COMM_WORLD,sf,ierr)) 126d8606c27SBarry Smith PetscCallA(PetscSFSetFromOptions(sf,ierr)) 127d8606c27SBarry Smith PetscCallA(PetscSFSetGraph(sf,nrootsalloc,nleaves,PETSC_NULL_INTEGER,PETSC_COPY_VALUES,remote,PETSC_COPY_VALUES,ierr)) 128d8606c27SBarry Smith PetscCallA(PetscSFSetUp(sf,ierr)) 1298dbb0df6SBarry Smith 1308dbb0df6SBarry Smith! View graph, mostly useful for debugging purposes. 131d8606c27SBarry Smith PetscCallA(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL,ierr)) 132d8606c27SBarry Smith PetscCallA(PetscSFView(sf,PETSC_VIEWER_STDOUT_WORLD,ierr)) 133d8606c27SBarry Smith PetscCallA(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD,ierr)) 1348dbb0df6SBarry Smith 135d8606c27SBarry Smith PetscCallA(PetscSFGetGraph(sf,gnroots,gnleaves,gmine,gremote,ierr)) 136*dcb3e689SBarry Smith PetscCheckA(loc(gmine) .eq. loc(PETSC_NULL_INTEGER),PETSC_COMM_WORLD,PETSC_ERR_PLIB,'Leaves from PetscSFGetGraph() not null as expected') 1378dbb0df6SBarry Smith deallocate(gremote) 138d8606c27SBarry Smith PetscCallA(PetscSFDestroy(sf,ierr)) 139d8606c27SBarry Smith PetscCallA(PetscFinalize(ierr)) 140c4762a1bSJed Brown end 141c4762a1bSJed Brown 142c4762a1bSJed Brown!/*TEST 143c4762a1bSJed Brown! build: 144dfd57a17SPierre Jolivet! requires: defined(PETSC_HAVE_FORTRAN_TYPE_STAR) 145c4762a1bSJed Brown! 146c4762a1bSJed Brown! test: 147c4762a1bSJed Brown! nsize: 3 148c4762a1bSJed Brown! 149c4762a1bSJed Brown!TEST*/ 150