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