147c6ae99SBarry Smith 247c6ae99SBarry Smith /* 347c6ae99SBarry Smith Code for manipulating distributed regular arrays in parallel. 447c6ae99SBarry Smith */ 547c6ae99SBarry Smith 6af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/ 7b2fff234SMatthew G. Knepley #include <petscbt.h> 80c312b8eSJed Brown #include <petscsf.h> 92764a2aaSMatthew G. Knepley #include <petscds.h> 109a800dd8SMatthew G. Knepley #include <petscfe.h> 1147c6ae99SBarry Smith 1247c6ae99SBarry Smith /* 13e3c5b3baSBarry Smith This allows the DMDA vectors to properly tell MATLAB their dimensions 1447c6ae99SBarry Smith */ 1547c6ae99SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 16c6db04a5SJed Brown #include <engine.h> /* MATLAB include file */ 17c6db04a5SJed Brown #include <mex.h> /* MATLAB include file */ 189371c9d4SSatish Balay static PetscErrorCode VecMatlabEnginePut_DA2d(PetscObject obj, void *mengine) { 1947c6ae99SBarry Smith PetscInt n, m; 2047c6ae99SBarry Smith Vec vec = (Vec)obj; 2147c6ae99SBarry Smith PetscScalar *array; 2247c6ae99SBarry Smith mxArray *mat; 239a42bb27SBarry Smith DM da; 2447c6ae99SBarry Smith 2547c6ae99SBarry Smith PetscFunctionBegin; 269566063dSJacob Faibussowitsch PetscCall(VecGetDM(vec, &da)); 277a8be351SBarry Smith PetscCheck(da, PetscObjectComm((PetscObject)vec), PETSC_ERR_ARG_WRONGSTATE, "Vector not associated with a DMDA"); 289566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, 0, 0, 0, &m, &n, 0)); 2947c6ae99SBarry Smith 309566063dSJacob Faibussowitsch PetscCall(VecGetArray(vec, &array)); 3147c6ae99SBarry Smith #if !defined(PETSC_USE_COMPLEX) 3247c6ae99SBarry Smith mat = mxCreateDoubleMatrix(m, n, mxREAL); 3347c6ae99SBarry Smith #else 3447c6ae99SBarry Smith mat = mxCreateDoubleMatrix(m, n, mxCOMPLEX); 3547c6ae99SBarry Smith #endif 369566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(mxGetPr(mat), array, n * m)); 379566063dSJacob Faibussowitsch PetscCall(PetscObjectName(obj)); 3847c6ae99SBarry Smith engPutVariable((Engine *)mengine, obj->name, mat); 3947c6ae99SBarry Smith 409566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(vec, &array)); 4147c6ae99SBarry Smith PetscFunctionReturn(0); 4247c6ae99SBarry Smith } 4347c6ae99SBarry Smith #endif 4447c6ae99SBarry Smith 459371c9d4SSatish Balay PetscErrorCode DMCreateLocalVector_DA(DM da, Vec *g) { 4647c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 4747c6ae99SBarry Smith 4847c6ae99SBarry Smith PetscFunctionBegin; 4947c6ae99SBarry Smith PetscValidHeaderSpecific(da, DM_CLASSID, 1); 5047c6ae99SBarry Smith PetscValidPointer(g, 2); 519566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, g)); 529566063dSJacob Faibussowitsch PetscCall(VecSetSizes(*g, dd->nlocal, PETSC_DETERMINE)); 539566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(*g, dd->w)); 549566063dSJacob Faibussowitsch PetscCall(VecSetType(*g, da->vectype)); 5574427ab1SRichard Tran Mills if (dd->nlocal < da->bind_below) { 569566063dSJacob Faibussowitsch PetscCall(VecSetBindingPropagates(*g, PETSC_TRUE)); 579566063dSJacob Faibussowitsch PetscCall(VecBindToCPU(*g, PETSC_TRUE)); 5874427ab1SRichard Tran Mills } 599566063dSJacob Faibussowitsch PetscCall(VecSetDM(*g, da)); 6047c6ae99SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE) 61*48a46eb9SPierre Jolivet if (dd->w == 1 && da->dim == 2) PetscCall(PetscObjectComposeFunction((PetscObject)*g, "PetscMatlabEnginePut_C", VecMatlabEnginePut_DA2d)); 6247c6ae99SBarry Smith #endif 6347c6ae99SBarry Smith PetscFunctionReturn(0); 6447c6ae99SBarry Smith } 6547c6ae99SBarry Smith 66939f6067SMatthew G. Knepley /*@ 67939f6067SMatthew G. Knepley DMDAGetNumCells - Get the number of cells in the local piece of the DMDA. This includes ghost cells. 68939f6067SMatthew G. Knepley 69939f6067SMatthew G. Knepley Input Parameter: 70939f6067SMatthew G. Knepley . dm - The DM object 71939f6067SMatthew G. Knepley 72939f6067SMatthew G. Knepley Output Parameters: 73939f6067SMatthew G. Knepley + numCellsX - The number of local cells in the x-direction 74939f6067SMatthew G. Knepley . numCellsY - The number of local cells in the y-direction 75939f6067SMatthew G. Knepley . numCellsZ - The number of local cells in the z-direction 76939f6067SMatthew G. Knepley - numCells - The number of local cells 77939f6067SMatthew G. Knepley 78939f6067SMatthew G. Knepley Level: developer 79939f6067SMatthew G. Knepley 80db781477SPatrick Sanan .seealso: `DMDAGetCellPoint()` 81939f6067SMatthew G. Knepley @*/ 829371c9d4SSatish Balay PetscErrorCode DMDAGetNumCells(DM dm, PetscInt *numCellsX, PetscInt *numCellsY, PetscInt *numCellsZ, PetscInt *numCells) { 8357459e9aSMatthew G Knepley DM_DA *da = (DM_DA *)dm->data; 84c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 8557459e9aSMatthew G Knepley const PetscInt mx = (da->Xe - da->Xs) / da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs; 8657459e9aSMatthew G Knepley const PetscInt nC = (mx) * (dim > 1 ? (my) * (dim > 2 ? (mz) : 1) : 1); 8757459e9aSMatthew G Knepley 8857459e9aSMatthew G Knepley PetscFunctionBegin; 89a9a02de4SBarry Smith PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMDA); 903582350cSMatthew G. Knepley if (numCellsX) { 913582350cSMatthew G. Knepley PetscValidIntPointer(numCellsX, 2); 923582350cSMatthew G. Knepley *numCellsX = mx; 933582350cSMatthew G. Knepley } 943582350cSMatthew G. Knepley if (numCellsY) { 958135c375SStefano Zampini PetscValidIntPointer(numCellsY, 3); 963582350cSMatthew G. Knepley *numCellsY = my; 973582350cSMatthew G. Knepley } 983582350cSMatthew G. Knepley if (numCellsZ) { 998135c375SStefano Zampini PetscValidIntPointer(numCellsZ, 4); 1003582350cSMatthew G. Knepley *numCellsZ = mz; 1013582350cSMatthew G. Knepley } 10257459e9aSMatthew G Knepley if (numCells) { 1033582350cSMatthew G. Knepley PetscValidIntPointer(numCells, 5); 10457459e9aSMatthew G Knepley *numCells = nC; 10557459e9aSMatthew G Knepley } 10657459e9aSMatthew G Knepley PetscFunctionReturn(0); 10757459e9aSMatthew G Knepley } 10857459e9aSMatthew G Knepley 109d6401b93SMatthew G. Knepley /*@ 110d6401b93SMatthew G. Knepley DMDAGetCellPoint - Get the DM point corresponding to the tuple (i, j, k) in the DMDA 111d6401b93SMatthew G. Knepley 112d6401b93SMatthew G. Knepley Input Parameters: 113d6401b93SMatthew G. Knepley + dm - The DM object 114d6401b93SMatthew G. Knepley - i,j,k - The global indices for the cell 115d6401b93SMatthew G. Knepley 116d6401b93SMatthew G. Knepley Output Parameters: 117d6401b93SMatthew G. Knepley . point - The local DM point 118d6401b93SMatthew G. Knepley 119d6401b93SMatthew G. Knepley Level: developer 120d6401b93SMatthew G. Knepley 121db781477SPatrick Sanan .seealso: `DMDAGetNumCells()` 122d6401b93SMatthew G. Knepley @*/ 1239371c9d4SSatish Balay PetscErrorCode DMDAGetCellPoint(DM dm, PetscInt i, PetscInt j, PetscInt k, PetscInt *point) { 124c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 125d6401b93SMatthew G. Knepley DMDALocalInfo info; 126d6401b93SMatthew G. Knepley 127d6401b93SMatthew G. Knepley PetscFunctionBegin; 128a9a02de4SBarry Smith PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMDA); 129d6401b93SMatthew G. Knepley PetscValidIntPointer(point, 5); 1309566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(dm, &info)); 13163a3b9bcSJacob Faibussowitsch if (dim > 0) PetscCheck(!(i < info.gxs) && !(i >= info.gxs + info.gxm), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "X index %" PetscInt_FMT " not in [%" PetscInt_FMT ", %" PetscInt_FMT ")", i, info.gxs, info.gxs + info.gxm); 13263a3b9bcSJacob Faibussowitsch if (dim > 1) PetscCheck(!(j < info.gys) && !(j >= info.gys + info.gym), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Y index %" PetscInt_FMT " not in [%" PetscInt_FMT ", %" PetscInt_FMT ")", j, info.gys, info.gys + info.gym); 13363a3b9bcSJacob Faibussowitsch if (dim > 2) PetscCheck(!(k < info.gzs) && !(k >= info.gzs + info.gzm), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Z index %" PetscInt_FMT PetscInt_FMT " not in [%" PetscInt_FMT ", %" PetscInt_FMT ")", k, info.gzs, info.gzs + info.gzm); 134d6401b93SMatthew G. Knepley *point = i + (dim > 1 ? (j + (dim > 2 ? k * info.gym : 0)) * info.gxm : 0); 135d6401b93SMatthew G. Knepley PetscFunctionReturn(0); 136d6401b93SMatthew G. Knepley } 137d6401b93SMatthew G. Knepley 1389371c9d4SSatish Balay PetscErrorCode DMDAGetNumVertices(DM dm, PetscInt *numVerticesX, PetscInt *numVerticesY, PetscInt *numVerticesZ, PetscInt *numVertices) { 13957459e9aSMatthew G Knepley DM_DA *da = (DM_DA *)dm->data; 140c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 14157459e9aSMatthew G Knepley const PetscInt mx = (da->Xe - da->Xs) / da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs; 14257459e9aSMatthew G Knepley const PetscInt nVx = mx + 1; 14357459e9aSMatthew G Knepley const PetscInt nVy = dim > 1 ? (my + 1) : 1; 14457459e9aSMatthew G Knepley const PetscInt nVz = dim > 2 ? (mz + 1) : 1; 14557459e9aSMatthew G Knepley const PetscInt nV = nVx * nVy * nVz; 14657459e9aSMatthew G Knepley 14757459e9aSMatthew G Knepley PetscFunctionBegin; 14857459e9aSMatthew G Knepley if (numVerticesX) { 14957459e9aSMatthew G Knepley PetscValidIntPointer(numVerticesX, 2); 15057459e9aSMatthew G Knepley *numVerticesX = nVx; 15157459e9aSMatthew G Knepley } 15257459e9aSMatthew G Knepley if (numVerticesY) { 15357459e9aSMatthew G Knepley PetscValidIntPointer(numVerticesY, 3); 15457459e9aSMatthew G Knepley *numVerticesY = nVy; 15557459e9aSMatthew G Knepley } 15657459e9aSMatthew G Knepley if (numVerticesZ) { 15757459e9aSMatthew G Knepley PetscValidIntPointer(numVerticesZ, 4); 15857459e9aSMatthew G Knepley *numVerticesZ = nVz; 15957459e9aSMatthew G Knepley } 16057459e9aSMatthew G Knepley if (numVertices) { 16157459e9aSMatthew G Knepley PetscValidIntPointer(numVertices, 5); 16257459e9aSMatthew G Knepley *numVertices = nV; 16357459e9aSMatthew G Knepley } 16457459e9aSMatthew G Knepley PetscFunctionReturn(0); 16557459e9aSMatthew G Knepley } 16657459e9aSMatthew G Knepley 1679371c9d4SSatish Balay PetscErrorCode DMDAGetNumFaces(DM dm, PetscInt *numXFacesX, PetscInt *numXFaces, PetscInt *numYFacesY, PetscInt *numYFaces, PetscInt *numZFacesZ, PetscInt *numZFaces) { 16857459e9aSMatthew G Knepley DM_DA *da = (DM_DA *)dm->data; 169c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 17057459e9aSMatthew G Knepley const PetscInt mx = (da->Xe - da->Xs) / da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs; 17157459e9aSMatthew G Knepley const PetscInt nxF = (dim > 1 ? (my) * (dim > 2 ? (mz) : 1) : 1); 17257459e9aSMatthew G Knepley const PetscInt nXF = (mx + 1) * nxF; 17357459e9aSMatthew G Knepley const PetscInt nyF = mx * (dim > 2 ? mz : 1); 17457459e9aSMatthew G Knepley const PetscInt nYF = dim > 1 ? (my + 1) * nyF : 0; 17557459e9aSMatthew G Knepley const PetscInt nzF = mx * (dim > 1 ? my : 0); 17657459e9aSMatthew G Knepley const PetscInt nZF = dim > 2 ? (mz + 1) * nzF : 0; 17757459e9aSMatthew G Knepley 17857459e9aSMatthew G Knepley PetscFunctionBegin; 17957459e9aSMatthew G Knepley if (numXFacesX) { 18057459e9aSMatthew G Knepley PetscValidIntPointer(numXFacesX, 2); 18157459e9aSMatthew G Knepley *numXFacesX = nxF; 18257459e9aSMatthew G Knepley } 18357459e9aSMatthew G Knepley if (numXFaces) { 18457459e9aSMatthew G Knepley PetscValidIntPointer(numXFaces, 3); 18557459e9aSMatthew G Knepley *numXFaces = nXF; 18657459e9aSMatthew G Knepley } 18757459e9aSMatthew G Knepley if (numYFacesY) { 18857459e9aSMatthew G Knepley PetscValidIntPointer(numYFacesY, 4); 18957459e9aSMatthew G Knepley *numYFacesY = nyF; 19057459e9aSMatthew G Knepley } 19157459e9aSMatthew G Knepley if (numYFaces) { 19257459e9aSMatthew G Knepley PetscValidIntPointer(numYFaces, 5); 19357459e9aSMatthew G Knepley *numYFaces = nYF; 19457459e9aSMatthew G Knepley } 19557459e9aSMatthew G Knepley if (numZFacesZ) { 19657459e9aSMatthew G Knepley PetscValidIntPointer(numZFacesZ, 6); 19757459e9aSMatthew G Knepley *numZFacesZ = nzF; 19857459e9aSMatthew G Knepley } 19957459e9aSMatthew G Knepley if (numZFaces) { 20057459e9aSMatthew G Knepley PetscValidIntPointer(numZFaces, 7); 20157459e9aSMatthew G Knepley *numZFaces = nZF; 20257459e9aSMatthew G Knepley } 20357459e9aSMatthew G Knepley PetscFunctionReturn(0); 20457459e9aSMatthew G Knepley } 20557459e9aSMatthew G Knepley 2069371c9d4SSatish Balay PetscErrorCode DMDAGetHeightStratum(DM dm, PetscInt height, PetscInt *pStart, PetscInt *pEnd) { 207c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 20857459e9aSMatthew G Knepley PetscInt nC, nV, nXF, nYF, nZF; 20957459e9aSMatthew G Knepley 21057459e9aSMatthew G Knepley PetscFunctionBegin; 2118865f1eaSKarl Rupp if (pStart) PetscValidIntPointer(pStart, 3); 2128865f1eaSKarl Rupp if (pEnd) PetscValidIntPointer(pEnd, 4); 2139566063dSJacob Faibussowitsch PetscCall(DMDAGetNumCells(dm, NULL, NULL, NULL, &nC)); 2149566063dSJacob Faibussowitsch PetscCall(DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV)); 2159566063dSJacob Faibussowitsch PetscCall(DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF)); 21657459e9aSMatthew G Knepley if (height == 0) { 21757459e9aSMatthew G Knepley /* Cells */ 2188865f1eaSKarl Rupp if (pStart) *pStart = 0; 2198865f1eaSKarl Rupp if (pEnd) *pEnd = nC; 22057459e9aSMatthew G Knepley } else if (height == 1) { 22157459e9aSMatthew G Knepley /* Faces */ 2228865f1eaSKarl Rupp if (pStart) *pStart = nC + nV; 2238865f1eaSKarl Rupp if (pEnd) *pEnd = nC + nV + nXF + nYF + nZF; 22457459e9aSMatthew G Knepley } else if (height == dim) { 22557459e9aSMatthew G Knepley /* Vertices */ 2268865f1eaSKarl Rupp if (pStart) *pStart = nC; 2278865f1eaSKarl Rupp if (pEnd) *pEnd = nC + nV; 22857459e9aSMatthew G Knepley } else if (height < 0) { 22957459e9aSMatthew G Knepley /* All points */ 2308865f1eaSKarl Rupp if (pStart) *pStart = 0; 2318865f1eaSKarl Rupp if (pEnd) *pEnd = nC + nV + nXF + nYF + nZF; 23263a3b9bcSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No points of height %" PetscInt_FMT " in the DA", height); 23357459e9aSMatthew G Knepley PetscFunctionReturn(0); 23457459e9aSMatthew G Knepley } 23557459e9aSMatthew G Knepley 2369371c9d4SSatish Balay PetscErrorCode DMDAGetDepthStratum(DM dm, PetscInt depth, PetscInt *pStart, PetscInt *pEnd) { 237c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 2383385ff7eSMatthew G. Knepley PetscInt nC, nV, nXF, nYF, nZF; 2393385ff7eSMatthew G. Knepley 2403385ff7eSMatthew G. Knepley PetscFunctionBegin; 2413385ff7eSMatthew G. Knepley if (pStart) PetscValidIntPointer(pStart, 3); 2423385ff7eSMatthew G. Knepley if (pEnd) PetscValidIntPointer(pEnd, 4); 2439566063dSJacob Faibussowitsch PetscCall(DMDAGetNumCells(dm, NULL, NULL, NULL, &nC)); 2449566063dSJacob Faibussowitsch PetscCall(DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV)); 2459566063dSJacob Faibussowitsch PetscCall(DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF)); 2463385ff7eSMatthew G. Knepley if (depth == dim) { 2473385ff7eSMatthew G. Knepley /* Cells */ 2483385ff7eSMatthew G. Knepley if (pStart) *pStart = 0; 2493385ff7eSMatthew G. Knepley if (pEnd) *pEnd = nC; 2503385ff7eSMatthew G. Knepley } else if (depth == dim - 1) { 2513385ff7eSMatthew G. Knepley /* Faces */ 2523385ff7eSMatthew G. Knepley if (pStart) *pStart = nC + nV; 2533385ff7eSMatthew G. Knepley if (pEnd) *pEnd = nC + nV + nXF + nYF + nZF; 2543385ff7eSMatthew G. Knepley } else if (depth == 0) { 2553385ff7eSMatthew G. Knepley /* Vertices */ 2563385ff7eSMatthew G. Knepley if (pStart) *pStart = nC; 2573385ff7eSMatthew G. Knepley if (pEnd) *pEnd = nC + nV; 2583385ff7eSMatthew G. Knepley } else if (depth < 0) { 2593385ff7eSMatthew G. Knepley /* All points */ 2603385ff7eSMatthew G. Knepley if (pStart) *pStart = 0; 2613385ff7eSMatthew G. Knepley if (pEnd) *pEnd = nC + nV + nXF + nYF + nZF; 26263a3b9bcSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No points of depth %" PetscInt_FMT " in the DA", depth); 2633385ff7eSMatthew G. Knepley PetscFunctionReturn(0); 2643385ff7eSMatthew G. Knepley } 2653385ff7eSMatthew G. Knepley 2669371c9d4SSatish Balay PetscErrorCode DMDAGetConeSize(DM dm, PetscInt p, PetscInt *coneSize) { 267c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 2683385ff7eSMatthew G. Knepley PetscInt nC, nV, nXF, nYF, nZF; 2693385ff7eSMatthew G. Knepley 2703385ff7eSMatthew G. Knepley PetscFunctionBegin; 2713385ff7eSMatthew G. Knepley *coneSize = 0; 2729566063dSJacob Faibussowitsch PetscCall(DMDAGetNumCells(dm, NULL, NULL, NULL, &nC)); 2739566063dSJacob Faibussowitsch PetscCall(DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV)); 2749566063dSJacob Faibussowitsch PetscCall(DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF)); 2753385ff7eSMatthew G. Knepley switch (dim) { 2763385ff7eSMatthew G. Knepley case 2: 2773385ff7eSMatthew G. Knepley if (p >= 0) { 2783385ff7eSMatthew G. Knepley if (p < nC) { 2793385ff7eSMatthew G. Knepley *coneSize = 4; 2803385ff7eSMatthew G. Knepley } else if (p < nC + nV) { 2813385ff7eSMatthew G. Knepley *coneSize = 0; 2823385ff7eSMatthew G. Knepley } else if (p < nC + nV + nXF + nYF + nZF) { 2833385ff7eSMatthew G. Knepley *coneSize = 2; 28463a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %" PetscInt_FMT " should be in [0, %" PetscInt_FMT ")", p, nC + nV + nXF + nYF + nZF); 28563a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative point %" PetscInt_FMT " is invalid", p); 2863385ff7eSMatthew G. Knepley break; 2879371c9d4SSatish Balay case 3: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do 3D"); 2883385ff7eSMatthew G. Knepley } 2893385ff7eSMatthew G. Knepley PetscFunctionReturn(0); 2903385ff7eSMatthew G. Knepley } 2913385ff7eSMatthew G. Knepley 2929371c9d4SSatish Balay PetscErrorCode DMDAGetCone(DM dm, PetscInt p, PetscInt *cone[]) { 293c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim; 2943385ff7eSMatthew G. Knepley PetscInt nCx, nCy, nCz, nC, nVx, nVy, nVz, nV, nxF, nyF, nzF, nXF, nYF, nZF; 2953385ff7eSMatthew G. Knepley 2963385ff7eSMatthew G. Knepley PetscFunctionBegin; 2979566063dSJacob Faibussowitsch if (!cone) PetscCall(DMGetWorkArray(dm, 6, MPIU_INT, cone)); 2989566063dSJacob Faibussowitsch PetscCall(DMDAGetNumCells(dm, &nCx, &nCy, &nCz, &nC)); 2999566063dSJacob Faibussowitsch PetscCall(DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV)); 3009566063dSJacob Faibussowitsch PetscCall(DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF)); 3013385ff7eSMatthew G. Knepley switch (dim) { 3023385ff7eSMatthew G. Knepley case 2: 3033385ff7eSMatthew G. Knepley if (p >= 0) { 3043385ff7eSMatthew G. Knepley if (p < nC) { 3053385ff7eSMatthew G. Knepley const PetscInt cy = p / nCx; 3063385ff7eSMatthew G. Knepley const PetscInt cx = p % nCx; 3073385ff7eSMatthew G. Knepley 3083385ff7eSMatthew G. Knepley (*cone)[0] = cy * nxF + cx + nC + nV; 3093385ff7eSMatthew G. Knepley (*cone)[1] = cx * nyF + cy + nyF + nC + nV + nXF; 3103385ff7eSMatthew G. Knepley (*cone)[2] = cy * nxF + cx + nxF + nC + nV; 3113385ff7eSMatthew G. Knepley (*cone)[3] = cx * nyF + cy + nC + nV + nXF; 3123385ff7eSMatthew G. Knepley SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do cell cones"); 3133385ff7eSMatthew G. Knepley } else if (p < nC + nV) { 3143385ff7eSMatthew G. Knepley } else if (p < nC + nV + nXF) { 3153385ff7eSMatthew G. Knepley const PetscInt fy = (p - nC + nV) / nxF; 3163385ff7eSMatthew G. Knepley const PetscInt fx = (p - nC + nV) % nxF; 3173385ff7eSMatthew G. Knepley 3183385ff7eSMatthew G. Knepley (*cone)[0] = fy * nVx + fx + nC; 3193385ff7eSMatthew G. Knepley (*cone)[1] = fy * nVx + fx + 1 + nC; 3203385ff7eSMatthew G. Knepley } else if (p < nC + nV + nXF + nYF) { 3213385ff7eSMatthew G. Knepley const PetscInt fx = (p - nC + nV + nXF) / nyF; 3223385ff7eSMatthew G. Knepley const PetscInt fy = (p - nC + nV + nXF) % nyF; 3233385ff7eSMatthew G. Knepley 3243385ff7eSMatthew G. Knepley (*cone)[0] = fy * nVx + fx + nC; 3253385ff7eSMatthew G. Knepley (*cone)[1] = fy * nVx + fx + nVx + nC; 32663a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %" PetscInt_FMT " should be in [0, %" PetscInt_FMT ")", p, nC + nV + nXF + nYF + nZF); 32763a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative point %" PetscInt_FMT " is invalid", p); 3283385ff7eSMatthew G. Knepley break; 3299371c9d4SSatish Balay case 3: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do 3D"); 3303385ff7eSMatthew G. Knepley } 3313385ff7eSMatthew G. Knepley PetscFunctionReturn(0); 3323385ff7eSMatthew G. Knepley } 3333385ff7eSMatthew G. Knepley 3349371c9d4SSatish Balay PetscErrorCode DMDARestoreCone(DM dm, PetscInt p, PetscInt *cone[]) { 3353385ff7eSMatthew G. Knepley PetscFunctionBegin; 3369566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(dm, 6, MPIU_INT, cone)); 3373385ff7eSMatthew G. Knepley PetscFunctionReturn(0); 3383385ff7eSMatthew G. Knepley } 3393385ff7eSMatthew G. Knepley 3409371c9d4SSatish Balay PetscErrorCode DMDASetVertexCoordinates(DM dm, PetscReal xl, PetscReal xu, PetscReal yl, PetscReal yu, PetscReal zl, PetscReal zu) { 3413385ff7eSMatthew G. Knepley DM_DA *da = (DM_DA *)dm->data; 3423385ff7eSMatthew G. Knepley Vec coordinates; 3433385ff7eSMatthew G. Knepley PetscSection section; 3443385ff7eSMatthew G. Knepley PetscScalar *coords; 3453385ff7eSMatthew G. Knepley PetscReal h[3]; 3463385ff7eSMatthew G. Knepley PetscInt dim, size, M, N, P, nVx, nVy, nVz, nV, vStart, vEnd, v, i, j, k; 3473385ff7eSMatthew G. Knepley 3483385ff7eSMatthew G. Knepley PetscFunctionBegin; 349a9a02de4SBarry Smith PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMDA); 3509566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(dm, &dim, &M, &N, &P, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 35108401ef6SPierre Jolivet PetscCheck(dim <= 3, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "The following code only works for dim <= 3"); 3523385ff7eSMatthew G. Knepley h[0] = (xu - xl) / M; 3533385ff7eSMatthew G. Knepley h[1] = (yu - yl) / N; 3543385ff7eSMatthew G. Knepley h[2] = (zu - zl) / P; 3559566063dSJacob Faibussowitsch PetscCall(DMDAGetDepthStratum(dm, 0, &vStart, &vEnd)); 3569566063dSJacob Faibussowitsch PetscCall(DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV)); 3579566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), §ion)); 3589566063dSJacob Faibussowitsch PetscCall(PetscSectionSetNumFields(section, 1)); 3599566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldComponents(section, 0, dim)); 3609566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(section, vStart, vEnd)); 361*48a46eb9SPierre Jolivet for (v = vStart; v < vEnd; ++v) PetscCall(PetscSectionSetDof(section, v, dim)); 3629566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(section)); 3639566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(section, &size)); 3649566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, size, &coordinates)); 3659566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates")); 3669566063dSJacob Faibussowitsch PetscCall(VecGetArray(coordinates, &coords)); 3673385ff7eSMatthew G. Knepley for (k = 0; k < nVz; ++k) { 368e4d69e08SBarry Smith PetscInt ind[3], d, off; 3693385ff7eSMatthew G. Knepley 370e4d69e08SBarry Smith ind[0] = 0; 371e4d69e08SBarry Smith ind[1] = 0; 372e4d69e08SBarry Smith ind[2] = k + da->zs; 3733385ff7eSMatthew G. Knepley for (j = 0; j < nVy; ++j) { 3743385ff7eSMatthew G. Knepley ind[1] = j + da->ys; 3753385ff7eSMatthew G. Knepley for (i = 0; i < nVx; ++i) { 3763385ff7eSMatthew G. Knepley const PetscInt vertex = (k * nVy + j) * nVx + i + vStart; 3773385ff7eSMatthew G. Knepley 3789566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, vertex, &off)); 3793385ff7eSMatthew G. Knepley ind[0] = i + da->xs; 3809371c9d4SSatish Balay for (d = 0; d < dim; ++d) { coords[off + d] = h[d] * ind[d]; } 3813385ff7eSMatthew G. Knepley } 3823385ff7eSMatthew G. Knepley } 3833385ff7eSMatthew G. Knepley } 3849566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coordinates, &coords)); 3859566063dSJacob Faibussowitsch PetscCall(DMSetCoordinateSection(dm, PETSC_DETERMINE, section)); 3869566063dSJacob Faibussowitsch PetscCall(DMSetCoordinatesLocal(dm, coordinates)); 3879566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(§ion)); 3889566063dSJacob Faibussowitsch PetscCall(VecDestroy(&coordinates)); 3893385ff7eSMatthew G. Knepley PetscFunctionReturn(0); 3903385ff7eSMatthew G. Knepley } 3919a800dd8SMatthew G. Knepley 39247c6ae99SBarry Smith /* ------------------------------------------------------------------- */ 39347c6ae99SBarry Smith 39447c6ae99SBarry Smith /*@C 395aa219208SBarry Smith DMDAGetArray - Gets a work array for a DMDA 39647c6ae99SBarry Smith 397d8d19677SJose E. Roman Input Parameters: 39847c6ae99SBarry Smith + da - information about my local patch 39947c6ae99SBarry Smith - ghosted - do you want arrays for the ghosted or nonghosted patch 40047c6ae99SBarry Smith 40147c6ae99SBarry Smith Output Parameters: 40247c6ae99SBarry Smith . vptr - array data structured 40347c6ae99SBarry Smith 40447c6ae99SBarry Smith Note: The vector values are NOT initialized and may have garbage in them, so you may need 40547c6ae99SBarry Smith to zero them. 40647c6ae99SBarry Smith 40747c6ae99SBarry Smith Level: advanced 40847c6ae99SBarry Smith 409db781477SPatrick Sanan .seealso: `DMDARestoreArray()` 41047c6ae99SBarry Smith 41147c6ae99SBarry Smith @*/ 4129371c9d4SSatish Balay PetscErrorCode DMDAGetArray(DM da, PetscBool ghosted, void *vptr) { 41347c6ae99SBarry Smith PetscInt j, i, xs, ys, xm, ym, zs, zm; 41447c6ae99SBarry Smith char *iarray_start; 41547c6ae99SBarry Smith void **iptr = (void **)vptr; 41647c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 41747c6ae99SBarry Smith 41847c6ae99SBarry Smith PetscFunctionBegin; 419a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 42047c6ae99SBarry Smith if (ghosted) { 421aa219208SBarry Smith for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) { 42247c6ae99SBarry Smith if (dd->arrayghostedin[i]) { 42347c6ae99SBarry Smith *iptr = dd->arrayghostedin[i]; 42447c6ae99SBarry Smith iarray_start = (char *)dd->startghostedin[i]; 4250298fd71SBarry Smith dd->arrayghostedin[i] = NULL; 4260298fd71SBarry Smith dd->startghostedin[i] = NULL; 42747c6ae99SBarry Smith 42847c6ae99SBarry Smith goto done; 42947c6ae99SBarry Smith } 43047c6ae99SBarry Smith } 43147c6ae99SBarry Smith xs = dd->Xs; 43247c6ae99SBarry Smith ys = dd->Ys; 43347c6ae99SBarry Smith zs = dd->Zs; 43447c6ae99SBarry Smith xm = dd->Xe - dd->Xs; 43547c6ae99SBarry Smith ym = dd->Ye - dd->Ys; 43647c6ae99SBarry Smith zm = dd->Ze - dd->Zs; 43747c6ae99SBarry Smith } else { 438aa219208SBarry Smith for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) { 43947c6ae99SBarry Smith if (dd->arrayin[i]) { 44047c6ae99SBarry Smith *iptr = dd->arrayin[i]; 44147c6ae99SBarry Smith iarray_start = (char *)dd->startin[i]; 4420298fd71SBarry Smith dd->arrayin[i] = NULL; 4430298fd71SBarry Smith dd->startin[i] = NULL; 44447c6ae99SBarry Smith 44547c6ae99SBarry Smith goto done; 44647c6ae99SBarry Smith } 44747c6ae99SBarry Smith } 44847c6ae99SBarry Smith xs = dd->xs; 44947c6ae99SBarry Smith ys = dd->ys; 45047c6ae99SBarry Smith zs = dd->zs; 45147c6ae99SBarry Smith xm = dd->xe - dd->xs; 45247c6ae99SBarry Smith ym = dd->ye - dd->ys; 45347c6ae99SBarry Smith zm = dd->ze - dd->zs; 45447c6ae99SBarry Smith } 45547c6ae99SBarry Smith 456c73cfb54SMatthew G. Knepley switch (da->dim) { 45747c6ae99SBarry Smith case 1: { 45847c6ae99SBarry Smith void *ptr; 45947c6ae99SBarry Smith 4609566063dSJacob Faibussowitsch PetscCall(PetscMalloc(xm * sizeof(PetscScalar), &iarray_start)); 46147c6ae99SBarry Smith 46247c6ae99SBarry Smith ptr = (void *)(iarray_start - xs * sizeof(PetscScalar)); 46347c6ae99SBarry Smith *iptr = (void *)ptr; 4648865f1eaSKarl Rupp break; 4658865f1eaSKarl Rupp } 46647c6ae99SBarry Smith case 2: { 46747c6ae99SBarry Smith void **ptr; 46847c6ae99SBarry Smith 4699566063dSJacob Faibussowitsch PetscCall(PetscMalloc((ym + 1) * sizeof(void *) + xm * ym * sizeof(PetscScalar), &iarray_start)); 47047c6ae99SBarry Smith 47147c6ae99SBarry Smith ptr = (void **)(iarray_start + xm * ym * sizeof(PetscScalar) - ys * sizeof(void *)); 4728865f1eaSKarl Rupp for (j = ys; j < ys + ym; j++) ptr[j] = iarray_start + sizeof(PetscScalar) * (xm * (j - ys) - xs); 47347c6ae99SBarry Smith *iptr = (void *)ptr; 4748865f1eaSKarl Rupp break; 4758865f1eaSKarl Rupp } 47647c6ae99SBarry Smith case 3: { 47747c6ae99SBarry Smith void ***ptr, **bptr; 47847c6ae99SBarry Smith 4799566063dSJacob Faibussowitsch PetscCall(PetscMalloc((zm + 1) * sizeof(void **) + (ym * zm + 1) * sizeof(void *) + xm * ym * zm * sizeof(PetscScalar), &iarray_start)); 48047c6ae99SBarry Smith 48147c6ae99SBarry Smith ptr = (void ***)(iarray_start + xm * ym * zm * sizeof(PetscScalar) - zs * sizeof(void *)); 48247c6ae99SBarry Smith bptr = (void **)(iarray_start + xm * ym * zm * sizeof(PetscScalar) + zm * sizeof(void **)); 4838865f1eaSKarl Rupp for (i = zs; i < zs + zm; i++) ptr[i] = bptr + ((i - zs) * ym - ys); 48447c6ae99SBarry Smith for (i = zs; i < zs + zm; i++) { 4859371c9d4SSatish Balay for (j = ys; j < ys + ym; j++) { ptr[i][j] = iarray_start + sizeof(PetscScalar) * (xm * ym * (i - zs) + xm * (j - ys) - xs); } 48647c6ae99SBarry Smith } 48747c6ae99SBarry Smith *iptr = (void *)ptr; 4888865f1eaSKarl Rupp break; 4898865f1eaSKarl Rupp } 4909371c9d4SSatish Balay default: SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Dimension %" PetscInt_FMT " not supported", da->dim); 49147c6ae99SBarry Smith } 49247c6ae99SBarry Smith 49347c6ae99SBarry Smith done: 49447c6ae99SBarry Smith /* add arrays to the checked out list */ 49547c6ae99SBarry Smith if (ghosted) { 496aa219208SBarry Smith for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) { 49747c6ae99SBarry Smith if (!dd->arrayghostedout[i]) { 49847c6ae99SBarry Smith dd->arrayghostedout[i] = *iptr; 49947c6ae99SBarry Smith dd->startghostedout[i] = iarray_start; 50047c6ae99SBarry Smith break; 50147c6ae99SBarry Smith } 50247c6ae99SBarry Smith } 50347c6ae99SBarry Smith } else { 504aa219208SBarry Smith for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) { 50547c6ae99SBarry Smith if (!dd->arrayout[i]) { 50647c6ae99SBarry Smith dd->arrayout[i] = *iptr; 50747c6ae99SBarry Smith dd->startout[i] = iarray_start; 50847c6ae99SBarry Smith break; 50947c6ae99SBarry Smith } 51047c6ae99SBarry Smith } 51147c6ae99SBarry Smith } 51247c6ae99SBarry Smith PetscFunctionReturn(0); 51347c6ae99SBarry Smith } 51447c6ae99SBarry Smith 51547c6ae99SBarry Smith /*@C 516aa219208SBarry Smith DMDARestoreArray - Restores an array of derivative types for a DMDA 51747c6ae99SBarry Smith 518d8d19677SJose E. Roman Input Parameters: 51947c6ae99SBarry Smith + da - information about my local patch 52047c6ae99SBarry Smith . ghosted - do you want arrays for the ghosted or nonghosted patch 521cb004a26SBarry Smith - vptr - array data structured 52247c6ae99SBarry Smith 52347c6ae99SBarry Smith Level: advanced 52447c6ae99SBarry Smith 525db781477SPatrick Sanan .seealso: `DMDAGetArray()` 52647c6ae99SBarry Smith 52747c6ae99SBarry Smith @*/ 5289371c9d4SSatish Balay PetscErrorCode DMDARestoreArray(DM da, PetscBool ghosted, void *vptr) { 52947c6ae99SBarry Smith PetscInt i; 530ea78f98cSLisandro Dalcin void **iptr = (void **)vptr, *iarray_start = NULL; 53147c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data; 53247c6ae99SBarry Smith 53347c6ae99SBarry Smith PetscFunctionBegin; 534a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 53547c6ae99SBarry Smith if (ghosted) { 536aa219208SBarry Smith for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) { 53747c6ae99SBarry Smith if (dd->arrayghostedout[i] == *iptr) { 53847c6ae99SBarry Smith iarray_start = dd->startghostedout[i]; 5390298fd71SBarry Smith dd->arrayghostedout[i] = NULL; 5400298fd71SBarry Smith dd->startghostedout[i] = NULL; 54147c6ae99SBarry Smith break; 54247c6ae99SBarry Smith } 54347c6ae99SBarry Smith } 544aa219208SBarry Smith for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) { 54547c6ae99SBarry Smith if (!dd->arrayghostedin[i]) { 54647c6ae99SBarry Smith dd->arrayghostedin[i] = *iptr; 54747c6ae99SBarry Smith dd->startghostedin[i] = iarray_start; 54847c6ae99SBarry Smith break; 54947c6ae99SBarry Smith } 55047c6ae99SBarry Smith } 55147c6ae99SBarry Smith } else { 552aa219208SBarry Smith for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) { 55347c6ae99SBarry Smith if (dd->arrayout[i] == *iptr) { 55447c6ae99SBarry Smith iarray_start = dd->startout[i]; 5550298fd71SBarry Smith dd->arrayout[i] = NULL; 5560298fd71SBarry Smith dd->startout[i] = NULL; 55747c6ae99SBarry Smith break; 55847c6ae99SBarry Smith } 55947c6ae99SBarry Smith } 560aa219208SBarry Smith for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) { 56147c6ae99SBarry Smith if (!dd->arrayin[i]) { 56247c6ae99SBarry Smith dd->arrayin[i] = *iptr; 56347c6ae99SBarry Smith dd->startin[i] = iarray_start; 56447c6ae99SBarry Smith break; 56547c6ae99SBarry Smith } 56647c6ae99SBarry Smith } 56747c6ae99SBarry Smith } 56847c6ae99SBarry Smith PetscFunctionReturn(0); 56947c6ae99SBarry Smith } 570