xref: /petsc/src/dm/impls/da/dalocal.c (revision dcca6d9d80ebd869fe6029bd05a3aa9faafef49e)
147c6ae99SBarry Smith 
247c6ae99SBarry Smith /*
347c6ae99SBarry Smith   Code for manipulating distributed regular arrays in parallel.
447c6ae99SBarry Smith */
547c6ae99SBarry Smith 
64035e84dSBarry Smith #include <petsc-private/dmdaimpl.h>    /*I   "petscdmda.h"   I*/
70c312b8eSJed Brown #include <petscsf.h>
847c6ae99SBarry Smith 
947c6ae99SBarry Smith /*
10e3c5b3baSBarry Smith    This allows the DMDA vectors to properly tell MATLAB their dimensions
1147c6ae99SBarry Smith */
1247c6ae99SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
13c6db04a5SJed Brown #include <engine.h>   /* MATLAB include file */
14c6db04a5SJed Brown #include <mex.h>      /* MATLAB include file */
1547c6ae99SBarry Smith #undef __FUNCT__
1647c6ae99SBarry Smith #define __FUNCT__ "VecMatlabEnginePut_DA2d"
171e6b0712SBarry Smith static 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);
28ce94432eSBarry Smith   if (!da) SETERRQ(PetscObjectComm((PetscObject)vec),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 #endif
4547c6ae99SBarry Smith 
4647c6ae99SBarry Smith 
4747c6ae99SBarry Smith #undef __FUNCT__
48564755cdSBarry Smith #define __FUNCT__ "DMCreateLocalVector_DA"
497087cfbeSBarry Smith PetscErrorCode  DMCreateLocalVector_DA(DM da,Vec *g)
5047c6ae99SBarry Smith {
5147c6ae99SBarry Smith   PetscErrorCode ierr;
5247c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
5347c6ae99SBarry Smith 
5447c6ae99SBarry Smith   PetscFunctionBegin;
5547c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
5647c6ae99SBarry Smith   PetscValidPointer(g,2);
5711689aebSJed Brown   if (da->defaultSection) {
5811689aebSJed Brown     ierr = DMCreateLocalVector_Section_Private(da,g);CHKERRQ(ierr);
5911689aebSJed Brown   } else {
6047c6ae99SBarry Smith     ierr = VecCreate(PETSC_COMM_SELF,g);CHKERRQ(ierr);
6147c6ae99SBarry Smith     ierr = VecSetSizes(*g,dd->nlocal,PETSC_DETERMINE);CHKERRQ(ierr);
6247c6ae99SBarry Smith     ierr = VecSetBlockSize(*g,dd->w);CHKERRQ(ierr);
63401ddaa8SBarry Smith     ierr = VecSetType(*g,da->vectype);CHKERRQ(ierr);
64c688c046SMatthew G Knepley     ierr = VecSetDM(*g, da);CHKERRQ(ierr);
6547c6ae99SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
6647c6ae99SBarry Smith     if (dd->w == 1  && dd->dim == 2) {
67bdf89e91SBarry Smith       ierr = PetscObjectComposeFunction((PetscObject)*g,"PetscMatlabEnginePut_C",VecMatlabEnginePut_DA2d);CHKERRQ(ierr);
6847c6ae99SBarry Smith     }
6947c6ae99SBarry Smith #endif
7011689aebSJed Brown   }
7147c6ae99SBarry Smith   PetscFunctionReturn(0);
7247c6ae99SBarry Smith }
7347c6ae99SBarry Smith 
74a66d4d66SMatthew G Knepley #undef __FUNCT__
7557459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetNumCells"
7657459e9aSMatthew G Knepley PetscErrorCode DMDAGetNumCells(DM dm, PetscInt *numCells)
7757459e9aSMatthew G Knepley {
7857459e9aSMatthew G Knepley   DM_DA          *da = (DM_DA*) dm->data;
7957459e9aSMatthew G Knepley   const PetscInt dim = da->dim;
8057459e9aSMatthew G Knepley   const PetscInt mx  = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs;
8157459e9aSMatthew G Knepley   const PetscInt nC  = (mx)*(dim > 1 ? (my)*(dim > 2 ? (mz) : 1) : 1);
8257459e9aSMatthew G Knepley 
8357459e9aSMatthew G Knepley   PetscFunctionBegin;
8457459e9aSMatthew G Knepley   if (numCells) {
8557459e9aSMatthew G Knepley     PetscValidIntPointer(numCells,2);
8657459e9aSMatthew G Knepley     *numCells = nC;
8757459e9aSMatthew G Knepley   }
8857459e9aSMatthew G Knepley   PetscFunctionReturn(0);
8957459e9aSMatthew G Knepley }
9057459e9aSMatthew G Knepley 
9157459e9aSMatthew G Knepley #undef __FUNCT__
9257459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetNumVertices"
9357459e9aSMatthew G Knepley PetscErrorCode DMDAGetNumVertices(DM dm, PetscInt *numVerticesX, PetscInt *numVerticesY, PetscInt *numVerticesZ, PetscInt *numVertices)
9457459e9aSMatthew G Knepley {
9557459e9aSMatthew G Knepley   DM_DA          *da = (DM_DA*) dm->data;
9657459e9aSMatthew G Knepley   const PetscInt dim = da->dim;
9757459e9aSMatthew G Knepley   const PetscInt mx  = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs;
9857459e9aSMatthew G Knepley   const PetscInt nVx = mx+1;
9957459e9aSMatthew G Knepley   const PetscInt nVy = dim > 1 ? (my+1) : 1;
10057459e9aSMatthew G Knepley   const PetscInt nVz = dim > 2 ? (mz+1) : 1;
10157459e9aSMatthew G Knepley   const PetscInt nV  = nVx*nVy*nVz;
10257459e9aSMatthew G Knepley 
10357459e9aSMatthew G Knepley   PetscFunctionBegin;
10457459e9aSMatthew G Knepley   if (numVerticesX) {
10557459e9aSMatthew G Knepley     PetscValidIntPointer(numVerticesX,2);
10657459e9aSMatthew G Knepley     *numVerticesX = nVx;
10757459e9aSMatthew G Knepley   }
10857459e9aSMatthew G Knepley   if (numVerticesY) {
10957459e9aSMatthew G Knepley     PetscValidIntPointer(numVerticesY,3);
11057459e9aSMatthew G Knepley     *numVerticesY = nVy;
11157459e9aSMatthew G Knepley   }
11257459e9aSMatthew G Knepley   if (numVerticesZ) {
11357459e9aSMatthew G Knepley     PetscValidIntPointer(numVerticesZ,4);
11457459e9aSMatthew G Knepley     *numVerticesZ = nVz;
11557459e9aSMatthew G Knepley   }
11657459e9aSMatthew G Knepley   if (numVertices) {
11757459e9aSMatthew G Knepley     PetscValidIntPointer(numVertices,5);
11857459e9aSMatthew G Knepley     *numVertices = nV;
11957459e9aSMatthew G Knepley   }
12057459e9aSMatthew G Knepley   PetscFunctionReturn(0);
12157459e9aSMatthew G Knepley }
12257459e9aSMatthew G Knepley 
12357459e9aSMatthew G Knepley #undef __FUNCT__
12457459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetNumFaces"
12557459e9aSMatthew G Knepley PetscErrorCode DMDAGetNumFaces(DM dm, PetscInt *numXFacesX, PetscInt *numXFaces, PetscInt *numYFacesY, PetscInt *numYFaces, PetscInt *numZFacesZ, PetscInt *numZFaces)
12657459e9aSMatthew G Knepley {
12757459e9aSMatthew G Knepley   DM_DA          *da = (DM_DA*) dm->data;
12857459e9aSMatthew G Knepley   const PetscInt dim = da->dim;
12957459e9aSMatthew G Knepley   const PetscInt mx  = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs;
13057459e9aSMatthew G Knepley   const PetscInt nxF = (dim > 1 ? (my)*(dim > 2 ? (mz) : 1) : 1);
13157459e9aSMatthew G Knepley   const PetscInt nXF = (mx+1)*nxF;
13257459e9aSMatthew G Knepley   const PetscInt nyF = mx*(dim > 2 ? mz : 1);
13357459e9aSMatthew G Knepley   const PetscInt nYF = dim > 1 ? (my+1)*nyF : 0;
13457459e9aSMatthew G Knepley   const PetscInt nzF = mx*(dim > 1 ? my : 0);
13557459e9aSMatthew G Knepley   const PetscInt nZF = dim > 2 ? (mz+1)*nzF : 0;
13657459e9aSMatthew G Knepley 
13757459e9aSMatthew G Knepley   PetscFunctionBegin;
13857459e9aSMatthew G Knepley   if (numXFacesX) {
13957459e9aSMatthew G Knepley     PetscValidIntPointer(numXFacesX,2);
14057459e9aSMatthew G Knepley     *numXFacesX = nxF;
14157459e9aSMatthew G Knepley   }
14257459e9aSMatthew G Knepley   if (numXFaces) {
14357459e9aSMatthew G Knepley     PetscValidIntPointer(numXFaces,3);
14457459e9aSMatthew G Knepley     *numXFaces = nXF;
14557459e9aSMatthew G Knepley   }
14657459e9aSMatthew G Knepley   if (numYFacesY) {
14757459e9aSMatthew G Knepley     PetscValidIntPointer(numYFacesY,4);
14857459e9aSMatthew G Knepley     *numYFacesY = nyF;
14957459e9aSMatthew G Knepley   }
15057459e9aSMatthew G Knepley   if (numYFaces) {
15157459e9aSMatthew G Knepley     PetscValidIntPointer(numYFaces,5);
15257459e9aSMatthew G Knepley     *numYFaces = nYF;
15357459e9aSMatthew G Knepley   }
15457459e9aSMatthew G Knepley   if (numZFacesZ) {
15557459e9aSMatthew G Knepley     PetscValidIntPointer(numZFacesZ,6);
15657459e9aSMatthew G Knepley     *numZFacesZ = nzF;
15757459e9aSMatthew G Knepley   }
15857459e9aSMatthew G Knepley   if (numZFaces) {
15957459e9aSMatthew G Knepley     PetscValidIntPointer(numZFaces,7);
16057459e9aSMatthew G Knepley     *numZFaces = nZF;
16157459e9aSMatthew G Knepley   }
16257459e9aSMatthew G Knepley   PetscFunctionReturn(0);
16357459e9aSMatthew G Knepley }
16457459e9aSMatthew G Knepley 
16557459e9aSMatthew G Knepley #undef __FUNCT__
16657459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetHeightStratum"
16757459e9aSMatthew G Knepley PetscErrorCode DMDAGetHeightStratum(DM dm, PetscInt height, PetscInt *pStart, PetscInt *pEnd)
16857459e9aSMatthew G Knepley {
16957459e9aSMatthew G Knepley   DM_DA          *da = (DM_DA*) dm->data;
17057459e9aSMatthew G Knepley   const PetscInt dim = da->dim;
17157459e9aSMatthew G Knepley   PetscInt       nC, nV, nXF, nYF, nZF;
17257459e9aSMatthew G Knepley   PetscErrorCode ierr;
17357459e9aSMatthew G Knepley 
17457459e9aSMatthew G Knepley   PetscFunctionBegin;
1758865f1eaSKarl Rupp   if (pStart) PetscValidIntPointer(pStart,3);
1768865f1eaSKarl Rupp   if (pEnd)   PetscValidIntPointer(pEnd,4);
17757459e9aSMatthew G Knepley   ierr = DMDAGetNumCells(dm, &nC);CHKERRQ(ierr);
1780298fd71SBarry Smith   ierr = DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV);CHKERRQ(ierr);
1790298fd71SBarry Smith   ierr = DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF);CHKERRQ(ierr);
18057459e9aSMatthew G Knepley   if (height == 0) {
18157459e9aSMatthew G Knepley     /* Cells */
1828865f1eaSKarl Rupp     if (pStart) *pStart = 0;
1838865f1eaSKarl Rupp     if (pEnd)   *pEnd   = nC;
18457459e9aSMatthew G Knepley   } else if (height == 1) {
18557459e9aSMatthew G Knepley     /* Faces */
1868865f1eaSKarl Rupp     if (pStart) *pStart = nC+nV;
1878865f1eaSKarl Rupp     if (pEnd)   *pEnd   = nC+nV+nXF+nYF+nZF;
18857459e9aSMatthew G Knepley   } else if (height == dim) {
18957459e9aSMatthew G Knepley     /* Vertices */
1908865f1eaSKarl Rupp     if (pStart) *pStart = nC;
1918865f1eaSKarl Rupp     if (pEnd)   *pEnd   = nC+nV;
19257459e9aSMatthew G Knepley   } else if (height < 0) {
19357459e9aSMatthew G Knepley     /* All points */
1948865f1eaSKarl Rupp     if (pStart) *pStart = 0;
1958865f1eaSKarl Rupp     if (pEnd)   *pEnd   = nC+nV+nXF+nYF+nZF;
19682f516ccSBarry Smith   } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No points of height %d in the DA", height);
19757459e9aSMatthew G Knepley   PetscFunctionReturn(0);
19857459e9aSMatthew G Knepley }
19957459e9aSMatthew G Knepley 
20057459e9aSMatthew G Knepley #undef __FUNCT__
201a66d4d66SMatthew G Knepley #define __FUNCT__ "DMDACreateSection"
202a66d4d66SMatthew G Knepley /*@C
203a66d4d66SMatthew G Knepley   DMDACreateSection - Create a PetscSection inside the DMDA that describes data layout. This allows multiple fields with
204a66d4d66SMatthew G Knepley   different numbers of dofs on vertices, cells, and faces in each direction.
205a66d4d66SMatthew G Knepley 
206a66d4d66SMatthew G Knepley   Input Parameters:
207a66d4d66SMatthew G Knepley + dm- The DMDA
208a66d4d66SMatthew G Knepley . numFields - The number of fields
2090298fd71SBarry Smith . numComp - The number of components in each field, or NULL for 1
2100298fd71SBarry Smith . numVertexDof - The number of dofs per vertex for each field, or NULL
2110298fd71SBarry Smith . numFaceDof - The number of dofs per face for each field and direction, or NULL
2120298fd71SBarry Smith - numCellDof - The number of dofs per cell for each field, or NULL
213a66d4d66SMatthew G Knepley 
214a66d4d66SMatthew G Knepley   Level: developer
215a66d4d66SMatthew G Knepley 
216a66d4d66SMatthew G Knepley   Note:
217a66d4d66SMatthew G Knepley   The default DMDA numbering is as follows:
218a66d4d66SMatthew G Knepley 
219a66d4d66SMatthew G Knepley     - Cells:    [0,             nC)
220a66d4d66SMatthew G Knepley     - Vertices: [nC,            nC+nV)
22188ed4aceSMatthew G Knepley     - X-Faces:  [nC+nV,         nC+nV+nXF)         normal is +- x-dir
22288ed4aceSMatthew G Knepley     - Y-Faces:  [nC+nV+nXF,     nC+nV+nXF+nYF)     normal is +- y-dir
22388ed4aceSMatthew G Knepley     - Z-Faces:  [nC+nV+nXF+nYF, nC+nV+nXF+nYF+nZF) normal is +- z-dir
224a66d4d66SMatthew G Knepley 
225a66d4d66SMatthew G Knepley   We interpret the default DMDA partition as a cell partition, and the data assignment as a cell assignment.
226a66d4d66SMatthew G Knepley @*/
22780800b1aSMatthew G Knepley PetscErrorCode DMDACreateSection(DM dm, PetscInt numComp[], PetscInt numVertexDof[], PetscInt numFaceDof[], PetscInt numCellDof[])
228a66d4d66SMatthew G Knepley {
229a66d4d66SMatthew G Knepley   DM_DA            *da  = (DM_DA*) dm->data;
230a4b60ecfSMatthew G. Knepley   PetscSection      section;
23188ed4aceSMatthew G Knepley   const PetscInt    dim = da->dim;
23280800b1aSMatthew G Knepley   PetscInt          numFields, numVertexTotDof = 0, numCellTotDof = 0, numFaceTotDof[3] = {0, 0, 0};
23388ed4aceSMatthew G Knepley   PetscSF           sf;
234feafbc5cSMatthew G Knepley   PetscMPIInt       rank;
23588ed4aceSMatthew G Knepley   const PetscMPIInt *neighbors;
23688ed4aceSMatthew G Knepley   PetscInt          *localPoints;
23788ed4aceSMatthew G Knepley   PetscSFNode       *remotePoints;
238f5eeb527SMatthew G Knepley   PetscInt          nleaves = 0,  nleavesCheck = 0, nL = 0;
23957459e9aSMatthew G Knepley   PetscInt          nC, nVx, nVy, nVz, nV, nxF, nXF, nyF, nYF, nzF, nZF;
24057459e9aSMatthew G Knepley   PetscInt          pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart, yfEnd, zfStart, zfEnd;
24188ed4aceSMatthew G Knepley   PetscInt          f, v, c, xf, yf, zf, xn, yn, zn;
242a66d4d66SMatthew G Knepley   PetscErrorCode    ierr;
243a66d4d66SMatthew G Knepley 
244a66d4d66SMatthew G Knepley   PetscFunctionBegin;
245a66d4d66SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
24682f516ccSBarry Smith   ierr    = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
24757459e9aSMatthew G Knepley   ierr    = DMDAGetNumCells(dm, &nC);CHKERRQ(ierr);
24857459e9aSMatthew G Knepley   ierr    = DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV);CHKERRQ(ierr);
24957459e9aSMatthew G Knepley   ierr    = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr);
25057459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, -1,  &pStart, &pEnd);CHKERRQ(ierr);
25157459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, 0,   &cStart, &cEnd);CHKERRQ(ierr);
25257459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, 1,   &fStart, &fEnd);CHKERRQ(ierr);
25357459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr);
25457459e9aSMatthew G Knepley   xfStart = vEnd;  xfEnd = xfStart+nXF;
25557459e9aSMatthew G Knepley   yfStart = xfEnd; yfEnd = yfStart+nYF;
25657459e9aSMatthew G Knepley   zfStart = yfEnd; zfEnd = zfStart+nZF;
25782f516ccSBarry Smith   if (zfEnd != fEnd) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Invalid face end %d, should be %d", zfEnd, fEnd);
25888ed4aceSMatthew G Knepley   /* Create local section */
25980800b1aSMatthew G Knepley   ierr = DMDAGetInfo(dm, 0,0,0,0,0,0,0, &numFields, 0,0,0,0,0);CHKERRQ(ierr);
260a66d4d66SMatthew G Knepley   for (f = 0; f < numFields; ++f) {
2618865f1eaSKarl Rupp     if (numVertexDof) numVertexTotDof  += numVertexDof[f];
2628865f1eaSKarl Rupp     if (numCellDof)   numCellTotDof    += numCellDof[f];
2638865f1eaSKarl Rupp     if (numFaceDof) {
2648865f1eaSKarl Rupp       numFaceTotDof[0] += numFaceDof[f*dim+0];
26588ed4aceSMatthew G Knepley       numFaceTotDof[1] += dim > 1 ? numFaceDof[f*dim+1] : 0;
2660ad7597dSKarl Rupp       numFaceTotDof[2] += dim > 2 ? numFaceDof[f*dim+2] : 0;
2670ad7597dSKarl Rupp     }
268a66d4d66SMatthew G Knepley   }
269a4b60ecfSMatthew G. Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &section);CHKERRQ(ierr);
270a66d4d66SMatthew G Knepley   if (numFields > 1) {
271a4b60ecfSMatthew G. Knepley     ierr = PetscSectionSetNumFields(section, numFields);CHKERRQ(ierr);
272a66d4d66SMatthew G Knepley     for (f = 0; f < numFields; ++f) {
27380800b1aSMatthew G Knepley       const char *name;
27480800b1aSMatthew G Knepley 
27580800b1aSMatthew G Knepley       ierr = DMDAGetFieldName(dm, f, &name);CHKERRQ(ierr);
276a4b60ecfSMatthew G. Knepley       ierr = PetscSectionSetFieldName(section, f, name);CHKERRQ(ierr);
27780800b1aSMatthew G Knepley       if (numComp) {
278a4b60ecfSMatthew G. Knepley         ierr = PetscSectionSetFieldComponents(section, f, numComp[f]);CHKERRQ(ierr);
279a66d4d66SMatthew G Knepley       }
280a66d4d66SMatthew G Knepley     }
281a66d4d66SMatthew G Knepley   } else {
282a66d4d66SMatthew G Knepley     numFields = 0;
283a66d4d66SMatthew G Knepley   }
284a4b60ecfSMatthew G. Knepley   ierr = PetscSectionSetChart(section, pStart, pEnd);CHKERRQ(ierr);
285a66d4d66SMatthew G Knepley   if (numVertexDof) {
286a66d4d66SMatthew G Knepley     for (v = vStart; v < vEnd; ++v) {
287a66d4d66SMatthew G Knepley       for (f = 0; f < numFields; ++f) {
288a4b60ecfSMatthew G. Knepley         ierr = PetscSectionSetFieldDof(section, v, f, numVertexDof[f]);CHKERRQ(ierr);
289a66d4d66SMatthew G Knepley       }
290a4b60ecfSMatthew G. Knepley       ierr = PetscSectionSetDof(section, v, numVertexTotDof);CHKERRQ(ierr);
291a66d4d66SMatthew G Knepley     }
292a66d4d66SMatthew G Knepley   }
293a66d4d66SMatthew G Knepley   if (numFaceDof) {
294a66d4d66SMatthew G Knepley     for (xf = xfStart; xf < xfEnd; ++xf) {
295a66d4d66SMatthew G Knepley       for (f = 0; f < numFields; ++f) {
296a4b60ecfSMatthew G. Knepley         ierr = PetscSectionSetFieldDof(section, xf, f, numFaceDof[f*dim+0]);CHKERRQ(ierr);
297a66d4d66SMatthew G Knepley       }
298a4b60ecfSMatthew G. Knepley       ierr = PetscSectionSetDof(section, xf, numFaceTotDof[0]);CHKERRQ(ierr);
299a66d4d66SMatthew G Knepley     }
300a66d4d66SMatthew G Knepley     for (yf = yfStart; yf < yfEnd; ++yf) {
301a66d4d66SMatthew G Knepley       for (f = 0; f < numFields; ++f) {
302a4b60ecfSMatthew G. Knepley         ierr = PetscSectionSetFieldDof(section, yf, f, numFaceDof[f*dim+1]);CHKERRQ(ierr);
303a66d4d66SMatthew G Knepley       }
304a4b60ecfSMatthew G. Knepley       ierr = PetscSectionSetDof(section, yf, numFaceTotDof[1]);CHKERRQ(ierr);
305a66d4d66SMatthew G Knepley     }
306a66d4d66SMatthew G Knepley     for (zf = zfStart; zf < zfEnd; ++zf) {
307a66d4d66SMatthew G Knepley       for (f = 0; f < numFields; ++f) {
308a4b60ecfSMatthew G. Knepley         ierr = PetscSectionSetFieldDof(section, zf, f, numFaceDof[f*dim+2]);CHKERRQ(ierr);
309a66d4d66SMatthew G Knepley       }
310a4b60ecfSMatthew G. Knepley       ierr = PetscSectionSetDof(section, zf, numFaceTotDof[2]);CHKERRQ(ierr);
311a66d4d66SMatthew G Knepley     }
312a66d4d66SMatthew G Knepley   }
313a66d4d66SMatthew G Knepley   if (numCellDof) {
314a66d4d66SMatthew G Knepley     for (c = cStart; c < cEnd; ++c) {
315a66d4d66SMatthew G Knepley       for (f = 0; f < numFields; ++f) {
316a4b60ecfSMatthew G. Knepley         ierr = PetscSectionSetFieldDof(section, c, f, numCellDof[f]);CHKERRQ(ierr);
317a66d4d66SMatthew G Knepley       }
318a4b60ecfSMatthew G. Knepley       ierr = PetscSectionSetDof(section, c, numCellTotDof);CHKERRQ(ierr);
319a66d4d66SMatthew G Knepley     }
320a66d4d66SMatthew G Knepley   }
321a4b60ecfSMatthew G. Knepley   ierr = PetscSectionSetUp(section);CHKERRQ(ierr);
32288ed4aceSMatthew G Knepley   /* Create mesh point SF */
32388ed4aceSMatthew G Knepley   ierr = DMDAGetNeighbors(dm, &neighbors);CHKERRQ(ierr);
32488ed4aceSMatthew G Knepley   for (zn = 0; zn < (dim > 2 ? 3 : 1); ++zn) {
32588ed4aceSMatthew G Knepley     for (yn = 0; yn < (dim > 1 ? 3 : 1); ++yn) {
32688ed4aceSMatthew G Knepley       for (xn = 0; xn < 3; ++xn) {
3277128ae9fSMatthew G Knepley         const PetscInt xp       = xn-1, yp = dim > 1 ? yn-1 : 0, zp = dim > 2 ? zn-1 : 0;
32888ed4aceSMatthew G Knepley         const PetscInt neighbor = neighbors[(zn*3+yn)*3+xn];
32988ed4aceSMatthew G Knepley 
3303814d604SMatthew G Knepley         if (neighbor >= 0 && neighbor < rank) {
331feafbc5cSMatthew G Knepley           nleaves += (!xp ? nVx : 1) * (!yp ? nVy : 1) * (!zp ? nVz : 1); /* vertices */
33288ed4aceSMatthew G Knepley           if (xp && !yp && !zp) {
33388ed4aceSMatthew G Knepley             nleaves += nxF; /* x faces */
33488ed4aceSMatthew G Knepley           } else if (yp && !zp && !xp) {
33588ed4aceSMatthew G Knepley             nleaves += nyF; /* y faces */
33688ed4aceSMatthew G Knepley           } else if (zp && !xp && !yp) {
33788ed4aceSMatthew G Knepley             nleaves += nzF; /* z faces */
33888ed4aceSMatthew G Knepley           }
33988ed4aceSMatthew G Knepley         }
34088ed4aceSMatthew G Knepley       }
34188ed4aceSMatthew G Knepley     }
34288ed4aceSMatthew G Knepley   }
343*dcca6d9dSJed Brown   ierr = PetscMalloc2(nleaves,&localPoints,nleaves,&remotePoints);CHKERRQ(ierr);
34488ed4aceSMatthew G Knepley   for (zn = 0; zn < (dim > 2 ? 3 : 1); ++zn) {
34588ed4aceSMatthew G Knepley     for (yn = 0; yn < (dim > 1 ? 3 : 1); ++yn) {
34688ed4aceSMatthew G Knepley       for (xn = 0; xn < 3; ++xn) {
3477128ae9fSMatthew G Knepley         const PetscInt xp       = xn-1, yp = dim > 1 ? yn-1 : 0, zp = dim > 2 ? zn-1 : 0;
34888ed4aceSMatthew G Knepley         const PetscInt neighbor = neighbors[(zn*3+yn)*3+xn];
349f5eeb527SMatthew G Knepley         PetscInt       xv, yv, zv;
35088ed4aceSMatthew G Knepley 
3513814d604SMatthew G Knepley         if (neighbor >= 0 && neighbor < rank) {
35288ed4aceSMatthew G Knepley           if (xp < 0) { /* left */
35388ed4aceSMatthew G Knepley             if (yp < 0) { /* bottom */
35488ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
355f5eeb527SMatthew G Knepley                 const PetscInt localVertex  = (      0*nVy +     0)*nVx +     0 + nC;
356f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = ((nVz-1)*nVy + nVy-1)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
357628e3b0dSSatish Balay                 nleavesCheck += 1; /* left bottom back vertex */
358f5eeb527SMatthew G Knepley 
359f5eeb527SMatthew G Knepley                 localPoints[nL]        = localVertex;
360f5eeb527SMatthew G Knepley                 remotePoints[nL].rank  = neighbor;
361f5eeb527SMatthew G Knepley                 remotePoints[nL].index = remoteVertex;
362f5eeb527SMatthew G Knepley                 ++nL;
36388ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
364f5eeb527SMatthew G Knepley                 const PetscInt localVertex  = ((nVz-1)*nVy +     0)*nVx +     0 + nC;
365f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = (      0*nVy + nVy-1)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
366628e3b0dSSatish Balay                 nleavesCheck += 1; /* left bottom front vertex */
367f5eeb527SMatthew G Knepley 
368f5eeb527SMatthew G Knepley                 localPoints[nL]        = localVertex;
369f5eeb527SMatthew G Knepley                 remotePoints[nL].rank  = neighbor;
370f5eeb527SMatthew G Knepley                 remotePoints[nL].index = remoteVertex;
371f5eeb527SMatthew G Knepley                 ++nL;
37288ed4aceSMatthew G Knepley               } else {
37388ed4aceSMatthew G Knepley                 nleavesCheck += nVz; /* left bottom vertices */
374f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv, ++nL) {
375f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = (zv*nVy +     0)*nVx +     0 + nC;
376f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (zv*nVy + nVy-1)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
377f5eeb527SMatthew G Knepley 
378f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
379f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
380f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
381f5eeb527SMatthew G Knepley                 }
38288ed4aceSMatthew G Knepley               }
38388ed4aceSMatthew G Knepley             } else if (yp > 0) { /* top */
38488ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
385f5eeb527SMatthew G Knepley                 const PetscInt localVertex  = (      0*nVy + nVy-1)*nVx +     0 + nC;
386f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = ((nVz-1)*nVy +     0)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
387628e3b0dSSatish Balay                 nleavesCheck += 1; /* left top back vertex */
388f5eeb527SMatthew G Knepley 
389f5eeb527SMatthew G Knepley                 localPoints[nL]        = localVertex;
390f5eeb527SMatthew G Knepley                 remotePoints[nL].rank  = neighbor;
391f5eeb527SMatthew G Knepley                 remotePoints[nL].index = remoteVertex;
392f5eeb527SMatthew G Knepley                 ++nL;
39388ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
394f5eeb527SMatthew G Knepley                 const PetscInt localVertex  = ((nVz-1)*nVy + nVy-1)*nVx +     0 + nC;
395f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = (      0*nVy +     0)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
396628e3b0dSSatish Balay                 nleavesCheck += 1; /* left top front vertex */
397f5eeb527SMatthew G Knepley 
398f5eeb527SMatthew G Knepley                 localPoints[nL]        = localVertex;
399f5eeb527SMatthew G Knepley                 remotePoints[nL].rank  = neighbor;
400f5eeb527SMatthew G Knepley                 remotePoints[nL].index = remoteVertex;
401f5eeb527SMatthew G Knepley                 ++nL;
40288ed4aceSMatthew G Knepley               } else {
40388ed4aceSMatthew G Knepley                 nleavesCheck += nVz; /* left top vertices */
404f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv, ++nL) {
405f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = (zv*nVy + nVy-1)*nVx +     0 + nC;
406f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (zv*nVy +     0)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
407f5eeb527SMatthew G Knepley 
408f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
409f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
410f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
411f5eeb527SMatthew G Knepley                 }
41288ed4aceSMatthew G Knepley               }
41388ed4aceSMatthew G Knepley             } else {
41488ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
41588ed4aceSMatthew G Knepley                 nleavesCheck += nVy; /* left back vertices */
416f5eeb527SMatthew G Knepley                 for (yv = 0; yv < nVy; ++yv, ++nL) {
417f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = (      0*nVy + yv)*nVx +     0 + nC;
418f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = ((nVz-1)*nVy + yv)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
419f5eeb527SMatthew G Knepley 
420f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
421f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
422f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
423f5eeb527SMatthew G Knepley                 }
42488ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
42588ed4aceSMatthew G Knepley                 nleavesCheck += nVy; /* left front vertices */
426f5eeb527SMatthew G Knepley                 for (yv = 0; yv < nVy; ++yv, ++nL) {
427f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = ((nVz-1)*nVy + yv)*nVx +     0 + nC;
428f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (      0*nVy + yv)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
429f5eeb527SMatthew G Knepley 
430f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
431f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
432f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
433f5eeb527SMatthew G Knepley                 }
43488ed4aceSMatthew G Knepley               } else {
43588ed4aceSMatthew G Knepley                 nleavesCheck += nVy*nVz; /* left vertices */
436f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv) {
437f5eeb527SMatthew G Knepley                   for (yv = 0; yv < nVy; ++yv, ++nL) {
438f5eeb527SMatthew G Knepley                     const PetscInt localVertex  = (zv*nVy + yv)*nVx +     0 + nC;
439f5eeb527SMatthew G Knepley                     const PetscInt remoteVertex = (zv*nVy + yv)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
440f5eeb527SMatthew G Knepley 
441f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
442f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
443f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
444f5eeb527SMatthew G Knepley                   }
445f5eeb527SMatthew G Knepley                 }
44688ed4aceSMatthew G Knepley                 nleavesCheck += nxF;     /* left faces */
447f5eeb527SMatthew G Knepley                 for (xf = 0; xf < nxF; ++xf, ++nL) {
448f5eeb527SMatthew G Knepley                   /* THIS IS WRONG */
449f5eeb527SMatthew G Knepley                   const PetscInt localFace  = 0 + nC+nV;
450f5eeb527SMatthew G Knepley                   const PetscInt remoteFace = 0 + nC+nV;
451f5eeb527SMatthew G Knepley 
452f5eeb527SMatthew G Knepley                   localPoints[nL]        = localFace;
453f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
454f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteFace;
455f5eeb527SMatthew G Knepley                 }
45688ed4aceSMatthew G Knepley               }
45788ed4aceSMatthew G Knepley             }
45888ed4aceSMatthew G Knepley           } else if (xp > 0) { /* right */
45988ed4aceSMatthew G Knepley             if (yp < 0) { /* bottom */
46088ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
461f5eeb527SMatthew G Knepley                 const PetscInt localVertex  = (      0*nVy +     0)*nVx + nVx-1 + nC;
462f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = ((nVz-1)*nVy + nVy-1)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
463628e3b0dSSatish Balay                 nleavesCheck += 1; /* right bottom back vertex */
464f5eeb527SMatthew G Knepley 
465f5eeb527SMatthew G Knepley                 localPoints[nL]        = localVertex;
466f5eeb527SMatthew G Knepley                 remotePoints[nL].rank  = neighbor;
467f5eeb527SMatthew G Knepley                 remotePoints[nL].index = remoteVertex;
468f5eeb527SMatthew G Knepley                 ++nL;
46988ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
470f5eeb527SMatthew G Knepley                 const PetscInt localVertex  = ((nVz-1)*nVy +     0)*nVx + nVx-1 + nC;
471f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = (      0*nVy + nVy-1)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
472628e3b0dSSatish Balay                 nleavesCheck += 1; /* right bottom front vertex */
473f5eeb527SMatthew G Knepley 
474f5eeb527SMatthew G Knepley                 localPoints[nL]        = localVertex;
475f5eeb527SMatthew G Knepley                 remotePoints[nL].rank  = neighbor;
476f5eeb527SMatthew G Knepley                 remotePoints[nL].index = remoteVertex;
477f5eeb527SMatthew G Knepley                 ++nL;
47888ed4aceSMatthew G Knepley               } else {
47988ed4aceSMatthew G Knepley                 nleavesCheck += nVz; /* right bottom vertices */
480f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv, ++nL) {
481f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = (zv*nVy +     0)*nVx + nVx-1 + nC;
482f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (zv*nVy + nVy-1)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
483f5eeb527SMatthew G Knepley 
484f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
485f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
486f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
487f5eeb527SMatthew G Knepley                 }
48888ed4aceSMatthew G Knepley               }
48988ed4aceSMatthew G Knepley             } else if (yp > 0) { /* top */
49088ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
491f5eeb527SMatthew G Knepley                 const PetscInt localVertex  = (      0*nVy + nVy-1)*nVx + nVx-1 + nC;
492f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = ((nVz-1)*nVy +     0)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
493628e3b0dSSatish Balay                 nleavesCheck += 1; /* right top back vertex */
494f5eeb527SMatthew G Knepley 
495f5eeb527SMatthew G Knepley                 localPoints[nL]        = localVertex;
496f5eeb527SMatthew G Knepley                 remotePoints[nL].rank  = neighbor;
497f5eeb527SMatthew G Knepley                 remotePoints[nL].index = remoteVertex;
498f5eeb527SMatthew G Knepley                 ++nL;
49988ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
500f5eeb527SMatthew G Knepley                 const PetscInt localVertex  = ((nVz-1)*nVy + nVy-1)*nVx + nVx-1 + nC;
501f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = (      0*nVy +     0)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
502628e3b0dSSatish Balay                 nleavesCheck += 1; /* right top front vertex */
503f5eeb527SMatthew G Knepley 
504f5eeb527SMatthew G Knepley                 localPoints[nL]        = localVertex;
505f5eeb527SMatthew G Knepley                 remotePoints[nL].rank  = neighbor;
506f5eeb527SMatthew G Knepley                 remotePoints[nL].index = remoteVertex;
507f5eeb527SMatthew G Knepley                 ++nL;
50888ed4aceSMatthew G Knepley               } else {
50988ed4aceSMatthew G Knepley                 nleavesCheck += nVz; /* right top vertices */
510f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv, ++nL) {
511f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = (zv*nVy + nVy-1)*nVx + nVx-1 + nC;
512f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (zv*nVy +     0)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
513f5eeb527SMatthew G Knepley 
514f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
515f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
516f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
517f5eeb527SMatthew G Knepley                 }
51888ed4aceSMatthew G Knepley               }
51988ed4aceSMatthew G Knepley             } else {
52088ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
52188ed4aceSMatthew G Knepley                 nleavesCheck += nVy; /* right back vertices */
522f5eeb527SMatthew G Knepley                 for (yv = 0; yv < nVy; ++yv, ++nL) {
523f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = (      0*nVy + yv)*nVx + nVx-1 + nC;
524f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = ((nVz-1)*nVy + yv)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
525f5eeb527SMatthew G Knepley 
526f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
527f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
528f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
529f5eeb527SMatthew G Knepley                 }
53088ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
53188ed4aceSMatthew G Knepley                 nleavesCheck += nVy; /* right front vertices */
532f5eeb527SMatthew G Knepley                 for (yv = 0; yv < nVy; ++yv, ++nL) {
533f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = ((nVz-1)*nVy + yv)*nVx + nVx-1 + nC;
534f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (      0*nVy + yv)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
535f5eeb527SMatthew G Knepley 
536f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
537f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
538f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
539f5eeb527SMatthew G Knepley                 }
54088ed4aceSMatthew G Knepley               } else {
54188ed4aceSMatthew G Knepley                 nleavesCheck += nVy*nVz; /* right vertices */
542f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv) {
543f5eeb527SMatthew G Knepley                   for (yv = 0; yv < nVy; ++yv, ++nL) {
544f5eeb527SMatthew G Knepley                     const PetscInt localVertex  = (zv*nVy + yv)*nVx + nVx-1 + nC;
545f5eeb527SMatthew G Knepley                     const PetscInt remoteVertex = (zv*nVy + yv)*nVx + 0     + nC; /* TODO: Correct this for neighbor sizes */
546f5eeb527SMatthew G Knepley 
547f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
548f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
549f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
550f5eeb527SMatthew G Knepley                   }
551f5eeb527SMatthew G Knepley                 }
55288ed4aceSMatthew G Knepley                 nleavesCheck += nxF;     /* right faces */
553f5eeb527SMatthew G Knepley                 for (xf = 0; xf < nxF; ++xf, ++nL) {
554f5eeb527SMatthew G Knepley                   /* THIS IS WRONG */
555f5eeb527SMatthew G Knepley                   const PetscInt localFace  = 0 + nC+nV;
556f5eeb527SMatthew G Knepley                   const PetscInt remoteFace = 0 + nC+nV;
557f5eeb527SMatthew G Knepley 
558f5eeb527SMatthew G Knepley                   localPoints[nL]        = localFace;
559f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
560f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteFace;
561f5eeb527SMatthew G Knepley                 }
56288ed4aceSMatthew G Knepley               }
56388ed4aceSMatthew G Knepley             }
56488ed4aceSMatthew G Knepley           } else {
56588ed4aceSMatthew G Knepley             if (yp < 0) { /* bottom */
56688ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
56788ed4aceSMatthew G Knepley                 nleavesCheck += nVx; /* bottom back vertices */
568f5eeb527SMatthew G Knepley                 for (xv = 0; xv < nVx; ++xv, ++nL) {
569f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = (      0*nVy +     0)*nVx + xv + nC;
570f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = ((nVz-1)*nVy + nVy-1)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
571f5eeb527SMatthew G Knepley 
572f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
573f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
574f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
575f5eeb527SMatthew G Knepley                 }
57688ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
57788ed4aceSMatthew G Knepley                 nleavesCheck += nVx; /* bottom front vertices */
578f5eeb527SMatthew G Knepley                 for (xv = 0; xv < nVx; ++xv, ++nL) {
579f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = ((nVz-1)*nVy +     0)*nVx + xv + nC;
580f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (      0*nVy + nVy-1)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
581f5eeb527SMatthew G Knepley 
582f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
583f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
584f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
585f5eeb527SMatthew G Knepley                 }
58688ed4aceSMatthew G Knepley               } else {
58788ed4aceSMatthew G Knepley                 nleavesCheck += nVx*nVz; /* bottom vertices */
588f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv) {
589f5eeb527SMatthew G Knepley                   for (xv = 0; xv < nVx; ++xv, ++nL) {
590f5eeb527SMatthew G Knepley                     const PetscInt localVertex  = (zv*nVy +     0)*nVx + xv + nC;
591f5eeb527SMatthew G Knepley                     const PetscInt remoteVertex = (zv*nVy + nVy-1)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
592f5eeb527SMatthew G Knepley 
593f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
594f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
595f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
596f5eeb527SMatthew G Knepley                   }
597f5eeb527SMatthew G Knepley                 }
59888ed4aceSMatthew G Knepley                 nleavesCheck += nyF;     /* bottom faces */
599f5eeb527SMatthew G Knepley                 for (yf = 0; yf < nyF; ++yf, ++nL) {
600f5eeb527SMatthew G Knepley                   /* THIS IS WRONG */
601f5eeb527SMatthew G Knepley                   const PetscInt localFace  = 0 + nC+nV;
602f5eeb527SMatthew G Knepley                   const PetscInt remoteFace = 0 + nC+nV;
603f5eeb527SMatthew G Knepley 
604f5eeb527SMatthew G Knepley                   localPoints[nL]        = localFace;
605f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
606f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteFace;
607f5eeb527SMatthew G Knepley                 }
60888ed4aceSMatthew G Knepley               }
60988ed4aceSMatthew G Knepley             } else if (yp > 0) { /* top */
61088ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
61188ed4aceSMatthew G Knepley                 nleavesCheck += nVx; /* top back vertices */
612f5eeb527SMatthew G Knepley                 for (xv = 0; xv < nVx; ++xv, ++nL) {
613f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = (      0*nVy + nVy-1)*nVx + xv + nC;
614f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = ((nVz-1)*nVy +     0)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
615f5eeb527SMatthew G Knepley 
616f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
617f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
618f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
619f5eeb527SMatthew G Knepley                 }
62088ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
62188ed4aceSMatthew G Knepley                 nleavesCheck += nVx; /* top front vertices */
622f5eeb527SMatthew G Knepley                 for (xv = 0; xv < nVx; ++xv, ++nL) {
623f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = ((nVz-1)*nVy + nVy-1)*nVx + xv + nC;
624f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (      0*nVy +     0)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
625f5eeb527SMatthew G Knepley 
626f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
627f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
628f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
629f5eeb527SMatthew G Knepley                 }
63088ed4aceSMatthew G Knepley               } else {
63188ed4aceSMatthew G Knepley                 nleavesCheck += nVx*nVz; /* top vertices */
632f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv) {
633f5eeb527SMatthew G Knepley                   for (xv = 0; xv < nVx; ++xv, ++nL) {
634f5eeb527SMatthew G Knepley                     const PetscInt localVertex  = (zv*nVy + nVy-1)*nVx + xv + nC;
635f5eeb527SMatthew G Knepley                     const PetscInt remoteVertex = (zv*nVy +     0)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
636f5eeb527SMatthew G Knepley 
637f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
638f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
639f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
640f5eeb527SMatthew G Knepley                   }
641f5eeb527SMatthew G Knepley                 }
64288ed4aceSMatthew G Knepley                 nleavesCheck += nyF;     /* top faces */
643f5eeb527SMatthew G Knepley                 for (yf = 0; yf < nyF; ++yf, ++nL) {
644f5eeb527SMatthew G Knepley                   /* THIS IS WRONG */
645f5eeb527SMatthew G Knepley                   const PetscInt localFace  = 0 + nC+nV;
646f5eeb527SMatthew G Knepley                   const PetscInt remoteFace = 0 + nC+nV;
647f5eeb527SMatthew G Knepley 
648f5eeb527SMatthew G Knepley                   localPoints[nL]        = localFace;
649f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
650f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteFace;
651f5eeb527SMatthew G Knepley                 }
65288ed4aceSMatthew G Knepley               }
65388ed4aceSMatthew G Knepley             } else {
65488ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
65588ed4aceSMatthew G Knepley                 nleavesCheck += nVx*nVy; /* back vertices */
656f5eeb527SMatthew G Knepley                 for (yv = 0; yv < nVy; ++yv) {
657f5eeb527SMatthew G Knepley                   for (xv = 0; xv < nVx; ++xv, ++nL) {
658f5eeb527SMatthew G Knepley                     const PetscInt localVertex  = (      0*nVy + yv)*nVx + xv + nC;
659f5eeb527SMatthew G Knepley                     const PetscInt remoteVertex = ((nVz-1)*nVy + yv)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
660f5eeb527SMatthew G Knepley 
661f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
662f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
663f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
664f5eeb527SMatthew G Knepley                   }
665f5eeb527SMatthew G Knepley                 }
66688ed4aceSMatthew G Knepley                 nleavesCheck += nzF;     /* back faces */
667f5eeb527SMatthew G Knepley                 for (zf = 0; zf < nzF; ++zf, ++nL) {
668f5eeb527SMatthew G Knepley                   /* THIS IS WRONG */
669f5eeb527SMatthew G Knepley                   const PetscInt localFace  = 0 + nC+nV;
670f5eeb527SMatthew G Knepley                   const PetscInt remoteFace = 0 + nC+nV;
671f5eeb527SMatthew G Knepley 
672f5eeb527SMatthew G Knepley                   localPoints[nL]        = localFace;
673f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
674f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteFace;
675f5eeb527SMatthew G Knepley                 }
67688ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
67788ed4aceSMatthew G Knepley                 nleavesCheck += nVx*nVy; /* front vertices */
678f5eeb527SMatthew G Knepley                 for (yv = 0; yv < nVy; ++yv) {
679f5eeb527SMatthew G Knepley                   for (xv = 0; xv < nVx; ++xv, ++nL) {
680f5eeb527SMatthew G Knepley                     const PetscInt localVertex  = ((nVz-1)*nVy + yv)*nVx + xv + nC;
681f5eeb527SMatthew G Knepley                     const PetscInt remoteVertex = (      0*nVy + yv)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
682f5eeb527SMatthew G Knepley 
683f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
684f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
685f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
686f5eeb527SMatthew G Knepley                   }
687f5eeb527SMatthew G Knepley                 }
68888ed4aceSMatthew G Knepley                 nleavesCheck += nzF;     /* front faces */
689f5eeb527SMatthew G Knepley                 for (zf = 0; zf < nzF; ++zf, ++nL) {
690f5eeb527SMatthew G Knepley                   /* THIS IS WRONG */
691f5eeb527SMatthew G Knepley                   const PetscInt localFace  = 0 + nC+nV;
692f5eeb527SMatthew G Knepley                   const PetscInt remoteFace = 0 + nC+nV;
693f5eeb527SMatthew G Knepley 
694f5eeb527SMatthew G Knepley                   localPoints[nL]        = localFace;
695f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
696f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteFace;
697f5eeb527SMatthew G Knepley                 }
69888ed4aceSMatthew G Knepley               } else {
69988ed4aceSMatthew G Knepley                 /* Nothing is shared from the interior */
70088ed4aceSMatthew G Knepley               }
70188ed4aceSMatthew G Knepley             }
70288ed4aceSMatthew G Knepley           }
70388ed4aceSMatthew G Knepley         }
70488ed4aceSMatthew G Knepley       }
70588ed4aceSMatthew G Knepley     }
70688ed4aceSMatthew G Knepley   }
7073814d604SMatthew G Knepley   /* TODO: Remove duplication in leaf determination */
70882f516ccSBarry Smith   if (nleaves != nleavesCheck) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "The number of leaves %d did not match the number of remote leaves %d", nleaves, nleavesCheck);
70982f516ccSBarry Smith   ierr = PetscSFCreate(PetscObjectComm((PetscObject)dm), &sf);CHKERRQ(ierr);
71088ed4aceSMatthew G Knepley   ierr = PetscSFSetGraph(sf, pEnd, nleaves, localPoints, PETSC_OWN_POINTER, remotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr);
711a4b60ecfSMatthew G. Knepley   ierr = DMSetPointSF(dm, sf);CHKERRQ(ierr);
71288ed4aceSMatthew G Knepley   ierr = PetscSFDestroy(&sf);CHKERRQ(ierr);
713a4b60ecfSMatthew G. Knepley   ierr = DMSetDefaultSection(dm, section);CHKERRQ(ierr);
714a4b60ecfSMatthew G. Knepley   ierr = PetscSectionDestroy(&section);CHKERRQ(ierr);
715a66d4d66SMatthew G Knepley   PetscFunctionReturn(0);
716a66d4d66SMatthew G Knepley }
717a66d4d66SMatthew G Knepley 
71847c6ae99SBarry Smith /* ------------------------------------------------------------------- */
71947c6ae99SBarry Smith 
72047c6ae99SBarry Smith #undef __FUNCT__
721aa219208SBarry Smith #define __FUNCT__ "DMDAGetArray"
72247c6ae99SBarry Smith /*@C
723aa219208SBarry Smith      DMDAGetArray - Gets a work array for a DMDA
72447c6ae99SBarry Smith 
72547c6ae99SBarry Smith     Input Parameter:
72647c6ae99SBarry Smith +    da - information about my local patch
72747c6ae99SBarry Smith -    ghosted - do you want arrays for the ghosted or nonghosted patch
72847c6ae99SBarry Smith 
72947c6ae99SBarry Smith     Output Parameters:
73047c6ae99SBarry Smith .    vptr - array data structured
73147c6ae99SBarry Smith 
73247c6ae99SBarry Smith     Note:  The vector values are NOT initialized and may have garbage in them, so you may need
73347c6ae99SBarry Smith            to zero them.
73447c6ae99SBarry Smith 
73547c6ae99SBarry Smith   Level: advanced
73647c6ae99SBarry Smith 
737bcaeba4dSBarry Smith .seealso: DMDARestoreArray()
73847c6ae99SBarry Smith 
73947c6ae99SBarry Smith @*/
7407087cfbeSBarry Smith PetscErrorCode  DMDAGetArray(DM da,PetscBool ghosted,void *vptr)
74147c6ae99SBarry Smith {
74247c6ae99SBarry Smith   PetscErrorCode ierr;
74347c6ae99SBarry Smith   PetscInt       j,i,xs,ys,xm,ym,zs,zm;
74447c6ae99SBarry Smith   char           *iarray_start;
74547c6ae99SBarry Smith   void           **iptr = (void**)vptr;
74647c6ae99SBarry Smith   DM_DA          *dd    = (DM_DA*)da->data;
74747c6ae99SBarry Smith 
74847c6ae99SBarry Smith   PetscFunctionBegin;
74947c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
75047c6ae99SBarry Smith   if (ghosted) {
751aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
75247c6ae99SBarry Smith       if (dd->arrayghostedin[i]) {
75347c6ae99SBarry Smith         *iptr                 = dd->arrayghostedin[i];
75447c6ae99SBarry Smith         iarray_start          = (char*)dd->startghostedin[i];
7550298fd71SBarry Smith         dd->arrayghostedin[i] = NULL;
7560298fd71SBarry Smith         dd->startghostedin[i] = NULL;
75747c6ae99SBarry Smith 
75847c6ae99SBarry Smith         goto done;
75947c6ae99SBarry Smith       }
76047c6ae99SBarry Smith     }
76147c6ae99SBarry Smith     xs = dd->Xs;
76247c6ae99SBarry Smith     ys = dd->Ys;
76347c6ae99SBarry Smith     zs = dd->Zs;
76447c6ae99SBarry Smith     xm = dd->Xe-dd->Xs;
76547c6ae99SBarry Smith     ym = dd->Ye-dd->Ys;
76647c6ae99SBarry Smith     zm = dd->Ze-dd->Zs;
76747c6ae99SBarry Smith   } else {
768aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
76947c6ae99SBarry Smith       if (dd->arrayin[i]) {
77047c6ae99SBarry Smith         *iptr          = dd->arrayin[i];
77147c6ae99SBarry Smith         iarray_start   = (char*)dd->startin[i];
7720298fd71SBarry Smith         dd->arrayin[i] = NULL;
7730298fd71SBarry Smith         dd->startin[i] = NULL;
77447c6ae99SBarry Smith 
77547c6ae99SBarry Smith         goto done;
77647c6ae99SBarry Smith       }
77747c6ae99SBarry Smith     }
77847c6ae99SBarry Smith     xs = dd->xs;
77947c6ae99SBarry Smith     ys = dd->ys;
78047c6ae99SBarry Smith     zs = dd->zs;
78147c6ae99SBarry Smith     xm = dd->xe-dd->xs;
78247c6ae99SBarry Smith     ym = dd->ye-dd->ys;
78347c6ae99SBarry Smith     zm = dd->ze-dd->zs;
78447c6ae99SBarry Smith   }
78547c6ae99SBarry Smith 
78647c6ae99SBarry Smith   switch (dd->dim) {
78747c6ae99SBarry Smith   case 1: {
78847c6ae99SBarry Smith     void *ptr;
78947c6ae99SBarry Smith 
79047c6ae99SBarry Smith     ierr = PetscMalloc(xm*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
79147c6ae99SBarry Smith 
79247c6ae99SBarry Smith     ptr   = (void*)(iarray_start - xs*sizeof(PetscScalar));
79347c6ae99SBarry Smith     *iptr = (void*)ptr;
7948865f1eaSKarl Rupp     break;
7958865f1eaSKarl Rupp   }
79647c6ae99SBarry Smith   case 2: {
79747c6ae99SBarry Smith     void **ptr;
79847c6ae99SBarry Smith 
79947c6ae99SBarry Smith     ierr = PetscMalloc((ym+1)*sizeof(void*)+xm*ym*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
80047c6ae99SBarry Smith 
80147c6ae99SBarry Smith     ptr = (void**)(iarray_start + xm*ym*sizeof(PetscScalar) - ys*sizeof(void*));
8028865f1eaSKarl Rupp     for (j=ys; j<ys+ym; j++) ptr[j] = iarray_start + sizeof(PetscScalar)*(xm*(j-ys) - xs);
80347c6ae99SBarry Smith     *iptr = (void*)ptr;
8048865f1eaSKarl Rupp     break;
8058865f1eaSKarl Rupp   }
80647c6ae99SBarry Smith   case 3: {
80747c6ae99SBarry Smith     void ***ptr,**bptr;
80847c6ae99SBarry Smith 
80947c6ae99SBarry Smith     ierr = PetscMalloc((zm+1)*sizeof(void**)+(ym*zm+1)*sizeof(void*)+xm*ym*zm*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
81047c6ae99SBarry Smith 
81147c6ae99SBarry Smith     ptr  = (void***)(iarray_start + xm*ym*zm*sizeof(PetscScalar) - zs*sizeof(void*));
81247c6ae99SBarry Smith     bptr = (void**)(iarray_start + xm*ym*zm*sizeof(PetscScalar) + zm*sizeof(void**));
8138865f1eaSKarl Rupp     for (i=zs; i<zs+zm; i++) ptr[i] = bptr + ((i-zs)*ym - ys);
81447c6ae99SBarry Smith     for (i=zs; i<zs+zm; i++) {
81547c6ae99SBarry Smith       for (j=ys; j<ys+ym; j++) {
81647c6ae99SBarry Smith         ptr[i][j] = iarray_start + sizeof(PetscScalar)*(xm*ym*(i-zs) + xm*(j-ys) - xs);
81747c6ae99SBarry Smith       }
81847c6ae99SBarry Smith     }
81947c6ae99SBarry Smith     *iptr = (void*)ptr;
8208865f1eaSKarl Rupp     break;
8218865f1eaSKarl Rupp   }
82247c6ae99SBarry Smith   default:
823ce94432eSBarry Smith     SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Dimension %D not supported",dd->dim);
82447c6ae99SBarry Smith   }
82547c6ae99SBarry Smith 
82647c6ae99SBarry Smith done:
82747c6ae99SBarry Smith   /* add arrays to the checked out list */
82847c6ae99SBarry Smith   if (ghosted) {
829aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
83047c6ae99SBarry Smith       if (!dd->arrayghostedout[i]) {
83147c6ae99SBarry Smith         dd->arrayghostedout[i] = *iptr;
83247c6ae99SBarry Smith         dd->startghostedout[i] = iarray_start;
83347c6ae99SBarry Smith         break;
83447c6ae99SBarry Smith       }
83547c6ae99SBarry Smith     }
83647c6ae99SBarry Smith   } else {
837aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
83847c6ae99SBarry Smith       if (!dd->arrayout[i]) {
83947c6ae99SBarry Smith         dd->arrayout[i] = *iptr;
84047c6ae99SBarry Smith         dd->startout[i] = iarray_start;
84147c6ae99SBarry Smith         break;
84247c6ae99SBarry Smith       }
84347c6ae99SBarry Smith     }
84447c6ae99SBarry Smith   }
84547c6ae99SBarry Smith   PetscFunctionReturn(0);
84647c6ae99SBarry Smith }
84747c6ae99SBarry Smith 
84847c6ae99SBarry Smith #undef __FUNCT__
849aa219208SBarry Smith #define __FUNCT__ "DMDARestoreArray"
85047c6ae99SBarry Smith /*@C
851aa219208SBarry Smith      DMDARestoreArray - Restores an array of derivative types for a DMDA
85247c6ae99SBarry Smith 
85347c6ae99SBarry Smith     Input Parameter:
85447c6ae99SBarry Smith +    da - information about my local patch
85547c6ae99SBarry Smith .    ghosted - do you want arrays for the ghosted or nonghosted patch
85647c6ae99SBarry Smith -    vptr - array data structured to be passed to ad_FormFunctionLocal()
85747c6ae99SBarry Smith 
85847c6ae99SBarry Smith      Level: advanced
85947c6ae99SBarry Smith 
860bcaeba4dSBarry Smith .seealso: DMDAGetArray()
86147c6ae99SBarry Smith 
86247c6ae99SBarry Smith @*/
8637087cfbeSBarry Smith PetscErrorCode  DMDARestoreArray(DM da,PetscBool ghosted,void *vptr)
86447c6ae99SBarry Smith {
86547c6ae99SBarry Smith   PetscInt i;
86647c6ae99SBarry Smith   void     **iptr = (void**)vptr,*iarray_start = 0;
86747c6ae99SBarry Smith   DM_DA    *dd    = (DM_DA*)da->data;
86847c6ae99SBarry Smith 
86947c6ae99SBarry Smith   PetscFunctionBegin;
87047c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
87147c6ae99SBarry Smith   if (ghosted) {
872aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
87347c6ae99SBarry Smith       if (dd->arrayghostedout[i] == *iptr) {
87447c6ae99SBarry Smith         iarray_start           = dd->startghostedout[i];
8750298fd71SBarry Smith         dd->arrayghostedout[i] = NULL;
8760298fd71SBarry Smith         dd->startghostedout[i] = NULL;
87747c6ae99SBarry Smith         break;
87847c6ae99SBarry Smith       }
87947c6ae99SBarry Smith     }
880aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
88147c6ae99SBarry Smith       if (!dd->arrayghostedin[i]) {
88247c6ae99SBarry Smith         dd->arrayghostedin[i] = *iptr;
88347c6ae99SBarry Smith         dd->startghostedin[i] = iarray_start;
88447c6ae99SBarry Smith         break;
88547c6ae99SBarry Smith       }
88647c6ae99SBarry Smith     }
88747c6ae99SBarry Smith   } else {
888aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
88947c6ae99SBarry Smith       if (dd->arrayout[i] == *iptr) {
89047c6ae99SBarry Smith         iarray_start    = dd->startout[i];
8910298fd71SBarry Smith         dd->arrayout[i] = NULL;
8920298fd71SBarry Smith         dd->startout[i] = NULL;
89347c6ae99SBarry Smith         break;
89447c6ae99SBarry Smith       }
89547c6ae99SBarry Smith     }
896aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
89747c6ae99SBarry Smith       if (!dd->arrayin[i]) {
89847c6ae99SBarry Smith         dd->arrayin[i] = *iptr;
89947c6ae99SBarry Smith         dd->startin[i] = iarray_start;
90047c6ae99SBarry Smith         break;
90147c6ae99SBarry Smith       }
90247c6ae99SBarry Smith     }
90347c6ae99SBarry Smith   }
90447c6ae99SBarry Smith   PetscFunctionReturn(0);
90547c6ae99SBarry Smith }
90647c6ae99SBarry Smith 
907