Actual source code: ex1.c

petsc-master 2019-02-16
Report Typos and Errors
  1: static const char help[] = "Test star forest communication (PetscSF)\n\n";

  3: /*T
  4:     Description: A star is a simple tree with one root and zero or more leaves.
  5:     A star forest is a union of disjoint stars.
  6:     Many common communication patterns can be expressed as updates of rootdata using leafdata and vice-versa.
  7:     This example creates a star forest, communicates values using the graph (see options for types of communication), views the graph, then destroys it.
  8: T*/

 10: /*
 11:   Include petscsf.h so we can use PetscSF objects. Note that this automatically
 12:   includes petscsys.h.
 13: */
 14:  #include <petscsf.h>
 15:  #include <petscviewer.h>

 17: /* like PetscSFView() but with alternative array of local indices */
 18: static PetscErrorCode PetscSFViewCustomLocals_Private(PetscSF sf,const PetscInt locals[],PetscViewer viewer)
 19: {
 20:   const PetscSFNode *iremote;
 21:   PetscInt          i,nroots,nleaves,nranks;
 22:   PetscMPIInt       rank;
 23:   PetscErrorCode    ierr;

 26:   MPI_Comm_rank(PetscObjectComm((PetscObject)sf),&rank);
 27:   PetscSFGetGraph(sf,&nroots,&nleaves,NULL,&iremote);
 28:   PetscSFGetRanks(sf,&nranks,NULL,NULL,NULL,NULL);
 29:   PetscViewerASCIIPushTab(viewer);
 30:   PetscViewerASCIIPushSynchronized(viewer);
 31:   PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Number of roots=%D, leaves=%D, remote ranks=%D\n",rank,nroots,nleaves,nranks);
 32:   for (i=0; i<nleaves; i++) {
 33:     PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %D <- (%D,%D)\n",rank,locals[i],iremote[i].rank,iremote[i].index);
 34:   }
 35:   PetscViewerFlush(viewer);
 36:   PetscViewerASCIIPopTab(viewer);
 37:   PetscViewerASCIIPopSynchronized(viewer);
 38:   return(0);
 39: }

 41: int main(int argc,char **argv)
 42: {
 44:   PetscInt       i,nroots,nrootsalloc,nleaves,nleavesalloc,*mine,stride;
 45:   PetscSFNode    *remote;
 46:   PetscMPIInt    rank,size;
 47:   PetscSF        sf;
 48:   PetscBool      test_bcast,test_reduce,test_degree,test_fetchandop,test_gather,test_scatter,test_embed,test_invert,test_sf_distribute;
 49:   MPI_Op         mop=MPI_OP_NULL; /* initialize to prevent compiler warnings with cxx_quad build */
 50:   char           opstring[256];
 51:   PetscBool      strflg;

 53:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 54:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 55:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

 57:   PetscOptionsBegin(PETSC_COMM_WORLD,"","PetscSF Test Options","none");
 58:   test_bcast      = PETSC_FALSE;
 59:   PetscOptionsBool("-test_bcast","Test broadcast","",test_bcast,&test_bcast,NULL);
 60:   test_reduce     = PETSC_FALSE;
 61:   PetscOptionsBool("-test_reduce","Test reduction","",test_reduce,&test_reduce,NULL);
 62:   mop             = MPI_SUM;
 63:   PetscStrcpy(opstring,"sum");
 64:   PetscOptionsString("-test_op","Designate which MPI_Op to use","",opstring,opstring,256,NULL);
 65:   PetscStrcmp("sum",opstring,&strflg);
 66:   if (strflg) {
 67:     mop = MPIU_SUM;
 68:   }
 69:   PetscStrcmp("prod",opstring,&strflg);
 70:   if (strflg) {
 71:     mop = MPI_PROD;
 72:   }
 73:   PetscStrcmp("max",opstring,&strflg);
 74:   if (strflg) {
 75:     mop = MPI_MAX;
 76:   }
 77:   PetscStrcmp("min",opstring,&strflg);
 78:   if (strflg) {
 79:     mop = MPI_MIN;
 80:   }
 81:   PetscStrcmp("land",opstring,&strflg);
 82:   if (strflg) {
 83:     mop = MPI_LAND;
 84:   }
 85:   PetscStrcmp("band",opstring,&strflg);
 86:   if (strflg) {
 87:     mop = MPI_BAND;
 88:   }
 89:   PetscStrcmp("lor",opstring,&strflg);
 90:   if (strflg) {
 91:     mop = MPI_LOR;
 92:   }
 93:   PetscStrcmp("bor",opstring,&strflg);
 94:   if (strflg) {
 95:     mop = MPI_BOR;
 96:   }
 97:   PetscStrcmp("lxor",opstring,&strflg);
 98:   if (strflg) {
 99:     mop = MPI_LXOR;
100:   }
101:   PetscStrcmp("bxor",opstring,&strflg);
102:   if (strflg) {
103:     mop = MPI_BXOR;
104:   }
105:   test_degree     = PETSC_FALSE;
106:   PetscOptionsBool("-test_degree","Test computation of vertex degree","",test_degree,&test_degree,NULL);
107:   test_fetchandop = PETSC_FALSE;
108:   PetscOptionsBool("-test_fetchandop","Test atomic Fetch-And-Op","",test_fetchandop,&test_fetchandop,NULL);
109:   test_gather     = PETSC_FALSE;
110:   PetscOptionsBool("-test_gather","Test point gather","",test_gather,&test_gather,NULL);
111:   test_scatter    = PETSC_FALSE;
112:   PetscOptionsBool("-test_scatter","Test point scatter","",test_scatter,&test_scatter,NULL);
113:   test_embed      = PETSC_FALSE;
114:   PetscOptionsBool("-test_embed","Test point embed","",test_embed,&test_embed,NULL);
115:   test_invert     = PETSC_FALSE;
116:   PetscOptionsBool("-test_invert","Test point invert","",test_invert,&test_invert,NULL);
117:   stride          = 1;
118:   PetscOptionsInt("-stride","Stride for leaf and root data","",stride,&stride,NULL);
119:   test_sf_distribute = PETSC_FALSE;
120:   PetscOptionsBool("-test_sf_distribute","Create an SF that 'distributes' to each process, like an alltoall","",test_sf_distribute,&test_sf_distribute,NULL);
121:   PetscOptionsEnd();

123:   if (test_sf_distribute) {
124:     nroots = size;
125:     nrootsalloc = size;
126:     nleaves = size;
127:     nleavesalloc = size;
128:     mine = NULL;
129:     PetscMalloc1(nleaves,&remote);
130:     for (i=0; i<size; i++) {
131:       remote[i].rank = i;
132:       remote[i].index = rank;
133:     }
134:   } else {
135:     nroots       = 2 + (PetscInt)(rank == 0);
136:     nrootsalloc  = nroots * stride;
137:     nleaves      = 2 + (PetscInt)(rank > 0);
138:     nleavesalloc = nleaves * stride;
139:     mine         = NULL;
140:     if (stride > 1) {
141:       PetscInt i;

143:       PetscMalloc1(nleaves,&mine);
144:       for (i = 0; i < nleaves; i++) {
145:         mine[i] = stride * i;
146:       }
147:     }
148:     PetscMalloc1(nleaves,&remote);
149:     /* Left periodic neighbor */
150:     remote[0].rank  = (rank+size-1)%size;
151:     remote[0].index = 1 * stride;
152:     /* Right periodic neighbor */
153:     remote[1].rank  = (rank+1)%size;
154:     remote[1].index = 0 * stride;
155:     if (rank > 0) {               /* All processes reference rank 0, index 1 */
156:       remote[2].rank  = 0;
157:       remote[2].index = 2 * stride;
158:     }
159:   }

161:   /* Create a star forest for communication. In this example, the leaf space is dense, so we pass NULL. */
162:   PetscSFCreate(PETSC_COMM_WORLD,&sf);
163:   PetscSFSetFromOptions(sf);
164:   PetscSFSetGraph(sf,nrootsalloc,nleaves,mine,PETSC_OWN_POINTER,remote,PETSC_OWN_POINTER);
165:   PetscSFSetUp(sf);

167:   /* View graph, mostly useful for debugging purposes. */
168:   PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
169:   PetscSFView(sf,PETSC_VIEWER_STDOUT_WORLD);
170:   PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);


