xref: /petsc/src/dm/impls/da/dalocal.c (revision d089267096544a113110f12e19de0de1dfbc8169)
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>
947c6ae99SBarry Smith 
1047c6ae99SBarry Smith /*
11e3c5b3baSBarry Smith    This allows the DMDA vectors to properly tell MATLAB their dimensions
1247c6ae99SBarry Smith */
1347c6ae99SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
14c6db04a5SJed Brown #include <engine.h>   /* MATLAB include file */
15c6db04a5SJed Brown #include <mex.h>      /* MATLAB include file */
1647c6ae99SBarry Smith #undef __FUNCT__
1747c6ae99SBarry Smith #define __FUNCT__ "VecMatlabEnginePut_DA2d"
181e6b0712SBarry Smith static PetscErrorCode  VecMatlabEnginePut_DA2d(PetscObject obj,void *mengine)
1947c6ae99SBarry Smith {
2047c6ae99SBarry Smith   PetscErrorCode ierr;
2147c6ae99SBarry Smith   PetscInt       n,m;
2247c6ae99SBarry Smith   Vec            vec = (Vec)obj;
2347c6ae99SBarry Smith   PetscScalar    *array;
2447c6ae99SBarry Smith   mxArray        *mat;
259a42bb27SBarry Smith   DM             da;
2647c6ae99SBarry Smith 
2747c6ae99SBarry Smith   PetscFunctionBegin;
28c688c046SMatthew G Knepley   ierr = VecGetDM(vec, &da);CHKERRQ(ierr);
29ce94432eSBarry Smith   if (!da) SETERRQ(PetscObjectComm((PetscObject)vec),PETSC_ERR_ARG_WRONGSTATE,"Vector not associated with a DMDA");
30aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,0,0,0,&m,&n,0);CHKERRQ(ierr);
3147c6ae99SBarry Smith 
3247c6ae99SBarry Smith   ierr = VecGetArray(vec,&array);CHKERRQ(ierr);
3347c6ae99SBarry Smith #if !defined(PETSC_USE_COMPLEX)
3447c6ae99SBarry Smith   mat = mxCreateDoubleMatrix(m,n,mxREAL);
3547c6ae99SBarry Smith #else
3647c6ae99SBarry Smith   mat = mxCreateDoubleMatrix(m,n,mxCOMPLEX);
3747c6ae99SBarry Smith #endif
3847c6ae99SBarry Smith   ierr = PetscMemcpy(mxGetPr(mat),array,n*m*sizeof(PetscScalar));CHKERRQ(ierr);
3947c6ae99SBarry Smith   ierr = PetscObjectName(obj);CHKERRQ(ierr);
4047c6ae99SBarry Smith   engPutVariable((Engine*)mengine,obj->name,mat);
4147c6ae99SBarry Smith 
4247c6ae99SBarry Smith   ierr = VecRestoreArray(vec,&array);CHKERRQ(ierr);
4347c6ae99SBarry Smith   PetscFunctionReturn(0);
4447c6ae99SBarry Smith }
4547c6ae99SBarry Smith #endif
4647c6ae99SBarry Smith 
4747c6ae99SBarry Smith 
4847c6ae99SBarry Smith #undef __FUNCT__
49564755cdSBarry Smith #define __FUNCT__ "DMCreateLocalVector_DA"
507087cfbeSBarry Smith PetscErrorCode  DMCreateLocalVector_DA(DM da,Vec *g)
5147c6ae99SBarry Smith {
5247c6ae99SBarry Smith   PetscErrorCode ierr;
5347c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
5447c6ae99SBarry Smith 
5547c6ae99SBarry Smith   PetscFunctionBegin;
5647c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
5747c6ae99SBarry Smith   PetscValidPointer(g,2);
5811689aebSJed Brown   if (da->defaultSection) {
5911689aebSJed Brown     ierr = DMCreateLocalVector_Section_Private(da,g);CHKERRQ(ierr);
6011689aebSJed Brown   } else {
6147c6ae99SBarry Smith     ierr = VecCreate(PETSC_COMM_SELF,g);CHKERRQ(ierr);
6247c6ae99SBarry Smith     ierr = VecSetSizes(*g,dd->nlocal,PETSC_DETERMINE);CHKERRQ(ierr);
6347c6ae99SBarry Smith     ierr = VecSetBlockSize(*g,dd->w);CHKERRQ(ierr);
64401ddaa8SBarry Smith     ierr = VecSetType(*g,da->vectype);CHKERRQ(ierr);
65c688c046SMatthew G Knepley     ierr = VecSetDM(*g, da);CHKERRQ(ierr);
6647c6ae99SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
6747c6ae99SBarry Smith     if (dd->w == 1  && dd->dim == 2) {
68bdf89e91SBarry Smith       ierr = PetscObjectComposeFunction((PetscObject)*g,"PetscMatlabEnginePut_C",VecMatlabEnginePut_DA2d);CHKERRQ(ierr);
6947c6ae99SBarry Smith     }
7047c6ae99SBarry Smith #endif
7111689aebSJed Brown   }
7247c6ae99SBarry Smith   PetscFunctionReturn(0);
7347c6ae99SBarry Smith }
7447c6ae99SBarry Smith 
75a66d4d66SMatthew G Knepley #undef __FUNCT__
7657459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetNumCells"
77e42e3c58SMatthew G. Knepley PetscErrorCode DMDAGetNumCells(DM dm, PetscInt *numCellsX, PetscInt *numCellsY, PetscInt *numCellsZ, PetscInt *numCells)
7857459e9aSMatthew G Knepley {
7957459e9aSMatthew G Knepley   DM_DA          *da = (DM_DA*) dm->data;
8057459e9aSMatthew G Knepley   const PetscInt dim = da->dim;
8157459e9aSMatthew G Knepley   const PetscInt mx  = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs;
8257459e9aSMatthew G Knepley   const PetscInt nC  = (mx)*(dim > 1 ? (my)*(dim > 2 ? (mz) : 1) : 1);
8357459e9aSMatthew G Knepley 
8457459e9aSMatthew G Knepley   PetscFunctionBegin;
85e42e3c58SMatthew G. Knepley   if (numCellsX) {
86e42e3c58SMatthew G. Knepley     PetscValidIntPointer(numCellsX,2);
87e42e3c58SMatthew G. Knepley     *numCellsX = mx;
88e42e3c58SMatthew G. Knepley   }
89e42e3c58SMatthew G. Knepley   if (numCellsY) {
90e42e3c58SMatthew G. Knepley     PetscValidIntPointer(numCellsX,3);
91e42e3c58SMatthew G. Knepley     *numCellsY = my;
92e42e3c58SMatthew G. Knepley   }
93e42e3c58SMatthew G. Knepley   if (numCellsZ) {
94e42e3c58SMatthew G. Knepley     PetscValidIntPointer(numCellsX,4);
95e42e3c58SMatthew G. Knepley     *numCellsZ = mz;
96e42e3c58SMatthew G. Knepley   }
9757459e9aSMatthew G Knepley   if (numCells) {
98e42e3c58SMatthew G. Knepley     PetscValidIntPointer(numCells,5);
9957459e9aSMatthew G Knepley     *numCells = nC;
10057459e9aSMatthew G Knepley   }
10157459e9aSMatthew G Knepley   PetscFunctionReturn(0);
10257459e9aSMatthew G Knepley }
10357459e9aSMatthew G Knepley 
10457459e9aSMatthew G Knepley #undef __FUNCT__
10557459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetNumVertices"
10657459e9aSMatthew G Knepley PetscErrorCode DMDAGetNumVertices(DM dm, PetscInt *numVerticesX, PetscInt *numVerticesY, PetscInt *numVerticesZ, PetscInt *numVertices)
10757459e9aSMatthew G Knepley {
10857459e9aSMatthew G Knepley   DM_DA          *da = (DM_DA*) dm->data;
10957459e9aSMatthew G Knepley   const PetscInt dim = da->dim;
11057459e9aSMatthew G Knepley   const PetscInt mx  = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs;
11157459e9aSMatthew G Knepley   const PetscInt nVx = mx+1;
11257459e9aSMatthew G Knepley   const PetscInt nVy = dim > 1 ? (my+1) : 1;
11357459e9aSMatthew G Knepley   const PetscInt nVz = dim > 2 ? (mz+1) : 1;
11457459e9aSMatthew G Knepley   const PetscInt nV  = nVx*nVy*nVz;
11557459e9aSMatthew G Knepley 
11657459e9aSMatthew G Knepley   PetscFunctionBegin;
11757459e9aSMatthew G Knepley   if (numVerticesX) {
11857459e9aSMatthew G Knepley     PetscValidIntPointer(numVerticesX,2);
11957459e9aSMatthew G Knepley     *numVerticesX = nVx;
12057459e9aSMatthew G Knepley   }
12157459e9aSMatthew G Knepley   if (numVerticesY) {
12257459e9aSMatthew G Knepley     PetscValidIntPointer(numVerticesY,3);
12357459e9aSMatthew G Knepley     *numVerticesY = nVy;
12457459e9aSMatthew G Knepley   }
12557459e9aSMatthew G Knepley   if (numVerticesZ) {
12657459e9aSMatthew G Knepley     PetscValidIntPointer(numVerticesZ,4);
12757459e9aSMatthew G Knepley     *numVerticesZ = nVz;
12857459e9aSMatthew G Knepley   }
12957459e9aSMatthew G Knepley   if (numVertices) {
13057459e9aSMatthew G Knepley     PetscValidIntPointer(numVertices,5);
13157459e9aSMatthew G Knepley     *numVertices = nV;
13257459e9aSMatthew G Knepley   }
13357459e9aSMatthew G Knepley   PetscFunctionReturn(0);
13457459e9aSMatthew G Knepley }
13557459e9aSMatthew G Knepley 
13657459e9aSMatthew G Knepley #undef __FUNCT__
13757459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetNumFaces"
13857459e9aSMatthew G Knepley PetscErrorCode DMDAGetNumFaces(DM dm, PetscInt *numXFacesX, PetscInt *numXFaces, PetscInt *numYFacesY, PetscInt *numYFaces, PetscInt *numZFacesZ, PetscInt *numZFaces)
13957459e9aSMatthew G Knepley {
14057459e9aSMatthew G Knepley   DM_DA          *da = (DM_DA*) dm->data;
14157459e9aSMatthew G Knepley   const PetscInt dim = da->dim;
14257459e9aSMatthew G Knepley   const PetscInt mx  = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs;
14357459e9aSMatthew G Knepley   const PetscInt nxF = (dim > 1 ? (my)*(dim > 2 ? (mz) : 1) : 1);
14457459e9aSMatthew G Knepley   const PetscInt nXF = (mx+1)*nxF;
14557459e9aSMatthew G Knepley   const PetscInt nyF = mx*(dim > 2 ? mz : 1);
14657459e9aSMatthew G Knepley   const PetscInt nYF = dim > 1 ? (my+1)*nyF : 0;
14757459e9aSMatthew G Knepley   const PetscInt nzF = mx*(dim > 1 ? my : 0);
14857459e9aSMatthew G Knepley   const PetscInt nZF = dim > 2 ? (mz+1)*nzF : 0;
14957459e9aSMatthew G Knepley 
15057459e9aSMatthew G Knepley   PetscFunctionBegin;
15157459e9aSMatthew G Knepley   if (numXFacesX) {
15257459e9aSMatthew G Knepley     PetscValidIntPointer(numXFacesX,2);
15357459e9aSMatthew G Knepley     *numXFacesX = nxF;
15457459e9aSMatthew G Knepley   }
15557459e9aSMatthew G Knepley   if (numXFaces) {
15657459e9aSMatthew G Knepley     PetscValidIntPointer(numXFaces,3);
15757459e9aSMatthew G Knepley     *numXFaces = nXF;
15857459e9aSMatthew G Knepley   }
15957459e9aSMatthew G Knepley   if (numYFacesY) {
16057459e9aSMatthew G Knepley     PetscValidIntPointer(numYFacesY,4);
16157459e9aSMatthew G Knepley     *numYFacesY = nyF;
16257459e9aSMatthew G Knepley   }
16357459e9aSMatthew G Knepley   if (numYFaces) {
16457459e9aSMatthew G Knepley     PetscValidIntPointer(numYFaces,5);
16557459e9aSMatthew G Knepley     *numYFaces = nYF;
16657459e9aSMatthew G Knepley   }
16757459e9aSMatthew G Knepley   if (numZFacesZ) {
16857459e9aSMatthew G Knepley     PetscValidIntPointer(numZFacesZ,6);
16957459e9aSMatthew G Knepley     *numZFacesZ = nzF;
17057459e9aSMatthew G Knepley   }
17157459e9aSMatthew G Knepley   if (numZFaces) {
17257459e9aSMatthew G Knepley     PetscValidIntPointer(numZFaces,7);
17357459e9aSMatthew G Knepley     *numZFaces = nZF;
17457459e9aSMatthew G Knepley   }
17557459e9aSMatthew G Knepley   PetscFunctionReturn(0);
17657459e9aSMatthew G Knepley }
17757459e9aSMatthew G Knepley 
17857459e9aSMatthew G Knepley #undef __FUNCT__
17957459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetHeightStratum"
18057459e9aSMatthew G Knepley PetscErrorCode DMDAGetHeightStratum(DM dm, PetscInt height, PetscInt *pStart, PetscInt *pEnd)
18157459e9aSMatthew G Knepley {
18257459e9aSMatthew G Knepley   DM_DA          *da = (DM_DA*) dm->data;
18357459e9aSMatthew G Knepley   const PetscInt dim = da->dim;
18457459e9aSMatthew G Knepley   PetscInt       nC, nV, nXF, nYF, nZF;
18557459e9aSMatthew G Knepley   PetscErrorCode ierr;
18657459e9aSMatthew G Knepley 
18757459e9aSMatthew G Knepley   PetscFunctionBegin;
1888865f1eaSKarl Rupp   if (pStart) PetscValidIntPointer(pStart,3);
1898865f1eaSKarl Rupp   if (pEnd)   PetscValidIntPointer(pEnd,4);
190e42e3c58SMatthew G. Knepley   ierr = DMDAGetNumCells(dm, NULL, NULL, NULL, &nC);CHKERRQ(ierr);
1910298fd71SBarry Smith   ierr = DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV);CHKERRQ(ierr);
1920298fd71SBarry Smith   ierr = DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF);CHKERRQ(ierr);
19357459e9aSMatthew G Knepley   if (height == 0) {
19457459e9aSMatthew G Knepley     /* Cells */
1958865f1eaSKarl Rupp     if (pStart) *pStart = 0;
1968865f1eaSKarl Rupp     if (pEnd)   *pEnd   = nC;
19757459e9aSMatthew G Knepley   } else if (height == 1) {
19857459e9aSMatthew G Knepley     /* Faces */
1998865f1eaSKarl Rupp     if (pStart) *pStart = nC+nV;
2008865f1eaSKarl Rupp     if (pEnd)   *pEnd   = nC+nV+nXF+nYF+nZF;
20157459e9aSMatthew G Knepley   } else if (height == dim) {
20257459e9aSMatthew G Knepley     /* Vertices */
2038865f1eaSKarl Rupp     if (pStart) *pStart = nC;
2048865f1eaSKarl Rupp     if (pEnd)   *pEnd   = nC+nV;
20557459e9aSMatthew G Knepley   } else if (height < 0) {
20657459e9aSMatthew G Knepley     /* All points */
2078865f1eaSKarl Rupp     if (pStart) *pStart = 0;
2088865f1eaSKarl Rupp     if (pEnd)   *pEnd   = nC+nV+nXF+nYF+nZF;
20982f516ccSBarry Smith   } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No points of height %d in the DA", height);
21057459e9aSMatthew G Knepley   PetscFunctionReturn(0);
21157459e9aSMatthew G Knepley }
21257459e9aSMatthew G Knepley 
21357459e9aSMatthew G Knepley #undef __FUNCT__
214*d0892670SMatthew G. Knepley #define __FUNCT__ "DMDAGetDepthStratum"
215*d0892670SMatthew G. Knepley PetscErrorCode DMDAGetDepthStratum(DM dm, PetscInt depth, PetscInt *pStart, PetscInt *pEnd)
216*d0892670SMatthew G. Knepley {
217*d0892670SMatthew G. Knepley   DM_DA         *da  = (DM_DA*) dm->data;
218*d0892670SMatthew G. Knepley   const PetscInt dim = da->dim;
219*d0892670SMatthew G. Knepley   PetscInt       nC, nV, nXF, nYF, nZF;
220*d0892670SMatthew G. Knepley   PetscErrorCode ierr;
221*d0892670SMatthew G. Knepley 
222*d0892670SMatthew G. Knepley   PetscFunctionBegin;
223*d0892670SMatthew G. Knepley   if (pStart) PetscValidIntPointer(pStart,3);
224*d0892670SMatthew G. Knepley   if (pEnd)   PetscValidIntPointer(pEnd,4);
225*d0892670SMatthew G. Knepley   ierr = DMDAGetNumCells(dm, NULL, NULL, NULL, &nC);CHKERRQ(ierr);
226*d0892670SMatthew G. Knepley   ierr = DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV);CHKERRQ(ierr);
227*d0892670SMatthew G. Knepley   ierr = DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF);CHKERRQ(ierr);
228*d0892670SMatthew G. Knepley   if (depth == dim) {
229*d0892670SMatthew G. Knepley     /* Cells */
230*d0892670SMatthew G. Knepley     if (pStart) *pStart = 0;
231*d0892670SMatthew G. Knepley     if (pEnd)   *pEnd   = nC;
232*d0892670SMatthew G. Knepley   } else if (depth == dim-1) {
233*d0892670SMatthew G. Knepley     /* Faces */
234*d0892670SMatthew G. Knepley     if (pStart) *pStart = nC+nV;
235*d0892670SMatthew G. Knepley     if (pEnd)   *pEnd   = nC+nV+nXF+nYF+nZF;
236*d0892670SMatthew G. Knepley   } else if (depth == 0) {
237*d0892670SMatthew G. Knepley     /* Vertices */
238*d0892670SMatthew G. Knepley     if (pStart) *pStart = nC;
239*d0892670SMatthew G. Knepley     if (pEnd)   *pEnd   = nC+nV;
240*d0892670SMatthew G. Knepley   } else if (depth < 0) {
241*d0892670SMatthew G. Knepley     /* All points */
242*d0892670SMatthew G. Knepley     if (pStart) *pStart = 0;
243*d0892670SMatthew G. Knepley     if (pEnd)   *pEnd   = nC+nV+nXF+nYF+nZF;
244*d0892670SMatthew G. Knepley   } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No points of depth %d in the DA", depth);
245*d0892670SMatthew G. Knepley   PetscFunctionReturn(0);
246*d0892670SMatthew G. Knepley }
247*d0892670SMatthew G. Knepley 
248*d0892670SMatthew G. Knepley #undef __FUNCT__
249*d0892670SMatthew G. Knepley #define __FUNCT__ "DMDAGetConeSize"
250*d0892670SMatthew G. Knepley PetscErrorCode DMDAGetConeSize(DM dm, PetscInt p, PetscInt *coneSize)
251*d0892670SMatthew G. Knepley {
252*d0892670SMatthew G. Knepley   DM_DA         *da  = (DM_DA*) dm->data;
253*d0892670SMatthew G. Knepley   const PetscInt dim = da->dim;
254*d0892670SMatthew G. Knepley   PetscInt       nC, nV, nXF, nYF, nZF;
255*d0892670SMatthew G. Knepley   PetscErrorCode ierr;
256*d0892670SMatthew G. Knepley 
257*d0892670SMatthew G. Knepley   PetscFunctionBegin;
258*d0892670SMatthew G. Knepley   *coneSize = 0;
259*d0892670SMatthew G. Knepley   ierr = DMDAGetNumCells(dm, NULL, NULL, NULL, &nC);CHKERRQ(ierr);
260*d0892670SMatthew G. Knepley   ierr = DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV);CHKERRQ(ierr);
261*d0892670SMatthew G. Knepley   ierr = DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF);CHKERRQ(ierr);
262*d0892670SMatthew G. Knepley   switch (dim) {
263*d0892670SMatthew G. Knepley   case 2:
264*d0892670SMatthew G. Knepley     if (p >= 0) {
265*d0892670SMatthew G. Knepley       if (p < nC) {
266*d0892670SMatthew G. Knepley         *coneSize = 4;
267*d0892670SMatthew G. Knepley       } else if (p < nC+nV) {
268*d0892670SMatthew G. Knepley         *coneSize = 0;
269*d0892670SMatthew G. Knepley       } else if (p < nC+nV+nXF+nYF+nZF) {
270*d0892670SMatthew G. Knepley         *coneSize = 2;
271*d0892670SMatthew G. Knepley       } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %d should be in [0, %d)", p, nC+nV+nXF+nYF+nZF);
272*d0892670SMatthew G. Knepley     } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative point %d is invalid", p);
273*d0892670SMatthew G. Knepley     break;
274*d0892670SMatthew G. Knepley   case 3:
275*d0892670SMatthew G. Knepley     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do 3D");
276*d0892670SMatthew G. Knepley     break;
277*d0892670SMatthew G. Knepley   }
278*d0892670SMatthew G. Knepley   PetscFunctionReturn(0);
279*d0892670SMatthew G. Knepley }
280*d0892670SMatthew G. Knepley 
281*d0892670SMatthew G. Knepley #undef __FUNCT__
282*d0892670SMatthew G. Knepley #define __FUNCT__ "DMDAGetCone"
283*d0892670SMatthew G. Knepley PetscErrorCode DMDAGetCone(DM dm, PetscInt p, PetscInt *cone[])
284*d0892670SMatthew G. Knepley {
285*d0892670SMatthew G. Knepley   DM_DA         *da  = (DM_DA*) dm->data;
286*d0892670SMatthew G. Knepley   const PetscInt dim = da->dim;
287*d0892670SMatthew G. Knepley   PetscInt       nCx, nCy, nCz, nC, nVx, nVy, nVz, nV, nxF, nyF, nzF, nXF, nYF, nZF;
288*d0892670SMatthew G. Knepley   PetscErrorCode ierr;
289*d0892670SMatthew G. Knepley 
290*d0892670SMatthew G. Knepley   PetscFunctionBegin;
291*d0892670SMatthew G. Knepley   if (!cone) {ierr = DMGetWorkArray(dm, 6, PETSC_INT, cone);CHKERRQ(ierr);}
292*d0892670SMatthew G. Knepley   ierr = DMDAGetNumCells(dm, &nCx, &nCy, &nCz, &nC);CHKERRQ(ierr);
293*d0892670SMatthew G. Knepley   ierr = DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV);CHKERRQ(ierr);
294*d0892670SMatthew G. Knepley   ierr = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr);
295*d0892670SMatthew G. Knepley   switch (dim) {
296*d0892670SMatthew G. Knepley   case 2:
297*d0892670SMatthew G. Knepley     if (p >= 0) {
298*d0892670SMatthew G. Knepley       if (p < nC) {
299*d0892670SMatthew G. Knepley         const PetscInt cy = p / nCx;
300*d0892670SMatthew G. Knepley         const PetscInt cx = p % nCx;
301*d0892670SMatthew G. Knepley 
302*d0892670SMatthew G. Knepley         (*cone)[0] = cy*nxF + cx + nC+nV;
303*d0892670SMatthew G. Knepley         (*cone)[1] = cx*nyF + cy + nyF + nC+nV+nXF;
304*d0892670SMatthew G. Knepley         (*cone)[2] = cy*nxF + cx + nxF + nC+nV;
305*d0892670SMatthew G. Knepley         (*cone)[3] = cx*nyF + cy + nC+nV+nXF;
306*d0892670SMatthew G. Knepley         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do cell cones");
307*d0892670SMatthew G. Knepley       } else if (p < nC+nV) {
308*d0892670SMatthew G. Knepley       } else if (p < nC+nV+nXF) {
309*d0892670SMatthew G. Knepley         const PetscInt fy = (p - nC+nV) / nxF;
310*d0892670SMatthew G. Knepley         const PetscInt fx = (p - nC+nV) % nxF;
311*d0892670SMatthew G. Knepley 
312*d0892670SMatthew G. Knepley         (*cone)[0] = fy*nVx + fx + nC;
313*d0892670SMatthew G. Knepley         (*cone)[1] = fy*nVx + fx + 1 + nC;
314*d0892670SMatthew G. Knepley       } else if (p < nC+nV+nXF+nYF) {
315*d0892670SMatthew G. Knepley         const PetscInt fx = (p - nC+nV+nXF) / nyF;
316*d0892670SMatthew G. Knepley         const PetscInt fy = (p - nC+nV+nXF) % nyF;
317*d0892670SMatthew G. Knepley 
318*d0892670SMatthew G. Knepley         (*cone)[0] = fy*nVx + fx + nC;
319*d0892670SMatthew G. Knepley         (*cone)[1] = fy*nVx + fx + nVx + nC;
320*d0892670SMatthew G. Knepley       } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %d should be in [0, %d)", p, nC+nV+nXF+nYF+nZF);
321*d0892670SMatthew G. Knepley     } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative point %d is invalid", p);
322*d0892670SMatthew G. Knepley     break;
323*d0892670SMatthew G. Knepley   case 3:
324*d0892670SMatthew G. Knepley     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do 3D");
325*d0892670SMatthew G. Knepley     break;
326*d0892670SMatthew G. Knepley   }
327*d0892670SMatthew G. Knepley   PetscFunctionReturn(0);
328*d0892670SMatthew G. Knepley }
329*d0892670SMatthew G. Knepley 
330*d0892670SMatthew G. Knepley #undef __FUNCT__
331*d0892670SMatthew G. Knepley #define __FUNCT__ "DMDARestoreCone"
332*d0892670SMatthew G. Knepley PetscErrorCode DMDARestoreCone(DM dm, PetscInt p, PetscInt *cone[])
333*d0892670SMatthew G. Knepley {
334*d0892670SMatthew G. Knepley   PetscErrorCode ierr;
335*d0892670SMatthew G. Knepley 
336*d0892670SMatthew G. Knepley   PetscFunctionBegin;
337*d0892670SMatthew G. Knepley   ierr = DMGetWorkArray(dm, 6, PETSC_INT, cone);CHKERRQ(ierr);
338*d0892670SMatthew G. Knepley   PetscFunctionReturn(0);
339*d0892670SMatthew G. Knepley }
340*d0892670SMatthew G. Knepley 
341*d0892670SMatthew G. Knepley #undef __FUNCT__
342a66d4d66SMatthew G Knepley #define __FUNCT__ "DMDACreateSection"
343a66d4d66SMatthew G Knepley /*@C
344a66d4d66SMatthew G Knepley   DMDACreateSection - Create a PetscSection inside the DMDA that describes data layout. This allows multiple fields with
345a66d4d66SMatthew G Knepley   different numbers of dofs on vertices, cells, and faces in each direction.
346a66d4d66SMatthew G Knepley 
347a66d4d66SMatthew G Knepley   Input Parameters:
348a66d4d66SMatthew G Knepley + dm- The DMDA
349a66d4d66SMatthew G Knepley . numFields - The number of fields
3500298fd71SBarry Smith . numComp - The number of components in each field, or NULL for 1
3510298fd71SBarry Smith . numVertexDof - The number of dofs per vertex for each field, or NULL
3520298fd71SBarry Smith . numFaceDof - The number of dofs per face for each field and direction, or NULL
3530298fd71SBarry Smith - numCellDof - The number of dofs per cell for each field, 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 @*/
36880800b1aSMatthew G Knepley PetscErrorCode DMDACreateSection(DM dm, PetscInt numComp[], PetscInt numVertexDof[], PetscInt numFaceDof[], PetscInt numCellDof[])
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) {
4048865f1eaSKarl Rupp     if (numVertexDof) numVertexTotDof  += numVertexDof[f];
4058865f1eaSKarl Rupp     if (numCellDof)   numCellTotDof    += numCellDof[f];
4068865f1eaSKarl Rupp     if (numFaceDof) {
4078865f1eaSKarl Rupp       numFaceTotDof[0] += numFaceDof[f*dim+0];
40888ed4aceSMatthew G Knepley       numFaceTotDof[1] += dim > 1 ? numFaceDof[f*dim+1] : 0;
4090ad7597dSKarl Rupp       numFaceTotDof[2] += dim > 2 ? numFaceDof[f*dim+2] : 0;
4100ad7597dSKarl Rupp     }
411a66d4d66SMatthew G Knepley   }
412a4b60ecfSMatthew G. Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &section);CHKERRQ(ierr);
413a66d4d66SMatthew G Knepley   if (numFields > 1) {
414a4b60ecfSMatthew G. Knepley     ierr = PetscSectionSetNumFields(section, numFields);CHKERRQ(ierr);
415a66d4d66SMatthew G Knepley     for (f = 0; f < numFields; ++f) {
41680800b1aSMatthew G Knepley       const char *name;
41780800b1aSMatthew G Knepley 
41880800b1aSMatthew G Knepley       ierr = DMDAGetFieldName(dm, f, &name);CHKERRQ(ierr);
419a4b60ecfSMatthew G. Knepley       ierr = PetscSectionSetFieldName(section, f, name);CHKERRQ(ierr);
42080800b1aSMatthew G Knepley       if (numComp) {
421a4b60ecfSMatthew G. Knepley         ierr = PetscSectionSetFieldComponents(section, f, numComp[f]);CHKERRQ(ierr);
422a66d4d66SMatthew G Knepley       }
423a66d4d66SMatthew G Knepley     }
424a66d4d66SMatthew G Knepley   } else {
425a66d4d66SMatthew G Knepley     numFields = 0;
426a66d4d66SMatthew G Knepley   }
427a4b60ecfSMatthew G. Knepley   ierr = PetscSectionSetChart(section, pStart, pEnd);CHKERRQ(ierr);
428a66d4d66SMatthew G Knepley   if (numVertexDof) {
429a66d4d66SMatthew G Knepley     for (v = vStart; v < vEnd; ++v) {
430a66d4d66SMatthew G Knepley       for (f = 0; f < numFields; ++f) {
431a4b60ecfSMatthew G. Knepley         ierr = PetscSectionSetFieldDof(section, v, f, numVertexDof[f]);CHKERRQ(ierr);
432a66d4d66SMatthew G Knepley       }
433a4b60ecfSMatthew G. Knepley       ierr = PetscSectionSetDof(section, v, numVertexTotDof);CHKERRQ(ierr);
434a66d4d66SMatthew G Knepley     }
435a66d4d66SMatthew G Knepley   }
436a66d4d66SMatthew G Knepley   if (numFaceDof) {
437a66d4d66SMatthew G Knepley     for (xf = xfStart; xf < xfEnd; ++xf) {
438a66d4d66SMatthew G Knepley       for (f = 0; f < numFields; ++f) {
439a4b60ecfSMatthew G. Knepley         ierr = PetscSectionSetFieldDof(section, xf, f, numFaceDof[f*dim+0]);CHKERRQ(ierr);
440a66d4d66SMatthew G Knepley       }
441a4b60ecfSMatthew G. Knepley       ierr = PetscSectionSetDof(section, xf, numFaceTotDof[0]);CHKERRQ(ierr);
442a66d4d66SMatthew G Knepley     }
443a66d4d66SMatthew G Knepley     for (yf = yfStart; yf < yfEnd; ++yf) {
444a66d4d66SMatthew G Knepley       for (f = 0; f < numFields; ++f) {
445a4b60ecfSMatthew G. Knepley         ierr = PetscSectionSetFieldDof(section, yf, f, numFaceDof[f*dim+1]);CHKERRQ(ierr);
446a66d4d66SMatthew G Knepley       }
447a4b60ecfSMatthew G. Knepley       ierr = PetscSectionSetDof(section, yf, numFaceTotDof[1]);CHKERRQ(ierr);
448a66d4d66SMatthew G Knepley     }
449a66d4d66SMatthew G Knepley     for (zf = zfStart; zf < zfEnd; ++zf) {
450a66d4d66SMatthew G Knepley       for (f = 0; f < numFields; ++f) {
451a4b60ecfSMatthew G. Knepley         ierr = PetscSectionSetFieldDof(section, zf, f, numFaceDof[f*dim+2]);CHKERRQ(ierr);
452a66d4d66SMatthew G Knepley       }
453a4b60ecfSMatthew G. Knepley       ierr = PetscSectionSetDof(section, zf, numFaceTotDof[2]);CHKERRQ(ierr);
454a66d4d66SMatthew G Knepley     }
455a66d4d66SMatthew G Knepley   }
456a66d4d66SMatthew G Knepley   if (numCellDof) {
457a66d4d66SMatthew G Knepley     for (c = cStart; c < cEnd; ++c) {
458a66d4d66SMatthew G Knepley       for (f = 0; f < numFields; ++f) {
459a4b60ecfSMatthew G. Knepley         ierr = PetscSectionSetFieldDof(section, c, f, numCellDof[f]);CHKERRQ(ierr);
460a66d4d66SMatthew G Knepley       }
461a4b60ecfSMatthew G. Knepley       ierr = PetscSectionSetDof(section, c, numCellTotDof);CHKERRQ(ierr);
462a66d4d66SMatthew G Knepley     }
463a66d4d66SMatthew G Knepley   }
464a4b60ecfSMatthew G. Knepley   ierr = PetscSectionSetUp(section);CHKERRQ(ierr);
46588ed4aceSMatthew G Knepley   /* Create mesh point SF */
466b2fff234SMatthew G. Knepley   ierr = PetscBTCreate(pEnd-pStart, &isLeaf);CHKERRQ(ierr);
46788ed4aceSMatthew G Knepley   ierr = DMDAGetNeighbors(dm, &neighbors);CHKERRQ(ierr);
46888ed4aceSMatthew G Knepley   for (zn = 0; zn < (dim > 2 ? 3 : 1); ++zn) {
46988ed4aceSMatthew G Knepley     for (yn = 0; yn < (dim > 1 ? 3 : 1); ++yn) {
47088ed4aceSMatthew G Knepley       for (xn = 0; xn < 3; ++xn) {
4717128ae9fSMatthew G Knepley         const PetscInt xp       = xn-1, yp = dim > 1 ? yn-1 : 0, zp = dim > 2 ? zn-1 : 0;
47288ed4aceSMatthew G Knepley         const PetscInt neighbor = neighbors[(zn*3+yn)*3+xn];
473b2fff234SMatthew G. Knepley         PetscInt       xv, yv, zv;
47488ed4aceSMatthew G Knepley 
4753814d604SMatthew G Knepley         if (neighbor >= 0 && neighbor < rank) {
476b2fff234SMatthew G. Knepley           if (xp < 0) { /* left */
477b2fff234SMatthew G. Knepley             if (yp < 0) { /* bottom */
478b2fff234SMatthew G. Knepley               if (zp < 0) { /* back */
479b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = (      0*nVy +     0)*nVx +     0 + nC; /* left bottom back vertex */
480b2fff234SMatthew G. Knepley                 if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
481b2fff234SMatthew G. Knepley               } else if (zp > 0) { /* front */
482b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = ((nVz-1)*nVy +     0)*nVx +     0 + nC; /* left bottom front vertex */
483b2fff234SMatthew G. Knepley                 if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
484b2fff234SMatthew G. Knepley               } else {
485b2fff234SMatthew G. Knepley                 for (zv = 0; zv < nVz; ++zv) {
486b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (zv*nVy +     0)*nVx +     0 + nC; /* left bottom vertices */
487b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
488b2fff234SMatthew G. Knepley                 }
489b2fff234SMatthew G. Knepley               }
490b2fff234SMatthew G. Knepley             } else if (yp > 0) { /* top */
491b2fff234SMatthew G. Knepley               if (zp < 0) { /* back */
492b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = (      0*nVy + nVy-1)*nVx +     0 + nC; /* left top back vertex */
493b2fff234SMatthew G. Knepley                 if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
494b2fff234SMatthew G. Knepley               } else if (zp > 0) { /* front */
495b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = ((nVz-1)*nVy + nVy-1)*nVx +     0 + nC; /* left top front vertex */
496b2fff234SMatthew G. Knepley                 if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
497b2fff234SMatthew G. Knepley               } else {
498b2fff234SMatthew G. Knepley                 for (zv = 0; zv < nVz; ++zv) {
499b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (zv*nVy + nVy-1)*nVx +     0 + nC; /* left top vertices */
500b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
501b2fff234SMatthew G. Knepley                 }
502b2fff234SMatthew G. Knepley               }
503b2fff234SMatthew G. Knepley             } else {
504b2fff234SMatthew G. Knepley               if (zp < 0) { /* back */
505b2fff234SMatthew G. Knepley                 for (yv = 0; yv < nVy; ++yv) {
506b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (      0*nVy + yv)*nVx +     0 + nC; /* left back vertices */
507b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
508b2fff234SMatthew G. Knepley                 }
509b2fff234SMatthew G. Knepley               } else if (zp > 0) { /* front */
510b2fff234SMatthew G. Knepley                 for (yv = 0; yv < nVy; ++yv) {
511b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = ((nVz-1)*nVy + yv)*nVx +     0 + nC; /* left front vertices */
512b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
513b2fff234SMatthew G. Knepley                 }
514b2fff234SMatthew G. Knepley               } else {
515b2fff234SMatthew G. Knepley                 for (zv = 0; zv < nVz; ++zv) {
516b2fff234SMatthew G. Knepley                   for (yv = 0; yv < nVy; ++yv) {
517b2fff234SMatthew G. Knepley                     const PetscInt localVertex  = (zv*nVy + yv)*nVx +     0 + nC; /* left vertices */
518b2fff234SMatthew G. Knepley                     if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
519b2fff234SMatthew G. Knepley                   }
520b2fff234SMatthew G. Knepley                 }
521b2fff234SMatthew G. Knepley #if 0
522b2fff234SMatthew G. Knepley                 for (xf = 0; xf < nxF; ++xf) {
523b2fff234SMatthew G. Knepley                   /* THIS IS WRONG */
524b2fff234SMatthew G. Knepley                   const PetscInt localFace  = 0 + nC+nV; /* left faces */
525b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localFace)) ++nleaves;
526b2fff234SMatthew G. Knepley                 }
527b2fff234SMatthew G. Knepley #endif
528b2fff234SMatthew G. Knepley               }
529b2fff234SMatthew G. Knepley             }
530b2fff234SMatthew G. Knepley           } else if (xp > 0) { /* right */
531b2fff234SMatthew G. Knepley             if (yp < 0) { /* bottom */
532b2fff234SMatthew G. Knepley               if (zp < 0) { /* back */
533b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = (      0*nVy +     0)*nVx + nVx-1 + nC; /* right bottom back vertex */
534b2fff234SMatthew G. Knepley                 if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
535b2fff234SMatthew G. Knepley               } else if (zp > 0) { /* front */
536b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = ((nVz-1)*nVy +     0)*nVx + nVx-1 + nC; /* right bottom front vertex */
537b2fff234SMatthew G. Knepley                 if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
538b2fff234SMatthew G. Knepley               } else {
539b2fff234SMatthew G. Knepley                 for (zv = 0; zv < nVz; ++zv) {
540b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (zv*nVy +     0)*nVx + nVx-1 + nC; /* right bottom vertices */
541b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
542b2fff234SMatthew G. Knepley                 }
543b2fff234SMatthew G. Knepley               }
544b2fff234SMatthew G. Knepley             } else if (yp > 0) { /* top */
545b2fff234SMatthew G. Knepley               if (zp < 0) { /* back */
546b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = (      0*nVy + nVy-1)*nVx + nVx-1 + nC; /* right top back vertex */
547b2fff234SMatthew G. Knepley                 if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
548b2fff234SMatthew G. Knepley               } else if (zp > 0) { /* front */
549b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = ((nVz-1)*nVy + nVy-1)*nVx + nVx-1 + nC; /* right top front vertex */
550b2fff234SMatthew G. Knepley                 if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
551b2fff234SMatthew G. Knepley               } else {
552b2fff234SMatthew G. Knepley                 for (zv = 0; zv < nVz; ++zv) {
553b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (zv*nVy + nVy-1)*nVx + nVx-1 + nC; /* right top vertices */
554b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
555b2fff234SMatthew G. Knepley                 }
556b2fff234SMatthew G. Knepley               }
557b2fff234SMatthew G. Knepley             } else {
558b2fff234SMatthew G. Knepley               if (zp < 0) { /* back */
559b2fff234SMatthew G. Knepley                 for (yv = 0; yv < nVy; ++yv) {
560b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (      0*nVy + yv)*nVx + nVx-1 + nC; /* right back vertices */
561b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
562b2fff234SMatthew G. Knepley                 }
563b2fff234SMatthew G. Knepley               } else if (zp > 0) { /* front */
564b2fff234SMatthew G. Knepley                 for (yv = 0; yv < nVy; ++yv) {
565b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = ((nVz-1)*nVy + yv)*nVx + nVx-1 + nC; /* right front vertices */
566b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
567b2fff234SMatthew G. Knepley                 }
568b2fff234SMatthew G. Knepley               } else {
569b2fff234SMatthew G. Knepley                 for (zv = 0; zv < nVz; ++zv) {
570b2fff234SMatthew G. Knepley                   for (yv = 0; yv < nVy; ++yv) {
571b2fff234SMatthew G. Knepley                     const PetscInt localVertex  = (zv*nVy + yv)*nVx + nVx-1 + nC; /* right vertices */
572b2fff234SMatthew G. Knepley                     if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
573b2fff234SMatthew G. Knepley                   }
574b2fff234SMatthew G. Knepley                 }
575b2fff234SMatthew G. Knepley #if 0
576b2fff234SMatthew G. Knepley                 for (xf = 0; xf < nxF; ++xf) {
577b2fff234SMatthew G. Knepley                   /* THIS IS WRONG */
578b2fff234SMatthew G. Knepley                   const PetscInt localFace  = 0 + nC+nV; /* right faces */
579b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localFace)) ++nleaves;
580b2fff234SMatthew G. Knepley                 }
581b2fff234SMatthew G. Knepley #endif
582b2fff234SMatthew G. Knepley               }
583b2fff234SMatthew G. Knepley             }
584b2fff234SMatthew G. Knepley           } else {
585b2fff234SMatthew G. Knepley             if (yp < 0) { /* bottom */
586b2fff234SMatthew G. Knepley               if (zp < 0) { /* back */
587b2fff234SMatthew G. Knepley                 for (xv = 0; xv < nVx; ++xv) {
588b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (      0*nVy +     0)*nVx + xv + nC; /* bottom back vertices */
589b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
590b2fff234SMatthew G. Knepley                 }
591b2fff234SMatthew G. Knepley               } else if (zp > 0) { /* front */
592b2fff234SMatthew G. Knepley                 for (xv = 0; xv < nVx; ++xv) {
593b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = ((nVz-1)*nVy +     0)*nVx + xv + nC; /* bottom front vertices */
594b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
595b2fff234SMatthew G. Knepley                 }
596b2fff234SMatthew G. Knepley               } else {
597b2fff234SMatthew G. Knepley                 for (zv = 0; zv < nVz; ++zv) {
598b2fff234SMatthew G. Knepley                   for (xv = 0; xv < nVx; ++xv) {
599b2fff234SMatthew G. Knepley                     const PetscInt localVertex  = (zv*nVy +     0)*nVx + xv + nC; /* bottom vertices */
600b2fff234SMatthew G. Knepley                     if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
601b2fff234SMatthew G. Knepley                   }
602b2fff234SMatthew G. Knepley                 }
603b2fff234SMatthew G. Knepley #if 0
604b2fff234SMatthew G. Knepley                 for (yf = 0; yf < nyF; ++yf) {
605b2fff234SMatthew G. Knepley                   /* THIS IS WRONG */
606b2fff234SMatthew G. Knepley                   const PetscInt localFace  = 0 + nC+nV; /* bottom faces */
607b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVFace)) ++nleaves;
608b2fff234SMatthew G. Knepley                 }
609b2fff234SMatthew G. Knepley #endif
610b2fff234SMatthew G. Knepley               }
611b2fff234SMatthew G. Knepley             } else if (yp > 0) { /* top */
612b2fff234SMatthew G. Knepley               if (zp < 0) { /* back */
613b2fff234SMatthew G. Knepley                 for (xv = 0; xv < nVx; ++xv) {
614b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (      0*nVy + nVy-1)*nVx + xv + nC; /* top back vertices */
615b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
616b2fff234SMatthew G. Knepley                 }
617b2fff234SMatthew G. Knepley               } else if (zp > 0) { /* front */
618b2fff234SMatthew G. Knepley                 for (xv = 0; xv < nVx; ++xv) {
619b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = ((nVz-1)*nVy + nVy-1)*nVx + xv + nC; /* top front vertices */
620b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
621b2fff234SMatthew G. Knepley                 }
622b2fff234SMatthew G. Knepley               } else {
623b2fff234SMatthew G. Knepley                 for (zv = 0; zv < nVz; ++zv) {
624b2fff234SMatthew G. Knepley                   for (xv = 0; xv < nVx; ++xv) {
625b2fff234SMatthew G. Knepley                     const PetscInt localVertex  = (zv*nVy + nVy-1)*nVx + xv + nC; /* top vertices */
626b2fff234SMatthew G. Knepley                     if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
627b2fff234SMatthew G. Knepley                   }
628b2fff234SMatthew G. Knepley                 }
629b2fff234SMatthew G. Knepley #if 0
630b2fff234SMatthew G. Knepley                 for (yf = 0; yf < nyF; ++yf) {
631b2fff234SMatthew G. Knepley                   /* THIS IS WRONG */
632b2fff234SMatthew G. Knepley                   const PetscInt localFace  = 0 + nC+nV; /* top faces */
633b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVFace)) ++nleaves;
634b2fff234SMatthew G. Knepley                 }
635b2fff234SMatthew G. Knepley #endif
636b2fff234SMatthew G. Knepley               }
637b2fff234SMatthew G. Knepley             } else {
638b2fff234SMatthew G. Knepley               if (zp < 0) { /* back */
639b2fff234SMatthew G. Knepley                 for (yv = 0; yv < nVy; ++yv) {
640b2fff234SMatthew G. Knepley                   for (xv = 0; xv < nVx; ++xv) {
641b2fff234SMatthew G. Knepley                     const PetscInt localVertex  = (      0*nVy + yv)*nVx + xv + nC; /* back vertices */
642b2fff234SMatthew G. Knepley                     if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
643b2fff234SMatthew G. Knepley                   }
644b2fff234SMatthew G. Knepley                 }
645b2fff234SMatthew G. Knepley #if 0
646b2fff234SMatthew G. Knepley                 for (zf = 0; zf < nzF; ++zf) {
647b2fff234SMatthew G. Knepley                   /* THIS IS WRONG */
648b2fff234SMatthew G. Knepley                   const PetscInt localFace  = 0 + nC+nV; /* back faces */
649b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localFace)) ++nleaves;
650b2fff234SMatthew G. Knepley                 }
651b2fff234SMatthew G. Knepley #endif
652b2fff234SMatthew G. Knepley               } else if (zp > 0) { /* front */
653b2fff234SMatthew G. Knepley                 for (yv = 0; yv < nVy; ++yv) {
654b2fff234SMatthew G. Knepley                   for (xv = 0; xv < nVx; ++xv) {
655b2fff234SMatthew G. Knepley                     const PetscInt localVertex  = ((nVz-1)*nVy + yv)*nVx + xv + nC; /* front vertices */
656b2fff234SMatthew G. Knepley                     if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
657b2fff234SMatthew G. Knepley                   }
658b2fff234SMatthew G. Knepley                 }
659b2fff234SMatthew G. Knepley #if 0
660b2fff234SMatthew G. Knepley                 for (zf = 0; zf < nzF; ++zf) {
661b2fff234SMatthew G. Knepley                   /* THIS IS WRONG */
662b2fff234SMatthew G. Knepley                   const PetscInt localFace  = 0 + nC+nV; /* front faces */
663b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localFace)) ++nleaves;
664b2fff234SMatthew G. Knepley                 }
665b2fff234SMatthew G. Knepley #endif
666b2fff234SMatthew G. Knepley               } else {
667b2fff234SMatthew G. Knepley                 /* Nothing is shared from the interior */
66888ed4aceSMatthew G Knepley               }
66988ed4aceSMatthew G Knepley             }
67088ed4aceSMatthew G Knepley           }
67188ed4aceSMatthew G Knepley         }
67288ed4aceSMatthew G Knepley       }
673b2fff234SMatthew G. Knepley     }
674b2fff234SMatthew G. Knepley   }
675b2fff234SMatthew G. Knepley   ierr = PetscBTMemzero(pEnd-pStart, isLeaf);CHKERRQ(ierr);
67688ed4aceSMatthew G Knepley   ierr = PetscMalloc2(nleaves,PetscInt,&localPoints,nleaves,PetscSFNode,&remotePoints);CHKERRQ(ierr);
67788ed4aceSMatthew G Knepley   for (zn = 0; zn < (dim > 2 ? 3 : 1); ++zn) {
67888ed4aceSMatthew G Knepley     for (yn = 0; yn < (dim > 1 ? 3 : 1); ++yn) {
67988ed4aceSMatthew G Knepley       for (xn = 0; xn < 3; ++xn) {
6807128ae9fSMatthew G Knepley         const PetscInt xp       = xn-1, yp = dim > 1 ? yn-1 : 0, zp = dim > 2 ? zn-1 : 0;
68188ed4aceSMatthew G Knepley         const PetscInt neighbor = neighbors[(zn*3+yn)*3+xn];
682f5eeb527SMatthew G Knepley         PetscInt       xv, yv, zv;
68388ed4aceSMatthew G Knepley 
6843814d604SMatthew G Knepley         if (neighbor >= 0 && neighbor < rank) {
68588ed4aceSMatthew G Knepley           if (xp < 0) { /* left */
68688ed4aceSMatthew G Knepley             if (yp < 0) { /* bottom */
68788ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
688b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = (      0*nVy +     0)*nVx +     0 + nC; /* left bottom back vertex */
689f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = ((nVz-1)*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 if (zp > 0) { /* front */
698b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = ((nVz-1)*nVy +     0)*nVx +     0 + nC; /* left bottom front vertex */
699f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = (      0*nVy + nVy-1)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
700f5eeb527SMatthew G Knepley 
701b2fff234SMatthew G. Knepley                 if (!PetscBTLookupSet(isLeaf, localVertex)) {
702f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
703f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
704f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
705f5eeb527SMatthew G Knepley                   ++nL;
706b2fff234SMatthew G. Knepley                 }
70788ed4aceSMatthew G Knepley               } else {
708b2fff234SMatthew G. Knepley                 for (zv = 0; zv < nVz; ++zv) {
709b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (zv*nVy +     0)*nVx +     0 + nC; /* left bottom vertices */
710f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (zv*nVy + nVy-1)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
711f5eeb527SMatthew G Knepley 
712b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) {
713f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
714f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
715f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
716b2fff234SMatthew G. Knepley                     ++nL;
717b2fff234SMatthew G. Knepley                   }
718f5eeb527SMatthew G Knepley                 }
71988ed4aceSMatthew G Knepley               }
72088ed4aceSMatthew G Knepley             } else if (yp > 0) { /* top */
72188ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
722b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = (      0*nVy + nVy-1)*nVx +     0 + nC; /* left top back vertex */
723f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = ((nVz-1)*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 if (zp > 0) { /* front */
732b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = ((nVz-1)*nVy + nVy-1)*nVx +     0 + nC; /* left top front vertex */
733f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = (      0*nVy +     0)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
734f5eeb527SMatthew G Knepley 
735b2fff234SMatthew G. Knepley                 if (!PetscBTLookupSet(isLeaf, localVertex)) {
736f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
737f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
738f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
739f5eeb527SMatthew G Knepley                   ++nL;
740b2fff234SMatthew G. Knepley                 }
74188ed4aceSMatthew G Knepley               } else {
742b2fff234SMatthew G. Knepley                 for (zv = 0; zv < nVz; ++zv) {
743b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (zv*nVy + nVy-1)*nVx +     0 + nC; /* left top vertices */
744f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (zv*nVy +     0)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
745f5eeb527SMatthew G Knepley 
746b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) {
747f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
748f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
749f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
750b2fff234SMatthew G. Knepley                     ++nL;
751b2fff234SMatthew G. Knepley                   }
752f5eeb527SMatthew G Knepley                 }
75388ed4aceSMatthew G Knepley               }
75488ed4aceSMatthew G Knepley             } else {
75588ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
756b2fff234SMatthew G. Knepley                 for (yv = 0; yv < nVy; ++yv) {
757b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (      0*nVy + yv)*nVx +     0 + nC; /* left back vertices */
758f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = ((nVz-1)*nVy + yv)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
759f5eeb527SMatthew G Knepley 
760b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) {
761f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
762f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
763f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
764b2fff234SMatthew G. Knepley                     ++nL;
765b2fff234SMatthew G. Knepley                   }
766f5eeb527SMatthew G Knepley                 }
76788ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
768b2fff234SMatthew G. Knepley                 for (yv = 0; yv < nVy; ++yv) {
769b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = ((nVz-1)*nVy + yv)*nVx +     0 + nC; /* left front vertices */
770f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (      0*nVy + yv)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
771f5eeb527SMatthew G Knepley 
772b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) {
773f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
774f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
775f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
776b2fff234SMatthew G. Knepley                     ++nL;
777b2fff234SMatthew G. Knepley                   }
778f5eeb527SMatthew G Knepley                 }
77988ed4aceSMatthew G Knepley               } else {
780f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv) {
781b2fff234SMatthew G. Knepley                   for (yv = 0; yv < nVy; ++yv) {
782b2fff234SMatthew G. Knepley                     const PetscInt localVertex  = (zv*nVy + yv)*nVx +     0 + nC; /* left vertices */
783f5eeb527SMatthew G Knepley                     const PetscInt remoteVertex = (zv*nVy + yv)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
784f5eeb527SMatthew G Knepley 
785b2fff234SMatthew G. Knepley                     if (!PetscBTLookupSet(isLeaf, localVertex)) {
786f5eeb527SMatthew G Knepley                       localPoints[nL]        = localVertex;
787f5eeb527SMatthew G Knepley                       remotePoints[nL].rank  = neighbor;
788f5eeb527SMatthew G Knepley                       remotePoints[nL].index = remoteVertex;
789b2fff234SMatthew G. Knepley                       ++nL;
790f5eeb527SMatthew G Knepley                     }
791f5eeb527SMatthew G Knepley                   }
792b2fff234SMatthew G. Knepley                 }
793b2fff234SMatthew G. Knepley #if 0
794b2fff234SMatthew G. Knepley                 for (xf = 0; xf < nxF; ++xf) {
795f5eeb527SMatthew G Knepley                   /* THIS IS WRONG */
796b2fff234SMatthew G. Knepley                   const PetscInt localFace  = 0 + nC+nV; /* left faces */
797f5eeb527SMatthew G Knepley                   const PetscInt remoteFace = 0 + nC+nV;
798f5eeb527SMatthew G Knepley 
799b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localFace)) {
800f5eeb527SMatthew G Knepley                     localPoints[nL]        = localFace;
801f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
802f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteFace;
803f5eeb527SMatthew G Knepley                   }
80488ed4aceSMatthew G Knepley                 }
805b2fff234SMatthew G. Knepley #endif
806b2fff234SMatthew G. Knepley               }
80788ed4aceSMatthew G Knepley             }
80888ed4aceSMatthew G Knepley           } else if (xp > 0) { /* right */
80988ed4aceSMatthew G Knepley             if (yp < 0) { /* bottom */
81088ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
811b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = (      0*nVy +     0)*nVx + nVx-1 + nC; /* right bottom back vertex */
812f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = ((nVz-1)*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 if (zp > 0) { /* front */
821b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = ((nVz-1)*nVy +     0)*nVx + nVx-1 + nC; /* right bottom front vertex */
822f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = (      0*nVy + nVy-1)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
823f5eeb527SMatthew G Knepley 
824b2fff234SMatthew G. Knepley                 if (!PetscBTLookupSet(isLeaf, localVertex)) {
825f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
826f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
827f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
828f5eeb527SMatthew G Knepley                   ++nL;
829b2fff234SMatthew G. Knepley                 }
83088ed4aceSMatthew G Knepley               } else {
831b2fff234SMatthew G. Knepley                 nleavesCheck += nVz;
832b2fff234SMatthew G. Knepley                 for (zv = 0; zv < nVz; ++zv) {
833b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (zv*nVy +     0)*nVx + nVx-1 + nC; /* right bottom vertices */
834f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (zv*nVy + nVy-1)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
835f5eeb527SMatthew G Knepley 
836b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) {
837f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
838f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
839f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
840b2fff234SMatthew G. Knepley                     ++nL;
841b2fff234SMatthew G. Knepley                   }
842f5eeb527SMatthew G Knepley                 }
84388ed4aceSMatthew G Knepley               }
84488ed4aceSMatthew G Knepley             } else if (yp > 0) { /* top */
84588ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
846b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = (      0*nVy + nVy-1)*nVx + nVx-1 + nC; /* right top back vertex */
847f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = ((nVz-1)*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 if (zp > 0) { /* front */
856b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = ((nVz-1)*nVy + nVy-1)*nVx + nVx-1 + nC; /* right top front vertex */
857f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = (      0*nVy +     0)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
858f5eeb527SMatthew G Knepley 
859b2fff234SMatthew G. Knepley                 if (!PetscBTLookupSet(isLeaf, localVertex)) {
860f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
861f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
862f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
863f5eeb527SMatthew G Knepley                   ++nL;
864b2fff234SMatthew G. Knepley                 }
86588ed4aceSMatthew G Knepley               } else {
866b2fff234SMatthew G. Knepley                 for (zv = 0; zv < nVz; ++zv) {
867b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (zv*nVy + nVy-1)*nVx + nVx-1 + nC; /* right top vertices */
868f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (zv*nVy +     0)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
869f5eeb527SMatthew G Knepley 
870b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) {
871f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
872f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
873f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
874b2fff234SMatthew G. Knepley                     ++nL;
875b2fff234SMatthew G. Knepley                   }
876f5eeb527SMatthew G Knepley                 }
87788ed4aceSMatthew G Knepley               }
87888ed4aceSMatthew G Knepley             } else {
87988ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
880b2fff234SMatthew G. Knepley                 for (yv = 0; yv < nVy; ++yv) {
881b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (      0*nVy + yv)*nVx + nVx-1 + nC; /* right back vertices */
882f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = ((nVz-1)*nVy + yv)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
883f5eeb527SMatthew G Knepley 
884b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) {
885f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
886f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
887f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
888b2fff234SMatthew G. Knepley                     ++nL;
889b2fff234SMatthew G. Knepley                   }
890f5eeb527SMatthew G Knepley                 }
89188ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
892b2fff234SMatthew G. Knepley                 for (yv = 0; yv < nVy; ++yv) {
893b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = ((nVz-1)*nVy + yv)*nVx + nVx-1 + nC; /* right front vertices */
894f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (      0*nVy + yv)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
895f5eeb527SMatthew G Knepley 
896b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) {
897f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
898f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
899f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
900b2fff234SMatthew G. Knepley                     ++nL;
901b2fff234SMatthew G. Knepley                   }
902f5eeb527SMatthew G Knepley                 }
90388ed4aceSMatthew G Knepley               } else {
904f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv) {
905b2fff234SMatthew G. Knepley                   for (yv = 0; yv < nVy; ++yv) {
906b2fff234SMatthew G. Knepley                     const PetscInt localVertex  = (zv*nVy + yv)*nVx + nVx-1 + nC; /* right vertices */
907f5eeb527SMatthew G Knepley                     const PetscInt remoteVertex = (zv*nVy + yv)*nVx + 0     + nC; /* TODO: Correct this for neighbor sizes */
908f5eeb527SMatthew G Knepley 
909b2fff234SMatthew G. Knepley                     if (!PetscBTLookupSet(isLeaf, localVertex)) {
910f5eeb527SMatthew G Knepley                       localPoints[nL]        = localVertex;
911f5eeb527SMatthew G Knepley                       remotePoints[nL].rank  = neighbor;
912f5eeb527SMatthew G Knepley                       remotePoints[nL].index = remoteVertex;
913b2fff234SMatthew G. Knepley                       ++nL;
914f5eeb527SMatthew G Knepley                     }
915f5eeb527SMatthew G Knepley                   }
916b2fff234SMatthew G. Knepley                 }
917b2fff234SMatthew G. Knepley #if 0
918b2fff234SMatthew G. Knepley                 for (xf = 0; xf < nxF; ++xf) {
919f5eeb527SMatthew G Knepley                   /* THIS IS WRONG */
920b2fff234SMatthew G. Knepley                   const PetscInt localFace  = 0 + nC+nV; /* right faces */
921f5eeb527SMatthew G Knepley                   const PetscInt remoteFace = 0 + nC+nV;
922f5eeb527SMatthew G Knepley 
923b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localFace)) {
924f5eeb527SMatthew G Knepley                     localPoints[nL]        = localFace;
925f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
926f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteFace;
927b2fff234SMatthew G. Knepley                     ++nL;
928f5eeb527SMatthew G Knepley                   }
92988ed4aceSMatthew G Knepley                 }
930b2fff234SMatthew G. Knepley #endif
931b2fff234SMatthew G. Knepley               }
93288ed4aceSMatthew G Knepley             }
93388ed4aceSMatthew G Knepley           } else {
93488ed4aceSMatthew G Knepley             if (yp < 0) { /* bottom */
93588ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
936b2fff234SMatthew G. Knepley                 for (xv = 0; xv < nVx; ++xv) {
937b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (      0*nVy +     0)*nVx + xv + nC; /* bottom back vertices */
938f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = ((nVz-1)*nVy + nVy-1)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
939f5eeb527SMatthew G Knepley 
940b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) {
941f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
942f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
943f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
944b2fff234SMatthew G. Knepley                     ++nL;
945b2fff234SMatthew G. Knepley                   }
946f5eeb527SMatthew G Knepley                 }
94788ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
948b2fff234SMatthew G. Knepley                 for (xv = 0; xv < nVx; ++xv) {
949b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = ((nVz-1)*nVy +     0)*nVx + xv + nC; /* bottom front vertices */
950f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (      0*nVy + nVy-1)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
951f5eeb527SMatthew G Knepley 
952b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) {
953f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
954f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
955f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
956b2fff234SMatthew G. Knepley                     ++nL;
957b2fff234SMatthew G. Knepley                   }
958f5eeb527SMatthew G Knepley                 }
95988ed4aceSMatthew G Knepley               } else {
960f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv) {
961b2fff234SMatthew G. Knepley                   for (xv = 0; xv < nVx; ++xv) {
962b2fff234SMatthew G. Knepley                     const PetscInt localVertex  = (zv*nVy +     0)*nVx + xv + nC; /* bottom vertices */
963f5eeb527SMatthew G Knepley                     const PetscInt remoteVertex = (zv*nVy + nVy-1)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
964f5eeb527SMatthew G Knepley 
965b2fff234SMatthew G. Knepley                     if (!PetscBTLookupSet(isLeaf, localVertex)) {
966f5eeb527SMatthew G Knepley                       localPoints[nL]        = localVertex;
967f5eeb527SMatthew G Knepley                       remotePoints[nL].rank  = neighbor;
968f5eeb527SMatthew G Knepley                       remotePoints[nL].index = remoteVertex;
969b2fff234SMatthew G. Knepley                       ++nL;
970f5eeb527SMatthew G Knepley                     }
971f5eeb527SMatthew G Knepley                   }
972b2fff234SMatthew G. Knepley                 }
973b2fff234SMatthew G. Knepley #if 0
974b2fff234SMatthew G. Knepley                 for (yf = 0; yf < nyF; ++yf) {
975f5eeb527SMatthew G Knepley                   /* THIS IS WRONG */
976b2fff234SMatthew G. Knepley                   const PetscInt localFace  = 0 + nC+nV; /* bottom faces */
977f5eeb527SMatthew G Knepley                   const PetscInt remoteFace = 0 + nC+nV;
978f5eeb527SMatthew G Knepley 
979b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localFace)) {
980f5eeb527SMatthew G Knepley                     localPoints[nL]        = localFace;
981f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
982f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteFace;
983b2fff234SMatthew G. Knepley                     ++nL;
984f5eeb527SMatthew G Knepley                   }
98588ed4aceSMatthew G Knepley                 }
986b2fff234SMatthew G. Knepley #endif
987b2fff234SMatthew G. Knepley               }
98888ed4aceSMatthew G Knepley             } else if (yp > 0) { /* top */
98988ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
990b2fff234SMatthew G. Knepley                 for (xv = 0; xv < nVx; ++xv) {
991b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (      0*nVy + nVy-1)*nVx + xv + nC; /* top back vertices */
992f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = ((nVz-1)*nVy +     0)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
993f5eeb527SMatthew G Knepley 
994b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) {
995f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
996f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
997f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
998b2fff234SMatthew G. Knepley                     ++nL;
999b2fff234SMatthew G. Knepley                   }
1000f5eeb527SMatthew G Knepley                 }
100188ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
1002b2fff234SMatthew G. Knepley                 for (xv = 0; xv < nVx; ++xv) {
1003b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = ((nVz-1)*nVy + nVy-1)*nVx + xv + nC; /* top front vertices */
1004f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (      0*nVy +     0)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
1005f5eeb527SMatthew G Knepley 
1006b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) {
1007f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
1008f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
1009f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
1010b2fff234SMatthew G. Knepley                     ++nL;
1011b2fff234SMatthew G. Knepley                   }
1012f5eeb527SMatthew G Knepley                 }
101388ed4aceSMatthew G Knepley               } else {
1014f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv) {
1015b2fff234SMatthew G. Knepley                   for (xv = 0; xv < nVx; ++xv) {
1016b2fff234SMatthew G. Knepley                     const PetscInt localVertex  = (zv*nVy + nVy-1)*nVx + xv + nC; /* top vertices */
1017f5eeb527SMatthew G Knepley                     const PetscInt remoteVertex = (zv*nVy +     0)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
1018f5eeb527SMatthew G Knepley 
1019b2fff234SMatthew G. Knepley                     if (!PetscBTLookupSet(isLeaf, localVertex)) {
1020f5eeb527SMatthew G Knepley                       localPoints[nL]        = localVertex;
1021f5eeb527SMatthew G Knepley                       remotePoints[nL].rank  = neighbor;
1022f5eeb527SMatthew G Knepley                       remotePoints[nL].index = remoteVertex;
1023b2fff234SMatthew G. Knepley                       ++nL;
1024f5eeb527SMatthew G Knepley                     }
1025f5eeb527SMatthew G Knepley                   }
1026b2fff234SMatthew G. Knepley                 }
1027b2fff234SMatthew G. Knepley #if 0
1028b2fff234SMatthew G. Knepley                 for (yf = 0; yf < nyF; ++yf) {
1029f5eeb527SMatthew G Knepley                   /* THIS IS WRONG */
1030b2fff234SMatthew G. Knepley                   const PetscInt localFace  = 0 + nC+nV; /* top faces */
1031f5eeb527SMatthew G Knepley                   const PetscInt remoteFace = 0 + nC+nV;
1032f5eeb527SMatthew G Knepley 
1033b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localFace)) {
1034f5eeb527SMatthew G Knepley                     localPoints[nL]        = localFace;
1035f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
1036f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteFace;
1037b2fff234SMatthew G. Knepley                     ++nL;
1038f5eeb527SMatthew G Knepley                   }
103988ed4aceSMatthew G Knepley                 }
1040b2fff234SMatthew G. Knepley #endif
1041b2fff234SMatthew G. Knepley               }
104288ed4aceSMatthew G Knepley             } else {
104388ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
1044f5eeb527SMatthew G Knepley                 for (yv = 0; yv < nVy; ++yv) {
1045b2fff234SMatthew G. Knepley                   for (xv = 0; xv < nVx; ++xv) {
1046b2fff234SMatthew G. Knepley                     const PetscInt localVertex  = (      0*nVy + yv)*nVx + xv + nC; /* back vertices */
1047f5eeb527SMatthew G Knepley                     const PetscInt remoteVertex = ((nVz-1)*nVy + yv)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
1048f5eeb527SMatthew G Knepley 
1049b2fff234SMatthew G. Knepley                     if (!PetscBTLookupSet(isLeaf, localVertex)) {
1050f5eeb527SMatthew G Knepley                       localPoints[nL]        = localVertex;
1051f5eeb527SMatthew G Knepley                       remotePoints[nL].rank  = neighbor;
1052f5eeb527SMatthew G Knepley                       remotePoints[nL].index = remoteVertex;
1053b2fff234SMatthew G. Knepley                       ++nL;
1054f5eeb527SMatthew G Knepley                     }
1055f5eeb527SMatthew G Knepley                   }
1056b2fff234SMatthew G. Knepley                 }
1057b2fff234SMatthew G. Knepley #if 0
1058b2fff234SMatthew G. Knepley                 for (zf = 0; zf < nzF; ++zf) {
1059f5eeb527SMatthew G Knepley                   /* THIS IS WRONG */
1060b2fff234SMatthew G. Knepley                   const PetscInt localFace  = 0 + nC+nV; /* back faces */
1061f5eeb527SMatthew G Knepley                   const PetscInt remoteFace = 0 + nC+nV;
1062f5eeb527SMatthew G Knepley 
1063b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localFace)) {
1064f5eeb527SMatthew G Knepley                     localPoints[nL]        = localFace;
1065f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
1066f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteFace;
1067b2fff234SMatthew G. Knepley                     ++nL;
1068f5eeb527SMatthew G Knepley                   }
1069b2fff234SMatthew G. Knepley                 }
1070b2fff234SMatthew G. Knepley #endif
107188ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
1072f5eeb527SMatthew G Knepley                 for (yv = 0; yv < nVy; ++yv) {
1073b2fff234SMatthew G. Knepley                   for (xv = 0; xv < nVx; ++xv) {
1074b2fff234SMatthew G. Knepley                     const PetscInt localVertex  = ((nVz-1)*nVy + yv)*nVx + xv + nC; /* front vertices */
1075f5eeb527SMatthew G Knepley                     const PetscInt remoteVertex = (      0*nVy + yv)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
1076f5eeb527SMatthew G Knepley 
1077b2fff234SMatthew G. Knepley                     if (!PetscBTLookupSet(isLeaf, localVertex)) {
1078f5eeb527SMatthew G Knepley                       localPoints[nL]        = localVertex;
1079f5eeb527SMatthew G Knepley                       remotePoints[nL].rank  = neighbor;
1080f5eeb527SMatthew G Knepley                       remotePoints[nL].index = remoteVertex;
1081b2fff234SMatthew G. Knepley                       ++nL;
1082f5eeb527SMatthew G Knepley                     }
1083f5eeb527SMatthew G Knepley                   }
1084b2fff234SMatthew G. Knepley                 }
1085b2fff234SMatthew G. Knepley #if 0
1086b2fff234SMatthew G. Knepley                 for (zf = 0; zf < nzF; ++zf) {
1087f5eeb527SMatthew G Knepley                   /* THIS IS WRONG */
1088b2fff234SMatthew G. Knepley                   const PetscInt localFace  = 0 + nC+nV; /* front faces */
1089f5eeb527SMatthew G Knepley                   const PetscInt remoteFace = 0 + nC+nV;
1090f5eeb527SMatthew G Knepley 
1091b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localFace)) {
1092f5eeb527SMatthew G Knepley                     localPoints[nL]        = localFace;
1093f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
1094f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteFace;
1095b2fff234SMatthew G. Knepley                     ++nL;
1096f5eeb527SMatthew G Knepley                   }
1097b2fff234SMatthew G. Knepley                 }
1098b2fff234SMatthew G. Knepley #endif
109988ed4aceSMatthew G Knepley               } else {
110088ed4aceSMatthew G Knepley                 /* Nothing is shared from the interior */
110188ed4aceSMatthew G Knepley               }
110288ed4aceSMatthew G Knepley             }
110388ed4aceSMatthew G Knepley           }
110488ed4aceSMatthew G Knepley         }
110588ed4aceSMatthew G Knepley       }
110688ed4aceSMatthew G Knepley     }
110788ed4aceSMatthew G Knepley   }
1108b2fff234SMatthew G. Knepley   ierr = PetscBTDestroy(&isLeaf);CHKERRQ(ierr);
1109b2fff234SMatthew G. Knepley   /* Remove duplication in leaf determination */
1110b2fff234SMatthew 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);
111182f516ccSBarry Smith   ierr = PetscSFCreate(PetscObjectComm((PetscObject)dm), &sf);CHKERRQ(ierr);
111288ed4aceSMatthew G Knepley   ierr = PetscSFSetGraph(sf, pEnd, nleaves, localPoints, PETSC_OWN_POINTER, remotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr);
1113a4b60ecfSMatthew G. Knepley   ierr = DMSetPointSF(dm, sf);CHKERRQ(ierr);
111488ed4aceSMatthew G Knepley   ierr = PetscSFDestroy(&sf);CHKERRQ(ierr);
1115a4b60ecfSMatthew G. Knepley   ierr = DMSetDefaultSection(dm, section);CHKERRQ(ierr);
1116*d0892670SMatthew G. Knepley #undef __FUNCT__
1117*d0892670SMatthew G. Knepley #define __FUNCT__ "DMDASetVertexCoordinates"
1118*d0892670SMatthew G. Knepley PetscErrorCode DMDASetVertexCoordinates(DM dm, PetscReal xl, PetscReal xu, PetscReal yl, PetscReal yu, PetscReal zl, PetscReal zu)
1119*d0892670SMatthew G. Knepley {
1120*d0892670SMatthew G. Knepley   DM_DA         *da = (DM_DA *) dm->data;
1121*d0892670SMatthew G. Knepley   Vec            coordinates;
1122*d0892670SMatthew G. Knepley   PetscSection   section;
1123*d0892670SMatthew G. Knepley   PetscScalar   *coords;
1124*d0892670SMatthew G. Knepley   PetscReal      h[3];
1125*d0892670SMatthew G. Knepley   PetscInt       dim, size, M, N, P, nVx, nVy, nVz, nV, vStart, vEnd, v, i, j, k;
1126*d0892670SMatthew G. Knepley   PetscErrorCode ierr;
1127*d0892670SMatthew G. Knepley 
1128*d0892670SMatthew G. Knepley   PetscFunctionBegin;
1129*d0892670SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1130*d0892670SMatthew G. Knepley   ierr = DMDAGetInfo(dm, &dim, &M, &N, &P, 0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
1131*d0892670SMatthew G. Knepley   h[0] = (xu - xl)/M;
1132*d0892670SMatthew G. Knepley   h[1] = (yu - yl)/N;
1133*d0892670SMatthew G. Knepley   h[2] = (zu - zl)/P;
1134*d0892670SMatthew G. Knepley   ierr = DMDAGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
1135*d0892670SMatthew G. Knepley   ierr = DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV);CHKERRQ(ierr);
1136*d0892670SMatthew G. Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &section);CHKERRQ(ierr);
1137*d0892670SMatthew G. Knepley   ierr = PetscSectionSetNumFields(section, 1);CHKERRQ(ierr);
1138*d0892670SMatthew G. Knepley   ierr = PetscSectionSetFieldComponents(section, 0, dim);CHKERRQ(ierr);
1139*d0892670SMatthew G. Knepley   ierr = PetscSectionSetChart(section, vStart, vEnd);CHKERRQ(ierr);
1140*d0892670SMatthew G. Knepley   for (v = vStart; v < vEnd; ++v) {
1141*d0892670SMatthew G. Knepley     ierr = PetscSectionSetDof(section, v, dim);CHKERRQ(ierr);
1142*d0892670SMatthew G. Knepley   }
1143*d0892670SMatthew G. Knepley   ierr = PetscSectionSetUp(section);CHKERRQ(ierr);
1144*d0892670SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(section, &size);CHKERRQ(ierr);
1145*d0892670SMatthew G. Knepley   ierr = VecCreateSeq(PETSC_COMM_SELF, size, &coordinates);CHKERRQ(ierr);
1146*d0892670SMatthew G. Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
1147*d0892670SMatthew G. Knepley   for (k = 0; k < nVz; ++k) {
1148*d0892670SMatthew G. Knepley     PetscInt ind[3] = {0, 0, k + da->zs}, d, off;
1149*d0892670SMatthew G. Knepley 
1150*d0892670SMatthew G. Knepley     for (j = 0; j < nVy; ++j) {
1151*d0892670SMatthew G. Knepley       ind[1] = j + da->ys;
1152*d0892670SMatthew G. Knepley       for (i = 0; i < nVx; ++i) {
1153*d0892670SMatthew G. Knepley         const PetscInt vertex = (k*nVy + j)*nVx + i + vStart;
1154*d0892670SMatthew G. Knepley 
1155*d0892670SMatthew G. Knepley         ierr = PetscSectionGetOffset(section, vertex, &off);CHKERRQ(ierr);
1156*d0892670SMatthew G. Knepley         ind[0] = i + da->xs;
1157*d0892670SMatthew G. Knepley         for (d = 0; d < dim; ++d) {
1158*d0892670SMatthew G. Knepley           coords[off+d] = h[d]*ind[d];
1159*d0892670SMatthew G. Knepley         }
1160*d0892670SMatthew G. Knepley       }
1161*d0892670SMatthew G. Knepley     }
1162*d0892670SMatthew G. Knepley   }
1163*d0892670SMatthew G. Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
1164*d0892670SMatthew G. Knepley   ierr = DMSetCoordinateSection(dm, section);CHKERRQ(ierr);
1165*d0892670SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
1166a4b60ecfSMatthew G. Knepley   ierr = PetscSectionDestroy(&section);CHKERRQ(ierr);
1167*d0892670SMatthew G. Knepley   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
1168*d0892670SMatthew G. Knepley   PetscFunctionReturn(0);
1169*d0892670SMatthew G. Knepley }
1170a66d4d66SMatthew G Knepley   PetscFunctionReturn(0);
1171a66d4d66SMatthew G Knepley }
1172a66d4d66SMatthew G Knepley 
117347c6ae99SBarry Smith /* ------------------------------------------------------------------- */
117447c6ae99SBarry Smith 
117547c6ae99SBarry Smith #undef __FUNCT__
1176aa219208SBarry Smith #define __FUNCT__ "DMDAGetArray"
117747c6ae99SBarry Smith /*@C
1178aa219208SBarry Smith      DMDAGetArray - Gets a work array for a DMDA
117947c6ae99SBarry Smith 
118047c6ae99SBarry Smith     Input Parameter:
118147c6ae99SBarry Smith +    da - information about my local patch
118247c6ae99SBarry Smith -    ghosted - do you want arrays for the ghosted or nonghosted patch
118347c6ae99SBarry Smith 
118447c6ae99SBarry Smith     Output Parameters:
118547c6ae99SBarry Smith .    vptr - array data structured
118647c6ae99SBarry Smith 
118747c6ae99SBarry Smith     Note:  The vector values are NOT initialized and may have garbage in them, so you may need
118847c6ae99SBarry Smith            to zero them.
118947c6ae99SBarry Smith 
119047c6ae99SBarry Smith   Level: advanced
119147c6ae99SBarry Smith 
1192bcaeba4dSBarry Smith .seealso: DMDARestoreArray()
119347c6ae99SBarry Smith 
119447c6ae99SBarry Smith @*/
11957087cfbeSBarry Smith PetscErrorCode  DMDAGetArray(DM da,PetscBool ghosted,void *vptr)
119647c6ae99SBarry Smith {
119747c6ae99SBarry Smith   PetscErrorCode ierr;
119847c6ae99SBarry Smith   PetscInt       j,i,xs,ys,xm,ym,zs,zm;
119947c6ae99SBarry Smith   char           *iarray_start;
120047c6ae99SBarry Smith   void           **iptr = (void**)vptr;
120147c6ae99SBarry Smith   DM_DA          *dd    = (DM_DA*)da->data;
120247c6ae99SBarry Smith 
120347c6ae99SBarry Smith   PetscFunctionBegin;
120447c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
120547c6ae99SBarry Smith   if (ghosted) {
1206aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
120747c6ae99SBarry Smith       if (dd->arrayghostedin[i]) {
120847c6ae99SBarry Smith         *iptr                 = dd->arrayghostedin[i];
120947c6ae99SBarry Smith         iarray_start          = (char*)dd->startghostedin[i];
12100298fd71SBarry Smith         dd->arrayghostedin[i] = NULL;
12110298fd71SBarry Smith         dd->startghostedin[i] = NULL;
121247c6ae99SBarry Smith 
121347c6ae99SBarry Smith         goto done;
121447c6ae99SBarry Smith       }
121547c6ae99SBarry Smith     }
121647c6ae99SBarry Smith     xs = dd->Xs;
121747c6ae99SBarry Smith     ys = dd->Ys;
121847c6ae99SBarry Smith     zs = dd->Zs;
121947c6ae99SBarry Smith     xm = dd->Xe-dd->Xs;
122047c6ae99SBarry Smith     ym = dd->Ye-dd->Ys;
122147c6ae99SBarry Smith     zm = dd->Ze-dd->Zs;
122247c6ae99SBarry Smith   } else {
1223aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
122447c6ae99SBarry Smith       if (dd->arrayin[i]) {
122547c6ae99SBarry Smith         *iptr          = dd->arrayin[i];
122647c6ae99SBarry Smith         iarray_start   = (char*)dd->startin[i];
12270298fd71SBarry Smith         dd->arrayin[i] = NULL;
12280298fd71SBarry Smith         dd->startin[i] = NULL;
122947c6ae99SBarry Smith 
123047c6ae99SBarry Smith         goto done;
123147c6ae99SBarry Smith       }
123247c6ae99SBarry Smith     }
123347c6ae99SBarry Smith     xs = dd->xs;
123447c6ae99SBarry Smith     ys = dd->ys;
123547c6ae99SBarry Smith     zs = dd->zs;
123647c6ae99SBarry Smith     xm = dd->xe-dd->xs;
123747c6ae99SBarry Smith     ym = dd->ye-dd->ys;
123847c6ae99SBarry Smith     zm = dd->ze-dd->zs;
123947c6ae99SBarry Smith   }
124047c6ae99SBarry Smith 
124147c6ae99SBarry Smith   switch (dd->dim) {
124247c6ae99SBarry Smith   case 1: {
124347c6ae99SBarry Smith     void *ptr;
124447c6ae99SBarry Smith 
124547c6ae99SBarry Smith     ierr = PetscMalloc(xm*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
124647c6ae99SBarry Smith 
124747c6ae99SBarry Smith     ptr   = (void*)(iarray_start - xs*sizeof(PetscScalar));
124847c6ae99SBarry Smith     *iptr = (void*)ptr;
12498865f1eaSKarl Rupp     break;
12508865f1eaSKarl Rupp   }
125147c6ae99SBarry Smith   case 2: {
125247c6ae99SBarry Smith     void **ptr;
125347c6ae99SBarry Smith 
125447c6ae99SBarry Smith     ierr = PetscMalloc((ym+1)*sizeof(void*)+xm*ym*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
125547c6ae99SBarry Smith 
125647c6ae99SBarry Smith     ptr = (void**)(iarray_start + xm*ym*sizeof(PetscScalar) - ys*sizeof(void*));
12578865f1eaSKarl Rupp     for (j=ys; j<ys+ym; j++) ptr[j] = iarray_start + sizeof(PetscScalar)*(xm*(j-ys) - xs);
125847c6ae99SBarry Smith     *iptr = (void*)ptr;
12598865f1eaSKarl Rupp     break;
12608865f1eaSKarl Rupp   }
126147c6ae99SBarry Smith   case 3: {
126247c6ae99SBarry Smith     void ***ptr,**bptr;
126347c6ae99SBarry Smith 
126447c6ae99SBarry Smith     ierr = PetscMalloc((zm+1)*sizeof(void**)+(ym*zm+1)*sizeof(void*)+xm*ym*zm*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
126547c6ae99SBarry Smith 
126647c6ae99SBarry Smith     ptr  = (void***)(iarray_start + xm*ym*zm*sizeof(PetscScalar) - zs*sizeof(void*));
126747c6ae99SBarry Smith     bptr = (void**)(iarray_start + xm*ym*zm*sizeof(PetscScalar) + zm*sizeof(void**));
12688865f1eaSKarl Rupp     for (i=zs; i<zs+zm; i++) ptr[i] = bptr + ((i-zs)*ym - ys);
126947c6ae99SBarry Smith     for (i=zs; i<zs+zm; i++) {
127047c6ae99SBarry Smith       for (j=ys; j<ys+ym; j++) {
127147c6ae99SBarry Smith         ptr[i][j] = iarray_start + sizeof(PetscScalar)*(xm*ym*(i-zs) + xm*(j-ys) - xs);
127247c6ae99SBarry Smith       }
127347c6ae99SBarry Smith     }
127447c6ae99SBarry Smith     *iptr = (void*)ptr;
12758865f1eaSKarl Rupp     break;
12768865f1eaSKarl Rupp   }
127747c6ae99SBarry Smith   default:
1278ce94432eSBarry Smith     SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Dimension %D not supported",dd->dim);
127947c6ae99SBarry Smith   }
128047c6ae99SBarry Smith 
128147c6ae99SBarry Smith done:
128247c6ae99SBarry Smith   /* add arrays to the checked out list */
128347c6ae99SBarry Smith   if (ghosted) {
1284aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
128547c6ae99SBarry Smith       if (!dd->arrayghostedout[i]) {
128647c6ae99SBarry Smith         dd->arrayghostedout[i] = *iptr;
128747c6ae99SBarry Smith         dd->startghostedout[i] = iarray_start;
128847c6ae99SBarry Smith         break;
128947c6ae99SBarry Smith       }
129047c6ae99SBarry Smith     }
129147c6ae99SBarry Smith   } else {
1292aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
129347c6ae99SBarry Smith       if (!dd->arrayout[i]) {
129447c6ae99SBarry Smith         dd->arrayout[i] = *iptr;
129547c6ae99SBarry Smith         dd->startout[i] = iarray_start;
129647c6ae99SBarry Smith         break;
129747c6ae99SBarry Smith       }
129847c6ae99SBarry Smith     }
129947c6ae99SBarry Smith   }
130047c6ae99SBarry Smith   PetscFunctionReturn(0);
130147c6ae99SBarry Smith }
130247c6ae99SBarry Smith 
130347c6ae99SBarry Smith #undef __FUNCT__
1304aa219208SBarry Smith #define __FUNCT__ "DMDARestoreArray"
130547c6ae99SBarry Smith /*@C
1306aa219208SBarry Smith      DMDARestoreArray - Restores an array of derivative types for a DMDA
130747c6ae99SBarry Smith 
130847c6ae99SBarry Smith     Input Parameter:
130947c6ae99SBarry Smith +    da - information about my local patch
131047c6ae99SBarry Smith .    ghosted - do you want arrays for the ghosted or nonghosted patch
131147c6ae99SBarry Smith -    vptr - array data structured to be passed to ad_FormFunctionLocal()
131247c6ae99SBarry Smith 
131347c6ae99SBarry Smith      Level: advanced
131447c6ae99SBarry Smith 
1315bcaeba4dSBarry Smith .seealso: DMDAGetArray()
131647c6ae99SBarry Smith 
131747c6ae99SBarry Smith @*/
13187087cfbeSBarry Smith PetscErrorCode  DMDARestoreArray(DM da,PetscBool ghosted,void *vptr)
131947c6ae99SBarry Smith {
132047c6ae99SBarry Smith   PetscInt i;
132147c6ae99SBarry Smith   void     **iptr = (void**)vptr,*iarray_start = 0;
132247c6ae99SBarry Smith   DM_DA    *dd    = (DM_DA*)da->data;
132347c6ae99SBarry Smith 
132447c6ae99SBarry Smith   PetscFunctionBegin;
132547c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
132647c6ae99SBarry Smith   if (ghosted) {
1327aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
132847c6ae99SBarry Smith       if (dd->arrayghostedout[i] == *iptr) {
132947c6ae99SBarry Smith         iarray_start           = dd->startghostedout[i];
13300298fd71SBarry Smith         dd->arrayghostedout[i] = NULL;
13310298fd71SBarry Smith         dd->startghostedout[i] = NULL;
133247c6ae99SBarry Smith         break;
133347c6ae99SBarry Smith       }
133447c6ae99SBarry Smith     }
1335aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
133647c6ae99SBarry Smith       if (!dd->arrayghostedin[i]) {
133747c6ae99SBarry Smith         dd->arrayghostedin[i] = *iptr;
133847c6ae99SBarry Smith         dd->startghostedin[i] = iarray_start;
133947c6ae99SBarry Smith         break;
134047c6ae99SBarry Smith       }
134147c6ae99SBarry Smith     }
134247c6ae99SBarry Smith   } else {
1343aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
134447c6ae99SBarry Smith       if (dd->arrayout[i] == *iptr) {
134547c6ae99SBarry Smith         iarray_start    = dd->startout[i];
13460298fd71SBarry Smith         dd->arrayout[i] = NULL;
13470298fd71SBarry Smith         dd->startout[i] = NULL;
134847c6ae99SBarry Smith         break;
134947c6ae99SBarry Smith       }
135047c6ae99SBarry Smith     }
1351aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
135247c6ae99SBarry Smith       if (!dd->arrayin[i]) {
135347c6ae99SBarry Smith         dd->arrayin[i] = *iptr;
135447c6ae99SBarry Smith         dd->startin[i] = iarray_start;
135547c6ae99SBarry Smith         break;
135647c6ae99SBarry Smith       }
135747c6ae99SBarry Smith     }
135847c6ae99SBarry Smith   }
135947c6ae99SBarry Smith   PetscFunctionReturn(0);
136047c6ae99SBarry Smith }
136147c6ae99SBarry Smith 
1362