xref: /petsc/src/vec/is/sf/tutorials/ex1f.F90 (revision f14a7c29b82d1117d8e3de344377442be395a55f)
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
8c4762a1bSJed Brown  use petscvec
9c4762a1bSJed Brown  implicit none
10c4762a1bSJed Brown
11c4762a1bSJed Brown  PetscErrorCode ierr
12c4762a1bSJed Brown  PetscInt i, nroots, nrootsalloc, nleaves, nleavesalloc, mine(6), stride
13ce78bad3SBarry Smith  PetscSFNode remote(6)
14c4762a1bSJed Brown  PetscMPIInt rank, size
15c4762a1bSJed Brown  PetscSF sf
16c4762a1bSJed Brown  PetscInt rootdata(6), leafdata(6)
17c4762a1bSJed Brown
18c4762a1bSJed Brown! used with PetscSFGetGraph()
19ce78bad3SBarry Smith  PetscSFNode, pointer ::       gremote(:)
20c4762a1bSJed Brown  PetscInt, pointer ::          gmine(:)
21ccfd86f1SBarry Smith  PetscInt gnroots, gnleaves
22c4762a1bSJed Brown
236497c311SBarry Smith  PetscMPIInt niranks, nranks
2494a885e8SJunchao Zhang  PetscMPIInt, pointer ::       iranks(:), ranks(:)
2594a885e8SJunchao Zhang  PetscInt, pointer ::          ioffset(:), irootloc(:), roffset(:), rmine(:), rremote(:)
2694a885e8SJunchao Zhang
27d8606c27SBarry Smith  PetscCallA(PetscInitialize(ierr))
28c4762a1bSJed Brown  stride = 2
29d8606c27SBarry Smith  PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD, rank, ierr))
30d8606c27SBarry Smith  PetscCallMPIA(MPI_Comm_size(PETSC_COMM_WORLD, size, ierr))
31c4762a1bSJed Brown
32c4762a1bSJed Brown  if (rank == 0) then
33c4762a1bSJed Brown    nroots = 3
34c4762a1bSJed Brown  else
35c4762a1bSJed Brown    nroots = 2
36c4762a1bSJed Brown  end if
37ccfd86f1SBarry Smith  nrootsalloc = nroots*stride
38c4762a1bSJed Brown  if (rank > 0) then
39c4762a1bSJed Brown    nleaves = 3
40c4762a1bSJed Brown  else
41c4762a1bSJed Brown    nleaves = 2
42c4762a1bSJed Brown  end if
43c4762a1bSJed Brown  nleavesalloc = nleaves*stride
44c4762a1bSJed Brown  if (stride > 1) then
45c4762a1bSJed Brown    do i = 1, nleaves
46c4762a1bSJed Brown      mine(i) = stride*(i - 1)
47c4762a1bSJed Brown    end do
48c4762a1bSJed Brown  end if
49c4762a1bSJed Brown
50c4762a1bSJed Brown! Left periodic neighbor
51c4762a1bSJed Brown  remote(1)%rank = modulo(rank + size - 1, size)
52c4762a1bSJed Brown  remote(1)%index = 1*stride
53c4762a1bSJed Brown! Right periodic neighbor
54c4762a1bSJed Brown  remote(2)%rank = modulo(rank + 1, size)
55c4762a1bSJed Brown  remote(2)%index = 0*stride
56c4762a1bSJed Brown  if (rank > 0) then !               All processes reference rank 0, index
57c4762a1bSJed Brown    remote(3)%rank = 0
58c4762a1bSJed Brown    remote(3)%index = 2*stride
59c4762a1bSJed Brown  end if
60c4762a1bSJed Brown
61c4762a1bSJed Brown!  Create a star forest for communication
62d8606c27SBarry Smith  PetscCallA(PetscSFCreate(PETSC_COMM_WORLD, sf, ierr))
63d8606c27SBarry Smith  PetscCallA(PetscSFSetFromOptions(sf, ierr))
64d8606c27SBarry Smith  PetscCallA(PetscSFSetGraph(sf, nrootsalloc, nleaves, mine, PETSC_COPY_VALUES, remote, PETSC_COPY_VALUES, ierr))
65d8606c27SBarry Smith  PetscCallA(PetscSFSetUp(sf, ierr))
66c4762a1bSJed Brown
67c4762a1bSJed Brown!   View graph, mostly useful for debugging purposes.
68d8606c27SBarry Smith  PetscCallA(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_INFO_DETAIL, ierr))
69d8606c27SBarry Smith  PetscCallA(PetscSFView(sf, PETSC_VIEWER_STDOUT_WORLD, ierr))
70d8606c27SBarry Smith  PetscCallA(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD, ierr))
71c4762a1bSJed Brown
722d4ee042Sprj-!   Allocate space for send and receive buffers. This example communicates PetscInt, but other types, including
73c4762a1bSJed Brown!     * user-defined structures, could also be used.
74c4762a1bSJed Brown!     Set rootdata buffer to be broadcast
75c4762a1bSJed Brown  do i = 1, nrootsalloc
76c4762a1bSJed Brown    rootdata(i) = -1
77c4762a1bSJed Brown  end do
78c4762a1bSJed Brown  do i = 1, nroots
79c4762a1bSJed Brown    rootdata(1 + (i - 1)*stride) = 100*(rank + 1) + i - 1
80c4762a1bSJed Brown  end do
81c4762a1bSJed Brown
82c4762a1bSJed Brown!     Initialize local buffer, these values are never used.
83c4762a1bSJed Brown  do i = 1, nleavesalloc
84c4762a1bSJed Brown    leafdata(i) = -1
85c4762a1bSJed Brown  end do
86c4762a1bSJed Brown
87c4762a1bSJed Brown!     Broadcast entries from rootdata to leafdata. Computation or other communication can be performed between the begin and end calls.
88d8606c27SBarry Smith  PetscCallA(PetscSFBcastBegin(sf, MPIU_INTEGER, rootdata, leafdata, MPI_REPLACE, ierr))
89d8606c27SBarry Smith  PetscCallA(PetscSFBcastEnd(sf, MPIU_INTEGER, rootdata, leafdata, MPI_REPLACE, ierr))
90dcb3e689SBarry Smith  PetscCallA(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, '## Bcast Rootdata\n', ierr))
91d8606c27SBarry Smith  PetscCallA(PetscIntView(nrootsalloc, rootdata, PETSC_VIEWER_STDOUT_WORLD, ierr))
92dcb3e689SBarry Smith  PetscCallA(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, '## Bcast Leafdata\n', ierr))
93d8606c27SBarry Smith  PetscCallA(PetscIntView(nleavesalloc, leafdata, PETSC_VIEWER_STDOUT_WORLD, ierr))
94c4762a1bSJed Brown
959037d788SNicholas Arnold-Medabalimi!     Reduce entries from leafdata to rootdata. Computation or other communication can be performed between the begin and end calls.
969037d788SNicholas Arnold-Medabalimi  PetscCallA(PetscSFReduceBegin(sf, MPIU_INTEGER, leafdata, rootdata, MPI_SUM, ierr))
979037d788SNicholas Arnold-Medabalimi  PetscCallA(PetscSFReduceEnd(sf, MPIU_INTEGER, leafdata, rootdata, MPI_SUM, ierr))
98dcb3e689SBarry Smith  PetscCallA(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, '## Reduce Leafdata\n', ierr))
999037d788SNicholas Arnold-Medabalimi  PetscCallA(PetscIntView(nleavesalloc, leafdata, PETSC_VIEWER_STDOUT_WORLD, ierr))
100dcb3e689SBarry Smith  PetscCallA(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, '## Reduce Rootdata\n', ierr))
1019037d788SNicholas Arnold-Medabalimi  PetscCallA(PetscIntView(nrootsalloc, rootdata, PETSC_VIEWER_STDOUT_WORLD, ierr))
1029037d788SNicholas Arnold-Medabalimi
103d8606c27SBarry Smith  PetscCallA(PetscSFGetGraph(sf, gnroots, gnleaves, gmine, gremote, ierr))
1044820e4eaSBarry Smith  PetscCheckA(gnleaves == nleaves, PETSC_COMM_WORLD, PETSC_ERR_PLIB, 'nleaves returned from PetscSFGetGraph() does not match that set with PetscSFSetGraph()')
105c4762a1bSJed Brown  do i = 1, nleaves
1064820e4eaSBarry Smith    PetscCheckA(gmine(i) == mine(i), PETSC_COMM_WORLD, PETSC_ERR_PLIB, 'Root from PetscSFGetGraph() does not match that set with PetscSFSetGraph()')
107c4762a1bSJed Brown  end do
108c4762a1bSJed Brown  do i = 1, nleaves
1094820e4eaSBarry Smith    PetscCheckA(gremote(i)%index == remote(i)%index, PETSC_COMM_WORLD, PETSC_ERR_PLIB, 'Leaf from PetscSFGetGraph() does not match that set with PetscSFSetGraph()')
110c4762a1bSJed Brown  end do
111ce78bad3SBarry Smith  PetscCallA(PetscSFRestoreGraph(sf, gnroots, gnleaves, gmine, gremote, ierr))
11294a885e8SJunchao Zhang
11394a885e8SJunchao Zhang! Test PetscSFGet{Leaf,Root}Ranks
11494a885e8SJunchao Zhang  PetscCallA(PetscSFGetLeafRanks(sf, niranks, iranks, ioffset, irootloc, ierr))
11594a885e8SJunchao Zhang  PetscCallA(PetscSFGetRootRanks(sf, nranks, ranks, roffset, rmine, rremote, ierr))
11694a885e8SJunchao Zhang
117c4762a1bSJed Brown!    Clean storage for star forest.
118d8606c27SBarry Smith  PetscCallA(PetscSFDestroy(sf, ierr))
1198dbb0df6SBarry Smith
120da81f932SPierre Jolivet!  Create a star forest with continuous leaves and hence no buffer
121d8606c27SBarry Smith  PetscCallA(PetscSFCreate(PETSC_COMM_WORLD, sf, ierr))
122d8606c27SBarry Smith  PetscCallA(PetscSFSetFromOptions(sf, ierr))
1235d83a8b1SBarry Smith  PetscCallA(PetscSFSetGraph(sf, nrootsalloc, nleaves, PETSC_NULL_INTEGER_ARRAY, PETSC_COPY_VALUES, remote, PETSC_COPY_VALUES, ierr))
124d8606c27SBarry Smith  PetscCallA(PetscSFSetUp(sf, ierr))
1258dbb0df6SBarry Smith
1268dbb0df6SBarry Smith!   View graph, mostly useful for debugging purposes.
127d8606c27SBarry Smith  PetscCallA(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_INFO_DETAIL, ierr))
128d8606c27SBarry Smith  PetscCallA(PetscSFView(sf, PETSC_VIEWER_STDOUT_WORLD, ierr))
129d8606c27SBarry Smith  PetscCallA(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD, ierr))
1308dbb0df6SBarry Smith
131d8606c27SBarry Smith  PetscCallA(PetscSFGetGraph(sf, gnroots, gnleaves, gmine, gremote, ierr))
1324820e4eaSBarry Smith  PetscCheckA(loc(gmine) == loc(PETSC_NULL_INTEGER), PETSC_COMM_WORLD, PETSC_ERR_PLIB, 'Leaves from PetscSFGetGraph() not null as expected')
133ce78bad3SBarry Smith  PetscCallA(PetscSFRestoreGraph(sf, gnroots, gnleaves, gmine, gremote, ierr))
134d8606c27SBarry Smith  PetscCallA(PetscSFDestroy(sf, ierr))
135d8606c27SBarry Smith  PetscCallA(PetscFinalize(ierr))
136c4762a1bSJed Brownend
137c4762a1bSJed Brown
138c4762a1bSJed Brown!/*TEST
139c4762a1bSJed Brown!
140c4762a1bSJed Brown!  test:
141c4762a1bSJed Brown!    nsize: 3
142c4762a1bSJed Brown!
143c4762a1bSJed Brown!TEST*/
144