173:   if (test_bcast) {             /* broadcast rootdata into leafdata */
174:     PetscInt *rootdata,*leafdata;
175:     /* Allocate space for send and recieve buffers. This example communicates PetscInt, but other types, including
176:      * user-defined structures, could also be used. */
177:     PetscMalloc2(nrootsalloc,&rootdata,nleavesalloc,&leafdata);
178:     /* Set rootdata buffer to be broadcast */
179:     for (i=0; i<nrootsalloc; i++) rootdata[i] = -1;
180:     for (i=0; i<nroots; i++) rootdata[i*stride] = 100*(rank+1) + i;
181:     /* Initialize local buffer, these values are never used. */
182:     for (i=0; i<nleavesalloc; i++) leafdata[i] = -1;
183:     /* Broadcast entries from rootdata to leafdata. Computation or other communication can be performed between the begin and end calls. */
184:     PetscSFBcastBegin(sf,MPIU_INT,rootdata,leafdata);
185:     PetscSFBcastEnd(sf,MPIU_INT,rootdata,leafdata);
186:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Bcast Rootdata\n");
187:     PetscIntView(nrootsalloc,rootdata,PETSC_VIEWER_STDOUT_WORLD);
188:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Bcast Leafdata\n");
189:     PetscIntView(nleavesalloc,leafdata,PETSC_VIEWER_STDOUT_WORLD);
190:     PetscFree2(rootdata,leafdata);
191:   }

