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