Actual source code: ex1.c

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

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

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

 16: int main(int argc,char **argv)
 17: {
 19:   PetscInt       i,nroots,nrootsalloc,nleaves,nleavesalloc,*mine,stride;
 20:   PetscSFNode    *remote;
 21:   PetscMPIInt    rank,size;
 22:   PetscSF        sf;
 23:   PetscBool      test_bcast,test_reduce,test_degree,test_fetchandop,test_gather,test_scatter,test_embed,test_invert;
 24:   MPI_Op         mop=MPI_OP_NULL; /* initialize to prevent compiler warnings with cxx_quad build */
 25:   char           opstring[256];
 26:   PetscBool      strflg;

 28:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 29:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 30:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

 32:   PetscOptionsBegin(PETSC_COMM_WORLD,"","PetscSF Test Options","none");
 33:   test_bcast      = PETSC_FALSE;
 34:   PetscOptionsBool("-test_bcast","Test broadcast","",test_bcast,&test_bcast,NULL);
 35:   test_reduce     = PETSC_FALSE;
 36:   PetscOptionsBool("-test_reduce","Test reduction","",test_reduce,&test_reduce,NULL);
 37:   mop             = MPI_SUM;
 38:   PetscStrcpy(opstring,"sum");
 39:   PetscOptionsString("-test_op","Designate which MPI_Op to use","",opstring,opstring,256,NULL);
 40:   PetscStrcmp("sum",opstring,&strflg);
 41:   if (strflg) {
 42:     mop = MPIU_SUM;
 43:   }
 44:   PetscStrcmp("prod",opstring,&strflg);
 45:   if (strflg) {
 46:     mop = MPI_PROD;
 47:   }
 48:   PetscStrcmp("max",opstring,&strflg);
 49:   if (strflg) {
 50:     mop = MPI_MAX;
 51:   }
 52:   PetscStrcmp("min",opstring,&strflg);
 53:   if (strflg) {
 54:     mop = MPI_MIN;
 55:   }
 56:   PetscStrcmp("land",opstring,&strflg);
 57:   if (strflg) {
 58:     mop = MPI_LAND;
 59:   }
 60:   PetscStrcmp("band",opstring,&strflg);
 61:   if (strflg) {
 62:     mop = MPI_BAND;
 63:   }
 64:   PetscStrcmp("lor",opstring,&strflg);
 65:   if (strflg) {
 66:     mop = MPI_LOR;
 67:   }
 68:   PetscStrcmp("bor",opstring,&strflg);
 69:   if (strflg) {
 70:     mop = MPI_BOR;
 71:   }
 72:   PetscStrcmp("lxor",opstring,&strflg);
 73:   if (strflg) {
 74:     mop = MPI_LXOR;
 75:   }
 76:   PetscStrcmp("bxor",opstring,&strflg);
 77:   if (strflg) {
 78:     mop = MPI_BXOR;
 79:   }
 80:   test_degree     = PETSC_FALSE;
 81:   PetscOptionsBool("-test_degree","Test computation of vertex degree","",test_degree,&test_degree,NULL);
 82:   test_fetchandop = PETSC_FALSE;
 83:   PetscOptionsBool("-test_fetchandop","Test atomic Fetch-And-Op","",test_fetchandop,&test_fetchandop,NULL);
 84:   test_gather     = PETSC_FALSE;
 85:   PetscOptionsBool("-test_gather","Test point gather","",test_gather,&test_gather,NULL);
 86:   test_scatter    = PETSC_FALSE;
 87:   PetscOptionsBool("-test_scatter","Test point scatter","",test_scatter,&test_scatter,NULL);
 88:   test_embed      = PETSC_FALSE;
 89:   PetscOptionsBool("-test_embed","Test point embed","",test_embed,&test_embed,NULL);
 90:   test_invert     = PETSC_FALSE;
 91:   PetscOptionsBool("-test_invert","Test point invert","",test_invert,&test_invert,NULL);
 92:   stride          = 1;
 93:   PetscOptionsInt("-stride","Stride for leaf and root data","",stride,&stride,NULL);
 94:   PetscOptionsEnd();

 96:   nroots       = 2 + (PetscInt)(rank == 0);
 97:   nrootsalloc  = nroots * stride;
 98:   nleaves      = 2 + (PetscInt)(rank > 0);
 99:   nleavesalloc = nleaves * stride;
100:   mine         = NULL;
101:   if (stride > 1) {
102:     PetscInt i;

104:     PetscMalloc1(nleaves,&mine);
105:     for (i = 0; i < nleaves; i++) {
106:       mine[i] = stride * i;
107:     }
108:   }
109:   PetscMalloc1(nleaves,&remote);
110:   /* Left periodic neighbor */
111:   remote[0].rank  = (rank+size-1)%size;
112:   remote[0].index = 1 * stride;
113:   /* Right periodic neighbor */
114:   remote[1].rank  = (rank+1)%size;
115:   remote[1].index = 0 * stride;
116:   if (rank > 0) {               /* All processes reference rank 0, index 1 */
117:     remote[2].rank  = 0;
118:     remote[2].index = 2 * stride;
119:   }

121:   /* Create a star forest for communication. In this example, the leaf space is dense, so we pass NULL. */
122:   PetscSFCreate(PETSC_COMM_WORLD,&sf);
123:   PetscSFSetFromOptions(sf);
124:   PetscSFSetGraph(sf,nrootsalloc,nleaves,mine,PETSC_OWN_POINTER,remote,PETSC_OWN_POINTER);

126:   /* View graph, mostly useful for debugging purposes. */
127:   PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
128:   PetscSFView(sf,PETSC_VIEWER_STDOUT_WORLD);
129:   PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);