193:   if (test_reduce) {            /* Reduce leafdata into rootdata */
194:     PetscInt *rootdata,*leafdata;
195:     PetscMalloc2(nrootsalloc,&rootdata,nleavesalloc,&leafdata);
196:     /* Initialize rootdata buffer in which the result of the reduction will appear. */
197:     for (i=0; i<nrootsalloc; i++) rootdata[i] = -1;
198:     for (i=0; i<nroots; i++) rootdata[i*stride] = 100*(rank+1) + i;
199:     /* Set leaf values to reduce. */
200:     for (i=0; i<nleavesalloc; i++) leafdata[i] = -1;
201:     for (i=0; i<nleaves; i++) leafdata[i*stride] = 1000*(rank+1) + 10*i;
202:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Pre-Reduce Rootdata\n");
203:     PetscIntView(nrootsalloc,rootdata,PETSC_VIEWER_STDOUT_WORLD);
204:     /* Perform reduction. Computation or other communication can be performed between the begin and end calls.
205:      * This example sums the values, but other MPI_Ops can be used (e.g MPI_MAX, MPI_PROD). */
206:     PetscSFReduceBegin(sf,MPIU_INT,leafdata,rootdata,mop);
207:     PetscSFReduceEnd(sf,MPIU_INT,leafdata,rootdata,mop);
208:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Reduce Leafdata\n");
209:     PetscIntView(nleavesalloc,leafdata,PETSC_VIEWER_STDOUT_WORLD);
210:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Reduce Rootdata\n");
211:     PetscIntView(nrootsalloc,rootdata,PETSC_VIEWER_STDOUT_WORLD);
212:     PetscFree2(rootdata,leafdata);
213:   }

215:   if (test_degree) {
216:     const PetscInt *degree;
217:     PetscSFComputeDegreeBegin(sf,&degree);
218:     PetscSFComputeDegreeEnd(sf,&degree);
219:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Root degrees\n");
220:     PetscIntView(nrootsalloc,degree,PETSC_VIEWER_STDOUT_WORLD);
221:   }

223:   if (test_fetchandop) {
224:     /* Cannot use text compare here because token ordering is not deterministic */
225:     PetscInt *leafdata,*leafupdate,*rootdata;
226:     PetscMalloc3(nleavesalloc,&leafdata,nleavesalloc,&leafupdate,nrootsalloc,&rootdata);
227:     for (i=0; i<nleavesalloc; i++) leafdata[i] = -1;
228:     for (i=0; i<nleaves; i++) leafdata[i*stride] = 1;
229:     for (i=0; i<nrootsalloc; i++) rootdata[i] = -1;
230:     for (i=0; i<nroots; i++) rootdata[i*stride] = 0;
231:     PetscSFFetchAndOpBegin(sf,MPIU_INT,rootdata,leafdata,leafupdate,mop);
232:     PetscSFFetchAndOpEnd(sf,MPIU_INT,rootdata,leafdata,leafupdate,mop);
233:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Rootdata (sum of 1 from each leaf)\n");
234:     PetscIntView(nrootsalloc,rootdata,PETSC_VIEWER_STDOUT_WORLD);
235:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Leafupdate (value at roots prior to my atomic update)\n");
236:     PetscIntView(nleavesalloc,leafupdate,PETSC_VIEWER_STDOUT_WORLD);
237:     PetscFree3(leafdata,leafupdate,rootdata);
238:   }

