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 16fd697ae2SVaclav Hapla static PetscErrorCode GetOptions(MPI_Comm comm, AppCtx *ctx) 17fd697ae2SVaclav Hapla { 18fd697ae2SVaclav Hapla PetscFunctionBegin; 19fd697ae2SVaclav Hapla ctx->comm = comm; 20fd697ae2SVaclav Hapla ctx->nLeavesPerRank = 4; 21fd697ae2SVaclav Hapla ctx->leaveStep = 1; 22fd697ae2SVaclav Hapla ctx->contiguousLeaves = PETSC_FALSE; 23fd697ae2SVaclav Hapla ctx->localmode = PETSC_OWN_POINTER; 24fd697ae2SVaclav Hapla ctx->remotemode = PETSC_OWN_POINTER; 25fd697ae2SVaclav Hapla ctx->ilocal = NULL; 26fd697ae2SVaclav Hapla ctx->iremote = NULL; 275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL, NULL, "-n_leaves_per_rank", &ctx->nLeavesPerRank, NULL)); 285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL, NULL, "-leave_step", &ctx->leaveStep, NULL)); 295f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetEnum(NULL, NULL, "-localmode", PetscCopyModes, (PetscEnum*) &ctx->localmode, NULL)); 305f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetEnum(NULL, NULL, "-remotemode", PetscCopyModes, (PetscEnum*) &ctx->remotemode, NULL)); 31fd697ae2SVaclav Hapla ctx->contiguousLeaves = (PetscBool) (ctx->leaveStep == 1); 325f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(comm, &ctx->size)); 335f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(comm, &ctx->rank)); 34fd697ae2SVaclav Hapla PetscFunctionReturn(0); 35fd697ae2SVaclav Hapla } 36fd697ae2SVaclav Hapla 37fd697ae2SVaclav Hapla static PetscErrorCode PetscSFCheckEqual_Private(PetscSF sf0, PetscSF sf1) 38fd697ae2SVaclav Hapla { 39fd697ae2SVaclav Hapla PetscInt nRoot, nLeave; 40fd697ae2SVaclav Hapla Vec vecRoot0, vecLeave0, vecRoot1, vecLeave1; 41fd697ae2SVaclav Hapla MPI_Comm comm; 42fd697ae2SVaclav Hapla PetscBool flg; 43fd697ae2SVaclav Hapla 44fd697ae2SVaclav Hapla PetscFunctionBegin; 455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetComm((PetscObject)sf0, &comm)); 465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFGetGraph(sf0, &nRoot, NULL, NULL, NULL)); 475f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFGetLeafRange(sf0, NULL, &nLeave)); 48fd697ae2SVaclav Hapla nLeave++; 495f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateMPI(comm, nRoot, PETSC_DECIDE, &vecRoot0)); 505f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateMPI(comm, nLeave, PETSC_DECIDE, &vecLeave0)); 515f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(vecRoot0, &vecRoot1)); 525f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(vecLeave0, &vecLeave1)); 53fd697ae2SVaclav Hapla { 54fd697ae2SVaclav Hapla PetscRandom rand; 55fd697ae2SVaclav Hapla 565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(comm, &rand)); 575f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetFromOptions(rand)); 585f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(vecRoot0, rand)); 595f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(vecLeave0, rand)); 605f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(vecRoot0, vecRoot1)); 615f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(vecLeave0, vecLeave1)); 625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rand)); 63fd697ae2SVaclav Hapla } 64fd697ae2SVaclav Hapla 655f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(sf0, vecRoot0, vecLeave0, ADD_VALUES, SCATTER_FORWARD)); 665f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd( sf0, vecRoot0, vecLeave0, ADD_VALUES, SCATTER_FORWARD)); 675f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(sf1, vecRoot1, vecLeave1, ADD_VALUES, SCATTER_FORWARD)); 685f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd( sf1, vecRoot1, vecLeave1, ADD_VALUES, SCATTER_FORWARD)); 695f80ce2aSJacob Faibussowitsch CHKERRQ(VecEqual(vecLeave0, vecLeave1, &flg)); 70fd697ae2SVaclav Hapla PetscCheck(flg, comm, PETSC_ERR_PLIB, "leave vectors differ"); 71fd697ae2SVaclav Hapla 725f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(sf0, vecLeave0, vecRoot0, ADD_VALUES, SCATTER_REVERSE)); 735f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd( sf0, vecLeave0, vecRoot0, ADD_VALUES, SCATTER_REVERSE)); 745f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(sf1, vecLeave1, vecRoot1, ADD_VALUES, SCATTER_REVERSE)); 755f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd( sf1, vecLeave1, vecRoot1, ADD_VALUES, SCATTER_REVERSE)); 765f80ce2aSJacob Faibussowitsch CHKERRQ(VecEqual(vecRoot0, vecRoot1, &flg)); 77fd697ae2SVaclav Hapla PetscCheck(flg, comm, PETSC_ERR_PLIB, "root vectors differ"); 78fd697ae2SVaclav Hapla 795f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&vecRoot0)); 805f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&vecRoot1)); 815f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&vecLeave0)); 825f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&vecLeave1)); 83fd697ae2SVaclav Hapla PetscFunctionReturn(0); 84fd697ae2SVaclav Hapla } 85fd697ae2SVaclav Hapla 86fd697ae2SVaclav Hapla PetscErrorCode CreateSF0(AppCtx *ctx, PetscSF *sf0) 87fd697ae2SVaclav Hapla { 88fd697ae2SVaclav Hapla PetscInt j, k, r; 89fd697ae2SVaclav Hapla PetscInt nLeaves = ctx->nLeavesPerRank * ctx->size; 90fd697ae2SVaclav Hapla PetscInt nroots = ctx->nLeavesPerRank; 91fd697ae2SVaclav Hapla PetscSF sf; 92fd697ae2SVaclav Hapla PetscInt *ilocal; 93fd697ae2SVaclav Hapla PetscSFNode *iremote; 94fd697ae2SVaclav Hapla 95fd697ae2SVaclav Hapla PetscFunctionBegin; 965f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nLeaves+1, &ctx->ilocal)); 975f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nLeaves, &ctx->iremote)); 98fd697ae2SVaclav Hapla ilocal = ctx->ilocal; 99fd697ae2SVaclav Hapla iremote = ctx->iremote; 100fd697ae2SVaclav Hapla ilocal[nLeaves] = -ctx->leaveStep; 1015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFCreate(ctx->comm, &sf)); 102fd697ae2SVaclav Hapla for (r=0, j=nLeaves-1; 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 } 1095f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFSetGraph(sf, nroots, nLeaves, ilocal, ctx->localmode, iremote, ctx->remotemode)); 110fd697ae2SVaclav Hapla { 111fd697ae2SVaclav Hapla const PetscInt *tlocal; 112fd697ae2SVaclav Hapla PetscBool sorted; 113fd697ae2SVaclav Hapla 1145f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFGetGraph(sf, NULL, NULL, &tlocal, NULL)); 115fd697ae2SVaclav Hapla PetscCheckFalse(ctx->contiguousLeaves && tlocal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"ilocal=NULL expected for contiguous case"); 116fd697ae2SVaclav Hapla if (tlocal) { 1175f80ce2aSJacob Faibussowitsch CHKERRQ(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; 122fd697ae2SVaclav Hapla PetscFunctionReturn(0); 123fd697ae2SVaclav Hapla } 124fd697ae2SVaclav Hapla 125fd697ae2SVaclav Hapla PetscErrorCode CreateSF1(AppCtx *ctx, PetscSF *sf1) 126fd697ae2SVaclav Hapla { 127fd697ae2SVaclav Hapla PetscInt j, k, r; 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; 136fd697ae2SVaclav Hapla if (!ctx->contiguousLeaves) { 1375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscCalloc1(nLeaves+1, &ilocal)); 138fd697ae2SVaclav Hapla } 1395f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nLeaves, &iremote)); 1405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFCreate(ctx->comm, &sf)); 141fd697ae2SVaclav Hapla for (r=0, j=0; r<ctx->size; r++) { 142fd697ae2SVaclav Hapla for (k=0; k<ctx->nLeavesPerRank; k++, j++) { 143fd697ae2SVaclav Hapla if (!ctx->contiguousLeaves) { 144fd697ae2SVaclav Hapla ilocal[j+1] = ilocal[j] + ctx->leaveStep; 145fd697ae2SVaclav Hapla } 146fd697ae2SVaclav Hapla iremote[j].rank = r; 147fd697ae2SVaclav Hapla iremote[j].index = k; 148fd697ae2SVaclav Hapla } 149fd697ae2SVaclav Hapla } 150fd697ae2SVaclav Hapla PetscCheck(j == nLeaves,PETSC_COMM_SELF,PETSC_ERR_PLIB,"j != nLeaves"); 1515f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFSetGraph(sf, nroots, nLeaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER)); 152fd697ae2SVaclav Hapla if (ctx->contiguousLeaves) { 153fd697ae2SVaclav Hapla const PetscInt *tlocal; 154fd697ae2SVaclav Hapla 1555f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFGetGraph(sf, NULL, NULL, &tlocal, NULL)); 156fd697ae2SVaclav Hapla PetscCheckFalse(tlocal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"ilocal=NULL expected for contiguous case"); 157fd697ae2SVaclav Hapla } 158fd697ae2SVaclav Hapla *sf1 = sf; 159fd697ae2SVaclav Hapla PetscFunctionReturn(0); 160fd697ae2SVaclav Hapla } 161fd697ae2SVaclav Hapla 162fd697ae2SVaclav Hapla int main(int argc, char **argv) 163fd697ae2SVaclav Hapla { 164fd697ae2SVaclav Hapla AppCtx ctx; 165fd697ae2SVaclav Hapla PetscSF sf0, sf1; 166fd697ae2SVaclav Hapla MPI_Comm comm; 167fd697ae2SVaclav Hapla 168*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,NULL,help)); 169fd697ae2SVaclav Hapla comm = PETSC_COMM_WORLD; 1705f80ce2aSJacob Faibussowitsch CHKERRQ(GetOptions(comm, &ctx)); 171fd697ae2SVaclav Hapla 1725f80ce2aSJacob Faibussowitsch CHKERRQ(CreateSF0(&ctx, &sf0)); 1735f80ce2aSJacob Faibussowitsch CHKERRQ(CreateSF1(&ctx, &sf1)); 174*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscSFViewFromOptions(sf0, NULL, "-sf0_view")); 175*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscSFViewFromOptions(sf1, NULL, "-sf1_view")); 1765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFCheckEqual_Private(sf0, sf1)); 177fd697ae2SVaclav Hapla 178*b122ec5aSJacob Faibussowitsch if (ctx.localmode != PETSC_OWN_POINTER) CHKERRQ(PetscFree(ctx.ilocal)); 179*b122ec5aSJacob Faibussowitsch if (ctx.remotemode != PETSC_OWN_POINTER) CHKERRQ(PetscFree(ctx.iremote)); 1805f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFDestroy(&sf0)); 1815f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSFDestroy(&sf1)); 182*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 183*b122ec5aSJacob Faibussowitsch return 0; 184fd697ae2SVaclav Hapla } 185fd697ae2SVaclav Hapla 186fd697ae2SVaclav Hapla /*TEST 187fd697ae2SVaclav Hapla testset: 188fd697ae2SVaclav Hapla suffix: 1 189fd697ae2SVaclav Hapla nsize: {{1 3}} 190fd697ae2SVaclav Hapla args: -n_leaves_per_rank {{0 5}} -leave_step {{1 3}} 191fd697ae2SVaclav Hapla test: 192fd697ae2SVaclav Hapla suffix: a 193fd697ae2SVaclav Hapla args: -localmode {{COPY_VALUES OWN_POINTER}} -remotemode {{COPY_VALUES OWN_POINTER}} 194fd697ae2SVaclav Hapla test: 195fd697ae2SVaclav Hapla suffix: b 196fd697ae2SVaclav Hapla args: -localmode USE_POINTER -remotemode {{COPY_VALUES OWN_POINTER USE_POINTER}} 197fd697ae2SVaclav Hapla TEST*/ 198