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