xref: /petsc/src/vec/is/sf/tutorials/ex1.c (revision d0609ced746bc51b019815ca91d747429db24893)
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 */
18c4762a1bSJed Brown static PetscErrorCode PetscSFViewCustomLocals_Private(PetscSF sf,const PetscInt locals[],PetscViewer viewer)
19c4762a1bSJed Brown {
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));
31c4762a1bSJed Brown   for (i=0; i<nleaves; i++) {
329566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %" PetscInt_FMT " <- (%" PetscInt_FMT ",%" PetscInt_FMT ")\n",rank,locals[i],iremote[i].rank,iremote[i].index));
33c4762a1bSJed Brown   }
349566063dSJacob Faibussowitsch   PetscCall(PetscViewerFlush(viewer));
359566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPopTab(viewer));
369566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPopSynchronized(viewer));
37c4762a1bSJed Brown   PetscFunctionReturn(0);
38c4762a1bSJed Brown }
39c4762a1bSJed Brown 
40c4762a1bSJed Brown int main(int argc,char **argv)
41c4762a1bSJed Brown {
42c4762a1bSJed Brown   PetscInt       i,nroots,nrootsalloc,nleaves,nleavesalloc,*mine,stride;
43c4762a1bSJed Brown   PetscSFNode    *remote;
44c4762a1bSJed Brown   PetscMPIInt    rank,size;
45c4762a1bSJed Brown   PetscSF        sf;
46c4762a1bSJed 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;
47c4762a1bSJed Brown   MPI_Op         mop=MPI_OP_NULL; /* initialize to prevent compiler warnings with cxx_quad build */
48c4762a1bSJed Brown   char           opstring[256];
49c4762a1bSJed Brown   PetscBool      strflg;
50c4762a1bSJed Brown 
519566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
529566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
539566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
54c4762a1bSJed Brown 
55*d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD,"","PetscSF Test Options","none");
56c4762a1bSJed Brown   test_all        = PETSC_FALSE;
579566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_all","Test all SF communications","",test_all,&test_all,NULL));
58c4762a1bSJed Brown   test_bcast      = test_all;
599566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_bcast","Test broadcast","",test_bcast,&test_bcast,NULL));
60c4762a1bSJed Brown   test_bcastop    = test_all;
619566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_bcastop","Test broadcast and reduce","",test_bcastop,&test_bcastop,NULL));
62c4762a1bSJed Brown   test_reduce     = test_all;
639566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_reduce","Test reduction","",test_reduce,&test_reduce,NULL));
64c4762a1bSJed Brown   test_char       = test_all;
659566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_char","Test signed char, unsigned char, and char","",test_char,&test_char,NULL));
66c4762a1bSJed Brown   mop             = MPI_SUM;
679566063dSJacob Faibussowitsch   PetscCall(PetscStrcpy(opstring,"sum"));
689566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-test_op","Designate which MPI_Op to use","",opstring,opstring,sizeof(opstring),NULL));
699566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp("sum",opstring,&strflg));
70c4762a1bSJed Brown   if (strflg) {
71c4762a1bSJed Brown     mop = MPIU_SUM;
72c4762a1bSJed Brown   }
739566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp("prod",opstring,&strflg));
74c4762a1bSJed Brown   if (strflg) {
75c4762a1bSJed Brown     mop = MPI_PROD;
76c4762a1bSJed Brown   }
779566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp("max",opstring,&strflg));
78c4762a1bSJed Brown   if (strflg) {
79c4762a1bSJed Brown     mop = MPI_MAX;
80c4762a1bSJed Brown   }
819566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp("min",opstring,&strflg));
82c4762a1bSJed Brown   if (strflg) {
83c4762a1bSJed Brown     mop = MPI_MIN;
84c4762a1bSJed Brown   }
859566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp("land",opstring,&strflg));
86c4762a1bSJed Brown   if (strflg) {
87c4762a1bSJed Brown     mop = MPI_LAND;
88c4762a1bSJed Brown   }
899566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp("band",opstring,&strflg));
90c4762a1bSJed Brown   if (strflg) {
91c4762a1bSJed Brown     mop = MPI_BAND;
92c4762a1bSJed Brown   }
939566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp("lor",opstring,&strflg));
94c4762a1bSJed Brown   if (strflg) {
95c4762a1bSJed Brown     mop = MPI_LOR;
96c4762a1bSJed Brown   }
979566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp("bor",opstring,&strflg));
98c4762a1bSJed Brown   if (strflg) {
99c4762a1bSJed Brown     mop = MPI_BOR;
100c4762a1bSJed Brown   }
1019566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp("lxor",opstring,&strflg));
102c4762a1bSJed Brown   if (strflg) {
103c4762a1bSJed Brown     mop = MPI_LXOR;
104c4762a1bSJed Brown   }
1059566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp("bxor",opstring,&strflg));
106c4762a1bSJed Brown   if (strflg) {
107c4762a1bSJed Brown     mop = MPI_BXOR;
108c4762a1bSJed Brown   }
109c4762a1bSJed Brown   test_degree     = test_all;
1109566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_degree","Test computation of vertex degree","",test_degree,&test_degree,NULL));
111c4762a1bSJed Brown   test_fetchandop = test_all;
1129566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_fetchandop","Test atomic Fetch-And-Op","",test_fetchandop,&test_fetchandop,NULL));
113c4762a1bSJed Brown   test_gather     = test_all;
1149566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_gather","Test point gather","",test_gather,&test_gather,NULL));
115c4762a1bSJed Brown   test_scatter    = test_all;
1169566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_scatter","Test point scatter","",test_scatter,&test_scatter,NULL));
117c4762a1bSJed Brown   test_embed      = test_all;
1189566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_embed","Test point embed","",test_embed,&test_embed,NULL));
119c4762a1bSJed Brown   test_invert     = test_all;
1209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_invert","Test point invert","",test_invert,&test_invert,NULL));
121c4762a1bSJed Brown   stride          = 1;
1229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-stride","Stride for leaf and root data","",stride,&stride,NULL));
123c4762a1bSJed Brown   test_sf_distribute = PETSC_FALSE;
1249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_sf_distribute","Create an SF that 'distributes' to each process, like an alltoall","",test_sf_distribute,&test_sf_distribute,NULL));
1259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-test_op","Designate which MPI_Op to use","",opstring,opstring,sizeof(opstring),NULL));
126*d0609cedSBarry Smith   PetscOptionsEnd();
127c4762a1bSJed Brown 
128c4762a1bSJed Brown   if (test_sf_distribute) {
129c4762a1bSJed Brown     nroots = size;
130c4762a1bSJed Brown     nrootsalloc = size;
131c4762a1bSJed Brown     nleaves = size;
132c4762a1bSJed Brown     nleavesalloc = size;
133c4762a1bSJed Brown     mine = NULL;
1349566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nleaves,&remote));
135c4762a1bSJed Brown     for (i=0; i<size; i++) {
136c4762a1bSJed Brown       remote[i].rank = i;
137c4762a1bSJed Brown       remote[i].index = rank;
138c4762a1bSJed Brown     }
139c4762a1bSJed Brown   } else {
140c4762a1bSJed Brown     nroots       = 2 + (PetscInt)(rank == 0);
141c4762a1bSJed Brown     nrootsalloc  = nroots * stride;
142c4762a1bSJed Brown     nleaves      = 2 + (PetscInt)(rank > 0);
143c4762a1bSJed Brown     nleavesalloc = nleaves * stride;
144c4762a1bSJed Brown     mine         = NULL;
145c4762a1bSJed Brown     if (stride > 1) {
146c4762a1bSJed Brown       PetscInt i;
147c4762a1bSJed Brown 
1489566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nleaves,&mine));
149c4762a1bSJed Brown       for (i = 0; i < nleaves; i++) {
150c4762a1bSJed Brown         mine[i] = stride * i;
151c4762a1bSJed Brown       }
152c4762a1bSJed Brown     }
1539566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nleaves,&remote));
154c4762a1bSJed Brown     /* Left periodic neighbor */
155c4762a1bSJed Brown     remote[0].rank  = (rank+size-1)%size;
156c4762a1bSJed Brown     remote[0].index = 1 * stride;
157c4762a1bSJed Brown     /* Right periodic neighbor */
158c4762a1bSJed Brown     remote[1].rank  = (rank+1)%size;
159c4762a1bSJed Brown     remote[1].index = 0 * stride;
160c4762a1bSJed Brown     if (rank > 0) {               /* All processes reference rank 0, index 1 */
161c4762a1bSJed Brown       remote[2].rank  = 0;
162c4762a1bSJed Brown       remote[2].index = 2 * stride;
163c4762a1bSJed Brown     }
164c4762a1bSJed Brown   }
165c4762a1bSJed Brown 
166c4762a1bSJed Brown   /* Create a star forest for communication. In this example, the leaf space is dense, so we pass NULL. */
1679566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(PETSC_COMM_WORLD,&sf));
1689566063dSJacob Faibussowitsch   PetscCall(PetscSFSetFromOptions(sf));
1699566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(sf,nrootsalloc,nleaves,mine,PETSC_OWN_POINTER,remote,PETSC_OWN_POINTER));
1709566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sf));
171c4762a1bSJed Brown 
172c4762a1bSJed Brown   /* View graph, mostly useful for debugging purposes. */
1739566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
1749566063dSJacob Faibussowitsch   PetscCall(PetscSFView(sf,PETSC_VIEWER_STDOUT_WORLD));
1759566063dSJacob Faibussowitsch   PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
176c4762a1bSJed Brown 
177c4762a1bSJed Brown   if (test_bcast) {             /* broadcast rootdata into leafdata */
178c4762a1bSJed Brown     PetscInt *rootdata,*leafdata;
1792d4ee042Sprj-     /* Allocate space for send and receive buffers. This example communicates PetscInt, but other types, including
180c4762a1bSJed Brown      * user-defined structures, could also be used. */
1819566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nrootsalloc,&rootdata,nleavesalloc,&leafdata));
182c4762a1bSJed Brown     /* Set rootdata buffer to be broadcast */
183c4762a1bSJed Brown     for (i=0; i<nrootsalloc; i++) rootdata[i] = -1;
184c4762a1bSJed Brown     for (i=0; i<nroots; i++) rootdata[i*stride] = 100*(rank+1) + i;
185c4762a1bSJed Brown     /* Initialize local buffer, these values are never used. */
186c4762a1bSJed Brown     for (i=0; i<nleavesalloc; i++) leafdata[i] = -1;
187c4762a1bSJed Brown     /* Broadcast entries from rootdata to leafdata. Computation or other communication can be performed between the begin and end calls. */
1889566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(sf,MPIU_INT,rootdata,leafdata,MPI_REPLACE));
1899566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(sf,MPIU_INT,rootdata,leafdata,MPI_REPLACE));
1909566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Bcast Rootdata\n"));
1919566063dSJacob Faibussowitsch     PetscCall(PetscIntView(nrootsalloc,rootdata,PETSC_VIEWER_STDOUT_WORLD));
1929566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Bcast Leafdata\n"));
1939566063dSJacob Faibussowitsch     PetscCall(PetscIntView(nleavesalloc,leafdata,PETSC_VIEWER_STDOUT_WORLD));
1949566063dSJacob Faibussowitsch     PetscCall(PetscFree2(rootdata,leafdata));
195c4762a1bSJed Brown   }
196c4762a1bSJed Brown 
197c4762a1bSJed Brown   if (test_bcast && test_char) { /* Bcast with char */
198c4762a1bSJed Brown     PetscInt len;
199c4762a1bSJed Brown     char buf[256];
200c4762a1bSJed Brown     char *rootdata,*leafdata;
2019566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nrootsalloc,&rootdata,nleavesalloc,&leafdata));
202c4762a1bSJed Brown     /* Set rootdata buffer to be broadcast */
203c4762a1bSJed Brown     for (i=0; i<nrootsalloc; i++) rootdata[i] = '*';
204c4762a1bSJed 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 */
205c4762a1bSJed Brown     /* Initialize local buffer, these values are never used. */
206c4762a1bSJed Brown     for (i=0; i<nleavesalloc; i++) leafdata[i] = '?';
207c4762a1bSJed Brown 
2089566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(sf,MPI_CHAR,rootdata,leafdata,MPI_REPLACE));
2099566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(sf,MPI_CHAR,rootdata,leafdata,MPI_REPLACE));
210c4762a1bSJed Brown 
2119566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Bcast Rootdata in type of char\n"));
2129566063dSJacob Faibussowitsch     len  = 0; PetscCall(PetscSNPrintf(buf,256,"%4d:",rank)); len += 5;
2139566063dSJacob Faibussowitsch     for (i=0; i<nrootsalloc; i++) {PetscCall(PetscSNPrintf(buf+len,256-len,"%5c",rootdata[i])); len += 5;}
2149566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"%s\n",buf));
2159566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT));
216c4762a1bSJed Brown 
2179566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Bcast Leafdata in type of char\n"));
2189566063dSJacob Faibussowitsch     len = 0; PetscCall(PetscSNPrintf(buf,256,"%4d:",rank)); len += 5;
2199566063dSJacob Faibussowitsch     for (i=0; i<nleavesalloc; i++) {PetscCall(PetscSNPrintf(buf+len,256-len,"%5c",leafdata[i])); len += 5;}
2209566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"%s\n",buf));
2219566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT));
222c4762a1bSJed Brown 
2239566063dSJacob Faibussowitsch     PetscCall(PetscFree2(rootdata,leafdata));
224c4762a1bSJed Brown   }
225c4762a1bSJed Brown 
226c4762a1bSJed Brown   if (test_bcastop) {         /* Reduce rootdata into leafdata */
227c4762a1bSJed Brown     PetscInt *rootdata,*leafdata;
228c4762a1bSJed Brown     /* Allocate space for send and receive buffers. This example communicates PetscInt, but other types, including
229c4762a1bSJed Brown      * user-defined structures, could also be used. */
2309566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nrootsalloc,&rootdata,nleavesalloc,&leafdata));
231c4762a1bSJed Brown     /* Set rootdata buffer to be broadcast */
232c4762a1bSJed Brown     for (i=0; i<nrootsalloc; i++) rootdata[i] = -1;
233c4762a1bSJed Brown     for (i=0; i<nroots; i++) rootdata[i*stride] = 100*(rank+1) + i;
234c4762a1bSJed Brown     /* Set leaf values to reduce with */
235c4762a1bSJed Brown     for (i=0; i<nleavesalloc; i++) leafdata[i] = -10*(rank+1) - i;
2369566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Pre-BcastAndOp Leafdata\n"));
2379566063dSJacob Faibussowitsch     PetscCall(PetscIntView(nleavesalloc,leafdata,PETSC_VIEWER_STDOUT_WORLD));
238c4762a1bSJed Brown     /* Broadcast entries from rootdata to leafdata. Computation or other communication can be performed between the begin and end calls. */
2399566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(sf,MPIU_INT,rootdata,leafdata,mop));
2409566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(sf,MPIU_INT,rootdata,leafdata,mop));
2419566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## BcastAndOp Rootdata\n"));
2429566063dSJacob Faibussowitsch     PetscCall(PetscIntView(nrootsalloc,rootdata,PETSC_VIEWER_STDOUT_WORLD));
2439566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## BcastAndOp Leafdata\n"));
2449566063dSJacob Faibussowitsch     PetscCall(PetscIntView(nleavesalloc,leafdata,PETSC_VIEWER_STDOUT_WORLD));
2459566063dSJacob Faibussowitsch     PetscCall(PetscFree2(rootdata,leafdata));
246c4762a1bSJed Brown   }
247c4762a1bSJed Brown 
248c4762a1bSJed Brown   if (test_reduce) {            /* Reduce leafdata into rootdata */
249c4762a1bSJed Brown     PetscInt *rootdata,*leafdata;
2509566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nrootsalloc,&rootdata,nleavesalloc,&leafdata));
251c4762a1bSJed Brown     /* Initialize rootdata buffer in which the result of the reduction will appear. */
252c4762a1bSJed Brown     for (i=0; i<nrootsalloc; i++) rootdata[i] = -1;
253c4762a1bSJed Brown     for (i=0; i<nroots; i++) rootdata[i*stride] = 100*(rank+1) + i;
254c4762a1bSJed Brown     /* Set leaf values to reduce. */
255c4762a1bSJed Brown     for (i=0; i<nleavesalloc; i++) leafdata[i] = -1;
256c4762a1bSJed Brown     for (i=0; i<nleaves; i++) leafdata[i*stride] = 1000*(rank+1) + 10*i;
2579566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Pre-Reduce Rootdata\n"));
2589566063dSJacob Faibussowitsch     PetscCall(PetscIntView(nrootsalloc,rootdata,PETSC_VIEWER_STDOUT_WORLD));
259c4762a1bSJed Brown     /* Perform reduction. Computation or other communication can be performed between the begin and end calls.
260c4762a1bSJed Brown      * This example sums the values, but other MPI_Ops can be used (e.g MPI_MAX, MPI_PROD). */
2619566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceBegin(sf,MPIU_INT,leafdata,rootdata,mop));
2629566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceEnd(sf,MPIU_INT,leafdata,rootdata,mop));
2639566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Reduce Leafdata\n"));
2649566063dSJacob Faibussowitsch     PetscCall(PetscIntView(nleavesalloc,leafdata,PETSC_VIEWER_STDOUT_WORLD));
2659566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Reduce Rootdata\n"));
2669566063dSJacob Faibussowitsch     PetscCall(PetscIntView(nrootsalloc,rootdata,PETSC_VIEWER_STDOUT_WORLD));
2679566063dSJacob Faibussowitsch     PetscCall(PetscFree2(rootdata,leafdata));
268c4762a1bSJed Brown   }
269c4762a1bSJed Brown 
270c4762a1bSJed Brown   if (test_reduce && test_char) { /* Reduce with signed char */
271c4762a1bSJed Brown     PetscInt len;
272c4762a1bSJed Brown     char buf[256];
273c4762a1bSJed Brown     signed char *rootdata,*leafdata;
2749566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nrootsalloc,&rootdata,nleavesalloc,&leafdata));
275c4762a1bSJed Brown     /* Initialize rootdata buffer in which the result of the reduction will appear. */
276c4762a1bSJed Brown     for (i=0; i<nrootsalloc; i++) rootdata[i] = -1;
277c4762a1bSJed Brown     for (i=0; i<nroots; i++) rootdata[i*stride] = 10*(rank+1) + i;
278c4762a1bSJed Brown     /* Set leaf values to reduce. */
279c4762a1bSJed Brown     for (i=0; i<nleavesalloc; i++) leafdata[i] = -1;
280c4762a1bSJed Brown     for (i=0; i<nleaves; i++) leafdata[i*stride] = 50*(rank+1) + 10*i;
2819566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Pre-Reduce Rootdata in type of signed char\n"));
282c4762a1bSJed Brown 
2839566063dSJacob Faibussowitsch     len = 0; PetscCall(PetscSNPrintf(buf,256,"%4d:",rank)); len += 5;
2849566063dSJacob Faibussowitsch     for (i=0; i<nrootsalloc; i++) {PetscCall(PetscSNPrintf(buf+len,256-len,"%5d",rootdata[i])); len += 5;}
2859566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"%s\n",buf));
2869566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT));
287c4762a1bSJed Brown 
288c4762a1bSJed Brown     /* Using MPI_CHAR should trigger an error since MPI standard does not support reduction on MPI_CHAR.
289c4762a1bSJed Brown        Testing with -test_op max, one can see the sign does take effect in MPI_MAX.
290c4762a1bSJed Brown      */
2919566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceBegin(sf,MPI_SIGNED_CHAR,leafdata,rootdata,mop));
2929566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceEnd(sf,MPI_SIGNED_CHAR,leafdata,rootdata,mop));
293c4762a1bSJed Brown 
2949566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Reduce Leafdata in type of signed char\n"));
2959566063dSJacob Faibussowitsch     len  = 0; PetscCall(PetscSNPrintf(buf,256,"%4d:",rank)); len += 5;
2969566063dSJacob Faibussowitsch     for (i=0; i<nleavesalloc; i++) {PetscCall(PetscSNPrintf(buf+len,256-len,"%5d",leafdata[i])); len += 5;}
2979566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"%s\n",buf));
2989566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT));
299c4762a1bSJed Brown 
3009566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Reduce Rootdata in type of signed char\n"));
3019566063dSJacob Faibussowitsch     len = 0; PetscCall(PetscSNPrintf(buf,256,"%4d:",rank)); len += 5;
3029566063dSJacob Faibussowitsch     for (i=0; i<nrootsalloc; i++) {PetscCall(PetscSNPrintf(buf+len,256-len,"%5d",rootdata[i])); len += 5;}
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 
3229566063dSJacob Faibussowitsch     len = 0; PetscCall(PetscSNPrintf(buf,256,"%4d:",rank)); len += 5;
3239566063dSJacob Faibussowitsch     for (i=0; i<nrootsalloc; i++) {PetscCall(PetscSNPrintf(buf+len,256-len,"%5u",rootdata[i])); len += 5;}
3249566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"%s\n",buf));
3259566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT));
326c4762a1bSJed Brown 
327c4762a1bSJed Brown     /* Using MPI_CHAR should trigger an error since MPI standard does not support reduction on MPI_CHAR.
328c4762a1bSJed Brown        Testing with -test_op max, one can see the sign does take effect in MPI_MAX.
329c4762a1bSJed Brown      */
3309566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceBegin(sf,MPI_UNSIGNED_CHAR,leafdata,rootdata,mop));
3319566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceEnd(sf,MPI_UNSIGNED_CHAR,leafdata,rootdata,mop));
332c4762a1bSJed Brown 
3339566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Reduce Leafdata in type of unsigned char\n"));
3349566063dSJacob Faibussowitsch     len  = 0; PetscCall(PetscSNPrintf(buf,256,"%4d:",rank)); len += 5;
3359566063dSJacob Faibussowitsch     for (i=0; i<nleavesalloc; i++) {PetscCall(PetscSNPrintf(buf+len,256-len,"%5u",leafdata[i])); len += 5;}
3369566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"%s\n",buf));
3379566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT));
338c4762a1bSJed Brown 
3399566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Reduce Rootdata in type of unsigned char\n"));
3409566063dSJacob Faibussowitsch     len = 0; PetscCall(PetscSNPrintf(buf,256,"%4d:",rank)); len += 5;
3419566063dSJacob Faibussowitsch     for (i=0; i<nrootsalloc; i++) {PetscCall(PetscSNPrintf(buf+len,256-len,"%5u",rootdata[i])); len += 5;}
3429566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"%s\n",buf));
3439566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT));
344c4762a1bSJed Brown 
3459566063dSJacob Faibussowitsch     PetscCall(PetscFree2(rootdata,leafdata));
346c4762a1bSJed Brown   }
347c4762a1bSJed Brown 
348c4762a1bSJed Brown   if (test_degree) {
349c4762a1bSJed Brown     const PetscInt *degree;
3509566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeBegin(sf,&degree));
3519566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeEnd(sf,&degree));
3529566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Root degrees\n"));
3539566063dSJacob Faibussowitsch     PetscCall(PetscIntView(nrootsalloc,degree,PETSC_VIEWER_STDOUT_WORLD));
354c4762a1bSJed Brown   }
355c4762a1bSJed Brown 
356c4762a1bSJed Brown   if (test_fetchandop) {
357c4762a1bSJed Brown     /* Cannot use text compare here because token ordering is not deterministic */
358c4762a1bSJed Brown     PetscInt *leafdata,*leafupdate,*rootdata;
3599566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(nleavesalloc,&leafdata,nleavesalloc,&leafupdate,nrootsalloc,&rootdata));
360c4762a1bSJed Brown     for (i=0; i<nleavesalloc; i++) leafdata[i] = -1;
361c4762a1bSJed Brown     for (i=0; i<nleaves; i++) leafdata[i*stride] = 1;
362c4762a1bSJed Brown     for (i=0; i<nrootsalloc; i++) rootdata[i] = -1;
363c4762a1bSJed Brown     for (i=0; i<nroots; i++) rootdata[i*stride] = 0;
3649566063dSJacob Faibussowitsch     PetscCall(PetscSFFetchAndOpBegin(sf,MPIU_INT,rootdata,leafdata,leafupdate,mop));
3659566063dSJacob Faibussowitsch     PetscCall(PetscSFFetchAndOpEnd(sf,MPIU_INT,rootdata,leafdata,leafupdate,mop));
3669566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Rootdata (sum of 1 from each leaf)\n"));
3679566063dSJacob Faibussowitsch     PetscCall(PetscIntView(nrootsalloc,rootdata,PETSC_VIEWER_STDOUT_WORLD));
3689566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Leafupdate (value at roots prior to my atomic update)\n"));
3699566063dSJacob Faibussowitsch     PetscCall(PetscIntView(nleavesalloc,leafupdate,PETSC_VIEWER_STDOUT_WORLD));
3709566063dSJacob Faibussowitsch     PetscCall(PetscFree3(leafdata,leafupdate,rootdata));
371c4762a1bSJed Brown   }
372c4762a1bSJed Brown 
373c4762a1bSJed Brown   if (test_gather) {
374c4762a1bSJed Brown     const PetscInt *degree;
375c4762a1bSJed Brown     PetscInt       inedges,*indata,*outdata;
3769566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeBegin(sf,&degree));
3779566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeEnd(sf,&degree));
378c4762a1bSJed Brown     for (i=0,inedges=0; i<nrootsalloc; i++) inedges += degree[i];
3799566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(inedges,&indata,nleavesalloc,&outdata));
380c4762a1bSJed Brown     for (i=0; i<nleavesalloc; i++) outdata[i] = -1;
381c4762a1bSJed Brown     for (i=0; i<nleaves; i++) outdata[i*stride] = 1000*(rank+1) + i;
3829566063dSJacob Faibussowitsch     PetscCall(PetscSFGatherBegin(sf,MPIU_INT,outdata,indata));
3839566063dSJacob Faibussowitsch     PetscCall(PetscSFGatherEnd(sf,MPIU_INT,outdata,indata));
3849566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Gathered data at multi-roots from leaves\n"));
3859566063dSJacob Faibussowitsch     PetscCall(PetscIntView(inedges,indata,PETSC_VIEWER_STDOUT_WORLD));
3869566063dSJacob Faibussowitsch     PetscCall(PetscFree2(indata,outdata));
387c4762a1bSJed Brown   }
388c4762a1bSJed Brown 
389c4762a1bSJed Brown   if (test_scatter) {
390c4762a1bSJed Brown     const PetscInt *degree;
391c4762a1bSJed Brown     PetscInt       j,count,inedges,*indata,*outdata;
3929566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeBegin(sf,&degree));
3939566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeEnd(sf,&degree));
394c4762a1bSJed Brown     for (i=0,inedges=0; i<nrootsalloc; i++) inedges += degree[i];
3959566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(inedges,&indata,nleavesalloc,&outdata));
396c4762a1bSJed Brown     for (i=0; i<nleavesalloc; i++) outdata[i] = -1;
397c4762a1bSJed Brown     for (i=0,count=0; i<nrootsalloc; i++) {
398c4762a1bSJed Brown       for (j=0; j<degree[i]; j++) indata[count++] = 1000*(rank+1) + 100*i + j;
399c4762a1bSJed Brown     }
4009566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Data at multi-roots, to scatter to leaves\n"));
4019566063dSJacob Faibussowitsch     PetscCall(PetscIntView(inedges,indata,PETSC_VIEWER_STDOUT_WORLD));
402c4762a1bSJed Brown 
4039566063dSJacob Faibussowitsch     PetscCall(PetscSFScatterBegin(sf,MPIU_INT,indata,outdata));
4049566063dSJacob Faibussowitsch     PetscCall(PetscSFScatterEnd(sf,MPIU_INT,indata,outdata));
4059566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Scattered data at leaves\n"));
4069566063dSJacob Faibussowitsch     PetscCall(PetscIntView(nleavesalloc,outdata,PETSC_VIEWER_STDOUT_WORLD));
4079566063dSJacob Faibussowitsch     PetscCall(PetscFree2(indata,outdata));
408c4762a1bSJed Brown   }
409c4762a1bSJed Brown 
410c4762a1bSJed Brown   if (test_embed) {
411dd400576SPatrick Sanan     const PetscInt nroots = 1 + (PetscInt) (rank == 0);
412c4762a1bSJed Brown     PetscInt       selected[2];
413c4762a1bSJed Brown     PetscSF        esf;
414c4762a1bSJed Brown 
415c4762a1bSJed Brown     selected[0] = stride;
416c4762a1bSJed Brown     selected[1] = 2*stride;
4179566063dSJacob Faibussowitsch     PetscCall(PetscSFCreateEmbeddedRootSF(sf,nroots,selected,&esf));
4189566063dSJacob Faibussowitsch     PetscCall(PetscSFSetUp(esf));
4199566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Embedded PetscSF\n"));
4209566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
4219566063dSJacob Faibussowitsch     PetscCall(PetscSFView(esf,PETSC_VIEWER_STDOUT_WORLD));
4229566063dSJacob Faibussowitsch     PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
4239566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&esf));
424c4762a1bSJed Brown   }
425c4762a1bSJed Brown 
426c4762a1bSJed Brown   if (test_invert) {
427c4762a1bSJed Brown     const PetscInt *degree;
428c4762a1bSJed Brown     PetscInt *mRootsOrigNumbering;
429c4762a1bSJed Brown     PetscInt inedges;
430c4762a1bSJed Brown     PetscSF msf,imsf;
431c4762a1bSJed Brown 
4329566063dSJacob Faibussowitsch     PetscCall(PetscSFGetMultiSF(sf,&msf));
4339566063dSJacob Faibussowitsch     PetscCall(PetscSFCreateInverseSF(msf,&imsf));
4349566063dSJacob Faibussowitsch     PetscCall(PetscSFSetUp(msf));
4359566063dSJacob Faibussowitsch     PetscCall(PetscSFSetUp(imsf));
4369566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Multi-SF\n"));
4379566063dSJacob Faibussowitsch     PetscCall(PetscSFView(msf,PETSC_VIEWER_STDOUT_WORLD));
4389566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Multi-SF roots indices in original SF roots numbering\n"));
4399566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeBegin(sf,&degree));
4409566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeEnd(sf,&degree));
4419566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeMultiRootOriginalNumbering(sf,degree,&inedges,&mRootsOrigNumbering));
4429566063dSJacob Faibussowitsch     PetscCall(PetscIntView(inedges,mRootsOrigNumbering,PETSC_VIEWER_STDOUT_WORLD));
4439566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Inverse of Multi-SF\n"));
4449566063dSJacob Faibussowitsch     PetscCall(PetscSFView(imsf,PETSC_VIEWER_STDOUT_WORLD));
4459566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Inverse of Multi-SF, original numbering\n"));
4469566063dSJacob Faibussowitsch     PetscCall(PetscSFViewCustomLocals_Private(imsf,mRootsOrigNumbering,PETSC_VIEWER_STDOUT_WORLD));
4479566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&imsf));
4489566063dSJacob Faibussowitsch     PetscCall(PetscFree(mRootsOrigNumbering));
449c4762a1bSJed Brown   }
450c4762a1bSJed Brown 
451c4762a1bSJed Brown   /* Clean storage for star forest. */
4529566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sf));
4539566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
454b122ec5aSJacob Faibussowitsch   return 0;
455c4762a1bSJed Brown }
456c4762a1bSJed Brown 
457c4762a1bSJed Brown /*TEST
458c4762a1bSJed Brown 
459c4762a1bSJed Brown    test:
460c4762a1bSJed Brown       nsize: 4
461c4762a1bSJed Brown       filter: grep -v "type" | grep -v "sort"
462c4762a1bSJed Brown       args: -test_bcast -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}}
463dfd57a17SPierre Jolivet       requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
464c4762a1bSJed Brown 
465c4762a1bSJed Brown    test:
466c4762a1bSJed Brown       suffix: 2
467c4762a1bSJed Brown       nsize: 4
468c4762a1bSJed Brown       filter: grep -v "type" | grep -v "sort"
469c4762a1bSJed Brown       args: -test_reduce -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}}
470dfd57a17SPierre Jolivet       requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
471c4762a1bSJed Brown 
472c4762a1bSJed Brown    test:
473c4762a1bSJed Brown       suffix: 2_basic
474c4762a1bSJed Brown       nsize: 4
475c4762a1bSJed Brown       args: -test_reduce -sf_type basic
476c4762a1bSJed Brown 
477c4762a1bSJed Brown    test:
478c4762a1bSJed Brown       suffix: 3
479c4762a1bSJed Brown       nsize: 4
480c4762a1bSJed Brown       filter: grep -v "type" | grep -v "sort"
481c4762a1bSJed Brown       args: -test_degree -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}}
482dfd57a17SPierre Jolivet       requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
483c4762a1bSJed Brown 
484c4762a1bSJed Brown    test:
485c4762a1bSJed Brown       suffix: 3_basic
486c4762a1bSJed Brown       nsize: 4
487c4762a1bSJed Brown       args: -test_degree -sf_type basic
488c4762a1bSJed Brown 
489c4762a1bSJed Brown    test:
490c4762a1bSJed Brown       suffix: 4
491c4762a1bSJed Brown       nsize: 4
492c4762a1bSJed Brown       filter: grep -v "type" | grep -v "sort"
493c4762a1bSJed Brown       args: -test_gather -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}}
494dfd57a17SPierre Jolivet       requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
495c4762a1bSJed Brown 
496c4762a1bSJed Brown    test:
497c4762a1bSJed Brown       suffix: 4_basic
498c4762a1bSJed Brown       nsize: 4
499c4762a1bSJed Brown       args: -test_gather -sf_type basic
500c4762a1bSJed Brown 
501c4762a1bSJed Brown    test:
502c4762a1bSJed Brown       suffix: 4_stride
503c4762a1bSJed Brown       nsize: 4
504c4762a1bSJed Brown       args: -test_gather -sf_type basic -stride 2
505c4762a1bSJed Brown 
506c4762a1bSJed Brown    test:
507c4762a1bSJed Brown       suffix: 5
508c4762a1bSJed Brown       nsize: 4
509c4762a1bSJed Brown       filter: grep -v "type" | grep -v "sort"
510c4762a1bSJed Brown       args: -test_scatter -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: 5_basic
515c4762a1bSJed Brown       nsize: 4
516c4762a1bSJed Brown       args: -test_scatter -sf_type basic
517c4762a1bSJed Brown 
518c4762a1bSJed Brown    test:
519c4762a1bSJed Brown       suffix: 5_stride
520c4762a1bSJed Brown       nsize: 4
521c4762a1bSJed Brown       args: -test_scatter -sf_type basic -stride 2
522c4762a1bSJed Brown 
523c4762a1bSJed Brown    test:
524c4762a1bSJed Brown       suffix: 6
525c4762a1bSJed Brown       nsize: 4
526c4762a1bSJed Brown       filter: grep -v "type" | grep -v "sort"
527c20d7725SJed Brown       # No -sf_window_flavor dynamic due to bug https://gitlab.com/petsc/petsc/issues/555
528c20d7725SJed Brown       args: -test_embed -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create allocate}}
529dfd57a17SPierre Jolivet       requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
530c4762a1bSJed Brown 
531c4762a1bSJed Brown    test:
532c4762a1bSJed Brown       suffix: 6_basic
533c4762a1bSJed Brown       nsize: 4
534c4762a1bSJed Brown       args: -test_embed -sf_type basic
535c4762a1bSJed Brown 
536c4762a1bSJed Brown    test:
537c4762a1bSJed Brown       suffix: 7
538c4762a1bSJed Brown       nsize: 4
539c4762a1bSJed Brown       filter: grep -v "type" | grep -v "sort"
540c4762a1bSJed Brown       args: -test_invert -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}}
541dfd57a17SPierre Jolivet       requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
542c4762a1bSJed Brown 
543c4762a1bSJed Brown    test:
544c4762a1bSJed Brown       suffix: 7_basic
545c4762a1bSJed Brown       nsize: 4
546c4762a1bSJed Brown       args: -test_invert -sf_type basic
547c4762a1bSJed Brown 
548c4762a1bSJed Brown    test:
549c4762a1bSJed Brown       suffix: basic
550c4762a1bSJed Brown       nsize: 4
551c4762a1bSJed Brown       args: -test_bcast -sf_type basic
552c4762a1bSJed Brown       output_file: output/ex1_1_basic.out
553c4762a1bSJed Brown 
554c4762a1bSJed Brown    test:
555c4762a1bSJed Brown       suffix: bcastop_basic
556c4762a1bSJed Brown       nsize: 4
557c4762a1bSJed Brown       args: -test_bcastop -sf_type basic
558c4762a1bSJed Brown       output_file: output/ex1_bcastop_basic.out
559c4762a1bSJed Brown 
560c4762a1bSJed Brown    test:
561c4762a1bSJed Brown       suffix: 8
562c4762a1bSJed Brown       nsize: 3
563c4762a1bSJed Brown       filter: grep -v "type" | grep -v "sort"
564c4762a1bSJed Brown       args: -test_bcast -test_sf_distribute -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create dynamic allocate}}
565dfd57a17SPierre Jolivet       requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
566c4762a1bSJed Brown 
567c4762a1bSJed Brown    test:
568c4762a1bSJed Brown       suffix: 8_basic
569c4762a1bSJed Brown       nsize: 3
570c4762a1bSJed Brown       args: -test_bcast -test_sf_distribute -sf_type basic
571c4762a1bSJed Brown 
572c4762a1bSJed Brown    test:
573c4762a1bSJed Brown       suffix: 9_char
574c4762a1bSJed Brown       nsize: 4
575c4762a1bSJed Brown       args: -sf_type basic -test_bcast -test_reduce -test_op max -test_char
576c4762a1bSJed Brown 
577c4762a1bSJed Brown    # Here we do not test -sf_window_flavor dynamic since it is designed for repeated SFs with few different rootdata pointers
578c4762a1bSJed Brown    test:
579c4762a1bSJed Brown       suffix: 10
580c4762a1bSJed Brown       filter: grep -v "type" | grep -v "sort"
581c4762a1bSJed Brown       nsize: 4
582c4762a1bSJed Brown       args: -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create allocate}} -test_all -test_bcastop 0 -test_fetchandop 0
583dfd57a17SPierre Jolivet       requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
584c4762a1bSJed Brown 
585c4762a1bSJed Brown    # The nightly test suite with MPICH uses ch3:sock, which is broken when winsize == 0 in some of the processes
586c4762a1bSJed Brown    test:
587c4762a1bSJed Brown       suffix: 10_shared
588c4762a1bSJed Brown       output_file: output/ex1_10.out
589c4762a1bSJed Brown       filter: grep -v "type" | grep -v "sort"
590c4762a1bSJed Brown       nsize: 4
591c4762a1bSJed Brown       args: -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor shared -test_all -test_bcastop 0 -test_fetchandop 0
592dfd57a17SPierre 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)
593c4762a1bSJed Brown 
594c4762a1bSJed Brown    test:
595c4762a1bSJed Brown       suffix: 10_basic
596c4762a1bSJed Brown       nsize: 4
597c4762a1bSJed Brown       args: -sf_type basic -test_all -test_bcastop 0 -test_fetchandop 0
598c4762a1bSJed Brown 
599c4762a1bSJed Brown TEST*/
600