xref: /petsc/src/sys/tests/ex8.c (revision 40badf4fbc550ac1f60bd080eaff6de6d55b946d)
1 static char help[] = "Demonstrates BuildTwoSided functions.\n";
2 
3 #include <petscsys.h>
4 
5 typedef struct {
6   PetscInt    rank;
7   PetscScalar value;
8   char        ok[3];
9 } Unit;
10 
11 static PetscErrorCode MakeDatatype(MPI_Datatype *dtype)
12 {
13   MPI_Datatype dtypes[3],tmptype;
14   PetscMPIInt  lengths[3];
15   MPI_Aint     displs[3];
16   Unit         dummy;
17 
18   PetscFunctionBegin;
19   dtypes[0] = MPIU_INT;
20   dtypes[1] = MPIU_SCALAR;
21   dtypes[2] = MPI_CHAR;
22   lengths[0] = 1;
23   lengths[1] = 1;
24   lengths[2] = 3;
25   /* Curse the evil beings that made std::complex a non-POD type. */
26   displs[0] = (char*)&dummy.rank - (char*)&dummy;  /* offsetof(Unit,rank); */
27   displs[1] = (char*)&dummy.value - (char*)&dummy; /* offsetof(Unit,value); */
28   displs[2] = (char*)&dummy.ok - (char*)&dummy;    /* offsetof(Unit,ok); */
29   CHKERRMPI(MPI_Type_create_struct(3,lengths,displs,dtypes,&tmptype));
30   CHKERRMPI(MPI_Type_commit(&tmptype));
31   CHKERRMPI(MPI_Type_create_resized(tmptype,0,sizeof(Unit),dtype));
32   CHKERRMPI(MPI_Type_commit(dtype));
33   CHKERRMPI(MPI_Type_free(&tmptype));
34   {
35     MPI_Aint lb,extent;
36     CHKERRMPI(MPI_Type_get_extent(*dtype,&lb,&extent));
37     PetscCheckFalse(extent != sizeof(Unit),PETSC_COMM_WORLD,PETSC_ERR_LIB,"New type has extent %d != sizeof(Unit) %d",(int)extent,(int)sizeof(Unit));
38   }
39   PetscFunctionReturn(0);
40 }
41 
42 struct FCtx {
43   PetscMPIInt rank;
44   PetscMPIInt nto;
45   PetscMPIInt *toranks;
46   Unit *todata;
47   PetscSegBuffer seg;
48 };
49 
50 static PetscErrorCode FSend(MPI_Comm comm,const PetscMPIInt tag[],PetscMPIInt tonum,PetscMPIInt rank,void *todata,MPI_Request req[],void *ctx)
51 {
52   struct FCtx *fctx = (struct FCtx*)ctx;
53 
54   PetscFunctionBegin;
55   PetscCheckFalse(rank != fctx->toranks[tonum],PETSC_COMM_SELF,PETSC_ERR_PLIB,"Rank %d does not match toranks[%d] %d",rank,tonum,fctx->toranks[tonum]);
56   PetscCheckFalse(fctx->rank != *(PetscMPIInt*)todata,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Todata %d does not match rank %d",*(PetscMPIInt*)todata,fctx->rank);
57   CHKERRMPI(MPI_Isend(&fctx->todata[tonum].rank,1,MPIU_INT,rank,tag[0],comm,&req[0]));
58   CHKERRMPI(MPI_Isend(&fctx->todata[tonum].value,1,MPIU_SCALAR,rank,tag[1],comm,&req[1]));
59   PetscFunctionReturn(0);
60 }
61 
62 static PetscErrorCode FRecv(MPI_Comm comm,const PetscMPIInt tag[],PetscMPIInt rank,void *fromdata,MPI_Request req[],void *ctx)
63 {
64   struct FCtx *fctx = (struct FCtx*)ctx;
65   Unit           *buf;
66 
67   PetscFunctionBegin;
68   PetscCheckFalse(*(PetscMPIInt*)fromdata != rank,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Dummy data %d from rank %d corrupt",*(PetscMPIInt*)fromdata,rank);
69   CHKERRQ(PetscSegBufferGet(fctx->seg,1,&buf));
70   CHKERRMPI(MPI_Irecv(&buf->rank,1,MPIU_INT,rank,tag[0],comm,&req[0]));
71   CHKERRMPI(MPI_Irecv(&buf->value,1,MPIU_SCALAR,rank,tag[1],comm,&req[1]));
72   buf->ok[0] = 'o';
73   buf->ok[1] = 'k';
74   buf->ok[2] = 0;
75   PetscFunctionReturn(0);
76 }
77 
78 int main(int argc,char **argv)
79 {
80   PetscErrorCode ierr;
81   PetscMPIInt    rank,size,*toranks,*fromranks,nto,nfrom;
82   PetscInt       i,n;
83   PetscBool      verbose,build_twosided_f;
84   Unit           *todata,*fromdata;
85   MPI_Datatype   dtype;
86 
87   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
88   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
89   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
90 
91   verbose = PETSC_FALSE;
92   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-verbose",&verbose,NULL));
93   build_twosided_f = PETSC_FALSE;
94   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-build_twosided_f",&build_twosided_f,NULL));
95 
96   for (i=1,nto=0; i<size; i*=2) nto++;
97   CHKERRQ(PetscMalloc2(nto,&todata,nto,&toranks));
98   for (n=0,i=1; i<size; n++,i*=2) {
99     toranks[n] = (rank+i) % size;
100     todata[n].rank  = (rank+i) % size;
101     todata[n].value = (PetscScalar)rank;
102     todata[n].ok[0] = 'o';
103     todata[n].ok[1] = 'k';
104     todata[n].ok[2] = 0;
105   }
106   if (verbose) {
107     for (i=0; i<nto; i++) {
108       CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d] TO %d: {%" PetscInt_FMT ", %g, \"%s\"}\n",rank,toranks[i],todata[i].rank,(double)PetscRealPart(todata[i].value),todata[i].ok));
109     }
110     CHKERRQ(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT));
111   }
112 
113   CHKERRQ(MakeDatatype(&dtype));
114 
115   if (build_twosided_f) {
116     struct FCtx fctx;
117     PetscMPIInt *todummy,*fromdummy;
118     fctx.rank    = rank;
119     fctx.nto     = nto;
120     fctx.toranks = toranks;
121     fctx.todata  = todata;
122     CHKERRQ(PetscSegBufferCreate(sizeof(Unit),1,&fctx.seg));
123     CHKERRQ(PetscMalloc1(nto,&todummy));
124     for (i=0; i<nto; i++) todummy[i] = rank;
125     CHKERRQ(PetscCommBuildTwoSidedF(PETSC_COMM_WORLD,1,MPI_INT,nto,toranks,todummy,&nfrom,&fromranks,&fromdummy,2,FSend,FRecv,&fctx));
126     CHKERRQ(PetscFree(todummy));
127     CHKERRQ(PetscFree(fromdummy));
128     CHKERRQ(PetscSegBufferExtractAlloc(fctx.seg,&fromdata));
129     CHKERRQ(PetscSegBufferDestroy(&fctx.seg));
130   } else {
131     CHKERRQ(PetscCommBuildTwoSided(PETSC_COMM_WORLD,1,dtype,nto,toranks,todata,&nfrom,&fromranks,&fromdata));
132   }
133   CHKERRMPI(MPI_Type_free(&dtype));
134 
135   if (verbose) {
136     PetscInt *iranks,*iperm;
137     CHKERRQ(PetscMalloc2(nfrom,&iranks,nfrom,&iperm));
138     for (i=0; i<nfrom; i++) {
139       iranks[i] = fromranks[i];
140       iperm[i] = i;
141     }
142     /* Receive ordering is non-deterministic in general, so sort to make verbose output deterministic. */
143     CHKERRQ(PetscSortIntWithPermutation(nfrom,iranks,iperm));
144     for (i=0; i<nfrom; i++) {
145       PetscInt ip = iperm[i];
146       CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d] FROM %d: {%" PetscInt_FMT ", %g, \"%s\"}\n",rank,fromranks[ip],fromdata[ip].rank,(double)PetscRealPart(fromdata[ip].value),fromdata[ip].ok));
147     }
148     CHKERRQ(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT));
149     CHKERRQ(PetscFree2(iranks,iperm));
150   }
151 
152   PetscCheckFalse(nto != nfrom,PETSC_COMM_SELF,PETSC_ERR_PLIB,"[%d] From ranks %d does not match To ranks %d",rank,nto,nfrom);
153   for (i=1; i<size; i*=2) {
154     PetscMPIInt expected_rank = (rank-i+size)%size;
155     PetscBool flg;
156     for (n=0; n<nfrom; n++) {
157       if (expected_rank == fromranks[n]) goto found;
158     }
159     SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"[%d] Could not find expected from rank %d",rank,expected_rank);
160     found:
161     PetscCheckFalse(PetscRealPart(fromdata[n].value) != expected_rank,PETSC_COMM_SELF,PETSC_ERR_PLIB,"[%d] Got data %g from rank %d",rank,(double)PetscRealPart(fromdata[n].value),expected_rank);
162     CHKERRQ(PetscStrcmp(fromdata[n].ok,"ok",&flg));
163     PetscCheckFalse(!flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"[%d] Got string %s from rank %d",rank,fromdata[n].ok,expected_rank);
164   }
165   CHKERRQ(PetscFree2(todata,toranks));
166   CHKERRQ(PetscFree(fromdata));
167   CHKERRQ(PetscFree(fromranks));
168   ierr = PetscFinalize();
169   return ierr;
170 }
171 
172 /*TEST
173 
174    test:
175       nsize: 4
176       args: -verbose -build_twosided allreduce
177 
178    test:
179       suffix: f
180       nsize: 4
181       args: -verbose -build_twosided_f -build_twosided allreduce
182       output_file: output/ex8_1.out
183 
184    test:
185       suffix: f_ibarrier
186       nsize: 4
187       args: -verbose -build_twosided_f -build_twosided ibarrier
188       output_file: output/ex8_1.out
189       requires: defined(PETSC_HAVE_MPI_NONBLOCKING_COLLECTIVES)
190 
191    test:
192       suffix: ibarrier
193       nsize: 4
194       args: -verbose -build_twosided ibarrier
195       output_file: output/ex8_1.out
196       requires: defined(PETSC_HAVE_MPI_NONBLOCKING_COLLECTIVES)
197 
198    test:
199       suffix: redscatter
200       requires: mpi_reduce_scatter_block
201       nsize: 4
202       args: -verbose -build_twosided redscatter
203       output_file: output/ex8_1.out
204 
205 TEST*/
206