1fd697ae2SVaclav Hapla static char help[] = "Test leaf sorting in PetscSFSetGraph()\n\n"; 2fd697ae2SVaclav Hapla 3fd697ae2SVaclav Hapla #include <petscsf.h> 4fd697ae2SVaclav Hapla 5fd697ae2SVaclav Hapla typedef struct { 6fd697ae2SVaclav Hapla MPI_Comm comm; 7fd697ae2SVaclav Hapla PetscMPIInt rank, size; 8fd697ae2SVaclav Hapla PetscInt leaveStep, nLeavesPerRank; 9fd697ae2SVaclav Hapla PetscBool contiguousLeaves; 10fd697ae2SVaclav Hapla PetscCopyMode localmode, remotemode; 11fd697ae2SVaclav Hapla PetscInt *ilocal; 12fd697ae2SVaclav Hapla PetscSFNode *iremote; 13fd697ae2SVaclav Hapla } AppCtx; 14fd697ae2SVaclav Hapla 15d71ae5a4SJacob Faibussowitsch static PetscErrorCode GetOptions(MPI_Comm comm, AppCtx *ctx) 16d71ae5a4SJacob Faibussowitsch { 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)); 333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 34fd697ae2SVaclav Hapla } 35fd697ae2SVaclav Hapla 36d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSFCheckEqual_Private(PetscSF sf0, PetscSF sf1) 37d71ae5a4SJacob Faibussowitsch { 38fd697ae2SVaclav Hapla PetscInt nRoot, nLeave; 39fd697ae2SVaclav Hapla Vec vecRoot0, vecLeave0, vecRoot1, vecLeave1; 40fd697ae2SVaclav Hapla MPI_Comm comm; 41fd697ae2SVaclav Hapla PetscBool flg; 42fd697ae2SVaclav Hapla 43fd697ae2SVaclav Hapla PetscFunctionBegin; 449566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)sf0, &comm)); 459566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf0, &nRoot, NULL, NULL, NULL)); 469566063dSJacob Faibussowitsch PetscCall(PetscSFGetLeafRange(sf0, NULL, &nLeave)); 47fd697ae2SVaclav Hapla nLeave++; 4877433607SBarry Smith PetscCall(VecCreateFromOptions(comm, NULL, 1, nRoot, PETSC_DECIDE, &vecRoot0)); 4977433607SBarry Smith PetscCall(VecCreateFromOptions(comm, NULL, 1, nLeave, PETSC_DECIDE, &vecLeave0)); 509566063dSJacob Faibussowitsch PetscCall(VecDuplicate(vecRoot0, &vecRoot1)); 519566063dSJacob Faibussowitsch PetscCall(VecDuplicate(vecLeave0, &vecLeave1)); 52fd697ae2SVaclav Hapla { 53fd697ae2SVaclav Hapla PetscRandom rand; 54fd697ae2SVaclav Hapla 559566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(comm, &rand)); 569566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rand)); 579566063dSJacob Faibussowitsch PetscCall(VecSetRandom(vecRoot0, rand)); 589566063dSJacob Faibussowitsch PetscCall(VecSetRandom(vecLeave0, rand)); 599566063dSJacob Faibussowitsch PetscCall(VecCopy(vecRoot0, vecRoot1)); 609566063dSJacob Faibussowitsch PetscCall(VecCopy(vecLeave0, vecLeave1)); 619566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rand)); 62fd697ae2SVaclav Hapla } 63fd697ae2SVaclav Hapla 649566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(sf0, vecRoot0, vecLeave0, ADD_VALUES, SCATTER_FORWARD)); 659566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(sf0, vecRoot0, vecLeave0, ADD_VALUES, SCATTER_FORWARD)); 669566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(sf1, vecRoot1, vecLeave1, ADD_VALUES, SCATTER_FORWARD)); 679566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(sf1, vecRoot1, vecLeave1, ADD_VALUES, SCATTER_FORWARD)); 689566063dSJacob Faibussowitsch PetscCall(VecEqual(vecLeave0, vecLeave1, &flg)); 69fd697ae2SVaclav Hapla PetscCheck(flg, comm, PETSC_ERR_PLIB, "leave vectors differ"); 70fd697ae2SVaclav Hapla 719566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(sf0, vecLeave0, vecRoot0, ADD_VALUES, SCATTER_REVERSE)); 729566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(sf0, vecLeave0, vecRoot0, ADD_VALUES, SCATTER_REVERSE)); 739566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(sf1, vecLeave1, vecRoot1, ADD_VALUES, SCATTER_REVERSE)); 749566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(sf1, vecLeave1, vecRoot1, ADD_VALUES, SCATTER_REVERSE)); 759566063dSJacob Faibussowitsch PetscCall(VecEqual(vecRoot0, vecRoot1, &flg)); 76fd697ae2SVaclav Hapla PetscCheck(flg, comm, PETSC_ERR_PLIB, "root vectors differ"); 77fd697ae2SVaclav Hapla 789566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vecRoot0)); 799566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vecRoot1)); 809566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vecLeave0)); 819566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vecLeave1)); 823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 83fd697ae2SVaclav Hapla } 84fd697ae2SVaclav Hapla 85d71ae5a4SJacob Faibussowitsch PetscErrorCode CreateSF0(AppCtx *ctx, PetscSF *sf0) 86d71ae5a4SJacob Faibussowitsch { 876497c311SBarry Smith PetscInt j, k; 88fd697ae2SVaclav Hapla PetscInt nLeaves = ctx->nLeavesPerRank * ctx->size; 89fd697ae2SVaclav Hapla PetscInt nroots = ctx->nLeavesPerRank; 90fd697ae2SVaclav Hapla PetscSF sf; 91fd697ae2SVaclav Hapla PetscInt *ilocal; 92fd697ae2SVaclav Hapla PetscSFNode *iremote; 93fd697ae2SVaclav Hapla 94fd697ae2SVaclav Hapla PetscFunctionBegin; 959566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nLeaves + 1, &ctx->ilocal)); 969566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nLeaves, &ctx->iremote)); 97fd697ae2SVaclav Hapla ilocal = ctx->ilocal; 98fd697ae2SVaclav Hapla iremote = ctx->iremote; 99fd697ae2SVaclav Hapla ilocal[nLeaves] = -ctx->leaveStep; 1009566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(ctx->comm, &sf)); 1016497c311SBarry Smith j = nLeaves - 1; 1026497c311SBarry Smith for (PetscMPIInt r = 0; r < ctx->size; r++) { 103fd697ae2SVaclav Hapla for (k = 0; k < ctx->nLeavesPerRank; k++, j--) { 104fd697ae2SVaclav Hapla ilocal[j] = ilocal[j + 1] + ctx->leaveStep; 105fd697ae2SVaclav Hapla iremote[j].rank = r; 106fd697ae2SVaclav Hapla iremote[j].index = k; 107fd697ae2SVaclav Hapla } 108fd697ae2SVaclav Hapla } 1099566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(sf, nroots, nLeaves, ilocal, ctx->localmode, iremote, ctx->remotemode)); 110fd697ae2SVaclav Hapla { 111fd697ae2SVaclav Hapla const PetscInt *tlocal; 112fd697ae2SVaclav Hapla PetscBool sorted; 113fd697ae2SVaclav Hapla 1149566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, NULL, NULL, &tlocal, NULL)); 11508401ef6SPierre Jolivet PetscCheck(!ctx->contiguousLeaves || !tlocal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ilocal=NULL expected for contiguous case"); 116fd697ae2SVaclav Hapla if (tlocal) { 1179566063dSJacob Faibussowitsch PetscCall(PetscSortedInt(nLeaves, tlocal, &sorted)); 118fd697ae2SVaclav Hapla PetscCheck(sorted, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ilocal expected to be sorted"); 119fd697ae2SVaclav Hapla } 120fd697ae2SVaclav Hapla } 121fd697ae2SVaclav Hapla *sf0 = sf; 1223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 123fd697ae2SVaclav Hapla } 124fd697ae2SVaclav Hapla 125d71ae5a4SJacob Faibussowitsch PetscErrorCode CreateSF1(AppCtx *ctx, PetscSF *sf1) 126d71ae5a4SJacob Faibussowitsch { 1276497c311SBarry Smith PetscInt j, k; 128fd697ae2SVaclav Hapla PetscInt *ilocal = NULL; 129fd697ae2SVaclav Hapla PetscSFNode *iremote; 130fd697ae2SVaclav Hapla PetscInt nLeaves = ctx->nLeavesPerRank * ctx->size; 131fd697ae2SVaclav Hapla PetscInt nroots = ctx->nLeavesPerRank; 132fd697ae2SVaclav Hapla PetscSF sf; 133fd697ae2SVaclav Hapla 134fd697ae2SVaclav Hapla PetscFunctionBegin; 135fd697ae2SVaclav Hapla ilocal = NULL; 13648a46eb9SPierre Jolivet if (!ctx->contiguousLeaves) PetscCall(PetscCalloc1(nLeaves + 1, &ilocal)); 1379566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nLeaves, &iremote)); 1389566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(ctx->comm, &sf)); 1396497c311SBarry Smith j = 0; 1406497c311SBarry Smith for (PetscMPIInt r = 0; r < ctx->size; r++) { 141fd697ae2SVaclav Hapla for (k = 0; k < ctx->nLeavesPerRank; k++, j++) { 142ad540459SPierre Jolivet if (!ctx->contiguousLeaves) ilocal[j + 1] = ilocal[j] + ctx->leaveStep; 143fd697ae2SVaclav Hapla iremote[j].rank = r; 144fd697ae2SVaclav Hapla iremote[j].index = k; 145fd697ae2SVaclav Hapla } 146fd697ae2SVaclav Hapla } 147fd697ae2SVaclav Hapla PetscCheck(j == nLeaves, PETSC_COMM_SELF, PETSC_ERR_PLIB, "j != nLeaves"); 1489566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(sf, nroots, nLeaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER)); 149fd697ae2SVaclav Hapla if (ctx->contiguousLeaves) { 150fd697ae2SVaclav Hapla const PetscInt *tlocal; 151fd697ae2SVaclav Hapla 1529566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, NULL, NULL, &tlocal, NULL)); 153c9cc58a2SBarry Smith PetscCheck(!tlocal, PETSC_COMM_SELF, PETSC_ERR_PLIB, "ilocal=NULL expected for contiguous case"); 154fd697ae2SVaclav Hapla } 155fd697ae2SVaclav Hapla *sf1 = sf; 1563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 157fd697ae2SVaclav Hapla } 158fd697ae2SVaclav Hapla 159d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 160d71ae5a4SJacob Faibussowitsch { 161fd697ae2SVaclav Hapla AppCtx ctx; 162fd697ae2SVaclav Hapla PetscSF sf0, sf1; 163fd697ae2SVaclav Hapla MPI_Comm comm; 164fd697ae2SVaclav Hapla 165327415f7SBarry Smith PetscFunctionBeginUser; 1669566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 167fd697ae2SVaclav Hapla comm = PETSC_COMM_WORLD; 1689566063dSJacob Faibussowitsch PetscCall(GetOptions(comm, &ctx)); 169fd697ae2SVaclav Hapla 1709566063dSJacob Faibussowitsch PetscCall(CreateSF0(&ctx, &sf0)); 1719566063dSJacob Faibussowitsch PetscCall(CreateSF1(&ctx, &sf1)); 1729566063dSJacob Faibussowitsch PetscCall(PetscSFViewFromOptions(sf0, NULL, "-sf0_view")); 1739566063dSJacob Faibussowitsch PetscCall(PetscSFViewFromOptions(sf1, NULL, "-sf1_view")); 1749566063dSJacob Faibussowitsch PetscCall(PetscSFCheckEqual_Private(sf0, sf1)); 175fd697ae2SVaclav Hapla 1769566063dSJacob Faibussowitsch if (ctx.localmode != PETSC_OWN_POINTER) PetscCall(PetscFree(ctx.ilocal)); 1779566063dSJacob Faibussowitsch if (ctx.remotemode != PETSC_OWN_POINTER) PetscCall(PetscFree(ctx.iremote)); 1789566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sf0)); 1799566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sf1)); 1809566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 181b122ec5aSJacob Faibussowitsch return 0; 182fd697ae2SVaclav Hapla } 183fd697ae2SVaclav Hapla 184fd697ae2SVaclav Hapla /*TEST 185fd697ae2SVaclav Hapla testset: 186fd697ae2SVaclav Hapla suffix: 1 187fd697ae2SVaclav Hapla nsize: {{1 3}} 188fd697ae2SVaclav Hapla args: -n_leaves_per_rank {{0 5}} -leave_step {{1 3}} 189*3886731fSPierre Jolivet output_file: output/empty.out 190fd697ae2SVaclav Hapla test: 191fd697ae2SVaclav Hapla suffix: a 192fd697ae2SVaclav Hapla args: -localmode {{COPY_VALUES OWN_POINTER}} -remotemode {{COPY_VALUES OWN_POINTER}} 193fd697ae2SVaclav Hapla test: 194fd697ae2SVaclav Hapla suffix: b 195fd697ae2SVaclav Hapla args: -localmode USE_POINTER -remotemode {{COPY_VALUES OWN_POINTER USE_POINTER}} 196fd697ae2SVaclav Hapla TEST*/ 197