1fd697ae2SVaclav Hapla 2fd697ae2SVaclav Hapla static char help[] = "Test leaf sorting in PetscSFSetGraph()\n\n"; 3fd697ae2SVaclav Hapla 4fd697ae2SVaclav Hapla #include <petscsf.h> 5fd697ae2SVaclav Hapla 6fd697ae2SVaclav Hapla typedef struct { 7fd697ae2SVaclav Hapla MPI_Comm comm; 8fd697ae2SVaclav Hapla PetscMPIInt rank, size; 9fd697ae2SVaclav Hapla PetscInt leaveStep, nLeavesPerRank; 10fd697ae2SVaclav Hapla PetscBool contiguousLeaves; 11fd697ae2SVaclav Hapla PetscCopyMode localmode, remotemode; 12fd697ae2SVaclav Hapla PetscInt *ilocal; 13fd697ae2SVaclav Hapla PetscSFNode *iremote; 14fd697ae2SVaclav Hapla } AppCtx; 15fd697ae2SVaclav Hapla 169371c9d4SSatish Balay static PetscErrorCode GetOptions(MPI_Comm comm, AppCtx *ctx) { 17fd697ae2SVaclav Hapla PetscFunctionBegin; 18fd697ae2SVaclav Hapla ctx->comm = comm; 19fd697ae2SVaclav Hapla ctx->nLeavesPerRank = 4; 20fd697ae2SVaclav Hapla ctx->leaveStep = 1; 21fd697ae2SVaclav Hapla ctx->contiguousLeaves = PETSC_FALSE; 22fd697ae2SVaclav Hapla ctx->localmode = PETSC_OWN_POINTER; 23fd697ae2SVaclav Hapla ctx->remotemode = PETSC_OWN_POINTER; 24fd697ae2SVaclav Hapla ctx->ilocal = NULL; 25fd697ae2SVaclav Hapla ctx->iremote = NULL; 269566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-n_leaves_per_rank", &ctx->nLeavesPerRank, NULL)); 279566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-leave_step", &ctx->leaveStep, NULL)); 289566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetEnum(NULL, NULL, "-localmode", PetscCopyModes, (PetscEnum *)&ctx->localmode, NULL)); 299566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetEnum(NULL, NULL, "-remotemode", PetscCopyModes, (PetscEnum *)&ctx->remotemode, NULL)); 30fd697ae2SVaclav Hapla ctx->contiguousLeaves = (PetscBool)(ctx->leaveStep == 1); 319566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &ctx->size)); 329566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &ctx->rank)); 33fd697ae2SVaclav Hapla PetscFunctionReturn(0); 34fd697ae2SVaclav Hapla } 35fd697ae2SVaclav Hapla 369371c9d4SSatish Balay static PetscErrorCode PetscSFCheckEqual_Private(PetscSF sf0, PetscSF sf1) { 37fd697ae2SVaclav Hapla PetscInt nRoot, nLeave; 38fd697ae2SVaclav Hapla Vec vecRoot0, vecLeave0, vecRoot1, vecLeave1; 39fd697ae2SVaclav Hapla MPI_Comm comm; 40fd697ae2SVaclav Hapla PetscBool flg; 41fd697ae2SVaclav Hapla 42fd697ae2SVaclav Hapla PetscFunctionBegin; 439566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)sf0, &comm)); 449566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf0, &nRoot, NULL, NULL, NULL)); 459566063dSJacob Faibussowitsch PetscCall(PetscSFGetLeafRange(sf0, NULL, &nLeave)); 46fd697ae2SVaclav Hapla nLeave++; 479566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(comm, nRoot, PETSC_DECIDE, &vecRoot0)); 489566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(comm, nLeave, PETSC_DECIDE, &vecLeave0)); 499566063dSJacob Faibussowitsch PetscCall(VecDuplicate(vecRoot0, &vecRoot1)); 509566063dSJacob Faibussowitsch PetscCall(VecDuplicate(vecLeave0, &vecLeave1)); 51fd697ae2SVaclav Hapla { 52fd697ae2SVaclav Hapla PetscRandom rand; 53fd697ae2SVaclav Hapla 549566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(comm, &rand)); 559566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rand)); 569566063dSJacob Faibussowitsch PetscCall(VecSetRandom(vecRoot0, rand)); 579566063dSJacob Faibussowitsch PetscCall(VecSetRandom(vecLeave0, rand)); 589566063dSJacob Faibussowitsch PetscCall(VecCopy(vecRoot0, vecRoot1)); 599566063dSJacob Faibussowitsch PetscCall(VecCopy(vecLeave0, vecLeave1)); 609566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rand)); 61fd697ae2SVaclav Hapla } 62fd697ae2SVaclav Hapla 639566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(sf0, vecRoot0, vecLeave0, ADD_VALUES, SCATTER_FORWARD)); 649566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(sf0, vecRoot0, vecLeave0, ADD_VALUES, SCATTER_FORWARD)); 659566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(sf1, vecRoot1, vecLeave1, ADD_VALUES, SCATTER_FORWARD)); 669566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(sf1, vecRoot1, vecLeave1, ADD_VALUES, SCATTER_FORWARD)); 679566063dSJacob Faibussowitsch PetscCall(VecEqual(vecLeave0, vecLeave1, &flg)); 68fd697ae2SVaclav Hapla PetscCheck(flg, comm, PETSC_ERR_PLIB, "leave vectors differ"); 69fd697ae2SVaclav Hapla 709566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(sf0, vecLeave0, vecRoot0, ADD_VALUES, SCATTER_REVERSE)); 719566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(sf0, vecLeave0, vecRoot0, ADD_VALUES, SCATTER_REVERSE)); 729566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(sf1, vecLeave1, vecRoot1, ADD_VALUES, SCATTER_REVERSE)); 739566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(sf1, vecLeave1, vecRoot1, ADD_VALUES, SCATTER_REVERSE)); 749566063dSJacob Faibussowitsch PetscCall(VecEqual(vecRoot0, vecRoot1, &flg)); 75fd697ae2SVaclav Hapla PetscCheck(flg, comm, PETSC_ERR_PLIB, "root vectors differ"); 76fd697ae2SVaclav Hapla 779566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vecRoot0)); 789566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vecRoot1)); 799566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vecLeave0)); 809566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vecLeave1)); 81fd697ae2SVaclav Hapla PetscFunctionReturn(0); 82fd697ae2SVaclav Hapla } 83fd697ae2SVaclav Hapla 849371c9d4SSatish Balay PetscErrorCode CreateSF0(AppCtx *ctx, PetscSF *sf0) { 85fd697ae2SVaclav Hapla PetscInt j, k, r; 86fd697ae2SVaclav Hapla PetscInt nLeaves = ctx->nLeavesPerRank * ctx->size; 87fd697ae2SVaclav Hapla PetscInt nroots = ctx->nLeavesPerRank; 88fd697ae2SVaclav Hapla PetscSF sf; 89fd697ae2SVaclav Hapla PetscInt *ilocal; 90fd697ae2SVaclav Hapla PetscSFNode *iremote; 91fd697ae2SVaclav Hapla 92fd697ae2SVaclav Hapla PetscFunctionBegin; 939566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nLeaves + 1, &ctx->ilocal)); 949566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nLeaves, &ctx->iremote)); 95fd697ae2SVaclav Hapla ilocal = ctx->ilocal; 96fd697ae2SVaclav Hapla iremote = ctx->iremote; 97fd697ae2SVaclav Hapla ilocal[nLeaves] = -ctx->leaveStep; 989566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(ctx->comm, &sf)); 99fd697ae2SVaclav Hapla for (r = 0, j = nLeaves - 1; r < ctx->size; r++) { 100fd697ae2SVaclav Hapla for (k = 0; k < ctx->nLeavesPerRank; k++, j--) { 101fd697ae2SVaclav Hapla ilocal[j] = ilocal[j + 1] + ctx->leaveStep; 102fd697ae2SVaclav Hapla iremote[j].rank = r; 103fd697ae2SVaclav Hapla iremote[j].index = k; 104fd697ae2SVaclav Hapla } 105fd697ae2SVaclav Hapla } 1069566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(sf, nroots, nLeaves, ilocal, ctx->localmode, iremote, ctx->remotemode)); 107fd697ae2SVaclav Hapla { 108fd697ae2SVaclav Hapla const PetscInt *tlocal; 109fd697ae2SVaclav Hapla PetscBool sorted; 110fd697ae2SVaclav Hapla 1119566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, NULL, NULL, &tlocal, NULL)); 11208401ef6SPierre Jolivet PetscCheck(!ctx->contiguousLeaves || !tlocal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ilocal=NULL expected for contiguous case"); 113fd697ae2SVaclav Hapla if (tlocal) { 1149566063dSJacob Faibussowitsch PetscCall(PetscSortedInt(nLeaves, tlocal, &sorted)); 115fd697ae2SVaclav Hapla PetscCheck(sorted, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ilocal expected to be sorted"); 116fd697ae2SVaclav Hapla } 117fd697ae2SVaclav Hapla } 118fd697ae2SVaclav Hapla *sf0 = sf; 119fd697ae2SVaclav Hapla PetscFunctionReturn(0); 120fd697ae2SVaclav Hapla } 121fd697ae2SVaclav Hapla 1229371c9d4SSatish Balay PetscErrorCode CreateSF1(AppCtx *ctx, PetscSF *sf1) { 123fd697ae2SVaclav Hapla PetscInt j, k, r; 124fd697ae2SVaclav Hapla PetscInt *ilocal = NULL; 125fd697ae2SVaclav Hapla PetscSFNode *iremote; 126fd697ae2SVaclav Hapla PetscInt nLeaves = ctx->nLeavesPerRank * ctx->size; 127fd697ae2SVaclav Hapla PetscInt nroots = ctx->nLeavesPerRank; 128fd697ae2SVaclav Hapla PetscSF sf; 129fd697ae2SVaclav Hapla 130fd697ae2SVaclav Hapla PetscFunctionBegin; 131fd697ae2SVaclav Hapla ilocal = NULL; 132*48a46eb9SPierre Jolivet if (!ctx->contiguousLeaves) PetscCall(PetscCalloc1(nLeaves + 1, &ilocal)); 1339566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nLeaves, &iremote)); 1349566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(ctx->comm, &sf)); 135fd697ae2SVaclav Hapla for (r = 0, j = 0; r < ctx->size; r++) { 136fd697ae2SVaclav Hapla for (k = 0; k < ctx->nLeavesPerRank; k++, j++) { 1379371c9d4SSatish Balay if (!ctx->contiguousLeaves) { ilocal[j + 1] = ilocal[j] + ctx->leaveStep; } 138fd697ae2SVaclav Hapla iremote[j].rank = r; 139fd697ae2SVaclav Hapla iremote[j].index = k; 140fd697ae2SVaclav Hapla } 141fd697ae2SVaclav Hapla } 142fd697ae2SVaclav Hapla PetscCheck(j == nLeaves, PETSC_COMM_SELF, PETSC_ERR_PLIB, "j != nLeaves"); 1439566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(sf, nroots, nLeaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER)); 144fd697ae2SVaclav Hapla if (ctx->contiguousLeaves) { 145fd697ae2SVaclav Hapla const PetscInt *tlocal; 146fd697ae2SVaclav Hapla 1479566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, NULL, NULL, &tlocal, NULL)); 148c9cc58a2SBarry Smith PetscCheck(!tlocal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ilocal=NULL expected for contiguous case"); 149fd697ae2SVaclav Hapla } 150fd697ae2SVaclav Hapla *sf1 = sf; 151fd697ae2SVaclav Hapla PetscFunctionReturn(0); 152fd697ae2SVaclav Hapla } 153fd697ae2SVaclav Hapla 1549371c9d4SSatish Balay int main(int argc, char **argv) { 155fd697ae2SVaclav Hapla AppCtx ctx; 156fd697ae2SVaclav Hapla PetscSF sf0, sf1; 157fd697ae2SVaclav Hapla MPI_Comm comm; 158fd697ae2SVaclav Hapla 159327415f7SBarry Smith PetscFunctionBeginUser; 1609566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 161fd697ae2SVaclav Hapla comm = PETSC_COMM_WORLD; 1629566063dSJacob Faibussowitsch PetscCall(GetOptions(comm, &ctx)); 163fd697ae2SVaclav Hapla 1649566063dSJacob Faibussowitsch PetscCall(CreateSF0(&ctx, &sf0)); 1659566063dSJacob Faibussowitsch PetscCall(CreateSF1(&ctx, &sf1)); 1669566063dSJacob Faibussowitsch PetscCall(PetscSFViewFromOptions(sf0, NULL, "-sf0_view")); 1679566063dSJacob Faibussowitsch PetscCall(PetscSFViewFromOptions(sf1, NULL, "-sf1_view")); 1689566063dSJacob Faibussowitsch PetscCall(PetscSFCheckEqual_Private(sf0, sf1)); 169fd697ae2SVaclav Hapla 1709566063dSJacob Faibussowitsch if (ctx.localmode != PETSC_OWN_POINTER) PetscCall(PetscFree(ctx.ilocal)); 1719566063dSJacob Faibussowitsch if (ctx.remotemode != PETSC_OWN_POINTER) PetscCall(PetscFree(ctx.iremote)); 1729566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sf0)); 1739566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sf1)); 1749566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 175b122ec5aSJacob Faibussowitsch return 0; 176fd697ae2SVaclav Hapla } 177fd697ae2SVaclav Hapla 178fd697ae2SVaclav Hapla /*TEST 179fd697ae2SVaclav Hapla testset: 180fd697ae2SVaclav Hapla suffix: 1 181fd697ae2SVaclav Hapla nsize: {{1 3}} 182fd697ae2SVaclav Hapla args: -n_leaves_per_rank {{0 5}} -leave_step {{1 3}} 183fd697ae2SVaclav Hapla test: 184fd697ae2SVaclav Hapla suffix: a 185fd697ae2SVaclav Hapla args: -localmode {{COPY_VALUES OWN_POINTER}} -remotemode {{COPY_VALUES OWN_POINTER}} 186fd697ae2SVaclav Hapla test: 187fd697ae2SVaclav Hapla suffix: b 188fd697ae2SVaclav Hapla args: -localmode USE_POINTER -remotemode {{COPY_VALUES OWN_POINTER USE_POINTER}} 189fd697ae2SVaclav Hapla TEST*/ 190