Actual source code: ex1.c

petsc-3.4.5 2014-06-29
  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;

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

 31:   PetscOptionsBegin(PETSC_COMM_WORLD,"","PetscSF Test Options","none");
 32:   test_bcast      = PETSC_FALSE;
 33:   PetscOptionsBool("-test_bcast","Test broadcast","",test_bcast,&test_bcast,NULL);
 34:   test_reduce     = PETSC_FALSE;
 35:   PetscOptionsBool("-test_reduce","Test reduction","",test_reduce,&test_reduce,NULL);
 36:   test_degree     = PETSC_FALSE;
 37:   PetscOptionsBool("-test_degree","Test computation of vertex degree","",test_degree,&test_degree,NULL);
 38:   test_fetchandop = PETSC_FALSE;
 39:   PetscOptionsBool("-test_fetchandop","Test atomic Fetch-And-Op","",test_fetchandop,&test_fetchandop,NULL);
 40:   test_gather     = PETSC_FALSE;
 41:   PetscOptionsBool("-test_gather","Test point gather","",test_gather,&test_gather,NULL);
 42:   test_scatter    = PETSC_FALSE;
 43:   PetscOptionsBool("-test_scatter","Test point scatter","",test_scatter,&test_scatter,NULL);
 44:   test_embed      = PETSC_FALSE;
 45:   PetscOptionsBool("-test_embed","Test point embed","",test_embed,&test_embed,NULL);
 46:   test_invert     = PETSC_FALSE;
 47:   PetscOptionsBool("-test_invert","Test point invert","",test_invert,&test_invert,NULL);
 48:   PetscOptionsEnd();

 50:   nroots  = 2 + (PetscInt)(rank == 0);
 51:   nleaves = 2 + (PetscInt)(rank > 0);
 52:   PetscMalloc(nleaves*sizeof(*remote),&remote);
 53:   /* Left periodic neighbor */
 54:   remote[0].rank  = (rank+size-1)%size;
 55:   remote[0].index = 1;
 56:   /* Right periodic neighbor */
 57:   remote[1].rank  = (rank+1)%size;
 58:   remote[1].index = 0;
 59:   if (rank > 0) {               /* All processes reference rank 0, index 1 */
 60:     remote[2].rank  = 0;
 61:     remote[2].index = 2;
 62:   }

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

 69:   /* View graph, mostly useful for debugging purposes. */
 70:   PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
 71:   PetscSFView(sf,PETSC_VIEWER_STDOUT_WORLD);
 72:   PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);


 75:   if (test_bcast) {             /* broadcast rootdata into leafdata */
 76:     PetscInt *rootdata,*leafdata;
 77:     /* Allocate space for send and recieve buffers. This example communicates PetscInt, but other types, including
 78:      * user-defined structures, could also be used. */
 79:     PetscMalloc2(nroots,PetscInt,&rootdata,nleaves,PetscInt,&leafdata);
 80:     /* Set rootdata buffer to be broadcast */
 81:     for (i=0; i<nroots; i++) rootdata[i] = 100*(rank+1) + i;
 82:     /* Initialize local buffer, these values are never used. */
 83:     for (i=0; i<nleaves; i++) leafdata[i] = -1;
 84:     /* Broadcast entries from rootdata to leafdata. Computation or other communication can be performed between the begin and end calls. */
 85:     PetscSFBcastBegin(sf,MPIU_INT,rootdata,leafdata);
 86:     PetscSFBcastEnd(sf,MPIU_INT,rootdata,leafdata);
 87:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Bcast Rootdata\n");
 88:     PetscIntView(nroots,rootdata,PETSC_VIEWER_STDOUT_WORLD);
 89:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Bcast Leafdata\n");
 90:     PetscIntView(nleaves,leafdata,PETSC_VIEWER_STDOUT_WORLD);
 91:     PetscFree2(rootdata,leafdata);
 92:   }

 94:   if (test_reduce) {            /* Reduce leafdata into rootdata */
 95:     PetscInt *rootdata,*leafdata;
 96:     PetscMalloc2(nroots,PetscInt,&rootdata,nleaves,PetscInt,&leafdata);
 97:     /* Initialize rootdata buffer in which the result of the reduction will appear. */
 98:     for (i=0; i<nroots; i++) rootdata[i] = 100*(rank+1) + i;
 99:     /* Set leaf values to reduce. */
100:     for (i=0; i<nleaves; i++) leafdata[i] = 1000*(rank+1) + 10*i;
101:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Pre-Reduce Rootdata\n");
102:     PetscIntView(nroots,rootdata,PETSC_VIEWER_STDOUT_WORLD);
103:     /* Perform reduction. Computation or other communication can be performed between the begin and end calls.
104:      * This example sums the values, but other MPI_Ops can be used (e.g MPI_MAX, MPI_PROD). */
105:     PetscSFReduceBegin(sf,MPIU_INT,leafdata,rootdata,MPIU_SUM);
106:     PetscSFReduceEnd(sf,MPIU_INT,leafdata,rootdata,MPIU_SUM);
107:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Reduce Leafdata\n");
108:     PetscIntView(nleaves,leafdata,PETSC_VIEWER_STDOUT_WORLD);
109:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Reduce Rootdata\n");
110:     PetscIntView(nroots,rootdata,PETSC_VIEWER_STDOUT_WORLD);
111:     PetscFree2(rootdata,leafdata);
112:   }

