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> 10c4762a1bSJed Brown use petscvec 11c4762a1bSJed Brown implicit none 12c4762a1bSJed Brown 13c4762a1bSJed Brown PetscErrorCode ierr 14c4762a1bSJed Brown PetscInt i,nroots,nrootsalloc,nleaves,nleavesalloc,mine(6),stride 15c4762a1bSJed Brown type(PetscSFNode) remote(6) 16c4762a1bSJed Brown PetscMPIInt rank,size 17c4762a1bSJed Brown PetscSF sf 18c4762a1bSJed Brown PetscInt rootdata(6),leafdata(6) 19c4762a1bSJed Brown 20c4762a1bSJed Brown! used with PetscSFGetGraph() 21c4762a1bSJed Brown type(PetscSFNode), pointer :: gremote(:) 22c4762a1bSJed Brown PetscInt, pointer :: gmine(:) 23c4762a1bSJed Brown PetscInt gnroots,gnleaves; 24c4762a1bSJed Brown 25c4762a1bSJed Brown call PetscInitialize(PETSC_NULL_CHARACTER,ierr) 26c4762a1bSJed Brown if (ierr .ne. 0) then 27c4762a1bSJed Brown print*,'Unable to initialize PETSc' 28c4762a1bSJed Brown stop 29c4762a1bSJed Brown endif 30c4762a1bSJed Brown stride = 2 31c4762a1bSJed Brown call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr);CHKERRA(ierr); 32c4762a1bSJed Brown call MPI_Comm_size(PETSC_COMM_WORLD,size,ierr);CHKERRA(ierr); 33c4762a1bSJed Brown 34c4762a1bSJed Brown if (rank == 0) then 35c4762a1bSJed Brown nroots = 3 36c4762a1bSJed Brown else 37c4762a1bSJed Brown nroots = 2 38c4762a1bSJed Brown endif 39c4762a1bSJed Brown nrootsalloc = nroots * stride; 40c4762a1bSJed Brown if (rank > 0) then 41c4762a1bSJed Brown nleaves = 3 42c4762a1bSJed Brown else 43c4762a1bSJed Brown nleaves = 2 44c4762a1bSJed Brown endif 45c4762a1bSJed Brown nleavesalloc = nleaves * stride 46c4762a1bSJed Brown if (stride > 1) then 47c4762a1bSJed Brown do i=1,nleaves 48c4762a1bSJed Brown mine(i) = stride * (i-1) 49c4762a1bSJed Brown enddo 50c4762a1bSJed Brown endif 51c4762a1bSJed Brown 52c4762a1bSJed Brown! Left periodic neighbor 53c4762a1bSJed Brown remote(1)%rank = modulo(rank+size-1,size) 54c4762a1bSJed Brown remote(1)%index = 1 * stride 55c4762a1bSJed Brown! Right periodic neighbor 56c4762a1bSJed Brown remote(2)%rank = modulo(rank+1,size) 57c4762a1bSJed Brown remote(2)%index = 0 * stride 58c4762a1bSJed Brown if (rank > 0) then ! All processes reference rank 0, index 59c4762a1bSJed Brown remote(3)%rank = 0 60c4762a1bSJed Brown remote(3)%index = 2 * stride 61c4762a1bSJed Brown endif 62c4762a1bSJed Brown 63c4762a1bSJed Brown! Create a star forest for communication 64c4762a1bSJed Brown call PetscSFCreate(PETSC_COMM_WORLD,sf,ierr);CHKERRA(ierr); 65c4762a1bSJed Brown call PetscSFSetFromOptions(sf,ierr);CHKERRA(ierr); 66c4762a1bSJed Brown call PetscSFSetGraph(sf,nrootsalloc,nleaves,mine,PETSC_COPY_VALUES,remote,PETSC_COPY_VALUES,ierr);CHKERRA(ierr); 67c4762a1bSJed Brown call PetscSFSetUp(sf,ierr);CHKERRA(ierr); 68c4762a1bSJed Brown 69c4762a1bSJed Brown! View graph, mostly useful for debugging purposes. 70c4762a1bSJed Brown call PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL,ierr);CHKERRA(ierr); 71c4762a1bSJed Brown call PetscSFView(sf,PETSC_VIEWER_STDOUT_WORLD,ierr);CHKERRA(ierr); 72c4762a1bSJed Brown call PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD,ierr);CHKERRA(ierr); 73c4762a1bSJed Brown 74*2d4ee042Sprj-! Allocate space for send and receive buffers. This example communicates PetscInt, but other types, including 75c4762a1bSJed Brown! * user-defined structures, could also be used. 76c4762a1bSJed Brown! Set rootdata buffer to be broadcast 77c4762a1bSJed Brown do i=1,nrootsalloc 78c4762a1bSJed Brown rootdata(i) = -1 79c4762a1bSJed Brown enddo 80c4762a1bSJed Brown do i=1,nroots 81c4762a1bSJed Brown rootdata(1 + (i-1)*stride) = 100*(rank+1) + i - 1 82c4762a1bSJed Brown enddo 83c4762a1bSJed Brown 84c4762a1bSJed Brown! Initialize local buffer, these values are never used. 85c4762a1bSJed Brown do i=1,nleavesalloc 86c4762a1bSJed Brown leafdata(i) = -1 87c4762a1bSJed Brown enddo 88c4762a1bSJed Brown 89c4762a1bSJed Brown! Broadcast entries from rootdata to leafdata. Computation or other communication can be performed between the begin and end calls. 90c4762a1bSJed Brown call PetscSFBcastBegin(sf,MPIU_INTEGER,rootdata,leafdata,ierr);CHKERRA(ierr); 91c4762a1bSJed Brown call PetscSFBcastEnd(sf,MPIU_INTEGER,rootdata,leafdata,ierr);CHKERRA(ierr); 92c4762a1bSJed Brown call PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Bcast Rootdata\n",ierr);CHKERRA(ierr); 93c4762a1bSJed Brown call PetscIntView(nrootsalloc,rootdata,PETSC_VIEWER_STDOUT_WORLD,ierr);CHKERRA(ierr); 94c4762a1bSJed Brown call PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Bcast Leafdata\n",ierr);CHKERRA(ierr); 95c4762a1bSJed Brown call PetscIntView(nleavesalloc,leafdata,PETSC_VIEWER_STDOUT_WORLD,ierr);CHKERRA(ierr); 96c4762a1bSJed Brown 97c4762a1bSJed Brown call PetscSFGetGraph(sf,gnroots,gnleaves,gmine,gremote,ierr);CHKERRA(ierr); 98c4762a1bSJed 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 99c4762a1bSJed Brown do i=1,nleaves 100c4762a1bSJed 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 101c4762a1bSJed Brown enddo 102c4762a1bSJed Brown do i=1,nleaves 103c4762a1bSJed 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 104c4762a1bSJed Brown enddo 105c4762a1bSJed Brown 106c4762a1bSJed Brown deallocate(gremote) 107c4762a1bSJed Brown! Clean storage for star forest. 108c4762a1bSJed Brown call PetscSFDestroy(sf,ierr);CHKERRA(ierr); 109c4762a1bSJed Brown call PetscFinalize(ierr); 110c4762a1bSJed Brown end 111c4762a1bSJed Brown 112c4762a1bSJed Brown!/*TEST 113c4762a1bSJed Brown! build: 114c4762a1bSJed Brown! requires: define(PETSC_HAVE_FORTRAN_TYPE_STAR) 115c4762a1bSJed Brown! 116c4762a1bSJed Brown! test: 117c4762a1bSJed Brown! nsize: 3 118c4762a1bSJed Brown! 119c4762a1bSJed Brown!TEST*/ 120