xref: /petsc/src/dm/field/impls/ds/dmfieldds.c (revision ccb4e88a40f0b86eaeca07ff64c64e4de2fae686)
1 #include <petsc/private/dmfieldimpl.h> /*I "petscdmfield.h" I*/
2 #include <petsc/private/petscfeimpl.h> /*I "petscdmfield.h" I*/
3 #include <petscfe.h>
4 #include <petscdmplex.h>
5 #include <petscds.h>
6 
7 typedef struct _n_DMField_DS
8 {
9   PetscInt    fieldNum;
10   Vec         vec;
11   PetscInt    height;
12   PetscObject *disc;
13   PetscBool   multifieldVec;
14 }
15 DMField_DS;
16 
17 static PetscErrorCode DMFieldDestroy_DS(DMField field)
18 {
19   DMField_DS     *dsfield;
20   PetscInt       i;
21   PetscErrorCode ierr;
22 
23   PetscFunctionBegin;
24   dsfield = (DMField_DS *) field->data;
25   ierr = VecDestroy(&dsfield->vec);CHKERRQ(ierr);
26   for (i = 0; i < dsfield->height; i++) {
27     ierr = PetscObjectDereference(dsfield->disc[i]);CHKERRQ(ierr);
28   }
29   ierr = PetscFree(dsfield->disc);CHKERRQ(ierr);
30   ierr = PetscFree(dsfield);CHKERRQ(ierr);
31   PetscFunctionReturn(0);
32 }
33 
34 static PetscErrorCode DMFieldView_DS(DMField field,PetscViewer viewer)
35 {
36   DMField_DS     *dsfield = (DMField_DS *) field->data;
37   PetscBool      iascii;
38   PetscObject    disc;
39   PetscErrorCode ierr;
40 
41   PetscFunctionBegin;
42   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
43   disc = dsfield->disc[0];
44   if (iascii) {
45     ierr = PetscViewerASCIIPrintf(viewer, "PetscDS field %D\n",dsfield->fieldNum);CHKERRQ(ierr);
46     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
47     ierr = PetscObjectView(disc,viewer);CHKERRQ(ierr);
48     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
49   }
50   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
51   if (dsfield->multifieldVec) SETERRQ(PetscObjectComm((PetscObject)field),PETSC_ERR_SUP,"View of subfield not implemented yet");
52   ierr = VecView(dsfield->vec,viewer);CHKERRQ(ierr);
53   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
54   PetscFunctionReturn(0);
55 }
56 
57 static PetscErrorCode DMFieldDSGetHeightDisc(DMField field, PetscInt height, PetscObject *disc)
58 {
59   DMField_DS     *dsfield = (DMField_DS *) field->data;
60   PetscErrorCode ierr;
61 
62   PetscFunctionBegin;
63   if (!dsfield->disc[height]) {
64     PetscClassId   id;
65 
66     ierr = PetscObjectGetClassId(dsfield->disc[0],&id);CHKERRQ(ierr);
67     if (id == PETSCFE_CLASSID) {
68       PetscFE fe = (PetscFE) dsfield->disc[0];
69 
70       ierr = PetscFECreateHeightTrace(fe,height,(PetscFE *)&dsfield->disc[height]);CHKERRQ(ierr);
71     }
72   }
73   *disc = dsfield->disc[height];
74   PetscFunctionReturn(0);
75 }
76 
77 /* y[m,c] = A[m,n,c] . b[n] */
78 #define DMFieldDSdot(y,A,b,m,n,c,cast)                                           \
79   do {                                                                           \
80     PetscInt _i, _j, _k;                                                         \
81     for (_i = 0; _i < (m); _i++) {                                               \
82       for (_k = 0; _k < (c); _k++) {                                             \
83         (y)[_i * (c) + _k] = 0.;                                                 \
84       }                                                                          \
85       for (_j = 0; _j < (n); _j++) {                                             \
86         for (_k = 0; _k < (c); _k++) {                                           \
87           (y)[_i * (c) + _k] += (A)[(_i * (n) + _j) * (c) + _k] * cast((b)[_j]); \
88         }                                                                        \
89       }                                                                          \
90     }                                                                            \
91   } while (0)
92 
93 /* TODO: Reorganize interface so that I can reuse a tabulation rather than mallocing each time */
94 static PetscErrorCode DMFieldEvaluateFE_DS(DMField field, IS pointIS, PetscQuadrature quad, PetscDataType type, void *B, void *D, void *H)
95 {
96   DMField_DS      *dsfield = (DMField_DS *) field->data;
97   DM              dm;
98   PetscObject     disc;
99   PetscClassId    classid;
100   PetscInt        nq, nc, dim, meshDim, numCells;
101   PetscSection    section;
102   const PetscReal *qpoints;
103   PetscBool       isStride;
104   const PetscInt  *points = NULL;
105   PetscInt        sfirst = -1, stride = -1;
106   PetscErrorCode  ierr;
107 
108   PetscFunctionBeginHot;
109   dm   = field->dm;
110   nc   = field->numComponents;
111   ierr = PetscQuadratureGetData(quad,&dim,NULL,&nq,&qpoints,NULL);CHKERRQ(ierr);
112   ierr = DMFieldDSGetHeightDisc(field,dsfield->height - 1 - dim,&disc);CHKERRQ(ierr);
113   ierr = DMGetDimension(dm,&meshDim);CHKERRQ(ierr);
114   ierr = DMGetLocalSection(dm,&section);CHKERRQ(ierr);
115   ierr = PetscSectionGetField(section,dsfield->fieldNum,&section);CHKERRQ(ierr);
116   ierr = PetscObjectGetClassId(disc,&classid);CHKERRQ(ierr);
117   /* TODO: batch */
118   ierr = PetscObjectTypeCompare((PetscObject)pointIS,ISSTRIDE,&isStride);CHKERRQ(ierr);
119   ierr = ISGetLocalSize(pointIS,&numCells);CHKERRQ(ierr);
120   if (isStride) {
121     ierr = ISStrideGetInfo(pointIS,&sfirst,&stride);CHKERRQ(ierr);
122   } else {
123     ierr = ISGetIndices(pointIS,&points);CHKERRQ(ierr);
124   }
125   if (classid == PETSCFE_CLASSID) {
126     PetscFE      fe = (PetscFE) disc;
127     PetscInt     feDim, i;
128     PetscInt          K = H ? 2 : (D ? 1 : (B ? 0 : -1));
129     PetscTabulation   T;
130 
131     ierr = PetscFEGetDimension(fe,&feDim);CHKERRQ(ierr);
132     ierr = PetscFECreateTabulation(fe,1,nq,qpoints,K,&T);CHKERRQ(ierr);
133     for (i = 0; i < numCells; i++) {
134       PetscInt     c = isStride ? (sfirst + i * stride) : points[i];
135       PetscInt     closureSize;
136       PetscScalar *elem = NULL;
137 
138       ierr = DMPlexVecGetClosure(dm,section,dsfield->vec,c,&closureSize,&elem);CHKERRQ(ierr);
139       if (B) {
140         /* field[c] = T[q,b,c] . coef[b], so v[c] = T[q,b,c] . coords[b] */
141         if (type == PETSC_SCALAR) {
142           PetscScalar *cB = &((PetscScalar *) B)[nc * nq * i];
143 
144           DMFieldDSdot(cB,T->T[0],elem,nq,feDim,nc,(PetscScalar));
145         } else {
146           PetscReal *cB = &((PetscReal *) B)[nc * nq * i];
147 
148           DMFieldDSdot(cB,T->T[0],elem,nq,feDim,nc,PetscRealPart);
149         }
150       }
151       if (D) {
152         if (type == PETSC_SCALAR) {
153           PetscScalar *cD = &((PetscScalar *) D)[nc * nq * dim * i];
154 
155           DMFieldDSdot(cD,T->T[1],elem,nq,feDim,(nc * dim),(PetscScalar));
156         } else {
157           PetscReal *cD = &((PetscReal *) D)[nc * nq * dim * i];
158 
159           DMFieldDSdot(cD,T->T[1],elem,nq,feDim,(nc * dim),PetscRealPart);
160         }
161       }
162       if (H) {
163         if (type == PETSC_SCALAR) {
164           PetscScalar *cH = &((PetscScalar *) H)[nc * nq * dim * dim * i];
165 
166           DMFieldDSdot(cH,T->T[2],elem,nq,feDim,(nc * dim * dim),(PetscScalar));
167         } else {
168           PetscReal *cH = &((PetscReal *) H)[nc * nq * dim * dim * i];
169 
170           DMFieldDSdot(cH,T->T[2],elem,nq,feDim,(nc * dim * dim),PetscRealPart);
171         }
172       }
173       ierr = DMPlexVecRestoreClosure(dm,section,dsfield->vec,c,&closureSize,&elem);CHKERRQ(ierr);
174     }
175     ierr = PetscTabulationDestroy(&T);CHKERRQ(ierr);
176   } else SETERRQ(PetscObjectComm((PetscObject)field),PETSC_ERR_SUP,"Not implemented");
177   if (!isStride) {
178     ierr = ISRestoreIndices(pointIS,&points);CHKERRQ(ierr);
179   }
180   PetscFunctionReturn(0);
181 }
182 
183 static PetscErrorCode DMFieldEvaluate_DS(DMField field, Vec points, PetscDataType datatype, void *B, void *D, void *H)
184 {
185   DMField_DS        *dsfield = (DMField_DS *) field->data;
186   PetscSF            cellSF = NULL;
187   const PetscSFNode *cells;
188   PetscInt           c, nFound, numCells, feDim, nc;
189   const PetscInt    *cellDegrees;
190   const PetscScalar *pointsArray;
191   PetscScalar       *cellPoints;
192   PetscInt           gatherSize, gatherMax;
193   PetscInt           dim, dimR, offset;
194   MPI_Datatype       pointType;
195   PetscObject        cellDisc;
196   PetscFE            cellFE;
197   PetscClassId       discID;
198   PetscReal         *coordsReal, *coordsRef;
199   PetscSection       section;
200   PetscScalar       *cellBs = NULL, *cellDs = NULL, *cellHs = NULL;
201   PetscReal         *cellBr = NULL, *cellDr = NULL, *cellHr = NULL;
202   PetscReal         *v, *J, *invJ, *detJ;
203   PetscErrorCode     ierr;
204 
205   PetscFunctionBegin;
206   nc   = field->numComponents;
207   ierr = DMGetLocalSection(field->dm,&section);CHKERRQ(ierr);
208   ierr = DMFieldDSGetHeightDisc(field,0,&cellDisc);CHKERRQ(ierr);
209   ierr = PetscObjectGetClassId(cellDisc, &discID);CHKERRQ(ierr);
210   if (discID != PETSCFE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB, "Discretization type not supported\n");
211   cellFE = (PetscFE) cellDisc;
212   ierr = PetscFEGetDimension(cellFE,&feDim);CHKERRQ(ierr);
213   ierr = DMGetCoordinateDim(field->dm, &dim);CHKERRQ(ierr);
214   ierr = DMGetDimension(field->dm, &dimR);CHKERRQ(ierr);
215   ierr = DMLocatePoints(field->dm, points, DM_POINTLOCATION_NONE, &cellSF);CHKERRQ(ierr);
216   ierr = PetscSFGetGraph(cellSF, &numCells, &nFound, NULL, &cells);CHKERRQ(ierr);
217   for (c = 0; c < nFound; c++) {
218     if (cells[c].index < 0) SETERRQ1(PetscObjectComm((PetscObject)points),PETSC_ERR_ARG_WRONG, "Point %D could not be located\n", c);
219   }
220   ierr = PetscSFComputeDegreeBegin(cellSF,&cellDegrees);CHKERRQ(ierr);
221   ierr = PetscSFComputeDegreeEnd(cellSF,&cellDegrees);CHKERRQ(ierr);
222   for (c = 0, gatherSize = 0, gatherMax = 0; c < numCells; c++) {
223     gatherMax = PetscMax(gatherMax,cellDegrees[c]);
224     gatherSize += cellDegrees[c];
225   }
226   ierr = PetscMalloc3(gatherSize*dim,&cellPoints,gatherMax*dim,&coordsReal,gatherMax*dimR,&coordsRef);CHKERRQ(ierr);
227   ierr = PetscMalloc4(gatherMax*dimR,&v,gatherMax*dimR*dimR,&J,gatherMax*dimR*dimR,&invJ,gatherMax,&detJ);CHKERRQ(ierr);
228   if (datatype == PETSC_SCALAR) {
229     ierr = PetscMalloc3(B ? nc * gatherSize : 0, &cellBs, D ? nc * dim * gatherSize : 0, &cellDs, H ? nc * dim * dim * gatherSize : 0, &cellHs);CHKERRQ(ierr);
230   } else {
231     ierr = PetscMalloc3(B ? nc * gatherSize : 0, &cellBr, D ? nc * dim * gatherSize : 0, &cellDr, H ? nc * dim * dim * gatherSize : 0, &cellHr);CHKERRQ(ierr);
232   }
233 
234   ierr = MPI_Type_contiguous(dim,MPIU_SCALAR,&pointType);CHKERRMPI(ierr);
235   ierr = MPI_Type_commit(&pointType);CHKERRMPI(ierr);
236   ierr = VecGetArrayRead(points,&pointsArray);CHKERRQ(ierr);
237   ierr = PetscSFGatherBegin(cellSF, pointType, pointsArray, cellPoints);CHKERRQ(ierr);
238   ierr = PetscSFGatherEnd(cellSF, pointType, pointsArray, cellPoints);CHKERRQ(ierr);
239   ierr = VecRestoreArrayRead(points,&pointsArray);CHKERRQ(ierr);
240   for (c = 0, offset = 0; c < numCells; c++) {
241     PetscInt nq = cellDegrees[c], p;
242 
243     if (nq) {
244       PetscInt          K = H ? 2 : (D ? 1 : (B ? 0 : -1));
245       PetscTabulation   T;
246       PetscInt     closureSize;
247       PetscScalar *elem = NULL;
248       PetscReal   *quadPoints;
249       PetscQuadrature quad;
250       PetscInt d, e, f, g;
251 
252       for (p = 0; p < dim * nq; p++) coordsReal[p] = PetscRealPart(cellPoints[dim * offset + p]);
253       ierr = DMPlexCoordinatesToReference(field->dm, c, nq, coordsReal, coordsRef);CHKERRQ(ierr);
254       ierr = PetscFECreateTabulation(cellFE,1,nq,coordsRef,K,&T);CHKERRQ(ierr);
255       ierr = PetscQuadratureCreate(PETSC_COMM_SELF, &quad);CHKERRQ(ierr);
256       ierr = PetscMalloc1(dimR * nq, &quadPoints);CHKERRQ(ierr);
257       for (p = 0; p < dimR * nq; p++) quadPoints[p] = coordsRef[p];
258       ierr = PetscQuadratureSetData(quad, dimR, 0, nq, quadPoints, NULL);CHKERRQ(ierr);
259       ierr = DMPlexComputeCellGeometryFEM(field->dm, c, quad, v, J, invJ, detJ);CHKERRQ(ierr);
260       ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr);
261       ierr = DMPlexVecGetClosure(field->dm,section,dsfield->vec,c,&closureSize,&elem);CHKERRQ(ierr);
262       if (B) {
263         if (datatype == PETSC_SCALAR) {
264           PetscScalar *cB = &cellBs[nc * offset];
265 
266           DMFieldDSdot(cB,T->T[0],elem,nq,feDim,nc,(PetscScalar));
267         } else {
268           PetscReal *cB = &cellBr[nc * offset];
269 
270           DMFieldDSdot(cB,T->T[0],elem,nq,feDim,nc,PetscRealPart);
271         }
272       }
273       if (D) {
274         if (datatype == PETSC_SCALAR) {
275           PetscScalar *cD = &cellDs[nc * dim * offset];
276 
277           DMFieldDSdot(cD,T->T[1],elem,nq,feDim,(nc * dim),(PetscScalar));
278           for (p = 0; p < nq; p++) {
279             for (g = 0; g < nc; g++) {
280               PetscScalar vs[3];
281 
282               for (d = 0; d < dimR; d++) {
283                 vs[d] = 0.;
284                 for (e = 0; e < dimR; e++) {
285                   vs[d] += invJ[dimR * dimR * p + e * dimR + d] * cD[(nc * p + g) * dimR + e];
286                 }
287               }
288               for (d = 0; d < dimR; d++) {
289                 cD[(nc * p + g) * dimR + d] = vs[d];
290               }
291             }
292           }
293         } else {
294           PetscReal *cD = &cellDr[nc * dim * offset];
295 
296           DMFieldDSdot(cD,T->T[1],elem,nq,feDim,(nc * dim),PetscRealPart);
297           for (p = 0; p < nq; p++) {
298             for (g = 0; g < nc; g++) {
299               for (d = 0; d < dimR; d++) {
300                 v[d] = 0.;
301                 for (e = 0; e < dimR; e++) {
302                   v[d] += invJ[dimR * dimR * p + e * dimR + d] * cD[(nc * p + g) * dimR + e];
303                 }
304               }
305               for (d = 0; d < dimR; d++) {
306                 cD[(nc * p + g) * dimR + d] = v[d];
307               }
308             }
309           }
310         }
311       }
312       if (H) {
313         if (datatype == PETSC_SCALAR) {
314           PetscScalar *cH = &cellHs[nc * dim * dim * offset];
315 
316           DMFieldDSdot(cH,T->T[2],elem,nq,feDim,(nc * dim * dim),(PetscScalar));
317           for (p = 0; p < nq; p++) {
318             for (g = 0; g < nc * dimR; g++) {
319               PetscScalar vs[3];
320 
321               for (d = 0; d < dimR; d++) {
322                 vs[d] = 0.;
323                 for (e = 0; e < dimR; e++) {
324                   vs[d] += invJ[dimR * dimR * p + e * dimR + d] * cH[(nc * dimR * p + g) * dimR + e];
325                 }
326               }
327               for (d = 0; d < dimR; d++) {
328                 cH[(nc * dimR * p + g) * dimR + d] = vs[d];
329               }
330             }
331             for (g = 0; g < nc; g++) {
332               for (f = 0; f < dimR; f++) {
333                 PetscScalar vs[3];
334 
335                 for (d = 0; d < dimR; d++) {
336                   vs[d] = 0.;
337                   for (e = 0; e < dimR; e++) {
338                     vs[d] += invJ[dimR * dimR * p + e * dimR + d] * cH[((nc * p + g) * dimR + e) * dimR + f];
339                   }
340                 }
341                 for (d = 0; d < dimR; d++) {
342                   cH[((nc * p + g) * dimR + d) * dimR + f] = vs[d];
343                 }
344               }
345             }
346           }
347         } else {
348           PetscReal *cH = &cellHr[nc * dim * dim * offset];
349 
350           DMFieldDSdot(cH,T->T[2],elem,nq,feDim,(nc * dim * dim),PetscRealPart);
351           for (p = 0; p < nq; p++) {
352             for (g = 0; g < nc * dimR; g++) {
353               for (d = 0; d < dimR; d++) {
354                 v[d] = 0.;
355                 for (e = 0; e < dimR; e++) {
356                   v[d] += invJ[dimR * dimR * p + e * dimR + d] * cH[(nc * dimR * p + g) * dimR + e];
357                 }
358               }
359               for (d = 0; d < dimR; d++) {
360                 cH[(nc * dimR * p + g) * dimR + d] = v[d];
361               }
362             }
363             for (g = 0; g < nc; g++) {
364               for (f = 0; f < dimR; f++) {
365                 for (d = 0; d < dimR; d++) {
366                   v[d] = 0.;
367                   for (e = 0; e < dimR; e++) {
368                     v[d] += invJ[dimR * dimR * p + e * dimR + d] * cH[((nc * p + g) * dimR + e) * dimR + f];
369                   }
370                 }
371                 for (d = 0; d < dimR; d++) {
372                   cH[((nc * p + g) * dimR + d) * dimR + f] = v[d];
373                 }
374               }
375             }
376           }
377         }
378       }
379       ierr = DMPlexVecRestoreClosure(field->dm,section,dsfield->vec,c,&closureSize,&elem);CHKERRQ(ierr);
380       ierr = PetscTabulationDestroy(&T);CHKERRQ(ierr);
381     }
382     offset += nq;
383   }
384   {
385     MPI_Datatype origtype;
386 
387     if (datatype == PETSC_SCALAR) {
388       origtype = MPIU_SCALAR;
389     } else {
390       origtype = MPIU_REAL;
391     }
392     if (B) {
393       MPI_Datatype Btype;
394 
395       ierr = MPI_Type_contiguous(nc, origtype, &Btype);CHKERRMPI(ierr);
396       ierr = MPI_Type_commit(&Btype);CHKERRMPI(ierr);
397       ierr = PetscSFScatterBegin(cellSF,Btype,(datatype == PETSC_SCALAR) ? (void *) cellBs : (void *) cellBr, B);CHKERRQ(ierr);
398       ierr = PetscSFScatterEnd(cellSF,Btype,(datatype == PETSC_SCALAR) ? (void *) cellBs : (void *) cellBr, B);CHKERRQ(ierr);
399       ierr = MPI_Type_free(&Btype);CHKERRMPI(ierr);
400     }
401     if (D) {
402       MPI_Datatype Dtype;
403 
404       ierr = MPI_Type_contiguous(nc * dim, origtype, &Dtype);CHKERRMPI(ierr);
405       ierr = MPI_Type_commit(&Dtype);CHKERRMPI(ierr);
406       ierr = PetscSFScatterBegin(cellSF,Dtype,(datatype == PETSC_SCALAR) ? (void *) cellDs : (void *) cellDr, D);CHKERRQ(ierr);
407       ierr = PetscSFScatterEnd(cellSF,Dtype,(datatype == PETSC_SCALAR) ? (void *) cellDs : (void *) cellDr, D);CHKERRQ(ierr);
408       ierr = MPI_Type_free(&Dtype);CHKERRMPI(ierr);
409     }
410     if (H) {
411       MPI_Datatype Htype;
412 
413       ierr = MPI_Type_contiguous(nc * dim * dim, origtype, &Htype);CHKERRMPI(ierr);
414       ierr = MPI_Type_commit(&Htype);CHKERRMPI(ierr);
415       ierr = PetscSFScatterBegin(cellSF,Htype,(datatype == PETSC_SCALAR) ? (void *) cellHs : (void *) cellHr, H);CHKERRQ(ierr);
416       ierr = PetscSFScatterEnd(cellSF,Htype,(datatype == PETSC_SCALAR) ? (void *) cellHs : (void *) cellHr, H);CHKERRQ(ierr);
417       ierr = MPI_Type_free(&Htype);CHKERRMPI(ierr);
418     }
419   }
420   ierr = PetscFree4(v,J,invJ,detJ);CHKERRQ(ierr);
421   ierr = PetscFree3(cellBr, cellDr, cellHr);CHKERRQ(ierr);
422   ierr = PetscFree3(cellBs, cellDs, cellHs);CHKERRQ(ierr);
423   ierr = PetscFree3(cellPoints,coordsReal,coordsRef);CHKERRQ(ierr);
424   ierr = MPI_Type_free(&pointType);CHKERRMPI(ierr);
425   ierr = PetscSFDestroy(&cellSF);CHKERRQ(ierr);
426   PetscFunctionReturn(0);
427 }
428 
429 static PetscErrorCode DMFieldEvaluateFV_DS(DMField field, IS pointIS, PetscDataType type, void *B, void *D, void *H)
430 {
431   DMField_DS      *dsfield = (DMField_DS *) field->data;
432   PetscInt         h, imin;
433   PetscInt         dim;
434   PetscClassId     id;
435   PetscQuadrature  quad = NULL;
436   PetscInt         maxDegree;
437   PetscFEGeom      *geom;
438   PetscInt         Nq, Nc, dimC, qNc, N;
439   PetscInt         numPoints;
440   void            *qB = NULL, *qD = NULL, *qH = NULL;
441   const PetscReal *weights;
442   MPI_Datatype     mpitype = type == PETSC_SCALAR ? MPIU_SCALAR : MPIU_REAL;
443   PetscObject      disc;
444   DMField          coordField;
445   PetscErrorCode   ierr;
446 
447   PetscFunctionBegin;
448   Nc = field->numComponents;
449   ierr = DMGetCoordinateDim(field->dm, &dimC);CHKERRQ(ierr);
450   ierr = DMGetDimension(field->dm, &dim);CHKERRQ(ierr);
451   ierr = ISGetLocalSize(pointIS, &numPoints);CHKERRQ(ierr);
452   ierr = ISGetMinMax(pointIS,&imin,NULL);CHKERRQ(ierr);
453   for (h = 0; h < dsfield->height; h++) {
454     PetscInt hEnd;
455 
456     ierr = DMPlexGetHeightStratum(field->dm,h,NULL,&hEnd);CHKERRQ(ierr);
457     if (imin < hEnd) break;
458   }
459   dim -= h;
460   ierr = DMFieldDSGetHeightDisc(field,h,&disc);CHKERRQ(ierr);
461   ierr = PetscObjectGetClassId(disc,&id);CHKERRQ(ierr);
462   if (id != PETSCFE_CLASSID) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Discretization not supported\n");
463   ierr = DMGetCoordinateField(field->dm, &coordField);CHKERRQ(ierr);
464   ierr = DMFieldGetDegree(coordField, pointIS, NULL, &maxDegree);CHKERRQ(ierr);
465   if (maxDegree <= 1) {
466     ierr = DMFieldCreateDefaultQuadrature(coordField, pointIS, &quad);CHKERRQ(ierr);
467   }
468   if (!quad) {ierr = DMFieldCreateDefaultQuadrature(field, pointIS, &quad);CHKERRQ(ierr);}
469   if (!quad) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not determine quadrature for cell averages\n");
470   ierr = DMFieldCreateFEGeom(coordField,pointIS,quad,PETSC_FALSE,&geom);CHKERRQ(ierr);
471   ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, NULL, &weights);CHKERRQ(ierr);
472   if (qNc != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Expected scalar quadrature components\n");
473   N = numPoints * Nq * Nc;
474   if (B) {ierr = DMGetWorkArray(field->dm, N, mpitype, &qB);CHKERRQ(ierr);}
475   if (D) {ierr = DMGetWorkArray(field->dm, N * dimC, mpitype, &qD);CHKERRQ(ierr);}
476   if (H) {ierr = DMGetWorkArray(field->dm, N * dimC * dimC, mpitype, &qH);CHKERRQ(ierr);}
477   ierr = DMFieldEvaluateFE(field,pointIS,quad,type,qB,qD,qH);CHKERRQ(ierr);
478   if (B) {
479     PetscInt i, j, k;
480 
481     if (type == PETSC_SCALAR) {
482       PetscScalar * sB  = (PetscScalar *) B;
483       PetscScalar * sqB = (PetscScalar *) qB;
484 
485       for (i = 0; i < numPoints; i++) {
486         PetscReal vol = 0.;
487 
488         for (j = 0; j < Nc; j++) {sB[i * Nc + j] = 0.;}
489         for (k = 0; k < Nq; k++) {
490           vol += geom->detJ[i * Nq + k] * weights[k];
491           for (j = 0; j < Nc; j++) {
492             sB[i * Nc + j] += geom->detJ[i * Nq + k] * weights[k] * sqB[ (i * Nq + k) * Nc + j];
493           }
494         }
495         for (k = 0; k < Nc; k++) sB[i * Nc + k] /= vol;
496       }
497     } else {
498       PetscReal * rB  = (PetscReal *) B;
499       PetscReal * rqB = (PetscReal *) qB;
500 
501       for (i = 0; i < numPoints; i++) {
502         PetscReal vol = 0.;
503 
504         for (j = 0; j < Nc; j++) {rB[i * Nc + j] = 0.;}
505         for (k = 0; k < Nq; k++) {
506           vol += geom->detJ[i * Nq + k] * weights[k];
507           for (j = 0; j < Nc; j++) {
508             rB[i * Nc + j] += weights[k] * rqB[ (i * Nq + k) * Nc + j];
509           }
510         }
511         for (k = 0; k < Nc; k++) rB[i * Nc + k] /= vol;
512       }
513     }
514   }
515   if (D) {
516     PetscInt i, j, k, l, m;
517 
518     if (type == PETSC_SCALAR) {
519       PetscScalar * sD  = (PetscScalar *) D;
520       PetscScalar * sqD = (PetscScalar *) qD;
521 
522       for (i = 0; i < numPoints; i++) {
523         PetscReal vol = 0.;
524 
525         for (j = 0; j < Nc * dimC; j++) {sD[i * Nc * dimC + j] = 0.;}
526         for (k = 0; k < Nq; k++) {
527           vol += geom->detJ[i * Nq + k] * weights[k];
528           for (j = 0; j < Nc; j++) {
529             PetscScalar pD[3] = {0.,0.,0.};
530 
531             for (l = 0; l < dimC; l++) {
532               for (m = 0; m < dim; m++) {
533                 pD[l] += geom->invJ[((i * Nq + k) * dimC + m) * dimC + l] * sqD[((i * Nq + k) * Nc + j) * dim + m];
534               }
535             }
536             for (l = 0; l < dimC; l++) {
537               sD[(i * Nc + j) * dimC + l] += geom->detJ[i * Nq + k] * weights[k] * pD[l];
538             }
539           }
540         }
541         for (k = 0; k < Nc * dimC; k++) sD[i * Nc * dimC + k] /= vol;
542       }
543     } else {
544       PetscReal * rD  = (PetscReal *) D;
545       PetscReal * rqD = (PetscReal *) qD;
546 
547       for (i = 0; i < numPoints; i++) {
548         PetscReal vol = 0.;
549 
550         for (j = 0; j < Nc * dimC; j++) {rD[i * Nc * dimC + j] = 0.;}
551         for (k = 0; k < Nq; k++) {
552           vol += geom->detJ[i * Nq + k] * weights[k];
553           for (j = 0; j < Nc; j++) {
554             PetscReal pD[3] = {0.,0.,0.};
555 
556             for (l = 0; l < dimC; l++) {
557               for (m = 0; m < dim; m++) {
558                 pD[l] += geom->invJ[((i * Nq + k) * dimC + m) * dimC + l] * rqD[((i * Nq + k) * Nc + j) * dim + m];
559               }
560             }
561             for (l = 0; l < dimC; l++) {
562               rD[(i * Nc + j) * dimC + l] += geom->detJ[i * Nq + k] * weights[k] * pD[l];
563             }
564           }
565         }
566         for (k = 0; k < Nc * dimC; k++) rD[i * Nc * dimC + k] /= vol;
567       }
568     }
569   }
570   if (H) {
571     PetscInt i, j, k, l, m, q, r;
572 
573     if (type == PETSC_SCALAR) {
574       PetscScalar * sH  = (PetscScalar *) H;
575       PetscScalar * sqH = (PetscScalar *) qH;
576 
577       for (i = 0; i < numPoints; i++) {
578         PetscReal vol = 0.;
579 
580         for (j = 0; j < Nc * dimC * dimC; j++) {sH[i * Nc * dimC * dimC + j] = 0.;}
581         for (k = 0; k < Nq; k++) {
582           const PetscReal *invJ = &geom->invJ[(i * Nq + k) * dimC * dimC];
583 
584           vol += geom->detJ[i * Nq + k] * weights[k];
585           for (j = 0; j < Nc; j++) {
586             PetscScalar pH[3][3] = {{0.,0.,0.},{0.,0.,0.},{0.,0.,0.}};
587             const PetscScalar *spH = &sqH[((i * Nq + k) * Nc + j) * dimC * dimC];
588 
589             for (l = 0; l < dimC; l++) {
590               for (m = 0; m < dimC; m++) {
591                 for (q = 0; q < dim; q++) {
592                   for (r = 0; r < dim; r++) {
593                     pH[l][m] += invJ[q * dimC + l] * invJ[r * dimC + m] * spH[q * dim + r];
594                   }
595                 }
596               }
597             }
598             for (l = 0; l < dimC; l++) {
599               for (m = 0; m < dimC; m++) {
600                 sH[(i * Nc + j) * dimC * dimC + l * dimC + m] += geom->detJ[i * Nq + k] * weights[k] * pH[l][m];
601               }
602             }
603           }
604         }
605         for (k = 0; k < Nc * dimC * dimC; k++) sH[i * Nc * dimC * dimC + k] /= vol;
606       }
607     } else {
608       PetscReal * rH  = (PetscReal *) H;
609       PetscReal * rqH = (PetscReal *) qH;
610 
611       for (i = 0; i < numPoints; i++) {
612         PetscReal vol = 0.;
613 
614         for (j = 0; j < Nc * dimC * dimC; j++) {rH[i * Nc * dimC * dimC + j] = 0.;}
615         for (k = 0; k < Nq; k++) {
616           const PetscReal *invJ = &geom->invJ[(i * Nq + k) * dimC * dimC];
617 
618           vol += geom->detJ[i * Nq + k] * weights[k];
619           for (j = 0; j < Nc; j++) {
620             PetscReal pH[3][3] = {{0.,0.,0.},{0.,0.,0.},{0.,0.,0.}};
621             const PetscReal *rpH = &rqH[((i * Nq + k) * Nc + j) * dimC * dimC];
622 
623             for (l = 0; l < dimC; l++) {
624               for (m = 0; m < dimC; m++) {
625                 for (q = 0; q < dim; q++) {
626                   for (r = 0; r < dim; r++) {
627                     pH[l][m] += invJ[q * dimC + l] * invJ[r * dimC + m] * rpH[q * dim + r];
628                   }
629                 }
630               }
631             }
632             for (l = 0; l < dimC; l++) {
633               for (m = 0; m < dimC; m++) {
634                 rH[(i * Nc + j) * dimC * dimC + l * dimC + m] += geom->detJ[i * Nq + k] * weights[k] * pH[l][m];
635               }
636             }
637           }
638         }
639         for (k = 0; k < Nc * dimC * dimC; k++) rH[i * Nc * dimC * dimC + k] /= vol;
640       }
641     }
642   }
643   if (B) {ierr = DMRestoreWorkArray(field->dm, N, mpitype, &qB);CHKERRQ(ierr);}
644   if (D) {ierr = DMRestoreWorkArray(field->dm, N * dimC, mpitype, &qD);CHKERRQ(ierr);}
645   if (H) {ierr = DMRestoreWorkArray(field->dm, N * dimC * dimC, mpitype, &qH);CHKERRQ(ierr);}
646   ierr = PetscFEGeomDestroy(&geom);CHKERRQ(ierr);
647   ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr);
648   PetscFunctionReturn(0);
649 }
650 
651 static PetscErrorCode DMFieldGetDegree_DS(DMField field, IS pointIS, PetscInt *minDegree, PetscInt *maxDegree)
652 {
653   DMField_DS     *dsfield;
654   PetscObject    disc;
655   PetscInt       h, imin, imax;
656   PetscClassId   id;
657   PetscErrorCode ierr;
658 
659   PetscFunctionBegin;
660   dsfield = (DMField_DS *) field->data;
661   ierr = ISGetMinMax(pointIS,&imin,&imax);CHKERRQ(ierr);
662   if (imin >= imax) {
663     h = 0;
664   } else {
665     for (h = 0; h < dsfield->height; h++) {
666       PetscInt hEnd;
667 
668       ierr = DMPlexGetHeightStratum(field->dm,h,NULL,&hEnd);CHKERRQ(ierr);
669       if (imin < hEnd) break;
670     }
671   }
672   ierr = DMFieldDSGetHeightDisc(field,h,&disc);CHKERRQ(ierr);
673   ierr = PetscObjectGetClassId(disc,&id);CHKERRQ(ierr);
674   if (id == PETSCFE_CLASSID) {
675     PetscFE    fe = (PetscFE) disc;
676     PetscSpace sp;
677 
678     ierr = PetscFEGetBasisSpace(fe, &sp);CHKERRQ(ierr);
679     ierr = PetscSpaceGetDegree(sp, minDegree, maxDegree);CHKERRQ(ierr);
680   }
681   PetscFunctionReturn(0);
682 }
683 
684 static PetscErrorCode DMFieldCreateDefaultQuadrature_DS(DMField field, IS pointIS, PetscQuadrature *quad)
685 {
686   PetscInt       h, dim, imax, imin, cellHeight;
687   DM             dm;
688   DMField_DS     *dsfield;
689   PetscObject    disc;
690   PetscFE        fe;
691   PetscClassId   id;
692   PetscErrorCode ierr;
693 
694   PetscFunctionBegin;
695   dm = field->dm;
696   dsfield = (DMField_DS *) field->data;
697   ierr = ISGetMinMax(pointIS,&imin,&imax);CHKERRQ(ierr);
698   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
699   for (h = 0; h <= dim; h++) {
700     PetscInt hStart, hEnd;
701 
702     ierr = DMPlexGetHeightStratum(dm,h,&hStart,&hEnd);CHKERRQ(ierr);
703     if (imax >= hStart && imin < hEnd) break;
704   }
705   ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr);
706   h -= cellHeight;
707   *quad = NULL;
708   if (h < dsfield->height) {
709     ierr = DMFieldDSGetHeightDisc(field,h,&disc);CHKERRQ(ierr);
710     ierr = PetscObjectGetClassId(disc,&id);CHKERRQ(ierr);
711     if (id != PETSCFE_CLASSID) PetscFunctionReturn(0);
712     fe = (PetscFE) disc;
713     ierr = PetscFEGetQuadrature(fe,quad);CHKERRQ(ierr);
714     ierr = PetscObjectReference((PetscObject)*quad);CHKERRQ(ierr);
715   }
716   PetscFunctionReturn(0);
717 }
718 
719 static PetscErrorCode DMFieldComputeFaceData_DS(DMField field, IS pointIS, PetscQuadrature quad, PetscFEGeom *geom)
720 {
721   const PetscInt *points;
722   PetscInt        p, dim, dE, numFaces, Nq;
723   PetscInt        maxDegree;
724   DMLabel         depthLabel;
725   IS              cellIS;
726   DM              dm = field->dm;
727   PetscErrorCode  ierr;
728 
729   PetscFunctionBegin;
730   dim = geom->dim;
731   dE  = geom->dimEmbed;
732   ierr = DMPlexGetDepthLabel(dm, &depthLabel);CHKERRQ(ierr);
733   ierr = DMLabelGetStratumIS(depthLabel, dim + 1, &cellIS);CHKERRQ(ierr);
734   ierr = DMFieldGetDegree(field,cellIS,NULL,&maxDegree);CHKERRQ(ierr);
735   ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
736   numFaces = geom->numCells;
737   Nq = geom->numPoints;
738   /* First, set local faces and flip normals so that they are outward for the first supporting cell */
739   for (p = 0; p < numFaces; p++) {
740     PetscInt        point = points[p];
741     PetscInt        suppSize, s, coneSize, c, numChildren;
742     const PetscInt *supp, *cone, *ornt;
743 
744     ierr = DMPlexGetTreeChildren(dm, point, &numChildren, NULL);CHKERRQ(ierr);
745     if (numChildren) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face data not valid for facets with children");
746     ierr = DMPlexGetSupportSize(dm, point, &suppSize);CHKERRQ(ierr);
747     if (suppSize > 2) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D has %D support, expected at most 2\n", point, suppSize);
748     if (!suppSize) continue;
749     ierr = DMPlexGetSupport(dm, point, &supp);CHKERRQ(ierr);
750     for (s = 0; s < suppSize; ++s) {
751       ierr = DMPlexGetConeSize(dm, supp[s], &coneSize);CHKERRQ(ierr);
752       ierr = DMPlexGetCone(dm, supp[s], &cone);CHKERRQ(ierr);
753       for (c = 0; c < coneSize; ++c) if (cone[c] == point) break;
754       if (c == coneSize) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid connectivity: point %D not found in cone of support point %D", point, supp[s]);
755       geom->face[p][s] = c;
756     }
757     ierr = DMPlexGetConeOrientation(dm, supp[0], &ornt);CHKERRQ(ierr);
758     if (ornt[geom->face[p][0]] < 0) {
759       PetscInt Np = geom->numPoints, q, dE = geom->dimEmbed, d;
760 
761       for (q = 0; q < Np; ++q) for (d = 0; d < dE; ++d) geom->n[(p*Np + q)*dE + d] = -geom->n[(p*Np + q)*dE + d];
762     }
763   }
764   if (maxDegree <= 1) {
765     PetscInt        numCells, offset, *cells;
766     PetscFEGeom     *cellGeom;
767     IS              suppIS;
768 
769     for (p = 0, numCells = 0; p < numFaces; p++) {
770       PetscInt        point = points[p];
771       PetscInt        numSupp, numChildren;
772 
773       ierr = DMPlexGetTreeChildren(dm, point, &numChildren, NULL);CHKERRQ(ierr);
774       if (numChildren) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face data not valid for facets with children");
775       ierr = DMPlexGetSupportSize(dm, point,&numSupp);CHKERRQ(ierr);
776       if (numSupp > 2) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D has %D support, expected at most 2\n", point, numSupp);
777       numCells += numSupp;
778     }
779     ierr = PetscMalloc1(numCells, &cells);CHKERRQ(ierr);
780     for (p = 0, offset = 0; p < numFaces; p++) {
781       PetscInt        point = points[p];
782       PetscInt        numSupp, s;
783       const PetscInt *supp;
784 
785       ierr = DMPlexGetSupportSize(dm, point,&numSupp);CHKERRQ(ierr);
786       ierr = DMPlexGetSupport(dm, point, &supp);CHKERRQ(ierr);
787       for (s = 0; s < numSupp; s++, offset++) {
788         cells[offset] = supp[s];
789       }
790     }
791     ierr = ISCreateGeneral(PETSC_COMM_SELF,numCells,cells,PETSC_USE_POINTER, &suppIS);CHKERRQ(ierr);
792     ierr = DMFieldCreateFEGeom(field,suppIS,quad,PETSC_FALSE,&cellGeom);CHKERRQ(ierr);
793     for (p = 0, offset = 0; p < numFaces; p++) {
794       PetscInt        point = points[p];
795       PetscInt        numSupp, s, q;
796       const PetscInt *supp;
797 
798       ierr = DMPlexGetSupportSize(dm, point,&numSupp);CHKERRQ(ierr);
799       ierr = DMPlexGetSupport(dm, point, &supp);CHKERRQ(ierr);
800       for (s = 0; s < numSupp; s++, offset++) {
801         for (q = 0; q < Nq * dE * dE; q++) {
802           geom->suppJ[s][p * Nq * dE * dE + q]    = cellGeom->J[offset * Nq * dE * dE + q];
803           geom->suppInvJ[s][p * Nq * dE * dE + q] = cellGeom->invJ[offset * Nq * dE * dE + q];
804         }
805         for (q = 0; q < Nq; q++) geom->suppDetJ[s][p * Nq + q] = cellGeom->detJ[offset * Nq + q];
806       }
807     }
808     ierr = PetscFEGeomDestroy(&cellGeom);CHKERRQ(ierr);
809     ierr = ISDestroy(&suppIS);CHKERRQ(ierr);
810     ierr = PetscFree(cells);CHKERRQ(ierr);
811   } else {
812     PetscObject          faceDisc, cellDisc;
813     PetscClassId         faceId, cellId;
814     PetscDualSpace       dsp;
815     DM                   K;
816     DMPolytopeType       ct;
817     PetscInt           (*co)[2][3];
818     PetscInt             coneSize;
819     PetscInt           **counts;
820     PetscInt             f, i, o, q, s;
821     const PetscInt      *coneK;
822     PetscInt             eStart, minOrient, maxOrient, numOrient;
823     PetscInt            *orients;
824     PetscReal          **orientPoints;
825     PetscReal           *cellPoints;
826     PetscReal           *dummyWeights;
827     PetscQuadrature      cellQuad = NULL;
828 
829     ierr = DMFieldDSGetHeightDisc(field, 1, &faceDisc);CHKERRQ(ierr);
830     ierr = DMFieldDSGetHeightDisc(field, 0, &cellDisc);CHKERRQ(ierr);
831     ierr = PetscObjectGetClassId(faceDisc,&faceId);CHKERRQ(ierr);
832     ierr = PetscObjectGetClassId(cellDisc,&cellId);CHKERRQ(ierr);
833     if (faceId != PETSCFE_CLASSID || cellId != PETSCFE_CLASSID) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Not supported\n");
834     ierr = PetscFEGetDualSpace((PetscFE)cellDisc, &dsp);CHKERRQ(ierr);
835     ierr = PetscDualSpaceGetDM(dsp, &K);CHKERRQ(ierr);
836     ierr = DMPlexGetHeightStratum(K, 1, &eStart, NULL);CHKERRQ(ierr);
837     ierr = DMPlexGetCellType(K, eStart, &ct);CHKERRQ(ierr);
838     ierr = DMPlexGetConeSize(K,0,&coneSize);CHKERRQ(ierr);
839     ierr = DMPlexGetCone(K,0,&coneK);CHKERRQ(ierr);
840     ierr = PetscMalloc2(numFaces, &co, coneSize, &counts);CHKERRQ(ierr);
841     ierr = PetscMalloc1(dE*Nq, &cellPoints);CHKERRQ(ierr);
842     ierr = PetscMalloc1(Nq, &dummyWeights);CHKERRQ(ierr);
843     ierr = PetscQuadratureCreate(PETSC_COMM_SELF, &cellQuad);CHKERRQ(ierr);
844     ierr = PetscQuadratureSetData(cellQuad, dE, 1, Nq, cellPoints, dummyWeights);CHKERRQ(ierr);
845     minOrient = PETSC_MAX_INT;
846     maxOrient = PETSC_MIN_INT;
847     for (p = 0; p < numFaces; p++) { /* record the orientation of the facet wrt the support cells */
848       PetscInt        point = points[p];
849       PetscInt        numSupp, numChildren;
850       const PetscInt *supp;
851 
852       ierr = DMPlexGetTreeChildren(dm, point, &numChildren, NULL);CHKERRQ(ierr);
853       if (numChildren) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face data not valid for facets with children");
854       ierr = DMPlexGetSupportSize(dm, point,&numSupp);CHKERRQ(ierr);
855       if (numSupp > 2) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D has %D support, expected at most 2\n", point, numSupp);
856       ierr = DMPlexGetSupport(dm, point, &supp);CHKERRQ(ierr);
857       for (s = 0; s < numSupp; s++) {
858         PetscInt        cell = supp[s];
859         PetscInt        numCone;
860         const PetscInt *cone, *orient;
861 
862         ierr = DMPlexGetConeSize(dm, cell, &numCone);CHKERRQ(ierr);
863         if (numCone != coneSize) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Support point does not match reference element");
864         ierr = DMPlexGetCone(dm, cell, &cone);CHKERRQ(ierr);
865         ierr = DMPlexGetConeOrientation(dm, cell, &orient);CHKERRQ(ierr);
866         for (f = 0; f < coneSize; f++) {
867           if (cone[f] == point) break;
868         }
869         co[p][s][0] = f;
870         co[p][s][1] = orient[f];
871         co[p][s][2] = cell;
872         minOrient = PetscMin(minOrient, orient[f]);
873         maxOrient = PetscMax(maxOrient, orient[f]);
874       }
875       for (; s < 2; s++) {
876         co[p][s][0] = -1;
877         co[p][s][1] = -1;
878         co[p][s][2] = -1;
879       }
880     }
881     numOrient = maxOrient + 1 - minOrient;
882     ierr = DMPlexGetCone(K,0,&coneK);CHKERRQ(ierr);
883     /* count all (face,orientation) doubles that appear */
884     ierr = PetscCalloc2(numOrient,&orients,numOrient,&orientPoints);CHKERRQ(ierr);
885     for (f = 0; f < coneSize; f++) {ierr = PetscCalloc1(numOrient+1, &counts[f]);CHKERRQ(ierr);}
886     for (p = 0; p < numFaces; p++) {
887       for (s = 0; s < 2; s++) {
888         if (co[p][s][0] >= 0) {
889           counts[co[p][s][0]][co[p][s][1] - minOrient]++;
890           orients[co[p][s][1] - minOrient]++;
891         }
892       }
893     }
894     for (o = 0; o < numOrient; o++) {
895       if (orients[o]) {
896         PetscInt orient = o + minOrient;
897         PetscInt q;
898 
899         ierr = PetscMalloc1(Nq * dim, &orientPoints[o]);CHKERRQ(ierr);
900         /* rotate the quadrature points appropriately */
901         switch (ct) {
902         case DM_POLYTOPE_POINT: break;
903         case DM_POLYTOPE_SEGMENT:
904           if (orient == -2 || orient == 1) {
905             for (q = 0; q < Nq; q++) {
906               orientPoints[o][q] = -geom->xi[q];
907             }
908           } else {
909             for (q = 0; q < Nq; q++) {
910               orientPoints[o][q] = geom->xi[q];
911             }
912           }
913           break;
914         case DM_POLYTOPE_TRIANGLE:
915           for (q = 0; q < Nq; q++) {
916             PetscReal lambda[3];
917             PetscReal lambdao[3];
918 
919             /* convert to barycentric */
920             lambda[0] = - (geom->xi[2 * q] + geom->xi[2 * q + 1]) / 2.;
921             lambda[1] = (geom->xi[2 * q] + 1.) / 2.;
922             lambda[2] = (geom->xi[2 * q + 1] + 1.) / 2.;
923             if (orient >= 0) {
924               for (i = 0; i < 3; i++) {
925                 lambdao[i] = lambda[(orient + i) % 3];
926               }
927             } else {
928               for (i = 0; i < 3; i++) {
929                 lambdao[i] = lambda[(-(orient + i) + 3) % 3];
930               }
931             }
932             /* convert to coordinates */
933             orientPoints[o][2 * q + 0] = -(lambdao[0] + lambdao[2]) + lambdao[1];
934             orientPoints[o][2 * q + 1] = -(lambdao[0] + lambdao[1]) + lambdao[2];
935           }
936           break;
937         case DM_POLYTOPE_QUADRILATERAL:
938           for (q = 0; q < Nq; q++) {
939             PetscReal xi[2], xio[2];
940             PetscInt oabs = (orient >= 0) ? orient : -(orient + 1);
941 
942             xi[0] = geom->xi[2 * q];
943             xi[1] = geom->xi[2 * q + 1];
944             switch (oabs) {
945             case 1:
946               xio[0] = xi[1];
947               xio[1] = -xi[0];
948               break;
949             case 2:
950               xio[0] = -xi[0];
951               xio[1] = -xi[1];
952             case 3:
953               xio[0] = -xi[1];
954               xio[1] = xi[0];
955             case 0:
956             default:
957               xio[0] = xi[0];
958               xio[1] = xi[1];
959               break;
960             }
961             if (orient < 0) {
962               xio[0] = -xio[0];
963             }
964             orientPoints[o][2 * q + 0] = xio[0];
965             orientPoints[o][2 * q + 1] = xio[1];
966           }
967           break;
968         default: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cell type %s not yet supported\n", DMPolytopeTypes[ct]);
969         }
970       }
971     }
972     for (f = 0; f < coneSize; f++) {
973       PetscInt face = coneK[f];
974       PetscReal v0[3];
975       PetscReal J[9], detJ;
976       PetscInt numCells, offset;
977       PetscInt *cells;
978       IS suppIS;
979 
980       ierr = DMPlexComputeCellGeometryFEM(K, face, NULL, v0, J, NULL, &detJ);CHKERRQ(ierr);
981       for (o = 0; o <= numOrient; o++) {
982         PetscFEGeom *cellGeom;
983 
984         if (!counts[f][o]) continue;
985         /* If this (face,orientation) double appears,
986          * convert the face quadrature points into volume quadrature points */
987         for (q = 0; q < Nq; q++) {
988           PetscReal xi0[3] = {-1., -1., -1.};
989 
990           CoordinatesRefToReal(dE, dim, xi0, v0, J, &orientPoints[o][dim * q + 0], &cellPoints[dE * q + 0]);
991         }
992         for (p = 0, numCells = 0; p < numFaces; p++) {
993           for (s = 0; s < 2; s++) {
994             if (co[p][s][0] == f && co[p][s][1] == o + minOrient) numCells++;
995           }
996         }
997         ierr = PetscMalloc1(numCells, &cells);CHKERRQ(ierr);
998         for (p = 0, offset = 0; p < numFaces; p++) {
999           for (s = 0; s < 2; s++) {
1000             if (co[p][s][0] == f && co[p][s][1] == o + minOrient) {
1001               cells[offset++] = co[p][s][2];
1002             }
1003           }
1004         }
1005         ierr = ISCreateGeneral(PETSC_COMM_SELF,numCells,cells,PETSC_USE_POINTER, &suppIS);CHKERRQ(ierr);
1006         ierr = DMFieldCreateFEGeom(field,suppIS,cellQuad,PETSC_FALSE,&cellGeom);CHKERRQ(ierr);
1007         for (p = 0, offset = 0; p < numFaces; p++) {
1008           for (s = 0; s < 2; s++) {
1009             if (co[p][s][0] == f && co[p][s][1] == o + minOrient) {
1010               for (q = 0; q < Nq * dE * dE; q++) {
1011                 geom->suppJ[s][p * Nq * dE * dE + q]    = cellGeom->J[offset * Nq * dE * dE + q];
1012                 geom->suppInvJ[s][p * Nq * dE * dE + q] = cellGeom->invJ[offset * Nq * dE * dE + q];
1013               }
1014               for (q = 0; q < Nq; q++) geom->suppDetJ[s][p * Nq + q] = cellGeom->detJ[offset * Nq + q];
1015               offset++;
1016             }
1017           }
1018         }
1019         ierr = PetscFEGeomDestroy(&cellGeom);CHKERRQ(ierr);
1020         ierr = ISDestroy(&suppIS);CHKERRQ(ierr);
1021         ierr = PetscFree(cells);CHKERRQ(ierr);
1022       }
1023     }
1024     for (o = 0; o < numOrient; o++) {
1025       if (orients[o]) {
1026         ierr = PetscFree(orientPoints[o]);CHKERRQ(ierr);
1027       }
1028     }
1029     ierr = PetscFree2(orients,orientPoints);CHKERRQ(ierr);
1030     ierr = PetscQuadratureDestroy(&cellQuad);CHKERRQ(ierr);
1031     for (f = 0; f < coneSize; f++) {ierr = PetscFree(counts[f]);CHKERRQ(ierr);}
1032     ierr = PetscFree2(co,counts);CHKERRQ(ierr);
1033   }
1034   ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
1035   ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
1036   PetscFunctionReturn(0);
1037 }
1038 
1039 static PetscErrorCode DMFieldInitialize_DS(DMField field)
1040 {
1041   PetscFunctionBegin;
1042   field->ops->destroy                 = DMFieldDestroy_DS;
1043   field->ops->evaluate                = DMFieldEvaluate_DS;
1044   field->ops->evaluateFE              = DMFieldEvaluateFE_DS;
1045   field->ops->evaluateFV              = DMFieldEvaluateFV_DS;
1046   field->ops->getDegree               = DMFieldGetDegree_DS;
1047   field->ops->createDefaultQuadrature = DMFieldCreateDefaultQuadrature_DS;
1048   field->ops->view                    = DMFieldView_DS;
1049   field->ops->computeFaceData         = DMFieldComputeFaceData_DS;
1050   PetscFunctionReturn(0);
1051 }
1052 
1053 PETSC_INTERN PetscErrorCode DMFieldCreate_DS(DMField field)
1054 {
1055   DMField_DS     *dsfield;
1056   PetscErrorCode ierr;
1057 
1058   PetscFunctionBegin;
1059   ierr = PetscNewLog(field,&dsfield);CHKERRQ(ierr);
1060   field->data = dsfield;
1061   ierr = DMFieldInitialize_DS(field);CHKERRQ(ierr);
1062   PetscFunctionReturn(0);
1063 }
1064 
1065 PetscErrorCode DMFieldCreateDS(DM dm, PetscInt fieldNum, Vec vec,DMField *field)
1066 {
1067   DMField        b;
1068   DMField_DS     *dsfield;
1069   PetscObject    disc = NULL;
1070   PetscBool      isContainer = PETSC_FALSE;
1071   PetscClassId   id = -1;
1072   PetscInt       numComponents = -1, dsNumFields;
1073   PetscSection   section;
1074   PetscErrorCode ierr;
1075 
1076   PetscFunctionBegin;
1077   ierr = DMGetLocalSection(dm,&section);CHKERRQ(ierr);
1078   ierr = PetscSectionGetFieldComponents(section,fieldNum,&numComponents);CHKERRQ(ierr);
1079   ierr = DMGetNumFields(dm,&dsNumFields);CHKERRQ(ierr);
1080   if (dsNumFields) {ierr = DMGetField(dm,fieldNum,NULL,&disc);CHKERRQ(ierr);}
1081   if (disc) {
1082     ierr = PetscObjectGetClassId(disc,&id);CHKERRQ(ierr);
1083     isContainer = (id == PETSC_CONTAINER_CLASSID) ? PETSC_TRUE : PETSC_FALSE;
1084   }
1085   if (!disc || isContainer) {
1086     MPI_Comm        comm = PetscObjectComm((PetscObject) dm);
1087     PetscInt        cStart, cEnd, dim, cellHeight;
1088     PetscInt        localConeSize = 0, coneSize;
1089     PetscFE         fe;
1090     PetscDualSpace  Q;
1091     PetscSpace      P;
1092     DM              K;
1093     PetscQuadrature quad, fquad;
1094     PetscBool       isSimplex;
1095 
1096     ierr = DMPlexGetVTKCellHeight(dm, &cellHeight);CHKERRQ(ierr);
1097     ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
1098     ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1099     if (cEnd > cStart) {
1100       ierr = DMPlexGetConeSize(dm, cStart, &localConeSize);CHKERRQ(ierr);
1101     }
1102     ierr = MPI_Allreduce(&localConeSize,&coneSize,1,MPIU_INT,MPI_MAX,comm);CHKERRMPI(ierr);
1103     isSimplex = (coneSize == (dim + 1)) ? PETSC_TRUE : PETSC_FALSE;
1104     ierr = PetscSpaceCreate(PETSC_COMM_SELF, &P);CHKERRQ(ierr);
1105     ierr = PetscSpaceSetType(P,PETSCSPACEPOLYNOMIAL);CHKERRQ(ierr);
1106     ierr = PetscSpaceSetDegree(P, 1, PETSC_DETERMINE);CHKERRQ(ierr);
1107     ierr = PetscSpaceSetNumComponents(P, numComponents);CHKERRQ(ierr);
1108     ierr = PetscSpaceSetNumVariables(P, dim);CHKERRQ(ierr);
1109     ierr = PetscSpacePolynomialSetTensor(P, isSimplex ? PETSC_FALSE : PETSC_TRUE);CHKERRQ(ierr);
1110     ierr = PetscSpaceSetUp(P);CHKERRQ(ierr);
1111     ierr = PetscDualSpaceCreate(PETSC_COMM_SELF, &Q);CHKERRQ(ierr);
1112     ierr = PetscDualSpaceSetType(Q,PETSCDUALSPACELAGRANGE);CHKERRQ(ierr);
1113     ierr = PetscDualSpaceCreateReferenceCell(Q, dim, isSimplex, &K);CHKERRQ(ierr);
1114     ierr = PetscDualSpaceSetDM(Q, K);CHKERRQ(ierr);
1115     ierr = DMDestroy(&K);CHKERRQ(ierr);
1116     ierr = PetscDualSpaceSetNumComponents(Q, numComponents);CHKERRQ(ierr);
1117     ierr = PetscDualSpaceSetOrder(Q, 1);CHKERRQ(ierr);
1118     ierr = PetscDualSpaceLagrangeSetTensor(Q, isSimplex ? PETSC_FALSE : PETSC_TRUE);CHKERRQ(ierr);
1119     ierr = PetscDualSpaceSetUp(Q);CHKERRQ(ierr);
1120     ierr = PetscFECreate(PETSC_COMM_SELF, &fe);CHKERRQ(ierr);
1121     ierr = PetscFESetType(fe,PETSCFEBASIC);CHKERRQ(ierr);
1122     ierr = PetscFESetBasisSpace(fe, P);CHKERRQ(ierr);
1123     ierr = PetscFESetDualSpace(fe, Q);CHKERRQ(ierr);
1124     ierr = PetscFESetNumComponents(fe, numComponents);CHKERRQ(ierr);
1125     ierr = PetscFESetUp(fe);CHKERRQ(ierr);
1126     ierr = PetscSpaceDestroy(&P);CHKERRQ(ierr);
1127     ierr = PetscDualSpaceDestroy(&Q);CHKERRQ(ierr);
1128     if (isSimplex) {
1129       ierr = PetscDTStroudConicalQuadrature(dim,   1, 1, -1.0, 1.0, &quad);CHKERRQ(ierr);
1130       ierr = PetscDTStroudConicalQuadrature(dim-1, 1, 1, -1.0, 1.0, &fquad);CHKERRQ(ierr);
1131     }
1132     else {
1133       ierr = PetscDTGaussTensorQuadrature(dim,   1, 1, -1.0, 1.0, &quad);CHKERRQ(ierr);
1134       ierr = PetscDTGaussTensorQuadrature(dim-1, 1, 1, -1.0, 1.0, &fquad);CHKERRQ(ierr);
1135     }
1136     ierr = PetscFESetQuadrature(fe, quad);CHKERRQ(ierr);
1137     ierr = PetscFESetFaceQuadrature(fe, fquad);CHKERRQ(ierr);
1138     ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr);
1139     ierr = PetscQuadratureDestroy(&fquad);CHKERRQ(ierr);
1140     disc = (PetscObject) fe;
1141   } else {
1142     ierr = PetscObjectReference(disc);CHKERRQ(ierr);
1143   }
1144   ierr = PetscObjectGetClassId(disc,&id);CHKERRQ(ierr);
1145   if (id == PETSCFE_CLASSID) {
1146     PetscFE fe = (PetscFE) disc;
1147 
1148     ierr = PetscFEGetNumComponents(fe,&numComponents);CHKERRQ(ierr);
1149   } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Not implemented");
1150   ierr = DMFieldCreate(dm,numComponents,DMFIELD_VERTEX,&b);CHKERRQ(ierr);
1151   ierr = DMFieldSetType(b,DMFIELDDS);CHKERRQ(ierr);
1152   dsfield = (DMField_DS *) b->data;
1153   dsfield->fieldNum = fieldNum;
1154   ierr = DMGetDimension(dm,&dsfield->height);CHKERRQ(ierr);
1155   dsfield->height++;
1156   ierr = PetscCalloc1(dsfield->height,&dsfield->disc);CHKERRQ(ierr);
1157   dsfield->disc[0] = disc;
1158   ierr = PetscObjectReference((PetscObject)vec);CHKERRQ(ierr);
1159   dsfield->vec = vec;
1160   *field = b;
1161 
1162   PetscFunctionReturn(0);
1163 }
1164