Actual source code: ex8.c

petsc-3.13.6 2020-09-29
Report Typos and Errors
  1: static char help[] = "Demonstrates BuildTwoSided functions.\n";

  3:  #include <petscsys.h>

  5: typedef struct {
  6:   PetscInt    rank;
  7:   PetscScalar value;
  8:   char        ok[3];
  9: } Unit;

 11: static PetscErrorCode MakeDatatype(MPI_Datatype *dtype)
 12: {
 14:   MPI_Datatype dtypes[3],tmptype;
 15:   PetscMPIInt  lengths[3];
 16:   MPI_Aint     displs[3];
 17:   Unit         dummy;

 20:   dtypes[0] = MPIU_INT;
 21:   dtypes[1] = MPIU_SCALAR;
 22:   dtypes[2] = MPI_CHAR;
 23:   lengths[0] = 1;
 24:   lengths[1] = 1;
 25:   lengths[2] = 3;
 26:   /* Curse the evil beings that made std::complex a non-POD type. */
 27:   displs[0] = (char*)&dummy.rank - (char*)&dummy;  /* offsetof(Unit,rank); */
 28:   displs[1] = (char*)&dummy.value - (char*)&dummy; /* offsetof(Unit,value); */
 29:   displs[2] = (char*)&dummy.ok - (char*)&dummy;    /* offsetof(Unit,ok); */
 30:   MPI_Type_create_struct(3,lengths,displs,dtypes,&tmptype);
 31:   MPI_Type_commit(&tmptype);
 32:   MPI_Type_create_resized(tmptype,0,sizeof(Unit),dtype);
 33:   MPI_Type_commit(dtype);
 34:   MPI_Type_free(&tmptype);
 35:   {
 36:     MPI_Aint lb,extent;
 37:     MPI_Type_get_extent(*dtype,&lb,&extent);
 38:     if (extent != sizeof(Unit)) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_LIB,"New type has extent %d != sizeof(Unit) %d",extent,(int)sizeof(Unit));
 39:   }
 40:   return(0);
 41: }

 43: struct FCtx {
 44:   PetscMPIInt rank;
 45:   PetscMPIInt nto;
 46:   PetscMPIInt *toranks;
 47:   Unit *todata;
 48:   PetscSegBuffer seg;
 49: };

 51: static PetscErrorCode FSend(MPI_Comm comm,const PetscMPIInt tag[],PetscMPIInt tonum,PetscMPIInt rank,void *todata,MPI_Request req[],void *ctx)
 52: {
 53:   struct FCtx *fctx = (struct FCtx*)ctx;

 57:   if (rank != fctx->toranks[tonum]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Rank %d does not match toranks[%d] %d",rank,tonum,fctx->toranks[tonum]);
 58:   if (fctx->rank != *(PetscMPIInt*)todata) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Todata %d does not match rank %d",*(PetscMPIInt*)todata,fctx->rank);
 59:   MPI_Isend(&fctx->todata[tonum].rank,1,MPIU_INT,rank,tag[0],comm,&req[0]);
 60:   MPI_Isend(&fctx->todata[tonum].value,1,MPIU_SCALAR,rank,tag[1],comm,&req[1]);
 61:   return(0);
 62: }

 64: static PetscErrorCode FRecv(MPI_Comm comm,const PetscMPIInt tag[],PetscMPIInt rank,void *fromdata,MPI_Request req[],void *ctx)
 65: {
 66:   struct FCtx *fctx = (struct FCtx*)ctx;
 68:   Unit           *buf;

 71:   if (*(PetscMPIInt*)fromdata != rank) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Dummy data %d from rank %d corrupt",*(PetscMPIInt*)fromdata,rank);
 72:   PetscSegBufferGet(fctx->seg,1,&buf);
 73:   MPI_Irecv(&buf->rank,1,MPIU_INT,rank,tag[0],comm,&req[0]);
 74:   MPI_Irecv(&buf->value,1,MPIU_SCALAR,rank,tag[1],comm,&req[1]);
 75:   buf->ok[0] = 'o';
 76:   buf->ok[1] = 'k';
 77:   buf->ok[2] = 0;
 78:   return(0);
 79: }

 81: int main(int argc,char **argv)
 82: {
 84:   PetscMPIInt    rank,size,*toranks,*fromranks,nto,nfrom;
 85:   PetscInt       i,n;
 86:   PetscBool      verbose,build_twosided_f;
 87:   Unit           *todata,*fromdata;
 88:   MPI_Datatype   dtype;

 90:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 91:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 92:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);

 94:   verbose = PETSC_FALSE;
 95:   PetscOptionsGetBool(NULL,NULL,"-verbose",&verbose,NULL);
 96:   build_twosided_f = PETSC_FALSE;
 97:   PetscOptionsGetBool(NULL,NULL,"-build_twosided_f",&build_twosided_f,NULL);

 99:   for (i=1,nto=0; i<size; i*=2) nto++;
100:   PetscMalloc2(nto,&todata,nto,&toranks);
101:   for (n=0,i=1; i<size; n++,i*=2) {
102:     toranks[n] = (rank+i) % size;
103:     todata[n].rank  = (rank+i) % size;
104:     todata[n].value = (PetscScalar)rank;
105:     todata[n].ok[0] = 'o';
106:     todata[n].ok[1] = 'k';
107:     todata[n].ok[2] = 0;
108:   }
109:   if (verbose) {
110:     for (i=0; i<nto; i++) {
111:       PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d] TO %d: {%D, %g, \"%s\"}\n",rank,toranks[i],todata[i].rank,(double)PetscRealPart(todata[i].value),todata[i].ok);
112:     }
113:     PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);
114:   }

