xref: /petsc/src/dm/impls/da/dageometry.c (revision 1a83f524638a4da8317a6bd80eb6d9a2936d8384)
1 #include <petsc-private/daimpl.h>     /*I  "petscdmda.h"   I*/
2 
3 #undef __FUNCT__
4 #define __FUNCT__ "FillClosureArray_Private"
5 PETSC_STATIC_INLINE PetscErrorCode FillClosureArray_Private(DM dm, PetscSection section, PetscInt nP, const PetscInt points[], PetscScalar *vArray, const PetscScalar **array)
6 {
7   PetscScalar   *a;
8   PetscInt       size = 0, dof, off, d, k, i;
9   PetscErrorCode ierr;
10 
11   PetscFunctionBegin;
12   for (i = 0; i < nP; ++i) {
13     ierr = PetscSectionGetDof(section, points[i], &dof);CHKERRQ(ierr);
14     size += dof;
15   }
16   ierr = DMGetWorkArray(dm, size, PETSC_SCALAR, &a);CHKERRQ(ierr);
17   for (i = 0, k = 0; i < nP; ++i) {
18     ierr = PetscSectionGetDof(section, points[i], &dof);CHKERRQ(ierr);
19     ierr = PetscSectionGetOffset(section, points[i], &off);CHKERRQ(ierr);
20 
21     for (d = 0; d < dof; ++d, ++k) {
22       a[k] = vArray[off+d];
23     }
24   }
25   *array = a;
26   PetscFunctionReturn(0);
27 }
28 
29 #undef __FUNCT__
30 #define __FUNCT__ "FillClosureVec_Private"
31 PETSC_STATIC_INLINE PetscErrorCode FillClosureVec_Private(DM dm, PetscSection section, PetscInt nP, const PetscInt points[], PetscScalar *vArray, const PetscScalar *array, InsertMode mode)
32 {
33   PetscInt       dof, off, d, k, i;
34   PetscErrorCode ierr;
35 
36   PetscFunctionBegin;
37   if ((mode == INSERT_VALUES) || (mode == INSERT_ALL_VALUES)) {
38     for (i = 0, k = 0; i < nP; ++i) {
39       ierr = PetscSectionGetDof(section, points[i], &dof);CHKERRQ(ierr);
40       ierr = PetscSectionGetOffset(section, points[i], &off);CHKERRQ(ierr);
41 
42       for (d = 0; d < dof; ++d, ++k) {
43         vArray[off+d] = array[k];
44       }
45     }
46   } else {
47     for (i = 0, k = 0; i < nP; ++i) {
48       ierr = PetscSectionGetDof(section, points[i], &dof);CHKERRQ(ierr);
49       ierr = PetscSectionGetOffset(section, points[i], &off);CHKERRQ(ierr);
50 
51       for (d = 0; d < dof; ++d, ++k) {
52         vArray[off+d] += array[k];
53       }
54     }
55   }
56   PetscFunctionReturn(0);
57 }
58 
59 #undef __FUNCT__
60 #define __FUNCT__ "GetPointArray_Private"
61 PETSC_STATIC_INLINE PetscErrorCode GetPointArray_Private(DM dm,PetscInt n,PetscInt *points,PetscInt *rn,const PetscInt **rpoints)
62 {
63   PetscErrorCode ierr;
64   PetscInt *work;
65 
66   PetscFunctionBegin;
67   if (rn) *rn = n;
68   if (rpoints) {
69     ierr = DMGetWorkArray(dm,n,PETSC_INT,&work);CHKERRQ(ierr);
70     ierr = PetscMemcpy(work,points,n*sizeof(PetscInt));CHKERRQ(ierr);
71     *rpoints = work;
72   }
73   PetscFunctionReturn(0);
74 }
75 
76 #undef __FUNCT__
77 #define __FUNCT__ "RestorePointArray_Private"
78 PETSC_STATIC_INLINE PetscErrorCode RestorePointArray_Private(DM dm,PetscInt *rn,const PetscInt **rpoints)
79 {
80   PetscErrorCode ierr;
81 
82   PetscFunctionBegin;
83   if (rn) *rn = 0;
84   if (rpoints) {ierr = DMRestoreWorkArray(dm,*rn,PETSC_INT,rpoints);CHKERRQ(ierr);}
85   PetscFunctionReturn(0);
86 }
87 
88 #undef __FUNCT__
89 #define __FUNCT__ "DMDAGetClosure"
90 PetscErrorCode DMDAGetClosure(DM dm, PetscSection section, PetscInt p,PetscInt *n,const PetscInt **closure)
91 {
92   DM_DA         *da  = (DM_DA *) dm->data;
93   PetscInt       dim = da->dim;
94   PetscInt       nVx, nVy, nxF, nXF, nyF, nYF, nzF, nZF;
95   PetscInt       pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart;
96   PetscErrorCode ierr;
97 
98   PetscFunctionBegin;
99   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
100   if (n) PetscValidIntPointer(n,4);
101   if (closure) PetscValidPointer(closure, 5);
102   if (!section) {ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);}
103   if (!section) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_ARG_WRONG, "This DM has not default PetscSection");
104   ierr = DMDAGetHeightStratum(dm, -1,  &pStart, &pEnd);CHKERRQ(ierr);
105   ierr = DMDAGetHeightStratum(dm, 0,   &cStart, &cEnd);CHKERRQ(ierr);
106   ierr = DMDAGetHeightStratum(dm, 1,   &fStart, &fEnd);CHKERRQ(ierr);
107   ierr = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr);
108   ierr = DMDAGetNumVertices(dm, &nVx, &nVy, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
109   ierr = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr);
110   xfStart = fStart; xfEnd = xfStart+nXF;
111   yfStart = xfEnd;
112   if ((p < pStart) || (p >= pEnd)) SETERRQ3(((PetscObject) dm)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid point %d should be in [%d, %d)", p, pStart, pEnd);
113   if ((p >= cStart) || (p < cEnd)) {
114     /* Cell */
115     if (dim == 1) {
116       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
117     } else if (dim == 2) {
118       /* 4 faces, 4 vertices
119          Bottom-left vertex follows same order as cells
120          Bottom y-face same order as cells
121          Left x-face follows same order as cells
122          We number the quad:
123 
124            8--3--7
125            |     |
126            4  0  2
127            |     |
128            5--1--6
129       */
130       PetscInt c         = p - cStart, cx = c % (nVx-1), cy = c / (nVx-1);
131       PetscInt v         = cy*nVx + cx +  vStart;
132       PetscInt xf        = cy*nxF + cx + xfStart;
133       PetscInt yf        = c + yfStart;
134       PetscInt points[9] = {p, yf, xf+1, yf+nyF, xf+0, v+0, v+1, v+nVx+1, v+nVx+0};
135 
136       ierr = GetPointArray_Private(dm,9,points,n,closure);CHKERRQ(ierr);
137     } else {
138       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
139     }
140   } else if ((p >= vStart) || (p < vEnd)) {
141     /* Vertex */
142     ierr = GetPointArray_Private(dm,1,&p,n,closure);CHKERRQ(ierr);
143   } else if ((p >= fStart) || (p < fStart + nXF)) {
144     /* X Face */
145     if (dim == 1) {
146       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D");
147     } else if (dim == 2) {
148       /* 2 vertices: The bottom vertex has the same numbering as the face */
149       PetscInt f         = p - xfStart;
150       PetscInt points[3] = {p, f, f+nVx};
151 
152       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Broken");
153       ierr = GetPointArray_Private(dm,3,points,n,closure);CHKERRQ(ierr);
154     } else if (dim == 3) {
155       /* 4 vertices */
156       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
157     }
158   } else if ((p >= fStart + nXF) || (p < fStart + nXF + nYF)) {
159     /* Y Face */
160     if (dim == 1) {
161       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D");
162     } else if (dim == 2) {
163       /* 2 vertices: The left vertex has the same numbering as the face */
164       PetscInt f         = p - yfStart;
165       PetscInt points[3] = {p, f, f+1};
166 
167       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Broken");
168       ierr = GetPointArray_Private(dm, 3, points, n, closure);CHKERRQ(ierr);
169     } else if (dim == 3) {
170       /* 4 vertices */
171       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
172     }
173   } else {
174     /* Z Face */
175     if (dim == 1) {
176       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D");
177     } else if (dim == 2) {
178       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no z-faces in 2D");
179     } else if (dim == 3) {
180       /* 4 vertices */
181       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
182     }
183   }
184   PetscFunctionReturn(0);
185 }
186 
187 #undef __FUNCT__
188 #define __FUNCT__ "DMDARestoreClosure"
189 PetscErrorCode DMDARestoreClosure(DM dm, PetscSection section, PetscInt p,PetscInt *n,const PetscInt **closure)
190 {
191   PetscErrorCode ierr;
192 
193   PetscFunctionBegin;
194   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
195   PetscValidIntPointer(n,4);
196   PetscValidPointer(closure, 5);
197   ierr = RestorePointArray_Private(dm,n,closure);CHKERRQ(ierr);
198   PetscFunctionReturn(0);
199 }
200 
201 #undef __FUNCT__
202 #define __FUNCT__ "DMDAGetClosureScalar"
203 /* If you did not pass PETSC_NULL for 'values', you must call DMDARestoreClosureScalar() */
204 PetscErrorCode DMDAGetClosureScalar(DM dm, PetscSection section,PetscInt p,PetscScalar *vArray,const PetscScalar **values)
205 {
206   DM_DA         *da  = (DM_DA *) dm->data;
207   PetscInt       dim = da->dim;
208   PetscInt       nVx, nVy, nxF, nXF, nyF, nYF, nzF, nZF;
209   PetscInt       pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart, yfEnd, zfStart;
210   PetscErrorCode ierr;
211 
212   PetscFunctionBegin;
213   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
214   PetscValidScalarPointer(vArray, 4);
215   PetscValidPointer(values, 5);
216   if (!section) {ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);}
217   if (!section) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_ARG_WRONG, "This DM has not default PetscSection");
218   ierr = DMDAGetHeightStratum(dm, -1,  &pStart, &pEnd);CHKERRQ(ierr);
219   ierr = DMDAGetHeightStratum(dm, 0,   &cStart, &cEnd);CHKERRQ(ierr);
220   ierr = DMDAGetHeightStratum(dm, 1,   &fStart, &fEnd);CHKERRQ(ierr);
221   ierr = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr);
222   ierr = DMDAGetNumVertices(dm, &nVx, &nVy, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
223   ierr = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr);
224   xfStart = fStart; xfEnd = xfStart+nXF;
225   yfStart = xfEnd;  yfEnd = yfStart+nYF;
226   zfStart = yfEnd;
227   if ((p < pStart) || (p >= pEnd)) SETERRQ3(((PetscObject) dm)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid point %d should be in [%d, %d)", p, pStart, pEnd);
228   if ((p >= cStart) || (p < cEnd)) {
229     /* Cell */
230     if (dim == 1) {
231       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
232     } else if (dim == 2) {
233       /* 4 faces, 4 vertices
234          Bottom-left vertex follows same order as cells
235          Bottom y-face same order as cells
236          Left x-face follows same order as cells
237          We number the quad:
238 
239            8--3--7
240            |     |
241            4  0  2
242            |     |
243            5--1--6
244       */
245       PetscInt c         = p - cStart, cx = c % (nVx-1), cy = c / (nVx-1);
246       PetscInt v         = cy*nVx + cx +  vStart;
247       PetscInt xf        = cy*nxF + cx + xfStart;
248       PetscInt yf        = c + yfStart;
249       PetscInt points[9] = {p, yf, xf+1, yf+nyF, xf+0, v+0, v+1, v+nVx+1, v+nVx+0};
250 
251       ierr = FillClosureArray_Private(dm, section, 9, points, vArray, values);CHKERRQ(ierr);
252     } else {
253       /* 6 faces, 8 vertices
254          Bottom-left-back vertex follows same order as cells
255          Back z-face follows same order as cells
256          Bottom y-face follows same order as cells
257          Left x-face follows same order as cells
258 
259               14-----13
260               /|    /|
261              / | 2 / |
262             / 5|  /  |
263           10-----9  4|
264            |  11-|---12
265            |6 /  |  /
266            | /1 3| /
267            |/    |/
268            7-----8
269       */
270       PetscInt c          = p - cStart;
271       PetscInt points[15] = {p, c+zfStart, c+zfStart+nzF, c+yfStart, c+xfStart+nxF, c+yfStart+nyF, c+xfStart,
272                              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};
273 
274       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Broken");
275       ierr = FillClosureArray_Private(dm, section, 15, points, vArray, values);CHKERRQ(ierr);
276     }
277   } else if ((p >= vStart) || (p < vEnd)) {
278     /* Vertex */
279     ierr = FillClosureArray_Private(dm, section, 1, &p, vArray, values);CHKERRQ(ierr);
280   } else if ((p >= fStart) || (p < fStart + nXF)) {
281     /* X Face */
282     if (dim == 1) {
283       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D");
284     } else if (dim == 2) {
285       /* 2 vertices: The bottom vertex has the same numbering as the face */
286       PetscInt f         = p - xfStart;
287       PetscInt points[3] = {p, f, f+nVx};
288 
289       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Broken");
290       ierr = FillClosureArray_Private(dm, section, 3, points, vArray, values);CHKERRQ(ierr);
291     } else if (dim == 3) {
292       /* 4 vertices */
293       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
294     }
295   } else if ((p >= fStart + nXF) || (p < fStart + nXF + nYF)) {
296     /* Y Face */
297     if (dim == 1) {
298       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D");
299     } else if (dim == 2) {
300       /* 2 vertices: The left vertex has the same numbering as the face */
301       PetscInt f         = p - yfStart;
302       PetscInt points[3] = {p, f, f+1};
303 
304       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Broken");
305       ierr = FillClosureArray_Private(dm, section, 3, points, vArray, values);CHKERRQ(ierr);
306     } else if (dim == 3) {
307       /* 4 vertices */
308       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
309     }
310   } else {
311     /* Z Face */
312     if (dim == 1) {
313       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D");
314     } else if (dim == 2) {
315       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no z-faces in 2D");
316     } else if (dim == 3) {
317       /* 4 vertices */
318       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
319     }
320   }
321   PetscFunctionReturn(0);
322 }
323 
324 #undef __FUNCT__
325 #define __FUNCT__ "DMDAVecGetClosure"
326 PetscErrorCode DMDAVecGetClosure(DM dm, PetscSection section, Vec v, PetscInt p, const PetscScalar **values)
327 {
328   PetscScalar *vArray;
329   PetscErrorCode ierr;
330 
331   PetscFunctionBegin;
332   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
333   PetscValidHeaderSpecific(v, VEC_CLASSID, 3);
334   PetscValidPointer(values, 5);
335   ierr = VecGetArray(v,&vArray);CHKERRQ(ierr);
336   ierr = DMDAGetClosureScalar(dm,section,p,vArray,values);CHKERRQ(ierr);
337   ierr = VecRestoreArray(v,&vArray);CHKERRQ(ierr);
338   PetscFunctionReturn(0);
339 }
340 
341 #undef __FUNCT__
342 #define __FUNCT__ "DMDARestoreClosureScalar"
343 PetscErrorCode DMDARestoreClosureScalar(DM dm, PetscSection section,PetscInt p,PetscScalar *vArray,const PetscScalar **values)
344 {
345   PetscErrorCode ierr;
346   PetscInt count;
347 
348   PetscFunctionBegin;
349   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
350   PetscValidPointer(values, 5);
351   count = 0;                    /* We are lying about the count */
352   ierr = DMRestoreWorkArray(dm,count,PETSC_SCALAR,values);CHKERRQ(ierr);
353   PetscFunctionReturn(0);
354 }
355 
356 #undef __FUNCT__
357 #define __FUNCT__ "DMDAVecRestoreClosure"
358 PetscErrorCode DMDAVecRestoreClosure(DM dm, PetscSection section, Vec v, PetscInt p, const PetscScalar **values)
359 {
360   PetscErrorCode ierr;
361 
362   PetscFunctionBegin;
363   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
364   PetscValidHeaderSpecific(v, VEC_CLASSID, 3);
365   PetscValidPointer(values, 5);
366   ierr = DMDARestoreClosureScalar(dm,section,p,PETSC_NULL,values);CHKERRQ(ierr);
367   PetscFunctionReturn(0);
368 }
369 
370 #undef __FUNCT__
371 #define __FUNCT__ "DMDASetClosureScalar"
372 PetscErrorCode DMDASetClosureScalar(DM dm, PetscSection section, PetscInt p,PetscScalar *vArray, const PetscScalar *values, InsertMode mode)
373 {
374   DM_DA         *da  = (DM_DA *) dm->data;
375   PetscInt       dim = da->dim;
376   PetscInt       nVx, nVy, nxF, nXF, nyF, nYF, nzF, nZF;
377   PetscInt       pStart, pEnd, cStart, cEnd, vStart, vEnd, fStart, fEnd, xfStart, xfEnd, yfStart, yfEnd, zfStart;
378   PetscErrorCode ierr;
379 
380   PetscFunctionBegin;
381   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
382   PetscValidScalarPointer(values, 4);
383   PetscValidPointer(values, 5);
384   if (!section) {ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);}
385   if (!section) SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_ARG_WRONG, "This DM has not default PetscSection");
386   ierr = DMDAGetHeightStratum(dm, -1,  &pStart, &pEnd);CHKERRQ(ierr);
387   ierr = DMDAGetHeightStratum(dm, 0,   &cStart, &cEnd);CHKERRQ(ierr);
388   ierr = DMDAGetHeightStratum(dm, 1,   &fStart, &fEnd);CHKERRQ(ierr);
389   ierr = DMDAGetHeightStratum(dm, dim, &vStart, &vEnd);CHKERRQ(ierr);
390   ierr = DMDAGetNumVertices(dm, &nVx, &nVy, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
391   ierr = DMDAGetNumFaces(dm, &nxF, &nXF, &nyF, &nYF, &nzF, &nZF);CHKERRQ(ierr);
392   xfStart = fStart; xfEnd = xfStart+nXF;
393   yfStart = xfEnd;  yfEnd = yfStart+nYF;
394   zfStart = yfEnd;
395   if ((p < pStart) || (p >= pEnd)) SETERRQ3(((PetscObject) dm)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid point %d should be in [%d, %d)", p, pStart, pEnd);
396   if ((p >= cStart) || (p < cEnd)) {
397     /* Cell */
398     if (dim == 1) {
399       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
400     } else if (dim == 2) {
401       /* 4 faces, 4 vertices
402          Bottom-left vertex follows same order as cells
403          Bottom y-face same order as cells
404          Left x-face follows same order as cells
405          We number the quad:
406 
407            8--3--7
408            |     |
409            4  0  2
410            |     |
411            5--1--6
412       */
413       PetscInt c         = p - cStart;
414       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};
415 
416       ierr = FillClosureVec_Private(dm, section, 9, points, vArray, values, mode);CHKERRQ(ierr);
417     } else {
418       /* 6 faces, 8 vertices
419          Bottom-left-back vertex follows same order as cells
420          Back z-face follows same order as cells
421          Bottom y-face follows same order as cells
422          Left x-face follows same order as cells
423 
424               14-----13
425               /|    /|
426              / | 2 / |
427             / 5|  /  |
428           10-----9  4|
429            |  11-|---12
430            |6 /  |  /
431            | /1 3| /
432            |/    |/
433            7-----8
434       */
435       PetscInt c          = p - cStart;
436       PetscInt points[15] = {p, c+zfStart, c+zfStart+nzF, c+yfStart, c+xfStart+nxF, c+yfStart+nyF, c+xfStart,
437                              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};
438 
439       ierr = FillClosureVec_Private(dm, section, 15, points, vArray, values, mode);CHKERRQ(ierr);
440     }
441   } else if ((p >= vStart) || (p < vEnd)) {
442     /* Vertex */
443     ierr = FillClosureVec_Private(dm, section, 1, &p, vArray, values, mode);CHKERRQ(ierr);
444   } else if ((p >= fStart) || (p < fStart + nXF)) {
445     /* X Face */
446     if (dim == 1) {
447       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D");
448     } else if (dim == 2) {
449       /* 2 vertices: The bottom vertex has the same numbering as the face */
450       PetscInt f         = p - xfStart;
451       PetscInt points[3] = {p, f, f+nVx};
452 
453       ierr = FillClosureVec_Private(dm, section, 3, points, vArray, values, mode);CHKERRQ(ierr);
454     } else if (dim == 3) {
455       /* 4 vertices */
456       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
457     }
458   } else if ((p >= fStart + nXF) || (p < fStart + nXF + nYF)) {
459     /* Y Face */
460     if (dim == 1) {
461       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D");
462     } else if (dim == 2) {
463       /* 2 vertices: The left vertex has the same numbering as the face */
464       PetscInt f         = p - yfStart;
465       PetscInt points[3] = {p, f, f+1};
466 
467       ierr = FillClosureVec_Private(dm, section, 3, points, vArray, values, mode);CHKERRQ(ierr);
468     } else if (dim == 3) {
469       /* 4 vertices */
470       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
471     }
472   } else {
473     /* Z Face */
474     if (dim == 1) {
475       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no faces in 1D");
476     } else if (dim == 2) {
477       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_PLIB, "There are no z-faces in 2D");
478     } else if (dim == 3) {
479       /* 4 vertices */
480       SETERRQ(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Not implemented");
481     }
482   }
483   PetscFunctionReturn(0);
484 }
485 
486 #undef __FUNCT__
487 #define __FUNCT__ "DMDAVecSetClosure"
488 PetscErrorCode DMDAVecSetClosure(DM dm, PetscSection section, Vec v, PetscInt p, const PetscScalar *values, InsertMode mode)
489 {
490   PetscScalar   *vArray;
491   PetscErrorCode ierr;
492 
493   PetscFunctionBegin;
494   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
495   PetscValidHeaderSpecific(v, VEC_CLASSID, 3);
496   PetscValidPointer(values, 5);
497   ierr = VecGetArray(v,&vArray);CHKERRQ(ierr);
498   ierr = DMDASetClosureScalar(dm,section,p,vArray,values,mode);CHKERRQ(ierr);
499   ierr = VecRestoreArray(v,&vArray);CHKERRQ(ierr);
500   PetscFunctionReturn(0);
501 }
502 
503 #undef __FUNCT__
504 #define __FUNCT__ "DMDACnvertToCell"
505 /*@
506   DMDAConvertToCell - Convert (i,j,k) to local cell number
507 
508   Not Collective
509 
510   Input Parameter:
511 + da - the distributed array
512 = s - A MatStencil giving (i,j,k)
513 
514   Output Parameter:
515 . cell - the local cell number
516 
517   Level: developer
518 
519 .seealso: DMDAVecGetClosure()
520 @*/
521 PetscErrorCode DMDAConvertToCell(DM dm, MatStencil s, PetscInt *cell)
522 {
523   DM_DA         *da  = (DM_DA *) dm->data;
524   const PetscInt dim = da->dim;
525   const PetscInt mx  = (da->Xe - da->Xs)/da->w, my = da->Ye - da->Ys/*, mz = da->Ze - da->Zs*/;
526   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;
527 
528   PetscFunctionBegin;
529   *cell = -1;
530   if ((s.i < da->Xs/da->w) || (s.i >= da->Xe/da->w))    SETERRQ3(((PetscObject) dm)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Stencil i %d should be in [%d, %d)", s.i, da->Xs/da->w, da->Xe/da->w);
531   if ((dim > 1) && ((s.j < da->Ys) || (s.j >= da->Ye))) SETERRQ3(((PetscObject) dm)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Stencil j %d should be in [%d, %d)", s.j, da->Ys, da->Ye);
532   if ((dim > 2) && ((s.k < da->Zs) || (s.k >= da->Ze))) SETERRQ3(((PetscObject) dm)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Stencil k %d should be in [%d, %d)", s.k, da->Zs, da->Ze);
533   *cell = (kl*my + jl)*mx + il;
534   PetscFunctionReturn(0);
535 }
536 
537 #undef __FUNCT__
538 #define __FUNCT__ "DMDAComputeCellGeometry_2D"
539 PetscErrorCode DMDAComputeCellGeometry_2D(DM dm, const PetscScalar vertices[], const PetscReal refPoint[], PetscReal J[], PetscReal invJ[], PetscReal *detJ)
540 {
541   const PetscScalar x0   = vertices[0];
542   const PetscScalar y0   = vertices[1];
543   const PetscScalar x1   = vertices[2];
544   const PetscScalar y1   = vertices[3];
545   const PetscScalar x2   = vertices[4];
546   const PetscScalar y2   = vertices[5];
547   const PetscScalar x3   = vertices[6];
548   const PetscScalar y3   = vertices[7];
549   const PetscScalar f_01 = x2 - x1 - x3 + x0;
550   const PetscScalar g_01 = y2 - y1 - y3 + y0;
551   const PetscScalar x    = refPoint[0];
552   const PetscScalar y    = refPoint[1];
553   PetscReal         invDet;
554   PetscErrorCode    ierr;
555 
556   PetscFunctionBegin;
557 #if defined(PETSC_USE_DEBUG)
558   ierr = PetscPrintf(PETSC_COMM_SELF, "Cell (%g,%g)--(%g,%g)--(%g,%g)--(%g,%g)\n",
559                      PetscRealPart(x0),PetscRealPart(y0),PetscRealPart(x1),PetscRealPart(y1),PetscRealPart(x2),PetscRealPart(y2),PetscRealPart(x3),PetscRealPart(y3));CHKERRQ(ierr);
560   ierr = PetscPrintf(PETSC_COMM_SELF, "Ref Point (%g,%g)\n", PetscRealPart(x), PetscRealPart(y));CHKERRQ(ierr);
561 #endif
562   J[0] = PetscRealPart(x1 - x0 + f_01*y) * 0.5; J[1] = PetscRealPart(x3 - x0 + f_01*x) * 0.5;
563   J[2] = PetscRealPart(y1 - y0 + g_01*y) * 0.5; J[3] = PetscRealPart(y3 - y0 + g_01*x) * 0.5;
564   *detJ   = J[0]*J[3] - J[1]*J[2];
565   invDet  = 1.0/(*detJ);
566   invJ[0] =  invDet*J[3]; invJ[1] = -invDet*J[1];
567   invJ[2] = -invDet*J[2]; invJ[3] =  invDet*J[0];
568   ierr = PetscLogFlops(30);CHKERRQ(ierr);
569   PetscFunctionReturn(0);
570 }
571 
572 #undef __FUNCT__
573 #define __FUNCT__ "DMDAComputeCellGeometry"
574 PetscErrorCode DMDAComputeCellGeometry(DM dm, PetscInt cell, PetscQuadrature *quad, PetscReal v0[], PetscReal J[], PetscReal invJ[], PetscReal detJ[])
575 {
576   DM                 cdm;
577   Vec                coordinates;
578   const PetscScalar *vertices;
579   PetscInt           dim, d, q;
580   PetscErrorCode     ierr;
581 
582   PetscFunctionBegin;
583   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
584   ierr = DMDAGetInfo(dm, &dim, 0,0,0,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
585   ierr = DMGetCoordinates(dm, &coordinates);CHKERRQ(ierr);
586   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
587   ierr = DMDAVecGetClosure(cdm, PETSC_NULL, coordinates, cell, &vertices);CHKERRQ(ierr);
588   for (d = 0; d < dim; ++d) {
589     v0[d] = PetscRealPart(vertices[d]);
590   }
591   switch(dim) {
592   case 2:
593     for (q = 0; q < quad->numQuadPoints; ++q) {
594       ierr = DMDAComputeCellGeometry_2D(dm, vertices, &quad->quadPoints[q*dim], J, invJ, detJ);CHKERRQ(ierr);
595     }
596     break;
597   default:
598     SETERRQ1(((PetscObject) dm)->comm, PETSC_ERR_SUP, "Dimension %d not supported", dim);
599   }
600   PetscFunctionReturn(0);
601 }
602