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