116:   MakeDatatype(&dtype);

118:   if (build_twosided_f) {
119:     struct FCtx fctx;
120:     PetscMPIInt *todummy,*fromdummy;
121:     fctx.rank    = rank;
122:     fctx.nto     = nto;
123:     fctx.toranks = toranks;
124:     fctx.todata  = todata;
125:     PetscSegBufferCreate(sizeof(Unit),1,&fctx.seg);
126:     PetscMalloc1(nto,&todummy);
127:     for (i=0; i<nto; i++) todummy[i] = rank;
128:     PetscCommBuildTwoSidedF(PETSC_COMM_WORLD,1,MPI_INT,nto,toranks,todummy,&nfrom,&fromranks,&fromdummy,2,FSend,FRecv,&fctx);
129:     PetscFree(todummy);
130:     PetscFree(fromdummy);
131:     PetscSegBufferExtractAlloc(fctx.seg,&fromdata);
132:     PetscSegBufferDestroy(&fctx.seg);
133:   } else {
134:     PetscCommBuildTwoSided(PETSC_COMM_WORLD,1,dtype,nto,toranks,todata,&nfrom,&fromranks,&fromdata);
135:   }
136:   MPI_Type_free(&dtype);

138:   if (verbose) {
139:     PetscInt *iranks,*iperm;
140:     PetscMalloc2(nfrom,&iranks,nfrom,&iperm);
141:     for (i=0; i<nfrom; i++) {
142:       iranks[i] = fromranks[i];
143:       iperm[i] = i;
144:     }
145:     /* Receive ordering is non-deterministic in general, so sort to make verbose output deterministic. */
146:     PetscSortIntWithPermutation(nfrom,iranks,iperm);
147:     for (i=0; i<nfrom; i++) {
148:       PetscInt ip = iperm[i];
149:       PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d] FROM %d: {%D, %g, \"%s\"}\n",rank,fromranks[ip],fromdata[ip].rank,(double)PetscRealPart(fromdata[ip].value),fromdata[ip].ok);
150:     }
151:     PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);
152:     PetscFree2(iranks,iperm);
153:   }

155:   if (nto != nfrom) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"[%d] From ranks %d does not match To ranks %d",rank,nto,nfrom);
156:   for (i=1; i<size; i*=2) {
157:     PetscMPIInt expected_rank = (rank-i+size)%size;
158:     PetscBool flg;
159:     for (n=0; n<nfrom; n++) {
160:       if (expected_rank == fromranks[n]) goto found;
161:     }
162:     SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"[%d] Could not find expected from rank %d",rank,expected_rank);
163:     found:
164:     if (PetscRealPart(fromdata[n].value) != expected_rank) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"[%d] Got data %g from rank %d",rank,(double)PetscRealPart(fromdata[n].value),expected_rank);
165:     PetscStrcmp(fromdata[n].ok,"ok",&flg);
166:     if (!flg) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"[%d] Got string %s from rank %d",rank,fromdata[n].ok,expected_rank);
167:   }
168:   PetscFree2(todata,toranks);
169:   PetscFree(fromdata);
170:   PetscFree(fromranks);
171:   PetscFinalize();
172:   return ierr;
173: }



177: /*TEST

179:    test:
180:       nsize: 4
181:       args: -verbose -build_twosided allreduce

183:    test:
184:       suffix: f
185:       nsize: 4
186:       args: -verbose -build_twosided_f -build_twosided allreduce
187:       output_file: output/ex8_1.out

189:    test:
190:       suffix: f_ibarrier
191:       nsize: 4
192:       args: -verbose -build_twosided_f -build_twosided ibarrier
193:       output_file: output/ex8_1.out

195:    test:
196:       suffix: ibarrier
197:       nsize: 4
198:       args: -verbose -build_twosided ibarrier
199:       output_file: output/ex8_1.out

201:    test:
202:       suffix: redscatter
203:       requires: mpi_reduce_scatter_block
204:       nsize: 4
205:       args: -verbose -build_twosided redscatter
206:       output_file: output/ex8_1.out

208: TEST*/