147c6ae99SBarry Smith /*
247c6ae99SBarry Smith Code for manipulating distributed regular arrays in parallel.
347c6ae99SBarry Smith */
447c6ae99SBarry Smith
5af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/
6b2fff234SMatthew G. Knepley #include <petscbt.h>
70c312b8eSJed Brown #include <petscsf.h>
82764a2aaSMatthew G. Knepley #include <petscds.h>
99a800dd8SMatthew G. Knepley #include <petscfe.h>
1047c6ae99SBarry Smith
1147c6ae99SBarry Smith /*
12e3c5b3baSBarry Smith This allows the DMDA vectors to properly tell MATLAB their dimensions
1347c6ae99SBarry Smith */
14d1e78c4fSBarry Smith #if defined(PETSC_HAVE_MATLAB)
15c6db04a5SJed Brown #include <engine.h> /* MATLAB include file */
16c6db04a5SJed Brown #include <mex.h> /* MATLAB include file */
VecMatlabEnginePut_DA2d(PetscObject obj,void * mengine)17d71ae5a4SJacob Faibussowitsch static PetscErrorCode VecMatlabEnginePut_DA2d(PetscObject obj, void *mengine)
18d71ae5a4SJacob Faibussowitsch {
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));
413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
4247c6ae99SBarry Smith }
4347c6ae99SBarry Smith #endif
4447c6ae99SBarry Smith
DMCreateLocalVector_DA(DM da,Vec * g)45d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateLocalVector_DA(DM da, Vec *g)
46d71ae5a4SJacob Faibussowitsch {
4747c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data;
4847c6ae99SBarry Smith
4947c6ae99SBarry Smith PetscFunctionBegin;
5047c6ae99SBarry Smith PetscValidHeaderSpecific(da, DM_CLASSID, 1);
514f572ea9SToby Isaac PetscAssertPointer(g, 2);
529566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, g));
539566063dSJacob Faibussowitsch PetscCall(VecSetSizes(*g, dd->nlocal, PETSC_DETERMINE));
549566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(*g, dd->w));
559566063dSJacob Faibussowitsch PetscCall(VecSetType(*g, da->vectype));
5674427ab1SRichard Tran Mills if (dd->nlocal < da->bind_below) {
579566063dSJacob Faibussowitsch PetscCall(VecSetBindingPropagates(*g, PETSC_TRUE));
589566063dSJacob Faibussowitsch PetscCall(VecBindToCPU(*g, PETSC_TRUE));
5974427ab1SRichard Tran Mills }
609566063dSJacob Faibussowitsch PetscCall(VecSetDM(*g, da));
61d1e78c4fSBarry Smith #if defined(PETSC_HAVE_MATLAB)
6248a46eb9SPierre Jolivet if (dd->w == 1 && da->dim == 2) PetscCall(PetscObjectComposeFunction((PetscObject)*g, "PetscMatlabEnginePut_C", VecMatlabEnginePut_DA2d));
6347c6ae99SBarry Smith #endif
643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
6547c6ae99SBarry Smith }
6647c6ae99SBarry Smith
67939f6067SMatthew G. Knepley /*@
6812b4a537SBarry Smith DMDAGetNumCells - Get the number of cells (or vertices) in the local piece of the `DMDA`. This includes ghost cells.
69939f6067SMatthew G. Knepley
70939f6067SMatthew G. Knepley Input Parameter:
71dce8aebaSBarry Smith . dm - The `DMDA` object
72939f6067SMatthew G. Knepley
73939f6067SMatthew G. Knepley Output Parameters:
74939f6067SMatthew G. Knepley + numCellsX - The number of local cells in the x-direction
75939f6067SMatthew G. Knepley . numCellsY - The number of local cells in the y-direction
76939f6067SMatthew G. Knepley . numCellsZ - The number of local cells in the z-direction
77939f6067SMatthew G. Knepley - numCells - The number of local cells
78939f6067SMatthew G. Knepley
79939f6067SMatthew G. Knepley Level: developer
80939f6067SMatthew G. Knepley
8112b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDAGetCellPoint()`
82939f6067SMatthew G. Knepley @*/
DMDAGetNumCells(DM dm,PeOp PetscInt * numCellsX,PeOp PetscInt * numCellsY,PeOp PetscInt * numCellsZ,PeOp PetscInt * numCells)83*ce78bad3SBarry Smith PetscErrorCode DMDAGetNumCells(DM dm, PeOp PetscInt *numCellsX, PeOp PetscInt *numCellsY, PeOp PetscInt *numCellsZ, PeOp PetscInt *numCells)
84d71ae5a4SJacob Faibussowitsch {
8557459e9aSMatthew G Knepley DM_DA *da = (DM_DA *)dm->data;
86c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim;
8757459e9aSMatthew G Knepley const PetscInt mx = (da->Xe - da->Xs) / da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs;
8857459e9aSMatthew G Knepley const PetscInt nC = (mx) * (dim > 1 ? (my) * (dim > 2 ? (mz) : 1) : 1);
8957459e9aSMatthew G Knepley
9057459e9aSMatthew G Knepley PetscFunctionBegin;
91a9a02de4SBarry Smith PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMDA);
923582350cSMatthew G. Knepley if (numCellsX) {
934f572ea9SToby Isaac PetscAssertPointer(numCellsX, 2);
943582350cSMatthew G. Knepley *numCellsX = mx;
953582350cSMatthew G. Knepley }
963582350cSMatthew G. Knepley if (numCellsY) {
974f572ea9SToby Isaac PetscAssertPointer(numCellsY, 3);
983582350cSMatthew G. Knepley *numCellsY = my;
993582350cSMatthew G. Knepley }
1003582350cSMatthew G. Knepley if (numCellsZ) {
1014f572ea9SToby Isaac PetscAssertPointer(numCellsZ, 4);
1023582350cSMatthew G. Knepley *numCellsZ = mz;
1033582350cSMatthew G. Knepley }
10457459e9aSMatthew G Knepley if (numCells) {
1054f572ea9SToby Isaac PetscAssertPointer(numCells, 5);
10657459e9aSMatthew G Knepley *numCells = nC;
10757459e9aSMatthew G Knepley }
1083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
10957459e9aSMatthew G Knepley }
11057459e9aSMatthew G Knepley
111d6401b93SMatthew G. Knepley /*@
11212b4a537SBarry Smith DMDAGetCellPoint - Get the `DM` point corresponding to the tuple (i, j, k) in the `DMDA`
113d6401b93SMatthew G. Knepley
114d6401b93SMatthew G. Knepley Input Parameters:
115dce8aebaSBarry Smith + dm - The `DMDA` object
116a4e35b19SJacob Faibussowitsch . i - The global x index for the cell
117a4e35b19SJacob Faibussowitsch . j - The global y index for the cell
118a4e35b19SJacob Faibussowitsch - k - The global z index for the cell
119d6401b93SMatthew G. Knepley
1202fe279fdSBarry Smith Output Parameter:
121dce8aebaSBarry Smith . point - The local `DM` point
122d6401b93SMatthew G. Knepley
123d6401b93SMatthew G. Knepley Level: developer
124d6401b93SMatthew G. Knepley
12512b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDAGetNumCells()`
126d6401b93SMatthew G. Knepley @*/
DMDAGetCellPoint(DM dm,PetscInt i,PetscInt j,PetscInt k,PetscInt * point)127d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetCellPoint(DM dm, PetscInt i, PetscInt j, PetscInt k, PetscInt *point)
128d71ae5a4SJacob Faibussowitsch {
129c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim;
130d6401b93SMatthew G. Knepley DMDALocalInfo info;
131d6401b93SMatthew G. Knepley
132d6401b93SMatthew G. Knepley PetscFunctionBegin;
133a9a02de4SBarry Smith PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMDA);
1344f572ea9SToby Isaac PetscAssertPointer(point, 5);
1359566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(dm, &info));
13663a3b9bcSJacob 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);
13763a3b9bcSJacob 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);
13815229ffcSPierre Jolivet if (dim > 2) PetscCheck(!(k < info.gzs) && !(k >= info.gzs + info.gzm), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Z index %" PetscInt_FMT " not in [%" PetscInt_FMT ", %" PetscInt_FMT ")", k, info.gzs, info.gzs + info.gzm);
139d6401b93SMatthew G. Knepley *point = i + (dim > 1 ? (j + (dim > 2 ? k * info.gym : 0)) * info.gxm : 0);
1403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
141d6401b93SMatthew G. Knepley }
142d6401b93SMatthew G. Knepley
DMDAGetNumVertices(DM dm,PetscInt * numVerticesX,PetscInt * numVerticesY,PetscInt * numVerticesZ,PetscInt * numVertices)143d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetNumVertices(DM dm, PetscInt *numVerticesX, PetscInt *numVerticesY, PetscInt *numVerticesZ, PetscInt *numVertices)
144d71ae5a4SJacob Faibussowitsch {
14557459e9aSMatthew G Knepley DM_DA *da = (DM_DA *)dm->data;
146c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim;
14757459e9aSMatthew G Knepley const PetscInt mx = (da->Xe - da->Xs) / da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs;
14857459e9aSMatthew G Knepley const PetscInt nVx = mx + 1;
14957459e9aSMatthew G Knepley const PetscInt nVy = dim > 1 ? (my + 1) : 1;
15057459e9aSMatthew G Knepley const PetscInt nVz = dim > 2 ? (mz + 1) : 1;
15157459e9aSMatthew G Knepley const PetscInt nV = nVx * nVy * nVz;
15257459e9aSMatthew G Knepley
15357459e9aSMatthew G Knepley PetscFunctionBegin;
15457459e9aSMatthew G Knepley if (numVerticesX) {
1554f572ea9SToby Isaac PetscAssertPointer(numVerticesX, 2);
15657459e9aSMatthew G Knepley *numVerticesX = nVx;
15757459e9aSMatthew G Knepley }
15857459e9aSMatthew G Knepley if (numVerticesY) {
1594f572ea9SToby Isaac PetscAssertPointer(numVerticesY, 3);
16057459e9aSMatthew G Knepley *numVerticesY = nVy;
16157459e9aSMatthew G Knepley }
16257459e9aSMatthew G Knepley if (numVerticesZ) {
1634f572ea9SToby Isaac PetscAssertPointer(numVerticesZ, 4);
16457459e9aSMatthew G Knepley *numVerticesZ = nVz;
16557459e9aSMatthew G Knepley }
16657459e9aSMatthew G Knepley if (numVertices) {
1674f572ea9SToby Isaac PetscAssertPointer(numVertices, 5);
16857459e9aSMatthew G Knepley *numVertices = nV;
16957459e9aSMatthew G Knepley }
1703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
17157459e9aSMatthew G Knepley }
17257459e9aSMatthew G Knepley
DMDAGetNumFaces(DM dm,PetscInt * numXFacesX,PetscInt * numXFaces,PetscInt * numYFacesY,PetscInt * numYFaces,PetscInt * numZFacesZ,PetscInt * numZFaces)173d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetNumFaces(DM dm, PetscInt *numXFacesX, PetscInt *numXFaces, PetscInt *numYFacesY, PetscInt *numYFaces, PetscInt *numZFacesZ, PetscInt *numZFaces)
174d71ae5a4SJacob Faibussowitsch {
17557459e9aSMatthew G Knepley DM_DA *da = (DM_DA *)dm->data;
176c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim;
17757459e9aSMatthew G Knepley const PetscInt mx = (da->Xe - da->Xs) / da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs;
17857459e9aSMatthew G Knepley const PetscInt nxF = (dim > 1 ? (my) * (dim > 2 ? (mz) : 1) : 1);
17957459e9aSMatthew G Knepley const PetscInt nXF = (mx + 1) * nxF;
18057459e9aSMatthew G Knepley const PetscInt nyF = mx * (dim > 2 ? mz : 1);
18157459e9aSMatthew G Knepley const PetscInt nYF = dim > 1 ? (my + 1) * nyF : 0;
18257459e9aSMatthew G Knepley const PetscInt nzF = mx * (dim > 1 ? my : 0);
18357459e9aSMatthew G Knepley const PetscInt nZF = dim > 2 ? (mz + 1) * nzF : 0;
18457459e9aSMatthew G Knepley
18557459e9aSMatthew G Knepley PetscFunctionBegin;
18657459e9aSMatthew G Knepley if (numXFacesX) {
1874f572ea9SToby Isaac PetscAssertPointer(numXFacesX, 2);
18857459e9aSMatthew G Knepley *numXFacesX = nxF;
18957459e9aSMatthew G Knepley }
19057459e9aSMatthew G Knepley if (numXFaces) {
1914f572ea9SToby Isaac PetscAssertPointer(numXFaces, 3);
19257459e9aSMatthew G Knepley *numXFaces = nXF;
19357459e9aSMatthew G Knepley }
19457459e9aSMatthew G Knepley if (numYFacesY) {
1954f572ea9SToby Isaac PetscAssertPointer(numYFacesY, 4);
19657459e9aSMatthew G Knepley *numYFacesY = nyF;
19757459e9aSMatthew G Knepley }
19857459e9aSMatthew G Knepley if (numYFaces) {
1994f572ea9SToby Isaac PetscAssertPointer(numYFaces, 5);
20057459e9aSMatthew G Knepley *numYFaces = nYF;
20157459e9aSMatthew G Knepley }
20257459e9aSMatthew G Knepley if (numZFacesZ) {
2034f572ea9SToby Isaac PetscAssertPointer(numZFacesZ, 6);
20457459e9aSMatthew G Knepley *numZFacesZ = nzF;
20557459e9aSMatthew G Knepley }
20657459e9aSMatthew G Knepley if (numZFaces) {
2074f572ea9SToby Isaac PetscAssertPointer(numZFaces, 7);
20857459e9aSMatthew G Knepley *numZFaces = nZF;
20957459e9aSMatthew G Knepley }
2103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
21157459e9aSMatthew G Knepley }
21257459e9aSMatthew G Knepley
213*ce78bad3SBarry Smith /*@
214*ce78bad3SBarry Smith DMDAGetHeightStratum - Get the bounds [`start`, `end`) for all points at a certain height.
215*ce78bad3SBarry Smith
216*ce78bad3SBarry Smith Not Collective
217*ce78bad3SBarry Smith
218*ce78bad3SBarry Smith Input Parameters:
219*ce78bad3SBarry Smith + dm - The `DMDA` object
220*ce78bad3SBarry Smith - height - The requested height
221*ce78bad3SBarry Smith
222*ce78bad3SBarry Smith Output Parameters:
223*ce78bad3SBarry Smith + pStart - The first point at this `height`
224*ce78bad3SBarry Smith - pEnd - One beyond the last point at this `height`
225*ce78bad3SBarry Smith
226*ce78bad3SBarry Smith Level: developer
227*ce78bad3SBarry Smith
228*ce78bad3SBarry Smith Note:
229*ce78bad3SBarry Smith See `DMPlexGetHeightStratum()` for the meaning of these values
230*ce78bad3SBarry Smith
231*ce78bad3SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMDA`, `DMPlexGetDepthStratum()`, `DMPlexGetHeightStratum()`, `DMPlexGetCellTypeStratum()`, `DMPlexGetDepth()`,
232*ce78bad3SBarry Smith `DMPlexGetDepthLabel()`, `DMPlexGetPointDepth()`, `DMPlexSymmetrize()`, `DMPlexInterpolate()`, `DMDAGetDepthStratum()`
233*ce78bad3SBarry Smith @*/
DMDAGetHeightStratum(DM dm,PetscInt height,PeOp PetscInt * pStart,PeOp PetscInt * pEnd)234*ce78bad3SBarry Smith PetscErrorCode DMDAGetHeightStratum(DM dm, PetscInt height, PeOp PetscInt *pStart, PeOp PetscInt *pEnd)
235d71ae5a4SJacob Faibussowitsch {
236c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim;
23757459e9aSMatthew G Knepley PetscInt nC, nV, nXF, nYF, nZF;
23857459e9aSMatthew G Knepley
23957459e9aSMatthew G Knepley PetscFunctionBegin;
2404f572ea9SToby Isaac if (pStart) PetscAssertPointer(pStart, 3);
2414f572ea9SToby Isaac if (pEnd) PetscAssertPointer(pEnd, 4);
2429566063dSJacob Faibussowitsch PetscCall(DMDAGetNumCells(dm, NULL, NULL, NULL, &nC));
2439566063dSJacob Faibussowitsch PetscCall(DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV));
2449566063dSJacob Faibussowitsch PetscCall(DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF));
24557459e9aSMatthew G Knepley if (height == 0) {
24657459e9aSMatthew G Knepley /* Cells */
2478865f1eaSKarl Rupp if (pStart) *pStart = 0;
2488865f1eaSKarl Rupp if (pEnd) *pEnd = nC;
24957459e9aSMatthew G Knepley } else if (height == 1) {
25057459e9aSMatthew G Knepley /* Faces */
2518865f1eaSKarl Rupp if (pStart) *pStart = nC + nV;
2528865f1eaSKarl Rupp if (pEnd) *pEnd = nC + nV + nXF + nYF + nZF;
25357459e9aSMatthew G Knepley } else if (height == dim) {
25457459e9aSMatthew G Knepley /* Vertices */
2558865f1eaSKarl Rupp if (pStart) *pStart = nC;
2568865f1eaSKarl Rupp if (pEnd) *pEnd = nC + nV;
25757459e9aSMatthew G Knepley } else if (height < 0) {
25857459e9aSMatthew G Knepley /* All points */
2598865f1eaSKarl Rupp if (pStart) *pStart = 0;
2608865f1eaSKarl Rupp if (pEnd) *pEnd = nC + nV + nXF + nYF + nZF;
26163a3b9bcSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No points of height %" PetscInt_FMT " in the DA", height);
2623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
26357459e9aSMatthew G Knepley }
26457459e9aSMatthew G Knepley
265*ce78bad3SBarry Smith /*@
266*ce78bad3SBarry Smith DMDAGetDepthStratum - Get the bounds [`start`, `end`) for all points at a certain depth.
267*ce78bad3SBarry Smith
268*ce78bad3SBarry Smith Not Collective
269*ce78bad3SBarry Smith
270*ce78bad3SBarry Smith Input Parameters:
271*ce78bad3SBarry Smith + dm - The `DMDA` object
272*ce78bad3SBarry Smith - depth - The requested depth
273*ce78bad3SBarry Smith
274*ce78bad3SBarry Smith Output Parameters:
275*ce78bad3SBarry Smith + pStart - The first point at this `depth`
276*ce78bad3SBarry Smith - pEnd - One beyond the last point at this `depth`
277*ce78bad3SBarry Smith
278*ce78bad3SBarry Smith Level: developer
279*ce78bad3SBarry Smith
280*ce78bad3SBarry Smith Note:
281*ce78bad3SBarry Smith See `DMPlexGetDepthStratum()` for the meaning of these values
282*ce78bad3SBarry Smith
283*ce78bad3SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMDA`, `DMPlexGetDepthStratum()`, `DMPlexGetHeightStratum()`, `DMPlexGetCellTypeStratum()`, `DMPlexGetDepth()`,
284*ce78bad3SBarry Smith `DMPlexGetDepthLabel()`, `DMPlexGetPointDepth()`, `DMPlexSymmetrize()`, `DMPlexInterpolate()`, `DMDAGetHeightStratum()`
285*ce78bad3SBarry Smith @*/
DMDAGetDepthStratum(DM dm,PetscInt depth,PeOp PetscInt * pStart,PeOp PetscInt * pEnd)286*ce78bad3SBarry Smith PetscErrorCode DMDAGetDepthStratum(DM dm, PetscInt depth, PeOp PetscInt *pStart, PeOp PetscInt *pEnd)
287d71ae5a4SJacob Faibussowitsch {
288c73cfb54SMatthew G. Knepley const PetscInt dim = dm->dim;
2893385ff7eSMatthew G. Knepley PetscInt nC, nV, nXF, nYF, nZF;
2903385ff7eSMatthew G. Knepley
2913385ff7eSMatthew G. Knepley PetscFunctionBegin;
2924f572ea9SToby Isaac if (pStart) PetscAssertPointer(pStart, 3);
2934f572ea9SToby Isaac if (pEnd) PetscAssertPointer(pEnd, 4);
2949566063dSJacob Faibussowitsch PetscCall(DMDAGetNumCells(dm, NULL, NULL, NULL, &nC));
2959566063dSJacob Faibussowitsch PetscCall(DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV));
2969566063dSJacob Faibussowitsch PetscCall(DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF));
2973385ff7eSMatthew G. Knepley if (depth == dim) {
2983385ff7eSMatthew G. Knepley /* Cells */
2993385ff7eSMatthew G. Knepley if (pStart) *pStart = 0;
3003385ff7eSMatthew G. Knepley if (pEnd) *pEnd = nC;
3013385ff7eSMatthew G. Knepley } else if (depth == dim - 1) {
3023385ff7eSMatthew G. Knepley /* Faces */
3033385ff7eSMatthew G. Knepley if (pStart) *pStart = nC + nV;
3043385ff7eSMatthew G. Knepley if (pEnd) *pEnd = nC + nV + nXF + nYF + nZF;
3053385ff7eSMatthew G. Knepley } else if (depth == 0) {
3063385ff7eSMatthew G. Knepley /* Vertices */
3073385ff7eSMatthew G. Knepley if (pStart) *pStart = nC;
3083385ff7eSMatthew G. Knepley if (pEnd) *pEnd = nC + nV;
3093385ff7eSMatthew G. Knepley } else if (depth < 0) {
3103385ff7eSMatthew G. Knepley /* All points */
3113385ff7eSMatthew G. Knepley if (pStart) *pStart = 0;
3123385ff7eSMatthew G. Knepley if (pEnd) *pEnd = nC + nV + nXF + nYF + nZF;
31363a3b9bcSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No points of depth %" PetscInt_FMT " in the DA", depth);
3143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3153385ff7eSMatthew G. Knepley }
3163385ff7eSMatthew G. Knepley
317*ce78bad3SBarry Smith /*@
318*ce78bad3SBarry Smith DMDASetVertexCoordinates - Sets the lower and upper coordinates for a `DMDA`
319*ce78bad3SBarry Smith
320*ce78bad3SBarry Smith Logically Collective
321*ce78bad3SBarry Smith
322*ce78bad3SBarry Smith Input Parameters:
323*ce78bad3SBarry Smith + dm - The `DMDA` object
324*ce78bad3SBarry Smith . xl - the lower x coordinate
325*ce78bad3SBarry Smith . xu - the upper x coordinate
326*ce78bad3SBarry Smith . yl - the lower y coordinate
327*ce78bad3SBarry Smith . yu - the upper y coordinate
328*ce78bad3SBarry Smith . zl - the lower z coordinate
329*ce78bad3SBarry Smith - zu - the upper z coordinate
330*ce78bad3SBarry Smith
331*ce78bad3SBarry Smith Level: intermediate
332*ce78bad3SBarry Smith
333*ce78bad3SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMDA`
334*ce78bad3SBarry Smith @*/
DMDASetVertexCoordinates(DM dm,PetscReal xl,PetscReal xu,PetscReal yl,PetscReal yu,PetscReal zl,PetscReal zu)335d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetVertexCoordinates(DM dm, PetscReal xl, PetscReal xu, PetscReal yl, PetscReal yu, PetscReal zl, PetscReal zu)
336d71ae5a4SJacob Faibussowitsch {
3373385ff7eSMatthew G. Knepley DM_DA *da = (DM_DA *)dm->data;
3383385ff7eSMatthew G. Knepley Vec coordinates;
3393385ff7eSMatthew G. Knepley PetscSection section;
3403385ff7eSMatthew G. Knepley PetscScalar *coords;
3413385ff7eSMatthew G. Knepley PetscReal h[3];
3423385ff7eSMatthew G. Knepley PetscInt dim, size, M, N, P, nVx, nVy, nVz, nV, vStart, vEnd, v, i, j, k;
3433385ff7eSMatthew G. Knepley
3443385ff7eSMatthew G. Knepley PetscFunctionBegin;
345a9a02de4SBarry Smith PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMDA);
3469566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(dm, &dim, &M, &N, &P, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
34708401ef6SPierre Jolivet PetscCheck(dim <= 3, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "The following code only works for dim <= 3");
3483385ff7eSMatthew G. Knepley h[0] = (xu - xl) / M;
3493385ff7eSMatthew G. Knepley h[1] = (yu - yl) / N;
3503385ff7eSMatthew G. Knepley h[2] = (zu - zl) / P;
3519566063dSJacob Faibussowitsch PetscCall(DMDAGetDepthStratum(dm, 0, &vStart, &vEnd));
3529566063dSJacob Faibussowitsch PetscCall(DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV));
3539566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), §ion));
3549566063dSJacob Faibussowitsch PetscCall(PetscSectionSetNumFields(section, 1));
3559566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldComponents(section, 0, dim));
3569566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(section, vStart, vEnd));
35748a46eb9SPierre Jolivet for (v = vStart; v < vEnd; ++v) PetscCall(PetscSectionSetDof(section, v, dim));
3589566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(section));
3599566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(section, &size));
3609566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, size, &coordinates));
3619566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates"));
3629566063dSJacob Faibussowitsch PetscCall(VecGetArray(coordinates, &coords));
3633385ff7eSMatthew G. Knepley for (k = 0; k < nVz; ++k) {
364e4d69e08SBarry Smith PetscInt ind[3], d, off;
3653385ff7eSMatthew G. Knepley
366e4d69e08SBarry Smith ind[0] = 0;
367e4d69e08SBarry Smith ind[1] = 0;
368e4d69e08SBarry Smith ind[2] = k + da->zs;
3693385ff7eSMatthew G. Knepley for (j = 0; j < nVy; ++j) {
3703385ff7eSMatthew G. Knepley ind[1] = j + da->ys;
3713385ff7eSMatthew G. Knepley for (i = 0; i < nVx; ++i) {
3723385ff7eSMatthew G. Knepley const PetscInt vertex = (k * nVy + j) * nVx + i + vStart;
3733385ff7eSMatthew G. Knepley
3749566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, vertex, &off));
3753385ff7eSMatthew G. Knepley ind[0] = i + da->xs;
376ad540459SPierre Jolivet for (d = 0; d < dim; ++d) coords[off + d] = h[d] * ind[d];
3773385ff7eSMatthew G. Knepley }
3783385ff7eSMatthew G. Knepley }
3793385ff7eSMatthew G. Knepley }
3809566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coordinates, &coords));
3819566063dSJacob Faibussowitsch PetscCall(DMSetCoordinateSection(dm, PETSC_DETERMINE, section));
3829566063dSJacob Faibussowitsch PetscCall(DMSetCoordinatesLocal(dm, coordinates));
3839566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(§ion));
3849566063dSJacob Faibussowitsch PetscCall(VecDestroy(&coordinates));
3853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3863385ff7eSMatthew G. Knepley }
3879a800dd8SMatthew G. Knepley
38847c6ae99SBarry Smith /*@C
389dce8aebaSBarry Smith DMDAGetArray - Gets a work array for a `DMDA`
39047c6ae99SBarry Smith
391d8d19677SJose E. Roman Input Parameters:
39212b4a537SBarry Smith + da - a `DMDA`
39347c6ae99SBarry Smith - ghosted - do you want arrays for the ghosted or nonghosted patch
39447c6ae99SBarry Smith
3952fe279fdSBarry Smith Output Parameter:
39647c6ae99SBarry Smith . vptr - array data structured
39747c6ae99SBarry Smith
39847c6ae99SBarry Smith Level: advanced
39947c6ae99SBarry Smith
40012b4a537SBarry Smith Notes:
401dce8aebaSBarry Smith The vector values are NOT initialized and may have garbage in them, so you may need
402dce8aebaSBarry Smith to zero them.
40347c6ae99SBarry Smith
40412b4a537SBarry Smith Use `DMDARestoreArray()` to return the array
40512b4a537SBarry Smith
40612b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDARestoreArray()`
40747c6ae99SBarry Smith @*/
DMDAGetArray(DM da,PetscBool ghosted,void * vptr)408d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetArray(DM da, PetscBool ghosted, void *vptr)
409d71ae5a4SJacob Faibussowitsch {
41047c6ae99SBarry Smith PetscInt j, i, xs, ys, xm, ym, zs, zm;
41147c6ae99SBarry Smith char *iarray_start;
41247c6ae99SBarry Smith void **iptr = (void **)vptr;
41347c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data;
41447c6ae99SBarry Smith
41547c6ae99SBarry Smith PetscFunctionBegin;
416a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
41747c6ae99SBarry Smith if (ghosted) {
418aa219208SBarry Smith for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) {
41947c6ae99SBarry Smith if (dd->arrayghostedin[i]) {
42047c6ae99SBarry Smith *iptr = dd->arrayghostedin[i];
42147c6ae99SBarry Smith iarray_start = (char *)dd->startghostedin[i];
4220298fd71SBarry Smith dd->arrayghostedin[i] = NULL;
4230298fd71SBarry Smith dd->startghostedin[i] = NULL;
42447c6ae99SBarry Smith
42547c6ae99SBarry Smith goto done;
42647c6ae99SBarry Smith }
42747c6ae99SBarry Smith }
42847c6ae99SBarry Smith xs = dd->Xs;
42947c6ae99SBarry Smith ys = dd->Ys;
43047c6ae99SBarry Smith zs = dd->Zs;
43147c6ae99SBarry Smith xm = dd->Xe - dd->Xs;
43247c6ae99SBarry Smith ym = dd->Ye - dd->Ys;
43347c6ae99SBarry Smith zm = dd->Ze - dd->Zs;
43447c6ae99SBarry Smith } else {
435aa219208SBarry Smith for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) {
43647c6ae99SBarry Smith if (dd->arrayin[i]) {
43747c6ae99SBarry Smith *iptr = dd->arrayin[i];
43847c6ae99SBarry Smith iarray_start = (char *)dd->startin[i];
4390298fd71SBarry Smith dd->arrayin[i] = NULL;
4400298fd71SBarry Smith dd->startin[i] = NULL;
44147c6ae99SBarry Smith
44247c6ae99SBarry Smith goto done;
44347c6ae99SBarry Smith }
44447c6ae99SBarry Smith }
44547c6ae99SBarry Smith xs = dd->xs;
44647c6ae99SBarry Smith ys = dd->ys;
44747c6ae99SBarry Smith zs = dd->zs;
44847c6ae99SBarry Smith xm = dd->xe - dd->xs;
44947c6ae99SBarry Smith ym = dd->ye - dd->ys;
45047c6ae99SBarry Smith zm = dd->ze - dd->zs;
45147c6ae99SBarry Smith }
45247c6ae99SBarry Smith
453c73cfb54SMatthew G. Knepley switch (da->dim) {
45447c6ae99SBarry Smith case 1: {
45547c6ae99SBarry Smith void *ptr;
45647c6ae99SBarry Smith
4579566063dSJacob Faibussowitsch PetscCall(PetscMalloc(xm * sizeof(PetscScalar), &iarray_start));
45847c6ae99SBarry Smith
459a680e639SStefano Zampini ptr = (void *)((PetscScalar *)iarray_start - xs);
460835f2295SStefano Zampini *iptr = ptr;
4618865f1eaSKarl Rupp break;
4628865f1eaSKarl Rupp }
46347c6ae99SBarry Smith case 2: {
46447c6ae99SBarry Smith void **ptr;
46547c6ae99SBarry Smith
4669566063dSJacob Faibussowitsch PetscCall(PetscMalloc((ym + 1) * sizeof(void *) + xm * ym * sizeof(PetscScalar), &iarray_start));
46747c6ae99SBarry Smith
46847c6ae99SBarry Smith ptr = (void **)(iarray_start + xm * ym * sizeof(PetscScalar) - ys * sizeof(void *));
4698865f1eaSKarl Rupp for (j = ys; j < ys + ym; j++) ptr[j] = iarray_start + sizeof(PetscScalar) * (xm * (j - ys) - xs);
47047c6ae99SBarry Smith *iptr = (void *)ptr;
4718865f1eaSKarl Rupp break;
4728865f1eaSKarl Rupp }
47347c6ae99SBarry Smith case 3: {
47447c6ae99SBarry Smith void ***ptr, **bptr;
47547c6ae99SBarry Smith
4769566063dSJacob Faibussowitsch PetscCall(PetscMalloc((zm + 1) * sizeof(void **) + (ym * zm + 1) * sizeof(void *) + xm * ym * zm * sizeof(PetscScalar), &iarray_start));
47747c6ae99SBarry Smith
47847c6ae99SBarry Smith ptr = (void ***)(iarray_start + xm * ym * zm * sizeof(PetscScalar) - zs * sizeof(void *));
47947c6ae99SBarry Smith bptr = (void **)(iarray_start + xm * ym * zm * sizeof(PetscScalar) + zm * sizeof(void **));
4808865f1eaSKarl Rupp for (i = zs; i < zs + zm; i++) ptr[i] = bptr + ((i - zs) * ym - ys);
48147c6ae99SBarry Smith for (i = zs; i < zs + zm; i++) {
482ad540459SPierre Jolivet for (j = ys; j < ys + ym; j++) ptr[i][j] = iarray_start + sizeof(PetscScalar) * (xm * ym * (i - zs) + xm * (j - ys) - xs);
48347c6ae99SBarry Smith }
48447c6ae99SBarry Smith *iptr = (void *)ptr;
4858865f1eaSKarl Rupp break;
4868865f1eaSKarl Rupp }
487d71ae5a4SJacob Faibussowitsch default:
488d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Dimension %" PetscInt_FMT " not supported", da->dim);
48947c6ae99SBarry Smith }
49047c6ae99SBarry Smith
49147c6ae99SBarry Smith done:
49247c6ae99SBarry Smith /* add arrays to the checked out list */
49347c6ae99SBarry Smith if (ghosted) {
494aa219208SBarry Smith for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) {
49547c6ae99SBarry Smith if (!dd->arrayghostedout[i]) {
49647c6ae99SBarry Smith dd->arrayghostedout[i] = *iptr;
49747c6ae99SBarry Smith dd->startghostedout[i] = iarray_start;
49847c6ae99SBarry Smith break;
49947c6ae99SBarry Smith }
50047c6ae99SBarry Smith }
50147c6ae99SBarry Smith } else {
502aa219208SBarry Smith for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) {
50347c6ae99SBarry Smith if (!dd->arrayout[i]) {
50447c6ae99SBarry Smith dd->arrayout[i] = *iptr;
50547c6ae99SBarry Smith dd->startout[i] = iarray_start;
50647c6ae99SBarry Smith break;
50747c6ae99SBarry Smith }
50847c6ae99SBarry Smith }
50947c6ae99SBarry Smith }
5103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
51147c6ae99SBarry Smith }
51247c6ae99SBarry Smith
51347c6ae99SBarry Smith /*@C
51412b4a537SBarry Smith DMDARestoreArray - Restores an array for a `DMDA` obtained with `DMDAGetArray()`
51547c6ae99SBarry Smith
516d8d19677SJose E. Roman Input Parameters:
51747c6ae99SBarry Smith + da - information about my local patch
51847c6ae99SBarry Smith . ghosted - do you want arrays for the ghosted or nonghosted patch
519cb004a26SBarry Smith - vptr - array data structured
52047c6ae99SBarry Smith
52147c6ae99SBarry Smith Level: advanced
52247c6ae99SBarry Smith
52312b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDAGetArray()`
52447c6ae99SBarry Smith @*/
DMDARestoreArray(DM da,PetscBool ghosted,void * vptr)525d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDARestoreArray(DM da, PetscBool ghosted, void *vptr)
526d71ae5a4SJacob Faibussowitsch {
52747c6ae99SBarry Smith PetscInt i;
528ea78f98cSLisandro Dalcin void **iptr = (void **)vptr, *iarray_start = NULL;
52947c6ae99SBarry Smith DM_DA *dd = (DM_DA *)da->data;
53047c6ae99SBarry Smith
53147c6ae99SBarry Smith PetscFunctionBegin;
532a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
53347c6ae99SBarry Smith if (ghosted) {
534aa219208SBarry Smith for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) {
53547c6ae99SBarry Smith if (dd->arrayghostedout[i] == *iptr) {
53647c6ae99SBarry Smith iarray_start = dd->startghostedout[i];
5370298fd71SBarry Smith dd->arrayghostedout[i] = NULL;
5380298fd71SBarry Smith dd->startghostedout[i] = NULL;
53947c6ae99SBarry Smith break;
54047c6ae99SBarry Smith }
54147c6ae99SBarry Smith }
542aa219208SBarry Smith for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) {
54347c6ae99SBarry Smith if (!dd->arrayghostedin[i]) {
54447c6ae99SBarry Smith dd->arrayghostedin[i] = *iptr;
54547c6ae99SBarry Smith dd->startghostedin[i] = iarray_start;
54647c6ae99SBarry Smith break;
54747c6ae99SBarry Smith }
54847c6ae99SBarry Smith }
54947c6ae99SBarry Smith } else {
550aa219208SBarry Smith for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) {
55147c6ae99SBarry Smith if (dd->arrayout[i] == *iptr) {
55247c6ae99SBarry Smith iarray_start = dd->startout[i];
5530298fd71SBarry Smith dd->arrayout[i] = NULL;
5540298fd71SBarry Smith dd->startout[i] = NULL;
55547c6ae99SBarry Smith break;
55647c6ae99SBarry Smith }
55747c6ae99SBarry Smith }
558aa219208SBarry Smith for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) {
55947c6ae99SBarry Smith if (!dd->arrayin[i]) {
56047c6ae99SBarry Smith dd->arrayin[i] = *iptr;
56147c6ae99SBarry Smith dd->startin[i] = iarray_start;
56247c6ae99SBarry Smith break;
56347c6ae99SBarry Smith }
56447c6ae99SBarry Smith }
56547c6ae99SBarry Smith }
5663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
56747c6ae99SBarry Smith }
568