xref: /petsc/src/dm/impls/da/dalocal.c (revision dce8aeba1c9b69b19f651c53d8a6b674bd7e9cbd)
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 */
15d1e78c4fSBarry Smith #if defined(PETSC_HAVE_MATLAB)
16c6db04a5SJed Brown   #include <engine.h> /* MATLAB include file */
17c6db04a5SJed Brown   #include <mex.h>    /* MATLAB include file */
18d71ae5a4SJacob Faibussowitsch static PetscErrorCode VecMatlabEnginePut_DA2d(PetscObject obj, void *mengine)
19d71ae5a4SJacob Faibussowitsch {
2047c6ae99SBarry Smith   PetscInt     n, m;
2147c6ae99SBarry Smith   Vec          vec = (Vec)obj;
2247c6ae99SBarry Smith   PetscScalar *array;
2347c6ae99SBarry Smith   mxArray     *mat;
249a42bb27SBarry Smith   DM           da;
2547c6ae99SBarry Smith 
2647c6ae99SBarry Smith   PetscFunctionBegin;
279566063dSJacob Faibussowitsch   PetscCall(VecGetDM(vec, &da));
287a8be351SBarry Smith   PetscCheck(da, PetscObjectComm((PetscObject)vec), PETSC_ERR_ARG_WRONGSTATE, "Vector not associated with a DMDA");
299566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, 0, 0, 0, &m, &n, 0));
3047c6ae99SBarry Smith 
319566063dSJacob Faibussowitsch   PetscCall(VecGetArray(vec, &array));
3247c6ae99SBarry Smith   #if !defined(PETSC_USE_COMPLEX)
3347c6ae99SBarry Smith   mat = mxCreateDoubleMatrix(m, n, mxREAL);
3447c6ae99SBarry Smith   #else
3547c6ae99SBarry Smith   mat = mxCreateDoubleMatrix(m, n, mxCOMPLEX);
3647c6ae99SBarry Smith   #endif
379566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(mxGetPr(mat), array, n * m));
389566063dSJacob Faibussowitsch   PetscCall(PetscObjectName(obj));
3947c6ae99SBarry Smith   engPutVariable((Engine *)mengine, obj->name, mat);
4047c6ae99SBarry Smith 
419566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(vec, &array));
4247c6ae99SBarry Smith   PetscFunctionReturn(0);
4347c6ae99SBarry Smith }
4447c6ae99SBarry Smith #endif
4547c6ae99SBarry Smith 
46d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateLocalVector_DA(DM da, Vec *g)
47d71ae5a4SJacob Faibussowitsch {
4847c6ae99SBarry Smith   DM_DA *dd = (DM_DA *)da->data;
4947c6ae99SBarry Smith 
5047c6ae99SBarry Smith   PetscFunctionBegin;
5147c6ae99SBarry Smith   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
5247c6ae99SBarry Smith   PetscValidPointer(g, 2);
539566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_SELF, g));
549566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(*g, dd->nlocal, PETSC_DETERMINE));
559566063dSJacob Faibussowitsch   PetscCall(VecSetBlockSize(*g, dd->w));
569566063dSJacob Faibussowitsch   PetscCall(VecSetType(*g, da->vectype));
5774427ab1SRichard Tran Mills   if (dd->nlocal < da->bind_below) {
589566063dSJacob Faibussowitsch     PetscCall(VecSetBindingPropagates(*g, PETSC_TRUE));
599566063dSJacob Faibussowitsch     PetscCall(VecBindToCPU(*g, PETSC_TRUE));
6074427ab1SRichard Tran Mills   }
619566063dSJacob Faibussowitsch   PetscCall(VecSetDM(*g, da));
62d1e78c4fSBarry Smith #if defined(PETSC_HAVE_MATLAB)
6348a46eb9SPierre Jolivet   if (dd->w == 1 && da->dim == 2) PetscCall(PetscObjectComposeFunction((PetscObject)*g, "PetscMatlabEnginePut_C", VecMatlabEnginePut_DA2d));
6447c6ae99SBarry Smith #endif
6547c6ae99SBarry Smith   PetscFunctionReturn(0);
6647c6ae99SBarry Smith }
6747c6ae99SBarry Smith 
68939f6067SMatthew G. Knepley /*@
69*dce8aebaSBarry Smith   DMDAGetNumCells - Get the number of cells in the local piece of the `DMDA`. This includes ghost cells.
70939f6067SMatthew G. Knepley 
71939f6067SMatthew G. Knepley   Input Parameter:
72*dce8aebaSBarry Smith . dm - The `DMDA` object
73939f6067SMatthew G. Knepley 
74939f6067SMatthew G. Knepley   Output Parameters:
75939f6067SMatthew G. Knepley + numCellsX - The number of local cells in the x-direction
76939f6067SMatthew G. Knepley . numCellsY - The number of local cells in the y-direction
77939f6067SMatthew G. Knepley . numCellsZ - The number of local cells in the z-direction
78939f6067SMatthew G. Knepley - numCells - The number of local cells
79939f6067SMatthew G. Knepley 
80939f6067SMatthew G. Knepley   Level: developer
81939f6067SMatthew G. Knepley 
82*dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetCellPoint()`
83939f6067SMatthew G. Knepley @*/
84d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetNumCells(DM dm, PetscInt *numCellsX, PetscInt *numCellsY, PetscInt *numCellsZ, PetscInt *numCells)
85d71ae5a4SJacob Faibussowitsch {
8657459e9aSMatthew G Knepley   DM_DA         *da  = (DM_DA *)dm->data;
87c73cfb54SMatthew G. Knepley   const PetscInt dim = dm->dim;
8857459e9aSMatthew G Knepley   const PetscInt mx = (da->Xe - da->Xs) / da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs;
8957459e9aSMatthew G Knepley   const PetscInt nC = (mx) * (dim > 1 ? (my) * (dim > 2 ? (mz) : 1) : 1);
9057459e9aSMatthew G Knepley 
9157459e9aSMatthew G Knepley   PetscFunctionBegin;
92a9a02de4SBarry Smith   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMDA);
933582350cSMatthew G. Knepley   if (numCellsX) {
943582350cSMatthew G. Knepley     PetscValidIntPointer(numCellsX, 2);
953582350cSMatthew G. Knepley     *numCellsX = mx;
963582350cSMatthew G. Knepley   }
973582350cSMatthew G. Knepley   if (numCellsY) {
988135c375SStefano Zampini     PetscValidIntPointer(numCellsY, 3);
993582350cSMatthew G. Knepley     *numCellsY = my;
1003582350cSMatthew G. Knepley   }
1013582350cSMatthew G. Knepley   if (numCellsZ) {
1028135c375SStefano Zampini     PetscValidIntPointer(numCellsZ, 4);
1033582350cSMatthew G. Knepley     *numCellsZ = mz;
1043582350cSMatthew G. Knepley   }
10557459e9aSMatthew G Knepley   if (numCells) {
1063582350cSMatthew G. Knepley     PetscValidIntPointer(numCells, 5);
10757459e9aSMatthew G Knepley     *numCells = nC;
10857459e9aSMatthew G Knepley   }
10957459e9aSMatthew G Knepley   PetscFunctionReturn(0);
11057459e9aSMatthew G Knepley }
11157459e9aSMatthew G Knepley 
112d6401b93SMatthew G. Knepley /*@
113*dce8aebaSBarry Smith   DMDAGetCellPoint - Get the DM point corresponding to the tuple (i, j, k) in the `DMDA`
114d6401b93SMatthew G. Knepley 
115d6401b93SMatthew G. Knepley   Input Parameters:
116*dce8aebaSBarry Smith + dm - The `DMDA` object
117d6401b93SMatthew G. Knepley - i,j,k - The global indices for the cell
118d6401b93SMatthew G. Knepley 
119d6401b93SMatthew G. Knepley   Output Parameters:
120*dce8aebaSBarry Smith . point - The local `DM` point
121d6401b93SMatthew G. Knepley 
122d6401b93SMatthew G. Knepley   Level: developer
123d6401b93SMatthew G. Knepley 
124*dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetNumCells()`
125d6401b93SMatthew G. Knepley @*/
126d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetCellPoint(DM dm, PetscInt i, PetscInt j, PetscInt k, PetscInt *point)
127d71ae5a4SJacob Faibussowitsch {
128c73cfb54SMatthew G. Knepley   const PetscInt dim = dm->dim;
129d6401b93SMatthew G. Knepley   DMDALocalInfo  info;
130d6401b93SMatthew G. Knepley 
131d6401b93SMatthew G. Knepley   PetscFunctionBegin;
132a9a02de4SBarry Smith   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMDA);
133d6401b93SMatthew G. Knepley   PetscValidIntPointer(point, 5);
1349566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(dm, &info));
13563a3b9bcSJacob 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);
13663a3b9bcSJacob 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);
13763a3b9bcSJacob 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);
138d6401b93SMatthew G. Knepley   *point = i + (dim > 1 ? (j + (dim > 2 ? k * info.gym : 0)) * info.gxm : 0);
139d6401b93SMatthew G. Knepley   PetscFunctionReturn(0);
140d6401b93SMatthew G. Knepley }
141d6401b93SMatthew G. Knepley 
142d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetNumVertices(DM dm, PetscInt *numVerticesX, PetscInt *numVerticesY, PetscInt *numVerticesZ, PetscInt *numVertices)
143d71ae5a4SJacob Faibussowitsch {
14457459e9aSMatthew G Knepley   DM_DA         *da  = (DM_DA *)dm->data;
145c73cfb54SMatthew G. Knepley   const PetscInt dim = dm->dim;
14657459e9aSMatthew G Knepley   const PetscInt mx = (da->Xe - da->Xs) / da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs;
14757459e9aSMatthew G Knepley   const PetscInt nVx = mx + 1;
14857459e9aSMatthew G Knepley   const PetscInt nVy = dim > 1 ? (my + 1) : 1;
14957459e9aSMatthew G Knepley   const PetscInt nVz = dim > 2 ? (mz + 1) : 1;
15057459e9aSMatthew G Knepley   const PetscInt nV  = nVx * nVy * nVz;
15157459e9aSMatthew G Knepley 
15257459e9aSMatthew G Knepley   PetscFunctionBegin;
15357459e9aSMatthew G Knepley   if (numVerticesX) {
15457459e9aSMatthew G Knepley     PetscValidIntPointer(numVerticesX, 2);
15557459e9aSMatthew G Knepley     *numVerticesX = nVx;
15657459e9aSMatthew G Knepley   }
15757459e9aSMatthew G Knepley   if (numVerticesY) {
15857459e9aSMatthew G Knepley     PetscValidIntPointer(numVerticesY, 3);
15957459e9aSMatthew G Knepley     *numVerticesY = nVy;
16057459e9aSMatthew G Knepley   }
16157459e9aSMatthew G Knepley   if (numVerticesZ) {
16257459e9aSMatthew G Knepley     PetscValidIntPointer(numVerticesZ, 4);
16357459e9aSMatthew G Knepley     *numVerticesZ = nVz;
16457459e9aSMatthew G Knepley   }
16557459e9aSMatthew G Knepley   if (numVertices) {
16657459e9aSMatthew G Knepley     PetscValidIntPointer(numVertices, 5);
16757459e9aSMatthew G Knepley     *numVertices = nV;
16857459e9aSMatthew G Knepley   }
16957459e9aSMatthew G Knepley   PetscFunctionReturn(0);
17057459e9aSMatthew G Knepley }
17157459e9aSMatthew G Knepley 
172d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetNumFaces(DM dm, PetscInt *numXFacesX, PetscInt *numXFaces, PetscInt *numYFacesY, PetscInt *numYFaces, PetscInt *numZFacesZ, PetscInt *numZFaces)
173d71ae5a4SJacob Faibussowitsch {
17457459e9aSMatthew G Knepley   DM_DA         *da  = (DM_DA *)dm->data;
175c73cfb54SMatthew G. Knepley   const PetscInt dim = dm->dim;
17657459e9aSMatthew G Knepley   const PetscInt mx = (da->Xe - da->Xs) / da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs;
17757459e9aSMatthew G Knepley   const PetscInt nxF = (dim > 1 ? (my) * (dim > 2 ? (mz) : 1) : 1);
17857459e9aSMatthew G Knepley   const PetscInt nXF = (mx + 1) * nxF;
17957459e9aSMatthew G Knepley   const PetscInt nyF = mx * (dim > 2 ? mz : 1);
18057459e9aSMatthew G Knepley   const PetscInt nYF = dim > 1 ? (my + 1) * nyF : 0;
18157459e9aSMatthew G Knepley   const PetscInt nzF = mx * (dim > 1 ? my : 0);
18257459e9aSMatthew G Knepley   const PetscInt nZF = dim > 2 ? (mz + 1) * nzF : 0;
18357459e9aSMatthew G Knepley 
18457459e9aSMatthew G Knepley   PetscFunctionBegin;
18557459e9aSMatthew G Knepley   if (numXFacesX) {
18657459e9aSMatthew G Knepley     PetscValidIntPointer(numXFacesX, 2);
18757459e9aSMatthew G Knepley     *numXFacesX = nxF;
18857459e9aSMatthew G Knepley   }
18957459e9aSMatthew G Knepley   if (numXFaces) {
19057459e9aSMatthew G Knepley     PetscValidIntPointer(numXFaces, 3);
19157459e9aSMatthew G Knepley     *numXFaces = nXF;
19257459e9aSMatthew G Knepley   }
19357459e9aSMatthew G Knepley   if (numYFacesY) {
19457459e9aSMatthew G Knepley     PetscValidIntPointer(numYFacesY, 4);
19557459e9aSMatthew G Knepley     *numYFacesY = nyF;
19657459e9aSMatthew G Knepley   }
19757459e9aSMatthew G Knepley   if (numYFaces) {
19857459e9aSMatthew G Knepley     PetscValidIntPointer(numYFaces, 5);
19957459e9aSMatthew G Knepley     *numYFaces = nYF;
20057459e9aSMatthew G Knepley   }
20157459e9aSMatthew G Knepley   if (numZFacesZ) {
20257459e9aSMatthew G Knepley     PetscValidIntPointer(numZFacesZ, 6);
20357459e9aSMatthew G Knepley     *numZFacesZ = nzF;
20457459e9aSMatthew G Knepley   }
20557459e9aSMatthew G Knepley   if (numZFaces) {
20657459e9aSMatthew G Knepley     PetscValidIntPointer(numZFaces, 7);
20757459e9aSMatthew G Knepley     *numZFaces = nZF;
20857459e9aSMatthew G Knepley   }
20957459e9aSMatthew G Knepley   PetscFunctionReturn(0);
21057459e9aSMatthew G Knepley }
21157459e9aSMatthew G Knepley 
212d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetHeightStratum(DM dm, PetscInt height, PetscInt *pStart, PetscInt *pEnd)
213d71ae5a4SJacob Faibussowitsch {
214c73cfb54SMatthew G. Knepley   const PetscInt dim = dm->dim;
21557459e9aSMatthew G Knepley   PetscInt       nC, nV, nXF, nYF, nZF;
21657459e9aSMatthew G Knepley 
21757459e9aSMatthew G Knepley   PetscFunctionBegin;
2188865f1eaSKarl Rupp   if (pStart) PetscValidIntPointer(pStart, 3);
2198865f1eaSKarl Rupp   if (pEnd) PetscValidIntPointer(pEnd, 4);
2209566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumCells(dm, NULL, NULL, NULL, &nC));
2219566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV));
2229566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF));
22357459e9aSMatthew G Knepley   if (height == 0) {
22457459e9aSMatthew G Knepley     /* Cells */
2258865f1eaSKarl Rupp     if (pStart) *pStart = 0;
2268865f1eaSKarl Rupp     if (pEnd) *pEnd = nC;
22757459e9aSMatthew G Knepley   } else if (height == 1) {
22857459e9aSMatthew G Knepley     /* Faces */
2298865f1eaSKarl Rupp     if (pStart) *pStart = nC + nV;
2308865f1eaSKarl Rupp     if (pEnd) *pEnd = nC + nV + nXF + nYF + nZF;
23157459e9aSMatthew G Knepley   } else if (height == dim) {
23257459e9aSMatthew G Knepley     /* Vertices */
2338865f1eaSKarl Rupp     if (pStart) *pStart = nC;
2348865f1eaSKarl Rupp     if (pEnd) *pEnd = nC + nV;
23557459e9aSMatthew G Knepley   } else if (height < 0) {
23657459e9aSMatthew G Knepley     /* All points */
2378865f1eaSKarl Rupp     if (pStart) *pStart = 0;
2388865f1eaSKarl Rupp     if (pEnd) *pEnd = nC + nV + nXF + nYF + nZF;
23963a3b9bcSJacob Faibussowitsch   } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No points of height %" PetscInt_FMT " in the DA", height);
24057459e9aSMatthew G Knepley   PetscFunctionReturn(0);
24157459e9aSMatthew G Knepley }
24257459e9aSMatthew G Knepley 
243d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetDepthStratum(DM dm, PetscInt depth, PetscInt *pStart, PetscInt *pEnd)
244d71ae5a4SJacob Faibussowitsch {
245c73cfb54SMatthew G. Knepley   const PetscInt dim = dm->dim;
2463385ff7eSMatthew G. Knepley   PetscInt       nC, nV, nXF, nYF, nZF;
2473385ff7eSMatthew G. Knepley 
2483385ff7eSMatthew G. Knepley   PetscFunctionBegin;
2493385ff7eSMatthew G. Knepley   if (pStart) PetscValidIntPointer(pStart, 3);
2503385ff7eSMatthew G. Knepley   if (pEnd) PetscValidIntPointer(pEnd, 4);
2519566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumCells(dm, NULL, NULL, NULL, &nC));
2529566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV));
2539566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF));
2543385ff7eSMatthew G. Knepley   if (depth == dim) {
2553385ff7eSMatthew G. Knepley     /* Cells */
2563385ff7eSMatthew G. Knepley     if (pStart) *pStart = 0;
2573385ff7eSMatthew G. Knepley     if (pEnd) *pEnd = nC;
2583385ff7eSMatthew G. Knepley   } else if (depth == dim - 1) {
2593385ff7eSMatthew G. Knepley     /* Faces */
2603385ff7eSMatthew G. Knepley     if (pStart) *pStart = nC + nV;
2613385ff7eSMatthew G. Knepley     if (pEnd) *pEnd = nC + nV + nXF + nYF + nZF;
2623385ff7eSMatthew G. Knepley   } else if (depth == 0) {
2633385ff7eSMatthew G. Knepley     /* Vertices */
2643385ff7eSMatthew G. Knepley     if (pStart) *pStart = nC;
2653385ff7eSMatthew G. Knepley     if (pEnd) *pEnd = nC + nV;
2663385ff7eSMatthew G. Knepley   } else if (depth < 0) {
2673385ff7eSMatthew G. Knepley     /* All points */
2683385ff7eSMatthew G. Knepley     if (pStart) *pStart = 0;
2693385ff7eSMatthew G. Knepley     if (pEnd) *pEnd = nC + nV + nXF + nYF + nZF;
27063a3b9bcSJacob Faibussowitsch   } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No points of depth %" PetscInt_FMT " in the DA", depth);
2713385ff7eSMatthew G. Knepley   PetscFunctionReturn(0);
2723385ff7eSMatthew G. Knepley }
2733385ff7eSMatthew G. Knepley 
274d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetConeSize(DM dm, PetscInt p, PetscInt *coneSize)
275d71ae5a4SJacob Faibussowitsch {
276c73cfb54SMatthew G. Knepley   const PetscInt dim = dm->dim;
2773385ff7eSMatthew G. Knepley   PetscInt       nC, nV, nXF, nYF, nZF;
2783385ff7eSMatthew G. Knepley 
2793385ff7eSMatthew G. Knepley   PetscFunctionBegin;
2803385ff7eSMatthew G. Knepley   *coneSize = 0;
2819566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumCells(dm, NULL, NULL, NULL, &nC));
2829566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV));
2839566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF));
2843385ff7eSMatthew G. Knepley   switch (dim) {
2853385ff7eSMatthew G. Knepley   case 2:
2863385ff7eSMatthew G. Knepley     if (p >= 0) {
2873385ff7eSMatthew G. Knepley       if (p < nC) {
2883385ff7eSMatthew G. Knepley         *coneSize = 4;
2893385ff7eSMatthew G. Knepley       } else if (p < nC + nV) {
2903385ff7eSMatthew G. Knepley         *coneSize = 0;
2913385ff7eSMatthew G. Knepley       } else if (p < nC + nV + nXF + nYF + nZF) {
2923385ff7eSMatthew G. Knepley         *coneSize = 2;
29363a3b9bcSJacob 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);
29463a3b9bcSJacob Faibussowitsch     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative point %" PetscInt_FMT " is invalid", p);
2953385ff7eSMatthew G. Knepley     break;
296d71ae5a4SJacob Faibussowitsch   case 3:
297d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do 3D");
2983385ff7eSMatthew G. Knepley   }
2993385ff7eSMatthew G. Knepley   PetscFunctionReturn(0);
3003385ff7eSMatthew G. Knepley }
3013385ff7eSMatthew G. Knepley 
302d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetCone(DM dm, PetscInt p, PetscInt *cone[])
303d71ae5a4SJacob Faibussowitsch {
304c73cfb54SMatthew G. Knepley   const PetscInt dim = dm->dim;
3053385ff7eSMatthew G. Knepley   PetscInt       nCx, nCy, nCz, nC, nVx, nVy, nVz, nV, nxF, nyF, nzF, nXF, nYF, nZF;
3063385ff7eSMatthew G. Knepley 
3073385ff7eSMatthew G. Knepley   PetscFunctionBegin;
3089566063dSJacob Faibussowitsch   if (!cone) PetscCall(DMGetWorkArray(dm, 6, MPIU_INT, cone));
3099566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumCells(dm, &nCx, &nCy, &nCz, &nC));
3109566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV));
3119566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF));
3123385ff7eSMatthew G. Knepley   switch (dim) {
3133385ff7eSMatthew G. Knepley   case 2:
3143385ff7eSMatthew G. Knepley     if (p >= 0) {
3153385ff7eSMatthew G. Knepley       if (p < nC) {
3163385ff7eSMatthew G. Knepley         const PetscInt cy = p / nCx;
3173385ff7eSMatthew G. Knepley         const PetscInt cx = p % nCx;
3183385ff7eSMatthew G. Knepley 
3193385ff7eSMatthew G. Knepley         (*cone)[0] = cy * nxF + cx + nC + nV;
3203385ff7eSMatthew G. Knepley         (*cone)[1] = cx * nyF + cy + nyF + nC + nV + nXF;
3213385ff7eSMatthew G. Knepley         (*cone)[2] = cy * nxF + cx + nxF + nC + nV;
3223385ff7eSMatthew G. Knepley         (*cone)[3] = cx * nyF + cy + nC + nV + nXF;
3233385ff7eSMatthew G. Knepley         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do cell cones");
3243385ff7eSMatthew G. Knepley       } else if (p < nC + nV) {
3253385ff7eSMatthew G. Knepley       } else if (p < nC + nV + nXF) {
3263385ff7eSMatthew G. Knepley         const PetscInt fy = (p - nC + nV) / nxF;
3273385ff7eSMatthew G. Knepley         const PetscInt fx = (p - nC + nV) % nxF;
3283385ff7eSMatthew G. Knepley 
3293385ff7eSMatthew G. Knepley         (*cone)[0] = fy * nVx + fx + nC;
3303385ff7eSMatthew G. Knepley         (*cone)[1] = fy * nVx + fx + 1 + nC;
3313385ff7eSMatthew G. Knepley       } else if (p < nC + nV + nXF + nYF) {
3323385ff7eSMatthew G. Knepley         const PetscInt fx = (p - nC + nV + nXF) / nyF;
3333385ff7eSMatthew G. Knepley         const PetscInt fy = (p - nC + nV + nXF) % nyF;
3343385ff7eSMatthew G. Knepley 
3353385ff7eSMatthew G. Knepley         (*cone)[0] = fy * nVx + fx + nC;
3363385ff7eSMatthew G. Knepley         (*cone)[1] = fy * nVx + fx + nVx + nC;
33763a3b9bcSJacob 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);
33863a3b9bcSJacob Faibussowitsch     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative point %" PetscInt_FMT " is invalid", p);
3393385ff7eSMatthew G. Knepley     break;
340d71ae5a4SJacob Faibussowitsch   case 3:
341d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do 3D");
3423385ff7eSMatthew G. Knepley   }
3433385ff7eSMatthew G. Knepley   PetscFunctionReturn(0);
3443385ff7eSMatthew G. Knepley }
3453385ff7eSMatthew G. Knepley 
346d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDARestoreCone(DM dm, PetscInt p, PetscInt *cone[])
347d71ae5a4SJacob Faibussowitsch {
3483385ff7eSMatthew G. Knepley   PetscFunctionBegin;
3499566063dSJacob Faibussowitsch   PetscCall(DMGetWorkArray(dm, 6, MPIU_INT, cone));
3503385ff7eSMatthew G. Knepley   PetscFunctionReturn(0);
3513385ff7eSMatthew G. Knepley }
3523385ff7eSMatthew G. Knepley 
353d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetVertexCoordinates(DM dm, PetscReal xl, PetscReal xu, PetscReal yl, PetscReal yu, PetscReal zl, PetscReal zu)
354d71ae5a4SJacob Faibussowitsch {
3553385ff7eSMatthew G. Knepley   DM_DA       *da = (DM_DA *)dm->data;
3563385ff7eSMatthew G. Knepley   Vec          coordinates;
3573385ff7eSMatthew G. Knepley   PetscSection section;
3583385ff7eSMatthew G. Knepley   PetscScalar *coords;
3593385ff7eSMatthew G. Knepley   PetscReal    h[3];
3603385ff7eSMatthew G. Knepley   PetscInt     dim, size, M, N, P, nVx, nVy, nVz, nV, vStart, vEnd, v, i, j, k;
3613385ff7eSMatthew G. Knepley 
3623385ff7eSMatthew G. Knepley   PetscFunctionBegin;
363a9a02de4SBarry Smith   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMDA);
3649566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dm, &dim, &M, &N, &P, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
36508401ef6SPierre Jolivet   PetscCheck(dim <= 3, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "The following code only works for dim <= 3");
3663385ff7eSMatthew G. Knepley   h[0] = (xu - xl) / M;
3673385ff7eSMatthew G. Knepley   h[1] = (yu - yl) / N;
3683385ff7eSMatthew G. Knepley   h[2] = (zu - zl) / P;
3699566063dSJacob Faibussowitsch   PetscCall(DMDAGetDepthStratum(dm, 0, &vStart, &vEnd));
3709566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV));
3719566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), &section));
3729566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetNumFields(section, 1));
3739566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetFieldComponents(section, 0, dim));
3749566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(section, vStart, vEnd));
37548a46eb9SPierre Jolivet   for (v = vStart; v < vEnd; ++v) PetscCall(PetscSectionSetDof(section, v, dim));
3769566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(section));
3779566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(section, &size));
3789566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, size, &coordinates));
3799566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates"));
3809566063dSJacob Faibussowitsch   PetscCall(VecGetArray(coordinates, &coords));
3813385ff7eSMatthew G. Knepley   for (k = 0; k < nVz; ++k) {
382e4d69e08SBarry Smith     PetscInt ind[3], d, off;
3833385ff7eSMatthew G. Knepley 
384e4d69e08SBarry Smith     ind[0] = 0;
385e4d69e08SBarry Smith     ind[1] = 0;
386e4d69e08SBarry Smith     ind[2] = k + da->zs;
3873385ff7eSMatthew G. Knepley     for (j = 0; j < nVy; ++j) {
3883385ff7eSMatthew G. Knepley       ind[1] = j + da->ys;
3893385ff7eSMatthew G. Knepley       for (i = 0; i < nVx; ++i) {
3903385ff7eSMatthew G. Knepley         const PetscInt vertex = (k * nVy + j) * nVx + i + vStart;
3913385ff7eSMatthew G. Knepley 
3929566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetOffset(section, vertex, &off));
3933385ff7eSMatthew G. Knepley         ind[0] = i + da->xs;
394ad540459SPierre Jolivet         for (d = 0; d < dim; ++d) coords[off + d] = h[d] * ind[d];
3953385ff7eSMatthew G. Knepley       }
3963385ff7eSMatthew G. Knepley     }
3973385ff7eSMatthew G. Knepley   }
3989566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(coordinates, &coords));
3999566063dSJacob Faibussowitsch   PetscCall(DMSetCoordinateSection(dm, PETSC_DETERMINE, section));
4009566063dSJacob Faibussowitsch   PetscCall(DMSetCoordinatesLocal(dm, coordinates));
4019566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&section));
4029566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&coordinates));
4033385ff7eSMatthew G. Knepley   PetscFunctionReturn(0);
4043385ff7eSMatthew G. Knepley }
4059a800dd8SMatthew G. Knepley 
40647c6ae99SBarry Smith /* ------------------------------------------------------------------- */
40747c6ae99SBarry Smith 
40847c6ae99SBarry Smith /*@C
409*dce8aebaSBarry Smith      DMDAGetArray - Gets a work array for a `DMDA`
41047c6ae99SBarry Smith 
411d8d19677SJose E. Roman     Input Parameters:
41247c6ae99SBarry Smith +    da - information about my local patch
41347c6ae99SBarry Smith -    ghosted - do you want arrays for the ghosted or nonghosted patch
41447c6ae99SBarry Smith 
41547c6ae99SBarry Smith     Output Parameters:
41647c6ae99SBarry Smith .    vptr - array data structured
41747c6ae99SBarry Smith 
41847c6ae99SBarry Smith   Level: advanced
41947c6ae99SBarry Smith 
420*dce8aebaSBarry Smith   Note:
421*dce8aebaSBarry Smith    The vector values are NOT initialized and may have garbage in them, so you may need
422*dce8aebaSBarry Smith    to zero them.
42347c6ae99SBarry Smith 
424*dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDARestoreArray()`
42547c6ae99SBarry Smith @*/
426d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetArray(DM da, PetscBool ghosted, void *vptr)
427d71ae5a4SJacob Faibussowitsch {
42847c6ae99SBarry Smith   PetscInt j, i, xs, ys, xm, ym, zs, zm;
42947c6ae99SBarry Smith   char    *iarray_start;
43047c6ae99SBarry Smith   void   **iptr = (void **)vptr;
43147c6ae99SBarry Smith   DM_DA   *dd   = (DM_DA *)da->data;
43247c6ae99SBarry Smith 
43347c6ae99SBarry Smith   PetscFunctionBegin;
434a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
43547c6ae99SBarry Smith   if (ghosted) {
436aa219208SBarry Smith     for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) {
43747c6ae99SBarry Smith       if (dd->arrayghostedin[i]) {
43847c6ae99SBarry Smith         *iptr                 = dd->arrayghostedin[i];
43947c6ae99SBarry Smith         iarray_start          = (char *)dd->startghostedin[i];
4400298fd71SBarry Smith         dd->arrayghostedin[i] = NULL;
4410298fd71SBarry Smith         dd->startghostedin[i] = NULL;
44247c6ae99SBarry Smith 
44347c6ae99SBarry Smith         goto done;
44447c6ae99SBarry Smith       }
44547c6ae99SBarry Smith     }
44647c6ae99SBarry Smith     xs = dd->Xs;
44747c6ae99SBarry Smith     ys = dd->Ys;
44847c6ae99SBarry Smith     zs = dd->Zs;
44947c6ae99SBarry Smith     xm = dd->Xe - dd->Xs;
45047c6ae99SBarry Smith     ym = dd->Ye - dd->Ys;
45147c6ae99SBarry Smith     zm = dd->Ze - dd->Zs;
45247c6ae99SBarry Smith   } else {
453aa219208SBarry Smith     for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) {
45447c6ae99SBarry Smith       if (dd->arrayin[i]) {
45547c6ae99SBarry Smith         *iptr          = dd->arrayin[i];
45647c6ae99SBarry Smith         iarray_start   = (char *)dd->startin[i];
4570298fd71SBarry Smith         dd->arrayin[i] = NULL;
4580298fd71SBarry Smith         dd->startin[i] = NULL;
45947c6ae99SBarry Smith 
46047c6ae99SBarry Smith         goto done;
46147c6ae99SBarry Smith       }
46247c6ae99SBarry Smith     }
46347c6ae99SBarry Smith     xs = dd->xs;
46447c6ae99SBarry Smith     ys = dd->ys;
46547c6ae99SBarry Smith     zs = dd->zs;
46647c6ae99SBarry Smith     xm = dd->xe - dd->xs;
46747c6ae99SBarry Smith     ym = dd->ye - dd->ys;
46847c6ae99SBarry Smith     zm = dd->ze - dd->zs;
46947c6ae99SBarry Smith   }
47047c6ae99SBarry Smith 
471c73cfb54SMatthew G. Knepley   switch (da->dim) {
47247c6ae99SBarry Smith   case 1: {
47347c6ae99SBarry Smith     void *ptr;
47447c6ae99SBarry Smith 
4759566063dSJacob Faibussowitsch     PetscCall(PetscMalloc(xm * sizeof(PetscScalar), &iarray_start));
47647c6ae99SBarry Smith 
47747c6ae99SBarry Smith     ptr   = (void *)(iarray_start - xs * sizeof(PetscScalar));
47847c6ae99SBarry Smith     *iptr = (void *)ptr;
4798865f1eaSKarl Rupp     break;
4808865f1eaSKarl Rupp   }
48147c6ae99SBarry Smith   case 2: {
48247c6ae99SBarry Smith     void **ptr;
48347c6ae99SBarry Smith 
4849566063dSJacob Faibussowitsch     PetscCall(PetscMalloc((ym + 1) * sizeof(void *) + xm * ym * sizeof(PetscScalar), &iarray_start));
48547c6ae99SBarry Smith 
48647c6ae99SBarry Smith     ptr = (void **)(iarray_start + xm * ym * sizeof(PetscScalar) - ys * sizeof(void *));
4878865f1eaSKarl Rupp     for (j = ys; j < ys + ym; j++) ptr[j] = iarray_start + sizeof(PetscScalar) * (xm * (j - ys) - xs);
48847c6ae99SBarry Smith     *iptr = (void *)ptr;
4898865f1eaSKarl Rupp     break;
4908865f1eaSKarl Rupp   }
49147c6ae99SBarry Smith   case 3: {
49247c6ae99SBarry Smith     void ***ptr, **bptr;
49347c6ae99SBarry Smith 
4949566063dSJacob Faibussowitsch     PetscCall(PetscMalloc((zm + 1) * sizeof(void **) + (ym * zm + 1) * sizeof(void *) + xm * ym * zm * sizeof(PetscScalar), &iarray_start));
49547c6ae99SBarry Smith 
49647c6ae99SBarry Smith     ptr  = (void ***)(iarray_start + xm * ym * zm * sizeof(PetscScalar) - zs * sizeof(void *));
49747c6ae99SBarry Smith     bptr = (void **)(iarray_start + xm * ym * zm * sizeof(PetscScalar) + zm * sizeof(void **));
4988865f1eaSKarl Rupp     for (i = zs; i < zs + zm; i++) ptr[i] = bptr + ((i - zs) * ym - ys);
49947c6ae99SBarry Smith     for (i = zs; i < zs + zm; i++) {
500ad540459SPierre Jolivet       for (j = ys; j < ys + ym; j++) ptr[i][j] = iarray_start + sizeof(PetscScalar) * (xm * ym * (i - zs) + xm * (j - ys) - xs);
50147c6ae99SBarry Smith     }
50247c6ae99SBarry Smith     *iptr = (void *)ptr;
5038865f1eaSKarl Rupp     break;
5048865f1eaSKarl Rupp   }
505d71ae5a4SJacob Faibussowitsch   default:
506d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Dimension %" PetscInt_FMT " not supported", da->dim);
50747c6ae99SBarry Smith   }
50847c6ae99SBarry Smith 
50947c6ae99SBarry Smith done:
51047c6ae99SBarry Smith   /* add arrays to the checked out list */
51147c6ae99SBarry Smith   if (ghosted) {
512aa219208SBarry Smith     for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) {
51347c6ae99SBarry Smith       if (!dd->arrayghostedout[i]) {
51447c6ae99SBarry Smith         dd->arrayghostedout[i] = *iptr;
51547c6ae99SBarry Smith         dd->startghostedout[i] = iarray_start;
51647c6ae99SBarry Smith         break;
51747c6ae99SBarry Smith       }
51847c6ae99SBarry Smith     }
51947c6ae99SBarry Smith   } else {
520aa219208SBarry Smith     for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) {
52147c6ae99SBarry Smith       if (!dd->arrayout[i]) {
52247c6ae99SBarry Smith         dd->arrayout[i] = *iptr;
52347c6ae99SBarry Smith         dd->startout[i] = iarray_start;
52447c6ae99SBarry Smith         break;
52547c6ae99SBarry Smith       }
52647c6ae99SBarry Smith     }
52747c6ae99SBarry Smith   }
52847c6ae99SBarry Smith   PetscFunctionReturn(0);
52947c6ae99SBarry Smith }
53047c6ae99SBarry Smith 
53147c6ae99SBarry Smith /*@C
532*dce8aebaSBarry Smith      DMDARestoreArray - Restores an array of derivative types for a `DMDA`
53347c6ae99SBarry Smith 
534d8d19677SJose E. Roman     Input Parameters:
53547c6ae99SBarry Smith +    da - information about my local patch
53647c6ae99SBarry Smith .    ghosted - do you want arrays for the ghosted or nonghosted patch
537cb004a26SBarry Smith -    vptr - array data structured
53847c6ae99SBarry Smith 
53947c6ae99SBarry Smith      Level: advanced
54047c6ae99SBarry Smith 
541*dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetArray()`
54247c6ae99SBarry Smith @*/
543d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDARestoreArray(DM da, PetscBool ghosted, void *vptr)
544d71ae5a4SJacob Faibussowitsch {
54547c6ae99SBarry Smith   PetscInt i;
546ea78f98cSLisandro Dalcin   void   **iptr = (void **)vptr, *iarray_start = NULL;
54747c6ae99SBarry Smith   DM_DA   *dd = (DM_DA *)da->data;
54847c6ae99SBarry Smith 
54947c6ae99SBarry Smith   PetscFunctionBegin;
550a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
55147c6ae99SBarry Smith   if (ghosted) {
552aa219208SBarry Smith     for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) {
55347c6ae99SBarry Smith       if (dd->arrayghostedout[i] == *iptr) {
55447c6ae99SBarry Smith         iarray_start           = dd->startghostedout[i];
5550298fd71SBarry Smith         dd->arrayghostedout[i] = NULL;
5560298fd71SBarry Smith         dd->startghostedout[i] = NULL;
55747c6ae99SBarry Smith         break;
55847c6ae99SBarry Smith       }
55947c6ae99SBarry Smith     }
560aa219208SBarry Smith     for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) {
56147c6ae99SBarry Smith       if (!dd->arrayghostedin[i]) {
56247c6ae99SBarry Smith         dd->arrayghostedin[i] = *iptr;
56347c6ae99SBarry Smith         dd->startghostedin[i] = iarray_start;
56447c6ae99SBarry Smith         break;
56547c6ae99SBarry Smith       }
56647c6ae99SBarry Smith     }
56747c6ae99SBarry Smith   } else {
568aa219208SBarry Smith     for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) {
56947c6ae99SBarry Smith       if (dd->arrayout[i] == *iptr) {
57047c6ae99SBarry Smith         iarray_start    = dd->startout[i];
5710298fd71SBarry Smith         dd->arrayout[i] = NULL;
5720298fd71SBarry Smith         dd->startout[i] = NULL;
57347c6ae99SBarry Smith         break;
57447c6ae99SBarry Smith       }
57547c6ae99SBarry Smith     }
576aa219208SBarry Smith     for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) {
57747c6ae99SBarry Smith       if (!dd->arrayin[i]) {
57847c6ae99SBarry Smith         dd->arrayin[i] = *iptr;
57947c6ae99SBarry Smith         dd->startin[i] = iarray_start;
58047c6ae99SBarry Smith         break;
58147c6ae99SBarry Smith       }
58247c6ae99SBarry Smith     }
58347c6ae99SBarry Smith   }
58447c6ae99SBarry Smith   PetscFunctionReturn(0);
58547c6ae99SBarry Smith }
586