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
MakeDatatype(MPI_Datatype * dtype)11d71ae5a4SJacob Faibussowitsch static PetscErrorCode MakeDatatype(MPI_Datatype *dtype)
12d71ae5a4SJacob Faibussowitsch {
13c4762a1bSJed Brown MPI_Datatype dtypes[3], tmptype;
14c4762a1bSJed Brown PetscMPIInt lengths[3];
15c4762a1bSJed Brown MPI_Aint displs[3];
16c4762a1bSJed Brown Unit dummy;
17c4762a1bSJed Brown
18c4762a1bSJed Brown PetscFunctionBegin;
19c4762a1bSJed Brown dtypes[0] = MPIU_INT;
20c4762a1bSJed Brown dtypes[1] = MPIU_SCALAR;
21c4762a1bSJed Brown dtypes[2] = MPI_CHAR;
22c4762a1bSJed Brown lengths[0] = 1;
23c4762a1bSJed Brown lengths[1] = 1;
24c4762a1bSJed Brown lengths[2] = 3;
25c4762a1bSJed Brown /* Curse the evil beings that made std::complex a non-POD type. */
26c4762a1bSJed Brown displs[0] = (char *)&dummy.rank - (char *)&dummy; /* offsetof(Unit,rank); */
27c4762a1bSJed Brown displs[1] = (char *)&dummy.value - (char *)&dummy; /* offsetof(Unit,value); */
28c4762a1bSJed Brown displs[2] = (char *)&dummy.ok - (char *)&dummy; /* offsetof(Unit,ok); */
299566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_create_struct(3, lengths, displs, dtypes, &tmptype));
309566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&tmptype));
319566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_create_resized(tmptype, 0, sizeof(Unit), dtype));
329566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(dtype));
339566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&tmptype));
34c4762a1bSJed Brown {
35c4762a1bSJed Brown MPI_Aint lb, extent;
369566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_get_extent(*dtype, &lb, &extent));
3708401ef6SPierre Jolivet PetscCheck(extent == sizeof(Unit), PETSC_COMM_WORLD, PETSC_ERR_LIB, "New type has extent %d != sizeof(Unit) %d", (int)extent, (int)sizeof(Unit));
38c4762a1bSJed Brown }
393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
40c4762a1bSJed Brown }
41c4762a1bSJed Brown
42c4762a1bSJed Brown struct FCtx {
43c4762a1bSJed Brown PetscMPIInt rank;
44c4762a1bSJed Brown PetscMPIInt nto;
45c4762a1bSJed Brown PetscMPIInt *toranks;
46c4762a1bSJed Brown Unit *todata;
47c4762a1bSJed Brown PetscSegBuffer seg;
48c4762a1bSJed Brown };
49c4762a1bSJed Brown
FSend(MPI_Comm comm,const PetscMPIInt tag[],PetscMPIInt tonum,PetscMPIInt rank,void * todata,MPI_Request req[],PetscCtx ctx)50*2a8381b2SBarry Smith static PetscErrorCode FSend(MPI_Comm comm, const PetscMPIInt tag[], PetscMPIInt tonum, PetscMPIInt rank, void *todata, MPI_Request req[], PetscCtx ctx)
51d71ae5a4SJacob Faibussowitsch {
52c4762a1bSJed Brown struct FCtx *fctx = (struct FCtx *)ctx;
53c4762a1bSJed Brown
54c4762a1bSJed Brown PetscFunctionBegin;
5508401ef6SPierre Jolivet PetscCheck(rank == fctx->toranks[tonum], PETSC_COMM_SELF, PETSC_ERR_PLIB, "Rank %d does not match toranks[%d] %d", rank, tonum, fctx->toranks[tonum]);
56cc73adaaSBarry Smith PetscCheck(fctx->rank == *(PetscMPIInt *)todata, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Todata %d does not match rank %d", *(PetscMPIInt *)todata, fctx->rank);
579566063dSJacob Faibussowitsch PetscCallMPI(MPI_Isend(&fctx->todata[tonum].rank, 1, MPIU_INT, rank, tag[0], comm, &req[0]));
589566063dSJacob Faibussowitsch PetscCallMPI(MPI_Isend(&fctx->todata[tonum].value, 1, MPIU_SCALAR, rank, tag[1], comm, &req[1]));
593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
60c4762a1bSJed Brown }
61c4762a1bSJed Brown
FRecv(MPI_Comm comm,const PetscMPIInt tag[],PetscMPIInt rank,void * fromdata,MPI_Request req[],PetscCtx ctx)62*2a8381b2SBarry Smith static PetscErrorCode FRecv(MPI_Comm comm, const PetscMPIInt tag[], PetscMPIInt rank, void *fromdata, MPI_Request req[], PetscCtx ctx)
63d71ae5a4SJacob Faibussowitsch {
64c4762a1bSJed Brown struct FCtx *fctx = (struct FCtx *)ctx;
65c4762a1bSJed Brown Unit *buf;
66c4762a1bSJed Brown
67c4762a1bSJed Brown PetscFunctionBegin;
68cc73adaaSBarry Smith PetscCheck(*(PetscMPIInt *)fromdata == rank, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Dummy data %d from rank %d corrupt", *(PetscMPIInt *)fromdata, rank);
699566063dSJacob Faibussowitsch PetscCall(PetscSegBufferGet(fctx->seg, 1, &buf));
709566063dSJacob Faibussowitsch PetscCallMPI(MPI_Irecv(&buf->rank, 1, MPIU_INT, rank, tag[0], comm, &req[0]));
719566063dSJacob Faibussowitsch PetscCallMPI(MPI_Irecv(&buf->value, 1, MPIU_SCALAR, rank, tag[1], comm, &req[1]));
72c4762a1bSJed Brown buf->ok[0] = 'o';
73c4762a1bSJed Brown buf->ok[1] = 'k';
74c4762a1bSJed Brown buf->ok[2] = 0;
753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
76c4762a1bSJed Brown }
77c4762a1bSJed Brown
main(int argc,char ** argv)78d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
79d71ae5a4SJacob Faibussowitsch {
80c4762a1bSJed Brown PetscMPIInt rank, size, *toranks, *fromranks, nto, nfrom;
81c4762a1bSJed Brown PetscInt i, n;
82c4762a1bSJed Brown PetscBool verbose, build_twosided_f;
83c4762a1bSJed Brown Unit *todata, *fromdata;
84c4762a1bSJed Brown MPI_Datatype dtype;
85c4762a1bSJed Brown
86327415f7SBarry Smith PetscFunctionBeginUser;
87c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
889566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
899566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
90c4762a1bSJed Brown
91c4762a1bSJed Brown verbose = PETSC_FALSE;
929566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-verbose", &verbose, NULL));
93c4762a1bSJed Brown build_twosided_f = PETSC_FALSE;
949566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-build_twosided_f", &build_twosided_f, NULL));
95c4762a1bSJed Brown
96c4762a1bSJed Brown for (i = 1, nto = 0; i < size; i *= 2) nto++;
979566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nto, &todata, nto, &toranks));
98c4762a1bSJed Brown for (n = 0, i = 1; i < size; n++, i *= 2) {
99c4762a1bSJed Brown toranks[n] = (rank + i) % size;
100c4762a1bSJed Brown todata[n].rank = (rank + i) % size;
101c4762a1bSJed Brown todata[n].value = (PetscScalar)rank;
102c4762a1bSJed Brown todata[n].ok[0] = 'o';
103c4762a1bSJed Brown todata[n].ok[1] = 'k';
104c4762a1bSJed Brown todata[n].ok[2] = 0;
105c4762a1bSJed Brown }
106c4762a1bSJed Brown if (verbose) {
10748a46eb9SPierre Jolivet for (i = 0; i < nto; i++) PetscCall(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));
1089566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT));
109c4762a1bSJed Brown }
110c4762a1bSJed Brown
1119566063dSJacob Faibussowitsch PetscCall(MakeDatatype(&dtype));
112c4762a1bSJed Brown
113c4762a1bSJed Brown if (build_twosided_f) {
114c4762a1bSJed Brown struct FCtx fctx;
115c4762a1bSJed Brown PetscMPIInt *todummy, *fromdummy;
116c4762a1bSJed Brown fctx.rank = rank;
117c4762a1bSJed Brown fctx.nto = nto;
118c4762a1bSJed Brown fctx.toranks = toranks;
119c4762a1bSJed Brown fctx.todata = todata;
1209566063dSJacob Faibussowitsch PetscCall(PetscSegBufferCreate(sizeof(Unit), 1, &fctx.seg));
1219566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nto, &todummy));
122c4762a1bSJed Brown for (i = 0; i < nto; i++) todummy[i] = rank;
1239566063dSJacob Faibussowitsch PetscCall(PetscCommBuildTwoSidedF(PETSC_COMM_WORLD, 1, MPI_INT, nto, toranks, todummy, &nfrom, &fromranks, &fromdummy, 2, FSend, FRecv, &fctx));
1249566063dSJacob Faibussowitsch PetscCall(PetscFree(todummy));
1259566063dSJacob Faibussowitsch PetscCall(PetscFree(fromdummy));
1269566063dSJacob Faibussowitsch PetscCall(PetscSegBufferExtractAlloc(fctx.seg, &fromdata));
1279566063dSJacob Faibussowitsch PetscCall(PetscSegBufferDestroy(&fctx.seg));
128c4762a1bSJed Brown } else {
1299566063dSJacob Faibussowitsch PetscCall(PetscCommBuildTwoSided(PETSC_COMM_WORLD, 1, dtype, nto, toranks, todata, &nfrom, &fromranks, &fromdata));
130c4762a1bSJed Brown }
1319566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&dtype));
132c4762a1bSJed Brown
133c4762a1bSJed Brown if (verbose) {
134c4762a1bSJed Brown PetscInt *iranks, *iperm;
1359566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nfrom, &iranks, nfrom, &iperm));
136c4762a1bSJed Brown for (i = 0; i < nfrom; i++) {
137c4762a1bSJed Brown iranks[i] = fromranks[i];
138c4762a1bSJed Brown iperm[i] = i;
139c4762a1bSJed Brown }
140c4762a1bSJed Brown /* Receive ordering is non-deterministic in general, so sort to make verbose output deterministic. */
1419566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithPermutation(nfrom, iranks, iperm));
142c4762a1bSJed Brown for (i = 0; i < nfrom; i++) {
143c4762a1bSJed Brown PetscInt ip = iperm[i];
1449566063dSJacob Faibussowitsch PetscCall(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));
145c4762a1bSJed Brown }
1469566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT));
1479566063dSJacob Faibussowitsch PetscCall(PetscFree2(iranks, iperm));
148c4762a1bSJed Brown }
149c4762a1bSJed Brown
15008401ef6SPierre Jolivet PetscCheck(nto == nfrom, PETSC_COMM_SELF, PETSC_ERR_PLIB, "[%d] From ranks %d does not match To ranks %d", rank, nto, nfrom);
151c4762a1bSJed Brown for (i = 1; i < size; i *= 2) {
152c4762a1bSJed Brown PetscMPIInt expected_rank = (rank - i + size) % size;
153c4762a1bSJed Brown PetscBool flg;
154c4762a1bSJed Brown for (n = 0; n < nfrom; n++) {
155c4762a1bSJed Brown if (expected_rank == fromranks[n]) goto found;
156c4762a1bSJed Brown }
15798921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "[%d] Could not find expected from rank %d", rank, expected_rank);
158c4762a1bSJed Brown found:
15908401ef6SPierre Jolivet PetscCheck(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);
1609566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(fromdata[n].ok, "ok", &flg));
16128b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_PLIB, "[%d] Got string %s from rank %d", rank, fromdata[n].ok, expected_rank);
162c4762a1bSJed Brown }
1639566063dSJacob Faibussowitsch PetscCall(PetscFree2(todata, toranks));
1649566063dSJacob Faibussowitsch PetscCall(PetscFree(fromdata));
1659566063dSJacob Faibussowitsch PetscCall(PetscFree(fromranks));
1669566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
167b122ec5aSJacob Faibussowitsch return 0;
168c4762a1bSJed Brown }
169c4762a1bSJed Brown
170c4762a1bSJed Brown /*TEST
171c4762a1bSJed Brown
172c4762a1bSJed Brown test:
173c4762a1bSJed Brown nsize: 4
174c4762a1bSJed Brown args: -verbose -build_twosided allreduce
175c4762a1bSJed Brown
176c4762a1bSJed Brown test:
177c4762a1bSJed Brown suffix: f
178c4762a1bSJed Brown nsize: 4
179c4762a1bSJed Brown args: -verbose -build_twosided_f -build_twosided allreduce
180c4762a1bSJed Brown output_file: output/ex8_1.out
181c4762a1bSJed Brown
182c4762a1bSJed Brown test:
183c4762a1bSJed Brown suffix: f_ibarrier
184c4762a1bSJed Brown nsize: 4
185c4762a1bSJed Brown args: -verbose -build_twosided_f -build_twosided ibarrier
186c4762a1bSJed Brown output_file: output/ex8_1.out
1878ee3e7ecSJunchao Zhang requires: defined(PETSC_HAVE_MPI_NONBLOCKING_COLLECTIVES)
188c4762a1bSJed Brown
189c4762a1bSJed Brown test:
190c4762a1bSJed Brown suffix: ibarrier
191c4762a1bSJed Brown nsize: 4
192c4762a1bSJed Brown args: -verbose -build_twosided ibarrier
193c4762a1bSJed Brown output_file: output/ex8_1.out
1948ee3e7ecSJunchao Zhang requires: defined(PETSC_HAVE_MPI_NONBLOCKING_COLLECTIVES)
195c4762a1bSJed Brown
196c4762a1bSJed Brown test:
197c4762a1bSJed Brown suffix: redscatter
198c4762a1bSJed Brown requires: mpi_reduce_scatter_block
199c4762a1bSJed Brown nsize: 4
200c4762a1bSJed Brown args: -verbose -build_twosided redscatter
201c4762a1bSJed Brown output_file: output/ex8_1.out
202c4762a1bSJed Brown
203c4762a1bSJed Brown TEST*/
204