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