xref: /petsc/src/vec/is/sf/tutorials/ex1.c (revision a69119a591a03a9d906b29c0a4e9802e4d7c9795)
1 static const char help[] = "Test star forest communication (PetscSF)\n\n";
2 
3 /*
4     Description: A star is a simple tree with one root and zero or more leaves.
5     A star forest is a union of disjoint stars.
6     Many common communication patterns can be expressed as updates of rootdata using leafdata and vice-versa.
7     This example creates a star forest, communicates values using the graph (see options for types of communication), views the graph, then destroys it.
8 */
9 
10 /*
11   Include petscsf.h so we can use PetscSF objects. Note that this automatically
12   includes petscsys.h.
13 */
14 #include <petscsf.h>
15 #include <petscviewer.h>
16 
17 /* like PetscSFView() but with alternative array of local indices */
18 static PetscErrorCode PetscSFViewCustomLocals_Private(PetscSF sf, const PetscInt locals[], PetscViewer viewer) {
19   const PetscSFNode *iremote;
20   PetscInt           i, nroots, nleaves, nranks;
21   PetscMPIInt        rank;
22 
23   PetscFunctionBeginUser;
24   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sf), &rank));
25   PetscCall(PetscSFGetGraph(sf, &nroots, &nleaves, NULL, &iremote));
26   PetscCall(PetscSFGetRootRanks(sf, &nranks, NULL, NULL, NULL, NULL));
27   PetscCall(PetscViewerASCIIPushTab(viewer));
28   PetscCall(PetscViewerASCIIPushSynchronized(viewer));
29   PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] Number of roots=%" PetscInt_FMT ", leaves=%" PetscInt_FMT ", remote ranks=%" PetscInt_FMT "\n", rank, nroots, nleaves, nranks));
30   for (i = 0; i < nleaves; i++) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] %" PetscInt_FMT " <- (%" PetscInt_FMT ",%" PetscInt_FMT ")\n", rank, locals[i], iremote[i].rank, iremote[i].index));
31   PetscCall(PetscViewerFlush(viewer));
32   PetscCall(PetscViewerASCIIPopTab(viewer));
33   PetscCall(PetscViewerASCIIPopSynchronized(viewer));
34   PetscFunctionReturn(0);
35 }
36 
37 int main(int argc, char **argv) {
38   PetscInt     i, nroots, nrootsalloc, nleaves, nleavesalloc, *mine, stride;
39   PetscSFNode *remote;
40   PetscMPIInt  rank, size;
41   PetscSF      sf;
42   PetscBool    test_all, test_bcast, test_bcastop, test_reduce, test_degree, test_fetchandop, test_gather, test_scatter, test_embed, test_invert, test_sf_distribute, test_char;
43   MPI_Op       mop = MPI_OP_NULL; /* initialize to prevent compiler warnings with cxx_quad build */
44   char         opstring[256];
45   PetscBool    strflg;
46 
47   PetscFunctionBeginUser;
48   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
49   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
50   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
51 
52   PetscOptionsBegin(PETSC_COMM_WORLD, "", "PetscSF Test Options", "none");
53   test_all = PETSC_FALSE;
54   PetscCall(PetscOptionsBool("-test_all", "Test all SF communications", "", test_all, &test_all, NULL));
55   test_bcast = test_all;
56   PetscCall(PetscOptionsBool("-test_bcast", "Test broadcast", "", test_bcast, &test_bcast, NULL));
57   test_bcastop = test_all;
58   PetscCall(PetscOptionsBool("-test_bcastop", "Test broadcast and reduce", "", test_bcastop, &test_bcastop, NULL));
59   test_reduce = test_all;
60   PetscCall(PetscOptionsBool("-test_reduce", "Test reduction", "", test_reduce, &test_reduce, NULL));
61   test_char = test_all;
62   PetscCall(PetscOptionsBool("-test_char", "Test signed char, unsigned char, and char", "", test_char, &test_char, NULL));
63   mop = MPI_SUM;
64   PetscCall(PetscStrcpy(opstring, "sum"));
65   PetscCall(PetscOptionsString("-test_op", "Designate which MPI_Op to use", "", opstring, opstring, sizeof(opstring), NULL));
66   PetscCall(PetscStrcmp("sum", opstring, &strflg));
67   if (strflg) mop = MPIU_SUM;
68   PetscCall(PetscStrcmp("prod", opstring, &strflg));
69   if (strflg) mop = MPI_PROD;
70   PetscCall(PetscStrcmp("max", opstring, &strflg));
71   if (strflg) mop = MPI_MAX;
72   PetscCall(PetscStrcmp("min", opstring, &strflg));
73   if (strflg) mop = MPI_MIN;
74   PetscCall(PetscStrcmp("land", opstring, &strflg));
75   if (strflg) mop = MPI_LAND;
76   PetscCall(PetscStrcmp("band", opstring, &strflg));
77   if (strflg) mop = MPI_BAND;
78   PetscCall(PetscStrcmp("lor", opstring, &strflg));
79   if (strflg) mop = MPI_LOR;
80   PetscCall(PetscStrcmp("bor", opstring, &strflg));
81   if (strflg) mop = MPI_BOR;
82   PetscCall(PetscStrcmp("lxor", opstring, &strflg));
83   if (strflg) mop = MPI_LXOR;
84   PetscCall(PetscStrcmp("bxor", opstring, &strflg));
85   if (strflg) mop = MPI_BXOR;
86   test_degree = test_all;
87   PetscCall(PetscOptionsBool("-test_degree", "Test computation of vertex degree", "", test_degree, &test_degree, NULL));
88   test_fetchandop = test_all;
89   PetscCall(PetscOptionsBool("-test_fetchandop", "Test atomic Fetch-And-Op", "", test_fetchandop, &test_fetchandop, NULL));
90   test_gather = test_all;
91   PetscCall(PetscOptionsBool("-test_gather", "Test point gather", "", test_gather, &test_gather, NULL));
92   test_scatter = test_all;
93   PetscCall(PetscOptionsBool("-test_scatter", "Test point scatter", "", test_scatter, &test_scatter, NULL));
94   test_embed = test_all;
95   PetscCall(PetscOptionsBool("-test_embed", "Test point embed", "", test_embed, &test_embed, NULL));
96   test_invert = test_all;
97   PetscCall(PetscOptionsBool("-test_invert", "Test point invert", "", test_invert, &test_invert, NULL));
98   stride = 1;
99   PetscCall(PetscOptionsInt("-stride", "Stride for leaf and root data", "", stride, &stride, NULL));
100   test_sf_distribute = PETSC_FALSE;
101   PetscCall(PetscOptionsBool("-test_sf_distribute", "Create an SF that 'distributes' to each process, like an alltoall", "", test_sf_distribute, &test_sf_distribute, NULL));
102   PetscCall(PetscOptionsString("-test_op", "Designate which MPI_Op to use", "", opstring, opstring, sizeof(opstring), NULL));
103   PetscOptionsEnd();
104 
105   if (test_sf_distribute) {
106     nroots       = size;
107     nrootsalloc  = size;
108     nleaves      = size;
109     nleavesalloc = size;
110     mine         = NULL;
111     PetscCall(PetscMalloc1(nleaves, &remote));
112     for (i = 0; i < size; i++) {
113       remote[i].rank  = i;
114       remote[i].index = rank;
115     }
116   } else {
117     nroots       = 2 + (PetscInt)(rank == 0);
118     nrootsalloc  = nroots * stride;
119     nleaves      = 2 + (PetscInt)(rank > 0);
120     nleavesalloc = nleaves * stride;
121     mine         = NULL;
122     if (stride > 1) {
123       PetscInt i;
124 
125       PetscCall(PetscMalloc1(nleaves, &mine));
126       for (i = 0; i < nleaves; i++) mine[i] = stride * i;
127     }
128     PetscCall(PetscMalloc1(nleaves, &remote));
129     /* Left periodic neighbor */
130     remote[0].rank  = (rank + size - 1) % size;
131     remote[0].index = 1 * stride;
132     /* Right periodic neighbor */
133     remote[1].rank  = (rank + 1) % size;
134     remote[1].index = 0 * stride;
135     if (rank > 0) { /* All processes reference rank 0, index 1 */
136       remote[2].rank  = 0;
137       remote[2].index = 2 * stride;
138     }
139   }
140 
141   /* Create a star forest for communication. In this example, the leaf space is dense, so we pass NULL. */
142   PetscCall(PetscSFCreate(PETSC_COMM_WORLD, &sf));
143   PetscCall(PetscSFSetFromOptions(sf));
144   PetscCall(PetscSFSetGraph(sf, nrootsalloc, nleaves, mine, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER));
145   PetscCall(PetscSFSetUp(sf));
146 
147   /* View graph, mostly useful for debugging purposes. */
148   PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_INFO_DETAIL));
149   PetscCall(PetscSFView(sf, PETSC_VIEWER_STDOUT_WORLD));
150   PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
151 
152   if (test_bcast) { /* broadcast rootdata into leafdata */
153     PetscInt *rootdata, *leafdata;
154     /* Allocate space for send and receive buffers. This example communicates PetscInt, but other types, including
155      * user-defined structures, could also be used. */
156     PetscCall(PetscMalloc2(nrootsalloc, &rootdata, nleavesalloc, &leafdata));
157     /* Set rootdata buffer to be broadcast */
158     for (i = 0; i < nrootsalloc; i++) rootdata[i] = -1;
159     for (i = 0; i < nroots; i++) rootdata[i * stride] = 100 * (rank + 1) + i;
160     /* Initialize local buffer, these values are never used. */
161     for (i = 0; i < nleavesalloc; i++) leafdata[i] = -1;
162     /* Broadcast entries from rootdata to leafdata. Computation or other communication can be performed between the begin and end calls. */
163     PetscCall(PetscSFBcastBegin(sf, MPIU_INT, rootdata, leafdata, MPI_REPLACE));
164     PetscCall(PetscSFBcastEnd(sf, MPIU_INT, rootdata, leafdata, MPI_REPLACE));
165     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Bcast Rootdata\n"));
166     PetscCall(PetscIntView(nrootsalloc, rootdata, PETSC_VIEWER_STDOUT_WORLD));
167     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Bcast Leafdata\n"));
168     PetscCall(PetscIntView(nleavesalloc, leafdata, PETSC_VIEWER_STDOUT_WORLD));
169     PetscCall(PetscFree2(rootdata, leafdata));
170   }
171 
172   if (test_bcast && test_char) { /* Bcast with char */
173     PetscInt len;
174     char     buf[256];
175     char    *rootdata, *leafdata;
176     PetscCall(PetscMalloc2(nrootsalloc, &rootdata, nleavesalloc, &leafdata));
177     /* Set rootdata buffer to be broadcast */
178     for (i = 0; i < nrootsalloc; i++) rootdata[i] = '*';
179     for (i = 0; i < nroots; i++) rootdata[i * stride] = 'A' + rank * 3 + i; /* rank is very small, so it is fine to compute a char */
180     /* Initialize local buffer, these values are never used. */
181     for (i = 0; i < nleavesalloc; i++) leafdata[i] = '?';
182 
183     PetscCall(PetscSFBcastBegin(sf, MPI_CHAR, rootdata, leafdata, MPI_REPLACE));
184     PetscCall(PetscSFBcastEnd(sf, MPI_CHAR, rootdata, leafdata, MPI_REPLACE));
185 
186     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Bcast Rootdata in type of char\n"));
187     len = 0;
188     PetscCall(PetscSNPrintf(buf, 256, "%4d:", rank));
189     len += 5;
190     for (i = 0; i < nrootsalloc; i++) {
191       PetscCall(PetscSNPrintf(buf + len, 256 - len, "%5c", rootdata[i]));
192       len += 5;
193     }
194     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s\n", buf));
195     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT));
196 
197     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Bcast Leafdata in type of char\n"));
198     len = 0;
199     PetscCall(PetscSNPrintf(buf, 256, "%4d:", rank));
200     len += 5;
201     for (i = 0; i < nleavesalloc; i++) {
202       PetscCall(PetscSNPrintf(buf + len, 256 - len, "%5c", leafdata[i]));
203       len += 5;
204     }
205     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s\n", buf));
206     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT));
207 
208     PetscCall(PetscFree2(rootdata, leafdata));
209   }
210 
211   if (test_bcastop) { /* Reduce rootdata into leafdata */
212     PetscInt *rootdata, *leafdata;
213     /* Allocate space for send and receive buffers. This example communicates PetscInt, but other types, including
214      * user-defined structures, could also be used. */
215     PetscCall(PetscMalloc2(nrootsalloc, &rootdata, nleavesalloc, &leafdata));
216     /* Set rootdata buffer to be broadcast */
217     for (i = 0; i < nrootsalloc; i++) rootdata[i] = -1;
218     for (i = 0; i < nroots; i++) rootdata[i * stride] = 100 * (rank + 1) + i;
219     /* Set leaf values to reduce with */
220     for (i = 0; i < nleavesalloc; i++) leafdata[i] = -10 * (rank + 1) - i;
221     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Pre-BcastAndOp Leafdata\n"));
222     PetscCall(PetscIntView(nleavesalloc, leafdata, PETSC_VIEWER_STDOUT_WORLD));
223     /* Broadcast entries from rootdata to leafdata. Computation or other communication can be performed between the begin and end calls. */
224     PetscCall(PetscSFBcastBegin(sf, MPIU_INT, rootdata, leafdata, mop));
225     PetscCall(PetscSFBcastEnd(sf, MPIU_INT, rootdata, leafdata, mop));
226     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## BcastAndOp Rootdata\n"));
227     PetscCall(PetscIntView(nrootsalloc, rootdata, PETSC_VIEWER_STDOUT_WORLD));
228     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## BcastAndOp Leafdata\n"));
229     PetscCall(PetscIntView(nleavesalloc, leafdata, PETSC_VIEWER_STDOUT_WORLD));
230     PetscCall(PetscFree2(rootdata, leafdata));
231   }
232 
233   if (test_reduce) { /* Reduce leafdata into rootdata */
234     PetscInt *rootdata, *leafdata;
235     PetscCall(PetscMalloc2(nrootsalloc, &rootdata, nleavesalloc, &leafdata));
236     /* Initialize rootdata buffer in which the result of the reduction will appear. */
237     for (i = 0; i < nrootsalloc; i++) rootdata[i] = -1;
238     for (i = 0; i < nroots; i++) rootdata[i * stride] = 100 * (rank + 1) + i;
239     /* Set leaf values to reduce. */
240     for (i = 0; i < nleavesalloc; i++) leafdata[i] = -1;
241     for (i = 0; i < nleaves; i++) leafdata[i * stride] = 1000 * (rank + 1) + 10 * i;
242     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Pre-Reduce Rootdata\n"));
243     PetscCall(PetscIntView(nrootsalloc, rootdata, PETSC_VIEWER_STDOUT_WORLD));
244     /* Perform reduction. Computation or other communication can be performed between the begin and end calls.
245      * This example sums the values, but other MPI_Ops can be used (e.g MPI_MAX, MPI_PROD). */
246     PetscCall(PetscSFReduceBegin(sf, MPIU_INT, leafdata, rootdata, mop));
247     PetscCall(PetscSFReduceEnd(sf, MPIU_INT, leafdata, rootdata, mop));
248     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Reduce Leafdata\n"));
249     PetscCall(PetscIntView(nleavesalloc, leafdata, PETSC_VIEWER_STDOUT_WORLD));
250     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Reduce Rootdata\n"));
251     PetscCall(PetscIntView(nrootsalloc, rootdata, PETSC_VIEWER_STDOUT_WORLD));
252     PetscCall(PetscFree2(rootdata, leafdata));
253   }
254 
255   if (test_reduce && test_char) { /* Reduce with signed char */
256     PetscInt     len;
257     char         buf[256];
258     signed char *rootdata, *leafdata;
259     PetscCall(PetscMalloc2(nrootsalloc, &rootdata, nleavesalloc, &leafdata));
260     /* Initialize rootdata buffer in which the result of the reduction will appear. */
261     for (i = 0; i < nrootsalloc; i++) rootdata[i] = -1;
262     for (i = 0; i < nroots; i++) rootdata[i * stride] = 10 * (rank + 1) + i;
263     /* Set leaf values to reduce. */
264     for (i = 0; i < nleavesalloc; i++) leafdata[i] = -1;
265     for (i = 0; i < nleaves; i++) leafdata[i * stride] = 50 * (rank + 1) + 10 * i;
266     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Pre-Reduce Rootdata in type of signed char\n"));
267 
268     len = 0;
269     PetscCall(PetscSNPrintf(buf, 256, "%4d:", rank));
270     len += 5;
271     for (i = 0; i < nrootsalloc; i++) {
272       PetscCall(PetscSNPrintf(buf + len, 256 - len, "%5d", rootdata[i]));
273       len += 5;
274     }
275     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s\n", buf));
276     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT));
277 
278     /* Using MPI_CHAR should trigger an error since MPI standard does not support reduction on MPI_CHAR.
279        Testing with -test_op max, one can see the sign does take effect in MPI_MAX.
280      */
281     PetscCall(PetscSFReduceBegin(sf, MPI_SIGNED_CHAR, leafdata, rootdata, mop));
282     PetscCall(PetscSFReduceEnd(sf, MPI_SIGNED_CHAR, leafdata, rootdata, mop));
283 
284     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Reduce Leafdata in type of signed char\n"));
285     len = 0;
286     PetscCall(PetscSNPrintf(buf, 256, "%4d:", rank));
287     len += 5;
288     for (i = 0; i < nleavesalloc; i++) {
289       PetscCall(PetscSNPrintf(buf + len, 256 - len, "%5d", leafdata[i]));
290       len += 5;
291     }
292     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s\n", buf));
293     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT));
294 
295     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Reduce Rootdata in type of signed char\n"));
296     len = 0;
297     PetscCall(PetscSNPrintf(buf, 256, "%4d:", rank));
298     len += 5;
299     for (i = 0; i < nrootsalloc; i++) {
300       PetscCall(PetscSNPrintf(buf + len, 256 - len, "%5d", rootdata[i]));
301       len += 5;
302     }
303     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s\n", buf));
304     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT));
305 
306     PetscCall(PetscFree2(rootdata, leafdata));
307   }
308 
309   if (test_reduce && test_char) { /* Reduce with unsigned char */
310     PetscInt       len;
311     char           buf[256];
312     unsigned char *rootdata, *leafdata;
313     PetscCall(PetscMalloc2(nrootsalloc, &rootdata, nleavesalloc, &leafdata));
314     /* Initialize rootdata buffer in which the result of the reduction will appear. */
315     for (i = 0; i < nrootsalloc; i++) rootdata[i] = 0;
316     for (i = 0; i < nroots; i++) rootdata[i * stride] = 10 * (rank + 1) + i;
317     /* Set leaf values to reduce. */
318     for (i = 0; i < nleavesalloc; i++) leafdata[i] = 0;
319     for (i = 0; i < nleaves; i++) leafdata[i * stride] = 50 * (rank + 1) + 10 * i;
320     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Pre-Reduce Rootdata in type of unsigned char\n"));
321 
322     len = 0;
323     PetscCall(PetscSNPrintf(buf, 256, "%4d:", rank));
324     len += 5;
325     for (i = 0; i < nrootsalloc; i++) {
326       PetscCall(PetscSNPrintf(buf + len, 256 - len, "%5u", rootdata[i]));
327       len += 5;
328     }
329     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s\n", buf));
330     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT));
331 
332     /* Using MPI_CHAR should trigger an error since MPI standard does not support reduction on MPI_CHAR.
333        Testing with -test_op max, one can see the sign does take effect in MPI_MAX.
334      */
335     PetscCall(PetscSFReduceBegin(sf, MPI_UNSIGNED_CHAR, leafdata, rootdata, mop));
336     PetscCall(PetscSFReduceEnd(sf, MPI_UNSIGNED_CHAR, leafdata, rootdata, mop));
337 
338     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Reduce Leafdata in type of unsigned char\n"));
339     len = 0;
340     PetscCall(PetscSNPrintf(buf, 256, "%4d:", rank));
341     len += 5;
342     for (i = 0; i < nleavesalloc; i++) {
343       PetscCall(PetscSNPrintf(buf + len, 256 - len, "%5u", leafdata[i]));
344       len += 5;
345     }
346     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s\n", buf));
347     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT));
348 
349     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Reduce Rootdata in type of unsigned char\n"));
350     len = 0;
351     PetscCall(PetscSNPrintf(buf, 256, "%4d:", rank));
352     len += 5;
353     for (i = 0; i < nrootsalloc; i++) {
354       PetscCall(PetscSNPrintf(buf + len, 256 - len, "%5u", rootdata[i]));
355       len += 5;
356     }
357     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s\n", buf));
358     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT));
359 
360     PetscCall(PetscFree2(rootdata, leafdata));
361   }
362 
363   if (test_degree) {
364     const PetscInt *degree;
365     PetscCall(PetscSFComputeDegreeBegin(sf, &degree));
366     PetscCall(PetscSFComputeDegreeEnd(sf, &degree));
367     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Root degrees\n"));
368     PetscCall(PetscIntView(nrootsalloc, degree, PETSC_VIEWER_STDOUT_WORLD));
369   }
370 
371   if (test_fetchandop) {
372     /* Cannot use text compare here because token ordering is not deterministic */
373     PetscInt *leafdata, *leafupdate, *rootdata;
374     PetscCall(PetscMalloc3(nleavesalloc, &leafdata, nleavesalloc, &leafupdate, nrootsalloc, &rootdata));
375     for (i = 0; i < nleavesalloc; i++) leafdata[i] = -1;
376     for (i = 0; i < nleaves; i++) leafdata[i * stride] = 1;
377     for (i = 0; i < nrootsalloc; i++) rootdata[i] = -1;
378     for (i = 0; i < nroots; i++) rootdata[i * stride] = 0;
379     PetscCall(PetscSFFetchAndOpBegin(sf, MPIU_INT, rootdata, leafdata, leafupdate, mop));
380     PetscCall(PetscSFFetchAndOpEnd(sf, MPIU_INT, rootdata, leafdata, leafupdate, mop));
381     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Rootdata (sum of 1 from each leaf)\n"));
382     PetscCall(PetscIntView(nrootsalloc, rootdata, PETSC_VIEWER_STDOUT_WORLD));
383     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Leafupdate (value at roots prior to my atomic update)\n"));
384     PetscCall(PetscIntView(nleavesalloc, leafupdate, PETSC_VIEWER_STDOUT_WORLD));
385     PetscCall(PetscFree3(leafdata, leafupdate, rootdata));
386   }
387 
388   if (test_gather) {
389     const PetscInt *degree;
390     PetscInt        inedges, *indata, *outdata;
391     PetscCall(PetscSFComputeDegreeBegin(sf, &degree));
392     PetscCall(PetscSFComputeDegreeEnd(sf, &degree));
393     for (i = 0, inedges = 0; i < nrootsalloc; i++) inedges += degree[i];
394     PetscCall(PetscMalloc2(inedges, &indata, nleavesalloc, &outdata));
395     for (i = 0; i < nleavesalloc; i++) outdata[i] = -1;
396     for (i = 0; i < nleaves; i++) outdata[i * stride] = 1000 * (rank + 1) + i;
397     PetscCall(PetscSFGatherBegin(sf, MPIU_INT, outdata, indata));
398     PetscCall(PetscSFGatherEnd(sf, MPIU_INT, outdata, indata));
399     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Gathered data at multi-roots from leaves\n"));
400     PetscCall(PetscIntView(inedges, indata, PETSC_VIEWER_STDOUT_WORLD));
401     PetscCall(PetscFree2(indata, outdata));
402   }
403 
404   if (test_scatter) {
405     const PetscInt *degree;
406     PetscInt        j, count, inedges, *indata, *outdata;
407     PetscCall(PetscSFComputeDegreeBegin(sf, &degree));
408     PetscCall(PetscSFComputeDegreeEnd(sf, &degree));
409     for (i = 0, inedges = 0; i < nrootsalloc; i++) inedges += degree[i];
410     PetscCall(PetscMalloc2(inedges, &indata, nleavesalloc, &outdata));
411     for (i = 0; i < nleavesalloc; i++) outdata[i] = -1;
412     for (i = 0, count = 0; i < nrootsalloc; i++) {
413       for (j = 0; j < degree[i]; j++) indata[count++] = 1000 * (rank + 1) + 100 * i + j;
414     }
415     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Data at multi-roots, to scatter to leaves\n"));
416     PetscCall(PetscIntView(inedges, indata, PETSC_VIEWER_STDOUT_WORLD));
417 
418     PetscCall(PetscSFScatterBegin(sf, MPIU_INT, indata, outdata));
419     PetscCall(PetscSFScatterEnd(sf, MPIU_INT, indata, outdata));
420     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Scattered data at leaves\n"));
421     PetscCall(PetscIntView(nleavesalloc, outdata, PETSC_VIEWER_STDOUT_WORLD));
422     PetscCall(PetscFree2(indata, outdata));
423   }
424 
425   if (test_embed) {
426     const PetscInt nroots = 1 + (PetscInt)(rank == 0);
427     PetscInt       selected[2];
428     PetscSF        esf;
429 
430     selected[0] = stride;
431     selected[1] = 2 * stride;
432     PetscCall(PetscSFCreateEmbeddedRootSF(sf, nroots, selected, &esf));
433     PetscCall(PetscSFSetUp(esf));
434     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Embedded PetscSF\n"));
435     PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_INFO_DETAIL));
436     PetscCall(PetscSFView(esf, PETSC_VIEWER_STDOUT_WORLD));
437     PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
438     PetscCall(PetscSFDestroy(&esf));
439   }
440 
441   if (test_invert) {
442     const PetscInt *degree;
443     PetscInt       *mRootsOrigNumbering;
444     PetscInt        inedges;
445     PetscSF         msf, imsf;
446 
447     PetscCall(PetscSFGetMultiSF(sf, &msf));
448     PetscCall(PetscSFCreateInverseSF(msf, &imsf));
449     PetscCall(PetscSFSetUp(msf));
450     PetscCall(PetscSFSetUp(imsf));
451     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Multi-SF\n"));
452     PetscCall(PetscSFView(msf, PETSC_VIEWER_STDOUT_WORLD));
453     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Multi-SF roots indices in original SF roots numbering\n"));
454     PetscCall(PetscSFComputeDegreeBegin(sf, &degree));
455     PetscCall(PetscSFComputeDegreeEnd(sf, &degree));
456     PetscCall(PetscSFComputeMultiRootOriginalNumbering(sf, degree, &inedges, &mRootsOrigNumbering));
457     PetscCall(PetscIntView(inedges, mRootsOrigNumbering, PETSC_VIEWER_STDOUT_WORLD));
458     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Inverse of Multi-SF\n"));
459     PetscCall(PetscSFView(imsf, PETSC_VIEWER_STDOUT_WORLD));
460     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Inverse of Multi-SF, original numbering\n"));
461     PetscCall(PetscSFViewCustomLocals_Private(imsf, mRootsOrigNumbering, PETSC_VIEWER_STDOUT_WORLD));
462     PetscCall(PetscSFDestroy(&imsf));
463     PetscCall(PetscFree(mRootsOrigNumbering));
464   }
465 
466   /* Clean storage for star forest. */
467   PetscCall(PetscSFDestroy(&sf));
468   PetscCall(PetscFinalize());
469   return 0;
470 }
471 
472 /*TEST
473 
474    test:
475       nsize: 4
476       filter: grep -v "type" | grep -v "sort"
477       args: -test_bcast -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}}
478       requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
479 
480    test:
481       suffix: 2
482       nsize: 4
483       filter: grep -v "type" | grep -v "sort"
484       args: -test_reduce -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}}
485       requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
486 
487    test:
488       suffix: 2_basic
489       nsize: 4
490       args: -test_reduce -sf_type basic
491 
492    test:
493       suffix: 3
494       nsize: 4
495       filter: grep -v "type" | grep -v "sort"
496       args: -test_degree -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}}
497       requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
498 
499    test:
500       suffix: 3_basic
501       nsize: 4
502       args: -test_degree -sf_type basic
503 
504    test:
505       suffix: 4
506       nsize: 4
507       filter: grep -v "type" | grep -v "sort"
508       args: -test_gather -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}}
509       requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
510 
511    test:
512       suffix: 4_basic
513       nsize: 4
514       args: -test_gather -sf_type basic
515 
516    test:
517       suffix: 4_stride
518       nsize: 4
519       args: -test_gather -sf_type basic -stride 2
520 
521    test:
522       suffix: 5
523       nsize: 4
524       filter: grep -v "type" | grep -v "sort"
525       args: -test_scatter -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}}
526       requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
527 
528    test:
529       suffix: 5_basic
530       nsize: 4
531       args: -test_scatter -sf_type basic
532 
533    test:
534       suffix: 5_stride
535       nsize: 4
536       args: -test_scatter -sf_type basic -stride 2
537 
538    test:
539       suffix: 6
540       nsize: 4
541       filter: grep -v "type" | grep -v "sort"
542       # No -sf_window_flavor dynamic due to bug https://gitlab.com/petsc/petsc/issues/555
543       args: -test_embed -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create allocate}}
544       requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
545 
546    test:
547       suffix: 6_basic
548       nsize: 4
549       args: -test_embed -sf_type basic
550 
551    test:
552       suffix: 7
553       nsize: 4
554       filter: grep -v "type" | grep -v "sort"
555       args: -test_invert -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}}
556       requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
557 
558    test:
559       suffix: 7_basic
560       nsize: 4
561       args: -test_invert -sf_type basic
562 
563    test:
564       suffix: basic
565       nsize: 4
566       args: -test_bcast -sf_type basic
567       output_file: output/ex1_1_basic.out
568 
569    test:
570       suffix: bcastop_basic
571       nsize: 4
572       args: -test_bcastop -sf_type basic
573       output_file: output/ex1_bcastop_basic.out
574 
575    test:
576       suffix: 8
577       nsize: 3
578       filter: grep -v "type" | grep -v "sort"
579       args: -test_bcast -test_sf_distribute -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}}
580       requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
581 
582    test:
583       suffix: 8_basic
584       nsize: 3
585       args: -test_bcast -test_sf_distribute -sf_type basic
586 
587    test:
588       suffix: 9_char
589       nsize: 4
590       args: -sf_type basic -test_bcast -test_reduce -test_op max -test_char
591 
592    # Here we do not test -sf_window_flavor dynamic since it is designed for repeated SFs with few different rootdata pointers
593    test:
594       suffix: 10
595       filter: grep -v "type" | grep -v "sort"
596       nsize: 4
597       args: -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create allocate}} -test_all -test_bcastop 0 -test_fetchandop 0
598       requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
599 
600    # The nightly test suite with MPICH uses ch3:sock, which is broken when winsize == 0 in some of the processes
601    test:
602       suffix: 10_shared
603       output_file: output/ex1_10.out
604       filter: grep -v "type" | grep -v "sort"
605       nsize: 4
606       args: -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor shared -test_all -test_bcastop 0 -test_fetchandop 0
607       requires: defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY) !defined(PETSC_HAVE_MPICH_NUMVERSION) defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
608 
609    test:
610       suffix: 10_basic
611       nsize: 4
612       args: -sf_type basic -test_all -test_bcastop 0 -test_fetchandop 0
613 
614 TEST*/
615