xref: /petsc/src/vec/is/sf/tutorials/ex1f.F90 (revision 6497c311e7b976d467be1503c1effce92a60525c)
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
7c4762a1bSJed Brown      program main
8c4762a1bSJed Brown#include <petsc/finclude/petscvec.h>
98c8af28eSPedro Ricardo C. Souza      use petscmpi  ! or mpi or mpi_f08
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
25*6497c311SBarry Smith      PetscMPIInt                   niranks,nranks
2694a885e8SJunchao Zhang      PetscMPIInt, pointer ::       iranks(:), ranks(:)
2794a885e8SJunchao Zhang      PetscInt, pointer ::          ioffset(:),irootloc(:),roffset(:),rmine(:),rremote(:)
2894a885e8SJunchao Zhang
29d8606c27SBarry Smith      PetscCallA(PetscInitialize(ierr))
30c4762a1bSJed Brown      stride = 2
31d8606c27SBarry Smith      PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr))
32d8606c27SBarry Smith      PetscCallMPIA(MPI_Comm_size(PETSC_COMM_WORLD,size,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
64d8606c27SBarry Smith      PetscCallA(PetscSFCreate(PETSC_COMM_WORLD,sf,ierr))
65d8606c27SBarry Smith      PetscCallA(PetscSFSetFromOptions(sf,ierr))
66d8606c27SBarry Smith      PetscCallA(PetscSFSetGraph(sf,nrootsalloc,nleaves,mine,PETSC_COPY_VALUES,remote,PETSC_COPY_VALUES,ierr))
67d8606c27SBarry Smith      PetscCallA(PetscSFSetUp(sf,ierr))
68c4762a1bSJed Brown
69c4762a1bSJed Brown!   View graph, mostly useful for debugging purposes.
70d8606c27SBarry Smith      PetscCallA(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL,ierr))
71d8606c27SBarry Smith      PetscCallA(PetscSFView(sf,PETSC_VIEWER_STDOUT_WORLD,ierr))
72d8606c27SBarry Smith      PetscCallA(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD,ierr))
73c4762a1bSJed Brown
742d4ee042Sprj-!   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.
90d8606c27SBarry Smith      PetscCallA(PetscSFBcastBegin(sf,MPIU_INTEGER,rootdata,leafdata,MPI_REPLACE,ierr))
91d8606c27SBarry Smith      PetscCallA(PetscSFBcastEnd(sf,MPIU_INTEGER,rootdata,leafdata,MPI_REPLACE,ierr))
92dcb3e689SBarry Smith      PetscCallA(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,'## Bcast Rootdata\n',ierr))
93d8606c27SBarry Smith      PetscCallA(PetscIntView(nrootsalloc,rootdata,PETSC_VIEWER_STDOUT_WORLD,ierr))
94dcb3e689SBarry Smith      PetscCallA(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,'## Bcast Leafdata\n',ierr))
95d8606c27SBarry Smith      PetscCallA(PetscIntView(nleavesalloc,leafdata,PETSC_VIEWER_STDOUT_WORLD,ierr))
96c4762a1bSJed Brown
979037d788SNicholas Arnold-Medabalimi!     Reduce entries from leafdata to rootdata. Computation or other communication can be performed between the begin and end calls.
989037d788SNicholas Arnold-Medabalimi      PetscCallA(PetscSFReduceBegin(sf,MPIU_INTEGER,leafdata,rootdata,MPI_SUM,ierr))
999037d788SNicholas Arnold-Medabalimi      PetscCallA(PetscSFReduceEnd(sf,MPIU_INTEGER,leafdata,rootdata,MPI_SUM,ierr))
100dcb3e689SBarry Smith      PetscCallA(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,'## Reduce Leafdata\n',ierr))
1019037d788SNicholas Arnold-Medabalimi      PetscCallA(PetscIntView(nleavesalloc,leafdata,PETSC_VIEWER_STDOUT_WORLD,ierr))
102dcb3e689SBarry Smith      PetscCallA(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,'## Reduce Rootdata\n',ierr))
1039037d788SNicholas Arnold-Medabalimi      PetscCallA(PetscIntView(nrootsalloc,rootdata,PETSC_VIEWER_STDOUT_WORLD,ierr))
1049037d788SNicholas Arnold-Medabalimi
105d8606c27SBarry Smith      PetscCallA(PetscSFGetGraph(sf,gnroots,gnleaves,gmine,gremote,ierr))
106dcb3e689SBarry Smith      PetscCheckA(gnleaves .eq. nleaves,PETSC_COMM_WORLD,PETSC_ERR_PLIB,'nleaves returned from PetscSFGetGraph() does not match that set with PetscSFSetGraph()')
107c4762a1bSJed Brown      do i=1,nleaves
108dcb3e689SBarry Smith        PetscCheckA(gmine(i) .eq. mine(i),PETSC_COMM_WORLD,PETSC_ERR_PLIB,'Root from PetscSFGetGraph() does not match that set with PetscSFSetGraph()')
109c4762a1bSJed Brown      enddo
110c4762a1bSJed Brown      do i=1,nleaves
111dcb3e689SBarry 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()')
112c4762a1bSJed Brown      enddo
113c4762a1bSJed Brown
114c4762a1bSJed Brown      deallocate(gremote)
11594a885e8SJunchao Zhang
11694a885e8SJunchao Zhang! Test PetscSFGet{Leaf,Root}Ranks
11794a885e8SJunchao Zhang      PetscCallA(PetscSFGetLeafRanks(sf,niranks,iranks,ioffset,irootloc,ierr))
11894a885e8SJunchao Zhang      PetscCallA(PetscSFGetRootRanks(sf,nranks,ranks,roffset,rmine,rremote,ierr))
11994a885e8SJunchao Zhang
120c4762a1bSJed Brown!    Clean storage for star forest.
121d8606c27SBarry Smith      PetscCallA(PetscSFDestroy(sf,ierr))
1228dbb0df6SBarry Smith
123da81f932SPierre Jolivet!  Create a star forest with continuous leaves and hence no buffer
124d8606c27SBarry Smith      PetscCallA(PetscSFCreate(PETSC_COMM_WORLD,sf,ierr))
125d8606c27SBarry Smith      PetscCallA(PetscSFSetFromOptions(sf,ierr))
1265d83a8b1SBarry Smith      PetscCallA(PetscSFSetGraph(sf,nrootsalloc,nleaves,PETSC_NULL_INTEGER_ARRAY,PETSC_COPY_VALUES,remote,PETSC_COPY_VALUES,ierr))
127d8606c27SBarry Smith      PetscCallA(PetscSFSetUp(sf,ierr))
1288dbb0df6SBarry Smith
1298dbb0df6SBarry Smith!   View graph, mostly useful for debugging purposes.
130d8606c27SBarry Smith      PetscCallA(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL,ierr))
131d8606c27SBarry Smith      PetscCallA(PetscSFView(sf,PETSC_VIEWER_STDOUT_WORLD,ierr))
132d8606c27SBarry Smith      PetscCallA(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD,ierr))
1338dbb0df6SBarry Smith
134d8606c27SBarry Smith      PetscCallA(PetscSFGetGraph(sf,gnroots,gnleaves,gmine,gremote,ierr))
135dcb3e689SBarry Smith      PetscCheckA(loc(gmine) .eq. loc(PETSC_NULL_INTEGER),PETSC_COMM_WORLD,PETSC_ERR_PLIB,'Leaves from PetscSFGetGraph() not null as expected')
1368dbb0df6SBarry Smith      deallocate(gremote)
137d8606c27SBarry Smith      PetscCallA(PetscSFDestroy(sf,ierr))
138d8606c27SBarry Smith      PetscCallA(PetscFinalize(ierr))
139c4762a1bSJed Brown  end
140c4762a1bSJed Brown
141c4762a1bSJed Brown!/*TEST
142c4762a1bSJed Brown!  build:
143dfd57a17SPierre Jolivet!    requires: defined(PETSC_HAVE_FORTRAN_TYPE_STAR)
144c4762a1bSJed Brown!
145c4762a1bSJed Brown!  test:
146c4762a1bSJed Brown!    nsize: 3
147c4762a1bSJed Brown!
148c4762a1bSJed Brown!TEST*/
149