xref: /petsc/src/vec/is/sf/tutorials/ex1f.F90 (revision 607e733f3db3ee7f6f605a13295c517df8dbb9c9)
1!    Description: A star forest is a simple tree with one root and zero or more leaves.
2!    Many common communication patterns can be expressed as updates of rootdata using leafdata and vice-versa.
3!     This example creates a star forest, communicates values using the graph  views the graph, then destroys it.
4!
5!     This is a copy of ex1.c but currently only tests the broadcast operation
6#include <petsc/finclude/petscvec.h>
7program main
8  use petscvec
9  implicit none
10
11  PetscErrorCode ierr
12  PetscInt i, nroots, nrootsalloc, nleaves, nleavesalloc, mine(6), stride
13  PetscSFNode remote(6)
14  PetscMPIInt rank, size
15  PetscSF sf
16  PetscInt rootdata(6), leafdata(6)
17
18! used with PetscSFGetGraph()
19  PetscSFNode, pointer ::       gremote(:)
20  PetscInt, pointer ::          gmine(:)
21  PetscInt gnroots, gnleaves
22
23  PetscMPIInt niranks, nranks
24  PetscMPIInt, pointer ::       iranks(:), ranks(:)
25  PetscInt, pointer ::          ioffset(:), irootloc(:), roffset(:), rmine(:), rremote(:)
26
27  PetscCallA(PetscInitialize(ierr))
28  stride = 2
29  PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD, rank, ierr))
30  PetscCallMPIA(MPI_Comm_size(PETSC_COMM_WORLD, size, ierr))
31
32  if (rank == 0) then
33    nroots = 3
34  else
35    nroots = 2
36  end if
37  nrootsalloc = nroots*stride
38  if (rank > 0) then
39    nleaves = 3
40  else
41    nleaves = 2
42  end if
43  nleavesalloc = nleaves*stride
44  if (stride > 1) then
45    do i = 1, nleaves
46      mine(i) = stride*(i - 1)
47    end do
48  end if
49
50! Left periodic neighbor
51  remote(1)%rank = modulo(rank + size - 1, size)
52  remote(1)%index = 1*stride
53! Right periodic neighbor
54  remote(2)%rank = modulo(rank + 1, size)
55  remote(2)%index = 0*stride
56  if (rank > 0) then !               All processes reference rank 0, index
57    remote(3)%rank = 0
58    remote(3)%index = 2*stride
59  end if
60
61!  Create a star forest for communication
62  PetscCallA(PetscSFCreate(PETSC_COMM_WORLD, sf, ierr))
63  PetscCallA(PetscSFSetFromOptions(sf, ierr))
64  PetscCallA(PetscSFSetGraph(sf, nrootsalloc, nleaves, mine, PETSC_COPY_VALUES, remote, PETSC_COPY_VALUES, ierr))
65  PetscCallA(PetscSFSetUp(sf, ierr))
66
67!   View graph, mostly useful for debugging purposes.
68  PetscCallA(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_INFO_DETAIL, ierr))
69  PetscCallA(PetscSFView(sf, PETSC_VIEWER_STDOUT_WORLD, ierr))
70  PetscCallA(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD, ierr))
71
72!   Allocate space for send and receive buffers. This example communicates PetscInt, but other types, including
73!     * user-defined structures, could also be used.
74!     Set rootdata buffer to be broadcast
75  do i = 1, nrootsalloc
76    rootdata(i) = -1
77  end do
78  do i = 1, nroots
79    rootdata(1 + (i - 1)*stride) = 100*(rank + 1) + i - 1
80  end do
81
82!     Initialize local buffer, these values are never used.
83  do i = 1, nleavesalloc
84    leafdata(i) = -1
85  end do
86
87!     Broadcast entries from rootdata to leafdata. Computation or other communication can be performed between the begin and end calls.
88  PetscCallA(PetscSFBcastBegin(sf, MPIU_INTEGER, rootdata, leafdata, MPI_REPLACE, ierr))
89  PetscCallA(PetscSFBcastEnd(sf, MPIU_INTEGER, rootdata, leafdata, MPI_REPLACE, ierr))
90  PetscCallA(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, '## Bcast Rootdata\n', ierr))
91  PetscCallA(PetscIntView(nrootsalloc, rootdata, PETSC_VIEWER_STDOUT_WORLD, ierr))
92  PetscCallA(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, '## Bcast Leafdata\n', ierr))
93  PetscCallA(PetscIntView(nleavesalloc, leafdata, PETSC_VIEWER_STDOUT_WORLD, ierr))
94
95!     Reduce entries from leafdata to rootdata. Computation or other communication can be performed between the begin and end calls.
96  PetscCallA(PetscSFReduceBegin(sf, MPIU_INTEGER, leafdata, rootdata, MPI_SUM, ierr))
97  PetscCallA(PetscSFReduceEnd(sf, MPIU_INTEGER, leafdata, rootdata, MPI_SUM, ierr))
98  PetscCallA(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, '## Reduce Leafdata\n', ierr))
99  PetscCallA(PetscIntView(nleavesalloc, leafdata, PETSC_VIEWER_STDOUT_WORLD, ierr))
100  PetscCallA(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, '## Reduce Rootdata\n', ierr))
101  PetscCallA(PetscIntView(nrootsalloc, rootdata, PETSC_VIEWER_STDOUT_WORLD, ierr))
102
103  PetscCallA(PetscSFGetGraph(sf, gnroots, gnleaves, gmine, gremote, ierr))
104  PetscCheckA(gnleaves == nleaves, PETSC_COMM_WORLD, PETSC_ERR_PLIB, 'nleaves returned from PetscSFGetGraph() does not match that set with PetscSFSetGraph()')
105  do i = 1, nleaves
106    PetscCheckA(gmine(i) == mine(i), PETSC_COMM_WORLD, PETSC_ERR_PLIB, 'Root from PetscSFGetGraph() does not match that set with PetscSFSetGraph()')
107  end do
108  do i = 1, nleaves
109    PetscCheckA(gremote(i)%index == remote(i)%index, PETSC_COMM_WORLD, PETSC_ERR_PLIB, 'Leaf from PetscSFGetGraph() does not match that set with PetscSFSetGraph()')
110  end do
111  PetscCallA(PetscSFRestoreGraph(sf, gnroots, gnleaves, gmine, gremote, ierr))
112
113! Test PetscSFGet{Leaf,Root}Ranks
114  PetscCallA(PetscSFGetLeafRanks(sf, niranks, iranks, ioffset, irootloc, ierr))
115  PetscCallA(PetscSFGetRootRanks(sf, nranks, ranks, roffset, rmine, rremote, ierr))
116
117!    Clean storage for star forest.
118  PetscCallA(PetscSFDestroy(sf, ierr))
119
120!  Create a star forest with continuous leaves and hence no buffer
121  PetscCallA(PetscSFCreate(PETSC_COMM_WORLD, sf, ierr))
122  PetscCallA(PetscSFSetFromOptions(sf, ierr))
123  PetscCallA(PetscSFSetGraph(sf, nrootsalloc, nleaves, PETSC_NULL_INTEGER_ARRAY, PETSC_COPY_VALUES, remote, PETSC_COPY_VALUES, ierr))
124  PetscCallA(PetscSFSetUp(sf, ierr))
125
126!   View graph, mostly useful for debugging purposes.
127  PetscCallA(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_INFO_DETAIL, ierr))
128  PetscCallA(PetscSFView(sf, PETSC_VIEWER_STDOUT_WORLD, ierr))
129  PetscCallA(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD, ierr))
130
131  PetscCallA(PetscSFGetGraph(sf, gnroots, gnleaves, gmine, gremote, ierr))
132  PetscCheckA(loc(gmine) == loc(PETSC_NULL_INTEGER), PETSC_COMM_WORLD, PETSC_ERR_PLIB, 'Leaves from PetscSFGetGraph() not null as expected')
133  PetscCallA(PetscSFRestoreGraph(sf, gnroots, gnleaves, gmine, gremote, ierr))
134  PetscCallA(PetscSFDestroy(sf, ierr))
135  PetscCallA(PetscFinalize(ierr))
136end
137
138!/*TEST
139!
140!  test:
141!    nsize: 3
142!
143!TEST*/
144