114:   if (test_degree) {
115:     const PetscInt *degree;
116:     PetscSFComputeDegreeBegin(sf,&degree);
117:     PetscSFComputeDegreeEnd(sf,&degree);
118:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Root degrees\n");
119:     PetscIntView(nroots,degree,PETSC_VIEWER_STDOUT_WORLD);
120:   }

122:   if (test_fetchandop) {
123:     /* Cannot use text compare here because token ordering is not deterministic */
124:     PetscInt *leafdata,*leafupdate,*rootdata;
125:     PetscMalloc3(nleaves,PetscInt,&leafdata,nleaves,PetscInt,&leafupdate,nroots,PetscInt,&rootdata);
126:     for (i=0; i<nleaves; i++) leafdata[i] = 1;
127:     for (i=0; i<nroots; i++) rootdata[i] = 0;
128:     PetscSFFetchAndOpBegin(sf,MPIU_INT,rootdata,leafdata,leafupdate,MPIU_SUM);
129:     PetscSFFetchAndOpEnd(sf,MPIU_INT,rootdata,leafdata,leafupdate,MPIU_SUM);
130:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Rootdata (sum of 1 from each leaf)\n");
131:     PetscIntView(nroots,rootdata,PETSC_VIEWER_STDOUT_WORLD);
132:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Leafupdate (value at roots prior to my atomic update)\n");
133:     PetscIntView(nleaves,leafupdate,PETSC_VIEWER_STDOUT_WORLD);
134:     PetscFree3(leafdata,leafupdate,rootdata);
135:   }

137:   if (test_gather) {
138:     const PetscInt *degree;
139:     PetscInt       inedges,*indata,*outdata;
140:     PetscSFComputeDegreeBegin(sf,&degree);
141:     PetscSFComputeDegreeEnd(sf,&degree);
142:     for (i=0,inedges=0; i<nroots; i++) inedges += degree[i];
143:     PetscMalloc2(inedges,PetscInt,&indata,nleaves,PetscInt,&outdata);
144:     for (i=0; i<nleaves; i++) outdata[i] = 1000*(rank+1) + i;
145:     PetscSFGatherBegin(sf,MPIU_INT,outdata,indata);
146:     PetscSFGatherEnd(sf,MPIU_INT,outdata,indata);
147:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Gathered data at multi-roots from leaves\n");
148:     PetscIntView(inedges,indata,PETSC_VIEWER_STDOUT_WORLD);
149:     PetscFree2(indata,outdata);
150:   }

152:   if (test_scatter) {
153:     const PetscInt *degree;
154:     PetscInt       j,count,inedges,*indata,*outdata;
155:     PetscSFComputeDegreeBegin(sf,&degree);
156:     PetscSFComputeDegreeEnd(sf,&degree);
157:     for (i=0,inedges=0; i<nroots; i++) inedges += degree[i];
158:     PetscMalloc2(inedges,PetscInt,&indata,nleaves,PetscInt,&outdata);
159:     for (i=0,count=0; i<nroots; i++) {
160:       for (j=0; j<degree[i]; j++) indata[count++] = 1000*(rank+1) + 100*i + j;
161:     }
162:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Data at multi-roots, to scatter to leaves\n");
163:     PetscIntView(inedges,indata,PETSC_VIEWER_STDOUT_WORLD);

165:     PetscSFScatterBegin(sf,MPIU_INT,indata,outdata);
166:     PetscSFScatterEnd(sf,MPIU_INT,indata,outdata);
167:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Scattered data at leaves\n");
168:     PetscIntView(nleaves,outdata,PETSC_VIEWER_STDOUT_WORLD);
169:     PetscFree2(indata,outdata);
170:   }

172:   if (test_embed) {
173:     const PetscInt nroots = 1 + (PetscInt) !rank,selected[] = {1,2};
174:     PetscSF        esf;
175:     PetscSFCreateEmbeddedSF(sf,nroots,selected,&esf);
176:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Embedded PetscSF\n");
177:     PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
178:     PetscSFView(esf,PETSC_VIEWER_STDOUT_WORLD);
179:     PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
180:     PetscSFDestroy(&esf);
181:   }

183:   if (test_invert) {
184:     PetscSF msf,imsf;
185:     PetscSFGetMultiSF(sf,&msf);
186:     PetscSFCreateInverseSF(msf,&imsf);
187:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Multi-SF\n");
188:     PetscSFView(msf,PETSC_VIEWER_STDOUT_WORLD);
189:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"## Inverse of Multi-SF\n");
190:     PetscSFView(imsf,PETSC_VIEWER_STDOUT_WORLD);
191:     PetscSFDestroy(&imsf);
192:   }

194:   /* Clean storage for star forest. */
195:   PetscSFDestroy(&sf);
196:   PetscFinalize();
197:   return 0;
198: }