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 119371c9d4SSatish Balay static PetscErrorCode MakeDatatype(MPI_Datatype *dtype) { 12c4762a1bSJed Brown MPI_Datatype dtypes[3], tmptype; 13c4762a1bSJed Brown PetscMPIInt lengths[3]; 14c4762a1bSJed Brown MPI_Aint displs[3]; 15c4762a1bSJed Brown Unit dummy; 16c4762a1bSJed Brown 17c4762a1bSJed Brown PetscFunctionBegin; 18c4762a1bSJed Brown dtypes[0] = MPIU_INT; 19c4762a1bSJed Brown dtypes[1] = MPIU_SCALAR; 20c4762a1bSJed Brown dtypes[2] = MPI_CHAR; 21c4762a1bSJed Brown lengths[0] = 1; 22c4762a1bSJed Brown lengths[1] = 1; 23c4762a1bSJed Brown lengths[2] = 3; 24c4762a1bSJed Brown /* Curse the evil beings that made std::complex a non-POD type. */ 25c4762a1bSJed Brown displs[0] = (char *)&dummy.rank - (char *)&dummy; /* offsetof(Unit,rank); */ 26c4762a1bSJed Brown displs[1] = (char *)&dummy.value - (char *)&dummy; /* offsetof(Unit,value); */ 27c4762a1bSJed Brown displs[2] = (char *)&dummy.ok - (char *)&dummy; /* offsetof(Unit,ok); */ 289566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_create_struct(3, lengths, displs, dtypes, &tmptype)); 299566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(&tmptype)); 309566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_create_resized(tmptype, 0, sizeof(Unit), dtype)); 319566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_commit(dtype)); 329566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&tmptype)); 33c4762a1bSJed Brown { 34c4762a1bSJed Brown MPI_Aint lb, extent; 359566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_get_extent(*dtype, &lb, &extent)); 3608401ef6SPierre Jolivet PetscCheck(extent == sizeof(Unit), PETSC_COMM_WORLD, PETSC_ERR_LIB, "New type has extent %d != sizeof(Unit) %d", (int)extent, (int)sizeof(Unit)); 37c4762a1bSJed Brown } 38c4762a1bSJed Brown PetscFunctionReturn(0); 39c4762a1bSJed Brown } 40c4762a1bSJed Brown 41c4762a1bSJed Brown struct FCtx { 42c4762a1bSJed Brown PetscMPIInt rank; 43c4762a1bSJed Brown PetscMPIInt nto; 44c4762a1bSJed Brown PetscMPIInt *toranks; 45c4762a1bSJed Brown Unit *todata; 46c4762a1bSJed Brown PetscSegBuffer seg; 47c4762a1bSJed Brown }; 48c4762a1bSJed Brown 499371c9d4SSatish Balay static PetscErrorCode FSend(MPI_Comm comm, const PetscMPIInt tag[], PetscMPIInt tonum, PetscMPIInt rank, void *todata, MPI_Request req[], void *ctx) { 50c4762a1bSJed Brown struct FCtx *fctx = (struct FCtx *)ctx; 51c4762a1bSJed Brown 52c4762a1bSJed Brown PetscFunctionBegin; 5308401ef6SPierre 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]); 54cc73adaaSBarry Smith PetscCheck(fctx->rank == *(PetscMPIInt *)todata, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Todata %d does not match rank %d", *(PetscMPIInt *)todata, fctx->rank); 559566063dSJacob Faibussowitsch PetscCallMPI(MPI_Isend(&fctx->todata[tonum].rank, 1, MPIU_INT, rank, tag[0], comm, &req[0])); 569566063dSJacob Faibussowitsch PetscCallMPI(MPI_Isend(&fctx->todata[tonum].value, 1, MPIU_SCALAR, rank, tag[1], comm, &req[1])); 57c4762a1bSJed Brown PetscFunctionReturn(0); 58c4762a1bSJed Brown } 59c4762a1bSJed Brown 609371c9d4SSatish Balay static PetscErrorCode FRecv(MPI_Comm comm, const PetscMPIInt tag[], PetscMPIInt rank, void *fromdata, MPI_Request req[], void *ctx) { 61c4762a1bSJed Brown struct FCtx *fctx = (struct FCtx *)ctx; 62c4762a1bSJed Brown Unit *buf; 63c4762a1bSJed Brown 64c4762a1bSJed Brown PetscFunctionBegin; 65cc73adaaSBarry Smith PetscCheck(*(PetscMPIInt *)fromdata == rank, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Dummy data %d from rank %d corrupt", *(PetscMPIInt *)fromdata, rank); 669566063dSJacob Faibussowitsch PetscCall(PetscSegBufferGet(fctx->seg, 1, &buf)); 679566063dSJacob Faibussowitsch PetscCallMPI(MPI_Irecv(&buf->rank, 1, MPIU_INT, rank, tag[0], comm, &req[0])); 689566063dSJacob Faibussowitsch PetscCallMPI(MPI_Irecv(&buf->value, 1, MPIU_SCALAR, rank, tag[1], comm, &req[1])); 69c4762a1bSJed Brown buf->ok[0] = 'o'; 70c4762a1bSJed Brown buf->ok[1] = 'k'; 71c4762a1bSJed Brown buf->ok[2] = 0; 72c4762a1bSJed Brown PetscFunctionReturn(0); 73c4762a1bSJed Brown } 74c4762a1bSJed Brown 759371c9d4SSatish Balay int main(int argc, char **argv) { 76c4762a1bSJed Brown PetscMPIInt rank, size, *toranks, *fromranks, nto, nfrom; 77c4762a1bSJed Brown PetscInt i, n; 78c4762a1bSJed Brown PetscBool verbose, build_twosided_f; 79c4762a1bSJed Brown Unit *todata, *fromdata; 80c4762a1bSJed Brown MPI_Datatype dtype; 81c4762a1bSJed Brown 82327415f7SBarry Smith PetscFunctionBeginUser; 839566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 849566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 859566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 86c4762a1bSJed Brown 87c4762a1bSJed Brown verbose = PETSC_FALSE; 889566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-verbose", &verbose, NULL)); 89c4762a1bSJed Brown build_twosided_f = PETSC_FALSE; 909566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-build_twosided_f", &build_twosided_f, NULL)); 91c4762a1bSJed Brown 92c4762a1bSJed Brown for (i = 1, nto = 0; i < size; i *= 2) nto++; 939566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nto, &todata, nto, &toranks)); 94c4762a1bSJed Brown for (n = 0, i = 1; i < size; n++, i *= 2) { 95c4762a1bSJed Brown toranks[n] = (rank + i) % size; 96c4762a1bSJed Brown todata[n].rank = (rank + i) % size; 97c4762a1bSJed Brown todata[n].value = (PetscScalar)rank; 98c4762a1bSJed Brown todata[n].ok[0] = 'o'; 99c4762a1bSJed Brown todata[n].ok[1] = 'k'; 100c4762a1bSJed Brown todata[n].ok[2] = 0; 101c4762a1bSJed Brown } 102c4762a1bSJed Brown if (verbose) { 103*48a46eb9SPierre 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)); 1049566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 105c4762a1bSJed Brown } 106c4762a1bSJed Brown 1079566063dSJacob Faibussowitsch PetscCall(MakeDatatype(&dtype)); 108c4762a1bSJed Brown 109c4762a1bSJed Brown if (build_twosided_f) { 110c4762a1bSJed Brown struct FCtx fctx; 111c4762a1bSJed Brown PetscMPIInt *todummy, *fromdummy; 112c4762a1bSJed Brown fctx.rank = rank; 113c4762a1bSJed Brown fctx.nto = nto; 114c4762a1bSJed Brown fctx.toranks = toranks; 115c4762a1bSJed Brown fctx.todata = todata; 1169566063dSJacob Faibussowitsch PetscCall(PetscSegBufferCreate(sizeof(Unit), 1, &fctx.seg)); 1179566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nto, &todummy)); 118c4762a1bSJed Brown for (i = 0; i < nto; i++) todummy[i] = rank; 1199566063dSJacob Faibussowitsch PetscCall(PetscCommBuildTwoSidedF(PETSC_COMM_WORLD, 1, MPI_INT, nto, toranks, todummy, &nfrom, &fromranks, &fromdummy, 2, FSend, FRecv, &fctx)); 1209566063dSJacob Faibussowitsch PetscCall(PetscFree(todummy)); 1219566063dSJacob Faibussowitsch PetscCall(PetscFree(fromdummy)); 1229566063dSJacob Faibussowitsch PetscCall(PetscSegBufferExtractAlloc(fctx.seg, &fromdata)); 1239566063dSJacob Faibussowitsch PetscCall(PetscSegBufferDestroy(&fctx.seg)); 124c4762a1bSJed Brown } else { 1259566063dSJacob Faibussowitsch PetscCall(PetscCommBuildTwoSided(PETSC_COMM_WORLD, 1, dtype, nto, toranks, todata, &nfrom, &fromranks, &fromdata)); 126c4762a1bSJed Brown } 1279566063dSJacob Faibussowitsch PetscCallMPI(MPI_Type_free(&dtype)); 128c4762a1bSJed Brown 129c4762a1bSJed Brown if (verbose) { 130c4762a1bSJed Brown PetscInt *iranks, *iperm; 1319566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nfrom, &iranks, nfrom, &iperm)); 132c4762a1bSJed Brown for (i = 0; i < nfrom; i++) { 133c4762a1bSJed Brown iranks[i] = fromranks[i]; 134c4762a1bSJed Brown iperm[i] = i; 135c4762a1bSJed Brown } 136c4762a1bSJed Brown /* Receive ordering is non-deterministic in general, so sort to make verbose output deterministic. */ 1379566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithPermutation(nfrom, iranks, iperm)); 138c4762a1bSJed Brown for (i = 0; i < nfrom; i++) { 139c4762a1bSJed Brown PetscInt ip = iperm[i]; 1409566063dSJacob 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)); 141c4762a1bSJed Brown } 1429566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 1439566063dSJacob Faibussowitsch PetscCall(PetscFree2(iranks, iperm)); 144c4762a1bSJed Brown } 145c4762a1bSJed Brown 14608401ef6SPierre Jolivet PetscCheck(nto == nfrom, PETSC_COMM_SELF, PETSC_ERR_PLIB, "[%d] From ranks %d does not match To ranks %d", rank, nto, nfrom); 147c4762a1bSJed Brown for (i = 1; i < size; i *= 2) { 148c4762a1bSJed Brown PetscMPIInt expected_rank = (rank - i + size) % size; 149c4762a1bSJed Brown PetscBool flg; 150c4762a1bSJed Brown for (n = 0; n < nfrom; n++) { 151c4762a1bSJed Brown if (expected_rank == fromranks[n]) goto found; 152c4762a1bSJed Brown } 15398921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "[%d] Could not find expected from rank %d", rank, expected_rank); 154c4762a1bSJed Brown found: 15508401ef6SPierre 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); 1569566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(fromdata[n].ok, "ok", &flg)); 15728b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_PLIB, "[%d] Got string %s from rank %d", rank, fromdata[n].ok, expected_rank); 158c4762a1bSJed Brown } 1599566063dSJacob Faibussowitsch PetscCall(PetscFree2(todata, toranks)); 1609566063dSJacob Faibussowitsch PetscCall(PetscFree(fromdata)); 1619566063dSJacob Faibussowitsch PetscCall(PetscFree(fromranks)); 1629566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 163b122ec5aSJacob Faibussowitsch return 0; 164c4762a1bSJed Brown } 165c4762a1bSJed Brown 166c4762a1bSJed Brown /*TEST 167c4762a1bSJed Brown 168c4762a1bSJed Brown test: 169c4762a1bSJed Brown nsize: 4 170c4762a1bSJed Brown args: -verbose -build_twosided allreduce 171c4762a1bSJed Brown 172c4762a1bSJed Brown test: 173c4762a1bSJed Brown suffix: f 174c4762a1bSJed Brown nsize: 4 175c4762a1bSJed Brown args: -verbose -build_twosided_f -build_twosided allreduce 176c4762a1bSJed Brown output_file: output/ex8_1.out 177c4762a1bSJed Brown 178c4762a1bSJed Brown test: 179c4762a1bSJed Brown suffix: f_ibarrier 180c4762a1bSJed Brown nsize: 4 181c4762a1bSJed Brown args: -verbose -build_twosided_f -build_twosided ibarrier 182c4762a1bSJed Brown output_file: output/ex8_1.out 1838ee3e7ecSJunchao Zhang requires: defined(PETSC_HAVE_MPI_NONBLOCKING_COLLECTIVES) 184c4762a1bSJed Brown 185c4762a1bSJed Brown test: 186c4762a1bSJed Brown suffix: ibarrier 187c4762a1bSJed Brown nsize: 4 188c4762a1bSJed Brown args: -verbose -build_twosided ibarrier 189c4762a1bSJed Brown output_file: output/ex8_1.out 1908ee3e7ecSJunchao Zhang requires: defined(PETSC_HAVE_MPI_NONBLOCKING_COLLECTIVES) 191c4762a1bSJed Brown 192c4762a1bSJed Brown test: 193c4762a1bSJed Brown suffix: redscatter 194c4762a1bSJed Brown requires: mpi_reduce_scatter_block 195c4762a1bSJed Brown nsize: 4 196c4762a1bSJed Brown args: -verbose -build_twosided redscatter 197c4762a1bSJed Brown output_file: output/ex8_1.out 198c4762a1bSJed Brown 199c4762a1bSJed Brown TEST*/ 200