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 26d8606c27SBarry Smith PetscCallA(PetscInitialize(ierr)) 27c4762a1bSJed Brown stride = 2 28d8606c27SBarry Smith PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)) 29d8606c27SBarry Smith PetscCallMPIA(MPI_Comm_size(PETSC_COMM_WORLD,size,ierr)) 30c4762a1bSJed Brown 31c4762a1bSJed Brown if (rank == 0) then 32c4762a1bSJed Brown nroots = 3 33c4762a1bSJed Brown else 34c4762a1bSJed Brown nroots = 2 35c4762a1bSJed Brown endif 36c4762a1bSJed Brown nrootsalloc = nroots * stride; 37c4762a1bSJed Brown if (rank > 0) then 38c4762a1bSJed Brown nleaves = 3 39c4762a1bSJed Brown else 40c4762a1bSJed Brown nleaves = 2 41c4762a1bSJed Brown endif 42c4762a1bSJed Brown nleavesalloc = nleaves * stride 43c4762a1bSJed Brown if (stride > 1) then 44c4762a1bSJed Brown do i=1,nleaves 45c4762a1bSJed Brown mine(i) = stride * (i-1) 46c4762a1bSJed Brown enddo 47c4762a1bSJed Brown endif 48c4762a1bSJed Brown 49c4762a1bSJed Brown! Left periodic neighbor 50c4762a1bSJed Brown remote(1)%rank = modulo(rank+size-1,size) 51c4762a1bSJed Brown remote(1)%index = 1 * stride 52c4762a1bSJed Brown! Right periodic neighbor 53c4762a1bSJed Brown remote(2)%rank = modulo(rank+1,size) 54c4762a1bSJed Brown remote(2)%index = 0 * stride 55c4762a1bSJed Brown if (rank > 0) then ! All processes reference rank 0, index 56c4762a1bSJed Brown remote(3)%rank = 0 57c4762a1bSJed Brown remote(3)%index = 2 * stride 58c4762a1bSJed Brown endif 59c4762a1bSJed Brown 60c4762a1bSJed Brown! Create a star forest for communication 61d8606c27SBarry Smith PetscCallA(PetscSFCreate(PETSC_COMM_WORLD,sf,ierr)) 62d8606c27SBarry Smith PetscCallA(PetscSFSetFromOptions(sf,ierr)) 63d8606c27SBarry Smith PetscCallA(PetscSFSetGraph(sf,nrootsalloc,nleaves,mine,PETSC_COPY_VALUES,remote,PETSC_COPY_VALUES,ierr)) 64d8606c27SBarry Smith PetscCallA(PetscSFSetUp(sf,ierr)) 65c4762a1bSJed Brown 66c4762a1bSJed Brown! View graph, mostly useful for debugging purposes. 67d8606c27SBarry Smith PetscCallA(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL,ierr)) 68d8606c27SBarry Smith PetscCallA(PetscSFView(sf,PETSC_VIEWER_STDOUT_WORLD,ierr)) 69d8606c27SBarry Smith PetscCallA(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD,ierr)) 70c4762a1bSJed Brown 712d4ee042Sprj-! Allocate space for send and receive buffers. This example communicates PetscInt, but other types, including 72c4762a1bSJed Brown! * user-defined structures, could also be used. 73c4762a1bSJed Brown! Set rootdata buffer to be broadcast 74c4762a1bSJed Brown do i=1,nrootsalloc 75c4762a1bSJed Brown rootdata(i) = -1 76c4762a1bSJed Brown enddo 77c4762a1bSJed Brown do i=1,nroots 78c4762a1bSJed Brown rootdata(1 + (i-1)*stride) = 100*(rank+1) + i - 1 79c4762a1bSJed Brown enddo 80c4762a1bSJed Brown 81c4762a1bSJed Brown! Initialize local buffer, these values are never used. 82c4762a1bSJed Brown do i=1,nleavesalloc 83c4762a1bSJed Brown leafdata(i) = -1 84c4762a1bSJed Brown enddo 85c4762a1bSJed Brown 86c4762a1bSJed Brown! Broadcast entries from rootdata to leafdata. Computation or other communication can be performed between the begin and end calls. 87d8606c27SBarry Smith PetscCallA(PetscSFBcastBegin(sf,MPIU_INTEGER,rootdata,leafdata,MPI_REPLACE,ierr)) 88d8606c27SBarry Smith PetscCallA(PetscSFBcastEnd(sf,MPIU_INTEGER,rootdata,leafdata,MPI_REPLACE,ierr)) 89d8606c27SBarry Smith PetscCallA(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Bcast Rootdata\n",ierr)) 90d8606c27SBarry Smith PetscCallA(PetscIntView(nrootsalloc,rootdata,PETSC_VIEWER_STDOUT_WORLD,ierr)) 91d8606c27SBarry Smith PetscCallA(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Bcast Leafdata\n",ierr)) 92d8606c27SBarry Smith PetscCallA(PetscIntView(nleavesalloc,leafdata,PETSC_VIEWER_STDOUT_WORLD,ierr)) 93c4762a1bSJed Brown 94d8606c27SBarry Smith PetscCallA(PetscSFGetGraph(sf,gnroots,gnleaves,gmine,gremote,ierr)) 95c4762a1bSJed Brown if (gnleaves .ne. nleaves) then; SETERRA(PETSC_COMM_WORLD,PETSC_ERR_PLIB,'nleaves returned from PetscSFGetGraph() does not match that set with PetscSFSetGraph()'); endif 96c4762a1bSJed Brown do i=1,nleaves 97c4762a1bSJed Brown if (gmine(i) .ne. mine(i)) then; SETERRA(PETSC_COMM_WORLD,PETSC_ERR_PLIB,'Root from PetscSFGetGraph() does not match that set with PetscSFSetGraph()'); endif 98c4762a1bSJed Brown enddo 99c4762a1bSJed Brown do i=1,nleaves 100c4762a1bSJed Brown if (gremote(i)%index .ne. remote(i)%index) then; SETERRA(PETSC_COMM_WORLD,PETSC_ERR_PLIB,'Leaf from PetscSFGetGraph() does not match that set with PetscSFSetGraph()'); endif 101c4762a1bSJed Brown enddo 102c4762a1bSJed Brown 103c4762a1bSJed Brown deallocate(gremote) 104c4762a1bSJed Brown! Clean storage for star forest. 105d8606c27SBarry Smith PetscCallA(PetscSFDestroy(sf,ierr)) 1068dbb0df6SBarry Smith 107*da81f932SPierre Jolivet! Create a star forest with continuous leaves and hence no buffer 108d8606c27SBarry Smith PetscCallA(PetscSFCreate(PETSC_COMM_WORLD,sf,ierr)) 109d8606c27SBarry Smith PetscCallA(PetscSFSetFromOptions(sf,ierr)) 110d8606c27SBarry Smith PetscCallA(PetscSFSetGraph(sf,nrootsalloc,nleaves,PETSC_NULL_INTEGER,PETSC_COPY_VALUES,remote,PETSC_COPY_VALUES,ierr)) 111d8606c27SBarry Smith PetscCallA(PetscSFSetUp(sf,ierr)) 1128dbb0df6SBarry Smith 1138dbb0df6SBarry Smith! View graph, mostly useful for debugging purposes. 114d8606c27SBarry Smith PetscCallA(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL,ierr)) 115d8606c27SBarry Smith PetscCallA(PetscSFView(sf,PETSC_VIEWER_STDOUT_WORLD,ierr)) 116d8606c27SBarry Smith PetscCallA(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD,ierr)) 1178dbb0df6SBarry Smith 118d8606c27SBarry Smith PetscCallA(PetscSFGetGraph(sf,gnroots,gnleaves,gmine,gremote,ierr)) 1198dbb0df6SBarry Smith if (loc(gmine) .ne. loc(PETSC_NULL_INTEGER)) then; SETERRA(PETSC_COMM_WORLD,PETSC_ERR_PLIB,'Leaves from PetscSFGetGraph() not null as expected'); endif 1208dbb0df6SBarry Smith deallocate(gremote) 121d8606c27SBarry Smith PetscCallA(PetscSFDestroy(sf,ierr)) 122d8606c27SBarry Smith PetscCallA(PetscFinalize(ierr)) 123c4762a1bSJed Brown end 124c4762a1bSJed Brown 125c4762a1bSJed Brown!/*TEST 126c4762a1bSJed Brown! build: 127dfd57a17SPierre Jolivet! requires: defined(PETSC_HAVE_FORTRAN_TYPE_STAR) 128c4762a1bSJed Brown! 129c4762a1bSJed Brown! test: 130c4762a1bSJed Brown! nsize: 3 131c4762a1bSJed Brown! 132c4762a1bSJed Brown!TEST*/ 133