Actual source code: ex8.c
petsc-3.7.7 2017-09-25
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;
13: static PetscErrorCode MakeDatatype(MPI_Datatype *dtype)
14: {
16: MPI_Datatype dtypes[3],tmptype;
17: PetscMPIInt lengths[3];
18: MPI_Aint displs[3];
19: Unit dummy;
22: dtypes[0] = MPIU_INT;
23: dtypes[1] = MPIU_SCALAR;
24: dtypes[2] = MPI_CHAR;
25: lengths[0] = 1;
26: lengths[1] = 1;
27: lengths[2] = 3;
28: /* Curse the evil beings that made std::complex a non-POD type. */
29: displs[0] = (char*)&dummy.rank - (char*)&dummy; /* offsetof(Unit,rank); */
30: displs[1] = (char*)&dummy.value - (char*)&dummy; /* offsetof(Unit,value); */
31: displs[2] = (char*)&dummy.ok - (char*)&dummy; /* offsetof(Unit,ok); */
32: MPI_Type_create_struct(3,lengths,displs,dtypes,&tmptype);
33: MPI_Type_commit(&tmptype);
34: MPI_Type_create_resized(tmptype,0,sizeof(Unit),dtype);
35: MPI_Type_commit(dtype);
36: MPI_Type_free(&tmptype);
37: {
38: MPI_Aint lb,extent;
39: MPI_Type_get_extent(*dtype,&lb,&extent);
40: if (extent != sizeof(Unit)) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_LIB,"New type has extent %d != sizeof(Unit) %d",extent,(int)sizeof(Unit));
41: }
42: return(0);
43: }
45: struct FCtx {
46: PetscMPIInt rank;
47: PetscMPIInt nto;
48: PetscMPIInt *toranks;
49: Unit *todata;
50: PetscSegBuffer seg;
51: };
55: static PetscErrorCode FSend(MPI_Comm comm,const PetscMPIInt tag[],PetscMPIInt tonum,PetscMPIInt rank,void *todata,MPI_Request req[],void *ctx)
56: {
57: struct FCtx *fctx = (struct FCtx*)ctx;
61: 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]);
62: if (fctx->rank != *(PetscMPIInt*)todata) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Todata %d does not match rank %d",*(PetscMPIInt*)todata,fctx->rank);
63: MPI_Isend(&fctx->todata[tonum].rank,1,MPIU_INT,rank,tag[0],comm,&req[0]);
64: MPI_Isend(&fctx->todata[tonum].value,1,MPIU_SCALAR,rank,tag[1],comm,&req[1]);
65: return(0);
66: }
70: static PetscErrorCode FRecv(MPI_Comm comm,const PetscMPIInt tag[],PetscMPIInt rank,void *fromdata,MPI_Request req[],void *ctx)
71: {
72: struct FCtx *fctx = (struct FCtx*)ctx;
74: Unit *buf;
77: if (*(PetscMPIInt*)fromdata != rank) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Dummy data %d from rank %d corrupt",*(PetscMPIInt*)fromdata,rank);
78: PetscSegBufferGet(fctx->seg,1,&buf);
79: MPI_Irecv(&buf->rank,1,MPIU_INT,rank,tag[0],comm,&req[0]);
80: MPI_Irecv(&buf->value,1,MPIU_SCALAR,rank,tag[1],comm,&req[1]);
81: buf->ok[0] = 'o';
82: buf->ok[1] = 'k';
83: buf->ok[2] = 0;
84: return(0);
85: }
89: int main(int argc,char **argv)
90: {
92: PetscMPIInt rank,size,*toranks,*fromranks,nto,nfrom;
93: PetscInt i,n;
94: PetscBool verbose,build_twosided_f;
95: Unit *todata,*fromdata;
96: MPI_Datatype dtype;
98: PetscInitialize(&argc,&argv,(char*)0,help);
99: MPI_Comm_size(PETSC_COMM_WORLD,&size);
100: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
102: verbose = PETSC_FALSE;
103: PetscOptionsGetBool(NULL,NULL,"-verbose",&verbose,NULL);
104: build_twosided_f = PETSC_FALSE;
105: PetscOptionsGetBool(NULL,NULL,"-build_twosided_f",&build_twosided_f,NULL);
107: for (i=1,nto=0; i<size; i*=2) nto++;
108: PetscMalloc2(nto,&todata,nto,&toranks);
109: for (n=0,i=1; i<size; n++,i*=2) {
110: toranks[n] = (rank+i) % size;
111: todata[n].rank = (rank+i) % size;
112: todata[n].value = (PetscScalar)rank;
113: todata[n].ok[0] = 'o';
114: todata[n].ok[1] = 'k';
115: todata[n].ok[2] = 0;
116: }
117: if (verbose) {
118: for (i=0; i<nto; i++) {
119: 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);
120: }
121: PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);
122: }
124: MakeDatatype(&dtype);
126: if (build_twosided_f) {
127: struct FCtx fctx;
128: PetscMPIInt *todummy,*fromdummy;
129: fctx.rank = rank;
130: fctx.nto = nto;
131: fctx.toranks = toranks;
132: fctx.todata = todata;
133: PetscSegBufferCreate(sizeof(Unit),1,&fctx.seg);
134: PetscMalloc1(nto,&todummy);
135: for (i=0; i<nto; i++) todummy[i] = rank;
136: PetscCommBuildTwoSidedF(PETSC_COMM_WORLD,1,MPI_INT,nto,toranks,todummy,&nfrom,&fromranks,&fromdummy,2,FSend,FRecv,&fctx);
137: PetscFree(todummy);
138: PetscFree(fromdummy);
139: PetscSegBufferExtractAlloc(fctx.seg,&fromdata);
140: PetscSegBufferDestroy(&fctx.seg);
141: } else {
142: PetscCommBuildTwoSided(PETSC_COMM_WORLD,1,dtype,nto,toranks,todata,&nfrom,&fromranks,&fromdata);
143: }
144: MPI_Type_free(&dtype);
146: if (verbose) {
147: PetscInt *iranks,*iperm;
148: PetscMalloc2(nfrom,&iranks,nfrom,&iperm);
149: for (i=0; i<nfrom; i++) {
150: iranks[i] = fromranks[i];
151: iperm[i] = i;
152: }
153: /* Receive ordering is non-deterministic in general, so sort to make verbose output deterministic. */
154: PetscSortIntWithPermutation(nfrom,iranks,iperm);
155: for (i=0; i<nfrom; i++) {
156: PetscInt ip = iperm[i];
157: 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);
158: }
159: PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);
160: PetscFree2(iranks,iperm);
161: }
163: if (nto != nfrom) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"[%d] From ranks %d does not match To ranks %d",rank,nto,nfrom);
164: for (i=1; i<size; i*=2) {
165: PetscMPIInt expected_rank = (rank-i+size)%size;
166: PetscBool flg;
167: for (n=0; n<nfrom; n++) {
168: if (expected_rank == fromranks[n]) goto found;
169: }
170: SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"[%d] Could not find expected from rank %d",rank,expected_rank);
171: found:
172: 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);
173: PetscStrcmp(fromdata[n].ok,"ok",&flg);
174: if (!flg) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"[%d] Got string %s from rank %d",rank,fromdata[n].ok,expected_rank);
175: }
176: PetscFree2(todata,toranks);
177: PetscFree(fromdata);
178: PetscFree(fromranks);
179: PetscFinalize();
180: return 0;
181: }