xref: /petsc/src/dm/impls/da/dageometry.c (revision 705a7f0eb0ef90644623ad86b4e8bfdbdabd856c)
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:
10*72fd0fbdSBarry Smith + dm - the `DMDA`
11*72fd0fbdSBarry 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 
18*72fd0fbdSBarry Smith   Note:
19*72fd0fbdSBarry Smith   The (i,j,k) are in the local numbering of the `DMDA`. That is they are non-negative offsets to the ghost corners returned by `DMDAGetGhostCorners()`
20*72fd0fbdSBarry Smith 
21*72fd0fbdSBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDAGetGhostCorners()`
22ed4e918cSMatthew G Knepley @*/
DMDAConvertToCell(DM dm,MatStencil s,PetscInt * cell)23d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAConvertToCell(DM dm, MatStencil s, PetscInt *cell)
24d71ae5a4SJacob Faibussowitsch {
25341c9bdaSMatthew G Knepley   DM_DA         *da  = (DM_DA *)dm->data;
26c73cfb54SMatthew G. Knepley   const PetscInt dim = dm->dim;
274e846693SMatthew G Knepley   const PetscInt mx = (da->Xe - da->Xs) / da->w, my = da->Ye - da->Ys /*, mz = da->Ze - da->Zs*/;
28ed4e918cSMatthew 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;
29341c9bdaSMatthew G Knepley 
30341c9bdaSMatthew G Knepley   PetscFunctionBegin;
31ed4e918cSMatthew G Knepley   *cell = -1;
3263a3b9bcSJacob 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);
331dca8a05SBarry 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);
341dca8a05SBarry 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);
35ed4e918cSMatthew G Knepley   *cell = (kl * my + jl) * mx + il;
363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
37341c9bdaSMatthew G Knepley }
38341c9bdaSMatthew G Knepley 
DMGetLocalBoundingBox_DA(DM da,PetscReal lmin[],PetscReal lmax[],PetscInt cs[],PetscInt ce[])396c6a6b79SMatthew G. Knepley PetscErrorCode DMGetLocalBoundingBox_DA(DM da, PetscReal lmin[], PetscReal lmax[], PetscInt cs[], PetscInt ce[])
40d71ae5a4SJacob Faibussowitsch {
416c6a6b79SMatthew G. Knepley   PetscInt           xs, xe, Xs, Xe;
426c6a6b79SMatthew G. Knepley   PetscInt           ys, ye, Ys, Ye;
436c6a6b79SMatthew G. Knepley   PetscInt           zs, ze, Zs, Ze;
446c6a6b79SMatthew G. Knepley   PetscInt           dim, M, N, P, c0, c1;
456c6a6b79SMatthew G. Knepley   PetscReal          gmax[3] = {0., 0., 0.};
466c6a6b79SMatthew G. Knepley   const PetscReal   *L, *Lstart;
476c6a6b79SMatthew G. Knepley   Vec                coordinates;
486c6a6b79SMatthew G. Knepley   const PetscScalar *coor;
496c6a6b79SMatthew G. Knepley   DMBoundaryType     bx, by, bz;
5005264a50SDave May 
51aab5bcd8SJed Brown   PetscFunctionBegin;
526c6a6b79SMatthew G. Knepley   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xe, &ye, &ze));
536c6a6b79SMatthew G. Knepley   PetscCall(DMDAGetGhostCorners(da, &Xs, &Ys, &Zs, &Xe, &Ye, &Ze));
546c6a6b79SMatthew G. Knepley   PetscCall(DMDAGetInfo(da, &dim, &M, &N, &P, NULL, NULL, NULL, NULL, NULL, &bx, &by, &bz, NULL));
556c6a6b79SMatthew G. Knepley   // Convert from widths to endpoints
569371c9d4SSatish Balay   xe += xs;
579371c9d4SSatish Balay   Xe += Xs;
589371c9d4SSatish Balay   ye += ys;
599371c9d4SSatish Balay   Ye += Ys;
606c6a6b79SMatthew G. Knepley   ze += zs;
616c6a6b79SMatthew G. Knepley   Ze += Zs;
626c6a6b79SMatthew G. Knepley   // What is this doing?
6351c42ae6SMatthew G. Knepley   if (xs != Xs && Xs >= 0) xs -= 1;
6451c42ae6SMatthew G. Knepley   if (ys != Ys && Ys >= 0) ys -= 1;
656c6a6b79SMatthew G. Knepley   if (zs != Zs && Zs >= 0) zs -= 1;
6605264a50SDave May 
676c6a6b79SMatthew G. Knepley   PetscCall(DMGetCoordinatesLocal(da, &coordinates));
686c6a6b79SMatthew G. Knepley   if (!coordinates) {
696c6a6b79SMatthew G. Knepley     PetscCall(DMGetLocalBoundingIndices_DMDA(da, lmin, lmax));
706c6a6b79SMatthew G. Knepley     for (PetscInt d = 0; d < dim; ++d) {
716c6a6b79SMatthew G. Knepley       if (cs) cs[d] = lmin[d];
726c6a6b79SMatthew G. Knepley       if (ce) ce[d] = lmax[d];
736c6a6b79SMatthew G. Knepley     }
746c6a6b79SMatthew G. Knepley     PetscFunctionReturn(PETSC_SUCCESS);
756c6a6b79SMatthew G. Knepley   }
766c6a6b79SMatthew G. Knepley   PetscCall(VecGetArrayRead(coordinates, &coor));
776c6a6b79SMatthew G. Knepley   switch (dim) {
786c6a6b79SMatthew G. Knepley   case 1:
796c6a6b79SMatthew G. Knepley     c0 = (xs - Xs);
806c6a6b79SMatthew G. Knepley     c1 = (xe - 2 - Xs + 1);
816c6a6b79SMatthew G. Knepley     break;
826c6a6b79SMatthew G. Knepley   case 2:
8305264a50SDave May     c0 = (xs - Xs) + (ys - Ys) * (Xe - Xs);
8405264a50SDave May     c1 = (xe - 2 - Xs + 1) + (ye - 2 - Ys + 1) * (Xe - Xs);
856c6a6b79SMatthew G. Knepley     break;
866c6a6b79SMatthew G. Knepley   case 3:
876c6a6b79SMatthew G. Knepley     c0 = (xs - Xs) + (ys - Ys) * (Xe - Xs) + (zs - Zs) * (Xe - Xs) * (Ye - Ys);
886c6a6b79SMatthew G. Knepley     c1 = (xe - 2 - Xs + 1) + (ye - 2 - Ys + 1) * (Xe - Xs) + (ze - 2 - Zs + 1) * (Xe - Xs) * (Ye - Ys);
896c6a6b79SMatthew G. Knepley     break;
906c6a6b79SMatthew G. Knepley   default:
916c6a6b79SMatthew G. Knepley     SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONG, "Invalid dimension %" PetscInt_FMT " for DMDA", dim);
926c6a6b79SMatthew G. Knepley   }
936c6a6b79SMatthew G. Knepley   for (PetscInt d = 0; d < dim; ++d) {
946c6a6b79SMatthew G. Knepley     lmin[d] = PetscRealPart(coor[c0 * dim + d]);
956c6a6b79SMatthew G. Knepley     lmax[d] = PetscRealPart(coor[c1 * dim + d]);
966c6a6b79SMatthew G. Knepley   }
976c6a6b79SMatthew G. Knepley   PetscCall(VecRestoreArrayRead(coordinates, &coor));
9805264a50SDave May 
996c6a6b79SMatthew G. Knepley   PetscCall(DMGetPeriodicity(da, NULL, &Lstart, &L));
1006c6a6b79SMatthew G. Knepley   if (L) {
1016c6a6b79SMatthew G. Knepley     for (PetscInt d = 0; d < dim; ++d)
1026c6a6b79SMatthew G. Knepley       if (L[d] > 0.0) gmax[d] = Lstart[d] + L[d];
1036c6a6b79SMatthew G. Knepley   }
1046c6a6b79SMatthew G. Knepley   // Must check for periodic boundary
1056c6a6b79SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC && xe == M) {
1066c6a6b79SMatthew G. Knepley     lmax[0] = gmax[0];
1076c6a6b79SMatthew G. Knepley     ++xe;
1086c6a6b79SMatthew G. Knepley   }
1096c6a6b79SMatthew G. Knepley   if (by == DM_BOUNDARY_PERIODIC && ye == N) {
1106c6a6b79SMatthew G. Knepley     lmax[1] = gmax[1];
1116c6a6b79SMatthew G. Knepley     ++ye;
1126c6a6b79SMatthew G. Knepley   }
1136c6a6b79SMatthew G. Knepley   if (bz == DM_BOUNDARY_PERIODIC && ze == P) {
1146c6a6b79SMatthew G. Knepley     lmax[2] = gmax[2];
1156c6a6b79SMatthew G. Knepley     ++ze;
1166c6a6b79SMatthew G. Knepley   }
1176c6a6b79SMatthew G. Knepley   if (cs) {
1186c6a6b79SMatthew G. Knepley     cs[0] = xs;
1196c6a6b79SMatthew G. Knepley     if (dim > 1) cs[1] = ys;
1206c6a6b79SMatthew G. Knepley     if (dim > 2) cs[2] = zs;
1216c6a6b79SMatthew G. Knepley   }
1226c6a6b79SMatthew G. Knepley   if (ce) {
1236c6a6b79SMatthew G. Knepley     ce[0] = xe;
1246c6a6b79SMatthew G. Knepley     if (dim > 1) ce[1] = ye;
1256c6a6b79SMatthew G. Knepley     if (dim > 2) ce[2] = ze;
1266c6a6b79SMatthew G. Knepley   }
1276c6a6b79SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
1286c6a6b79SMatthew G. Knepley }
12905264a50SDave May 
private_DMDALocatePointsIS_2D_Regular(DM dmregular,Vec pos,IS * iscell)1306c6a6b79SMatthew G. Knepley static PetscErrorCode private_DMDALocatePointsIS_2D_Regular(DM dmregular, Vec pos, IS *iscell)
1316c6a6b79SMatthew G. Knepley {
1326c6a6b79SMatthew G. Knepley   PetscInt           n, bs, npoints;
1336c6a6b79SMatthew G. Knepley   PetscInt           cs[2], ce[2];
1346c6a6b79SMatthew G. Knepley   PetscInt           xs, xe, mxlocal;
1356c6a6b79SMatthew G. Knepley   PetscInt           ys, ye, mylocal;
1366c6a6b79SMatthew G. Knepley   PetscReal          gmin_l[2], gmax_l[2], dx[2];
1376c6a6b79SMatthew G. Knepley   PetscReal          gmin[2], gmax[2];
1386c6a6b79SMatthew G. Knepley   PetscInt          *cellidx;
1396c6a6b79SMatthew G. Knepley   const PetscScalar *coor;
1406c6a6b79SMatthew G. Knepley 
1416c6a6b79SMatthew G. Knepley   PetscFunctionBegin;
1426c6a6b79SMatthew G. Knepley   PetscCall(DMGetLocalBoundingBox_DA(dmregular, gmin_l, gmax_l, cs, ce));
1436c6a6b79SMatthew G. Knepley   xs = cs[0];
1446c6a6b79SMatthew G. Knepley   ys = cs[1];
1456c6a6b79SMatthew G. Knepley   xe = ce[0];
1466c6a6b79SMatthew G. Knepley   ye = ce[1];
1476c6a6b79SMatthew G. Knepley   PetscCall(DMGetBoundingBox(dmregular, gmin, gmax));
14851c42ae6SMatthew G. Knepley 
14951c42ae6SMatthew G. Knepley   mxlocal = xe - xs - 1;
15051c42ae6SMatthew G. Knepley   mylocal = ye - ys - 1;
15105264a50SDave May 
15205264a50SDave May   dx[0] = (gmax_l[0] - gmin_l[0]) / ((PetscReal)mxlocal);
15305264a50SDave May   dx[1] = (gmax_l[1] - gmin_l[1]) / ((PetscReal)mylocal);
15405264a50SDave May 
1559566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(pos, &n));
1569566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(pos, &bs));
15705264a50SDave May   npoints = n / bs;
15805264a50SDave May 
1599566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(npoints, &cellidx));
1606c6a6b79SMatthew G. Knepley   PetscCall(VecGetArrayRead(pos, &coor));
1616c6a6b79SMatthew G. Knepley   for (PetscInt p = 0; p < npoints; p++) {
1620ca99263SDave May     PetscReal coor_p[2];
16305264a50SDave May     PetscInt  mi[2];
16405264a50SDave May 
1656c6a6b79SMatthew G. Knepley     coor_p[0] = PetscRealPart(coor[2 * p]);
1666c6a6b79SMatthew G. Knepley     coor_p[1] = PetscRealPart(coor[2 * p + 1]);
16705264a50SDave May 
16805264a50SDave May     cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND;
16905264a50SDave May 
170ad540459SPierre Jolivet     if (coor_p[0] < gmin_l[0]) continue;
171ad540459SPierre Jolivet     if (coor_p[0] > gmax_l[0]) continue;
172ad540459SPierre Jolivet     if (coor_p[1] < gmin_l[1]) continue;
173ad540459SPierre Jolivet     if (coor_p[1] > gmax_l[1]) continue;
17405264a50SDave May 
1756c6a6b79SMatthew G. Knepley     for (PetscInt d = 0; d < 2; d++) mi[d] = (PetscInt)((coor_p[d] - gmin[d]) / dx[d]);
17605264a50SDave May 
177ad540459SPierre Jolivet     if (mi[0] < xs) continue;
178ad540459SPierre Jolivet     if (mi[0] > (xe - 1)) continue;
179ad540459SPierre Jolivet     if (mi[1] < ys) continue;
180ad540459SPierre Jolivet     if (mi[1] > (ye - 1)) continue;
18105264a50SDave May 
182ad540459SPierre Jolivet     if (mi[0] == (xe - 1)) mi[0]--;
183ad540459SPierre Jolivet     if (mi[1] == (ye - 1)) mi[1]--;
18405264a50SDave May 
18505264a50SDave May     cellidx[p] = (mi[0] - xs) + (mi[1] - ys) * mxlocal;
18605264a50SDave May   }
1876c6a6b79SMatthew G. Knepley   PetscCall(VecRestoreArrayRead(pos, &coor));
1889566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, npoints, cellidx, PETSC_OWN_POINTER, iscell));
1893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19005264a50SDave May }
19105264a50SDave May 
private_DMDALocatePointsIS_3D_Regular(DM dmregular,Vec pos,IS * iscell)19266976f2fSJacob Faibussowitsch static PetscErrorCode private_DMDALocatePointsIS_3D_Regular(DM dmregular, Vec pos, IS *iscell)
193d71ae5a4SJacob Faibussowitsch {
1946c6a6b79SMatthew G. Knepley   PetscInt           n, bs, npoints;
1956c6a6b79SMatthew G. Knepley   PetscInt           cs[3], ce[3];
1966c6a6b79SMatthew G. Knepley   PetscInt           xs, xe, mxlocal;
1976c6a6b79SMatthew G. Knepley   PetscInt           ys, ye, mylocal;
1986c6a6b79SMatthew G. Knepley   PetscInt           zs, ze, mzlocal;
19905264a50SDave May   PetscReal          gmin_l[3], gmax_l[3], dx[3];
20005264a50SDave May   PetscReal          gmin[3], gmax[3];
20105264a50SDave May   PetscInt          *cellidx;
2026c6a6b79SMatthew G. Knepley   const PetscScalar *coor;
20305264a50SDave May 
204aab5bcd8SJed Brown   PetscFunctionBegin;
2056c6a6b79SMatthew G. Knepley   PetscCall(DMGetLocalBoundingBox_DA(dmregular, gmin_l, gmax_l, cs, ce));
2066c6a6b79SMatthew G. Knepley   xs = cs[0];
2076c6a6b79SMatthew G. Knepley   ys = cs[1];
2086c6a6b79SMatthew G. Knepley   zs = cs[2];
2096c6a6b79SMatthew G. Knepley   xe = ce[0];
2106c6a6b79SMatthew G. Knepley   ye = ce[1];
2116c6a6b79SMatthew G. Knepley   ze = ce[2];
2126c6a6b79SMatthew G. Knepley   PetscCall(DMGetBoundingBox(dmregular, gmin, gmax));
21351c42ae6SMatthew G. Knepley 
21451c42ae6SMatthew G. Knepley   mxlocal = xe - xs - 1;
21551c42ae6SMatthew G. Knepley   mylocal = ye - ys - 1;
21651c42ae6SMatthew G. Knepley   mzlocal = ze - zs - 1;
21705264a50SDave May 
21805264a50SDave May   dx[0] = (gmax_l[0] - gmin_l[0]) / ((PetscReal)mxlocal);
21905264a50SDave May   dx[1] = (gmax_l[1] - gmin_l[1]) / ((PetscReal)mylocal);
22005264a50SDave May   dx[2] = (gmax_l[2] - gmin_l[2]) / ((PetscReal)mzlocal);
22105264a50SDave May 
2229566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(pos, &n));
2239566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(pos, &bs));
22405264a50SDave May   npoints = n / bs;
22505264a50SDave May 
2269566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(npoints, &cellidx));
2276c6a6b79SMatthew G. Knepley   PetscCall(VecGetArrayRead(pos, &coor));
2286c6a6b79SMatthew G. Knepley   for (PetscInt p = 0; p < npoints; p++) {
2290ca99263SDave May     PetscReal coor_p[3];
23005264a50SDave May     PetscInt  mi[3];
23105264a50SDave May 
2326c6a6b79SMatthew G. Knepley     coor_p[0] = PetscRealPart(coor[3 * p]);
2336c6a6b79SMatthew G. Knepley     coor_p[1] = PetscRealPart(coor[3 * p + 1]);
2346c6a6b79SMatthew G. Knepley     coor_p[2] = PetscRealPart(coor[3 * p + 2]);
23505264a50SDave May 
23605264a50SDave May     cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND;
23705264a50SDave May 
238ad540459SPierre Jolivet     if (coor_p[0] < gmin_l[0]) continue;
239ad540459SPierre Jolivet     if (coor_p[0] > gmax_l[0]) continue;
240ad540459SPierre Jolivet     if (coor_p[1] < gmin_l[1]) continue;
241ad540459SPierre Jolivet     if (coor_p[1] > gmax_l[1]) continue;
242ad540459SPierre Jolivet     if (coor_p[2] < gmin_l[2]) continue;
243ad540459SPierre Jolivet     if (coor_p[2] > gmax_l[2]) continue;
24405264a50SDave May 
2456c6a6b79SMatthew G. Knepley     for (PetscInt d = 0; d < 3; d++) mi[d] = (PetscInt)((coor_p[d] - gmin[d]) / dx[d]);
24605264a50SDave May 
2476c6a6b79SMatthew G. Knepley     // TODO: Check for periodicity here
248ad540459SPierre Jolivet     if (mi[0] < xs) continue;
249ad540459SPierre Jolivet     if (mi[0] > (xe - 1)) continue;
250ad540459SPierre Jolivet     if (mi[1] < ys) continue;
251ad540459SPierre Jolivet     if (mi[1] > (ye - 1)) continue;
252ad540459SPierre Jolivet     if (mi[2] < zs) continue;
253ad540459SPierre Jolivet     if (mi[2] > (ze - 1)) continue;
25405264a50SDave May 
255ad540459SPierre Jolivet     if (mi[0] == (xe - 1)) mi[0]--;
256ad540459SPierre Jolivet     if (mi[1] == (ye - 1)) mi[1]--;
257ad540459SPierre Jolivet     if (mi[2] == (ze - 1)) mi[2]--;
25805264a50SDave May 
25905264a50SDave May     cellidx[p] = (mi[0] - xs) + (mi[1] - ys) * mxlocal + (mi[2] - zs) * mxlocal * mylocal;
26005264a50SDave May   }
2616c6a6b79SMatthew G. Knepley   PetscCall(VecRestoreArrayRead(pos, &coor));
2629566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, npoints, cellidx, PETSC_OWN_POINTER, iscell));
2633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
26405264a50SDave May }
26505264a50SDave May 
DMLocatePoints_DA_Regular(DM dm,Vec pos,DMPointLocationType ltype,PetscSF cellSF)266d71ae5a4SJacob Faibussowitsch PetscErrorCode DMLocatePoints_DA_Regular(DM dm, Vec pos, DMPointLocationType ltype, PetscSF cellSF)
267d71ae5a4SJacob Faibussowitsch {
26805264a50SDave May   IS              iscell;
26905264a50SDave May   PetscSFNode    *cells;
27005264a50SDave May   PetscInt        p, bs, dim, npoints, nfound;
27105264a50SDave May   const PetscInt *boxCells;
27205264a50SDave May 
27305264a50SDave May   PetscFunctionBegin;
2749566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(pos, &dim));
27505264a50SDave May   switch (dim) {
276d71ae5a4SJacob Faibussowitsch   case 1:
277d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Support not provided for 1D");
278d71ae5a4SJacob Faibussowitsch   case 2:
279d71ae5a4SJacob Faibussowitsch     PetscCall(private_DMDALocatePointsIS_2D_Regular(dm, pos, &iscell));
280d71ae5a4SJacob Faibussowitsch     break;
281d71ae5a4SJacob Faibussowitsch   case 3:
282d71ae5a4SJacob Faibussowitsch     PetscCall(private_DMDALocatePointsIS_3D_Regular(dm, pos, &iscell));
283d71ae5a4SJacob Faibussowitsch     break;
284d71ae5a4SJacob Faibussowitsch   default:
285da81f932SPierre Jolivet     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported spatial dimension");
28605264a50SDave May   }
28705264a50SDave May 
2889566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(pos, &npoints));
2899566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(pos, &bs));
29005264a50SDave May   npoints = npoints / bs;
29105264a50SDave May 
2929566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(npoints, &cells));
2939566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(iscell, &boxCells));
29405264a50SDave May 
29505264a50SDave May   for (p = 0; p < npoints; p++) {
29605264a50SDave May     cells[p].rank  = 0;
29705264a50SDave May     cells[p].index = boxCells[p];
29805264a50SDave May   }
2999566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(iscell, &boxCells));
30005264a50SDave May 
30105264a50SDave May   nfound = npoints;
3029566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(cellSF, npoints, nfound, NULL, PETSC_OWN_POINTER, cells, PETSC_OWN_POINTER));
3039566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iscell));
3043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
30505264a50SDave May }
306