105264a50SDave May #include <petscsf.h> 2af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/ 357459e9aSMatthew G Knepley 4ed4e918cSMatthew G Knepley /*@ 512b4a537SBarry Smith DMDAConvertToCell - Convert a (i,j,k) location in a `DMDA` to its local cell or vertex number 6341c9bdaSMatthew G Knepley 7ed4e918cSMatthew G Knepley Not Collective 8ed4e918cSMatthew G Knepley 9d8d19677SJose E. Roman Input Parameters: 1060225df5SJacob Faibussowitsch + dm - the distributed array 1112b4a537SBarry Smith - s - A `MatStencil` that provides (i,j,k) 12ed4e918cSMatthew G Knepley 13ed4e918cSMatthew G Knepley Output Parameter: 1412b4a537SBarry Smith . cell - the local cell or vertext number 15ed4e918cSMatthew G Knepley 16ed4e918cSMatthew G Knepley Level: developer 17dce8aebaSBarry Smith 1812b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA` 19ed4e918cSMatthew G Knepley @*/ 20d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAConvertToCell(DM dm, MatStencil s, PetscInt *cell) 21d71ae5a4SJacob Faibussowitsch { 22341c9bdaSMatthew G Knepley DM_DA *da = (DM_DA *)dm->data; 23c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 244e846693SMatthew G Knepley const PetscInt mx = (da->Xe - da->Xs) / da->w, my = da->Ye - da->Ys /*, mz = da->Ze - da->Zs*/; 25ed4e918cSMatthew G Knepley const PetscInt il = s.i - da->Xs / da->w, jl = dim > 1 ? s.j - da->Ys : 0, kl = dim > 2 ? s.k - da->Zs : 0; 26341c9bdaSMatthew G Knepley 27341c9bdaSMatthew G Knepley PetscFunctionBegin; 28ed4e918cSMatthew G Knepley *cell = -1; 2963a3b9bcSJacob Faibussowitsch PetscCheck(!(s.i < da->Xs / da->w) && !(s.i >= da->Xe / da->w), PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Stencil i %" PetscInt_FMT " should be in [%" PetscInt_FMT ", %" PetscInt_FMT ")", s.i, da->Xs / da->w, da->Xe / da->w); 301dca8a05SBarry Smith PetscCheck(dim <= 1 || (s.j >= da->Ys && s.j < da->Ye), PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Stencil j %" PetscInt_FMT " should be in [%" PetscInt_FMT ", %" PetscInt_FMT ")", s.j, da->Ys, da->Ye); 311dca8a05SBarry Smith PetscCheck(dim <= 2 || (s.k >= da->Zs && s.k < da->Ze), PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Stencil k %" PetscInt_FMT " should be in [%" PetscInt_FMT ", %" PetscInt_FMT ")", s.k, da->Zs, da->Ze); 32ed4e918cSMatthew G Knepley *cell = (kl * my + jl) * mx + il; 333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 34341c9bdaSMatthew G Knepley } 35341c9bdaSMatthew G Knepley 36*6c6a6b79SMatthew G. Knepley PetscErrorCode DMGetLocalBoundingBox_DA(DM da, PetscReal lmin[], PetscReal lmax[], PetscInt cs[], PetscInt ce[]) 37d71ae5a4SJacob Faibussowitsch { 38*6c6a6b79SMatthew G. Knepley PetscInt xs, xe, Xs, Xe; 39*6c6a6b79SMatthew G. Knepley PetscInt ys, ye, Ys, Ye; 40*6c6a6b79SMatthew G. Knepley PetscInt zs, ze, Zs, Ze; 41*6c6a6b79SMatthew G. Knepley PetscInt dim, M, N, P, c0, c1; 42*6c6a6b79SMatthew G. Knepley PetscReal gmax[3] = {0., 0., 0.}; 43*6c6a6b79SMatthew G. Knepley const PetscReal *L, *Lstart; 44*6c6a6b79SMatthew G. Knepley Vec coordinates; 45*6c6a6b79SMatthew G. Knepley const PetscScalar *coor; 46*6c6a6b79SMatthew G. Knepley DMBoundaryType bx, by, bz; 4705264a50SDave May 48aab5bcd8SJed Brown PetscFunctionBegin; 49*6c6a6b79SMatthew G. Knepley PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xe, &ye, &ze)); 50*6c6a6b79SMatthew G. Knepley PetscCall(DMDAGetGhostCorners(da, &Xs, &Ys, &Zs, &Xe, &Ye, &Ze)); 51*6c6a6b79SMatthew G. Knepley PetscCall(DMDAGetInfo(da, &dim, &M, &N, &P, NULL, NULL, NULL, NULL, NULL, &bx, &by, &bz, NULL)); 52*6c6a6b79SMatthew G. Knepley // Convert from widths to endpoints 539371c9d4SSatish Balay xe += xs; 549371c9d4SSatish Balay Xe += Xs; 559371c9d4SSatish Balay ye += ys; 569371c9d4SSatish Balay Ye += Ys; 57*6c6a6b79SMatthew G. Knepley ze += zs; 58*6c6a6b79SMatthew G. Knepley Ze += Zs; 59*6c6a6b79SMatthew G. Knepley // What is this doing? 6051c42ae6SMatthew G. Knepley if (xs != Xs && Xs >= 0) xs -= 1; 6151c42ae6SMatthew G. Knepley if (ys != Ys && Ys >= 0) ys -= 1; 62*6c6a6b79SMatthew G. Knepley if (zs != Zs && Zs >= 0) zs -= 1; 6305264a50SDave May 64*6c6a6b79SMatthew G. Knepley PetscCall(DMGetCoordinatesLocal(da, &coordinates)); 65*6c6a6b79SMatthew G. Knepley if (!coordinates) { 66*6c6a6b79SMatthew G. Knepley PetscCall(DMGetLocalBoundingIndices_DMDA(da, lmin, lmax)); 67*6c6a6b79SMatthew G. Knepley for (PetscInt d = 0; d < dim; ++d) { 68*6c6a6b79SMatthew G. Knepley if (cs) cs[d] = lmin[d]; 69*6c6a6b79SMatthew G. Knepley if (ce) ce[d] = lmax[d]; 70*6c6a6b79SMatthew G. Knepley } 71*6c6a6b79SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 72*6c6a6b79SMatthew G. Knepley } 73*6c6a6b79SMatthew G. Knepley PetscCall(VecGetArrayRead(coordinates, &coor)); 74*6c6a6b79SMatthew G. Knepley switch (dim) { 75*6c6a6b79SMatthew G. Knepley case 1: 76*6c6a6b79SMatthew G. Knepley c0 = (xs - Xs); 77*6c6a6b79SMatthew G. Knepley c1 = (xe - 2 - Xs + 1); 78*6c6a6b79SMatthew G. Knepley break; 79*6c6a6b79SMatthew G. Knepley case 2: 8005264a50SDave May c0 = (xs - Xs) + (ys - Ys) * (Xe - Xs); 8105264a50SDave May c1 = (xe - 2 - Xs + 1) + (ye - 2 - Ys + 1) * (Xe - Xs); 82*6c6a6b79SMatthew G. Knepley break; 83*6c6a6b79SMatthew G. Knepley case 3: 84*6c6a6b79SMatthew G. Knepley c0 = (xs - Xs) + (ys - Ys) * (Xe - Xs) + (zs - Zs) * (Xe - Xs) * (Ye - Ys); 85*6c6a6b79SMatthew G. Knepley c1 = (xe - 2 - Xs + 1) + (ye - 2 - Ys + 1) * (Xe - Xs) + (ze - 2 - Zs + 1) * (Xe - Xs) * (Ye - Ys); 86*6c6a6b79SMatthew G. Knepley break; 87*6c6a6b79SMatthew G. Knepley default: 88*6c6a6b79SMatthew G. Knepley SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONG, "Invalid dimension %" PetscInt_FMT " for DMDA", dim); 89*6c6a6b79SMatthew G. Knepley } 90*6c6a6b79SMatthew G. Knepley for (PetscInt d = 0; d < dim; ++d) { 91*6c6a6b79SMatthew G. Knepley lmin[d] = PetscRealPart(coor[c0 * dim + d]); 92*6c6a6b79SMatthew G. Knepley lmax[d] = PetscRealPart(coor[c1 * dim + d]); 93*6c6a6b79SMatthew G. Knepley } 94*6c6a6b79SMatthew G. Knepley PetscCall(VecRestoreArrayRead(coordinates, &coor)); 9505264a50SDave May 96*6c6a6b79SMatthew G. Knepley PetscCall(DMGetPeriodicity(da, NULL, &Lstart, &L)); 97*6c6a6b79SMatthew G. Knepley if (L) { 98*6c6a6b79SMatthew G. Knepley for (PetscInt d = 0; d < dim; ++d) 99*6c6a6b79SMatthew G. Knepley if (L[d] > 0.0) gmax[d] = Lstart[d] + L[d]; 100*6c6a6b79SMatthew G. Knepley } 101*6c6a6b79SMatthew G. Knepley // Must check for periodic boundary 102*6c6a6b79SMatthew G. Knepley if (bx == DM_BOUNDARY_PERIODIC && xe == M) { 103*6c6a6b79SMatthew G. Knepley lmax[0] = gmax[0]; 104*6c6a6b79SMatthew G. Knepley ++xe; 105*6c6a6b79SMatthew G. Knepley } 106*6c6a6b79SMatthew G. Knepley if (by == DM_BOUNDARY_PERIODIC && ye == N) { 107*6c6a6b79SMatthew G. Knepley lmax[1] = gmax[1]; 108*6c6a6b79SMatthew G. Knepley ++ye; 109*6c6a6b79SMatthew G. Knepley } 110*6c6a6b79SMatthew G. Knepley if (bz == DM_BOUNDARY_PERIODIC && ze == P) { 111*6c6a6b79SMatthew G. Knepley lmax[2] = gmax[2]; 112*6c6a6b79SMatthew G. Knepley ++ze; 113*6c6a6b79SMatthew G. Knepley } 114*6c6a6b79SMatthew G. Knepley if (cs) { 115*6c6a6b79SMatthew G. Knepley cs[0] = xs; 116*6c6a6b79SMatthew G. Knepley if (dim > 1) cs[1] = ys; 117*6c6a6b79SMatthew G. Knepley if (dim > 2) cs[2] = zs; 118*6c6a6b79SMatthew G. Knepley } 119*6c6a6b79SMatthew G. Knepley if (ce) { 120*6c6a6b79SMatthew G. Knepley ce[0] = xe; 121*6c6a6b79SMatthew G. Knepley if (dim > 1) ce[1] = ye; 122*6c6a6b79SMatthew G. Knepley if (dim > 2) ce[2] = ze; 123*6c6a6b79SMatthew G. Knepley } 124*6c6a6b79SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 125*6c6a6b79SMatthew G. Knepley } 12605264a50SDave May 127*6c6a6b79SMatthew G. Knepley static PetscErrorCode private_DMDALocatePointsIS_2D_Regular(DM dmregular, Vec pos, IS *iscell) 128*6c6a6b79SMatthew G. Knepley { 129*6c6a6b79SMatthew G. Knepley PetscInt n, bs, npoints; 130*6c6a6b79SMatthew G. Knepley PetscInt cs[2], ce[2]; 131*6c6a6b79SMatthew G. Knepley PetscInt xs, xe, mxlocal; 132*6c6a6b79SMatthew G. Knepley PetscInt ys, ye, mylocal; 133*6c6a6b79SMatthew G. Knepley PetscReal gmin_l[2], gmax_l[2], dx[2]; 134*6c6a6b79SMatthew G. Knepley PetscReal gmin[2], gmax[2]; 135*6c6a6b79SMatthew G. Knepley PetscInt *cellidx; 136*6c6a6b79SMatthew G. Knepley const PetscScalar *coor; 137*6c6a6b79SMatthew G. Knepley 138*6c6a6b79SMatthew G. Knepley PetscFunctionBegin; 139*6c6a6b79SMatthew G. Knepley PetscCall(DMGetLocalBoundingBox_DA(dmregular, gmin_l, gmax_l, cs, ce)); 140*6c6a6b79SMatthew G. Knepley xs = cs[0]; 141*6c6a6b79SMatthew G. Knepley ys = cs[1]; 142*6c6a6b79SMatthew G. Knepley xe = ce[0]; 143*6c6a6b79SMatthew G. Knepley ye = ce[1]; 144*6c6a6b79SMatthew G. Knepley PetscCall(DMGetBoundingBox(dmregular, gmin, gmax)); 14551c42ae6SMatthew G. Knepley 14651c42ae6SMatthew G. Knepley mxlocal = xe - xs - 1; 14751c42ae6SMatthew G. Knepley mylocal = ye - ys - 1; 14805264a50SDave May 14905264a50SDave May dx[0] = (gmax_l[0] - gmin_l[0]) / ((PetscReal)mxlocal); 15005264a50SDave May dx[1] = (gmax_l[1] - gmin_l[1]) / ((PetscReal)mylocal); 15105264a50SDave May 1529566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(pos, &n)); 1539566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(pos, &bs)); 15405264a50SDave May npoints = n / bs; 15505264a50SDave May 1569566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(npoints, &cellidx)); 157*6c6a6b79SMatthew G. Knepley PetscCall(VecGetArrayRead(pos, &coor)); 158*6c6a6b79SMatthew G. Knepley for (PetscInt p = 0; p < npoints; p++) { 1590ca99263SDave May PetscReal coor_p[2]; 16005264a50SDave May PetscInt mi[2]; 16105264a50SDave May 162*6c6a6b79SMatthew G. Knepley coor_p[0] = PetscRealPart(coor[2 * p]); 163*6c6a6b79SMatthew G. Knepley coor_p[1] = PetscRealPart(coor[2 * p + 1]); 16405264a50SDave May 16505264a50SDave May cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND; 16605264a50SDave May 167ad540459SPierre Jolivet if (coor_p[0] < gmin_l[0]) continue; 168ad540459SPierre Jolivet if (coor_p[0] > gmax_l[0]) continue; 169ad540459SPierre Jolivet if (coor_p[1] < gmin_l[1]) continue; 170ad540459SPierre Jolivet if (coor_p[1] > gmax_l[1]) continue; 17105264a50SDave May 172*6c6a6b79SMatthew G. Knepley for (PetscInt d = 0; d < 2; d++) mi[d] = (PetscInt)((coor_p[d] - gmin[d]) / dx[d]); 17305264a50SDave May 174ad540459SPierre Jolivet if (mi[0] < xs) continue; 175ad540459SPierre Jolivet if (mi[0] > (xe - 1)) continue; 176ad540459SPierre Jolivet if (mi[1] < ys) continue; 177ad540459SPierre Jolivet if (mi[1] > (ye - 1)) continue; 17805264a50SDave May 179ad540459SPierre Jolivet if (mi[0] == (xe - 1)) mi[0]--; 180ad540459SPierre Jolivet if (mi[1] == (ye - 1)) mi[1]--; 18105264a50SDave May 18205264a50SDave May cellidx[p] = (mi[0] - xs) + (mi[1] - ys) * mxlocal; 18305264a50SDave May } 184*6c6a6b79SMatthew G. Knepley PetscCall(VecRestoreArrayRead(pos, &coor)); 1859566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, npoints, cellidx, PETSC_OWN_POINTER, iscell)); 1863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 18705264a50SDave May } 18805264a50SDave May 18966976f2fSJacob Faibussowitsch static PetscErrorCode private_DMDALocatePointsIS_3D_Regular(DM dmregular, Vec pos, IS *iscell) 190d71ae5a4SJacob Faibussowitsch { 191*6c6a6b79SMatthew G. Knepley PetscInt n, bs, npoints; 192*6c6a6b79SMatthew G. Knepley PetscInt cs[3], ce[3]; 193*6c6a6b79SMatthew G. Knepley PetscInt xs, xe, mxlocal; 194*6c6a6b79SMatthew G. Knepley PetscInt ys, ye, mylocal; 195*6c6a6b79SMatthew G. Knepley PetscInt zs, ze, mzlocal; 19605264a50SDave May PetscReal gmin_l[3], gmax_l[3], dx[3]; 19705264a50SDave May PetscReal gmin[3], gmax[3]; 19805264a50SDave May PetscInt *cellidx; 199*6c6a6b79SMatthew G. Knepley const PetscScalar *coor; 20005264a50SDave May 201aab5bcd8SJed Brown PetscFunctionBegin; 202*6c6a6b79SMatthew G. Knepley PetscCall(DMGetLocalBoundingBox_DA(dmregular, gmin_l, gmax_l, cs, ce)); 203*6c6a6b79SMatthew G. Knepley xs = cs[0]; 204*6c6a6b79SMatthew G. Knepley ys = cs[1]; 205*6c6a6b79SMatthew G. Knepley zs = cs[2]; 206*6c6a6b79SMatthew G. Knepley xe = ce[0]; 207*6c6a6b79SMatthew G. Knepley ye = ce[1]; 208*6c6a6b79SMatthew G. Knepley ze = ce[2]; 209*6c6a6b79SMatthew G. Knepley PetscCall(DMGetBoundingBox(dmregular, gmin, gmax)); 21051c42ae6SMatthew G. Knepley 21151c42ae6SMatthew G. Knepley mxlocal = xe - xs - 1; 21251c42ae6SMatthew G. Knepley mylocal = ye - ys - 1; 21351c42ae6SMatthew G. Knepley mzlocal = ze - zs - 1; 21405264a50SDave May 21505264a50SDave May dx[0] = (gmax_l[0] - gmin_l[0]) / ((PetscReal)mxlocal); 21605264a50SDave May dx[1] = (gmax_l[1] - gmin_l[1]) / ((PetscReal)mylocal); 21705264a50SDave May dx[2] = (gmax_l[2] - gmin_l[2]) / ((PetscReal)mzlocal); 21805264a50SDave May 2199566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(pos, &n)); 2209566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(pos, &bs)); 22105264a50SDave May npoints = n / bs; 22205264a50SDave May 2239566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(npoints, &cellidx)); 224*6c6a6b79SMatthew G. Knepley PetscCall(VecGetArrayRead(pos, &coor)); 225*6c6a6b79SMatthew G. Knepley for (PetscInt p = 0; p < npoints; p++) { 2260ca99263SDave May PetscReal coor_p[3]; 22705264a50SDave May PetscInt mi[3]; 22805264a50SDave May 229*6c6a6b79SMatthew G. Knepley coor_p[0] = PetscRealPart(coor[3 * p]); 230*6c6a6b79SMatthew G. Knepley coor_p[1] = PetscRealPart(coor[3 * p + 1]); 231*6c6a6b79SMatthew G. Knepley coor_p[2] = PetscRealPart(coor[3 * p + 2]); 23205264a50SDave May 23305264a50SDave May cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND; 23405264a50SDave May 235ad540459SPierre Jolivet if (coor_p[0] < gmin_l[0]) continue; 236ad540459SPierre Jolivet if (coor_p[0] > gmax_l[0]) continue; 237ad540459SPierre Jolivet if (coor_p[1] < gmin_l[1]) continue; 238ad540459SPierre Jolivet if (coor_p[1] > gmax_l[1]) continue; 239ad540459SPierre Jolivet if (coor_p[2] < gmin_l[2]) continue; 240ad540459SPierre Jolivet if (coor_p[2] > gmax_l[2]) continue; 24105264a50SDave May 242*6c6a6b79SMatthew G. Knepley for (PetscInt d = 0; d < 3; d++) mi[d] = (PetscInt)((coor_p[d] - gmin[d]) / dx[d]); 24305264a50SDave May 244*6c6a6b79SMatthew G. Knepley // TODO: Check for periodicity here 245ad540459SPierre Jolivet if (mi[0] < xs) continue; 246ad540459SPierre Jolivet if (mi[0] > (xe - 1)) continue; 247ad540459SPierre Jolivet if (mi[1] < ys) continue; 248ad540459SPierre Jolivet if (mi[1] > (ye - 1)) continue; 249ad540459SPierre Jolivet if (mi[2] < zs) continue; 250ad540459SPierre Jolivet if (mi[2] > (ze - 1)) continue; 25105264a50SDave May 252ad540459SPierre Jolivet if (mi[0] == (xe - 1)) mi[0]--; 253ad540459SPierre Jolivet if (mi[1] == (ye - 1)) mi[1]--; 254ad540459SPierre Jolivet if (mi[2] == (ze - 1)) mi[2]--; 25505264a50SDave May 25605264a50SDave May cellidx[p] = (mi[0] - xs) + (mi[1] - ys) * mxlocal + (mi[2] - zs) * mxlocal * mylocal; 25705264a50SDave May } 258*6c6a6b79SMatthew G. Knepley PetscCall(VecRestoreArrayRead(pos, &coor)); 2599566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, npoints, cellidx, PETSC_OWN_POINTER, iscell)); 2603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 26105264a50SDave May } 26205264a50SDave May 263d71ae5a4SJacob Faibussowitsch PetscErrorCode DMLocatePoints_DA_Regular(DM dm, Vec pos, DMPointLocationType ltype, PetscSF cellSF) 264d71ae5a4SJacob Faibussowitsch { 26505264a50SDave May IS iscell; 26605264a50SDave May PetscSFNode *cells; 26705264a50SDave May PetscInt p, bs, dim, npoints, nfound; 26805264a50SDave May const PetscInt *boxCells; 26905264a50SDave May 27005264a50SDave May PetscFunctionBegin; 2719566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(pos, &dim)); 27205264a50SDave May switch (dim) { 273d71ae5a4SJacob Faibussowitsch case 1: 274d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Support not provided for 1D"); 275d71ae5a4SJacob Faibussowitsch case 2: 276d71ae5a4SJacob Faibussowitsch PetscCall(private_DMDALocatePointsIS_2D_Regular(dm, pos, &iscell)); 277d71ae5a4SJacob Faibussowitsch break; 278d71ae5a4SJacob Faibussowitsch case 3: 279d71ae5a4SJacob Faibussowitsch PetscCall(private_DMDALocatePointsIS_3D_Regular(dm, pos, &iscell)); 280d71ae5a4SJacob Faibussowitsch break; 281d71ae5a4SJacob Faibussowitsch default: 282da81f932SPierre Jolivet SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported spatial dimension"); 28305264a50SDave May } 28405264a50SDave May 2859566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(pos, &npoints)); 2869566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(pos, &bs)); 28705264a50SDave May npoints = npoints / bs; 28805264a50SDave May 2899566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(npoints, &cells)); 2909566063dSJacob Faibussowitsch PetscCall(ISGetIndices(iscell, &boxCells)); 29105264a50SDave May 29205264a50SDave May for (p = 0; p < npoints; p++) { 29305264a50SDave May cells[p].rank = 0; 29405264a50SDave May cells[p].index = boxCells[p]; 29505264a50SDave May } 2969566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(iscell, &boxCells)); 29705264a50SDave May 29805264a50SDave May nfound = npoints; 2999566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(cellSF, npoints, nfound, NULL, PETSC_OWN_POINTER, cells, PETSC_OWN_POINTER)); 3009566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscell)); 3013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 30205264a50SDave May } 303