xref: /petsc/src/dm/impls/da/dalocal.c (revision bcaeba4d41d6ca6c6dc4189db20683073a9959ce)
147c6ae99SBarry Smith 
247c6ae99SBarry Smith /*
347c6ae99SBarry Smith   Code for manipulating distributed regular arrays in parallel.
447c6ae99SBarry Smith */
547c6ae99SBarry Smith 
6b45d2f2cSJed Brown #include <petsc-private/daimpl.h>    /*I   "petscdmda.h"   I*/
747c6ae99SBarry Smith 
847c6ae99SBarry Smith /*
9e3c5b3baSBarry Smith    This allows the DMDA vectors to properly tell MATLAB their dimensions
1047c6ae99SBarry Smith */
1147c6ae99SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
12c6db04a5SJed Brown #include <engine.h>   /* MATLAB include file */
13c6db04a5SJed Brown #include <mex.h>      /* MATLAB include file */
1447c6ae99SBarry Smith EXTERN_C_BEGIN
1547c6ae99SBarry Smith #undef __FUNCT__
1647c6ae99SBarry Smith #define __FUNCT__ "VecMatlabEnginePut_DA2d"
177087cfbeSBarry Smith PetscErrorCode  VecMatlabEnginePut_DA2d(PetscObject obj,void *mengine)
1847c6ae99SBarry Smith {
1947c6ae99SBarry Smith   PetscErrorCode ierr;
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;
27c688c046SMatthew G Knepley   ierr = VecGetDM(vec, &da);CHKERRQ(ierr);
28aa219208SBarry Smith   if (!da) SETERRQ(((PetscObject)vec)->comm,PETSC_ERR_ARG_WRONGSTATE,"Vector not associated with a DMDA");
29aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,0,0,0,&m,&n,0);CHKERRQ(ierr);
3047c6ae99SBarry Smith 
3147c6ae99SBarry Smith   ierr = VecGetArray(vec,&array);CHKERRQ(ierr);
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
3747c6ae99SBarry Smith   ierr = PetscMemcpy(mxGetPr(mat),array,n*m*sizeof(PetscScalar));CHKERRQ(ierr);
3847c6ae99SBarry Smith   ierr = PetscObjectName(obj);CHKERRQ(ierr);
3947c6ae99SBarry Smith   engPutVariable((Engine *)mengine,obj->name,mat);
4047c6ae99SBarry Smith 
4147c6ae99SBarry Smith   ierr = VecRestoreArray(vec,&array);CHKERRQ(ierr);
4247c6ae99SBarry Smith   PetscFunctionReturn(0);
4347c6ae99SBarry Smith }
4447c6ae99SBarry Smith EXTERN_C_END
4547c6ae99SBarry Smith #endif
4647c6ae99SBarry Smith 
4747c6ae99SBarry Smith 
4847c6ae99SBarry Smith #undef __FUNCT__
49564755cdSBarry Smith #define __FUNCT__ "DMCreateLocalVector_DA"
507087cfbeSBarry Smith PetscErrorCode  DMCreateLocalVector_DA(DM da,Vec* g)
5147c6ae99SBarry Smith {
5247c6ae99SBarry Smith   PetscErrorCode ierr;
5347c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
5447c6ae99SBarry Smith 
5547c6ae99SBarry Smith   PetscFunctionBegin;
5647c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
5747c6ae99SBarry Smith   PetscValidPointer(g,2);
5847c6ae99SBarry Smith   ierr = VecCreate(PETSC_COMM_SELF,g);CHKERRQ(ierr);
5947c6ae99SBarry Smith   ierr = VecSetSizes(*g,dd->nlocal,PETSC_DETERMINE);CHKERRQ(ierr);
6047c6ae99SBarry Smith   ierr = VecSetBlockSize(*g,dd->w);CHKERRQ(ierr);
61401ddaa8SBarry Smith   ierr = VecSetType(*g,da->vectype);CHKERRQ(ierr);
62c688c046SMatthew G Knepley   ierr = VecSetDM(*g, da);CHKERRQ(ierr);
6347c6ae99SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
6447c6ae99SBarry Smith   if (dd->w == 1  && dd->dim == 2) {
6547c6ae99SBarry Smith     ierr = PetscObjectComposeFunctionDynamic((PetscObject)*g,"PetscMatlabEnginePut_C","VecMatlabEnginePut_DA2d",VecMatlabEnginePut_DA2d);CHKERRQ(ierr);
6647c6ae99SBarry Smith   }
6747c6ae99SBarry Smith #endif
6847c6ae99SBarry Smith   PetscFunctionReturn(0);
6947c6ae99SBarry Smith }
7047c6ae99SBarry Smith 
71a66d4d66SMatthew G Knepley #undef __FUNCT__
7257459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetNumCells"
7357459e9aSMatthew G Knepley PetscErrorCode DMDAGetNumCells(DM dm, PetscInt *numCells)
7457459e9aSMatthew G Knepley {
7557459e9aSMatthew G Knepley   DM_DA         *da  = (DM_DA *) dm->data;
7657459e9aSMatthew G Knepley   const PetscInt dim = da->dim;
7757459e9aSMatthew G Knepley   const PetscInt mx  = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs;
7857459e9aSMatthew G Knepley   const PetscInt nC  = (mx  )*(dim > 1 ? (my  )*(dim > 2 ? (mz  ) : 1) : 1);
7957459e9aSMatthew G Knepley 
8057459e9aSMatthew G Knepley   PetscFunctionBegin;
8157459e9aSMatthew G Knepley   if (numCells) {
8257459e9aSMatthew G Knepley     PetscValidIntPointer(numCells,2);
8357459e9aSMatthew G Knepley     *numCells = nC;
8457459e9aSMatthew G Knepley   }
8557459e9aSMatthew G Knepley   PetscFunctionReturn(0);
8657459e9aSMatthew G Knepley }
8757459e9aSMatthew G Knepley 
8857459e9aSMatthew G Knepley #undef __FUNCT__
8957459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetNumVertices"
9057459e9aSMatthew G Knepley PetscErrorCode DMDAGetNumVertices(DM dm, PetscInt *numVerticesX, PetscInt *numVerticesY, PetscInt *numVerticesZ, PetscInt *numVertices)
9157459e9aSMatthew G Knepley {
9257459e9aSMatthew G Knepley   DM_DA         *da  = (DM_DA *) dm->data;
9357459e9aSMatthew G Knepley   const PetscInt dim = da->dim;
9457459e9aSMatthew G Knepley   const PetscInt mx  = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs;
9557459e9aSMatthew G Knepley   const PetscInt nVx = mx+1;
9657459e9aSMatthew G Knepley   const PetscInt nVy = dim > 1 ? (my+1) : 1;
9757459e9aSMatthew G Knepley   const PetscInt nVz = dim > 2 ? (mz+1) : 1;
9857459e9aSMatthew G Knepley   const PetscInt nV  = nVx*nVy*nVz;
9957459e9aSMatthew G Knepley 
10057459e9aSMatthew G Knepley   PetscFunctionBegin;
10157459e9aSMatthew G Knepley   if (numVerticesX) {
10257459e9aSMatthew G Knepley     PetscValidIntPointer(numVerticesX,2);
10357459e9aSMatthew G Knepley     *numVerticesX = nVx;
10457459e9aSMatthew G Knepley   }
10557459e9aSMatthew G Knepley   if (numVerticesY) {
10657459e9aSMatthew G Knepley     PetscValidIntPointer(numVerticesY,3);
10757459e9aSMatthew G Knepley     *numVerticesY = nVy;
10857459e9aSMatthew G Knepley   }
10957459e9aSMatthew G Knepley   if (numVerticesZ) {
11057459e9aSMatthew G Knepley     PetscValidIntPointer(numVerticesZ,4);
11157459e9aSMatthew G Knepley     *numVerticesZ = nVz;
11257459e9aSMatthew G Knepley   }
11357459e9aSMatthew G Knepley   if (numVertices) {
11457459e9aSMatthew G Knepley     PetscValidIntPointer(numVertices,5);
11557459e9aSMatthew G Knepley     *numVertices = nV;
11657459e9aSMatthew G Knepley   }
11757459e9aSMatthew G Knepley   PetscFunctionReturn(0);
11857459e9aSMatthew G Knepley }
11957459e9aSMatthew G Knepley 
12057459e9aSMatthew G Knepley #undef __FUNCT__
12157459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetNumFaces"
12257459e9aSMatthew G Knepley PetscErrorCode DMDAGetNumFaces(DM dm, PetscInt *numXFacesX, PetscInt *numXFaces, PetscInt *numYFacesY, PetscInt *numYFaces, PetscInt *numZFacesZ, PetscInt *numZFaces)
12357459e9aSMatthew G Knepley {
12457459e9aSMatthew G Knepley   DM_DA         *da  = (DM_DA *) dm->data;
12557459e9aSMatthew G Knepley   const PetscInt dim = da->dim;
12657459e9aSMatthew G Knepley   const PetscInt mx  = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs;
12757459e9aSMatthew G Knepley   const PetscInt nxF = (dim > 1 ? (my  )*(dim > 2 ? (mz  ) : 1) : 1);
12857459e9aSMatthew G Knepley   const PetscInt nXF = (mx+1)*nxF;
12957459e9aSMatthew G Knepley   const PetscInt nyF = mx*(dim > 2 ? mz : 1);
13057459e9aSMatthew G Knepley   const PetscInt nYF = dim > 1 ? (my+1)*nyF : 0;
13157459e9aSMatthew G Knepley   const PetscInt nzF = mx*(dim > 1 ? my : 0);
13257459e9aSMatthew G Knepley   const PetscInt nZF = dim > 2 ? (mz+1)*nzF : 0;
13357459e9aSMatthew G Knepley 
13457459e9aSMatthew G Knepley   PetscFunctionBegin;
13557459e9aSMatthew G Knepley   if (numXFacesX) {
13657459e9aSMatthew G Knepley     PetscValidIntPointer(numXFacesX,2);
13757459e9aSMatthew G Knepley     *numXFacesX = nxF;
13857459e9aSMatthew G Knepley   }
13957459e9aSMatthew G Knepley   if (numXFaces) {
14057459e9aSMatthew G Knepley     PetscValidIntPointer(numXFaces,3);
14157459e9aSMatthew G Knepley     *numXFaces = nXF;
14257459e9aSMatthew G Knepley   }
14357459e9aSMatthew G Knepley   if (numYFacesY) {
14457459e9aSMatthew G Knepley     PetscValidIntPointer(numYFacesY,4);
14557459e9aSMatthew G Knepley     *numYFacesY = nyF;
14657459e9aSMatthew G Knepley   }
14757459e9aSMatthew G Knepley   if (numYFaces) {
14857459e9aSMatthew G Knepley     PetscValidIntPointer(numYFaces,5);
14957459e9aSMatthew G Knepley     *numYFaces = nYF;
15057459e9aSMatthew G Knepley   }
15157459e9aSMatthew G Knepley   if (numZFacesZ) {
15257459e9aSMatthew G Knepley     PetscValidIntPointer(numZFacesZ,6);
15357459e9aSMatthew G Knepley     *numZFacesZ = nzF;
15457459e9aSMatthew G Knepley   }
15557459e9aSMatthew G Knepley   if (numZFaces) {
15657459e9aSMatthew G Knepley     PetscValidIntPointer(numZFaces,7);
15757459e9aSMatthew G Knepley     *numZFaces = nZF;
15857459e9aSMatthew G Knepley   }
15957459e9aSMatthew G Knepley   PetscFunctionReturn(0);
16057459e9aSMatthew G Knepley }
16157459e9aSMatthew G Knepley 
16257459e9aSMatthew G Knepley #undef __FUNCT__
16357459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetHeightStratum"
16457459e9aSMatthew G Knepley PetscErrorCode DMDAGetHeightStratum(DM dm, PetscInt height, PetscInt *pStart, PetscInt *pEnd)
16557459e9aSMatthew G Knepley {
16657459e9aSMatthew G Knepley   DM_DA         *da  = (DM_DA *) dm->data;
16757459e9aSMatthew G Knepley   const PetscInt dim = da->dim;
16857459e9aSMatthew G Knepley   PetscInt       nC, nV, nXF, nYF, nZF;
16957459e9aSMatthew G Knepley   PetscErrorCode ierr;
17057459e9aSMatthew G Knepley 
17157459e9aSMatthew G Knepley   PetscFunctionBegin;
17257459e9aSMatthew G Knepley   if (pStart) {PetscValidIntPointer(pStart,3);}
17357459e9aSMatthew G Knepley   if (pEnd)   {PetscValidIntPointer(pEnd,4);}
17457459e9aSMatthew G Knepley   ierr = DMDAGetNumCells(dm, &nC);CHKERRQ(ierr);
17557459e9aSMatthew G Knepley   ierr = DMDAGetNumVertices(dm, PETSC_NULL, PETSC_NULL, PETSC_NULL, &nV);CHKERRQ(ierr);
17657459e9aSMatthew G Knepley   ierr = DMDAGetNumFaces(dm, PETSC_NULL, &nXF, PETSC_NULL, &nYF, PETSC_NULL, &nZF);CHKERRQ(ierr);
17757459e9aSMatthew G Knepley   if (height == 0) {
17857459e9aSMatthew G Knepley     /* Cells */
17957459e9aSMatthew G Knepley     if (pStart) {*pStart = 0;}
18057459e9aSMatthew G Knepley     if (pEnd)   {*pEnd   = nC;}
18157459e9aSMatthew G Knepley   } else if (height == 1) {
18257459e9aSMatthew G Knepley     /* Faces */
18357459e9aSMatthew G Knepley     if (pStart) {*pStart = nC+nV;}
18457459e9aSMatthew G Knepley     if (pEnd)   {*pEnd   = nC+nV+nXF+nYF+nZF;}
18557459e9aSMatthew G Knepley   } else if (height == dim) {
18657459e9aSMatthew G Knepley     /* Vertices */
18757459e9aSMatthew G Knepley     if (pStart) {*pStart = nC;}
18857459e9aSMatthew G Knepley     if (pEnd)   {*pEnd   = nC+nV;}
18957459e9aSMatthew G Knepley   } else if (height < 0) {
19057459e9aSMatthew G Knepley     /* All points */
19157459e9aSMatthew G Knepley     if (pStart) {*pStart = 0;}
19257459e9aSMatthew G Knepley     if (pEnd)   {*pEnd   = nC+nV+nXF+nYF+nZF;}
19357459e9aSMatthew G Knepley   } else {
19457459e9aSMatthew G Knepley     SETERRQ1(((PetscObject) dm)->comm, PETSC_ERR_ARG_OUTOFRANGE, "No points of height %d in the DA", height);
19557459e9aSMatthew G Knepley   }
19657459e9aSMatthew G Knepley   PetscFunctionReturn(0);
19757459e9aSMatthew G Knepley }
19857459e9aSMatthew G Knepley 
19957459e9aSMatthew G Knepley #undef __FUNCT__
200a66d4d66SMatthew G Knepley #define __FUNCT__ "DMDACreateSection"
201a66d4d66SMatthew G Knepley /*@C
202a66d4d66SMatthew G Knepley   DMDACreateSection - Create a PetscSection inside the DMDA that describes data layout. This allows multiple fields with
203a66d4d66SMatthew G Knepley   different numbers of dofs on vertices, cells, and faces in each direction.
204a66d4d66SMatthew G Knepley 
205a66d4d66SMatthew G Knepley   Input Parameters:
206a66d4d66SMatthew G Knepley + dm- The DMDA
207a66d4d66SMatthew G Knepley . numFields - The number of fields
208a66d4d66SMatthew G Knepley . numComp - The number of components in each field, or PETSC_NULL for 1
209a66d4d66SMatthew G Knepley . numVertexDof - The number of dofs per vertex for each field, or PETSC_NULL
210a66d4d66SMatthew G Knepley . numFaceDof - The number of dofs per face for each field and direction, or PETSC_NULL
211a66d4d66SMatthew G Knepley - numCellDof - The number of dofs per cell for each field, or PETSC_NULL
212a66d4d66SMatthew G Knepley 
213a66d4d66SMatthew G Knepley   Level: developer
214a66d4d66SMatthew G Knepley 
215a66d4d66SMatthew G Knepley   Note:
216a66d4d66SMatthew G Knepley   The default DMDA numbering is as follows:
217a66d4d66SMatthew G Knepley 
218a66d4d66SMatthew G Knepley     - Cells:    [0,             nC)
219a66d4d66SMatthew G Knepley     - Vertices: [nC,            nC+nV)
22088ed4aceSMatthew G Knepley     - X-Faces:  [nC+nV,         nC+nV+nXF)         normal is +- x-dir
22188ed4aceSMatthew G Knepley     - Y-Faces:  [nC+nV+nXF,     nC+nV+nXF+nYF)     normal is +- y-dir
22288ed4aceSMatthew G Knepley     - Z-Faces:  [nC+nV+nXF+nYF, nC+nV+nXF+nYF+nZF) normal is +- z-dir
223a66d4d66SMatthew G Knepley 
224a66d4d66SMatthew G Knepley   We interpret the default DMDA partition as a cell partition, and the data assignment as a cell assignment.
225a66d4d66SMatthew G Knepley @*/
22680800b1aSMatthew G Knepley PetscErrorCode DMDACreateSection(DM dm, PetscInt numComp[], PetscInt numVertexDof[], PetscInt numFaceDof[], PetscInt numCellDof[])
227a66d4d66SMatthew G Knepley {
228a66d4d66SMatthew G Knepley   DM_DA         *da  = (DM_DA *) dm->data;
22988ed4aceSMatthew G Knepley   const PetscInt dim = da->dim;
23080800b1aSMatthew G Knepley   PetscInt       numFields, numVertexTotDof = 0, numCellTotDof = 0, numFaceTotDof[3] = {0, 0, 0};
23188ed4aceSMatthew G Knepley   PetscSF        sf;
232feafbc5cSMatthew G Knepley   PetscMPIInt    rank;
23388ed4aceSMatthew G Knepley   const PetscMPIInt *neighbors;
23488ed4aceSMatthew G Knepley   PetscInt      *localPoints;
23588ed4aceSMatthew G Knepley   PetscSFNode   *remotePoints;
236f5eeb527SMatthew G Knepley   PetscInt       nleaves = 0,  nleavesCheck = 0, nL = 0;
23757459e9aSMatthew G Knepley   PetscInt       nC, nVx, nVy, nVz, nV, nxF, nXF, nyF, nYF, nzF, nZF;
23857459e9aSMatthew G Knepley   PetscInt       pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart, yfEnd, zfStart, zfEnd;
23988ed4aceSMatthew G Knepley   PetscInt       f, v, c, xf, yf, zf, xn, yn, zn;
240a66d4d66SMatthew G Knepley   PetscErrorCode ierr;
241a66d4d66SMatthew G Knepley 
242a66d4d66SMatthew G Knepley   PetscFunctionBegin;
243a66d4d66SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
244feafbc5cSMatthew G Knepley   ierr = MPI_Comm_rank(((PetscObject) dm)->comm, &rank);CHKERRQ(ierr);
24557459e9aSMatthew G Knepley   ierr = DMDAGetNumCells(dm, &nC);CHKERRQ(ierr);
24657459e9aSMatthew G Knepley   ierr = DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV);CHKERRQ(ierr);
24757459e9aSMatthew G Knepley   ierr = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr);
24857459e9aSMatthew G Knepley   ierr = DMDAGetHeightStratum(dm, -1,  &pStart, &pEnd);CHKERRQ(ierr);
24957459e9aSMatthew G Knepley   ierr = DMDAGetHeightStratum(dm, 0,   &cStart, &cEnd);CHKERRQ(ierr);
25057459e9aSMatthew G Knepley   ierr = DMDAGetHeightStratum(dm, 1,   &fStart, &fEnd);CHKERRQ(ierr);
25157459e9aSMatthew G Knepley   ierr = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr);
25257459e9aSMatthew G Knepley   xfStart = vEnd;  xfEnd = xfStart+nXF;
25357459e9aSMatthew G Knepley   yfStart = xfEnd; yfEnd = yfStart+nYF;
25457459e9aSMatthew G Knepley   zfStart = yfEnd; zfEnd = zfStart+nZF;
25557459e9aSMatthew G Knepley   if (zfEnd != fEnd) SETERRQ2(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "Invalid face end %d, should be %d", zfEnd, fEnd);
25688ed4aceSMatthew G Knepley   /* Create local section */
25780800b1aSMatthew G Knepley   ierr = DMDAGetInfo(dm, 0,0,0,0,0,0,0, &numFields, 0,0,0,0,0);CHKERRQ(ierr);
258a66d4d66SMatthew G Knepley   for (f = 0; f < numFields; ++f) {
259a66d4d66SMatthew G Knepley     if (numVertexDof) {numVertexTotDof  += numVertexDof[f];}
260a66d4d66SMatthew G Knepley     if (numCellDof)   {numCellTotDof    += numCellDof[f];}
26188ed4aceSMatthew G Knepley     if (numFaceDof)   {numFaceTotDof[0] += numFaceDof[f*dim+0];
26288ed4aceSMatthew G Knepley                        numFaceTotDof[1] += dim > 1 ? numFaceDof[f*dim+1] : 0;
26388ed4aceSMatthew G Knepley                        numFaceTotDof[2] += dim > 2 ? numFaceDof[f*dim+2] : 0;}
264a66d4d66SMatthew G Knepley   }
26588ed4aceSMatthew G Knepley   ierr = PetscSectionCreate(((PetscObject) dm)->comm, &dm->defaultSection);CHKERRQ(ierr);
266a66d4d66SMatthew G Knepley   if (numFields > 1) {
26788ed4aceSMatthew G Knepley     ierr = PetscSectionSetNumFields(dm->defaultSection, numFields);CHKERRQ(ierr);
268a66d4d66SMatthew G Knepley     for (f = 0; f < numFields; ++f) {
26980800b1aSMatthew G Knepley       const char *name;
27080800b1aSMatthew G Knepley 
27180800b1aSMatthew G Knepley       ierr = DMDAGetFieldName(dm, f, &name);CHKERRQ(ierr);
27288ed4aceSMatthew G Knepley       ierr = PetscSectionSetFieldName(dm->defaultSection, f, name);CHKERRQ(ierr);
27380800b1aSMatthew G Knepley       if (numComp) {
27488ed4aceSMatthew G Knepley         ierr = PetscSectionSetFieldComponents(dm->defaultSection, f, numComp[f]);CHKERRQ(ierr);
275a66d4d66SMatthew G Knepley       }
276a66d4d66SMatthew G Knepley     }
277a66d4d66SMatthew G Knepley   } else {
278a66d4d66SMatthew G Knepley     numFields = 0;
279a66d4d66SMatthew G Knepley   }
28088ed4aceSMatthew G Knepley   ierr = PetscSectionSetChart(dm->defaultSection, pStart, pEnd);CHKERRQ(ierr);
281a66d4d66SMatthew G Knepley   if (numVertexDof) {
282a66d4d66SMatthew G Knepley     for (v = vStart; v < vEnd; ++v) {
283a66d4d66SMatthew G Knepley       for (f = 0; f < numFields; ++f) {
28488ed4aceSMatthew G Knepley         ierr = PetscSectionSetFieldDof(dm->defaultSection, v, f, numVertexDof[f]);CHKERRQ(ierr);
285a66d4d66SMatthew G Knepley       }
28688ed4aceSMatthew G Knepley       ierr = PetscSectionSetDof(dm->defaultSection, v, numVertexTotDof);CHKERRQ(ierr);
287a66d4d66SMatthew G Knepley     }
288a66d4d66SMatthew G Knepley   }
289a66d4d66SMatthew G Knepley   if (numFaceDof) {
290a66d4d66SMatthew G Knepley     for (xf = xfStart; xf < xfEnd; ++xf) {
291a66d4d66SMatthew G Knepley       for (f = 0; f < numFields; ++f) {
29288ed4aceSMatthew G Knepley         ierr = PetscSectionSetFieldDof(dm->defaultSection, xf, f, numFaceDof[f*dim+0]);CHKERRQ(ierr);
293a66d4d66SMatthew G Knepley       }
29488ed4aceSMatthew G Knepley       ierr = PetscSectionSetDof(dm->defaultSection, xf, numFaceTotDof[0]);CHKERRQ(ierr);
295a66d4d66SMatthew G Knepley     }
296a66d4d66SMatthew G Knepley     for (yf = yfStart; yf < yfEnd; ++yf) {
297a66d4d66SMatthew G Knepley       for (f = 0; f < numFields; ++f) {
29888ed4aceSMatthew G Knepley         ierr = PetscSectionSetFieldDof(dm->defaultSection, yf, f, numFaceDof[f*dim+1]);CHKERRQ(ierr);
299a66d4d66SMatthew G Knepley       }
30088ed4aceSMatthew G Knepley       ierr = PetscSectionSetDof(dm->defaultSection, yf, numFaceTotDof[1]);CHKERRQ(ierr);
301a66d4d66SMatthew G Knepley     }
302a66d4d66SMatthew G Knepley     for (zf = zfStart; zf < zfEnd; ++zf) {
303a66d4d66SMatthew G Knepley       for (f = 0; f < numFields; ++f) {
30488ed4aceSMatthew G Knepley         ierr = PetscSectionSetFieldDof(dm->defaultSection, zf, f, numFaceDof[f*dim+2]);CHKERRQ(ierr);
305a66d4d66SMatthew G Knepley       }
30688ed4aceSMatthew G Knepley       ierr = PetscSectionSetDof(dm->defaultSection, zf, numFaceTotDof[2]);CHKERRQ(ierr);
307a66d4d66SMatthew G Knepley     }
308a66d4d66SMatthew G Knepley   }
309a66d4d66SMatthew G Knepley   if (numCellDof) {
310a66d4d66SMatthew G Knepley     for (c = cStart; c < cEnd; ++c) {
311a66d4d66SMatthew G Knepley       for (f = 0; f < numFields; ++f) {
31288ed4aceSMatthew G Knepley         ierr = PetscSectionSetFieldDof(dm->defaultSection, c, f, numCellDof[f]);CHKERRQ(ierr);
313a66d4d66SMatthew G Knepley       }
31488ed4aceSMatthew G Knepley       ierr = PetscSectionSetDof(dm->defaultSection, c, numCellTotDof);CHKERRQ(ierr);
315a66d4d66SMatthew G Knepley     }
316a66d4d66SMatthew G Knepley   }
31788ed4aceSMatthew G Knepley   ierr = PetscSectionSetUp(dm->defaultSection);CHKERRQ(ierr);
31888ed4aceSMatthew G Knepley   /* Create mesh point SF */
31988ed4aceSMatthew G Knepley   ierr = DMDAGetNeighbors(dm, &neighbors);CHKERRQ(ierr);
32088ed4aceSMatthew G Knepley   for (zn = 0; zn < (dim > 2 ? 3 : 1); ++zn) {
32188ed4aceSMatthew G Knepley     for (yn = 0; yn < (dim > 1 ? 3 : 1); ++yn) {
32288ed4aceSMatthew G Knepley       for (xn = 0; xn < 3; ++xn) {
3237128ae9fSMatthew G Knepley         const PetscInt xp = xn-1, yp = dim > 1 ? yn-1 : 0, zp = dim > 2 ? zn-1 : 0;
32488ed4aceSMatthew G Knepley         const PetscInt neighbor = neighbors[(zn*3+yn)*3+xn];
32588ed4aceSMatthew G Knepley 
3263814d604SMatthew G Knepley         if (neighbor >= 0 && neighbor < rank) {
327feafbc5cSMatthew G Knepley           nleaves += (!xp ? nVx : 1) * (!yp ? nVy : 1) * (!zp ? nVz : 1); /* vertices */
32888ed4aceSMatthew G Knepley           if (xp && !yp && !zp) {
32988ed4aceSMatthew G Knepley             nleaves += nxF; /* x faces */
33088ed4aceSMatthew G Knepley           } else if (yp && !zp && !xp) {
33188ed4aceSMatthew G Knepley             nleaves += nyF; /* y faces */
33288ed4aceSMatthew G Knepley           } else if (zp && !xp && !yp) {
33388ed4aceSMatthew G Knepley             nleaves += nzF; /* z faces */
33488ed4aceSMatthew G Knepley           }
33588ed4aceSMatthew G Knepley         }
33688ed4aceSMatthew G Knepley       }
33788ed4aceSMatthew G Knepley     }
33888ed4aceSMatthew G Knepley   }
33988ed4aceSMatthew G Knepley   ierr = PetscMalloc2(nleaves,PetscInt,&localPoints,nleaves,PetscSFNode,&remotePoints);CHKERRQ(ierr);
34088ed4aceSMatthew G Knepley   for (zn = 0; zn < (dim > 2 ? 3 : 1); ++zn) {
34188ed4aceSMatthew G Knepley     for (yn = 0; yn < (dim > 1 ? 3 : 1); ++yn) {
34288ed4aceSMatthew G Knepley       for (xn = 0; xn < 3; ++xn) {
3437128ae9fSMatthew G Knepley         const PetscInt xp = xn-1, yp = dim > 1 ? yn-1 : 0, zp = dim > 2 ? zn-1 : 0;
34488ed4aceSMatthew G Knepley         const PetscInt neighbor = neighbors[(zn*3+yn)*3+xn];
345f5eeb527SMatthew G Knepley         PetscInt       xv, yv, zv;
34688ed4aceSMatthew G Knepley 
3473814d604SMatthew G Knepley         if (neighbor >= 0 && neighbor < rank) {
34888ed4aceSMatthew G Knepley           if (xp < 0) { /* left */
34988ed4aceSMatthew G Knepley             if (yp < 0) { /* bottom */
35088ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
351f5eeb527SMatthew G Knepley                 const PetscInt localVertex  = (      0*nVy +     0)*nVx +     0 + nC;
352f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = ((nVz-1)*nVy + nVy-1)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
353628e3b0dSSatish Balay                 nleavesCheck += 1; /* left bottom back vertex */
354f5eeb527SMatthew G Knepley 
355f5eeb527SMatthew G Knepley                 localPoints[nL]        = localVertex;
356f5eeb527SMatthew G Knepley                 remotePoints[nL].rank  = neighbor;
357f5eeb527SMatthew G Knepley                 remotePoints[nL].index = remoteVertex;
358f5eeb527SMatthew G Knepley                 ++nL;
35988ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
360f5eeb527SMatthew G Knepley                 const PetscInt localVertex  = ((nVz-1)*nVy +     0)*nVx +     0 + nC;
361f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = (      0*nVy + nVy-1)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
362628e3b0dSSatish Balay                 nleavesCheck += 1; /* left bottom front vertex */
363f5eeb527SMatthew G Knepley 
364f5eeb527SMatthew G Knepley                 localPoints[nL]        = localVertex;
365f5eeb527SMatthew G Knepley                 remotePoints[nL].rank  = neighbor;
366f5eeb527SMatthew G Knepley                 remotePoints[nL].index = remoteVertex;
367f5eeb527SMatthew G Knepley                 ++nL;
36888ed4aceSMatthew G Knepley               } else {
36988ed4aceSMatthew G Knepley                 nleavesCheck += nVz; /* left bottom vertices */
370f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv, ++nL) {
371f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = (zv*nVy +     0)*nVx +     0 + nC;
372f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (zv*nVy + nVy-1)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
373f5eeb527SMatthew G Knepley 
374f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
375f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
376f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
377f5eeb527SMatthew G Knepley                 }
37888ed4aceSMatthew G Knepley               }
37988ed4aceSMatthew G Knepley             } else if (yp > 0) { /* top */
38088ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
381f5eeb527SMatthew G Knepley                 const PetscInt localVertex  = (      0*nVy + nVy-1)*nVx +     0 + nC;
382f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = ((nVz-1)*nVy +     0)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
383628e3b0dSSatish Balay                 nleavesCheck += 1; /* left top back vertex */
384f5eeb527SMatthew G Knepley 
385f5eeb527SMatthew G Knepley                 localPoints[nL]        = localVertex;
386f5eeb527SMatthew G Knepley                 remotePoints[nL].rank  = neighbor;
387f5eeb527SMatthew G Knepley                 remotePoints[nL].index = remoteVertex;
388f5eeb527SMatthew G Knepley                 ++nL;
38988ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
390f5eeb527SMatthew G Knepley                 const PetscInt localVertex  = ((nVz-1)*nVy + nVy-1)*nVx +     0 + nC;
391f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = (      0*nVy +     0)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
392628e3b0dSSatish Balay                 nleavesCheck += 1; /* left top front vertex */
393f5eeb527SMatthew G Knepley 
394f5eeb527SMatthew G Knepley                 localPoints[nL]        = localVertex;
395f5eeb527SMatthew G Knepley                 remotePoints[nL].rank  = neighbor;
396f5eeb527SMatthew G Knepley                 remotePoints[nL].index = remoteVertex;
397f5eeb527SMatthew G Knepley                 ++nL;
39888ed4aceSMatthew G Knepley               } else {
39988ed4aceSMatthew G Knepley                 nleavesCheck += nVz; /* left top vertices */
400f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv, ++nL) {
401f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = (zv*nVy + nVy-1)*nVx +     0 + nC;
402f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (zv*nVy +     0)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
403f5eeb527SMatthew G Knepley 
404f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
405f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
406f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
407f5eeb527SMatthew G Knepley                 }
40888ed4aceSMatthew G Knepley               }
40988ed4aceSMatthew G Knepley             } else {
41088ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
41188ed4aceSMatthew G Knepley                 nleavesCheck += nVy; /* left back vertices */
412f5eeb527SMatthew G Knepley                 for (yv = 0; yv < nVy; ++yv, ++nL) {
413f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = (      0*nVy + yv)*nVx +     0 + nC;
414f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = ((nVz-1)*nVy + yv)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
415f5eeb527SMatthew G Knepley 
416f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
417f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
418f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
419f5eeb527SMatthew G Knepley                 }
42088ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
42188ed4aceSMatthew G Knepley                 nleavesCheck += nVy; /* left front vertices */
422f5eeb527SMatthew G Knepley                 for (yv = 0; yv < nVy; ++yv, ++nL) {
423f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = ((nVz-1)*nVy + yv)*nVx +     0 + nC;
424f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (      0*nVy + yv)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
425f5eeb527SMatthew G Knepley 
426f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
427f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
428f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
429f5eeb527SMatthew G Knepley                 }
43088ed4aceSMatthew G Knepley               } else {
43188ed4aceSMatthew G Knepley                 nleavesCheck += nVy*nVz; /* left vertices */
432f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv) {
433f5eeb527SMatthew G Knepley                   for (yv = 0; yv < nVy; ++yv, ++nL) {
434f5eeb527SMatthew G Knepley                     const PetscInt localVertex  = (zv*nVy + yv)*nVx +     0 + nC;
435f5eeb527SMatthew G Knepley                     const PetscInt remoteVertex = (zv*nVy + yv)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
436f5eeb527SMatthew G Knepley 
437f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
438f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
439f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
440f5eeb527SMatthew G Knepley                   }
441f5eeb527SMatthew G Knepley                 }
44288ed4aceSMatthew G Knepley                 nleavesCheck += nxF;     /* left faces */
443f5eeb527SMatthew G Knepley                 for (xf = 0; xf < nxF; ++xf, ++nL) {
444f5eeb527SMatthew G Knepley                   /* THIS IS WRONG */
445f5eeb527SMatthew G Knepley                   const PetscInt localFace  = 0 + nC+nV;
446f5eeb527SMatthew G Knepley                   const PetscInt remoteFace = 0 + nC+nV;
447f5eeb527SMatthew G Knepley 
448f5eeb527SMatthew G Knepley                   localPoints[nL]        = localFace;
449f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
450f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteFace;
451f5eeb527SMatthew G Knepley                 }
45288ed4aceSMatthew G Knepley               }
45388ed4aceSMatthew G Knepley             }
45488ed4aceSMatthew G Knepley           } else if (xp > 0) { /* right */
45588ed4aceSMatthew G Knepley             if (yp < 0) { /* bottom */
45688ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
457f5eeb527SMatthew G Knepley                 const PetscInt localVertex  = (      0*nVy +     0)*nVx + nVx-1 + nC;
458f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = ((nVz-1)*nVy + nVy-1)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
459628e3b0dSSatish Balay                 nleavesCheck += 1; /* right bottom back vertex */
460f5eeb527SMatthew G Knepley 
461f5eeb527SMatthew G Knepley                 localPoints[nL]        = localVertex;
462f5eeb527SMatthew G Knepley                 remotePoints[nL].rank  = neighbor;
463f5eeb527SMatthew G Knepley                 remotePoints[nL].index = remoteVertex;
464f5eeb527SMatthew G Knepley                 ++nL;
46588ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
466f5eeb527SMatthew G Knepley                 const PetscInt localVertex  = ((nVz-1)*nVy +     0)*nVx + nVx-1 + nC;
467f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = (      0*nVy + nVy-1)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
468628e3b0dSSatish Balay                 nleavesCheck += 1; /* right bottom front vertex */
469f5eeb527SMatthew G Knepley 
470f5eeb527SMatthew G Knepley                 localPoints[nL]        = localVertex;
471f5eeb527SMatthew G Knepley                 remotePoints[nL].rank  = neighbor;
472f5eeb527SMatthew G Knepley                 remotePoints[nL].index = remoteVertex;
473f5eeb527SMatthew G Knepley                 ++nL;
47488ed4aceSMatthew G Knepley               } else {
47588ed4aceSMatthew G Knepley                 nleavesCheck += nVz; /* right bottom vertices */
476f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv, ++nL) {
477f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = (zv*nVy +     0)*nVx + nVx-1 + nC;
478f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (zv*nVy + nVy-1)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
479f5eeb527SMatthew G Knepley 
480f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
481f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
482f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
483f5eeb527SMatthew G Knepley                 }
48488ed4aceSMatthew G Knepley               }
48588ed4aceSMatthew G Knepley             } else if (yp > 0) { /* top */
48688ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
487f5eeb527SMatthew G Knepley                 const PetscInt localVertex  = (      0*nVy + nVy-1)*nVx + nVx-1 + nC;
488f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = ((nVz-1)*nVy +     0)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
489628e3b0dSSatish Balay                 nleavesCheck += 1; /* right top back vertex */
490f5eeb527SMatthew G Knepley 
491f5eeb527SMatthew G Knepley                 localPoints[nL]        = localVertex;
492f5eeb527SMatthew G Knepley                 remotePoints[nL].rank  = neighbor;
493f5eeb527SMatthew G Knepley                 remotePoints[nL].index = remoteVertex;
494f5eeb527SMatthew G Knepley                 ++nL;
49588ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
496f5eeb527SMatthew G Knepley                 const PetscInt localVertex  = ((nVz-1)*nVy + nVy-1)*nVx + nVx-1 + nC;
497f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = (      0*nVy +     0)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
498628e3b0dSSatish Balay                 nleavesCheck += 1; /* right top front vertex */
499f5eeb527SMatthew G Knepley 
500f5eeb527SMatthew G Knepley                 localPoints[nL]        = localVertex;
501f5eeb527SMatthew G Knepley                 remotePoints[nL].rank  = neighbor;
502f5eeb527SMatthew G Knepley                 remotePoints[nL].index = remoteVertex;
503f5eeb527SMatthew G Knepley                 ++nL;
50488ed4aceSMatthew G Knepley               } else {
50588ed4aceSMatthew G Knepley                 nleavesCheck += nVz; /* right top vertices */
506f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv, ++nL) {
507f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = (zv*nVy + nVy-1)*nVx + nVx-1 + nC;
508f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (zv*nVy +     0)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
509f5eeb527SMatthew G Knepley 
510f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
511f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
512f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
513f5eeb527SMatthew G Knepley                 }
51488ed4aceSMatthew G Knepley               }
51588ed4aceSMatthew G Knepley             } else {
51688ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
51788ed4aceSMatthew G Knepley                 nleavesCheck += nVy; /* right back vertices */
518f5eeb527SMatthew G Knepley                 for (yv = 0; yv < nVy; ++yv, ++nL) {
519f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = (      0*nVy + yv)*nVx + nVx-1 + nC;
520f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = ((nVz-1)*nVy + yv)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
521f5eeb527SMatthew G Knepley 
522f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
523f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
524f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
525f5eeb527SMatthew G Knepley                 }
52688ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
52788ed4aceSMatthew G Knepley                 nleavesCheck += nVy; /* right front vertices */
528f5eeb527SMatthew G Knepley                 for (yv = 0; yv < nVy; ++yv, ++nL) {
529f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = ((nVz-1)*nVy + yv)*nVx + nVx-1 + nC;
530f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (      0*nVy + yv)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
531f5eeb527SMatthew G Knepley 
532f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
533f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
534f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
535f5eeb527SMatthew G Knepley                 }
53688ed4aceSMatthew G Knepley               } else {
53788ed4aceSMatthew G Knepley                 nleavesCheck += nVy*nVz; /* right vertices */
538f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv) {
539f5eeb527SMatthew G Knepley                   for (yv = 0; yv < nVy; ++yv, ++nL) {
540f5eeb527SMatthew G Knepley                     const PetscInt localVertex  = (zv*nVy + yv)*nVx + nVx-1 + nC;
541f5eeb527SMatthew G Knepley                     const PetscInt remoteVertex = (zv*nVy + yv)*nVx + 0     + nC; /* TODO: Correct this for neighbor sizes */
542f5eeb527SMatthew G Knepley 
543f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
544f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
545f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
546f5eeb527SMatthew G Knepley                   }
547f5eeb527SMatthew G Knepley                 }
54888ed4aceSMatthew G Knepley                 nleavesCheck += nxF;     /* right faces */
549f5eeb527SMatthew G Knepley                 for (xf = 0; xf < nxF; ++xf, ++nL) {
550f5eeb527SMatthew G Knepley                   /* THIS IS WRONG */
551f5eeb527SMatthew G Knepley                   const PetscInt localFace  = 0 + nC+nV;
552f5eeb527SMatthew G Knepley                   const PetscInt remoteFace = 0 + nC+nV;
553f5eeb527SMatthew G Knepley 
554f5eeb527SMatthew G Knepley                   localPoints[nL]        = localFace;
555f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
556f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteFace;
557f5eeb527SMatthew G Knepley                 }
55888ed4aceSMatthew G Knepley               }
55988ed4aceSMatthew G Knepley             }
56088ed4aceSMatthew G Knepley           } else {
56188ed4aceSMatthew G Knepley             if (yp < 0) { /* bottom */
56288ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
56388ed4aceSMatthew G Knepley                 nleavesCheck += nVx; /* bottom back vertices */
564f5eeb527SMatthew G Knepley                 for (xv = 0; xv < nVx; ++xv, ++nL) {
565f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = (      0*nVy +     0)*nVx + xv + nC;
566f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = ((nVz-1)*nVy + nVy-1)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
567f5eeb527SMatthew G Knepley 
568f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
569f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
570f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
571f5eeb527SMatthew G Knepley                 }
57288ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
57388ed4aceSMatthew G Knepley                 nleavesCheck += nVx; /* bottom front vertices */
574f5eeb527SMatthew G Knepley                 for (xv = 0; xv < nVx; ++xv, ++nL) {
575f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = ((nVz-1)*nVy +     0)*nVx + xv + nC;
576f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (      0*nVy + nVy-1)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
577f5eeb527SMatthew G Knepley 
578f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
579f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
580f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
581f5eeb527SMatthew G Knepley                 }
58288ed4aceSMatthew G Knepley               } else {
58388ed4aceSMatthew G Knepley                 nleavesCheck += nVx*nVz; /* bottom vertices */
584f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv) {
585f5eeb527SMatthew G Knepley                   for (xv = 0; xv < nVx; ++xv, ++nL) {
586f5eeb527SMatthew G Knepley                     const PetscInt localVertex  = (zv*nVy +     0)*nVx + xv + nC;
587f5eeb527SMatthew G Knepley                     const PetscInt remoteVertex = (zv*nVy + nVy-1)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
588f5eeb527SMatthew G Knepley 
589f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
590f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
591f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
592f5eeb527SMatthew G Knepley                   }
593f5eeb527SMatthew G Knepley                 }
59488ed4aceSMatthew G Knepley                 nleavesCheck += nyF;     /* bottom faces */
595f5eeb527SMatthew G Knepley                 for (yf = 0; yf < nyF; ++yf, ++nL) {
596f5eeb527SMatthew G Knepley                   /* THIS IS WRONG */
597f5eeb527SMatthew G Knepley                   const PetscInt localFace  = 0 + nC+nV;
598f5eeb527SMatthew G Knepley                   const PetscInt remoteFace = 0 + nC+nV;
599f5eeb527SMatthew G Knepley 
600f5eeb527SMatthew G Knepley                   localPoints[nL]        = localFace;
601f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
602f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteFace;
603f5eeb527SMatthew G Knepley                 }
60488ed4aceSMatthew G Knepley               }
60588ed4aceSMatthew G Knepley             } else if (yp > 0) { /* top */
60688ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
60788ed4aceSMatthew G Knepley                 nleavesCheck += nVx; /* top back vertices */
608f5eeb527SMatthew G Knepley                 for (xv = 0; xv < nVx; ++xv, ++nL) {
609f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = (      0*nVy + nVy-1)*nVx + xv + nC;
610f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = ((nVz-1)*nVy +     0)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
611f5eeb527SMatthew G Knepley 
612f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
613f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
614f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
615f5eeb527SMatthew G Knepley                 }
61688ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
61788ed4aceSMatthew G Knepley                 nleavesCheck += nVx; /* top front vertices */
618f5eeb527SMatthew G Knepley                 for (xv = 0; xv < nVx; ++xv, ++nL) {
619f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = ((nVz-1)*nVy + nVy-1)*nVx + xv + nC;
620f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (      0*nVy +     0)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
621f5eeb527SMatthew G Knepley 
622f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
623f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
624f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
625f5eeb527SMatthew G Knepley                 }
62688ed4aceSMatthew G Knepley               } else {
62788ed4aceSMatthew G Knepley                 nleavesCheck += nVx*nVz; /* top vertices */
628f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv) {
629f5eeb527SMatthew G Knepley                   for (xv = 0; xv < nVx; ++xv, ++nL) {
630f5eeb527SMatthew G Knepley                     const PetscInt localVertex  = (zv*nVy + nVy-1)*nVx + xv + nC;
631f5eeb527SMatthew G Knepley                     const PetscInt remoteVertex = (zv*nVy +     0)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
632f5eeb527SMatthew G Knepley 
633f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
634f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
635f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
636f5eeb527SMatthew G Knepley                   }
637f5eeb527SMatthew G Knepley                 }
63888ed4aceSMatthew G Knepley                 nleavesCheck += nyF;     /* top faces */
639f5eeb527SMatthew G Knepley                 for (yf = 0; yf < nyF; ++yf, ++nL) {
640f5eeb527SMatthew G Knepley                   /* THIS IS WRONG */
641f5eeb527SMatthew G Knepley                   const PetscInt localFace  = 0 + nC+nV;
642f5eeb527SMatthew G Knepley                   const PetscInt remoteFace = 0 + nC+nV;
643f5eeb527SMatthew G Knepley 
644f5eeb527SMatthew G Knepley                   localPoints[nL]        = localFace;
645f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
646f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteFace;
647f5eeb527SMatthew G Knepley                 }
64888ed4aceSMatthew G Knepley               }
64988ed4aceSMatthew G Knepley             } else {
65088ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
65188ed4aceSMatthew G Knepley                 nleavesCheck += nVx*nVy; /* back vertices */
652f5eeb527SMatthew G Knepley                 for (yv = 0; yv < nVy; ++yv) {
653f5eeb527SMatthew G Knepley                   for (xv = 0; xv < nVx; ++xv, ++nL) {
654f5eeb527SMatthew G Knepley                     const PetscInt localVertex  = (      0*nVy + yv)*nVx + xv + nC;
655f5eeb527SMatthew G Knepley                     const PetscInt remoteVertex = ((nVz-1)*nVy + yv)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
656f5eeb527SMatthew G Knepley 
657f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
658f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
659f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
660f5eeb527SMatthew G Knepley                   }
661f5eeb527SMatthew G Knepley                 }
66288ed4aceSMatthew G Knepley                 nleavesCheck += nzF;     /* back faces */
663f5eeb527SMatthew G Knepley                 for (zf = 0; zf < nzF; ++zf, ++nL) {
664f5eeb527SMatthew G Knepley                   /* THIS IS WRONG */
665f5eeb527SMatthew G Knepley                   const PetscInt localFace  = 0 + nC+nV;
666f5eeb527SMatthew G Knepley                   const PetscInt remoteFace = 0 + nC+nV;
667f5eeb527SMatthew G Knepley 
668f5eeb527SMatthew G Knepley                   localPoints[nL]        = localFace;
669f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
670f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteFace;
671f5eeb527SMatthew G Knepley                 }
67288ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
67388ed4aceSMatthew G Knepley                 nleavesCheck += nVx*nVy; /* front vertices */
674f5eeb527SMatthew G Knepley                 for (yv = 0; yv < nVy; ++yv) {
675f5eeb527SMatthew G Knepley                   for (xv = 0; xv < nVx; ++xv, ++nL) {
676f5eeb527SMatthew G Knepley                     const PetscInt localVertex  = ((nVz-1)*nVy + yv)*nVx + xv + nC;
677f5eeb527SMatthew G Knepley                     const PetscInt remoteVertex = (      0*nVy + yv)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
678f5eeb527SMatthew G Knepley 
679f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
680f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
681f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
682f5eeb527SMatthew G Knepley                   }
683f5eeb527SMatthew G Knepley                 }
68488ed4aceSMatthew G Knepley                 nleavesCheck += nzF;     /* front faces */
685f5eeb527SMatthew G Knepley                 for (zf = 0; zf < nzF; ++zf, ++nL) {
686f5eeb527SMatthew G Knepley                   /* THIS IS WRONG */
687f5eeb527SMatthew G Knepley                   const PetscInt localFace  = 0 + nC+nV;
688f5eeb527SMatthew G Knepley                   const PetscInt remoteFace = 0 + nC+nV;
689f5eeb527SMatthew G Knepley 
690f5eeb527SMatthew G Knepley                   localPoints[nL]        = localFace;
691f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
692f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteFace;
693f5eeb527SMatthew G Knepley                 }
69488ed4aceSMatthew G Knepley               } else {
69588ed4aceSMatthew G Knepley                 /* Nothing is shared from the interior */
69688ed4aceSMatthew G Knepley               }
69788ed4aceSMatthew G Knepley             }
69888ed4aceSMatthew G Knepley           }
69988ed4aceSMatthew G Knepley         }
70088ed4aceSMatthew G Knepley       }
70188ed4aceSMatthew G Knepley     }
70288ed4aceSMatthew G Knepley   }
7033814d604SMatthew G Knepley   /* TODO: Remove duplication in leaf determination */
70488ed4aceSMatthew G Knepley   if (nleaves != nleavesCheck) SETERRQ2(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "The number of leaves %d did not match the number of remote leaves %d", nleaves, nleavesCheck);
70588ed4aceSMatthew G Knepley   ierr = PetscSFCreate(((PetscObject) dm)->comm, &sf);CHKERRQ(ierr);
70688ed4aceSMatthew G Knepley   ierr = PetscSFSetGraph(sf, pEnd, nleaves, localPoints, PETSC_OWN_POINTER, remotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr);
70788ed4aceSMatthew G Knepley   /* Create global section */
708eaf8d80aSMatthew G Knepley   ierr = PetscSectionCreateGlobalSection(dm->defaultSection, sf, PETSC_FALSE, &dm->defaultGlobalSection);CHKERRQ(ierr);
70988ed4aceSMatthew G Knepley   ierr = PetscSFDestroy(&sf);CHKERRQ(ierr);
71088ed4aceSMatthew G Knepley   /* Create default SF */
71188ed4aceSMatthew G Knepley   ierr = DMCreateDefaultSF(dm, dm->defaultSection, dm->defaultGlobalSection);CHKERRQ(ierr);
712a66d4d66SMatthew G Knepley   PetscFunctionReturn(0);
713a66d4d66SMatthew G Knepley }
714a66d4d66SMatthew G Knepley 
71547c6ae99SBarry Smith /* ------------------------------------------------------------------- */
71647c6ae99SBarry Smith 
71747c6ae99SBarry Smith #undef __FUNCT__
718aa219208SBarry Smith #define __FUNCT__ "DMDAGetArray"
71947c6ae99SBarry Smith /*@C
720aa219208SBarry Smith      DMDAGetArray - Gets a work array for a DMDA
72147c6ae99SBarry Smith 
72247c6ae99SBarry Smith     Input Parameter:
72347c6ae99SBarry Smith +    da - information about my local patch
72447c6ae99SBarry Smith -    ghosted - do you want arrays for the ghosted or nonghosted patch
72547c6ae99SBarry Smith 
72647c6ae99SBarry Smith     Output Parameters:
72747c6ae99SBarry Smith .    vptr - array data structured
72847c6ae99SBarry Smith 
72947c6ae99SBarry Smith     Note:  The vector values are NOT initialized and may have garbage in them, so you may need
73047c6ae99SBarry Smith            to zero them.
73147c6ae99SBarry Smith 
73247c6ae99SBarry Smith   Level: advanced
73347c6ae99SBarry Smith 
734*bcaeba4dSBarry Smith .seealso: DMDARestoreArray()
73547c6ae99SBarry Smith 
73647c6ae99SBarry Smith @*/
7377087cfbeSBarry Smith PetscErrorCode  DMDAGetArray(DM da,PetscBool  ghosted,void *vptr)
73847c6ae99SBarry Smith {
73947c6ae99SBarry Smith   PetscErrorCode ierr;
74047c6ae99SBarry Smith   PetscInt       j,i,xs,ys,xm,ym,zs,zm;
74147c6ae99SBarry Smith   char           *iarray_start;
74247c6ae99SBarry Smith   void           **iptr = (void**)vptr;
74347c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
74447c6ae99SBarry Smith 
74547c6ae99SBarry Smith   PetscFunctionBegin;
74647c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
74747c6ae99SBarry Smith   if (ghosted) {
748aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
74947c6ae99SBarry Smith       if (dd->arrayghostedin[i]) {
75047c6ae99SBarry Smith         *iptr                 = dd->arrayghostedin[i];
75147c6ae99SBarry Smith         iarray_start          = (char*)dd->startghostedin[i];
75247c6ae99SBarry Smith         dd->arrayghostedin[i] = PETSC_NULL;
75347c6ae99SBarry Smith         dd->startghostedin[i] = PETSC_NULL;
75447c6ae99SBarry Smith 
75547c6ae99SBarry Smith         goto done;
75647c6ae99SBarry Smith       }
75747c6ae99SBarry Smith     }
75847c6ae99SBarry Smith     xs = dd->Xs;
75947c6ae99SBarry Smith     ys = dd->Ys;
76047c6ae99SBarry Smith     zs = dd->Zs;
76147c6ae99SBarry Smith     xm = dd->Xe-dd->Xs;
76247c6ae99SBarry Smith     ym = dd->Ye-dd->Ys;
76347c6ae99SBarry Smith     zm = dd->Ze-dd->Zs;
76447c6ae99SBarry Smith   } else {
765aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
76647c6ae99SBarry Smith       if (dd->arrayin[i]) {
76747c6ae99SBarry Smith         *iptr          = dd->arrayin[i];
76847c6ae99SBarry Smith         iarray_start   = (char*)dd->startin[i];
76947c6ae99SBarry Smith         dd->arrayin[i] = PETSC_NULL;
77047c6ae99SBarry Smith         dd->startin[i] = PETSC_NULL;
77147c6ae99SBarry Smith 
77247c6ae99SBarry Smith         goto done;
77347c6ae99SBarry Smith       }
77447c6ae99SBarry Smith     }
77547c6ae99SBarry Smith     xs = dd->xs;
77647c6ae99SBarry Smith     ys = dd->ys;
77747c6ae99SBarry Smith     zs = dd->zs;
77847c6ae99SBarry Smith     xm = dd->xe-dd->xs;
77947c6ae99SBarry Smith     ym = dd->ye-dd->ys;
78047c6ae99SBarry Smith     zm = dd->ze-dd->zs;
78147c6ae99SBarry Smith   }
78247c6ae99SBarry Smith 
78347c6ae99SBarry Smith   switch (dd->dim) {
78447c6ae99SBarry Smith     case 1: {
78547c6ae99SBarry Smith       void *ptr;
78647c6ae99SBarry Smith 
78747c6ae99SBarry Smith       ierr  = PetscMalloc(xm*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
78847c6ae99SBarry Smith 
78947c6ae99SBarry Smith       ptr   = (void*)(iarray_start - xs*sizeof(PetscScalar));
79047c6ae99SBarry Smith       *iptr = (void*)ptr;
79147c6ae99SBarry Smith       break;}
79247c6ae99SBarry Smith     case 2: {
79347c6ae99SBarry Smith       void **ptr;
79447c6ae99SBarry Smith 
79547c6ae99SBarry Smith       ierr  = PetscMalloc((ym+1)*sizeof(void*)+xm*ym*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
79647c6ae99SBarry Smith 
79747c6ae99SBarry Smith       ptr  = (void**)(iarray_start + xm*ym*sizeof(PetscScalar) - ys*sizeof(void*));
79847c6ae99SBarry Smith       for (j=ys;j<ys+ym;j++) {
79947c6ae99SBarry Smith         ptr[j] = iarray_start + sizeof(PetscScalar)*(xm*(j-ys) - xs);
80047c6ae99SBarry Smith       }
80147c6ae99SBarry Smith       *iptr = (void*)ptr;
80247c6ae99SBarry Smith       break;}
80347c6ae99SBarry Smith     case 3: {
80447c6ae99SBarry Smith       void ***ptr,**bptr;
80547c6ae99SBarry Smith 
80647c6ae99SBarry Smith       ierr  = PetscMalloc((zm+1)*sizeof(void **)+(ym*zm+1)*sizeof(void*)+xm*ym*zm*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
80747c6ae99SBarry Smith 
80847c6ae99SBarry Smith       ptr  = (void***)(iarray_start + xm*ym*zm*sizeof(PetscScalar) - zs*sizeof(void*));
80947c6ae99SBarry Smith       bptr = (void**)(iarray_start + xm*ym*zm*sizeof(PetscScalar) + zm*sizeof(void**));
81047c6ae99SBarry Smith       for (i=zs;i<zs+zm;i++) {
81147c6ae99SBarry Smith         ptr[i] = bptr + ((i-zs)*ym - ys);
81247c6ae99SBarry Smith       }
81347c6ae99SBarry Smith       for (i=zs; i<zs+zm; i++) {
81447c6ae99SBarry Smith         for (j=ys; j<ys+ym; j++) {
81547c6ae99SBarry Smith           ptr[i][j] = iarray_start + sizeof(PetscScalar)*(xm*ym*(i-zs) + xm*(j-ys) - xs);
81647c6ae99SBarry Smith         }
81747c6ae99SBarry Smith       }
81847c6ae99SBarry Smith 
81947c6ae99SBarry Smith       *iptr = (void*)ptr;
82047c6ae99SBarry Smith       break;}
82147c6ae99SBarry Smith     default:
82247c6ae99SBarry Smith       SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Dimension %D not supported",dd->dim);
82347c6ae99SBarry Smith   }
82447c6ae99SBarry Smith 
82547c6ae99SBarry Smith   done:
82647c6ae99SBarry Smith   /* add arrays to the checked out list */
82747c6ae99SBarry Smith   if (ghosted) {
828aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
82947c6ae99SBarry Smith       if (!dd->arrayghostedout[i]) {
83047c6ae99SBarry Smith         dd->arrayghostedout[i] = *iptr ;
83147c6ae99SBarry Smith         dd->startghostedout[i] = iarray_start;
83247c6ae99SBarry Smith         break;
83347c6ae99SBarry Smith       }
83447c6ae99SBarry Smith     }
83547c6ae99SBarry Smith   } else {
836aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
83747c6ae99SBarry Smith       if (!dd->arrayout[i]) {
83847c6ae99SBarry Smith         dd->arrayout[i] = *iptr ;
83947c6ae99SBarry Smith         dd->startout[i] = iarray_start;
84047c6ae99SBarry Smith         break;
84147c6ae99SBarry Smith       }
84247c6ae99SBarry Smith     }
84347c6ae99SBarry Smith   }
84447c6ae99SBarry Smith   PetscFunctionReturn(0);
84547c6ae99SBarry Smith }
84647c6ae99SBarry Smith 
84747c6ae99SBarry Smith #undef __FUNCT__
848aa219208SBarry Smith #define __FUNCT__ "DMDARestoreArray"
84947c6ae99SBarry Smith /*@C
850aa219208SBarry Smith      DMDARestoreArray - Restores an array of derivative types for a DMDA
85147c6ae99SBarry Smith 
85247c6ae99SBarry Smith     Input Parameter:
85347c6ae99SBarry Smith +    da - information about my local patch
85447c6ae99SBarry Smith .    ghosted - do you want arrays for the ghosted or nonghosted patch
85547c6ae99SBarry Smith -    vptr - array data structured to be passed to ad_FormFunctionLocal()
85647c6ae99SBarry Smith 
85747c6ae99SBarry Smith      Level: advanced
85847c6ae99SBarry Smith 
859*bcaeba4dSBarry Smith .seealso: DMDAGetArray()
86047c6ae99SBarry Smith 
86147c6ae99SBarry Smith @*/
8627087cfbeSBarry Smith PetscErrorCode  DMDARestoreArray(DM da,PetscBool  ghosted,void *vptr)
86347c6ae99SBarry Smith {
86447c6ae99SBarry Smith   PetscInt  i;
86547c6ae99SBarry Smith   void      **iptr = (void**)vptr,*iarray_start = 0;
86647c6ae99SBarry Smith   DM_DA     *dd = (DM_DA*)da->data;
86747c6ae99SBarry Smith 
86847c6ae99SBarry Smith   PetscFunctionBegin;
86947c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
87047c6ae99SBarry Smith   if (ghosted) {
871aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
87247c6ae99SBarry Smith       if (dd->arrayghostedout[i] == *iptr) {
87347c6ae99SBarry Smith         iarray_start           = dd->startghostedout[i];
87447c6ae99SBarry Smith         dd->arrayghostedout[i] = PETSC_NULL;
87547c6ae99SBarry Smith         dd->startghostedout[i] = PETSC_NULL;
87647c6ae99SBarry Smith         break;
87747c6ae99SBarry Smith       }
87847c6ae99SBarry Smith     }
879aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
88047c6ae99SBarry Smith       if (!dd->arrayghostedin[i]){
88147c6ae99SBarry Smith         dd->arrayghostedin[i] = *iptr;
88247c6ae99SBarry Smith         dd->startghostedin[i] = iarray_start;
88347c6ae99SBarry Smith         break;
88447c6ae99SBarry Smith       }
88547c6ae99SBarry Smith     }
88647c6ae99SBarry Smith   } else {
887aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
88847c6ae99SBarry Smith       if (dd->arrayout[i] == *iptr) {
88947c6ae99SBarry Smith         iarray_start    = dd->startout[i];
89047c6ae99SBarry Smith         dd->arrayout[i] = PETSC_NULL;
89147c6ae99SBarry Smith         dd->startout[i] = PETSC_NULL;
89247c6ae99SBarry Smith         break;
89347c6ae99SBarry Smith       }
89447c6ae99SBarry Smith     }
895aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
89647c6ae99SBarry Smith       if (!dd->arrayin[i]){
89747c6ae99SBarry Smith         dd->arrayin[i]  = *iptr;
89847c6ae99SBarry Smith         dd->startin[i]  = iarray_start;
89947c6ae99SBarry Smith         break;
90047c6ae99SBarry Smith       }
90147c6ae99SBarry Smith     }
90247c6ae99SBarry Smith   }
90347c6ae99SBarry Smith   PetscFunctionReturn(0);
90447c6ae99SBarry Smith }
90547c6ae99SBarry Smith 
906