1*f210f596SVaclav Hapla 2*f210f596SVaclav Hapla static char help[]= "Test PetscSFConcatenate()\n\n"; 3*f210f596SVaclav Hapla 4*f210f596SVaclav Hapla #include <petscsf.h> 5*f210f596SVaclav Hapla 6*f210f596SVaclav Hapla typedef struct { 7*f210f596SVaclav Hapla MPI_Comm comm; 8*f210f596SVaclav Hapla PetscMPIInt rank, size; 9*f210f596SVaclav Hapla PetscInt leaveStep, nsfs, nLeavesPerRank; 10*f210f596SVaclav Hapla PetscBool shareRoots, sparseLeaves; 11*f210f596SVaclav Hapla } AppCtx; 12*f210f596SVaclav Hapla 13*f210f596SVaclav Hapla static PetscErrorCode GetOptions(MPI_Comm comm, AppCtx *ctx) 14*f210f596SVaclav Hapla { 15*f210f596SVaclav Hapla PetscErrorCode ierr; 16*f210f596SVaclav Hapla 17*f210f596SVaclav Hapla PetscFunctionBegin; 18*f210f596SVaclav Hapla ctx->comm = comm; 19*f210f596SVaclav Hapla ctx->nsfs = 3; 20*f210f596SVaclav Hapla ctx->nLeavesPerRank = 4; 21*f210f596SVaclav Hapla ctx->leaveStep = 1; 22*f210f596SVaclav Hapla ctx->shareRoots = PETSC_FALSE; 23*f210f596SVaclav Hapla ctx->sparseLeaves = PETSC_FALSE; 24*f210f596SVaclav Hapla ierr = PetscOptionsGetInt(NULL, NULL, "-nsfs", &ctx->nsfs, NULL);CHKERRQ(ierr); 25*f210f596SVaclav Hapla ierr = PetscOptionsGetInt(NULL, NULL, "-n_leaves_per_rank", &ctx->nLeavesPerRank, NULL);CHKERRQ(ierr); 26*f210f596SVaclav Hapla ierr = PetscOptionsGetInt(NULL, NULL, "-leave_step", &ctx->leaveStep, NULL);CHKERRQ(ierr); 27*f210f596SVaclav Hapla ierr = PetscOptionsGetBool(NULL, NULL, "-share_roots", &ctx->shareRoots, NULL);CHKERRQ(ierr); 28*f210f596SVaclav Hapla ctx->sparseLeaves = (PetscBool) (ctx->leaveStep != 1); 29*f210f596SVaclav Hapla ierr = MPI_Comm_size(comm, &ctx->size);CHKERRMPI(ierr); 30*f210f596SVaclav Hapla ierr = MPI_Comm_rank(comm, &ctx->rank);CHKERRMPI(ierr); 31*f210f596SVaclav Hapla PetscFunctionReturn(0); 32*f210f596SVaclav Hapla } 33*f210f596SVaclav Hapla 34*f210f596SVaclav Hapla static PetscErrorCode PetscSFCheckEqual_Private(PetscSF sf0, PetscSF sf1) 35*f210f596SVaclav Hapla { 36*f210f596SVaclav Hapla PetscInt nRoot, nLeave; 37*f210f596SVaclav Hapla Vec vecRoot0, vecLeave0, vecRoot1, vecLeave1; 38*f210f596SVaclav Hapla MPI_Comm comm; 39*f210f596SVaclav Hapla PetscBool flg; 40*f210f596SVaclav Hapla PetscErrorCode ierr; 41*f210f596SVaclav Hapla 42*f210f596SVaclav Hapla PetscFunctionBegin; 43*f210f596SVaclav Hapla ierr = PetscObjectGetComm((PetscObject)sf0, &comm);CHKERRQ(ierr); 44*f210f596SVaclav Hapla ierr = PetscSFGetGraph(sf0, &nRoot, NULL, NULL, NULL);CHKERRQ(ierr); 45*f210f596SVaclav Hapla ierr = PetscSFGetLeafRange(sf0, NULL, &nLeave);CHKERRQ(ierr); 46*f210f596SVaclav Hapla nLeave++; 47*f210f596SVaclav Hapla ierr = VecCreateMPI(comm, nRoot, PETSC_DECIDE, &vecRoot0);CHKERRQ(ierr); 48*f210f596SVaclav Hapla ierr = VecCreateMPI(comm, nLeave, PETSC_DECIDE, &vecLeave0);CHKERRQ(ierr); 49*f210f596SVaclav Hapla ierr = VecDuplicate(vecRoot0, &vecRoot1);CHKERRQ(ierr); 50*f210f596SVaclav Hapla ierr = VecDuplicate(vecLeave0, &vecLeave1);CHKERRQ(ierr); 51*f210f596SVaclav Hapla { 52*f210f596SVaclav Hapla PetscRandom rand; 53*f210f596SVaclav Hapla 54*f210f596SVaclav Hapla ierr = PetscRandomCreate(comm, &rand);CHKERRQ(ierr); 55*f210f596SVaclav Hapla ierr = PetscRandomSetFromOptions(rand);CHKERRQ(ierr); 56*f210f596SVaclav Hapla ierr = VecSetRandom(vecRoot0, rand);CHKERRQ(ierr); 57*f210f596SVaclav Hapla ierr = VecSetRandom(vecLeave0, rand);CHKERRQ(ierr); 58*f210f596SVaclav Hapla ierr = VecCopy(vecRoot0, vecRoot1);CHKERRQ(ierr); 59*f210f596SVaclav Hapla ierr = VecCopy(vecLeave0, vecLeave1);CHKERRQ(ierr); 60*f210f596SVaclav Hapla ierr = PetscRandomDestroy(&rand);CHKERRQ(ierr); 61*f210f596SVaclav Hapla } 62*f210f596SVaclav Hapla 63*f210f596SVaclav Hapla ierr = VecScatterBegin(sf0, vecRoot0, vecLeave0, ADD_VALUES, SCATTER_FORWARD);CHKERRQ(ierr); 64*f210f596SVaclav Hapla ierr = VecScatterEnd( sf0, vecRoot0, vecLeave0, ADD_VALUES, SCATTER_FORWARD);CHKERRQ(ierr); 65*f210f596SVaclav Hapla ierr = VecScatterBegin(sf1, vecRoot1, vecLeave1, ADD_VALUES, SCATTER_FORWARD);CHKERRQ(ierr); 66*f210f596SVaclav Hapla ierr = VecScatterEnd( sf1, vecRoot1, vecLeave1, ADD_VALUES, SCATTER_FORWARD);CHKERRQ(ierr); 67*f210f596SVaclav Hapla ierr = VecEqual(vecLeave0, vecLeave1, &flg);CHKERRQ(ierr); 68*f210f596SVaclav Hapla PetscCheck(flg, comm, PETSC_ERR_PLIB, "leave vectors differ"); 69*f210f596SVaclav Hapla 70*f210f596SVaclav Hapla ierr = VecScatterBegin(sf0, vecLeave0, vecRoot0, ADD_VALUES, SCATTER_REVERSE);CHKERRQ(ierr); 71*f210f596SVaclav Hapla ierr = VecScatterEnd( sf0, vecLeave0, vecRoot0, ADD_VALUES, SCATTER_REVERSE);CHKERRQ(ierr); 72*f210f596SVaclav Hapla ierr = VecScatterBegin(sf1, vecLeave1, vecRoot1, ADD_VALUES, SCATTER_REVERSE);CHKERRQ(ierr); 73*f210f596SVaclav Hapla ierr = VecScatterEnd( sf1, vecLeave1, vecRoot1, ADD_VALUES, SCATTER_REVERSE);CHKERRQ(ierr); 74*f210f596SVaclav Hapla ierr = VecEqual(vecRoot0, vecRoot1, &flg);CHKERRQ(ierr); 75*f210f596SVaclav Hapla PetscCheck(flg, comm, PETSC_ERR_PLIB, "root vectors differ"); 76*f210f596SVaclav Hapla 77*f210f596SVaclav Hapla ierr = VecDestroy(&vecRoot0);CHKERRQ(ierr); 78*f210f596SVaclav Hapla ierr = VecDestroy(&vecRoot1);CHKERRQ(ierr); 79*f210f596SVaclav Hapla ierr = VecDestroy(&vecLeave0);CHKERRQ(ierr); 80*f210f596SVaclav Hapla ierr = VecDestroy(&vecLeave1);CHKERRQ(ierr); 81*f210f596SVaclav Hapla PetscFunctionReturn(0); 82*f210f596SVaclav Hapla } 83*f210f596SVaclav Hapla 84*f210f596SVaclav Hapla PetscErrorCode CreateReferenceSF(AppCtx *ctx, PetscSF *refSF) 85*f210f596SVaclav Hapla { 86*f210f596SVaclav Hapla PetscInt i, j, k, r; 87*f210f596SVaclav Hapla PetscInt *ilocal = NULL; 88*f210f596SVaclav Hapla PetscSFNode *iremote; 89*f210f596SVaclav Hapla PetscInt nLeaves = ctx->nsfs * ctx->nLeavesPerRank * ctx->size; 90*f210f596SVaclav Hapla PetscInt nroots = ctx->nLeavesPerRank * ctx->nsfs; 91*f210f596SVaclav Hapla PetscSF sf; 92*f210f596SVaclav Hapla PetscErrorCode ierr; 93*f210f596SVaclav Hapla 94*f210f596SVaclav Hapla PetscFunctionBegin; 95*f210f596SVaclav Hapla ilocal = NULL; 96*f210f596SVaclav Hapla if (ctx->sparseLeaves) { 97*f210f596SVaclav Hapla ierr = PetscCalloc1(nLeaves+1, &ilocal);CHKERRQ(ierr); 98*f210f596SVaclav Hapla } 99*f210f596SVaclav Hapla ierr = PetscMalloc1(nLeaves, &iremote);CHKERRQ(ierr); 100*f210f596SVaclav Hapla ierr = PetscSFCreate(ctx->comm, &sf);CHKERRQ(ierr); 101*f210f596SVaclav Hapla for (i=0, j=0; i<ctx->nsfs; i++) { 102*f210f596SVaclav Hapla for (r=0; r<ctx->size; r++) { 103*f210f596SVaclav Hapla for (k=0; k<ctx->nLeavesPerRank; k++, j++) { 104*f210f596SVaclav Hapla if (ctx->sparseLeaves) { 105*f210f596SVaclav Hapla ilocal[j+1] = ilocal[j] + ctx->leaveStep; 106*f210f596SVaclav Hapla } 107*f210f596SVaclav Hapla iremote[j].rank = r; 108*f210f596SVaclav Hapla iremote[j].index = k + i * ctx->nLeavesPerRank; 109*f210f596SVaclav Hapla } 110*f210f596SVaclav Hapla } 111*f210f596SVaclav Hapla } 112*f210f596SVaclav Hapla ierr = PetscSFSetGraph(sf, nroots, nLeaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER);CHKERRQ(ierr); 113*f210f596SVaclav Hapla *refSF = sf; 114*f210f596SVaclav Hapla PetscFunctionReturn(0); 115*f210f596SVaclav Hapla } 116*f210f596SVaclav Hapla 117*f210f596SVaclav Hapla PetscErrorCode CreateSFs(AppCtx *ctx, PetscSF *newSFs[], PetscInt *leafOffsets[]) 118*f210f596SVaclav Hapla { 119*f210f596SVaclav Hapla PetscInt i; 120*f210f596SVaclav Hapla PetscInt *lOffsets = NULL; 121*f210f596SVaclav Hapla PetscSF *sfs; 122*f210f596SVaclav Hapla PetscInt nLeaves = ctx->nLeavesPerRank * ctx->size; 123*f210f596SVaclav Hapla PetscInt nroots = ctx->shareRoots ? ctx->nLeavesPerRank * ctx->nsfs : ctx->nLeavesPerRank; 124*f210f596SVaclav Hapla PetscErrorCode ierr; 125*f210f596SVaclav Hapla 126*f210f596SVaclav Hapla PetscFunctionBegin; 127*f210f596SVaclav Hapla if (ctx->sparseLeaves) { 128*f210f596SVaclav Hapla ierr = PetscCalloc1(ctx->nsfs+1, &lOffsets);CHKERRQ(ierr); 129*f210f596SVaclav Hapla } 130*f210f596SVaclav Hapla ierr = PetscMalloc1(ctx->nsfs, &sfs);CHKERRQ(ierr); 131*f210f596SVaclav Hapla for (i=0; i<ctx->nsfs; i++) { 132*f210f596SVaclav Hapla PetscInt j, k; 133*f210f596SVaclav Hapla PetscMPIInt r; 134*f210f596SVaclav Hapla PetscInt *ilocal = NULL; 135*f210f596SVaclav Hapla PetscSFNode *iremote; 136*f210f596SVaclav Hapla 137*f210f596SVaclav Hapla if (ctx->sparseLeaves) { 138*f210f596SVaclav Hapla ierr = PetscCalloc1(nLeaves+1, &ilocal);CHKERRQ(ierr); 139*f210f596SVaclav Hapla } 140*f210f596SVaclav Hapla ierr = PetscMalloc1(nLeaves, &iremote);CHKERRQ(ierr); 141*f210f596SVaclav Hapla for (r=0, j=0; r<ctx->size; r++) { 142*f210f596SVaclav Hapla for (k=0; k<ctx->nLeavesPerRank; k++, j++) { 143*f210f596SVaclav Hapla if (ctx->sparseLeaves) { 144*f210f596SVaclav Hapla ilocal[j+1] = ilocal[j] + ctx->leaveStep; 145*f210f596SVaclav Hapla } 146*f210f596SVaclav Hapla iremote[j].rank = r; 147*f210f596SVaclav Hapla iremote[j].index = ctx->shareRoots ? k + i * ctx->nLeavesPerRank : k; 148*f210f596SVaclav Hapla } 149*f210f596SVaclav Hapla } 150*f210f596SVaclav Hapla if (ctx->sparseLeaves) lOffsets[i+1] = lOffsets[i] + ilocal[j]; 151*f210f596SVaclav Hapla 152*f210f596SVaclav Hapla ierr = PetscSFCreate(ctx->comm, &sfs[i]);CHKERRQ(ierr); 153*f210f596SVaclav Hapla ierr = PetscSFSetGraph(sfs[i], nroots, nLeaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER);CHKERRQ(ierr); 154*f210f596SVaclav Hapla } 155*f210f596SVaclav Hapla *newSFs = sfs; 156*f210f596SVaclav Hapla *leafOffsets = lOffsets; 157*f210f596SVaclav Hapla PetscFunctionReturn(0); 158*f210f596SVaclav Hapla } 159*f210f596SVaclav Hapla 160*f210f596SVaclav Hapla PetscErrorCode DestroySFs(AppCtx *ctx, PetscSF *sfs[]) 161*f210f596SVaclav Hapla { 162*f210f596SVaclav Hapla PetscInt i; 163*f210f596SVaclav Hapla PetscErrorCode ierr; 164*f210f596SVaclav Hapla 165*f210f596SVaclav Hapla PetscFunctionBegin; 166*f210f596SVaclav Hapla for (i=0; i<ctx->nsfs; i++) { 167*f210f596SVaclav Hapla ierr = PetscSFDestroy(&(*sfs)[i]);CHKERRQ(ierr); 168*f210f596SVaclav Hapla } 169*f210f596SVaclav Hapla ierr = PetscFree(*sfs);CHKERRQ(ierr); 170*f210f596SVaclav Hapla PetscFunctionReturn(0); 171*f210f596SVaclav Hapla } 172*f210f596SVaclav Hapla 173*f210f596SVaclav Hapla int main(int argc, char **argv) 174*f210f596SVaclav Hapla { 175*f210f596SVaclav Hapla AppCtx ctx; 176*f210f596SVaclav Hapla PetscSF sf, sfRef; 177*f210f596SVaclav Hapla PetscSF *sfs; 178*f210f596SVaclav Hapla PetscInt *leafOffsets; 179*f210f596SVaclav Hapla MPI_Comm comm; 180*f210f596SVaclav Hapla PetscErrorCode ierr; 181*f210f596SVaclav Hapla 182*f210f596SVaclav Hapla ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 183*f210f596SVaclav Hapla comm = PETSC_COMM_WORLD; 184*f210f596SVaclav Hapla ierr = GetOptions(comm, &ctx);CHKERRQ(ierr); 185*f210f596SVaclav Hapla 186*f210f596SVaclav Hapla ierr = CreateSFs(&ctx, &sfs, &leafOffsets);CHKERRQ(ierr); 187*f210f596SVaclav Hapla ierr = PetscSFConcatenate(comm, ctx.nsfs, sfs, ctx.shareRoots, leafOffsets, &sf);CHKERRQ(ierr); 188*f210f596SVaclav Hapla ierr = CreateReferenceSF(&ctx, &sfRef);CHKERRQ(ierr); 189*f210f596SVaclav Hapla ierr = PetscSFCheckEqual_Private(sf, sfRef);CHKERRQ(ierr); 190*f210f596SVaclav Hapla 191*f210f596SVaclav Hapla ierr = DestroySFs(&ctx, &sfs);CHKERRQ(ierr); 192*f210f596SVaclav Hapla ierr = PetscFree(leafOffsets);CHKERRQ(ierr); 193*f210f596SVaclav Hapla ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 194*f210f596SVaclav Hapla ierr = PetscSFDestroy(&sfRef);CHKERRQ(ierr); 195*f210f596SVaclav Hapla ierr = PetscFinalize(); 196*f210f596SVaclav Hapla return ierr; 197*f210f596SVaclav Hapla } 198*f210f596SVaclav Hapla 199*f210f596SVaclav Hapla /*TEST 200*f210f596SVaclav Hapla test: 201*f210f596SVaclav Hapla nsize: {{1 3}} 202*f210f596SVaclav Hapla args: -nsfs {{1 3}} -n_leaves_per_rank {{0 1 5}} -leave_step {{1 3}} -share_roots {{true false}} 203*f210f596SVaclav Hapla TEST*/ 204