xref: /petsc/src/dm/impls/da/dalocal.c (revision 580bdb303e1ee3b1222b2042810b4c26340259c6)
147c6ae99SBarry Smith 
247c6ae99SBarry Smith /*
347c6ae99SBarry Smith   Code for manipulating distributed regular arrays in parallel.
447c6ae99SBarry Smith */
547c6ae99SBarry Smith 
6af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h>    /*I   "petscdmda.h"   I*/
7b2fff234SMatthew G. Knepley #include <petscbt.h>
80c312b8eSJed Brown #include <petscsf.h>
92764a2aaSMatthew G. Knepley #include <petscds.h>
109a800dd8SMatthew G. Knepley #include <petscfe.h>
1147c6ae99SBarry Smith 
1247c6ae99SBarry Smith /*
13e3c5b3baSBarry Smith    This allows the DMDA vectors to properly tell MATLAB their dimensions
1447c6ae99SBarry Smith */
1547c6ae99SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
16c6db04a5SJed Brown #include <engine.h>   /* MATLAB include file */
17c6db04a5SJed Brown #include <mex.h>      /* MATLAB include file */
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
38*580bdb30SBarry Smith   ierr = PetscArraycpy(mxGetPr(mat),array,n*m);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 
477087cfbeSBarry Smith PetscErrorCode  DMCreateLocalVector_DA(DM da,Vec *g)
4847c6ae99SBarry Smith {
4947c6ae99SBarry Smith   PetscErrorCode ierr;
5047c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
5147c6ae99SBarry Smith 
5247c6ae99SBarry Smith   PetscFunctionBegin;
5347c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
5447c6ae99SBarry Smith   PetscValidPointer(g,2);
5511689aebSJed Brown   if (da->defaultSection) {
5611689aebSJed Brown     ierr = DMCreateLocalVector_Section_Private(da,g);CHKERRQ(ierr);
5711689aebSJed Brown   } else {
5847c6ae99SBarry Smith     ierr = VecCreate(PETSC_COMM_SELF,g);CHKERRQ(ierr);
5947c6ae99SBarry Smith     ierr = VecSetSizes(*g,dd->nlocal,PETSC_DETERMINE);CHKERRQ(ierr);
6047c6ae99SBarry Smith     ierr = VecSetBlockSize(*g,dd->w);CHKERRQ(ierr);
61401ddaa8SBarry Smith     ierr = VecSetType(*g,da->vectype);CHKERRQ(ierr);
62c688c046SMatthew G Knepley     ierr = VecSetDM(*g, da);CHKERRQ(ierr);
6347c6ae99SBarry Smith #if defined(PETSC_HAVE_MATLAB_ENGINE)
647bc9226cSBarry Smith     if (dd->w == 1  && da->dim == 2) {
65bdf89e91SBarry Smith       ierr = PetscObjectComposeFunction((PetscObject)*g,"PetscMatlabEnginePut_C",VecMatlabEnginePut_DA2d);CHKERRQ(ierr);
6647c6ae99SBarry Smith     }
6747c6ae99SBarry Smith #endif
6811689aebSJed Brown   }
6947c6ae99SBarry Smith   PetscFunctionReturn(0);
7047c6ae99SBarry Smith }
7147c6ae99SBarry Smith 
72939f6067SMatthew G. Knepley /*@
73939f6067SMatthew G. Knepley   DMDAGetNumCells - Get the number of cells in the local piece of the DMDA. This includes ghost cells.
74939f6067SMatthew G. Knepley 
75939f6067SMatthew G. Knepley   Input Parameter:
76939f6067SMatthew G. Knepley . dm - The DM object
77939f6067SMatthew G. Knepley 
78939f6067SMatthew G. Knepley   Output Parameters:
79939f6067SMatthew G. Knepley + numCellsX - The number of local cells in the x-direction
80939f6067SMatthew G. Knepley . numCellsY - The number of local cells in the y-direction
81939f6067SMatthew G. Knepley . numCellsZ - The number of local cells in the z-direction
82939f6067SMatthew G. Knepley - numCells - The number of local cells
83939f6067SMatthew G. Knepley 
84939f6067SMatthew G. Knepley   Level: developer
85939f6067SMatthew G. Knepley 
86939f6067SMatthew G. Knepley .seealso: DMDAGetCellPoint()
87939f6067SMatthew G. Knepley @*/
883582350cSMatthew G. Knepley PetscErrorCode DMDAGetNumCells(DM dm, PetscInt *numCellsX, PetscInt *numCellsY, PetscInt *numCellsZ, PetscInt *numCells)
8957459e9aSMatthew G Knepley {
9057459e9aSMatthew G Knepley   DM_DA         *da  = (DM_DA*) dm->data;
91c73cfb54SMatthew G. Knepley   const PetscInt dim = dm->dim;
9257459e9aSMatthew G Knepley   const PetscInt mx  = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs;
9357459e9aSMatthew G Knepley   const PetscInt nC  = (mx)*(dim > 1 ? (my)*(dim > 2 ? (mz) : 1) : 1);
9457459e9aSMatthew G Knepley 
9557459e9aSMatthew G Knepley   PetscFunctionBegin;
96a9a02de4SBarry Smith   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1,DMDA);
973582350cSMatthew G. Knepley   if (numCellsX) {
983582350cSMatthew G. Knepley     PetscValidIntPointer(numCellsX,2);
993582350cSMatthew G. Knepley     *numCellsX = mx;
1003582350cSMatthew G. Knepley   }
1013582350cSMatthew G. Knepley   if (numCellsY) {
1028135c375SStefano Zampini     PetscValidIntPointer(numCellsY,3);
1033582350cSMatthew G. Knepley     *numCellsY = my;
1043582350cSMatthew G. Knepley   }
1053582350cSMatthew G. Knepley   if (numCellsZ) {
1068135c375SStefano Zampini     PetscValidIntPointer(numCellsZ,4);
1073582350cSMatthew G. Knepley     *numCellsZ = mz;
1083582350cSMatthew G. Knepley   }
10957459e9aSMatthew G Knepley   if (numCells) {
1103582350cSMatthew G. Knepley     PetscValidIntPointer(numCells,5);
11157459e9aSMatthew G Knepley     *numCells = nC;
11257459e9aSMatthew G Knepley   }
11357459e9aSMatthew G Knepley   PetscFunctionReturn(0);
11457459e9aSMatthew G Knepley }
11557459e9aSMatthew G Knepley 
116d6401b93SMatthew G. Knepley /*@
117d6401b93SMatthew G. Knepley   DMDAGetCellPoint - Get the DM point corresponding to the tuple (i, j, k) in the DMDA
118d6401b93SMatthew G. Knepley 
119d6401b93SMatthew G. Knepley   Input Parameters:
120d6401b93SMatthew G. Knepley + dm - The DM object
121d6401b93SMatthew G. Knepley - i,j,k - The global indices for the cell
122d6401b93SMatthew G. Knepley 
123d6401b93SMatthew G. Knepley   Output Parameters:
124d6401b93SMatthew G. Knepley . point - The local DM point
125d6401b93SMatthew G. Knepley 
126d6401b93SMatthew G. Knepley   Level: developer
127d6401b93SMatthew G. Knepley 
128d6401b93SMatthew G. Knepley .seealso: DMDAGetNumCells()
129d6401b93SMatthew G. Knepley @*/
130d6401b93SMatthew G. Knepley PetscErrorCode DMDAGetCellPoint(DM dm, PetscInt i, PetscInt j, PetscInt k, PetscInt *point)
131d6401b93SMatthew G. Knepley {
132c73cfb54SMatthew G. Knepley   const PetscInt dim = dm->dim;
133d6401b93SMatthew G. Knepley   DMDALocalInfo  info;
134d6401b93SMatthew G. Knepley   PetscErrorCode ierr;
135d6401b93SMatthew G. Knepley 
136d6401b93SMatthew G. Knepley   PetscFunctionBegin;
137a9a02de4SBarry Smith   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1,DMDA);
138d6401b93SMatthew G. Knepley   PetscValidIntPointer(point,5);
139d6401b93SMatthew G. Knepley   ierr = DMDAGetLocalInfo(dm, &info);CHKERRQ(ierr);
140d6401b93SMatthew G. Knepley   if (dim > 0) {if ((i < info.gxs) || (i >= info.gxs+info.gxm)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "X index %d not in [%d, %d)", i, info.gxs, info.gxs+info.gxm);}
14163b6bffdSJed Brown   if (dim > 1) {if ((j < info.gys) || (j >= info.gys+info.gym)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Y index %d not in [%d, %d)", j, info.gys, info.gys+info.gym);}
14263b6bffdSJed Brown   if (dim > 2) {if ((k < info.gzs) || (k >= info.gzs+info.gzm)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Z index %d not in [%d, %d)", k, info.gzs, info.gzs+info.gzm);}
143d6401b93SMatthew G. Knepley   *point = i + (dim > 1 ? (j + (dim > 2 ? k*info.gym : 0))*info.gxm : 0);
144d6401b93SMatthew G. Knepley   PetscFunctionReturn(0);
145d6401b93SMatthew G. Knepley }
146d6401b93SMatthew G. Knepley 
14757459e9aSMatthew G Knepley PetscErrorCode DMDAGetNumVertices(DM dm, PetscInt *numVerticesX, PetscInt *numVerticesY, PetscInt *numVerticesZ, PetscInt *numVertices)
14857459e9aSMatthew G Knepley {
14957459e9aSMatthew G Knepley   DM_DA          *da = (DM_DA*) dm->data;
150c73cfb54SMatthew G. Knepley   const PetscInt dim = dm->dim;
15157459e9aSMatthew G Knepley   const PetscInt mx  = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs;
15257459e9aSMatthew G Knepley   const PetscInt nVx = mx+1;
15357459e9aSMatthew G Knepley   const PetscInt nVy = dim > 1 ? (my+1) : 1;
15457459e9aSMatthew G Knepley   const PetscInt nVz = dim > 2 ? (mz+1) : 1;
15557459e9aSMatthew G Knepley   const PetscInt nV  = nVx*nVy*nVz;
15657459e9aSMatthew G Knepley 
15757459e9aSMatthew G Knepley   PetscFunctionBegin;
15857459e9aSMatthew G Knepley   if (numVerticesX) {
15957459e9aSMatthew G Knepley     PetscValidIntPointer(numVerticesX,2);
16057459e9aSMatthew G Knepley     *numVerticesX = nVx;
16157459e9aSMatthew G Knepley   }
16257459e9aSMatthew G Knepley   if (numVerticesY) {
16357459e9aSMatthew G Knepley     PetscValidIntPointer(numVerticesY,3);
16457459e9aSMatthew G Knepley     *numVerticesY = nVy;
16557459e9aSMatthew G Knepley   }
16657459e9aSMatthew G Knepley   if (numVerticesZ) {
16757459e9aSMatthew G Knepley     PetscValidIntPointer(numVerticesZ,4);
16857459e9aSMatthew G Knepley     *numVerticesZ = nVz;
16957459e9aSMatthew G Knepley   }
17057459e9aSMatthew G Knepley   if (numVertices) {
17157459e9aSMatthew G Knepley     PetscValidIntPointer(numVertices,5);
17257459e9aSMatthew G Knepley     *numVertices = nV;
17357459e9aSMatthew G Knepley   }
17457459e9aSMatthew G Knepley   PetscFunctionReturn(0);
17557459e9aSMatthew G Knepley }
17657459e9aSMatthew G Knepley 
17757459e9aSMatthew G Knepley PetscErrorCode DMDAGetNumFaces(DM dm, PetscInt *numXFacesX, PetscInt *numXFaces, PetscInt *numYFacesY, PetscInt *numYFaces, PetscInt *numZFacesZ, PetscInt *numZFaces)
17857459e9aSMatthew G Knepley {
17957459e9aSMatthew G Knepley   DM_DA          *da = (DM_DA*) dm->data;
180c73cfb54SMatthew G. Knepley   const PetscInt dim = dm->dim;
18157459e9aSMatthew G Knepley   const PetscInt mx  = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys, mz = da->Ze - da->Zs;
18257459e9aSMatthew G Knepley   const PetscInt nxF = (dim > 1 ? (my)*(dim > 2 ? (mz) : 1) : 1);
18357459e9aSMatthew G Knepley   const PetscInt nXF = (mx+1)*nxF;
18457459e9aSMatthew G Knepley   const PetscInt nyF = mx*(dim > 2 ? mz : 1);
18557459e9aSMatthew G Knepley   const PetscInt nYF = dim > 1 ? (my+1)*nyF : 0;
18657459e9aSMatthew G Knepley   const PetscInt nzF = mx*(dim > 1 ? my : 0);
18757459e9aSMatthew G Knepley   const PetscInt nZF = dim > 2 ? (mz+1)*nzF : 0;
18857459e9aSMatthew G Knepley 
18957459e9aSMatthew G Knepley   PetscFunctionBegin;
19057459e9aSMatthew G Knepley   if (numXFacesX) {
19157459e9aSMatthew G Knepley     PetscValidIntPointer(numXFacesX,2);
19257459e9aSMatthew G Knepley     *numXFacesX = nxF;
19357459e9aSMatthew G Knepley   }
19457459e9aSMatthew G Knepley   if (numXFaces) {
19557459e9aSMatthew G Knepley     PetscValidIntPointer(numXFaces,3);
19657459e9aSMatthew G Knepley     *numXFaces = nXF;
19757459e9aSMatthew G Knepley   }
19857459e9aSMatthew G Knepley   if (numYFacesY) {
19957459e9aSMatthew G Knepley     PetscValidIntPointer(numYFacesY,4);
20057459e9aSMatthew G Knepley     *numYFacesY = nyF;
20157459e9aSMatthew G Knepley   }
20257459e9aSMatthew G Knepley   if (numYFaces) {
20357459e9aSMatthew G Knepley     PetscValidIntPointer(numYFaces,5);
20457459e9aSMatthew G Knepley     *numYFaces = nYF;
20557459e9aSMatthew G Knepley   }
20657459e9aSMatthew G Knepley   if (numZFacesZ) {
20757459e9aSMatthew G Knepley     PetscValidIntPointer(numZFacesZ,6);
20857459e9aSMatthew G Knepley     *numZFacesZ = nzF;
20957459e9aSMatthew G Knepley   }
21057459e9aSMatthew G Knepley   if (numZFaces) {
21157459e9aSMatthew G Knepley     PetscValidIntPointer(numZFaces,7);
21257459e9aSMatthew G Knepley     *numZFaces = nZF;
21357459e9aSMatthew G Knepley   }
21457459e9aSMatthew G Knepley   PetscFunctionReturn(0);
21557459e9aSMatthew G Knepley }
21657459e9aSMatthew G Knepley 
21757459e9aSMatthew G Knepley PetscErrorCode DMDAGetHeightStratum(DM dm, PetscInt height, PetscInt *pStart, PetscInt *pEnd)
21857459e9aSMatthew G Knepley {
219c73cfb54SMatthew G. Knepley   const PetscInt dim = dm->dim;
22057459e9aSMatthew G Knepley   PetscInt       nC, nV, nXF, nYF, nZF;
22157459e9aSMatthew G Knepley   PetscErrorCode ierr;
22257459e9aSMatthew G Knepley 
22357459e9aSMatthew G Knepley   PetscFunctionBegin;
2248865f1eaSKarl Rupp   if (pStart) PetscValidIntPointer(pStart,3);
2258865f1eaSKarl Rupp   if (pEnd)   PetscValidIntPointer(pEnd,4);
2263582350cSMatthew G. Knepley   ierr = DMDAGetNumCells(dm, NULL, NULL, NULL, &nC);CHKERRQ(ierr);
2270298fd71SBarry Smith   ierr = DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV);CHKERRQ(ierr);
2280298fd71SBarry Smith   ierr = DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF);CHKERRQ(ierr);
22957459e9aSMatthew G Knepley   if (height == 0) {
23057459e9aSMatthew G Knepley     /* Cells */
2318865f1eaSKarl Rupp     if (pStart) *pStart = 0;
2328865f1eaSKarl Rupp     if (pEnd)   *pEnd   = nC;
23357459e9aSMatthew G Knepley   } else if (height == 1) {
23457459e9aSMatthew G Knepley     /* Faces */
2358865f1eaSKarl Rupp     if (pStart) *pStart = nC+nV;
2368865f1eaSKarl Rupp     if (pEnd)   *pEnd   = nC+nV+nXF+nYF+nZF;
23757459e9aSMatthew G Knepley   } else if (height == dim) {
23857459e9aSMatthew G Knepley     /* Vertices */
2398865f1eaSKarl Rupp     if (pStart) *pStart = nC;
2408865f1eaSKarl Rupp     if (pEnd)   *pEnd   = nC+nV;
24157459e9aSMatthew G Knepley   } else if (height < 0) {
24257459e9aSMatthew G Knepley     /* All points */
2438865f1eaSKarl Rupp     if (pStart) *pStart = 0;
2448865f1eaSKarl Rupp     if (pEnd)   *pEnd   = nC+nV+nXF+nYF+nZF;
24582f516ccSBarry Smith   } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No points of height %d in the DA", height);
24657459e9aSMatthew G Knepley   PetscFunctionReturn(0);
24757459e9aSMatthew G Knepley }
24857459e9aSMatthew G Knepley 
2493385ff7eSMatthew G. Knepley PetscErrorCode DMDAGetDepthStratum(DM dm, PetscInt depth, PetscInt *pStart, PetscInt *pEnd)
2503385ff7eSMatthew G. Knepley {
251c73cfb54SMatthew G. Knepley   const PetscInt dim = dm->dim;
2523385ff7eSMatthew G. Knepley   PetscInt       nC, nV, nXF, nYF, nZF;
2533385ff7eSMatthew G. Knepley   PetscErrorCode ierr;
2543385ff7eSMatthew G. Knepley 
2553385ff7eSMatthew G. Knepley   PetscFunctionBegin;
2563385ff7eSMatthew G. Knepley   if (pStart) PetscValidIntPointer(pStart,3);
2573385ff7eSMatthew G. Knepley   if (pEnd)   PetscValidIntPointer(pEnd,4);
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   if (depth == dim) {
2623385ff7eSMatthew G. Knepley     /* Cells */
2633385ff7eSMatthew G. Knepley     if (pStart) *pStart = 0;
2643385ff7eSMatthew G. Knepley     if (pEnd)   *pEnd   = nC;
2653385ff7eSMatthew G. Knepley   } else if (depth == dim-1) {
2663385ff7eSMatthew G. Knepley     /* Faces */
2673385ff7eSMatthew G. Knepley     if (pStart) *pStart = nC+nV;
2683385ff7eSMatthew G. Knepley     if (pEnd)   *pEnd   = nC+nV+nXF+nYF+nZF;
2693385ff7eSMatthew G. Knepley   } else if (depth == 0) {
2703385ff7eSMatthew G. Knepley     /* Vertices */
2713385ff7eSMatthew G. Knepley     if (pStart) *pStart = nC;
2723385ff7eSMatthew G. Knepley     if (pEnd)   *pEnd   = nC+nV;
2733385ff7eSMatthew G. Knepley   } else if (depth < 0) {
2743385ff7eSMatthew G. Knepley     /* All points */
2753385ff7eSMatthew G. Knepley     if (pStart) *pStart = 0;
2763385ff7eSMatthew G. Knepley     if (pEnd)   *pEnd   = nC+nV+nXF+nYF+nZF;
2773385ff7eSMatthew G. Knepley   } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No points of depth %d in the DA", depth);
2783385ff7eSMatthew G. Knepley   PetscFunctionReturn(0);
2793385ff7eSMatthew G. Knepley }
2803385ff7eSMatthew G. Knepley 
2813385ff7eSMatthew G. Knepley PetscErrorCode DMDAGetConeSize(DM dm, PetscInt p, PetscInt *coneSize)
2823385ff7eSMatthew G. Knepley {
283c73cfb54SMatthew G. Knepley   const PetscInt dim = dm->dim;
2843385ff7eSMatthew G. Knepley   PetscInt       nC, nV, nXF, nYF, nZF;
2853385ff7eSMatthew G. Knepley   PetscErrorCode ierr;
2863385ff7eSMatthew G. Knepley 
2873385ff7eSMatthew G. Knepley   PetscFunctionBegin;
2883385ff7eSMatthew G. Knepley   *coneSize = 0;
2893385ff7eSMatthew G. Knepley   ierr = DMDAGetNumCells(dm, NULL, NULL, NULL, &nC);CHKERRQ(ierr);
2903385ff7eSMatthew G. Knepley   ierr = DMDAGetNumVertices(dm, NULL, NULL, NULL, &nV);CHKERRQ(ierr);
2913385ff7eSMatthew G. Knepley   ierr = DMDAGetNumFaces(dm, NULL, &nXF, NULL, &nYF, NULL, &nZF);CHKERRQ(ierr);
2923385ff7eSMatthew G. Knepley   switch (dim) {
2933385ff7eSMatthew G. Knepley   case 2:
2943385ff7eSMatthew G. Knepley     if (p >= 0) {
2953385ff7eSMatthew G. Knepley       if (p < nC) {
2963385ff7eSMatthew G. Knepley         *coneSize = 4;
2973385ff7eSMatthew G. Knepley       } else if (p < nC+nV) {
2983385ff7eSMatthew G. Knepley         *coneSize = 0;
2993385ff7eSMatthew G. Knepley       } else if (p < nC+nV+nXF+nYF+nZF) {
3003385ff7eSMatthew G. Knepley         *coneSize = 2;
3013385ff7eSMatthew G. Knepley       } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %d should be in [0, %d)", p, nC+nV+nXF+nYF+nZF);
3023385ff7eSMatthew G. Knepley     } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative point %d is invalid", p);
3033385ff7eSMatthew G. Knepley     break;
3043385ff7eSMatthew G. Knepley   case 3:
3053385ff7eSMatthew G. Knepley     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do 3D");
3063385ff7eSMatthew G. Knepley     break;
3073385ff7eSMatthew G. Knepley   }
3083385ff7eSMatthew G. Knepley   PetscFunctionReturn(0);
3093385ff7eSMatthew G. Knepley }
3103385ff7eSMatthew G. Knepley 
3113385ff7eSMatthew G. Knepley PetscErrorCode DMDAGetCone(DM dm, PetscInt p, PetscInt *cone[])
3123385ff7eSMatthew G. Knepley {
313c73cfb54SMatthew G. Knepley   const PetscInt dim = dm->dim;
3143385ff7eSMatthew G. Knepley   PetscInt       nCx, nCy, nCz, nC, nVx, nVy, nVz, nV, nxF, nyF, nzF, nXF, nYF, nZF;
3153385ff7eSMatthew G. Knepley   PetscErrorCode ierr;
3163385ff7eSMatthew G. Knepley 
3173385ff7eSMatthew G. Knepley   PetscFunctionBegin;
31869291d52SBarry Smith   if (!cone) {ierr = DMGetWorkArray(dm, 6, MPIU_INT, cone);CHKERRQ(ierr);}
3193385ff7eSMatthew G. Knepley   ierr = DMDAGetNumCells(dm, &nCx, &nCy, &nCz, &nC);CHKERRQ(ierr);
3203385ff7eSMatthew G. Knepley   ierr = DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV);CHKERRQ(ierr);
3213385ff7eSMatthew G. Knepley   ierr = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr);
3223385ff7eSMatthew G. Knepley   switch (dim) {
3233385ff7eSMatthew G. Knepley   case 2:
3243385ff7eSMatthew G. Knepley     if (p >= 0) {
3253385ff7eSMatthew G. Knepley       if (p < nC) {
3263385ff7eSMatthew G. Knepley         const PetscInt cy = p / nCx;
3273385ff7eSMatthew G. Knepley         const PetscInt cx = p % nCx;
3283385ff7eSMatthew G. Knepley 
3293385ff7eSMatthew G. Knepley         (*cone)[0] = cy*nxF + cx + nC+nV;
3303385ff7eSMatthew G. Knepley         (*cone)[1] = cx*nyF + cy + nyF + nC+nV+nXF;
3313385ff7eSMatthew G. Knepley         (*cone)[2] = cy*nxF + cx + nxF + nC+nV;
3323385ff7eSMatthew G. Knepley         (*cone)[3] = cx*nyF + cy + nC+nV+nXF;
3333385ff7eSMatthew G. Knepley         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do cell cones");
3343385ff7eSMatthew G. Knepley       } else if (p < nC+nV) {
3353385ff7eSMatthew G. Knepley       } else if (p < nC+nV+nXF) {
3363385ff7eSMatthew G. Knepley         const PetscInt fy = (p - nC+nV) / nxF;
3373385ff7eSMatthew G. Knepley         const PetscInt fx = (p - nC+nV) % nxF;
3383385ff7eSMatthew G. Knepley 
3393385ff7eSMatthew G. Knepley         (*cone)[0] = fy*nVx + fx + nC;
3403385ff7eSMatthew G. Knepley         (*cone)[1] = fy*nVx + fx + 1 + nC;
3413385ff7eSMatthew G. Knepley       } else if (p < nC+nV+nXF+nYF) {
3423385ff7eSMatthew G. Knepley         const PetscInt fx = (p - nC+nV+nXF) / nyF;
3433385ff7eSMatthew G. Knepley         const PetscInt fy = (p - nC+nV+nXF) % nyF;
3443385ff7eSMatthew G. Knepley 
3453385ff7eSMatthew G. Knepley         (*cone)[0] = fy*nVx + fx + nC;
3463385ff7eSMatthew G. Knepley         (*cone)[1] = fy*nVx + fx + nVx + nC;
3473385ff7eSMatthew G. Knepley       } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %d should be in [0, %d)", p, nC+nV+nXF+nYF+nZF);
3483385ff7eSMatthew G. Knepley     } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Negative point %d is invalid", p);
3493385ff7eSMatthew G. Knepley     break;
3503385ff7eSMatthew G. Knepley   case 3:
3513385ff7eSMatthew G. Knepley     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to do 3D");
3523385ff7eSMatthew G. Knepley     break;
3533385ff7eSMatthew G. Knepley   }
3543385ff7eSMatthew G. Knepley   PetscFunctionReturn(0);
3553385ff7eSMatthew G. Knepley }
3563385ff7eSMatthew G. Knepley 
3573385ff7eSMatthew G. Knepley PetscErrorCode DMDARestoreCone(DM dm, PetscInt p, PetscInt *cone[])
3583385ff7eSMatthew G. Knepley {
3593385ff7eSMatthew G. Knepley   PetscErrorCode ierr;
3603385ff7eSMatthew G. Knepley 
3613385ff7eSMatthew G. Knepley   PetscFunctionBegin;
36269291d52SBarry Smith   ierr = DMGetWorkArray(dm, 6, MPIU_INT, cone);CHKERRQ(ierr);
3633385ff7eSMatthew G. Knepley   PetscFunctionReturn(0);
3643385ff7eSMatthew G. Knepley }
3653385ff7eSMatthew G. Knepley 
3663385ff7eSMatthew G. Knepley PetscErrorCode DMDASetVertexCoordinates(DM dm, PetscReal xl, PetscReal xu, PetscReal yl, PetscReal yu, PetscReal zl, PetscReal zu)
3673385ff7eSMatthew G. Knepley {
3683385ff7eSMatthew G. Knepley   DM_DA         *da = (DM_DA *) dm->data;
3693385ff7eSMatthew G. Knepley   Vec            coordinates;
3703385ff7eSMatthew G. Knepley   PetscSection   section;
3713385ff7eSMatthew G. Knepley   PetscScalar   *coords;
3723385ff7eSMatthew G. Knepley   PetscReal      h[3];
3733385ff7eSMatthew G. Knepley   PetscInt       dim, size, M, N, P, nVx, nVy, nVz, nV, vStart, vEnd, v, i, j, k;
3743385ff7eSMatthew G. Knepley   PetscErrorCode ierr;
3753385ff7eSMatthew G. Knepley 
3763385ff7eSMatthew G. Knepley   PetscFunctionBegin;
377a9a02de4SBarry Smith   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1,DMDA);
3783385ff7eSMatthew G. Knepley   ierr = DMDAGetInfo(dm, &dim, &M, &N, &P, 0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
3794a2f8832SBarry Smith   if (dim > 3) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"The following code only works for dim <= 3");
3803385ff7eSMatthew G. Knepley   h[0] = (xu - xl)/M;
3813385ff7eSMatthew G. Knepley   h[1] = (yu - yl)/N;
3823385ff7eSMatthew G. Knepley   h[2] = (zu - zl)/P;
3833385ff7eSMatthew G. Knepley   ierr = DMDAGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
3843385ff7eSMatthew G. Knepley   ierr = DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, &nV);CHKERRQ(ierr);
3853385ff7eSMatthew G. Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &section);CHKERRQ(ierr);
3863385ff7eSMatthew G. Knepley   ierr = PetscSectionSetNumFields(section, 1);CHKERRQ(ierr);
3873385ff7eSMatthew G. Knepley   ierr = PetscSectionSetFieldComponents(section, 0, dim);CHKERRQ(ierr);
3883385ff7eSMatthew G. Knepley   ierr = PetscSectionSetChart(section, vStart, vEnd);CHKERRQ(ierr);
3893385ff7eSMatthew G. Knepley   for (v = vStart; v < vEnd; ++v) {
3903385ff7eSMatthew G. Knepley     ierr = PetscSectionSetDof(section, v, dim);CHKERRQ(ierr);
3913385ff7eSMatthew G. Knepley   }
3923385ff7eSMatthew G. Knepley   ierr = PetscSectionSetUp(section);CHKERRQ(ierr);
3933385ff7eSMatthew G. Knepley   ierr = PetscSectionGetStorageSize(section, &size);CHKERRQ(ierr);
3943385ff7eSMatthew G. Knepley   ierr = VecCreateSeq(PETSC_COMM_SELF, size, &coordinates);CHKERRQ(ierr);
39524640c55SToby Isaac   ierr = PetscObjectSetName((PetscObject)coordinates,"coordinates");CHKERRQ(ierr);
3963385ff7eSMatthew G. Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
3973385ff7eSMatthew G. Knepley   for (k = 0; k < nVz; ++k) {
398e4d69e08SBarry Smith     PetscInt ind[3], d, off;
3993385ff7eSMatthew G. Knepley 
400e4d69e08SBarry Smith     ind[0] = 0;
401e4d69e08SBarry Smith     ind[1] = 0;
402e4d69e08SBarry Smith     ind[2] = k + da->zs;
4033385ff7eSMatthew G. Knepley     for (j = 0; j < nVy; ++j) {
4043385ff7eSMatthew G. Knepley       ind[1] = j + da->ys;
4053385ff7eSMatthew G. Knepley       for (i = 0; i < nVx; ++i) {
4063385ff7eSMatthew G. Knepley         const PetscInt vertex = (k*nVy + j)*nVx + i + vStart;
4073385ff7eSMatthew G. Knepley 
4083385ff7eSMatthew G. Knepley         ierr = PetscSectionGetOffset(section, vertex, &off);CHKERRQ(ierr);
4093385ff7eSMatthew G. Knepley         ind[0] = i + da->xs;
4103385ff7eSMatthew G. Knepley         for (d = 0; d < dim; ++d) {
4113385ff7eSMatthew G. Knepley           coords[off+d] = h[d]*ind[d];
4123385ff7eSMatthew G. Knepley         }
4133385ff7eSMatthew G. Knepley       }
4143385ff7eSMatthew G. Knepley     }
4153385ff7eSMatthew G. Knepley   }
4163385ff7eSMatthew G. Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
41746e270d4SMatthew G. Knepley   ierr = DMSetCoordinateSection(dm, PETSC_DETERMINE, section);CHKERRQ(ierr);
4183385ff7eSMatthew G. Knepley   ierr = DMSetCoordinatesLocal(dm, coordinates);CHKERRQ(ierr);
419a4b60ecfSMatthew G. Knepley   ierr = PetscSectionDestroy(&section);CHKERRQ(ierr);
4203385ff7eSMatthew G. Knepley   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
4213385ff7eSMatthew G. Knepley   PetscFunctionReturn(0);
4223385ff7eSMatthew G. Knepley }
4239a800dd8SMatthew G. Knepley 
424c970b36aSToby Isaac PetscErrorCode DMProjectFunctionLocal_DA(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
4259a800dd8SMatthew G. Knepley {
4262764a2aaSMatthew G. Knepley   PetscDS         prob;
427b7a4d31fSMatthew G. Knepley   PetscFE         fe;
428b7a4d31fSMatthew G. Knepley   PetscDualSpace  sp;
4299a800dd8SMatthew G. Knepley   PetscQuadrature q;
4309a800dd8SMatthew G. Knepley   PetscSection    section;
4319a800dd8SMatthew G. Knepley   PetscScalar    *values;
432c330f8ffSToby Isaac   PetscInt        numFields, Nc, dim, dimEmbed, spDim, totDim, numValues, cStart, cEnd, f, c, v, d;
433c330f8ffSToby Isaac   PetscFEGeom    *geom;
4349a800dd8SMatthew G. Knepley   PetscErrorCode  ierr;
4359a800dd8SMatthew G. Knepley 
4369a800dd8SMatthew G. Knepley   PetscFunctionBegin;
4372764a2aaSMatthew G. Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
4382764a2aaSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
4392764a2aaSMatthew G. Knepley   ierr = PetscDSGetDiscretization(prob, 0, (PetscObject *) &fe);CHKERRQ(ierr);
440b7a4d31fSMatthew G. Knepley   ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
441e87a4003SBarry Smith   ierr = DMGetSection(dm, &section);CHKERRQ(ierr);
4429a800dd8SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
4439a800dd8SMatthew G. Knepley   ierr = DMDAGetInfo(dm, &dim,0,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
4447a1a1bd4SToby Isaac   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
4459a800dd8SMatthew G. Knepley   ierr = DMDAGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
4469a800dd8SMatthew G. Knepley   ierr = DMDAVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr);
4479a800dd8SMatthew 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);
44869291d52SBarry Smith   ierr = DMGetWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr);
449c330f8ffSToby Isaac   ierr = PetscFEGeomCreate(q,1,dimEmbed,PETSC_FALSE,&geom);CHKERRQ(ierr);
4509a800dd8SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
451c330f8ffSToby Isaac     ierr          = DMDAComputeCellGeometryFEM(dm, c, q, geom->v, geom->J, NULL, geom->detJ);CHKERRQ(ierr);
4529a800dd8SMatthew G. Knepley     for (f = 0, v = 0; f < numFields; ++f) {
453c110b1eeSGeoffrey Irving       void * const ctx = ctxs ? ctxs[f] : NULL;
454b7a4d31fSMatthew G. Knepley 
4552764a2aaSMatthew G. Knepley       ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
4569c3cf19fSMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
457b7a4d31fSMatthew G. Knepley       ierr = PetscFEGetDualSpace(fe, &sp);CHKERRQ(ierr);
458b7a4d31fSMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp, &spDim);CHKERRQ(ierr);
4599a800dd8SMatthew G. Knepley       for (d = 0; d < spDim; ++d) {
460c330f8ffSToby Isaac         ierr = PetscDualSpaceApply(sp, d, time, geom, Nc, funcs[f], ctx, &values[v]);CHKERRQ(ierr);
4619a800dd8SMatthew G. Knepley       }
4629a800dd8SMatthew G. Knepley     }
4639a800dd8SMatthew G. Knepley     ierr = DMDAVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr);
4649a800dd8SMatthew G. Knepley   }
465c330f8ffSToby Isaac   ierr = PetscFEGeomDestroy(&geom);CHKERRQ(ierr);
46669291d52SBarry Smith   ierr = DMRestoreWorkArray(dm, numValues, MPIU_SCALAR, &values);CHKERRQ(ierr);
4679a800dd8SMatthew G. Knepley   PetscFunctionReturn(0);
4689a800dd8SMatthew G. Knepley }
4699a800dd8SMatthew G. Knepley 
470c970b36aSToby Isaac PetscErrorCode DMComputeL2Diff_DA(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)
4719a800dd8SMatthew G. Knepley {
4729a800dd8SMatthew G. Knepley   const PetscInt  debug = 0;
4732764a2aaSMatthew G. Knepley   PetscDS         prob;
474b7a4d31fSMatthew G. Knepley   PetscFE         fe;
4759a800dd8SMatthew G. Knepley   PetscQuadrature quad;
476b7a4d31fSMatthew G. Knepley   PetscSection    section;
4779a800dd8SMatthew G. Knepley   Vec             localX;
4789a800dd8SMatthew G. Knepley   PetscScalar    *funcVal;
4799a800dd8SMatthew G. Knepley   PetscReal      *coords, *v0, *J, *invJ, detJ;
4809a800dd8SMatthew G. Knepley   PetscReal       localDiff = 0.0;
481b7a4d31fSMatthew G. Knepley   PetscInt        dim, numFields, Nc, cStart, cEnd, c, field, fieldOffset, comp;
4829a800dd8SMatthew G. Knepley   PetscErrorCode  ierr;
4839a800dd8SMatthew G. Knepley 
4849a800dd8SMatthew G. Knepley   PetscFunctionBegin;
4859a800dd8SMatthew G. Knepley   ierr = DMDAGetInfo(dm, &dim,0,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
4862764a2aaSMatthew G. Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
4872764a2aaSMatthew G. Knepley   ierr = PetscDSGetTotalComponents(prob, &Nc);CHKERRQ(ierr);
4882764a2aaSMatthew G. Knepley   ierr = PetscDSGetDiscretization(prob, 0, (PetscObject *) &fe);CHKERRQ(ierr);
489b7a4d31fSMatthew G. Knepley   ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
490e87a4003SBarry Smith   ierr = DMGetSection(dm, &section);CHKERRQ(ierr);
4919a800dd8SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
4929a800dd8SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
4939a800dd8SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
4949a800dd8SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
4959a800dd8SMatthew G. Knepley   /* There are no BC values in DAs right now: ierr = DMDAProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */
496b7a4d31fSMatthew G. Knepley   ierr = PetscMalloc5(Nc,&funcVal,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
4979a800dd8SMatthew G. Knepley   ierr = DMDAGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
4989a800dd8SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
4999a800dd8SMatthew G. Knepley     PetscScalar *x = NULL;
5009a800dd8SMatthew G. Knepley     PetscReal    elemDiff = 0.0;
5019a800dd8SMatthew G. Knepley 
5028e0841e0SMatthew G. Knepley     ierr = DMDAComputeCellGeometryFEM(dm, c, quad, v0, J, invJ, &detJ);CHKERRQ(ierr);
5039a800dd8SMatthew G. Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
5049a800dd8SMatthew G. Knepley     ierr = DMDAVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
5059a800dd8SMatthew G. Knepley 
5069a800dd8SMatthew G. Knepley     for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) {
507c110b1eeSGeoffrey Irving       void * const ctx = ctxs ? ctxs[field] : NULL;
50821454ff5SMatthew G. Knepley       const PetscReal *quadPoints, *quadWeights;
5099a800dd8SMatthew G. Knepley       PetscReal       *basis;
5109c3cf19fSMatthew G. Knepley       PetscInt         qNc, Nq, Nb, Nc, q, d, e, fc, f;
5119a800dd8SMatthew G. Knepley 
5122764a2aaSMatthew G. Knepley       ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);CHKERRQ(ierr);
5139c3cf19fSMatthew G. Knepley       ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr);
5149c3cf19fSMatthew G. Knepley       ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
5159c3cf19fSMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
5169c3cf19fSMatthew G. Knepley       if (qNc != Nc) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Quadrature components %D != %D field components\n", qNc, Nc);
517b7a4d31fSMatthew G. Knepley       ierr = PetscFEGetDefaultTabulation(fe, &basis, NULL, NULL);CHKERRQ(ierr);
5189a800dd8SMatthew G. Knepley       if (debug) {
5199a800dd8SMatthew G. Knepley         char title[1024];
5209a800dd8SMatthew G. Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
5219c3cf19fSMatthew G. Knepley         ierr = DMPrintCellVector(c, title, Nb, &x[fieldOffset]);CHKERRQ(ierr);
5229a800dd8SMatthew G. Knepley       }
5239c3cf19fSMatthew G. Knepley       for (q = 0; q < Nq; ++q) {
5249a800dd8SMatthew G. Knepley         for (d = 0; d < dim; d++) {
5259a800dd8SMatthew G. Knepley           coords[d] = v0[d];
5269a800dd8SMatthew G. Knepley           for (e = 0; e < dim; e++) {
5279a800dd8SMatthew G. Knepley             coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0);
5289a800dd8SMatthew G. Knepley           }
5299a800dd8SMatthew G. Knepley         }
530a98e6df5SMatthew G. Knepley         ierr = (*funcs[field])(dim, time, coords, numFields, funcVal, ctx);CHKERRQ(ierr);
5319c3cf19fSMatthew G. Knepley         for (fc = 0; fc < Nc; ++fc) {
5329a800dd8SMatthew G. Knepley           PetscScalar interpolant = 0.0;
5339a800dd8SMatthew G. Knepley 
5349c3cf19fSMatthew G. Knepley           for (f = 0; f < Nb; ++f) {
5359c3cf19fSMatthew G. Knepley             interpolant += x[fieldOffset+f]*basis[(q*Nb+f)*Nc+fc];
5369a800dd8SMatthew G. Knepley           }
5379c3cf19fSMatthew 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*Nc+fc]*detJ);CHKERRQ(ierr);}
5389c3cf19fSMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q*Nc+c]*detJ;
5399a800dd8SMatthew G. Knepley         }
5409a800dd8SMatthew G. Knepley       }
5419c3cf19fSMatthew G. Knepley       comp        += Nc;
5429c3cf19fSMatthew G. Knepley       fieldOffset += Nb;
5439a800dd8SMatthew G. Knepley     }
5449a800dd8SMatthew G. Knepley     ierr = DMDAVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
5459a800dd8SMatthew G. Knepley     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
5469a800dd8SMatthew G. Knepley     localDiff += elemDiff;
5479a800dd8SMatthew G. Knepley   }
5489a800dd8SMatthew G. Knepley   ierr  = PetscFree5(funcVal,coords,v0,J,invJ);CHKERRQ(ierr);
5499a800dd8SMatthew G. Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
550b2566f29SBarry Smith   ierr  = MPIU_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
5519a800dd8SMatthew G. Knepley   *diff = PetscSqrtReal(*diff);
552a66d4d66SMatthew G Knepley   PetscFunctionReturn(0);
553a66d4d66SMatthew G Knepley }
554a66d4d66SMatthew G Knepley 
555b698f381SToby Isaac PetscErrorCode DMComputeL2GradientDiff_DA(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, const PetscReal n[], PetscReal *diff)
55623d86601SMatthew G. Knepley {
55723d86601SMatthew G. Knepley   const PetscInt  debug = 0;
5582764a2aaSMatthew G. Knepley   PetscDS         prob;
559b7a4d31fSMatthew G. Knepley   PetscFE         fe;
56023d86601SMatthew G. Knepley   PetscQuadrature quad;
561b7a4d31fSMatthew G. Knepley   PetscSection    section;
56223d86601SMatthew G. Knepley   Vec             localX;
56323d86601SMatthew G. Knepley   PetscScalar    *funcVal, *interpolantVec;
56423d86601SMatthew G. Knepley   PetscReal      *coords, *realSpaceDer, *v0, *J, *invJ, detJ;
56523d86601SMatthew G. Knepley   PetscReal       localDiff = 0.0;
566b7a4d31fSMatthew G. Knepley   PetscInt        dim, numFields, Nc, cStart, cEnd, c, field, fieldOffset, comp;
56723d86601SMatthew G. Knepley   PetscErrorCode  ierr;
56823d86601SMatthew G. Knepley 
56923d86601SMatthew G. Knepley   PetscFunctionBegin;
57023d86601SMatthew G. Knepley   ierr = DMDAGetInfo(dm, &dim,0,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
5712764a2aaSMatthew G. Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
5722764a2aaSMatthew G. Knepley   ierr = PetscDSGetTotalComponents(prob, &Nc);CHKERRQ(ierr);
5732764a2aaSMatthew G. Knepley   ierr = PetscDSGetDiscretization(prob, 0, (PetscObject *) &fe);CHKERRQ(ierr);
574b7a4d31fSMatthew G. Knepley   ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
575e87a4003SBarry Smith   ierr = DMGetSection(dm, &section);CHKERRQ(ierr);
57623d86601SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
57723d86601SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
57823d86601SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
57923d86601SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
58023d86601SMatthew G. Knepley   /* There are no BC values in DAs right now: ierr = DMDAProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */
581b7a4d31fSMatthew G. Knepley   ierr = PetscMalloc7(Nc,&funcVal,dim,&coords,dim,&realSpaceDer,dim,&v0,dim*dim,&J,dim*dim,&invJ,dim,&interpolantVec);CHKERRQ(ierr);
58223d86601SMatthew G. Knepley   ierr = DMDAGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
58323d86601SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
58423d86601SMatthew G. Knepley     PetscScalar *x = NULL;
58523d86601SMatthew G. Knepley     PetscReal    elemDiff = 0.0;
58623d86601SMatthew G. Knepley 
5878e0841e0SMatthew G. Knepley     ierr = DMDAComputeCellGeometryFEM(dm, c, quad, v0, J, invJ, &detJ);CHKERRQ(ierr);
58823d86601SMatthew G. Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
58923d86601SMatthew G. Knepley     ierr = DMDAVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
59023d86601SMatthew G. Knepley 
59123d86601SMatthew G. Knepley     for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) {
59223d86601SMatthew G. Knepley       void * const ctx = ctxs ? ctxs[field] : NULL;
59323d86601SMatthew G. Knepley       const PetscReal *quadPoints, *quadWeights;
59423d86601SMatthew G. Knepley       PetscReal       *basisDer;
5959c3cf19fSMatthew G. Knepley       PetscInt         qNc, Nq, Nb, Nc, q, d, e, fc, f, g;
59623d86601SMatthew G. Knepley 
5972764a2aaSMatthew G. Knepley       ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);CHKERRQ(ierr);
5989c3cf19fSMatthew G. Knepley       ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr);
599b7a4d31fSMatthew G. Knepley       ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
600b7a4d31fSMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
6019c3cf19fSMatthew G. Knepley       if (qNc != Nc) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Quadrature components %D != %D field components\n", qNc, Nc);
602b7a4d31fSMatthew G. Knepley       ierr = PetscFEGetDefaultTabulation(fe, NULL, &basisDer, NULL);CHKERRQ(ierr);
60323d86601SMatthew G. Knepley       if (debug) {
60423d86601SMatthew G. Knepley         char title[1024];
60523d86601SMatthew G. Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
60623d86601SMatthew G. Knepley         ierr = DMPrintCellVector(c, title, Nb*Nc, &x[fieldOffset]);CHKERRQ(ierr);
60723d86601SMatthew G. Knepley       }
60823d86601SMatthew G. Knepley       for (q = 0; q < Nq; ++q) {
60923d86601SMatthew G. Knepley         for (d = 0; d < dim; d++) {
61023d86601SMatthew G. Knepley           coords[d] = v0[d];
61123d86601SMatthew G. Knepley           for (e = 0; e < dim; e++) {
61223d86601SMatthew G. Knepley             coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0);
61323d86601SMatthew G. Knepley           }
61423d86601SMatthew G. Knepley         }
615a98e6df5SMatthew G. Knepley         ierr = (*funcs[field])(dim, time, coords, n, numFields, funcVal, ctx);CHKERRQ(ierr);
61623d86601SMatthew G. Knepley         for (fc = 0; fc < Nc; ++fc) {
61723d86601SMatthew G. Knepley           PetscScalar interpolant = 0.0;
61823d86601SMatthew G. Knepley 
61923d86601SMatthew G. Knepley           for (d = 0; d < dim; ++d) interpolantVec[d] = 0.0;
62023d86601SMatthew G. Knepley           for (f = 0; f < Nb; ++f) {
62123d86601SMatthew G. Knepley             for (d = 0; d < dim; ++d) {
62223d86601SMatthew G. Knepley               realSpaceDer[d] = 0.0;
62323d86601SMatthew G. Knepley               for (g = 0; g < dim; ++g) {
6249c3cf19fSMatthew G. Knepley                 realSpaceDer[d] += invJ[g*dim+d]*basisDer[(q*Nb+f)*dim+g];
62523d86601SMatthew G. Knepley               }
6269c3cf19fSMatthew G. Knepley               interpolantVec[d] += x[fieldOffset+f]*realSpaceDer[d];
62723d86601SMatthew G. Knepley             }
62823d86601SMatthew G. Knepley           }
62923d86601SMatthew G. Knepley           for (d = 0; d < dim; ++d) interpolant += interpolantVec[d]*n[d];
6309c3cf19fSMatthew G. Knepley           if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "    elem %d fieldDer %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q*Nc+c]*detJ);CHKERRQ(ierr);}
6319c3cf19fSMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q*Nc+c]*detJ;
63223d86601SMatthew G. Knepley         }
63323d86601SMatthew G. Knepley       }
63423d86601SMatthew G. Knepley       comp        += Nc;
63523d86601SMatthew G. Knepley       fieldOffset += Nb*Nc;
63623d86601SMatthew G. Knepley     }
63723d86601SMatthew G. Knepley     ierr = DMDAVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
63823d86601SMatthew G. Knepley     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
63923d86601SMatthew G. Knepley     localDiff += elemDiff;
64023d86601SMatthew G. Knepley   }
64123d86601SMatthew G. Knepley   ierr  = PetscFree7(funcVal,coords,realSpaceDer,v0,J,invJ,interpolantVec);CHKERRQ(ierr);
64223d86601SMatthew G. Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
643b2566f29SBarry Smith   ierr  = MPIU_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
64423d86601SMatthew G. Knepley   *diff = PetscSqrtReal(*diff);
64523d86601SMatthew G. Knepley   PetscFunctionReturn(0);
64623d86601SMatthew G. Knepley }
64723d86601SMatthew G. Knepley 
64847c6ae99SBarry Smith /* ------------------------------------------------------------------- */
64947c6ae99SBarry Smith 
65047c6ae99SBarry Smith /*@C
651aa219208SBarry Smith      DMDAGetArray - Gets a work array for a DMDA
65247c6ae99SBarry Smith 
65347c6ae99SBarry Smith     Input Parameter:
65447c6ae99SBarry Smith +    da - information about my local patch
65547c6ae99SBarry Smith -    ghosted - do you want arrays for the ghosted or nonghosted patch
65647c6ae99SBarry Smith 
65747c6ae99SBarry Smith     Output Parameters:
65847c6ae99SBarry Smith .    vptr - array data structured
65947c6ae99SBarry Smith 
66047c6ae99SBarry Smith     Note:  The vector values are NOT initialized and may have garbage in them, so you may need
66147c6ae99SBarry Smith            to zero them.
66247c6ae99SBarry Smith 
66347c6ae99SBarry Smith   Level: advanced
66447c6ae99SBarry Smith 
665bcaeba4dSBarry Smith .seealso: DMDARestoreArray()
66647c6ae99SBarry Smith 
66747c6ae99SBarry Smith @*/
6687087cfbeSBarry Smith PetscErrorCode  DMDAGetArray(DM da,PetscBool ghosted,void *vptr)
66947c6ae99SBarry Smith {
67047c6ae99SBarry Smith   PetscErrorCode ierr;
67147c6ae99SBarry Smith   PetscInt       j,i,xs,ys,xm,ym,zs,zm;
67247c6ae99SBarry Smith   char           *iarray_start;
67347c6ae99SBarry Smith   void           **iptr = (void**)vptr;
67447c6ae99SBarry Smith   DM_DA          *dd    = (DM_DA*)da->data;
67547c6ae99SBarry Smith 
67647c6ae99SBarry Smith   PetscFunctionBegin;
677a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
67847c6ae99SBarry Smith   if (ghosted) {
679aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
68047c6ae99SBarry Smith       if (dd->arrayghostedin[i]) {
68147c6ae99SBarry Smith         *iptr                 = dd->arrayghostedin[i];
68247c6ae99SBarry Smith         iarray_start          = (char*)dd->startghostedin[i];
6830298fd71SBarry Smith         dd->arrayghostedin[i] = NULL;
6840298fd71SBarry Smith         dd->startghostedin[i] = NULL;
68547c6ae99SBarry Smith 
68647c6ae99SBarry Smith         goto done;
68747c6ae99SBarry Smith       }
68847c6ae99SBarry Smith     }
68947c6ae99SBarry Smith     xs = dd->Xs;
69047c6ae99SBarry Smith     ys = dd->Ys;
69147c6ae99SBarry Smith     zs = dd->Zs;
69247c6ae99SBarry Smith     xm = dd->Xe-dd->Xs;
69347c6ae99SBarry Smith     ym = dd->Ye-dd->Ys;
69447c6ae99SBarry Smith     zm = dd->Ze-dd->Zs;
69547c6ae99SBarry Smith   } else {
696aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
69747c6ae99SBarry Smith       if (dd->arrayin[i]) {
69847c6ae99SBarry Smith         *iptr          = dd->arrayin[i];
69947c6ae99SBarry Smith         iarray_start   = (char*)dd->startin[i];
7000298fd71SBarry Smith         dd->arrayin[i] = NULL;
7010298fd71SBarry Smith         dd->startin[i] = NULL;
70247c6ae99SBarry Smith 
70347c6ae99SBarry Smith         goto done;
70447c6ae99SBarry Smith       }
70547c6ae99SBarry Smith     }
70647c6ae99SBarry Smith     xs = dd->xs;
70747c6ae99SBarry Smith     ys = dd->ys;
70847c6ae99SBarry Smith     zs = dd->zs;
70947c6ae99SBarry Smith     xm = dd->xe-dd->xs;
71047c6ae99SBarry Smith     ym = dd->ye-dd->ys;
71147c6ae99SBarry Smith     zm = dd->ze-dd->zs;
71247c6ae99SBarry Smith   }
71347c6ae99SBarry Smith 
714c73cfb54SMatthew G. Knepley   switch (da->dim) {
71547c6ae99SBarry Smith   case 1: {
71647c6ae99SBarry Smith     void *ptr;
71747c6ae99SBarry Smith 
718901f1932SJed Brown     ierr = PetscMalloc(xm*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
71947c6ae99SBarry Smith 
72047c6ae99SBarry Smith     ptr   = (void*)(iarray_start - xs*sizeof(PetscScalar));
72147c6ae99SBarry Smith     *iptr = (void*)ptr;
7228865f1eaSKarl Rupp     break;
7238865f1eaSKarl Rupp   }
72447c6ae99SBarry Smith   case 2: {
72547c6ae99SBarry Smith     void **ptr;
72647c6ae99SBarry Smith 
727901f1932SJed Brown     ierr = PetscMalloc((ym+1)*sizeof(void*)+xm*ym*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
72847c6ae99SBarry Smith 
72947c6ae99SBarry Smith     ptr = (void**)(iarray_start + xm*ym*sizeof(PetscScalar) - ys*sizeof(void*));
7308865f1eaSKarl Rupp     for (j=ys; j<ys+ym; j++) ptr[j] = iarray_start + sizeof(PetscScalar)*(xm*(j-ys) - xs);
73147c6ae99SBarry Smith     *iptr = (void*)ptr;
7328865f1eaSKarl Rupp     break;
7338865f1eaSKarl Rupp   }
73447c6ae99SBarry Smith   case 3: {
73547c6ae99SBarry Smith     void ***ptr,**bptr;
73647c6ae99SBarry Smith 
737901f1932SJed Brown     ierr = PetscMalloc((zm+1)*sizeof(void**)+(ym*zm+1)*sizeof(void*)+xm*ym*zm*sizeof(PetscScalar),&iarray_start);CHKERRQ(ierr);
73847c6ae99SBarry Smith 
73947c6ae99SBarry Smith     ptr  = (void***)(iarray_start + xm*ym*zm*sizeof(PetscScalar) - zs*sizeof(void*));
74047c6ae99SBarry Smith     bptr = (void**)(iarray_start + xm*ym*zm*sizeof(PetscScalar) + zm*sizeof(void**));
7418865f1eaSKarl Rupp     for (i=zs; i<zs+zm; i++) ptr[i] = bptr + ((i-zs)*ym - ys);
74247c6ae99SBarry Smith     for (i=zs; i<zs+zm; i++) {
74347c6ae99SBarry Smith       for (j=ys; j<ys+ym; j++) {
74447c6ae99SBarry Smith         ptr[i][j] = iarray_start + sizeof(PetscScalar)*(xm*ym*(i-zs) + xm*(j-ys) - xs);
74547c6ae99SBarry Smith       }
74647c6ae99SBarry Smith     }
74747c6ae99SBarry Smith     *iptr = (void*)ptr;
7488865f1eaSKarl Rupp     break;
7498865f1eaSKarl Rupp   }
75047c6ae99SBarry Smith   default:
751c73cfb54SMatthew G. Knepley     SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Dimension %D not supported",da->dim);
75247c6ae99SBarry Smith   }
75347c6ae99SBarry Smith 
75447c6ae99SBarry Smith done:
75547c6ae99SBarry Smith   /* add arrays to the checked out list */
75647c6ae99SBarry Smith   if (ghosted) {
757aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
75847c6ae99SBarry Smith       if (!dd->arrayghostedout[i]) {
75947c6ae99SBarry Smith         dd->arrayghostedout[i] = *iptr;
76047c6ae99SBarry Smith         dd->startghostedout[i] = iarray_start;
76147c6ae99SBarry Smith         break;
76247c6ae99SBarry Smith       }
76347c6ae99SBarry Smith     }
76447c6ae99SBarry Smith   } else {
765aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
76647c6ae99SBarry Smith       if (!dd->arrayout[i]) {
76747c6ae99SBarry Smith         dd->arrayout[i] = *iptr;
76847c6ae99SBarry Smith         dd->startout[i] = iarray_start;
76947c6ae99SBarry Smith         break;
77047c6ae99SBarry Smith       }
77147c6ae99SBarry Smith     }
77247c6ae99SBarry Smith   }
77347c6ae99SBarry Smith   PetscFunctionReturn(0);
77447c6ae99SBarry Smith }
77547c6ae99SBarry Smith 
77647c6ae99SBarry Smith /*@C
777aa219208SBarry Smith      DMDARestoreArray - Restores an array of derivative types for a DMDA
77847c6ae99SBarry Smith 
77947c6ae99SBarry Smith     Input Parameter:
78047c6ae99SBarry Smith +    da - information about my local patch
78147c6ae99SBarry Smith .    ghosted - do you want arrays for the ghosted or nonghosted patch
78247c6ae99SBarry Smith -    vptr - array data structured to be passed to ad_FormFunctionLocal()
78347c6ae99SBarry Smith 
78447c6ae99SBarry Smith      Level: advanced
78547c6ae99SBarry Smith 
786bcaeba4dSBarry Smith .seealso: DMDAGetArray()
78747c6ae99SBarry Smith 
78847c6ae99SBarry Smith @*/
7897087cfbeSBarry Smith PetscErrorCode  DMDARestoreArray(DM da,PetscBool ghosted,void *vptr)
79047c6ae99SBarry Smith {
79147c6ae99SBarry Smith   PetscInt i;
79247c6ae99SBarry Smith   void     **iptr = (void**)vptr,*iarray_start = 0;
79347c6ae99SBarry Smith   DM_DA    *dd    = (DM_DA*)da->data;
79447c6ae99SBarry Smith 
79547c6ae99SBarry Smith   PetscFunctionBegin;
796a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da,DM_CLASSID,1,DMDA);
79747c6ae99SBarry Smith   if (ghosted) {
798aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
79947c6ae99SBarry Smith       if (dd->arrayghostedout[i] == *iptr) {
80047c6ae99SBarry Smith         iarray_start           = dd->startghostedout[i];
8010298fd71SBarry Smith         dd->arrayghostedout[i] = NULL;
8020298fd71SBarry Smith         dd->startghostedout[i] = NULL;
80347c6ae99SBarry Smith         break;
80447c6ae99SBarry Smith       }
80547c6ae99SBarry Smith     }
806aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
80747c6ae99SBarry Smith       if (!dd->arrayghostedin[i]) {
80847c6ae99SBarry Smith         dd->arrayghostedin[i] = *iptr;
80947c6ae99SBarry Smith         dd->startghostedin[i] = iarray_start;
81047c6ae99SBarry Smith         break;
81147c6ae99SBarry Smith       }
81247c6ae99SBarry Smith     }
81347c6ae99SBarry Smith   } else {
814aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
81547c6ae99SBarry Smith       if (dd->arrayout[i] == *iptr) {
81647c6ae99SBarry Smith         iarray_start    = dd->startout[i];
8170298fd71SBarry Smith         dd->arrayout[i] = NULL;
8180298fd71SBarry Smith         dd->startout[i] = NULL;
81947c6ae99SBarry Smith         break;
82047c6ae99SBarry Smith       }
82147c6ae99SBarry Smith     }
822aa219208SBarry Smith     for (i=0; i<DMDA_MAX_WORK_ARRAYS; i++) {
82347c6ae99SBarry Smith       if (!dd->arrayin[i]) {
82447c6ae99SBarry Smith         dd->arrayin[i] = *iptr;
82547c6ae99SBarry Smith         dd->startin[i] = iarray_start;
82647c6ae99SBarry Smith         break;
82747c6ae99SBarry Smith       }
82847c6ae99SBarry Smith     }
82947c6ae99SBarry Smith   }
83047c6ae99SBarry Smith   PetscFunctionReturn(0);
83147c6ae99SBarry Smith }
83247c6ae99SBarry Smith 
833