132:   if (test_bcast) {             /* broadcast rootdata into leafdata */
133:     PetscInt *rootdata,*leafdata;
134:     /* Allocate space for send and recieve buffers. This example communicates PetscInt, but other types, including
135:      * user-defined structures, could also be used. */
136:     PetscMalloc2(nrootsalloc,&rootdata,nleavesalloc,&leafdata);
137:     /* Set rootdata buffer to be broadcast */
138:     for (i=0; i<nrootsalloc; i++) rootdata[i] = -1;
139:     for (i=0; i<nroots; i++) rootdata[i*stride] = 100*(rank+1) + i;
140:     /* Initialize local buffer, these values are never used. */
141:     for (i=0; i<nleavesalloc; i++) leafdata[i] = -1;
142:     /* Broadcast entries from rootdata to leafdata. Computation or other communication can be performed between the begin and end calls. */
143:     PetscSFBcastBegin(sf,MPIU_INT,rootdata,leafdata);
144:     PetscSFBcastEnd(sf,MPIU_INT,rootdata,leafdata);
145:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Bcast Rootdata\n");
146:     PetscIntView(nrootsalloc,rootdata,PETSC_VIEWER_STDOUT_WORLD);
147:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Bcast Leafdata\n");
148:     PetscIntView(nleavesalloc,leafdata,PETSC_VIEWER_STDOUT_WORLD);
149:     PetscFree2(rootdata,leafdata);
150:   }

152:   if (test_reduce) {            /* Reduce leafdata into rootdata */
153:     PetscInt *rootdata,*leafdata;
154:     PetscMalloc2(nrootsalloc,&rootdata,nleavesalloc,&leafdata);
155:     /* Initialize rootdata buffer in which the result of the reduction will appear. */
156:     for (i=0; i<nrootsalloc; i++) rootdata[i] = -1;
157:     for (i=0; i<nroots; i++) rootdata[i*stride] = 100*(rank+1) + i;
158:     /* Set leaf values to reduce. */
159:     for (i=0; i<nleavesalloc; i++) leafdata[i] = -1;
160:     for (i=0; i<nleaves; i++) leafdata[i*stride] = 1000*(rank+1) + 10*i;
161:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Pre-Reduce Rootdata\n");
162:     PetscIntView(nrootsalloc,rootdata,PETSC_VIEWER_STDOUT_WORLD);
163:     /* Perform reduction. Computation or other communication can be performed between the begin and end calls.
164:      * This example sums the values, but other MPI_Ops can be used (e.g MPI_MAX, MPI_PROD). */
165:     PetscSFReduceBegin(sf,MPIU_INT,leafdata,rootdata,mop);
166:     PetscSFReduceEnd(sf,MPIU_INT,leafdata,rootdata,mop);
167:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Reduce Leafdata\n");
168:     PetscIntView(nleavesalloc,leafdata,PETSC_VIEWER_STDOUT_WORLD);
169:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Reduce Rootdata\n");
170:     PetscIntView(nrootsalloc,rootdata,PETSC_VIEWER_STDOUT_WORLD);
171:     PetscFree2(rootdata,leafdata);
172:   }

174:   if (test_degree) {
175:     const PetscInt *degree;
176:     PetscSFComputeDegreeBegin(sf,&degree);
177:     PetscSFComputeDegreeEnd(sf,&degree);
178:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Root degrees\n");
179:     PetscIntView(nrootsalloc,degree,PETSC_VIEWER_STDOUT_WORLD);
180:   }

