xref: /petsc/src/vec/is/sf/tutorials/ex1.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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 */
189371c9d4SSatish Balay static PetscErrorCode PetscSFViewCustomLocals_Private(PetscSF sf, const PetscInt locals[], PetscViewer viewer) {
19c4762a1bSJed Brown   const PetscSFNode *iremote;
20c4762a1bSJed Brown   PetscInt           i, nroots, nleaves, nranks;
21c4762a1bSJed Brown   PetscMPIInt        rank;
22c4762a1bSJed Brown 
23c4762a1bSJed Brown   PetscFunctionBeginUser;
249566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)sf), &rank));
259566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sf, &nroots, &nleaves, NULL, &iremote));
269566063dSJacob Faibussowitsch   PetscCall(PetscSFGetRootRanks(sf, &nranks, NULL, NULL, NULL, NULL));
279566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPushTab(viewer));
289566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPushSynchronized(viewer));
299566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] Number of roots=%" PetscInt_FMT ", leaves=%" PetscInt_FMT ", remote ranks=%" PetscInt_FMT "\n", rank, nroots, nleaves, nranks));
30*48a46eb9SPierre 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));
319566063dSJacob Faibussowitsch   PetscCall(PetscViewerFlush(viewer));
329566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPopTab(viewer));
339566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPopSynchronized(viewer));
34c4762a1bSJed Brown   PetscFunctionReturn(0);
35c4762a1bSJed Brown }
36c4762a1bSJed Brown 
379371c9d4SSatish Balay int main(int argc, char **argv) {
38c4762a1bSJed Brown   PetscInt     i, nroots, nrootsalloc, nleaves, nleavesalloc, *mine, stride;
39c4762a1bSJed Brown   PetscSFNode *remote;
40c4762a1bSJed Brown   PetscMPIInt  rank, size;
41c4762a1bSJed Brown   PetscSF      sf;
42c4762a1bSJed 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;
43c4762a1bSJed Brown   MPI_Op       mop = MPI_OP_NULL; /* initialize to prevent compiler warnings with cxx_quad build */
44c4762a1bSJed Brown   char         opstring[256];
45c4762a1bSJed Brown   PetscBool    strflg;
46c4762a1bSJed Brown 
47327415f7SBarry Smith   PetscFunctionBeginUser;
489566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
499566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
509566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
51c4762a1bSJed Brown 
52d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, "", "PetscSF Test Options", "none");
53c4762a1bSJed Brown   test_all = PETSC_FALSE;
549566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_all", "Test all SF communications", "", test_all, &test_all, NULL));
55c4762a1bSJed Brown   test_bcast = test_all;
569566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_bcast", "Test broadcast", "", test_bcast, &test_bcast, NULL));
57c4762a1bSJed Brown   test_bcastop = test_all;
589566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_bcastop", "Test broadcast and reduce", "", test_bcastop, &test_bcastop, NULL));
59c4762a1bSJed Brown   test_reduce = test_all;
609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_reduce", "Test reduction", "", test_reduce, &test_reduce, NULL));
61c4762a1bSJed Brown   test_char = test_all;
629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_char", "Test signed char, unsigned char, and char", "", test_char, &test_char, NULL));
63c4762a1bSJed Brown   mop = MPI_SUM;
649566063dSJacob Faibussowitsch   PetscCall(PetscStrcpy(opstring, "sum"));
659566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-test_op", "Designate which MPI_Op to use", "", opstring, opstring, sizeof(opstring), NULL));
669566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp("sum", opstring, &strflg));
679371c9d4SSatish Balay   if (strflg) { mop = MPIU_SUM; }
689566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp("prod", opstring, &strflg));
699371c9d4SSatish Balay   if (strflg) { mop = MPI_PROD; }
709566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp("max", opstring, &strflg));
719371c9d4SSatish Balay   if (strflg) { mop = MPI_MAX; }
729566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp("min", opstring, &strflg));
739371c9d4SSatish Balay   if (strflg) { mop = MPI_MIN; }
749566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp("land", opstring, &strflg));
759371c9d4SSatish Balay   if (strflg) { mop = MPI_LAND; }
769566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp("band", opstring, &strflg));
779371c9d4SSatish Balay   if (strflg) { mop = MPI_BAND; }
789566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp("lor", opstring, &strflg));
799371c9d4SSatish Balay   if (strflg) { mop = MPI_LOR; }
809566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp("bor", opstring, &strflg));
819371c9d4SSatish Balay   if (strflg) { mop = MPI_BOR; }
829566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp("lxor", opstring, &strflg));
839371c9d4SSatish Balay   if (strflg) { mop = MPI_LXOR; }
849566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp("bxor", opstring, &strflg));
859371c9d4SSatish Balay   if (strflg) { mop = MPI_BXOR; }
86c4762a1bSJed Brown   test_degree = test_all;
879566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_degree", "Test computation of vertex degree", "", test_degree, &test_degree, NULL));
88c4762a1bSJed Brown   test_fetchandop = test_all;
899566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_fetchandop", "Test atomic Fetch-And-Op", "", test_fetchandop, &test_fetchandop, NULL));
90c4762a1bSJed Brown   test_gather = test_all;
919566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_gather", "Test point gather", "", test_gather, &test_gather, NULL));
92c4762a1bSJed Brown   test_scatter = test_all;
939566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_scatter", "Test point scatter", "", test_scatter, &test_scatter, NULL));
94c4762a1bSJed Brown   test_embed = test_all;
959566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_embed", "Test point embed", "", test_embed, &test_embed, NULL));
96c4762a1bSJed Brown   test_invert = test_all;
979566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_invert", "Test point invert", "", test_invert, &test_invert, NULL));
98c4762a1bSJed Brown   stride = 1;
999566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-stride", "Stride for leaf and root data", "", stride, &stride, NULL));
100c4762a1bSJed Brown   test_sf_distribute = PETSC_FALSE;
1019566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_sf_distribute", "Create an SF that 'distributes' to each process, like an alltoall", "", test_sf_distribute, &test_sf_distribute, NULL));
1029566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-test_op", "Designate which MPI_Op to use", "", opstring, opstring, sizeof(opstring), NULL));
103d0609cedSBarry Smith   PetscOptionsEnd();
104c4762a1bSJed Brown 
105c4762a1bSJed Brown   if (test_sf_distribute) {
106c4762a1bSJed Brown     nroots       = size;
107c4762a1bSJed Brown     nrootsalloc  = size;
108c4762a1bSJed Brown     nleaves      = size;
109c4762a1bSJed Brown     nleavesalloc = size;
110c4762a1bSJed Brown     mine         = NULL;
1119566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nleaves, &remote));
112c4762a1bSJed Brown     for (i = 0; i < size; i++) {
113c4762a1bSJed Brown       remote[i].rank  = i;
114c4762a1bSJed Brown       remote[i].index = rank;
115c4762a1bSJed Brown     }
116c4762a1bSJed Brown   } else {
117c4762a1bSJed Brown     nroots       = 2 + (PetscInt)(rank == 0);
118c4762a1bSJed Brown     nrootsalloc  = nroots * stride;
119c4762a1bSJed Brown     nleaves      = 2 + (PetscInt)(rank > 0);
120c4762a1bSJed Brown     nleavesalloc = nleaves * stride;
121c4762a1bSJed Brown     mine         = NULL;
122c4762a1bSJed Brown     if (stride > 1) {
123c4762a1bSJed Brown       PetscInt i;
124c4762a1bSJed Brown 
1259566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nleaves, &mine));
1269371c9d4SSatish Balay       for (i = 0; i < nleaves; i++) { mine[i] = stride * i; }
127c4762a1bSJed Brown     }
1289566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nleaves, &remote));
129c4762a1bSJed Brown     /* Left periodic neighbor */
130c4762a1bSJed Brown     remote[0].rank  = (rank + size - 1) % size;
131c4762a1bSJed Brown     remote[0].index = 1 * stride;
132c4762a1bSJed Brown     /* Right periodic neighbor */
133c4762a1bSJed Brown     remote[1].rank  = (rank + 1) % size;
134c4762a1bSJed Brown     remote[1].index = 0 * stride;
135c4762a1bSJed Brown     if (rank > 0) { /* All processes reference rank 0, index 1 */
136c4762a1bSJed Brown       remote[2].rank  = 0;
137c4762a1bSJed Brown       remote[2].index = 2 * stride;
138c4762a1bSJed Brown     }
139c4762a1bSJed Brown   }
140c4762a1bSJed Brown 
141c4762a1bSJed Brown   /* Create a star forest for communication. In this example, the leaf space is dense, so we pass NULL. */
1429566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(PETSC_COMM_WORLD, &sf));
1439566063dSJacob Faibussowitsch   PetscCall(PetscSFSetFromOptions(sf));
1449566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(sf, nrootsalloc, nleaves, mine, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER));
1459566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sf));
146c4762a1bSJed Brown 
147c4762a1bSJed Brown   /* View graph, mostly useful for debugging purposes. */
1489566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_INFO_DETAIL));
1499566063dSJacob Faibussowitsch   PetscCall(PetscSFView(sf, PETSC_VIEWER_STDOUT_WORLD));
1509566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
151c4762a1bSJed Brown 
152c4762a1bSJed Brown   if (test_bcast) { /* broadcast rootdata into leafdata */
153c4762a1bSJed Brown     PetscInt *rootdata, *leafdata;
1542d4ee042Sprj-     /* Allocate space for send and receive buffers. This example communicates PetscInt, but other types, including
155c4762a1bSJed Brown      * user-defined structures, could also be used. */
1569566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nrootsalloc, &rootdata, nleavesalloc, &leafdata));
157c4762a1bSJed Brown     /* Set rootdata buffer to be broadcast */
158c4762a1bSJed Brown     for (i = 0; i < nrootsalloc; i++) rootdata[i] = -1;
159c4762a1bSJed Brown     for (i = 0; i < nroots; i++) rootdata[i * stride] = 100 * (rank + 1) + i;
160c4762a1bSJed Brown     /* Initialize local buffer, these values are never used. */
161c4762a1bSJed Brown     for (i = 0; i < nleavesalloc; i++) leafdata[i] = -1;
162c4762a1bSJed Brown     /* Broadcast entries from rootdata to leafdata. Computation or other communication can be performed between the begin and end calls. */
1639566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(sf, MPIU_INT, rootdata, leafdata, MPI_REPLACE));
1649566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(sf, MPIU_INT, rootdata, leafdata, MPI_REPLACE));
1659566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Bcast Rootdata\n"));
1669566063dSJacob Faibussowitsch     PetscCall(PetscIntView(nrootsalloc, rootdata, PETSC_VIEWER_STDOUT_WORLD));
1679566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Bcast Leafdata\n"));
1689566063dSJacob Faibussowitsch     PetscCall(PetscIntView(nleavesalloc, leafdata, PETSC_VIEWER_STDOUT_WORLD));
1699566063dSJacob Faibussowitsch     PetscCall(PetscFree2(rootdata, leafdata));
170c4762a1bSJed Brown   }
171c4762a1bSJed Brown 
172c4762a1bSJed Brown   if (test_bcast && test_char) { /* Bcast with char */
173c4762a1bSJed Brown     PetscInt len;
174c4762a1bSJed Brown     char     buf[256];
175c4762a1bSJed Brown     char    *rootdata, *leafdata;
1769566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nrootsalloc, &rootdata, nleavesalloc, &leafdata));
177c4762a1bSJed Brown     /* Set rootdata buffer to be broadcast */
178c4762a1bSJed Brown     for (i = 0; i < nrootsalloc; i++) rootdata[i] = '*';
179c4762a1bSJed 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 */
180c4762a1bSJed Brown     /* Initialize local buffer, these values are never used. */
181c4762a1bSJed Brown     for (i = 0; i < nleavesalloc; i++) leafdata[i] = '?';
182c4762a1bSJed Brown 
1839566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(sf, MPI_CHAR, rootdata, leafdata, MPI_REPLACE));
1849566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(sf, MPI_CHAR, rootdata, leafdata, MPI_REPLACE));
185c4762a1bSJed Brown 
1869566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Bcast Rootdata in type of char\n"));
1879371c9d4SSatish Balay     len = 0;
1889371c9d4SSatish Balay     PetscCall(PetscSNPrintf(buf, 256, "%4d:", rank));
1899371c9d4SSatish Balay     len += 5;
1909371c9d4SSatish Balay     for (i = 0; i < nrootsalloc; i++) {
1919371c9d4SSatish Balay       PetscCall(PetscSNPrintf(buf + len, 256 - len, "%5c", rootdata[i]));
1929371c9d4SSatish Balay       len += 5;
1939371c9d4SSatish Balay     }
1949566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s\n", buf));
1959566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT));
196c4762a1bSJed Brown 
1979566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Bcast Leafdata in type of char\n"));
1989371c9d4SSatish Balay     len = 0;
1999371c9d4SSatish Balay     PetscCall(PetscSNPrintf(buf, 256, "%4d:", rank));
2009371c9d4SSatish Balay     len += 5;
2019371c9d4SSatish Balay     for (i = 0; i < nleavesalloc; i++) {
2029371c9d4SSatish Balay       PetscCall(PetscSNPrintf(buf + len, 256 - len, "%5c", leafdata[i]));
2039371c9d4SSatish Balay       len += 5;
2049371c9d4SSatish Balay     }
2059566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s\n", buf));
2069566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT));
207c4762a1bSJed Brown 
2089566063dSJacob Faibussowitsch     PetscCall(PetscFree2(rootdata, leafdata));
209c4762a1bSJed Brown   }
210c4762a1bSJed Brown 
211c4762a1bSJed Brown   if (test_bcastop) { /* Reduce rootdata into leafdata */
212c4762a1bSJed Brown     PetscInt *rootdata, *leafdata;
213c4762a1bSJed Brown     /* Allocate space for send and receive buffers. This example communicates PetscInt, but other types, including
214c4762a1bSJed Brown      * user-defined structures, could also be used. */
2159566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nrootsalloc, &rootdata, nleavesalloc, &leafdata));
216c4762a1bSJed Brown     /* Set rootdata buffer to be broadcast */
217c4762a1bSJed Brown     for (i = 0; i < nrootsalloc; i++) rootdata[i] = -1;
218c4762a1bSJed Brown     for (i = 0; i < nroots; i++) rootdata[i * stride] = 100 * (rank + 1) + i;
219c4762a1bSJed Brown     /* Set leaf values to reduce with */
220c4762a1bSJed Brown     for (i = 0; i < nleavesalloc; i++) leafdata[i] = -10 * (rank + 1) - i;
2219566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Pre-BcastAndOp Leafdata\n"));
2229566063dSJacob Faibussowitsch     PetscCall(PetscIntView(nleavesalloc, leafdata, PETSC_VIEWER_STDOUT_WORLD));
223c4762a1bSJed Brown     /* Broadcast entries from rootdata to leafdata. Computation or other communication can be performed between the begin and end calls. */
2249566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(sf, MPIU_INT, rootdata, leafdata, mop));
2259566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(sf, MPIU_INT, rootdata, leafdata, mop));
2269566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## BcastAndOp Rootdata\n"));
2279566063dSJacob Faibussowitsch     PetscCall(PetscIntView(nrootsalloc, rootdata, PETSC_VIEWER_STDOUT_WORLD));
2289566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## BcastAndOp Leafdata\n"));
2299566063dSJacob Faibussowitsch     PetscCall(PetscIntView(nleavesalloc, leafdata, PETSC_VIEWER_STDOUT_WORLD));
2309566063dSJacob Faibussowitsch     PetscCall(PetscFree2(rootdata, leafdata));
231c4762a1bSJed Brown   }
232c4762a1bSJed Brown 
233c4762a1bSJed Brown   if (test_reduce) { /* Reduce leafdata into rootdata */
234c4762a1bSJed Brown     PetscInt *rootdata, *leafdata;
2359566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nrootsalloc, &rootdata, nleavesalloc, &leafdata));
236c4762a1bSJed Brown     /* Initialize rootdata buffer in which the result of the reduction will appear. */
237c4762a1bSJed Brown     for (i = 0; i < nrootsalloc; i++) rootdata[i] = -1;
238c4762a1bSJed Brown     for (i = 0; i < nroots; i++) rootdata[i * stride] = 100 * (rank + 1) + i;
239c4762a1bSJed Brown     /* Set leaf values to reduce. */
240c4762a1bSJed Brown     for (i = 0; i < nleavesalloc; i++) leafdata[i] = -1;
241c4762a1bSJed Brown     for (i = 0; i < nleaves; i++) leafdata[i * stride] = 1000 * (rank + 1) + 10 * i;
2429566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Pre-Reduce Rootdata\n"));
2439566063dSJacob Faibussowitsch     PetscCall(PetscIntView(nrootsalloc, rootdata, PETSC_VIEWER_STDOUT_WORLD));
244c4762a1bSJed Brown     /* Perform reduction. Computation or other communication can be performed between the begin and end calls.
245c4762a1bSJed Brown      * This example sums the values, but other MPI_Ops can be used (e.g MPI_MAX, MPI_PROD). */
2469566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceBegin(sf, MPIU_INT, leafdata, rootdata, mop));
2479566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceEnd(sf, MPIU_INT, leafdata, rootdata, mop));
2489566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Reduce Leafdata\n"));
2499566063dSJacob Faibussowitsch     PetscCall(PetscIntView(nleavesalloc, leafdata, PETSC_VIEWER_STDOUT_WORLD));
2509566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Reduce Rootdata\n"));
2519566063dSJacob Faibussowitsch     PetscCall(PetscIntView(nrootsalloc, rootdata, PETSC_VIEWER_STDOUT_WORLD));
2529566063dSJacob Faibussowitsch     PetscCall(PetscFree2(rootdata, leafdata));
253c4762a1bSJed Brown   }
254c4762a1bSJed Brown 
255c4762a1bSJed Brown   if (test_reduce && test_char) { /* Reduce with signed char */
256c4762a1bSJed Brown     PetscInt     len;
257c4762a1bSJed Brown     char         buf[256];
258c4762a1bSJed Brown     signed char *rootdata, *leafdata;
2599566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nrootsalloc, &rootdata, nleavesalloc, &leafdata));
260c4762a1bSJed Brown     /* Initialize rootdata buffer in which the result of the reduction will appear. */
261c4762a1bSJed Brown     for (i = 0; i < nrootsalloc; i++) rootdata[i] = -1;
262c4762a1bSJed Brown     for (i = 0; i < nroots; i++) rootdata[i * stride] = 10 * (rank + 1) + i;
263c4762a1bSJed Brown     /* Set leaf values to reduce. */
264c4762a1bSJed Brown     for (i = 0; i < nleavesalloc; i++) leafdata[i] = -1;
265c4762a1bSJed Brown     for (i = 0; i < nleaves; i++) leafdata[i * stride] = 50 * (rank + 1) + 10 * i;
2669566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Pre-Reduce Rootdata in type of signed char\n"));
267c4762a1bSJed Brown 
2689371c9d4SSatish Balay     len = 0;
2699371c9d4SSatish Balay     PetscCall(PetscSNPrintf(buf, 256, "%4d:", rank));
2709371c9d4SSatish Balay     len += 5;
2719371c9d4SSatish Balay     for (i = 0; i < nrootsalloc; i++) {
2729371c9d4SSatish Balay       PetscCall(PetscSNPrintf(buf + len, 256 - len, "%5d", rootdata[i]));
2739371c9d4SSatish Balay       len += 5;
2749371c9d4SSatish Balay     }
2759566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s\n", buf));
2769566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT));
277c4762a1bSJed Brown 
278c4762a1bSJed Brown     /* Using MPI_CHAR should trigger an error since MPI standard does not support reduction on MPI_CHAR.
279c4762a1bSJed Brown        Testing with -test_op max, one can see the sign does take effect in MPI_MAX.
280c4762a1bSJed Brown      */
2819566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceBegin(sf, MPI_SIGNED_CHAR, leafdata, rootdata, mop));
2829566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceEnd(sf, MPI_SIGNED_CHAR, leafdata, rootdata, mop));
283c4762a1bSJed Brown 
2849566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Reduce Leafdata in type of signed char\n"));
2859371c9d4SSatish Balay     len = 0;
2869371c9d4SSatish Balay     PetscCall(PetscSNPrintf(buf, 256, "%4d:", rank));
2879371c9d4SSatish Balay     len += 5;
2889371c9d4SSatish Balay     for (i = 0; i < nleavesalloc; i++) {
2899371c9d4SSatish Balay       PetscCall(PetscSNPrintf(buf + len, 256 - len, "%5d", leafdata[i]));
2909371c9d4SSatish Balay       len += 5;
2919371c9d4SSatish Balay     }
2929566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s\n", buf));
2939566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT));
294c4762a1bSJed Brown 
2959566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Reduce Rootdata in type of signed char\n"));
2969371c9d4SSatish Balay     len = 0;
2979371c9d4SSatish Balay     PetscCall(PetscSNPrintf(buf, 256, "%4d:", rank));
2989371c9d4SSatish Balay     len += 5;
2999371c9d4SSatish Balay     for (i = 0; i < nrootsalloc; i++) {
3009371c9d4SSatish Balay       PetscCall(PetscSNPrintf(buf + len, 256 - len, "%5d", rootdata[i]));
3019371c9d4SSatish Balay       len += 5;
3029371c9d4SSatish Balay     }
3039566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s\n", buf));
3049566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT));
305c4762a1bSJed Brown 
3069566063dSJacob Faibussowitsch     PetscCall(PetscFree2(rootdata, leafdata));
307c4762a1bSJed Brown   }
308c4762a1bSJed Brown 
309c4762a1bSJed Brown   if (test_reduce && test_char) { /* Reduce with unsigned char */
310c4762a1bSJed Brown     PetscInt       len;
311c4762a1bSJed Brown     char           buf[256];
312c4762a1bSJed Brown     unsigned char *rootdata, *leafdata;
3139566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nrootsalloc, &rootdata, nleavesalloc, &leafdata));
314c4762a1bSJed Brown     /* Initialize rootdata buffer in which the result of the reduction will appear. */
315b9d03b0cSStefano Zampini     for (i = 0; i < nrootsalloc; i++) rootdata[i] = 0;
316c4762a1bSJed Brown     for (i = 0; i < nroots; i++) rootdata[i * stride] = 10 * (rank + 1) + i;
317c4762a1bSJed Brown     /* Set leaf values to reduce. */
318b9d03b0cSStefano Zampini     for (i = 0; i < nleavesalloc; i++) leafdata[i] = 0;
319c4762a1bSJed Brown     for (i = 0; i < nleaves; i++) leafdata[i * stride] = 50 * (rank + 1) + 10 * i;
3209566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Pre-Reduce Rootdata in type of unsigned char\n"));
321c4762a1bSJed Brown 
3229371c9d4SSatish Balay     len = 0;
3239371c9d4SSatish Balay     PetscCall(PetscSNPrintf(buf, 256, "%4d:", rank));
3249371c9d4SSatish Balay     len += 5;
3259371c9d4SSatish Balay     for (i = 0; i < nrootsalloc; i++) {
3269371c9d4SSatish Balay       PetscCall(PetscSNPrintf(buf + len, 256 - len, "%5u", rootdata[i]));
3279371c9d4SSatish Balay       len += 5;
3289371c9d4SSatish Balay     }
3299566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s\n", buf));
3309566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT));
331c4762a1bSJed Brown 
332c4762a1bSJed Brown     /* Using MPI_CHAR should trigger an error since MPI standard does not support reduction on MPI_CHAR.
333c4762a1bSJed Brown        Testing with -test_op max, one can see the sign does take effect in MPI_MAX.
334c4762a1bSJed Brown      */
3359566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceBegin(sf, MPI_UNSIGNED_CHAR, leafdata, rootdata, mop));
3369566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceEnd(sf, MPI_UNSIGNED_CHAR, leafdata, rootdata, mop));
337c4762a1bSJed Brown 
3389566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Reduce Leafdata in type of unsigned char\n"));
3399371c9d4SSatish Balay     len = 0;
3409371c9d4SSatish Balay     PetscCall(PetscSNPrintf(buf, 256, "%4d:", rank));
3419371c9d4SSatish Balay     len += 5;
3429371c9d4SSatish Balay     for (i = 0; i < nleavesalloc; i++) {
3439371c9d4SSatish Balay       PetscCall(PetscSNPrintf(buf + len, 256 - len, "%5u", leafdata[i]));
3449371c9d4SSatish Balay       len += 5;
3459371c9d4SSatish Balay     }
3469566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s\n", buf));
3479566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT));
348c4762a1bSJed Brown 
3499566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Reduce Rootdata in type of unsigned char\n"));
3509371c9d4SSatish Balay     len = 0;
3519371c9d4SSatish Balay     PetscCall(PetscSNPrintf(buf, 256, "%4d:", rank));
3529371c9d4SSatish Balay     len += 5;
3539371c9d4SSatish Balay     for (i = 0; i < nrootsalloc; i++) {
3549371c9d4SSatish Balay       PetscCall(PetscSNPrintf(buf + len, 256 - len, "%5u", rootdata[i]));
3559371c9d4SSatish Balay       len += 5;
3569371c9d4SSatish Balay     }
3579566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%s\n", buf));
3589566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT));
359c4762a1bSJed Brown 
3609566063dSJacob Faibussowitsch     PetscCall(PetscFree2(rootdata, leafdata));
361c4762a1bSJed Brown   }
362c4762a1bSJed Brown 
363c4762a1bSJed Brown   if (test_degree) {
364c4762a1bSJed Brown     const PetscInt *degree;
3659566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeBegin(sf, &degree));
3669566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeEnd(sf, &degree));
3679566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Root degrees\n"));
3689566063dSJacob Faibussowitsch     PetscCall(PetscIntView(nrootsalloc, degree, PETSC_VIEWER_STDOUT_WORLD));
369c4762a1bSJed Brown   }
370c4762a1bSJed Brown 
371c4762a1bSJed Brown   if (test_fetchandop) {
372c4762a1bSJed Brown     /* Cannot use text compare here because token ordering is not deterministic */
373c4762a1bSJed Brown     PetscInt *leafdata, *leafupdate, *rootdata;
3749566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(nleavesalloc, &leafdata, nleavesalloc, &leafupdate, nrootsalloc, &rootdata));
375c4762a1bSJed Brown     for (i = 0; i < nleavesalloc; i++) leafdata[i] = -1;
376c4762a1bSJed Brown     for (i = 0; i < nleaves; i++) leafdata[i * stride] = 1;
377c4762a1bSJed Brown     for (i = 0; i < nrootsalloc; i++) rootdata[i] = -1;
378c4762a1bSJed Brown     for (i = 0; i < nroots; i++) rootdata[i * stride] = 0;
3799566063dSJacob Faibussowitsch     PetscCall(PetscSFFetchAndOpBegin(sf, MPIU_INT, rootdata, leafdata, leafupdate, mop));
3809566063dSJacob Faibussowitsch     PetscCall(PetscSFFetchAndOpEnd(sf, MPIU_INT, rootdata, leafdata, leafupdate, mop));
3819566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Rootdata (sum of 1 from each leaf)\n"));
3829566063dSJacob Faibussowitsch     PetscCall(PetscIntView(nrootsalloc, rootdata, PETSC_VIEWER_STDOUT_WORLD));
3839566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Leafupdate (value at roots prior to my atomic update)\n"));
3849566063dSJacob Faibussowitsch     PetscCall(PetscIntView(nleavesalloc, leafupdate, PETSC_VIEWER_STDOUT_WORLD));
3859566063dSJacob Faibussowitsch     PetscCall(PetscFree3(leafdata, leafupdate, rootdata));
386c4762a1bSJed Brown   }
387c4762a1bSJed Brown 
388c4762a1bSJed Brown   if (test_gather) {
389c4762a1bSJed Brown     const PetscInt *degree;
390c4762a1bSJed Brown     PetscInt        inedges, *indata, *outdata;
3919566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeBegin(sf, &degree));
3929566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeEnd(sf, &degree));
393c4762a1bSJed Brown     for (i = 0, inedges = 0; i < nrootsalloc; i++) inedges += degree[i];
3949566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(inedges, &indata, nleavesalloc, &outdata));
395c4762a1bSJed Brown     for (i = 0; i < nleavesalloc; i++) outdata[i] = -1;
396c4762a1bSJed Brown     for (i = 0; i < nleaves; i++) outdata[i * stride] = 1000 * (rank + 1) + i;
3979566063dSJacob Faibussowitsch     PetscCall(PetscSFGatherBegin(sf, MPIU_INT, outdata, indata));
3989566063dSJacob Faibussowitsch     PetscCall(PetscSFGatherEnd(sf, MPIU_INT, outdata, indata));
3999566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Gathered data at multi-roots from leaves\n"));
4009566063dSJacob Faibussowitsch     PetscCall(PetscIntView(inedges, indata, PETSC_VIEWER_STDOUT_WORLD));
4019566063dSJacob Faibussowitsch     PetscCall(PetscFree2(indata, outdata));
402c4762a1bSJed Brown   }
403c4762a1bSJed Brown 
404c4762a1bSJed Brown   if (test_scatter) {
405c4762a1bSJed Brown     const PetscInt *degree;
406c4762a1bSJed Brown     PetscInt        j, count, inedges, *indata, *outdata;
4079566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeBegin(sf, &degree));
4089566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeEnd(sf, &degree));
409c4762a1bSJed Brown     for (i = 0, inedges = 0; i < nrootsalloc; i++) inedges += degree[i];
4109566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(inedges, &indata, nleavesalloc, &outdata));
411c4762a1bSJed Brown     for (i = 0; i < nleavesalloc; i++) outdata[i] = -1;
412c4762a1bSJed Brown     for (i = 0, count = 0; i < nrootsalloc; i++) {
413c4762a1bSJed Brown       for (j = 0; j < degree[i]; j++) indata[count++] = 1000 * (rank + 1) + 100 * i + j;
414c4762a1bSJed Brown     }
4159566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Data at multi-roots, to scatter to leaves\n"));
4169566063dSJacob Faibussowitsch     PetscCall(PetscIntView(inedges, indata, PETSC_VIEWER_STDOUT_WORLD));
417c4762a1bSJed Brown 
4189566063dSJacob Faibussowitsch     PetscCall(PetscSFScatterBegin(sf, MPIU_INT, indata, outdata));
4199566063dSJacob Faibussowitsch     PetscCall(PetscSFScatterEnd(sf, MPIU_INT, indata, outdata));
4209566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Scattered data at leaves\n"));
4219566063dSJacob Faibussowitsch     PetscCall(PetscIntView(nleavesalloc, outdata, PETSC_VIEWER_STDOUT_WORLD));
4229566063dSJacob Faibussowitsch     PetscCall(PetscFree2(indata, outdata));
423c4762a1bSJed Brown   }
424c4762a1bSJed Brown 
425c4762a1bSJed Brown   if (test_embed) {
426dd400576SPatrick Sanan     const PetscInt nroots = 1 + (PetscInt)(rank == 0);
427c4762a1bSJed Brown     PetscInt       selected[2];
428c4762a1bSJed Brown     PetscSF        esf;
429c4762a1bSJed Brown 
430c4762a1bSJed Brown     selected[0] = stride;
431c4762a1bSJed Brown     selected[1] = 2 * stride;
4329566063dSJacob Faibussowitsch     PetscCall(PetscSFCreateEmbeddedRootSF(sf, nroots, selected, &esf));
4339566063dSJacob Faibussowitsch     PetscCall(PetscSFSetUp(esf));
4349566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Embedded PetscSF\n"));
4359566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_INFO_DETAIL));
4369566063dSJacob Faibussowitsch     PetscCall(PetscSFView(esf, PETSC_VIEWER_STDOUT_WORLD));
4379566063dSJacob Faibussowitsch     PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
4389566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&esf));
439c4762a1bSJed Brown   }
440c4762a1bSJed Brown 
441c4762a1bSJed Brown   if (test_invert) {
442c4762a1bSJed Brown     const PetscInt *degree;
443c4762a1bSJed Brown     PetscInt       *mRootsOrigNumbering;
444c4762a1bSJed Brown     PetscInt        inedges;
445c4762a1bSJed Brown     PetscSF         msf, imsf;
446c4762a1bSJed Brown 
4479566063dSJacob Faibussowitsch     PetscCall(PetscSFGetMultiSF(sf, &msf));
4489566063dSJacob Faibussowitsch     PetscCall(PetscSFCreateInverseSF(msf, &imsf));
4499566063dSJacob Faibussowitsch     PetscCall(PetscSFSetUp(msf));
4509566063dSJacob Faibussowitsch     PetscCall(PetscSFSetUp(imsf));
4519566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Multi-SF\n"));
4529566063dSJacob Faibussowitsch     PetscCall(PetscSFView(msf, PETSC_VIEWER_STDOUT_WORLD));
4539566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Multi-SF roots indices in original SF roots numbering\n"));
4549566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeBegin(sf, &degree));
4559566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeEnd(sf, &degree));
4569566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeMultiRootOriginalNumbering(sf, degree, &inedges, &mRootsOrigNumbering));
4579566063dSJacob Faibussowitsch     PetscCall(PetscIntView(inedges, mRootsOrigNumbering, PETSC_VIEWER_STDOUT_WORLD));
4589566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Inverse of Multi-SF\n"));
4599566063dSJacob Faibussowitsch     PetscCall(PetscSFView(imsf, PETSC_VIEWER_STDOUT_WORLD));
4609566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "## Inverse of Multi-SF, original numbering\n"));
4619566063dSJacob Faibussowitsch     PetscCall(PetscSFViewCustomLocals_Private(imsf, mRootsOrigNumbering, PETSC_VIEWER_STDOUT_WORLD));
4629566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&imsf));
4639566063dSJacob Faibussowitsch     PetscCall(PetscFree(mRootsOrigNumbering));
464c4762a1bSJed Brown   }
465c4762a1bSJed Brown 
466c4762a1bSJed Brown   /* Clean storage for star forest. */
4679566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sf));
4689566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
469b122ec5aSJacob Faibussowitsch   return 0;
470c4762a1bSJed Brown }
471c4762a1bSJed Brown 
472c4762a1bSJed Brown /*TEST
473c4762a1bSJed Brown 
474c4762a1bSJed Brown    test:
475c4762a1bSJed Brown       nsize: 4
476c4762a1bSJed Brown       filter: grep -v "type" | grep -v "sort"
477c4762a1bSJed Brown       args: -test_bcast -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}}
478dfd57a17SPierre Jolivet       requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
479c4762a1bSJed Brown 
480c4762a1bSJed Brown    test:
481c4762a1bSJed Brown       suffix: 2
482c4762a1bSJed Brown       nsize: 4
483c4762a1bSJed Brown       filter: grep -v "type" | grep -v "sort"
484c4762a1bSJed Brown       args: -test_reduce -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}}
485dfd57a17SPierre Jolivet       requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
486c4762a1bSJed Brown 
487c4762a1bSJed Brown    test:
488c4762a1bSJed Brown       suffix: 2_basic
489c4762a1bSJed Brown       nsize: 4
490c4762a1bSJed Brown       args: -test_reduce -sf_type basic
491c4762a1bSJed Brown 
492c4762a1bSJed Brown    test:
493c4762a1bSJed Brown       suffix: 3
494c4762a1bSJed Brown       nsize: 4
495c4762a1bSJed Brown       filter: grep -v "type" | grep -v "sort"
496c4762a1bSJed Brown       args: -test_degree -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}}
497dfd57a17SPierre Jolivet       requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
498c4762a1bSJed Brown 
499c4762a1bSJed Brown    test:
500c4762a1bSJed Brown       suffix: 3_basic
501c4762a1bSJed Brown       nsize: 4
502c4762a1bSJed Brown       args: -test_degree -sf_type basic
503c4762a1bSJed Brown 
504c4762a1bSJed Brown    test:
505c4762a1bSJed Brown       suffix: 4
506c4762a1bSJed Brown       nsize: 4
507c4762a1bSJed Brown       filter: grep -v "type" | grep -v "sort"
508c4762a1bSJed Brown       args: -test_gather -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}}
509dfd57a17SPierre Jolivet       requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
510c4762a1bSJed Brown 
511c4762a1bSJed Brown    test:
512c4762a1bSJed Brown       suffix: 4_basic
513c4762a1bSJed Brown       nsize: 4
514c4762a1bSJed Brown       args: -test_gather -sf_type basic
515c4762a1bSJed Brown 
516c4762a1bSJed Brown    test:
517c4762a1bSJed Brown       suffix: 4_stride
518c4762a1bSJed Brown       nsize: 4
519c4762a1bSJed Brown       args: -test_gather -sf_type basic -stride 2
520c4762a1bSJed Brown 
521c4762a1bSJed Brown    test:
522c4762a1bSJed Brown       suffix: 5
523c4762a1bSJed Brown       nsize: 4
524c4762a1bSJed Brown       filter: grep -v "type" | grep -v "sort"
525c4762a1bSJed Brown       args: -test_scatter -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}}
526dfd57a17SPierre Jolivet       requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
527c4762a1bSJed Brown 
528c4762a1bSJed Brown    test:
529c4762a1bSJed Brown       suffix: 5_basic
530c4762a1bSJed Brown       nsize: 4
531c4762a1bSJed Brown       args: -test_scatter -sf_type basic
532c4762a1bSJed Brown 
533c4762a1bSJed Brown    test:
534c4762a1bSJed Brown       suffix: 5_stride
535c4762a1bSJed Brown       nsize: 4
536c4762a1bSJed Brown       args: -test_scatter -sf_type basic -stride 2
537c4762a1bSJed Brown 
538c4762a1bSJed Brown    test:
539c4762a1bSJed Brown       suffix: 6
540c4762a1bSJed Brown       nsize: 4
541c4762a1bSJed Brown       filter: grep -v "type" | grep -v "sort"
542c20d7725SJed Brown       # No -sf_window_flavor dynamic due to bug https://gitlab.com/petsc/petsc/issues/555
543c20d7725SJed Brown       args: -test_embed -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create allocate}}
544dfd57a17SPierre Jolivet       requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
545c4762a1bSJed Brown 
546c4762a1bSJed Brown    test:
547c4762a1bSJed Brown       suffix: 6_basic
548c4762a1bSJed Brown       nsize: 4
549c4762a1bSJed Brown       args: -test_embed -sf_type basic
550c4762a1bSJed Brown 
551c4762a1bSJed Brown    test:
552c4762a1bSJed Brown       suffix: 7
553c4762a1bSJed Brown       nsize: 4
554c4762a1bSJed Brown       filter: grep -v "type" | grep -v "sort"
555c4762a1bSJed Brown       args: -test_invert -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}}
556dfd57a17SPierre Jolivet       requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
557c4762a1bSJed Brown 
558c4762a1bSJed Brown    test:
559c4762a1bSJed Brown       suffix: 7_basic
560c4762a1bSJed Brown       nsize: 4
561c4762a1bSJed Brown       args: -test_invert -sf_type basic
562c4762a1bSJed Brown 
563c4762a1bSJed Brown    test:
564c4762a1bSJed Brown       suffix: basic
565c4762a1bSJed Brown       nsize: 4
566c4762a1bSJed Brown       args: -test_bcast -sf_type basic
567c4762a1bSJed Brown       output_file: output/ex1_1_basic.out
568c4762a1bSJed Brown 
569c4762a1bSJed Brown    test:
570c4762a1bSJed Brown       suffix: bcastop_basic
571c4762a1bSJed Brown       nsize: 4
572c4762a1bSJed Brown       args: -test_bcastop -sf_type basic
573c4762a1bSJed Brown       output_file: output/ex1_bcastop_basic.out
574c4762a1bSJed Brown 
575c4762a1bSJed Brown    test:
576c4762a1bSJed Brown       suffix: 8
577c4762a1bSJed Brown       nsize: 3
578c4762a1bSJed Brown       filter: grep -v "type" | grep -v "sort"
579c4762a1bSJed Brown       args: -test_bcast -test_sf_distribute -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}}
580dfd57a17SPierre Jolivet       requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
581c4762a1bSJed Brown 
582c4762a1bSJed Brown    test:
583c4762a1bSJed Brown       suffix: 8_basic
584c4762a1bSJed Brown       nsize: 3
585c4762a1bSJed Brown       args: -test_bcast -test_sf_distribute -sf_type basic
586c4762a1bSJed Brown 
587c4762a1bSJed Brown    test:
588c4762a1bSJed Brown       suffix: 9_char
589c4762a1bSJed Brown       nsize: 4
590c4762a1bSJed Brown       args: -sf_type basic -test_bcast -test_reduce -test_op max -test_char
591c4762a1bSJed Brown 
592c4762a1bSJed Brown    # Here we do not test -sf_window_flavor dynamic since it is designed for repeated SFs with few different rootdata pointers
593c4762a1bSJed Brown    test:
594c4762a1bSJed Brown       suffix: 10
595c4762a1bSJed Brown       filter: grep -v "type" | grep -v "sort"
596c4762a1bSJed Brown       nsize: 4
597c4762a1bSJed Brown       args: -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create allocate}} -test_all -test_bcastop 0 -test_fetchandop 0
598dfd57a17SPierre Jolivet       requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
599c4762a1bSJed Brown 
600c4762a1bSJed Brown    # The nightly test suite with MPICH uses ch3:sock, which is broken when winsize == 0 in some of the processes
601c4762a1bSJed Brown    test:
602c4762a1bSJed Brown       suffix: 10_shared
603c4762a1bSJed Brown       output_file: output/ex1_10.out
604c4762a1bSJed Brown       filter: grep -v "type" | grep -v "sort"
605c4762a1bSJed Brown       nsize: 4
606c4762a1bSJed Brown       args: -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor shared -test_all -test_bcastop 0 -test_fetchandop 0
607dfd57a17SPierre 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)
608c4762a1bSJed Brown 
609c4762a1bSJed Brown    test:
610c4762a1bSJed Brown       suffix: 10_basic
611c4762a1bSJed Brown       nsize: 4
612c4762a1bSJed Brown       args: -sf_type basic -test_all -test_bcastop 0 -test_fetchandop 0
613c4762a1bSJed Brown 
614c4762a1bSJed Brown TEST*/
615