240:   if (test_gather) {
241:     const PetscInt *degree;
242:     PetscInt       inedges,*indata,*outdata;
243:     PetscSFComputeDegreeBegin(sf,&degree);
244:     PetscSFComputeDegreeEnd(sf,&degree);
245:     for (i=0,inedges=0; i<nrootsalloc; i++) inedges += degree[i];
246:     PetscMalloc2(inedges,&indata,nleavesalloc,&outdata);
247:     for (i=0; i<nleavesalloc; i++) outdata[i] = -1;
248:     for (i=0; i<nleaves; i++) outdata[i*stride] = 1000*(rank+1) + i;
249:     PetscSFGatherBegin(sf,MPIU_INT,outdata,indata);
250:     PetscSFGatherEnd(sf,MPIU_INT,outdata,indata);
251:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Gathered data at multi-roots from leaves\n");
252:     PetscIntView(inedges,indata,PETSC_VIEWER_STDOUT_WORLD);
253:     PetscFree2(indata,outdata);
254:   }

256:   if (test_scatter) {
257:     const PetscInt *degree;
258:     PetscInt       j,count,inedges,*indata,*outdata;
259:     PetscSFComputeDegreeBegin(sf,&degree);
260:     PetscSFComputeDegreeEnd(sf,&degree);
261:     for (i=0,inedges=0; i<nrootsalloc; i++) inedges += degree[i];
262:     PetscMalloc2(inedges,&indata,nleavesalloc,&outdata);
263:     for (i=0; i<nleavesalloc; i++) outdata[i] = -1;
264:     for (i=0,count=0; i<nrootsalloc; i++) {
265:       for (j=0; j<degree[i]; j++) indata[count++] = 1000*(rank+1) + 100*i + j;
266:     }
267:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Data at multi-roots, to scatter to leaves\n");
268:     PetscIntView(inedges,indata,PETSC_VIEWER_STDOUT_WORLD);

270:     PetscSFScatterBegin(sf,MPIU_INT,indata,outdata);
271:     PetscSFScatterEnd(sf,MPIU_INT,indata,outdata);
272:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Scattered data at leaves\n");
273:     PetscIntView(nleavesalloc,outdata,PETSC_VIEWER_STDOUT_WORLD);
274:     PetscFree2(indata,outdata);
275:   }

277:   if (test_embed) {
278:     const PetscInt nroots = 1 + (PetscInt) !rank;
279:     PetscInt       selected[2];
280:     PetscSF        esf;

282:     selected[0] = stride;
283:     selected[1] = 2*stride;
284:     PetscSFCreateEmbeddedSF(sf,nroots,selected,&esf);
285:     PetscSFSetUp(esf);
286:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Embedded PetscSF\n");
287:     PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
288:     PetscSFView(esf,PETSC_VIEWER_STDOUT_WORLD);
289:     PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
290:     PetscSFDestroy(&esf);
291:   }

293:   if (test_invert) {
294:     const PetscInt *degree;
295:     PetscInt *mRootsOrigNumbering;
296:     PetscInt inedges;
297:     PetscSF msf,imsf;

299:     PetscSFGetMultiSF(sf,&msf);
300:     PetscSFCreateInverseSF(msf,&imsf);
301:     PetscSFSetUp(msf);
302:     PetscSFSetUp(imsf);
303:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Multi-SF\n");
304:     PetscSFView(msf,PETSC_VIEWER_STDOUT_WORLD);
305:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Multi-SF roots indices in original SF roots numbering\n");
306:     PetscSFGetGraph(msf,&inedges,NULL,NULL,NULL);
307:     PetscSFComputeDegreeBegin(sf,&degree);
308:     PetscSFComputeDegreeEnd(sf,&degree);
309:     PetscSFComputeMultiRootOriginalNumbering(sf,degree,&mRootsOrigNumbering);
310:     PetscIntView(inedges,mRootsOrigNumbering,PETSC_VIEWER_STDOUT_WORLD);
311:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Inverse of Multi-SF\n");
312:     PetscSFView(imsf,PETSC_VIEWER_STDOUT_WORLD);
313:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Inverse of Multi-SF, original numbering\n");
314:     PetscSFComputeDegreeBegin(sf,&degree);
315:     PetscSFComputeDegreeEnd(sf,&degree);
316:     PetscSFViewCustomLocals_Private(imsf,mRootsOrigNumbering,PETSC_VIEWER_STDOUT_WORLD);
317:     PetscSFDestroy(&imsf);
318:     PetscFree(mRootsOrigNumbering);
319:   }

