xref: /petsc/src/dm/impls/da/dalocal.c (revision df66330e08ea809e16b0b335d346f2d47352060b)
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*/
7b2fff234SMatthew G. Knepley #include <petscbt.h>
80c312b8eSJed Brown #include <petscsf.h>
937b16e26SMatthew G. Knepley #include <petscfe.h>
1047c6ae99SBarry Smith 
1147c6ae99SBarry Smith /*
12e3c5b3baSBarry Smith    This allows the DMDA vectors to properly tell MATLAB their dimensions
1347c6ae99SBarry Smith */
1447c6ae99SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
15c6db04a5SJed Brown #include <engine.h>   /* MATLAB include file */
16c6db04a5SJed Brown #include <mex.h>      /* MATLAB include file */
1747c6ae99SBarry Smith #undef __FUNCT__
1847c6ae99SBarry Smith #define __FUNCT__ "VecMatlabEnginePut_DA2d"
191e6b0712SBarry Smith static PetscErrorCode  VecMatlabEnginePut_DA2d(PetscObject obj,void *mengine)
2047c6ae99SBarry Smith {
2147c6ae99SBarry Smith   PetscErrorCode ierr;
2247c6ae99SBarry Smith   PetscInt       n,m;
2347c6ae99SBarry Smith   Vec            vec = (Vec)obj;
2447c6ae99SBarry Smith   PetscScalar    *array;
2547c6ae99SBarry Smith   mxArray        *mat;
269a42bb27SBarry Smith   DM             da;
2747c6ae99SBarry Smith 
2847c6ae99SBarry Smith   PetscFunctionBegin;
29c688c046SMatthew G Knepley   ierr = VecGetDM(vec, &da);CHKERRQ(ierr);
30ce94432eSBarry Smith   if (!da) SETERRQ(PetscObjectComm((PetscObject)vec),PETSC_ERR_ARG_WRONGSTATE,"Vector not associated with a DMDA");
31aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,0,0,0,&m,&n,0);CHKERRQ(ierr);
3247c6ae99SBarry Smith 
3347c6ae99SBarry Smith   ierr = VecGetArray(vec,&array);CHKERRQ(ierr);
3447c6ae99SBarry Smith #if !defined(PETSC_USE_COMPLEX)
3547c6ae99SBarry Smith   mat = mxCreateDoubleMatrix(m,n,mxREAL);
3647c6ae99SBarry Smith #else
3747c6ae99SBarry Smith   mat = mxCreateDoubleMatrix(m,n,mxCOMPLEX);
3847c6ae99SBarry Smith #endif
3947c6ae99SBarry Smith   ierr = PetscMemcpy(mxGetPr(mat),array,n*m*sizeof(PetscScalar));CHKERRQ(ierr);
4047c6ae99SBarry Smith   ierr = PetscObjectName(obj);CHKERRQ(ierr);
4147c6ae99SBarry Smith   engPutVariable((Engine*)mengine,obj->name,mat);
4247c6ae99SBarry Smith 
4347c6ae99SBarry Smith   ierr = VecRestoreArray(vec,&array);CHKERRQ(ierr);
4447c6ae99SBarry Smith   PetscFunctionReturn(0);
4547c6ae99SBarry Smith }
4647c6ae99SBarry Smith #endif
4747c6ae99SBarry Smith 
4847c6ae99SBarry Smith 
4947c6ae99SBarry Smith #undef __FUNCT__
50564755cdSBarry Smith #define __FUNCT__ "DMCreateLocalVector_DA"
517087cfbeSBarry Smith PetscErrorCode  DMCreateLocalVector_DA(DM da,Vec *g)
5247c6ae99SBarry Smith {
5347c6ae99SBarry Smith   PetscErrorCode ierr;
5447c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
5547c6ae99SBarry Smith 
5647c6ae99SBarry Smith   PetscFunctionBegin;
5747c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
5847c6ae99SBarry Smith   PetscValidPointer(g,2);
5911689aebSJed Brown   if (da->defaultSection) {
6011689aebSJed Brown     ierr = DMCreateLocalVector_Section_Private(da,g);CHKERRQ(ierr);
6111689aebSJed Brown   } else {
6247c6ae99SBarry Smith     ierr = VecCreate(PETSC_COMM_SELF,g);CHKERRQ(ierr);
6347c6ae99SBarry Smith     ierr = VecSetSizes(*g,dd->nlocal,PETSC_DETERMINE);CHKERRQ(ierr);
6447c6ae99SBarry Smith     ierr = VecSetBlockSize(*g,dd->w);CHKERRQ(ierr);
65401ddaa8SBarry Smith     ierr = VecSetType(*g,da->vectype);CHKERRQ(ierr);
66c688c046SMatthew G Knepley     ierr = VecSetDM(*g, da);CHKERRQ(ierr);
6747c6ae99SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
6847c6ae99SBarry Smith     if (dd->w == 1  && dd->dim == 2) {
69bdf89e91SBarry Smith       ierr = PetscObjectComposeFunction((PetscObject)*g,"PetscMatlabEnginePut_C",VecMatlabEnginePut_DA2d);CHKERRQ(ierr);
7047c6ae99SBarry Smith     }
7147c6ae99SBarry Smith #endif
7211689aebSJed Brown   }
7347c6ae99SBarry Smith   PetscFunctionReturn(0);
7447c6ae99SBarry Smith }
7547c6ae99SBarry Smith 
76a66d4d66SMatthew G Knepley #undef __FUNCT__
7757459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetNumCells"
78e42e3c58SMatthew G. Knepley PetscErrorCode DMDAGetNumCells(DM dm, PetscInt *numCellsX, PetscInt *numCellsY, PetscInt *numCellsZ, PetscInt *numCells)
7957459e9aSMatthew G Knepley {
8057459e9aSMatthew G Knepley   DM_DA          *da = (DM_DA*) dm->data;
8157459e9aSMatthew G Knepley   const PetscInt dim = da->dim;
8257459e9aSMatthew G Knepley   const PetscInt mx  = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs;
8357459e9aSMatthew G Knepley   const PetscInt nC  = (mx)*(dim > 1 ? (my)*(dim > 2 ? (mz) : 1) : 1);
8457459e9aSMatthew G Knepley 
8557459e9aSMatthew G Knepley   PetscFunctionBegin;
86e42e3c58SMatthew G. Knepley   if (numCellsX) {
87e42e3c58SMatthew G. Knepley     PetscValidIntPointer(numCellsX,2);
88e42e3c58SMatthew G. Knepley     *numCellsX = mx;
89e42e3c58SMatthew G. Knepley   }
90e42e3c58SMatthew G. Knepley   if (numCellsY) {
91e42e3c58SMatthew G. Knepley     PetscValidIntPointer(numCellsX,3);
92e42e3c58SMatthew G. Knepley     *numCellsY = my;
93e42e3c58SMatthew G. Knepley   }
94e42e3c58SMatthew G. Knepley   if (numCellsZ) {
95e42e3c58SMatthew G. Knepley     PetscValidIntPointer(numCellsX,4);
96e42e3c58SMatthew G. Knepley     *numCellsZ = mz;
97e42e3c58SMatthew G. Knepley   }
9857459e9aSMatthew G Knepley   if (numCells) {
99e42e3c58SMatthew G. Knepley     PetscValidIntPointer(numCells,5);
10057459e9aSMatthew G Knepley     *numCells = nC;
10157459e9aSMatthew G Knepley   }
10257459e9aSMatthew G Knepley   PetscFunctionReturn(0);
10357459e9aSMatthew G Knepley }
10457459e9aSMatthew G Knepley 
10557459e9aSMatthew G Knepley #undef __FUNCT__
10657459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetNumVertices"
10757459e9aSMatthew G Knepley PetscErrorCode DMDAGetNumVertices(DM dm, PetscInt *numVerticesX, PetscInt *numVerticesY, PetscInt *numVerticesZ, PetscInt *numVertices)
10857459e9aSMatthew G Knepley {
10957459e9aSMatthew G Knepley   DM_DA          *da = (DM_DA*) dm->data;
11057459e9aSMatthew G Knepley   const PetscInt dim = da->dim;
11157459e9aSMatthew G Knepley   const PetscInt mx  = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs;
11257459e9aSMatthew G Knepley   const PetscInt nVx = mx+1;
11357459e9aSMatthew G Knepley   const PetscInt nVy = dim > 1 ? (my+1) : 1;
11457459e9aSMatthew G Knepley   const PetscInt nVz = dim > 2 ? (mz+1) : 1;
11557459e9aSMatthew G Knepley   const PetscInt nV  = nVx*nVy*nVz;
11657459e9aSMatthew G Knepley 
11757459e9aSMatthew G Knepley   PetscFunctionBegin;
11857459e9aSMatthew G Knepley   if (numVerticesX) {
11957459e9aSMatthew G Knepley     PetscValidIntPointer(numVerticesX,2);
12057459e9aSMatthew G Knepley     *numVerticesX = nVx;
12157459e9aSMatthew G Knepley   }
12257459e9aSMatthew G Knepley   if (numVerticesY) {
12357459e9aSMatthew G Knepley     PetscValidIntPointer(numVerticesY,3);
12457459e9aSMatthew G Knepley     *numVerticesY = nVy;
12557459e9aSMatthew G Knepley   }
12657459e9aSMatthew G Knepley   if (numVerticesZ) {
12757459e9aSMatthew G Knepley     PetscValidIntPointer(numVerticesZ,4);
12857459e9aSMatthew G Knepley     *numVerticesZ = nVz;
12957459e9aSMatthew G Knepley   }
13057459e9aSMatthew G Knepley   if (numVertices) {
13157459e9aSMatthew G Knepley     PetscValidIntPointer(numVertices,5);
13257459e9aSMatthew G Knepley     *numVertices = nV;
13357459e9aSMatthew G Knepley   }
13457459e9aSMatthew G Knepley   PetscFunctionReturn(0);
13557459e9aSMatthew G Knepley }
13657459e9aSMatthew G Knepley 
13757459e9aSMatthew G Knepley #undef __FUNCT__
13857459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetNumFaces"
13957459e9aSMatthew G Knepley PetscErrorCode DMDAGetNumFaces(DM dm, PetscInt *numXFacesX, PetscInt *numXFaces, PetscInt *numYFacesY, PetscInt *numYFaces, PetscInt *numZFacesZ, PetscInt *numZFaces)
14057459e9aSMatthew G Knepley {
14157459e9aSMatthew G Knepley   DM_DA          *da = (DM_DA*) dm->data;
14257459e9aSMatthew G Knepley   const PetscInt dim = da->dim;
14357459e9aSMatthew G Knepley   const PetscInt mx  = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs;
14457459e9aSMatthew G Knepley   const PetscInt nxF = (dim > 1 ? (my)*(dim > 2 ? (mz) : 1) : 1);
14557459e9aSMatthew G Knepley   const PetscInt nXF = (mx+1)*nxF;
14657459e9aSMatthew G Knepley   const PetscInt nyF = mx*(dim > 2 ? mz : 1);
14757459e9aSMatthew G Knepley   const PetscInt nYF = dim > 1 ? (my+1)*nyF : 0;
14857459e9aSMatthew G Knepley   const PetscInt nzF = mx*(dim > 1 ? my : 0);
14957459e9aSMatthew G Knepley   const PetscInt nZF = dim > 2 ? (mz+1)*nzF : 0;
15057459e9aSMatthew G Knepley 
15157459e9aSMatthew G Knepley   PetscFunctionBegin;
15257459e9aSMatthew G Knepley   if (numXFacesX) {
15357459e9aSMatthew G Knepley     PetscValidIntPointer(numXFacesX,2);
15457459e9aSMatthew G Knepley     *numXFacesX = nxF;
15557459e9aSMatthew G Knepley   }
15657459e9aSMatthew G Knepley   if (numXFaces) {
15757459e9aSMatthew G Knepley     PetscValidIntPointer(numXFaces,3);
15857459e9aSMatthew G Knepley     *numXFaces = nXF;
15957459e9aSMatthew G Knepley   }
16057459e9aSMatthew G Knepley   if (numYFacesY) {
16157459e9aSMatthew G Knepley     PetscValidIntPointer(numYFacesY,4);
16257459e9aSMatthew G Knepley     *numYFacesY = nyF;
16357459e9aSMatthew G Knepley   }
16457459e9aSMatthew G Knepley   if (numYFaces) {
16557459e9aSMatthew G Knepley     PetscValidIntPointer(numYFaces,5);
16657459e9aSMatthew G Knepley     *numYFaces = nYF;
16757459e9aSMatthew G Knepley   }
16857459e9aSMatthew G Knepley   if (numZFacesZ) {
16957459e9aSMatthew G Knepley     PetscValidIntPointer(numZFacesZ,6);
17057459e9aSMatthew G Knepley     *numZFacesZ = nzF;
17157459e9aSMatthew G Knepley   }
17257459e9aSMatthew G Knepley   if (numZFaces) {
17357459e9aSMatthew G Knepley     PetscValidIntPointer(numZFaces,7);
17457459e9aSMatthew G Knepley     *numZFaces = nZF;
17557459e9aSMatthew G Knepley   }
17657459e9aSMatthew G Knepley   PetscFunctionReturn(0);
17757459e9aSMatthew G Knepley }
17857459e9aSMatthew G Knepley 
17957459e9aSMatthew G Knepley #undef __FUNCT__
18057459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetHeightStratum"
18157459e9aSMatthew G Knepley PetscErrorCode DMDAGetHeightStratum(DM dm, PetscInt height, PetscInt *pStart, PetscInt *pEnd)
18257459e9aSMatthew G Knepley {
18357459e9aSMatthew G Knepley   DM_DA          *da = (DM_DA*) dm->data;
18457459e9aSMatthew G Knepley   const PetscInt dim = da->dim;
18557459e9aSMatthew G Knepley   PetscInt       nC, nV, nXF, nYF, nZF;
18657459e9aSMatthew G Knepley   PetscErrorCode ierr;
18757459e9aSMatthew G Knepley 
18857459e9aSMatthew G Knepley   PetscFunctionBegin;
1898865f1eaSKarl Rupp   if (pStart) PetscValidIntPointer(pStart,3);
1908865f1eaSKarl Rupp   if (pEnd)   PetscValidIntPointer(pEnd,4);
191e42e3c58SMatthew G. Knepley   ierr = DMDAGetNumCells(dm, NULL, NULL, NULL, &nC);CHKERRQ(ierr);
1920298fd71SBarry Smith   ierr = DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV);CHKERRQ(ierr);
1930298fd71SBarry Smith   ierr = DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF);CHKERRQ(ierr);
19457459e9aSMatthew G Knepley   if (height == 0) {
19557459e9aSMatthew G Knepley     /* Cells */
1968865f1eaSKarl Rupp     if (pStart) *pStart = 0;
1978865f1eaSKarl Rupp     if (pEnd)   *pEnd   = nC;
19857459e9aSMatthew G Knepley   } else if (height == 1) {
19957459e9aSMatthew G Knepley     /* Faces */
2008865f1eaSKarl Rupp     if (pStart) *pStart = nC+nV;
2018865f1eaSKarl Rupp     if (pEnd)   *pEnd   = nC+nV+nXF+nYF+nZF;
20257459e9aSMatthew G Knepley   } else if (height == dim) {
20357459e9aSMatthew G Knepley     /* Vertices */
2048865f1eaSKarl Rupp     if (pStart) *pStart = nC;
2058865f1eaSKarl Rupp     if (pEnd)   *pEnd   = nC+nV;
20657459e9aSMatthew G Knepley   } else if (height < 0) {
20757459e9aSMatthew G Knepley     /* All points */
2088865f1eaSKarl Rupp     if (pStart) *pStart = 0;
2098865f1eaSKarl Rupp     if (pEnd)   *pEnd   = nC+nV+nXF+nYF+nZF;
21082f516ccSBarry Smith   } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No points of height %d in the DA", height);
21157459e9aSMatthew G Knepley   PetscFunctionReturn(0);
21257459e9aSMatthew G Knepley }
21357459e9aSMatthew G Knepley 
21457459e9aSMatthew G Knepley #undef __FUNCT__
215d0892670SMatthew G. Knepley #define __FUNCT__ "DMDAGetDepthStratum"
216d0892670SMatthew G. Knepley PetscErrorCode DMDAGetDepthStratum(DM dm, PetscInt depth, PetscInt *pStart, PetscInt *pEnd)
217d0892670SMatthew G. Knepley {
218d0892670SMatthew G. Knepley   DM_DA         *da  = (DM_DA*) dm->data;
219d0892670SMatthew G. Knepley   const PetscInt dim = da->dim;
220d0892670SMatthew G. Knepley   PetscInt       nC, nV, nXF, nYF, nZF;
221d0892670SMatthew G. Knepley   PetscErrorCode ierr;
222d0892670SMatthew G. Knepley 
223d0892670SMatthew G. Knepley   PetscFunctionBegin;
224d0892670SMatthew G. Knepley   if (pStart) PetscValidIntPointer(pStart,3);
225d0892670SMatthew G. Knepley   if (pEnd)   PetscValidIntPointer(pEnd,4);
226d0892670SMatthew G. Knepley   ierr = DMDAGetNumCells(dm, NULL, NULL, NULL, &nC);CHKERRQ(ierr);
227d0892670SMatthew G. Knepley   ierr = DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV);CHKERRQ(ierr);
228d0892670SMatthew G. Knepley   ierr = DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF);CHKERRQ(ierr);
229d0892670SMatthew G. Knepley   if (depth == dim) {
230d0892670SMatthew G. Knepley     /* Cells */
231d0892670SMatthew G. Knepley     if (pStart) *pStart = 0;
232d0892670SMatthew G. Knepley     if (pEnd)   *pEnd   = nC;
233d0892670SMatthew G. Knepley   } else if (depth == dim-1) {
234d0892670SMatthew G. Knepley     /* Faces */
235d0892670SMatthew G. Knepley     if (pStart) *pStart = nC+nV;
236d0892670SMatthew G. Knepley     if (pEnd)   *pEnd   = nC+nV+nXF+nYF+nZF;
237d0892670SMatthew G. Knepley   } else if (depth == 0) {
238d0892670SMatthew G. Knepley     /* Vertices */
239d0892670SMatthew G. Knepley     if (pStart) *pStart = nC;
240d0892670SMatthew G. Knepley     if (pEnd)   *pEnd   = nC+nV;
241d0892670SMatthew G. Knepley   } else if (depth < 0) {
242d0892670SMatthew G. Knepley     /* All points */
243d0892670SMatthew G. Knepley     if (pStart) *pStart = 0;
244d0892670SMatthew G. Knepley     if (pEnd)   *pEnd   = nC+nV+nXF+nYF+nZF;
245d0892670SMatthew G. Knepley   } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No points of depth %d in the DA", depth);
246d0892670SMatthew G. Knepley   PetscFunctionReturn(0);
247d0892670SMatthew G. Knepley }
248d0892670SMatthew G. Knepley 
249d0892670SMatthew G. Knepley #undef __FUNCT__
250d0892670SMatthew G. Knepley #define __FUNCT__ "DMDAGetConeSize"
251d0892670SMatthew G. Knepley PetscErrorCode DMDAGetConeSize(DM dm, PetscInt p, PetscInt *coneSize)
252d0892670SMatthew G. Knepley {
253d0892670SMatthew G. Knepley   DM_DA         *da  = (DM_DA*) dm->data;
254d0892670SMatthew G. Knepley   const PetscInt dim = da->dim;
255d0892670SMatthew G. Knepley   PetscInt       nC, nV, nXF, nYF, nZF;
256d0892670SMatthew G. Knepley   PetscErrorCode ierr;
257d0892670SMatthew G. Knepley 
258d0892670SMatthew G. Knepley   PetscFunctionBegin;
259d0892670SMatthew G. Knepley   *coneSize = 0;
260d0892670SMatthew G. Knepley   ierr = DMDAGetNumCells(dm, NULL, NULL, NULL, &nC);CHKERRQ(ierr);
261d0892670SMatthew G. Knepley   ierr = DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV);CHKERRQ(ierr);
262d0892670SMatthew G. Knepley   ierr = DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF);CHKERRQ(ierr);
263d0892670SMatthew G. Knepley   switch (dim) {
264d0892670SMatthew G. Knepley   case 2:
265d0892670SMatthew G. Knepley     if (p >= 0) {
266d0892670SMatthew G. Knepley       if (p < nC) {
267d0892670SMatthew G. Knepley         *coneSize = 4;
268d0892670SMatthew G. Knepley       } else if (p < nC+nV) {
269d0892670SMatthew G. Knepley         *coneSize = 0;
270d0892670SMatthew G. Knepley       } else if (p < nC+nV+nXF+nYF+nZF) {
271d0892670SMatthew G. Knepley         *coneSize = 2;
272d0892670SMatthew G. Knepley       } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %d should be in [0, %d)", p, nC+nV+nXF+nYF+nZF);
273d0892670SMatthew G. Knepley     } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative point %d is invalid", p);
274d0892670SMatthew G. Knepley     break;
275d0892670SMatthew G. Knepley   case 3:
276d0892670SMatthew G. Knepley     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do 3D");
277d0892670SMatthew G. Knepley     break;
278d0892670SMatthew G. Knepley   }
279d0892670SMatthew G. Knepley   PetscFunctionReturn(0);
280d0892670SMatthew G. Knepley }
281d0892670SMatthew G. Knepley 
282d0892670SMatthew G. Knepley #undef __FUNCT__
283d0892670SMatthew G. Knepley #define __FUNCT__ "DMDAGetCone"
284d0892670SMatthew G. Knepley PetscErrorCode DMDAGetCone(DM dm, PetscInt p, PetscInt *cone[])
285d0892670SMatthew G. Knepley {
286d0892670SMatthew G. Knepley   DM_DA         *da  = (DM_DA*) dm->data;
287d0892670SMatthew G. Knepley   const PetscInt dim = da->dim;
288d0892670SMatthew G. Knepley   PetscInt       nCx, nCy, nCz, nC, nVx, nVy, nVz, nV, nxF, nyF, nzF, nXF, nYF, nZF;
289d0892670SMatthew G. Knepley   PetscErrorCode ierr;
290d0892670SMatthew G. Knepley 
291d0892670SMatthew G. Knepley   PetscFunctionBegin;
292d0892670SMatthew G. Knepley   if (!cone) {ierr = DMGetWorkArray(dm, 6, PETSC_INT, cone);CHKERRQ(ierr);}
293d0892670SMatthew G. Knepley   ierr = DMDAGetNumCells(dm, &nCx, &nCy, &nCz, &nC);CHKERRQ(ierr);
294d0892670SMatthew G. Knepley   ierr = DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV);CHKERRQ(ierr);
295d0892670SMatthew G. Knepley   ierr = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr);
296d0892670SMatthew G. Knepley   switch (dim) {
297d0892670SMatthew G. Knepley   case 2:
298d0892670SMatthew G. Knepley     if (p >= 0) {
299d0892670SMatthew G. Knepley       if (p < nC) {
300d0892670SMatthew G. Knepley         const PetscInt cy = p / nCx;
301d0892670SMatthew G. Knepley         const PetscInt cx = p % nCx;
302d0892670SMatthew G. Knepley 
303d0892670SMatthew G. Knepley         (*cone)[0] = cy*nxF + cx + nC+nV;
304d0892670SMatthew G. Knepley         (*cone)[1] = cx*nyF + cy + nyF + nC+nV+nXF;
305d0892670SMatthew G. Knepley         (*cone)[2] = cy*nxF + cx + nxF + nC+nV;
306d0892670SMatthew G. Knepley         (*cone)[3] = cx*nyF + cy + nC+nV+nXF;
307d0892670SMatthew G. Knepley         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do cell cones");
308d0892670SMatthew G. Knepley       } else if (p < nC+nV) {
309d0892670SMatthew G. Knepley       } else if (p < nC+nV+nXF) {
310d0892670SMatthew G. Knepley         const PetscInt fy = (p - nC+nV) / nxF;
311d0892670SMatthew G. Knepley         const PetscInt fx = (p - nC+nV) % nxF;
312d0892670SMatthew G. Knepley 
313d0892670SMatthew G. Knepley         (*cone)[0] = fy*nVx + fx + nC;
314d0892670SMatthew G. Knepley         (*cone)[1] = fy*nVx + fx + 1 + nC;
315d0892670SMatthew G. Knepley       } else if (p < nC+nV+nXF+nYF) {
316d0892670SMatthew G. Knepley         const PetscInt fx = (p - nC+nV+nXF) / nyF;
317d0892670SMatthew G. Knepley         const PetscInt fy = (p - nC+nV+nXF) % nyF;
318d0892670SMatthew G. Knepley 
319d0892670SMatthew G. Knepley         (*cone)[0] = fy*nVx + fx + nC;
320d0892670SMatthew G. Knepley         (*cone)[1] = fy*nVx + fx + nVx + nC;
321d0892670SMatthew G. Knepley       } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %d should be in [0, %d)", p, nC+nV+nXF+nYF+nZF);
322d0892670SMatthew G. Knepley     } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative point %d is invalid", p);
323d0892670SMatthew G. Knepley     break;
324d0892670SMatthew G. Knepley   case 3:
325d0892670SMatthew G. Knepley     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do 3D");
326d0892670SMatthew G. Knepley     break;
327d0892670SMatthew G. Knepley   }
328d0892670SMatthew G. Knepley   PetscFunctionReturn(0);
329d0892670SMatthew G. Knepley }
330d0892670SMatthew G. Knepley 
331d0892670SMatthew G. Knepley #undef __FUNCT__
332d0892670SMatthew G. Knepley #define __FUNCT__ "DMDARestoreCone"
333d0892670SMatthew G. Knepley PetscErrorCode DMDARestoreCone(DM dm, PetscInt p, PetscInt *cone[])
334d0892670SMatthew G. Knepley {
335d0892670SMatthew G. Knepley   PetscErrorCode ierr;
336d0892670SMatthew G. Knepley 
337d0892670SMatthew G. Knepley   PetscFunctionBegin;
338d0892670SMatthew G. Knepley   ierr = DMGetWorkArray(dm, 6, PETSC_INT, cone);CHKERRQ(ierr);
339d0892670SMatthew G. Knepley   PetscFunctionReturn(0);
340d0892670SMatthew G. Knepley }
341d0892670SMatthew G. Knepley 
342d0892670SMatthew G. Knepley #undef __FUNCT__
343a66d4d66SMatthew G Knepley #define __FUNCT__ "DMDACreateSection"
344a66d4d66SMatthew G Knepley /*@C
345a66d4d66SMatthew G Knepley   DMDACreateSection - Create a PetscSection inside the DMDA that describes data layout. This allows multiple fields with
346a66d4d66SMatthew G Knepley   different numbers of dofs on vertices, cells, and faces in each direction.
347a66d4d66SMatthew G Knepley 
348a66d4d66SMatthew G Knepley   Input Parameters:
349a66d4d66SMatthew G Knepley + dm- The DMDA
350a66d4d66SMatthew G Knepley . numFields - The number of fields
351a4294c51SMatthew G. Knepley . numComp - The number of components in each field
352a4294c51SMatthew G. Knepley . numDof - The number of dofs per dimension for each field
3530298fd71SBarry Smith . numFaceDof - The number of dofs per face for each field and direction, or NULL
354a66d4d66SMatthew G Knepley 
355a66d4d66SMatthew G Knepley   Level: developer
356a66d4d66SMatthew G Knepley 
357a66d4d66SMatthew G Knepley   Note:
358a66d4d66SMatthew G Knepley   The default DMDA numbering is as follows:
359a66d4d66SMatthew G Knepley 
360a66d4d66SMatthew G Knepley     - Cells:    [0,             nC)
361a66d4d66SMatthew G Knepley     - Vertices: [nC,            nC+nV)
36288ed4aceSMatthew G Knepley     - X-Faces:  [nC+nV,         nC+nV+nXF)         normal is +- x-dir
36388ed4aceSMatthew G Knepley     - Y-Faces:  [nC+nV+nXF,     nC+nV+nXF+nYF)     normal is +- y-dir
36488ed4aceSMatthew G Knepley     - Z-Faces:  [nC+nV+nXF+nYF, nC+nV+nXF+nYF+nZF) normal is +- z-dir
365a66d4d66SMatthew G Knepley 
366a66d4d66SMatthew G Knepley   We interpret the default DMDA partition as a cell partition, and the data assignment as a cell assignment.
367a66d4d66SMatthew G Knepley @*/
368a4294c51SMatthew G. Knepley PetscErrorCode DMDACreateSection(DM dm, PetscInt numComp[], PetscInt numDof[], PetscInt numFaceDof[], PetscSection *s)
369a66d4d66SMatthew G Knepley {
370a66d4d66SMatthew G Knepley   DM_DA            *da  = (DM_DA*) dm->data;
371a4b60ecfSMatthew G. Knepley   PetscSection      section;
37288ed4aceSMatthew G Knepley   const PetscInt    dim = da->dim;
37380800b1aSMatthew G Knepley   PetscInt          numFields, numVertexTotDof = 0, numCellTotDof = 0, numFaceTotDof[3] = {0, 0, 0};
374b2fff234SMatthew G. Knepley   PetscBT           isLeaf;
37588ed4aceSMatthew G Knepley   PetscSF           sf;
376feafbc5cSMatthew G Knepley   PetscMPIInt       rank;
37788ed4aceSMatthew G Knepley   const PetscMPIInt *neighbors;
37888ed4aceSMatthew G Knepley   PetscInt          *localPoints;
37988ed4aceSMatthew G Knepley   PetscSFNode       *remotePoints;
380f5eeb527SMatthew G Knepley   PetscInt          nleaves = 0,  nleavesCheck = 0, nL = 0;
38157459e9aSMatthew G Knepley   PetscInt          nC, nVx, nVy, nVz, nV, nxF, nXF, nyF, nYF, nzF, nZF;
38257459e9aSMatthew G Knepley   PetscInt          pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart, yfEnd, zfStart, zfEnd;
38388ed4aceSMatthew G Knepley   PetscInt          f, v, c, xf, yf, zf, xn, yn, zn;
384a66d4d66SMatthew G Knepley   PetscErrorCode    ierr;
385a66d4d66SMatthew G Knepley 
386a66d4d66SMatthew G Knepley   PetscFunctionBegin;
387a66d4d66SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
388e42e3c58SMatthew G. Knepley   PetscValidPointer(s, 4);
38982f516ccSBarry Smith   ierr    = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
390e42e3c58SMatthew G. Knepley   ierr    = DMDAGetNumCells(dm, NULL, NULL, NULL, &nC);CHKERRQ(ierr);
39157459e9aSMatthew G Knepley   ierr    = DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV);CHKERRQ(ierr);
39257459e9aSMatthew G Knepley   ierr    = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr);
39357459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, -1,  &pStart, &pEnd);CHKERRQ(ierr);
39457459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, 0,   &cStart, &cEnd);CHKERRQ(ierr);
39557459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, 1,   &fStart, &fEnd);CHKERRQ(ierr);
39657459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr);
39757459e9aSMatthew G Knepley   xfStart = vEnd;  xfEnd = xfStart+nXF;
39857459e9aSMatthew G Knepley   yfStart = xfEnd; yfEnd = yfStart+nYF;
39957459e9aSMatthew G Knepley   zfStart = yfEnd; zfEnd = zfStart+nZF;
40082f516ccSBarry Smith   if (zfEnd != fEnd) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Invalid face end %d, should be %d", zfEnd, fEnd);
40188ed4aceSMatthew G Knepley   /* Create local section */
40280800b1aSMatthew G Knepley   ierr = DMDAGetInfo(dm, 0,0,0,0,0,0,0, &numFields, 0,0,0,0,0);CHKERRQ(ierr);
403a66d4d66SMatthew G Knepley   for (f = 0; f < numFields; ++f) {
404a4294c51SMatthew G. Knepley     numVertexTotDof  += numDof[f*(dim+1)+0];
405a4294c51SMatthew G. Knepley     numCellTotDof    += numDof[f*(dim+1)+dim];
406a4294c51SMatthew G. Knepley     numFaceTotDof[0] += dim > 0 ? (numFaceDof ? numFaceDof[f*dim+0] : numDof[f*(dim+1)+dim-1]) : 0;
407a4294c51SMatthew G. Knepley     numFaceTotDof[1] += dim > 1 ? (numFaceDof ? numFaceDof[f*dim+1] : numDof[f*(dim+1)+dim-1]) : 0;
408a4294c51SMatthew G. Knepley     numFaceTotDof[2] += dim > 2 ? (numFaceDof ? numFaceDof[f*dim+2] : numDof[f*(dim+1)+dim-1]) : 0;
409a66d4d66SMatthew G Knepley   }
410a4b60ecfSMatthew G. Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &section);CHKERRQ(ierr);
411a4294c51SMatthew G. Knepley   if (numFields > 0) {
412a4b60ecfSMatthew G. Knepley     ierr = PetscSectionSetNumFields(section, numFields);CHKERRQ(ierr);
413a66d4d66SMatthew G Knepley     for (f = 0; f < numFields; ++f) {
41480800b1aSMatthew G Knepley       const char *name;
41580800b1aSMatthew G Knepley 
41680800b1aSMatthew G Knepley       ierr = DMDAGetFieldName(dm, f, &name);CHKERRQ(ierr);
417a4294c51SMatthew G. Knepley       ierr = PetscSectionSetFieldName(section, f, name ? name : "Field");CHKERRQ(ierr);
41880800b1aSMatthew G Knepley       if (numComp) {
419a4b60ecfSMatthew G. Knepley         ierr = PetscSectionSetFieldComponents(section, f, numComp[f]);CHKERRQ(ierr);
420a66d4d66SMatthew G Knepley       }
421a66d4d66SMatthew G Knepley     }
422a66d4d66SMatthew G Knepley   }
423a4b60ecfSMatthew G. Knepley   ierr = PetscSectionSetChart(section, pStart, pEnd);CHKERRQ(ierr);
424a66d4d66SMatthew G Knepley   for (v = vStart; v < vEnd; ++v) {
425a66d4d66SMatthew G Knepley     for (f = 0; f < numFields; ++f) {
426a4294c51SMatthew G. Knepley       ierr = PetscSectionSetFieldDof(section, v, f, numDof[f*(dim+1)+0]);CHKERRQ(ierr);
427a66d4d66SMatthew G Knepley     }
428a4b60ecfSMatthew G. Knepley     ierr = PetscSectionSetDof(section, v, numVertexTotDof);CHKERRQ(ierr);
429a66d4d66SMatthew G Knepley   }
430a66d4d66SMatthew G Knepley   for (xf = xfStart; xf < xfEnd; ++xf) {
431a66d4d66SMatthew G Knepley     for (f = 0; f < numFields; ++f) {
432a4294c51SMatthew G. Knepley       ierr = PetscSectionSetFieldDof(section, xf, f, numFaceDof ? numFaceDof[f*dim+0] : numDof[f*(dim+1)+dim-1]);CHKERRQ(ierr);
433a66d4d66SMatthew G Knepley     }
434a4b60ecfSMatthew G. Knepley     ierr = PetscSectionSetDof(section, xf, numFaceTotDof[0]);CHKERRQ(ierr);
435a66d4d66SMatthew G Knepley   }
436a66d4d66SMatthew G Knepley   for (yf = yfStart; yf < yfEnd; ++yf) {
437a66d4d66SMatthew G Knepley     for (f = 0; f < numFields; ++f) {
438a4294c51SMatthew G. Knepley       ierr = PetscSectionSetFieldDof(section, yf, f, numFaceDof ? numFaceDof[f*dim+1] : numDof[f*(dim+1)+dim-1]);CHKERRQ(ierr);
439a66d4d66SMatthew G Knepley     }
440a4b60ecfSMatthew G. Knepley     ierr = PetscSectionSetDof(section, yf, numFaceTotDof[1]);CHKERRQ(ierr);
441a66d4d66SMatthew G Knepley   }
442a66d4d66SMatthew G Knepley   for (zf = zfStart; zf < zfEnd; ++zf) {
443a66d4d66SMatthew G Knepley     for (f = 0; f < numFields; ++f) {
444a4294c51SMatthew G. Knepley       ierr = PetscSectionSetFieldDof(section, zf, f, numFaceDof ? numFaceDof[f*dim+2] : numDof[f*(dim+1)+dim-1]);CHKERRQ(ierr);
445a66d4d66SMatthew G Knepley     }
446a4b60ecfSMatthew G. Knepley     ierr = PetscSectionSetDof(section, zf, numFaceTotDof[2]);CHKERRQ(ierr);
447a66d4d66SMatthew G Knepley   }
448a66d4d66SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
449a66d4d66SMatthew G Knepley     for (f = 0; f < numFields; ++f) {
450a4294c51SMatthew G. Knepley       ierr = PetscSectionSetFieldDof(section, c, f, numDof[f*(dim+1)+dim]);CHKERRQ(ierr);
451a66d4d66SMatthew G Knepley     }
452a4b60ecfSMatthew G. Knepley     ierr = PetscSectionSetDof(section, c, numCellTotDof);CHKERRQ(ierr);
453a66d4d66SMatthew G Knepley   }
454a4b60ecfSMatthew G. Knepley   ierr = PetscSectionSetUp(section);CHKERRQ(ierr);
45588ed4aceSMatthew G Knepley   /* Create mesh point SF */
456b2fff234SMatthew G. Knepley   ierr = PetscBTCreate(pEnd-pStart, &isLeaf);CHKERRQ(ierr);
45788ed4aceSMatthew G Knepley   ierr = DMDAGetNeighbors(dm, &neighbors);CHKERRQ(ierr);
45888ed4aceSMatthew G Knepley   for (zn = 0; zn < (dim > 2 ? 3 : 1); ++zn) {
45988ed4aceSMatthew G Knepley     for (yn = 0; yn < (dim > 1 ? 3 : 1); ++yn) {
46088ed4aceSMatthew G Knepley       for (xn = 0; xn < 3; ++xn) {
4617128ae9fSMatthew G Knepley         const PetscInt xp       = xn-1, yp = dim > 1 ? yn-1 : 0, zp = dim > 2 ? zn-1 : 0;
46288ed4aceSMatthew G Knepley         const PetscInt neighbor = neighbors[(zn*3+yn)*3+xn];
463b2fff234SMatthew G. Knepley         PetscInt       xv, yv, zv;
46488ed4aceSMatthew G Knepley 
4653814d604SMatthew G Knepley         if (neighbor >= 0 && neighbor < rank) {
466b2fff234SMatthew G. Knepley           if (xp < 0) { /* left */
467b2fff234SMatthew G. Knepley             if (yp < 0) { /* bottom */
468b2fff234SMatthew G. Knepley               if (zp < 0) { /* back */
469b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = (      0*nVy +     0)*nVx +     0 + nC; /* left bottom back vertex */
470b2fff234SMatthew G. Knepley                 if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
471b2fff234SMatthew G. Knepley               } else if (zp > 0) { /* front */
472b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = ((nVz-1)*nVy +     0)*nVx +     0 + nC; /* left bottom front vertex */
473b2fff234SMatthew G. Knepley                 if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
474b2fff234SMatthew G. Knepley               } else {
475b2fff234SMatthew G. Knepley                 for (zv = 0; zv < nVz; ++zv) {
476b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (zv*nVy +     0)*nVx +     0 + nC; /* left bottom vertices */
477b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
478b2fff234SMatthew G. Knepley                 }
479b2fff234SMatthew G. Knepley               }
480b2fff234SMatthew G. Knepley             } else if (yp > 0) { /* top */
481b2fff234SMatthew G. Knepley               if (zp < 0) { /* back */
482b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = (      0*nVy + nVy-1)*nVx +     0 + nC; /* left top back vertex */
483b2fff234SMatthew G. Knepley                 if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
484b2fff234SMatthew G. Knepley               } else if (zp > 0) { /* front */
485b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = ((nVz-1)*nVy + nVy-1)*nVx +     0 + nC; /* left top front vertex */
486b2fff234SMatthew G. Knepley                 if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
487b2fff234SMatthew G. Knepley               } else {
488b2fff234SMatthew G. Knepley                 for (zv = 0; zv < nVz; ++zv) {
489b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (zv*nVy + nVy-1)*nVx +     0 + nC; /* left top vertices */
490b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
491b2fff234SMatthew G. Knepley                 }
492b2fff234SMatthew G. Knepley               }
493b2fff234SMatthew G. Knepley             } else {
494b2fff234SMatthew G. Knepley               if (zp < 0) { /* back */
495b2fff234SMatthew G. Knepley                 for (yv = 0; yv < nVy; ++yv) {
496b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (      0*nVy + yv)*nVx +     0 + nC; /* left back vertices */
497b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
498b2fff234SMatthew G. Knepley                 }
499b2fff234SMatthew G. Knepley               } else if (zp > 0) { /* front */
500b2fff234SMatthew G. Knepley                 for (yv = 0; yv < nVy; ++yv) {
501b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = ((nVz-1)*nVy + yv)*nVx +     0 + nC; /* left front vertices */
502b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
503b2fff234SMatthew G. Knepley                 }
504b2fff234SMatthew G. Knepley               } else {
505b2fff234SMatthew G. Knepley                 for (zv = 0; zv < nVz; ++zv) {
506b2fff234SMatthew G. Knepley                   for (yv = 0; yv < nVy; ++yv) {
507b2fff234SMatthew G. Knepley                     const PetscInt localVertex  = (zv*nVy + yv)*nVx +     0 + nC; /* left vertices */
508b2fff234SMatthew G. Knepley                     if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
509b2fff234SMatthew G. Knepley                   }
510b2fff234SMatthew G. Knepley                 }
511b2fff234SMatthew G. Knepley #if 0
512b2fff234SMatthew G. Knepley                 for (xf = 0; xf < nxF; ++xf) {
513b2fff234SMatthew G. Knepley                   /* THIS IS WRONG */
514b2fff234SMatthew G. Knepley                   const PetscInt localFace  = 0 + nC+nV; /* left faces */
515b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localFace)) ++nleaves;
516b2fff234SMatthew G. Knepley                 }
517b2fff234SMatthew G. Knepley #endif
518b2fff234SMatthew G. Knepley               }
519b2fff234SMatthew G. Knepley             }
520b2fff234SMatthew G. Knepley           } else if (xp > 0) { /* right */
521b2fff234SMatthew G. Knepley             if (yp < 0) { /* bottom */
522b2fff234SMatthew G. Knepley               if (zp < 0) { /* back */
523b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = (      0*nVy +     0)*nVx + nVx-1 + nC; /* right bottom back vertex */
524b2fff234SMatthew G. Knepley                 if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
525b2fff234SMatthew G. Knepley               } else if (zp > 0) { /* front */
526b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = ((nVz-1)*nVy +     0)*nVx + nVx-1 + nC; /* right bottom front vertex */
527b2fff234SMatthew G. Knepley                 if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
528b2fff234SMatthew G. Knepley               } else {
529b2fff234SMatthew G. Knepley                 for (zv = 0; zv < nVz; ++zv) {
530b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (zv*nVy +     0)*nVx + nVx-1 + nC; /* right bottom vertices */
531b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
532b2fff234SMatthew G. Knepley                 }
533b2fff234SMatthew G. Knepley               }
534b2fff234SMatthew G. Knepley             } else if (yp > 0) { /* top */
535b2fff234SMatthew G. Knepley               if (zp < 0) { /* back */
536b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = (      0*nVy + nVy-1)*nVx + nVx-1 + nC; /* right top back vertex */
537b2fff234SMatthew G. Knepley                 if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
538b2fff234SMatthew G. Knepley               } else if (zp > 0) { /* front */
539b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = ((nVz-1)*nVy + nVy-1)*nVx + nVx-1 + nC; /* right top front vertex */
540b2fff234SMatthew G. Knepley                 if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
541b2fff234SMatthew G. Knepley               } else {
542b2fff234SMatthew G. Knepley                 for (zv = 0; zv < nVz; ++zv) {
543b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (zv*nVy + nVy-1)*nVx + nVx-1 + nC; /* right top vertices */
544b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
545b2fff234SMatthew G. Knepley                 }
546b2fff234SMatthew G. Knepley               }
547b2fff234SMatthew G. Knepley             } else {
548b2fff234SMatthew G. Knepley               if (zp < 0) { /* back */
549b2fff234SMatthew G. Knepley                 for (yv = 0; yv < nVy; ++yv) {
550b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (      0*nVy + yv)*nVx + nVx-1 + nC; /* right back vertices */
551b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
552b2fff234SMatthew G. Knepley                 }
553b2fff234SMatthew G. Knepley               } else if (zp > 0) { /* front */
554b2fff234SMatthew G. Knepley                 for (yv = 0; yv < nVy; ++yv) {
555b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = ((nVz-1)*nVy + yv)*nVx + nVx-1 + nC; /* right front vertices */
556b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
557b2fff234SMatthew G. Knepley                 }
558b2fff234SMatthew G. Knepley               } else {
559b2fff234SMatthew G. Knepley                 for (zv = 0; zv < nVz; ++zv) {
560b2fff234SMatthew G. Knepley                   for (yv = 0; yv < nVy; ++yv) {
561b2fff234SMatthew G. Knepley                     const PetscInt localVertex  = (zv*nVy + yv)*nVx + nVx-1 + nC; /* right vertices */
562b2fff234SMatthew G. Knepley                     if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
563b2fff234SMatthew G. Knepley                   }
564b2fff234SMatthew G. Knepley                 }
565b2fff234SMatthew G. Knepley #if 0
566b2fff234SMatthew G. Knepley                 for (xf = 0; xf < nxF; ++xf) {
567b2fff234SMatthew G. Knepley                   /* THIS IS WRONG */
568b2fff234SMatthew G. Knepley                   const PetscInt localFace  = 0 + nC+nV; /* right faces */
569b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localFace)) ++nleaves;
570b2fff234SMatthew G. Knepley                 }
571b2fff234SMatthew G. Knepley #endif
572b2fff234SMatthew G. Knepley               }
573b2fff234SMatthew G. Knepley             }
574b2fff234SMatthew G. Knepley           } else {
575b2fff234SMatthew G. Knepley             if (yp < 0) { /* bottom */
576b2fff234SMatthew G. Knepley               if (zp < 0) { /* back */
577b2fff234SMatthew G. Knepley                 for (xv = 0; xv < nVx; ++xv) {
578b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (      0*nVy +     0)*nVx + xv + nC; /* bottom back vertices */
579b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
580b2fff234SMatthew G. Knepley                 }
581b2fff234SMatthew G. Knepley               } else if (zp > 0) { /* front */
582b2fff234SMatthew G. Knepley                 for (xv = 0; xv < nVx; ++xv) {
583b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = ((nVz-1)*nVy +     0)*nVx + xv + nC; /* bottom front vertices */
584b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
585b2fff234SMatthew G. Knepley                 }
586b2fff234SMatthew G. Knepley               } else {
587b2fff234SMatthew G. Knepley                 for (zv = 0; zv < nVz; ++zv) {
588b2fff234SMatthew G. Knepley                   for (xv = 0; xv < nVx; ++xv) {
589b2fff234SMatthew G. Knepley                     const PetscInt localVertex  = (zv*nVy +     0)*nVx + xv + nC; /* bottom vertices */
590b2fff234SMatthew G. Knepley                     if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
591b2fff234SMatthew G. Knepley                   }
592b2fff234SMatthew G. Knepley                 }
593b2fff234SMatthew G. Knepley #if 0
594b2fff234SMatthew G. Knepley                 for (yf = 0; yf < nyF; ++yf) {
595b2fff234SMatthew G. Knepley                   /* THIS IS WRONG */
596b2fff234SMatthew G. Knepley                   const PetscInt localFace  = 0 + nC+nV; /* bottom faces */
597b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVFace)) ++nleaves;
598b2fff234SMatthew G. Knepley                 }
599b2fff234SMatthew G. Knepley #endif
600b2fff234SMatthew G. Knepley               }
601b2fff234SMatthew G. Knepley             } else if (yp > 0) { /* top */
602b2fff234SMatthew G. Knepley               if (zp < 0) { /* back */
603b2fff234SMatthew G. Knepley                 for (xv = 0; xv < nVx; ++xv) {
604b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (      0*nVy + nVy-1)*nVx + xv + nC; /* top back vertices */
605b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
606b2fff234SMatthew G. Knepley                 }
607b2fff234SMatthew G. Knepley               } else if (zp > 0) { /* front */
608b2fff234SMatthew G. Knepley                 for (xv = 0; xv < nVx; ++xv) {
609b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = ((nVz-1)*nVy + nVy-1)*nVx + xv + nC; /* top front vertices */
610b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
611b2fff234SMatthew G. Knepley                 }
612b2fff234SMatthew G. Knepley               } else {
613b2fff234SMatthew G. Knepley                 for (zv = 0; zv < nVz; ++zv) {
614b2fff234SMatthew G. Knepley                   for (xv = 0; xv < nVx; ++xv) {
615b2fff234SMatthew G. Knepley                     const PetscInt localVertex  = (zv*nVy + nVy-1)*nVx + xv + nC; /* top vertices */
616b2fff234SMatthew G. Knepley                     if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
617b2fff234SMatthew G. Knepley                   }
618b2fff234SMatthew G. Knepley                 }
619b2fff234SMatthew G. Knepley #if 0
620b2fff234SMatthew G. Knepley                 for (yf = 0; yf < nyF; ++yf) {
621b2fff234SMatthew G. Knepley                   /* THIS IS WRONG */
622b2fff234SMatthew G. Knepley                   const PetscInt localFace  = 0 + nC+nV; /* top faces */
623b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVFace)) ++nleaves;
624b2fff234SMatthew G. Knepley                 }
625b2fff234SMatthew G. Knepley #endif
626b2fff234SMatthew G. Knepley               }
627b2fff234SMatthew G. Knepley             } else {
628b2fff234SMatthew G. Knepley               if (zp < 0) { /* back */
629b2fff234SMatthew G. Knepley                 for (yv = 0; yv < nVy; ++yv) {
630b2fff234SMatthew G. Knepley                   for (xv = 0; xv < nVx; ++xv) {
631b2fff234SMatthew G. Knepley                     const PetscInt localVertex  = (      0*nVy + yv)*nVx + xv + nC; /* back vertices */
632b2fff234SMatthew G. Knepley                     if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
633b2fff234SMatthew G. Knepley                   }
634b2fff234SMatthew G. Knepley                 }
635b2fff234SMatthew G. Knepley #if 0
636b2fff234SMatthew G. Knepley                 for (zf = 0; zf < nzF; ++zf) {
637b2fff234SMatthew G. Knepley                   /* THIS IS WRONG */
638b2fff234SMatthew G. Knepley                   const PetscInt localFace  = 0 + nC+nV; /* back faces */
639b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localFace)) ++nleaves;
640b2fff234SMatthew G. Knepley                 }
641b2fff234SMatthew G. Knepley #endif
642b2fff234SMatthew G. Knepley               } else if (zp > 0) { /* front */
643b2fff234SMatthew G. Knepley                 for (yv = 0; yv < nVy; ++yv) {
644b2fff234SMatthew G. Knepley                   for (xv = 0; xv < nVx; ++xv) {
645b2fff234SMatthew G. Knepley                     const PetscInt localVertex  = ((nVz-1)*nVy + yv)*nVx + xv + nC; /* front vertices */
646b2fff234SMatthew G. Knepley                     if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
647b2fff234SMatthew G. Knepley                   }
648b2fff234SMatthew G. Knepley                 }
649b2fff234SMatthew G. Knepley #if 0
650b2fff234SMatthew G. Knepley                 for (zf = 0; zf < nzF; ++zf) {
651b2fff234SMatthew G. Knepley                   /* THIS IS WRONG */
652b2fff234SMatthew G. Knepley                   const PetscInt localFace  = 0 + nC+nV; /* front faces */
653b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localFace)) ++nleaves;
654b2fff234SMatthew G. Knepley                 }
655b2fff234SMatthew G. Knepley #endif
656b2fff234SMatthew G. Knepley               } else {
657b2fff234SMatthew G. Knepley                 /* Nothing is shared from the interior */
65888ed4aceSMatthew G Knepley               }
65988ed4aceSMatthew G Knepley             }
66088ed4aceSMatthew G Knepley           }
66188ed4aceSMatthew G Knepley         }
66288ed4aceSMatthew G Knepley       }
663b2fff234SMatthew G. Knepley     }
664b2fff234SMatthew G. Knepley   }
665b2fff234SMatthew G. Knepley   ierr = PetscBTMemzero(pEnd-pStart, isLeaf);CHKERRQ(ierr);
66688ed4aceSMatthew G Knepley   ierr = PetscMalloc2(nleaves,PetscInt,&localPoints,nleaves,PetscSFNode,&remotePoints);CHKERRQ(ierr);
66788ed4aceSMatthew G Knepley   for (zn = 0; zn < (dim > 2 ? 3 : 1); ++zn) {
66888ed4aceSMatthew G Knepley     for (yn = 0; yn < (dim > 1 ? 3 : 1); ++yn) {
66988ed4aceSMatthew G Knepley       for (xn = 0; xn < 3; ++xn) {
6707128ae9fSMatthew G Knepley         const PetscInt xp       = xn-1, yp = dim > 1 ? yn-1 : 0, zp = dim > 2 ? zn-1 : 0;
67188ed4aceSMatthew G Knepley         const PetscInt neighbor = neighbors[(zn*3+yn)*3+xn];
672f5eeb527SMatthew G Knepley         PetscInt       xv, yv, zv;
67388ed4aceSMatthew G Knepley 
6743814d604SMatthew G Knepley         if (neighbor >= 0 && neighbor < rank) {
67588ed4aceSMatthew G Knepley           if (xp < 0) { /* left */
67688ed4aceSMatthew G Knepley             if (yp < 0) { /* bottom */
67788ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
678b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = (      0*nVy +     0)*nVx +     0 + nC; /* left bottom back vertex */
679f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = ((nVz-1)*nVy + nVy-1)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
680f5eeb527SMatthew G Knepley 
681b2fff234SMatthew G. Knepley                 if (!PetscBTLookupSet(isLeaf, localVertex)) {
682f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
683f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
684f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
685f5eeb527SMatthew G Knepley                   ++nL;
686b2fff234SMatthew G. Knepley                 }
68788ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
688b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = ((nVz-1)*nVy +     0)*nVx +     0 + nC; /* left bottom front vertex */
689f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = (      0*nVy + nVy-1)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
690f5eeb527SMatthew G Knepley 
691b2fff234SMatthew G. Knepley                 if (!PetscBTLookupSet(isLeaf, localVertex)) {
692f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
693f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
694f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
695f5eeb527SMatthew G Knepley                   ++nL;
696b2fff234SMatthew G. Knepley                 }
69788ed4aceSMatthew G Knepley               } else {
698b2fff234SMatthew G. Knepley                 for (zv = 0; zv < nVz; ++zv) {
699b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (zv*nVy +     0)*nVx +     0 + nC; /* left bottom vertices */
700f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (zv*nVy + nVy-1)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
701f5eeb527SMatthew G Knepley 
702b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) {
703f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
704f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
705f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
706b2fff234SMatthew G. Knepley                     ++nL;
707b2fff234SMatthew G. Knepley                   }
708f5eeb527SMatthew G Knepley                 }
70988ed4aceSMatthew G Knepley               }
71088ed4aceSMatthew G Knepley             } else if (yp > 0) { /* top */
71188ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
712b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = (      0*nVy + nVy-1)*nVx +     0 + nC; /* left top back vertex */
713f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = ((nVz-1)*nVy +     0)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
714f5eeb527SMatthew G Knepley 
715b2fff234SMatthew G. Knepley                 if (!PetscBTLookupSet(isLeaf, localVertex)) {
716f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
717f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
718f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
719f5eeb527SMatthew G Knepley                   ++nL;
720b2fff234SMatthew G. Knepley                 }
72188ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
722b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = ((nVz-1)*nVy + nVy-1)*nVx +     0 + nC; /* left top front vertex */
723f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = (      0*nVy +     0)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
724f5eeb527SMatthew G Knepley 
725b2fff234SMatthew G. Knepley                 if (!PetscBTLookupSet(isLeaf, localVertex)) {
726f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
727f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
728f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
729f5eeb527SMatthew G Knepley                   ++nL;
730b2fff234SMatthew G. Knepley                 }
73188ed4aceSMatthew G Knepley               } else {
732b2fff234SMatthew G. Knepley                 for (zv = 0; zv < nVz; ++zv) {
733b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (zv*nVy + nVy-1)*nVx +     0 + nC; /* left top vertices */
734f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (zv*nVy +     0)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
735f5eeb527SMatthew G Knepley 
736b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) {
737f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
738f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
739f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
740b2fff234SMatthew G. Knepley                     ++nL;
741b2fff234SMatthew G. Knepley                   }
742f5eeb527SMatthew G Knepley                 }
74388ed4aceSMatthew G Knepley               }
74488ed4aceSMatthew G Knepley             } else {
74588ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
746b2fff234SMatthew G. Knepley                 for (yv = 0; yv < nVy; ++yv) {
747b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (      0*nVy + yv)*nVx +     0 + nC; /* left back vertices */
748f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = ((nVz-1)*nVy + yv)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
749f5eeb527SMatthew G Knepley 
750b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) {
751f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
752f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
753f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
754b2fff234SMatthew G. Knepley                     ++nL;
755b2fff234SMatthew G. Knepley                   }
756f5eeb527SMatthew G Knepley                 }
75788ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
758b2fff234SMatthew G. Knepley                 for (yv = 0; yv < nVy; ++yv) {
759b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = ((nVz-1)*nVy + yv)*nVx +     0 + nC; /* left front vertices */
760f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (      0*nVy + yv)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
761f5eeb527SMatthew G Knepley 
762b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) {
763f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
764f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
765f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
766b2fff234SMatthew G. Knepley                     ++nL;
767b2fff234SMatthew G. Knepley                   }
768f5eeb527SMatthew G Knepley                 }
76988ed4aceSMatthew G Knepley               } else {
770f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv) {
771b2fff234SMatthew G. Knepley                   for (yv = 0; yv < nVy; ++yv) {
772b2fff234SMatthew G. Knepley                     const PetscInt localVertex  = (zv*nVy + yv)*nVx +     0 + nC; /* left vertices */
773f5eeb527SMatthew G Knepley                     const PetscInt remoteVertex = (zv*nVy + yv)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
774f5eeb527SMatthew G Knepley 
775b2fff234SMatthew G. Knepley                     if (!PetscBTLookupSet(isLeaf, localVertex)) {
776f5eeb527SMatthew G Knepley                       localPoints[nL]        = localVertex;
777f5eeb527SMatthew G Knepley                       remotePoints[nL].rank  = neighbor;
778f5eeb527SMatthew G Knepley                       remotePoints[nL].index = remoteVertex;
779b2fff234SMatthew G. Knepley                       ++nL;
780f5eeb527SMatthew G Knepley                     }
781f5eeb527SMatthew G Knepley                   }
782b2fff234SMatthew G. Knepley                 }
783b2fff234SMatthew G. Knepley #if 0
784b2fff234SMatthew G. Knepley                 for (xf = 0; xf < nxF; ++xf) {
785f5eeb527SMatthew G Knepley                   /* THIS IS WRONG */
786b2fff234SMatthew G. Knepley                   const PetscInt localFace  = 0 + nC+nV; /* left faces */
787f5eeb527SMatthew G Knepley                   const PetscInt remoteFace = 0 + nC+nV;
788f5eeb527SMatthew G Knepley 
789b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localFace)) {
790f5eeb527SMatthew G Knepley                     localPoints[nL]        = localFace;
791f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
792f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteFace;
793f5eeb527SMatthew G Knepley                   }
79488ed4aceSMatthew G Knepley                 }
795b2fff234SMatthew G. Knepley #endif
796b2fff234SMatthew G. Knepley               }
79788ed4aceSMatthew G Knepley             }
79888ed4aceSMatthew G Knepley           } else if (xp > 0) { /* right */
79988ed4aceSMatthew G Knepley             if (yp < 0) { /* bottom */
80088ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
801b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = (      0*nVy +     0)*nVx + nVx-1 + nC; /* right bottom back vertex */
802f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = ((nVz-1)*nVy + nVy-1)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
803f5eeb527SMatthew G Knepley 
804b2fff234SMatthew G. Knepley                 if (!PetscBTLookupSet(isLeaf, localVertex)) {
805f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
806f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
807f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
808f5eeb527SMatthew G Knepley                   ++nL;
809b2fff234SMatthew G. Knepley                 }
81088ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
811b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = ((nVz-1)*nVy +     0)*nVx + nVx-1 + nC; /* right bottom front vertex */
812f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = (      0*nVy + nVy-1)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
813f5eeb527SMatthew G Knepley 
814b2fff234SMatthew G. Knepley                 if (!PetscBTLookupSet(isLeaf, localVertex)) {
815f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
816f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
817f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
818f5eeb527SMatthew G Knepley                   ++nL;
819b2fff234SMatthew G. Knepley                 }
82088ed4aceSMatthew G Knepley               } else {
821b2fff234SMatthew G. Knepley                 nleavesCheck += nVz;
822b2fff234SMatthew G. Knepley                 for (zv = 0; zv < nVz; ++zv) {
823b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (zv*nVy +     0)*nVx + nVx-1 + nC; /* right bottom vertices */
824f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (zv*nVy + nVy-1)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
825f5eeb527SMatthew G Knepley 
826b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) {
827f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
828f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
829f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
830b2fff234SMatthew G. Knepley                     ++nL;
831b2fff234SMatthew G. Knepley                   }
832f5eeb527SMatthew G Knepley                 }
83388ed4aceSMatthew G Knepley               }
83488ed4aceSMatthew G Knepley             } else if (yp > 0) { /* top */
83588ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
836b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = (      0*nVy + nVy-1)*nVx + nVx-1 + nC; /* right top back vertex */
837f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = ((nVz-1)*nVy +     0)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
838f5eeb527SMatthew G Knepley 
839b2fff234SMatthew G. Knepley                 if (!PetscBTLookupSet(isLeaf, localVertex)) {
840f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
841f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
842f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
843f5eeb527SMatthew G Knepley                   ++nL;
844b2fff234SMatthew G. Knepley                 }
84588ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
846b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = ((nVz-1)*nVy + nVy-1)*nVx + nVx-1 + nC; /* right top front vertex */
847f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = (      0*nVy +     0)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
848f5eeb527SMatthew G Knepley 
849b2fff234SMatthew G. Knepley                 if (!PetscBTLookupSet(isLeaf, localVertex)) {
850f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
851f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
852f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
853f5eeb527SMatthew G Knepley                   ++nL;
854b2fff234SMatthew G. Knepley                 }
85588ed4aceSMatthew G Knepley               } else {
856b2fff234SMatthew G. Knepley                 for (zv = 0; zv < nVz; ++zv) {
857b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (zv*nVy + nVy-1)*nVx + nVx-1 + nC; /* right top vertices */
858f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (zv*nVy +     0)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
859f5eeb527SMatthew G Knepley 
860b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) {
861f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
862f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
863f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
864b2fff234SMatthew G. Knepley                     ++nL;
865b2fff234SMatthew G. Knepley                   }
866f5eeb527SMatthew G Knepley                 }
86788ed4aceSMatthew G Knepley               }
86888ed4aceSMatthew G Knepley             } else {
86988ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
870b2fff234SMatthew G. Knepley                 for (yv = 0; yv < nVy; ++yv) {
871b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (      0*nVy + yv)*nVx + nVx-1 + nC; /* right back vertices */
872f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = ((nVz-1)*nVy + yv)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
873f5eeb527SMatthew G Knepley 
874b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) {
875f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
876f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
877f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
878b2fff234SMatthew G. Knepley                     ++nL;
879b2fff234SMatthew G. Knepley                   }
880f5eeb527SMatthew G Knepley                 }
88188ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
882b2fff234SMatthew G. Knepley                 for (yv = 0; yv < nVy; ++yv) {
883b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = ((nVz-1)*nVy + yv)*nVx + nVx-1 + nC; /* right front vertices */
884f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (      0*nVy + yv)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
885f5eeb527SMatthew G Knepley 
886b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) {
887f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
888f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
889f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
890b2fff234SMatthew G. Knepley                     ++nL;
891b2fff234SMatthew G. Knepley                   }
892f5eeb527SMatthew G Knepley                 }
89388ed4aceSMatthew G Knepley               } else {
894f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv) {
895b2fff234SMatthew G. Knepley                   for (yv = 0; yv < nVy; ++yv) {
896b2fff234SMatthew G. Knepley                     const PetscInt localVertex  = (zv*nVy + yv)*nVx + nVx-1 + nC; /* right vertices */
897f5eeb527SMatthew G Knepley                     const PetscInt remoteVertex = (zv*nVy + yv)*nVx + 0     + nC; /* TODO: Correct this for neighbor sizes */
898f5eeb527SMatthew G Knepley 
899b2fff234SMatthew G. Knepley                     if (!PetscBTLookupSet(isLeaf, localVertex)) {
900f5eeb527SMatthew G Knepley                       localPoints[nL]        = localVertex;
901f5eeb527SMatthew G Knepley                       remotePoints[nL].rank  = neighbor;
902f5eeb527SMatthew G Knepley                       remotePoints[nL].index = remoteVertex;
903b2fff234SMatthew G. Knepley                       ++nL;
904f5eeb527SMatthew G Knepley                     }
905f5eeb527SMatthew G Knepley                   }
906b2fff234SMatthew G. Knepley                 }
907b2fff234SMatthew G. Knepley #if 0
908b2fff234SMatthew G. Knepley                 for (xf = 0; xf < nxF; ++xf) {
909f5eeb527SMatthew G Knepley                   /* THIS IS WRONG */
910b2fff234SMatthew G. Knepley                   const PetscInt localFace  = 0 + nC+nV; /* right faces */
911f5eeb527SMatthew G Knepley                   const PetscInt remoteFace = 0 + nC+nV;
912f5eeb527SMatthew G Knepley 
913b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localFace)) {
914f5eeb527SMatthew G Knepley                     localPoints[nL]        = localFace;
915f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
916f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteFace;
917b2fff234SMatthew G. Knepley                     ++nL;
918f5eeb527SMatthew G Knepley                   }
91988ed4aceSMatthew G Knepley                 }
920b2fff234SMatthew G. Knepley #endif
921b2fff234SMatthew G. Knepley               }
92288ed4aceSMatthew G Knepley             }
92388ed4aceSMatthew G Knepley           } else {
92488ed4aceSMatthew G Knepley             if (yp < 0) { /* bottom */
92588ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
926b2fff234SMatthew G. Knepley                 for (xv = 0; xv < nVx; ++xv) {
927b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (      0*nVy +     0)*nVx + xv + nC; /* bottom back vertices */
928f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = ((nVz-1)*nVy + nVy-1)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
929f5eeb527SMatthew G Knepley 
930b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) {
931f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
932f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
933f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
934b2fff234SMatthew G. Knepley                     ++nL;
935b2fff234SMatthew G. Knepley                   }
936f5eeb527SMatthew G Knepley                 }
93788ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
938b2fff234SMatthew G. Knepley                 for (xv = 0; xv < nVx; ++xv) {
939b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = ((nVz-1)*nVy +     0)*nVx + xv + nC; /* bottom front vertices */
940f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (      0*nVy + nVy-1)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
941f5eeb527SMatthew G Knepley 
942b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) {
943f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
944f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
945f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
946b2fff234SMatthew G. Knepley                     ++nL;
947b2fff234SMatthew G. Knepley                   }
948f5eeb527SMatthew G Knepley                 }
94988ed4aceSMatthew G Knepley               } else {
950f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv) {
951b2fff234SMatthew G. Knepley                   for (xv = 0; xv < nVx; ++xv) {
952b2fff234SMatthew G. Knepley                     const PetscInt localVertex  = (zv*nVy +     0)*nVx + xv + nC; /* bottom vertices */
953f5eeb527SMatthew G Knepley                     const PetscInt remoteVertex = (zv*nVy + nVy-1)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
954f5eeb527SMatthew G Knepley 
955b2fff234SMatthew G. Knepley                     if (!PetscBTLookupSet(isLeaf, localVertex)) {
956f5eeb527SMatthew G Knepley                       localPoints[nL]        = localVertex;
957f5eeb527SMatthew G Knepley                       remotePoints[nL].rank  = neighbor;
958f5eeb527SMatthew G Knepley                       remotePoints[nL].index = remoteVertex;
959b2fff234SMatthew G. Knepley                       ++nL;
960f5eeb527SMatthew G Knepley                     }
961f5eeb527SMatthew G Knepley                   }
962b2fff234SMatthew G. Knepley                 }
963b2fff234SMatthew G. Knepley #if 0
964b2fff234SMatthew G. Knepley                 for (yf = 0; yf < nyF; ++yf) {
965f5eeb527SMatthew G Knepley                   /* THIS IS WRONG */
966b2fff234SMatthew G. Knepley                   const PetscInt localFace  = 0 + nC+nV; /* bottom faces */
967f5eeb527SMatthew G Knepley                   const PetscInt remoteFace = 0 + nC+nV;
968f5eeb527SMatthew G Knepley 
969b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localFace)) {
970f5eeb527SMatthew G Knepley                     localPoints[nL]        = localFace;
971f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
972f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteFace;
973b2fff234SMatthew G. Knepley                     ++nL;
974f5eeb527SMatthew G Knepley                   }
97588ed4aceSMatthew G Knepley                 }
976b2fff234SMatthew G. Knepley #endif
977b2fff234SMatthew G. Knepley               }
97888ed4aceSMatthew G Knepley             } else if (yp > 0) { /* top */
97988ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
980b2fff234SMatthew G. Knepley                 for (xv = 0; xv < nVx; ++xv) {
981b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (      0*nVy + nVy-1)*nVx + xv + nC; /* top back vertices */
982f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = ((nVz-1)*nVy +     0)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
983f5eeb527SMatthew G Knepley 
984b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) {
985f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
986f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
987f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
988b2fff234SMatthew G. Knepley                     ++nL;
989b2fff234SMatthew G. Knepley                   }
990f5eeb527SMatthew G Knepley                 }
99188ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
992b2fff234SMatthew G. Knepley                 for (xv = 0; xv < nVx; ++xv) {
993b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = ((nVz-1)*nVy + nVy-1)*nVx + xv + nC; /* top front vertices */
994f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (      0*nVy +     0)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
995f5eeb527SMatthew G Knepley 
996b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) {
997f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
998f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
999f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
1000b2fff234SMatthew G. Knepley                     ++nL;
1001b2fff234SMatthew G. Knepley                   }
1002f5eeb527SMatthew G Knepley                 }
100388ed4aceSMatthew G Knepley               } else {
1004f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv) {
1005b2fff234SMatthew G. Knepley                   for (xv = 0; xv < nVx; ++xv) {
1006b2fff234SMatthew G. Knepley                     const PetscInt localVertex  = (zv*nVy + nVy-1)*nVx + xv + nC; /* top vertices */
1007f5eeb527SMatthew G Knepley                     const PetscInt remoteVertex = (zv*nVy +     0)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
1008f5eeb527SMatthew G Knepley 
1009b2fff234SMatthew G. Knepley                     if (!PetscBTLookupSet(isLeaf, localVertex)) {
1010f5eeb527SMatthew G Knepley                       localPoints[nL]        = localVertex;
1011f5eeb527SMatthew G Knepley                       remotePoints[nL].rank  = neighbor;
1012f5eeb527SMatthew G Knepley                       remotePoints[nL].index = remoteVertex;
1013b2fff234SMatthew G. Knepley                       ++nL;
1014f5eeb527SMatthew G Knepley                     }
1015f5eeb527SMatthew G Knepley                   }
1016b2fff234SMatthew G. Knepley                 }
1017b2fff234SMatthew G. Knepley #if 0
1018b2fff234SMatthew G. Knepley                 for (yf = 0; yf < nyF; ++yf) {
1019f5eeb527SMatthew G Knepley                   /* THIS IS WRONG */
1020b2fff234SMatthew G. Knepley                   const PetscInt localFace  = 0 + nC+nV; /* top faces */
1021f5eeb527SMatthew G Knepley                   const PetscInt remoteFace = 0 + nC+nV;
1022f5eeb527SMatthew G Knepley 
1023b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localFace)) {
1024f5eeb527SMatthew G Knepley                     localPoints[nL]        = localFace;
1025f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
1026f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteFace;
1027b2fff234SMatthew G. Knepley                     ++nL;
1028f5eeb527SMatthew G Knepley                   }
102988ed4aceSMatthew G Knepley                 }
1030b2fff234SMatthew G. Knepley #endif
1031b2fff234SMatthew G. Knepley               }
103288ed4aceSMatthew G Knepley             } else {
103388ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
1034f5eeb527SMatthew G Knepley                 for (yv = 0; yv < nVy; ++yv) {
1035b2fff234SMatthew G. Knepley                   for (xv = 0; xv < nVx; ++xv) {
1036b2fff234SMatthew G. Knepley                     const PetscInt localVertex  = (      0*nVy + yv)*nVx + xv + nC; /* back vertices */
1037f5eeb527SMatthew G Knepley                     const PetscInt remoteVertex = ((nVz-1)*nVy + yv)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
1038f5eeb527SMatthew G Knepley 
1039b2fff234SMatthew G. Knepley                     if (!PetscBTLookupSet(isLeaf, localVertex)) {
1040f5eeb527SMatthew G Knepley                       localPoints[nL]        = localVertex;
1041f5eeb527SMatthew G Knepley                       remotePoints[nL].rank  = neighbor;
1042f5eeb527SMatthew G Knepley                       remotePoints[nL].index = remoteVertex;
1043b2fff234SMatthew G. Knepley                       ++nL;
1044f5eeb527SMatthew G Knepley                     }
1045f5eeb527SMatthew G Knepley                   }
1046b2fff234SMatthew G. Knepley                 }
1047b2fff234SMatthew G. Knepley #if 0
1048b2fff234SMatthew G. Knepley                 for (zf = 0; zf < nzF; ++zf) {
1049f5eeb527SMatthew G Knepley                   /* THIS IS WRONG */
1050b2fff234SMatthew G. Knepley                   const PetscInt localFace  = 0 + nC+nV; /* back faces */
1051f5eeb527SMatthew G Knepley                   const PetscInt remoteFace = 0 + nC+nV;
1052f5eeb527SMatthew G Knepley 
1053b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localFace)) {
1054f5eeb527SMatthew G Knepley                     localPoints[nL]        = localFace;
1055f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
1056f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteFace;
1057b2fff234SMatthew G. Knepley                     ++nL;
1058f5eeb527SMatthew G Knepley                   }
1059b2fff234SMatthew G. Knepley                 }
1060b2fff234SMatthew G. Knepley #endif
106188ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
1062f5eeb527SMatthew G Knepley                 for (yv = 0; yv < nVy; ++yv) {
1063b2fff234SMatthew G. Knepley                   for (xv = 0; xv < nVx; ++xv) {
1064b2fff234SMatthew G. Knepley                     const PetscInt localVertex  = ((nVz-1)*nVy + yv)*nVx + xv + nC; /* front vertices */
1065f5eeb527SMatthew G Knepley                     const PetscInt remoteVertex = (      0*nVy + yv)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
1066f5eeb527SMatthew G Knepley 
1067b2fff234SMatthew G. Knepley                     if (!PetscBTLookupSet(isLeaf, localVertex)) {
1068f5eeb527SMatthew G Knepley                       localPoints[nL]        = localVertex;
1069f5eeb527SMatthew G Knepley                       remotePoints[nL].rank  = neighbor;
1070f5eeb527SMatthew G Knepley                       remotePoints[nL].index = remoteVertex;
1071b2fff234SMatthew G. Knepley                       ++nL;
1072f5eeb527SMatthew G Knepley                     }
1073f5eeb527SMatthew G Knepley                   }
1074b2fff234SMatthew G. Knepley                 }
1075b2fff234SMatthew G. Knepley #if 0
1076b2fff234SMatthew G. Knepley                 for (zf = 0; zf < nzF; ++zf) {
1077f5eeb527SMatthew G Knepley                   /* THIS IS WRONG */
1078b2fff234SMatthew G. Knepley                   const PetscInt localFace  = 0 + nC+nV; /* front faces */
1079f5eeb527SMatthew G Knepley                   const PetscInt remoteFace = 0 + nC+nV;
1080f5eeb527SMatthew G Knepley 
1081b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localFace)) {
1082f5eeb527SMatthew G Knepley                     localPoints[nL]        = localFace;
1083f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
1084f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteFace;
1085b2fff234SMatthew G. Knepley                     ++nL;
1086f5eeb527SMatthew G Knepley                   }
1087b2fff234SMatthew G. Knepley                 }
1088b2fff234SMatthew G. Knepley #endif
108988ed4aceSMatthew G Knepley               } else {
109088ed4aceSMatthew G Knepley                 /* Nothing is shared from the interior */
109188ed4aceSMatthew G Knepley               }
109288ed4aceSMatthew G Knepley             }
109388ed4aceSMatthew G Knepley           }
109488ed4aceSMatthew G Knepley         }
109588ed4aceSMatthew G Knepley       }
109688ed4aceSMatthew G Knepley     }
109788ed4aceSMatthew G Knepley   }
1098b2fff234SMatthew G. Knepley   ierr = PetscBTDestroy(&isLeaf);CHKERRQ(ierr);
1099b2fff234SMatthew G. Knepley   /* Remove duplication in leaf determination */
1100b2fff234SMatthew G. Knepley   if (nleaves != nL) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "The number of leaves %d did not match the number of remote leaves %d", nleaves, nleavesCheck);
110182f516ccSBarry Smith   ierr = PetscSFCreate(PetscObjectComm((PetscObject)dm), &sf);CHKERRQ(ierr);
110288ed4aceSMatthew G Knepley   ierr = PetscSFSetGraph(sf, pEnd, nleaves, localPoints, PETSC_OWN_POINTER, remotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr);
1103a4b60ecfSMatthew G. Knepley   ierr = DMSetPointSF(dm, sf);CHKERRQ(ierr);
110488ed4aceSMatthew G Knepley   ierr = PetscSFDestroy(&sf);CHKERRQ(ierr);
1105a4294c51SMatthew G. Knepley   *s = section;
1106a4294c51SMatthew G. Knepley   PetscFunctionReturn(0);
1107a4294c51SMatthew G. Knepley }
1108a4294c51SMatthew G. Knepley 
1109d0892670SMatthew G. Knepley #undef __FUNCT__
1110d0892670SMatthew G. Knepley #define __FUNCT__ "DMDASetVertexCoordinates"
1111d0892670SMatthew G. Knepley PetscErrorCode DMDASetVertexCoordinates(DM dm, PetscReal xl, PetscReal xu, PetscReal yl, PetscReal yu, PetscReal zl, PetscReal zu)
1112d0892670SMatthew G. Knepley {
1113d0892670SMatthew G. Knepley   DM_DA         *da = (DM_DA *) dm->data;
1114d0892670SMatthew G. Knepley   Vec            coordinates;
1115d0892670SMatthew G. Knepley   PetscSection   section;
1116d0892670SMatthew G. Knepley   PetscScalar   *coords;
1117d0892670SMatthew G. Knepley   PetscReal      h[3];
1118d0892670SMatthew G. Knepley   PetscInt       dim, size, M, N, P, nVx, nVy, nVz, nV, vStart, vEnd, v, i, j, k;
1119d0892670SMatthew G. Knepley   PetscErrorCode ierr;
1120d0892670SMatthew G. Knepley 
1121d0892670SMatthew G. Knepley   PetscFunctionBegin;
1122d0892670SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1123d0892670SMatthew G. Knepley   ierr = DMDAGetInfo(dm, &dim, &M, &N, &P, 0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
1124d0892670SMatthew G. Knepley   h[0] = (xu - xl)/M;
1125d0892670SMatthew G. Knepley   h[1] = (yu - yl)/N;
1126d0892670SMatthew G. Knepley   h[2] = (zu - zl)/P;
1127d0892670SMatthew G. Knepley   ierr = DMDAGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
1128d0892670SMatthew G. Knepley   ierr = DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV);CHKERRQ(ierr);
1129d0892670SMatthew G. Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &section);CHKERRQ(ierr);
1130d0892670SMatthew G. Knepley   ierr = PetscSectionSetNumFields(section, 1);CHKERRQ(ierr);
1131d0892670SMatthew G. Knepley   ierr = PetscSectionSetFieldComponents(section, 0, dim);CHKERRQ(ierr);
1132d0892670SMatthew G. Knepley   ierr = PetscSectionSetChart(section, vStart, vEnd);CHKERRQ(ierr);
1133d0892670SMatthew G. Knepley   for (v = vStart; v < vEnd; ++v) {
1134d0892670SMatthew G. Knepley     ierr = PetscSectionSetDof(section, v, dim);CHKERRQ(ierr);
1135d0892670SMatthew G. Knepley   }
1136d0892670SMatthew G. Knepley   ierr = PetscSectionSetUp(section);CHKERRQ(ierr);
1137d0892670SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(section, &size);CHKERRQ(ierr);
1138d0892670SMatthew G. Knepley   ierr = VecCreateSeq(PETSC_COMM_SELF, size, &coordinates);CHKERRQ(ierr);
1139d0892670SMatthew G. Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
1140d0892670SMatthew G. Knepley   for (k = 0; k < nVz; ++k) {
1141d0892670SMatthew G. Knepley     PetscInt ind[3] = {0, 0, k + da->zs}, d, off;
1142d0892670SMatthew G. Knepley 
1143d0892670SMatthew G. Knepley     for (j = 0; j < nVy; ++j) {
1144d0892670SMatthew G. Knepley       ind[1] = j + da->ys;
1145d0892670SMatthew G. Knepley       for (i = 0; i < nVx; ++i) {
1146d0892670SMatthew G. Knepley         const PetscInt vertex = (k*nVy + j)*nVx + i + vStart;
1147d0892670SMatthew G. Knepley 
1148d0892670SMatthew G. Knepley         ierr = PetscSectionGetOffset(section, vertex, &off);CHKERRQ(ierr);
1149d0892670SMatthew G. Knepley         ind[0] = i + da->xs;
1150d0892670SMatthew G. Knepley         for (d = 0; d < dim; ++d) {
1151d0892670SMatthew G. Knepley           coords[off+d] = h[d]*ind[d];
1152d0892670SMatthew G. Knepley         }
1153d0892670SMatthew G. Knepley       }
1154d0892670SMatthew G. Knepley     }
1155d0892670SMatthew G. Knepley   }
1156d0892670SMatthew G. Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
1157d0892670SMatthew G. Knepley   ierr = DMSetCoordinateSection(dm, section);CHKERRQ(ierr);
1158d0892670SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
1159a4b60ecfSMatthew G. Knepley   ierr = PetscSectionDestroy(&section);CHKERRQ(ierr);
1160d0892670SMatthew G. Knepley   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
1161d0892670SMatthew G. Knepley   PetscFunctionReturn(0);
1162d0892670SMatthew G. Knepley }
116337b16e26SMatthew G. Knepley 
116437b16e26SMatthew G. Knepley #undef __FUNCT__
116537b16e26SMatthew G. Knepley #define __FUNCT__ "DMDAProjectFunctionLocal"
116637b16e26SMatthew G. Knepley PetscErrorCode DMDAProjectFunctionLocal(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *), InsertMode mode, Vec localX)
116737b16e26SMatthew G. Knepley {
116837b16e26SMatthew G. Knepley   PetscDualSpace *sp;
116937b16e26SMatthew G. Knepley   PetscQuadrature q;
117037b16e26SMatthew G. Knepley   PetscSection    section;
117137b16e26SMatthew G. Knepley   PetscScalar    *values;
117237b16e26SMatthew G. Knepley   PetscReal      *v0, *J, *detJ;
117337b16e26SMatthew G. Knepley   PetscInt        numFields, numComp, dim, spDim, totDim = 0, numValues, cStart, cEnd, f, c, v, d;
117437b16e26SMatthew G. Knepley   PetscErrorCode  ierr;
117537b16e26SMatthew G. Knepley 
117637b16e26SMatthew G. Knepley   PetscFunctionBegin;
117737b16e26SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
117837b16e26SMatthew G. Knepley   ierr = PetscFEGetQuadrature(fe[0], &q);CHKERRQ(ierr);
117937b16e26SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
118037b16e26SMatthew G. Knepley   ierr = PetscMalloc(numFields * sizeof(PetscDualSpace), &sp);CHKERRQ(ierr);
118137b16e26SMatthew G. Knepley   for (f = 0; f < numFields; ++f) {
118237b16e26SMatthew G. Knepley     ierr = PetscFEGetDualSpace(fe[f], &sp[f]);CHKERRQ(ierr);
118337b16e26SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr);
118437b16e26SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
118537b16e26SMatthew G. Knepley     totDim += spDim*numComp;
118637b16e26SMatthew G. Knepley   }
118737b16e26SMatthew G. Knepley   ierr = DMDAGetInfo(dm, &dim,0,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
118837b16e26SMatthew G. Knepley   ierr = DMDAGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
118937b16e26SMatthew G. Knepley   ierr = DMDAVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr);
119037b16e26SMatthew G. Knepley   if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section cell closure size %d != dual space dimension %d", numValues, totDim);
119137b16e26SMatthew G. Knepley   ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
119237b16e26SMatthew G. Knepley   ierr = PetscMalloc3(dim*q.numPoints,PetscReal,&v0,dim*dim*q.numPoints,PetscReal,&J,q.numPoints,PetscReal,&detJ);CHKERRQ(ierr);
119337b16e26SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
119437b16e26SMatthew G. Knepley     PetscCellGeometry geom;
119537b16e26SMatthew G. Knepley 
119637b16e26SMatthew G. Knepley     ierr = DMDAComputeCellGeometry(dm, c, &q, v0, J, NULL, detJ);CHKERRQ(ierr);
119737b16e26SMatthew G. Knepley     geom.v0   = v0;
119837b16e26SMatthew G. Knepley     geom.J    = J;
119937b16e26SMatthew G. Knepley     geom.detJ = detJ;
120037b16e26SMatthew G. Knepley     for (f = 0, v = 0; f < numFields; ++f) {
120137b16e26SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr);
120237b16e26SMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
120337b16e26SMatthew G. Knepley       for (d = 0; d < spDim; ++d) {
120437b16e26SMatthew G. Knepley         ierr = PetscDualSpaceApply(sp[f], d, geom, numComp, funcs[f], &values[v]);CHKERRQ(ierr);
120537b16e26SMatthew G. Knepley         v += numComp;
120637b16e26SMatthew G. Knepley       }
120737b16e26SMatthew G. Knepley     }
120837b16e26SMatthew G. Knepley     ierr = DMDAVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr);
120937b16e26SMatthew G. Knepley   }
121037b16e26SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
121137b16e26SMatthew G. Knepley   ierr = PetscFree3(v0,J,detJ);CHKERRQ(ierr);
121237b16e26SMatthew G. Knepley   ierr = PetscFree(sp);CHKERRQ(ierr);
121337b16e26SMatthew G. Knepley   PetscFunctionReturn(0);
121437b16e26SMatthew G. Knepley }
121537b16e26SMatthew G. Knepley 
121637b16e26SMatthew G. Knepley #undef __FUNCT__
121737b16e26SMatthew G. Knepley #define __FUNCT__ "DMDAProjectFunction"
121837b16e26SMatthew G. Knepley /*@C
121937b16e26SMatthew G. Knepley   DMDAProjectFunction - This projects the given function into the function space provided.
122037b16e26SMatthew G. Knepley 
122137b16e26SMatthew G. Knepley   Input Parameters:
122237b16e26SMatthew G. Knepley + dm      - The DM
122337b16e26SMatthew G. Knepley . fe      - The PetscFE associated with the field
122437b16e26SMatthew G. Knepley . funcs   - The coordinate functions to evaluate
122537b16e26SMatthew G. Knepley - mode    - The insertion mode for values
122637b16e26SMatthew G. Knepley 
122737b16e26SMatthew G. Knepley   Output Parameter:
122837b16e26SMatthew G. Knepley . X - vector
122937b16e26SMatthew G. Knepley 
123037b16e26SMatthew G. Knepley   Level: developer
123137b16e26SMatthew G. Knepley 
123237b16e26SMatthew G. Knepley .seealso: DMDAComputeL2Diff()
123337b16e26SMatthew G. Knepley @*/
123437b16e26SMatthew G. Knepley PetscErrorCode DMDAProjectFunction(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *), InsertMode mode, Vec X)
123537b16e26SMatthew G. Knepley {
123637b16e26SMatthew G. Knepley   Vec            localX;
123737b16e26SMatthew G. Knepley   PetscErrorCode ierr;
123837b16e26SMatthew G. Knepley 
123937b16e26SMatthew G. Knepley   PetscFunctionBegin;
124037b16e26SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
124137b16e26SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
124237b16e26SMatthew G. Knepley   ierr = DMDAProjectFunctionLocal(dm, fe, funcs, mode, localX);CHKERRQ(ierr);
124337b16e26SMatthew G. Knepley   ierr = DMLocalToGlobalBegin(dm, localX, mode, X);CHKERRQ(ierr);
124437b16e26SMatthew G. Knepley   ierr = DMLocalToGlobalEnd(dm, localX, mode, X);CHKERRQ(ierr);
124537b16e26SMatthew G. Knepley   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
124637b16e26SMatthew G. Knepley   PetscFunctionReturn(0);
124737b16e26SMatthew G. Knepley }
124837b16e26SMatthew G. Knepley 
124937b16e26SMatthew G. Knepley #undef __FUNCT__
125037b16e26SMatthew G. Knepley #define __FUNCT__ "DMDAComputeL2Diff"
125137b16e26SMatthew G. Knepley /*@C
125237b16e26SMatthew G. Knepley   DMDAComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h.
125337b16e26SMatthew G. Knepley 
125437b16e26SMatthew G. Knepley   Input Parameters:
125537b16e26SMatthew G. Knepley + dm    - The DM
125637b16e26SMatthew G. Knepley . fe    - The PetscFE object for each field
125737b16e26SMatthew G. Knepley . funcs - The functions to evaluate for each field component
125837b16e26SMatthew G. Knepley - X     - The coefficient vector u_h
125937b16e26SMatthew G. Knepley 
126037b16e26SMatthew G. Knepley   Output Parameter:
126137b16e26SMatthew G. Knepley . diff - The diff ||u - u_h||_2
126237b16e26SMatthew G. Knepley 
126337b16e26SMatthew G. Knepley   Level: developer
126437b16e26SMatthew G. Knepley 
126537b16e26SMatthew G. Knepley .seealso: DMDAProjectFunction()
126637b16e26SMatthew G. Knepley @*/
126737b16e26SMatthew G. Knepley PetscErrorCode DMDAComputeL2Diff(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *), Vec X, PetscReal *diff)
126837b16e26SMatthew G. Knepley {
126937b16e26SMatthew G. Knepley   const PetscInt  debug = 0;
127037b16e26SMatthew G. Knepley   PetscSection    section;
127137b16e26SMatthew G. Knepley   PetscQuadrature quad;
127237b16e26SMatthew G. Knepley   Vec             localX;
127337b16e26SMatthew G. Knepley   PetscScalar    *funcVal;
127437b16e26SMatthew G. Knepley   PetscReal      *coords, *v0, *J, *invJ, detJ;
127537b16e26SMatthew G. Knepley   PetscReal       localDiff = 0.0;
127637b16e26SMatthew G. Knepley   PetscInt        dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp;
127737b16e26SMatthew G. Knepley   PetscErrorCode  ierr;
127837b16e26SMatthew G. Knepley 
127937b16e26SMatthew G. Knepley   PetscFunctionBegin;
128037b16e26SMatthew G. Knepley   ierr = DMDAGetInfo(dm, &dim,0,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
128137b16e26SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
128237b16e26SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
128337b16e26SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
128437b16e26SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
128537b16e26SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
128637b16e26SMatthew G. Knepley   for (field = 0; field < numFields; ++field) {
128737b16e26SMatthew G. Knepley     PetscInt Nc;
128837b16e26SMatthew G. Knepley 
128937b16e26SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr);
129037b16e26SMatthew G. Knepley     numComponents += Nc;
129137b16e26SMatthew G. Knepley   }
129237b16e26SMatthew G. Knepley   /* There are no BC values in DAs right now: ierr = DMDAProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */
129337b16e26SMatthew G. Knepley   ierr = PetscMalloc5(numComponents,PetscScalar,&funcVal,dim,PetscReal,&coords,dim,PetscReal,&v0,dim*dim,PetscReal,&J,dim*dim,PetscReal,&invJ);CHKERRQ(ierr);
129437b16e26SMatthew G. Knepley   ierr = DMDAGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
129537b16e26SMatthew G. Knepley   ierr = PetscFEGetQuadrature(fe[0], &quad);CHKERRQ(ierr);
129637b16e26SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
1297*df66330eSJed Brown     const PetscScalar *x = NULL;
129837b16e26SMatthew G. Knepley     PetscReal    elemDiff = 0.0;
129937b16e26SMatthew G. Knepley 
130037b16e26SMatthew G. Knepley     ierr = DMDAComputeCellGeometry(dm, c, &quad, v0, J, invJ, &detJ);CHKERRQ(ierr);
130137b16e26SMatthew G. Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
130237b16e26SMatthew G. Knepley     ierr = DMDAVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
130337b16e26SMatthew G. Knepley 
130437b16e26SMatthew G. Knepley     for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) {
130537b16e26SMatthew G. Knepley       const PetscInt   numQuadPoints = quad.numPoints;
130637b16e26SMatthew G. Knepley       const PetscReal *quadPoints    = quad.points;
130737b16e26SMatthew G. Knepley       const PetscReal *quadWeights   = quad.weights;
130837b16e26SMatthew G. Knepley       PetscReal       *basis;
130937b16e26SMatthew G. Knepley       PetscInt         numBasisFuncs, numBasisComps, q, d, e, fc, f;
131037b16e26SMatthew G. Knepley 
131137b16e26SMatthew G. Knepley       ierr = PetscFEGetDimension(fe[field], &numBasisFuncs);CHKERRQ(ierr);
131237b16e26SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe[field], &numBasisComps);CHKERRQ(ierr);
131337b16e26SMatthew G. Knepley       ierr = PetscFEGetDefaultTabulation(fe[field], &basis, NULL, NULL);CHKERRQ(ierr);
131437b16e26SMatthew G. Knepley       if (debug) {
131537b16e26SMatthew G. Knepley         char title[1024];
131637b16e26SMatthew G. Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
131737b16e26SMatthew G. Knepley         ierr = DMPrintCellVector(c, title, numBasisFuncs*numBasisComps, &x[fieldOffset]);CHKERRQ(ierr);
131837b16e26SMatthew G. Knepley       }
131937b16e26SMatthew G. Knepley       for (q = 0; q < numQuadPoints; ++q) {
132037b16e26SMatthew G. Knepley         for (d = 0; d < dim; d++) {
132137b16e26SMatthew G. Knepley           coords[d] = v0[d];
132237b16e26SMatthew G. Knepley           for (e = 0; e < dim; e++) {
132337b16e26SMatthew G. Knepley             coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0);
132437b16e26SMatthew G. Knepley           }
132537b16e26SMatthew G. Knepley         }
132637b16e26SMatthew G. Knepley         (*funcs[field])(coords, funcVal);
132737b16e26SMatthew G. Knepley         for (fc = 0; fc < numBasisComps; ++fc) {
132837b16e26SMatthew G. Knepley           PetscScalar interpolant = 0.0;
132937b16e26SMatthew G. Knepley 
133037b16e26SMatthew G. Knepley           for (f = 0; f < numBasisFuncs; ++f) {
133137b16e26SMatthew G. Knepley             const PetscInt fidx = f*numBasisComps+fc;
133237b16e26SMatthew G. Knepley             interpolant += x[fieldOffset+fidx]*basis[q*numBasisFuncs*numBasisComps+fidx];
133337b16e26SMatthew G. Knepley           }
133437b16e26SMatthew G. Knepley           if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "    elem %d field %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ);CHKERRQ(ierr);}
133537b16e26SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ;
133637b16e26SMatthew G. Knepley         }
133737b16e26SMatthew G. Knepley       }
133837b16e26SMatthew G. Knepley       comp        += numBasisComps;
133937b16e26SMatthew G. Knepley       fieldOffset += numBasisFuncs*numBasisComps;
134037b16e26SMatthew G. Knepley     }
134137b16e26SMatthew G. Knepley     ierr = DMDAVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
134237b16e26SMatthew G. Knepley     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
134337b16e26SMatthew G. Knepley     localDiff += elemDiff;
134437b16e26SMatthew G. Knepley   }
134537b16e26SMatthew G. Knepley   ierr  = PetscFree5(funcVal,coords,v0,J,invJ);CHKERRQ(ierr);
134637b16e26SMatthew G. Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
134737b16e26SMatthew G. Knepley   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
134837b16e26SMatthew G. Knepley   *diff = PetscSqrtReal(*diff);
1349a66d4d66SMatthew G Knepley   PetscFunctionReturn(0);
1350a66d4d66SMatthew G Knepley }
1351a66d4d66SMatthew G Knepley 
135247c6ae99SBarry Smith /* ------------------------------------------------------------------- */
135347c6ae99SBarry Smith 
135447c6ae99SBarry Smith #undef __FUNCT__
1355aa219208SBarry Smith #define __FUNCT__ "DMDAGetArray"
135647c6ae99SBarry Smith /*@C
1357aa219208SBarry Smith      DMDAGetArray - Gets a work array for a DMDA
135847c6ae99SBarry Smith 
135947c6ae99SBarry Smith     Input Parameter:
136047c6ae99SBarry Smith +    da - information about my local patch
136147c6ae99SBarry Smith -    ghosted - do you want arrays for the ghosted or nonghosted patch
136247c6ae99SBarry Smith 
136347c6ae99SBarry Smith     Output Parameters:
136447c6ae99SBarry Smith .    vptr - array data structured
136547c6ae99SBarry Smith 
136647c6ae99SBarry Smith     Note:  The vector values are NOT initialized and may have garbage in them, so you may need
136747c6ae99SBarry Smith            to zero them.
136847c6ae99SBarry Smith 
136947c6ae99SBarry Smith   Level: advanced
137047c6ae99SBarry Smith 
1371bcaeba4dSBarry Smith .seealso: DMDARestoreArray()
137247c6ae99SBarry Smith 
137347c6ae99SBarry Smith @*/
13747087cfbeSBarry Smith PetscErrorCode  DMDAGetArray(DM da,PetscBool ghosted,void *vptr)
137547c6ae99SBarry Smith {
137647c6ae99SBarry Smith   PetscErrorCode ierr;
137747c6ae99SBarry Smith   PetscInt       j,i,xs,ys,xm,ym,zs,zm;
137847c6ae99SBarry Smith   char           *iarray_start;
137947c6ae99SBarry Smith   void           **iptr = (void**)vptr;
138047c6ae99SBarry Smith   DM_DA          *dd    = (DM_DA*)da->data;
138147c6ae99SBarry Smith 
138247c6ae99SBarry Smith   PetscFunctionBegin;
138347c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
138447c6ae99SBarry Smith   if (ghosted) {
1385aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
138647c6ae99SBarry Smith       if (dd->arrayghostedin[i]) {
138747c6ae99SBarry Smith         *iptr                 = dd->arrayghostedin[i];
138847c6ae99SBarry Smith         iarray_start          = (char*)dd->startghostedin[i];
13890298fd71SBarry Smith         dd->arrayghostedin[i] = NULL;
13900298fd71SBarry Smith         dd->startghostedin[i] = NULL;
139147c6ae99SBarry Smith 
139247c6ae99SBarry Smith         goto done;
139347c6ae99SBarry Smith       }
139447c6ae99SBarry Smith     }
139547c6ae99SBarry Smith     xs = dd->Xs;
139647c6ae99SBarry Smith     ys = dd->Ys;
139747c6ae99SBarry Smith     zs = dd->Zs;
139847c6ae99SBarry Smith     xm = dd->Xe-dd->Xs;
139947c6ae99SBarry Smith     ym = dd->Ye-dd->Ys;
140047c6ae99SBarry Smith     zm = dd->Ze-dd->Zs;
140147c6ae99SBarry Smith   } else {
1402aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
140347c6ae99SBarry Smith       if (dd->arrayin[i]) {
140447c6ae99SBarry Smith         *iptr          = dd->arrayin[i];
140547c6ae99SBarry Smith         iarray_start   = (char*)dd->startin[i];
14060298fd71SBarry Smith         dd->arrayin[i] = NULL;
14070298fd71SBarry Smith         dd->startin[i] = NULL;
140847c6ae99SBarry Smith 
140947c6ae99SBarry Smith         goto done;
141047c6ae99SBarry Smith       }
141147c6ae99SBarry Smith     }
141247c6ae99SBarry Smith     xs = dd->xs;
141347c6ae99SBarry Smith     ys = dd->ys;
141447c6ae99SBarry Smith     zs = dd->zs;
141547c6ae99SBarry Smith     xm = dd->xe-dd->xs;
141647c6ae99SBarry Smith     ym = dd->ye-dd->ys;
141747c6ae99SBarry Smith     zm = dd->ze-dd->zs;
141847c6ae99SBarry Smith   }
141947c6ae99SBarry Smith 
142047c6ae99SBarry Smith   switch (dd->dim) {
142147c6ae99SBarry Smith   case 1: {
142247c6ae99SBarry Smith     void *ptr;
142347c6ae99SBarry Smith 
142447c6ae99SBarry Smith     ierr = PetscMalloc(xm*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
142547c6ae99SBarry Smith 
142647c6ae99SBarry Smith     ptr   = (void*)(iarray_start - xs*sizeof(PetscScalar));
142747c6ae99SBarry Smith     *iptr = (void*)ptr;
14288865f1eaSKarl Rupp     break;
14298865f1eaSKarl Rupp   }
143047c6ae99SBarry Smith   case 2: {
143147c6ae99SBarry Smith     void **ptr;
143247c6ae99SBarry Smith 
143347c6ae99SBarry Smith     ierr = PetscMalloc((ym+1)*sizeof(void*)+xm*ym*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
143447c6ae99SBarry Smith 
143547c6ae99SBarry Smith     ptr = (void**)(iarray_start + xm*ym*sizeof(PetscScalar) - ys*sizeof(void*));
14368865f1eaSKarl Rupp     for (j=ys; j<ys+ym; j++) ptr[j] = iarray_start + sizeof(PetscScalar)*(xm*(j-ys) - xs);
143747c6ae99SBarry Smith     *iptr = (void*)ptr;
14388865f1eaSKarl Rupp     break;
14398865f1eaSKarl Rupp   }
144047c6ae99SBarry Smith   case 3: {
144147c6ae99SBarry Smith     void ***ptr,**bptr;
144247c6ae99SBarry Smith 
144347c6ae99SBarry Smith     ierr = PetscMalloc((zm+1)*sizeof(void**)+(ym*zm+1)*sizeof(void*)+xm*ym*zm*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
144447c6ae99SBarry Smith 
144547c6ae99SBarry Smith     ptr  = (void***)(iarray_start + xm*ym*zm*sizeof(PetscScalar) - zs*sizeof(void*));
144647c6ae99SBarry Smith     bptr = (void**)(iarray_start + xm*ym*zm*sizeof(PetscScalar) + zm*sizeof(void**));
14478865f1eaSKarl Rupp     for (i=zs; i<zs+zm; i++) ptr[i] = bptr + ((i-zs)*ym - ys);
144847c6ae99SBarry Smith     for (i=zs; i<zs+zm; i++) {
144947c6ae99SBarry Smith       for (j=ys; j<ys+ym; j++) {
145047c6ae99SBarry Smith         ptr[i][j] = iarray_start + sizeof(PetscScalar)*(xm*ym*(i-zs) + xm*(j-ys) - xs);
145147c6ae99SBarry Smith       }
145247c6ae99SBarry Smith     }
145347c6ae99SBarry Smith     *iptr = (void*)ptr;
14548865f1eaSKarl Rupp     break;
14558865f1eaSKarl Rupp   }
145647c6ae99SBarry Smith   default:
1457ce94432eSBarry Smith     SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Dimension %D not supported",dd->dim);
145847c6ae99SBarry Smith   }
145947c6ae99SBarry Smith 
146047c6ae99SBarry Smith done:
146147c6ae99SBarry Smith   /* add arrays to the checked out list */
146247c6ae99SBarry Smith   if (ghosted) {
1463aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
146447c6ae99SBarry Smith       if (!dd->arrayghostedout[i]) {
146547c6ae99SBarry Smith         dd->arrayghostedout[i] = *iptr;
146647c6ae99SBarry Smith         dd->startghostedout[i] = iarray_start;
146747c6ae99SBarry Smith         break;
146847c6ae99SBarry Smith       }
146947c6ae99SBarry Smith     }
147047c6ae99SBarry Smith   } else {
1471aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
147247c6ae99SBarry Smith       if (!dd->arrayout[i]) {
147347c6ae99SBarry Smith         dd->arrayout[i] = *iptr;
147447c6ae99SBarry Smith         dd->startout[i] = iarray_start;
147547c6ae99SBarry Smith         break;
147647c6ae99SBarry Smith       }
147747c6ae99SBarry Smith     }
147847c6ae99SBarry Smith   }
147947c6ae99SBarry Smith   PetscFunctionReturn(0);
148047c6ae99SBarry Smith }
148147c6ae99SBarry Smith 
148247c6ae99SBarry Smith #undef __FUNCT__
1483aa219208SBarry Smith #define __FUNCT__ "DMDARestoreArray"
148447c6ae99SBarry Smith /*@C
1485aa219208SBarry Smith      DMDARestoreArray - Restores an array of derivative types for a DMDA
148647c6ae99SBarry Smith 
148747c6ae99SBarry Smith     Input Parameter:
148847c6ae99SBarry Smith +    da - information about my local patch
148947c6ae99SBarry Smith .    ghosted - do you want arrays for the ghosted or nonghosted patch
149047c6ae99SBarry Smith -    vptr - array data structured to be passed to ad_FormFunctionLocal()
149147c6ae99SBarry Smith 
149247c6ae99SBarry Smith      Level: advanced
149347c6ae99SBarry Smith 
1494bcaeba4dSBarry Smith .seealso: DMDAGetArray()
149547c6ae99SBarry Smith 
149647c6ae99SBarry Smith @*/
14977087cfbeSBarry Smith PetscErrorCode  DMDARestoreArray(DM da,PetscBool ghosted,void *vptr)
149847c6ae99SBarry Smith {
149947c6ae99SBarry Smith   PetscInt i;
150047c6ae99SBarry Smith   void     **iptr = (void**)vptr,*iarray_start = 0;
150147c6ae99SBarry Smith   DM_DA    *dd    = (DM_DA*)da->data;
150247c6ae99SBarry Smith 
150347c6ae99SBarry Smith   PetscFunctionBegin;
150447c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
150547c6ae99SBarry Smith   if (ghosted) {
1506aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
150747c6ae99SBarry Smith       if (dd->arrayghostedout[i] == *iptr) {
150847c6ae99SBarry Smith         iarray_start           = dd->startghostedout[i];
15090298fd71SBarry Smith         dd->arrayghostedout[i] = NULL;
15100298fd71SBarry Smith         dd->startghostedout[i] = NULL;
151147c6ae99SBarry Smith         break;
151247c6ae99SBarry Smith       }
151347c6ae99SBarry Smith     }
1514aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
151547c6ae99SBarry Smith       if (!dd->arrayghostedin[i]) {
151647c6ae99SBarry Smith         dd->arrayghostedin[i] = *iptr;
151747c6ae99SBarry Smith         dd->startghostedin[i] = iarray_start;
151847c6ae99SBarry Smith         break;
151947c6ae99SBarry Smith       }
152047c6ae99SBarry Smith     }
152147c6ae99SBarry Smith   } else {
1522aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
152347c6ae99SBarry Smith       if (dd->arrayout[i] == *iptr) {
152447c6ae99SBarry Smith         iarray_start    = dd->startout[i];
15250298fd71SBarry Smith         dd->arrayout[i] = NULL;
15260298fd71SBarry Smith         dd->startout[i] = NULL;
152747c6ae99SBarry Smith         break;
152847c6ae99SBarry Smith       }
152947c6ae99SBarry Smith     }
1530aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
153147c6ae99SBarry Smith       if (!dd->arrayin[i]) {
153247c6ae99SBarry Smith         dd->arrayin[i] = *iptr;
153347c6ae99SBarry Smith         dd->startin[i] = iarray_start;
153447c6ae99SBarry Smith         break;
153547c6ae99SBarry Smith       }
153647c6ae99SBarry Smith     }
153747c6ae99SBarry Smith   }
153847c6ae99SBarry Smith   PetscFunctionReturn(0);
153947c6ae99SBarry Smith }
154047c6ae99SBarry Smith 
1541