xref: /petsc/src/dm/impls/da/dalocal.c (revision 44e1f9ab22ea38f6440f1f67290e6cbbffb8f0be)
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>
99a800dd8SMatthew G. Knepley #include <petscfe.h>
1047c6ae99SBarry Smith 
1147c6ae99SBarry Smith /*
12e3c5b3baSBarry Smith    This allows the DMDA vectors to properly tell MATLAB their dimensions
1347c6ae99SBarry Smith */
1447c6ae99SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
15c6db04a5SJed Brown #include <engine.h>   /* MATLAB include file */
16c6db04a5SJed Brown #include <mex.h>      /* MATLAB include file */
1747c6ae99SBarry Smith #undef __FUNCT__
1847c6ae99SBarry Smith #define __FUNCT__ "VecMatlabEnginePut_DA2d"
191e6b0712SBarry Smith static PetscErrorCode  VecMatlabEnginePut_DA2d(PetscObject obj,void *mengine)
2047c6ae99SBarry Smith {
2147c6ae99SBarry Smith   PetscErrorCode ierr;
2247c6ae99SBarry Smith   PetscInt       n,m;
2347c6ae99SBarry Smith   Vec            vec = (Vec)obj;
2447c6ae99SBarry Smith   PetscScalar    *array;
2547c6ae99SBarry Smith   mxArray        *mat;
269a42bb27SBarry Smith   DM             da;
2747c6ae99SBarry Smith 
2847c6ae99SBarry Smith   PetscFunctionBegin;
29c688c046SMatthew G Knepley   ierr = VecGetDM(vec, &da);CHKERRQ(ierr);
30ce94432eSBarry Smith   if (!da) SETERRQ(PetscObjectComm((PetscObject)vec),PETSC_ERR_ARG_WRONGSTATE,"Vector not associated with a DMDA");
31aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,0,0,0,&m,&n,0);CHKERRQ(ierr);
3247c6ae99SBarry Smith 
3347c6ae99SBarry Smith   ierr = VecGetArray(vec,&array);CHKERRQ(ierr);
3447c6ae99SBarry Smith #if !defined(PETSC_USE_COMPLEX)
3547c6ae99SBarry Smith   mat = mxCreateDoubleMatrix(m,n,mxREAL);
3647c6ae99SBarry Smith #else
3747c6ae99SBarry Smith   mat = mxCreateDoubleMatrix(m,n,mxCOMPLEX);
3847c6ae99SBarry Smith #endif
3947c6ae99SBarry Smith   ierr = PetscMemcpy(mxGetPr(mat),array,n*m*sizeof(PetscScalar));CHKERRQ(ierr);
4047c6ae99SBarry Smith   ierr = PetscObjectName(obj);CHKERRQ(ierr);
4147c6ae99SBarry Smith   engPutVariable((Engine*)mengine,obj->name,mat);
4247c6ae99SBarry Smith 
4347c6ae99SBarry Smith   ierr = VecRestoreArray(vec,&array);CHKERRQ(ierr);
4447c6ae99SBarry Smith   PetscFunctionReturn(0);
4547c6ae99SBarry Smith }
4647c6ae99SBarry Smith #endif
4747c6ae99SBarry Smith 
4847c6ae99SBarry Smith 
4947c6ae99SBarry Smith #undef __FUNCT__
50564755cdSBarry Smith #define __FUNCT__ "DMCreateLocalVector_DA"
517087cfbeSBarry Smith PetscErrorCode  DMCreateLocalVector_DA(DM da,Vec *g)
5247c6ae99SBarry Smith {
5347c6ae99SBarry Smith   PetscErrorCode ierr;
5447c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
5547c6ae99SBarry Smith 
5647c6ae99SBarry Smith   PetscFunctionBegin;
5747c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
5847c6ae99SBarry Smith   PetscValidPointer(g,2);
5911689aebSJed Brown   if (da->defaultSection) {
6011689aebSJed Brown     ierr = DMCreateLocalVector_Section_Private(da,g);CHKERRQ(ierr);
6111689aebSJed Brown   } else {
6247c6ae99SBarry Smith     ierr = VecCreate(PETSC_COMM_SELF,g);CHKERRQ(ierr);
6347c6ae99SBarry Smith     ierr = VecSetSizes(*g,dd->nlocal,PETSC_DETERMINE);CHKERRQ(ierr);
6447c6ae99SBarry Smith     ierr = VecSetBlockSize(*g,dd->w);CHKERRQ(ierr);
65401ddaa8SBarry Smith     ierr = VecSetType(*g,da->vectype);CHKERRQ(ierr);
66c688c046SMatthew G Knepley     ierr = VecSetDM(*g, da);CHKERRQ(ierr);
6747c6ae99SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
6847c6ae99SBarry Smith     if (dd->w == 1  && dd->dim == 2) {
69bdf89e91SBarry Smith       ierr = PetscObjectComposeFunction((PetscObject)*g,"PetscMatlabEnginePut_C",VecMatlabEnginePut_DA2d);CHKERRQ(ierr);
7047c6ae99SBarry Smith     }
7147c6ae99SBarry Smith #endif
7211689aebSJed Brown   }
7347c6ae99SBarry Smith   PetscFunctionReturn(0);
7447c6ae99SBarry Smith }
7547c6ae99SBarry Smith 
76a66d4d66SMatthew G Knepley #undef __FUNCT__
7757459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetNumCells"
78939f6067SMatthew G. Knepley /*@
79939f6067SMatthew G. Knepley   DMDAGetNumCells - Get the number of cells in the local piece of the DMDA. This includes ghost cells.
80939f6067SMatthew G. Knepley 
81939f6067SMatthew G. Knepley   Input Parameter:
82939f6067SMatthew G. Knepley . dm - The DM object
83939f6067SMatthew G. Knepley 
84939f6067SMatthew G. Knepley   Output Parameters:
85939f6067SMatthew G. Knepley + numCellsX - The number of local cells in the x-direction
86939f6067SMatthew G. Knepley . numCellsY - The number of local cells in the y-direction
87939f6067SMatthew G. Knepley . numCellsZ - The number of local cells in the z-direction
88939f6067SMatthew G. Knepley - numCells - The number of local cells
89939f6067SMatthew G. Knepley 
90939f6067SMatthew G. Knepley   Level: developer
91939f6067SMatthew G. Knepley 
92939f6067SMatthew G. Knepley .seealso: DMDAGetCellPoint()
93939f6067SMatthew G. Knepley @*/
943582350cSMatthew G. Knepley PetscErrorCode DMDAGetNumCells(DM dm, PetscInt *numCellsX, PetscInt *numCellsY, PetscInt *numCellsZ, PetscInt *numCells)
9557459e9aSMatthew G Knepley {
9657459e9aSMatthew G Knepley   DM_DA         *da  = (DM_DA*) dm->data;
9757459e9aSMatthew G Knepley   const PetscInt dim = da->dim;
9857459e9aSMatthew G Knepley   const PetscInt mx  = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs;
9957459e9aSMatthew G Knepley   const PetscInt nC  = (mx)*(dim > 1 ? (my)*(dim > 2 ? (mz) : 1) : 1);
10057459e9aSMatthew G Knepley 
10157459e9aSMatthew G Knepley   PetscFunctionBegin;
102d6401b93SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1033582350cSMatthew G. Knepley   if (numCellsX) {
1043582350cSMatthew G. Knepley     PetscValidIntPointer(numCellsX,2);
1053582350cSMatthew G. Knepley     *numCellsX = mx;
1063582350cSMatthew G. Knepley   }
1073582350cSMatthew G. Knepley   if (numCellsY) {
1083582350cSMatthew G. Knepley     PetscValidIntPointer(numCellsX,3);
1093582350cSMatthew G. Knepley     *numCellsY = my;
1103582350cSMatthew G. Knepley   }
1113582350cSMatthew G. Knepley   if (numCellsZ) {
1123582350cSMatthew G. Knepley     PetscValidIntPointer(numCellsX,4);
1133582350cSMatthew G. Knepley     *numCellsZ = mz;
1143582350cSMatthew G. Knepley   }
11557459e9aSMatthew G Knepley   if (numCells) {
1163582350cSMatthew G. Knepley     PetscValidIntPointer(numCells,5);
11757459e9aSMatthew G Knepley     *numCells = nC;
11857459e9aSMatthew G Knepley   }
11957459e9aSMatthew G Knepley   PetscFunctionReturn(0);
12057459e9aSMatthew G Knepley }
12157459e9aSMatthew G Knepley 
12257459e9aSMatthew G Knepley #undef __FUNCT__
123d6401b93SMatthew G. Knepley #define __FUNCT__ "DMDAGetCellPoint"
124d6401b93SMatthew G. Knepley /*@
125d6401b93SMatthew G. Knepley   DMDAGetCellPoint - Get the DM point corresponding to the tuple (i, j, k) in the DMDA
126d6401b93SMatthew G. Knepley 
127d6401b93SMatthew G. Knepley   Input Parameters:
128d6401b93SMatthew G. Knepley + dm - The DM object
129d6401b93SMatthew G. Knepley - i,j,k - The global indices for the cell
130d6401b93SMatthew G. Knepley 
131d6401b93SMatthew G. Knepley   Output Parameters:
132d6401b93SMatthew G. Knepley . point - The local DM point
133d6401b93SMatthew G. Knepley 
134d6401b93SMatthew G. Knepley   Level: developer
135d6401b93SMatthew G. Knepley 
136d6401b93SMatthew G. Knepley .seealso: DMDAGetNumCells()
137d6401b93SMatthew G. Knepley @*/
138d6401b93SMatthew G. Knepley PetscErrorCode DMDAGetCellPoint(DM dm, PetscInt i, PetscInt j, PetscInt k, PetscInt *point)
139d6401b93SMatthew G. Knepley {
140d6401b93SMatthew G. Knepley   DM_DA         *da  = (DM_DA*) dm->data;
141d6401b93SMatthew G. Knepley   const PetscInt dim = da->dim;
142d6401b93SMatthew G. Knepley   DMDALocalInfo  info;
143d6401b93SMatthew G. Knepley   PetscErrorCode ierr;
144d6401b93SMatthew G. Knepley 
145d6401b93SMatthew G. Knepley   PetscFunctionBegin;
146d6401b93SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
147d6401b93SMatthew G. Knepley   PetscValidIntPointer(point,5);
148d6401b93SMatthew G. Knepley   ierr = DMDAGetLocalInfo(dm, &info);CHKERRQ(ierr);
149d6401b93SMatthew G. Knepley   if (dim > 0) {if ((i < info.gxs) || (i >= info.gxs+info.gxm)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "X index %d not in [%d, %d)", i, info.gxs, info.gxs+info.gxm);}
150d6401b93SMatthew G. Knepley   if (dim > 1) {if ((i < info.gys) || (i >= info.gys+info.gym)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Y index %d not in [%d, %d)", i, info.gys, info.gys+info.gym);}
151d6401b93SMatthew G. Knepley   if (dim > 2) {if ((i < info.gzs) || (i >= info.gzs+info.gzm)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Z index %d not in [%d, %d)", i, info.gzs, info.gzs+info.gzm);}
152d6401b93SMatthew G. Knepley   *point = i + (dim > 1 ? (j + (dim > 2 ? k*info.gym : 0))*info.gxm : 0);
153d6401b93SMatthew G. Knepley   PetscFunctionReturn(0);
154d6401b93SMatthew G. Knepley }
155d6401b93SMatthew G. Knepley 
156d6401b93SMatthew G. Knepley #undef __FUNCT__
15757459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetNumVertices"
15857459e9aSMatthew G Knepley PetscErrorCode DMDAGetNumVertices(DM dm, PetscInt *numVerticesX, PetscInt *numVerticesY, PetscInt *numVerticesZ, PetscInt *numVertices)
15957459e9aSMatthew G Knepley {
16057459e9aSMatthew G Knepley   DM_DA          *da = (DM_DA*) dm->data;
16157459e9aSMatthew G Knepley   const PetscInt dim = da->dim;
16257459e9aSMatthew G Knepley   const PetscInt mx  = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs;
16357459e9aSMatthew G Knepley   const PetscInt nVx = mx+1;
16457459e9aSMatthew G Knepley   const PetscInt nVy = dim > 1 ? (my+1) : 1;
16557459e9aSMatthew G Knepley   const PetscInt nVz = dim > 2 ? (mz+1) : 1;
16657459e9aSMatthew G Knepley   const PetscInt nV  = nVx*nVy*nVz;
16757459e9aSMatthew G Knepley 
16857459e9aSMatthew G Knepley   PetscFunctionBegin;
16957459e9aSMatthew G Knepley   if (numVerticesX) {
17057459e9aSMatthew G Knepley     PetscValidIntPointer(numVerticesX,2);
17157459e9aSMatthew G Knepley     *numVerticesX = nVx;
17257459e9aSMatthew G Knepley   }
17357459e9aSMatthew G Knepley   if (numVerticesY) {
17457459e9aSMatthew G Knepley     PetscValidIntPointer(numVerticesY,3);
17557459e9aSMatthew G Knepley     *numVerticesY = nVy;
17657459e9aSMatthew G Knepley   }
17757459e9aSMatthew G Knepley   if (numVerticesZ) {
17857459e9aSMatthew G Knepley     PetscValidIntPointer(numVerticesZ,4);
17957459e9aSMatthew G Knepley     *numVerticesZ = nVz;
18057459e9aSMatthew G Knepley   }
18157459e9aSMatthew G Knepley   if (numVertices) {
18257459e9aSMatthew G Knepley     PetscValidIntPointer(numVertices,5);
18357459e9aSMatthew G Knepley     *numVertices = nV;
18457459e9aSMatthew G Knepley   }
18557459e9aSMatthew G Knepley   PetscFunctionReturn(0);
18657459e9aSMatthew G Knepley }
18757459e9aSMatthew G Knepley 
18857459e9aSMatthew G Knepley #undef __FUNCT__
18957459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetNumFaces"
19057459e9aSMatthew G Knepley PetscErrorCode DMDAGetNumFaces(DM dm, PetscInt *numXFacesX, PetscInt *numXFaces, PetscInt *numYFacesY, PetscInt *numYFaces, PetscInt *numZFacesZ, PetscInt *numZFaces)
19157459e9aSMatthew G Knepley {
19257459e9aSMatthew G Knepley   DM_DA          *da = (DM_DA*) dm->data;
19357459e9aSMatthew G Knepley   const PetscInt dim = da->dim;
19457459e9aSMatthew G Knepley   const PetscInt mx  = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs;
19557459e9aSMatthew G Knepley   const PetscInt nxF = (dim > 1 ? (my)*(dim > 2 ? (mz) : 1) : 1);
19657459e9aSMatthew G Knepley   const PetscInt nXF = (mx+1)*nxF;
19757459e9aSMatthew G Knepley   const PetscInt nyF = mx*(dim > 2 ? mz : 1);
19857459e9aSMatthew G Knepley   const PetscInt nYF = dim > 1 ? (my+1)*nyF : 0;
19957459e9aSMatthew G Knepley   const PetscInt nzF = mx*(dim > 1 ? my : 0);
20057459e9aSMatthew G Knepley   const PetscInt nZF = dim > 2 ? (mz+1)*nzF : 0;
20157459e9aSMatthew G Knepley 
20257459e9aSMatthew G Knepley   PetscFunctionBegin;
20357459e9aSMatthew G Knepley   if (numXFacesX) {
20457459e9aSMatthew G Knepley     PetscValidIntPointer(numXFacesX,2);
20557459e9aSMatthew G Knepley     *numXFacesX = nxF;
20657459e9aSMatthew G Knepley   }
20757459e9aSMatthew G Knepley   if (numXFaces) {
20857459e9aSMatthew G Knepley     PetscValidIntPointer(numXFaces,3);
20957459e9aSMatthew G Knepley     *numXFaces = nXF;
21057459e9aSMatthew G Knepley   }
21157459e9aSMatthew G Knepley   if (numYFacesY) {
21257459e9aSMatthew G Knepley     PetscValidIntPointer(numYFacesY,4);
21357459e9aSMatthew G Knepley     *numYFacesY = nyF;
21457459e9aSMatthew G Knepley   }
21557459e9aSMatthew G Knepley   if (numYFaces) {
21657459e9aSMatthew G Knepley     PetscValidIntPointer(numYFaces,5);
21757459e9aSMatthew G Knepley     *numYFaces = nYF;
21857459e9aSMatthew G Knepley   }
21957459e9aSMatthew G Knepley   if (numZFacesZ) {
22057459e9aSMatthew G Knepley     PetscValidIntPointer(numZFacesZ,6);
22157459e9aSMatthew G Knepley     *numZFacesZ = nzF;
22257459e9aSMatthew G Knepley   }
22357459e9aSMatthew G Knepley   if (numZFaces) {
22457459e9aSMatthew G Knepley     PetscValidIntPointer(numZFaces,7);
22557459e9aSMatthew G Knepley     *numZFaces = nZF;
22657459e9aSMatthew G Knepley   }
22757459e9aSMatthew G Knepley   PetscFunctionReturn(0);
22857459e9aSMatthew G Knepley }
22957459e9aSMatthew G Knepley 
23057459e9aSMatthew G Knepley #undef __FUNCT__
23157459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetHeightStratum"
23257459e9aSMatthew G Knepley PetscErrorCode DMDAGetHeightStratum(DM dm, PetscInt height, PetscInt *pStart, PetscInt *pEnd)
23357459e9aSMatthew G Knepley {
23457459e9aSMatthew G Knepley   DM_DA          *da = (DM_DA*) dm->data;
23557459e9aSMatthew G Knepley   const PetscInt dim = da->dim;
23657459e9aSMatthew G Knepley   PetscInt       nC, nV, nXF, nYF, nZF;
23757459e9aSMatthew G Knepley   PetscErrorCode ierr;
23857459e9aSMatthew G Knepley 
23957459e9aSMatthew G Knepley   PetscFunctionBegin;
2408865f1eaSKarl Rupp   if (pStart) PetscValidIntPointer(pStart,3);
2418865f1eaSKarl Rupp   if (pEnd)   PetscValidIntPointer(pEnd,4);
2423582350cSMatthew G. Knepley   ierr = DMDAGetNumCells(dm, NULL, NULL, NULL, &nC);CHKERRQ(ierr);
2430298fd71SBarry Smith   ierr = DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV);CHKERRQ(ierr);
2440298fd71SBarry Smith   ierr = DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF);CHKERRQ(ierr);
24557459e9aSMatthew G Knepley   if (height == 0) {
24657459e9aSMatthew G Knepley     /* Cells */
2478865f1eaSKarl Rupp     if (pStart) *pStart = 0;
2488865f1eaSKarl Rupp     if (pEnd)   *pEnd   = nC;
24957459e9aSMatthew G Knepley   } else if (height == 1) {
25057459e9aSMatthew G Knepley     /* Faces */
2518865f1eaSKarl Rupp     if (pStart) *pStart = nC+nV;
2528865f1eaSKarl Rupp     if (pEnd)   *pEnd   = nC+nV+nXF+nYF+nZF;
25357459e9aSMatthew G Knepley   } else if (height == dim) {
25457459e9aSMatthew G Knepley     /* Vertices */
2558865f1eaSKarl Rupp     if (pStart) *pStart = nC;
2568865f1eaSKarl Rupp     if (pEnd)   *pEnd   = nC+nV;
25757459e9aSMatthew G Knepley   } else if (height < 0) {
25857459e9aSMatthew G Knepley     /* All points */
2598865f1eaSKarl Rupp     if (pStart) *pStart = 0;
2608865f1eaSKarl Rupp     if (pEnd)   *pEnd   = nC+nV+nXF+nYF+nZF;
26182f516ccSBarry Smith   } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No points of height %d in the DA", height);
26257459e9aSMatthew G Knepley   PetscFunctionReturn(0);
26357459e9aSMatthew G Knepley }
26457459e9aSMatthew G Knepley 
26557459e9aSMatthew G Knepley #undef __FUNCT__
2663385ff7eSMatthew G. Knepley #define __FUNCT__ "DMDAGetDepthStratum"
2673385ff7eSMatthew G. Knepley PetscErrorCode DMDAGetDepthStratum(DM dm, PetscInt depth, PetscInt *pStart, PetscInt *pEnd)
2683385ff7eSMatthew G. Knepley {
2693385ff7eSMatthew G. Knepley   DM_DA         *da  = (DM_DA*) dm->data;
2703385ff7eSMatthew G. Knepley   const PetscInt dim = da->dim;
2713385ff7eSMatthew G. Knepley   PetscInt       nC, nV, nXF, nYF, nZF;
2723385ff7eSMatthew G. Knepley   PetscErrorCode ierr;
2733385ff7eSMatthew G. Knepley 
2743385ff7eSMatthew G. Knepley   PetscFunctionBegin;
2753385ff7eSMatthew G. Knepley   if (pStart) PetscValidIntPointer(pStart,3);
2763385ff7eSMatthew G. Knepley   if (pEnd)   PetscValidIntPointer(pEnd,4);
2773385ff7eSMatthew G. Knepley   ierr = DMDAGetNumCells(dm, NULL, NULL, NULL, &nC);CHKERRQ(ierr);
2783385ff7eSMatthew G. Knepley   ierr = DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV);CHKERRQ(ierr);
2793385ff7eSMatthew G. Knepley   ierr = DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF);CHKERRQ(ierr);
2803385ff7eSMatthew G. Knepley   if (depth == dim) {
2813385ff7eSMatthew G. Knepley     /* Cells */
2823385ff7eSMatthew G. Knepley     if (pStart) *pStart = 0;
2833385ff7eSMatthew G. Knepley     if (pEnd)   *pEnd   = nC;
2843385ff7eSMatthew G. Knepley   } else if (depth == dim-1) {
2853385ff7eSMatthew G. Knepley     /* Faces */
2863385ff7eSMatthew G. Knepley     if (pStart) *pStart = nC+nV;
2873385ff7eSMatthew G. Knepley     if (pEnd)   *pEnd   = nC+nV+nXF+nYF+nZF;
2883385ff7eSMatthew G. Knepley   } else if (depth == 0) {
2893385ff7eSMatthew G. Knepley     /* Vertices */
2903385ff7eSMatthew G. Knepley     if (pStart) *pStart = nC;
2913385ff7eSMatthew G. Knepley     if (pEnd)   *pEnd   = nC+nV;
2923385ff7eSMatthew G. Knepley   } else if (depth < 0) {
2933385ff7eSMatthew G. Knepley     /* All points */
2943385ff7eSMatthew G. Knepley     if (pStart) *pStart = 0;
2953385ff7eSMatthew G. Knepley     if (pEnd)   *pEnd   = nC+nV+nXF+nYF+nZF;
2963385ff7eSMatthew G. Knepley   } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No points of depth %d in the DA", depth);
2973385ff7eSMatthew G. Knepley   PetscFunctionReturn(0);
2983385ff7eSMatthew G. Knepley }
2993385ff7eSMatthew G. Knepley 
3003385ff7eSMatthew G. Knepley #undef __FUNCT__
3013385ff7eSMatthew G. Knepley #define __FUNCT__ "DMDAGetConeSize"
3023385ff7eSMatthew G. Knepley PetscErrorCode DMDAGetConeSize(DM dm, PetscInt p, PetscInt *coneSize)
3033385ff7eSMatthew G. Knepley {
3043385ff7eSMatthew G. Knepley   DM_DA         *da  = (DM_DA*) dm->data;
3053385ff7eSMatthew G. Knepley   const PetscInt dim = da->dim;
3063385ff7eSMatthew G. Knepley   PetscInt       nC, nV, nXF, nYF, nZF;
3073385ff7eSMatthew G. Knepley   PetscErrorCode ierr;
3083385ff7eSMatthew G. Knepley 
3093385ff7eSMatthew G. Knepley   PetscFunctionBegin;
3103385ff7eSMatthew G. Knepley   *coneSize = 0;
3113385ff7eSMatthew G. Knepley   ierr = DMDAGetNumCells(dm, NULL, NULL, NULL, &nC);CHKERRQ(ierr);
3123385ff7eSMatthew G. Knepley   ierr = DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV);CHKERRQ(ierr);
3133385ff7eSMatthew G. Knepley   ierr = DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF);CHKERRQ(ierr);
3143385ff7eSMatthew G. Knepley   switch (dim) {
3153385ff7eSMatthew G. Knepley   case 2:
3163385ff7eSMatthew G. Knepley     if (p >= 0) {
3173385ff7eSMatthew G. Knepley       if (p < nC) {
3183385ff7eSMatthew G. Knepley         *coneSize = 4;
3193385ff7eSMatthew G. Knepley       } else if (p < nC+nV) {
3203385ff7eSMatthew G. Knepley         *coneSize = 0;
3213385ff7eSMatthew G. Knepley       } else if (p < nC+nV+nXF+nYF+nZF) {
3223385ff7eSMatthew G. Knepley         *coneSize = 2;
3233385ff7eSMatthew G. Knepley       } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %d should be in [0, %d)", p, nC+nV+nXF+nYF+nZF);
3243385ff7eSMatthew G. Knepley     } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative point %d is invalid", p);
3253385ff7eSMatthew G. Knepley     break;
3263385ff7eSMatthew G. Knepley   case 3:
3273385ff7eSMatthew G. Knepley     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do 3D");
3283385ff7eSMatthew G. Knepley     break;
3293385ff7eSMatthew G. Knepley   }
3303385ff7eSMatthew G. Knepley   PetscFunctionReturn(0);
3313385ff7eSMatthew G. Knepley }
3323385ff7eSMatthew G. Knepley 
3333385ff7eSMatthew G. Knepley #undef __FUNCT__
3343385ff7eSMatthew G. Knepley #define __FUNCT__ "DMDAGetCone"
3353385ff7eSMatthew G. Knepley PetscErrorCode DMDAGetCone(DM dm, PetscInt p, PetscInt *cone[])
3363385ff7eSMatthew G. Knepley {
3373385ff7eSMatthew G. Knepley   DM_DA         *da  = (DM_DA*) dm->data;
3383385ff7eSMatthew G. Knepley   const PetscInt dim = da->dim;
3393385ff7eSMatthew G. Knepley   PetscInt       nCx, nCy, nCz, nC, nVx, nVy, nVz, nV, nxF, nyF, nzF, nXF, nYF, nZF;
3403385ff7eSMatthew G. Knepley   PetscErrorCode ierr;
3413385ff7eSMatthew G. Knepley 
3423385ff7eSMatthew G. Knepley   PetscFunctionBegin;
3433385ff7eSMatthew G. Knepley   if (!cone) {ierr = DMGetWorkArray(dm, 6, PETSC_INT, cone);CHKERRQ(ierr);}
3443385ff7eSMatthew G. Knepley   ierr = DMDAGetNumCells(dm, &nCx, &nCy, &nCz, &nC);CHKERRQ(ierr);
3453385ff7eSMatthew G. Knepley   ierr = DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV);CHKERRQ(ierr);
3463385ff7eSMatthew G. Knepley   ierr = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr);
3473385ff7eSMatthew G. Knepley   switch (dim) {
3483385ff7eSMatthew G. Knepley   case 2:
3493385ff7eSMatthew G. Knepley     if (p >= 0) {
3503385ff7eSMatthew G. Knepley       if (p < nC) {
3513385ff7eSMatthew G. Knepley         const PetscInt cy = p / nCx;
3523385ff7eSMatthew G. Knepley         const PetscInt cx = p % nCx;
3533385ff7eSMatthew G. Knepley 
3543385ff7eSMatthew G. Knepley         (*cone)[0] = cy*nxF + cx + nC+nV;
3553385ff7eSMatthew G. Knepley         (*cone)[1] = cx*nyF + cy + nyF + nC+nV+nXF;
3563385ff7eSMatthew G. Knepley         (*cone)[2] = cy*nxF + cx + nxF + nC+nV;
3573385ff7eSMatthew G. Knepley         (*cone)[3] = cx*nyF + cy + nC+nV+nXF;
3583385ff7eSMatthew G. Knepley         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do cell cones");
3593385ff7eSMatthew G. Knepley       } else if (p < nC+nV) {
3603385ff7eSMatthew G. Knepley       } else if (p < nC+nV+nXF) {
3613385ff7eSMatthew G. Knepley         const PetscInt fy = (p - nC+nV) / nxF;
3623385ff7eSMatthew G. Knepley         const PetscInt fx = (p - nC+nV) % nxF;
3633385ff7eSMatthew G. Knepley 
3643385ff7eSMatthew G. Knepley         (*cone)[0] = fy*nVx + fx + nC;
3653385ff7eSMatthew G. Knepley         (*cone)[1] = fy*nVx + fx + 1 + nC;
3663385ff7eSMatthew G. Knepley       } else if (p < nC+nV+nXF+nYF) {
3673385ff7eSMatthew G. Knepley         const PetscInt fx = (p - nC+nV+nXF) / nyF;
3683385ff7eSMatthew G. Knepley         const PetscInt fy = (p - nC+nV+nXF) % nyF;
3693385ff7eSMatthew G. Knepley 
3703385ff7eSMatthew G. Knepley         (*cone)[0] = fy*nVx + fx + nC;
3713385ff7eSMatthew G. Knepley         (*cone)[1] = fy*nVx + fx + nVx + nC;
3723385ff7eSMatthew G. Knepley       } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %d should be in [0, %d)", p, nC+nV+nXF+nYF+nZF);
3733385ff7eSMatthew G. Knepley     } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative point %d is invalid", p);
3743385ff7eSMatthew G. Knepley     break;
3753385ff7eSMatthew G. Knepley   case 3:
3763385ff7eSMatthew G. Knepley     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do 3D");
3773385ff7eSMatthew G. Knepley     break;
3783385ff7eSMatthew G. Knepley   }
3793385ff7eSMatthew G. Knepley   PetscFunctionReturn(0);
3803385ff7eSMatthew G. Knepley }
3813385ff7eSMatthew G. Knepley 
3823385ff7eSMatthew G. Knepley #undef __FUNCT__
3833385ff7eSMatthew G. Knepley #define __FUNCT__ "DMDARestoreCone"
3843385ff7eSMatthew G. Knepley PetscErrorCode DMDARestoreCone(DM dm, PetscInt p, PetscInt *cone[])
3853385ff7eSMatthew G. Knepley {
3863385ff7eSMatthew G. Knepley   PetscErrorCode ierr;
3873385ff7eSMatthew G. Knepley 
3883385ff7eSMatthew G. Knepley   PetscFunctionBegin;
3893385ff7eSMatthew G. Knepley   ierr = DMGetWorkArray(dm, 6, PETSC_INT, cone);CHKERRQ(ierr);
3903385ff7eSMatthew G. Knepley   PetscFunctionReturn(0);
3913385ff7eSMatthew G. Knepley }
3923385ff7eSMatthew G. Knepley 
3933385ff7eSMatthew G. Knepley #undef __FUNCT__
394a66d4d66SMatthew G Knepley #define __FUNCT__ "DMDACreateSection"
395a66d4d66SMatthew G Knepley /*@C
396a66d4d66SMatthew G Knepley   DMDACreateSection - Create a PetscSection inside the DMDA that describes data layout. This allows multiple fields with
397a66d4d66SMatthew G Knepley   different numbers of dofs on vertices, cells, and faces in each direction.
398a66d4d66SMatthew G Knepley 
399a66d4d66SMatthew G Knepley   Input Parameters:
400a66d4d66SMatthew G Knepley + dm- The DMDA
401a66d4d66SMatthew G Knepley . numFields - The number of fields
402affa55c7SMatthew G. Knepley . numComp - The number of components in each field
403affa55c7SMatthew G. Knepley . numDof - The number of dofs per dimension for each field
4040298fd71SBarry Smith . numFaceDof - The number of dofs per face for each field and direction, or NULL
405a66d4d66SMatthew G Knepley 
406a66d4d66SMatthew G Knepley   Level: developer
407a66d4d66SMatthew G Knepley 
408a66d4d66SMatthew G Knepley   Note:
409a66d4d66SMatthew G Knepley   The default DMDA numbering is as follows:
410a66d4d66SMatthew G Knepley 
411a66d4d66SMatthew G Knepley     - Cells:    [0,             nC)
412a66d4d66SMatthew G Knepley     - Vertices: [nC,            nC+nV)
41388ed4aceSMatthew G Knepley     - X-Faces:  [nC+nV,         nC+nV+nXF)         normal is +- x-dir
41488ed4aceSMatthew G Knepley     - Y-Faces:  [nC+nV+nXF,     nC+nV+nXF+nYF)     normal is +- y-dir
41588ed4aceSMatthew G Knepley     - Z-Faces:  [nC+nV+nXF+nYF, nC+nV+nXF+nYF+nZF) normal is +- z-dir
416a66d4d66SMatthew G Knepley 
417a66d4d66SMatthew G Knepley   We interpret the default DMDA partition as a cell partition, and the data assignment as a cell assignment.
418a66d4d66SMatthew G Knepley @*/
419*44e1f9abSMatthew G. Knepley PetscErrorCode DMDACreateSection(DM dm, const PetscInt numComp[], const PetscInt numDof[], const PetscInt numFaceDof[], PetscSection *s)
420a66d4d66SMatthew G Knepley {
421a66d4d66SMatthew G Knepley   DM_DA            *da  = (DM_DA*) dm->data;
422a4b60ecfSMatthew G. Knepley   PetscSection      section;
42388ed4aceSMatthew G Knepley   const PetscInt    dim = da->dim;
42480800b1aSMatthew G Knepley   PetscInt          numFields, numVertexTotDof = 0, numCellTotDof = 0, numFaceTotDof[3] = {0, 0, 0};
425b2fff234SMatthew G. Knepley   PetscBT           isLeaf;
42688ed4aceSMatthew G Knepley   PetscSF           sf;
427feafbc5cSMatthew G Knepley   PetscMPIInt       rank;
42888ed4aceSMatthew G Knepley   const PetscMPIInt *neighbors;
42988ed4aceSMatthew G Knepley   PetscInt          *localPoints;
43088ed4aceSMatthew G Knepley   PetscSFNode       *remotePoints;
431f5eeb527SMatthew G Knepley   PetscInt          nleaves = 0,  nleavesCheck = 0, nL = 0;
43257459e9aSMatthew G Knepley   PetscInt          nC, nVx, nVy, nVz, nV, nxF, nXF, nyF, nYF, nzF, nZF;
43357459e9aSMatthew G Knepley   PetscInt          pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart, yfEnd, zfStart, zfEnd;
43488ed4aceSMatthew G Knepley   PetscInt          f, v, c, xf, yf, zf, xn, yn, zn;
435a66d4d66SMatthew G Knepley   PetscErrorCode    ierr;
436a66d4d66SMatthew G Knepley 
437a66d4d66SMatthew G Knepley   PetscFunctionBegin;
438a66d4d66SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4393582350cSMatthew G. Knepley   PetscValidPointer(s, 4);
44082f516ccSBarry Smith   ierr    = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
4413582350cSMatthew G. Knepley   ierr    = DMDAGetNumCells(dm, NULL, NULL, NULL, &nC);CHKERRQ(ierr);
44257459e9aSMatthew G Knepley   ierr    = DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV);CHKERRQ(ierr);
44357459e9aSMatthew G Knepley   ierr    = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr);
44457459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, -1,  &pStart, &pEnd);CHKERRQ(ierr);
44557459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, 0,   &cStart, &cEnd);CHKERRQ(ierr);
44657459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, 1,   &fStart, &fEnd);CHKERRQ(ierr);
44757459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr);
44857459e9aSMatthew G Knepley   xfStart = vEnd;  xfEnd = xfStart+nXF;
44957459e9aSMatthew G Knepley   yfStart = xfEnd; yfEnd = yfStart+nYF;
45057459e9aSMatthew G Knepley   zfStart = yfEnd; zfEnd = zfStart+nZF;
45182f516ccSBarry Smith   if (zfEnd != fEnd) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Invalid face end %d, should be %d", zfEnd, fEnd);
45288ed4aceSMatthew G Knepley   /* Create local section */
45380800b1aSMatthew G Knepley   ierr = DMDAGetInfo(dm, 0,0,0,0,0,0,0, &numFields, 0,0,0,0,0);CHKERRQ(ierr);
454a66d4d66SMatthew G Knepley   for (f = 0; f < numFields; ++f) {
455affa55c7SMatthew G. Knepley     numVertexTotDof  += numDof[f*(dim+1)+0];
456affa55c7SMatthew G. Knepley     numCellTotDof    += numDof[f*(dim+1)+dim];
457affa55c7SMatthew G. Knepley     numFaceTotDof[0] += dim > 0 ? (numFaceDof ? numFaceDof[f*dim+0] : numDof[f*(dim+1)+dim-1]) : 0;
458affa55c7SMatthew G. Knepley     numFaceTotDof[1] += dim > 1 ? (numFaceDof ? numFaceDof[f*dim+1] : numDof[f*(dim+1)+dim-1]) : 0;
459affa55c7SMatthew G. Knepley     numFaceTotDof[2] += dim > 2 ? (numFaceDof ? numFaceDof[f*dim+2] : numDof[f*(dim+1)+dim-1]) : 0;
460a66d4d66SMatthew G Knepley   }
461a4b60ecfSMatthew G. Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &section);CHKERRQ(ierr);
462affa55c7SMatthew G. Knepley   if (numFields > 0) {
463a4b60ecfSMatthew G. Knepley     ierr = PetscSectionSetNumFields(section, numFields);CHKERRQ(ierr);
464a66d4d66SMatthew G Knepley     for (f = 0; f < numFields; ++f) {
46580800b1aSMatthew G Knepley       const char *name;
46680800b1aSMatthew G Knepley 
46780800b1aSMatthew G Knepley       ierr = DMDAGetFieldName(dm, f, &name);CHKERRQ(ierr);
468affa55c7SMatthew G. Knepley       ierr = PetscSectionSetFieldName(section, f, name ? name : "Field");CHKERRQ(ierr);
46980800b1aSMatthew G Knepley       if (numComp) {
470a4b60ecfSMatthew G. Knepley         ierr = PetscSectionSetFieldComponents(section, f, numComp[f]);CHKERRQ(ierr);
471a66d4d66SMatthew G Knepley       }
472a66d4d66SMatthew G Knepley     }
473a66d4d66SMatthew G Knepley   }
474a4b60ecfSMatthew G. Knepley   ierr = PetscSectionSetChart(section, pStart, pEnd);CHKERRQ(ierr);
475a66d4d66SMatthew G Knepley   for (v = vStart; v < vEnd; ++v) {
476a66d4d66SMatthew G Knepley     for (f = 0; f < numFields; ++f) {
477affa55c7SMatthew G. Knepley       ierr = PetscSectionSetFieldDof(section, v, f, numDof[f*(dim+1)+0]);CHKERRQ(ierr);
478a66d4d66SMatthew G Knepley     }
479a4b60ecfSMatthew G. Knepley     ierr = PetscSectionSetDof(section, v, numVertexTotDof);CHKERRQ(ierr);
480a66d4d66SMatthew G Knepley   }
481a66d4d66SMatthew G Knepley   for (xf = xfStart; xf < xfEnd; ++xf) {
482a66d4d66SMatthew G Knepley     for (f = 0; f < numFields; ++f) {
483affa55c7SMatthew G. Knepley       ierr = PetscSectionSetFieldDof(section, xf, f, numFaceDof ? numFaceDof[f*dim+0] : numDof[f*(dim+1)+dim-1]);CHKERRQ(ierr);
484a66d4d66SMatthew G Knepley     }
485a4b60ecfSMatthew G. Knepley     ierr = PetscSectionSetDof(section, xf, numFaceTotDof[0]);CHKERRQ(ierr);
486a66d4d66SMatthew G Knepley   }
487a66d4d66SMatthew G Knepley   for (yf = yfStart; yf < yfEnd; ++yf) {
488a66d4d66SMatthew G Knepley     for (f = 0; f < numFields; ++f) {
489affa55c7SMatthew G. Knepley       ierr = PetscSectionSetFieldDof(section, yf, f, numFaceDof ? numFaceDof[f*dim+1] : numDof[f*(dim+1)+dim-1]);CHKERRQ(ierr);
490a66d4d66SMatthew G Knepley     }
491a4b60ecfSMatthew G. Knepley     ierr = PetscSectionSetDof(section, yf, numFaceTotDof[1]);CHKERRQ(ierr);
492a66d4d66SMatthew G Knepley   }
493a66d4d66SMatthew G Knepley   for (zf = zfStart; zf < zfEnd; ++zf) {
494a66d4d66SMatthew G Knepley     for (f = 0; f < numFields; ++f) {
495affa55c7SMatthew G. Knepley       ierr = PetscSectionSetFieldDof(section, zf, f, numFaceDof ? numFaceDof[f*dim+2] : numDof[f*(dim+1)+dim-1]);CHKERRQ(ierr);
496a66d4d66SMatthew G Knepley     }
497a4b60ecfSMatthew G. Knepley     ierr = PetscSectionSetDof(section, zf, numFaceTotDof[2]);CHKERRQ(ierr);
498a66d4d66SMatthew G Knepley   }
499a66d4d66SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
500a66d4d66SMatthew G Knepley     for (f = 0; f < numFields; ++f) {
501affa55c7SMatthew G. Knepley       ierr = PetscSectionSetFieldDof(section, c, f, numDof[f*(dim+1)+dim]);CHKERRQ(ierr);
502a66d4d66SMatthew G Knepley     }
503a4b60ecfSMatthew G. Knepley     ierr = PetscSectionSetDof(section, c, numCellTotDof);CHKERRQ(ierr);
504a66d4d66SMatthew G Knepley   }
505a4b60ecfSMatthew G. Knepley   ierr = PetscSectionSetUp(section);CHKERRQ(ierr);
50688ed4aceSMatthew G Knepley   /* Create mesh point SF */
507b2fff234SMatthew G. Knepley   ierr = PetscBTCreate(pEnd-pStart, &isLeaf);CHKERRQ(ierr);
50888ed4aceSMatthew G Knepley   ierr = DMDAGetNeighbors(dm, &neighbors);CHKERRQ(ierr);
50988ed4aceSMatthew G Knepley   for (zn = 0; zn < (dim > 2 ? 3 : 1); ++zn) {
51088ed4aceSMatthew G Knepley     for (yn = 0; yn < (dim > 1 ? 3 : 1); ++yn) {
51188ed4aceSMatthew G Knepley       for (xn = 0; xn < 3; ++xn) {
5127128ae9fSMatthew G Knepley         const PetscInt xp       = xn-1, yp = dim > 1 ? yn-1 : 0, zp = dim > 2 ? zn-1 : 0;
51388ed4aceSMatthew G Knepley         const PetscInt neighbor = neighbors[(zn*3+yn)*3+xn];
514b2fff234SMatthew G. Knepley         PetscInt       xv, yv, zv;
51588ed4aceSMatthew G Knepley 
5163814d604SMatthew G Knepley         if (neighbor >= 0 && neighbor < rank) {
517b2fff234SMatthew G. Knepley           if (xp < 0) { /* left */
518b2fff234SMatthew G. Knepley             if (yp < 0) { /* bottom */
519b2fff234SMatthew G. Knepley               if (zp < 0) { /* back */
520b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = (      0*nVy +     0)*nVx +     0 + nC; /* left bottom back vertex */
521b2fff234SMatthew G. Knepley                 if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
522b2fff234SMatthew G. Knepley               } else if (zp > 0) { /* front */
523b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = ((nVz-1)*nVy +     0)*nVx +     0 + nC; /* left bottom front vertex */
524b2fff234SMatthew G. Knepley                 if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
525b2fff234SMatthew G. Knepley               } else {
526b2fff234SMatthew G. Knepley                 for (zv = 0; zv < nVz; ++zv) {
527b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (zv*nVy +     0)*nVx +     0 + nC; /* left bottom vertices */
528b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
529b2fff234SMatthew G. Knepley                 }
530b2fff234SMatthew G. Knepley               }
531b2fff234SMatthew G. Knepley             } else if (yp > 0) { /* top */
532b2fff234SMatthew G. Knepley               if (zp < 0) { /* back */
533b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = (      0*nVy + nVy-1)*nVx +     0 + nC; /* left top 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 + nVy-1)*nVx +     0 + nC; /* left top 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 + nVy-1)*nVx +     0 + nC; /* left top vertices */
541b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
542b2fff234SMatthew G. Knepley                 }
543b2fff234SMatthew G. Knepley               }
544b2fff234SMatthew G. Knepley             } else {
545b2fff234SMatthew G. Knepley               if (zp < 0) { /* back */
546b2fff234SMatthew G. Knepley                 for (yv = 0; yv < nVy; ++yv) {
547b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (      0*nVy + yv)*nVx +     0 + nC; /* left back vertices */
548b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
549b2fff234SMatthew G. Knepley                 }
550b2fff234SMatthew G. Knepley               } else if (zp > 0) { /* front */
551b2fff234SMatthew G. Knepley                 for (yv = 0; yv < nVy; ++yv) {
552b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = ((nVz-1)*nVy + yv)*nVx +     0 + nC; /* left front vertices */
553b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
554b2fff234SMatthew G. Knepley                 }
555b2fff234SMatthew G. Knepley               } else {
556b2fff234SMatthew G. Knepley                 for (zv = 0; zv < nVz; ++zv) {
557b2fff234SMatthew G. Knepley                   for (yv = 0; yv < nVy; ++yv) {
558b2fff234SMatthew G. Knepley                     const PetscInt localVertex  = (zv*nVy + yv)*nVx +     0 + nC; /* left vertices */
559b2fff234SMatthew G. Knepley                     if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
560b2fff234SMatthew G. Knepley                   }
561b2fff234SMatthew G. Knepley                 }
562b2fff234SMatthew G. Knepley #if 0
563b2fff234SMatthew G. Knepley                 for (xf = 0; xf < nxF; ++xf) {
564b2fff234SMatthew G. Knepley                   /* THIS IS WRONG */
565b2fff234SMatthew G. Knepley                   const PetscInt localFace  = 0 + nC+nV; /* left faces */
566b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localFace)) ++nleaves;
567b2fff234SMatthew G. Knepley                 }
568b2fff234SMatthew G. Knepley #endif
569b2fff234SMatthew G. Knepley               }
570b2fff234SMatthew G. Knepley             }
571b2fff234SMatthew G. Knepley           } else if (xp > 0) { /* right */
572b2fff234SMatthew G. Knepley             if (yp < 0) { /* bottom */
573b2fff234SMatthew G. Knepley               if (zp < 0) { /* back */
574b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = (      0*nVy +     0)*nVx + nVx-1 + nC; /* right bottom back vertex */
575b2fff234SMatthew G. Knepley                 if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
576b2fff234SMatthew G. Knepley               } else if (zp > 0) { /* front */
577b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = ((nVz-1)*nVy +     0)*nVx + nVx-1 + nC; /* right bottom front vertex */
578b2fff234SMatthew G. Knepley                 if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
579b2fff234SMatthew G. Knepley               } else {
580b2fff234SMatthew G. Knepley                 for (zv = 0; zv < nVz; ++zv) {
581b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (zv*nVy +     0)*nVx + nVx-1 + nC; /* right bottom vertices */
582b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
583b2fff234SMatthew G. Knepley                 }
584b2fff234SMatthew G. Knepley               }
585b2fff234SMatthew G. Knepley             } else if (yp > 0) { /* top */
586b2fff234SMatthew G. Knepley               if (zp < 0) { /* back */
587b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = (      0*nVy + nVy-1)*nVx + nVx-1 + nC; /* right top back vertex */
588b2fff234SMatthew G. Knepley                 if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
589b2fff234SMatthew G. Knepley               } else if (zp > 0) { /* front */
590b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = ((nVz-1)*nVy + nVy-1)*nVx + nVx-1 + nC; /* right top front vertex */
591b2fff234SMatthew G. Knepley                 if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
592b2fff234SMatthew G. Knepley               } else {
593b2fff234SMatthew G. Knepley                 for (zv = 0; zv < nVz; ++zv) {
594b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (zv*nVy + nVy-1)*nVx + nVx-1 + nC; /* right top vertices */
595b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
596b2fff234SMatthew G. Knepley                 }
597b2fff234SMatthew G. Knepley               }
598b2fff234SMatthew G. Knepley             } else {
599b2fff234SMatthew G. Knepley               if (zp < 0) { /* back */
600b2fff234SMatthew G. Knepley                 for (yv = 0; yv < nVy; ++yv) {
601b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (      0*nVy + yv)*nVx + nVx-1 + nC; /* right back vertices */
602b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
603b2fff234SMatthew G. Knepley                 }
604b2fff234SMatthew G. Knepley               } else if (zp > 0) { /* front */
605b2fff234SMatthew G. Knepley                 for (yv = 0; yv < nVy; ++yv) {
606b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = ((nVz-1)*nVy + yv)*nVx + nVx-1 + nC; /* right front vertices */
607b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
608b2fff234SMatthew G. Knepley                 }
609b2fff234SMatthew G. Knepley               } else {
610b2fff234SMatthew G. Knepley                 for (zv = 0; zv < nVz; ++zv) {
611b2fff234SMatthew G. Knepley                   for (yv = 0; yv < nVy; ++yv) {
612b2fff234SMatthew G. Knepley                     const PetscInt localVertex  = (zv*nVy + yv)*nVx + nVx-1 + nC; /* right vertices */
613b2fff234SMatthew G. Knepley                     if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
614b2fff234SMatthew G. Knepley                   }
615b2fff234SMatthew G. Knepley                 }
616b2fff234SMatthew G. Knepley #if 0
617b2fff234SMatthew G. Knepley                 for (xf = 0; xf < nxF; ++xf) {
618b2fff234SMatthew G. Knepley                   /* THIS IS WRONG */
619b2fff234SMatthew G. Knepley                   const PetscInt localFace  = 0 + nC+nV; /* right faces */
620b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localFace)) ++nleaves;
621b2fff234SMatthew G. Knepley                 }
622b2fff234SMatthew G. Knepley #endif
623b2fff234SMatthew G. Knepley               }
624b2fff234SMatthew G. Knepley             }
625b2fff234SMatthew G. Knepley           } else {
626b2fff234SMatthew G. Knepley             if (yp < 0) { /* bottom */
627b2fff234SMatthew G. Knepley               if (zp < 0) { /* back */
628b2fff234SMatthew G. Knepley                 for (xv = 0; xv < nVx; ++xv) {
629b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (      0*nVy +     0)*nVx + xv + nC; /* bottom back vertices */
630b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
631b2fff234SMatthew G. Knepley                 }
632b2fff234SMatthew G. Knepley               } else if (zp > 0) { /* front */
633b2fff234SMatthew G. Knepley                 for (xv = 0; xv < nVx; ++xv) {
634b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = ((nVz-1)*nVy +     0)*nVx + xv + nC; /* bottom front vertices */
635b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
636b2fff234SMatthew G. Knepley                 }
637b2fff234SMatthew G. Knepley               } else {
638b2fff234SMatthew G. Knepley                 for (zv = 0; zv < nVz; ++zv) {
639b2fff234SMatthew G. Knepley                   for (xv = 0; xv < nVx; ++xv) {
640b2fff234SMatthew G. Knepley                     const PetscInt localVertex  = (zv*nVy +     0)*nVx + xv + nC; /* bottom vertices */
641b2fff234SMatthew G. Knepley                     if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
642b2fff234SMatthew G. Knepley                   }
643b2fff234SMatthew G. Knepley                 }
644b2fff234SMatthew G. Knepley #if 0
645b2fff234SMatthew G. Knepley                 for (yf = 0; yf < nyF; ++yf) {
646b2fff234SMatthew G. Knepley                   /* THIS IS WRONG */
647b2fff234SMatthew G. Knepley                   const PetscInt localFace  = 0 + nC+nV; /* bottom faces */
648b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVFace)) ++nleaves;
649b2fff234SMatthew G. Knepley                 }
650b2fff234SMatthew G. Knepley #endif
651b2fff234SMatthew G. Knepley               }
652b2fff234SMatthew G. Knepley             } else if (yp > 0) { /* top */
653b2fff234SMatthew G. Knepley               if (zp < 0) { /* back */
654b2fff234SMatthew G. Knepley                 for (xv = 0; xv < nVx; ++xv) {
655b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (      0*nVy + nVy-1)*nVx + xv + nC; /* top back vertices */
656b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
657b2fff234SMatthew G. Knepley                 }
658b2fff234SMatthew G. Knepley               } else if (zp > 0) { /* front */
659b2fff234SMatthew G. Knepley                 for (xv = 0; xv < nVx; ++xv) {
660b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = ((nVz-1)*nVy + nVy-1)*nVx + xv + nC; /* top front vertices */
661b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
662b2fff234SMatthew G. Knepley                 }
663b2fff234SMatthew G. Knepley               } else {
664b2fff234SMatthew G. Knepley                 for (zv = 0; zv < nVz; ++zv) {
665b2fff234SMatthew G. Knepley                   for (xv = 0; xv < nVx; ++xv) {
666b2fff234SMatthew G. Knepley                     const PetscInt localVertex  = (zv*nVy + nVy-1)*nVx + xv + nC; /* top vertices */
667b2fff234SMatthew G. Knepley                     if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
668b2fff234SMatthew G. Knepley                   }
669b2fff234SMatthew G. Knepley                 }
670b2fff234SMatthew G. Knepley #if 0
671b2fff234SMatthew G. Knepley                 for (yf = 0; yf < nyF; ++yf) {
672b2fff234SMatthew G. Knepley                   /* THIS IS WRONG */
673b2fff234SMatthew G. Knepley                   const PetscInt localFace  = 0 + nC+nV; /* top faces */
674b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVFace)) ++nleaves;
675b2fff234SMatthew G. Knepley                 }
676b2fff234SMatthew G. Knepley #endif
677b2fff234SMatthew G. Knepley               }
678b2fff234SMatthew G. Knepley             } else {
679b2fff234SMatthew G. Knepley               if (zp < 0) { /* back */
680b2fff234SMatthew G. Knepley                 for (yv = 0; yv < nVy; ++yv) {
681b2fff234SMatthew G. Knepley                   for (xv = 0; xv < nVx; ++xv) {
682b2fff234SMatthew G. Knepley                     const PetscInt localVertex  = (      0*nVy + yv)*nVx + xv + nC; /* back vertices */
683b2fff234SMatthew G. Knepley                     if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
684b2fff234SMatthew G. Knepley                   }
685b2fff234SMatthew G. Knepley                 }
686b2fff234SMatthew G. Knepley #if 0
687b2fff234SMatthew G. Knepley                 for (zf = 0; zf < nzF; ++zf) {
688b2fff234SMatthew G. Knepley                   /* THIS IS WRONG */
689b2fff234SMatthew G. Knepley                   const PetscInt localFace  = 0 + nC+nV; /* back faces */
690b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localFace)) ++nleaves;
691b2fff234SMatthew G. Knepley                 }
692b2fff234SMatthew G. Knepley #endif
693b2fff234SMatthew G. Knepley               } else if (zp > 0) { /* front */
694b2fff234SMatthew G. Knepley                 for (yv = 0; yv < nVy; ++yv) {
695b2fff234SMatthew G. Knepley                   for (xv = 0; xv < nVx; ++xv) {
696b2fff234SMatthew G. Knepley                     const PetscInt localVertex  = ((nVz-1)*nVy + yv)*nVx + xv + nC; /* front vertices */
697b2fff234SMatthew G. Knepley                     if (!PetscBTLookupSet(isLeaf, localVertex)) ++nleaves;
698b2fff234SMatthew G. Knepley                   }
699b2fff234SMatthew G. Knepley                 }
700b2fff234SMatthew G. Knepley #if 0
701b2fff234SMatthew G. Knepley                 for (zf = 0; zf < nzF; ++zf) {
702b2fff234SMatthew G. Knepley                   /* THIS IS WRONG */
703b2fff234SMatthew G. Knepley                   const PetscInt localFace  = 0 + nC+nV; /* front faces */
704b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localFace)) ++nleaves;
705b2fff234SMatthew G. Knepley                 }
706b2fff234SMatthew G. Knepley #endif
707b2fff234SMatthew G. Knepley               } else {
708b2fff234SMatthew G. Knepley                 /* Nothing is shared from the interior */
70988ed4aceSMatthew G Knepley               }
71088ed4aceSMatthew G Knepley             }
71188ed4aceSMatthew G Knepley           }
71288ed4aceSMatthew G Knepley         }
71388ed4aceSMatthew G Knepley       }
714b2fff234SMatthew G. Knepley     }
715b2fff234SMatthew G. Knepley   }
716b2fff234SMatthew G. Knepley   ierr = PetscBTMemzero(pEnd-pStart, isLeaf);CHKERRQ(ierr);
717dcca6d9dSJed Brown   ierr = PetscMalloc2(nleaves,&localPoints,nleaves,&remotePoints);CHKERRQ(ierr);
71888ed4aceSMatthew G Knepley   for (zn = 0; zn < (dim > 2 ? 3 : 1); ++zn) {
71988ed4aceSMatthew G Knepley     for (yn = 0; yn < (dim > 1 ? 3 : 1); ++yn) {
72088ed4aceSMatthew G Knepley       for (xn = 0; xn < 3; ++xn) {
7217128ae9fSMatthew G Knepley         const PetscInt xp       = xn-1, yp = dim > 1 ? yn-1 : 0, zp = dim > 2 ? zn-1 : 0;
72288ed4aceSMatthew G Knepley         const PetscInt neighbor = neighbors[(zn*3+yn)*3+xn];
723f5eeb527SMatthew G Knepley         PetscInt       xv, yv, zv;
72488ed4aceSMatthew G Knepley 
7253814d604SMatthew G Knepley         if (neighbor >= 0 && neighbor < rank) {
72688ed4aceSMatthew G Knepley           if (xp < 0) { /* left */
72788ed4aceSMatthew G Knepley             if (yp < 0) { /* bottom */
72888ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
729b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = (      0*nVy +     0)*nVx +     0 + nC; /* left bottom back vertex */
730f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = ((nVz-1)*nVy + nVy-1)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
731f5eeb527SMatthew G Knepley 
732b2fff234SMatthew G. Knepley                 if (!PetscBTLookupSet(isLeaf, localVertex)) {
733f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
734f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
735f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
736f5eeb527SMatthew G Knepley                   ++nL;
737b2fff234SMatthew G. Knepley                 }
73888ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
739b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = ((nVz-1)*nVy +     0)*nVx +     0 + nC; /* left bottom front vertex */
740f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = (      0*nVy + nVy-1)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
741f5eeb527SMatthew G Knepley 
742b2fff234SMatthew G. Knepley                 if (!PetscBTLookupSet(isLeaf, localVertex)) {
743f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
744f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
745f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
746f5eeb527SMatthew G Knepley                   ++nL;
747b2fff234SMatthew G. Knepley                 }
74888ed4aceSMatthew G Knepley               } else {
749b2fff234SMatthew G. Knepley                 for (zv = 0; zv < nVz; ++zv) {
750b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (zv*nVy +     0)*nVx +     0 + nC; /* left bottom vertices */
751f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (zv*nVy + nVy-1)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
752f5eeb527SMatthew G Knepley 
753b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) {
754f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
755f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
756f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
757b2fff234SMatthew G. Knepley                     ++nL;
758b2fff234SMatthew G. Knepley                   }
759f5eeb527SMatthew G Knepley                 }
76088ed4aceSMatthew G Knepley               }
76188ed4aceSMatthew G Knepley             } else if (yp > 0) { /* top */
76288ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
763b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = (      0*nVy + nVy-1)*nVx +     0 + nC; /* left top back vertex */
764f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = ((nVz-1)*nVy +     0)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
765f5eeb527SMatthew G Knepley 
766b2fff234SMatthew G. Knepley                 if (!PetscBTLookupSet(isLeaf, localVertex)) {
767f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
768f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
769f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
770f5eeb527SMatthew G Knepley                   ++nL;
771b2fff234SMatthew G. Knepley                 }
77288ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
773b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = ((nVz-1)*nVy + nVy-1)*nVx +     0 + nC; /* left top front vertex */
774f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = (      0*nVy +     0)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
775f5eeb527SMatthew G Knepley 
776b2fff234SMatthew G. Knepley                 if (!PetscBTLookupSet(isLeaf, localVertex)) {
777f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
778f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
779f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
780f5eeb527SMatthew G Knepley                   ++nL;
781b2fff234SMatthew G. Knepley                 }
78288ed4aceSMatthew G Knepley               } else {
783b2fff234SMatthew G. Knepley                 for (zv = 0; zv < nVz; ++zv) {
784b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (zv*nVy + nVy-1)*nVx +     0 + nC; /* left top vertices */
785f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (zv*nVy +     0)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
786f5eeb527SMatthew G Knepley 
787b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) {
788f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
789f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
790f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
791b2fff234SMatthew G. Knepley                     ++nL;
792b2fff234SMatthew G. Knepley                   }
793f5eeb527SMatthew G Knepley                 }
79488ed4aceSMatthew G Knepley               }
79588ed4aceSMatthew G Knepley             } else {
79688ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
797b2fff234SMatthew G. Knepley                 for (yv = 0; yv < nVy; ++yv) {
798b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (      0*nVy + yv)*nVx +     0 + nC; /* left back vertices */
799f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = ((nVz-1)*nVy + yv)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
800f5eeb527SMatthew G Knepley 
801b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) {
802f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
803f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
804f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
805b2fff234SMatthew G. Knepley                     ++nL;
806b2fff234SMatthew G. Knepley                   }
807f5eeb527SMatthew G Knepley                 }
80888ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
809b2fff234SMatthew G. Knepley                 for (yv = 0; yv < nVy; ++yv) {
810b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = ((nVz-1)*nVy + yv)*nVx +     0 + nC; /* left front vertices */
811f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (      0*nVy + yv)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
812f5eeb527SMatthew G Knepley 
813b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) {
814f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
815f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
816f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
817b2fff234SMatthew G. Knepley                     ++nL;
818b2fff234SMatthew G. Knepley                   }
819f5eeb527SMatthew G Knepley                 }
82088ed4aceSMatthew G Knepley               } else {
821f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv) {
822b2fff234SMatthew G. Knepley                   for (yv = 0; yv < nVy; ++yv) {
823b2fff234SMatthew G. Knepley                     const PetscInt localVertex  = (zv*nVy + yv)*nVx +     0 + nC; /* left vertices */
824f5eeb527SMatthew G Knepley                     const PetscInt remoteVertex = (zv*nVy + yv)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
825f5eeb527SMatthew G Knepley 
826b2fff234SMatthew G. Knepley                     if (!PetscBTLookupSet(isLeaf, localVertex)) {
827f5eeb527SMatthew G Knepley                       localPoints[nL]        = localVertex;
828f5eeb527SMatthew G Knepley                       remotePoints[nL].rank  = neighbor;
829f5eeb527SMatthew G Knepley                       remotePoints[nL].index = remoteVertex;
830b2fff234SMatthew G. Knepley                       ++nL;
831f5eeb527SMatthew G Knepley                     }
832f5eeb527SMatthew G Knepley                   }
833b2fff234SMatthew G. Knepley                 }
834b2fff234SMatthew G. Knepley #if 0
835b2fff234SMatthew G. Knepley                 for (xf = 0; xf < nxF; ++xf) {
836f5eeb527SMatthew G Knepley                   /* THIS IS WRONG */
837b2fff234SMatthew G. Knepley                   const PetscInt localFace  = 0 + nC+nV; /* left faces */
838f5eeb527SMatthew G Knepley                   const PetscInt remoteFace = 0 + nC+nV;
839f5eeb527SMatthew G Knepley 
840b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localFace)) {
841f5eeb527SMatthew G Knepley                     localPoints[nL]        = localFace;
842f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
843f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteFace;
844f5eeb527SMatthew G Knepley                   }
84588ed4aceSMatthew G Knepley                 }
846b2fff234SMatthew G. Knepley #endif
847b2fff234SMatthew G. Knepley               }
84888ed4aceSMatthew G Knepley             }
84988ed4aceSMatthew G Knepley           } else if (xp > 0) { /* right */
85088ed4aceSMatthew G Knepley             if (yp < 0) { /* bottom */
85188ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
852b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = (      0*nVy +     0)*nVx + nVx-1 + nC; /* right bottom back vertex */
853f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = ((nVz-1)*nVy + nVy-1)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
854f5eeb527SMatthew G Knepley 
855b2fff234SMatthew G. Knepley                 if (!PetscBTLookupSet(isLeaf, localVertex)) {
856f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
857f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
858f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
859f5eeb527SMatthew G Knepley                   ++nL;
860b2fff234SMatthew G. Knepley                 }
86188ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
862b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = ((nVz-1)*nVy +     0)*nVx + nVx-1 + nC; /* right bottom front vertex */
863f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = (      0*nVy + nVy-1)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
864f5eeb527SMatthew G Knepley 
865b2fff234SMatthew G. Knepley                 if (!PetscBTLookupSet(isLeaf, localVertex)) {
866f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
867f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
868f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
869f5eeb527SMatthew G Knepley                   ++nL;
870b2fff234SMatthew G. Knepley                 }
87188ed4aceSMatthew G Knepley               } else {
872b2fff234SMatthew G. Knepley                 nleavesCheck += nVz;
873b2fff234SMatthew G. Knepley                 for (zv = 0; zv < nVz; ++zv) {
874b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (zv*nVy +     0)*nVx + nVx-1 + nC; /* right bottom vertices */
875f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (zv*nVy + nVy-1)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
876f5eeb527SMatthew G Knepley 
877b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) {
878f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
879f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
880f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
881b2fff234SMatthew G. Knepley                     ++nL;
882b2fff234SMatthew G. Knepley                   }
883f5eeb527SMatthew G Knepley                 }
88488ed4aceSMatthew G Knepley               }
88588ed4aceSMatthew G Knepley             } else if (yp > 0) { /* top */
88688ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
887b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = (      0*nVy + nVy-1)*nVx + nVx-1 + nC; /* right top back vertex */
888f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = ((nVz-1)*nVy +     0)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
889f5eeb527SMatthew G Knepley 
890b2fff234SMatthew G. Knepley                 if (!PetscBTLookupSet(isLeaf, localVertex)) {
891f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
892f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
893f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
894f5eeb527SMatthew G Knepley                   ++nL;
895b2fff234SMatthew G. Knepley                 }
89688ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
897b2fff234SMatthew G. Knepley                 const PetscInt localVertex  = ((nVz-1)*nVy + nVy-1)*nVx + nVx-1 + nC; /* right top front vertex */
898f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = (      0*nVy +     0)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
899f5eeb527SMatthew G Knepley 
900b2fff234SMatthew G. Knepley                 if (!PetscBTLookupSet(isLeaf, localVertex)) {
901f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
902f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
903f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
904f5eeb527SMatthew G Knepley                   ++nL;
905b2fff234SMatthew G. Knepley                 }
90688ed4aceSMatthew G Knepley               } else {
907b2fff234SMatthew G. Knepley                 for (zv = 0; zv < nVz; ++zv) {
908b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (zv*nVy + nVy-1)*nVx + nVx-1 + nC; /* right top vertices */
909f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (zv*nVy +     0)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
910f5eeb527SMatthew G Knepley 
911b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) {
912f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
913f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
914f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
915b2fff234SMatthew G. Knepley                     ++nL;
916b2fff234SMatthew G. Knepley                   }
917f5eeb527SMatthew G Knepley                 }
91888ed4aceSMatthew G Knepley               }
91988ed4aceSMatthew G Knepley             } else {
92088ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
921b2fff234SMatthew G. Knepley                 for (yv = 0; yv < nVy; ++yv) {
922b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (      0*nVy + yv)*nVx + nVx-1 + nC; /* right back vertices */
923f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = ((nVz-1)*nVy + yv)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
924f5eeb527SMatthew G Knepley 
925b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) {
926f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
927f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
928f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
929b2fff234SMatthew G. Knepley                     ++nL;
930b2fff234SMatthew G. Knepley                   }
931f5eeb527SMatthew G Knepley                 }
93288ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
933b2fff234SMatthew G. Knepley                 for (yv = 0; yv < nVy; ++yv) {
934b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = ((nVz-1)*nVy + yv)*nVx + nVx-1 + nC; /* right front vertices */
935f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (      0*nVy + yv)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
936f5eeb527SMatthew G Knepley 
937b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) {
938f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
939f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
940f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
941b2fff234SMatthew G. Knepley                     ++nL;
942b2fff234SMatthew G. Knepley                   }
943f5eeb527SMatthew G Knepley                 }
94488ed4aceSMatthew G Knepley               } else {
945f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv) {
946b2fff234SMatthew G. Knepley                   for (yv = 0; yv < nVy; ++yv) {
947b2fff234SMatthew G. Knepley                     const PetscInt localVertex  = (zv*nVy + yv)*nVx + nVx-1 + nC; /* right vertices */
948f5eeb527SMatthew G Knepley                     const PetscInt remoteVertex = (zv*nVy + yv)*nVx + 0     + nC; /* TODO: Correct this for neighbor sizes */
949f5eeb527SMatthew G Knepley 
950b2fff234SMatthew G. Knepley                     if (!PetscBTLookupSet(isLeaf, localVertex)) {
951f5eeb527SMatthew G Knepley                       localPoints[nL]        = localVertex;
952f5eeb527SMatthew G Knepley                       remotePoints[nL].rank  = neighbor;
953f5eeb527SMatthew G Knepley                       remotePoints[nL].index = remoteVertex;
954b2fff234SMatthew G. Knepley                       ++nL;
955f5eeb527SMatthew G Knepley                     }
956f5eeb527SMatthew G Knepley                   }
957b2fff234SMatthew G. Knepley                 }
958b2fff234SMatthew G. Knepley #if 0
959b2fff234SMatthew G. Knepley                 for (xf = 0; xf < nxF; ++xf) {
960f5eeb527SMatthew G Knepley                   /* THIS IS WRONG */
961b2fff234SMatthew G. Knepley                   const PetscInt localFace  = 0 + nC+nV; /* right faces */
962f5eeb527SMatthew G Knepley                   const PetscInt remoteFace = 0 + nC+nV;
963f5eeb527SMatthew G Knepley 
964b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localFace)) {
965f5eeb527SMatthew G Knepley                     localPoints[nL]        = localFace;
966f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
967f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteFace;
968b2fff234SMatthew G. Knepley                     ++nL;
969f5eeb527SMatthew G Knepley                   }
97088ed4aceSMatthew G Knepley                 }
971b2fff234SMatthew G. Knepley #endif
972b2fff234SMatthew G. Knepley               }
97388ed4aceSMatthew G Knepley             }
97488ed4aceSMatthew G Knepley           } else {
97588ed4aceSMatthew G Knepley             if (yp < 0) { /* bottom */
97688ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
977b2fff234SMatthew G. Knepley                 for (xv = 0; xv < nVx; ++xv) {
978b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (      0*nVy +     0)*nVx + xv + nC; /* bottom back vertices */
979f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = ((nVz-1)*nVy + nVy-1)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
980f5eeb527SMatthew G Knepley 
981b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) {
982f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
983f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
984f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
985b2fff234SMatthew G. Knepley                     ++nL;
986b2fff234SMatthew G. Knepley                   }
987f5eeb527SMatthew G Knepley                 }
98888ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
989b2fff234SMatthew G. Knepley                 for (xv = 0; xv < nVx; ++xv) {
990b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = ((nVz-1)*nVy +     0)*nVx + xv + nC; /* bottom front vertices */
991f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (      0*nVy + nVy-1)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
992f5eeb527SMatthew G Knepley 
993b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) {
994f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
995f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
996f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
997b2fff234SMatthew G. Knepley                     ++nL;
998b2fff234SMatthew G. Knepley                   }
999f5eeb527SMatthew G Knepley                 }
100088ed4aceSMatthew G Knepley               } else {
1001f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv) {
1002b2fff234SMatthew G. Knepley                   for (xv = 0; xv < nVx; ++xv) {
1003b2fff234SMatthew G. Knepley                     const PetscInt localVertex  = (zv*nVy +     0)*nVx + xv + nC; /* bottom vertices */
1004f5eeb527SMatthew G Knepley                     const PetscInt remoteVertex = (zv*nVy + nVy-1)*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;
1011f5eeb527SMatthew G Knepley                     }
1012f5eeb527SMatthew G Knepley                   }
1013b2fff234SMatthew G. Knepley                 }
1014b2fff234SMatthew G. Knepley #if 0
1015b2fff234SMatthew G. Knepley                 for (yf = 0; yf < nyF; ++yf) {
1016f5eeb527SMatthew G Knepley                   /* THIS IS WRONG */
1017b2fff234SMatthew G. Knepley                   const PetscInt localFace  = 0 + nC+nV; /* bottom faces */
1018f5eeb527SMatthew G Knepley                   const PetscInt remoteFace = 0 + nC+nV;
1019f5eeb527SMatthew G Knepley 
1020b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localFace)) {
1021f5eeb527SMatthew G Knepley                     localPoints[nL]        = localFace;
1022f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
1023f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteFace;
1024b2fff234SMatthew G. Knepley                     ++nL;
1025f5eeb527SMatthew G Knepley                   }
102688ed4aceSMatthew G Knepley                 }
1027b2fff234SMatthew G. Knepley #endif
1028b2fff234SMatthew G. Knepley               }
102988ed4aceSMatthew G Knepley             } else if (yp > 0) { /* top */
103088ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
1031b2fff234SMatthew G. Knepley                 for (xv = 0; xv < nVx; ++xv) {
1032b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = (      0*nVy + nVy-1)*nVx + xv + nC; /* top back vertices */
1033f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = ((nVz-1)*nVy +     0)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
1034f5eeb527SMatthew G Knepley 
1035b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) {
1036f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
1037f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
1038f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
1039b2fff234SMatthew G. Knepley                     ++nL;
1040b2fff234SMatthew G. Knepley                   }
1041f5eeb527SMatthew G Knepley                 }
104288ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
1043b2fff234SMatthew G. Knepley                 for (xv = 0; xv < nVx; ++xv) {
1044b2fff234SMatthew G. Knepley                   const PetscInt localVertex  = ((nVz-1)*nVy + nVy-1)*nVx + xv + nC; /* top front vertices */
1045f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (      0*nVy +     0)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
1046f5eeb527SMatthew G Knepley 
1047b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localVertex)) {
1048f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
1049f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
1050f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
1051b2fff234SMatthew G. Knepley                     ++nL;
1052b2fff234SMatthew G. Knepley                   }
1053f5eeb527SMatthew G Knepley                 }
105488ed4aceSMatthew G Knepley               } else {
1055f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv) {
1056b2fff234SMatthew G. Knepley                   for (xv = 0; xv < nVx; ++xv) {
1057b2fff234SMatthew G. Knepley                     const PetscInt localVertex  = (zv*nVy + nVy-1)*nVx + xv + nC; /* top vertices */
1058f5eeb527SMatthew G Knepley                     const PetscInt remoteVertex = (zv*nVy +     0)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
1059f5eeb527SMatthew G Knepley 
1060b2fff234SMatthew G. Knepley                     if (!PetscBTLookupSet(isLeaf, localVertex)) {
1061f5eeb527SMatthew G Knepley                       localPoints[nL]        = localVertex;
1062f5eeb527SMatthew G Knepley                       remotePoints[nL].rank  = neighbor;
1063f5eeb527SMatthew G Knepley                       remotePoints[nL].index = remoteVertex;
1064b2fff234SMatthew G. Knepley                       ++nL;
1065f5eeb527SMatthew G Knepley                     }
1066f5eeb527SMatthew G Knepley                   }
1067b2fff234SMatthew G. Knepley                 }
1068b2fff234SMatthew G. Knepley #if 0
1069b2fff234SMatthew G. Knepley                 for (yf = 0; yf < nyF; ++yf) {
1070f5eeb527SMatthew G Knepley                   /* THIS IS WRONG */
1071b2fff234SMatthew G. Knepley                   const PetscInt localFace  = 0 + nC+nV; /* top faces */
1072f5eeb527SMatthew G Knepley                   const PetscInt remoteFace = 0 + nC+nV;
1073f5eeb527SMatthew G Knepley 
1074b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localFace)) {
1075f5eeb527SMatthew G Knepley                     localPoints[nL]        = localFace;
1076f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
1077f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteFace;
1078b2fff234SMatthew G. Knepley                     ++nL;
1079f5eeb527SMatthew G Knepley                   }
108088ed4aceSMatthew G Knepley                 }
1081b2fff234SMatthew G. Knepley #endif
1082b2fff234SMatthew G. Knepley               }
108388ed4aceSMatthew G Knepley             } else {
108488ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
1085f5eeb527SMatthew G Knepley                 for (yv = 0; yv < nVy; ++yv) {
1086b2fff234SMatthew G. Knepley                   for (xv = 0; xv < nVx; ++xv) {
1087b2fff234SMatthew G. Knepley                     const PetscInt localVertex  = (      0*nVy + yv)*nVx + xv + nC; /* back vertices */
1088f5eeb527SMatthew G Knepley                     const PetscInt remoteVertex = ((nVz-1)*nVy + yv)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
1089f5eeb527SMatthew G Knepley 
1090b2fff234SMatthew G. Knepley                     if (!PetscBTLookupSet(isLeaf, localVertex)) {
1091f5eeb527SMatthew G Knepley                       localPoints[nL]        = localVertex;
1092f5eeb527SMatthew G Knepley                       remotePoints[nL].rank  = neighbor;
1093f5eeb527SMatthew G Knepley                       remotePoints[nL].index = remoteVertex;
1094b2fff234SMatthew G. Knepley                       ++nL;
1095f5eeb527SMatthew G Knepley                     }
1096f5eeb527SMatthew G Knepley                   }
1097b2fff234SMatthew G. Knepley                 }
1098b2fff234SMatthew G. Knepley #if 0
1099b2fff234SMatthew G. Knepley                 for (zf = 0; zf < nzF; ++zf) {
1100f5eeb527SMatthew G Knepley                   /* THIS IS WRONG */
1101b2fff234SMatthew G. Knepley                   const PetscInt localFace  = 0 + nC+nV; /* back faces */
1102f5eeb527SMatthew G Knepley                   const PetscInt remoteFace = 0 + nC+nV;
1103f5eeb527SMatthew G Knepley 
1104b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localFace)) {
1105f5eeb527SMatthew G Knepley                     localPoints[nL]        = localFace;
1106f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
1107f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteFace;
1108b2fff234SMatthew G. Knepley                     ++nL;
1109f5eeb527SMatthew G Knepley                   }
1110b2fff234SMatthew G. Knepley                 }
1111b2fff234SMatthew G. Knepley #endif
111288ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
1113f5eeb527SMatthew G Knepley                 for (yv = 0; yv < nVy; ++yv) {
1114b2fff234SMatthew G. Knepley                   for (xv = 0; xv < nVx; ++xv) {
1115b2fff234SMatthew G. Knepley                     const PetscInt localVertex  = ((nVz-1)*nVy + yv)*nVx + xv + nC; /* front vertices */
1116f5eeb527SMatthew G Knepley                     const PetscInt remoteVertex = (      0*nVy + yv)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
1117f5eeb527SMatthew G Knepley 
1118b2fff234SMatthew G. Knepley                     if (!PetscBTLookupSet(isLeaf, localVertex)) {
1119f5eeb527SMatthew G Knepley                       localPoints[nL]        = localVertex;
1120f5eeb527SMatthew G Knepley                       remotePoints[nL].rank  = neighbor;
1121f5eeb527SMatthew G Knepley                       remotePoints[nL].index = remoteVertex;
1122b2fff234SMatthew G. Knepley                       ++nL;
1123f5eeb527SMatthew G Knepley                     }
1124f5eeb527SMatthew G Knepley                   }
1125b2fff234SMatthew G. Knepley                 }
1126b2fff234SMatthew G. Knepley #if 0
1127b2fff234SMatthew G. Knepley                 for (zf = 0; zf < nzF; ++zf) {
1128f5eeb527SMatthew G Knepley                   /* THIS IS WRONG */
1129b2fff234SMatthew G. Knepley                   const PetscInt localFace  = 0 + nC+nV; /* front faces */
1130f5eeb527SMatthew G Knepley                   const PetscInt remoteFace = 0 + nC+nV;
1131f5eeb527SMatthew G Knepley 
1132b2fff234SMatthew G. Knepley                   if (!PetscBTLookupSet(isLeaf, localFace)) {
1133f5eeb527SMatthew G Knepley                     localPoints[nL]        = localFace;
1134f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
1135f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteFace;
1136b2fff234SMatthew G. Knepley                     ++nL;
1137f5eeb527SMatthew G Knepley                   }
1138b2fff234SMatthew G. Knepley                 }
1139b2fff234SMatthew G. Knepley #endif
114088ed4aceSMatthew G Knepley               } else {
114188ed4aceSMatthew G Knepley                 /* Nothing is shared from the interior */
114288ed4aceSMatthew G Knepley               }
114388ed4aceSMatthew G Knepley             }
114488ed4aceSMatthew G Knepley           }
114588ed4aceSMatthew G Knepley         }
114688ed4aceSMatthew G Knepley       }
114788ed4aceSMatthew G Knepley     }
114888ed4aceSMatthew G Knepley   }
1149b2fff234SMatthew G. Knepley   ierr = PetscBTDestroy(&isLeaf);CHKERRQ(ierr);
1150b2fff234SMatthew G. Knepley   /* Remove duplication in leaf determination */
1151b2fff234SMatthew 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);
115282f516ccSBarry Smith   ierr = PetscSFCreate(PetscObjectComm((PetscObject)dm), &sf);CHKERRQ(ierr);
115388ed4aceSMatthew G Knepley   ierr = PetscSFSetGraph(sf, pEnd, nleaves, localPoints, PETSC_OWN_POINTER, remotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr);
1154a4b60ecfSMatthew G. Knepley   ierr = DMSetPointSF(dm, sf);CHKERRQ(ierr);
115588ed4aceSMatthew G Knepley   ierr = PetscSFDestroy(&sf);CHKERRQ(ierr);
1156affa55c7SMatthew G. Knepley   *s = section;
1157affa55c7SMatthew G. Knepley   PetscFunctionReturn(0);
1158affa55c7SMatthew G. Knepley }
1159affa55c7SMatthew G. Knepley 
11603385ff7eSMatthew G. Knepley #undef __FUNCT__
11613385ff7eSMatthew G. Knepley #define __FUNCT__ "DMDASetVertexCoordinates"
11623385ff7eSMatthew G. Knepley PetscErrorCode DMDASetVertexCoordinates(DM dm, PetscReal xl, PetscReal xu, PetscReal yl, PetscReal yu, PetscReal zl, PetscReal zu)
11633385ff7eSMatthew G. Knepley {
11643385ff7eSMatthew G. Knepley   DM_DA         *da = (DM_DA *) dm->data;
11653385ff7eSMatthew G. Knepley   Vec            coordinates;
11663385ff7eSMatthew G. Knepley   PetscSection   section;
11673385ff7eSMatthew G. Knepley   PetscScalar   *coords;
11683385ff7eSMatthew G. Knepley   PetscReal      h[3];
11693385ff7eSMatthew G. Knepley   PetscInt       dim, size, M, N, P, nVx, nVy, nVz, nV, vStart, vEnd, v, i, j, k;
11703385ff7eSMatthew G. Knepley   PetscErrorCode ierr;
11713385ff7eSMatthew G. Knepley 
11723385ff7eSMatthew G. Knepley   PetscFunctionBegin;
11733385ff7eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
11743385ff7eSMatthew G. Knepley   ierr = DMDAGetInfo(dm, &dim, &M, &N, &P, 0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
11753385ff7eSMatthew G. Knepley   h[0] = (xu - xl)/M;
11763385ff7eSMatthew G. Knepley   h[1] = (yu - yl)/N;
11773385ff7eSMatthew G. Knepley   h[2] = (zu - zl)/P;
11783385ff7eSMatthew G. Knepley   ierr = DMDAGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
11793385ff7eSMatthew G. Knepley   ierr = DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV);CHKERRQ(ierr);
11803385ff7eSMatthew G. Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &section);CHKERRQ(ierr);
11813385ff7eSMatthew G. Knepley   ierr = PetscSectionSetNumFields(section, 1);CHKERRQ(ierr);
11823385ff7eSMatthew G. Knepley   ierr = PetscSectionSetFieldComponents(section, 0, dim);CHKERRQ(ierr);
11833385ff7eSMatthew G. Knepley   ierr = PetscSectionSetChart(section, vStart, vEnd);CHKERRQ(ierr);
11843385ff7eSMatthew G. Knepley   for (v = vStart; v < vEnd; ++v) {
11853385ff7eSMatthew G. Knepley     ierr = PetscSectionSetDof(section, v, dim);CHKERRQ(ierr);
11863385ff7eSMatthew G. Knepley   }
11873385ff7eSMatthew G. Knepley   ierr = PetscSectionSetUp(section);CHKERRQ(ierr);
11883385ff7eSMatthew G. Knepley   ierr = PetscSectionGetStorageSize(section, &size);CHKERRQ(ierr);
11893385ff7eSMatthew G. Knepley   ierr = VecCreateSeq(PETSC_COMM_SELF, size, &coordinates);CHKERRQ(ierr);
11903385ff7eSMatthew G. Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
11913385ff7eSMatthew G. Knepley   for (k = 0; k < nVz; ++k) {
1192e4d69e08SBarry Smith     PetscInt ind[3], d, off;
11933385ff7eSMatthew G. Knepley 
1194e4d69e08SBarry Smith     ind[0] = 0;
1195e4d69e08SBarry Smith     ind[1] = 0;
1196e4d69e08SBarry Smith     ind[2] = k + da->zs;
11973385ff7eSMatthew G. Knepley     for (j = 0; j < nVy; ++j) {
11983385ff7eSMatthew G. Knepley       ind[1] = j + da->ys;
11993385ff7eSMatthew G. Knepley       for (i = 0; i < nVx; ++i) {
12003385ff7eSMatthew G. Knepley         const PetscInt vertex = (k*nVy + j)*nVx + i + vStart;
12013385ff7eSMatthew G. Knepley 
12023385ff7eSMatthew G. Knepley         ierr = PetscSectionGetOffset(section, vertex, &off);CHKERRQ(ierr);
12033385ff7eSMatthew G. Knepley         ind[0] = i + da->xs;
12043385ff7eSMatthew G. Knepley         for (d = 0; d < dim; ++d) {
12053385ff7eSMatthew G. Knepley           coords[off+d] = h[d]*ind[d];
12063385ff7eSMatthew G. Knepley         }
12073385ff7eSMatthew G. Knepley       }
12083385ff7eSMatthew G. Knepley     }
12093385ff7eSMatthew G. Knepley   }
12103385ff7eSMatthew G. Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
12113385ff7eSMatthew G. Knepley   ierr = DMSetCoordinateSection(dm, section);CHKERRQ(ierr);
12123385ff7eSMatthew G. Knepley   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
1213a4b60ecfSMatthew G. Knepley   ierr = PetscSectionDestroy(&section);CHKERRQ(ierr);
12143385ff7eSMatthew G. Knepley   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
12153385ff7eSMatthew G. Knepley   PetscFunctionReturn(0);
12163385ff7eSMatthew G. Knepley }
12179a800dd8SMatthew G. Knepley 
12189a800dd8SMatthew G. Knepley #undef __FUNCT__
12199a800dd8SMatthew G. Knepley #define __FUNCT__ "DMDAProjectFunctionLocal"
1220c110b1eeSGeoffrey Irving PetscErrorCode DMDAProjectFunctionLocal(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
12219a800dd8SMatthew G. Knepley {
12229a800dd8SMatthew G. Knepley   PetscDualSpace *sp;
12239a800dd8SMatthew G. Knepley   PetscQuadrature q;
12249a800dd8SMatthew G. Knepley   PetscSection    section;
12259a800dd8SMatthew G. Knepley   PetscScalar    *values;
12269a800dd8SMatthew G. Knepley   PetscReal      *v0, *J, *detJ;
122721454ff5SMatthew G. Knepley   PetscInt        numFields, numComp, numPoints, dim, spDim, totDim = 0, numValues, cStart, cEnd, f, c, v, d;
12289a800dd8SMatthew G. Knepley   PetscErrorCode  ierr;
12299a800dd8SMatthew G. Knepley 
12309a800dd8SMatthew G. Knepley   PetscFunctionBegin;
12319a800dd8SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
12329a800dd8SMatthew G. Knepley   ierr = PetscFEGetQuadrature(fe[0], &q);CHKERRQ(ierr);
12339a800dd8SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
12349a800dd8SMatthew G. Knepley   ierr = PetscMalloc(numFields * sizeof(PetscDualSpace), &sp);CHKERRQ(ierr);
12359a800dd8SMatthew G. Knepley   for (f = 0; f < numFields; ++f) {
12369a800dd8SMatthew G. Knepley     ierr = PetscFEGetDualSpace(fe[f], &sp[f]);CHKERRQ(ierr);
12379a800dd8SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr);
12389a800dd8SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
12399a800dd8SMatthew G. Knepley     totDim += spDim*numComp;
12409a800dd8SMatthew G. Knepley   }
12419a800dd8SMatthew G. Knepley   ierr = DMDAGetInfo(dm, &dim,0,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
12429a800dd8SMatthew G. Knepley   ierr = DMDAGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
12439a800dd8SMatthew G. Knepley   ierr = DMDAVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr);
12449a800dd8SMatthew G. Knepley   if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section cell closure size %d != dual space dimension %d", numValues, totDim);
12459a800dd8SMatthew G. Knepley   ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
124621454ff5SMatthew G. Knepley   ierr = PetscQuadratureGetData(q, NULL, &numPoints, NULL, NULL);CHKERRQ(ierr);
124721454ff5SMatthew G. Knepley   ierr = PetscMalloc3(dim*numPoints,&v0,dim*dim*numPoints,&J,numPoints,&detJ);CHKERRQ(ierr);
12489a800dd8SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
12499a800dd8SMatthew G. Knepley     PetscCellGeometry geom;
12509a800dd8SMatthew G. Knepley 
125121454ff5SMatthew G. Knepley     ierr = DMDAComputeCellGeometry(dm, c, q, v0, J, NULL, detJ);CHKERRQ(ierr);
12529a800dd8SMatthew G. Knepley     geom.v0   = v0;
12539a800dd8SMatthew G. Knepley     geom.J    = J;
12549a800dd8SMatthew G. Knepley     geom.detJ = detJ;
12559a800dd8SMatthew G. Knepley     for (f = 0, v = 0; f < numFields; ++f) {
1256c110b1eeSGeoffrey Irving       void * const ctx = ctxs ? ctxs[f] : NULL;
12579a800dd8SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr);
12589a800dd8SMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
12599a800dd8SMatthew G. Knepley       for (d = 0; d < spDim; ++d) {
1260c110b1eeSGeoffrey Irving         ierr = PetscDualSpaceApply(sp[f], d, geom, numComp, funcs[f], ctx, &values[v]);CHKERRQ(ierr);
12619a800dd8SMatthew G. Knepley         v += numComp;
12629a800dd8SMatthew G. Knepley       }
12639a800dd8SMatthew G. Knepley     }
12649a800dd8SMatthew G. Knepley     ierr = DMDAVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr);
12659a800dd8SMatthew G. Knepley   }
12669a800dd8SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
12679a800dd8SMatthew G. Knepley   ierr = PetscFree3(v0,J,detJ);CHKERRQ(ierr);
12689a800dd8SMatthew G. Knepley   ierr = PetscFree(sp);CHKERRQ(ierr);
12699a800dd8SMatthew G. Knepley   PetscFunctionReturn(0);
12709a800dd8SMatthew G. Knepley }
12719a800dd8SMatthew G. Knepley 
12729a800dd8SMatthew G. Knepley #undef __FUNCT__
12739a800dd8SMatthew G. Knepley #define __FUNCT__ "DMDAProjectFunction"
12749a800dd8SMatthew G. Knepley /*@C
12759a800dd8SMatthew G. Knepley   DMDAProjectFunction - This projects the given function into the function space provided.
12769a800dd8SMatthew G. Knepley 
12779a800dd8SMatthew G. Knepley   Input Parameters:
12789a800dd8SMatthew G. Knepley + dm      - The DM
12799a800dd8SMatthew G. Knepley . fe      - The PetscFE associated with the field
12809a800dd8SMatthew G. Knepley . funcs   - The coordinate functions to evaluate
1281c110b1eeSGeoffrey Irving . ctxs    - Optional array of contexts to pass to each coordinate function.  ctxs itself may be null.
12829a800dd8SMatthew G. Knepley - mode    - The insertion mode for values
12839a800dd8SMatthew G. Knepley 
12849a800dd8SMatthew G. Knepley   Output Parameter:
12859a800dd8SMatthew G. Knepley . X - vector
12869a800dd8SMatthew G. Knepley 
12879a800dd8SMatthew G. Knepley   Level: developer
12889a800dd8SMatthew G. Knepley 
12899a800dd8SMatthew G. Knepley .seealso: DMDAComputeL2Diff()
12909a800dd8SMatthew G. Knepley @*/
1291c110b1eeSGeoffrey Irving PetscErrorCode DMDAProjectFunction(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec X)
12929a800dd8SMatthew G. Knepley {
12939a800dd8SMatthew G. Knepley   Vec            localX;
12949a800dd8SMatthew G. Knepley   PetscErrorCode ierr;
12959a800dd8SMatthew G. Knepley 
12969a800dd8SMatthew G. Knepley   PetscFunctionBegin;
12979a800dd8SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
12989a800dd8SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
1299c110b1eeSGeoffrey Irving   ierr = DMDAProjectFunctionLocal(dm, fe, funcs, ctxs, mode, localX);CHKERRQ(ierr);
13009a800dd8SMatthew G. Knepley   ierr = DMLocalToGlobalBegin(dm, localX, mode, X);CHKERRQ(ierr);
13019a800dd8SMatthew G. Knepley   ierr = DMLocalToGlobalEnd(dm, localX, mode, X);CHKERRQ(ierr);
13029a800dd8SMatthew G. Knepley   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
13039a800dd8SMatthew G. Knepley   PetscFunctionReturn(0);
13049a800dd8SMatthew G. Knepley }
13059a800dd8SMatthew G. Knepley 
13069a800dd8SMatthew G. Knepley #undef __FUNCT__
13079a800dd8SMatthew G. Knepley #define __FUNCT__ "DMDAComputeL2Diff"
13089a800dd8SMatthew G. Knepley /*@C
13099a800dd8SMatthew G. Knepley   DMDAComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h.
13109a800dd8SMatthew G. Knepley 
13119a800dd8SMatthew G. Knepley   Input Parameters:
13129a800dd8SMatthew G. Knepley + dm    - The DM
13139a800dd8SMatthew G. Knepley . fe    - The PetscFE object for each field
13149a800dd8SMatthew G. Knepley . funcs - The functions to evaluate for each field component
1315c110b1eeSGeoffrey Irving . ctxs  - Optional array of contexts to pass to each coordinate function.  ctxs itself may be null.
13169a800dd8SMatthew G. Knepley - X     - The coefficient vector u_h
13179a800dd8SMatthew G. Knepley 
13189a800dd8SMatthew G. Knepley   Output Parameter:
13199a800dd8SMatthew G. Knepley . diff - The diff ||u - u_h||_2
13209a800dd8SMatthew G. Knepley 
13219a800dd8SMatthew G. Knepley   Level: developer
13229a800dd8SMatthew G. Knepley 
13239a800dd8SMatthew G. Knepley .seealso: DMDAProjectFunction()
13249a800dd8SMatthew G. Knepley @*/
1325c110b1eeSGeoffrey Irving PetscErrorCode DMDAComputeL2Diff(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)
13269a800dd8SMatthew G. Knepley {
13279a800dd8SMatthew G. Knepley   const PetscInt  debug = 0;
13289a800dd8SMatthew G. Knepley   PetscSection    section;
13299a800dd8SMatthew G. Knepley   PetscQuadrature quad;
13309a800dd8SMatthew G. Knepley   Vec             localX;
13319a800dd8SMatthew G. Knepley   PetscScalar    *funcVal;
13329a800dd8SMatthew G. Knepley   PetscReal      *coords, *v0, *J, *invJ, detJ;
13339a800dd8SMatthew G. Knepley   PetscReal       localDiff = 0.0;
13349a800dd8SMatthew G. Knepley   PetscInt        dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp;
13359a800dd8SMatthew G. Knepley   PetscErrorCode  ierr;
13369a800dd8SMatthew G. Knepley 
13379a800dd8SMatthew G. Knepley   PetscFunctionBegin;
13389a800dd8SMatthew G. Knepley   ierr = DMDAGetInfo(dm, &dim,0,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
13399a800dd8SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
13409a800dd8SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
13419a800dd8SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
13429a800dd8SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
13439a800dd8SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
13449a800dd8SMatthew G. Knepley   for (field = 0; field < numFields; ++field) {
13459a800dd8SMatthew G. Knepley     PetscInt Nc;
13469a800dd8SMatthew G. Knepley 
13479a800dd8SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr);
13489a800dd8SMatthew G. Knepley     numComponents += Nc;
13499a800dd8SMatthew G. Knepley   }
13509a800dd8SMatthew G. Knepley   /* There are no BC values in DAs right now: ierr = DMDAProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */
1351d7a12f93SMatthew G. Knepley   ierr = PetscMalloc5(numComponents,&funcVal,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
13529a800dd8SMatthew G. Knepley   ierr = DMDAGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
13539a800dd8SMatthew G. Knepley   ierr = PetscFEGetQuadrature(fe[0], &quad);CHKERRQ(ierr);
13549a800dd8SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
13559a800dd8SMatthew G. Knepley     PetscScalar *x = NULL;
13569a800dd8SMatthew G. Knepley     PetscReal    elemDiff = 0.0;
13579a800dd8SMatthew G. Knepley 
135821454ff5SMatthew G. Knepley     ierr = DMDAComputeCellGeometry(dm, c, quad, v0, J, invJ, &detJ);CHKERRQ(ierr);
13599a800dd8SMatthew G. Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
13609a800dd8SMatthew G. Knepley     ierr = DMDAVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
13619a800dd8SMatthew G. Knepley 
13629a800dd8SMatthew G. Knepley     for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) {
1363c110b1eeSGeoffrey Irving       void * const ctx = ctxs ? ctxs[field] : NULL;
136421454ff5SMatthew G. Knepley       const PetscReal *quadPoints, *quadWeights;
13659a800dd8SMatthew G. Knepley       PetscReal       *basis;
136621454ff5SMatthew G. Knepley       PetscInt         numQuadPoints, numBasisFuncs, numBasisComps, q, d, e, fc, f;
13679a800dd8SMatthew G. Knepley 
136821454ff5SMatthew G. Knepley       ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
13699a800dd8SMatthew G. Knepley       ierr = PetscFEGetDimension(fe[field], &numBasisFuncs);CHKERRQ(ierr);
13709a800dd8SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe[field], &numBasisComps);CHKERRQ(ierr);
13719a800dd8SMatthew G. Knepley       ierr = PetscFEGetDefaultTabulation(fe[field], &basis, NULL, NULL);CHKERRQ(ierr);
13729a800dd8SMatthew G. Knepley       if (debug) {
13739a800dd8SMatthew G. Knepley         char title[1024];
13749a800dd8SMatthew G. Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
13759a800dd8SMatthew G. Knepley         ierr = DMPrintCellVector(c, title, numBasisFuncs*numBasisComps, &x[fieldOffset]);CHKERRQ(ierr);
13769a800dd8SMatthew G. Knepley       }
13779a800dd8SMatthew G. Knepley       for (q = 0; q < numQuadPoints; ++q) {
13789a800dd8SMatthew G. Knepley         for (d = 0; d < dim; d++) {
13799a800dd8SMatthew G. Knepley           coords[d] = v0[d];
13809a800dd8SMatthew G. Knepley           for (e = 0; e < dim; e++) {
13819a800dd8SMatthew G. Knepley             coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0);
13829a800dd8SMatthew G. Knepley           }
13839a800dd8SMatthew G. Knepley         }
1384c110b1eeSGeoffrey Irving         (*funcs[field])(coords, funcVal, ctx);
13859a800dd8SMatthew G. Knepley         for (fc = 0; fc < numBasisComps; ++fc) {
13869a800dd8SMatthew G. Knepley           PetscScalar interpolant = 0.0;
13879a800dd8SMatthew G. Knepley 
13889a800dd8SMatthew G. Knepley           for (f = 0; f < numBasisFuncs; ++f) {
13899a800dd8SMatthew G. Knepley             const PetscInt fidx = f*numBasisComps+fc;
13909a800dd8SMatthew G. Knepley             interpolant += x[fieldOffset+fidx]*basis[q*numBasisFuncs*numBasisComps+fidx];
13919a800dd8SMatthew G. Knepley           }
13929a800dd8SMatthew G. Knepley           if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "    elem %d field %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ);CHKERRQ(ierr);}
13939a800dd8SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ;
13949a800dd8SMatthew G. Knepley         }
13959a800dd8SMatthew G. Knepley       }
13969a800dd8SMatthew G. Knepley       comp        += numBasisComps;
13979a800dd8SMatthew G. Knepley       fieldOffset += numBasisFuncs*numBasisComps;
13989a800dd8SMatthew G. Knepley     }
13999a800dd8SMatthew G. Knepley     ierr = DMDAVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
14009a800dd8SMatthew G. Knepley     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
14019a800dd8SMatthew G. Knepley     localDiff += elemDiff;
14029a800dd8SMatthew G. Knepley   }
14039a800dd8SMatthew G. Knepley   ierr  = PetscFree5(funcVal,coords,v0,J,invJ);CHKERRQ(ierr);
14049a800dd8SMatthew G. Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
14059a800dd8SMatthew G. Knepley   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
14069a800dd8SMatthew G. Knepley   *diff = PetscSqrtReal(*diff);
1407a66d4d66SMatthew G Knepley   PetscFunctionReturn(0);
1408a66d4d66SMatthew G Knepley }
1409a66d4d66SMatthew G Knepley 
141047c6ae99SBarry Smith /* ------------------------------------------------------------------- */
141147c6ae99SBarry Smith 
141247c6ae99SBarry Smith #undef __FUNCT__
1413aa219208SBarry Smith #define __FUNCT__ "DMDAGetArray"
141447c6ae99SBarry Smith /*@C
1415aa219208SBarry Smith      DMDAGetArray - Gets a work array for a DMDA
141647c6ae99SBarry Smith 
141747c6ae99SBarry Smith     Input Parameter:
141847c6ae99SBarry Smith +    da - information about my local patch
141947c6ae99SBarry Smith -    ghosted - do you want arrays for the ghosted or nonghosted patch
142047c6ae99SBarry Smith 
142147c6ae99SBarry Smith     Output Parameters:
142247c6ae99SBarry Smith .    vptr - array data structured
142347c6ae99SBarry Smith 
142447c6ae99SBarry Smith     Note:  The vector values are NOT initialized and may have garbage in them, so you may need
142547c6ae99SBarry Smith            to zero them.
142647c6ae99SBarry Smith 
142747c6ae99SBarry Smith   Level: advanced
142847c6ae99SBarry Smith 
1429bcaeba4dSBarry Smith .seealso: DMDARestoreArray()
143047c6ae99SBarry Smith 
143147c6ae99SBarry Smith @*/
14327087cfbeSBarry Smith PetscErrorCode  DMDAGetArray(DM da,PetscBool ghosted,void *vptr)
143347c6ae99SBarry Smith {
143447c6ae99SBarry Smith   PetscErrorCode ierr;
143547c6ae99SBarry Smith   PetscInt       j,i,xs,ys,xm,ym,zs,zm;
143647c6ae99SBarry Smith   char           *iarray_start;
143747c6ae99SBarry Smith   void           **iptr = (void**)vptr;
143847c6ae99SBarry Smith   DM_DA          *dd    = (DM_DA*)da->data;
143947c6ae99SBarry Smith 
144047c6ae99SBarry Smith   PetscFunctionBegin;
144147c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
144247c6ae99SBarry Smith   if (ghosted) {
1443aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
144447c6ae99SBarry Smith       if (dd->arrayghostedin[i]) {
144547c6ae99SBarry Smith         *iptr                 = dd->arrayghostedin[i];
144647c6ae99SBarry Smith         iarray_start          = (char*)dd->startghostedin[i];
14470298fd71SBarry Smith         dd->arrayghostedin[i] = NULL;
14480298fd71SBarry Smith         dd->startghostedin[i] = NULL;
144947c6ae99SBarry Smith 
145047c6ae99SBarry Smith         goto done;
145147c6ae99SBarry Smith       }
145247c6ae99SBarry Smith     }
145347c6ae99SBarry Smith     xs = dd->Xs;
145447c6ae99SBarry Smith     ys = dd->Ys;
145547c6ae99SBarry Smith     zs = dd->Zs;
145647c6ae99SBarry Smith     xm = dd->Xe-dd->Xs;
145747c6ae99SBarry Smith     ym = dd->Ye-dd->Ys;
145847c6ae99SBarry Smith     zm = dd->Ze-dd->Zs;
145947c6ae99SBarry Smith   } else {
1460aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
146147c6ae99SBarry Smith       if (dd->arrayin[i]) {
146247c6ae99SBarry Smith         *iptr          = dd->arrayin[i];
146347c6ae99SBarry Smith         iarray_start   = (char*)dd->startin[i];
14640298fd71SBarry Smith         dd->arrayin[i] = NULL;
14650298fd71SBarry Smith         dd->startin[i] = NULL;
146647c6ae99SBarry Smith 
146747c6ae99SBarry Smith         goto done;
146847c6ae99SBarry Smith       }
146947c6ae99SBarry Smith     }
147047c6ae99SBarry Smith     xs = dd->xs;
147147c6ae99SBarry Smith     ys = dd->ys;
147247c6ae99SBarry Smith     zs = dd->zs;
147347c6ae99SBarry Smith     xm = dd->xe-dd->xs;
147447c6ae99SBarry Smith     ym = dd->ye-dd->ys;
147547c6ae99SBarry Smith     zm = dd->ze-dd->zs;
147647c6ae99SBarry Smith   }
147747c6ae99SBarry Smith 
147847c6ae99SBarry Smith   switch (dd->dim) {
147947c6ae99SBarry Smith   case 1: {
148047c6ae99SBarry Smith     void *ptr;
148147c6ae99SBarry Smith 
1482901f1932SJed Brown     ierr = PetscMalloc(xm*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
148347c6ae99SBarry Smith 
148447c6ae99SBarry Smith     ptr   = (void*)(iarray_start - xs*sizeof(PetscScalar));
148547c6ae99SBarry Smith     *iptr = (void*)ptr;
14868865f1eaSKarl Rupp     break;
14878865f1eaSKarl Rupp   }
148847c6ae99SBarry Smith   case 2: {
148947c6ae99SBarry Smith     void **ptr;
149047c6ae99SBarry Smith 
1491901f1932SJed Brown     ierr = PetscMalloc((ym+1)*sizeof(void*)+xm*ym*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
149247c6ae99SBarry Smith 
149347c6ae99SBarry Smith     ptr = (void**)(iarray_start + xm*ym*sizeof(PetscScalar) - ys*sizeof(void*));
14948865f1eaSKarl Rupp     for (j=ys; j<ys+ym; j++) ptr[j] = iarray_start + sizeof(PetscScalar)*(xm*(j-ys) - xs);
149547c6ae99SBarry Smith     *iptr = (void*)ptr;
14968865f1eaSKarl Rupp     break;
14978865f1eaSKarl Rupp   }
149847c6ae99SBarry Smith   case 3: {
149947c6ae99SBarry Smith     void ***ptr,**bptr;
150047c6ae99SBarry Smith 
1501901f1932SJed Brown     ierr = PetscMalloc((zm+1)*sizeof(void**)+(ym*zm+1)*sizeof(void*)+xm*ym*zm*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
150247c6ae99SBarry Smith 
150347c6ae99SBarry Smith     ptr  = (void***)(iarray_start + xm*ym*zm*sizeof(PetscScalar) - zs*sizeof(void*));
150447c6ae99SBarry Smith     bptr = (void**)(iarray_start + xm*ym*zm*sizeof(PetscScalar) + zm*sizeof(void**));
15058865f1eaSKarl Rupp     for (i=zs; i<zs+zm; i++) ptr[i] = bptr + ((i-zs)*ym - ys);
150647c6ae99SBarry Smith     for (i=zs; i<zs+zm; i++) {
150747c6ae99SBarry Smith       for (j=ys; j<ys+ym; j++) {
150847c6ae99SBarry Smith         ptr[i][j] = iarray_start + sizeof(PetscScalar)*(xm*ym*(i-zs) + xm*(j-ys) - xs);
150947c6ae99SBarry Smith       }
151047c6ae99SBarry Smith     }
151147c6ae99SBarry Smith     *iptr = (void*)ptr;
15128865f1eaSKarl Rupp     break;
15138865f1eaSKarl Rupp   }
151447c6ae99SBarry Smith   default:
1515ce94432eSBarry Smith     SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Dimension %D not supported",dd->dim);
151647c6ae99SBarry Smith   }
151747c6ae99SBarry Smith 
151847c6ae99SBarry Smith done:
151947c6ae99SBarry Smith   /* add arrays to the checked out list */
152047c6ae99SBarry Smith   if (ghosted) {
1521aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
152247c6ae99SBarry Smith       if (!dd->arrayghostedout[i]) {
152347c6ae99SBarry Smith         dd->arrayghostedout[i] = *iptr;
152447c6ae99SBarry Smith         dd->startghostedout[i] = iarray_start;
152547c6ae99SBarry Smith         break;
152647c6ae99SBarry Smith       }
152747c6ae99SBarry Smith     }
152847c6ae99SBarry Smith   } else {
1529aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
153047c6ae99SBarry Smith       if (!dd->arrayout[i]) {
153147c6ae99SBarry Smith         dd->arrayout[i] = *iptr;
153247c6ae99SBarry Smith         dd->startout[i] = iarray_start;
153347c6ae99SBarry Smith         break;
153447c6ae99SBarry Smith       }
153547c6ae99SBarry Smith     }
153647c6ae99SBarry Smith   }
153747c6ae99SBarry Smith   PetscFunctionReturn(0);
153847c6ae99SBarry Smith }
153947c6ae99SBarry Smith 
154047c6ae99SBarry Smith #undef __FUNCT__
1541aa219208SBarry Smith #define __FUNCT__ "DMDARestoreArray"
154247c6ae99SBarry Smith /*@C
1543aa219208SBarry Smith      DMDARestoreArray - Restores an array of derivative types for a DMDA
154447c6ae99SBarry Smith 
154547c6ae99SBarry Smith     Input Parameter:
154647c6ae99SBarry Smith +    da - information about my local patch
154747c6ae99SBarry Smith .    ghosted - do you want arrays for the ghosted or nonghosted patch
154847c6ae99SBarry Smith -    vptr - array data structured to be passed to ad_FormFunctionLocal()
154947c6ae99SBarry Smith 
155047c6ae99SBarry Smith      Level: advanced
155147c6ae99SBarry Smith 
1552bcaeba4dSBarry Smith .seealso: DMDAGetArray()
155347c6ae99SBarry Smith 
155447c6ae99SBarry Smith @*/
15557087cfbeSBarry Smith PetscErrorCode  DMDARestoreArray(DM da,PetscBool ghosted,void *vptr)
155647c6ae99SBarry Smith {
155747c6ae99SBarry Smith   PetscInt i;
155847c6ae99SBarry Smith   void     **iptr = (void**)vptr,*iarray_start = 0;
155947c6ae99SBarry Smith   DM_DA    *dd    = (DM_DA*)da->data;
156047c6ae99SBarry Smith 
156147c6ae99SBarry Smith   PetscFunctionBegin;
156247c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
156347c6ae99SBarry Smith   if (ghosted) {
1564aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
156547c6ae99SBarry Smith       if (dd->arrayghostedout[i] == *iptr) {
156647c6ae99SBarry Smith         iarray_start           = dd->startghostedout[i];
15670298fd71SBarry Smith         dd->arrayghostedout[i] = NULL;
15680298fd71SBarry Smith         dd->startghostedout[i] = NULL;
156947c6ae99SBarry Smith         break;
157047c6ae99SBarry Smith       }
157147c6ae99SBarry Smith     }
1572aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
157347c6ae99SBarry Smith       if (!dd->arrayghostedin[i]) {
157447c6ae99SBarry Smith         dd->arrayghostedin[i] = *iptr;
157547c6ae99SBarry Smith         dd->startghostedin[i] = iarray_start;
157647c6ae99SBarry Smith         break;
157747c6ae99SBarry Smith       }
157847c6ae99SBarry Smith     }
157947c6ae99SBarry Smith   } else {
1580aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
158147c6ae99SBarry Smith       if (dd->arrayout[i] == *iptr) {
158247c6ae99SBarry Smith         iarray_start    = dd->startout[i];
15830298fd71SBarry Smith         dd->arrayout[i] = NULL;
15840298fd71SBarry Smith         dd->startout[i] = NULL;
158547c6ae99SBarry Smith         break;
158647c6ae99SBarry Smith       }
158747c6ae99SBarry Smith     }
1588aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
158947c6ae99SBarry Smith       if (!dd->arrayin[i]) {
159047c6ae99SBarry Smith         dd->arrayin[i] = *iptr;
159147c6ae99SBarry Smith         dd->startin[i] = iarray_start;
159247c6ae99SBarry Smith         break;
159347c6ae99SBarry Smith       }
159447c6ae99SBarry Smith     }
159547c6ae99SBarry Smith   }
159647c6ae99SBarry Smith   PetscFunctionReturn(0);
159747c6ae99SBarry Smith }
159847c6ae99SBarry Smith 
1599