xref: /petsc/src/dm/impls/da/dalocal.c (revision a4e35b1925eceef64945ea472b84f2bf06a67b5e)
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));
423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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);
524f572ea9SToby Isaac   PetscAssertPointer(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
653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6647c6ae99SBarry Smith }
6747c6ae99SBarry Smith 
68939f6067SMatthew G. Knepley /*@
69dce8aebaSBarry 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:
72dce8aebaSBarry 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 
82dce8aebaSBarry 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) {
944f572ea9SToby Isaac     PetscAssertPointer(numCellsX, 2);
953582350cSMatthew G. Knepley     *numCellsX = mx;
963582350cSMatthew G. Knepley   }
973582350cSMatthew G. Knepley   if (numCellsY) {
984f572ea9SToby Isaac     PetscAssertPointer(numCellsY, 3);
993582350cSMatthew G. Knepley     *numCellsY = my;
1003582350cSMatthew G. Knepley   }
1013582350cSMatthew G. Knepley   if (numCellsZ) {
1024f572ea9SToby Isaac     PetscAssertPointer(numCellsZ, 4);
1033582350cSMatthew G. Knepley     *numCellsZ = mz;
1043582350cSMatthew G. Knepley   }
10557459e9aSMatthew G Knepley   if (numCells) {
1064f572ea9SToby Isaac     PetscAssertPointer(numCells, 5);
10757459e9aSMatthew G Knepley     *numCells = nC;
10857459e9aSMatthew G Knepley   }
1093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11057459e9aSMatthew G Knepley }
11157459e9aSMatthew G Knepley 
112d6401b93SMatthew G. Knepley /*@
113dce8aebaSBarry Smith   DMDAGetCellPoint - Get the DM point corresponding to the tuple (i, j, k) in the `DMDA`
114d6401b93SMatthew G. Knepley 
115d6401b93SMatthew G. Knepley   Input Parameters:
116dce8aebaSBarry Smith + dm - The `DMDA` object
117*a4e35b19SJacob Faibussowitsch . i  - The global x index for the cell
118*a4e35b19SJacob Faibussowitsch . j  - The global y index for the cell
119*a4e35b19SJacob Faibussowitsch - k  - The global z index for the cell
120d6401b93SMatthew G. Knepley 
1212fe279fdSBarry Smith   Output Parameter:
122dce8aebaSBarry Smith . point - The local `DM` point
123d6401b93SMatthew G. Knepley 
124d6401b93SMatthew G. Knepley   Level: developer
125d6401b93SMatthew G. Knepley 
126dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetNumCells()`
127d6401b93SMatthew G. Knepley @*/
128d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetCellPoint(DM dm, PetscInt i, PetscInt j, PetscInt k, PetscInt *point)
129d71ae5a4SJacob Faibussowitsch {
130c73cfb54SMatthew G. Knepley   const PetscInt dim = dm->dim;
131d6401b93SMatthew G. Knepley   DMDALocalInfo  info;
132d6401b93SMatthew G. Knepley 
133d6401b93SMatthew G. Knepley   PetscFunctionBegin;
134a9a02de4SBarry Smith   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMDA);
1354f572ea9SToby Isaac   PetscAssertPointer(point, 5);
1369566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(dm, &info));
13763a3b9bcSJacob 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);
13863a3b9bcSJacob 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);
13963a3b9bcSJacob 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);
140d6401b93SMatthew G. Knepley   *point = i + (dim > 1 ? (j + (dim > 2 ? k * info.gym : 0)) * info.gxm : 0);
1413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
142d6401b93SMatthew G. Knepley }
143d6401b93SMatthew G. Knepley 
144d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetNumVertices(DM dm, PetscInt *numVerticesX, PetscInt *numVerticesY, PetscInt *numVerticesZ, PetscInt *numVertices)
145d71ae5a4SJacob Faibussowitsch {
14657459e9aSMatthew G Knepley   DM_DA         *da  = (DM_DA *)dm->data;
147c73cfb54SMatthew G. Knepley   const PetscInt dim = dm->dim;
14857459e9aSMatthew G Knepley   const PetscInt mx = (da->Xe - da->Xs) / da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs;
14957459e9aSMatthew G Knepley   const PetscInt nVx = mx + 1;
15057459e9aSMatthew G Knepley   const PetscInt nVy = dim > 1 ? (my + 1) : 1;
15157459e9aSMatthew G Knepley   const PetscInt nVz = dim > 2 ? (mz + 1) : 1;
15257459e9aSMatthew G Knepley   const PetscInt nV  = nVx * nVy * nVz;
15357459e9aSMatthew G Knepley 
15457459e9aSMatthew G Knepley   PetscFunctionBegin;
15557459e9aSMatthew G Knepley   if (numVerticesX) {
1564f572ea9SToby Isaac     PetscAssertPointer(numVerticesX, 2);
15757459e9aSMatthew G Knepley     *numVerticesX = nVx;
15857459e9aSMatthew G Knepley   }
15957459e9aSMatthew G Knepley   if (numVerticesY) {
1604f572ea9SToby Isaac     PetscAssertPointer(numVerticesY, 3);
16157459e9aSMatthew G Knepley     *numVerticesY = nVy;
16257459e9aSMatthew G Knepley   }
16357459e9aSMatthew G Knepley   if (numVerticesZ) {
1644f572ea9SToby Isaac     PetscAssertPointer(numVerticesZ, 4);
16557459e9aSMatthew G Knepley     *numVerticesZ = nVz;
16657459e9aSMatthew G Knepley   }
16757459e9aSMatthew G Knepley   if (numVertices) {
1684f572ea9SToby Isaac     PetscAssertPointer(numVertices, 5);
16957459e9aSMatthew G Knepley     *numVertices = nV;
17057459e9aSMatthew G Knepley   }
1713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17257459e9aSMatthew G Knepley }
17357459e9aSMatthew G Knepley 
174d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetNumFaces(DM dm, PetscInt *numXFacesX, PetscInt *numXFaces, PetscInt *numYFacesY, PetscInt *numYFaces, PetscInt *numZFacesZ, PetscInt *numZFaces)
175d71ae5a4SJacob Faibussowitsch {
17657459e9aSMatthew G Knepley   DM_DA         *da  = (DM_DA *)dm->data;
177c73cfb54SMatthew G. Knepley   const PetscInt dim = dm->dim;
17857459e9aSMatthew G Knepley   const PetscInt mx = (da->Xe - da->Xs) / da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs;
17957459e9aSMatthew G Knepley   const PetscInt nxF = (dim > 1 ? (my) * (dim > 2 ? (mz) : 1) : 1);
18057459e9aSMatthew G Knepley   const PetscInt nXF = (mx + 1) * nxF;
18157459e9aSMatthew G Knepley   const PetscInt nyF = mx * (dim > 2 ? mz : 1);
18257459e9aSMatthew G Knepley   const PetscInt nYF = dim > 1 ? (my + 1) * nyF : 0;
18357459e9aSMatthew G Knepley   const PetscInt nzF = mx * (dim > 1 ? my : 0);
18457459e9aSMatthew G Knepley   const PetscInt nZF = dim > 2 ? (mz + 1) * nzF : 0;
18557459e9aSMatthew G Knepley 
18657459e9aSMatthew G Knepley   PetscFunctionBegin;
18757459e9aSMatthew G Knepley   if (numXFacesX) {
1884f572ea9SToby Isaac     PetscAssertPointer(numXFacesX, 2);
18957459e9aSMatthew G Knepley     *numXFacesX = nxF;
19057459e9aSMatthew G Knepley   }
19157459e9aSMatthew G Knepley   if (numXFaces) {
1924f572ea9SToby Isaac     PetscAssertPointer(numXFaces, 3);
19357459e9aSMatthew G Knepley     *numXFaces = nXF;
19457459e9aSMatthew G Knepley   }
19557459e9aSMatthew G Knepley   if (numYFacesY) {
1964f572ea9SToby Isaac     PetscAssertPointer(numYFacesY, 4);
19757459e9aSMatthew G Knepley     *numYFacesY = nyF;
19857459e9aSMatthew G Knepley   }
19957459e9aSMatthew G Knepley   if (numYFaces) {
2004f572ea9SToby Isaac     PetscAssertPointer(numYFaces, 5);
20157459e9aSMatthew G Knepley     *numYFaces = nYF;
20257459e9aSMatthew G Knepley   }
20357459e9aSMatthew G Knepley   if (numZFacesZ) {
2044f572ea9SToby Isaac     PetscAssertPointer(numZFacesZ, 6);
20557459e9aSMatthew G Knepley     *numZFacesZ = nzF;
20657459e9aSMatthew G Knepley   }
20757459e9aSMatthew G Knepley   if (numZFaces) {
2084f572ea9SToby Isaac     PetscAssertPointer(numZFaces, 7);
20957459e9aSMatthew G Knepley     *numZFaces = nZF;
21057459e9aSMatthew G Knepley   }
2113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
21257459e9aSMatthew G Knepley }
21357459e9aSMatthew G Knepley 
214d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetHeightStratum(DM dm, PetscInt height, PetscInt *pStart, PetscInt *pEnd)
215d71ae5a4SJacob Faibussowitsch {
216c73cfb54SMatthew G. Knepley   const PetscInt dim = dm->dim;
21757459e9aSMatthew G Knepley   PetscInt       nC, nV, nXF, nYF, nZF;
21857459e9aSMatthew G Knepley 
21957459e9aSMatthew G Knepley   PetscFunctionBegin;
2204f572ea9SToby Isaac   if (pStart) PetscAssertPointer(pStart, 3);
2214f572ea9SToby Isaac   if (pEnd) PetscAssertPointer(pEnd, 4);
2229566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumCells(dm, NULL, NULL, NULL, &nC));
2239566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV));
2249566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF));
22557459e9aSMatthew G Knepley   if (height == 0) {
22657459e9aSMatthew G Knepley     /* Cells */
2278865f1eaSKarl Rupp     if (pStart) *pStart = 0;
2288865f1eaSKarl Rupp     if (pEnd) *pEnd = nC;
22957459e9aSMatthew G Knepley   } else if (height == 1) {
23057459e9aSMatthew G Knepley     /* Faces */
2318865f1eaSKarl Rupp     if (pStart) *pStart = nC + nV;
2328865f1eaSKarl Rupp     if (pEnd) *pEnd = nC + nV + nXF + nYF + nZF;
23357459e9aSMatthew G Knepley   } else if (height == dim) {
23457459e9aSMatthew G Knepley     /* Vertices */
2358865f1eaSKarl Rupp     if (pStart) *pStart = nC;
2368865f1eaSKarl Rupp     if (pEnd) *pEnd = nC + nV;
23757459e9aSMatthew G Knepley   } else if (height < 0) {
23857459e9aSMatthew G Knepley     /* All points */
2398865f1eaSKarl Rupp     if (pStart) *pStart = 0;
2408865f1eaSKarl Rupp     if (pEnd) *pEnd = nC + nV + nXF + nYF + nZF;
24163a3b9bcSJacob Faibussowitsch   } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No points of height %" PetscInt_FMT " in the DA", height);
2423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
24357459e9aSMatthew G Knepley }
24457459e9aSMatthew G Knepley 
245d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetDepthStratum(DM dm, PetscInt depth, PetscInt *pStart, PetscInt *pEnd)
246d71ae5a4SJacob Faibussowitsch {
247c73cfb54SMatthew G. Knepley   const PetscInt dim = dm->dim;
2483385ff7eSMatthew G. Knepley   PetscInt       nC, nV, nXF, nYF, nZF;
2493385ff7eSMatthew G. Knepley 
2503385ff7eSMatthew G. Knepley   PetscFunctionBegin;
2514f572ea9SToby Isaac   if (pStart) PetscAssertPointer(pStart, 3);
2524f572ea9SToby Isaac   if (pEnd) PetscAssertPointer(pEnd, 4);
2539566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumCells(dm, NULL, NULL, NULL, &nC));
2549566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV));
2559566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF));
2563385ff7eSMatthew G. Knepley   if (depth == dim) {
2573385ff7eSMatthew G. Knepley     /* Cells */
2583385ff7eSMatthew G. Knepley     if (pStart) *pStart = 0;
2593385ff7eSMatthew G. Knepley     if (pEnd) *pEnd = nC;
2603385ff7eSMatthew G. Knepley   } else if (depth == dim - 1) {
2613385ff7eSMatthew G. Knepley     /* Faces */
2623385ff7eSMatthew G. Knepley     if (pStart) *pStart = nC + nV;
2633385ff7eSMatthew G. Knepley     if (pEnd) *pEnd = nC + nV + nXF + nYF + nZF;
2643385ff7eSMatthew G. Knepley   } else if (depth == 0) {
2653385ff7eSMatthew G. Knepley     /* Vertices */
2663385ff7eSMatthew G. Knepley     if (pStart) *pStart = nC;
2673385ff7eSMatthew G. Knepley     if (pEnd) *pEnd = nC + nV;
2683385ff7eSMatthew G. Knepley   } else if (depth < 0) {
2693385ff7eSMatthew G. Knepley     /* All points */
2703385ff7eSMatthew G. Knepley     if (pStart) *pStart = 0;
2713385ff7eSMatthew G. Knepley     if (pEnd) *pEnd = nC + nV + nXF + nYF + nZF;
27263a3b9bcSJacob Faibussowitsch   } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No points of depth %" PetscInt_FMT " in the DA", depth);
2733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2743385ff7eSMatthew G. Knepley }
2753385ff7eSMatthew G. Knepley 
276d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetConeSize(DM dm, PetscInt p, PetscInt *coneSize)
277d71ae5a4SJacob Faibussowitsch {
278c73cfb54SMatthew G. Knepley   const PetscInt dim = dm->dim;
2793385ff7eSMatthew G. Knepley   PetscInt       nC, nV, nXF, nYF, nZF;
2803385ff7eSMatthew G. Knepley 
2813385ff7eSMatthew G. Knepley   PetscFunctionBegin;
2823385ff7eSMatthew G. Knepley   *coneSize = 0;
2839566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumCells(dm, NULL, NULL, NULL, &nC));
2849566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV));
2859566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF));
2863385ff7eSMatthew G. Knepley   switch (dim) {
2873385ff7eSMatthew G. Knepley   case 2:
2883385ff7eSMatthew G. Knepley     if (p >= 0) {
2893385ff7eSMatthew G. Knepley       if (p < nC) {
2903385ff7eSMatthew G. Knepley         *coneSize = 4;
2913385ff7eSMatthew G. Knepley       } else if (p < nC + nV) {
2923385ff7eSMatthew G. Knepley         *coneSize = 0;
2933385ff7eSMatthew G. Knepley       } else if (p < nC + nV + nXF + nYF + nZF) {
2943385ff7eSMatthew G. Knepley         *coneSize = 2;
29563a3b9bcSJacob 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);
29663a3b9bcSJacob Faibussowitsch     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative point %" PetscInt_FMT " is invalid", p);
2973385ff7eSMatthew G. Knepley     break;
298d71ae5a4SJacob Faibussowitsch   case 3:
299d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do 3D");
3003385ff7eSMatthew G. Knepley   }
3013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3023385ff7eSMatthew G. Knepley }
3033385ff7eSMatthew G. Knepley 
304d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetCone(DM dm, PetscInt p, PetscInt *cone[])
305d71ae5a4SJacob Faibussowitsch {
306c73cfb54SMatthew G. Knepley   const PetscInt dim = dm->dim;
3073385ff7eSMatthew G. Knepley   PetscInt       nCx, nCy, nCz, nC, nVx, nVy, nVz, nV, nxF, nyF, nzF, nXF, nYF, nZF;
3083385ff7eSMatthew G. Knepley 
3093385ff7eSMatthew G. Knepley   PetscFunctionBegin;
3109566063dSJacob Faibussowitsch   if (!cone) PetscCall(DMGetWorkArray(dm, 6, MPIU_INT, cone));
3119566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumCells(dm, &nCx, &nCy, &nCz, &nC));
3129566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV));
3139566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF));
3143385ff7eSMatthew G. Knepley   switch (dim) {
3153385ff7eSMatthew G. Knepley   case 2:
3163385ff7eSMatthew G. Knepley     if (p >= 0) {
3173385ff7eSMatthew G. Knepley       if (p < nC) {
3183385ff7eSMatthew G. Knepley         const PetscInt cy = p / nCx;
3193385ff7eSMatthew G. Knepley         const PetscInt cx = p % nCx;
3203385ff7eSMatthew G. Knepley 
3213385ff7eSMatthew G. Knepley         (*cone)[0] = cy * nxF + cx + nC + nV;
3223385ff7eSMatthew G. Knepley         (*cone)[1] = cx * nyF + cy + nyF + nC + nV + nXF;
3233385ff7eSMatthew G. Knepley         (*cone)[2] = cy * nxF + cx + nxF + nC + nV;
3243385ff7eSMatthew G. Knepley         (*cone)[3] = cx * nyF + cy + nC + nV + nXF;
3253385ff7eSMatthew G. Knepley         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do cell cones");
3263385ff7eSMatthew G. Knepley       } else if (p < nC + nV) {
3273385ff7eSMatthew G. Knepley       } else if (p < nC + nV + nXF) {
3283385ff7eSMatthew G. Knepley         const PetscInt fy = (p - nC + nV) / nxF;
3293385ff7eSMatthew G. Knepley         const PetscInt fx = (p - nC + nV) % nxF;
3303385ff7eSMatthew G. Knepley 
3313385ff7eSMatthew G. Knepley         (*cone)[0] = fy * nVx + fx + nC;
3323385ff7eSMatthew G. Knepley         (*cone)[1] = fy * nVx + fx + 1 + nC;
3333385ff7eSMatthew G. Knepley       } else if (p < nC + nV + nXF + nYF) {
3343385ff7eSMatthew G. Knepley         const PetscInt fx = (p - nC + nV + nXF) / nyF;
3353385ff7eSMatthew G. Knepley         const PetscInt fy = (p - nC + nV + nXF) % nyF;
3363385ff7eSMatthew G. Knepley 
3373385ff7eSMatthew G. Knepley         (*cone)[0] = fy * nVx + fx + nC;
3383385ff7eSMatthew G. Knepley         (*cone)[1] = fy * nVx + fx + nVx + nC;
33963a3b9bcSJacob 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);
34063a3b9bcSJacob Faibussowitsch     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative point %" PetscInt_FMT " is invalid", p);
3413385ff7eSMatthew G. Knepley     break;
342d71ae5a4SJacob Faibussowitsch   case 3:
343d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do 3D");
3443385ff7eSMatthew G. Knepley   }
3453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3463385ff7eSMatthew G. Knepley }
3473385ff7eSMatthew G. Knepley 
348d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDARestoreCone(DM dm, PetscInt p, PetscInt *cone[])
349d71ae5a4SJacob Faibussowitsch {
3503385ff7eSMatthew G. Knepley   PetscFunctionBegin;
3519566063dSJacob Faibussowitsch   PetscCall(DMGetWorkArray(dm, 6, MPIU_INT, cone));
3523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3533385ff7eSMatthew G. Knepley }
3543385ff7eSMatthew G. Knepley 
355d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDASetVertexCoordinates(DM dm, PetscReal xl, PetscReal xu, PetscReal yl, PetscReal yu, PetscReal zl, PetscReal zu)
356d71ae5a4SJacob Faibussowitsch {
3573385ff7eSMatthew G. Knepley   DM_DA       *da = (DM_DA *)dm->data;
3583385ff7eSMatthew G. Knepley   Vec          coordinates;
3593385ff7eSMatthew G. Knepley   PetscSection section;
3603385ff7eSMatthew G. Knepley   PetscScalar *coords;
3613385ff7eSMatthew G. Knepley   PetscReal    h[3];
3623385ff7eSMatthew G. Knepley   PetscInt     dim, size, M, N, P, nVx, nVy, nVz, nV, vStart, vEnd, v, i, j, k;
3633385ff7eSMatthew G. Knepley 
3643385ff7eSMatthew G. Knepley   PetscFunctionBegin;
365a9a02de4SBarry Smith   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMDA);
3669566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dm, &dim, &M, &N, &P, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
36708401ef6SPierre Jolivet   PetscCheck(dim <= 3, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "The following code only works for dim <= 3");
3683385ff7eSMatthew G. Knepley   h[0] = (xu - xl) / M;
3693385ff7eSMatthew G. Knepley   h[1] = (yu - yl) / N;
3703385ff7eSMatthew G. Knepley   h[2] = (zu - zl) / P;
3719566063dSJacob Faibussowitsch   PetscCall(DMDAGetDepthStratum(dm, 0, &vStart, &vEnd));
3729566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV));
3739566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), &section));
3749566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetNumFields(section, 1));
3759566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetFieldComponents(section, 0, dim));
3769566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(section, vStart, vEnd));
37748a46eb9SPierre Jolivet   for (v = vStart; v < vEnd; ++v) PetscCall(PetscSectionSetDof(section, v, dim));
3789566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(section));
3799566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(section, &size));
3809566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, size, &coordinates));
3819566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates"));
3829566063dSJacob Faibussowitsch   PetscCall(VecGetArray(coordinates, &coords));
3833385ff7eSMatthew G. Knepley   for (k = 0; k < nVz; ++k) {
384e4d69e08SBarry Smith     PetscInt ind[3], d, off;
3853385ff7eSMatthew G. Knepley 
386e4d69e08SBarry Smith     ind[0] = 0;
387e4d69e08SBarry Smith     ind[1] = 0;
388e4d69e08SBarry Smith     ind[2] = k + da->zs;
3893385ff7eSMatthew G. Knepley     for (j = 0; j < nVy; ++j) {
3903385ff7eSMatthew G. Knepley       ind[1] = j + da->ys;
3913385ff7eSMatthew G. Knepley       for (i = 0; i < nVx; ++i) {
3923385ff7eSMatthew G. Knepley         const PetscInt vertex = (k * nVy + j) * nVx + i + vStart;
3933385ff7eSMatthew G. Knepley 
3949566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetOffset(section, vertex, &off));
3953385ff7eSMatthew G. Knepley         ind[0] = i + da->xs;
396ad540459SPierre Jolivet         for (d = 0; d < dim; ++d) coords[off + d] = h[d] * ind[d];
3973385ff7eSMatthew G. Knepley       }
3983385ff7eSMatthew G. Knepley     }
3993385ff7eSMatthew G. Knepley   }
4009566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(coordinates, &coords));
4019566063dSJacob Faibussowitsch   PetscCall(DMSetCoordinateSection(dm, PETSC_DETERMINE, section));
4029566063dSJacob Faibussowitsch   PetscCall(DMSetCoordinatesLocal(dm, coordinates));
4039566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&section));
4049566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&coordinates));
4053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4063385ff7eSMatthew G. Knepley }
4079a800dd8SMatthew G. Knepley 
40847c6ae99SBarry Smith /* ------------------------------------------------------------------- */
40947c6ae99SBarry Smith 
41047c6ae99SBarry Smith /*@C
411dce8aebaSBarry Smith   DMDAGetArray - Gets a work array for a `DMDA`
41247c6ae99SBarry Smith 
413d8d19677SJose E. Roman   Input Parameters:
41447c6ae99SBarry Smith + da      - information about my local patch
41547c6ae99SBarry Smith - ghosted - do you want arrays for the ghosted or nonghosted patch
41647c6ae99SBarry Smith 
4172fe279fdSBarry Smith   Output Parameter:
41847c6ae99SBarry Smith . vptr - array data structured
41947c6ae99SBarry Smith 
42047c6ae99SBarry Smith   Level: advanced
42147c6ae99SBarry Smith 
422dce8aebaSBarry Smith   Note:
423dce8aebaSBarry Smith   The vector values are NOT initialized and may have garbage in them, so you may need
424dce8aebaSBarry Smith   to zero them.
42547c6ae99SBarry Smith 
426dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDARestoreArray()`
42747c6ae99SBarry Smith @*/
428d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAGetArray(DM da, PetscBool ghosted, void *vptr)
429d71ae5a4SJacob Faibussowitsch {
43047c6ae99SBarry Smith   PetscInt j, i, xs, ys, xm, ym, zs, zm;
43147c6ae99SBarry Smith   char    *iarray_start;
43247c6ae99SBarry Smith   void   **iptr = (void **)vptr;
43347c6ae99SBarry Smith   DM_DA   *dd   = (DM_DA *)da->data;
43447c6ae99SBarry Smith 
43547c6ae99SBarry Smith   PetscFunctionBegin;
436a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
43747c6ae99SBarry Smith   if (ghosted) {
438aa219208SBarry Smith     for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) {
43947c6ae99SBarry Smith       if (dd->arrayghostedin[i]) {
44047c6ae99SBarry Smith         *iptr                 = dd->arrayghostedin[i];
44147c6ae99SBarry Smith         iarray_start          = (char *)dd->startghostedin[i];
4420298fd71SBarry Smith         dd->arrayghostedin[i] = NULL;
4430298fd71SBarry Smith         dd->startghostedin[i] = NULL;
44447c6ae99SBarry Smith 
44547c6ae99SBarry Smith         goto done;
44647c6ae99SBarry Smith       }
44747c6ae99SBarry Smith     }
44847c6ae99SBarry Smith     xs = dd->Xs;
44947c6ae99SBarry Smith     ys = dd->Ys;
45047c6ae99SBarry Smith     zs = dd->Zs;
45147c6ae99SBarry Smith     xm = dd->Xe - dd->Xs;
45247c6ae99SBarry Smith     ym = dd->Ye - dd->Ys;
45347c6ae99SBarry Smith     zm = dd->Ze - dd->Zs;
45447c6ae99SBarry Smith   } else {
455aa219208SBarry Smith     for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) {
45647c6ae99SBarry Smith       if (dd->arrayin[i]) {
45747c6ae99SBarry Smith         *iptr          = dd->arrayin[i];
45847c6ae99SBarry Smith         iarray_start   = (char *)dd->startin[i];
4590298fd71SBarry Smith         dd->arrayin[i] = NULL;
4600298fd71SBarry Smith         dd->startin[i] = NULL;
46147c6ae99SBarry Smith 
46247c6ae99SBarry Smith         goto done;
46347c6ae99SBarry Smith       }
46447c6ae99SBarry Smith     }
46547c6ae99SBarry Smith     xs = dd->xs;
46647c6ae99SBarry Smith     ys = dd->ys;
46747c6ae99SBarry Smith     zs = dd->zs;
46847c6ae99SBarry Smith     xm = dd->xe - dd->xs;
46947c6ae99SBarry Smith     ym = dd->ye - dd->ys;
47047c6ae99SBarry Smith     zm = dd->ze - dd->zs;
47147c6ae99SBarry Smith   }
47247c6ae99SBarry Smith 
473c73cfb54SMatthew G. Knepley   switch (da->dim) {
47447c6ae99SBarry Smith   case 1: {
47547c6ae99SBarry Smith     void *ptr;
47647c6ae99SBarry Smith 
4779566063dSJacob Faibussowitsch     PetscCall(PetscMalloc(xm * sizeof(PetscScalar), &iarray_start));
47847c6ae99SBarry Smith 
47947c6ae99SBarry Smith     ptr   = (void *)(iarray_start - xs * sizeof(PetscScalar));
48047c6ae99SBarry Smith     *iptr = (void *)ptr;
4818865f1eaSKarl Rupp     break;
4828865f1eaSKarl Rupp   }
48347c6ae99SBarry Smith   case 2: {
48447c6ae99SBarry Smith     void **ptr;
48547c6ae99SBarry Smith 
4869566063dSJacob Faibussowitsch     PetscCall(PetscMalloc((ym + 1) * sizeof(void *) + xm * ym * sizeof(PetscScalar), &iarray_start));
48747c6ae99SBarry Smith 
48847c6ae99SBarry Smith     ptr = (void **)(iarray_start + xm * ym * sizeof(PetscScalar) - ys * sizeof(void *));
4898865f1eaSKarl Rupp     for (j = ys; j < ys + ym; j++) ptr[j] = iarray_start + sizeof(PetscScalar) * (xm * (j - ys) - xs);
49047c6ae99SBarry Smith     *iptr = (void *)ptr;
4918865f1eaSKarl Rupp     break;
4928865f1eaSKarl Rupp   }
49347c6ae99SBarry Smith   case 3: {
49447c6ae99SBarry Smith     void ***ptr, **bptr;
49547c6ae99SBarry Smith 
4969566063dSJacob Faibussowitsch     PetscCall(PetscMalloc((zm + 1) * sizeof(void **) + (ym * zm + 1) * sizeof(void *) + xm * ym * zm * sizeof(PetscScalar), &iarray_start));
49747c6ae99SBarry Smith 
49847c6ae99SBarry Smith     ptr  = (void ***)(iarray_start + xm * ym * zm * sizeof(PetscScalar) - zs * sizeof(void *));
49947c6ae99SBarry Smith     bptr = (void **)(iarray_start + xm * ym * zm * sizeof(PetscScalar) + zm * sizeof(void **));
5008865f1eaSKarl Rupp     for (i = zs; i < zs + zm; i++) ptr[i] = bptr + ((i - zs) * ym - ys);
50147c6ae99SBarry Smith     for (i = zs; i < zs + zm; i++) {
502ad540459SPierre Jolivet       for (j = ys; j < ys + ym; j++) ptr[i][j] = iarray_start + sizeof(PetscScalar) * (xm * ym * (i - zs) + xm * (j - ys) - xs);
50347c6ae99SBarry Smith     }
50447c6ae99SBarry Smith     *iptr = (void *)ptr;
5058865f1eaSKarl Rupp     break;
5068865f1eaSKarl Rupp   }
507d71ae5a4SJacob Faibussowitsch   default:
508d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Dimension %" PetscInt_FMT " not supported", da->dim);
50947c6ae99SBarry Smith   }
51047c6ae99SBarry Smith 
51147c6ae99SBarry Smith done:
51247c6ae99SBarry Smith   /* add arrays to the checked out list */
51347c6ae99SBarry Smith   if (ghosted) {
514aa219208SBarry Smith     for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) {
51547c6ae99SBarry Smith       if (!dd->arrayghostedout[i]) {
51647c6ae99SBarry Smith         dd->arrayghostedout[i] = *iptr;
51747c6ae99SBarry Smith         dd->startghostedout[i] = iarray_start;
51847c6ae99SBarry Smith         break;
51947c6ae99SBarry Smith       }
52047c6ae99SBarry Smith     }
52147c6ae99SBarry Smith   } else {
522aa219208SBarry Smith     for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) {
52347c6ae99SBarry Smith       if (!dd->arrayout[i]) {
52447c6ae99SBarry Smith         dd->arrayout[i] = *iptr;
52547c6ae99SBarry Smith         dd->startout[i] = iarray_start;
52647c6ae99SBarry Smith         break;
52747c6ae99SBarry Smith       }
52847c6ae99SBarry Smith     }
52947c6ae99SBarry Smith   }
5303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
53147c6ae99SBarry Smith }
53247c6ae99SBarry Smith 
53347c6ae99SBarry Smith /*@C
534dce8aebaSBarry Smith   DMDARestoreArray - Restores an array of derivative types for a `DMDA`
53547c6ae99SBarry Smith 
536d8d19677SJose E. Roman   Input Parameters:
53747c6ae99SBarry Smith + da      - information about my local patch
53847c6ae99SBarry Smith . ghosted - do you want arrays for the ghosted or nonghosted patch
539cb004a26SBarry Smith - vptr    - array data structured
54047c6ae99SBarry Smith 
54147c6ae99SBarry Smith   Level: advanced
54247c6ae99SBarry Smith 
543dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetArray()`
54447c6ae99SBarry Smith @*/
545d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDARestoreArray(DM da, PetscBool ghosted, void *vptr)
546d71ae5a4SJacob Faibussowitsch {
54747c6ae99SBarry Smith   PetscInt i;
548ea78f98cSLisandro Dalcin   void   **iptr = (void **)vptr, *iarray_start = NULL;
54947c6ae99SBarry Smith   DM_DA   *dd = (DM_DA *)da->data;
55047c6ae99SBarry Smith 
55147c6ae99SBarry Smith   PetscFunctionBegin;
552a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
55347c6ae99SBarry Smith   if (ghosted) {
554aa219208SBarry Smith     for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) {
55547c6ae99SBarry Smith       if (dd->arrayghostedout[i] == *iptr) {
55647c6ae99SBarry Smith         iarray_start           = dd->startghostedout[i];
5570298fd71SBarry Smith         dd->arrayghostedout[i] = NULL;
5580298fd71SBarry Smith         dd->startghostedout[i] = NULL;
55947c6ae99SBarry Smith         break;
56047c6ae99SBarry Smith       }
56147c6ae99SBarry Smith     }
562aa219208SBarry Smith     for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) {
56347c6ae99SBarry Smith       if (!dd->arrayghostedin[i]) {
56447c6ae99SBarry Smith         dd->arrayghostedin[i] = *iptr;
56547c6ae99SBarry Smith         dd->startghostedin[i] = iarray_start;
56647c6ae99SBarry Smith         break;
56747c6ae99SBarry Smith       }
56847c6ae99SBarry Smith     }
56947c6ae99SBarry Smith   } else {
570aa219208SBarry Smith     for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) {
57147c6ae99SBarry Smith       if (dd->arrayout[i] == *iptr) {
57247c6ae99SBarry Smith         iarray_start    = dd->startout[i];
5730298fd71SBarry Smith         dd->arrayout[i] = NULL;
5740298fd71SBarry Smith         dd->startout[i] = NULL;
57547c6ae99SBarry Smith         break;
57647c6ae99SBarry Smith       }
57747c6ae99SBarry Smith     }
578aa219208SBarry Smith     for (i = 0; i < DMDA_MAX_WORK_ARRAYS; i++) {
57947c6ae99SBarry Smith       if (!dd->arrayin[i]) {
58047c6ae99SBarry Smith         dd->arrayin[i] = *iptr;
58147c6ae99SBarry Smith         dd->startin[i] = iarray_start;
58247c6ae99SBarry Smith         break;
58347c6ae99SBarry Smith       }
58447c6ae99SBarry Smith     }
58547c6ae99SBarry Smith   }
5863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
58747c6ae99SBarry Smith }
588