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