xref: /petsc/src/dm/impls/da/dageometry.c (revision 6693a7312c227fe970795483f14c4ef91ae5294e)
14035e84dSBarry Smith #include <petsc-private/dmdaimpl.h>     /*I  "petscdmda.h"   I*/
257459e9aSMatthew G Knepley 
357459e9aSMatthew G Knepley #undef __FUNCT__
40a3ada39SMatthew G. Knepley #define __FUNCT__ "FillClosureArray_Static"
50a3ada39SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode FillClosureArray_Static(DM dm, PetscSection section, PetscInt nP, const PetscInt points[], PetscScalar *vArray, PetscInt *csize, const PetscScalar **array)
657459e9aSMatthew G Knepley {
757459e9aSMatthew G Knepley   PetscScalar    *a;
80a3ada39SMatthew G. Knepley   PetscInt       pStart, pEnd, size = 0, dof, off, d, k, i;
957459e9aSMatthew G Knepley   PetscErrorCode ierr;
1057459e9aSMatthew G Knepley 
1157459e9aSMatthew G Knepley   PetscFunctionBegin;
120a3ada39SMatthew G. Knepley   ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
1357459e9aSMatthew G Knepley   for (i = 0; i < nP; ++i) {
140a3ada39SMatthew G. Knepley     const PetscInt p = points[i];
150a3ada39SMatthew G. Knepley 
160a3ada39SMatthew G. Knepley     if ((p < pStart) || (p >= pEnd)) continue;
170a3ada39SMatthew G. Knepley     ierr  = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
1857459e9aSMatthew G Knepley     size += dof;
1957459e9aSMatthew G Knepley   }
200a3ada39SMatthew G. Knepley   if (csize) *csize = size;
210a3ada39SMatthew G. Knepley   if (array) {
22aa1993deSMatthew G Knepley     ierr = DMGetWorkArray(dm, size, PETSC_SCALAR, &a);CHKERRQ(ierr);
2357459e9aSMatthew G Knepley     for (i = 0, k = 0; i < nP; ++i) {
240a3ada39SMatthew G. Knepley       const PetscInt p = points[i];
250a3ada39SMatthew G. Knepley 
260a3ada39SMatthew G. Knepley       if ((p < pStart) || (p >= pEnd)) continue;
270a3ada39SMatthew G. Knepley       ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
280a3ada39SMatthew G. Knepley       ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
2957459e9aSMatthew G Knepley 
308865f1eaSKarl Rupp       for (d = 0; d < dof; ++d, ++k) a[k] = vArray[off+d];
3157459e9aSMatthew G Knepley     }
3257459e9aSMatthew G Knepley     *array = a;
330a3ada39SMatthew G. Knepley   }
3457459e9aSMatthew G Knepley   PetscFunctionReturn(0);
3557459e9aSMatthew G Knepley }
3657459e9aSMatthew G Knepley 
3757459e9aSMatthew G Knepley #undef __FUNCT__
3857459e9aSMatthew G Knepley #define __FUNCT__ "FillClosureVec_Private"
3957459e9aSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode FillClosureVec_Private(DM dm, PetscSection section, PetscInt nP, const PetscInt points[], PetscScalar *vArray, const PetscScalar *array, InsertMode mode)
4057459e9aSMatthew G Knepley {
4157459e9aSMatthew G Knepley   PetscInt       dof, off, d, k, i;
4257459e9aSMatthew G Knepley   PetscErrorCode ierr;
4357459e9aSMatthew G Knepley 
4457459e9aSMatthew G Knepley   PetscFunctionBegin;
4557459e9aSMatthew G Knepley   if ((mode == INSERT_VALUES) || (mode == INSERT_ALL_VALUES)) {
4657459e9aSMatthew G Knepley     for (i = 0, k = 0; i < nP; ++i) {
4757459e9aSMatthew G Knepley       ierr = PetscSectionGetDof(section, points[i], &dof);CHKERRQ(ierr);
4857459e9aSMatthew G Knepley       ierr = PetscSectionGetOffset(section, points[i], &off);CHKERRQ(ierr);
4957459e9aSMatthew G Knepley 
508865f1eaSKarl Rupp       for (d = 0; d < dof; ++d, ++k) vArray[off+d] = array[k];
5157459e9aSMatthew G Knepley     }
5257459e9aSMatthew G Knepley   } else {
5357459e9aSMatthew G Knepley     for (i = 0, k = 0; i < nP; ++i) {
5457459e9aSMatthew G Knepley       ierr = PetscSectionGetDof(section, points[i], &dof);CHKERRQ(ierr);
5557459e9aSMatthew G Knepley       ierr = PetscSectionGetOffset(section, points[i], &off);CHKERRQ(ierr);
5657459e9aSMatthew G Knepley 
578865f1eaSKarl Rupp       for (d = 0; d < dof; ++d, ++k) vArray[off+d] += array[k];
5857459e9aSMatthew G Knepley     }
5957459e9aSMatthew G Knepley   }
6057459e9aSMatthew G Knepley   PetscFunctionReturn(0);
6157459e9aSMatthew G Knepley }
6257459e9aSMatthew G Knepley 
6357459e9aSMatthew G Knepley #undef __FUNCT__
64aa1993deSMatthew G Knepley #define __FUNCT__ "GetPointArray_Private"
65aa1993deSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode GetPointArray_Private(DM dm,PetscInt n,PetscInt *points,PetscInt *rn,const PetscInt **rpoints)
66aa1993deSMatthew G Knepley {
67aa1993deSMatthew G Knepley   PetscErrorCode ierr;
68aa1993deSMatthew G Knepley   PetscInt       *work;
69aa1993deSMatthew G Knepley 
70aa1993deSMatthew G Knepley   PetscFunctionBegin;
71aa1993deSMatthew G Knepley   if (rn) *rn = n;
72aa1993deSMatthew G Knepley   if (rpoints) {
73aa1993deSMatthew G Knepley     ierr     = DMGetWorkArray(dm,n,PETSC_INT,&work);CHKERRQ(ierr);
74aa1993deSMatthew G Knepley     ierr     = PetscMemcpy(work,points,n*sizeof(PetscInt));CHKERRQ(ierr);
75aa1993deSMatthew G Knepley     *rpoints = work;
76aa1993deSMatthew G Knepley   }
77aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
78aa1993deSMatthew G Knepley }
79aa1993deSMatthew G Knepley 
80aa1993deSMatthew G Knepley #undef __FUNCT__
81aa1993deSMatthew G Knepley #define __FUNCT__ "RestorePointArray_Private"
82aa1993deSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode RestorePointArray_Private(DM dm,PetscInt *rn,const PetscInt **rpoints)
83aa1993deSMatthew G Knepley {
84aa1993deSMatthew G Knepley   PetscErrorCode ierr;
85aa1993deSMatthew G Knepley 
86aa1993deSMatthew G Knepley   PetscFunctionBegin;
87aa1993deSMatthew G Knepley   if (rn) *rn = 0;
888865f1eaSKarl Rupp   if (rpoints) {
898865f1eaSKarl Rupp     ierr = DMRestoreWorkArray(dm,*rn, PETSC_INT, (void*) rpoints);CHKERRQ(ierr);
908865f1eaSKarl Rupp   }
91aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
92aa1993deSMatthew G Knepley }
93aa1993deSMatthew G Knepley 
94aa1993deSMatthew G Knepley #undef __FUNCT__
95*6693a731SMatthew G. Knepley #define __FUNCT__ "DMDAGetTransitiveClosure"
96*6693a731SMatthew G. Knepley PetscErrorCode DMDAGetTransitiveClosure(DM dm, PetscInt p, PetscBool useClosure, PetscInt *closureSize, PetscInt **closure)
97*6693a731SMatthew G. Knepley {
98*6693a731SMatthew G. Knepley   DM_DA         *da = (DM_DA *) dm->data;
99*6693a731SMatthew G. Knepley   PetscInt       dim = da->dim;
100*6693a731SMatthew G. Knepley   PetscInt       nVx, nVy, nVz, nxF, nXF, nyF, nYF, nzF, nZF;
101*6693a731SMatthew G. Knepley   PetscInt       pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart;
102*6693a731SMatthew G. Knepley   PetscErrorCode ierr;
103*6693a731SMatthew G. Knepley 
104*6693a731SMatthew G. Knepley   PetscFunctionBegin;
105*6693a731SMatthew G. Knepley   if (!useClosure) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Star operation is not yet supported");
106*6693a731SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
107*6693a731SMatthew G. Knepley   ierr = DMDAGetHeightStratum(dm, -1,  &pStart, &pEnd);CHKERRQ(ierr);
108*6693a731SMatthew G. Knepley   ierr = DMDAGetHeightStratum(dm,  0,  &cStart, &cEnd);CHKERRQ(ierr);
109*6693a731SMatthew G. Knepley   ierr = DMDAGetHeightStratum(dm,  1,  &fStart, &fEnd);CHKERRQ(ierr);
110*6693a731SMatthew G. Knepley   ierr = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr);
111*6693a731SMatthew G. Knepley   ierr = DMDAGetNumVertices(dm, &nVx, &nVy, &nVz, NULL);CHKERRQ(ierr);
112*6693a731SMatthew G. Knepley   ierr = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr);
113*6693a731SMatthew G. Knepley   xfStart = fStart; xfEnd = xfStart+nXF;
114*6693a731SMatthew G. Knepley   yfStart = xfEnd;
115*6693a731SMatthew G. Knepley   if ((p < pStart) || (p >= pEnd)) SETERRQ3(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid point %d should be in [%d, %d)", p, pStart, pEnd);
116*6693a731SMatthew G. Knepley   if ((p >= cStart) || (p < cEnd)) {
117*6693a731SMatthew G. Knepley     /* Cell */
118*6693a731SMatthew G. Knepley     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
119*6693a731SMatthew G. Knepley     else if (dim == 2) {
120*6693a731SMatthew G. Knepley       /* 4 faces, 4 vertices
121*6693a731SMatthew G. Knepley          Bottom-left vertex follows same order as cells
122*6693a731SMatthew G. Knepley          Bottom y-face same order as cells
123*6693a731SMatthew G. Knepley          Left x-face follows same order as cells
124*6693a731SMatthew G. Knepley          We number the quad:
125*6693a731SMatthew G. Knepley 
126*6693a731SMatthew G. Knepley            8--3--7
127*6693a731SMatthew G. Knepley            |     |
128*6693a731SMatthew G. Knepley            4  0  2
129*6693a731SMatthew G. Knepley            |     |
130*6693a731SMatthew G. Knepley            5--1--6
131*6693a731SMatthew G. Knepley       */
132*6693a731SMatthew G. Knepley       PetscInt c  = p - cStart, cx = c % (nVx-1), cy = c / (nVx-1);
133*6693a731SMatthew G. Knepley       PetscInt v  = cy*nVx + cx +  vStart;
134*6693a731SMatthew G. Knepley       PetscInt xf = cy*nxF + cx + xfStart;
135*6693a731SMatthew G. Knepley       PetscInt yf = c + yfStart;
136*6693a731SMatthew G. Knepley 
137*6693a731SMatthew G. Knepley       if (closureSize) {PetscValidPointer(closureSize, 4); *closureSize = 9;}
138*6693a731SMatthew G. Knepley       if (!(*closure)) {ierr = DMGetWorkArray(dm, *closureSize, PETSC_INT, closure);CHKERRQ(ierr);}
139*6693a731SMatthew G. Knepley       (*closure)[0] = p; (*closure)[1] = yf; (*closure)[2] = xf+1; (*closure)[3] = yf+nyF; (*closure)[4] = xf+0; (*closure)[5] = v+0; (*closure)[6]= v+1; (*closure)[7] = v+nVx+1; (*closure)[8] = v+nVx+0;
140*6693a731SMatthew G. Knepley     } else {
141*6693a731SMatthew G. Knepley       /* 6 faces, 12 edges, 8 vertices
142*6693a731SMatthew G. Knepley          Bottom-left vertex follows same order as cells
143*6693a731SMatthew G. Knepley          Bottom y-face same order as cells
144*6693a731SMatthew G. Knepley          Left x-face follows same order as cells
145*6693a731SMatthew G. Knepley          We number the quad:
146*6693a731SMatthew G. Knepley 
147*6693a731SMatthew G. Knepley            8--3--7
148*6693a731SMatthew G. Knepley            |     |
149*6693a731SMatthew G. Knepley            4  0  2
150*6693a731SMatthew G. Knepley            |     |
151*6693a731SMatthew G. Knepley            5--1--6
152*6693a731SMatthew G. Knepley       */
153*6693a731SMatthew G. Knepley       PetscInt c  = p - cStart, cx = c % (nVx-1), cy = c / (nVx-1), cz = c / ((nVx-1)*(nVy-1));
154*6693a731SMatthew G. Knepley       PetscInt v  = cy*nVx + cx +  vStart;
155*6693a731SMatthew G. Knepley       PetscInt xf = cy*nxF + cx + xfStart;
156*6693a731SMatthew G. Knepley       PetscInt yf = c + yfStart;
157*6693a731SMatthew G. Knepley 
158*6693a731SMatthew G. Knepley       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
159*6693a731SMatthew G. Knepley       if (closureSize) {PetscValidPointer(closureSize, 4); *closureSize = 26;}
160*6693a731SMatthew G. Knepley       if (!(*closure)) {ierr = DMGetWorkArray(dm, *closureSize, PETSC_INT, closure);CHKERRQ(ierr);}
161*6693a731SMatthew G. Knepley     }
162*6693a731SMatthew G. Knepley   } else if ((p >= vStart) || (p < vEnd)) {
163*6693a731SMatthew G. Knepley     /* Vertex */
164*6693a731SMatthew G. Knepley     if (closureSize) {PetscValidPointer(closureSize, 4); *closureSize = 1;}
165*6693a731SMatthew G. Knepley     if (!(*closure)) {ierr = DMGetWorkArray(dm, *closureSize, PETSC_INT, closure);CHKERRQ(ierr);}
166*6693a731SMatthew G. Knepley     (*closure)[0] = p;
167*6693a731SMatthew G. Knepley   } else if ((p >= fStart) || (p < fStart + nXF)) {
168*6693a731SMatthew G. Knepley     /* X Face */
169*6693a731SMatthew G. Knepley     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D");
170*6693a731SMatthew G. Knepley     else if (dim == 2) {
171*6693a731SMatthew G. Knepley       /* 2 vertices: The bottom vertex has the same numbering as the face */
172*6693a731SMatthew G. Knepley       PetscInt f = p - xfStart;
173*6693a731SMatthew G. Knepley 
174*6693a731SMatthew G. Knepley       if (closureSize) {PetscValidPointer(closureSize, 4); *closureSize = 3;}
175*6693a731SMatthew G. Knepley       if (!(*closure)) {ierr = DMGetWorkArray(dm, *closureSize, PETSC_INT, closure);CHKERRQ(ierr);}
176*6693a731SMatthew G. Knepley       (*closure)[0] = p; (*closure)[1] = f; (*closure)[2] = f+nVx;
177*6693a731SMatthew G. Knepley     } else if (dim == 3) {
178*6693a731SMatthew G. Knepley       /* 4 vertices */
179*6693a731SMatthew G. Knepley       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
180*6693a731SMatthew G. Knepley     }
181*6693a731SMatthew G. Knepley   } else if ((p >= fStart + nXF) || (p < fStart + nXF + nYF)) {
182*6693a731SMatthew G. Knepley     /* Y Face */
183*6693a731SMatthew G. Knepley     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D");
184*6693a731SMatthew G. Knepley     else if (dim == 2) {
185*6693a731SMatthew G. Knepley       /* 2 vertices: The left vertex has the same numbering as the face */
186*6693a731SMatthew G. Knepley       PetscInt f = p - yfStart;
187*6693a731SMatthew G. Knepley 
188*6693a731SMatthew G. Knepley       if (closureSize) {PetscValidPointer(closureSize, 4); *closureSize = 3;}
189*6693a731SMatthew G. Knepley       if (!(*closure)) {ierr = DMGetWorkArray(dm, *closureSize, PETSC_INT, closure);CHKERRQ(ierr);}
190*6693a731SMatthew G. Knepley       (*closure)[0] = p; (*closure)[1] = f; (*closure)[2]= f+1;
191*6693a731SMatthew G. Knepley     } else if (dim == 3) {
192*6693a731SMatthew G. Knepley       /* 4 vertices */
193*6693a731SMatthew G. Knepley       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
194*6693a731SMatthew G. Knepley     }
195*6693a731SMatthew G. Knepley   } else {
196*6693a731SMatthew G. Knepley     /* Z Face */
197*6693a731SMatthew G. Knepley     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D");
198*6693a731SMatthew G. Knepley     else if (dim == 2) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no z-faces in 2D");
199*6693a731SMatthew G. Knepley     else if (dim == 3) {
200*6693a731SMatthew G. Knepley       /* 4 vertices */
201*6693a731SMatthew G. Knepley       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
202*6693a731SMatthew G. Knepley     }
203*6693a731SMatthew G. Knepley   }
204*6693a731SMatthew G. Knepley   PetscFunctionReturn(0);
205*6693a731SMatthew G. Knepley }
206*6693a731SMatthew G. Knepley 
207*6693a731SMatthew G. Knepley #undef __FUNCT__
208*6693a731SMatthew G. Knepley #define __FUNCT__ "DMDARestoreTransitiveClosure"
209*6693a731SMatthew G. Knepley PetscErrorCode DMDARestoreTransitiveClosure(DM dm, PetscInt p, PetscBool useClosure, PetscInt *closureSize, PetscInt **closure)
210*6693a731SMatthew G. Knepley {
211*6693a731SMatthew G. Knepley   PetscErrorCode ierr;
212*6693a731SMatthew G. Knepley 
213*6693a731SMatthew G. Knepley   PetscFunctionBegin;
214*6693a731SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, 0, PETSC_INT, closure);CHKERRQ(ierr);
215*6693a731SMatthew G. Knepley   PetscFunctionReturn(0);
216*6693a731SMatthew G. Knepley }
217*6693a731SMatthew G. Knepley 
218*6693a731SMatthew G. Knepley #undef __FUNCT__
219aa1993deSMatthew G Knepley #define __FUNCT__ "DMDAGetClosure"
220aa1993deSMatthew G Knepley PetscErrorCode DMDAGetClosure(DM dm, PetscSection section, PetscInt p,PetscInt *n,const PetscInt **closure)
22157459e9aSMatthew G Knepley {
22257459e9aSMatthew G Knepley   DM_DA          *da = (DM_DA*) dm->data;
22357459e9aSMatthew G Knepley   PetscInt       dim = da->dim;
224aa1993deSMatthew G Knepley   PetscInt       nVx, nVy, nxF, nXF, nyF, nYF, nzF, nZF;
225aa1993deSMatthew G Knepley   PetscInt       pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart;
226aa1993deSMatthew G Knepley   PetscErrorCode ierr;
227aa1993deSMatthew G Knepley 
228aa1993deSMatthew G Knepley   PetscFunctionBegin;
229aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
230aa1993deSMatthew G Knepley   if (n) PetscValidIntPointer(n,4);
231aa1993deSMatthew G Knepley   if (closure) PetscValidPointer(closure, 5);
232aa1993deSMatthew G Knepley   if (!section) {ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);}
23382f516ccSBarry Smith   if (!section) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "This DM has not default PetscSection");
234aa1993deSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, -1,  &pStart, &pEnd);CHKERRQ(ierr);
235aa1993deSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, 0,   &cStart, &cEnd);CHKERRQ(ierr);
236aa1993deSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, 1,   &fStart, &fEnd);CHKERRQ(ierr);
237aa1993deSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr);
2380298fd71SBarry Smith   ierr    = DMDAGetNumVertices(dm, &nVx, &nVy, NULL, NULL);CHKERRQ(ierr);
239aa1993deSMatthew G Knepley   ierr    = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr);
240aa1993deSMatthew G Knepley   xfStart = fStart; xfEnd = xfStart+nXF;
241aa1993deSMatthew G Knepley   yfStart = xfEnd;
24282f516ccSBarry Smith   if ((p < pStart) || (p >= pEnd)) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid point %d should be in [%d, %d)", p, pStart, pEnd);
243aa1993deSMatthew G Knepley   if ((p >= cStart) || (p < cEnd)) {
244aa1993deSMatthew G Knepley     /* Cell */
24582f516ccSBarry Smith     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
246201c01f8SBarry Smith     else if (dim == 2) {
247aa1993deSMatthew G Knepley       /* 4 faces, 4 vertices
248aa1993deSMatthew G Knepley          Bottom-left vertex follows same order as cells
249aa1993deSMatthew G Knepley          Bottom y-face same order as cells
250aa1993deSMatthew G Knepley          Left x-face follows same order as cells
251aa1993deSMatthew G Knepley          We number the quad:
252aa1993deSMatthew G Knepley 
253aa1993deSMatthew G Knepley            8--3--7
254aa1993deSMatthew G Knepley            |     |
255aa1993deSMatthew G Knepley            4  0  2
256aa1993deSMatthew G Knepley            |     |
257aa1993deSMatthew G Knepley            5--1--6
258aa1993deSMatthew G Knepley       */
259aa1993deSMatthew G Knepley       PetscInt c  = p - cStart, cx = c % (nVx-1), cy = c / (nVx-1);
260aa1993deSMatthew G Knepley       PetscInt v  = cy*nVx + cx +  vStart;
261aa1993deSMatthew G Knepley       PetscInt xf = cy*nxF + cx + xfStart;
262aa1993deSMatthew G Knepley       PetscInt yf = c + yfStart;
263da80777bSKarl Rupp       PetscInt points[9];
264da80777bSKarl Rupp 
265da80777bSKarl Rupp       /* Note: initializer list is not computable at compile time, hence we fill the array manually */
266da80777bSKarl Rupp       points[0] = p; points[1] = yf; points[2] = xf+1; points[3] = yf+nyF; points[4] = xf+0; points[5] = v+0; points[6]= v+1; points[7] = v+nVx+1; points[8] = v+nVx+0;
267aa1993deSMatthew G Knepley 
268aa1993deSMatthew G Knepley       ierr = GetPointArray_Private(dm,9,points,n,closure);CHKERRQ(ierr);
26982f516ccSBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
270aa1993deSMatthew G Knepley   } else if ((p >= vStart) || (p < vEnd)) {
271aa1993deSMatthew G Knepley     /* Vertex */
272aa1993deSMatthew G Knepley     ierr = GetPointArray_Private(dm,1,&p,n,closure);CHKERRQ(ierr);
273aa1993deSMatthew G Knepley   } else if ((p >= fStart) || (p < fStart + nXF)) {
274aa1993deSMatthew G Knepley     /* X Face */
27582f516ccSBarry Smith     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D");
276201c01f8SBarry Smith     else if (dim == 2) {
277aa1993deSMatthew G Knepley       /* 2 vertices: The bottom vertex has the same numbering as the face */
278aa1993deSMatthew G Knepley       PetscInt f = p - xfStart;
279da80777bSKarl Rupp       PetscInt points[3];
280aa1993deSMatthew G Knepley 
281da80777bSKarl Rupp       points[0] = p; points[1] = f; points[2] = f+nVx;
28282f516ccSBarry Smith       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Broken");
283aa1993deSMatthew G Knepley       ierr = GetPointArray_Private(dm,3,points,n,closure);CHKERRQ(ierr);
284aa1993deSMatthew G Knepley     } else if (dim == 3) {
285aa1993deSMatthew G Knepley       /* 4 vertices */
28682f516ccSBarry Smith       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
287aa1993deSMatthew G Knepley     }
288aa1993deSMatthew G Knepley   } else if ((p >= fStart + nXF) || (p < fStart + nXF + nYF)) {
289aa1993deSMatthew G Knepley     /* Y Face */
29082f516ccSBarry Smith     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D");
291201c01f8SBarry Smith     else if (dim == 2) {
292aa1993deSMatthew G Knepley       /* 2 vertices: The left vertex has the same numbering as the face */
293aa1993deSMatthew G Knepley       PetscInt f = p - yfStart;
294da80777bSKarl Rupp       PetscInt points[3];
295aa1993deSMatthew G Knepley 
296da80777bSKarl Rupp       points[0] = p; points[1] = f; points[2]= f+1;
29782f516ccSBarry Smith       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Broken");
298aa1993deSMatthew G Knepley       ierr = GetPointArray_Private(dm, 3, points, n, closure);CHKERRQ(ierr);
299aa1993deSMatthew G Knepley     } else if (dim == 3) {
300aa1993deSMatthew G Knepley       /* 4 vertices */
30182f516ccSBarry Smith       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
302aa1993deSMatthew G Knepley     }
303aa1993deSMatthew G Knepley   } else {
304aa1993deSMatthew G Knepley     /* Z Face */
30582f516ccSBarry Smith     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D");
30682f516ccSBarry Smith     else if (dim == 2) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no z-faces in 2D");
307201c01f8SBarry Smith     else if (dim == 3) {
308aa1993deSMatthew G Knepley       /* 4 vertices */
30982f516ccSBarry Smith       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
310aa1993deSMatthew G Knepley     }
311aa1993deSMatthew G Knepley   }
312aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
313aa1993deSMatthew G Knepley }
314aa1993deSMatthew G Knepley 
315aa1993deSMatthew G Knepley #undef __FUNCT__
316aa1993deSMatthew G Knepley #define __FUNCT__ "DMDARestoreClosure"
317e5c84f05SJed Brown /* Zeros n and closure. */
318aa1993deSMatthew G Knepley PetscErrorCode DMDARestoreClosure(DM dm, PetscSection section, PetscInt p,PetscInt *n,const PetscInt **closure)
319aa1993deSMatthew G Knepley {
320aa1993deSMatthew G Knepley   PetscErrorCode ierr;
321aa1993deSMatthew G Knepley 
322aa1993deSMatthew G Knepley   PetscFunctionBegin;
323aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
324e5c84f05SJed Brown   if (n) PetscValidIntPointer(n,4);
325e5c84f05SJed Brown   if (closure) PetscValidPointer(closure, 5);
326aa1993deSMatthew G Knepley   ierr = RestorePointArray_Private(dm,n,closure);CHKERRQ(ierr);
327aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
328aa1993deSMatthew G Knepley }
329aa1993deSMatthew G Knepley 
330aa1993deSMatthew G Knepley #undef __FUNCT__
331aa1993deSMatthew G Knepley #define __FUNCT__ "DMDAGetClosureScalar"
3320a3ada39SMatthew G. Knepley PetscErrorCode DMDAGetClosureScalar(DM dm, PetscSection section, PetscInt p, PetscScalar *vArray, PetscInt *csize, const PetscScalar **values)
333aa1993deSMatthew G Knepley {
334aa1993deSMatthew G Knepley   DM_DA          *da = (DM_DA*) dm->data;
335aa1993deSMatthew G Knepley   PetscInt       dim = da->dim;
33657459e9aSMatthew G Knepley   PetscInt       nVx, nVy, nxF, nXF, nyF, nYF, nzF, nZF;
33795babc02SJed Brown   PetscInt       pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart, yfEnd, zfStart;
33857459e9aSMatthew G Knepley   PetscErrorCode ierr;
33957459e9aSMatthew G Knepley 
34057459e9aSMatthew G Knepley   PetscFunctionBegin;
34157459e9aSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
342aa1993deSMatthew G Knepley   PetscValidScalarPointer(vArray, 4);
3430a3ada39SMatthew G. Knepley   if (values) PetscValidPointer(values, 6);
34457459e9aSMatthew G Knepley   if (!section) {ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);}
34582f516ccSBarry Smith   if (!section) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "This DM has not default PetscSection");
34657459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, -1,  &pStart, &pEnd);CHKERRQ(ierr);
34757459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, 0,   &cStart, &cEnd);CHKERRQ(ierr);
34857459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, 1,   &fStart, &fEnd);CHKERRQ(ierr);
34957459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr);
3500298fd71SBarry Smith   ierr    = DMDAGetNumVertices(dm, &nVx, &nVy, NULL, NULL);CHKERRQ(ierr);
35157459e9aSMatthew G Knepley   ierr    = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr);
35257459e9aSMatthew G Knepley   xfStart = fStart; xfEnd = xfStart+nXF;
35357459e9aSMatthew G Knepley   yfStart = xfEnd;  yfEnd = yfStart+nYF;
35495babc02SJed Brown   zfStart = yfEnd;
35582f516ccSBarry Smith   if ((p < pStart) || (p >= pEnd)) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid point %d should be in [%d, %d)", p, pStart, pEnd);
35657459e9aSMatthew G Knepley   if ((p >= cStart) || (p < cEnd)) {
35757459e9aSMatthew G Knepley     /* Cell */
35882f516ccSBarry Smith     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
359201c01f8SBarry Smith     else if (dim == 2) {
36057459e9aSMatthew G Knepley       /* 4 faces, 4 vertices
36157459e9aSMatthew G Knepley          Bottom-left vertex follows same order as cells
36257459e9aSMatthew G Knepley          Bottom y-face same order as cells
36357459e9aSMatthew G Knepley          Left x-face follows same order as cells
36457459e9aSMatthew G Knepley          We number the quad:
36557459e9aSMatthew G Knepley 
36657459e9aSMatthew G Knepley            8--3--7
36757459e9aSMatthew G Knepley            |     |
36857459e9aSMatthew G Knepley            4  0  2
36957459e9aSMatthew G Knepley            |     |
37057459e9aSMatthew G Knepley            5--1--6
37157459e9aSMatthew G Knepley       */
37257459e9aSMatthew G Knepley       PetscInt c  = p - cStart, cx = c % (nVx-1), cy = c / (nVx-1);
37357459e9aSMatthew G Knepley       PetscInt v  = cy*nVx + cx +  vStart;
37457459e9aSMatthew G Knepley       PetscInt xf = cy*nxF + cx + xfStart;
37557459e9aSMatthew G Knepley       PetscInt yf = c + yfStart;
376da80777bSKarl Rupp       PetscInt points[9];
37757459e9aSMatthew G Knepley 
378da80777bSKarl Rupp       points[0] = p; points[1] = yf; points[2] = xf+1; points[3] = yf+nyF; points[4] = xf+0; points[5] = v+0; points[6] = v+1; points[7] = v+nVx+1; points[8] = v+nVx+0;
3790a3ada39SMatthew G. Knepley       ierr = FillClosureArray_Static(dm, section, 9, points, vArray, csize, values);CHKERRQ(ierr);
38057459e9aSMatthew G Knepley     } else {
38157459e9aSMatthew G Knepley       /* 6 faces, 8 vertices
38257459e9aSMatthew G Knepley          Bottom-left-back vertex follows same order as cells
38357459e9aSMatthew G Knepley          Back z-face follows same order as cells
38457459e9aSMatthew G Knepley          Bottom y-face follows same order as cells
38557459e9aSMatthew G Knepley          Left x-face follows same order as cells
38657459e9aSMatthew G Knepley 
38757459e9aSMatthew G Knepley               14-----13
38857459e9aSMatthew G Knepley               /|    /|
38957459e9aSMatthew G Knepley              / | 2 / |
39057459e9aSMatthew G Knepley             / 5|  /  |
39157459e9aSMatthew G Knepley           10-----9  4|
39257459e9aSMatthew G Knepley            |  11-|---12
39357459e9aSMatthew G Knepley            |6 /  |  /
39457459e9aSMatthew G Knepley            | /1 3| /
39557459e9aSMatthew G Knepley            |/    |/
39657459e9aSMatthew G Knepley            7-----8
39757459e9aSMatthew G Knepley       */
39857459e9aSMatthew G Knepley       PetscInt c = p - cStart;
399da80777bSKarl Rupp       PetscInt points[15];
400da80777bSKarl Rupp 
401da80777bSKarl Rupp       points[0]  = p; points[1] = c+zfStart; points[2] = c+zfStart+nzF; points[3] = c+yfStart; points[4] = c+xfStart+nxF; points[5] = c+yfStart+nyF; points[6] = c+xfStart;
402da80777bSKarl Rupp       points[7]  = c+vStart+0; points[8] = c+vStart+1; points[9] = c+vStart+nVx+1; points[10] = c+vStart+nVx+0; points[11] = c+vStart+nVx*nVy+0;
403da80777bSKarl Rupp       points[12] = c+vStart+nVx*nVy+1; points[13] = c+vStart+nVx*nVy+nVx+1; points[14] = c+vStart+nVx*nVy+nVx+0;
40457459e9aSMatthew G Knepley 
40582f516ccSBarry Smith       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Broken");
4060a3ada39SMatthew G. Knepley       ierr = FillClosureArray_Static(dm, section, 15, points, vArray, csize, values);CHKERRQ(ierr);
40757459e9aSMatthew G Knepley     }
40857459e9aSMatthew G Knepley   } else if ((p >= vStart) || (p < vEnd)) {
40957459e9aSMatthew G Knepley     /* Vertex */
4100a3ada39SMatthew G. Knepley     ierr = FillClosureArray_Static(dm, section, 1, &p, vArray, csize, values);CHKERRQ(ierr);
41157459e9aSMatthew G Knepley   } else if ((p >= fStart) || (p < fStart + nXF)) {
41257459e9aSMatthew G Knepley     /* X Face */
41382f516ccSBarry Smith     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D");
414201c01f8SBarry Smith     else if (dim == 2) {
41557459e9aSMatthew G Knepley       /* 2 vertices: The bottom vertex has the same numbering as the face */
41657459e9aSMatthew G Knepley       PetscInt f = p - xfStart;
417da80777bSKarl Rupp       PetscInt points[3];
41857459e9aSMatthew G Knepley 
419da80777bSKarl Rupp       points[0] = p; points[1] = f; points[2] = f+nVx;
42082f516ccSBarry Smith       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Broken");
4210a3ada39SMatthew G. Knepley       ierr = FillClosureArray_Static(dm, section, 3, points, vArray, csize, values);CHKERRQ(ierr);
42282f516ccSBarry Smith     } else if (dim == 3) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
42357459e9aSMatthew G Knepley   } else if ((p >= fStart + nXF) || (p < fStart + nXF + nYF)) {
42457459e9aSMatthew G Knepley     /* Y Face */
42582f516ccSBarry Smith     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D");
426201c01f8SBarry Smith     else if (dim == 2) {
42757459e9aSMatthew G Knepley       /* 2 vertices: The left vertex has the same numbering as the face */
42857459e9aSMatthew G Knepley       PetscInt f = p - yfStart;
429da80777bSKarl Rupp       PetscInt points[3];
43057459e9aSMatthew G Knepley 
431da80777bSKarl Rupp       points[0] = p; points[1] = f; points[2] = f+1;
43282f516ccSBarry Smith       SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Broken");
4330a3ada39SMatthew G. Knepley       ierr = FillClosureArray_Static(dm, section, 3, points, vArray, csize, values);CHKERRQ(ierr);
43482f516ccSBarry Smith     } else if (dim == 3) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
43557459e9aSMatthew G Knepley   } else {
43657459e9aSMatthew G Knepley     /* Z Face */
43782f516ccSBarry Smith     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D");
43882f516ccSBarry Smith     else if (dim == 2) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no z-faces in 2D");
43982f516ccSBarry Smith     else if (dim == 3) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
44057459e9aSMatthew G Knepley   }
441aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
442aa1993deSMatthew G Knepley }
443aa1993deSMatthew G Knepley 
444aa1993deSMatthew G Knepley #undef __FUNCT__
445aa1993deSMatthew G Knepley #define __FUNCT__ "DMDAVecGetClosure"
4460a3ada39SMatthew G. Knepley PetscErrorCode DMDAVecGetClosure(DM dm, PetscSection section, Vec v, PetscInt p, PetscInt *csize, const PetscScalar **values)
447aa1993deSMatthew G Knepley {
448aa1993deSMatthew G Knepley   PetscScalar    *vArray;
449aa1993deSMatthew G Knepley   PetscErrorCode ierr;
450aa1993deSMatthew G Knepley 
451aa1993deSMatthew G Knepley   PetscFunctionBegin;
452aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
453aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(v, VEC_CLASSID, 3);
4540a3ada39SMatthew G. Knepley   if (values) PetscValidPointer(values, 6);
455aa1993deSMatthew G Knepley   ierr = VecGetArray(v, &vArray);CHKERRQ(ierr);
4560a3ada39SMatthew G. Knepley   ierr = DMDAGetClosureScalar(dm, section, p, vArray, csize, values);CHKERRQ(ierr);
45757459e9aSMatthew G Knepley   ierr = VecRestoreArray(v, &vArray);CHKERRQ(ierr);
45857459e9aSMatthew G Knepley   PetscFunctionReturn(0);
45957459e9aSMatthew G Knepley }
46057459e9aSMatthew G Knepley 
46157459e9aSMatthew G Knepley #undef __FUNCT__
462aa1993deSMatthew G Knepley #define __FUNCT__ "DMDARestoreClosureScalar"
4630a3ada39SMatthew G. Knepley PetscErrorCode DMDARestoreClosureScalar(DM dm, PetscSection section, PetscInt p, PetscScalar *vArray, PetscInt *csize, const PetscScalar **values)
464aa1993deSMatthew G Knepley {
465aa1993deSMatthew G Knepley   PetscErrorCode ierr;
466aa1993deSMatthew G Knepley 
467aa1993deSMatthew G Knepley   PetscFunctionBegin;
468aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4690a3ada39SMatthew G. Knepley   PetscValidPointer(values, 6);
4700a3ada39SMatthew G. Knepley   ierr  = DMRestoreWorkArray(dm, csize ? *csize : 0, PETSC_SCALAR, (void*) values);CHKERRQ(ierr);
471aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
472aa1993deSMatthew G Knepley }
473aa1993deSMatthew G Knepley 
474aa1993deSMatthew G Knepley #undef __FUNCT__
475aa1993deSMatthew G Knepley #define __FUNCT__ "DMDAVecRestoreClosure"
4760a3ada39SMatthew G. Knepley PetscErrorCode DMDAVecRestoreClosure(DM dm, PetscSection section, Vec v, PetscInt p, PetscInt *csize, const PetscScalar **values)
477aa1993deSMatthew G Knepley {
478aa1993deSMatthew G Knepley   PetscErrorCode ierr;
479aa1993deSMatthew G Knepley 
480aa1993deSMatthew G Knepley   PetscFunctionBegin;
481aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
482aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(v, VEC_CLASSID, 3);
483aa1993deSMatthew G Knepley   PetscValidPointer(values, 5);
4840a3ada39SMatthew G. Knepley   ierr = DMDARestoreClosureScalar(dm, section, p, NULL, csize, values);CHKERRQ(ierr);
485aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
486aa1993deSMatthew G Knepley }
487aa1993deSMatthew G Knepley 
488aa1993deSMatthew G Knepley #undef __FUNCT__
489aa1993deSMatthew G Knepley #define __FUNCT__ "DMDASetClosureScalar"
490aa1993deSMatthew G Knepley PetscErrorCode DMDASetClosureScalar(DM dm, PetscSection section, PetscInt p,PetscScalar *vArray, const PetscScalar *values, InsertMode mode)
49157459e9aSMatthew G Knepley {
49257459e9aSMatthew G Knepley   DM_DA          *da = (DM_DA*) dm->data;
49357459e9aSMatthew G Knepley   PetscInt       dim = da->dim;
4946a7fb054SMatthew G. Knepley   PetscInt       nVx, nVy, nxF, nXF, nyF, nYF, nzF, nZF, nCx, nCy;
49595babc02SJed Brown   PetscInt       pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart, yfEnd, zfStart;
49657459e9aSMatthew G Knepley   PetscErrorCode ierr;
49757459e9aSMatthew G Knepley 
49857459e9aSMatthew G Knepley   PetscFunctionBegin;
49957459e9aSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
500aa1993deSMatthew G Knepley   PetscValidScalarPointer(values, 4);
50157459e9aSMatthew G Knepley   PetscValidPointer(values, 5);
50257459e9aSMatthew G Knepley   if (!section) {ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);}
50382f516ccSBarry Smith   if (!section) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "This DM has not default PetscSection");
50457459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, -1,  &pStart, &pEnd);CHKERRQ(ierr);
50557459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, 0,   &cStart, &cEnd);CHKERRQ(ierr);
50657459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, 1,   &fStart, &fEnd);CHKERRQ(ierr);
50757459e9aSMatthew G Knepley   ierr    = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr);
5086a7fb054SMatthew G. Knepley   ierr    = DMDAGetNumCells(dm, &nCx, &nCy, NULL, NULL);CHKERRQ(ierr);
5090298fd71SBarry Smith   ierr    = DMDAGetNumVertices(dm, &nVx, &nVy, NULL, NULL);CHKERRQ(ierr);
51057459e9aSMatthew G Knepley   ierr    = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr);
51157459e9aSMatthew G Knepley   xfStart = fStart; xfEnd = xfStart+nXF;
51257459e9aSMatthew G Knepley   yfStart = xfEnd;  yfEnd = yfStart+nYF;
51395babc02SJed Brown   zfStart = yfEnd;
51482f516ccSBarry Smith   if ((p < pStart) || (p >= pEnd)) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid point %d should be in [%d, %d)", p, pStart, pEnd);
51557459e9aSMatthew G Knepley   if ((p >= cStart) || (p < cEnd)) {
51657459e9aSMatthew G Knepley     /* Cell */
51782f516ccSBarry Smith     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
518201c01f8SBarry Smith     else if (dim == 2) {
51957459e9aSMatthew G Knepley       /* 4 faces, 4 vertices
52057459e9aSMatthew G Knepley          Bottom-left vertex follows same order as cells
52157459e9aSMatthew G Knepley          Bottom y-face same order as cells
52257459e9aSMatthew G Knepley          Left x-face follows same order as cells
52357459e9aSMatthew G Knepley          We number the quad:
52457459e9aSMatthew G Knepley 
52557459e9aSMatthew G Knepley            8--3--7
52657459e9aSMatthew G Knepley            |     |
52757459e9aSMatthew G Knepley            4  0  2
52857459e9aSMatthew G Knepley            |     |
52957459e9aSMatthew G Knepley            5--1--6
53057459e9aSMatthew G Knepley       */
5316a7fb054SMatthew G. Knepley       PetscInt c = p - cStart, l = c/nCx;
532da80777bSKarl Rupp       PetscInt points[9];
53357459e9aSMatthew G Knepley 
5346a7fb054SMatthew G. Knepley       points[0] = p;
5356a7fb054SMatthew G. Knepley       points[1] = yfStart+c+0;  points[2] = xfStart+l+c+1; points[3] = yfStart+nyF+c+0;  points[4] = xfStart+l+c+0;
5366a7fb054SMatthew G. Knepley       points[5] = vStart+l+c+0; points[6] = vStart+l+c+1;  points[7] = vStart+nVx+l+c+1; points[8] = vStart+nVx+l+c+0;
53757459e9aSMatthew G Knepley       ierr      = FillClosureVec_Private(dm, section, 9, points, vArray, values, mode);CHKERRQ(ierr);
53857459e9aSMatthew G Knepley     } else {
53957459e9aSMatthew G Knepley       /* 6 faces, 8 vertices
54057459e9aSMatthew G Knepley          Bottom-left-back vertex follows same order as cells
54157459e9aSMatthew G Knepley          Back z-face follows same order as cells
54257459e9aSMatthew G Knepley          Bottom y-face follows same order as cells
54357459e9aSMatthew G Knepley          Left x-face follows same order as cells
54457459e9aSMatthew G Knepley 
54557459e9aSMatthew G Knepley               14-----13
54657459e9aSMatthew G Knepley               /|    /|
54757459e9aSMatthew G Knepley              / | 2 / |
54857459e9aSMatthew G Knepley             / 5|  /  |
54957459e9aSMatthew G Knepley           10-----9  4|
55057459e9aSMatthew G Knepley            |  11-|---12
55157459e9aSMatthew G Knepley            |6 /  |  /
55257459e9aSMatthew G Knepley            | /1 3| /
55357459e9aSMatthew G Knepley            |/    |/
55457459e9aSMatthew G Knepley            7-----8
55557459e9aSMatthew G Knepley       */
55657459e9aSMatthew G Knepley       PetscInt c = p - cStart;
557da80777bSKarl Rupp       PetscInt points[15];
55857459e9aSMatthew G Knepley 
559da80777bSKarl Rupp       points[0]  = p; points[1] = c+zfStart; points[2] = c+zfStart+nzF; points[3] = c+yfStart; points[4] = c+xfStart+nxF; points[5] = c+yfStart+nyF; points[6] = c+xfStart;
560da80777bSKarl Rupp       points[7]  = c+vStart+0; points[8] = c+vStart+1; points[9] = c+vStart+nVx+1; points[10] = c+vStart+nVx+0; points[11] = c+vStart+nVx*nVy+0; points[12] = c+vStart+nVx*nVy+1;
561da80777bSKarl Rupp       points[13] = c+vStart+nVx*nVy+nVx+1; points[14] = c+vStart+nVx*nVy+nVx+0;
56257459e9aSMatthew G Knepley       ierr       = FillClosureVec_Private(dm, section, 15, points, vArray, values, mode);CHKERRQ(ierr);
56357459e9aSMatthew G Knepley     }
56457459e9aSMatthew G Knepley   } else if ((p >= vStart) || (p < vEnd)) {
56557459e9aSMatthew G Knepley     /* Vertex */
56657459e9aSMatthew G Knepley     ierr = FillClosureVec_Private(dm, section, 1, &p, vArray, values, mode);CHKERRQ(ierr);
56757459e9aSMatthew G Knepley   } else if ((p >= fStart) || (p < fStart + nXF)) {
56857459e9aSMatthew G Knepley     /* X Face */
56982f516ccSBarry Smith     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D");
570201c01f8SBarry Smith     else if (dim == 2) {
57157459e9aSMatthew G Knepley       /* 2 vertices: The bottom vertex has the same numbering as the face */
57257459e9aSMatthew G Knepley       PetscInt f = p - xfStart;
573da80777bSKarl Rupp       PetscInt points[3];
57457459e9aSMatthew G Knepley 
575da80777bSKarl Rupp       points[0] = p; points[1] = f; points[2] = f+nVx;
57657459e9aSMatthew G Knepley       ierr      = FillClosureVec_Private(dm, section, 3, points, vArray, values, mode);CHKERRQ(ierr);
57782f516ccSBarry Smith     } else if (dim == 3) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
57857459e9aSMatthew G Knepley   } else if ((p >= fStart + nXF) || (p < fStart + nXF + nYF)) {
57957459e9aSMatthew G Knepley     /* Y Face */
58082f516ccSBarry Smith     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D");
581201c01f8SBarry Smith     else if (dim == 2) {
58257459e9aSMatthew G Knepley       /* 2 vertices: The left vertex has the same numbering as the face */
58357459e9aSMatthew G Knepley       PetscInt f = p - yfStart;
584da80777bSKarl Rupp       PetscInt points[3];
58557459e9aSMatthew G Knepley 
586da80777bSKarl Rupp       points[0] = p; points[1] = f; points[2] = f+1;
58757459e9aSMatthew G Knepley       ierr      = FillClosureVec_Private(dm, section, 3, points, vArray, values, mode);CHKERRQ(ierr);
58882f516ccSBarry Smith     } else if (dim == 3) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
58957459e9aSMatthew G Knepley   } else {
59057459e9aSMatthew G Knepley     /* Z Face */
59182f516ccSBarry Smith     if (dim == 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no faces in 1D");
59282f516ccSBarry Smith     else if (dim == 2) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "There are no z-faces in 2D");
59382f516ccSBarry Smith     else if (dim == 3) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Not implemented");
59457459e9aSMatthew G Knepley   }
595aa1993deSMatthew G Knepley   PetscFunctionReturn(0);
596aa1993deSMatthew G Knepley }
597aa1993deSMatthew G Knepley 
598aa1993deSMatthew G Knepley #undef __FUNCT__
599aa1993deSMatthew G Knepley #define __FUNCT__ "DMDAVecSetClosure"
600aa1993deSMatthew G Knepley PetscErrorCode DMDAVecSetClosure(DM dm, PetscSection section, Vec v, PetscInt p, const PetscScalar *values, InsertMode mode)
601aa1993deSMatthew G Knepley {
602aa1993deSMatthew G Knepley   PetscScalar    *vArray;
603aa1993deSMatthew G Knepley   PetscErrorCode ierr;
604aa1993deSMatthew G Knepley 
605aa1993deSMatthew G Knepley   PetscFunctionBegin;
606aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
607aa1993deSMatthew G Knepley   PetscValidHeaderSpecific(v, VEC_CLASSID, 3);
608aa1993deSMatthew G Knepley   PetscValidPointer(values, 5);
609aa1993deSMatthew G Knepley   ierr = VecGetArray(v,&vArray);CHKERRQ(ierr);
610aa1993deSMatthew G Knepley   ierr = DMDASetClosureScalar(dm,section,p,vArray,values,mode);CHKERRQ(ierr);
61157459e9aSMatthew G Knepley   ierr = VecRestoreArray(v,&vArray);CHKERRQ(ierr);
61257459e9aSMatthew G Knepley   PetscFunctionReturn(0);
61357459e9aSMatthew G Knepley }
61457459e9aSMatthew G Knepley 
615ed4e918cSMatthew G Knepley #undef __FUNCT__
616ed4e918cSMatthew G Knepley #define __FUNCT__ "DMDACnvertToCell"
617ed4e918cSMatthew G Knepley /*@
618f0e3914cSSatish Balay   DMDAConvertToCell - Convert (i,j,k) to local cell number
619341c9bdaSMatthew G Knepley 
620ed4e918cSMatthew G Knepley   Not Collective
621ed4e918cSMatthew G Knepley 
622ed4e918cSMatthew G Knepley   Input Parameter:
623ed4e918cSMatthew G Knepley + da - the distributed array
624ed4e918cSMatthew G Knepley = s - A MatStencil giving (i,j,k)
625ed4e918cSMatthew G Knepley 
626ed4e918cSMatthew G Knepley   Output Parameter:
627ed4e918cSMatthew G Knepley . cell - the local cell number
628ed4e918cSMatthew G Knepley 
629ed4e918cSMatthew G Knepley   Level: developer
630ed4e918cSMatthew G Knepley 
631ed4e918cSMatthew G Knepley .seealso: DMDAVecGetClosure()
632ed4e918cSMatthew G Knepley @*/
633ed4e918cSMatthew G Knepley PetscErrorCode DMDAConvertToCell(DM dm, MatStencil s, PetscInt *cell)
634341c9bdaSMatthew G Knepley {
635341c9bdaSMatthew G Knepley   DM_DA          *da = (DM_DA*) dm->data;
636341c9bdaSMatthew G Knepley   const PetscInt dim = da->dim;
6374e846693SMatthew G Knepley   const PetscInt mx  = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys /*, mz = da->Ze - da->Zs*/;
638ed4e918cSMatthew G Knepley   const PetscInt il  = s.i - da->Xs/da->w, jl = dim > 1 ? s.j - da->Ys : 0, kl = dim > 2 ? s.k - da->Zs : 0;
639341c9bdaSMatthew G Knepley 
640341c9bdaSMatthew G Knepley   PetscFunctionBegin;
641ed4e918cSMatthew G Knepley   *cell = -1;
64282f516ccSBarry Smith   if ((s.i < da->Xs/da->w) || (s.i >= da->Xe/da->w))    SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Stencil i %d should be in [%d, %d)", s.i, da->Xs/da->w, da->Xe/da->w);
64382f516ccSBarry Smith   if ((dim > 1) && ((s.j < da->Ys) || (s.j >= da->Ye))) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Stencil j %d should be in [%d, %d)", s.j, da->Ys, da->Ye);
64482f516ccSBarry Smith   if ((dim > 2) && ((s.k < da->Zs) || (s.k >= da->Ze))) SETERRQ3(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Stencil k %d should be in [%d, %d)", s.k, da->Zs, da->Ze);
645ed4e918cSMatthew G Knepley   *cell = (kl*my + jl)*mx + il;
646ed4e918cSMatthew G Knepley   PetscFunctionReturn(0);
647341c9bdaSMatthew G Knepley }
648341c9bdaSMatthew G Knepley 
64957459e9aSMatthew G Knepley #undef __FUNCT__
65057459e9aSMatthew G Knepley #define __FUNCT__ "DMDAComputeCellGeometry_2D"
65115d2e57cSJed Brown PetscErrorCode DMDAComputeCellGeometry_2D(DM dm, const PetscScalar vertices[], const PetscReal refPoint[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
65257459e9aSMatthew G Knepley {
65357459e9aSMatthew G Knepley   const PetscScalar x0   = vertices[0];
65457459e9aSMatthew G Knepley   const PetscScalar y0   = vertices[1];
65557459e9aSMatthew G Knepley   const PetscScalar x1   = vertices[2];
65657459e9aSMatthew G Knepley   const PetscScalar y1   = vertices[3];
65757459e9aSMatthew G Knepley   const PetscScalar x2   = vertices[4];
65857459e9aSMatthew G Knepley   const PetscScalar y2   = vertices[5];
65957459e9aSMatthew G Knepley   const PetscScalar x3   = vertices[6];
66057459e9aSMatthew G Knepley   const PetscScalar y3   = vertices[7];
66157459e9aSMatthew G Knepley   const PetscScalar f_01 = x2 - x1 - x3 + x0;
66257459e9aSMatthew G Knepley   const PetscScalar g_01 = y2 - y1 - y3 + y0;
66357459e9aSMatthew G Knepley   const PetscScalar x    = refPoint[0];
66457459e9aSMatthew G Knepley   const PetscScalar y    = refPoint[1];
66557459e9aSMatthew G Knepley   PetscReal         invDet;
66657459e9aSMatthew G Knepley   PetscErrorCode    ierr;
66757459e9aSMatthew G Knepley 
66857459e9aSMatthew G Knepley   PetscFunctionBegin;
669e477d84eSMatthew G. Knepley #if 0
67015d2e57cSJed Brown   ierr = PetscPrintf(PETSC_COMM_SELF, "Cell (%g,%g)--(%g,%g)--(%g,%g)--(%g,%g)\n",
6714c06c320SMatthew G Knepley                      PetscRealPart(x0),PetscRealPart(y0),PetscRealPart(x1),PetscRealPart(y1),PetscRealPart(x2),PetscRealPart(y2),PetscRealPart(x3),PetscRealPart(y3));CHKERRQ(ierr);
67215d2e57cSJed Brown   ierr = PetscPrintf(PETSC_COMM_SELF, "Ref Point (%g,%g)\n", PetscRealPart(x), PetscRealPart(y));CHKERRQ(ierr);
67315d2e57cSJed Brown #endif
67415d2e57cSJed Brown   J[0]    = PetscRealPart(x1 - x0 + f_01*y) * 0.5; J[1] = PetscRealPart(x3 - x0 + f_01*x) * 0.5;
67515d2e57cSJed Brown   J[2]    = PetscRealPart(y1 - y0 + g_01*y) * 0.5; J[3] = PetscRealPart(y3 - y0 + g_01*x) * 0.5;
67657459e9aSMatthew G Knepley   *detJ   = J[0]*J[3] - J[1]*J[2];
67757459e9aSMatthew G Knepley   invDet  = 1.0/(*detJ);
678e477d84eSMatthew G. Knepley   if (invJ) {
67957459e9aSMatthew G Knepley     invJ[0] =  invDet*J[3]; invJ[1] = -invDet*J[1];
68057459e9aSMatthew G Knepley     invJ[2] = -invDet*J[2]; invJ[3] =  invDet*J[0];
681e477d84eSMatthew G. Knepley   }
68257459e9aSMatthew G Knepley   ierr    = PetscLogFlops(30);CHKERRQ(ierr);
68357459e9aSMatthew G Knepley   PetscFunctionReturn(0);
68457459e9aSMatthew G Knepley }
68557459e9aSMatthew G Knepley 
68657459e9aSMatthew G Knepley #undef __FUNCT__
68757459e9aSMatthew G Knepley #define __FUNCT__ "DMDAComputeCellGeometry"
68857459e9aSMatthew G Knepley PetscErrorCode DMDAComputeCellGeometry(DM dm, PetscInt cell, PetscQuadrature *quad, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal detJ[])
68957459e9aSMatthew G Knepley {
69057459e9aSMatthew G Knepley   DM                 cdm;
69157459e9aSMatthew G Knepley   Vec                coordinates;
692a1e44745SMatthew G. Knepley   const PetscScalar *vertices = NULL;
6930a3ada39SMatthew G. Knepley   PetscInt           csize, dim, d, q;
69457459e9aSMatthew G Knepley   PetscErrorCode     ierr;
69557459e9aSMatthew G Knepley 
69657459e9aSMatthew G Knepley   PetscFunctionBegin;
69757459e9aSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
69857459e9aSMatthew G Knepley   ierr = DMDAGetInfo(dm, &dim, 0,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
699e477d84eSMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
7006636e97aSMatthew G Knepley   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
7010a3ada39SMatthew G. Knepley   ierr = DMDAVecGetClosure(cdm, NULL, coordinates, cell, &csize, &vertices);CHKERRQ(ierr);
7028865f1eaSKarl Rupp   for (d = 0; d < dim; ++d) v0[d] = PetscRealPart(vertices[d]);
70357459e9aSMatthew G Knepley   switch (dim) {
70457459e9aSMatthew G Knepley   case 2:
705f9fd7fdbSMatthew G. Knepley     for (q = 0; q < quad->numPoints; ++q) {
706f9fd7fdbSMatthew G. Knepley       ierr = DMDAComputeCellGeometry_2D(dm, vertices, &quad->points[q*dim], J, invJ, detJ);CHKERRQ(ierr);
70757459e9aSMatthew G Knepley     }
70857459e9aSMatthew G Knepley     break;
70957459e9aSMatthew G Knepley   default:
71082f516ccSBarry Smith     SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Dimension %d not supported", dim);
71157459e9aSMatthew G Knepley   }
7120a3ada39SMatthew G. Knepley   ierr = DMDAVecRestoreClosure(cdm, NULL, coordinates, cell, &csize, &vertices);CHKERRQ(ierr);
71357459e9aSMatthew G Knepley   PetscFunctionReturn(0);
71457459e9aSMatthew G Knepley }
715