xref: /petsc/src/dm/impls/da/dalocal.c (revision affa55c76964d8ba770af915a7f2fed557cdab72)
147c6ae99SBarry Smith 
247c6ae99SBarry Smith /*
347c6ae99SBarry Smith   Code for manipulating distributed regular arrays in parallel.
447c6ae99SBarry Smith */
547c6ae99SBarry Smith 
64035e84dSBarry Smith #include <petsc-private/dmdaimpl.h>    /*I   "petscdmda.h"   I*/
70c312b8eSJed Brown #include <petscsf.h>
847c6ae99SBarry Smith 
947c6ae99SBarry Smith /*
10e3c5b3baSBarry Smith    This allows the DMDA vectors to properly tell MATLAB their dimensions
1147c6ae99SBarry Smith */
1247c6ae99SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
13c6db04a5SJed Brown #include <engine.h>   /* MATLAB include file */
14c6db04a5SJed Brown #include <mex.h>      /* MATLAB include file */
1547c6ae99SBarry Smith #undef __FUNCT__
1647c6ae99SBarry Smith #define __FUNCT__ "VecMatlabEnginePut_DA2d"
171e6b0712SBarry Smith static PetscErrorCode  VecMatlabEnginePut_DA2d(PetscObject obj,void *mengine)
1847c6ae99SBarry Smith {
1947c6ae99SBarry Smith   PetscErrorCode ierr;
2047c6ae99SBarry Smith   PetscInt       n,m;
2147c6ae99SBarry Smith   Vec            vec = (Vec)obj;
2247c6ae99SBarry Smith   PetscScalar    *array;
2347c6ae99SBarry Smith   mxArray        *mat;
249a42bb27SBarry Smith   DM             da;
2547c6ae99SBarry Smith 
2647c6ae99SBarry Smith   PetscFunctionBegin;
27c688c046SMatthew G Knepley   ierr = VecGetDM(vec, &da);CHKERRQ(ierr);
28ce94432eSBarry Smith   if (!da) SETERRQ(PetscObjectComm((PetscObject)vec),PETSC_ERR_ARG_WRONGSTATE,"Vector not associated with a DMDA");
29aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,0,0,0,&m,&n,0);CHKERRQ(ierr);
3047c6ae99SBarry Smith 
3147c6ae99SBarry Smith   ierr = VecGetArray(vec,&array);CHKERRQ(ierr);
3247c6ae99SBarry Smith #if !defined(PETSC_USE_COMPLEX)
3347c6ae99SBarry Smith   mat = mxCreateDoubleMatrix(m,n,mxREAL);
3447c6ae99SBarry Smith #else
3547c6ae99SBarry Smith   mat = mxCreateDoubleMatrix(m,n,mxCOMPLEX);
3647c6ae99SBarry Smith #endif
3747c6ae99SBarry Smith   ierr = PetscMemcpy(mxGetPr(mat),array,n*m*sizeof(PetscScalar));CHKERRQ(ierr);
3847c6ae99SBarry Smith   ierr = PetscObjectName(obj);CHKERRQ(ierr);
3947c6ae99SBarry Smith   engPutVariable((Engine*)mengine,obj->name,mat);
4047c6ae99SBarry Smith 
4147c6ae99SBarry Smith   ierr = VecRestoreArray(vec,&array);CHKERRQ(ierr);
4247c6ae99SBarry Smith   PetscFunctionReturn(0);
4347c6ae99SBarry Smith }
4447c6ae99SBarry Smith #endif
4547c6ae99SBarry Smith 
4647c6ae99SBarry Smith 
4747c6ae99SBarry Smith #undef __FUNCT__
48564755cdSBarry Smith #define __FUNCT__ "DMCreateLocalVector_DA"
497087cfbeSBarry Smith PetscErrorCode  DMCreateLocalVector_DA(DM da,Vec *g)
5047c6ae99SBarry Smith {
5147c6ae99SBarry Smith   PetscErrorCode ierr;
5247c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
5347c6ae99SBarry Smith 
5447c6ae99SBarry Smith   PetscFunctionBegin;
5547c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
5647c6ae99SBarry Smith   PetscValidPointer(g,2);
5711689aebSJed Brown   if (da->defaultSection) {
5811689aebSJed Brown     ierr = DMCreateLocalVector_Section_Private(da,g);CHKERRQ(ierr);
5911689aebSJed Brown   } else {
6047c6ae99SBarry Smith     ierr = VecCreate(PETSC_COMM_SELF,g);CHKERRQ(ierr);
6147c6ae99SBarry Smith     ierr = VecSetSizes(*g,dd->nlocal,PETSC_DETERMINE);CHKERRQ(ierr);
6247c6ae99SBarry Smith     ierr = VecSetBlockSize(*g,dd->w);CHKERRQ(ierr);
63401ddaa8SBarry Smith     ierr = VecSetType(*g,da->vectype);CHKERRQ(ierr);
64c688c046SMatthew G Knepley     ierr = VecSetDM(*g, da);CHKERRQ(ierr);
6547c6ae99SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
6647c6ae99SBarry Smith     if (dd->w == 1  && dd->dim == 2) {
67bdf89e91SBarry Smith       ierr = PetscObjectComposeFunction((PetscObject)*g,"PetscMatlabEnginePut_C",VecMatlabEnginePut_DA2d);CHKERRQ(ierr);
6847c6ae99SBarry Smith     }
6947c6ae99SBarry Smith #endif
7011689aebSJed Brown   }
7147c6ae99SBarry Smith   PetscFunctionReturn(0);
7247c6ae99SBarry Smith }
7347c6ae99SBarry Smith 
74a66d4d66SMatthew G Knepley #undef __FUNCT__
7557459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetNumCells"
763582350cSMatthew G. Knepley PetscErrorCode DMDAGetNumCells(DM dm, PetscInt *numCellsX, PetscInt *numCellsY, PetscInt *numCellsZ, PetscInt *numCells)
7757459e9aSMatthew G Knepley {
7857459e9aSMatthew G Knepley   DM_DA          *da = (DM_DA*) dm->data;
7957459e9aSMatthew G Knepley   const PetscInt dim = da->dim;
8057459e9aSMatthew G Knepley   const PetscInt mx  = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs;
8157459e9aSMatthew G Knepley   const PetscInt nC  = (mx)*(dim > 1 ? (my)*(dim > 2 ? (mz) : 1) : 1);
8257459e9aSMatthew G Knepley 
8357459e9aSMatthew G Knepley   PetscFunctionBegin;
843582350cSMatthew G. Knepley   if (numCellsX) {
853582350cSMatthew G. Knepley     PetscValidIntPointer(numCellsX,2);
863582350cSMatthew G. Knepley     *numCellsX = mx;
873582350cSMatthew G. Knepley   }
883582350cSMatthew G. Knepley   if (numCellsY) {
893582350cSMatthew G. Knepley     PetscValidIntPointer(numCellsX,3);
903582350cSMatthew G. Knepley     *numCellsY = my;
913582350cSMatthew G. Knepley   }
923582350cSMatthew G. Knepley   if (numCellsZ) {
933582350cSMatthew G. Knepley     PetscValidIntPointer(numCellsX,4);
943582350cSMatthew G. Knepley     *numCellsZ = mz;
953582350cSMatthew G. Knepley   }
9657459e9aSMatthew G Knepley   if (numCells) {
973582350cSMatthew G. Knepley     PetscValidIntPointer(numCells,5);
9857459e9aSMatthew G Knepley     *numCells = nC;
9957459e9aSMatthew G Knepley   }
10057459e9aSMatthew G Knepley   PetscFunctionReturn(0);
10157459e9aSMatthew G Knepley }
10257459e9aSMatthew G Knepley 
10357459e9aSMatthew G Knepley #undef __FUNCT__
10457459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetNumVertices"
10557459e9aSMatthew G Knepley PetscErrorCode DMDAGetNumVertices(DM dm, PetscInt *numVerticesX, PetscInt *numVerticesY, PetscInt *numVerticesZ, PetscInt *numVertices)
10657459e9aSMatthew G Knepley {
10757459e9aSMatthew G Knepley   DM_DA          *da = (DM_DA*) dm->data;
10857459e9aSMatthew G Knepley   const PetscInt dim = da->dim;
10957459e9aSMatthew G Knepley   const PetscInt mx  = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs;
11057459e9aSMatthew G Knepley   const PetscInt nVx = mx+1;
11157459e9aSMatthew G Knepley   const PetscInt nVy = dim > 1 ? (my+1) : 1;
11257459e9aSMatthew G Knepley   const PetscInt nVz = dim > 2 ? (mz+1) : 1;
11357459e9aSMatthew G Knepley   const PetscInt nV  = nVx*nVy*nVz;
11457459e9aSMatthew G Knepley 
11557459e9aSMatthew G Knepley   PetscFunctionBegin;
11657459e9aSMatthew G Knepley   if (numVerticesX) {
11757459e9aSMatthew G Knepley     PetscValidIntPointer(numVerticesX,2);
11857459e9aSMatthew G Knepley     *numVerticesX = nVx;
11957459e9aSMatthew G Knepley   }
12057459e9aSMatthew G Knepley   if (numVerticesY) {
12157459e9aSMatthew G Knepley     PetscValidIntPointer(numVerticesY,3);
12257459e9aSMatthew G Knepley     *numVerticesY = nVy;
12357459e9aSMatthew G Knepley   }
12457459e9aSMatthew G Knepley   if (numVerticesZ) {
12557459e9aSMatthew G Knepley     PetscValidIntPointer(numVerticesZ,4);
12657459e9aSMatthew G Knepley     *numVerticesZ = nVz;
12757459e9aSMatthew G Knepley   }
12857459e9aSMatthew G Knepley   if (numVertices) {
12957459e9aSMatthew G Knepley     PetscValidIntPointer(numVertices,5);
13057459e9aSMatthew G Knepley     *numVertices = nV;
13157459e9aSMatthew G Knepley   }
13257459e9aSMatthew G Knepley   PetscFunctionReturn(0);
13357459e9aSMatthew G Knepley }
13457459e9aSMatthew G Knepley 
13557459e9aSMatthew G Knepley #undef __FUNCT__
13657459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetNumFaces"
13757459e9aSMatthew G Knepley PetscErrorCode DMDAGetNumFaces(DM dm, PetscInt *numXFacesX, PetscInt *numXFaces, PetscInt *numYFacesY, PetscInt *numYFaces, PetscInt *numZFacesZ, PetscInt *numZFaces)
13857459e9aSMatthew G Knepley {
13957459e9aSMatthew G Knepley   DM_DA          *da = (DM_DA*) dm->data;
14057459e9aSMatthew G Knepley   const PetscInt dim = da->dim;
14157459e9aSMatthew G Knepley   const PetscInt mx  = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs;
14257459e9aSMatthew G Knepley   const PetscInt nxF = (dim > 1 ? (my)*(dim > 2 ? (mz) : 1) : 1);
14357459e9aSMatthew G Knepley   const PetscInt nXF = (mx+1)*nxF;
14457459e9aSMatthew G Knepley   const PetscInt nyF = mx*(dim > 2 ? mz : 1);
14557459e9aSMatthew G Knepley   const PetscInt nYF = dim > 1 ? (my+1)*nyF : 0;
14657459e9aSMatthew G Knepley   const PetscInt nzF = mx*(dim > 1 ? my : 0);
14757459e9aSMatthew G Knepley   const PetscInt nZF = dim > 2 ? (mz+1)*nzF : 0;
14857459e9aSMatthew G Knepley 
14957459e9aSMatthew G Knepley   PetscFunctionBegin;
15057459e9aSMatthew G Knepley   if (numXFacesX) {
15157459e9aSMatthew G Knepley     PetscValidIntPointer(numXFacesX,2);
15257459e9aSMatthew G Knepley     *numXFacesX = nxF;
15357459e9aSMatthew G Knepley   }
15457459e9aSMatthew G Knepley   if (numXFaces) {
15557459e9aSMatthew G Knepley     PetscValidIntPointer(numXFaces,3);
15657459e9aSMatthew G Knepley     *numXFaces = nXF;
15757459e9aSMatthew G Knepley   }
15857459e9aSMatthew G Knepley   if (numYFacesY) {
15957459e9aSMatthew G Knepley     PetscValidIntPointer(numYFacesY,4);
16057459e9aSMatthew G Knepley     *numYFacesY = nyF;
16157459e9aSMatthew G Knepley   }
16257459e9aSMatthew G Knepley   if (numYFaces) {
16357459e9aSMatthew G Knepley     PetscValidIntPointer(numYFaces,5);
16457459e9aSMatthew G Knepley     *numYFaces = nYF;
16557459e9aSMatthew G Knepley   }
16657459e9aSMatthew G Knepley   if (numZFacesZ) {
16757459e9aSMatthew G Knepley     PetscValidIntPointer(numZFacesZ,6);
16857459e9aSMatthew G Knepley     *numZFacesZ = nzF;
16957459e9aSMatthew G Knepley   }
17057459e9aSMatthew G Knepley   if (numZFaces) {
17157459e9aSMatthew G Knepley     PetscValidIntPointer(numZFaces,7);
17257459e9aSMatthew G Knepley     *numZFaces = nZF;
17357459e9aSMatthew G Knepley   }
17457459e9aSMatthew G Knepley   PetscFunctionReturn(0);
17557459e9aSMatthew G Knepley }
17657459e9aSMatthew G Knepley 
17757459e9aSMatthew G Knepley #undef __FUNCT__
17857459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetHeightStratum"
17957459e9aSMatthew G Knepley PetscErrorCode DMDAGetHeightStratum(DM dm, PetscInt height, PetscInt *pStart, PetscInt *pEnd)
18057459e9aSMatthew G Knepley {
18157459e9aSMatthew G Knepley   DM_DA          *da = (DM_DA*) dm->data;
18257459e9aSMatthew G Knepley   const PetscInt dim = da->dim;
18357459e9aSMatthew G Knepley   PetscInt       nC, nV, nXF, nYF, nZF;
18457459e9aSMatthew G Knepley   PetscErrorCode ierr;
18557459e9aSMatthew G Knepley 
18657459e9aSMatthew G Knepley   PetscFunctionBegin;
1878865f1eaSKarl Rupp   if (pStart) PetscValidIntPointer(pStart,3);
1888865f1eaSKarl Rupp   if (pEnd)   PetscValidIntPointer(pEnd,4);
1893582350cSMatthew G. Knepley   ierr = DMDAGetNumCells(dm, NULL, NULL, NULL, &nC);CHKERRQ(ierr);
1900298fd71SBarry Smith   ierr = DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV);CHKERRQ(ierr);
1910298fd71SBarry Smith   ierr = DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF);CHKERRQ(ierr);
19257459e9aSMatthew G Knepley   if (height == 0) {
19357459e9aSMatthew G Knepley     /* Cells */
1948865f1eaSKarl Rupp     if (pStart) *pStart = 0;
1958865f1eaSKarl Rupp     if (pEnd)   *pEnd   = nC;
19657459e9aSMatthew G Knepley   } else if (height == 1) {
19757459e9aSMatthew G Knepley     /* Faces */
1988865f1eaSKarl Rupp     if (pStart) *pStart = nC+nV;
1998865f1eaSKarl Rupp     if (pEnd)   *pEnd   = nC+nV+nXF+nYF+nZF;
20057459e9aSMatthew G Knepley   } else if (height == dim) {
20157459e9aSMatthew G Knepley     /* Vertices */
2028865f1eaSKarl Rupp     if (pStart) *pStart = nC;
2038865f1eaSKarl Rupp     if (pEnd)   *pEnd   = nC+nV;
20457459e9aSMatthew G Knepley   } else if (height < 0) {
20557459e9aSMatthew G Knepley     /* All points */
2068865f1eaSKarl Rupp     if (pStart) *pStart = 0;
2078865f1eaSKarl Rupp     if (pEnd)   *pEnd   = nC+nV+nXF+nYF+nZF;
20882f516ccSBarry Smith   } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No points of height %d in the DA", height);
20957459e9aSMatthew G Knepley   PetscFunctionReturn(0);
21057459e9aSMatthew G Knepley }
21157459e9aSMatthew G Knepley 
21257459e9aSMatthew G Knepley #undef __FUNCT__
2133385ff7eSMatthew G. Knepley #define __FUNCT__ "DMDAGetDepthStratum"
2143385ff7eSMatthew G. Knepley PetscErrorCode DMDAGetDepthStratum(DM dm, PetscInt depth, PetscInt *pStart, PetscInt *pEnd)
2153385ff7eSMatthew G. Knepley {
2163385ff7eSMatthew G. Knepley   DM_DA         *da  = (DM_DA*) dm->data;
2173385ff7eSMatthew G. Knepley   const PetscInt dim = da->dim;
2183385ff7eSMatthew G. Knepley   PetscInt       nC, nV, nXF, nYF, nZF;
2193385ff7eSMatthew G. Knepley   PetscErrorCode ierr;
2203385ff7eSMatthew G. Knepley 
2213385ff7eSMatthew G. Knepley   PetscFunctionBegin;
2223385ff7eSMatthew G. Knepley   if (pStart) PetscValidIntPointer(pStart,3);
2233385ff7eSMatthew G. Knepley   if (pEnd)   PetscValidIntPointer(pEnd,4);
2243385ff7eSMatthew G. Knepley   ierr = DMDAGetNumCells(dm, NULL, NULL, NULL, &nC);CHKERRQ(ierr);
2253385ff7eSMatthew G. Knepley   ierr = DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV);CHKERRQ(ierr);
2263385ff7eSMatthew G. Knepley   ierr = DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF);CHKERRQ(ierr);
2273385ff7eSMatthew G. Knepley   if (depth == dim) {
2283385ff7eSMatthew G. Knepley     /* Cells */
2293385ff7eSMatthew G. Knepley     if (pStart) *pStart = 0;
2303385ff7eSMatthew G. Knepley     if (pEnd)   *pEnd   = nC;
2313385ff7eSMatthew G. Knepley   } else if (depth == dim-1) {
2323385ff7eSMatthew G. Knepley     /* Faces */
2333385ff7eSMatthew G. Knepley     if (pStart) *pStart = nC+nV;
2343385ff7eSMatthew G. Knepley     if (pEnd)   *pEnd   = nC+nV+nXF+nYF+nZF;
2353385ff7eSMatthew G. Knepley   } else if (depth == 0) {
2363385ff7eSMatthew G. Knepley     /* Vertices */
2373385ff7eSMatthew G. Knepley     if (pStart) *pStart = nC;
2383385ff7eSMatthew G. Knepley     if (pEnd)   *pEnd   = nC+nV;
2393385ff7eSMatthew G. Knepley   } else if (depth < 0) {
2403385ff7eSMatthew G. Knepley     /* All points */
2413385ff7eSMatthew G. Knepley     if (pStart) *pStart = 0;
2423385ff7eSMatthew G. Knepley     if (pEnd)   *pEnd   = nC+nV+nXF+nYF+nZF;
2433385ff7eSMatthew G. Knepley   } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No points of depth %d in the DA", depth);
2443385ff7eSMatthew G. Knepley   PetscFunctionReturn(0);
2453385ff7eSMatthew G. Knepley }
2463385ff7eSMatthew G. Knepley 
2473385ff7eSMatthew G. Knepley #undef __FUNCT__
2483385ff7eSMatthew G. Knepley #define __FUNCT__ "DMDAGetConeSize"
2493385ff7eSMatthew G. Knepley PetscErrorCode DMDAGetConeSize(DM dm, PetscInt p, PetscInt *coneSize)
2503385ff7eSMatthew G. Knepley {
2513385ff7eSMatthew G. Knepley   DM_DA         *da  = (DM_DA*) dm->data;
2523385ff7eSMatthew G. Knepley   const PetscInt dim = da->dim;
2533385ff7eSMatthew G. Knepley   PetscInt       nC, nV, nXF, nYF, nZF;
2543385ff7eSMatthew G. Knepley   PetscErrorCode ierr;
2553385ff7eSMatthew G. Knepley 
2563385ff7eSMatthew G. Knepley   PetscFunctionBegin;
2573385ff7eSMatthew G. Knepley   *coneSize = 0;
2583385ff7eSMatthew G. Knepley   ierr = DMDAGetNumCells(dm, NULL, NULL, NULL, &nC);CHKERRQ(ierr);
2593385ff7eSMatthew G. Knepley   ierr = DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV);CHKERRQ(ierr);
2603385ff7eSMatthew G. Knepley   ierr = DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF);CHKERRQ(ierr);
2613385ff7eSMatthew G. Knepley   switch (dim) {
2623385ff7eSMatthew G. Knepley   case 2:
2633385ff7eSMatthew G. Knepley     if (p >= 0) {
2643385ff7eSMatthew G. Knepley       if (p < nC) {
2653385ff7eSMatthew G. Knepley         *coneSize = 4;
2663385ff7eSMatthew G. Knepley       } else if (p < nC+nV) {
2673385ff7eSMatthew G. Knepley         *coneSize = 0;
2683385ff7eSMatthew G. Knepley       } else if (p < nC+nV+nXF+nYF+nZF) {
2693385ff7eSMatthew G. Knepley         *coneSize = 2;
2703385ff7eSMatthew G. Knepley       } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %d should be in [0, %d)", p, nC+nV+nXF+nYF+nZF);
2713385ff7eSMatthew G. Knepley     } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative point %d is invalid", p);
2723385ff7eSMatthew G. Knepley     break;
2733385ff7eSMatthew G. Knepley   case 3:
2743385ff7eSMatthew G. Knepley     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do 3D");
2753385ff7eSMatthew G. Knepley     break;
2763385ff7eSMatthew G. Knepley   }
2773385ff7eSMatthew G. Knepley   PetscFunctionReturn(0);
2783385ff7eSMatthew G. Knepley }
2793385ff7eSMatthew G. Knepley 
2803385ff7eSMatthew G. Knepley #undef __FUNCT__
2813385ff7eSMatthew G. Knepley #define __FUNCT__ "DMDAGetCone"
2823385ff7eSMatthew G. Knepley PetscErrorCode DMDAGetCone(DM dm, PetscInt p, PetscInt *cone[])
2833385ff7eSMatthew G. Knepley {
2843385ff7eSMatthew G. Knepley   DM_DA         *da  = (DM_DA*) dm->data;
2853385ff7eSMatthew G. Knepley   const PetscInt dim = da->dim;
2863385ff7eSMatthew G. Knepley   PetscInt       nCx, nCy, nCz, nC, nVx, nVy, nVz, nV, nxF, nyF, nzF, nXF, nYF, nZF;
2873385ff7eSMatthew G. Knepley   PetscErrorCode ierr;
2883385ff7eSMatthew G. Knepley 
2893385ff7eSMatthew G. Knepley   PetscFunctionBegin;
2903385ff7eSMatthew G. Knepley   if (!cone) {ierr = DMGetWorkArray(dm, 6, PETSC_INT, cone);CHKERRQ(ierr);}
2913385ff7eSMatthew G. Knepley   ierr = DMDAGetNumCells(dm, &nCx, &nCy, &nCz, &nC);CHKERRQ(ierr);
2923385ff7eSMatthew G. Knepley   ierr = DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV);CHKERRQ(ierr);
2933385ff7eSMatthew G. Knepley   ierr = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr);
2943385ff7eSMatthew G. Knepley   switch (dim) {
2953385ff7eSMatthew G. Knepley   case 2:
2963385ff7eSMatthew G. Knepley     if (p >= 0) {
2973385ff7eSMatthew G. Knepley       if (p < nC) {
2983385ff7eSMatthew G. Knepley         const PetscInt cy = p / nCx;
2993385ff7eSMatthew G. Knepley         const PetscInt cx = p % nCx;
3003385ff7eSMatthew G. Knepley 
3013385ff7eSMatthew G. Knepley         (*cone)[0] = cy*nxF + cx + nC+nV;
3023385ff7eSMatthew G. Knepley         (*cone)[1] = cx*nyF + cy + nyF + nC+nV+nXF;
3033385ff7eSMatthew G. Knepley         (*cone)[2] = cy*nxF + cx + nxF + nC+nV;
3043385ff7eSMatthew G. Knepley         (*cone)[3] = cx*nyF + cy + nC+nV+nXF;
3053385ff7eSMatthew G. Knepley         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do cell cones");
3063385ff7eSMatthew G. Knepley       } else if (p < nC+nV) {
3073385ff7eSMatthew G. Knepley       } else if (p < nC+nV+nXF) {
3083385ff7eSMatthew G. Knepley         const PetscInt fy = (p - nC+nV) / nxF;
3093385ff7eSMatthew G. Knepley         const PetscInt fx = (p - nC+nV) % nxF;
3103385ff7eSMatthew G. Knepley 
3113385ff7eSMatthew G. Knepley         (*cone)[0] = fy*nVx + fx + nC;
3123385ff7eSMatthew G. Knepley         (*cone)[1] = fy*nVx + fx + 1 + nC;
3133385ff7eSMatthew G. Knepley       } else if (p < nC+nV+nXF+nYF) {
3143385ff7eSMatthew G. Knepley         const PetscInt fx = (p - nC+nV+nXF) / nyF;
3153385ff7eSMatthew G. Knepley         const PetscInt fy = (p - nC+nV+nXF) % nyF;
3163385ff7eSMatthew G. Knepley 
3173385ff7eSMatthew G. Knepley         (*cone)[0] = fy*nVx + fx + nC;
3183385ff7eSMatthew G. Knepley         (*cone)[1] = fy*nVx + fx + nVx + nC;
3193385ff7eSMatthew G. Knepley       } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %d should be in [0, %d)", p, nC+nV+nXF+nYF+nZF);
3203385ff7eSMatthew G. Knepley     } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative point %d is invalid", p);
3213385ff7eSMatthew G. Knepley     break;
3223385ff7eSMatthew G. Knepley   case 3:
3233385ff7eSMatthew G. Knepley     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do 3D");
3243385ff7eSMatthew G. Knepley     break;
3253385ff7eSMatthew G. Knepley   }
3263385ff7eSMatthew G. Knepley   PetscFunctionReturn(0);
3273385ff7eSMatthew G. Knepley }
3283385ff7eSMatthew G. Knepley 
3293385ff7eSMatthew G. Knepley #undef __FUNCT__
3303385ff7eSMatthew G. Knepley #define __FUNCT__ "DMDARestoreCone"
3313385ff7eSMatthew G. Knepley PetscErrorCode DMDARestoreCone(DM dm, PetscInt p, PetscInt *cone[])
3323385ff7eSMatthew G. Knepley {
3333385ff7eSMatthew G. Knepley   PetscErrorCode ierr;
3343385ff7eSMatthew G. Knepley 
3353385ff7eSMatthew G. Knepley   PetscFunctionBegin;
3363385ff7eSMatthew G. Knepley   ierr = DMGetWorkArray(dm, 6, PETSC_INT, cone);CHKERRQ(ierr);
3373385ff7eSMatthew G. Knepley   PetscFunctionReturn(0);
3383385ff7eSMatthew G. Knepley }
3393385ff7eSMatthew G. Knepley 
3403385ff7eSMatthew G. Knepley #undef __FUNCT__
341a66d4d66SMatthew G Knepley #define __FUNCT__ "DMDACreateSection"
342a66d4d66SMatthew G Knepley /*@C
343a66d4d66SMatthew G Knepley   DMDACreateSection - Create a PetscSection inside the DMDA that describes data layout. This allows multiple fields with
344a66d4d66SMatthew G Knepley   different numbers of dofs on vertices, cells, and faces in each direction.
345a66d4d66SMatthew G Knepley 
346a66d4d66SMatthew G Knepley   Input Parameters:
347a66d4d66SMatthew G Knepley + dm- The DMDA
348a66d4d66SMatthew G Knepley . numFields - The number of fields
349*affa55c7SMatthew G. Knepley . numComp - The number of components in each field
350*affa55c7SMatthew G. Knepley . numDof - The number of dofs per dimension for each field
3510298fd71SBarry Smith . numFaceDof - The number of dofs per face for each field and direction, or NULL
352a66d4d66SMatthew G Knepley 
353a66d4d66SMatthew G Knepley   Level: developer
354a66d4d66SMatthew G Knepley 
355a66d4d66SMatthew G Knepley   Note:
356a66d4d66SMatthew G Knepley   The default DMDA numbering is as follows:
357a66d4d66SMatthew G Knepley 
358a66d4d66SMatthew G Knepley     - Cells:    [0,             nC)
359a66d4d66SMatthew G Knepley     - Vertices: [nC,            nC+nV)
36088ed4aceSMatthew G Knepley     - X-Faces:  [nC+nV,         nC+nV+nXF)         normal is +- x-dir
36188ed4aceSMatthew G Knepley     - Y-Faces:  [nC+nV+nXF,     nC+nV+nXF+nYF)     normal is +- y-dir
36288ed4aceSMatthew G Knepley     - Z-Faces:  [nC+nV+nXF+nYF, nC+nV+nXF+nYF+nZF) normal is +- z-dir
363a66d4d66SMatthew G Knepley 
364a66d4d66SMatthew G Knepley   We interpret the default DMDA partition as a cell partition, and the data assignment as a cell assignment.
365a66d4d66SMatthew G Knepley @*/
366*affa55c7SMatthew G. Knepley PetscErrorCode DMDACreateSection(DM dm, PetscInt numComp[], PetscInt numDof[], PetscInt numFaceDof[], PetscSection *s)
367a66d4d66SMatthew G Knepley {
368a66d4d66SMatthew G Knepley   DM_DA            *da  = (DM_DA*) dm->data;
369a4b60ecfSMatthew G. Knepley   PetscSection      section;
37088ed4aceSMatthew G Knepley   const PetscInt    dim = da->dim;
37180800b1aSMatthew G Knepley   PetscInt          numFields, numVertexTotDof = 0, numCellTotDof = 0, numFaceTotDof[3] = {0, 0, 0};
37288ed4aceSMatthew G Knepley   PetscSF           sf;
373feafbc5cSMatthew G Knepley   PetscMPIInt       rank;
37488ed4aceSMatthew G Knepley   const PetscMPIInt *neighbors;
37588ed4aceSMatthew G Knepley   PetscInt          *localPoints;
37688ed4aceSMatthew G Knepley   PetscSFNode       *remotePoints;
377f5eeb527SMatthew G Knepley   PetscInt          nleaves = 0,  nleavesCheck = 0, nL = 0;
37857459e9aSMatthew G Knepley   PetscInt          nC, nVx, nVy, nVz, nV, nxF, nXF, nyF, nYF, nzF, nZF;
37957459e9aSMatthew G Knepley   PetscInt          pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart, yfEnd, zfStart, zfEnd;
38088ed4aceSMatthew G Knepley   PetscInt          f, v, c, xf, yf, zf, xn, yn, zn;
381a66d4d66SMatthew G Knepley   PetscErrorCode    ierr;
382a66d4d66SMatthew G Knepley 
383a66d4d66SMatthew G Knepley   PetscFunctionBegin;
384a66d4d66SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3853582350cSMatthew G. Knepley   PetscValidPointer(s, 4);
38682f516ccSBarry Smith   ierr    = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
3873582350cSMatthew G. Knepley   ierr    = DMDAGetNumCells(dm, NULL, NULL, NULL, &nC);CHKERRQ(ierr);
38857459e9aSMatthew G Knepley   ierr    = DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV);CHKERRQ(ierr);
38957459e9aSMatthew G Knepley   ierr    = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr);
39057459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, -1,  &pStart, &pEnd);CHKERRQ(ierr);
39157459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, 0,   &cStart, &cEnd);CHKERRQ(ierr);
39257459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, 1,   &fStart, &fEnd);CHKERRQ(ierr);
39357459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr);
39457459e9aSMatthew G Knepley   xfStart = vEnd;  xfEnd = xfStart+nXF;
39557459e9aSMatthew G Knepley   yfStart = xfEnd; yfEnd = yfStart+nYF;
39657459e9aSMatthew G Knepley   zfStart = yfEnd; zfEnd = zfStart+nZF;
39782f516ccSBarry Smith   if (zfEnd != fEnd) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Invalid face end %d, should be %d", zfEnd, fEnd);
39888ed4aceSMatthew G Knepley   /* Create local section */
39980800b1aSMatthew G Knepley   ierr = DMDAGetInfo(dm, 0,0,0,0,0,0,0, &numFields, 0,0,0,0,0);CHKERRQ(ierr);
400a66d4d66SMatthew G Knepley   for (f = 0; f < numFields; ++f) {
401*affa55c7SMatthew G. Knepley     numVertexTotDof  += numDof[f*(dim+1)+0];
402*affa55c7SMatthew G. Knepley     numCellTotDof    += numDof[f*(dim+1)+dim];
403*affa55c7SMatthew G. Knepley     numFaceTotDof[0] += dim > 0 ? (numFaceDof ? numFaceDof[f*dim+0] : numDof[f*(dim+1)+dim-1]) : 0;
404*affa55c7SMatthew G. Knepley     numFaceTotDof[1] += dim > 1 ? (numFaceDof ? numFaceDof[f*dim+1] : numDof[f*(dim+1)+dim-1]) : 0;
405*affa55c7SMatthew G. Knepley     numFaceTotDof[2] += dim > 2 ? (numFaceDof ? numFaceDof[f*dim+2] : numDof[f*(dim+1)+dim-1]) : 0;
406a66d4d66SMatthew G Knepley   }
407a4b60ecfSMatthew G. Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &section);CHKERRQ(ierr);
408*affa55c7SMatthew G. Knepley   if (numFields > 0) {
409a4b60ecfSMatthew G. Knepley     ierr = PetscSectionSetNumFields(section, numFields);CHKERRQ(ierr);
410a66d4d66SMatthew G Knepley     for (f = 0; f < numFields; ++f) {
41180800b1aSMatthew G Knepley       const char *name;
41280800b1aSMatthew G Knepley 
41380800b1aSMatthew G Knepley       ierr = DMDAGetFieldName(dm, f, &name);CHKERRQ(ierr);
414*affa55c7SMatthew G. Knepley       ierr = PetscSectionSetFieldName(section, f, name ? name : "Field");CHKERRQ(ierr);
41580800b1aSMatthew G Knepley       if (numComp) {
416a4b60ecfSMatthew G. Knepley         ierr = PetscSectionSetFieldComponents(section, f, numComp[f]);CHKERRQ(ierr);
417a66d4d66SMatthew G Knepley       }
418a66d4d66SMatthew G Knepley     }
419a66d4d66SMatthew G Knepley   }
420a4b60ecfSMatthew G. Knepley   ierr = PetscSectionSetChart(section, pStart, pEnd);CHKERRQ(ierr);
421a66d4d66SMatthew G Knepley   for (v = vStart; v < vEnd; ++v) {
422a66d4d66SMatthew G Knepley     for (f = 0; f < numFields; ++f) {
423*affa55c7SMatthew G. Knepley       ierr = PetscSectionSetFieldDof(section, v, f, numDof[f*(dim+1)+0]);CHKERRQ(ierr);
424a66d4d66SMatthew G Knepley     }
425a4b60ecfSMatthew G. Knepley     ierr = PetscSectionSetDof(section, v, numVertexTotDof);CHKERRQ(ierr);
426a66d4d66SMatthew G Knepley   }
427a66d4d66SMatthew G Knepley   for (xf = xfStart; xf < xfEnd; ++xf) {
428a66d4d66SMatthew G Knepley     for (f = 0; f < numFields; ++f) {
429*affa55c7SMatthew G. Knepley       ierr = PetscSectionSetFieldDof(section, xf, f, numFaceDof ? numFaceDof[f*dim+0] : numDof[f*(dim+1)+dim-1]);CHKERRQ(ierr);
430a66d4d66SMatthew G Knepley     }
431a4b60ecfSMatthew G. Knepley     ierr = PetscSectionSetDof(section, xf, numFaceTotDof[0]);CHKERRQ(ierr);
432a66d4d66SMatthew G Knepley   }
433a66d4d66SMatthew G Knepley   for (yf = yfStart; yf < yfEnd; ++yf) {
434a66d4d66SMatthew G Knepley     for (f = 0; f < numFields; ++f) {
435*affa55c7SMatthew G. Knepley       ierr = PetscSectionSetFieldDof(section, yf, f, numFaceDof ? numFaceDof[f*dim+1] : numDof[f*(dim+1)+dim-1]);CHKERRQ(ierr);
436a66d4d66SMatthew G Knepley     }
437a4b60ecfSMatthew G. Knepley     ierr = PetscSectionSetDof(section, yf, numFaceTotDof[1]);CHKERRQ(ierr);
438a66d4d66SMatthew G Knepley   }
439a66d4d66SMatthew G Knepley   for (zf = zfStart; zf < zfEnd; ++zf) {
440a66d4d66SMatthew G Knepley     for (f = 0; f < numFields; ++f) {
441*affa55c7SMatthew G. Knepley       ierr = PetscSectionSetFieldDof(section, zf, f, numFaceDof ? numFaceDof[f*dim+2] : numDof[f*(dim+1)+dim-1]);CHKERRQ(ierr);
442a66d4d66SMatthew G Knepley     }
443a4b60ecfSMatthew G. Knepley     ierr = PetscSectionSetDof(section, zf, numFaceTotDof[2]);CHKERRQ(ierr);
444a66d4d66SMatthew G Knepley   }
445a66d4d66SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
446a66d4d66SMatthew G Knepley     for (f = 0; f < numFields; ++f) {
447*affa55c7SMatthew G. Knepley       ierr = PetscSectionSetFieldDof(section, c, f, numDof[f*(dim+1)+dim]);CHKERRQ(ierr);
448a66d4d66SMatthew G Knepley     }
449a4b60ecfSMatthew G. Knepley     ierr = PetscSectionSetDof(section, c, numCellTotDof);CHKERRQ(ierr);
450a66d4d66SMatthew G Knepley   }
451a4b60ecfSMatthew G. Knepley   ierr = PetscSectionSetUp(section);CHKERRQ(ierr);
45288ed4aceSMatthew G Knepley   /* Create mesh point SF */
45388ed4aceSMatthew G Knepley   ierr = DMDAGetNeighbors(dm, &neighbors);CHKERRQ(ierr);
45488ed4aceSMatthew G Knepley   for (zn = 0; zn < (dim > 2 ? 3 : 1); ++zn) {
45588ed4aceSMatthew G Knepley     for (yn = 0; yn < (dim > 1 ? 3 : 1); ++yn) {
45688ed4aceSMatthew G Knepley       for (xn = 0; xn < 3; ++xn) {
4577128ae9fSMatthew G Knepley         const PetscInt xp       = xn-1, yp = dim > 1 ? yn-1 : 0, zp = dim > 2 ? zn-1 : 0;
45888ed4aceSMatthew G Knepley         const PetscInt neighbor = neighbors[(zn*3+yn)*3+xn];
45988ed4aceSMatthew G Knepley 
4603814d604SMatthew G Knepley         if (neighbor >= 0 && neighbor < rank) {
461feafbc5cSMatthew G Knepley           nleaves += (!xp ? nVx : 1) * (!yp ? nVy : 1) * (!zp ? nVz : 1); /* vertices */
46288ed4aceSMatthew G Knepley           if (xp && !yp && !zp) {
46388ed4aceSMatthew G Knepley             nleaves += nxF; /* x faces */
46488ed4aceSMatthew G Knepley           } else if (yp && !zp && !xp) {
46588ed4aceSMatthew G Knepley             nleaves += nyF; /* y faces */
46688ed4aceSMatthew G Knepley           } else if (zp && !xp && !yp) {
46788ed4aceSMatthew G Knepley             nleaves += nzF; /* z faces */
46888ed4aceSMatthew G Knepley           }
46988ed4aceSMatthew G Knepley         }
47088ed4aceSMatthew G Knepley       }
47188ed4aceSMatthew G Knepley     }
47288ed4aceSMatthew G Knepley   }
473dcca6d9dSJed Brown   ierr = PetscMalloc2(nleaves,&localPoints,nleaves,&remotePoints);CHKERRQ(ierr);
47488ed4aceSMatthew G Knepley   for (zn = 0; zn < (dim > 2 ? 3 : 1); ++zn) {
47588ed4aceSMatthew G Knepley     for (yn = 0; yn < (dim > 1 ? 3 : 1); ++yn) {
47688ed4aceSMatthew G Knepley       for (xn = 0; xn < 3; ++xn) {
4777128ae9fSMatthew G Knepley         const PetscInt xp       = xn-1, yp = dim > 1 ? yn-1 : 0, zp = dim > 2 ? zn-1 : 0;
47888ed4aceSMatthew G Knepley         const PetscInt neighbor = neighbors[(zn*3+yn)*3+xn];
479f5eeb527SMatthew G Knepley         PetscInt       xv, yv, zv;
48088ed4aceSMatthew G Knepley 
4813814d604SMatthew G Knepley         if (neighbor >= 0 && neighbor < rank) {
48288ed4aceSMatthew G Knepley           if (xp < 0) { /* left */
48388ed4aceSMatthew G Knepley             if (yp < 0) { /* bottom */
48488ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
485f5eeb527SMatthew G Knepley                 const PetscInt localVertex  = (      0*nVy +     0)*nVx +     0 + nC;
486f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = ((nVz-1)*nVy + nVy-1)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
487628e3b0dSSatish Balay                 nleavesCheck += 1; /* left bottom back vertex */
488f5eeb527SMatthew G Knepley 
489f5eeb527SMatthew G Knepley                 localPoints[nL]        = localVertex;
490f5eeb527SMatthew G Knepley                 remotePoints[nL].rank  = neighbor;
491f5eeb527SMatthew G Knepley                 remotePoints[nL].index = remoteVertex;
492f5eeb527SMatthew G Knepley                 ++nL;
49388ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
494f5eeb527SMatthew G Knepley                 const PetscInt localVertex  = ((nVz-1)*nVy +     0)*nVx +     0 + nC;
495f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = (      0*nVy + nVy-1)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
496628e3b0dSSatish Balay                 nleavesCheck += 1; /* left bottom front vertex */
497f5eeb527SMatthew G Knepley 
498f5eeb527SMatthew G Knepley                 localPoints[nL]        = localVertex;
499f5eeb527SMatthew G Knepley                 remotePoints[nL].rank  = neighbor;
500f5eeb527SMatthew G Knepley                 remotePoints[nL].index = remoteVertex;
501f5eeb527SMatthew G Knepley                 ++nL;
50288ed4aceSMatthew G Knepley               } else {
50388ed4aceSMatthew G Knepley                 nleavesCheck += nVz; /* left bottom vertices */
504f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv, ++nL) {
505f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = (zv*nVy +     0)*nVx +     0 + nC;
506f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (zv*nVy + nVy-1)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
507f5eeb527SMatthew G Knepley 
508f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
509f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
510f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
511f5eeb527SMatthew G Knepley                 }
51288ed4aceSMatthew G Knepley               }
51388ed4aceSMatthew G Knepley             } else if (yp > 0) { /* top */
51488ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
515f5eeb527SMatthew G Knepley                 const PetscInt localVertex  = (      0*nVy + nVy-1)*nVx +     0 + nC;
516f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = ((nVz-1)*nVy +     0)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
517628e3b0dSSatish Balay                 nleavesCheck += 1; /* left top back vertex */
518f5eeb527SMatthew G Knepley 
519f5eeb527SMatthew G Knepley                 localPoints[nL]        = localVertex;
520f5eeb527SMatthew G Knepley                 remotePoints[nL].rank  = neighbor;
521f5eeb527SMatthew G Knepley                 remotePoints[nL].index = remoteVertex;
522f5eeb527SMatthew G Knepley                 ++nL;
52388ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
524f5eeb527SMatthew G Knepley                 const PetscInt localVertex  = ((nVz-1)*nVy + nVy-1)*nVx +     0 + nC;
525f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = (      0*nVy +     0)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
526628e3b0dSSatish Balay                 nleavesCheck += 1; /* left top front vertex */
527f5eeb527SMatthew G Knepley 
528f5eeb527SMatthew G Knepley                 localPoints[nL]        = localVertex;
529f5eeb527SMatthew G Knepley                 remotePoints[nL].rank  = neighbor;
530f5eeb527SMatthew G Knepley                 remotePoints[nL].index = remoteVertex;
531f5eeb527SMatthew G Knepley                 ++nL;
53288ed4aceSMatthew G Knepley               } else {
53388ed4aceSMatthew G Knepley                 nleavesCheck += nVz; /* left top vertices */
534f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv, ++nL) {
535f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = (zv*nVy + nVy-1)*nVx +     0 + nC;
536f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (zv*nVy +     0)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
537f5eeb527SMatthew G Knepley 
538f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
539f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
540f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
541f5eeb527SMatthew G Knepley                 }
54288ed4aceSMatthew G Knepley               }
54388ed4aceSMatthew G Knepley             } else {
54488ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
54588ed4aceSMatthew G Knepley                 nleavesCheck += nVy; /* left back vertices */
546f5eeb527SMatthew G Knepley                 for (yv = 0; yv < nVy; ++yv, ++nL) {
547f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = (      0*nVy + yv)*nVx +     0 + nC;
548f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = ((nVz-1)*nVy + yv)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
549f5eeb527SMatthew G Knepley 
550f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
551f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
552f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
553f5eeb527SMatthew G Knepley                 }
55488ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
55588ed4aceSMatthew G Knepley                 nleavesCheck += nVy; /* left front vertices */
556f5eeb527SMatthew G Knepley                 for (yv = 0; yv < nVy; ++yv, ++nL) {
557f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = ((nVz-1)*nVy + yv)*nVx +     0 + nC;
558f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (      0*nVy + yv)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
559f5eeb527SMatthew G Knepley 
560f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
561f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
562f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
563f5eeb527SMatthew G Knepley                 }
56488ed4aceSMatthew G Knepley               } else {
56588ed4aceSMatthew G Knepley                 nleavesCheck += nVy*nVz; /* left vertices */
566f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv) {
567f5eeb527SMatthew G Knepley                   for (yv = 0; yv < nVy; ++yv, ++nL) {
568f5eeb527SMatthew G Knepley                     const PetscInt localVertex  = (zv*nVy + yv)*nVx +     0 + nC;
569f5eeb527SMatthew G Knepley                     const PetscInt remoteVertex = (zv*nVy + yv)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
570f5eeb527SMatthew G Knepley 
571f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
572f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
573f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
574f5eeb527SMatthew G Knepley                   }
575f5eeb527SMatthew G Knepley                 }
57688ed4aceSMatthew G Knepley                 nleavesCheck += nxF;     /* left faces */
577f5eeb527SMatthew G Knepley                 for (xf = 0; xf < nxF; ++xf, ++nL) {
578f5eeb527SMatthew G Knepley                   /* THIS IS WRONG */
579f5eeb527SMatthew G Knepley                   const PetscInt localFace  = 0 + nC+nV;
580f5eeb527SMatthew G Knepley                   const PetscInt remoteFace = 0 + nC+nV;
581f5eeb527SMatthew G Knepley 
582f5eeb527SMatthew G Knepley                   localPoints[nL]        = localFace;
583f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
584f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteFace;
585f5eeb527SMatthew G Knepley                 }
58688ed4aceSMatthew G Knepley               }
58788ed4aceSMatthew G Knepley             }
58888ed4aceSMatthew G Knepley           } else if (xp > 0) { /* right */
58988ed4aceSMatthew G Knepley             if (yp < 0) { /* bottom */
59088ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
591f5eeb527SMatthew G Knepley                 const PetscInt localVertex  = (      0*nVy +     0)*nVx + nVx-1 + nC;
592f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = ((nVz-1)*nVy + nVy-1)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
593628e3b0dSSatish Balay                 nleavesCheck += 1; /* right bottom back vertex */
594f5eeb527SMatthew G Knepley 
595f5eeb527SMatthew G Knepley                 localPoints[nL]        = localVertex;
596f5eeb527SMatthew G Knepley                 remotePoints[nL].rank  = neighbor;
597f5eeb527SMatthew G Knepley                 remotePoints[nL].index = remoteVertex;
598f5eeb527SMatthew G Knepley                 ++nL;
59988ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
600f5eeb527SMatthew G Knepley                 const PetscInt localVertex  = ((nVz-1)*nVy +     0)*nVx + nVx-1 + nC;
601f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = (      0*nVy + nVy-1)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
602628e3b0dSSatish Balay                 nleavesCheck += 1; /* right bottom front vertex */
603f5eeb527SMatthew G Knepley 
604f5eeb527SMatthew G Knepley                 localPoints[nL]        = localVertex;
605f5eeb527SMatthew G Knepley                 remotePoints[nL].rank  = neighbor;
606f5eeb527SMatthew G Knepley                 remotePoints[nL].index = remoteVertex;
607f5eeb527SMatthew G Knepley                 ++nL;
60888ed4aceSMatthew G Knepley               } else {
60988ed4aceSMatthew G Knepley                 nleavesCheck += nVz; /* right bottom vertices */
610f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv, ++nL) {
611f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = (zv*nVy +     0)*nVx + nVx-1 + nC;
612f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (zv*nVy + nVy-1)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
613f5eeb527SMatthew G Knepley 
614f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
615f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
616f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
617f5eeb527SMatthew G Knepley                 }
61888ed4aceSMatthew G Knepley               }
61988ed4aceSMatthew G Knepley             } else if (yp > 0) { /* top */
62088ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
621f5eeb527SMatthew G Knepley                 const PetscInt localVertex  = (      0*nVy + nVy-1)*nVx + nVx-1 + nC;
622f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = ((nVz-1)*nVy +     0)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
623628e3b0dSSatish Balay                 nleavesCheck += 1; /* right top back vertex */
624f5eeb527SMatthew G Knepley 
625f5eeb527SMatthew G Knepley                 localPoints[nL]        = localVertex;
626f5eeb527SMatthew G Knepley                 remotePoints[nL].rank  = neighbor;
627f5eeb527SMatthew G Knepley                 remotePoints[nL].index = remoteVertex;
628f5eeb527SMatthew G Knepley                 ++nL;
62988ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
630f5eeb527SMatthew G Knepley                 const PetscInt localVertex  = ((nVz-1)*nVy + nVy-1)*nVx + nVx-1 + nC;
631f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = (      0*nVy +     0)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
632628e3b0dSSatish Balay                 nleavesCheck += 1; /* right top front vertex */
633f5eeb527SMatthew G Knepley 
634f5eeb527SMatthew G Knepley                 localPoints[nL]        = localVertex;
635f5eeb527SMatthew G Knepley                 remotePoints[nL].rank  = neighbor;
636f5eeb527SMatthew G Knepley                 remotePoints[nL].index = remoteVertex;
637f5eeb527SMatthew G Knepley                 ++nL;
63888ed4aceSMatthew G Knepley               } else {
63988ed4aceSMatthew G Knepley                 nleavesCheck += nVz; /* right top vertices */
640f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv, ++nL) {
641f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = (zv*nVy + nVy-1)*nVx + nVx-1 + nC;
642f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (zv*nVy +     0)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
643f5eeb527SMatthew G Knepley 
644f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
645f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
646f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
647f5eeb527SMatthew G Knepley                 }
64888ed4aceSMatthew G Knepley               }
64988ed4aceSMatthew G Knepley             } else {
65088ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
65188ed4aceSMatthew G Knepley                 nleavesCheck += nVy; /* right back vertices */
652f5eeb527SMatthew G Knepley                 for (yv = 0; yv < nVy; ++yv, ++nL) {
653f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = (      0*nVy + yv)*nVx + nVx-1 + nC;
654f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = ((nVz-1)*nVy + yv)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
655f5eeb527SMatthew G Knepley 
656f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
657f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
658f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
659f5eeb527SMatthew G Knepley                 }
66088ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
66188ed4aceSMatthew G Knepley                 nleavesCheck += nVy; /* right front vertices */
662f5eeb527SMatthew G Knepley                 for (yv = 0; yv < nVy; ++yv, ++nL) {
663f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = ((nVz-1)*nVy + yv)*nVx + nVx-1 + nC;
664f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (      0*nVy + yv)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
665f5eeb527SMatthew G Knepley 
666f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
667f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
668f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
669f5eeb527SMatthew G Knepley                 }
67088ed4aceSMatthew G Knepley               } else {
67188ed4aceSMatthew G Knepley                 nleavesCheck += nVy*nVz; /* right vertices */
672f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv) {
673f5eeb527SMatthew G Knepley                   for (yv = 0; yv < nVy; ++yv, ++nL) {
674f5eeb527SMatthew G Knepley                     const PetscInt localVertex  = (zv*nVy + yv)*nVx + nVx-1 + nC;
675f5eeb527SMatthew G Knepley                     const PetscInt remoteVertex = (zv*nVy + yv)*nVx + 0     + nC; /* TODO: Correct this for neighbor sizes */
676f5eeb527SMatthew G Knepley 
677f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
678f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
679f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
680f5eeb527SMatthew G Knepley                   }
681f5eeb527SMatthew G Knepley                 }
68288ed4aceSMatthew G Knepley                 nleavesCheck += nxF;     /* right faces */
683f5eeb527SMatthew G Knepley                 for (xf = 0; xf < nxF; ++xf, ++nL) {
684f5eeb527SMatthew G Knepley                   /* THIS IS WRONG */
685f5eeb527SMatthew G Knepley                   const PetscInt localFace  = 0 + nC+nV;
686f5eeb527SMatthew G Knepley                   const PetscInt remoteFace = 0 + nC+nV;
687f5eeb527SMatthew G Knepley 
688f5eeb527SMatthew G Knepley                   localPoints[nL]        = localFace;
689f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
690f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteFace;
691f5eeb527SMatthew G Knepley                 }
69288ed4aceSMatthew G Knepley               }
69388ed4aceSMatthew G Knepley             }
69488ed4aceSMatthew G Knepley           } else {
69588ed4aceSMatthew G Knepley             if (yp < 0) { /* bottom */
69688ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
69788ed4aceSMatthew G Knepley                 nleavesCheck += nVx; /* bottom back vertices */
698f5eeb527SMatthew G Knepley                 for (xv = 0; xv < nVx; ++xv, ++nL) {
699f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = (      0*nVy +     0)*nVx + xv + nC;
700f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = ((nVz-1)*nVy + nVy-1)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
701f5eeb527SMatthew G Knepley 
702f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
703f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
704f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
705f5eeb527SMatthew G Knepley                 }
70688ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
70788ed4aceSMatthew G Knepley                 nleavesCheck += nVx; /* bottom front vertices */
708f5eeb527SMatthew G Knepley                 for (xv = 0; xv < nVx; ++xv, ++nL) {
709f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = ((nVz-1)*nVy +     0)*nVx + xv + nC;
710f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (      0*nVy + nVy-1)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
711f5eeb527SMatthew G Knepley 
712f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
713f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
714f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
715f5eeb527SMatthew G Knepley                 }
71688ed4aceSMatthew G Knepley               } else {
71788ed4aceSMatthew G Knepley                 nleavesCheck += nVx*nVz; /* bottom vertices */
718f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv) {
719f5eeb527SMatthew G Knepley                   for (xv = 0; xv < nVx; ++xv, ++nL) {
720f5eeb527SMatthew G Knepley                     const PetscInt localVertex  = (zv*nVy +     0)*nVx + xv + nC;
721f5eeb527SMatthew G Knepley                     const PetscInt remoteVertex = (zv*nVy + nVy-1)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
722f5eeb527SMatthew G Knepley 
723f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
724f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
725f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
726f5eeb527SMatthew G Knepley                   }
727f5eeb527SMatthew G Knepley                 }
72888ed4aceSMatthew G Knepley                 nleavesCheck += nyF;     /* bottom faces */
729f5eeb527SMatthew G Knepley                 for (yf = 0; yf < nyF; ++yf, ++nL) {
730f5eeb527SMatthew G Knepley                   /* THIS IS WRONG */
731f5eeb527SMatthew G Knepley                   const PetscInt localFace  = 0 + nC+nV;
732f5eeb527SMatthew G Knepley                   const PetscInt remoteFace = 0 + nC+nV;
733f5eeb527SMatthew G Knepley 
734f5eeb527SMatthew G Knepley                   localPoints[nL]        = localFace;
735f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
736f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteFace;
737f5eeb527SMatthew G Knepley                 }
73888ed4aceSMatthew G Knepley               }
73988ed4aceSMatthew G Knepley             } else if (yp > 0) { /* top */
74088ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
74188ed4aceSMatthew G Knepley                 nleavesCheck += nVx; /* top back vertices */
742f5eeb527SMatthew G Knepley                 for (xv = 0; xv < nVx; ++xv, ++nL) {
743f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = (      0*nVy + nVy-1)*nVx + xv + nC;
744f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = ((nVz-1)*nVy +     0)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
745f5eeb527SMatthew G Knepley 
746f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
747f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
748f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
749f5eeb527SMatthew G Knepley                 }
75088ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
75188ed4aceSMatthew G Knepley                 nleavesCheck += nVx; /* top front vertices */
752f5eeb527SMatthew G Knepley                 for (xv = 0; xv < nVx; ++xv, ++nL) {
753f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = ((nVz-1)*nVy + nVy-1)*nVx + xv + nC;
754f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (      0*nVy +     0)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
755f5eeb527SMatthew G Knepley 
756f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
757f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
758f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
759f5eeb527SMatthew G Knepley                 }
76088ed4aceSMatthew G Knepley               } else {
76188ed4aceSMatthew G Knepley                 nleavesCheck += nVx*nVz; /* top vertices */
762f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv) {
763f5eeb527SMatthew G Knepley                   for (xv = 0; xv < nVx; ++xv, ++nL) {
764f5eeb527SMatthew G Knepley                     const PetscInt localVertex  = (zv*nVy + nVy-1)*nVx + xv + nC;
765f5eeb527SMatthew G Knepley                     const PetscInt remoteVertex = (zv*nVy +     0)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
766f5eeb527SMatthew G Knepley 
767f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
768f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
769f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
770f5eeb527SMatthew G Knepley                   }
771f5eeb527SMatthew G Knepley                 }
77288ed4aceSMatthew G Knepley                 nleavesCheck += nyF;     /* top faces */
773f5eeb527SMatthew G Knepley                 for (yf = 0; yf < nyF; ++yf, ++nL) {
774f5eeb527SMatthew G Knepley                   /* THIS IS WRONG */
775f5eeb527SMatthew G Knepley                   const PetscInt localFace  = 0 + nC+nV;
776f5eeb527SMatthew G Knepley                   const PetscInt remoteFace = 0 + nC+nV;
777f5eeb527SMatthew G Knepley 
778f5eeb527SMatthew G Knepley                   localPoints[nL]        = localFace;
779f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
780f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteFace;
781f5eeb527SMatthew G Knepley                 }
78288ed4aceSMatthew G Knepley               }
78388ed4aceSMatthew G Knepley             } else {
78488ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
78588ed4aceSMatthew G Knepley                 nleavesCheck += nVx*nVy; /* back vertices */
786f5eeb527SMatthew G Knepley                 for (yv = 0; yv < nVy; ++yv) {
787f5eeb527SMatthew G Knepley                   for (xv = 0; xv < nVx; ++xv, ++nL) {
788f5eeb527SMatthew G Knepley                     const PetscInt localVertex  = (      0*nVy + yv)*nVx + xv + nC;
789f5eeb527SMatthew G Knepley                     const PetscInt remoteVertex = ((nVz-1)*nVy + yv)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
790f5eeb527SMatthew G Knepley 
791f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
792f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
793f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
794f5eeb527SMatthew G Knepley                   }
795f5eeb527SMatthew G Knepley                 }
79688ed4aceSMatthew G Knepley                 nleavesCheck += nzF;     /* back faces */
797f5eeb527SMatthew G Knepley                 for (zf = 0; zf < nzF; ++zf, ++nL) {
798f5eeb527SMatthew G Knepley                   /* THIS IS WRONG */
799f5eeb527SMatthew G Knepley                   const PetscInt localFace  = 0 + nC+nV;
800f5eeb527SMatthew G Knepley                   const PetscInt remoteFace = 0 + nC+nV;
801f5eeb527SMatthew G Knepley 
802f5eeb527SMatthew G Knepley                   localPoints[nL]        = localFace;
803f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
804f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteFace;
805f5eeb527SMatthew G Knepley                 }
80688ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
80788ed4aceSMatthew G Knepley                 nleavesCheck += nVx*nVy; /* front vertices */
808f5eeb527SMatthew G Knepley                 for (yv = 0; yv < nVy; ++yv) {
809f5eeb527SMatthew G Knepley                   for (xv = 0; xv < nVx; ++xv, ++nL) {
810f5eeb527SMatthew G Knepley                     const PetscInt localVertex  = ((nVz-1)*nVy + yv)*nVx + xv + nC;
811f5eeb527SMatthew G Knepley                     const PetscInt remoteVertex = (      0*nVy + yv)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
812f5eeb527SMatthew G Knepley 
813f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
814f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
815f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
816f5eeb527SMatthew G Knepley                   }
817f5eeb527SMatthew G Knepley                 }
81888ed4aceSMatthew G Knepley                 nleavesCheck += nzF;     /* front faces */
819f5eeb527SMatthew G Knepley                 for (zf = 0; zf < nzF; ++zf, ++nL) {
820f5eeb527SMatthew G Knepley                   /* THIS IS WRONG */
821f5eeb527SMatthew G Knepley                   const PetscInt localFace  = 0 + nC+nV;
822f5eeb527SMatthew G Knepley                   const PetscInt remoteFace = 0 + nC+nV;
823f5eeb527SMatthew G Knepley 
824f5eeb527SMatthew G Knepley                   localPoints[nL]        = localFace;
825f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
826f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteFace;
827f5eeb527SMatthew G Knepley                 }
82888ed4aceSMatthew G Knepley               } else {
82988ed4aceSMatthew G Knepley                 /* Nothing is shared from the interior */
83088ed4aceSMatthew G Knepley               }
83188ed4aceSMatthew G Knepley             }
83288ed4aceSMatthew G Knepley           }
83388ed4aceSMatthew G Knepley         }
83488ed4aceSMatthew G Knepley       }
83588ed4aceSMatthew G Knepley     }
83688ed4aceSMatthew G Knepley   }
8373814d604SMatthew G Knepley   /* TODO: Remove duplication in leaf determination */
83882f516ccSBarry Smith   if (nleaves != nleavesCheck) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "The number of leaves %d did not match the number of remote leaves %d", nleaves, nleavesCheck);
83982f516ccSBarry Smith   ierr = PetscSFCreate(PetscObjectComm((PetscObject)dm), &sf);CHKERRQ(ierr);
84088ed4aceSMatthew G Knepley   ierr = PetscSFSetGraph(sf, pEnd, nleaves, localPoints, PETSC_OWN_POINTER, remotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr);
841a4b60ecfSMatthew G. Knepley   ierr = DMSetPointSF(dm, sf);CHKERRQ(ierr);
84288ed4aceSMatthew G Knepley   ierr = PetscSFDestroy(&sf);CHKERRQ(ierr);
843*affa55c7SMatthew G. Knepley   *s = section;
844*affa55c7SMatthew G. Knepley   PetscFunctionReturn(0);
845*affa55c7SMatthew G. Knepley }
846*affa55c7SMatthew G. Knepley 
8473385ff7eSMatthew G. Knepley #undef __FUNCT__
8483385ff7eSMatthew G. Knepley #define __FUNCT__ "DMDASetVertexCoordinates"
8493385ff7eSMatthew G. Knepley PetscErrorCode DMDASetVertexCoordinates(DM dm, PetscReal xl, PetscReal xu, PetscReal yl, PetscReal yu, PetscReal zl, PetscReal zu)
8503385ff7eSMatthew G. Knepley {
8513385ff7eSMatthew G. Knepley   DM_DA         *da = (DM_DA *) dm->data;
8523385ff7eSMatthew G. Knepley   Vec            coordinates;
8533385ff7eSMatthew G. Knepley   PetscSection   section;
8543385ff7eSMatthew G. Knepley   PetscScalar   *coords;
8553385ff7eSMatthew G. Knepley   PetscReal      h[3];
8563385ff7eSMatthew G. Knepley   PetscInt       dim, size, M, N, P, nVx, nVy, nVz, nV, vStart, vEnd, v, i, j, k;
8573385ff7eSMatthew G. Knepley   PetscErrorCode ierr;
8583385ff7eSMatthew G. Knepley 
8593385ff7eSMatthew G. Knepley   PetscFunctionBegin;
8603385ff7eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8613385ff7eSMatthew G. Knepley   ierr = DMDAGetInfo(dm, &dim, &M, &N, &P, 0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
8623385ff7eSMatthew G. Knepley   h[0] = (xu - xl)/M;
8633385ff7eSMatthew G. Knepley   h[1] = (yu - yl)/N;
8643385ff7eSMatthew G. Knepley   h[2] = (zu - zl)/P;
8653385ff7eSMatthew G. Knepley   ierr = DMDAGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
8663385ff7eSMatthew G. Knepley   ierr = DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV);CHKERRQ(ierr);
8673385ff7eSMatthew G. Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &section);CHKERRQ(ierr);
8683385ff7eSMatthew G. Knepley   ierr = PetscSectionSetNumFields(section, 1);CHKERRQ(ierr);
8693385ff7eSMatthew G. Knepley   ierr = PetscSectionSetFieldComponents(section, 0, dim);CHKERRQ(ierr);
8703385ff7eSMatthew G. Knepley   ierr = PetscSectionSetChart(section, vStart, vEnd);CHKERRQ(ierr);
8713385ff7eSMatthew G. Knepley   for (v = vStart; v < vEnd; ++v) {
8723385ff7eSMatthew G. Knepley     ierr = PetscSectionSetDof(section, v, dim);CHKERRQ(ierr);
8733385ff7eSMatthew G. Knepley   }
8743385ff7eSMatthew G. Knepley   ierr = PetscSectionSetUp(section);CHKERRQ(ierr);
8753385ff7eSMatthew G. Knepley   ierr = PetscSectionGetStorageSize(section, &size);CHKERRQ(ierr);
8763385ff7eSMatthew G. Knepley   ierr = VecCreateSeq(PETSC_COMM_SELF, size, &coordinates);CHKERRQ(ierr);
8773385ff7eSMatthew G. Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
8783385ff7eSMatthew G. Knepley   for (k = 0; k < nVz; ++k) {
8793385ff7eSMatthew G. Knepley     PetscInt ind[3] = {0, 0, k + da->zs}, d, off;
8803385ff7eSMatthew G. Knepley 
8813385ff7eSMatthew G. Knepley     for (j = 0; j < nVy; ++j) {
8823385ff7eSMatthew G. Knepley       ind[1] = j + da->ys;
8833385ff7eSMatthew G. Knepley       for (i = 0; i < nVx; ++i) {
8843385ff7eSMatthew G. Knepley         const PetscInt vertex = (k*nVy + j)*nVx + i + vStart;
8853385ff7eSMatthew G. Knepley 
8863385ff7eSMatthew G. Knepley         ierr = PetscSectionGetOffset(section, vertex, &off);CHKERRQ(ierr);
8873385ff7eSMatthew G. Knepley         ind[0] = i + da->xs;
8883385ff7eSMatthew G. Knepley         for (d = 0; d < dim; ++d) {
8893385ff7eSMatthew G. Knepley           coords[off+d] = h[d]*ind[d];
8903385ff7eSMatthew G. Knepley         }
8913385ff7eSMatthew G. Knepley       }
8923385ff7eSMatthew G. Knepley     }
8933385ff7eSMatthew G. Knepley   }
8943385ff7eSMatthew G. Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
8953385ff7eSMatthew G. Knepley   ierr = DMSetCoordinateSection(dm, section);CHKERRQ(ierr);
8963385ff7eSMatthew G. Knepley   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
897a4b60ecfSMatthew G. Knepley   ierr = PetscSectionDestroy(&section);CHKERRQ(ierr);
8983385ff7eSMatthew G. Knepley   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
8993385ff7eSMatthew G. Knepley   PetscFunctionReturn(0);
9003385ff7eSMatthew G. Knepley }
901a66d4d66SMatthew G Knepley   PetscFunctionReturn(0);
902a66d4d66SMatthew G Knepley }
903a66d4d66SMatthew G Knepley 
90447c6ae99SBarry Smith /* ------------------------------------------------------------------- */
90547c6ae99SBarry Smith 
90647c6ae99SBarry Smith #undef __FUNCT__
907aa219208SBarry Smith #define __FUNCT__ "DMDAGetArray"
90847c6ae99SBarry Smith /*@C
909aa219208SBarry Smith      DMDAGetArray - Gets a work array for a DMDA
91047c6ae99SBarry Smith 
91147c6ae99SBarry Smith     Input Parameter:
91247c6ae99SBarry Smith +    da - information about my local patch
91347c6ae99SBarry Smith -    ghosted - do you want arrays for the ghosted or nonghosted patch
91447c6ae99SBarry Smith 
91547c6ae99SBarry Smith     Output Parameters:
91647c6ae99SBarry Smith .    vptr - array data structured
91747c6ae99SBarry Smith 
91847c6ae99SBarry Smith     Note:  The vector values are NOT initialized and may have garbage in them, so you may need
91947c6ae99SBarry Smith            to zero them.
92047c6ae99SBarry Smith 
92147c6ae99SBarry Smith   Level: advanced
92247c6ae99SBarry Smith 
923bcaeba4dSBarry Smith .seealso: DMDARestoreArray()
92447c6ae99SBarry Smith 
92547c6ae99SBarry Smith @*/
9267087cfbeSBarry Smith PetscErrorCode  DMDAGetArray(DM da,PetscBool ghosted,void *vptr)
92747c6ae99SBarry Smith {
92847c6ae99SBarry Smith   PetscErrorCode ierr;
92947c6ae99SBarry Smith   PetscInt       j,i,xs,ys,xm,ym,zs,zm;
93047c6ae99SBarry Smith   char           *iarray_start;
93147c6ae99SBarry Smith   void           **iptr = (void**)vptr;
93247c6ae99SBarry Smith   DM_DA          *dd    = (DM_DA*)da->data;
93347c6ae99SBarry Smith 
93447c6ae99SBarry Smith   PetscFunctionBegin;
93547c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
93647c6ae99SBarry Smith   if (ghosted) {
937aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
93847c6ae99SBarry Smith       if (dd->arrayghostedin[i]) {
93947c6ae99SBarry Smith         *iptr                 = dd->arrayghostedin[i];
94047c6ae99SBarry Smith         iarray_start          = (char*)dd->startghostedin[i];
9410298fd71SBarry Smith         dd->arrayghostedin[i] = NULL;
9420298fd71SBarry Smith         dd->startghostedin[i] = NULL;
94347c6ae99SBarry Smith 
94447c6ae99SBarry Smith         goto done;
94547c6ae99SBarry Smith       }
94647c6ae99SBarry Smith     }
94747c6ae99SBarry Smith     xs = dd->Xs;
94847c6ae99SBarry Smith     ys = dd->Ys;
94947c6ae99SBarry Smith     zs = dd->Zs;
95047c6ae99SBarry Smith     xm = dd->Xe-dd->Xs;
95147c6ae99SBarry Smith     ym = dd->Ye-dd->Ys;
95247c6ae99SBarry Smith     zm = dd->Ze-dd->Zs;
95347c6ae99SBarry Smith   } else {
954aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
95547c6ae99SBarry Smith       if (dd->arrayin[i]) {
95647c6ae99SBarry Smith         *iptr          = dd->arrayin[i];
95747c6ae99SBarry Smith         iarray_start   = (char*)dd->startin[i];
9580298fd71SBarry Smith         dd->arrayin[i] = NULL;
9590298fd71SBarry Smith         dd->startin[i] = NULL;
96047c6ae99SBarry Smith 
96147c6ae99SBarry Smith         goto done;
96247c6ae99SBarry Smith       }
96347c6ae99SBarry Smith     }
96447c6ae99SBarry Smith     xs = dd->xs;
96547c6ae99SBarry Smith     ys = dd->ys;
96647c6ae99SBarry Smith     zs = dd->zs;
96747c6ae99SBarry Smith     xm = dd->xe-dd->xs;
96847c6ae99SBarry Smith     ym = dd->ye-dd->ys;
96947c6ae99SBarry Smith     zm = dd->ze-dd->zs;
97047c6ae99SBarry Smith   }
97147c6ae99SBarry Smith 
97247c6ae99SBarry Smith   switch (dd->dim) {
97347c6ae99SBarry Smith   case 1: {
97447c6ae99SBarry Smith     void *ptr;
97547c6ae99SBarry Smith 
976901f1932SJed Brown     ierr = PetscMalloc(xm*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
97747c6ae99SBarry Smith 
97847c6ae99SBarry Smith     ptr   = (void*)(iarray_start - xs*sizeof(PetscScalar));
97947c6ae99SBarry Smith     *iptr = (void*)ptr;
9808865f1eaSKarl Rupp     break;
9818865f1eaSKarl Rupp   }
98247c6ae99SBarry Smith   case 2: {
98347c6ae99SBarry Smith     void **ptr;
98447c6ae99SBarry Smith 
985901f1932SJed Brown     ierr = PetscMalloc((ym+1)*sizeof(void*)+xm*ym*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
98647c6ae99SBarry Smith 
98747c6ae99SBarry Smith     ptr = (void**)(iarray_start + xm*ym*sizeof(PetscScalar) - ys*sizeof(void*));
9888865f1eaSKarl Rupp     for (j=ys; j<ys+ym; j++) ptr[j] = iarray_start + sizeof(PetscScalar)*(xm*(j-ys) - xs);
98947c6ae99SBarry Smith     *iptr = (void*)ptr;
9908865f1eaSKarl Rupp     break;
9918865f1eaSKarl Rupp   }
99247c6ae99SBarry Smith   case 3: {
99347c6ae99SBarry Smith     void ***ptr,**bptr;
99447c6ae99SBarry Smith 
995901f1932SJed Brown     ierr = PetscMalloc((zm+1)*sizeof(void**)+(ym*zm+1)*sizeof(void*)+xm*ym*zm*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
99647c6ae99SBarry Smith 
99747c6ae99SBarry Smith     ptr  = (void***)(iarray_start + xm*ym*zm*sizeof(PetscScalar) - zs*sizeof(void*));
99847c6ae99SBarry Smith     bptr = (void**)(iarray_start + xm*ym*zm*sizeof(PetscScalar) + zm*sizeof(void**));
9998865f1eaSKarl Rupp     for (i=zs; i<zs+zm; i++) ptr[i] = bptr + ((i-zs)*ym - ys);
100047c6ae99SBarry Smith     for (i=zs; i<zs+zm; i++) {
100147c6ae99SBarry Smith       for (j=ys; j<ys+ym; j++) {
100247c6ae99SBarry Smith         ptr[i][j] = iarray_start + sizeof(PetscScalar)*(xm*ym*(i-zs) + xm*(j-ys) - xs);
100347c6ae99SBarry Smith       }
100447c6ae99SBarry Smith     }
100547c6ae99SBarry Smith     *iptr = (void*)ptr;
10068865f1eaSKarl Rupp     break;
10078865f1eaSKarl Rupp   }
100847c6ae99SBarry Smith   default:
1009ce94432eSBarry Smith     SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Dimension %D not supported",dd->dim);
101047c6ae99SBarry Smith   }
101147c6ae99SBarry Smith 
101247c6ae99SBarry Smith done:
101347c6ae99SBarry Smith   /* add arrays to the checked out list */
101447c6ae99SBarry Smith   if (ghosted) {
1015aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
101647c6ae99SBarry Smith       if (!dd->arrayghostedout[i]) {
101747c6ae99SBarry Smith         dd->arrayghostedout[i] = *iptr;
101847c6ae99SBarry Smith         dd->startghostedout[i] = iarray_start;
101947c6ae99SBarry Smith         break;
102047c6ae99SBarry Smith       }
102147c6ae99SBarry Smith     }
102247c6ae99SBarry Smith   } else {
1023aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
102447c6ae99SBarry Smith       if (!dd->arrayout[i]) {
102547c6ae99SBarry Smith         dd->arrayout[i] = *iptr;
102647c6ae99SBarry Smith         dd->startout[i] = iarray_start;
102747c6ae99SBarry Smith         break;
102847c6ae99SBarry Smith       }
102947c6ae99SBarry Smith     }
103047c6ae99SBarry Smith   }
103147c6ae99SBarry Smith   PetscFunctionReturn(0);
103247c6ae99SBarry Smith }
103347c6ae99SBarry Smith 
103447c6ae99SBarry Smith #undef __FUNCT__
1035aa219208SBarry Smith #define __FUNCT__ "DMDARestoreArray"
103647c6ae99SBarry Smith /*@C
1037aa219208SBarry Smith      DMDARestoreArray - Restores an array of derivative types for a DMDA
103847c6ae99SBarry Smith 
103947c6ae99SBarry Smith     Input Parameter:
104047c6ae99SBarry Smith +    da - information about my local patch
104147c6ae99SBarry Smith .    ghosted - do you want arrays for the ghosted or nonghosted patch
104247c6ae99SBarry Smith -    vptr - array data structured to be passed to ad_FormFunctionLocal()
104347c6ae99SBarry Smith 
104447c6ae99SBarry Smith      Level: advanced
104547c6ae99SBarry Smith 
1046bcaeba4dSBarry Smith .seealso: DMDAGetArray()
104747c6ae99SBarry Smith 
104847c6ae99SBarry Smith @*/
10497087cfbeSBarry Smith PetscErrorCode  DMDARestoreArray(DM da,PetscBool ghosted,void *vptr)
105047c6ae99SBarry Smith {
105147c6ae99SBarry Smith   PetscInt i;
105247c6ae99SBarry Smith   void     **iptr = (void**)vptr,*iarray_start = 0;
105347c6ae99SBarry Smith   DM_DA    *dd    = (DM_DA*)da->data;
105447c6ae99SBarry Smith 
105547c6ae99SBarry Smith   PetscFunctionBegin;
105647c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
105747c6ae99SBarry Smith   if (ghosted) {
1058aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
105947c6ae99SBarry Smith       if (dd->arrayghostedout[i] == *iptr) {
106047c6ae99SBarry Smith         iarray_start           = dd->startghostedout[i];
10610298fd71SBarry Smith         dd->arrayghostedout[i] = NULL;
10620298fd71SBarry Smith         dd->startghostedout[i] = NULL;
106347c6ae99SBarry Smith         break;
106447c6ae99SBarry Smith       }
106547c6ae99SBarry Smith     }
1066aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
106747c6ae99SBarry Smith       if (!dd->arrayghostedin[i]) {
106847c6ae99SBarry Smith         dd->arrayghostedin[i] = *iptr;
106947c6ae99SBarry Smith         dd->startghostedin[i] = iarray_start;
107047c6ae99SBarry Smith         break;
107147c6ae99SBarry Smith       }
107247c6ae99SBarry Smith     }
107347c6ae99SBarry Smith   } else {
1074aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
107547c6ae99SBarry Smith       if (dd->arrayout[i] == *iptr) {
107647c6ae99SBarry Smith         iarray_start    = dd->startout[i];
10770298fd71SBarry Smith         dd->arrayout[i] = NULL;
10780298fd71SBarry Smith         dd->startout[i] = NULL;
107947c6ae99SBarry Smith         break;
108047c6ae99SBarry Smith       }
108147c6ae99SBarry Smith     }
1082aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
108347c6ae99SBarry Smith       if (!dd->arrayin[i]) {
108447c6ae99SBarry Smith         dd->arrayin[i] = *iptr;
108547c6ae99SBarry Smith         dd->startin[i] = iarray_start;
108647c6ae99SBarry Smith         break;
108747c6ae99SBarry Smith       }
108847c6ae99SBarry Smith     }
108947c6ae99SBarry Smith   }
109047c6ae99SBarry Smith   PetscFunctionReturn(0);
109147c6ae99SBarry Smith }
109247c6ae99SBarry Smith 
1093