xref: /petsc/src/dm/impls/da/dalocal.c (revision 15229ffc342989b2bf0590a733d92c152a3348fc)
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 */
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 
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 @*/
83d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetNumCells(DM dm, PetscInt *numCellsX, PetscInt *numCellsY, PetscInt *numCellsZ, 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 @*/
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);
138*15229ffcSPierre 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 
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 
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 
213d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetHeightStratum(DM dm, PetscInt height, PetscInt *pStart, PetscInt *pEnd)
214d71ae5a4SJacob Faibussowitsch {
215c73cfb54SMatthew G. Knepley   const PetscInt dim = dm->dim;
21657459e9aSMatthew G Knepley   PetscInt       nC, nV, nXF, nYF, nZF;
21757459e9aSMatthew G Knepley 
21857459e9aSMatthew G Knepley   PetscFunctionBegin;
2194f572ea9SToby Isaac   if (pStart) PetscAssertPointer(pStart, 3);
2204f572ea9SToby Isaac   if (pEnd) PetscAssertPointer(pEnd, 4);
2219566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumCells(dm, NULL, NULL, NULL, &nC));
2229566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV));
2239566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF));
22457459e9aSMatthew G Knepley   if (height == 0) {
22557459e9aSMatthew G Knepley     /* Cells */
2268865f1eaSKarl Rupp     if (pStart) *pStart = 0;
2278865f1eaSKarl Rupp     if (pEnd) *pEnd = nC;
22857459e9aSMatthew G Knepley   } else if (height == 1) {
22957459e9aSMatthew G Knepley     /* Faces */
2308865f1eaSKarl Rupp     if (pStart) *pStart = nC + nV;
2318865f1eaSKarl Rupp     if (pEnd) *pEnd = nC + nV + nXF + nYF + nZF;
23257459e9aSMatthew G Knepley   } else if (height == dim) {
23357459e9aSMatthew G Knepley     /* Vertices */
2348865f1eaSKarl Rupp     if (pStart) *pStart = nC;
2358865f1eaSKarl Rupp     if (pEnd) *pEnd = nC + nV;
23657459e9aSMatthew G Knepley   } else if (height < 0) {
23757459e9aSMatthew G Knepley     /* All points */
2388865f1eaSKarl Rupp     if (pStart) *pStart = 0;
2398865f1eaSKarl Rupp     if (pEnd) *pEnd = nC + nV + nXF + nYF + nZF;
24063a3b9bcSJacob Faibussowitsch   } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No points of height %" PetscInt_FMT " in the DA", height);
2413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
24257459e9aSMatthew G Knepley }
24357459e9aSMatthew G Knepley 
244d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetDepthStratum(DM dm, PetscInt depth, PetscInt *pStart, PetscInt *pEnd)
245d71ae5a4SJacob Faibussowitsch {
246c73cfb54SMatthew G. Knepley   const PetscInt dim = dm->dim;
2473385ff7eSMatthew G. Knepley   PetscInt       nC, nV, nXF, nYF, nZF;
2483385ff7eSMatthew G. Knepley 
2493385ff7eSMatthew G. Knepley   PetscFunctionBegin;
2504f572ea9SToby Isaac   if (pStart) PetscAssertPointer(pStart, 3);
2514f572ea9SToby Isaac   if (pEnd) PetscAssertPointer(pEnd, 4);
2529566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumCells(dm, NULL, NULL, NULL, &nC));
2539566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV));
2549566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF));
2553385ff7eSMatthew G. Knepley   if (depth == dim) {
2563385ff7eSMatthew G. Knepley     /* Cells */
2573385ff7eSMatthew G. Knepley     if (pStart) *pStart = 0;
2583385ff7eSMatthew G. Knepley     if (pEnd) *pEnd = nC;
2593385ff7eSMatthew G. Knepley   } else if (depth == dim - 1) {
2603385ff7eSMatthew G. Knepley     /* Faces */
2613385ff7eSMatthew G. Knepley     if (pStart) *pStart = nC + nV;
2623385ff7eSMatthew G. Knepley     if (pEnd) *pEnd = nC + nV + nXF + nYF + nZF;
2633385ff7eSMatthew G. Knepley   } else if (depth == 0) {
2643385ff7eSMatthew G. Knepley     /* Vertices */
2653385ff7eSMatthew G. Knepley     if (pStart) *pStart = nC;
2663385ff7eSMatthew G. Knepley     if (pEnd) *pEnd = nC + nV;
2673385ff7eSMatthew G. Knepley   } else if (depth < 0) {
2683385ff7eSMatthew G. Knepley     /* All points */
2693385ff7eSMatthew G. Knepley     if (pStart) *pStart = 0;
2703385ff7eSMatthew G. Knepley     if (pEnd) *pEnd = nC + nV + nXF + nYF + nZF;
27163a3b9bcSJacob Faibussowitsch   } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No points of depth %" PetscInt_FMT " in the DA", depth);
2723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2733385ff7eSMatthew G. Knepley }
2743385ff7eSMatthew G. Knepley 
275d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetVertexCoordinates(DM dm, PetscReal xl, PetscReal xu, PetscReal yl, PetscReal yu, PetscReal zl, PetscReal zu)
276d71ae5a4SJacob Faibussowitsch {
2773385ff7eSMatthew G. Knepley   DM_DA       *da = (DM_DA *)dm->data;
2783385ff7eSMatthew G. Knepley   Vec          coordinates;
2793385ff7eSMatthew G. Knepley   PetscSection section;
2803385ff7eSMatthew G. Knepley   PetscScalar *coords;
2813385ff7eSMatthew G. Knepley   PetscReal    h[3];
2823385ff7eSMatthew G. Knepley   PetscInt     dim, size, M, N, P, nVx, nVy, nVz, nV, vStart, vEnd, v, i, j, k;
2833385ff7eSMatthew G. Knepley 
2843385ff7eSMatthew G. Knepley   PetscFunctionBegin;
285a9a02de4SBarry Smith   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMDA);
2869566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dm, &dim, &M, &N, &P, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
28708401ef6SPierre Jolivet   PetscCheck(dim <= 3, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "The following code only works for dim <= 3");
2883385ff7eSMatthew G. Knepley   h[0] = (xu - xl) / M;
2893385ff7eSMatthew G. Knepley   h[1] = (yu - yl) / N;
2903385ff7eSMatthew G. Knepley   h[2] = (zu - zl) / P;
2919566063dSJacob Faibussowitsch   PetscCall(DMDAGetDepthStratum(dm, 0, &vStart, &vEnd));
2929566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV));
2939566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), &section));
2949566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetNumFields(section, 1));
2959566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetFieldComponents(section, 0, dim));
2969566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(section, vStart, vEnd));
29748a46eb9SPierre Jolivet   for (v = vStart; v < vEnd; ++v) PetscCall(PetscSectionSetDof(section, v, dim));
2989566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(section));
2999566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(section, &size));
3009566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, size, &coordinates));
3019566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates"));
3029566063dSJacob Faibussowitsch   PetscCall(VecGetArray(coordinates, &coords));
3033385ff7eSMatthew G. Knepley   for (k = 0; k < nVz; ++k) {
304e4d69e08SBarry Smith     PetscInt ind[3], d, off;
3053385ff7eSMatthew G. Knepley 
306e4d69e08SBarry Smith     ind[0] = 0;
307e4d69e08SBarry Smith     ind[1] = 0;
308e4d69e08SBarry Smith     ind[2] = k + da->zs;
3093385ff7eSMatthew G. Knepley     for (j = 0; j < nVy; ++j) {
3103385ff7eSMatthew G. Knepley       ind[1] = j + da->ys;
3113385ff7eSMatthew G. Knepley       for (i = 0; i < nVx; ++i) {
3123385ff7eSMatthew G. Knepley         const PetscInt vertex = (k * nVy + j) * nVx + i + vStart;
3133385ff7eSMatthew G. Knepley 
3149566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetOffset(section, vertex, &off));
3153385ff7eSMatthew G. Knepley         ind[0] = i + da->xs;
316ad540459SPierre Jolivet         for (d = 0; d < dim; ++d) coords[off + d] = h[d] * ind[d];
3173385ff7eSMatthew G. Knepley       }
3183385ff7eSMatthew G. Knepley     }
3193385ff7eSMatthew G. Knepley   }
3209566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(coordinates, &coords));
3219566063dSJacob Faibussowitsch   PetscCall(DMSetCoordinateSection(dm, PETSC_DETERMINE, section));
3229566063dSJacob Faibussowitsch   PetscCall(DMSetCoordinatesLocal(dm, coordinates));
3239566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&section));
3249566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&coordinates));
3253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3263385ff7eSMatthew G. Knepley }
3279a800dd8SMatthew G. Knepley 
32847c6ae99SBarry Smith /* ------------------------------------------------------------------- */
32947c6ae99SBarry Smith 
33047c6ae99SBarry Smith /*@C
331dce8aebaSBarry Smith   DMDAGetArray - Gets a work array for a `DMDA`
33247c6ae99SBarry Smith 
333d8d19677SJose E. Roman   Input Parameters:
33412b4a537SBarry Smith + da      - a `DMDA`
33547c6ae99SBarry Smith - ghosted - do you want arrays for the ghosted or nonghosted patch
33647c6ae99SBarry Smith 
3372fe279fdSBarry Smith   Output Parameter:
33847c6ae99SBarry Smith . vptr - array data structured
33947c6ae99SBarry Smith 
34047c6ae99SBarry Smith   Level: advanced
34147c6ae99SBarry Smith 
34212b4a537SBarry Smith   Notes:
343dce8aebaSBarry Smith   The vector values are NOT initialized and may have garbage in them, so you may need
344dce8aebaSBarry Smith   to zero them.
34547c6ae99SBarry Smith 
34612b4a537SBarry Smith   Use `DMDARestoreArray()` to return the array
34712b4a537SBarry Smith 
34812b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDARestoreArray()`
34947c6ae99SBarry Smith @*/
350d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetArray(DM da, PetscBool ghosted, void *vptr)
351d71ae5a4SJacob Faibussowitsch {
35247c6ae99SBarry Smith   PetscInt j, i, xs, ys, xm, ym, zs, zm;
35347c6ae99SBarry Smith   char    *iarray_start;
35447c6ae99SBarry Smith   void   **iptr = (void **)vptr;
35547c6ae99SBarry Smith   DM_DA   *dd   = (DM_DA *)da->data;
35647c6ae99SBarry Smith 
35747c6ae99SBarry Smith   PetscFunctionBegin;
358a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
35947c6ae99SBarry Smith   if (ghosted) {
360aa219208SBarry Smith     for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) {
36147c6ae99SBarry Smith       if (dd->arrayghostedin[i]) {
36247c6ae99SBarry Smith         *iptr                 = dd->arrayghostedin[i];
36347c6ae99SBarry Smith         iarray_start          = (char *)dd->startghostedin[i];
3640298fd71SBarry Smith         dd->arrayghostedin[i] = NULL;
3650298fd71SBarry Smith         dd->startghostedin[i] = NULL;
36647c6ae99SBarry Smith 
36747c6ae99SBarry Smith         goto done;
36847c6ae99SBarry Smith       }
36947c6ae99SBarry Smith     }
37047c6ae99SBarry Smith     xs = dd->Xs;
37147c6ae99SBarry Smith     ys = dd->Ys;
37247c6ae99SBarry Smith     zs = dd->Zs;
37347c6ae99SBarry Smith     xm = dd->Xe - dd->Xs;
37447c6ae99SBarry Smith     ym = dd->Ye - dd->Ys;
37547c6ae99SBarry Smith     zm = dd->Ze - dd->Zs;
37647c6ae99SBarry Smith   } else {
377aa219208SBarry Smith     for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) {
37847c6ae99SBarry Smith       if (dd->arrayin[i]) {
37947c6ae99SBarry Smith         *iptr          = dd->arrayin[i];
38047c6ae99SBarry Smith         iarray_start   = (char *)dd->startin[i];
3810298fd71SBarry Smith         dd->arrayin[i] = NULL;
3820298fd71SBarry Smith         dd->startin[i] = NULL;
38347c6ae99SBarry Smith 
38447c6ae99SBarry Smith         goto done;
38547c6ae99SBarry Smith       }
38647c6ae99SBarry Smith     }
38747c6ae99SBarry Smith     xs = dd->xs;
38847c6ae99SBarry Smith     ys = dd->ys;
38947c6ae99SBarry Smith     zs = dd->zs;
39047c6ae99SBarry Smith     xm = dd->xe - dd->xs;
39147c6ae99SBarry Smith     ym = dd->ye - dd->ys;
39247c6ae99SBarry Smith     zm = dd->ze - dd->zs;
39347c6ae99SBarry Smith   }
39447c6ae99SBarry Smith 
395c73cfb54SMatthew G. Knepley   switch (da->dim) {
39647c6ae99SBarry Smith   case 1: {
39747c6ae99SBarry Smith     void *ptr;
39847c6ae99SBarry Smith 
3999566063dSJacob Faibussowitsch     PetscCall(PetscMalloc(xm * sizeof(PetscScalar), &iarray_start));
40047c6ae99SBarry Smith 
40147c6ae99SBarry Smith     ptr   = (void *)(iarray_start - xs * sizeof(PetscScalar));
40247c6ae99SBarry Smith     *iptr = (void *)ptr;
4038865f1eaSKarl Rupp     break;
4048865f1eaSKarl Rupp   }
40547c6ae99SBarry Smith   case 2: {
40647c6ae99SBarry Smith     void **ptr;
40747c6ae99SBarry Smith 
4089566063dSJacob Faibussowitsch     PetscCall(PetscMalloc((ym + 1) * sizeof(void *) + xm * ym * sizeof(PetscScalar), &iarray_start));
40947c6ae99SBarry Smith 
41047c6ae99SBarry Smith     ptr = (void **)(iarray_start + xm * ym * sizeof(PetscScalar) - ys * sizeof(void *));
4118865f1eaSKarl Rupp     for (j = ys; j < ys + ym; j++) ptr[j] = iarray_start + sizeof(PetscScalar) * (xm * (j - ys) - xs);
41247c6ae99SBarry Smith     *iptr = (void *)ptr;
4138865f1eaSKarl Rupp     break;
4148865f1eaSKarl Rupp   }
41547c6ae99SBarry Smith   case 3: {
41647c6ae99SBarry Smith     void ***ptr, **bptr;
41747c6ae99SBarry Smith 
4189566063dSJacob Faibussowitsch     PetscCall(PetscMalloc((zm + 1) * sizeof(void **) + (ym * zm + 1) * sizeof(void *) + xm * ym * zm * sizeof(PetscScalar), &iarray_start));
41947c6ae99SBarry Smith 
42047c6ae99SBarry Smith     ptr  = (void ***)(iarray_start + xm * ym * zm * sizeof(PetscScalar) - zs * sizeof(void *));
42147c6ae99SBarry Smith     bptr = (void **)(iarray_start + xm * ym * zm * sizeof(PetscScalar) + zm * sizeof(void **));
4228865f1eaSKarl Rupp     for (i = zs; i < zs + zm; i++) ptr[i] = bptr + ((i - zs) * ym - ys);
42347c6ae99SBarry Smith     for (i = zs; i < zs + zm; i++) {
424ad540459SPierre Jolivet       for (j = ys; j < ys + ym; j++) ptr[i][j] = iarray_start + sizeof(PetscScalar) * (xm * ym * (i - zs) + xm * (j - ys) - xs);
42547c6ae99SBarry Smith     }
42647c6ae99SBarry Smith     *iptr = (void *)ptr;
4278865f1eaSKarl Rupp     break;
4288865f1eaSKarl Rupp   }
429d71ae5a4SJacob Faibussowitsch   default:
430d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Dimension %" PetscInt_FMT " not supported", da->dim);
43147c6ae99SBarry Smith   }
43247c6ae99SBarry Smith 
43347c6ae99SBarry Smith done:
43447c6ae99SBarry Smith   /* add arrays to the checked out list */
43547c6ae99SBarry Smith   if (ghosted) {
436aa219208SBarry Smith     for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) {
43747c6ae99SBarry Smith       if (!dd->arrayghostedout[i]) {
43847c6ae99SBarry Smith         dd->arrayghostedout[i] = *iptr;
43947c6ae99SBarry Smith         dd->startghostedout[i] = iarray_start;
44047c6ae99SBarry Smith         break;
44147c6ae99SBarry Smith       }
44247c6ae99SBarry Smith     }
44347c6ae99SBarry Smith   } else {
444aa219208SBarry Smith     for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) {
44547c6ae99SBarry Smith       if (!dd->arrayout[i]) {
44647c6ae99SBarry Smith         dd->arrayout[i] = *iptr;
44747c6ae99SBarry Smith         dd->startout[i] = iarray_start;
44847c6ae99SBarry Smith         break;
44947c6ae99SBarry Smith       }
45047c6ae99SBarry Smith     }
45147c6ae99SBarry Smith   }
4523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
45347c6ae99SBarry Smith }
45447c6ae99SBarry Smith 
45547c6ae99SBarry Smith /*@C
45612b4a537SBarry Smith   DMDARestoreArray - Restores an array for a `DMDA` obtained with  `DMDAGetArray()`
45747c6ae99SBarry Smith 
458d8d19677SJose E. Roman   Input Parameters:
45947c6ae99SBarry Smith + da      - information about my local patch
46047c6ae99SBarry Smith . ghosted - do you want arrays for the ghosted or nonghosted patch
461cb004a26SBarry Smith - vptr    - array data structured
46247c6ae99SBarry Smith 
46347c6ae99SBarry Smith   Level: advanced
46447c6ae99SBarry Smith 
46512b4a537SBarry Smith .seealso: [](sec_struct), `DM`, `DMDA`, `DMDAGetArray()`
46647c6ae99SBarry Smith @*/
467d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDARestoreArray(DM da, PetscBool ghosted, void *vptr)
468d71ae5a4SJacob Faibussowitsch {
46947c6ae99SBarry Smith   PetscInt i;
470ea78f98cSLisandro Dalcin   void   **iptr = (void **)vptr, *iarray_start = NULL;
47147c6ae99SBarry Smith   DM_DA   *dd = (DM_DA *)da->data;
47247c6ae99SBarry Smith 
47347c6ae99SBarry Smith   PetscFunctionBegin;
474a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
47547c6ae99SBarry Smith   if (ghosted) {
476aa219208SBarry Smith     for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) {
47747c6ae99SBarry Smith       if (dd->arrayghostedout[i] == *iptr) {
47847c6ae99SBarry Smith         iarray_start           = dd->startghostedout[i];
4790298fd71SBarry Smith         dd->arrayghostedout[i] = NULL;
4800298fd71SBarry Smith         dd->startghostedout[i] = NULL;
48147c6ae99SBarry Smith         break;
48247c6ae99SBarry Smith       }
48347c6ae99SBarry Smith     }
484aa219208SBarry Smith     for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) {
48547c6ae99SBarry Smith       if (!dd->arrayghostedin[i]) {
48647c6ae99SBarry Smith         dd->arrayghostedin[i] = *iptr;
48747c6ae99SBarry Smith         dd->startghostedin[i] = iarray_start;
48847c6ae99SBarry Smith         break;
48947c6ae99SBarry Smith       }
49047c6ae99SBarry Smith     }
49147c6ae99SBarry Smith   } else {
492aa219208SBarry Smith     for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) {
49347c6ae99SBarry Smith       if (dd->arrayout[i] == *iptr) {
49447c6ae99SBarry Smith         iarray_start    = dd->startout[i];
4950298fd71SBarry Smith         dd->arrayout[i] = NULL;
4960298fd71SBarry Smith         dd->startout[i] = NULL;
49747c6ae99SBarry Smith         break;
49847c6ae99SBarry Smith       }
49947c6ae99SBarry Smith     }
500aa219208SBarry Smith     for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) {
50147c6ae99SBarry Smith       if (!dd->arrayin[i]) {
50247c6ae99SBarry Smith         dd->arrayin[i] = *iptr;
50347c6ae99SBarry Smith         dd->startin[i] = iarray_start;
50447c6ae99SBarry Smith         break;
50547c6ae99SBarry Smith       }
50647c6ae99SBarry Smith     }
50747c6ae99SBarry Smith   }
5083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
50947c6ae99SBarry Smith }
510