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