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*9371c9d4SSatish Balay PetscErrorCode DMDAConvertToCell(DM dm, MatStencil s, PetscInt *cell) { 19341c9bdaSMatthew G Knepley DM_DA *da = (DM_DA *)dm->data; 20c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 214e846693SMatthew G Knepley const PetscInt mx = (da->Xe - da->Xs) / da->w, my = da->Ye - da->Ys /*, mz = da->Ze - da->Zs*/; 22ed4e918cSMatthew 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; 23341c9bdaSMatthew G Knepley 24341c9bdaSMatthew G Knepley PetscFunctionBegin; 25ed4e918cSMatthew G Knepley *cell = -1; 2663a3b9bcSJacob 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); 271dca8a05SBarry 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); 281dca8a05SBarry 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); 29ed4e918cSMatthew G Knepley *cell = (kl * my + jl) * mx + il; 30ed4e918cSMatthew G Knepley PetscFunctionReturn(0); 31341c9bdaSMatthew G Knepley } 32341c9bdaSMatthew G Knepley 33*9371c9d4SSatish Balay PetscErrorCode private_DMDALocatePointsIS_2D_Regular(DM dmregular, Vec pos, IS *iscell) { 3405264a50SDave May PetscInt n, bs, p, npoints; 3505264a50SDave May PetscInt xs, xe, Xs, Xe, mxlocal; 3605264a50SDave May PetscInt ys, ye, Ys, Ye, mylocal; 37aab5bcd8SJed Brown PetscInt d, c0, c1; 3805264a50SDave May PetscReal gmin_l[2], gmax_l[2], dx[2]; 3905264a50SDave May PetscReal gmin[2], gmax[2]; 4005264a50SDave May PetscInt *cellidx; 4105264a50SDave May Vec coor; 4205264a50SDave May const PetscScalar *_coor; 4305264a50SDave May 44aab5bcd8SJed Brown PetscFunctionBegin; 459566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dmregular, &xs, &ys, NULL, &xe, &ye, NULL)); 469566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(dmregular, &Xs, &Ys, NULL, &Xe, &Ye, NULL)); 47*9371c9d4SSatish Balay xe += xs; 48*9371c9d4SSatish Balay Xe += Xs; 49*9371c9d4SSatish Balay if (xs != Xs) xs -= 1; 50*9371c9d4SSatish Balay ye += ys; 51*9371c9d4SSatish Balay Ye += Ys; 52*9371c9d4SSatish Balay if (ys != Ys) ys -= 1; 5305264a50SDave May 5405264a50SDave May mxlocal = xe - xs - 1; 5505264a50SDave May mylocal = ye - ys - 1; 5605264a50SDave May 579566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dmregular, &coor)); 589566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(coor, &_coor)); 5905264a50SDave May c0 = (xs - Xs) + (ys - Ys) * (Xe - Xs); 6005264a50SDave May c1 = (xe - 2 - Xs + 1) + (ye - 2 - Ys + 1) * (Xe - Xs); 6105264a50SDave May 62a5f152d1SDave May gmin_l[0] = PetscRealPart(_coor[2 * c0 + 0]); 63a5f152d1SDave May gmin_l[1] = PetscRealPart(_coor[2 * c0 + 1]); 6405264a50SDave May 65a5f152d1SDave May gmax_l[0] = PetscRealPart(_coor[2 * c1 + 0]); 66a5f152d1SDave May gmax_l[1] = PetscRealPart(_coor[2 * c1 + 1]); 6705264a50SDave May 6805264a50SDave May dx[0] = (gmax_l[0] - gmin_l[0]) / ((PetscReal)mxlocal); 6905264a50SDave May dx[1] = (gmax_l[1] - gmin_l[1]) / ((PetscReal)mylocal); 7005264a50SDave May 719566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(coor, &_coor)); 7205264a50SDave May 739566063dSJacob Faibussowitsch PetscCall(DMGetBoundingBox(dmregular, gmin, gmax)); 7405264a50SDave May 759566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(pos, &n)); 769566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(pos, &bs)); 7705264a50SDave May npoints = n / bs; 7805264a50SDave May 799566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(npoints, &cellidx)); 809566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pos, &_coor)); 8105264a50SDave May for (p = 0; p < npoints; p++) { 820ca99263SDave May PetscReal coor_p[2]; 8305264a50SDave May PetscInt mi[2]; 8405264a50SDave May 850ca99263SDave May coor_p[0] = PetscRealPart(_coor[2 * p]); 860ca99263SDave May coor_p[1] = PetscRealPart(_coor[2 * p + 1]); 8705264a50SDave May 8805264a50SDave May cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND; 8905264a50SDave May 900ca99263SDave May if (coor_p[0] < gmin_l[0]) { continue; } 910ca99263SDave May if (coor_p[0] > gmax_l[0]) { continue; } 920ca99263SDave May if (coor_p[1] < gmin_l[1]) { continue; } 930ca99263SDave May if (coor_p[1] > gmax_l[1]) { continue; } 9405264a50SDave May 95*9371c9d4SSatish Balay for (d = 0; d < 2; d++) { mi[d] = (PetscInt)((coor_p[d] - gmin[d]) / dx[d]); } 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 112*9371c9d4SSatish Balay PetscErrorCode private_DMDALocatePointsIS_3D_Regular(DM dmregular, Vec pos, IS *iscell) { 11305264a50SDave May PetscInt n, bs, p, npoints; 11405264a50SDave May PetscInt xs, xe, Xs, Xe, mxlocal; 11505264a50SDave May PetscInt ys, ye, Ys, Ye, mylocal; 11605264a50SDave May PetscInt zs, ze, Zs, Ze, mzlocal; 117aab5bcd8SJed Brown PetscInt d, c0, c1; 11805264a50SDave May PetscReal gmin_l[3], gmax_l[3], dx[3]; 11905264a50SDave May PetscReal gmin[3], gmax[3]; 12005264a50SDave May PetscInt *cellidx; 12105264a50SDave May Vec coor; 12205264a50SDave May const PetscScalar *_coor; 12305264a50SDave May 124aab5bcd8SJed Brown PetscFunctionBegin; 1259566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dmregular, &xs, &ys, &zs, &xe, &ye, &ze)); 1269566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(dmregular, &Xs, &Ys, &Zs, &Xe, &Ye, &Ze)); 127*9371c9d4SSatish Balay xe += xs; 128*9371c9d4SSatish Balay Xe += Xs; 129*9371c9d4SSatish Balay if (xs != Xs) xs -= 1; 130*9371c9d4SSatish Balay ye += ys; 131*9371c9d4SSatish Balay Ye += Ys; 132*9371c9d4SSatish Balay if (ys != Ys) ys -= 1; 133*9371c9d4SSatish Balay ze += zs; 134*9371c9d4SSatish Balay Ze += Zs; 135*9371c9d4SSatish Balay if (zs != Zs) zs -= 1; 13605264a50SDave May 13705264a50SDave May mxlocal = xe - xs - 1; 13805264a50SDave May mylocal = ye - ys - 1; 13905264a50SDave May mzlocal = ze - zs - 1; 14005264a50SDave May 1419566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dmregular, &coor)); 1429566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(coor, &_coor)); 14305264a50SDave May c0 = (xs - Xs) + (ys - Ys) * (Xe - Xs) + (zs - Zs) * (Xe - Xs) * (Ye - Ys); 14405264a50SDave May c1 = (xe - 2 - Xs + 1) + (ye - 2 - Ys + 1) * (Xe - Xs) + (ze - 2 - Zs + 1) * (Xe - Xs) * (Ye - Ys); 14505264a50SDave May 146a5f152d1SDave May gmin_l[0] = PetscRealPart(_coor[3 * c0 + 0]); 147a5f152d1SDave May gmin_l[1] = PetscRealPart(_coor[3 * c0 + 1]); 148a5f152d1SDave May gmin_l[2] = PetscRealPart(_coor[3 * c0 + 2]); 14905264a50SDave May 150a5f152d1SDave May gmax_l[0] = PetscRealPart(_coor[3 * c1 + 0]); 151a5f152d1SDave May gmax_l[1] = PetscRealPart(_coor[3 * c1 + 1]); 152a5f152d1SDave May gmax_l[2] = PetscRealPart(_coor[3 * c1 + 2]); 15305264a50SDave May 15405264a50SDave May dx[0] = (gmax_l[0] - gmin_l[0]) / ((PetscReal)mxlocal); 15505264a50SDave May dx[1] = (gmax_l[1] - gmin_l[1]) / ((PetscReal)mylocal); 15605264a50SDave May dx[2] = (gmax_l[2] - gmin_l[2]) / ((PetscReal)mzlocal); 15705264a50SDave May 1589566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(coor, &_coor)); 15905264a50SDave May 1609566063dSJacob Faibussowitsch PetscCall(DMGetBoundingBox(dmregular, gmin, gmax)); 16105264a50SDave May 1629566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(pos, &n)); 1639566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(pos, &bs)); 16405264a50SDave May npoints = n / bs; 16505264a50SDave May 1669566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(npoints, &cellidx)); 1679566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pos, &_coor)); 16805264a50SDave May for (p = 0; p < npoints; p++) { 1690ca99263SDave May PetscReal coor_p[3]; 17005264a50SDave May PetscInt mi[3]; 17105264a50SDave May 1720ca99263SDave May coor_p[0] = PetscRealPart(_coor[3 * p]); 1730ca99263SDave May coor_p[1] = PetscRealPart(_coor[3 * p + 1]); 1740ca99263SDave May coor_p[2] = PetscRealPart(_coor[3 * p + 2]); 17505264a50SDave May 17605264a50SDave May cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND; 17705264a50SDave May 1780ca99263SDave May if (coor_p[0] < gmin_l[0]) { continue; } 1790ca99263SDave May if (coor_p[0] > gmax_l[0]) { continue; } 1800ca99263SDave May if (coor_p[1] < gmin_l[1]) { continue; } 1810ca99263SDave May if (coor_p[1] > gmax_l[1]) { continue; } 1820ca99263SDave May if (coor_p[2] < gmin_l[2]) { continue; } 1830ca99263SDave May if (coor_p[2] > gmax_l[2]) { continue; } 18405264a50SDave May 185*9371c9d4SSatish Balay for (d = 0; d < 3; d++) { mi[d] = (PetscInt)((coor_p[d] - gmin[d]) / dx[d]); } 18605264a50SDave May 18705264a50SDave May if (mi[0] < xs) { continue; } 18805264a50SDave May if (mi[0] > (xe - 1)) { continue; } 18905264a50SDave May if (mi[1] < ys) { continue; } 19005264a50SDave May if (mi[1] > (ye - 1)) { continue; } 19105264a50SDave May if (mi[2] < zs) { continue; } 19205264a50SDave May if (mi[2] > (ze - 1)) { continue; } 19305264a50SDave May 19405264a50SDave May if (mi[0] == (xe - 1)) { mi[0]--; } 19505264a50SDave May if (mi[1] == (ye - 1)) { mi[1]--; } 19605264a50SDave May if (mi[2] == (ze - 1)) { mi[2]--; } 19705264a50SDave May 19805264a50SDave May cellidx[p] = (mi[0] - xs) + (mi[1] - ys) * mxlocal + (mi[2] - zs) * mxlocal * mylocal; 19905264a50SDave May } 2009566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pos, &_coor)); 2019566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, npoints, cellidx, PETSC_OWN_POINTER, iscell)); 20205264a50SDave May PetscFunctionReturn(0); 20305264a50SDave May } 20405264a50SDave May 205*9371c9d4SSatish Balay PetscErrorCode DMLocatePoints_DA_Regular(DM dm, Vec pos, DMPointLocationType ltype, PetscSF cellSF) { 20605264a50SDave May IS iscell; 20705264a50SDave May PetscSFNode *cells; 20805264a50SDave May PetscInt p, bs, dim, npoints, nfound; 20905264a50SDave May const PetscInt *boxCells; 21005264a50SDave May 21105264a50SDave May PetscFunctionBegin; 2129566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(pos, &dim)); 21305264a50SDave May switch (dim) { 214*9371c9d4SSatish Balay case 1: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Support not provided for 1D"); 215*9371c9d4SSatish Balay case 2: PetscCall(private_DMDALocatePointsIS_2D_Regular(dm, pos, &iscell)); break; 216*9371c9d4SSatish Balay case 3: PetscCall(private_DMDALocatePointsIS_3D_Regular(dm, pos, &iscell)); break; 217*9371c9d4SSatish Balay default: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupport spatial dimension"); 21805264a50SDave May } 21905264a50SDave May 2209566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(pos, &npoints)); 2219566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(pos, &bs)); 22205264a50SDave May npoints = npoints / bs; 22305264a50SDave May 2249566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(npoints, &cells)); 2259566063dSJacob Faibussowitsch PetscCall(ISGetIndices(iscell, &boxCells)); 22605264a50SDave May 22705264a50SDave May for (p = 0; p < npoints; p++) { 22805264a50SDave May cells[p].rank = 0; 22905264a50SDave May cells[p].index = boxCells[p]; 23005264a50SDave May } 2319566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(iscell, &boxCells)); 23205264a50SDave May 23305264a50SDave May nfound = npoints; 2349566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(cellSF, npoints, nfound, NULL, PETSC_OWN_POINTER, cells, PETSC_OWN_POINTER)); 2359566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscell)); 23605264a50SDave May PetscFunctionReturn(0); 23705264a50SDave May } 238