xref: /petsc/src/dm/impls/da/dageometry.c (revision d8d19677bbccf95218448bee62e6b87f4513e133)
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 
9*d8d19677SJose 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 @*/
18ed4e918cSMatthew G Knepley PetscErrorCode DMDAConvertToCell(DM dm, MatStencil s, PetscInt *cell)
19341c9bdaSMatthew G Knepley {
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;
27cb004a26SBarry Smith   if ((s.i < da->Xs/da->w) || (s.i >= da->Xe/da->w))    SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Stencil i %D should be in [%D, %D)", s.i, da->Xs/da->w, da->Xe/da->w);
28cb004a26SBarry Smith   if ((dim > 1) && ((s.j < da->Ys) || (s.j >= da->Ye))) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Stencil j %D should be in [%D, %D)", s.j, da->Ys, da->Ye);
29cb004a26SBarry Smith   if ((dim > 2) && ((s.k < da->Zs) || (s.k >= da->Ze))) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Stencil k %D should be in [%D, %D)", 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 
3405264a50SDave May PetscErrorCode private_DMDALocatePointsIS_2D_Regular(DM dmregular,Vec pos,IS *iscell)
3505264a50SDave May {
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   PetscErrorCode    ierr;
4605264a50SDave May 
47aab5bcd8SJed Brown   PetscFunctionBegin;
48ea78f98cSLisandro Dalcin   ierr = DMDAGetCorners(dmregular,&xs,&ys,NULL,&xe,&ye,NULL);CHKERRQ(ierr);
49ea78f98cSLisandro Dalcin   ierr = DMDAGetGhostCorners(dmregular,&Xs,&Ys,NULL,&Xe,&Ye,NULL);CHKERRQ(ierr);
5005264a50SDave May   xe += xs; Xe += Xs; if (xs != Xs) xs -= 1;
5105264a50SDave May   ye += ys; Ye += Ys; if (ys != Ys) ys -= 1;
5205264a50SDave May 
5305264a50SDave May   mxlocal = xe - xs - 1;
5405264a50SDave May   mylocal = ye - ys - 1;
5505264a50SDave May 
5605264a50SDave May   ierr = DMGetCoordinatesLocal(dmregular,&coor);CHKERRQ(ierr);
5705264a50SDave May   ierr = VecGetArrayRead(coor,&_coor);CHKERRQ(ierr);
5805264a50SDave May   c0 = (xs-Xs) + (ys-Ys)*(Xe-Xs);
5905264a50SDave May   c1 = (xe-2-Xs+1) + (ye-2-Ys+1)*(Xe-Xs);
6005264a50SDave May 
61a5f152d1SDave May   gmin_l[0] = PetscRealPart(_coor[2*c0+0]);
62a5f152d1SDave May   gmin_l[1] = PetscRealPart(_coor[2*c0+1]);
6305264a50SDave May 
64a5f152d1SDave May   gmax_l[0] = PetscRealPart(_coor[2*c1+0]);
65a5f152d1SDave May   gmax_l[1] = PetscRealPart(_coor[2*c1+1]);
6605264a50SDave May 
6705264a50SDave May   dx[0] = (gmax_l[0]-gmin_l[0])/((PetscReal)mxlocal);
6805264a50SDave May   dx[1] = (gmax_l[1]-gmin_l[1])/((PetscReal)mylocal);
6905264a50SDave May 
7005264a50SDave May   ierr = VecRestoreArrayRead(coor,&_coor);CHKERRQ(ierr);
7105264a50SDave May 
72c5daea06SMatthew G. Knepley   ierr = DMGetBoundingBox(dmregular,gmin,gmax);CHKERRQ(ierr);
7305264a50SDave May 
7405264a50SDave May   ierr = VecGetLocalSize(pos,&n);CHKERRQ(ierr);
7505264a50SDave May   ierr = VecGetBlockSize(pos,&bs);CHKERRQ(ierr);
7605264a50SDave May   npoints = n/bs;
7705264a50SDave May 
7805264a50SDave May   ierr = PetscMalloc1(npoints,&cellidx);CHKERRQ(ierr);
7905264a50SDave May   ierr = VecGetArrayRead(pos,&_coor);CHKERRQ(ierr);
8005264a50SDave May   for (p=0; p<npoints; p++) {
810ca99263SDave May     PetscReal coor_p[2];
8205264a50SDave May     PetscInt  mi[2];
8305264a50SDave May 
840ca99263SDave May     coor_p[0] = PetscRealPart(_coor[2*p]);
850ca99263SDave May     coor_p[1] = PetscRealPart(_coor[2*p+1]);
8605264a50SDave May 
8705264a50SDave May     cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND;
8805264a50SDave May 
890ca99263SDave May     if (coor_p[0] < gmin_l[0]) { continue; }
900ca99263SDave May     if (coor_p[0] > gmax_l[0]) { continue; }
910ca99263SDave May     if (coor_p[1] < gmin_l[1]) { continue; }
920ca99263SDave May     if (coor_p[1] > gmax_l[1]) { continue; }
9305264a50SDave May 
94aab5bcd8SJed Brown     for (d=0; d<2; d++) {
9505264a50SDave May       mi[d] = (PetscInt)((coor_p[d] - gmin[d])/dx[d]);
9605264a50SDave May     }
9705264a50SDave May 
9805264a50SDave May     if (mi[0] < xs)     { continue; }
9905264a50SDave May     if (mi[0] > (xe-1)) { continue; }
10005264a50SDave May     if (mi[1] < ys)     { continue; }
10105264a50SDave May     if (mi[1] > (ye-1)) { continue; }
10205264a50SDave May 
10305264a50SDave May     if (mi[0] == (xe-1)) { mi[0]--; }
10405264a50SDave May     if (mi[1] == (ye-1)) { mi[1]--; }
10505264a50SDave May 
10605264a50SDave May     cellidx[p] = (mi[0]-xs) + (mi[1]-ys) * mxlocal;
10705264a50SDave May   }
10805264a50SDave May   ierr = VecRestoreArrayRead(pos,&_coor);CHKERRQ(ierr);
10905264a50SDave May   ierr = ISCreateGeneral(PETSC_COMM_SELF,npoints,cellidx,PETSC_OWN_POINTER,iscell);CHKERRQ(ierr);
11005264a50SDave May   PetscFunctionReturn(0);
11105264a50SDave May }
11205264a50SDave May 
11305264a50SDave May PetscErrorCode private_DMDALocatePointsIS_3D_Regular(DM dmregular,Vec pos,IS *iscell)
11405264a50SDave May {
11505264a50SDave May   PetscInt          n,bs,p,npoints;
11605264a50SDave May   PetscInt          xs,xe,Xs,Xe,mxlocal;
11705264a50SDave May   PetscInt          ys,ye,Ys,Ye,mylocal;
11805264a50SDave May   PetscInt          zs,ze,Zs,Ze,mzlocal;
119aab5bcd8SJed Brown   PetscInt          d,c0,c1;
12005264a50SDave May   PetscReal         gmin_l[3],gmax_l[3],dx[3];
12105264a50SDave May   PetscReal         gmin[3],gmax[3];
12205264a50SDave May   PetscInt          *cellidx;
12305264a50SDave May   Vec               coor;
12405264a50SDave May   const PetscScalar *_coor;
12505264a50SDave May   PetscErrorCode    ierr;
12605264a50SDave May 
127aab5bcd8SJed Brown   PetscFunctionBegin;
12805264a50SDave May   ierr = DMDAGetCorners(dmregular,&xs,&ys,&zs,&xe,&ye,&ze);CHKERRQ(ierr);
12905264a50SDave May   ierr = DMDAGetGhostCorners(dmregular,&Xs,&Ys,&Zs,&Xe,&Ye,&Ze);CHKERRQ(ierr);
13005264a50SDave May   xe += xs; Xe += Xs; if (xs != Xs) xs -= 1;
13105264a50SDave May   ye += ys; Ye += Ys; if (ys != Ys) ys -= 1;
13205264a50SDave May   ze += zs; Ze += Zs; if (zs != Zs) zs -= 1;
13305264a50SDave May 
13405264a50SDave May   mxlocal = xe - xs - 1;
13505264a50SDave May   mylocal = ye - ys - 1;
13605264a50SDave May   mzlocal = ze - zs - 1;
13705264a50SDave May 
13805264a50SDave May   ierr = DMGetCoordinatesLocal(dmregular,&coor);CHKERRQ(ierr);
13905264a50SDave May   ierr = VecGetArrayRead(coor,&_coor);CHKERRQ(ierr);
14005264a50SDave May   c0 = (xs-Xs)     + (ys-Ys)    *(Xe-Xs) + (zs-Zs)    *(Xe-Xs)*(Ye-Ys);
14105264a50SDave May   c1 = (xe-2-Xs+1) + (ye-2-Ys+1)*(Xe-Xs) + (ze-2-Zs+1)*(Xe-Xs)*(Ye-Ys);
14205264a50SDave May 
143a5f152d1SDave May   gmin_l[0] = PetscRealPart(_coor[3*c0+0]);
144a5f152d1SDave May   gmin_l[1] = PetscRealPart(_coor[3*c0+1]);
145a5f152d1SDave May   gmin_l[2] = PetscRealPart(_coor[3*c0+2]);
14605264a50SDave May 
147a5f152d1SDave May   gmax_l[0] = PetscRealPart(_coor[3*c1+0]);
148a5f152d1SDave May   gmax_l[1] = PetscRealPart(_coor[3*c1+1]);
149a5f152d1SDave May   gmax_l[2] = PetscRealPart(_coor[3*c1+2]);
15005264a50SDave May 
15105264a50SDave May   dx[0] = (gmax_l[0]-gmin_l[0])/((PetscReal)mxlocal);
15205264a50SDave May   dx[1] = (gmax_l[1]-gmin_l[1])/((PetscReal)mylocal);
15305264a50SDave May   dx[2] = (gmax_l[2]-gmin_l[2])/((PetscReal)mzlocal);
15405264a50SDave May 
15505264a50SDave May   ierr = VecRestoreArrayRead(coor,&_coor);CHKERRQ(ierr);
15605264a50SDave May 
157c5daea06SMatthew G. Knepley   ierr = DMGetBoundingBox(dmregular,gmin,gmax);CHKERRQ(ierr);
15805264a50SDave May 
15905264a50SDave May   ierr = VecGetLocalSize(pos,&n);CHKERRQ(ierr);
16005264a50SDave May   ierr = VecGetBlockSize(pos,&bs);CHKERRQ(ierr);
16105264a50SDave May   npoints = n/bs;
16205264a50SDave May 
16305264a50SDave May   ierr = PetscMalloc1(npoints,&cellidx);CHKERRQ(ierr);
16405264a50SDave May   ierr = VecGetArrayRead(pos,&_coor);CHKERRQ(ierr);
16505264a50SDave May   for (p=0; p<npoints; p++) {
1660ca99263SDave May     PetscReal coor_p[3];
16705264a50SDave May     PetscInt  mi[3];
16805264a50SDave May 
1690ca99263SDave May     coor_p[0] = PetscRealPart(_coor[3*p]);
1700ca99263SDave May     coor_p[1] = PetscRealPart(_coor[3*p+1]);
1710ca99263SDave May     coor_p[2] = PetscRealPart(_coor[3*p+2]);
17205264a50SDave May 
17305264a50SDave May     cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND;
17405264a50SDave May 
1750ca99263SDave May     if (coor_p[0] < gmin_l[0]) { continue; }
1760ca99263SDave May     if (coor_p[0] > gmax_l[0]) { continue; }
1770ca99263SDave May     if (coor_p[1] < gmin_l[1]) { continue; }
1780ca99263SDave May     if (coor_p[1] > gmax_l[1]) { continue; }
1790ca99263SDave May     if (coor_p[2] < gmin_l[2]) { continue; }
1800ca99263SDave May     if (coor_p[2] > gmax_l[2]) { continue; }
18105264a50SDave May 
182aab5bcd8SJed Brown     for (d=0; d<3; d++) {
18305264a50SDave May       mi[d] = (PetscInt)((coor_p[d] - gmin[d])/dx[d]);
18405264a50SDave May     }
18505264a50SDave May 
18605264a50SDave May     if (mi[0] < xs)     { continue; }
18705264a50SDave May     if (mi[0] > (xe-1)) { continue; }
18805264a50SDave May     if (mi[1] < ys)     { continue; }
18905264a50SDave May     if (mi[1] > (ye-1)) { continue; }
19005264a50SDave May     if (mi[2] < zs)     { continue; }
19105264a50SDave May     if (mi[2] > (ze-1)) { continue; }
19205264a50SDave May 
19305264a50SDave May     if (mi[0] == (xe-1)) { mi[0]--; }
19405264a50SDave May     if (mi[1] == (ye-1)) { mi[1]--; }
19505264a50SDave May     if (mi[2] == (ze-1)) { mi[2]--; }
19605264a50SDave May 
19705264a50SDave May     cellidx[p] = (mi[0]-xs) + (mi[1]-ys) * mxlocal + (mi[2]-zs) * mxlocal * mylocal;
19805264a50SDave May   }
19905264a50SDave May   ierr = VecRestoreArrayRead(pos,&_coor);CHKERRQ(ierr);
20005264a50SDave May   ierr = ISCreateGeneral(PETSC_COMM_SELF,npoints,cellidx,PETSC_OWN_POINTER,iscell);CHKERRQ(ierr);
20105264a50SDave May   PetscFunctionReturn(0);
20205264a50SDave May }
20305264a50SDave May 
20405264a50SDave May PetscErrorCode DMLocatePoints_DA_Regular(DM dm,Vec pos,DMPointLocationType ltype,PetscSF cellSF)
20505264a50SDave May {
20605264a50SDave May   IS             iscell;
20705264a50SDave May   PetscSFNode    *cells;
20805264a50SDave May   PetscInt       p,bs,dim,npoints,nfound;
20905264a50SDave May   const PetscInt *boxCells;
21005264a50SDave May   PetscErrorCode ierr;
21105264a50SDave May 
21205264a50SDave May   PetscFunctionBegin;
21305264a50SDave May   ierr = VecGetBlockSize(pos,&dim);CHKERRQ(ierr);
21405264a50SDave May   switch (dim) {
21505264a50SDave May     case 1:
21605264a50SDave May       SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Support not provided for 1D");
21705264a50SDave May     case 2:
21805264a50SDave May       ierr = private_DMDALocatePointsIS_2D_Regular(dm,pos,&iscell);CHKERRQ(ierr);
21905264a50SDave May       break;
22005264a50SDave May     case 3:
22105264a50SDave May       ierr = private_DMDALocatePointsIS_3D_Regular(dm,pos,&iscell);CHKERRQ(ierr);
22205264a50SDave May       break;
22305264a50SDave May     default:
22405264a50SDave May       SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unsupport spatial dimension");
22505264a50SDave May   }
22605264a50SDave May 
22705264a50SDave May   ierr = VecGetLocalSize(pos,&npoints);CHKERRQ(ierr);
22805264a50SDave May   ierr = VecGetBlockSize(pos,&bs);CHKERRQ(ierr);
22905264a50SDave May   npoints = npoints / bs;
23005264a50SDave May 
23105264a50SDave May   ierr = PetscMalloc1(npoints, &cells);CHKERRQ(ierr);
23205264a50SDave May   ierr = ISGetIndices(iscell, &boxCells);CHKERRQ(ierr);
23305264a50SDave May 
23405264a50SDave May   for (p=0; p<npoints; p++) {
23505264a50SDave May     cells[p].rank  = 0;
23605264a50SDave May     cells[p].index = boxCells[p];
23705264a50SDave May   }
23805264a50SDave May   ierr = ISRestoreIndices(iscell, &boxCells);CHKERRQ(ierr);
23905264a50SDave May 
24005264a50SDave May   nfound = npoints;
24105264a50SDave May   ierr = PetscSFSetGraph(cellSF, npoints, nfound, NULL, PETSC_OWN_POINTER, cells, PETSC_OWN_POINTER);CHKERRQ(ierr);
24205264a50SDave May   ierr = ISDestroy(&iscell);CHKERRQ(ierr);
24305264a50SDave May   PetscFunctionReturn(0);
24405264a50SDave May }
245