xref: /petsc/src/dm/impls/da/dalocal.c (revision e4d69e081bcf8bbabfb8a3c65677035c2cd07561)
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>
89a800dd8SMatthew G. Knepley #include <petscfe.h>
947c6ae99SBarry Smith 
1047c6ae99SBarry Smith /*
11e3c5b3baSBarry Smith    This allows the DMDA vectors to properly tell MATLAB their dimensions
1247c6ae99SBarry Smith */
1347c6ae99SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
14c6db04a5SJed Brown #include <engine.h>   /* MATLAB include file */
15c6db04a5SJed Brown #include <mex.h>      /* MATLAB include file */
1647c6ae99SBarry Smith #undef __FUNCT__
1747c6ae99SBarry Smith #define __FUNCT__ "VecMatlabEnginePut_DA2d"
181e6b0712SBarry Smith static PetscErrorCode  VecMatlabEnginePut_DA2d(PetscObject obj,void *mengine)
1947c6ae99SBarry Smith {
2047c6ae99SBarry Smith   PetscErrorCode ierr;
2147c6ae99SBarry Smith   PetscInt       n,m;
2247c6ae99SBarry Smith   Vec            vec = (Vec)obj;
2347c6ae99SBarry Smith   PetscScalar    *array;
2447c6ae99SBarry Smith   mxArray        *mat;
259a42bb27SBarry Smith   DM             da;
2647c6ae99SBarry Smith 
2747c6ae99SBarry Smith   PetscFunctionBegin;
28c688c046SMatthew G Knepley   ierr = VecGetDM(vec, &da);CHKERRQ(ierr);
29ce94432eSBarry Smith   if (!da) SETERRQ(PetscObjectComm((PetscObject)vec),PETSC_ERR_ARG_WRONGSTATE,"Vector not associated with a DMDA");
30aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,0,0,0,&m,&n,0);CHKERRQ(ierr);
3147c6ae99SBarry Smith 
3247c6ae99SBarry Smith   ierr = VecGetArray(vec,&array);CHKERRQ(ierr);
3347c6ae99SBarry Smith #if !defined(PETSC_USE_COMPLEX)
3447c6ae99SBarry Smith   mat = mxCreateDoubleMatrix(m,n,mxREAL);
3547c6ae99SBarry Smith #else
3647c6ae99SBarry Smith   mat = mxCreateDoubleMatrix(m,n,mxCOMPLEX);
3747c6ae99SBarry Smith #endif
3847c6ae99SBarry Smith   ierr = PetscMemcpy(mxGetPr(mat),array,n*m*sizeof(PetscScalar));CHKERRQ(ierr);
3947c6ae99SBarry Smith   ierr = PetscObjectName(obj);CHKERRQ(ierr);
4047c6ae99SBarry Smith   engPutVariable((Engine*)mengine,obj->name,mat);
4147c6ae99SBarry Smith 
4247c6ae99SBarry Smith   ierr = VecRestoreArray(vec,&array);CHKERRQ(ierr);
4347c6ae99SBarry Smith   PetscFunctionReturn(0);
4447c6ae99SBarry Smith }
4547c6ae99SBarry Smith #endif
4647c6ae99SBarry Smith 
4747c6ae99SBarry Smith 
4847c6ae99SBarry Smith #undef __FUNCT__
49564755cdSBarry Smith #define __FUNCT__ "DMCreateLocalVector_DA"
507087cfbeSBarry Smith PetscErrorCode  DMCreateLocalVector_DA(DM da,Vec *g)
5147c6ae99SBarry Smith {
5247c6ae99SBarry Smith   PetscErrorCode ierr;
5347c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
5447c6ae99SBarry Smith 
5547c6ae99SBarry Smith   PetscFunctionBegin;
5647c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
5747c6ae99SBarry Smith   PetscValidPointer(g,2);
5811689aebSJed Brown   if (da->defaultSection) {
5911689aebSJed Brown     ierr = DMCreateLocalVector_Section_Private(da,g);CHKERRQ(ierr);
6011689aebSJed Brown   } else {
6147c6ae99SBarry Smith     ierr = VecCreate(PETSC_COMM_SELF,g);CHKERRQ(ierr);
6247c6ae99SBarry Smith     ierr = VecSetSizes(*g,dd->nlocal,PETSC_DETERMINE);CHKERRQ(ierr);
6347c6ae99SBarry Smith     ierr = VecSetBlockSize(*g,dd->w);CHKERRQ(ierr);
64401ddaa8SBarry Smith     ierr = VecSetType(*g,da->vectype);CHKERRQ(ierr);
65c688c046SMatthew G Knepley     ierr = VecSetDM(*g, da);CHKERRQ(ierr);
6647c6ae99SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
6747c6ae99SBarry Smith     if (dd->w == 1  && dd->dim == 2) {
68bdf89e91SBarry Smith       ierr = PetscObjectComposeFunction((PetscObject)*g,"PetscMatlabEnginePut_C",VecMatlabEnginePut_DA2d);CHKERRQ(ierr);
6947c6ae99SBarry Smith     }
7047c6ae99SBarry Smith #endif
7111689aebSJed Brown   }
7247c6ae99SBarry Smith   PetscFunctionReturn(0);
7347c6ae99SBarry Smith }
7447c6ae99SBarry Smith 
75a66d4d66SMatthew G Knepley #undef __FUNCT__
7657459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetNumCells"
773582350cSMatthew G. Knepley PetscErrorCode DMDAGetNumCells(DM dm, PetscInt *numCellsX, PetscInt *numCellsY, PetscInt *numCellsZ, PetscInt *numCells)
7857459e9aSMatthew G Knepley {
7957459e9aSMatthew G Knepley   DM_DA          *da = (DM_DA*) dm->data;
8057459e9aSMatthew G Knepley   const PetscInt dim = da->dim;
8157459e9aSMatthew G Knepley   const PetscInt mx  = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs;
8257459e9aSMatthew G Knepley   const PetscInt nC  = (mx)*(dim > 1 ? (my)*(dim > 2 ? (mz) : 1) : 1);
8357459e9aSMatthew G Knepley 
8457459e9aSMatthew G Knepley   PetscFunctionBegin;
853582350cSMatthew G. Knepley   if (numCellsX) {
863582350cSMatthew G. Knepley     PetscValidIntPointer(numCellsX,2);
873582350cSMatthew G. Knepley     *numCellsX = mx;
883582350cSMatthew G. Knepley   }
893582350cSMatthew G. Knepley   if (numCellsY) {
903582350cSMatthew G. Knepley     PetscValidIntPointer(numCellsX,3);
913582350cSMatthew G. Knepley     *numCellsY = my;
923582350cSMatthew G. Knepley   }
933582350cSMatthew G. Knepley   if (numCellsZ) {
943582350cSMatthew G. Knepley     PetscValidIntPointer(numCellsX,4);
953582350cSMatthew G. Knepley     *numCellsZ = mz;
963582350cSMatthew G. Knepley   }
9757459e9aSMatthew G Knepley   if (numCells) {
983582350cSMatthew G. Knepley     PetscValidIntPointer(numCells,5);
9957459e9aSMatthew G Knepley     *numCells = nC;
10057459e9aSMatthew G Knepley   }
10157459e9aSMatthew G Knepley   PetscFunctionReturn(0);
10257459e9aSMatthew G Knepley }
10357459e9aSMatthew G Knepley 
10457459e9aSMatthew G Knepley #undef __FUNCT__
10557459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetNumVertices"
10657459e9aSMatthew G Knepley PetscErrorCode DMDAGetNumVertices(DM dm, PetscInt *numVerticesX, PetscInt *numVerticesY, PetscInt *numVerticesZ, PetscInt *numVertices)
10757459e9aSMatthew G Knepley {
10857459e9aSMatthew G Knepley   DM_DA          *da = (DM_DA*) dm->data;
10957459e9aSMatthew G Knepley   const PetscInt dim = da->dim;
11057459e9aSMatthew G Knepley   const PetscInt mx  = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs;
11157459e9aSMatthew G Knepley   const PetscInt nVx = mx+1;
11257459e9aSMatthew G Knepley   const PetscInt nVy = dim > 1 ? (my+1) : 1;
11357459e9aSMatthew G Knepley   const PetscInt nVz = dim > 2 ? (mz+1) : 1;
11457459e9aSMatthew G Knepley   const PetscInt nV  = nVx*nVy*nVz;
11557459e9aSMatthew G Knepley 
11657459e9aSMatthew G Knepley   PetscFunctionBegin;
11757459e9aSMatthew G Knepley   if (numVerticesX) {
11857459e9aSMatthew G Knepley     PetscValidIntPointer(numVerticesX,2);
11957459e9aSMatthew G Knepley     *numVerticesX = nVx;
12057459e9aSMatthew G Knepley   }
12157459e9aSMatthew G Knepley   if (numVerticesY) {
12257459e9aSMatthew G Knepley     PetscValidIntPointer(numVerticesY,3);
12357459e9aSMatthew G Knepley     *numVerticesY = nVy;
12457459e9aSMatthew G Knepley   }
12557459e9aSMatthew G Knepley   if (numVerticesZ) {
12657459e9aSMatthew G Knepley     PetscValidIntPointer(numVerticesZ,4);
12757459e9aSMatthew G Knepley     *numVerticesZ = nVz;
12857459e9aSMatthew G Knepley   }
12957459e9aSMatthew G Knepley   if (numVertices) {
13057459e9aSMatthew G Knepley     PetscValidIntPointer(numVertices,5);
13157459e9aSMatthew G Knepley     *numVertices = nV;
13257459e9aSMatthew G Knepley   }
13357459e9aSMatthew G Knepley   PetscFunctionReturn(0);
13457459e9aSMatthew G Knepley }
13557459e9aSMatthew G Knepley 
13657459e9aSMatthew G Knepley #undef __FUNCT__
13757459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetNumFaces"
13857459e9aSMatthew G Knepley PetscErrorCode DMDAGetNumFaces(DM dm, PetscInt *numXFacesX, PetscInt *numXFaces, PetscInt *numYFacesY, PetscInt *numYFaces, PetscInt *numZFacesZ, PetscInt *numZFaces)
13957459e9aSMatthew G Knepley {
14057459e9aSMatthew G Knepley   DM_DA          *da = (DM_DA*) dm->data;
14157459e9aSMatthew G Knepley   const PetscInt dim = da->dim;
14257459e9aSMatthew G Knepley   const PetscInt mx  = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs;
14357459e9aSMatthew G Knepley   const PetscInt nxF = (dim > 1 ? (my)*(dim > 2 ? (mz) : 1) : 1);
14457459e9aSMatthew G Knepley   const PetscInt nXF = (mx+1)*nxF;
14557459e9aSMatthew G Knepley   const PetscInt nyF = mx*(dim > 2 ? mz : 1);
14657459e9aSMatthew G Knepley   const PetscInt nYF = dim > 1 ? (my+1)*nyF : 0;
14757459e9aSMatthew G Knepley   const PetscInt nzF = mx*(dim > 1 ? my : 0);
14857459e9aSMatthew G Knepley   const PetscInt nZF = dim > 2 ? (mz+1)*nzF : 0;
14957459e9aSMatthew G Knepley 
15057459e9aSMatthew G Knepley   PetscFunctionBegin;
15157459e9aSMatthew G Knepley   if (numXFacesX) {
15257459e9aSMatthew G Knepley     PetscValidIntPointer(numXFacesX,2);
15357459e9aSMatthew G Knepley     *numXFacesX = nxF;
15457459e9aSMatthew G Knepley   }
15557459e9aSMatthew G Knepley   if (numXFaces) {
15657459e9aSMatthew G Knepley     PetscValidIntPointer(numXFaces,3);
15757459e9aSMatthew G Knepley     *numXFaces = nXF;
15857459e9aSMatthew G Knepley   }
15957459e9aSMatthew G Knepley   if (numYFacesY) {
16057459e9aSMatthew G Knepley     PetscValidIntPointer(numYFacesY,4);
16157459e9aSMatthew G Knepley     *numYFacesY = nyF;
16257459e9aSMatthew G Knepley   }
16357459e9aSMatthew G Knepley   if (numYFaces) {
16457459e9aSMatthew G Knepley     PetscValidIntPointer(numYFaces,5);
16557459e9aSMatthew G Knepley     *numYFaces = nYF;
16657459e9aSMatthew G Knepley   }
16757459e9aSMatthew G Knepley   if (numZFacesZ) {
16857459e9aSMatthew G Knepley     PetscValidIntPointer(numZFacesZ,6);
16957459e9aSMatthew G Knepley     *numZFacesZ = nzF;
17057459e9aSMatthew G Knepley   }
17157459e9aSMatthew G Knepley   if (numZFaces) {
17257459e9aSMatthew G Knepley     PetscValidIntPointer(numZFaces,7);
17357459e9aSMatthew G Knepley     *numZFaces = nZF;
17457459e9aSMatthew G Knepley   }
17557459e9aSMatthew G Knepley   PetscFunctionReturn(0);
17657459e9aSMatthew G Knepley }
17757459e9aSMatthew G Knepley 
17857459e9aSMatthew G Knepley #undef __FUNCT__
17957459e9aSMatthew G Knepley #define __FUNCT__ "DMDAGetHeightStratum"
18057459e9aSMatthew G Knepley PetscErrorCode DMDAGetHeightStratum(DM dm, PetscInt height, PetscInt *pStart, PetscInt *pEnd)
18157459e9aSMatthew G Knepley {
18257459e9aSMatthew G Knepley   DM_DA          *da = (DM_DA*) dm->data;
18357459e9aSMatthew G Knepley   const PetscInt dim = da->dim;
18457459e9aSMatthew G Knepley   PetscInt       nC, nV, nXF, nYF, nZF;
18557459e9aSMatthew G Knepley   PetscErrorCode ierr;
18657459e9aSMatthew G Knepley 
18757459e9aSMatthew G Knepley   PetscFunctionBegin;
1888865f1eaSKarl Rupp   if (pStart) PetscValidIntPointer(pStart,3);
1898865f1eaSKarl Rupp   if (pEnd)   PetscValidIntPointer(pEnd,4);
1903582350cSMatthew G. Knepley   ierr = DMDAGetNumCells(dm, NULL, NULL, NULL, &nC);CHKERRQ(ierr);
1910298fd71SBarry Smith   ierr = DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV);CHKERRQ(ierr);
1920298fd71SBarry Smith   ierr = DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF);CHKERRQ(ierr);
19357459e9aSMatthew G Knepley   if (height == 0) {
19457459e9aSMatthew G Knepley     /* Cells */
1958865f1eaSKarl Rupp     if (pStart) *pStart = 0;
1968865f1eaSKarl Rupp     if (pEnd)   *pEnd   = nC;
19757459e9aSMatthew G Knepley   } else if (height == 1) {
19857459e9aSMatthew G Knepley     /* Faces */
1998865f1eaSKarl Rupp     if (pStart) *pStart = nC+nV;
2008865f1eaSKarl Rupp     if (pEnd)   *pEnd   = nC+nV+nXF+nYF+nZF;
20157459e9aSMatthew G Knepley   } else if (height == dim) {
20257459e9aSMatthew G Knepley     /* Vertices */
2038865f1eaSKarl Rupp     if (pStart) *pStart = nC;
2048865f1eaSKarl Rupp     if (pEnd)   *pEnd   = nC+nV;
20557459e9aSMatthew G Knepley   } else if (height < 0) {
20657459e9aSMatthew G Knepley     /* All points */
2078865f1eaSKarl Rupp     if (pStart) *pStart = 0;
2088865f1eaSKarl Rupp     if (pEnd)   *pEnd   = nC+nV+nXF+nYF+nZF;
20982f516ccSBarry Smith   } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No points of height %d in the DA", height);
21057459e9aSMatthew G Knepley   PetscFunctionReturn(0);
21157459e9aSMatthew G Knepley }
21257459e9aSMatthew G Knepley 
21357459e9aSMatthew G Knepley #undef __FUNCT__
2143385ff7eSMatthew G. Knepley #define __FUNCT__ "DMDAGetDepthStratum"
2153385ff7eSMatthew G. Knepley PetscErrorCode DMDAGetDepthStratum(DM dm, PetscInt depth, PetscInt *pStart, PetscInt *pEnd)
2163385ff7eSMatthew G. Knepley {
2173385ff7eSMatthew G. Knepley   DM_DA         *da  = (DM_DA*) dm->data;
2183385ff7eSMatthew G. Knepley   const PetscInt dim = da->dim;
2193385ff7eSMatthew G. Knepley   PetscInt       nC, nV, nXF, nYF, nZF;
2203385ff7eSMatthew G. Knepley   PetscErrorCode ierr;
2213385ff7eSMatthew G. Knepley 
2223385ff7eSMatthew G. Knepley   PetscFunctionBegin;
2233385ff7eSMatthew G. Knepley   if (pStart) PetscValidIntPointer(pStart,3);
2243385ff7eSMatthew G. Knepley   if (pEnd)   PetscValidIntPointer(pEnd,4);
2253385ff7eSMatthew G. Knepley   ierr = DMDAGetNumCells(dm, NULL, NULL, NULL, &nC);CHKERRQ(ierr);
2263385ff7eSMatthew G. Knepley   ierr = DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV);CHKERRQ(ierr);
2273385ff7eSMatthew G. Knepley   ierr = DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF);CHKERRQ(ierr);
2283385ff7eSMatthew G. Knepley   if (depth == dim) {
2293385ff7eSMatthew G. Knepley     /* Cells */
2303385ff7eSMatthew G. Knepley     if (pStart) *pStart = 0;
2313385ff7eSMatthew G. Knepley     if (pEnd)   *pEnd   = nC;
2323385ff7eSMatthew G. Knepley   } else if (depth == dim-1) {
2333385ff7eSMatthew G. Knepley     /* Faces */
2343385ff7eSMatthew G. Knepley     if (pStart) *pStart = nC+nV;
2353385ff7eSMatthew G. Knepley     if (pEnd)   *pEnd   = nC+nV+nXF+nYF+nZF;
2363385ff7eSMatthew G. Knepley   } else if (depth == 0) {
2373385ff7eSMatthew G. Knepley     /* Vertices */
2383385ff7eSMatthew G. Knepley     if (pStart) *pStart = nC;
2393385ff7eSMatthew G. Knepley     if (pEnd)   *pEnd   = nC+nV;
2403385ff7eSMatthew G. Knepley   } else if (depth < 0) {
2413385ff7eSMatthew G. Knepley     /* All points */
2423385ff7eSMatthew G. Knepley     if (pStart) *pStart = 0;
2433385ff7eSMatthew G. Knepley     if (pEnd)   *pEnd   = nC+nV+nXF+nYF+nZF;
2443385ff7eSMatthew G. Knepley   } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No points of depth %d in the DA", depth);
2453385ff7eSMatthew G. Knepley   PetscFunctionReturn(0);
2463385ff7eSMatthew G. Knepley }
2473385ff7eSMatthew G. Knepley 
2483385ff7eSMatthew G. Knepley #undef __FUNCT__
2493385ff7eSMatthew G. Knepley #define __FUNCT__ "DMDAGetConeSize"
2503385ff7eSMatthew G. Knepley PetscErrorCode DMDAGetConeSize(DM dm, PetscInt p, PetscInt *coneSize)
2513385ff7eSMatthew G. Knepley {
2523385ff7eSMatthew G. Knepley   DM_DA         *da  = (DM_DA*) dm->data;
2533385ff7eSMatthew G. Knepley   const PetscInt dim = da->dim;
2543385ff7eSMatthew G. Knepley   PetscInt       nC, nV, nXF, nYF, nZF;
2553385ff7eSMatthew G. Knepley   PetscErrorCode ierr;
2563385ff7eSMatthew G. Knepley 
2573385ff7eSMatthew G. Knepley   PetscFunctionBegin;
2583385ff7eSMatthew G. Knepley   *coneSize = 0;
2593385ff7eSMatthew G. Knepley   ierr = DMDAGetNumCells(dm, NULL, NULL, NULL, &nC);CHKERRQ(ierr);
2603385ff7eSMatthew G. Knepley   ierr = DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV);CHKERRQ(ierr);
2613385ff7eSMatthew G. Knepley   ierr = DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF);CHKERRQ(ierr);
2623385ff7eSMatthew G. Knepley   switch (dim) {
2633385ff7eSMatthew G. Knepley   case 2:
2643385ff7eSMatthew G. Knepley     if (p >= 0) {
2653385ff7eSMatthew G. Knepley       if (p < nC) {
2663385ff7eSMatthew G. Knepley         *coneSize = 4;
2673385ff7eSMatthew G. Knepley       } else if (p < nC+nV) {
2683385ff7eSMatthew G. Knepley         *coneSize = 0;
2693385ff7eSMatthew G. Knepley       } else if (p < nC+nV+nXF+nYF+nZF) {
2703385ff7eSMatthew G. Knepley         *coneSize = 2;
2713385ff7eSMatthew G. Knepley       } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %d should be in [0, %d)", p, nC+nV+nXF+nYF+nZF);
2723385ff7eSMatthew G. Knepley     } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative point %d is invalid", p);
2733385ff7eSMatthew G. Knepley     break;
2743385ff7eSMatthew G. Knepley   case 3:
2753385ff7eSMatthew G. Knepley     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do 3D");
2763385ff7eSMatthew G. Knepley     break;
2773385ff7eSMatthew G. Knepley   }
2783385ff7eSMatthew G. Knepley   PetscFunctionReturn(0);
2793385ff7eSMatthew G. Knepley }
2803385ff7eSMatthew G. Knepley 
2813385ff7eSMatthew G. Knepley #undef __FUNCT__
2823385ff7eSMatthew G. Knepley #define __FUNCT__ "DMDAGetCone"
2833385ff7eSMatthew G. Knepley PetscErrorCode DMDAGetCone(DM dm, PetscInt p, PetscInt *cone[])
2843385ff7eSMatthew G. Knepley {
2853385ff7eSMatthew G. Knepley   DM_DA         *da  = (DM_DA*) dm->data;
2863385ff7eSMatthew G. Knepley   const PetscInt dim = da->dim;
2873385ff7eSMatthew G. Knepley   PetscInt       nCx, nCy, nCz, nC, nVx, nVy, nVz, nV, nxF, nyF, nzF, nXF, nYF, nZF;
2883385ff7eSMatthew G. Knepley   PetscErrorCode ierr;
2893385ff7eSMatthew G. Knepley 
2903385ff7eSMatthew G. Knepley   PetscFunctionBegin;
2913385ff7eSMatthew G. Knepley   if (!cone) {ierr = DMGetWorkArray(dm, 6, PETSC_INT, cone);CHKERRQ(ierr);}
2923385ff7eSMatthew G. Knepley   ierr = DMDAGetNumCells(dm, &nCx, &nCy, &nCz, &nC);CHKERRQ(ierr);
2933385ff7eSMatthew G. Knepley   ierr = DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV);CHKERRQ(ierr);
2943385ff7eSMatthew G. Knepley   ierr = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr);
2953385ff7eSMatthew G. Knepley   switch (dim) {
2963385ff7eSMatthew G. Knepley   case 2:
2973385ff7eSMatthew G. Knepley     if (p >= 0) {
2983385ff7eSMatthew G. Knepley       if (p < nC) {
2993385ff7eSMatthew G. Knepley         const PetscInt cy = p / nCx;
3003385ff7eSMatthew G. Knepley         const PetscInt cx = p % nCx;
3013385ff7eSMatthew G. Knepley 
3023385ff7eSMatthew G. Knepley         (*cone)[0] = cy*nxF + cx + nC+nV;
3033385ff7eSMatthew G. Knepley         (*cone)[1] = cx*nyF + cy + nyF + nC+nV+nXF;
3043385ff7eSMatthew G. Knepley         (*cone)[2] = cy*nxF + cx + nxF + nC+nV;
3053385ff7eSMatthew G. Knepley         (*cone)[3] = cx*nyF + cy + nC+nV+nXF;
3063385ff7eSMatthew G. Knepley         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do cell cones");
3073385ff7eSMatthew G. Knepley       } else if (p < nC+nV) {
3083385ff7eSMatthew G. Knepley       } else if (p < nC+nV+nXF) {
3093385ff7eSMatthew G. Knepley         const PetscInt fy = (p - nC+nV) / nxF;
3103385ff7eSMatthew G. Knepley         const PetscInt fx = (p - nC+nV) % nxF;
3113385ff7eSMatthew G. Knepley 
3123385ff7eSMatthew G. Knepley         (*cone)[0] = fy*nVx + fx + nC;
3133385ff7eSMatthew G. Knepley         (*cone)[1] = fy*nVx + fx + 1 + nC;
3143385ff7eSMatthew G. Knepley       } else if (p < nC+nV+nXF+nYF) {
3153385ff7eSMatthew G. Knepley         const PetscInt fx = (p - nC+nV+nXF) / nyF;
3163385ff7eSMatthew G. Knepley         const PetscInt fy = (p - nC+nV+nXF) % nyF;
3173385ff7eSMatthew G. Knepley 
3183385ff7eSMatthew G. Knepley         (*cone)[0] = fy*nVx + fx + nC;
3193385ff7eSMatthew G. Knepley         (*cone)[1] = fy*nVx + fx + nVx + nC;
3203385ff7eSMatthew G. Knepley       } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %d should be in [0, %d)", p, nC+nV+nXF+nYF+nZF);
3213385ff7eSMatthew G. Knepley     } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative point %d is invalid", p);
3223385ff7eSMatthew G. Knepley     break;
3233385ff7eSMatthew G. Knepley   case 3:
3243385ff7eSMatthew G. Knepley     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do 3D");
3253385ff7eSMatthew G. Knepley     break;
3263385ff7eSMatthew G. Knepley   }
3273385ff7eSMatthew G. Knepley   PetscFunctionReturn(0);
3283385ff7eSMatthew G. Knepley }
3293385ff7eSMatthew G. Knepley 
3303385ff7eSMatthew G. Knepley #undef __FUNCT__
3313385ff7eSMatthew G. Knepley #define __FUNCT__ "DMDARestoreCone"
3323385ff7eSMatthew G. Knepley PetscErrorCode DMDARestoreCone(DM dm, PetscInt p, PetscInt *cone[])
3333385ff7eSMatthew G. Knepley {
3343385ff7eSMatthew G. Knepley   PetscErrorCode ierr;
3353385ff7eSMatthew G. Knepley 
3363385ff7eSMatthew G. Knepley   PetscFunctionBegin;
3373385ff7eSMatthew G. Knepley   ierr = DMGetWorkArray(dm, 6, PETSC_INT, cone);CHKERRQ(ierr);
3383385ff7eSMatthew G. Knepley   PetscFunctionReturn(0);
3393385ff7eSMatthew G. Knepley }
3403385ff7eSMatthew G. Knepley 
3413385ff7eSMatthew G. Knepley #undef __FUNCT__
342a66d4d66SMatthew G Knepley #define __FUNCT__ "DMDACreateSection"
343a66d4d66SMatthew G Knepley /*@C
344a66d4d66SMatthew G Knepley   DMDACreateSection - Create a PetscSection inside the DMDA that describes data layout. This allows multiple fields with
345a66d4d66SMatthew G Knepley   different numbers of dofs on vertices, cells, and faces in each direction.
346a66d4d66SMatthew G Knepley 
347a66d4d66SMatthew G Knepley   Input Parameters:
348a66d4d66SMatthew G Knepley + dm- The DMDA
349a66d4d66SMatthew G Knepley . numFields - The number of fields
350affa55c7SMatthew G. Knepley . numComp - The number of components in each field
351affa55c7SMatthew G. Knepley . numDof - The number of dofs per dimension for each field
3520298fd71SBarry Smith . numFaceDof - The number of dofs per face for each field and direction, or NULL
353a66d4d66SMatthew G Knepley 
354a66d4d66SMatthew G Knepley   Level: developer
355a66d4d66SMatthew G Knepley 
356a66d4d66SMatthew G Knepley   Note:
357a66d4d66SMatthew G Knepley   The default DMDA numbering is as follows:
358a66d4d66SMatthew G Knepley 
359a66d4d66SMatthew G Knepley     - Cells:    [0,             nC)
360a66d4d66SMatthew G Knepley     - Vertices: [nC,            nC+nV)
36188ed4aceSMatthew G Knepley     - X-Faces:  [nC+nV,         nC+nV+nXF)         normal is +- x-dir
36288ed4aceSMatthew G Knepley     - Y-Faces:  [nC+nV+nXF,     nC+nV+nXF+nYF)     normal is +- y-dir
36388ed4aceSMatthew G Knepley     - Z-Faces:  [nC+nV+nXF+nYF, nC+nV+nXF+nYF+nZF) normal is +- z-dir
364a66d4d66SMatthew G Knepley 
365a66d4d66SMatthew G Knepley   We interpret the default DMDA partition as a cell partition, and the data assignment as a cell assignment.
366a66d4d66SMatthew G Knepley @*/
367affa55c7SMatthew G. Knepley PetscErrorCode DMDACreateSection(DM dm, PetscInt numComp[], PetscInt numDof[], PetscInt numFaceDof[], PetscSection *s)
368a66d4d66SMatthew G Knepley {
369a66d4d66SMatthew G Knepley   DM_DA            *da  = (DM_DA*) dm->data;
370a4b60ecfSMatthew G. Knepley   PetscSection      section;
37188ed4aceSMatthew G Knepley   const PetscInt    dim = da->dim;
37280800b1aSMatthew G Knepley   PetscInt          numFields, numVertexTotDof = 0, numCellTotDof = 0, numFaceTotDof[3] = {0, 0, 0};
37388ed4aceSMatthew G Knepley   PetscSF           sf;
374feafbc5cSMatthew G Knepley   PetscMPIInt       rank;
37588ed4aceSMatthew G Knepley   const PetscMPIInt *neighbors;
37688ed4aceSMatthew G Knepley   PetscInt          *localPoints;
37788ed4aceSMatthew G Knepley   PetscSFNode       *remotePoints;
378f5eeb527SMatthew G Knepley   PetscInt          nleaves = 0,  nleavesCheck = 0, nL = 0;
37957459e9aSMatthew G Knepley   PetscInt          nC, nVx, nVy, nVz, nV, nxF, nXF, nyF, nYF, nzF, nZF;
38057459e9aSMatthew G Knepley   PetscInt          pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart, yfEnd, zfStart, zfEnd;
38188ed4aceSMatthew G Knepley   PetscInt          f, v, c, xf, yf, zf, xn, yn, zn;
382a66d4d66SMatthew G Knepley   PetscErrorCode    ierr;
383a66d4d66SMatthew G Knepley 
384a66d4d66SMatthew G Knepley   PetscFunctionBegin;
385a66d4d66SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3863582350cSMatthew G. Knepley   PetscValidPointer(s, 4);
38782f516ccSBarry Smith   ierr    = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr);
3883582350cSMatthew G. Knepley   ierr    = DMDAGetNumCells(dm, NULL, NULL, NULL, &nC);CHKERRQ(ierr);
38957459e9aSMatthew G Knepley   ierr    = DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV);CHKERRQ(ierr);
39057459e9aSMatthew G Knepley   ierr    = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr);
39157459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, -1,  &pStart, &pEnd);CHKERRQ(ierr);
39257459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, 0,   &cStart, &cEnd);CHKERRQ(ierr);
39357459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, 1,   &fStart, &fEnd);CHKERRQ(ierr);
39457459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr);
39557459e9aSMatthew G Knepley   xfStart = vEnd;  xfEnd = xfStart+nXF;
39657459e9aSMatthew G Knepley   yfStart = xfEnd; yfEnd = yfStart+nYF;
39757459e9aSMatthew G Knepley   zfStart = yfEnd; zfEnd = zfStart+nZF;
39882f516ccSBarry Smith   if (zfEnd != fEnd) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Invalid face end %d, should be %d", zfEnd, fEnd);
39988ed4aceSMatthew G Knepley   /* Create local section */
40080800b1aSMatthew G Knepley   ierr = DMDAGetInfo(dm, 0,0,0,0,0,0,0, &numFields, 0,0,0,0,0);CHKERRQ(ierr);
401a66d4d66SMatthew G Knepley   for (f = 0; f < numFields; ++f) {
402affa55c7SMatthew G. Knepley     numVertexTotDof  += numDof[f*(dim+1)+0];
403affa55c7SMatthew G. Knepley     numCellTotDof    += numDof[f*(dim+1)+dim];
404affa55c7SMatthew G. Knepley     numFaceTotDof[0] += dim > 0 ? (numFaceDof ? numFaceDof[f*dim+0] : numDof[f*(dim+1)+dim-1]) : 0;
405affa55c7SMatthew G. Knepley     numFaceTotDof[1] += dim > 1 ? (numFaceDof ? numFaceDof[f*dim+1] : numDof[f*(dim+1)+dim-1]) : 0;
406affa55c7SMatthew G. Knepley     numFaceTotDof[2] += dim > 2 ? (numFaceDof ? numFaceDof[f*dim+2] : numDof[f*(dim+1)+dim-1]) : 0;
407a66d4d66SMatthew G Knepley   }
408a4b60ecfSMatthew G. Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), &section);CHKERRQ(ierr);
409affa55c7SMatthew G. Knepley   if (numFields > 0) {
410a4b60ecfSMatthew G. Knepley     ierr = PetscSectionSetNumFields(section, numFields);CHKERRQ(ierr);
411a66d4d66SMatthew G Knepley     for (f = 0; f < numFields; ++f) {
41280800b1aSMatthew G Knepley       const char *name;
41380800b1aSMatthew G Knepley 
41480800b1aSMatthew G Knepley       ierr = DMDAGetFieldName(dm, f, &name);CHKERRQ(ierr);
415affa55c7SMatthew G. Knepley       ierr = PetscSectionSetFieldName(section, f, name ? name : "Field");CHKERRQ(ierr);
41680800b1aSMatthew G Knepley       if (numComp) {
417a4b60ecfSMatthew G. Knepley         ierr = PetscSectionSetFieldComponents(section, f, numComp[f]);CHKERRQ(ierr);
418a66d4d66SMatthew G Knepley       }
419a66d4d66SMatthew G Knepley     }
420a66d4d66SMatthew G Knepley   }
421a4b60ecfSMatthew G. Knepley   ierr = PetscSectionSetChart(section, pStart, pEnd);CHKERRQ(ierr);
422a66d4d66SMatthew G Knepley   for (v = vStart; v < vEnd; ++v) {
423a66d4d66SMatthew G Knepley     for (f = 0; f < numFields; ++f) {
424affa55c7SMatthew G. Knepley       ierr = PetscSectionSetFieldDof(section, v, f, numDof[f*(dim+1)+0]);CHKERRQ(ierr);
425a66d4d66SMatthew G Knepley     }
426a4b60ecfSMatthew G. Knepley     ierr = PetscSectionSetDof(section, v, numVertexTotDof);CHKERRQ(ierr);
427a66d4d66SMatthew G Knepley   }
428a66d4d66SMatthew G Knepley   for (xf = xfStart; xf < xfEnd; ++xf) {
429a66d4d66SMatthew G Knepley     for (f = 0; f < numFields; ++f) {
430affa55c7SMatthew G. Knepley       ierr = PetscSectionSetFieldDof(section, xf, f, numFaceDof ? numFaceDof[f*dim+0] : numDof[f*(dim+1)+dim-1]);CHKERRQ(ierr);
431a66d4d66SMatthew G Knepley     }
432a4b60ecfSMatthew G. Knepley     ierr = PetscSectionSetDof(section, xf, numFaceTotDof[0]);CHKERRQ(ierr);
433a66d4d66SMatthew G Knepley   }
434a66d4d66SMatthew G Knepley   for (yf = yfStart; yf < yfEnd; ++yf) {
435a66d4d66SMatthew G Knepley     for (f = 0; f < numFields; ++f) {
436affa55c7SMatthew G. Knepley       ierr = PetscSectionSetFieldDof(section, yf, f, numFaceDof ? numFaceDof[f*dim+1] : numDof[f*(dim+1)+dim-1]);CHKERRQ(ierr);
437a66d4d66SMatthew G Knepley     }
438a4b60ecfSMatthew G. Knepley     ierr = PetscSectionSetDof(section, yf, numFaceTotDof[1]);CHKERRQ(ierr);
439a66d4d66SMatthew G Knepley   }
440a66d4d66SMatthew G Knepley   for (zf = zfStart; zf < zfEnd; ++zf) {
441a66d4d66SMatthew G Knepley     for (f = 0; f < numFields; ++f) {
442affa55c7SMatthew G. Knepley       ierr = PetscSectionSetFieldDof(section, zf, f, numFaceDof ? numFaceDof[f*dim+2] : numDof[f*(dim+1)+dim-1]);CHKERRQ(ierr);
443a66d4d66SMatthew G Knepley     }
444a4b60ecfSMatthew G. Knepley     ierr = PetscSectionSetDof(section, zf, numFaceTotDof[2]);CHKERRQ(ierr);
445a66d4d66SMatthew G Knepley   }
446a66d4d66SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
447a66d4d66SMatthew G Knepley     for (f = 0; f < numFields; ++f) {
448affa55c7SMatthew G. Knepley       ierr = PetscSectionSetFieldDof(section, c, f, numDof[f*(dim+1)+dim]);CHKERRQ(ierr);
449a66d4d66SMatthew G Knepley     }
450a4b60ecfSMatthew G. Knepley     ierr = PetscSectionSetDof(section, c, numCellTotDof);CHKERRQ(ierr);
451a66d4d66SMatthew G Knepley   }
452a4b60ecfSMatthew G. Knepley   ierr = PetscSectionSetUp(section);CHKERRQ(ierr);
45388ed4aceSMatthew G Knepley   /* Create mesh point SF */
45488ed4aceSMatthew G Knepley   ierr = DMDAGetNeighbors(dm, &neighbors);CHKERRQ(ierr);
45588ed4aceSMatthew G Knepley   for (zn = 0; zn < (dim > 2 ? 3 : 1); ++zn) {
45688ed4aceSMatthew G Knepley     for (yn = 0; yn < (dim > 1 ? 3 : 1); ++yn) {
45788ed4aceSMatthew G Knepley       for (xn = 0; xn < 3; ++xn) {
4587128ae9fSMatthew G Knepley         const PetscInt xp       = xn-1, yp = dim > 1 ? yn-1 : 0, zp = dim > 2 ? zn-1 : 0;
45988ed4aceSMatthew G Knepley         const PetscInt neighbor = neighbors[(zn*3+yn)*3+xn];
46088ed4aceSMatthew G Knepley 
4613814d604SMatthew G Knepley         if (neighbor >= 0 && neighbor < rank) {
462feafbc5cSMatthew G Knepley           nleaves += (!xp ? nVx : 1) * (!yp ? nVy : 1) * (!zp ? nVz : 1); /* vertices */
46388ed4aceSMatthew G Knepley           if (xp && !yp && !zp) {
46488ed4aceSMatthew G Knepley             nleaves += nxF; /* x faces */
46588ed4aceSMatthew G Knepley           } else if (yp && !zp && !xp) {
46688ed4aceSMatthew G Knepley             nleaves += nyF; /* y faces */
46788ed4aceSMatthew G Knepley           } else if (zp && !xp && !yp) {
46888ed4aceSMatthew G Knepley             nleaves += nzF; /* z faces */
46988ed4aceSMatthew G Knepley           }
47088ed4aceSMatthew G Knepley         }
47188ed4aceSMatthew G Knepley       }
47288ed4aceSMatthew G Knepley     }
47388ed4aceSMatthew G Knepley   }
474dcca6d9dSJed Brown   ierr = PetscMalloc2(nleaves,&localPoints,nleaves,&remotePoints);CHKERRQ(ierr);
47588ed4aceSMatthew G Knepley   for (zn = 0; zn < (dim > 2 ? 3 : 1); ++zn) {
47688ed4aceSMatthew G Knepley     for (yn = 0; yn < (dim > 1 ? 3 : 1); ++yn) {
47788ed4aceSMatthew G Knepley       for (xn = 0; xn < 3; ++xn) {
4787128ae9fSMatthew G Knepley         const PetscInt xp       = xn-1, yp = dim > 1 ? yn-1 : 0, zp = dim > 2 ? zn-1 : 0;
47988ed4aceSMatthew G Knepley         const PetscInt neighbor = neighbors[(zn*3+yn)*3+xn];
480f5eeb527SMatthew G Knepley         PetscInt       xv, yv, zv;
48188ed4aceSMatthew G Knepley 
4823814d604SMatthew G Knepley         if (neighbor >= 0 && neighbor < rank) {
48388ed4aceSMatthew G Knepley           if (xp < 0) { /* left */
48488ed4aceSMatthew G Knepley             if (yp < 0) { /* bottom */
48588ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
486f5eeb527SMatthew G Knepley                 const PetscInt localVertex  = (      0*nVy +     0)*nVx +     0 + nC;
487f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = ((nVz-1)*nVy + nVy-1)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
488628e3b0dSSatish Balay                 nleavesCheck += 1; /* left bottom back vertex */
489f5eeb527SMatthew G Knepley 
490f5eeb527SMatthew G Knepley                 localPoints[nL]        = localVertex;
491f5eeb527SMatthew G Knepley                 remotePoints[nL].rank  = neighbor;
492f5eeb527SMatthew G Knepley                 remotePoints[nL].index = remoteVertex;
493f5eeb527SMatthew G Knepley                 ++nL;
49488ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
495f5eeb527SMatthew G Knepley                 const PetscInt localVertex  = ((nVz-1)*nVy +     0)*nVx +     0 + nC;
496f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = (      0*nVy + nVy-1)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
497628e3b0dSSatish Balay                 nleavesCheck += 1; /* left bottom front vertex */
498f5eeb527SMatthew G Knepley 
499f5eeb527SMatthew G Knepley                 localPoints[nL]        = localVertex;
500f5eeb527SMatthew G Knepley                 remotePoints[nL].rank  = neighbor;
501f5eeb527SMatthew G Knepley                 remotePoints[nL].index = remoteVertex;
502f5eeb527SMatthew G Knepley                 ++nL;
50388ed4aceSMatthew G Knepley               } else {
50488ed4aceSMatthew G Knepley                 nleavesCheck += nVz; /* left bottom vertices */
505f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv, ++nL) {
506f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = (zv*nVy +     0)*nVx +     0 + nC;
507f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (zv*nVy + nVy-1)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
508f5eeb527SMatthew G Knepley 
509f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
510f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
511f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
512f5eeb527SMatthew G Knepley                 }
51388ed4aceSMatthew G Knepley               }
51488ed4aceSMatthew G Knepley             } else if (yp > 0) { /* top */
51588ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
516f5eeb527SMatthew G Knepley                 const PetscInt localVertex  = (      0*nVy + nVy-1)*nVx +     0 + nC;
517f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = ((nVz-1)*nVy +     0)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
518628e3b0dSSatish Balay                 nleavesCheck += 1; /* left top back vertex */
519f5eeb527SMatthew G Knepley 
520f5eeb527SMatthew G Knepley                 localPoints[nL]        = localVertex;
521f5eeb527SMatthew G Knepley                 remotePoints[nL].rank  = neighbor;
522f5eeb527SMatthew G Knepley                 remotePoints[nL].index = remoteVertex;
523f5eeb527SMatthew G Knepley                 ++nL;
52488ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
525f5eeb527SMatthew G Knepley                 const PetscInt localVertex  = ((nVz-1)*nVy + nVy-1)*nVx +     0 + nC;
526f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = (      0*nVy +     0)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
527628e3b0dSSatish Balay                 nleavesCheck += 1; /* left top front vertex */
528f5eeb527SMatthew G Knepley 
529f5eeb527SMatthew G Knepley                 localPoints[nL]        = localVertex;
530f5eeb527SMatthew G Knepley                 remotePoints[nL].rank  = neighbor;
531f5eeb527SMatthew G Knepley                 remotePoints[nL].index = remoteVertex;
532f5eeb527SMatthew G Knepley                 ++nL;
53388ed4aceSMatthew G Knepley               } else {
53488ed4aceSMatthew G Knepley                 nleavesCheck += nVz; /* left top vertices */
535f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv, ++nL) {
536f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = (zv*nVy + nVy-1)*nVx +     0 + nC;
537f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (zv*nVy +     0)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
538f5eeb527SMatthew G Knepley 
539f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
540f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
541f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
542f5eeb527SMatthew G Knepley                 }
54388ed4aceSMatthew G Knepley               }
54488ed4aceSMatthew G Knepley             } else {
54588ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
54688ed4aceSMatthew G Knepley                 nleavesCheck += nVy; /* left back vertices */
547f5eeb527SMatthew G Knepley                 for (yv = 0; yv < nVy; ++yv, ++nL) {
548f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = (      0*nVy + yv)*nVx +     0 + nC;
549f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = ((nVz-1)*nVy + yv)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
550f5eeb527SMatthew G Knepley 
551f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
552f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
553f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
554f5eeb527SMatthew G Knepley                 }
55588ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
55688ed4aceSMatthew G Knepley                 nleavesCheck += nVy; /* left front vertices */
557f5eeb527SMatthew G Knepley                 for (yv = 0; yv < nVy; ++yv, ++nL) {
558f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = ((nVz-1)*nVy + yv)*nVx +     0 + nC;
559f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (      0*nVy + yv)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
560f5eeb527SMatthew G Knepley 
561f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
562f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
563f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
564f5eeb527SMatthew G Knepley                 }
56588ed4aceSMatthew G Knepley               } else {
56688ed4aceSMatthew G Knepley                 nleavesCheck += nVy*nVz; /* left vertices */
567f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv) {
568f5eeb527SMatthew G Knepley                   for (yv = 0; yv < nVy; ++yv, ++nL) {
569f5eeb527SMatthew G Knepley                     const PetscInt localVertex  = (zv*nVy + yv)*nVx +     0 + nC;
570f5eeb527SMatthew G Knepley                     const PetscInt remoteVertex = (zv*nVy + yv)*nVx + nVx-1 + nC; /* TODO: Correct this for neighbor sizes */
571f5eeb527SMatthew G Knepley 
572f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
573f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
574f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
575f5eeb527SMatthew G Knepley                   }
576f5eeb527SMatthew G Knepley                 }
57788ed4aceSMatthew G Knepley                 nleavesCheck += nxF;     /* left faces */
578f5eeb527SMatthew G Knepley                 for (xf = 0; xf < nxF; ++xf, ++nL) {
579f5eeb527SMatthew G Knepley                   /* THIS IS WRONG */
580f5eeb527SMatthew G Knepley                   const PetscInt localFace  = 0 + nC+nV;
581f5eeb527SMatthew G Knepley                   const PetscInt remoteFace = 0 + nC+nV;
582f5eeb527SMatthew G Knepley 
583f5eeb527SMatthew G Knepley                   localPoints[nL]        = localFace;
584f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
585f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteFace;
586f5eeb527SMatthew G Knepley                 }
58788ed4aceSMatthew G Knepley               }
58888ed4aceSMatthew G Knepley             }
58988ed4aceSMatthew G Knepley           } else if (xp > 0) { /* right */
59088ed4aceSMatthew G Knepley             if (yp < 0) { /* bottom */
59188ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
592f5eeb527SMatthew G Knepley                 const PetscInt localVertex  = (      0*nVy +     0)*nVx + nVx-1 + nC;
593f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = ((nVz-1)*nVy + nVy-1)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
594628e3b0dSSatish Balay                 nleavesCheck += 1; /* right bottom back vertex */
595f5eeb527SMatthew G Knepley 
596f5eeb527SMatthew G Knepley                 localPoints[nL]        = localVertex;
597f5eeb527SMatthew G Knepley                 remotePoints[nL].rank  = neighbor;
598f5eeb527SMatthew G Knepley                 remotePoints[nL].index = remoteVertex;
599f5eeb527SMatthew G Knepley                 ++nL;
60088ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
601f5eeb527SMatthew G Knepley                 const PetscInt localVertex  = ((nVz-1)*nVy +     0)*nVx + nVx-1 + nC;
602f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = (      0*nVy + nVy-1)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
603628e3b0dSSatish Balay                 nleavesCheck += 1; /* right bottom front vertex */
604f5eeb527SMatthew G Knepley 
605f5eeb527SMatthew G Knepley                 localPoints[nL]        = localVertex;
606f5eeb527SMatthew G Knepley                 remotePoints[nL].rank  = neighbor;
607f5eeb527SMatthew G Knepley                 remotePoints[nL].index = remoteVertex;
608f5eeb527SMatthew G Knepley                 ++nL;
60988ed4aceSMatthew G Knepley               } else {
61088ed4aceSMatthew G Knepley                 nleavesCheck += nVz; /* right bottom vertices */
611f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv, ++nL) {
612f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = (zv*nVy +     0)*nVx + nVx-1 + nC;
613f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (zv*nVy + nVy-1)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
614f5eeb527SMatthew G Knepley 
615f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
616f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
617f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
618f5eeb527SMatthew G Knepley                 }
61988ed4aceSMatthew G Knepley               }
62088ed4aceSMatthew G Knepley             } else if (yp > 0) { /* top */
62188ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
622f5eeb527SMatthew G Knepley                 const PetscInt localVertex  = (      0*nVy + nVy-1)*nVx + nVx-1 + nC;
623f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = ((nVz-1)*nVy +     0)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
624628e3b0dSSatish Balay                 nleavesCheck += 1; /* right top back vertex */
625f5eeb527SMatthew G Knepley 
626f5eeb527SMatthew G Knepley                 localPoints[nL]        = localVertex;
627f5eeb527SMatthew G Knepley                 remotePoints[nL].rank  = neighbor;
628f5eeb527SMatthew G Knepley                 remotePoints[nL].index = remoteVertex;
629f5eeb527SMatthew G Knepley                 ++nL;
63088ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
631f5eeb527SMatthew G Knepley                 const PetscInt localVertex  = ((nVz-1)*nVy + nVy-1)*nVx + nVx-1 + nC;
632f5eeb527SMatthew G Knepley                 const PetscInt remoteVertex = (      0*nVy +     0)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
633628e3b0dSSatish Balay                 nleavesCheck += 1; /* right top front vertex */
634f5eeb527SMatthew G Knepley 
635f5eeb527SMatthew G Knepley                 localPoints[nL]        = localVertex;
636f5eeb527SMatthew G Knepley                 remotePoints[nL].rank  = neighbor;
637f5eeb527SMatthew G Knepley                 remotePoints[nL].index = remoteVertex;
638f5eeb527SMatthew G Knepley                 ++nL;
63988ed4aceSMatthew G Knepley               } else {
64088ed4aceSMatthew G Knepley                 nleavesCheck += nVz; /* right top vertices */
641f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv, ++nL) {
642f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = (zv*nVy + nVy-1)*nVx + nVx-1 + nC;
643f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (zv*nVy +     0)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
644f5eeb527SMatthew G Knepley 
645f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
646f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
647f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
648f5eeb527SMatthew G Knepley                 }
64988ed4aceSMatthew G Knepley               }
65088ed4aceSMatthew G Knepley             } else {
65188ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
65288ed4aceSMatthew G Knepley                 nleavesCheck += nVy; /* right back vertices */
653f5eeb527SMatthew G Knepley                 for (yv = 0; yv < nVy; ++yv, ++nL) {
654f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = (      0*nVy + yv)*nVx + nVx-1 + nC;
655f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = ((nVz-1)*nVy + yv)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
656f5eeb527SMatthew G Knepley 
657f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
658f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
659f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
660f5eeb527SMatthew G Knepley                 }
66188ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
66288ed4aceSMatthew G Knepley                 nleavesCheck += nVy; /* right front vertices */
663f5eeb527SMatthew G Knepley                 for (yv = 0; yv < nVy; ++yv, ++nL) {
664f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = ((nVz-1)*nVy + yv)*nVx + nVx-1 + nC;
665f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (      0*nVy + yv)*nVx +     0 + nC; /* TODO: Correct this for neighbor sizes */
666f5eeb527SMatthew G Knepley 
667f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
668f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
669f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
670f5eeb527SMatthew G Knepley                 }
67188ed4aceSMatthew G Knepley               } else {
67288ed4aceSMatthew G Knepley                 nleavesCheck += nVy*nVz; /* right vertices */
673f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv) {
674f5eeb527SMatthew G Knepley                   for (yv = 0; yv < nVy; ++yv, ++nL) {
675f5eeb527SMatthew G Knepley                     const PetscInt localVertex  = (zv*nVy + yv)*nVx + nVx-1 + nC;
676f5eeb527SMatthew G Knepley                     const PetscInt remoteVertex = (zv*nVy + yv)*nVx + 0     + nC; /* TODO: Correct this for neighbor sizes */
677f5eeb527SMatthew G Knepley 
678f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
679f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
680f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
681f5eeb527SMatthew G Knepley                   }
682f5eeb527SMatthew G Knepley                 }
68388ed4aceSMatthew G Knepley                 nleavesCheck += nxF;     /* right faces */
684f5eeb527SMatthew G Knepley                 for (xf = 0; xf < nxF; ++xf, ++nL) {
685f5eeb527SMatthew G Knepley                   /* THIS IS WRONG */
686f5eeb527SMatthew G Knepley                   const PetscInt localFace  = 0 + nC+nV;
687f5eeb527SMatthew G Knepley                   const PetscInt remoteFace = 0 + nC+nV;
688f5eeb527SMatthew G Knepley 
689f5eeb527SMatthew G Knepley                   localPoints[nL]        = localFace;
690f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
691f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteFace;
692f5eeb527SMatthew G Knepley                 }
69388ed4aceSMatthew G Knepley               }
69488ed4aceSMatthew G Knepley             }
69588ed4aceSMatthew G Knepley           } else {
69688ed4aceSMatthew G Knepley             if (yp < 0) { /* bottom */
69788ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
69888ed4aceSMatthew G Knepley                 nleavesCheck += nVx; /* bottom back vertices */
699f5eeb527SMatthew G Knepley                 for (xv = 0; xv < nVx; ++xv, ++nL) {
700f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = (      0*nVy +     0)*nVx + xv + nC;
701f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = ((nVz-1)*nVy + nVy-1)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
702f5eeb527SMatthew G Knepley 
703f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
704f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
705f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
706f5eeb527SMatthew G Knepley                 }
70788ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
70888ed4aceSMatthew G Knepley                 nleavesCheck += nVx; /* bottom front vertices */
709f5eeb527SMatthew G Knepley                 for (xv = 0; xv < nVx; ++xv, ++nL) {
710f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = ((nVz-1)*nVy +     0)*nVx + xv + nC;
711f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (      0*nVy + nVy-1)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
712f5eeb527SMatthew G Knepley 
713f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
714f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
715f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
716f5eeb527SMatthew G Knepley                 }
71788ed4aceSMatthew G Knepley               } else {
71888ed4aceSMatthew G Knepley                 nleavesCheck += nVx*nVz; /* bottom vertices */
719f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv) {
720f5eeb527SMatthew G Knepley                   for (xv = 0; xv < nVx; ++xv, ++nL) {
721f5eeb527SMatthew G Knepley                     const PetscInt localVertex  = (zv*nVy +     0)*nVx + xv + nC;
722f5eeb527SMatthew G Knepley                     const PetscInt remoteVertex = (zv*nVy + nVy-1)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
723f5eeb527SMatthew G Knepley 
724f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
725f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
726f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
727f5eeb527SMatthew G Knepley                   }
728f5eeb527SMatthew G Knepley                 }
72988ed4aceSMatthew G Knepley                 nleavesCheck += nyF;     /* bottom faces */
730f5eeb527SMatthew G Knepley                 for (yf = 0; yf < nyF; ++yf, ++nL) {
731f5eeb527SMatthew G Knepley                   /* THIS IS WRONG */
732f5eeb527SMatthew G Knepley                   const PetscInt localFace  = 0 + nC+nV;
733f5eeb527SMatthew G Knepley                   const PetscInt remoteFace = 0 + nC+nV;
734f5eeb527SMatthew G Knepley 
735f5eeb527SMatthew G Knepley                   localPoints[nL]        = localFace;
736f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
737f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteFace;
738f5eeb527SMatthew G Knepley                 }
73988ed4aceSMatthew G Knepley               }
74088ed4aceSMatthew G Knepley             } else if (yp > 0) { /* top */
74188ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
74288ed4aceSMatthew G Knepley                 nleavesCheck += nVx; /* top back vertices */
743f5eeb527SMatthew G Knepley                 for (xv = 0; xv < nVx; ++xv, ++nL) {
744f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = (      0*nVy + nVy-1)*nVx + xv + nC;
745f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = ((nVz-1)*nVy +     0)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
746f5eeb527SMatthew G Knepley 
747f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
748f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
749f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
750f5eeb527SMatthew G Knepley                 }
75188ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
75288ed4aceSMatthew G Knepley                 nleavesCheck += nVx; /* top front vertices */
753f5eeb527SMatthew G Knepley                 for (xv = 0; xv < nVx; ++xv, ++nL) {
754f5eeb527SMatthew G Knepley                   const PetscInt localVertex  = ((nVz-1)*nVy + nVy-1)*nVx + xv + nC;
755f5eeb527SMatthew G Knepley                   const PetscInt remoteVertex = (      0*nVy +     0)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
756f5eeb527SMatthew G Knepley 
757f5eeb527SMatthew G Knepley                   localPoints[nL]        = localVertex;
758f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
759f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteVertex;
760f5eeb527SMatthew G Knepley                 }
76188ed4aceSMatthew G Knepley               } else {
76288ed4aceSMatthew G Knepley                 nleavesCheck += nVx*nVz; /* top vertices */
763f5eeb527SMatthew G Knepley                 for (zv = 0; zv < nVz; ++zv) {
764f5eeb527SMatthew G Knepley                   for (xv = 0; xv < nVx; ++xv, ++nL) {
765f5eeb527SMatthew G Knepley                     const PetscInt localVertex  = (zv*nVy + nVy-1)*nVx + xv + nC;
766f5eeb527SMatthew G Knepley                     const PetscInt remoteVertex = (zv*nVy +     0)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
767f5eeb527SMatthew G Knepley 
768f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
769f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
770f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
771f5eeb527SMatthew G Knepley                   }
772f5eeb527SMatthew G Knepley                 }
77388ed4aceSMatthew G Knepley                 nleavesCheck += nyF;     /* top faces */
774f5eeb527SMatthew G Knepley                 for (yf = 0; yf < nyF; ++yf, ++nL) {
775f5eeb527SMatthew G Knepley                   /* THIS IS WRONG */
776f5eeb527SMatthew G Knepley                   const PetscInt localFace  = 0 + nC+nV;
777f5eeb527SMatthew G Knepley                   const PetscInt remoteFace = 0 + nC+nV;
778f5eeb527SMatthew G Knepley 
779f5eeb527SMatthew G Knepley                   localPoints[nL]        = localFace;
780f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
781f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteFace;
782f5eeb527SMatthew G Knepley                 }
78388ed4aceSMatthew G Knepley               }
78488ed4aceSMatthew G Knepley             } else {
78588ed4aceSMatthew G Knepley               if (zp < 0) { /* back */
78688ed4aceSMatthew G Knepley                 nleavesCheck += nVx*nVy; /* back vertices */
787f5eeb527SMatthew G Knepley                 for (yv = 0; yv < nVy; ++yv) {
788f5eeb527SMatthew G Knepley                   for (xv = 0; xv < nVx; ++xv, ++nL) {
789f5eeb527SMatthew G Knepley                     const PetscInt localVertex  = (      0*nVy + yv)*nVx + xv + nC;
790f5eeb527SMatthew G Knepley                     const PetscInt remoteVertex = ((nVz-1)*nVy + yv)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
791f5eeb527SMatthew G Knepley 
792f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
793f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
794f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
795f5eeb527SMatthew G Knepley                   }
796f5eeb527SMatthew G Knepley                 }
79788ed4aceSMatthew G Knepley                 nleavesCheck += nzF;     /* back faces */
798f5eeb527SMatthew G Knepley                 for (zf = 0; zf < nzF; ++zf, ++nL) {
799f5eeb527SMatthew G Knepley                   /* THIS IS WRONG */
800f5eeb527SMatthew G Knepley                   const PetscInt localFace  = 0 + nC+nV;
801f5eeb527SMatthew G Knepley                   const PetscInt remoteFace = 0 + nC+nV;
802f5eeb527SMatthew G Knepley 
803f5eeb527SMatthew G Knepley                   localPoints[nL]        = localFace;
804f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
805f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteFace;
806f5eeb527SMatthew G Knepley                 }
80788ed4aceSMatthew G Knepley               } else if (zp > 0) { /* front */
80888ed4aceSMatthew G Knepley                 nleavesCheck += nVx*nVy; /* front vertices */
809f5eeb527SMatthew G Knepley                 for (yv = 0; yv < nVy; ++yv) {
810f5eeb527SMatthew G Knepley                   for (xv = 0; xv < nVx; ++xv, ++nL) {
811f5eeb527SMatthew G Knepley                     const PetscInt localVertex  = ((nVz-1)*nVy + yv)*nVx + xv + nC;
812f5eeb527SMatthew G Knepley                     const PetscInt remoteVertex = (      0*nVy + yv)*nVx + xv + nC; /* TODO: Correct this for neighbor sizes */
813f5eeb527SMatthew G Knepley 
814f5eeb527SMatthew G Knepley                     localPoints[nL]        = localVertex;
815f5eeb527SMatthew G Knepley                     remotePoints[nL].rank  = neighbor;
816f5eeb527SMatthew G Knepley                     remotePoints[nL].index = remoteVertex;
817f5eeb527SMatthew G Knepley                   }
818f5eeb527SMatthew G Knepley                 }
81988ed4aceSMatthew G Knepley                 nleavesCheck += nzF;     /* front faces */
820f5eeb527SMatthew G Knepley                 for (zf = 0; zf < nzF; ++zf, ++nL) {
821f5eeb527SMatthew G Knepley                   /* THIS IS WRONG */
822f5eeb527SMatthew G Knepley                   const PetscInt localFace  = 0 + nC+nV;
823f5eeb527SMatthew G Knepley                   const PetscInt remoteFace = 0 + nC+nV;
824f5eeb527SMatthew G Knepley 
825f5eeb527SMatthew G Knepley                   localPoints[nL]        = localFace;
826f5eeb527SMatthew G Knepley                   remotePoints[nL].rank  = neighbor;
827f5eeb527SMatthew G Knepley                   remotePoints[nL].index = remoteFace;
828f5eeb527SMatthew G Knepley                 }
82988ed4aceSMatthew G Knepley               } else {
83088ed4aceSMatthew G Knepley                 /* Nothing is shared from the interior */
83188ed4aceSMatthew G Knepley               }
83288ed4aceSMatthew G Knepley             }
83388ed4aceSMatthew G Knepley           }
83488ed4aceSMatthew G Knepley         }
83588ed4aceSMatthew G Knepley       }
83688ed4aceSMatthew G Knepley     }
83788ed4aceSMatthew G Knepley   }
8383814d604SMatthew G Knepley   /* TODO: Remove duplication in leaf determination */
83982f516ccSBarry 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);
84082f516ccSBarry Smith   ierr = PetscSFCreate(PetscObjectComm((PetscObject)dm), &sf);CHKERRQ(ierr);
84188ed4aceSMatthew G Knepley   ierr = PetscSFSetGraph(sf, pEnd, nleaves, localPoints, PETSC_OWN_POINTER, remotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr);
842a4b60ecfSMatthew G. Knepley   ierr = DMSetPointSF(dm, sf);CHKERRQ(ierr);
84388ed4aceSMatthew G Knepley   ierr = PetscSFDestroy(&sf);CHKERRQ(ierr);
844affa55c7SMatthew G. Knepley   *s = section;
845affa55c7SMatthew G. Knepley   PetscFunctionReturn(0);
846affa55c7SMatthew G. Knepley }
847affa55c7SMatthew G. Knepley 
8483385ff7eSMatthew G. Knepley #undef __FUNCT__
8493385ff7eSMatthew G. Knepley #define __FUNCT__ "DMDASetVertexCoordinates"
8503385ff7eSMatthew G. Knepley PetscErrorCode DMDASetVertexCoordinates(DM dm, PetscReal xl, PetscReal xu, PetscReal yl, PetscReal yu, PetscReal zl, PetscReal zu)
8513385ff7eSMatthew G. Knepley {
8523385ff7eSMatthew G. Knepley   DM_DA         *da = (DM_DA *) dm->data;
8533385ff7eSMatthew G. Knepley   Vec            coordinates;
8543385ff7eSMatthew G. Knepley   PetscSection   section;
8553385ff7eSMatthew G. Knepley   PetscScalar   *coords;
8563385ff7eSMatthew G. Knepley   PetscReal      h[3];
8573385ff7eSMatthew G. Knepley   PetscInt       dim, size, M, N, P, nVx, nVy, nVz, nV, vStart, vEnd, v, i, j, k;
8583385ff7eSMatthew G. Knepley   PetscErrorCode ierr;
8593385ff7eSMatthew G. Knepley 
8603385ff7eSMatthew G. Knepley   PetscFunctionBegin;
8613385ff7eSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8623385ff7eSMatthew G. Knepley   ierr = DMDAGetInfo(dm, &dim, &M, &N, &P, 0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
8633385ff7eSMatthew G. Knepley   h[0] = (xu - xl)/M;
8643385ff7eSMatthew G. Knepley   h[1] = (yu - yl)/N;
8653385ff7eSMatthew G. Knepley   h[2] = (zu - zl)/P;
8663385ff7eSMatthew G. Knepley   ierr = DMDAGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
8673385ff7eSMatthew G. Knepley   ierr = DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV);CHKERRQ(ierr);
8683385ff7eSMatthew G. Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &section);CHKERRQ(ierr);
8693385ff7eSMatthew G. Knepley   ierr = PetscSectionSetNumFields(section, 1);CHKERRQ(ierr);
8703385ff7eSMatthew G. Knepley   ierr = PetscSectionSetFieldComponents(section, 0, dim);CHKERRQ(ierr);
8713385ff7eSMatthew G. Knepley   ierr = PetscSectionSetChart(section, vStart, vEnd);CHKERRQ(ierr);
8723385ff7eSMatthew G. Knepley   for (v = vStart; v < vEnd; ++v) {
8733385ff7eSMatthew G. Knepley     ierr = PetscSectionSetDof(section, v, dim);CHKERRQ(ierr);
8743385ff7eSMatthew G. Knepley   }
8753385ff7eSMatthew G. Knepley   ierr = PetscSectionSetUp(section);CHKERRQ(ierr);
8763385ff7eSMatthew G. Knepley   ierr = PetscSectionGetStorageSize(section, &size);CHKERRQ(ierr);
8773385ff7eSMatthew G. Knepley   ierr = VecCreateSeq(PETSC_COMM_SELF, size, &coordinates);CHKERRQ(ierr);
8783385ff7eSMatthew G. Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
8793385ff7eSMatthew G. Knepley   for (k = 0; k < nVz; ++k) {
880*e4d69e08SBarry Smith     PetscInt ind[3], d, off;
8813385ff7eSMatthew G. Knepley 
882*e4d69e08SBarry Smith     ind[0] = 0;
883*e4d69e08SBarry Smith     ind[1] = 0;
884*e4d69e08SBarry Smith     ind[2] = k + da->zs;
8853385ff7eSMatthew G. Knepley     for (j = 0; j < nVy; ++j) {
8863385ff7eSMatthew G. Knepley       ind[1] = j + da->ys;
8873385ff7eSMatthew G. Knepley       for (i = 0; i < nVx; ++i) {
8883385ff7eSMatthew G. Knepley         const PetscInt vertex = (k*nVy + j)*nVx + i + vStart;
8893385ff7eSMatthew G. Knepley 
8903385ff7eSMatthew G. Knepley         ierr = PetscSectionGetOffset(section, vertex, &off);CHKERRQ(ierr);
8913385ff7eSMatthew G. Knepley         ind[0] = i + da->xs;
8923385ff7eSMatthew G. Knepley         for (d = 0; d < dim; ++d) {
8933385ff7eSMatthew G. Knepley           coords[off+d] = h[d]*ind[d];
8943385ff7eSMatthew G. Knepley         }
8953385ff7eSMatthew G. Knepley       }
8963385ff7eSMatthew G. Knepley     }
8973385ff7eSMatthew G. Knepley   }
8983385ff7eSMatthew G. Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
8993385ff7eSMatthew G. Knepley   ierr = DMSetCoordinateSection(dm, section);CHKERRQ(ierr);
9003385ff7eSMatthew G. Knepley   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
901a4b60ecfSMatthew G. Knepley   ierr = PetscSectionDestroy(&section);CHKERRQ(ierr);
9023385ff7eSMatthew G. Knepley   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
9033385ff7eSMatthew G. Knepley   PetscFunctionReturn(0);
9043385ff7eSMatthew G. Knepley }
9059a800dd8SMatthew G. Knepley 
9069a800dd8SMatthew G. Knepley #undef __FUNCT__
9079a800dd8SMatthew G. Knepley #define __FUNCT__ "DMDAProjectFunctionLocal"
9089a800dd8SMatthew G. Knepley PetscErrorCode DMDAProjectFunctionLocal(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *), InsertMode mode, Vec localX)
9099a800dd8SMatthew G. Knepley {
9109a800dd8SMatthew G. Knepley   PetscDualSpace *sp;
9119a800dd8SMatthew G. Knepley   PetscQuadrature q;
9129a800dd8SMatthew G. Knepley   PetscSection    section;
9139a800dd8SMatthew G. Knepley   PetscScalar    *values;
9149a800dd8SMatthew G. Knepley   PetscReal      *v0, *J, *detJ;
9159a800dd8SMatthew G. Knepley   PetscInt        numFields, numComp, dim, spDim, totDim = 0, numValues, cStart, cEnd, f, c, v, d;
9169a800dd8SMatthew G. Knepley   PetscErrorCode  ierr;
9179a800dd8SMatthew G. Knepley 
9189a800dd8SMatthew G. Knepley   PetscFunctionBegin;
9199a800dd8SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
9209a800dd8SMatthew G. Knepley   ierr = PetscFEGetQuadrature(fe[0], &q);CHKERRQ(ierr);
9219a800dd8SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
9229a800dd8SMatthew G. Knepley   ierr = PetscMalloc(numFields * sizeof(PetscDualSpace), &sp);CHKERRQ(ierr);
9239a800dd8SMatthew G. Knepley   for (f = 0; f < numFields; ++f) {
9249a800dd8SMatthew G. Knepley     ierr = PetscFEGetDualSpace(fe[f], &sp[f]);CHKERRQ(ierr);
9259a800dd8SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr);
9269a800dd8SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
9279a800dd8SMatthew G. Knepley     totDim += spDim*numComp;
9289a800dd8SMatthew G. Knepley   }
9299a800dd8SMatthew G. Knepley   ierr = DMDAGetInfo(dm, &dim,0,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
9309a800dd8SMatthew G. Knepley   ierr = DMDAGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
9319a800dd8SMatthew G. Knepley   ierr = DMDAVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr);
9329a800dd8SMatthew 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);
9339a800dd8SMatthew G. Knepley   ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
934d7a12f93SMatthew G. Knepley   ierr = PetscMalloc3(dim*q.numPoints,&v0,dim*dim*q.numPoints,&J,q.numPoints,&detJ);CHKERRQ(ierr);
9359a800dd8SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
9369a800dd8SMatthew G. Knepley     PetscCellGeometry geom;
9379a800dd8SMatthew G. Knepley 
9389a800dd8SMatthew G. Knepley     ierr = DMDAComputeCellGeometry(dm, c, &q, v0, J, NULL, detJ);CHKERRQ(ierr);
9399a800dd8SMatthew G. Knepley     geom.v0   = v0;
9409a800dd8SMatthew G. Knepley     geom.J    = J;
9419a800dd8SMatthew G. Knepley     geom.detJ = detJ;
9429a800dd8SMatthew G. Knepley     for (f = 0, v = 0; f < numFields; ++f) {
9439a800dd8SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr);
9449a800dd8SMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
9459a800dd8SMatthew G. Knepley       for (d = 0; d < spDim; ++d) {
9469a800dd8SMatthew G. Knepley         ierr = PetscDualSpaceApply(sp[f], d, geom, numComp, funcs[f], &values[v]);CHKERRQ(ierr);
9479a800dd8SMatthew G. Knepley         v += numComp;
9489a800dd8SMatthew G. Knepley       }
9499a800dd8SMatthew G. Knepley     }
9509a800dd8SMatthew G. Knepley     ierr = DMDAVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr);
9519a800dd8SMatthew G. Knepley   }
9529a800dd8SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
9539a800dd8SMatthew G. Knepley   ierr = PetscFree3(v0,J,detJ);CHKERRQ(ierr);
9549a800dd8SMatthew G. Knepley   ierr = PetscFree(sp);CHKERRQ(ierr);
9559a800dd8SMatthew G. Knepley   PetscFunctionReturn(0);
9569a800dd8SMatthew G. Knepley }
9579a800dd8SMatthew G. Knepley 
9589a800dd8SMatthew G. Knepley #undef __FUNCT__
9599a800dd8SMatthew G. Knepley #define __FUNCT__ "DMDAProjectFunction"
9609a800dd8SMatthew G. Knepley /*@C
9619a800dd8SMatthew G. Knepley   DMDAProjectFunction - This projects the given function into the function space provided.
9629a800dd8SMatthew G. Knepley 
9639a800dd8SMatthew G. Knepley   Input Parameters:
9649a800dd8SMatthew G. Knepley + dm      - The DM
9659a800dd8SMatthew G. Knepley . fe      - The PetscFE associated with the field
9669a800dd8SMatthew G. Knepley . funcs   - The coordinate functions to evaluate
9679a800dd8SMatthew G. Knepley - mode    - The insertion mode for values
9689a800dd8SMatthew G. Knepley 
9699a800dd8SMatthew G. Knepley   Output Parameter:
9709a800dd8SMatthew G. Knepley . X - vector
9719a800dd8SMatthew G. Knepley 
9729a800dd8SMatthew G. Knepley   Level: developer
9739a800dd8SMatthew G. Knepley 
9749a800dd8SMatthew G. Knepley .seealso: DMDAComputeL2Diff()
9759a800dd8SMatthew G. Knepley @*/
9769a800dd8SMatthew G. Knepley PetscErrorCode DMDAProjectFunction(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *), InsertMode mode, Vec X)
9779a800dd8SMatthew G. Knepley {
9789a800dd8SMatthew G. Knepley   Vec            localX;
9799a800dd8SMatthew G. Knepley   PetscErrorCode ierr;
9809a800dd8SMatthew G. Knepley 
9819a800dd8SMatthew G. Knepley   PetscFunctionBegin;
9829a800dd8SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
9839a800dd8SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
9849a800dd8SMatthew G. Knepley   ierr = DMDAProjectFunctionLocal(dm, fe, funcs, mode, localX);CHKERRQ(ierr);
9859a800dd8SMatthew G. Knepley   ierr = DMLocalToGlobalBegin(dm, localX, mode, X);CHKERRQ(ierr);
9869a800dd8SMatthew G. Knepley   ierr = DMLocalToGlobalEnd(dm, localX, mode, X);CHKERRQ(ierr);
9879a800dd8SMatthew G. Knepley   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
9889a800dd8SMatthew G. Knepley   PetscFunctionReturn(0);
9899a800dd8SMatthew G. Knepley }
9909a800dd8SMatthew G. Knepley 
9919a800dd8SMatthew G. Knepley #undef __FUNCT__
9929a800dd8SMatthew G. Knepley #define __FUNCT__ "DMDAComputeL2Diff"
9939a800dd8SMatthew G. Knepley /*@C
9949a800dd8SMatthew G. Knepley   DMDAComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h.
9959a800dd8SMatthew G. Knepley 
9969a800dd8SMatthew G. Knepley   Input Parameters:
9979a800dd8SMatthew G. Knepley + dm    - The DM
9989a800dd8SMatthew G. Knepley . fe    - The PetscFE object for each field
9999a800dd8SMatthew G. Knepley . funcs - The functions to evaluate for each field component
10009a800dd8SMatthew G. Knepley - X     - The coefficient vector u_h
10019a800dd8SMatthew G. Knepley 
10029a800dd8SMatthew G. Knepley   Output Parameter:
10039a800dd8SMatthew G. Knepley . diff - The diff ||u - u_h||_2
10049a800dd8SMatthew G. Knepley 
10059a800dd8SMatthew G. Knepley   Level: developer
10069a800dd8SMatthew G. Knepley 
10079a800dd8SMatthew G. Knepley .seealso: DMDAProjectFunction()
10089a800dd8SMatthew G. Knepley @*/
10099a800dd8SMatthew G. Knepley PetscErrorCode DMDAComputeL2Diff(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *), Vec X, PetscReal *diff)
10109a800dd8SMatthew G. Knepley {
10119a800dd8SMatthew G. Knepley   const PetscInt  debug = 0;
10129a800dd8SMatthew G. Knepley   PetscSection    section;
10139a800dd8SMatthew G. Knepley   PetscQuadrature quad;
10149a800dd8SMatthew G. Knepley   Vec             localX;
10159a800dd8SMatthew G. Knepley   PetscScalar    *funcVal;
10169a800dd8SMatthew G. Knepley   PetscReal      *coords, *v0, *J, *invJ, detJ;
10179a800dd8SMatthew G. Knepley   PetscReal       localDiff = 0.0;
10189a800dd8SMatthew G. Knepley   PetscInt        dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp;
10199a800dd8SMatthew G. Knepley   PetscErrorCode  ierr;
10209a800dd8SMatthew G. Knepley 
10219a800dd8SMatthew G. Knepley   PetscFunctionBegin;
10229a800dd8SMatthew G. Knepley   ierr = DMDAGetInfo(dm, &dim,0,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
10239a800dd8SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
10249a800dd8SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
10259a800dd8SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
10269a800dd8SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
10279a800dd8SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
10289a800dd8SMatthew G. Knepley   for (field = 0; field < numFields; ++field) {
10299a800dd8SMatthew G. Knepley     PetscInt Nc;
10309a800dd8SMatthew G. Knepley 
10319a800dd8SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr);
10329a800dd8SMatthew G. Knepley     numComponents += Nc;
10339a800dd8SMatthew G. Knepley   }
10349a800dd8SMatthew G. Knepley   /* There are no BC values in DAs right now: ierr = DMDAProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */
1035d7a12f93SMatthew G. Knepley   ierr = PetscMalloc5(numComponents,&funcVal,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
10369a800dd8SMatthew G. Knepley   ierr = DMDAGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
10379a800dd8SMatthew G. Knepley   ierr = PetscFEGetQuadrature(fe[0], &quad);CHKERRQ(ierr);
10389a800dd8SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
10399a800dd8SMatthew G. Knepley     PetscScalar *x = NULL;
10409a800dd8SMatthew G. Knepley     PetscReal    elemDiff = 0.0;
10419a800dd8SMatthew G. Knepley 
10429a800dd8SMatthew G. Knepley     ierr = DMDAComputeCellGeometry(dm, c, &quad, v0, J, invJ, &detJ);CHKERRQ(ierr);
10439a800dd8SMatthew G. Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
10449a800dd8SMatthew G. Knepley     ierr = DMDAVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
10459a800dd8SMatthew G. Knepley 
10469a800dd8SMatthew G. Knepley     for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) {
10479a800dd8SMatthew G. Knepley       const PetscInt   numQuadPoints = quad.numPoints;
10489a800dd8SMatthew G. Knepley       const PetscReal *quadPoints    = quad.points;
10499a800dd8SMatthew G. Knepley       const PetscReal *quadWeights   = quad.weights;
10509a800dd8SMatthew G. Knepley       PetscReal       *basis;
10519a800dd8SMatthew G. Knepley       PetscInt         numBasisFuncs, numBasisComps, q, d, e, fc, f;
10529a800dd8SMatthew G. Knepley 
10539a800dd8SMatthew G. Knepley       ierr = PetscFEGetDimension(fe[field], &numBasisFuncs);CHKERRQ(ierr);
10549a800dd8SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe[field], &numBasisComps);CHKERRQ(ierr);
10559a800dd8SMatthew G. Knepley       ierr = PetscFEGetDefaultTabulation(fe[field], &basis, NULL, NULL);CHKERRQ(ierr);
10569a800dd8SMatthew G. Knepley       if (debug) {
10579a800dd8SMatthew G. Knepley         char title[1024];
10589a800dd8SMatthew G. Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
10599a800dd8SMatthew G. Knepley         ierr = DMPrintCellVector(c, title, numBasisFuncs*numBasisComps, &x[fieldOffset]);CHKERRQ(ierr);
10609a800dd8SMatthew G. Knepley       }
10619a800dd8SMatthew G. Knepley       for (q = 0; q < numQuadPoints; ++q) {
10629a800dd8SMatthew G. Knepley         for (d = 0; d < dim; d++) {
10639a800dd8SMatthew G. Knepley           coords[d] = v0[d];
10649a800dd8SMatthew G. Knepley           for (e = 0; e < dim; e++) {
10659a800dd8SMatthew G. Knepley             coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0);
10669a800dd8SMatthew G. Knepley           }
10679a800dd8SMatthew G. Knepley         }
10689a800dd8SMatthew G. Knepley         (*funcs[field])(coords, funcVal);
10699a800dd8SMatthew G. Knepley         for (fc = 0; fc < numBasisComps; ++fc) {
10709a800dd8SMatthew G. Knepley           PetscScalar interpolant = 0.0;
10719a800dd8SMatthew G. Knepley 
10729a800dd8SMatthew G. Knepley           for (f = 0; f < numBasisFuncs; ++f) {
10739a800dd8SMatthew G. Knepley             const PetscInt fidx = f*numBasisComps+fc;
10749a800dd8SMatthew G. Knepley             interpolant += x[fieldOffset+fidx]*basis[q*numBasisFuncs*numBasisComps+fidx];
10759a800dd8SMatthew G. Knepley           }
10769a800dd8SMatthew 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);}
10779a800dd8SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ;
10789a800dd8SMatthew G. Knepley         }
10799a800dd8SMatthew G. Knepley       }
10809a800dd8SMatthew G. Knepley       comp        += numBasisComps;
10819a800dd8SMatthew G. Knepley       fieldOffset += numBasisFuncs*numBasisComps;
10829a800dd8SMatthew G. Knepley     }
10839a800dd8SMatthew G. Knepley     ierr = DMDAVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
10849a800dd8SMatthew G. Knepley     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
10859a800dd8SMatthew G. Knepley     localDiff += elemDiff;
10869a800dd8SMatthew G. Knepley   }
10879a800dd8SMatthew G. Knepley   ierr  = PetscFree5(funcVal,coords,v0,J,invJ);CHKERRQ(ierr);
10889a800dd8SMatthew G. Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
10899a800dd8SMatthew G. Knepley   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
10909a800dd8SMatthew G. Knepley   *diff = PetscSqrtReal(*diff);
1091a66d4d66SMatthew G Knepley   PetscFunctionReturn(0);
1092a66d4d66SMatthew G Knepley }
1093a66d4d66SMatthew G Knepley 
109447c6ae99SBarry Smith /* ------------------------------------------------------------------- */
109547c6ae99SBarry Smith 
109647c6ae99SBarry Smith #undef __FUNCT__
1097aa219208SBarry Smith #define __FUNCT__ "DMDAGetArray"
109847c6ae99SBarry Smith /*@C
1099aa219208SBarry Smith      DMDAGetArray - Gets a work array for a DMDA
110047c6ae99SBarry Smith 
110147c6ae99SBarry Smith     Input Parameter:
110247c6ae99SBarry Smith +    da - information about my local patch
110347c6ae99SBarry Smith -    ghosted - do you want arrays for the ghosted or nonghosted patch
110447c6ae99SBarry Smith 
110547c6ae99SBarry Smith     Output Parameters:
110647c6ae99SBarry Smith .    vptr - array data structured
110747c6ae99SBarry Smith 
110847c6ae99SBarry Smith     Note:  The vector values are NOT initialized and may have garbage in them, so you may need
110947c6ae99SBarry Smith            to zero them.
111047c6ae99SBarry Smith 
111147c6ae99SBarry Smith   Level: advanced
111247c6ae99SBarry Smith 
1113bcaeba4dSBarry Smith .seealso: DMDARestoreArray()
111447c6ae99SBarry Smith 
111547c6ae99SBarry Smith @*/
11167087cfbeSBarry Smith PetscErrorCode  DMDAGetArray(DM da,PetscBool ghosted,void *vptr)
111747c6ae99SBarry Smith {
111847c6ae99SBarry Smith   PetscErrorCode ierr;
111947c6ae99SBarry Smith   PetscInt       j,i,xs,ys,xm,ym,zs,zm;
112047c6ae99SBarry Smith   char           *iarray_start;
112147c6ae99SBarry Smith   void           **iptr = (void**)vptr;
112247c6ae99SBarry Smith   DM_DA          *dd    = (DM_DA*)da->data;
112347c6ae99SBarry Smith 
112447c6ae99SBarry Smith   PetscFunctionBegin;
112547c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
112647c6ae99SBarry Smith   if (ghosted) {
1127aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
112847c6ae99SBarry Smith       if (dd->arrayghostedin[i]) {
112947c6ae99SBarry Smith         *iptr                 = dd->arrayghostedin[i];
113047c6ae99SBarry Smith         iarray_start          = (char*)dd->startghostedin[i];
11310298fd71SBarry Smith         dd->arrayghostedin[i] = NULL;
11320298fd71SBarry Smith         dd->startghostedin[i] = NULL;
113347c6ae99SBarry Smith 
113447c6ae99SBarry Smith         goto done;
113547c6ae99SBarry Smith       }
113647c6ae99SBarry Smith     }
113747c6ae99SBarry Smith     xs = dd->Xs;
113847c6ae99SBarry Smith     ys = dd->Ys;
113947c6ae99SBarry Smith     zs = dd->Zs;
114047c6ae99SBarry Smith     xm = dd->Xe-dd->Xs;
114147c6ae99SBarry Smith     ym = dd->Ye-dd->Ys;
114247c6ae99SBarry Smith     zm = dd->Ze-dd->Zs;
114347c6ae99SBarry Smith   } else {
1144aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
114547c6ae99SBarry Smith       if (dd->arrayin[i]) {
114647c6ae99SBarry Smith         *iptr          = dd->arrayin[i];
114747c6ae99SBarry Smith         iarray_start   = (char*)dd->startin[i];
11480298fd71SBarry Smith         dd->arrayin[i] = NULL;
11490298fd71SBarry Smith         dd->startin[i] = NULL;
115047c6ae99SBarry Smith 
115147c6ae99SBarry Smith         goto done;
115247c6ae99SBarry Smith       }
115347c6ae99SBarry Smith     }
115447c6ae99SBarry Smith     xs = dd->xs;
115547c6ae99SBarry Smith     ys = dd->ys;
115647c6ae99SBarry Smith     zs = dd->zs;
115747c6ae99SBarry Smith     xm = dd->xe-dd->xs;
115847c6ae99SBarry Smith     ym = dd->ye-dd->ys;
115947c6ae99SBarry Smith     zm = dd->ze-dd->zs;
116047c6ae99SBarry Smith   }
116147c6ae99SBarry Smith 
116247c6ae99SBarry Smith   switch (dd->dim) {
116347c6ae99SBarry Smith   case 1: {
116447c6ae99SBarry Smith     void *ptr;
116547c6ae99SBarry Smith 
1166901f1932SJed Brown     ierr = PetscMalloc(xm*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
116747c6ae99SBarry Smith 
116847c6ae99SBarry Smith     ptr   = (void*)(iarray_start - xs*sizeof(PetscScalar));
116947c6ae99SBarry Smith     *iptr = (void*)ptr;
11708865f1eaSKarl Rupp     break;
11718865f1eaSKarl Rupp   }
117247c6ae99SBarry Smith   case 2: {
117347c6ae99SBarry Smith     void **ptr;
117447c6ae99SBarry Smith 
1175901f1932SJed Brown     ierr = PetscMalloc((ym+1)*sizeof(void*)+xm*ym*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
117647c6ae99SBarry Smith 
117747c6ae99SBarry Smith     ptr = (void**)(iarray_start + xm*ym*sizeof(PetscScalar) - ys*sizeof(void*));
11788865f1eaSKarl Rupp     for (j=ys; j<ys+ym; j++) ptr[j] = iarray_start + sizeof(PetscScalar)*(xm*(j-ys) - xs);
117947c6ae99SBarry Smith     *iptr = (void*)ptr;
11808865f1eaSKarl Rupp     break;
11818865f1eaSKarl Rupp   }
118247c6ae99SBarry Smith   case 3: {
118347c6ae99SBarry Smith     void ***ptr,**bptr;
118447c6ae99SBarry Smith 
1185901f1932SJed Brown     ierr = PetscMalloc((zm+1)*sizeof(void**)+(ym*zm+1)*sizeof(void*)+xm*ym*zm*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
118647c6ae99SBarry Smith 
118747c6ae99SBarry Smith     ptr  = (void***)(iarray_start + xm*ym*zm*sizeof(PetscScalar) - zs*sizeof(void*));
118847c6ae99SBarry Smith     bptr = (void**)(iarray_start + xm*ym*zm*sizeof(PetscScalar) + zm*sizeof(void**));
11898865f1eaSKarl Rupp     for (i=zs; i<zs+zm; i++) ptr[i] = bptr + ((i-zs)*ym - ys);
119047c6ae99SBarry Smith     for (i=zs; i<zs+zm; i++) {
119147c6ae99SBarry Smith       for (j=ys; j<ys+ym; j++) {
119247c6ae99SBarry Smith         ptr[i][j] = iarray_start + sizeof(PetscScalar)*(xm*ym*(i-zs) + xm*(j-ys) - xs);
119347c6ae99SBarry Smith       }
119447c6ae99SBarry Smith     }
119547c6ae99SBarry Smith     *iptr = (void*)ptr;
11968865f1eaSKarl Rupp     break;
11978865f1eaSKarl Rupp   }
119847c6ae99SBarry Smith   default:
1199ce94432eSBarry Smith     SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Dimension %D not supported",dd->dim);
120047c6ae99SBarry Smith   }
120147c6ae99SBarry Smith 
120247c6ae99SBarry Smith done:
120347c6ae99SBarry Smith   /* add arrays to the checked out list */
120447c6ae99SBarry Smith   if (ghosted) {
1205aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
120647c6ae99SBarry Smith       if (!dd->arrayghostedout[i]) {
120747c6ae99SBarry Smith         dd->arrayghostedout[i] = *iptr;
120847c6ae99SBarry Smith         dd->startghostedout[i] = iarray_start;
120947c6ae99SBarry Smith         break;
121047c6ae99SBarry Smith       }
121147c6ae99SBarry Smith     }
121247c6ae99SBarry Smith   } else {
1213aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
121447c6ae99SBarry Smith       if (!dd->arrayout[i]) {
121547c6ae99SBarry Smith         dd->arrayout[i] = *iptr;
121647c6ae99SBarry Smith         dd->startout[i] = iarray_start;
121747c6ae99SBarry Smith         break;
121847c6ae99SBarry Smith       }
121947c6ae99SBarry Smith     }
122047c6ae99SBarry Smith   }
122147c6ae99SBarry Smith   PetscFunctionReturn(0);
122247c6ae99SBarry Smith }
122347c6ae99SBarry Smith 
122447c6ae99SBarry Smith #undef __FUNCT__
1225aa219208SBarry Smith #define __FUNCT__ "DMDARestoreArray"
122647c6ae99SBarry Smith /*@C
1227aa219208SBarry Smith      DMDARestoreArray - Restores an array of derivative types for a DMDA
122847c6ae99SBarry Smith 
122947c6ae99SBarry Smith     Input Parameter:
123047c6ae99SBarry Smith +    da - information about my local patch
123147c6ae99SBarry Smith .    ghosted - do you want arrays for the ghosted or nonghosted patch
123247c6ae99SBarry Smith -    vptr - array data structured to be passed to ad_FormFunctionLocal()
123347c6ae99SBarry Smith 
123447c6ae99SBarry Smith      Level: advanced
123547c6ae99SBarry Smith 
1236bcaeba4dSBarry Smith .seealso: DMDAGetArray()
123747c6ae99SBarry Smith 
123847c6ae99SBarry Smith @*/
12397087cfbeSBarry Smith PetscErrorCode  DMDARestoreArray(DM da,PetscBool ghosted,void *vptr)
124047c6ae99SBarry Smith {
124147c6ae99SBarry Smith   PetscInt i;
124247c6ae99SBarry Smith   void     **iptr = (void**)vptr,*iarray_start = 0;
124347c6ae99SBarry Smith   DM_DA    *dd    = (DM_DA*)da->data;
124447c6ae99SBarry Smith 
124547c6ae99SBarry Smith   PetscFunctionBegin;
124647c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
124747c6ae99SBarry Smith   if (ghosted) {
1248aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
124947c6ae99SBarry Smith       if (dd->arrayghostedout[i] == *iptr) {
125047c6ae99SBarry Smith         iarray_start           = dd->startghostedout[i];
12510298fd71SBarry Smith         dd->arrayghostedout[i] = NULL;
12520298fd71SBarry Smith         dd->startghostedout[i] = NULL;
125347c6ae99SBarry Smith         break;
125447c6ae99SBarry Smith       }
125547c6ae99SBarry Smith     }
1256aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
125747c6ae99SBarry Smith       if (!dd->arrayghostedin[i]) {
125847c6ae99SBarry Smith         dd->arrayghostedin[i] = *iptr;
125947c6ae99SBarry Smith         dd->startghostedin[i] = iarray_start;
126047c6ae99SBarry Smith         break;
126147c6ae99SBarry Smith       }
126247c6ae99SBarry Smith     }
126347c6ae99SBarry Smith   } else {
1264aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
126547c6ae99SBarry Smith       if (dd->arrayout[i] == *iptr) {
126647c6ae99SBarry Smith         iarray_start    = dd->startout[i];
12670298fd71SBarry Smith         dd->arrayout[i] = NULL;
12680298fd71SBarry Smith         dd->startout[i] = NULL;
126947c6ae99SBarry Smith         break;
127047c6ae99SBarry Smith       }
127147c6ae99SBarry Smith     }
1272aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
127347c6ae99SBarry Smith       if (!dd->arrayin[i]) {
127447c6ae99SBarry Smith         dd->arrayin[i] = *iptr;
127547c6ae99SBarry Smith         dd->startin[i] = iarray_start;
127647c6ae99SBarry Smith         break;
127747c6ae99SBarry Smith       }
127847c6ae99SBarry Smith     }
127947c6ae99SBarry Smith   }
128047c6ae99SBarry Smith   PetscFunctionReturn(0);
128147c6ae99SBarry Smith }
128247c6ae99SBarry Smith 
1283