321:   /* Clean storage for star forest. */
322:   PetscSFDestroy(&sf);
323:   PetscFinalize();
324:   return ierr;
325: }

327: /*TEST

329:    test:
330:       nsize: 4
331:       args: -test_bcast -sf_type window
332:       requires: define(PETSC_HAVE_MPI_WIN_CREATE) define(PETSC_HAVE_MPICH_NUMVERSION)


335:    test:
336:       suffix: 2
337:       nsize: 4
338:       args: -test_reduce -sf_type window
339:       requires: define(PETSC_HAVE_MPI_WIN_CREATE) define(PETSC_HAVE_MPICH_NUMVERSION)

341:    test:
342:       suffix: 2_basic
343:       nsize: 4
344:       args: -test_reduce -sf_type basic

346:    test:
347:       suffix: 3
348:       nsize: 4
349:       args: -test_degree -sf_type window
350:       requires: define(PETSC_HAVE_MPI_WIN_CREATE) define(PETSC_HAVE_MPICH_NUMVERSION)

352:    test:
353:       suffix: 3_basic
354:       nsize: 4
355:       args: -test_degree -sf_type basic

357:    test:
358:       suffix: 4
359:       nsize: 4
360:       args: -test_gather -sf_type window
361:       requires: define(PETSC_HAVE_MPI_WIN_CREATE) define(PETSC_HAVE_MPICH_NUMVERSION)

363:    test:
364:       suffix: 4_basic
365:       nsize: 4
366:       args: -test_gather -sf_type basic

368:    test:
369:       suffix: 4_stride
370:       nsize: 4
371:       args: -test_gather -sf_type basic -stride 2

373:    test:
374:       suffix: 5
375:       nsize: 4
376:       args: -test_scatter -sf_type window
377:       requires: define(PETSC_HAVE_MPI_WIN_CREATE) define(PETSC_HAVE_MPICH_NUMVERSION)

379:    test:
380:       suffix: 5_basic
381:       nsize: 4
382:       args: -test_scatter -sf_type basic

384:    test:
385:       suffix: 5_stride
386:       nsize: 4
387:       args: -test_scatter -sf_type basic -stride 2

389:    test:
390:       suffix: 6
391:       nsize: 4
392:       args: -test_embed -sf_type window
393:       requires: define(PETSC_HAVE_MPI_WIN_CREATE) define(PETSC_HAVE_MPICH_NUMVERSION)

395:    test:
396:       suffix: 6_basic
397:       nsize: 4
398:       args: -test_embed -sf_type basic

400:    test:
401:       suffix: 7
402:       nsize: 4
403:       args: -test_invert -sf_type window
404:       requires: define(PETSC_HAVE_MPI_WIN_CREATE) define(PETSC_HAVE_MPICH_NUMVERSION)

406:    test:
407:       suffix: 7_basic
408:       nsize: 4
409:       args: -test_invert -sf_type basic

411:    test:
412:       suffix: basic
413:       nsize: 4
414:       args: -test_bcast -sf_type basic
415:       output_file: output/ex1_1_basic.out

417:    test:
418:       suffix: 8
419:       nsize: 3
420:       args: -test_bcast -test_sf_distribute -sf_type window
421:       requires: define(PETSC_HAVE_MPI_WIN_CREATE)

423:    test:
424:       suffix: 8_basic
425:       nsize: 3
426:       args: -test_bcast -test_sf_distribute -sf_type basic

428: TEST*/