Actual source code: ex1.c

petsc-3.6.1 2015-08-06
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>

 18: int main(int argc,char **argv)
 19: {
 21:   PetscInt       i,nroots,nleaves;
 22:   PetscSFNode    *remote;
 23:   PetscMPIInt    rank,size;
 24:   PetscSF        sf;
 25:   PetscBool      test_bcast,test_reduce,test_degree,test_fetchandop,test_gather,test_scatter,test_embed,test_invert;
 26:   MPI_Op         mop;
 27:   char           opstring[256];
 28:   PetscBool      strflg;

 30:   PetscInitialize(&argc,&argv,(char*)0,help);
 31:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 32:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

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

 96:   nroots  = 2 + (PetscInt)(rank == 0);
 97:   nleaves = 2 + (PetscInt)(rank > 0);
 98:   PetscMalloc1(nleaves,&remote);
 99:   /* Left periodic neighbor */
100:   remote[0].rank  = (rank+size-1)%size;
101:   remote[0].index = 1;
102:   /* Right periodic neighbor */
103:   remote[1].rank  = (rank+1)%size;
104:   remote[1].index = 0;
105:   if (rank > 0) {               /* All processes reference rank 0, index 1 */
106:     remote[2].rank  = 0;
107:     remote[2].index = 2;
108:   }

110:   /* Create a star forest for communication. In this example, the leaf space is dense, so we pass NULL. */
111:   PetscSFCreate(PETSC_COMM_WORLD,&sf);
112:   PetscSFSetFromOptions(sf);
113:   PetscSFSetGraph(sf,nroots,nleaves,NULL,PETSC_COPY_VALUES,remote,PETSC_OWN_POINTER);

115:   /* View graph, mostly useful for debugging purposes. */
116:   PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
117:   PetscSFView(sf,PETSC_VIEWER_STDOUT_WORLD);
118:   PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);


121:   if (test_bcast) {             /* broadcast rootdata into leafdata */
122:     PetscInt *rootdata,*leafdata;
123:     /* Allocate space for send and recieve buffers. This example communicates PetscInt, but other types, including
124:      * user-defined structures, could also be used. */
125:     PetscMalloc2(nroots,&rootdata,nleaves,&leafdata);
126:     /* Set rootdata buffer to be broadcast */
127:     for (i=0; i<nroots; i++) rootdata[i] = 100*(rank+1) + i;
128:     /* Initialize local buffer, these values are never used. */
129:     for (i=0; i<nleaves; i++) leafdata[i] = -1;
130:     /* Broadcast entries from rootdata to leafdata. Computation or other communication can be performed between the begin and end calls. */
131:     PetscSFBcastBegin(sf,MPIU_INT,rootdata,leafdata);
132:     PetscSFBcastEnd(sf,MPIU_INT,rootdata,leafdata);
133:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Bcast Rootdata\n");
134:     PetscIntView(nroots,rootdata,PETSC_VIEWER_STDOUT_WORLD);
135:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Bcast Leafdata\n");
136:     PetscIntView(nleaves,leafdata,PETSC_VIEWER_STDOUT_WORLD);
137:     PetscFree2(rootdata,leafdata);
138:   }

140:   if (test_reduce) {            /* Reduce leafdata into rootdata */
141:     PetscInt *rootdata,*leafdata;
142:     PetscMalloc2(nroots,&rootdata,nleaves,&leafdata);
143:     /* Initialize rootdata buffer in which the result of the reduction will appear. */
144:     for (i=0; i<nroots; i++) rootdata[i] = 100*(rank+1) + i;
145:     /* Set leaf values to reduce. */
146:     for (i=0; i<nleaves; i++) leafdata[i] = 1000*(rank+1) + 10*i;
147:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Pre-Reduce Rootdata\n");
148:     PetscIntView(nroots,rootdata,PETSC_VIEWER_STDOUT_WORLD);
149:     /* Perform reduction. Computation or other communication can be performed between the begin and end calls.
150:      * This example sums the values, but other MPI_Ops can be used (e.g MPI_MAX, MPI_PROD). */
151:     PetscSFReduceBegin(sf,MPIU_INT,leafdata,rootdata,mop);
152:     PetscSFReduceEnd(sf,MPIU_INT,leafdata,rootdata,mop);
153:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Reduce Leafdata\n");
154:     PetscIntView(nleaves,leafdata,PETSC_VIEWER_STDOUT_WORLD);
155:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Reduce Rootdata\n");
156:     PetscIntView(nroots,rootdata,PETSC_VIEWER_STDOUT_WORLD);
157:     PetscFree2(rootdata,leafdata);
158:   }