182:   if (test_fetchandop) {
183:     /* Cannot use text compare here because token ordering is not deterministic */
184:     PetscInt *leafdata,*leafupdate,*rootdata;
185:     PetscMalloc3(nleavesalloc,&leafdata,nleavesalloc,&leafupdate,nrootsalloc,&rootdata);
186:     for (i=0; i<nleavesalloc; i++) leafdata[i] = -1;
187:     for (i=0; i<nleaves; i++) leafdata[i*stride] = 1;
188:     for (i=0; i<nrootsalloc; i++) rootdata[i] = -1;
189:     for (i=0; i<nroots; i++) rootdata[i*stride] = 0;
190:     PetscSFFetchAndOpBegin(sf,MPIU_INT,rootdata,leafdata,leafupdate,mop);
191:     PetscSFFetchAndOpEnd(sf,MPIU_INT,rootdata,leafdata,leafupdate,mop);
192:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Rootdata (sum of 1 from each leaf)\n");
193:     PetscIntView(nrootsalloc,rootdata,PETSC_VIEWER_STDOUT_WORLD);
194:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Leafupdate (value at roots prior to my atomic update)\n");
195:     PetscIntView(nleavesalloc,leafupdate,PETSC_VIEWER_STDOUT_WORLD);
196:     PetscFree3(leafdata,leafupdate,rootdata);
197:   }

199:   if (test_gather) {
200:     const PetscInt *degree;
201:     PetscInt       inedges,*indata,*outdata;
202:     PetscSFComputeDegreeBegin(sf,&degree);
203:     PetscSFComputeDegreeEnd(sf,&degree);
204:     for (i=0,inedges=0; i<nrootsalloc; i++) inedges += degree[i];
205:     PetscMalloc2(inedges,&indata,nleavesalloc,&outdata);
206:     for (i=0; i<nleavesalloc; i++) outdata[i] = -1;
207:     for (i=0; i<nleaves; i++) outdata[i*stride] = 1000*(rank+1) + i;
208:     PetscSFGatherBegin(sf,MPIU_INT,outdata,indata);
209:     PetscSFGatherEnd(sf,MPIU_INT,outdata,indata);
210:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Gathered data at multi-roots from leaves\n");
211:     PetscIntView(inedges,indata,PETSC_VIEWER_STDOUT_WORLD);
212:     PetscFree2(indata,outdata);
213:   }

215:   if (test_scatter) {
216:     const PetscInt *degree;
217:     PetscInt       j,count,inedges,*indata,*outdata;
218:     PetscSFComputeDegreeBegin(sf,&degree);
219:     PetscSFComputeDegreeEnd(sf,&degree);
220:     for (i=0,inedges=0; i<nrootsalloc; i++) inedges += degree[i];
221:     PetscMalloc2(inedges,&indata,nleavesalloc,&outdata);
222:     for (i=0; i<nleavesalloc; i++) outdata[i] = -1;
223:     for (i=0,count=0; i<nrootsalloc; i++) {
224:       for (j=0; j<degree[i]; j++) indata[count++] = 1000*(rank+1) + 100*i + j;
225:     }
226:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Data at multi-roots, to scatter to leaves\n");
227:     PetscIntView(inedges,indata,PETSC_VIEWER_STDOUT_WORLD);

229:     PetscSFScatterBegin(sf,MPIU_INT,indata,outdata);
230:     PetscSFScatterEnd(sf,MPIU_INT,indata,outdata);
231:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Scattered data at leaves\n");
232:     PetscIntView(nleavesalloc,outdata,PETSC_VIEWER_STDOUT_WORLD);
233:     PetscFree2(indata,outdata);
234:   }

236:   if (test_embed) {
237:     const PetscInt nroots = 1 + (PetscInt) !rank,selected[] = {1*stride,2*stride};
238:     PetscSF        esf;
239:     PetscSFCreateEmbeddedSF(sf,nroots,selected,&esf);
240:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Embedded PetscSF\n");
241:     PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
242:     PetscSFView(esf,PETSC_VIEWER_STDOUT_WORLD);
243:     PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
244:     PetscSFDestroy(&esf);
245:   }

247:   if (test_invert) {
248:     PetscSF msf,imsf;
249:     PetscSFGetMultiSF(sf,&msf);
250:     PetscSFCreateInverseSF(msf,&imsf);
251:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Multi-SF\n");
252:     PetscSFView(msf,PETSC_VIEWER_STDOUT_WORLD);
253:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Inverse of Multi-SF\n");
254:     PetscSFView(imsf,PETSC_VIEWER_STDOUT_WORLD);
255:     PetscSFDestroy(&imsf);
256:   }

258:   /* Clean storage for star forest. */
259:   PetscSFDestroy(&sf);
260:   PetscFinalize();
261:   return ierr;
262: }