xref: /petsc/src/vec/is/sf/tutorials/ex1.c (revision c6a7a37075f8bf8d34d92c4910d42445b7a3482d)
1c4762a1bSJed Brown static const char help[] = "Test star forest communication (PetscSF)\n\n";
2c4762a1bSJed Brown 
31bb3edfdSPatrick Sanan /*
4c4762a1bSJed Brown     Description: A star is a simple tree with one root and zero or more leaves.
5c4762a1bSJed Brown     A star forest is a union of disjoint stars.
6c4762a1bSJed Brown     Many common communication patterns can be expressed as updates of rootdata using leafdata and vice-versa.
7c4762a1bSJed Brown     This example creates a star forest, communicates values using the graph (see options for types of communication), views the graph, then destroys it.
81bb3edfdSPatrick Sanan */
9c4762a1bSJed Brown 
10c4762a1bSJed Brown /*
11c4762a1bSJed Brown   Include petscsf.h so we can use PetscSF objects. Note that this automatically
12c4762a1bSJed Brown   includes petscsys.h.
13c4762a1bSJed Brown */
14c4762a1bSJed Brown #include <petscsf.h>
15c4762a1bSJed Brown #include <petscviewer.h>
16c4762a1bSJed Brown 
17c4762a1bSJed Brown /* like PetscSFView() but with alternative array of local indices */
18d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSFViewCustomLocals_Private(PetscSF sf, const PetscInt locals[], PetscViewer viewer)
19d71ae5a4SJacob Faibussowitsch {
20c4762a1bSJed Brown   const PetscSFNode *iremote;
21c4762a1bSJed Brown   PetscInt           i, nroots, nleaves, nranks;
22c4762a1bSJed Brown   PetscMPIInt        rank;
23c4762a1bSJed Brown 
24c4762a1bSJed Brown   PetscFunctionBeginUser;
259566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sf), &rank));
269566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sf, &nroots, &nleaves, NULL, &iremote));
279566063dSJacob Faibussowitsch   PetscCall(PetscSFGetRootRanks(sf, &nranks, NULL, NULL, NULL, NULL));
289566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPushTab(viewer));
299566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPushSynchronized(viewer));
309566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] Number of roots=%" PetscInt_FMT ", leaves=%" PetscInt_FMT ", remote ranks=%" PetscInt_FMT "\n", rank, nroots, nleaves, nranks));
3148a46eb9SPierre Jolivet   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));
329566063dSJacob Faibussowitsch   PetscCall(PetscViewerFlush(viewer));
339566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPopTab(viewer));
349566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPopSynchronized(viewer));
353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
36c4762a1bSJed Brown }
37c4762a1bSJed Brown 
38d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
39d71ae5a4SJacob Faibussowitsch {
40c4762a1bSJed Brown   PetscInt     i, nroots, nrootsalloc, nleaves, nleavesalloc, *mine, stride;
41c4762a1bSJed Brown   PetscSFNode *remote;
42c4762a1bSJed Brown   PetscMPIInt  rank, size;
43c4762a1bSJed Brown   PetscSF      sf;
44c4762a1bSJed Brown   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;
45c4762a1bSJed Brown   MPI_Op       mop = MPI_OP_NULL; /* initialize to prevent compiler warnings with cxx_quad build */
46c4762a1bSJed Brown   char         opstring[256];
47c4762a1bSJed Brown   PetscBool    strflg;
48c4762a1bSJed Brown 
49327415f7SBarry Smith   PetscFunctionBeginUser;
509566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
519566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
529566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
53c4762a1bSJed Brown 
54d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, "", "PetscSF Test Options", "none");
55c4762a1bSJed Brown   test_all = PETSC_FALSE;
569566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_all", "Test all SF communications", "", test_all, &test_all, NULL));
57c4762a1bSJed Brown   test_bcast = test_all;
589566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_bcast", "Test broadcast", "", test_bcast, &test_bcast, NULL));
59c4762a1bSJed Brown   test_bcastop = test_all;
609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_bcastop", "Test broadcast and reduce", "", test_bcastop, &test_bcastop, NULL));
61c4762a1bSJed Brown   test_reduce = test_all;
629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_reduce", "Test reduction", "", test_reduce, &test_reduce, NULL));
63c4762a1bSJed Brown   test_char = test_all;
649566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_char", "Test signed char, unsigned char, and char", "", test_char, &test_char, NULL));
65c4762a1bSJed Brown   mop = MPI_SUM;
66*c6a7a370SJeremy L Thompson   PetscCall(PetscStrncpy(opstring, "sum", sizeof(opstring)));
679566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-test_op", "Designate which MPI_Op to use", "", opstring, opstring, sizeof(opstring), NULL));
689566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp("sum", opstring, &strflg));
69ad540459SPierre Jolivet   if (strflg) mop = MPIU_SUM;
709566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp("prod", opstring, &strflg));
71ad540459SPierre Jolivet   if (strflg) mop = MPI_PROD;
729566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp("max", opstring, &strflg));
73ad540459SPierre Jolivet   if (strflg) mop = MPI_MAX;
749566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp("min", opstring, &strflg));
75ad540459SPierre Jolivet   if (strflg) mop = MPI_MIN;
769566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp("land", opstring, &strflg));
77ad540459SPierre Jolivet   if (strflg) mop = MPI_LAND;
789566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp("band", opstring, &strflg));
79ad540459SPierre Jolivet   if (strflg) mop = MPI_BAND;
809566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp("lor", opstring, &strflg));
81ad540459SPierre Jolivet   if (strflg) mop = MPI_LOR;
829566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp("bor", opstring, &strflg));
83ad540459SPierre Jolivet   if (strflg) mop = MPI_BOR;
849566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp("lxor", opstring, &strflg));
85ad540459SPierre Jolivet   if (strflg) mop = MPI_LXOR;
869566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp("bxor", opstring, &strflg));
87ad540459SPierre Jolivet   if (strflg) mop = MPI_BXOR;
88c4762a1bSJed Brown   test_degree = test_all;
899566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_degree", "Test computation of vertex degree", "", test_degree, &test_degree, NULL));
90c4762a1bSJed Brown   test_fetchandop = test_all;
919566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_fetchandop", "Test atomic Fetch-And-Op", "", test_fetchandop, &test_fetchandop, NULL));
92c4762a1bSJed Brown   test_gather = test_all;
939566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_gather", "Test point gather", "", test_gather, &test_gather, NULL));
94c4762a1bSJed Brown   test_scatter = test_all;
959566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_scatter", "Test point scatter", "", test_scatter, &test_scatter, NULL));
96c4762a1bSJed Brown   test_embed = test_all;
979566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_embed", "Test point embed", "", test_embed, &test_embed, NULL));
98c4762a1bSJed Brown   test_invert = test_all;
999566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_invert", "Test point invert", "", test_invert, &test_invert, NULL));
100c4762a1bSJed Brown   stride = 1;
1019566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-stride", "Stride for leaf and root data", "", stride, &stride, NULL));
102c4762a1bSJed Brown   test_sf_distribute = PETSC_FALSE;
1039566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_sf_distribute", "Create an SF that 'distributes' to each process, like an alltoall", "", test_sf_distribute, &test_sf_distribute, NULL));
1049566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-test_op", "Designate which MPI_Op to use", "", opstring, opstring, sizeof(opstring), NULL));
105d0609cedSBarry Smith   PetscOptionsEnd();
106c4762a1bSJed Brown 
107c4762a1bSJed Brown   if (test_sf_distribute) {
108c4762a1bSJed Brown     nroots       = size;
109c4762a1bSJed Brown     nrootsalloc  = size;
110c4762a1bSJed Brown     nleaves      = size;
111c4762a1bSJed Brown     nleavesalloc = size;
112c4762a1bSJed Brown     mine         = NULL;
1139566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nleaves, &remote));
114c4762a1bSJed Brown     for (i = 0; i < size; i++) {
115c4762a1bSJed Brown       remote[i].rank  = i;
116c4762a1bSJed Brown       remote[i].index = rank;
117c4762a1bSJed Brown     }
118c4762a1bSJed Brown   } else {
119c4762a1bSJed Brown     nroots       = 2 + (PetscInt)(rank == 0);
120c4762a1bSJed Brown     nrootsalloc  = nroots * stride;
121c4762a1bSJed Brown     nleaves      = 2 + (PetscInt)(rank > 0);
122c4762a1bSJed Brown     nleavesalloc = nleaves * stride;
123c4762a1bSJed Brown     mine         = NULL;
124c4762a1bSJed Brown     if (stride > 1) {
125c4762a1bSJed Brown       PetscInt i;
126c4762a1bSJed Brown 
1279566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nleaves, &mine));
128ad540459SPierre Jolivet       for (i = 0; i < nleaves; i++) mine[i] = stride * i;
129c4762a1bSJed Brown     }
1309566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nleaves, &remote));
131c4762a1bSJed Brown     /* Left periodic neighbor */
132c4762a1bSJed Brown     remote[0].rank  = (rank + size - 1) % size;
133c4762a1bSJed Brown     remote[0].index = 1 * stride;
134c4762a1bSJed Brown     /* Right periodic neighbor */
135c4762a1bSJed Brown     remote[1].rank  = (rank + 1) % size;
136c4762a1bSJed Brown     remote[1].index = 0 * stride;
137c4762a1bSJed Brown     if (rank > 0) { /* All processes reference rank 0, index 1 */
138c4762a1bSJed Brown       remote[2].rank  = 0;
139c4762a1bSJed Brown       remote[2].index = 2 * stride;
140c4762a1bSJed Brown     }
141c4762a1bSJed Brown   }
142c4762a1bSJed Brown 
143c4762a1bSJed Brown   /* Create a star forest for communication. In this example, the leaf space is dense, so we pass NULL. */
1449566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(PETSC_COMM_WORLD, &sf));
1459566063dSJacob Faibussowitsch   PetscCall(PetscSFSetFromOptions(sf));
1469566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(sf, nrootsalloc, nleaves, mine, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER));
1479566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sf));
148c4762a1bSJed Brown 
149c4762a1bSJed Brown   /* View graph, mostly useful for debugging purposes. */
1509566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_INFO_DETAIL));
1519566063dSJacob Faibussowitsch   PetscCall(PetscSFView(sf, PETSC_VIEWER_STDOUT_WORLD));
1529566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
153c4762a1bSJed Brown 
154c4762a1bSJed Brown   if (test_bcast) { /* broadcast rootdata into leafdata */
155c4762a1bSJed Brown     PetscInt *rootdata, *leafdata;
1562d4ee042Sprj-     /* Allocate space for send and receive buffers. This example communicates PetscInt, but other types, including
157c4762a1bSJed Brown      * user-defined structures, could also be used. */
1589566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nrootsalloc, &rootdata, nleavesalloc, &leafdata));
159c4762a1bSJed Brown     /* Set rootdata buffer to be broadcast */
160c4762a1bSJed Brown     for (i = 0; i < nrootsalloc; i++) rootdata[i] = -1;
161c4762a1bSJed Brown     for (i = 0; i < nroots; i++) rootdata[i * stride] = 100 * (rank + 1) + i;
162c4762a1bSJed Brown     /* Initialize local buffer, these values are never used. */
163c4762a1bSJed Brown     for (i = 0; i < nleavesalloc; i++) leafdata[i] = -1;
164c4762a1bSJed Brown     /* Broadcast entries from rootdata to leafdata. Computation or other communication can be performed between the begin and end calls. */
1659566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(sf, MPIU_INT, rootdata, leafdata, MPI_REPLACE));
1669566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(sf, MPIU_INT, rootdata, leafdata, MPI_REPLACE));
1679566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Bcast Rootdata\n"));
1689566063dSJacob Faibussowitsch     PetscCall(PetscIntView(nrootsalloc, rootdata, PETSC_VIEWER_STDOUT_WORLD));
1699566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Bcast Leafdata\n"));
1709566063dSJacob Faibussowitsch     PetscCall(PetscIntView(nleavesalloc, leafdata, PETSC_VIEWER_STDOUT_WORLD));
1719566063dSJacob Faibussowitsch     PetscCall(PetscFree2(rootdata, leafdata));
172c4762a1bSJed Brown   }
173c4762a1bSJed Brown 
174c4762a1bSJed Brown   if (test_bcast && test_char) { /* Bcast with char */
175c4762a1bSJed Brown     PetscInt len;
176c4762a1bSJed Brown     char     buf[256];
177c4762a1bSJed Brown     char    *rootdata, *leafdata;
1789566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nrootsalloc, &rootdata, nleavesalloc, &leafdata));
179c4762a1bSJed Brown     /* Set rootdata buffer to be broadcast */
180c4762a1bSJed Brown     for (i = 0; i < nrootsalloc; i++) rootdata[i] = '*';
181c4762a1bSJed Brown     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 */
182c4762a1bSJed Brown     /* Initialize local buffer, these values are never used. */
183c4762a1bSJed Brown     for (i = 0; i < nleavesalloc; i++) leafdata[i] = '?';
184c4762a1bSJed Brown 
1859566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(sf, MPI_CHAR, rootdata, leafdata, MPI_REPLACE));
1869566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(sf, MPI_CHAR, rootdata, leafdata, MPI_REPLACE));
187c4762a1bSJed Brown 
1889566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Bcast Rootdata in type of char\n"));
1899371c9d4SSatish Balay     len = 0;
1909371c9d4SSatish Balay     PetscCall(PetscSNPrintf(buf, 256, "%4d:", rank));
1919371c9d4SSatish Balay     len += 5;
1929371c9d4SSatish Balay     for (i = 0; i < nrootsalloc; i++) {
1939371c9d4SSatish Balay       PetscCall(PetscSNPrintf(buf + len, 256 - len, "%5c", rootdata[i]));
1949371c9d4SSatish Balay       len += 5;
1959371c9d4SSatish Balay     }
1969566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s\n", buf));
1979566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT));
198c4762a1bSJed Brown 
1999566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Bcast Leafdata in type of char\n"));
2009371c9d4SSatish Balay     len = 0;
2019371c9d4SSatish Balay     PetscCall(PetscSNPrintf(buf, 256, "%4d:", rank));
2029371c9d4SSatish Balay     len += 5;
2039371c9d4SSatish Balay     for (i = 0; i < nleavesalloc; i++) {
2049371c9d4SSatish Balay       PetscCall(PetscSNPrintf(buf + len, 256 - len, "%5c", leafdata[i]));
2059371c9d4SSatish Balay       len += 5;
2069371c9d4SSatish Balay     }
2079566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s\n", buf));
2089566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT));
209c4762a1bSJed Brown 
2109566063dSJacob Faibussowitsch     PetscCall(PetscFree2(rootdata, leafdata));
211c4762a1bSJed Brown   }
212c4762a1bSJed Brown 
213c4762a1bSJed Brown   if (test_bcastop) { /* Reduce rootdata into leafdata */
214c4762a1bSJed Brown     PetscInt *rootdata, *leafdata;
215c4762a1bSJed Brown     /* Allocate space for send and receive buffers. This example communicates PetscInt, but other types, including
216c4762a1bSJed Brown      * user-defined structures, could also be used. */
2179566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nrootsalloc, &rootdata, nleavesalloc, &leafdata));
218c4762a1bSJed Brown     /* Set rootdata buffer to be broadcast */
219c4762a1bSJed Brown     for (i = 0; i < nrootsalloc; i++) rootdata[i] = -1;
220c4762a1bSJed Brown     for (i = 0; i < nroots; i++) rootdata[i * stride] = 100 * (rank + 1) + i;
221c4762a1bSJed Brown     /* Set leaf values to reduce with */
222c4762a1bSJed Brown     for (i = 0; i < nleavesalloc; i++) leafdata[i] = -10 * (rank + 1) - i;
2239566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Pre-BcastAndOp Leafdata\n"));
2249566063dSJacob Faibussowitsch     PetscCall(PetscIntView(nleavesalloc, leafdata, PETSC_VIEWER_STDOUT_WORLD));
225c4762a1bSJed Brown     /* Broadcast entries from rootdata to leafdata. Computation or other communication can be performed between the begin and end calls. */
2269566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(sf, MPIU_INT, rootdata, leafdata, mop));
2279566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(sf, MPIU_INT, rootdata, leafdata, mop));
2289566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## BcastAndOp Rootdata\n"));
2299566063dSJacob Faibussowitsch     PetscCall(PetscIntView(nrootsalloc, rootdata, PETSC_VIEWER_STDOUT_WORLD));
2309566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## BcastAndOp Leafdata\n"));
2319566063dSJacob Faibussowitsch     PetscCall(PetscIntView(nleavesalloc, leafdata, PETSC_VIEWER_STDOUT_WORLD));
2329566063dSJacob Faibussowitsch     PetscCall(PetscFree2(rootdata, leafdata));
233c4762a1bSJed Brown   }
234c4762a1bSJed Brown 
235c4762a1bSJed Brown   if (test_reduce) { /* Reduce leafdata into rootdata */
236c4762a1bSJed Brown     PetscInt *rootdata, *leafdata;
2379566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nrootsalloc, &rootdata, nleavesalloc, &leafdata));
238c4762a1bSJed Brown     /* Initialize rootdata buffer in which the result of the reduction will appear. */
239c4762a1bSJed Brown     for (i = 0; i < nrootsalloc; i++) rootdata[i] = -1;
240c4762a1bSJed Brown     for (i = 0; i < nroots; i++) rootdata[i * stride] = 100 * (rank + 1) + i;
241c4762a1bSJed Brown     /* Set leaf values to reduce. */
242c4762a1bSJed Brown     for (i = 0; i < nleavesalloc; i++) leafdata[i] = -1;
243c4762a1bSJed Brown     for (i = 0; i < nleaves; i++) leafdata[i * stride] = 1000 * (rank + 1) + 10 * i;
2449566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Pre-Reduce Rootdata\n"));
2459566063dSJacob Faibussowitsch     PetscCall(PetscIntView(nrootsalloc, rootdata, PETSC_VIEWER_STDOUT_WORLD));
246c4762a1bSJed Brown     /* Perform reduction. Computation or other communication can be performed between the begin and end calls.
247c4762a1bSJed Brown      * This example sums the values, but other MPI_Ops can be used (e.g MPI_MAX, MPI_PROD). */
2489566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceBegin(sf, MPIU_INT, leafdata, rootdata, mop));
2499566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceEnd(sf, MPIU_INT, leafdata, rootdata, mop));
2509566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Reduce Leafdata\n"));
2519566063dSJacob Faibussowitsch     PetscCall(PetscIntView(nleavesalloc, leafdata, PETSC_VIEWER_STDOUT_WORLD));
2529566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Reduce Rootdata\n"));
2539566063dSJacob Faibussowitsch     PetscCall(PetscIntView(nrootsalloc, rootdata, PETSC_VIEWER_STDOUT_WORLD));
2549566063dSJacob Faibussowitsch     PetscCall(PetscFree2(rootdata, leafdata));
255c4762a1bSJed Brown   }
256c4762a1bSJed Brown 
257c4762a1bSJed Brown   if (test_reduce && test_char) { /* Reduce with signed char */
258c4762a1bSJed Brown     PetscInt     len;
259c4762a1bSJed Brown     char         buf[256];
260c4762a1bSJed Brown     signed char *rootdata, *leafdata;
2619566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nrootsalloc, &rootdata, nleavesalloc, &leafdata));
262c4762a1bSJed Brown     /* Initialize rootdata buffer in which the result of the reduction will appear. */
263c4762a1bSJed Brown     for (i = 0; i < nrootsalloc; i++) rootdata[i] = -1;
264c4762a1bSJed Brown     for (i = 0; i < nroots; i++) rootdata[i * stride] = 10 * (rank + 1) + i;
265c4762a1bSJed Brown     /* Set leaf values to reduce. */
266c4762a1bSJed Brown     for (i = 0; i < nleavesalloc; i++) leafdata[i] = -1;
267c4762a1bSJed Brown     for (i = 0; i < nleaves; i++) leafdata[i * stride] = 50 * (rank + 1) + 10 * i;
2689566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Pre-Reduce Rootdata in type of signed char\n"));
269c4762a1bSJed Brown 
2709371c9d4SSatish Balay     len = 0;
2719371c9d4SSatish Balay     PetscCall(PetscSNPrintf(buf, 256, "%4d:", rank));
2729371c9d4SSatish Balay     len += 5;
2739371c9d4SSatish Balay     for (i = 0; i < nrootsalloc; i++) {
2749371c9d4SSatish Balay       PetscCall(PetscSNPrintf(buf + len, 256 - len, "%5d", rootdata[i]));
2759371c9d4SSatish Balay       len += 5;
2769371c9d4SSatish Balay     }
2779566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s\n", buf));
2789566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT));
279c4762a1bSJed Brown 
280c4762a1bSJed Brown     /* Using MPI_CHAR should trigger an error since MPI standard does not support reduction on MPI_CHAR.
281c4762a1bSJed Brown        Testing with -test_op max, one can see the sign does take effect in MPI_MAX.
282c4762a1bSJed Brown      */
2839566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceBegin(sf, MPI_SIGNED_CHAR, leafdata, rootdata, mop));
2849566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceEnd(sf, MPI_SIGNED_CHAR, leafdata, rootdata, mop));
285c4762a1bSJed Brown 
2869566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Reduce Leafdata in type of signed char\n"));
2879371c9d4SSatish Balay     len = 0;
2889371c9d4SSatish Balay     PetscCall(PetscSNPrintf(buf, 256, "%4d:", rank));
2899371c9d4SSatish Balay     len += 5;
2909371c9d4SSatish Balay     for (i = 0; i < nleavesalloc; i++) {
2919371c9d4SSatish Balay       PetscCall(PetscSNPrintf(buf + len, 256 - len, "%5d", leafdata[i]));
2929371c9d4SSatish Balay       len += 5;
2939371c9d4SSatish Balay     }
2949566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s\n", buf));
2959566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT));
296c4762a1bSJed Brown 
2979566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Reduce Rootdata in type of signed char\n"));
2989371c9d4SSatish Balay     len = 0;
2999371c9d4SSatish Balay     PetscCall(PetscSNPrintf(buf, 256, "%4d:", rank));
3009371c9d4SSatish Balay     len += 5;
3019371c9d4SSatish Balay     for (i = 0; i < nrootsalloc; i++) {
3029371c9d4SSatish Balay       PetscCall(PetscSNPrintf(buf + len, 256 - len, "%5d", rootdata[i]));
3039371c9d4SSatish Balay       len += 5;
3049371c9d4SSatish Balay     }
3059566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s\n", buf));
3069566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT));
307c4762a1bSJed Brown 
3089566063dSJacob Faibussowitsch     PetscCall(PetscFree2(rootdata, leafdata));
309c4762a1bSJed Brown   }
310c4762a1bSJed Brown 
311c4762a1bSJed Brown   if (test_reduce && test_char) { /* Reduce with unsigned char */
312c4762a1bSJed Brown     PetscInt       len;
313c4762a1bSJed Brown     char           buf[256];
314c4762a1bSJed Brown     unsigned char *rootdata, *leafdata;
3159566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nrootsalloc, &rootdata, nleavesalloc, &leafdata));
316c4762a1bSJed Brown     /* Initialize rootdata buffer in which the result of the reduction will appear. */
317b9d03b0cSStefano Zampini     for (i = 0; i < nrootsalloc; i++) rootdata[i] = 0;
318c4762a1bSJed Brown     for (i = 0; i < nroots; i++) rootdata[i * stride] = 10 * (rank + 1) + i;
319c4762a1bSJed Brown     /* Set leaf values to reduce. */
320b9d03b0cSStefano Zampini     for (i = 0; i < nleavesalloc; i++) leafdata[i] = 0;
321c4762a1bSJed Brown     for (i = 0; i < nleaves; i++) leafdata[i * stride] = 50 * (rank + 1) + 10 * i;
3229566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Pre-Reduce Rootdata in type of unsigned char\n"));
323c4762a1bSJed Brown 
3249371c9d4SSatish Balay     len = 0;
3259371c9d4SSatish Balay     PetscCall(PetscSNPrintf(buf, 256, "%4d:", rank));
3269371c9d4SSatish Balay     len += 5;
3279371c9d4SSatish Balay     for (i = 0; i < nrootsalloc; i++) {
3289371c9d4SSatish Balay       PetscCall(PetscSNPrintf(buf + len, 256 - len, "%5u", rootdata[i]));
3299371c9d4SSatish Balay       len += 5;
3309371c9d4SSatish Balay     }
3319566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s\n", buf));
3329566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT));
333c4762a1bSJed Brown 
334c4762a1bSJed Brown     /* Using MPI_CHAR should trigger an error since MPI standard does not support reduction on MPI_CHAR.
335c4762a1bSJed Brown        Testing with -test_op max, one can see the sign does take effect in MPI_MAX.
336c4762a1bSJed Brown      */
3379566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceBegin(sf, MPI_UNSIGNED_CHAR, leafdata, rootdata, mop));
3389566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceEnd(sf, MPI_UNSIGNED_CHAR, leafdata, rootdata, mop));
339c4762a1bSJed Brown 
3409566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Reduce Leafdata in type of unsigned char\n"));
3419371c9d4SSatish Balay     len = 0;
3429371c9d4SSatish Balay     PetscCall(PetscSNPrintf(buf, 256, "%4d:", rank));
3439371c9d4SSatish Balay     len += 5;
3449371c9d4SSatish Balay     for (i = 0; i < nleavesalloc; i++) {
3459371c9d4SSatish Balay       PetscCall(PetscSNPrintf(buf + len, 256 - len, "%5u", leafdata[i]));
3469371c9d4SSatish Balay       len += 5;
3479371c9d4SSatish Balay     }
3489566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s\n", buf));
3499566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT));
350c4762a1bSJed Brown 
3519566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Reduce Rootdata in type of unsigned char\n"));
3529371c9d4SSatish Balay     len = 0;
3539371c9d4SSatish Balay     PetscCall(PetscSNPrintf(buf, 256, "%4d:", rank));
3549371c9d4SSatish Balay     len += 5;
3559371c9d4SSatish Balay     for (i = 0; i < nrootsalloc; i++) {
3569371c9d4SSatish Balay       PetscCall(PetscSNPrintf(buf + len, 256 - len, "%5u", rootdata[i]));
3579371c9d4SSatish Balay       len += 5;
3589371c9d4SSatish Balay     }
3599566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s\n", buf));
3609566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT));
361c4762a1bSJed Brown 
3629566063dSJacob Faibussowitsch     PetscCall(PetscFree2(rootdata, leafdata));
363c4762a1bSJed Brown   }
364c4762a1bSJed Brown 
365c4762a1bSJed Brown   if (test_degree) {
366c4762a1bSJed Brown     const PetscInt *degree;
3679566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeBegin(sf, &degree));
3689566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeEnd(sf, &degree));
3699566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Root degrees\n"));
3709566063dSJacob Faibussowitsch     PetscCall(PetscIntView(nrootsalloc, degree, PETSC_VIEWER_STDOUT_WORLD));
371c4762a1bSJed Brown   }
372c4762a1bSJed Brown 
373c4762a1bSJed Brown   if (test_fetchandop) {
374c4762a1bSJed Brown     /* Cannot use text compare here because token ordering is not deterministic */
375c4762a1bSJed Brown     PetscInt *leafdata, *leafupdate, *rootdata;
3769566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(nleavesalloc, &leafdata, nleavesalloc, &leafupdate, nrootsalloc, &rootdata));
377c4762a1bSJed Brown     for (i = 0; i < nleavesalloc; i++) leafdata[i] = -1;
378c4762a1bSJed Brown     for (i = 0; i < nleaves; i++) leafdata[i * stride] = 1;
379c4762a1bSJed Brown     for (i = 0; i < nrootsalloc; i++) rootdata[i] = -1;
380c4762a1bSJed Brown     for (i = 0; i < nroots; i++) rootdata[i * stride] = 0;
3819566063dSJacob Faibussowitsch     PetscCall(PetscSFFetchAndOpBegin(sf, MPIU_INT, rootdata, leafdata, leafupdate, mop));
3829566063dSJacob Faibussowitsch     PetscCall(PetscSFFetchAndOpEnd(sf, MPIU_INT, rootdata, leafdata, leafupdate, mop));
3839566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Rootdata (sum of 1 from each leaf)\n"));
3849566063dSJacob Faibussowitsch     PetscCall(PetscIntView(nrootsalloc, rootdata, PETSC_VIEWER_STDOUT_WORLD));
3859566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Leafupdate (value at roots prior to my atomic update)\n"));
3869566063dSJacob Faibussowitsch     PetscCall(PetscIntView(nleavesalloc, leafupdate, PETSC_VIEWER_STDOUT_WORLD));
3879566063dSJacob Faibussowitsch     PetscCall(PetscFree3(leafdata, leafupdate, rootdata));
388c4762a1bSJed Brown   }
389c4762a1bSJed Brown 
390c4762a1bSJed Brown   if (test_gather) {
391c4762a1bSJed Brown     const PetscInt *degree;
392c4762a1bSJed Brown     PetscInt        inedges, *indata, *outdata;
3939566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeBegin(sf, &degree));
3949566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeEnd(sf, &degree));
395c4762a1bSJed Brown     for (i = 0, inedges = 0; i < nrootsalloc; i++) inedges += degree[i];
3969566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(inedges, &indata, nleavesalloc, &outdata));
397c4762a1bSJed Brown     for (i = 0; i < nleavesalloc; i++) outdata[i] = -1;
398c4762a1bSJed Brown     for (i = 0; i < nleaves; i++) outdata[i * stride] = 1000 * (rank + 1) + i;
3999566063dSJacob Faibussowitsch     PetscCall(PetscSFGatherBegin(sf, MPIU_INT, outdata, indata));
4009566063dSJacob Faibussowitsch     PetscCall(PetscSFGatherEnd(sf, MPIU_INT, outdata, indata));
4019566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Gathered data at multi-roots from leaves\n"));
4029566063dSJacob Faibussowitsch     PetscCall(PetscIntView(inedges, indata, PETSC_VIEWER_STDOUT_WORLD));
4039566063dSJacob Faibussowitsch     PetscCall(PetscFree2(indata, outdata));
404c4762a1bSJed Brown   }
405c4762a1bSJed Brown 
406c4762a1bSJed Brown   if (test_scatter) {
407c4762a1bSJed Brown     const PetscInt *degree;
408c4762a1bSJed Brown     PetscInt        j, count, inedges, *indata, *outdata;
4099566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeBegin(sf, &degree));
4109566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeEnd(sf, &degree));
411c4762a1bSJed Brown     for (i = 0, inedges = 0; i < nrootsalloc; i++) inedges += degree[i];
4129566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(inedges, &indata, nleavesalloc, &outdata));
413c4762a1bSJed Brown     for (i = 0; i < nleavesalloc; i++) outdata[i] = -1;
414c4762a1bSJed Brown     for (i = 0, count = 0; i < nrootsalloc; i++) {
415c4762a1bSJed Brown       for (j = 0; j < degree[i]; j++) indata[count++] = 1000 * (rank + 1) + 100 * i + j;
416c4762a1bSJed Brown     }
4179566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Data at multi-roots, to scatter to leaves\n"));
4189566063dSJacob Faibussowitsch     PetscCall(PetscIntView(inedges, indata, PETSC_VIEWER_STDOUT_WORLD));
419c4762a1bSJed Brown 
4209566063dSJacob Faibussowitsch     PetscCall(PetscSFScatterBegin(sf, MPIU_INT, indata, outdata));
4219566063dSJacob Faibussowitsch     PetscCall(PetscSFScatterEnd(sf, MPIU_INT, indata, outdata));
4229566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Scattered data at leaves\n"));
4239566063dSJacob Faibussowitsch     PetscCall(PetscIntView(nleavesalloc, outdata, PETSC_VIEWER_STDOUT_WORLD));
4249566063dSJacob Faibussowitsch     PetscCall(PetscFree2(indata, outdata));
425c4762a1bSJed Brown   }
426c4762a1bSJed Brown 
427c4762a1bSJed Brown   if (test_embed) {
428dd400576SPatrick Sanan     const PetscInt nroots = 1 + (PetscInt)(rank == 0);
429c4762a1bSJed Brown     PetscInt       selected[2];
430c4762a1bSJed Brown     PetscSF        esf;
431c4762a1bSJed Brown 
432c4762a1bSJed Brown     selected[0] = stride;
433c4762a1bSJed Brown     selected[1] = 2 * stride;
4349566063dSJacob Faibussowitsch     PetscCall(PetscSFCreateEmbeddedRootSF(sf, nroots, selected, &esf));
4359566063dSJacob Faibussowitsch     PetscCall(PetscSFSetUp(esf));
4369566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Embedded PetscSF\n"));
4379566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_INFO_DETAIL));
4389566063dSJacob Faibussowitsch     PetscCall(PetscSFView(esf, PETSC_VIEWER_STDOUT_WORLD));
4399566063dSJacob Faibussowitsch     PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
4409566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&esf));
441c4762a1bSJed Brown   }
442c4762a1bSJed Brown 
443c4762a1bSJed Brown   if (test_invert) {
444c4762a1bSJed Brown     const PetscInt *degree;
445c4762a1bSJed Brown     PetscInt       *mRootsOrigNumbering;
446c4762a1bSJed Brown     PetscInt        inedges;
447c4762a1bSJed Brown     PetscSF         msf, imsf;
448c4762a1bSJed Brown 
4499566063dSJacob Faibussowitsch     PetscCall(PetscSFGetMultiSF(sf, &msf));
4509566063dSJacob Faibussowitsch     PetscCall(PetscSFCreateInverseSF(msf, &imsf));
4519566063dSJacob Faibussowitsch     PetscCall(PetscSFSetUp(msf));
4529566063dSJacob Faibussowitsch     PetscCall(PetscSFSetUp(imsf));
4539566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Multi-SF\n"));
4549566063dSJacob Faibussowitsch     PetscCall(PetscSFView(msf, PETSC_VIEWER_STDOUT_WORLD));
4559566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Multi-SF roots indices in original SF roots numbering\n"));
4569566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeBegin(sf, &degree));
4579566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeEnd(sf, &degree));
4589566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeMultiRootOriginalNumbering(sf, degree, &inedges, &mRootsOrigNumbering));
4599566063dSJacob Faibussowitsch     PetscCall(PetscIntView(inedges, mRootsOrigNumbering, PETSC_VIEWER_STDOUT_WORLD));
4609566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Inverse of Multi-SF\n"));
4619566063dSJacob Faibussowitsch     PetscCall(PetscSFView(imsf, PETSC_VIEWER_STDOUT_WORLD));
4629566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Inverse of Multi-SF, original numbering\n"));
4639566063dSJacob Faibussowitsch     PetscCall(PetscSFViewCustomLocals_Private(imsf, mRootsOrigNumbering, PETSC_VIEWER_STDOUT_WORLD));
4649566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&imsf));
4659566063dSJacob Faibussowitsch     PetscCall(PetscFree(mRootsOrigNumbering));
466c4762a1bSJed Brown   }
467c4762a1bSJed Brown 
468c4762a1bSJed Brown   /* Clean storage for star forest. */
4699566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sf));
4709566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
471b122ec5aSJacob Faibussowitsch   return 0;
472c4762a1bSJed Brown }
473c4762a1bSJed Brown 
474c4762a1bSJed Brown /*TEST
475c4762a1bSJed Brown 
476c4762a1bSJed Brown    test:
477c4762a1bSJed Brown       nsize: 4
478c4762a1bSJed Brown       filter: grep -v "type" | grep -v "sort"
479c4762a1bSJed Brown       args: -test_bcast -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}}
480dfd57a17SPierre Jolivet       requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
481c4762a1bSJed Brown 
482c4762a1bSJed Brown    test:
483c4762a1bSJed Brown       suffix: 2
484c4762a1bSJed Brown       nsize: 4
485c4762a1bSJed Brown       filter: grep -v "type" | grep -v "sort"
486c4762a1bSJed Brown       args: -test_reduce -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}}
487dfd57a17SPierre Jolivet       requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
488c4762a1bSJed Brown 
489c4762a1bSJed Brown    test:
490c4762a1bSJed Brown       suffix: 2_basic
491c4762a1bSJed Brown       nsize: 4
492c4762a1bSJed Brown       args: -test_reduce -sf_type basic
493c4762a1bSJed Brown 
494c4762a1bSJed Brown    test:
495c4762a1bSJed Brown       suffix: 3
496c4762a1bSJed Brown       nsize: 4
497c4762a1bSJed Brown       filter: grep -v "type" | grep -v "sort"
498c4762a1bSJed Brown       args: -test_degree -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}}
499dfd57a17SPierre Jolivet       requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
500c4762a1bSJed Brown 
501c4762a1bSJed Brown    test:
502c4762a1bSJed Brown       suffix: 3_basic
503c4762a1bSJed Brown       nsize: 4
504c4762a1bSJed Brown       args: -test_degree -sf_type basic
505c4762a1bSJed Brown 
506c4762a1bSJed Brown    test:
507c4762a1bSJed Brown       suffix: 4
508c4762a1bSJed Brown       nsize: 4
509c4762a1bSJed Brown       filter: grep -v "type" | grep -v "sort"
510c4762a1bSJed Brown       args: -test_gather -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}}
511dfd57a17SPierre Jolivet       requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
512c4762a1bSJed Brown 
513c4762a1bSJed Brown    test:
514c4762a1bSJed Brown       suffix: 4_basic
515c4762a1bSJed Brown       nsize: 4
516c4762a1bSJed Brown       args: -test_gather -sf_type basic
517c4762a1bSJed Brown 
518c4762a1bSJed Brown    test:
519c4762a1bSJed Brown       suffix: 4_stride
520c4762a1bSJed Brown       nsize: 4
521c4762a1bSJed Brown       args: -test_gather -sf_type basic -stride 2
522c4762a1bSJed Brown 
523c4762a1bSJed Brown    test:
524c4762a1bSJed Brown       suffix: 5
525c4762a1bSJed Brown       nsize: 4
526c4762a1bSJed Brown       filter: grep -v "type" | grep -v "sort"
527c4762a1bSJed Brown       args: -test_scatter -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}}
528dfd57a17SPierre Jolivet       requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
529c4762a1bSJed Brown 
530c4762a1bSJed Brown    test:
531c4762a1bSJed Brown       suffix: 5_basic
532c4762a1bSJed Brown       nsize: 4
533c4762a1bSJed Brown       args: -test_scatter -sf_type basic
534c4762a1bSJed Brown 
535c4762a1bSJed Brown    test:
536c4762a1bSJed Brown       suffix: 5_stride
537c4762a1bSJed Brown       nsize: 4
538c4762a1bSJed Brown       args: -test_scatter -sf_type basic -stride 2
539c4762a1bSJed Brown 
540c4762a1bSJed Brown    test:
541c4762a1bSJed Brown       suffix: 6
542c4762a1bSJed Brown       nsize: 4
543c4762a1bSJed Brown       filter: grep -v "type" | grep -v "sort"
544c20d7725SJed Brown       # No -sf_window_flavor dynamic due to bug https://gitlab.com/petsc/petsc/issues/555
545c20d7725SJed Brown       args: -test_embed -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create allocate}}
546dfd57a17SPierre Jolivet       requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
547c4762a1bSJed Brown 
548c4762a1bSJed Brown    test:
549c4762a1bSJed Brown       suffix: 6_basic
550c4762a1bSJed Brown       nsize: 4
551c4762a1bSJed Brown       args: -test_embed -sf_type basic
552c4762a1bSJed Brown 
553c4762a1bSJed Brown    test:
554c4762a1bSJed Brown       suffix: 7
555c4762a1bSJed Brown       nsize: 4
556c4762a1bSJed Brown       filter: grep -v "type" | grep -v "sort"
557c4762a1bSJed Brown       args: -test_invert -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}}
558dfd57a17SPierre Jolivet       requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
559c4762a1bSJed Brown 
560c4762a1bSJed Brown    test:
561c4762a1bSJed Brown       suffix: 7_basic
562c4762a1bSJed Brown       nsize: 4
563c4762a1bSJed Brown       args: -test_invert -sf_type basic
564c4762a1bSJed Brown 
565c4762a1bSJed Brown    test:
566c4762a1bSJed Brown       suffix: basic
567c4762a1bSJed Brown       nsize: 4
568c4762a1bSJed Brown       args: -test_bcast -sf_type basic
569c4762a1bSJed Brown       output_file: output/ex1_1_basic.out
570c4762a1bSJed Brown 
571c4762a1bSJed Brown    test:
572c4762a1bSJed Brown       suffix: bcastop_basic
573c4762a1bSJed Brown       nsize: 4
574c4762a1bSJed Brown       args: -test_bcastop -sf_type basic
575c4762a1bSJed Brown       output_file: output/ex1_bcastop_basic.out
576c4762a1bSJed Brown 
577c4762a1bSJed Brown    test:
578c4762a1bSJed Brown       suffix: 8
579c4762a1bSJed Brown       nsize: 3
580c4762a1bSJed Brown       filter: grep -v "type" | grep -v "sort"
581c4762a1bSJed Brown       args: -test_bcast -test_sf_distribute -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}}
582dfd57a17SPierre Jolivet       requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
583c4762a1bSJed Brown 
584c4762a1bSJed Brown    test:
585c4762a1bSJed Brown       suffix: 8_basic
586c4762a1bSJed Brown       nsize: 3
587c4762a1bSJed Brown       args: -test_bcast -test_sf_distribute -sf_type basic
588c4762a1bSJed Brown 
589c4762a1bSJed Brown    test:
590c4762a1bSJed Brown       suffix: 9_char
591c4762a1bSJed Brown       nsize: 4
592c4762a1bSJed Brown       args: -sf_type basic -test_bcast -test_reduce -test_op max -test_char
593c4762a1bSJed Brown 
594c4762a1bSJed Brown    # Here we do not test -sf_window_flavor dynamic since it is designed for repeated SFs with few different rootdata pointers
595c4762a1bSJed Brown    test:
596c4762a1bSJed Brown       suffix: 10
597c4762a1bSJed Brown       filter: grep -v "type" | grep -v "sort"
598c4762a1bSJed Brown       nsize: 4
599c4762a1bSJed Brown       args: -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create allocate}} -test_all -test_bcastop 0 -test_fetchandop 0
600dfd57a17SPierre Jolivet       requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
601c4762a1bSJed Brown 
602c4762a1bSJed Brown    # The nightly test suite with MPICH uses ch3:sock, which is broken when winsize == 0 in some of the processes
603c4762a1bSJed Brown    test:
604c4762a1bSJed Brown       suffix: 10_shared
605c4762a1bSJed Brown       output_file: output/ex1_10.out
606c4762a1bSJed Brown       filter: grep -v "type" | grep -v "sort"
607c4762a1bSJed Brown       nsize: 4
608c4762a1bSJed Brown       args: -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor shared -test_all -test_bcastop 0 -test_fetchandop 0
609dfd57a17SPierre Jolivet       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)
610c4762a1bSJed Brown 
611c4762a1bSJed Brown    test:
612c4762a1bSJed Brown       suffix: 10_basic
613c4762a1bSJed Brown       nsize: 4
614c4762a1bSJed Brown       args: -sf_type basic -test_all -test_bcastop 0 -test_fetchandop 0
615c4762a1bSJed Brown 
616c4762a1bSJed Brown TEST*/
617