160:   if (test_degree) {
161:     const PetscInt *degree;
162:     PetscSFComputeDegreeBegin(sf,&degree);
163:     PetscSFComputeDegreeEnd(sf,&degree);
164:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Root degrees\n");
165:     PetscIntView(nroots,degree,PETSC_VIEWER_STDOUT_WORLD);
166:   }

168:   if (test_fetchandop) {
169:     /* Cannot use text compare here because token ordering is not deterministic */
170:     PetscInt *leafdata,*leafupdate,*rootdata;
171:     PetscMalloc3(nleaves,&leafdata,nleaves,&leafupdate,nroots,&rootdata);
172:     for (i=0; i<nleaves; i++) leafdata[i] = 1;
173:     for (i=0; i<nroots; i++) rootdata[i] = 0;
174:     PetscSFFetchAndOpBegin(sf,MPIU_INT,rootdata,leafdata,leafupdate,mop);
175:     PetscSFFetchAndOpEnd(sf,MPIU_INT,rootdata,leafdata,leafupdate,mop);
176:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Rootdata (sum of 1 from each leaf)\n");
177:     PetscIntView(nroots,rootdata,PETSC_VIEWER_STDOUT_WORLD);
178:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Leafupdate (value at roots prior to my atomic update)\n");
179:     PetscIntView(nleaves,leafupdate,PETSC_VIEWER_STDOUT_WORLD);
180:     PetscFree3(leafdata,leafupdate,rootdata);
181:   }

183:   if (test_gather) {
184:     const PetscInt *degree;
185:     PetscInt       inedges,*indata,*outdata;
186:     PetscSFComputeDegreeBegin(sf,&degree);
187:     PetscSFComputeDegreeEnd(sf,&degree);
188:     for (i=0,inedges=0; i<nroots; i++) inedges += degree[i];
189:     PetscMalloc2(inedges,&indata,nleaves,&outdata);
190:     for (i=0; i<nleaves; i++) outdata[i] = 1000*(rank+1) + i;
191:     PetscSFGatherBegin(sf,MPIU_INT,outdata,indata);
192:     PetscSFGatherEnd(sf,MPIU_INT,outdata,indata);
193:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Gathered data at multi-roots from leaves\n");
194:     PetscIntView(inedges,indata,PETSC_VIEWER_STDOUT_WORLD);
195:     PetscFree2(indata,outdata);
196:   }

198:   if (test_scatter) {
199:     const PetscInt *degree;
200:     PetscInt       j,count,inedges,*indata,*outdata;
201:     PetscSFComputeDegreeBegin(sf,&degree);
202:     PetscSFComputeDegreeEnd(sf,&degree);
203:     for (i=0,inedges=0; i<nroots; i++) inedges += degree[i];
204:     PetscMalloc2(inedges,&indata,nleaves,&outdata);
205:     for (i=0,count=0; i<nroots; i++) {
206:       for (j=0; j<degree[i]; j++) indata[count++] = 1000*(rank+1) + 100*i + j;
207:     }
208:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Data at multi-roots, to scatter to leaves\n");
209:     PetscIntView(inedges,indata,PETSC_VIEWER_STDOUT_WORLD);

211:     PetscSFScatterBegin(sf,MPIU_INT,indata,outdata);
212:     PetscSFScatterEnd(sf,MPIU_INT,indata,outdata);
213:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Scattered data at leaves\n");
214:     PetscIntView(nleaves,outdata,PETSC_VIEWER_STDOUT_WORLD);
215:     PetscFree2(indata,outdata);
216:   }

218:   if (test_embed) {
219:     const PetscInt nroots = 1 + (PetscInt) !rank,selected[] = {1,2};
220:     PetscSF        esf;
221:     PetscSFCreateEmbeddedSF(sf,nroots,selected,&esf);
222:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Embedded PetscSF\n");
223:     PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
224:     PetscSFView(esf,PETSC_VIEWER_STDOUT_WORLD);
225:     PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
226:     PetscSFDestroy(&esf);
227:   }

229:   if (test_invert) {
230:     PetscSF msf,imsf;
231:     PetscSFGetMultiSF(sf,&msf);
232:     PetscSFCreateInverseSF(msf,&imsf);
233:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Multi-SF\n");
234:     PetscSFView(msf,PETSC_VIEWER_STDOUT_WORLD);
235:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Inverse of Multi-SF\n");
236:     PetscSFView(imsf,PETSC_VIEWER_STDOUT_WORLD);
237:     PetscSFDestroy(&imsf);
238:   }

240:   /* Clean storage for star forest. */
241:   PetscSFDestroy(&sf);
242:   PetscFinalize();
243:   return 0;
244: }