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 @*/ 18d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAConvertToCell(DM dm, MatStencil s, PetscInt *cell) 19d71ae5a4SJacob 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 34d71ae5a4SJacob Faibussowitsch PetscErrorCode private_DMDALocatePointsIS_2D_Regular(DM dmregular, Vec pos, IS *iscell) 35d71ae5a4SJacob 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 ye += ys; 529371c9d4SSatish Balay Ye += Ys; 53*51c42ae6SMatthew G. Knepley if (xs != Xs && Xs >= 0) xs -= 1; 54*51c42ae6SMatthew G. Knepley if (ys != Ys && Ys >= 0) ys -= 1; 5505264a50SDave May 569566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dmregular, &coor)); 579566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(coor, &_coor)); 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]); 66*51c42ae6SMatthew G. Knepley PetscCall(VecRestoreArrayRead(coor, &_coor)); 67*51c42ae6SMatthew G. Knepley 68*51c42ae6SMatthew G. Knepley mxlocal = xe - xs - 1; 69*51c42ae6SMatthew G. Knepley mylocal = ye - ys - 1; 7005264a50SDave May 7105264a50SDave May dx[0] = (gmax_l[0] - gmin_l[0]) / ((PetscReal)mxlocal); 7205264a50SDave May dx[1] = (gmax_l[1] - gmin_l[1]) / ((PetscReal)mylocal); 7305264a50SDave May 749566063dSJacob Faibussowitsch PetscCall(DMGetBoundingBox(dmregular, gmin, gmax)); 7505264a50SDave May 769566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(pos, &n)); 779566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(pos, &bs)); 7805264a50SDave May npoints = n / bs; 7905264a50SDave May 809566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(npoints, &cellidx)); 819566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pos, &_coor)); 8205264a50SDave May for (p = 0; p < npoints; p++) { 830ca99263SDave May PetscReal coor_p[2]; 8405264a50SDave May PetscInt mi[2]; 8505264a50SDave May 860ca99263SDave May coor_p[0] = PetscRealPart(_coor[2 * p]); 870ca99263SDave May coor_p[1] = PetscRealPart(_coor[2 * p + 1]); 8805264a50SDave May 8905264a50SDave May cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND; 9005264a50SDave May 91ad540459SPierre Jolivet if (coor_p[0] < gmin_l[0]) continue; 92ad540459SPierre Jolivet if (coor_p[0] > gmax_l[0]) continue; 93ad540459SPierre Jolivet if (coor_p[1] < gmin_l[1]) continue; 94ad540459SPierre Jolivet if (coor_p[1] > gmax_l[1]) continue; 9505264a50SDave May 96ad540459SPierre Jolivet for (d = 0; d < 2; d++) mi[d] = (PetscInt)((coor_p[d] - gmin[d]) / dx[d]); 9705264a50SDave May 98ad540459SPierre Jolivet if (mi[0] < xs) continue; 99ad540459SPierre Jolivet if (mi[0] > (xe - 1)) continue; 100ad540459SPierre Jolivet if (mi[1] < ys) continue; 101ad540459SPierre Jolivet if (mi[1] > (ye - 1)) continue; 10205264a50SDave May 103ad540459SPierre Jolivet if (mi[0] == (xe - 1)) mi[0]--; 104ad540459SPierre Jolivet if (mi[1] == (ye - 1)) mi[1]--; 10505264a50SDave May 10605264a50SDave May cellidx[p] = (mi[0] - xs) + (mi[1] - ys) * mxlocal; 10705264a50SDave May } 1089566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pos, &_coor)); 1099566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, npoints, cellidx, PETSC_OWN_POINTER, iscell)); 11005264a50SDave May PetscFunctionReturn(0); 11105264a50SDave May } 11205264a50SDave May 113d71ae5a4SJacob Faibussowitsch PetscErrorCode private_DMDALocatePointsIS_3D_Regular(DM dmregular, Vec pos, IS *iscell) 114d71ae5a4SJacob Faibussowitsch { 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 126aab5bcd8SJed Brown PetscFunctionBegin; 1279566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(dmregular, &xs, &ys, &zs, &xe, &ye, &ze)); 1289566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(dmregular, &Xs, &Ys, &Zs, &Xe, &Ye, &Ze)); 1299371c9d4SSatish Balay xe += xs; 1309371c9d4SSatish Balay Xe += Xs; 1319371c9d4SSatish Balay ye += ys; 1329371c9d4SSatish Balay Ye += Ys; 1339371c9d4SSatish Balay ze += zs; 1349371c9d4SSatish Balay Ze += Zs; 135*51c42ae6SMatthew G. Knepley if (xs != Xs && Xs >= 0) xs -= 1; 136*51c42ae6SMatthew G. Knepley if (ys != Ys && Ys >= 0) ys -= 1; 137*51c42ae6SMatthew G. Knepley if (zs != Zs && Zs >= 0) zs -= 1; 13805264a50SDave May 1399566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dmregular, &coor)); 1409566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(coor, &_coor)); 14105264a50SDave May c0 = (xs - Xs) + (ys - Ys) * (Xe - Xs) + (zs - Zs) * (Xe - Xs) * (Ye - Ys); 14205264a50SDave May c1 = (xe - 2 - Xs + 1) + (ye - 2 - Ys + 1) * (Xe - Xs) + (ze - 2 - Zs + 1) * (Xe - Xs) * (Ye - Ys); 14305264a50SDave May 144a5f152d1SDave May gmin_l[0] = PetscRealPart(_coor[3 * c0 + 0]); 145a5f152d1SDave May gmin_l[1] = PetscRealPart(_coor[3 * c0 + 1]); 146a5f152d1SDave May gmin_l[2] = PetscRealPart(_coor[3 * c0 + 2]); 14705264a50SDave May 148a5f152d1SDave May gmax_l[0] = PetscRealPart(_coor[3 * c1 + 0]); 149a5f152d1SDave May gmax_l[1] = PetscRealPart(_coor[3 * c1 + 1]); 150a5f152d1SDave May gmax_l[2] = PetscRealPart(_coor[3 * c1 + 2]); 151*51c42ae6SMatthew G. Knepley PetscCall(VecRestoreArrayRead(coor, &_coor)); 152*51c42ae6SMatthew G. Knepley 153*51c42ae6SMatthew G. Knepley if (xs != Xs) xs -= 1; 154*51c42ae6SMatthew G. Knepley if (ys != Ys) ys -= 1; 155*51c42ae6SMatthew G. Knepley if (zs != Zs) zs -= 1; 156*51c42ae6SMatthew G. Knepley 157*51c42ae6SMatthew G. Knepley mxlocal = xe - xs - 1; 158*51c42ae6SMatthew G. Knepley mylocal = ye - ys - 1; 159*51c42ae6SMatthew G. Knepley mzlocal = ze - zs - 1; 16005264a50SDave May 16105264a50SDave May dx[0] = (gmax_l[0] - gmin_l[0]) / ((PetscReal)mxlocal); 16205264a50SDave May dx[1] = (gmax_l[1] - gmin_l[1]) / ((PetscReal)mylocal); 16305264a50SDave May dx[2] = (gmax_l[2] - gmin_l[2]) / ((PetscReal)mzlocal); 16405264a50SDave May 1659566063dSJacob Faibussowitsch PetscCall(DMGetBoundingBox(dmregular, gmin, gmax)); 16605264a50SDave May 1679566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(pos, &n)); 1689566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(pos, &bs)); 16905264a50SDave May npoints = n / bs; 17005264a50SDave May 1719566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(npoints, &cellidx)); 1729566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(pos, &_coor)); 17305264a50SDave May for (p = 0; p < npoints; p++) { 1740ca99263SDave May PetscReal coor_p[3]; 17505264a50SDave May PetscInt mi[3]; 17605264a50SDave May 1770ca99263SDave May coor_p[0] = PetscRealPart(_coor[3 * p]); 1780ca99263SDave May coor_p[1] = PetscRealPart(_coor[3 * p + 1]); 1790ca99263SDave May coor_p[2] = PetscRealPart(_coor[3 * p + 2]); 18005264a50SDave May 18105264a50SDave May cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND; 18205264a50SDave May 183ad540459SPierre Jolivet if (coor_p[0] < gmin_l[0]) continue; 184ad540459SPierre Jolivet if (coor_p[0] > gmax_l[0]) continue; 185ad540459SPierre Jolivet if (coor_p[1] < gmin_l[1]) continue; 186ad540459SPierre Jolivet if (coor_p[1] > gmax_l[1]) continue; 187ad540459SPierre Jolivet if (coor_p[2] < gmin_l[2]) continue; 188ad540459SPierre Jolivet if (coor_p[2] > gmax_l[2]) continue; 18905264a50SDave May 190ad540459SPierre Jolivet for (d = 0; d < 3; d++) mi[d] = (PetscInt)((coor_p[d] - gmin[d]) / dx[d]); 19105264a50SDave May 192ad540459SPierre Jolivet if (mi[0] < xs) continue; 193ad540459SPierre Jolivet if (mi[0] > (xe - 1)) continue; 194ad540459SPierre Jolivet if (mi[1] < ys) continue; 195ad540459SPierre Jolivet if (mi[1] > (ye - 1)) continue; 196ad540459SPierre Jolivet if (mi[2] < zs) continue; 197ad540459SPierre Jolivet if (mi[2] > (ze - 1)) continue; 19805264a50SDave May 199ad540459SPierre Jolivet if (mi[0] == (xe - 1)) mi[0]--; 200ad540459SPierre Jolivet if (mi[1] == (ye - 1)) mi[1]--; 201ad540459SPierre Jolivet if (mi[2] == (ze - 1)) mi[2]--; 20205264a50SDave May 20305264a50SDave May cellidx[p] = (mi[0] - xs) + (mi[1] - ys) * mxlocal + (mi[2] - zs) * mxlocal * mylocal; 20405264a50SDave May } 2059566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(pos, &_coor)); 2069566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, npoints, cellidx, PETSC_OWN_POINTER, iscell)); 20705264a50SDave May PetscFunctionReturn(0); 20805264a50SDave May } 20905264a50SDave May 210d71ae5a4SJacob Faibussowitsch PetscErrorCode DMLocatePoints_DA_Regular(DM dm, Vec pos, DMPointLocationType ltype, PetscSF cellSF) 211d71ae5a4SJacob Faibussowitsch { 21205264a50SDave May IS iscell; 21305264a50SDave May PetscSFNode *cells; 21405264a50SDave May PetscInt p, bs, dim, npoints, nfound; 21505264a50SDave May const PetscInt *boxCells; 21605264a50SDave May 21705264a50SDave May PetscFunctionBegin; 2189566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(pos, &dim)); 21905264a50SDave May switch (dim) { 220d71ae5a4SJacob Faibussowitsch case 1: 221d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Support not provided for 1D"); 222d71ae5a4SJacob Faibussowitsch case 2: 223d71ae5a4SJacob Faibussowitsch PetscCall(private_DMDALocatePointsIS_2D_Regular(dm, pos, &iscell)); 224d71ae5a4SJacob Faibussowitsch break; 225d71ae5a4SJacob Faibussowitsch case 3: 226d71ae5a4SJacob Faibussowitsch PetscCall(private_DMDALocatePointsIS_3D_Regular(dm, pos, &iscell)); 227d71ae5a4SJacob Faibussowitsch break; 228d71ae5a4SJacob Faibussowitsch default: 229d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupport spatial dimension"); 23005264a50SDave May } 23105264a50SDave May 2329566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(pos, &npoints)); 2339566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(pos, &bs)); 23405264a50SDave May npoints = npoints / bs; 23505264a50SDave May 2369566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(npoints, &cells)); 2379566063dSJacob Faibussowitsch PetscCall(ISGetIndices(iscell, &boxCells)); 23805264a50SDave May 23905264a50SDave May for (p = 0; p < npoints; p++) { 24005264a50SDave May cells[p].rank = 0; 24105264a50SDave May cells[p].index = boxCells[p]; 24205264a50SDave May } 2439566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(iscell, &boxCells)); 24405264a50SDave May 24505264a50SDave May nfound = npoints; 2469566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(cellSF, npoints, nfound, NULL, PETSC_OWN_POINTER, cells, PETSC_OWN_POINTER)); 2479566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscell)); 24805264a50SDave May PetscFunctionReturn(0); 24905264a50SDave May } 250