xref: /petsc/src/dm/impls/da/dageometry.c (revision 57459e9a2fccad9a1157d35efda8c8494065cac3)
1*57459e9aSMatthew G Knepley #include <petsc-private/daimpl.h>     /*I  "petscdmda.h"   I*/
2*57459e9aSMatthew G Knepley 
3*57459e9aSMatthew G Knepley #undef __FUNCT__
4*57459e9aSMatthew G Knepley #define __FUNCT__ "FillClosureArray_Private"
5*57459e9aSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode FillClosureArray_Private(DM dm, PetscSection section, PetscInt nP, const PetscInt points[], PetscScalar *vArray, const PetscScalar **array)
6*57459e9aSMatthew G Knepley {
7*57459e9aSMatthew G Knepley   PetscScalar   *a;
8*57459e9aSMatthew G Knepley   PetscInt       size = 0, dof, off, d, k, i;
9*57459e9aSMatthew G Knepley   PetscErrorCode ierr;
10*57459e9aSMatthew G Knepley 
11*57459e9aSMatthew G Knepley   PetscFunctionBegin;
12*57459e9aSMatthew G Knepley   for(i = 0; i < nP; ++i) {
13*57459e9aSMatthew G Knepley     ierr = PetscSectionGetDof(section, points[i], &dof);CHKERRQ(ierr);
14*57459e9aSMatthew G Knepley     size += dof;
15*57459e9aSMatthew G Knepley   }
16*57459e9aSMatthew G Knepley   ierr = DMGetWorkArray(dm, 2*size, &a);CHKERRQ(ierr);
17*57459e9aSMatthew G Knepley   for(i = 0, k = 0; i < nP; ++i) {
18*57459e9aSMatthew G Knepley     ierr = PetscSectionGetDof(section, points[i], &dof);CHKERRQ(ierr);
19*57459e9aSMatthew G Knepley     ierr = PetscSectionGetOffset(section, points[i], &off);CHKERRQ(ierr);
20*57459e9aSMatthew G Knepley 
21*57459e9aSMatthew G Knepley     for(d = 0; d < dof; ++d, ++k) {
22*57459e9aSMatthew G Knepley       a[k] = vArray[off+d];
23*57459e9aSMatthew G Knepley     }
24*57459e9aSMatthew G Knepley   }
25*57459e9aSMatthew G Knepley   *array = a;
26*57459e9aSMatthew G Knepley   PetscFunctionReturn(0);
27*57459e9aSMatthew G Knepley }
28*57459e9aSMatthew G Knepley 
29*57459e9aSMatthew G Knepley #undef __FUNCT__
30*57459e9aSMatthew G Knepley #define __FUNCT__ "FillClosureVec_Private"
31*57459e9aSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode FillClosureVec_Private(DM dm, PetscSection section, PetscInt nP, const PetscInt points[], PetscScalar *vArray, const PetscScalar *array, InsertMode mode)
32*57459e9aSMatthew G Knepley {
33*57459e9aSMatthew G Knepley   PetscInt       dof, off, d, k, i;
34*57459e9aSMatthew G Knepley   PetscErrorCode ierr;
35*57459e9aSMatthew G Knepley 
36*57459e9aSMatthew G Knepley   PetscFunctionBegin;
37*57459e9aSMatthew G Knepley   if ((mode == INSERT_VALUES) || (mode == INSERT_ALL_VALUES)) {
38*57459e9aSMatthew G Knepley     for(i = 0, k = 0; i < nP; ++i) {
39*57459e9aSMatthew G Knepley       ierr = PetscSectionGetDof(section, points[i], &dof);CHKERRQ(ierr);
40*57459e9aSMatthew G Knepley       ierr = PetscSectionGetOffset(section, points[i], &off);CHKERRQ(ierr);
41*57459e9aSMatthew G Knepley 
42*57459e9aSMatthew G Knepley       for(d = 0; d < dof; ++d, ++k) {
43*57459e9aSMatthew G Knepley         vArray[off+d] = array[k];
44*57459e9aSMatthew G Knepley       }
45*57459e9aSMatthew G Knepley     }
46*57459e9aSMatthew G Knepley   } else {
47*57459e9aSMatthew G Knepley     for(i = 0, k = 0; i < nP; ++i) {
48*57459e9aSMatthew G Knepley       ierr = PetscSectionGetDof(section, points[i], &dof);CHKERRQ(ierr);
49*57459e9aSMatthew G Knepley       ierr = PetscSectionGetOffset(section, points[i], &off);CHKERRQ(ierr);
50*57459e9aSMatthew G Knepley 
51*57459e9aSMatthew G Knepley       for(d = 0; d < dof; ++d, ++k) {
52*57459e9aSMatthew G Knepley         vArray[off+d] += array[k];
53*57459e9aSMatthew G Knepley       }
54*57459e9aSMatthew G Knepley     }
55*57459e9aSMatthew G Knepley   }
56*57459e9aSMatthew G Knepley   PetscFunctionReturn(0);
57*57459e9aSMatthew G Knepley }
58*57459e9aSMatthew G Knepley 
59*57459e9aSMatthew G Knepley #undef __FUNCT__
60*57459e9aSMatthew G Knepley #define __FUNCT__ "DMDAVecGetClosure"
61*57459e9aSMatthew G Knepley PetscErrorCode DMDAVecGetClosure(DM dm, PetscSection section, Vec v, PetscInt p, const PetscScalar **values)
62*57459e9aSMatthew G Knepley {
63*57459e9aSMatthew G Knepley   DM_DA         *da  = (DM_DA *) dm->data;
64*57459e9aSMatthew G Knepley   PetscInt       dim = da->dim;
65*57459e9aSMatthew G Knepley   PetscScalar   *vArray;
66*57459e9aSMatthew G Knepley   PetscInt       nVx, nVy, nxF, nXF, nyF, nYF, nzF, nZF;
67*57459e9aSMatthew G Knepley   PetscInt       pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart, yfEnd, zfStart, zfEnd;
68*57459e9aSMatthew G Knepley   PetscErrorCode ierr;
69*57459e9aSMatthew G Knepley 
70*57459e9aSMatthew G Knepley   PetscFunctionBegin;
71*57459e9aSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
72*57459e9aSMatthew G Knepley   PetscValidHeaderSpecific(v, VEC_CLASSID, 3);
73*57459e9aSMatthew G Knepley   PetscValidPointer(values, 5);
74*57459e9aSMatthew G Knepley   if (!section) {ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);}
75*57459e9aSMatthew G Knepley   if (!section) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_ARG_WRONG, "This DM has not default PetscSection");
76*57459e9aSMatthew G Knepley   ierr = DMDAGetHeightStratum(dm, -1,  &pStart, &pEnd);CHKERRQ(ierr);
77*57459e9aSMatthew G Knepley   ierr = DMDAGetHeightStratum(dm, 0,   &cStart, &cEnd);CHKERRQ(ierr);
78*57459e9aSMatthew G Knepley   ierr = DMDAGetHeightStratum(dm, 1,   &fStart, &fEnd);CHKERRQ(ierr);
79*57459e9aSMatthew G Knepley   ierr = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr);
80*57459e9aSMatthew G Knepley   ierr = DMDAGetNumVertices(dm, &nVx, &nVy, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
81*57459e9aSMatthew G Knepley   ierr = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr);
82*57459e9aSMatthew G Knepley   xfStart = fStart; xfEnd = xfStart+nXF;
83*57459e9aSMatthew G Knepley   yfStart = xfEnd;  yfEnd = yfStart+nYF;
84*57459e9aSMatthew G Knepley   zfStart = yfEnd;  zfEnd = zfStart+nZF;
85*57459e9aSMatthew G Knepley   if ((p < pStart) || (p >= pEnd)) SETERRQ3(((PetscObject) dm)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid point %d should be in [%d, %d)", p, pStart, pEnd);
86*57459e9aSMatthew G Knepley   ierr = VecGetArray(v, &vArray);CHKERRQ(ierr);
87*57459e9aSMatthew G Knepley   if ((p >= cStart) || (p < cEnd)) {
88*57459e9aSMatthew G Knepley     /* Cell */
89*57459e9aSMatthew G Knepley     if (dim == 1) {
90*57459e9aSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
91*57459e9aSMatthew G Knepley     } else if (dim == 2) {
92*57459e9aSMatthew G Knepley       /* 4 faces, 4 vertices
93*57459e9aSMatthew G Knepley          Bottom-left vertex follows same order as cells
94*57459e9aSMatthew G Knepley          Bottom y-face same order as cells
95*57459e9aSMatthew G Knepley          Left x-face follows same order as cells
96*57459e9aSMatthew G Knepley          We number the quad:
97*57459e9aSMatthew G Knepley 
98*57459e9aSMatthew G Knepley            8--3--7
99*57459e9aSMatthew G Knepley            |     |
100*57459e9aSMatthew G Knepley            4  0  2
101*57459e9aSMatthew G Knepley            |     |
102*57459e9aSMatthew G Knepley            5--1--6
103*57459e9aSMatthew G Knepley       */
104*57459e9aSMatthew G Knepley       PetscInt c         = p - cStart, cx = c % (nVx-1), cy = c / (nVx-1);
105*57459e9aSMatthew G Knepley       PetscInt v         = cy*nVx + cx +  vStart;
106*57459e9aSMatthew G Knepley       PetscInt xf        = cy*nxF + cx + xfStart;
107*57459e9aSMatthew G Knepley       PetscInt yf        = c + yfStart;
108*57459e9aSMatthew G Knepley       PetscInt points[9] = {p, yf, xf+1, yf+nyF, xf+0, v+0, v+1, v+nVx+1, v+nVx+0};
109*57459e9aSMatthew G Knepley 
110*57459e9aSMatthew G Knepley       ierr = FillClosureArray_Private(dm, section, 9, points, vArray, values);CHKERRQ(ierr);
111*57459e9aSMatthew G Knepley     } else {
112*57459e9aSMatthew G Knepley       /* 6 faces, 8 vertices
113*57459e9aSMatthew G Knepley          Bottom-left-back vertex follows same order as cells
114*57459e9aSMatthew G Knepley          Back z-face follows same order as cells
115*57459e9aSMatthew G Knepley          Bottom y-face follows same order as cells
116*57459e9aSMatthew G Knepley          Left x-face follows same order as cells
117*57459e9aSMatthew G Knepley 
118*57459e9aSMatthew G Knepley               14-----13
119*57459e9aSMatthew G Knepley               /|    /|
120*57459e9aSMatthew G Knepley              / | 2 / |
121*57459e9aSMatthew G Knepley             / 5|  /  |
122*57459e9aSMatthew G Knepley           10-----9  4|
123*57459e9aSMatthew G Knepley            |  11-|---12
124*57459e9aSMatthew G Knepley            |6 /  |  /
125*57459e9aSMatthew G Knepley            | /1 3| /
126*57459e9aSMatthew G Knepley            |/    |/
127*57459e9aSMatthew G Knepley            7-----8
128*57459e9aSMatthew G Knepley       */
129*57459e9aSMatthew G Knepley       PetscInt c          = p - cStart;
130*57459e9aSMatthew G Knepley       PetscInt points[15] = {p, c+zfStart, c+zfStart+nzF, c+yfStart, c+xfStart+nxF, c+yfStart+nyF, c+xfStart,
131*57459e9aSMatthew G Knepley                              c+vStart+0, c+vStart+1, c+vStart+nVx+1, c+vStart+nVx+0, c+vStart+nVx*nVy+0, c+vStart+nVx*nVy+1, c+vStart+nVx*nVy+nVx+1, c+vStart+nVx*nVy+nVx+0};
132*57459e9aSMatthew G Knepley 
133*57459e9aSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Broken");
134*57459e9aSMatthew G Knepley       ierr = FillClosureArray_Private(dm, section, 15, points, vArray, values);CHKERRQ(ierr);
135*57459e9aSMatthew G Knepley     }
136*57459e9aSMatthew G Knepley   } else if ((p >= vStart) || (p < vEnd)) {
137*57459e9aSMatthew G Knepley     /* Vertex */
138*57459e9aSMatthew G Knepley     ierr = FillClosureArray_Private(dm, section, 1, &p, vArray, values);CHKERRQ(ierr);
139*57459e9aSMatthew G Knepley   } else if ((p >= fStart) || (p < fStart + nXF)) {
140*57459e9aSMatthew G Knepley     /* X Face */
141*57459e9aSMatthew G Knepley     if (dim == 1) {
142*57459e9aSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D");
143*57459e9aSMatthew G Knepley     } else if (dim == 2) {
144*57459e9aSMatthew G Knepley       /* 2 vertices: The bottom vertex has the same numbering as the face */
145*57459e9aSMatthew G Knepley       PetscInt f         = p - xfStart;
146*57459e9aSMatthew G Knepley       PetscInt points[3] = {p, f, f+nVx};
147*57459e9aSMatthew G Knepley 
148*57459e9aSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Broken");
149*57459e9aSMatthew G Knepley       ierr = FillClosureArray_Private(dm, section, 3, points, vArray, values);CHKERRQ(ierr);
150*57459e9aSMatthew G Knepley     } else if (dim == 3) {
151*57459e9aSMatthew G Knepley       /* 4 vertices */
152*57459e9aSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
153*57459e9aSMatthew G Knepley     }
154*57459e9aSMatthew G Knepley   } else if ((p >= fStart + nXF) || (p < fStart + nXF + nYF)) {
155*57459e9aSMatthew G Knepley     /* Y Face */
156*57459e9aSMatthew G Knepley     if (dim == 1) {
157*57459e9aSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D");
158*57459e9aSMatthew G Knepley     } else if (dim == 2) {
159*57459e9aSMatthew G Knepley       /* 2 vertices: The left vertex has the same numbering as the face */
160*57459e9aSMatthew G Knepley       PetscInt f         = p - yfStart;
161*57459e9aSMatthew G Knepley       PetscInt points[3] = {p, f, f+1};
162*57459e9aSMatthew G Knepley 
163*57459e9aSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Broken");
164*57459e9aSMatthew G Knepley       ierr = FillClosureArray_Private(dm, section, 3, points, vArray, values);CHKERRQ(ierr);
165*57459e9aSMatthew G Knepley     } else if (dim == 3) {
166*57459e9aSMatthew G Knepley       /* 4 vertices */
167*57459e9aSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
168*57459e9aSMatthew G Knepley     }
169*57459e9aSMatthew G Knepley   } else {
170*57459e9aSMatthew G Knepley     /* Z Face */
171*57459e9aSMatthew G Knepley     if (dim == 1) {
172*57459e9aSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D");
173*57459e9aSMatthew G Knepley     } else if (dim == 2) {
174*57459e9aSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no z-faces in 2D");
175*57459e9aSMatthew G Knepley     } else if (dim == 3) {
176*57459e9aSMatthew G Knepley       /* 4 vertices */
177*57459e9aSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
178*57459e9aSMatthew G Knepley     }
179*57459e9aSMatthew G Knepley   }
180*57459e9aSMatthew G Knepley   ierr = VecRestoreArray(v, &vArray);CHKERRQ(ierr);
181*57459e9aSMatthew G Knepley   PetscFunctionReturn(0);
182*57459e9aSMatthew G Knepley }
183*57459e9aSMatthew G Knepley 
184*57459e9aSMatthew G Knepley #undef __FUNCT__
185*57459e9aSMatthew G Knepley #define __FUNCT__ "DMDAVecSetClosure"
186*57459e9aSMatthew G Knepley PetscErrorCode DMDAVecSetClosure(DM dm, PetscSection section, Vec v, PetscInt p, const PetscScalar *values, InsertMode mode)
187*57459e9aSMatthew G Knepley {
188*57459e9aSMatthew G Knepley   DM_DA         *da  = (DM_DA *) dm->data;
189*57459e9aSMatthew G Knepley   PetscInt       dim = da->dim;
190*57459e9aSMatthew G Knepley   PetscScalar   *vArray;
191*57459e9aSMatthew G Knepley   PetscInt       nVx, nVy, nxF, nXF, nyF, nYF, nzF, nZF;
192*57459e9aSMatthew G Knepley   PetscInt       pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart, yfEnd, zfStart, zfEnd;
193*57459e9aSMatthew G Knepley   PetscErrorCode ierr;
194*57459e9aSMatthew G Knepley 
195*57459e9aSMatthew G Knepley   PetscFunctionBegin;
196*57459e9aSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
197*57459e9aSMatthew G Knepley   PetscValidHeaderSpecific(v, VEC_CLASSID, 3);
198*57459e9aSMatthew G Knepley   PetscValidPointer(values, 5);
199*57459e9aSMatthew G Knepley   if (!section) {ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);}
200*57459e9aSMatthew G Knepley   if (!section) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_ARG_WRONG, "This DM has not default PetscSection");
201*57459e9aSMatthew G Knepley   ierr = DMDAGetHeightStratum(dm, -1,  &pStart, &pEnd);CHKERRQ(ierr);
202*57459e9aSMatthew G Knepley   ierr = DMDAGetHeightStratum(dm, 0,   &cStart, &cEnd);CHKERRQ(ierr);
203*57459e9aSMatthew G Knepley   ierr = DMDAGetHeightStratum(dm, 1,   &fStart, &fEnd);CHKERRQ(ierr);
204*57459e9aSMatthew G Knepley   ierr = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr);
205*57459e9aSMatthew G Knepley   ierr = DMDAGetNumVertices(dm, &nVx, &nVy, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
206*57459e9aSMatthew G Knepley   ierr = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr);
207*57459e9aSMatthew G Knepley   xfStart = fStart; xfEnd = xfStart+nXF;
208*57459e9aSMatthew G Knepley   yfStart = xfEnd;  yfEnd = yfStart+nYF;
209*57459e9aSMatthew G Knepley   zfStart = yfEnd;  zfEnd = zfStart+nZF;
210*57459e9aSMatthew G Knepley   if ((p < pStart) || (p >= pEnd)) SETERRQ3(((PetscObject) dm)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid point %d should be in [%d, %d)", p, pStart, pEnd);
211*57459e9aSMatthew G Knepley   ierr = VecGetArray(v, &vArray);CHKERRQ(ierr);
212*57459e9aSMatthew G Knepley   if ((p >= cStart) || (p < cEnd)) {
213*57459e9aSMatthew G Knepley     /* Cell */
214*57459e9aSMatthew G Knepley     if (dim == 1) {
215*57459e9aSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
216*57459e9aSMatthew G Knepley     } else if (dim == 2) {
217*57459e9aSMatthew G Knepley       /* 4 faces, 4 vertices
218*57459e9aSMatthew G Knepley          Bottom-left vertex follows same order as cells
219*57459e9aSMatthew G Knepley          Bottom y-face same order as cells
220*57459e9aSMatthew G Knepley          Left x-face follows same order as cells
221*57459e9aSMatthew G Knepley          We number the quad:
222*57459e9aSMatthew G Knepley 
223*57459e9aSMatthew G Knepley            8--3--7
224*57459e9aSMatthew G Knepley            |     |
225*57459e9aSMatthew G Knepley            4  0  2
226*57459e9aSMatthew G Knepley            |     |
227*57459e9aSMatthew G Knepley            5--1--6
228*57459e9aSMatthew G Knepley       */
229*57459e9aSMatthew G Knepley       PetscInt c         = p - cStart;
230*57459e9aSMatthew G Knepley       PetscInt points[9] = {p, c+yfStart, c+xfStart+1, c+yfStart+nyF, c+xfStart+0, c+vStart+0, c+vStart+1, c+vStart+nVx+1, c+vStart+nVx+0};
231*57459e9aSMatthew G Knepley 
232*57459e9aSMatthew G Knepley       ierr = FillClosureVec_Private(dm, section, 9, points, vArray, values, mode);CHKERRQ(ierr);
233*57459e9aSMatthew G Knepley     } else {
234*57459e9aSMatthew G Knepley       /* 6 faces, 8 vertices
235*57459e9aSMatthew G Knepley          Bottom-left-back vertex follows same order as cells
236*57459e9aSMatthew G Knepley          Back z-face follows same order as cells
237*57459e9aSMatthew G Knepley          Bottom y-face follows same order as cells
238*57459e9aSMatthew G Knepley          Left x-face follows same order as cells
239*57459e9aSMatthew G Knepley 
240*57459e9aSMatthew G Knepley               14-----13
241*57459e9aSMatthew G Knepley               /|    /|
242*57459e9aSMatthew G Knepley              / | 2 / |
243*57459e9aSMatthew G Knepley             / 5|  /  |
244*57459e9aSMatthew G Knepley           10-----9  4|
245*57459e9aSMatthew G Knepley            |  11-|---12
246*57459e9aSMatthew G Knepley            |6 /  |  /
247*57459e9aSMatthew G Knepley            | /1 3| /
248*57459e9aSMatthew G Knepley            |/    |/
249*57459e9aSMatthew G Knepley            7-----8
250*57459e9aSMatthew G Knepley       */
251*57459e9aSMatthew G Knepley       PetscInt c          = p - cStart;
252*57459e9aSMatthew G Knepley       PetscInt points[15] = {p, c+zfStart, c+zfStart+nzF, c+yfStart, c+xfStart+nxF, c+yfStart+nyF, c+xfStart,
253*57459e9aSMatthew G Knepley                              c+vStart+0, c+vStart+1, c+vStart+nVx+1, c+vStart+nVx+0, c+vStart+nVx*nVy+0, c+vStart+nVx*nVy+1, c+vStart+nVx*nVy+nVx+1, c+vStart+nVx*nVy+nVx+0};
254*57459e9aSMatthew G Knepley 
255*57459e9aSMatthew G Knepley       ierr = FillClosureVec_Private(dm, section, 15, points, vArray, values, mode);CHKERRQ(ierr);
256*57459e9aSMatthew G Knepley     }
257*57459e9aSMatthew G Knepley   } else if ((p >= vStart) || (p < vEnd)) {
258*57459e9aSMatthew G Knepley     /* Vertex */
259*57459e9aSMatthew G Knepley     ierr = FillClosureVec_Private(dm, section, 1, &p, vArray, values, mode);CHKERRQ(ierr);
260*57459e9aSMatthew G Knepley   } else if ((p >= fStart) || (p < fStart + nXF)) {
261*57459e9aSMatthew G Knepley     /* X Face */
262*57459e9aSMatthew G Knepley     if (dim == 1) {
263*57459e9aSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D");
264*57459e9aSMatthew G Knepley     } else if (dim == 2) {
265*57459e9aSMatthew G Knepley       /* 2 vertices: The bottom vertex has the same numbering as the face */
266*57459e9aSMatthew G Knepley       PetscInt f         = p - xfStart;
267*57459e9aSMatthew G Knepley       PetscInt points[3] = {p, f, f+nVx};
268*57459e9aSMatthew G Knepley 
269*57459e9aSMatthew G Knepley       ierr = FillClosureVec_Private(dm, section, 3, points, vArray, values, mode);CHKERRQ(ierr);
270*57459e9aSMatthew G Knepley     } else if (dim == 3) {
271*57459e9aSMatthew G Knepley       /* 4 vertices */
272*57459e9aSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
273*57459e9aSMatthew G Knepley     }
274*57459e9aSMatthew G Knepley   } else if ((p >= fStart + nXF) || (p < fStart + nXF + nYF)) {
275*57459e9aSMatthew G Knepley     /* Y Face */
276*57459e9aSMatthew G Knepley     if (dim == 1) {
277*57459e9aSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D");
278*57459e9aSMatthew G Knepley     } else if (dim == 2) {
279*57459e9aSMatthew G Knepley       /* 2 vertices: The left vertex has the same numbering as the face */
280*57459e9aSMatthew G Knepley       PetscInt f         = p - yfStart;
281*57459e9aSMatthew G Knepley       PetscInt points[3] = {p, f, f+1};
282*57459e9aSMatthew G Knepley 
283*57459e9aSMatthew G Knepley       ierr = FillClosureVec_Private(dm, section, 3, points, vArray, values, mode);CHKERRQ(ierr);
284*57459e9aSMatthew G Knepley     } else if (dim == 3) {
285*57459e9aSMatthew G Knepley       /* 4 vertices */
286*57459e9aSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
287*57459e9aSMatthew G Knepley     }
288*57459e9aSMatthew G Knepley   } else {
289*57459e9aSMatthew G Knepley     /* Z Face */
290*57459e9aSMatthew G Knepley     if (dim == 1) {
291*57459e9aSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D");
292*57459e9aSMatthew G Knepley     } else if (dim == 2) {
293*57459e9aSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no z-faces in 2D");
294*57459e9aSMatthew G Knepley     } else if (dim == 3) {
295*57459e9aSMatthew G Knepley       /* 4 vertices */
296*57459e9aSMatthew G Knepley       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
297*57459e9aSMatthew G Knepley     }
298*57459e9aSMatthew G Knepley   }
299*57459e9aSMatthew G Knepley   ierr = VecRestoreArray(v, &vArray);CHKERRQ(ierr);
300*57459e9aSMatthew G Knepley   PetscFunctionReturn(0);
301*57459e9aSMatthew G Knepley }
302*57459e9aSMatthew G Knepley 
303*57459e9aSMatthew G Knepley #undef __FUNCT__
304*57459e9aSMatthew G Knepley #define __FUNCT__ "DMDAComputeCellGeometry_2D"
305*57459e9aSMatthew G Knepley PetscErrorCode DMDAComputeCellGeometry_2D(DM dm, const PetscScalar vertices[], PetscReal refPoint[], PetscReal J[], PetscReal invJ[], PetscScalar *detJ)
306*57459e9aSMatthew G Knepley {
307*57459e9aSMatthew G Knepley   const PetscScalar x0   = vertices[0];
308*57459e9aSMatthew G Knepley   const PetscScalar y0   = vertices[1];
309*57459e9aSMatthew G Knepley   const PetscScalar x1   = vertices[2];
310*57459e9aSMatthew G Knepley   const PetscScalar y1   = vertices[3];
311*57459e9aSMatthew G Knepley   const PetscScalar x2   = vertices[4];
312*57459e9aSMatthew G Knepley   const PetscScalar y2   = vertices[5];
313*57459e9aSMatthew G Knepley   const PetscScalar x3   = vertices[6];
314*57459e9aSMatthew G Knepley   const PetscScalar y3   = vertices[7];
315*57459e9aSMatthew G Knepley   const PetscScalar f_01 = x2 - x1 - x3 + x0;
316*57459e9aSMatthew G Knepley   const PetscScalar g_01 = y2 - y1 - y3 + y0;
317*57459e9aSMatthew G Knepley   const PetscScalar x    = refPoint[0];
318*57459e9aSMatthew G Knepley   const PetscScalar y    = refPoint[1];
319*57459e9aSMatthew G Knepley   PetscReal         invDet;
320*57459e9aSMatthew G Knepley   PetscErrorCode    ierr;
321*57459e9aSMatthew G Knepley 
322*57459e9aSMatthew G Knepley   PetscFunctionBegin;
323*57459e9aSMatthew G Knepley   ierr = PetscPrintf(PETSC_COMM_SELF, "Cell (%g,%g)--(%g,%g)--(%g,%g)--(%g,%g)\n", x0, y0, x1, y1, x2, y2, x3, y3);CHKERRQ(ierr);
324*57459e9aSMatthew G Knepley   ierr = PetscPrintf(PETSC_COMM_SELF, "Ref Point (%g,%g)\n", x, y);CHKERRQ(ierr);
325*57459e9aSMatthew G Knepley   J[0] = (x1 - x0 + f_01*y) * 0.5; J[1] = (x3 - x0 + f_01*x) * 0.5;
326*57459e9aSMatthew G Knepley   J[2] = (y1 - y0 + g_01*y) * 0.5; J[3] = (y3 - y0 + g_01*x) * 0.5;
327*57459e9aSMatthew G Knepley   *detJ   = J[0]*J[3] - J[1]*J[2];
328*57459e9aSMatthew G Knepley   invDet  = 1.0/(*detJ);
329*57459e9aSMatthew G Knepley   invJ[0] =  invDet*J[3]; invJ[1] = -invDet*J[1];
330*57459e9aSMatthew G Knepley   invJ[2] = -invDet*J[2]; invJ[3] =  invDet*J[0];
331*57459e9aSMatthew G Knepley   ierr = PetscLogFlops(30);CHKERRQ(ierr);
332*57459e9aSMatthew G Knepley   PetscFunctionReturn(0);
333*57459e9aSMatthew G Knepley }
334*57459e9aSMatthew G Knepley 
335*57459e9aSMatthew G Knepley #undef __FUNCT__
336*57459e9aSMatthew G Knepley #define __FUNCT__ "DMDAComputeCellGeometry"
337*57459e9aSMatthew G Knepley PetscErrorCode DMDAComputeCellGeometry(DM dm, PetscInt cell, PetscQuadrature *quad, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal detJ[])
338*57459e9aSMatthew G Knepley {
339*57459e9aSMatthew G Knepley   DM                 cdm;
340*57459e9aSMatthew G Knepley   Vec                coordinates;
341*57459e9aSMatthew G Knepley   const PetscScalar *vertices;
342*57459e9aSMatthew G Knepley   PetscInt           dim, d, q;
343*57459e9aSMatthew G Knepley   PetscErrorCode     ierr;
344*57459e9aSMatthew G Knepley 
345*57459e9aSMatthew G Knepley   PetscFunctionBegin;
346*57459e9aSMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
347*57459e9aSMatthew G Knepley   ierr = DMDAGetInfo(dm, &dim, 0,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
348*57459e9aSMatthew G Knepley   ierr = DMDAGetCoordinates(dm, &coordinates);CHKERRQ(ierr);
349*57459e9aSMatthew G Knepley   ierr = DMDAGetCoordinateDA(dm, &cdm);CHKERRQ(ierr);
350*57459e9aSMatthew G Knepley   ierr = DMDAVecGetClosure(cdm, PETSC_NULL, coordinates, cell, &vertices);CHKERRQ(ierr);
351*57459e9aSMatthew G Knepley   for(d = 0; d < dim; ++d) {
352*57459e9aSMatthew G Knepley     v0[d] = vertices[d];
353*57459e9aSMatthew G Knepley   }
354*57459e9aSMatthew G Knepley   switch(dim) {
355*57459e9aSMatthew G Knepley   case 2:
356*57459e9aSMatthew G Knepley     for(q = 0; q < quad->numQuadPoints; ++q) {
357*57459e9aSMatthew G Knepley       ierr = DMDAComputeCellGeometry_2D(dm, vertices, &quad->quadPoints[q*dim], J, invJ, detJ);CHKERRQ(ierr);
358*57459e9aSMatthew G Knepley     }
359*57459e9aSMatthew G Knepley     break;
360*57459e9aSMatthew G Knepley   default:
361*57459e9aSMatthew G Knepley     SETERRQ1(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Dimension %d not supported", dim);
362*57459e9aSMatthew G Knepley   }
363*57459e9aSMatthew G Knepley   PetscFunctionReturn(0);
364*57459e9aSMatthew G Knepley }
365