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 @*/ 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; 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); 28*1dca8a05SBarry 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); 29*1dca8a05SBarry 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 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 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)); 4905264a50SDave May xe += xs; Xe += Xs; if (xs != Xs) xs -= 1; 5005264a50SDave May ye += ys; Ye += Ys; if (ys != Ys) ys -= 1; 5105264a50SDave May 5205264a50SDave May mxlocal = xe - xs - 1; 5305264a50SDave May mylocal = ye - ys - 1; 5405264a50SDave May 559566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dmregular,&coor)); 569566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(coor,&_coor)); 5705264a50SDave May c0 = (xs-Xs) + (ys-Ys)*(Xe-Xs); 5805264a50SDave May c1 = (xe-2-Xs+1) + (ye-2-Ys+1)*(Xe-Xs); 5905264a50SDave May 60a5f152d1SDave May gmin_l[0] = PetscRealPart(_coor[2*c0+0]); 61a5f152d1SDave May gmin_l[1] = PetscRealPart(_coor[2*c0+1]); 6205264a50SDave May 63a5f152d1SDave May gmax_l[0] = PetscRealPart(_coor[2*c1+0]); 64a5f152d1SDave May gmax_l[1] = PetscRealPart(_coor[2*c1+1]); 6505264a50SDave May 6605264a50SDave May dx[0] = (gmax_l[0]-gmin_l[0])/((PetscReal)mxlocal); 6705264a50SDave May dx[1] = (gmax_l[1]-gmin_l[1])/((PetscReal)mylocal); 6805264a50SDave May 699566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(coor,&_coor)); 7005264a50SDave May 719566063dSJacob Faibussowitsch PetscCall(DMGetBoundingBox(dmregular,gmin,gmax)); 7205264a50SDave May 739566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(pos,&n)); 749566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(pos,&bs)); 7505264a50SDave May npoints = n/bs; 7605264a50SDave May 779566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(npoints,&cellidx)); 789566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pos,&_coor)); 7905264a50SDave May for (p=0; p<npoints; p++) { 800ca99263SDave May PetscReal coor_p[2]; 8105264a50SDave May PetscInt mi[2]; 8205264a50SDave May 830ca99263SDave May coor_p[0] = PetscRealPart(_coor[2*p]); 840ca99263SDave May coor_p[1] = PetscRealPart(_coor[2*p+1]); 8505264a50SDave May 8605264a50SDave May cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND; 8705264a50SDave May 880ca99263SDave May if (coor_p[0] < gmin_l[0]) { continue; } 890ca99263SDave May if (coor_p[0] > gmax_l[0]) { continue; } 900ca99263SDave May if (coor_p[1] < gmin_l[1]) { continue; } 910ca99263SDave May if (coor_p[1] > gmax_l[1]) { continue; } 9205264a50SDave May 93aab5bcd8SJed Brown for (d=0; d<2; d++) { 9405264a50SDave May mi[d] = (PetscInt)((coor_p[d] - gmin[d])/dx[d]); 9505264a50SDave May } 9605264a50SDave May 9705264a50SDave May if (mi[0] < xs) { continue; } 9805264a50SDave May if (mi[0] > (xe-1)) { continue; } 9905264a50SDave May if (mi[1] < ys) { continue; } 10005264a50SDave May if (mi[1] > (ye-1)) { continue; } 10105264a50SDave May 10205264a50SDave May if (mi[0] == (xe-1)) { mi[0]--; } 10305264a50SDave May if (mi[1] == (ye-1)) { mi[1]--; } 10405264a50SDave May 10505264a50SDave May cellidx[p] = (mi[0]-xs) + (mi[1]-ys) * mxlocal; 10605264a50SDave May } 1079566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pos,&_coor)); 1089566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,npoints,cellidx,PETSC_OWN_POINTER,iscell)); 10905264a50SDave May PetscFunctionReturn(0); 11005264a50SDave May } 11105264a50SDave May 11205264a50SDave May PetscErrorCode private_DMDALocatePointsIS_3D_Regular(DM dmregular,Vec pos,IS *iscell) 11305264a50SDave May { 11405264a50SDave May PetscInt n,bs,p,npoints; 11505264a50SDave May PetscInt xs,xe,Xs,Xe,mxlocal; 11605264a50SDave May PetscInt ys,ye,Ys,Ye,mylocal; 11705264a50SDave May PetscInt zs,ze,Zs,Ze,mzlocal; 118aab5bcd8SJed Brown PetscInt d,c0,c1; 11905264a50SDave May PetscReal gmin_l[3],gmax_l[3],dx[3]; 12005264a50SDave May PetscReal gmin[3],gmax[3]; 12105264a50SDave May PetscInt *cellidx; 12205264a50SDave May Vec coor; 12305264a50SDave May const PetscScalar *_coor; 12405264a50SDave May 125aab5bcd8SJed Brown PetscFunctionBegin; 1269566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dmregular,&xs,&ys,&zs,&xe,&ye,&ze)); 1279566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(dmregular,&Xs,&Ys,&Zs,&Xe,&Ye,&Ze)); 12805264a50SDave May xe += xs; Xe += Xs; if (xs != Xs) xs -= 1; 12905264a50SDave May ye += ys; Ye += Ys; if (ys != Ys) ys -= 1; 13005264a50SDave May ze += zs; Ze += Zs; if (zs != Zs) zs -= 1; 13105264a50SDave May 13205264a50SDave May mxlocal = xe - xs - 1; 13305264a50SDave May mylocal = ye - ys - 1; 13405264a50SDave May mzlocal = ze - zs - 1; 13505264a50SDave May 1369566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dmregular,&coor)); 1379566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(coor,&_coor)); 13805264a50SDave May c0 = (xs-Xs) + (ys-Ys) *(Xe-Xs) + (zs-Zs) *(Xe-Xs)*(Ye-Ys); 13905264a50SDave May c1 = (xe-2-Xs+1) + (ye-2-Ys+1)*(Xe-Xs) + (ze-2-Zs+1)*(Xe-Xs)*(Ye-Ys); 14005264a50SDave May 141a5f152d1SDave May gmin_l[0] = PetscRealPart(_coor[3*c0+0]); 142a5f152d1SDave May gmin_l[1] = PetscRealPart(_coor[3*c0+1]); 143a5f152d1SDave May gmin_l[2] = PetscRealPart(_coor[3*c0+2]); 14405264a50SDave May 145a5f152d1SDave May gmax_l[0] = PetscRealPart(_coor[3*c1+0]); 146a5f152d1SDave May gmax_l[1] = PetscRealPart(_coor[3*c1+1]); 147a5f152d1SDave May gmax_l[2] = PetscRealPart(_coor[3*c1+2]); 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 dx[2] = (gmax_l[2]-gmin_l[2])/((PetscReal)mzlocal); 15205264a50SDave May 1539566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(coor,&_coor)); 15405264a50SDave May 1559566063dSJacob Faibussowitsch PetscCall(DMGetBoundingBox(dmregular,gmin,gmax)); 15605264a50SDave May 1579566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(pos,&n)); 1589566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(pos,&bs)); 15905264a50SDave May npoints = n/bs; 16005264a50SDave May 1619566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(npoints,&cellidx)); 1629566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pos,&_coor)); 16305264a50SDave May for (p=0; p<npoints; p++) { 1640ca99263SDave May PetscReal coor_p[3]; 16505264a50SDave May PetscInt mi[3]; 16605264a50SDave May 1670ca99263SDave May coor_p[0] = PetscRealPart(_coor[3*p]); 1680ca99263SDave May coor_p[1] = PetscRealPart(_coor[3*p+1]); 1690ca99263SDave May coor_p[2] = PetscRealPart(_coor[3*p+2]); 17005264a50SDave May 17105264a50SDave May cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND; 17205264a50SDave May 1730ca99263SDave May if (coor_p[0] < gmin_l[0]) { continue; } 1740ca99263SDave May if (coor_p[0] > gmax_l[0]) { continue; } 1750ca99263SDave May if (coor_p[1] < gmin_l[1]) { continue; } 1760ca99263SDave May if (coor_p[1] > gmax_l[1]) { continue; } 1770ca99263SDave May if (coor_p[2] < gmin_l[2]) { continue; } 1780ca99263SDave May if (coor_p[2] > gmax_l[2]) { continue; } 17905264a50SDave May 180aab5bcd8SJed Brown for (d=0; d<3; d++) { 18105264a50SDave May mi[d] = (PetscInt)((coor_p[d] - gmin[d])/dx[d]); 18205264a50SDave May } 18305264a50SDave May 18405264a50SDave May if (mi[0] < xs) { continue; } 18505264a50SDave May if (mi[0] > (xe-1)) { continue; } 18605264a50SDave May if (mi[1] < ys) { continue; } 18705264a50SDave May if (mi[1] > (ye-1)) { continue; } 18805264a50SDave May if (mi[2] < zs) { continue; } 18905264a50SDave May if (mi[2] > (ze-1)) { continue; } 19005264a50SDave May 19105264a50SDave May if (mi[0] == (xe-1)) { mi[0]--; } 19205264a50SDave May if (mi[1] == (ye-1)) { mi[1]--; } 19305264a50SDave May if (mi[2] == (ze-1)) { mi[2]--; } 19405264a50SDave May 19505264a50SDave May cellidx[p] = (mi[0]-xs) + (mi[1]-ys) * mxlocal + (mi[2]-zs) * mxlocal * mylocal; 19605264a50SDave May } 1979566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pos,&_coor)); 1989566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,npoints,cellidx,PETSC_OWN_POINTER,iscell)); 19905264a50SDave May PetscFunctionReturn(0); 20005264a50SDave May } 20105264a50SDave May 20205264a50SDave May PetscErrorCode DMLocatePoints_DA_Regular(DM dm,Vec pos,DMPointLocationType ltype,PetscSF cellSF) 20305264a50SDave May { 20405264a50SDave May IS iscell; 20505264a50SDave May PetscSFNode *cells; 20605264a50SDave May PetscInt p,bs,dim,npoints,nfound; 20705264a50SDave May const PetscInt *boxCells; 20805264a50SDave May 20905264a50SDave May PetscFunctionBegin; 2109566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(pos,&dim)); 21105264a50SDave May switch (dim) { 21205264a50SDave May case 1: 21305264a50SDave May SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Support not provided for 1D"); 21405264a50SDave May case 2: 2159566063dSJacob Faibussowitsch PetscCall(private_DMDALocatePointsIS_2D_Regular(dm,pos,&iscell)); 21605264a50SDave May break; 21705264a50SDave May case 3: 2189566063dSJacob Faibussowitsch PetscCall(private_DMDALocatePointsIS_3D_Regular(dm,pos,&iscell)); 21905264a50SDave May break; 22005264a50SDave May default: 22105264a50SDave May SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Unsupport spatial dimension"); 22205264a50SDave May } 22305264a50SDave May 2249566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(pos,&npoints)); 2259566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(pos,&bs)); 22605264a50SDave May npoints = npoints / bs; 22705264a50SDave May 2289566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(npoints, &cells)); 2299566063dSJacob Faibussowitsch PetscCall(ISGetIndices(iscell, &boxCells)); 23005264a50SDave May 23105264a50SDave May for (p=0; p<npoints; p++) { 23205264a50SDave May cells[p].rank = 0; 23305264a50SDave May cells[p].index = boxCells[p]; 23405264a50SDave May } 2359566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(iscell, &boxCells)); 23605264a50SDave May 23705264a50SDave May nfound = npoints; 2389566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(cellSF, npoints, nfound, NULL, PETSC_OWN_POINTER, cells, PETSC_OWN_POINTER)); 2399566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscell)); 24005264a50SDave May PetscFunctionReturn(0); 24105264a50SDave May } 242