xref: /petsc/src/dm/impls/da/dageometry.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
105264a50SDave May #include <petscsf.h>
2af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I  "petscdmda.h"   I*/
357459e9aSMatthew G Knepley 
4ed4e918cSMatthew G Knepley /*@
5f0e3914cSSatish Balay   DMDAConvertToCell - Convert (i,j,k) to local cell number
6341c9bdaSMatthew G Knepley 
7ed4e918cSMatthew G Knepley   Not Collective
8ed4e918cSMatthew G Knepley 
9d8d19677SJose E. Roman   Input Parameters:
10ed4e918cSMatthew G Knepley + da - the distributed array
11907376e6SBarry Smith - s - A MatStencil giving (i,j,k)
12ed4e918cSMatthew G Knepley 
13ed4e918cSMatthew G Knepley   Output Parameter:
14ed4e918cSMatthew G Knepley . cell - the local cell number
15ed4e918cSMatthew G Knepley 
16ed4e918cSMatthew G Knepley   Level: developer
17ed4e918cSMatthew G Knepley @*/
18*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAConvertToCell(DM dm, MatStencil s, PetscInt *cell)
19*d71ae5a4SJacob Faibussowitsch {
20341c9bdaSMatthew G Knepley   DM_DA         *da  = (DM_DA *)dm->data;
21c73cfb54SMatthew G. Knepley   const PetscInt dim = dm->dim;
224e846693SMatthew G Knepley   const PetscInt mx = (da->Xe - da->Xs) / da->w, my = da->Ye - da->Ys /*, mz = da->Ze - da->Zs*/;
23ed4e918cSMatthew 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;
24341c9bdaSMatthew G Knepley 
25341c9bdaSMatthew G Knepley   PetscFunctionBegin;
26ed4e918cSMatthew G Knepley   *cell = -1;
2763a3b9bcSJacob 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);
281dca8a05SBarry 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);
291dca8a05SBarry 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);
30ed4e918cSMatthew G Knepley   *cell = (kl * my + jl) * mx + il;
31ed4e918cSMatthew G Knepley   PetscFunctionReturn(0);
32341c9bdaSMatthew G Knepley }
33341c9bdaSMatthew G Knepley 
34*d71ae5a4SJacob Faibussowitsch PetscErrorCode private_DMDALocatePointsIS_2D_Regular(DM dmregular, Vec pos, IS *iscell)
35*d71ae5a4SJacob Faibussowitsch {
3605264a50SDave May   PetscInt           n, bs, p, npoints;
3705264a50SDave May   PetscInt           xs, xe, Xs, Xe, mxlocal;
3805264a50SDave May   PetscInt           ys, ye, Ys, Ye, mylocal;
39aab5bcd8SJed Brown   PetscInt           d, c0, c1;
4005264a50SDave May   PetscReal          gmin_l[2], gmax_l[2], dx[2];
4105264a50SDave May   PetscReal          gmin[2], gmax[2];
4205264a50SDave May   PetscInt          *cellidx;
4305264a50SDave May   Vec                coor;
4405264a50SDave May   const PetscScalar *_coor;
4505264a50SDave May 
46aab5bcd8SJed Brown   PetscFunctionBegin;
479566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dmregular, &xs, &ys, NULL, &xe, &ye, NULL));
489566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(dmregular, &Xs, &Ys, NULL, &Xe, &Ye, NULL));
499371c9d4SSatish Balay   xe += xs;
509371c9d4SSatish Balay   Xe += Xs;
519371c9d4SSatish Balay   if (xs != Xs) xs -= 1;
529371c9d4SSatish Balay   ye += ys;
539371c9d4SSatish Balay   Ye += Ys;
549371c9d4SSatish Balay   if (ys != Ys) ys -= 1;
5505264a50SDave May 
5605264a50SDave May   mxlocal = xe - xs - 1;
5705264a50SDave May   mylocal = ye - ys - 1;
5805264a50SDave May 
599566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dmregular, &coor));
609566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(coor, &_coor));
6105264a50SDave May   c0 = (xs - Xs) + (ys - Ys) * (Xe - Xs);
6205264a50SDave May   c1 = (xe - 2 - Xs + 1) + (ye - 2 - Ys + 1) * (Xe - Xs);
6305264a50SDave May 
64a5f152d1SDave May   gmin_l[0] = PetscRealPart(_coor[2 * c0 + 0]);
65a5f152d1SDave May   gmin_l[1] = PetscRealPart(_coor[2 * c0 + 1]);
6605264a50SDave May 
67a5f152d1SDave May   gmax_l[0] = PetscRealPart(_coor[2 * c1 + 0]);
68a5f152d1SDave May   gmax_l[1] = PetscRealPart(_coor[2 * c1 + 1]);
6905264a50SDave May 
7005264a50SDave May   dx[0] = (gmax_l[0] - gmin_l[0]) / ((PetscReal)mxlocal);
7105264a50SDave May   dx[1] = (gmax_l[1] - gmin_l[1]) / ((PetscReal)mylocal);
7205264a50SDave May 
739566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(coor, &_coor));
7405264a50SDave May 
759566063dSJacob Faibussowitsch   PetscCall(DMGetBoundingBox(dmregular, gmin, gmax));
7605264a50SDave May 
779566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(pos, &n));
789566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(pos, &bs));
7905264a50SDave May   npoints = n / bs;
8005264a50SDave May 
819566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(npoints, &cellidx));
829566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(pos, &_coor));
8305264a50SDave May   for (p = 0; p < npoints; p++) {
840ca99263SDave May     PetscReal coor_p[2];
8505264a50SDave May     PetscInt  mi[2];
8605264a50SDave May 
870ca99263SDave May     coor_p[0] = PetscRealPart(_coor[2 * p]);
880ca99263SDave May     coor_p[1] = PetscRealPart(_coor[2 * p + 1]);
8905264a50SDave May 
9005264a50SDave May     cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND;
9105264a50SDave May 
92ad540459SPierre Jolivet     if (coor_p[0] < gmin_l[0]) continue;
93ad540459SPierre Jolivet     if (coor_p[0] > gmax_l[0]) continue;
94ad540459SPierre Jolivet     if (coor_p[1] < gmin_l[1]) continue;
95ad540459SPierre Jolivet     if (coor_p[1] > gmax_l[1]) continue;
9605264a50SDave May 
97ad540459SPierre Jolivet     for (d = 0; d < 2; d++) mi[d] = (PetscInt)((coor_p[d] - gmin[d]) / dx[d]);
9805264a50SDave May 
99ad540459SPierre Jolivet     if (mi[0] < xs) continue;
100ad540459SPierre Jolivet     if (mi[0] > (xe - 1)) continue;
101ad540459SPierre Jolivet     if (mi[1] < ys) continue;
102ad540459SPierre Jolivet     if (mi[1] > (ye - 1)) continue;
10305264a50SDave May 
104ad540459SPierre Jolivet     if (mi[0] == (xe - 1)) mi[0]--;
105ad540459SPierre Jolivet     if (mi[1] == (ye - 1)) mi[1]--;
10605264a50SDave May 
10705264a50SDave May     cellidx[p] = (mi[0] - xs) + (mi[1] - ys) * mxlocal;
10805264a50SDave May   }
1099566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(pos, &_coor));
1109566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, npoints, cellidx, PETSC_OWN_POINTER, iscell));
11105264a50SDave May   PetscFunctionReturn(0);
11205264a50SDave May }
11305264a50SDave May 
114*d71ae5a4SJacob Faibussowitsch PetscErrorCode private_DMDALocatePointsIS_3D_Regular(DM dmregular, Vec pos, IS *iscell)
115*d71ae5a4SJacob Faibussowitsch {
11605264a50SDave May   PetscInt           n, bs, p, npoints;
11705264a50SDave May   PetscInt           xs, xe, Xs, Xe, mxlocal;
11805264a50SDave May   PetscInt           ys, ye, Ys, Ye, mylocal;
11905264a50SDave May   PetscInt           zs, ze, Zs, Ze, mzlocal;
120aab5bcd8SJed Brown   PetscInt           d, c0, c1;
12105264a50SDave May   PetscReal          gmin_l[3], gmax_l[3], dx[3];
12205264a50SDave May   PetscReal          gmin[3], gmax[3];
12305264a50SDave May   PetscInt          *cellidx;
12405264a50SDave May   Vec                coor;
12505264a50SDave May   const PetscScalar *_coor;
12605264a50SDave May 
127aab5bcd8SJed Brown   PetscFunctionBegin;
1289566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dmregular, &xs, &ys, &zs, &xe, &ye, &ze));
1299566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(dmregular, &Xs, &Ys, &Zs, &Xe, &Ye, &Ze));
1309371c9d4SSatish Balay   xe += xs;
1319371c9d4SSatish Balay   Xe += Xs;
1329371c9d4SSatish Balay   if (xs != Xs) xs -= 1;
1339371c9d4SSatish Balay   ye += ys;
1349371c9d4SSatish Balay   Ye += Ys;
1359371c9d4SSatish Balay   if (ys != Ys) ys -= 1;
1369371c9d4SSatish Balay   ze += zs;
1379371c9d4SSatish Balay   Ze += Zs;
1389371c9d4SSatish Balay   if (zs != Zs) zs -= 1;
13905264a50SDave May 
14005264a50SDave May   mxlocal = xe - xs - 1;
14105264a50SDave May   mylocal = ye - ys - 1;
14205264a50SDave May   mzlocal = ze - zs - 1;
14305264a50SDave May 
1449566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dmregular, &coor));
1459566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(coor, &_coor));
14605264a50SDave May   c0 = (xs - Xs) + (ys - Ys) * (Xe - Xs) + (zs - Zs) * (Xe - Xs) * (Ye - Ys);
14705264a50SDave May   c1 = (xe - 2 - Xs + 1) + (ye - 2 - Ys + 1) * (Xe - Xs) + (ze - 2 - Zs + 1) * (Xe - Xs) * (Ye - Ys);
14805264a50SDave May 
149a5f152d1SDave May   gmin_l[0] = PetscRealPart(_coor[3 * c0 + 0]);
150a5f152d1SDave May   gmin_l[1] = PetscRealPart(_coor[3 * c0 + 1]);
151a5f152d1SDave May   gmin_l[2] = PetscRealPart(_coor[3 * c0 + 2]);
15205264a50SDave May 
153a5f152d1SDave May   gmax_l[0] = PetscRealPart(_coor[3 * c1 + 0]);
154a5f152d1SDave May   gmax_l[1] = PetscRealPart(_coor[3 * c1 + 1]);
155a5f152d1SDave May   gmax_l[2] = PetscRealPart(_coor[3 * c1 + 2]);
15605264a50SDave May 
15705264a50SDave May   dx[0] = (gmax_l[0] - gmin_l[0]) / ((PetscReal)mxlocal);
15805264a50SDave May   dx[1] = (gmax_l[1] - gmin_l[1]) / ((PetscReal)mylocal);
15905264a50SDave May   dx[2] = (gmax_l[2] - gmin_l[2]) / ((PetscReal)mzlocal);
16005264a50SDave May 
1619566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(coor, &_coor));
16205264a50SDave May 
1639566063dSJacob Faibussowitsch   PetscCall(DMGetBoundingBox(dmregular, gmin, gmax));
16405264a50SDave May 
1659566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(pos, &n));
1669566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(pos, &bs));
16705264a50SDave May   npoints = n / bs;
16805264a50SDave May 
1699566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(npoints, &cellidx));
1709566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(pos, &_coor));
17105264a50SDave May   for (p = 0; p < npoints; p++) {
1720ca99263SDave May     PetscReal coor_p[3];
17305264a50SDave May     PetscInt  mi[3];
17405264a50SDave May 
1750ca99263SDave May     coor_p[0] = PetscRealPart(_coor[3 * p]);
1760ca99263SDave May     coor_p[1] = PetscRealPart(_coor[3 * p + 1]);
1770ca99263SDave May     coor_p[2] = PetscRealPart(_coor[3 * p + 2]);
17805264a50SDave May 
17905264a50SDave May     cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND;
18005264a50SDave May 
181ad540459SPierre Jolivet     if (coor_p[0] < gmin_l[0]) continue;
182ad540459SPierre Jolivet     if (coor_p[0] > gmax_l[0]) continue;
183ad540459SPierre Jolivet     if (coor_p[1] < gmin_l[1]) continue;
184ad540459SPierre Jolivet     if (coor_p[1] > gmax_l[1]) continue;
185ad540459SPierre Jolivet     if (coor_p[2] < gmin_l[2]) continue;
186ad540459SPierre Jolivet     if (coor_p[2] > gmax_l[2]) continue;
18705264a50SDave May 
188ad540459SPierre Jolivet     for (d = 0; d < 3; d++) mi[d] = (PetscInt)((coor_p[d] - gmin[d]) / dx[d]);
18905264a50SDave May 
190ad540459SPierre Jolivet     if (mi[0] < xs) continue;
191ad540459SPierre Jolivet     if (mi[0] > (xe - 1)) continue;
192ad540459SPierre Jolivet     if (mi[1] < ys) continue;
193ad540459SPierre Jolivet     if (mi[1] > (ye - 1)) continue;
194ad540459SPierre Jolivet     if (mi[2] < zs) continue;
195ad540459SPierre Jolivet     if (mi[2] > (ze - 1)) continue;
19605264a50SDave May 
197ad540459SPierre Jolivet     if (mi[0] == (xe - 1)) mi[0]--;
198ad540459SPierre Jolivet     if (mi[1] == (ye - 1)) mi[1]--;
199ad540459SPierre Jolivet     if (mi[2] == (ze - 1)) mi[2]--;
20005264a50SDave May 
20105264a50SDave May     cellidx[p] = (mi[0] - xs) + (mi[1] - ys) * mxlocal + (mi[2] - zs) * mxlocal * mylocal;
20205264a50SDave May   }
2039566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(pos, &_coor));
2049566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, npoints, cellidx, PETSC_OWN_POINTER, iscell));
20505264a50SDave May   PetscFunctionReturn(0);
20605264a50SDave May }
20705264a50SDave May 
208*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMLocatePoints_DA_Regular(DM dm, Vec pos, DMPointLocationType ltype, PetscSF cellSF)
209*d71ae5a4SJacob Faibussowitsch {
21005264a50SDave May   IS              iscell;
21105264a50SDave May   PetscSFNode    *cells;
21205264a50SDave May   PetscInt        p, bs, dim, npoints, nfound;
21305264a50SDave May   const PetscInt *boxCells;
21405264a50SDave May 
21505264a50SDave May   PetscFunctionBegin;
2169566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(pos, &dim));
21705264a50SDave May   switch (dim) {
218*d71ae5a4SJacob Faibussowitsch   case 1:
219*d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Support not provided for 1D");
220*d71ae5a4SJacob Faibussowitsch   case 2:
221*d71ae5a4SJacob Faibussowitsch     PetscCall(private_DMDALocatePointsIS_2D_Regular(dm, pos, &iscell));
222*d71ae5a4SJacob Faibussowitsch     break;
223*d71ae5a4SJacob Faibussowitsch   case 3:
224*d71ae5a4SJacob Faibussowitsch     PetscCall(private_DMDALocatePointsIS_3D_Regular(dm, pos, &iscell));
225*d71ae5a4SJacob Faibussowitsch     break;
226*d71ae5a4SJacob Faibussowitsch   default:
227*d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupport spatial dimension");
22805264a50SDave May   }
22905264a50SDave May 
2309566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(pos, &npoints));
2319566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(pos, &bs));
23205264a50SDave May   npoints = npoints / bs;
23305264a50SDave May 
2349566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(npoints, &cells));
2359566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(iscell, &boxCells));
23605264a50SDave May 
23705264a50SDave May   for (p = 0; p < npoints; p++) {
23805264a50SDave May     cells[p].rank  = 0;
23905264a50SDave May     cells[p].index = boxCells[p];
24005264a50SDave May   }
2419566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(iscell, &boxCells));
24205264a50SDave May 
24305264a50SDave May   nfound = npoints;
2449566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(cellSF, npoints, nfound, NULL, PETSC_OWN_POINTER, cells, PETSC_OWN_POINTER));
2459566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iscell));
24605264a50SDave May   PetscFunctionReturn(0);
24705264a50SDave May }
248