xref: /petsc/src/dm/dt/dualspace/impls/lagrange/dspacelagrange.c (revision 69cc43acf904673348835f1e68d7a5ddf92b44cb)
1 #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
2 #include <petscdmplex.h>
3 #include <petscblaslapack.h>
4 
5 PetscErrorCode DMPlexGetTransitiveClosure_Internal(DM, PetscInt, PetscInt, PetscBool, PetscInt *, PetscInt *[]);
6 
7 struct _n_Petsc1DNodeFamily
8 {
9   PetscInt         refct;
10   PetscDTNodeType  nodeFamily;
11   PetscReal        gaussJacobiExp;
12   PetscInt         nComputed;
13   PetscReal      **nodesets;
14   PetscBool        endpoints;
15 };
16 
17 /* users set node families for PETSCDUALSPACELAGRANGE with just the inputs to this function, but internally we create
18  * an object that can cache the computations across multiple dual spaces */
19 static PetscErrorCode Petsc1DNodeFamilyCreate(PetscDTNodeType family, PetscReal gaussJacobiExp, PetscBool endpoints, Petsc1DNodeFamily *nf)
20 {
21   Petsc1DNodeFamily f;
22   PetscErrorCode ierr;
23 
24   PetscFunctionBegin;
25   ierr = PetscNew(&f);CHKERRQ(ierr);
26   switch (family) {
27   case PETSCDTNODES_GAUSSJACOBI:
28   case PETSCDTNODES_EQUISPACED:
29     f->nodeFamily = family;
30     break;
31   default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unknown 1D node family");
32   }
33   f->endpoints = endpoints;
34   f->gaussJacobiExp = 0.;
35   if (family == PETSCDTNODES_GAUSSJACOBI) {
36     PetscCheckFalse(gaussJacobiExp <= -1.,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Gauss-Jacobi exponent must be > -1.");
37     f->gaussJacobiExp = gaussJacobiExp;
38   }
39   f->refct = 1;
40   *nf = f;
41   PetscFunctionReturn(0);
42 }
43 
44 static PetscErrorCode Petsc1DNodeFamilyReference(Petsc1DNodeFamily nf)
45 {
46   PetscFunctionBegin;
47   if (nf) nf->refct++;
48   PetscFunctionReturn(0);
49 }
50 
51 static PetscErrorCode Petsc1DNodeFamilyDestroy(Petsc1DNodeFamily *nf)
52 {
53   PetscInt       i, nc;
54   PetscErrorCode ierr;
55 
56   PetscFunctionBegin;
57   if (!(*nf)) PetscFunctionReturn(0);
58   if (--(*nf)->refct > 0) {
59     *nf = NULL;
60     PetscFunctionReturn(0);
61   }
62   nc = (*nf)->nComputed;
63   for (i = 0; i < nc; i++) {
64     ierr = PetscFree((*nf)->nodesets[i]);CHKERRQ(ierr);
65   }
66   ierr = PetscFree((*nf)->nodesets);CHKERRQ(ierr);
67   ierr = PetscFree(*nf);CHKERRQ(ierr);
68   *nf = NULL;
69   PetscFunctionReturn(0);
70 }
71 
72 static PetscErrorCode Petsc1DNodeFamilyGetNodeSets(Petsc1DNodeFamily f, PetscInt degree, PetscReal ***nodesets)
73 {
74   PetscInt       nc;
75   PetscErrorCode ierr;
76 
77   PetscFunctionBegin;
78   nc = f->nComputed;
79   if (degree >= nc) {
80     PetscInt    i, j;
81     PetscReal **new_nodesets;
82     PetscReal  *w;
83 
84     ierr = PetscMalloc1(degree + 1, &new_nodesets);CHKERRQ(ierr);
85     ierr = PetscArraycpy(new_nodesets, f->nodesets, nc);CHKERRQ(ierr);
86     ierr = PetscFree(f->nodesets);CHKERRQ(ierr);
87     f->nodesets = new_nodesets;
88     ierr = PetscMalloc1(degree + 1, &w);CHKERRQ(ierr);
89     for (i = nc; i < degree + 1; i++) {
90       ierr = PetscMalloc1(i + 1, &(f->nodesets[i]));CHKERRQ(ierr);
91       if (!i) {
92         f->nodesets[i][0] = 0.5;
93       } else {
94         switch (f->nodeFamily) {
95         case PETSCDTNODES_EQUISPACED:
96           if (f->endpoints) {
97             for (j = 0; j <= i; j++) f->nodesets[i][j] = (PetscReal) j / (PetscReal) i;
98           } else {
99             /* these nodes are at the centroids of the small simplices created by the equispaced nodes that include
100              * the endpoints */
101             for (j = 0; j <= i; j++) f->nodesets[i][j] = ((PetscReal) j + 0.5) / ((PetscReal) i + 1.);
102           }
103           break;
104         case PETSCDTNODES_GAUSSJACOBI:
105           if (f->endpoints) {
106             ierr = PetscDTGaussLobattoJacobiQuadrature(i + 1, 0., 1., f->gaussJacobiExp, f->gaussJacobiExp, f->nodesets[i], w);CHKERRQ(ierr);
107           } else {
108             ierr = PetscDTGaussJacobiQuadrature(i + 1, 0., 1., f->gaussJacobiExp, f->gaussJacobiExp, f->nodesets[i], w);CHKERRQ(ierr);
109           }
110           break;
111         default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unknown 1D node family");
112         }
113       }
114     }
115     ierr = PetscFree(w);CHKERRQ(ierr);
116     f->nComputed = degree + 1;
117   }
118   *nodesets = f->nodesets;
119   PetscFunctionReturn(0);
120 }
121 
122 /* http://arxiv.org/abs/2002.09421 for details */
123 static PetscErrorCode PetscNodeRecursive_Internal(PetscInt dim, PetscInt degree, PetscReal **nodesets, PetscInt tup[], PetscReal node[])
124 {
125   PetscReal w;
126   PetscInt i, j;
127   PetscErrorCode ierr;
128 
129   PetscFunctionBeginHot;
130   w = 0.;
131   if (dim == 1) {
132     node[0] = nodesets[degree][tup[0]];
133     node[1] = nodesets[degree][tup[1]];
134   } else {
135     for (i = 0; i < dim + 1; i++) node[i] = 0.;
136     for (i = 0; i < dim + 1; i++) {
137       PetscReal wi = nodesets[degree][degree-tup[i]];
138 
139       for (j = 0; j < dim+1; j++) tup[dim+1+j] = tup[j+(j>=i)];
140       ierr = PetscNodeRecursive_Internal(dim-1,degree-tup[i],nodesets,&tup[dim+1],&node[dim+1]);CHKERRQ(ierr);
141       for (j = 0; j < dim+1; j++) node[j+(j>=i)] += wi * node[dim+1+j];
142       w += wi;
143     }
144     for (i = 0; i < dim+1; i++) node[i] /= w;
145   }
146   PetscFunctionReturn(0);
147 }
148 
149 /* compute simplex nodes for the biunit simplex from the 1D node family */
150 static PetscErrorCode Petsc1DNodeFamilyComputeSimplexNodes(Petsc1DNodeFamily f, PetscInt dim, PetscInt degree, PetscReal points[])
151 {
152   PetscInt      *tup;
153   PetscInt       k;
154   PetscInt       npoints;
155   PetscReal    **nodesets = NULL;
156   PetscInt       worksize;
157   PetscReal     *nodework;
158   PetscInt      *tupwork;
159   PetscErrorCode ierr;
160 
161   PetscFunctionBegin;
162   PetscCheckFalse(dim < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must have non-negative dimension");
163   PetscCheckFalse(degree < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must have non-negative degree");
164   if (!dim) PetscFunctionReturn(0);
165   ierr = PetscCalloc1(dim+2, &tup);CHKERRQ(ierr);
166   k = 0;
167   ierr = PetscDTBinomialInt(degree + dim, dim, &npoints);CHKERRQ(ierr);
168   ierr = Petsc1DNodeFamilyGetNodeSets(f, degree, &nodesets);CHKERRQ(ierr);
169   worksize = ((dim + 2) * (dim + 3)) / 2;
170   ierr = PetscMalloc2(worksize, &nodework, worksize, &tupwork);CHKERRQ(ierr);
171   /* loop over the tuples of length dim with sum at most degree */
172   for (k = 0; k < npoints; k++) {
173     PetscInt i;
174 
175     /* turn thm into tuples of length dim + 1 with sum equal to degree (barycentric indice) */
176     tup[0] = degree;
177     for (i = 0; i < dim; i++) {
178       tup[0] -= tup[i+1];
179     }
180     switch(f->nodeFamily) {
181     case PETSCDTNODES_EQUISPACED:
182       /* compute equispaces nodes on the unit reference triangle */
183       if (f->endpoints) {
184         for (i = 0; i < dim; i++) {
185           points[dim*k + i] = (PetscReal) tup[i+1] / (PetscReal) degree;
186         }
187       } else {
188         for (i = 0; i < dim; i++) {
189           /* these nodes are at the centroids of the small simplices created by the equispaced nodes that include
190            * the endpoints */
191           points[dim*k + i] = ((PetscReal) tup[i+1] + 1./(dim+1.)) / (PetscReal) (degree + 1.);
192         }
193       }
194       break;
195     default:
196       /* compute equispaces nodes on the barycentric reference triangle (the trace on the first dim dimensions are the
197        * unit reference triangle nodes */
198       for (i = 0; i < dim + 1; i++) tupwork[i] = tup[i];
199       ierr = PetscNodeRecursive_Internal(dim, degree, nodesets, tupwork, nodework);CHKERRQ(ierr);
200       for (i = 0; i < dim; i++) points[dim*k + i] = nodework[i + 1];
201       break;
202     }
203     ierr = PetscDualSpaceLatticePointLexicographic_Internal(dim, degree, &tup[1]);CHKERRQ(ierr);
204   }
205   /* map from unit simplex to biunit simplex */
206   for (k = 0; k < npoints * dim; k++) points[k] = points[k] * 2. - 1.;
207   ierr = PetscFree2(nodework, tupwork);CHKERRQ(ierr);
208   ierr = PetscFree(tup);CHKERRQ(ierr);
209   PetscFunctionReturn(0);
210 }
211 
212 /* If we need to get the dofs from a mesh point, or add values into dofs at a mesh point, and there is more than one dof
213  * on that mesh point, we have to be careful about getting/adding everything in the right place.
214  *
215  * With nodal dofs like PETSCDUALSPACELAGRANGE makes, the general approach to calculate the value of dofs associate
216  * with a node A is
217  * - transform the node locations x(A) by the map that takes the mesh point to its reorientation, x' = phi(x(A))
218  * - figure out which node was originally at the location of the transformed point, A' = idx(x')
219  * - if the dofs are not scalars, figure out how to represent the transformed dofs in terms of the basis
220  *   of dofs at A' (using pushforward/pullback rules)
221  *
222  * The one sticky point with this approach is the "A' = idx(x')" step: trying to go from real valued coordinates
223  * back to indices.  I don't want to rely on floating point tolerances.  Additionally, PETSCDUALSPACELAGRANGE may
224  * eventually support quasi-Lagrangian dofs, which could involve quadrature at multiple points, so the location "x(A)"
225  * would be ambiguous.
226  *
227  * So each dof gets an integer value coordinate (nodeIdx in the structure below).  The choice of integer coordinates
228  * is somewhat arbitrary, as long as all of the relevant symmetries of the mesh point correspond to *permutations* of
229  * the integer coordinates, which do not depend on numerical precision.
230  *
231  * So
232  *
233  * - DMPlexGetTransitiveClosure_Internal() tells me how an orientation turns into a permutation of the vertices of a
234  *   mesh point
235  * - The permutation of the vertices, and the nodeIdx values assigned to them, tells what permutation in index space
236  *   is associated with the orientation
237  * - I uses that permutation to get xi' = phi(xi(A)), the integer coordinate of the transformed dof
238  * - I can without numerical issues compute A' = idx(xi')
239  *
240  * Here are some examples of how the process works
241  *
242  * - With a triangle:
243  *
244  *   The triangle has the following integer coordinates for vertices, taken from the barycentric triangle
245  *
246  *     closure order 2
247  *     nodeIdx (0,0,1)
248  *      \
249  *       +
250  *       |\
251  *       | \
252  *       |  \
253  *       |   \    closure order 1
254  *       |    \ / nodeIdx (0,1,0)
255  *       +-----+
256  *        \
257  *      closure order 0
258  *      nodeIdx (1,0,0)
259  *
260  *   If I do DMPlexGetTransitiveClosure_Internal() with orientation 1, the vertices would appear
261  *   in the order (1, 2, 0)
262  *
263  *   If I list the nodeIdx of each vertex in closure order for orientation 0 (0, 1, 2) and orientation 1 (1, 2, 0), I
264  *   see
265  *
266  *   orientation 0  | orientation 1
267  *
268  *   [0] (1,0,0)      [1] (0,1,0)
269  *   [1] (0,1,0)      [2] (0,0,1)
270  *   [2] (0,0,1)      [0] (1,0,0)
271  *          A                B
272  *
273  *   In other words, B is the result of a row permutation of A.  But, there is also
274  *   a column permutation that accomplishes the same result, (2,0,1).
275  *
276  *   So if a dof has nodeIdx coordinate (a,b,c), after the transformation its nodeIdx coordinate
277  *   is (c,a,b), and the transformed degree of freedom will be a linear combination of dofs
278  *   that originally had coordinate (c,a,b).
279  *
280  * - With a quadrilateral:
281  *
282  *   The quadrilateral has the following integer coordinates for vertices, taken from concatenating barycentric
283  *   coordinates for two segments:
284  *
285  *     closure order 3      closure order 2
286  *     nodeIdx (1,0,0,1)    nodeIdx (0,1,0,1)
287  *                   \      /
288  *                    +----+
289  *                    |    |
290  *                    |    |
291  *                    +----+
292  *                   /      \
293  *     closure order 0      closure order 1
294  *     nodeIdx (1,0,1,0)    nodeIdx (0,1,1,0)
295  *
296  *   If I do DMPlexGetTransitiveClosure_Internal() with orientation 1, the vertices would appear
297  *   in the order (1, 2, 3, 0)
298  *
299  *   If I list the nodeIdx of each vertex in closure order for orientation 0 (0, 1, 2, 3) and
300  *   orientation 1 (1, 2, 3, 0), I see
301  *
302  *   orientation 0  | orientation 1
303  *
304  *   [0] (1,0,1,0)    [1] (0,1,1,0)
305  *   [1] (0,1,1,0)    [2] (0,1,0,1)
306  *   [2] (0,1,0,1)    [3] (1,0,0,1)
307  *   [3] (1,0,0,1)    [0] (1,0,1,0)
308  *          A                B
309  *
310  *   The column permutation that accomplishes the same result is (3,2,0,1).
311  *
312  *   So if a dof has nodeIdx coordinate (a,b,c,d), after the transformation its nodeIdx coordinate
313  *   is (d,c,a,b), and the transformed degree of freedom will be a linear combination of dofs
314  *   that originally had coordinate (d,c,a,b).
315  *
316  * Previously PETSCDUALSPACELAGRANGE had hardcoded symmetries for the triangle and quadrilateral,
317  * but this approach will work for any polytope, such as the wedge (triangular prism).
318  */
319 struct _n_PetscLagNodeIndices
320 {
321   PetscInt   refct;
322   PetscInt   nodeIdxDim;
323   PetscInt   nodeVecDim;
324   PetscInt   nNodes;
325   PetscInt  *nodeIdx;      /* for each node an index of size nodeIdxDim */
326   PetscReal *nodeVec;      /* for each node a vector of size nodeVecDim */
327   PetscInt  *perm;         /* if these are vertices, perm takes DMPlex point index to closure order;
328                               if these are nodes, perm lists nodes in index revlex order */
329 };
330 
331 /* this is just here so I can access the values in tests/ex1.c outside the library */
332 PetscErrorCode PetscLagNodeIndicesGetData_Internal(PetscLagNodeIndices ni, PetscInt *nodeIdxDim, PetscInt *nodeVecDim, PetscInt *nNodes, const PetscInt *nodeIdx[], const PetscReal *nodeVec[])
333 {
334   PetscFunctionBegin;
335   *nodeIdxDim = ni->nodeIdxDim;
336   *nodeVecDim = ni->nodeVecDim;
337   *nNodes = ni->nNodes;
338   *nodeIdx = ni->nodeIdx;
339   *nodeVec = ni->nodeVec;
340   PetscFunctionReturn(0);
341 }
342 
343 static PetscErrorCode PetscLagNodeIndicesReference(PetscLagNodeIndices ni)
344 {
345   PetscFunctionBegin;
346   if (ni) ni->refct++;
347   PetscFunctionReturn(0);
348 }
349 
350 static PetscErrorCode PetscLagNodeIndicesDuplicate(PetscLagNodeIndices ni, PetscLagNodeIndices *niNew)
351 {
352   PetscErrorCode ierr;
353 
354   PetscFunctionBegin;
355   ierr = PetscNew(niNew);CHKERRQ(ierr);
356   (*niNew)->refct = 1;
357   (*niNew)->nodeIdxDim = ni->nodeIdxDim;
358   (*niNew)->nodeVecDim = ni->nodeVecDim;
359   (*niNew)->nNodes = ni->nNodes;
360   ierr = PetscMalloc1(ni->nNodes * ni->nodeIdxDim, &((*niNew)->nodeIdx));CHKERRQ(ierr);
361   ierr = PetscArraycpy((*niNew)->nodeIdx, ni->nodeIdx, ni->nNodes * ni->nodeIdxDim);CHKERRQ(ierr);
362   ierr = PetscMalloc1(ni->nNodes * ni->nodeVecDim, &((*niNew)->nodeVec));CHKERRQ(ierr);
363   ierr = PetscArraycpy((*niNew)->nodeVec, ni->nodeVec, ni->nNodes * ni->nodeVecDim);CHKERRQ(ierr);
364   (*niNew)->perm = NULL;
365   PetscFunctionReturn(0);
366 }
367 
368 static PetscErrorCode PetscLagNodeIndicesDestroy(PetscLagNodeIndices *ni)
369 {
370   PetscErrorCode ierr;
371 
372   PetscFunctionBegin;
373   if (!(*ni)) PetscFunctionReturn(0);
374   if (--(*ni)->refct > 0) {
375     *ni = NULL;
376     PetscFunctionReturn(0);
377   }
378   ierr = PetscFree((*ni)->nodeIdx);CHKERRQ(ierr);
379   ierr = PetscFree((*ni)->nodeVec);CHKERRQ(ierr);
380   ierr = PetscFree((*ni)->perm);CHKERRQ(ierr);
381   ierr = PetscFree(*ni);CHKERRQ(ierr);
382   *ni = NULL;
383   PetscFunctionReturn(0);
384 }
385 
386 /* The vertices are given nodeIdx coordinates (e.g. the corners of the barycentric triangle).  Those coordinates are
387  * in some other order, and to understand the effect of different symmetries, we need them to be in closure order.
388  *
389  * If sortIdx is PETSC_FALSE, the coordinates are already in revlex order, otherwise we must sort them
390  * to that order before we do the real work of this function, which is
391  *
392  * - mark the vertices in closure order
393  * - sort them in revlex order
394  * - use the resulting permutation to list the vertex coordinates in closure order
395  */
396 static PetscErrorCode PetscLagNodeIndicesComputeVertexOrder(DM dm, PetscLagNodeIndices ni, PetscBool sortIdx)
397 {
398   PetscInt        v, w, vStart, vEnd, c, d;
399   PetscInt        nVerts;
400   PetscInt        closureSize = 0;
401   PetscInt       *closure = NULL;
402   PetscInt       *closureOrder;
403   PetscInt       *invClosureOrder;
404   PetscInt       *revlexOrder;
405   PetscInt       *newNodeIdx;
406   PetscInt        dim;
407   Vec             coordVec;
408   const PetscScalar *coords;
409   PetscErrorCode  ierr;
410 
411   PetscFunctionBegin;
412   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
413   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
414   nVerts = vEnd - vStart;
415   ierr = PetscMalloc1(nVerts, &closureOrder);CHKERRQ(ierr);
416   ierr = PetscMalloc1(nVerts, &invClosureOrder);CHKERRQ(ierr);
417   ierr = PetscMalloc1(nVerts, &revlexOrder);CHKERRQ(ierr);
418   if (sortIdx) { /* bubble sort nodeIdx into revlex order */
419     PetscInt nodeIdxDim = ni->nodeIdxDim;
420     PetscInt *idxOrder;
421 
422     ierr = PetscMalloc1(nVerts * nodeIdxDim, &newNodeIdx);CHKERRQ(ierr);
423     ierr = PetscMalloc1(nVerts, &idxOrder);CHKERRQ(ierr);
424     for (v = 0; v < nVerts; v++) idxOrder[v] = v;
425     for (v = 0; v < nVerts; v++) {
426       for (w = v + 1; w < nVerts; w++) {
427         const PetscInt *iv = &(ni->nodeIdx[idxOrder[v] * nodeIdxDim]);
428         const PetscInt *iw = &(ni->nodeIdx[idxOrder[w] * nodeIdxDim]);
429         PetscInt diff = 0;
430 
431         for (d = nodeIdxDim - 1; d >= 0; d--) if ((diff = (iv[d] - iw[d]))) break;
432         if (diff > 0) {
433           PetscInt swap = idxOrder[v];
434 
435           idxOrder[v] = idxOrder[w];
436           idxOrder[w] = swap;
437         }
438       }
439     }
440     for (v = 0; v < nVerts; v++) {
441       for (d = 0; d < nodeIdxDim; d++) {
442         newNodeIdx[v * ni->nodeIdxDim + d] = ni->nodeIdx[idxOrder[v] * nodeIdxDim + d];
443       }
444     }
445     ierr = PetscFree(ni->nodeIdx);CHKERRQ(ierr);
446     ni->nodeIdx = newNodeIdx;
447     newNodeIdx = NULL;
448     ierr = PetscFree(idxOrder);CHKERRQ(ierr);
449   }
450   ierr = DMPlexGetTransitiveClosure(dm, 0, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
451   c = closureSize - nVerts;
452   for (v = 0; v < nVerts; v++) closureOrder[v] = closure[2 * (c + v)] - vStart;
453   for (v = 0; v < nVerts; v++) invClosureOrder[closureOrder[v]] = v;
454   ierr = DMPlexRestoreTransitiveClosure(dm, 0, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
455   ierr = DMGetCoordinatesLocal(dm, &coordVec);CHKERRQ(ierr);
456   ierr = VecGetArrayRead(coordVec, &coords);CHKERRQ(ierr);
457   /* bubble sort closure vertices by coordinates in revlex order */
458   for (v = 0; v < nVerts; v++) revlexOrder[v] = v;
459   for (v = 0; v < nVerts; v++) {
460     for (w = v + 1; w < nVerts; w++) {
461       const PetscScalar *cv = &coords[closureOrder[revlexOrder[v]] * dim];
462       const PetscScalar *cw = &coords[closureOrder[revlexOrder[w]] * dim];
463       PetscReal diff = 0;
464 
465       for (d = dim - 1; d >= 0; d--) if ((diff = PetscRealPart(cv[d] - cw[d])) != 0.) break;
466       if (diff > 0.) {
467         PetscInt swap = revlexOrder[v];
468 
469         revlexOrder[v] = revlexOrder[w];
470         revlexOrder[w] = swap;
471       }
472     }
473   }
474   ierr = VecRestoreArrayRead(coordVec, &coords);CHKERRQ(ierr);
475   ierr = PetscMalloc1(ni->nodeIdxDim * nVerts, &newNodeIdx);CHKERRQ(ierr);
476   /* reorder nodeIdx to be in closure order */
477   for (v = 0; v < nVerts; v++) {
478     for (d = 0; d < ni->nodeIdxDim; d++) {
479       newNodeIdx[revlexOrder[v] * ni->nodeIdxDim + d] = ni->nodeIdx[v * ni->nodeIdxDim + d];
480     }
481   }
482   ierr = PetscFree(ni->nodeIdx);CHKERRQ(ierr);
483   ni->nodeIdx = newNodeIdx;
484   ni->perm = invClosureOrder;
485   ierr = PetscFree(revlexOrder);CHKERRQ(ierr);
486   ierr = PetscFree(closureOrder);CHKERRQ(ierr);
487   PetscFunctionReturn(0);
488 }
489 
490 /* the coordinates of the simplex vertices are the corners of the barycentric simplex.
491  * When we stack them on top of each other in revlex order, they look like the identity matrix */
492 static PetscErrorCode PetscLagNodeIndicesCreateSimplexVertices(DM dm, PetscLagNodeIndices *nodeIndices)
493 {
494   PetscLagNodeIndices ni;
495   PetscInt       dim, d;
496 
497   PetscErrorCode ierr;
498 
499   PetscFunctionBegin;
500   ierr = PetscNew(&ni);CHKERRQ(ierr);
501   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
502   ni->nodeIdxDim = dim + 1;
503   ni->nodeVecDim = 0;
504   ni->nNodes = dim + 1;
505   ni->refct = 1;
506   ierr = PetscCalloc1((dim + 1)*(dim + 1), &(ni->nodeIdx));CHKERRQ(ierr);
507   for (d = 0; d < dim + 1; d++) ni->nodeIdx[d*(dim + 2)] = 1;
508   ierr = PetscLagNodeIndicesComputeVertexOrder(dm, ni, PETSC_FALSE);CHKERRQ(ierr);
509   *nodeIndices = ni;
510   PetscFunctionReturn(0);
511 }
512 
513 /* A polytope that is a tensor product of a facet and a segment.
514  * We take whatever coordinate system was being used for the facet
515  * and we concatenate the barycentric coordinates for the vertices
516  * at the end of the segment, (1,0) and (0,1), to get a coordinate
517  * system for the tensor product element */
518 static PetscErrorCode PetscLagNodeIndicesCreateTensorVertices(DM dm, PetscLagNodeIndices facetni, PetscLagNodeIndices *nodeIndices)
519 {
520   PetscLagNodeIndices ni;
521   PetscInt       nodeIdxDim, subNodeIdxDim = facetni->nodeIdxDim;
522   PetscInt       nVerts, nSubVerts = facetni->nNodes;
523   PetscInt       dim, d, e, f, g;
524 
525   PetscErrorCode ierr;
526 
527   PetscFunctionBegin;
528   ierr = PetscNew(&ni);CHKERRQ(ierr);
529   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
530   ni->nodeIdxDim = nodeIdxDim = subNodeIdxDim + 2;
531   ni->nodeVecDim = 0;
532   ni->nNodes = nVerts = 2 * nSubVerts;
533   ni->refct = 1;
534   ierr = PetscCalloc1(nodeIdxDim * nVerts, &(ni->nodeIdx));CHKERRQ(ierr);
535   for (f = 0, d = 0; d < 2; d++) {
536     for (e = 0; e < nSubVerts; e++, f++) {
537       for (g = 0; g < subNodeIdxDim; g++) {
538         ni->nodeIdx[f * nodeIdxDim + g] = facetni->nodeIdx[e * subNodeIdxDim + g];
539       }
540       ni->nodeIdx[f * nodeIdxDim + subNodeIdxDim] = (1 - d);
541       ni->nodeIdx[f * nodeIdxDim + subNodeIdxDim + 1] = d;
542     }
543   }
544   ierr = PetscLagNodeIndicesComputeVertexOrder(dm, ni, PETSC_TRUE);CHKERRQ(ierr);
545   *nodeIndices = ni;
546   PetscFunctionReturn(0);
547 }
548 
549 /* This helps us compute symmetries, and it also helps us compute coordinates for dofs that are being pushed
550  * forward from a boundary mesh point.
551  *
552  * Input:
553  *
554  * dm - the target reference cell where we want new coordinates and dof directions to be valid
555  * vert - the vertex coordinate system for the target reference cell
556  * p - the point in the target reference cell that the dofs are coming from
557  * vertp - the vertex coordinate system for p's reference cell
558  * ornt - the resulting coordinates and dof vectors will be for p under this orientation
559  * nodep - the node coordinates and dof vectors in p's reference cell
560  * formDegree - the form degree that the dofs transform as
561  *
562  * Output:
563  *
564  * pfNodeIdx - the node coordinates for p's dofs, in the dm reference cell, from the ornt perspective
565  * pfNodeVec - the node dof vectors for p's dofs, in the dm reference cell, from the ornt perspective
566  */
567 static PetscErrorCode PetscLagNodeIndicesPushForward(DM dm, PetscLagNodeIndices vert, PetscInt p, PetscLagNodeIndices vertp, PetscLagNodeIndices nodep, PetscInt ornt, PetscInt formDegree, PetscInt pfNodeIdx[], PetscReal pfNodeVec[])
568 {
569   PetscInt       *closureVerts;
570   PetscInt        closureSize = 0;
571   PetscInt       *closure = NULL;
572   PetscInt        dim, pdim, c, i, j, k, n, v, vStart, vEnd;
573   PetscInt        nSubVert = vertp->nNodes;
574   PetscInt        nodeIdxDim = vert->nodeIdxDim;
575   PetscInt        subNodeIdxDim = vertp->nodeIdxDim;
576   PetscInt        nNodes = nodep->nNodes;
577   const PetscInt  *vertIdx = vert->nodeIdx;
578   const PetscInt  *subVertIdx = vertp->nodeIdx;
579   const PetscInt  *nodeIdx = nodep->nodeIdx;
580   const PetscReal *nodeVec = nodep->nodeVec;
581   PetscReal       *J, *Jstar;
582   PetscReal       detJ;
583   PetscInt        depth, pdepth, Nk, pNk;
584   Vec             coordVec;
585   PetscScalar      *newCoords = NULL;
586   const PetscScalar *oldCoords = NULL;
587   PetscErrorCode  ierr;
588 
589   PetscFunctionBegin;
590   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
591   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
592   ierr = DMGetCoordinatesLocal(dm, &coordVec);CHKERRQ(ierr);
593   ierr = DMPlexGetPointDepth(dm, p, &pdepth);CHKERRQ(ierr);
594   pdim = pdepth != depth ? pdepth != 0 ? pdepth : 0 : dim;
595   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
596   ierr = DMGetWorkArray(dm, nSubVert, MPIU_INT, &closureVerts);CHKERRQ(ierr);
597   ierr = DMPlexGetTransitiveClosure_Internal(dm, p, ornt, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
598   c = closureSize - nSubVert;
599   /* we want which cell closure indices the closure of this point corresponds to */
600   for (v = 0; v < nSubVert; v++) closureVerts[v] = vert->perm[closure[2 * (c + v)] - vStart];
601   ierr = DMPlexRestoreTransitiveClosure(dm, p, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
602   /* push forward indices */
603   for (i = 0; i < nodeIdxDim; i++) { /* for every component of the target index space */
604     /* check if this is a component that all vertices around this point have in common */
605     for (j = 1; j < nSubVert; j++) {
606       if (vertIdx[closureVerts[j] * nodeIdxDim + i] != vertIdx[closureVerts[0] * nodeIdxDim + i]) break;
607     }
608     if (j == nSubVert) { /* all vertices have this component in common, directly copy to output */
609       PetscInt val = vertIdx[closureVerts[0] * nodeIdxDim + i];
610       for (n = 0; n < nNodes; n++) pfNodeIdx[n * nodeIdxDim + i] = val;
611     } else {
612       PetscInt subi = -1;
613       /* there must be a component in vertp that looks the same */
614       for (k = 0; k < subNodeIdxDim; k++) {
615         for (j = 0; j < nSubVert; j++) {
616           if (vertIdx[closureVerts[j] * nodeIdxDim + i] != subVertIdx[j * subNodeIdxDim + k]) break;
617         }
618         if (j == nSubVert) {
619           subi = k;
620           break;
621         }
622       }
623       PetscCheckFalse(subi < 0,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Did not find matching coordinate");
624       /* that component in the vertp system becomes component i in the vert system for each dof */
625       for (n = 0; n < nNodes; n++) pfNodeIdx[n * nodeIdxDim + i] = nodeIdx[n * subNodeIdxDim + subi];
626     }
627   }
628   /* push forward vectors */
629   ierr = DMGetWorkArray(dm, dim * dim, MPIU_REAL, &J);CHKERRQ(ierr);
630   if (ornt != 0) { /* temporarily change the coordinate vector so
631                       DMPlexComputeCellGeometryAffineFEM gives us the Jacobian we want */
632     PetscInt        closureSize2 = 0;
633     PetscInt       *closure2 = NULL;
634 
635     ierr = DMPlexGetTransitiveClosure_Internal(dm, p, 0, PETSC_TRUE, &closureSize2, &closure2);CHKERRQ(ierr);
636     ierr = PetscMalloc1(dim * nSubVert, &newCoords);CHKERRQ(ierr);
637     ierr = VecGetArrayRead(coordVec, &oldCoords);CHKERRQ(ierr);
638     for (v = 0; v < nSubVert; v++) {
639       PetscInt d;
640       for (d = 0; d < dim; d++) {
641         newCoords[(closure2[2 * (c + v)] - vStart) * dim + d] = oldCoords[closureVerts[v] * dim + d];
642       }
643     }
644     ierr = VecRestoreArrayRead(coordVec, &oldCoords);CHKERRQ(ierr);
645     ierr = DMPlexRestoreTransitiveClosure(dm, p, PETSC_TRUE, &closureSize2, &closure2);CHKERRQ(ierr);
646     ierr = VecPlaceArray(coordVec, newCoords);CHKERRQ(ierr);
647   }
648   ierr = DMPlexComputeCellGeometryAffineFEM(dm, p, NULL, J, NULL, &detJ);CHKERRQ(ierr);
649   if (ornt != 0) {
650     ierr = VecResetArray(coordVec);CHKERRQ(ierr);
651     ierr = PetscFree(newCoords);CHKERRQ(ierr);
652   }
653   ierr = DMRestoreWorkArray(dm, nSubVert, MPIU_INT, &closureVerts);CHKERRQ(ierr);
654   /* compactify */
655   for (i = 0; i < dim; i++) for (j = 0; j < pdim; j++) J[i * pdim + j] = J[i * dim + j];
656   /* We have the Jacobian mapping the point's reference cell to this reference cell:
657    * pulling back a function to the point and applying the dof is what we want,
658    * so we get the pullback matrix and multiply the dof by that matrix on the right */
659   ierr = PetscDTBinomialInt(dim, PetscAbsInt(formDegree), &Nk);CHKERRQ(ierr);
660   ierr = PetscDTBinomialInt(pdim, PetscAbsInt(formDegree), &pNk);CHKERRQ(ierr);
661   ierr = DMGetWorkArray(dm, pNk * Nk, MPIU_REAL, &Jstar);CHKERRQ(ierr);
662   ierr = PetscDTAltVPullbackMatrix(pdim, dim, J, formDegree, Jstar);CHKERRQ(ierr);
663   for (n = 0; n < nNodes; n++) {
664     for (i = 0; i < Nk; i++) {
665       PetscReal val = 0.;
666       for (j = 0; j < pNk; j++) val += nodeVec[n * pNk + j] * Jstar[j * Nk + i];
667       pfNodeVec[n * Nk + i] = val;
668     }
669   }
670   ierr = DMRestoreWorkArray(dm, pNk * Nk, MPIU_REAL, &Jstar);CHKERRQ(ierr);
671   ierr = DMRestoreWorkArray(dm, dim * dim, MPIU_REAL, &J);CHKERRQ(ierr);
672   PetscFunctionReturn(0);
673 }
674 
675 /* given to sets of nodes, take the tensor product, where the product of the dof indices is concatenation and the
676  * product of the dof vectors is the wedge product */
677 static PetscErrorCode PetscLagNodeIndicesTensor(PetscLagNodeIndices tracei, PetscInt dimT, PetscInt kT, PetscLagNodeIndices fiberi, PetscInt dimF, PetscInt kF, PetscLagNodeIndices *nodeIndices)
678 {
679   PetscInt       dim = dimT + dimF;
680   PetscInt       nodeIdxDim, nNodes;
681   PetscInt       formDegree = kT + kF;
682   PetscInt       Nk, NkT, NkF;
683   PetscInt       MkT, MkF;
684   PetscLagNodeIndices ni;
685   PetscInt       i, j, l;
686   PetscReal      *projF, *projT;
687   PetscReal      *projFstar, *projTstar;
688   PetscReal      *workF, *workF2, *workT, *workT2, *work, *work2;
689   PetscReal      *wedgeMat;
690   PetscReal      sign;
691   PetscErrorCode ierr;
692 
693   PetscFunctionBegin;
694   ierr = PetscDTBinomialInt(dim, PetscAbsInt(formDegree), &Nk);CHKERRQ(ierr);
695   ierr = PetscDTBinomialInt(dimT, PetscAbsInt(kT), &NkT);CHKERRQ(ierr);
696   ierr = PetscDTBinomialInt(dimF, PetscAbsInt(kF), &NkF);CHKERRQ(ierr);
697   ierr = PetscDTBinomialInt(dim, PetscAbsInt(kT), &MkT);CHKERRQ(ierr);
698   ierr = PetscDTBinomialInt(dim, PetscAbsInt(kF), &MkF);CHKERRQ(ierr);
699   ierr = PetscNew(&ni);CHKERRQ(ierr);
700   ni->nodeIdxDim = nodeIdxDim = tracei->nodeIdxDim + fiberi->nodeIdxDim;
701   ni->nodeVecDim = Nk;
702   ni->nNodes = nNodes = tracei->nNodes * fiberi->nNodes;
703   ni->refct = 1;
704   ierr = PetscMalloc1(nNodes * nodeIdxDim, &(ni->nodeIdx));CHKERRQ(ierr);
705   /* first concatenate the indices */
706   for (l = 0, j = 0; j < fiberi->nNodes; j++) {
707     for (i = 0; i < tracei->nNodes; i++, l++) {
708       PetscInt m, n = 0;
709 
710       for (m = 0; m < tracei->nodeIdxDim; m++) ni->nodeIdx[l * nodeIdxDim + n++] = tracei->nodeIdx[i * tracei->nodeIdxDim + m];
711       for (m = 0; m < fiberi->nodeIdxDim; m++) ni->nodeIdx[l * nodeIdxDim + n++] = fiberi->nodeIdx[j * fiberi->nodeIdxDim + m];
712     }
713   }
714 
715   /* now wedge together the push-forward vectors */
716   ierr = PetscMalloc1(nNodes * Nk, &(ni->nodeVec));CHKERRQ(ierr);
717   ierr = PetscCalloc2(dimT*dim, &projT, dimF*dim, &projF);CHKERRQ(ierr);
718   for (i = 0; i < dimT; i++) projT[i * (dim + 1)] = 1.;
719   for (i = 0; i < dimF; i++) projF[i * (dim + dimT + 1) + dimT] = 1.;
720   ierr = PetscMalloc2(MkT*NkT, &projTstar, MkF*NkF, &projFstar);CHKERRQ(ierr);
721   ierr = PetscDTAltVPullbackMatrix(dim, dimT, projT, kT, projTstar);CHKERRQ(ierr);
722   ierr = PetscDTAltVPullbackMatrix(dim, dimF, projF, kF, projFstar);CHKERRQ(ierr);
723   ierr = PetscMalloc6(MkT, &workT, MkT, &workT2, MkF, &workF, MkF, &workF2, Nk, &work, Nk, &work2);CHKERRQ(ierr);
724   ierr = PetscMalloc1(Nk * MkT, &wedgeMat);CHKERRQ(ierr);
725   sign = (PetscAbsInt(kT * kF) & 1) ? -1. : 1.;
726   for (l = 0, j = 0; j < fiberi->nNodes; j++) {
727     PetscInt d, e;
728 
729     /* push forward fiber k-form */
730     for (d = 0; d < MkF; d++) {
731       PetscReal val = 0.;
732       for (e = 0; e < NkF; e++) val += projFstar[d * NkF + e] * fiberi->nodeVec[j * NkF + e];
733       workF[d] = val;
734     }
735     /* Hodge star to proper form if necessary */
736     if (kF < 0) {
737       for (d = 0; d < MkF; d++) workF2[d] = workF[d];
738       ierr = PetscDTAltVStar(dim, PetscAbsInt(kF), 1, workF2, workF);CHKERRQ(ierr);
739     }
740     /* Compute the matrix that wedges this form with one of the trace k-form */
741     ierr = PetscDTAltVWedgeMatrix(dim, PetscAbsInt(kF), PetscAbsInt(kT), workF, wedgeMat);CHKERRQ(ierr);
742     for (i = 0; i < tracei->nNodes; i++, l++) {
743       /* push forward trace k-form */
744       for (d = 0; d < MkT; d++) {
745         PetscReal val = 0.;
746         for (e = 0; e < NkT; e++) val += projTstar[d * NkT + e] * tracei->nodeVec[i * NkT + e];
747         workT[d] = val;
748       }
749       /* Hodge star to proper form if necessary */
750       if (kT < 0) {
751         for (d = 0; d < MkT; d++) workT2[d] = workT[d];
752         ierr = PetscDTAltVStar(dim, PetscAbsInt(kT), 1, workT2, workT);CHKERRQ(ierr);
753       }
754       /* compute the wedge product of the push-forward trace form and firer forms */
755       for (d = 0; d < Nk; d++) {
756         PetscReal val = 0.;
757         for (e = 0; e < MkT; e++) val += wedgeMat[d * MkT + e] * workT[e];
758         work[d] = val;
759       }
760       /* inverse Hodge star from proper form if necessary */
761       if (formDegree < 0) {
762         for (d = 0; d < Nk; d++) work2[d] = work[d];
763         ierr = PetscDTAltVStar(dim, PetscAbsInt(formDegree), -1, work2, work);CHKERRQ(ierr);
764       }
765       /* insert into the array (adjusting for sign) */
766       for (d = 0; d < Nk; d++) ni->nodeVec[l * Nk + d] = sign * work[d];
767     }
768   }
769   ierr = PetscFree(wedgeMat);CHKERRQ(ierr);
770   ierr = PetscFree6(workT, workT2, workF, workF2, work, work2);CHKERRQ(ierr);
771   ierr = PetscFree2(projTstar, projFstar);CHKERRQ(ierr);
772   ierr = PetscFree2(projT, projF);CHKERRQ(ierr);
773   *nodeIndices = ni;
774   PetscFunctionReturn(0);
775 }
776 
777 /* simple union of two sets of nodes */
778 static PetscErrorCode PetscLagNodeIndicesMerge(PetscLagNodeIndices niA, PetscLagNodeIndices niB, PetscLagNodeIndices *nodeIndices)
779 {
780   PetscLagNodeIndices ni;
781   PetscInt            nodeIdxDim, nodeVecDim, nNodes;
782   PetscErrorCode      ierr;
783 
784   PetscFunctionBegin;
785   ierr = PetscNew(&ni);CHKERRQ(ierr);
786   ni->nodeIdxDim = nodeIdxDim = niA->nodeIdxDim;
787   PetscCheckFalse(niB->nodeIdxDim != nodeIdxDim,PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Cannot merge PetscLagNodeIndices with different nodeIdxDim");
788   ni->nodeVecDim = nodeVecDim = niA->nodeVecDim;
789   PetscCheckFalse(niB->nodeVecDim != nodeVecDim,PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Cannot merge PetscLagNodeIndices with different nodeVecDim");
790   ni->nNodes = nNodes = niA->nNodes + niB->nNodes;
791   ni->refct = 1;
792   ierr = PetscMalloc1(nNodes * nodeIdxDim, &(ni->nodeIdx));CHKERRQ(ierr);
793   ierr = PetscMalloc1(nNodes * nodeVecDim, &(ni->nodeVec));CHKERRQ(ierr);
794   ierr = PetscArraycpy(ni->nodeIdx, niA->nodeIdx, niA->nNodes * nodeIdxDim);CHKERRQ(ierr);
795   ierr = PetscArraycpy(ni->nodeVec, niA->nodeVec, niA->nNodes * nodeVecDim);CHKERRQ(ierr);
796   ierr = PetscArraycpy(&(ni->nodeIdx[niA->nNodes * nodeIdxDim]), niB->nodeIdx, niB->nNodes * nodeIdxDim);CHKERRQ(ierr);
797   ierr = PetscArraycpy(&(ni->nodeVec[niA->nNodes * nodeVecDim]), niB->nodeVec, niB->nNodes * nodeVecDim);CHKERRQ(ierr);
798   *nodeIndices = ni;
799   PetscFunctionReturn(0);
800 }
801 
802 #define PETSCTUPINTCOMPREVLEX(N)                                   \
803 static int PetscConcat_(PetscTupIntCompRevlex_,N)(const void *a, const void *b) \
804 {                                                                  \
805   const PetscInt *A = (const PetscInt *) a;                        \
806   const PetscInt *B = (const PetscInt *) b;                        \
807   int i;                                                           \
808   PetscInt diff = 0;                                               \
809   for (i = 0; i < N; i++) {                                        \
810     diff = A[N - i] - B[N - i];                                    \
811     if (diff) break;                                               \
812   }                                                                \
813   return (diff <= 0) ? (diff < 0) ? -1 : 0 : 1;                    \
814 }
815 
816 PETSCTUPINTCOMPREVLEX(3)
817 PETSCTUPINTCOMPREVLEX(4)
818 PETSCTUPINTCOMPREVLEX(5)
819 PETSCTUPINTCOMPREVLEX(6)
820 PETSCTUPINTCOMPREVLEX(7)
821 
822 static int PetscTupIntCompRevlex_N(const void *a, const void *b)
823 {
824   const PetscInt *A = (const PetscInt *) a;
825   const PetscInt *B = (const PetscInt *) b;
826   int i;
827   int N = A[0];
828   PetscInt diff = 0;
829   for (i = 0; i < N; i++) {
830     diff = A[N - i] - B[N - i];
831     if (diff) break;
832   }
833   return (diff <= 0) ? (diff < 0) ? -1 : 0 : 1;
834 }
835 
836 /* The nodes are not necessarily in revlex order wrt nodeIdx: get the permutation
837  * that puts them in that order */
838 static PetscErrorCode PetscLagNodeIndicesGetPermutation(PetscLagNodeIndices ni, PetscInt *perm[])
839 {
840   PetscErrorCode ierr;
841 
842   PetscFunctionBegin;
843   if (!(ni->perm)) {
844     PetscInt *sorter;
845     PetscInt m = ni->nNodes;
846     PetscInt nodeIdxDim = ni->nodeIdxDim;
847     PetscInt i, j, k, l;
848     PetscInt *prm;
849     int (*comp) (const void *, const void *);
850 
851     ierr = PetscMalloc1((nodeIdxDim + 2) * m, &sorter);CHKERRQ(ierr);
852     for (k = 0, l = 0, i = 0; i < m; i++) {
853       sorter[k++] = nodeIdxDim + 1;
854       sorter[k++] = i;
855       for (j = 0; j < nodeIdxDim; j++) sorter[k++] = ni->nodeIdx[l++];
856     }
857     switch (nodeIdxDim) {
858     case 2:
859       comp = PetscTupIntCompRevlex_3;
860       break;
861     case 3:
862       comp = PetscTupIntCompRevlex_4;
863       break;
864     case 4:
865       comp = PetscTupIntCompRevlex_5;
866       break;
867     case 5:
868       comp = PetscTupIntCompRevlex_6;
869       break;
870     case 6:
871       comp = PetscTupIntCompRevlex_7;
872       break;
873     default:
874       comp = PetscTupIntCompRevlex_N;
875       break;
876     }
877     qsort(sorter, m, (nodeIdxDim + 2) * sizeof(PetscInt), comp);
878     ierr = PetscMalloc1(m, &prm);CHKERRQ(ierr);
879     for (i = 0; i < m; i++) prm[i] = sorter[(nodeIdxDim + 2) * i + 1];
880     ni->perm = prm;
881     ierr = PetscFree(sorter);CHKERRQ(ierr);
882   }
883   *perm = ni->perm;
884   PetscFunctionReturn(0);
885 }
886 
887 static PetscErrorCode PetscDualSpaceDestroy_Lagrange(PetscDualSpace sp)
888 {
889   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data;
890   PetscErrorCode      ierr;
891 
892   PetscFunctionBegin;
893   if (lag->symperms) {
894     PetscInt **selfSyms = lag->symperms[0];
895 
896     if (selfSyms) {
897       PetscInt i, **allocated = &selfSyms[-lag->selfSymOff];
898 
899       for (i = 0; i < lag->numSelfSym; i++) {
900         ierr = PetscFree(allocated[i]);CHKERRQ(ierr);
901       }
902       ierr = PetscFree(allocated);CHKERRQ(ierr);
903     }
904     ierr = PetscFree(lag->symperms);CHKERRQ(ierr);
905   }
906   if (lag->symflips) {
907     PetscScalar **selfSyms = lag->symflips[0];
908 
909     if (selfSyms) {
910       PetscInt i;
911       PetscScalar **allocated = &selfSyms[-lag->selfSymOff];
912 
913       for (i = 0; i < lag->numSelfSym; i++) {
914         ierr = PetscFree(allocated[i]);CHKERRQ(ierr);
915       }
916       ierr = PetscFree(allocated);CHKERRQ(ierr);
917     }
918     ierr = PetscFree(lag->symflips);CHKERRQ(ierr);
919   }
920   ierr = Petsc1DNodeFamilyDestroy(&(lag->nodeFamily));CHKERRQ(ierr);
921   ierr = PetscLagNodeIndicesDestroy(&(lag->vertIndices));CHKERRQ(ierr);
922   ierr = PetscLagNodeIndicesDestroy(&(lag->intNodeIndices));CHKERRQ(ierr);
923   ierr = PetscLagNodeIndicesDestroy(&(lag->allNodeIndices));CHKERRQ(ierr);
924   ierr = PetscFree(lag);CHKERRQ(ierr);
925   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetContinuity_C", NULL);CHKERRQ(ierr);
926   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetContinuity_C", NULL);CHKERRQ(ierr);
927   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetTensor_C", NULL);CHKERRQ(ierr);
928   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetTensor_C", NULL);CHKERRQ(ierr);
929   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetTrimmed_C", NULL);CHKERRQ(ierr);
930   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetTrimmed_C", NULL);CHKERRQ(ierr);
931   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetNodeType_C", NULL);CHKERRQ(ierr);
932   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetNodeType_C", NULL);CHKERRQ(ierr);
933   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetUseMoments_C", NULL);CHKERRQ(ierr);
934   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetUseMoments_C", NULL);CHKERRQ(ierr);
935   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetMomentOrder_C", NULL);CHKERRQ(ierr);
936   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetMomentOrder_C", NULL);CHKERRQ(ierr);
937   PetscFunctionReturn(0);
938 }
939 
940 static PetscErrorCode PetscDualSpaceLagrangeView_Ascii(PetscDualSpace sp, PetscViewer viewer)
941 {
942   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data;
943   PetscErrorCode      ierr;
944 
945   PetscFunctionBegin;
946   ierr = PetscViewerASCIIPrintf(viewer, "%s %s%sLagrange dual space\n", lag->continuous ? "Continuous" : "Discontinuous", lag->tensorSpace ? "tensor " : "", lag->trimmed ? "trimmed " : "");CHKERRQ(ierr);
947   PetscFunctionReturn(0);
948 }
949 
950 static PetscErrorCode PetscDualSpaceView_Lagrange(PetscDualSpace sp, PetscViewer viewer)
951 {
952   PetscBool      iascii;
953   PetscErrorCode ierr;
954 
955   PetscFunctionBegin;
956   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
957   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
958   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
959   if (iascii) {ierr = PetscDualSpaceLagrangeView_Ascii(sp, viewer);CHKERRQ(ierr);}
960   PetscFunctionReturn(0);
961 }
962 
963 static PetscErrorCode PetscDualSpaceSetFromOptions_Lagrange(PetscOptionItems *PetscOptionsObject,PetscDualSpace sp)
964 {
965   PetscBool      continuous, tensor, trimmed, flg, flg2, flg3;
966   PetscDTNodeType nodeType;
967   PetscReal      nodeExponent;
968   PetscInt       momentOrder;
969   PetscBool      nodeEndpoints, useMoments;
970   PetscErrorCode ierr;
971 
972   PetscFunctionBegin;
973   ierr = PetscDualSpaceLagrangeGetContinuity(sp, &continuous);CHKERRQ(ierr);
974   ierr = PetscDualSpaceLagrangeGetTensor(sp, &tensor);CHKERRQ(ierr);
975   ierr = PetscDualSpaceLagrangeGetTrimmed(sp, &trimmed);CHKERRQ(ierr);
976   ierr = PetscDualSpaceLagrangeGetNodeType(sp, &nodeType, &nodeEndpoints, &nodeExponent);CHKERRQ(ierr);
977   if (nodeType == PETSCDTNODES_DEFAULT) nodeType = PETSCDTNODES_GAUSSJACOBI;
978   ierr = PetscDualSpaceLagrangeGetUseMoments(sp, &useMoments);CHKERRQ(ierr);
979   ierr = PetscDualSpaceLagrangeGetMomentOrder(sp, &momentOrder);CHKERRQ(ierr);
980   ierr = PetscOptionsHead(PetscOptionsObject,"PetscDualSpace Lagrange Options");CHKERRQ(ierr);
981   ierr = PetscOptionsBool("-petscdualspace_lagrange_continuity", "Flag for continuous element", "PetscDualSpaceLagrangeSetContinuity", continuous, &continuous, &flg);CHKERRQ(ierr);
982   if (flg) {ierr = PetscDualSpaceLagrangeSetContinuity(sp, continuous);CHKERRQ(ierr);}
983   ierr = PetscOptionsBool("-petscdualspace_lagrange_tensor", "Flag for tensor dual space", "PetscDualSpaceLagrangeSetTensor", tensor, &tensor, &flg);CHKERRQ(ierr);
984   if (flg) {ierr = PetscDualSpaceLagrangeSetTensor(sp, tensor);CHKERRQ(ierr);}
985   ierr = PetscOptionsBool("-petscdualspace_lagrange_trimmed", "Flag for trimmed dual space", "PetscDualSpaceLagrangeSetTrimmed", trimmed, &trimmed, &flg);CHKERRQ(ierr);
986   if (flg) {ierr = PetscDualSpaceLagrangeSetTrimmed(sp, trimmed);CHKERRQ(ierr);}
987   ierr = PetscOptionsEnum("-petscdualspace_lagrange_node_type", "Lagrange node location type", "PetscDualSpaceLagrangeSetNodeType", PetscDTNodeTypes, (PetscEnum)nodeType, (PetscEnum *)&nodeType, &flg);CHKERRQ(ierr);
988   ierr = PetscOptionsBool("-petscdualspace_lagrange_node_endpoints", "Flag for nodes that include endpoints", "PetscDualSpaceLagrangeSetNodeType", nodeEndpoints, &nodeEndpoints, &flg2);CHKERRQ(ierr);
989   flg3 = PETSC_FALSE;
990   if (nodeType == PETSCDTNODES_GAUSSJACOBI) {
991     ierr = PetscOptionsReal("-petscdualspace_lagrange_node_exponent", "Gauss-Jacobi weight function exponent", "PetscDualSpaceLagrangeSetNodeType", nodeExponent, &nodeExponent, &flg3);CHKERRQ(ierr);
992   }
993   if (flg || flg2 || flg3) {ierr = PetscDualSpaceLagrangeSetNodeType(sp, nodeType, nodeEndpoints, nodeExponent);CHKERRQ(ierr);}
994   ierr = PetscOptionsBool("-petscdualspace_lagrange_use_moments", "Use moments (where appropriate) for functionals", "PetscDualSpaceLagrangeSetUseMoments", useMoments, &useMoments, &flg);CHKERRQ(ierr);
995   if (flg) {ierr = PetscDualSpaceLagrangeSetUseMoments(sp, useMoments);CHKERRQ(ierr);}
996   ierr = PetscOptionsInt("-petscdualspace_lagrange_moment_order", "Quadrature order for moment functionals", "PetscDualSpaceLagrangeSetMomentOrder", momentOrder, &momentOrder, &flg);CHKERRQ(ierr);
997   if (flg) {ierr = PetscDualSpaceLagrangeSetMomentOrder(sp, momentOrder);CHKERRQ(ierr);}
998   ierr = PetscOptionsTail();CHKERRQ(ierr);
999   PetscFunctionReturn(0);
1000 }
1001 
1002 static PetscErrorCode PetscDualSpaceDuplicate_Lagrange(PetscDualSpace sp, PetscDualSpace spNew)
1003 {
1004   PetscBool           cont, tensor, trimmed, boundary;
1005   PetscDTNodeType     nodeType;
1006   PetscReal           exponent;
1007   PetscDualSpace_Lag *lag    = (PetscDualSpace_Lag *) sp->data;
1008   PetscErrorCode      ierr;
1009 
1010   PetscFunctionBegin;
1011   ierr = PetscDualSpaceLagrangeGetContinuity(sp, &cont);CHKERRQ(ierr);
1012   ierr = PetscDualSpaceLagrangeSetContinuity(spNew, cont);CHKERRQ(ierr);
1013   ierr = PetscDualSpaceLagrangeGetTensor(sp, &tensor);CHKERRQ(ierr);
1014   ierr = PetscDualSpaceLagrangeSetTensor(spNew, tensor);CHKERRQ(ierr);
1015   ierr = PetscDualSpaceLagrangeGetTrimmed(sp, &trimmed);CHKERRQ(ierr);
1016   ierr = PetscDualSpaceLagrangeSetTrimmed(spNew, trimmed);CHKERRQ(ierr);
1017   ierr = PetscDualSpaceLagrangeGetNodeType(sp, &nodeType, &boundary, &exponent);CHKERRQ(ierr);
1018   ierr = PetscDualSpaceLagrangeSetNodeType(spNew, nodeType, boundary, exponent);CHKERRQ(ierr);
1019   if (lag->nodeFamily) {
1020     PetscDualSpace_Lag *lagnew = (PetscDualSpace_Lag *) spNew->data;
1021 
1022     ierr = Petsc1DNodeFamilyReference(lag->nodeFamily);CHKERRQ(ierr);
1023     lagnew->nodeFamily = lag->nodeFamily;
1024   }
1025   PetscFunctionReturn(0);
1026 }
1027 
1028 /* for making tensor product spaces: take a dual space and product a segment space that has all the same
1029  * specifications (trimmed, continuous, order, node set), except for the form degree */
1030 static PetscErrorCode PetscDualSpaceCreateEdgeSubspace_Lagrange(PetscDualSpace sp, PetscInt order, PetscInt k, PetscInt Nc, PetscBool interiorOnly, PetscDualSpace *bdsp)
1031 {
1032   DM                 K;
1033   PetscDualSpace_Lag *newlag;
1034   PetscErrorCode     ierr;
1035 
1036   PetscFunctionBegin;
1037   ierr = PetscDualSpaceDuplicate(sp,bdsp);CHKERRQ(ierr);
1038   ierr = PetscDualSpaceSetFormDegree(*bdsp, k);CHKERRQ(ierr);
1039   ierr = DMPlexCreateReferenceCell(PETSC_COMM_SELF, DMPolytopeTypeSimpleShape(1, PETSC_TRUE), &K);CHKERRQ(ierr);
1040   ierr = PetscDualSpaceSetDM(*bdsp, K);CHKERRQ(ierr);
1041   ierr = DMDestroy(&K);CHKERRQ(ierr);
1042   ierr = PetscDualSpaceSetOrder(*bdsp, order);CHKERRQ(ierr);
1043   ierr = PetscDualSpaceSetNumComponents(*bdsp, Nc);CHKERRQ(ierr);
1044   newlag = (PetscDualSpace_Lag *) (*bdsp)->data;
1045   newlag->interiorOnly = interiorOnly;
1046   ierr = PetscDualSpaceSetUp(*bdsp);CHKERRQ(ierr);
1047   PetscFunctionReturn(0);
1048 }
1049 
1050 /* just the points, weights aren't handled */
1051 static PetscErrorCode PetscQuadratureCreateTensor(PetscQuadrature trace, PetscQuadrature fiber, PetscQuadrature *product)
1052 {
1053   PetscInt         dimTrace, dimFiber;
1054   PetscInt         numPointsTrace, numPointsFiber;
1055   PetscInt         dim, numPoints;
1056   const PetscReal *pointsTrace;
1057   const PetscReal *pointsFiber;
1058   PetscReal       *points;
1059   PetscInt         i, j, k, p;
1060   PetscErrorCode   ierr;
1061 
1062   PetscFunctionBegin;
1063   ierr = PetscQuadratureGetData(trace, &dimTrace, NULL, &numPointsTrace, &pointsTrace, NULL);CHKERRQ(ierr);
1064   ierr = PetscQuadratureGetData(fiber, &dimFiber, NULL, &numPointsFiber, &pointsFiber, NULL);CHKERRQ(ierr);
1065   dim = dimTrace + dimFiber;
1066   numPoints = numPointsFiber * numPointsTrace;
1067   ierr = PetscMalloc1(numPoints * dim, &points);CHKERRQ(ierr);
1068   for (p = 0, j = 0; j < numPointsFiber; j++) {
1069     for (i = 0; i < numPointsTrace; i++, p++) {
1070       for (k = 0; k < dimTrace; k++) points[p * dim +            k] = pointsTrace[i * dimTrace + k];
1071       for (k = 0; k < dimFiber; k++) points[p * dim + dimTrace + k] = pointsFiber[j * dimFiber + k];
1072     }
1073   }
1074   ierr = PetscQuadratureCreate(PETSC_COMM_SELF, product);CHKERRQ(ierr);
1075   ierr = PetscQuadratureSetData(*product, dim, 0, numPoints, points, NULL);CHKERRQ(ierr);
1076   PetscFunctionReturn(0);
1077 }
1078 
1079 /* Kronecker tensor product where matrix is considered a matrix of k-forms, so that
1080  * the entries in the product matrix are wedge products of the entries in the original matrices */
1081 static PetscErrorCode MatTensorAltV(Mat trace, Mat fiber, PetscInt dimTrace, PetscInt kTrace, PetscInt dimFiber, PetscInt kFiber, Mat *product)
1082 {
1083   PetscInt mTrace, nTrace, mFiber, nFiber, m, n, k, i, j, l;
1084   PetscInt dim, NkTrace, NkFiber, Nk;
1085   PetscInt dT, dF;
1086   PetscInt *nnzTrace, *nnzFiber, *nnz;
1087   PetscInt iT, iF, jT, jF, il, jl;
1088   PetscReal *workT, *workT2, *workF, *workF2, *work, *workstar;
1089   PetscReal *projT, *projF;
1090   PetscReal *projTstar, *projFstar;
1091   PetscReal *wedgeMat;
1092   PetscReal sign;
1093   PetscScalar *workS;
1094   Mat prod;
1095   /* this produces dof groups that look like the identity */
1096   PetscErrorCode ierr;
1097 
1098   PetscFunctionBegin;
1099   ierr = MatGetSize(trace, &mTrace, &nTrace);CHKERRQ(ierr);
1100   ierr = PetscDTBinomialInt(dimTrace, PetscAbsInt(kTrace), &NkTrace);CHKERRQ(ierr);
1101   PetscCheckFalse(nTrace % NkTrace,PETSC_COMM_SELF, PETSC_ERR_PLIB, "point value space of trace matrix is not a multiple of k-form size");
1102   ierr = MatGetSize(fiber, &mFiber, &nFiber);CHKERRQ(ierr);
1103   ierr = PetscDTBinomialInt(dimFiber, PetscAbsInt(kFiber), &NkFiber);CHKERRQ(ierr);
1104   PetscCheckFalse(nFiber % NkFiber,PETSC_COMM_SELF, PETSC_ERR_PLIB, "point value space of fiber matrix is not a multiple of k-form size");
1105   ierr = PetscMalloc2(mTrace, &nnzTrace, mFiber, &nnzFiber);CHKERRQ(ierr);
1106   for (i = 0; i < mTrace; i++) {
1107     ierr = MatGetRow(trace, i, &(nnzTrace[i]), NULL, NULL);CHKERRQ(ierr);
1108     PetscCheckFalse(nnzTrace[i] % NkTrace,PETSC_COMM_SELF, PETSC_ERR_PLIB, "nonzeros in trace matrix are not in k-form size blocks");
1109   }
1110   for (i = 0; i < mFiber; i++) {
1111     ierr = MatGetRow(fiber, i, &(nnzFiber[i]), NULL, NULL);CHKERRQ(ierr);
1112     PetscCheckFalse(nnzFiber[i] % NkFiber,PETSC_COMM_SELF, PETSC_ERR_PLIB, "nonzeros in fiber matrix are not in k-form size blocks");
1113   }
1114   dim = dimTrace + dimFiber;
1115   k = kFiber + kTrace;
1116   ierr = PetscDTBinomialInt(dim, PetscAbsInt(k), &Nk);CHKERRQ(ierr);
1117   m = mTrace * mFiber;
1118   ierr = PetscMalloc1(m, &nnz);CHKERRQ(ierr);
1119   for (l = 0, j = 0; j < mFiber; j++) for (i = 0; i < mTrace; i++, l++) nnz[l] = (nnzTrace[i] / NkTrace) * (nnzFiber[j] / NkFiber) * Nk;
1120   n = (nTrace / NkTrace) * (nFiber / NkFiber) * Nk;
1121   ierr = MatCreateSeqAIJ(PETSC_COMM_SELF, m, n, 0, nnz, &prod);CHKERRQ(ierr);
1122   ierr = PetscFree(nnz);CHKERRQ(ierr);
1123   ierr = PetscFree2(nnzTrace,nnzFiber);CHKERRQ(ierr);
1124   /* reasoning about which points each dof needs depends on having zeros computed at points preserved */
1125   ierr = MatSetOption(prod, MAT_IGNORE_ZERO_ENTRIES, PETSC_FALSE);CHKERRQ(ierr);
1126   /* compute pullbacks */
1127   ierr = PetscDTBinomialInt(dim, PetscAbsInt(kTrace), &dT);CHKERRQ(ierr);
1128   ierr = PetscDTBinomialInt(dim, PetscAbsInt(kFiber), &dF);CHKERRQ(ierr);
1129   ierr = PetscMalloc4(dimTrace * dim, &projT, dimFiber * dim, &projF, dT * NkTrace, &projTstar, dF * NkFiber, &projFstar);CHKERRQ(ierr);
1130   ierr = PetscArrayzero(projT, dimTrace * dim);CHKERRQ(ierr);
1131   for (i = 0; i < dimTrace; i++) projT[i * (dim + 1)] = 1.;
1132   ierr = PetscArrayzero(projF, dimFiber * dim);CHKERRQ(ierr);
1133   for (i = 0; i < dimFiber; i++) projF[i * (dim + 1) + dimTrace] = 1.;
1134   ierr = PetscDTAltVPullbackMatrix(dim, dimTrace, projT, kTrace, projTstar);CHKERRQ(ierr);
1135   ierr = PetscDTAltVPullbackMatrix(dim, dimFiber, projF, kFiber, projFstar);CHKERRQ(ierr);
1136   ierr = PetscMalloc5(dT, &workT, dF, &workF, Nk, &work, Nk, &workstar, Nk, &workS);CHKERRQ(ierr);
1137   ierr = PetscMalloc2(dT, &workT2, dF, &workF2);CHKERRQ(ierr);
1138   ierr = PetscMalloc1(Nk * dT, &wedgeMat);CHKERRQ(ierr);
1139   sign = (PetscAbsInt(kTrace * kFiber) & 1) ? -1. : 1.;
1140   for (i = 0, iF = 0; iF < mFiber; iF++) {
1141     PetscInt           ncolsF, nformsF;
1142     const PetscInt    *colsF;
1143     const PetscScalar *valsF;
1144 
1145     ierr = MatGetRow(fiber, iF, &ncolsF, &colsF, &valsF);CHKERRQ(ierr);
1146     nformsF = ncolsF / NkFiber;
1147     for (iT = 0; iT < mTrace; iT++, i++) {
1148       PetscInt           ncolsT, nformsT;
1149       const PetscInt    *colsT;
1150       const PetscScalar *valsT;
1151 
1152       ierr = MatGetRow(trace, iT, &ncolsT, &colsT, &valsT);CHKERRQ(ierr);
1153       nformsT = ncolsT / NkTrace;
1154       for (j = 0, jF = 0; jF < nformsF; jF++) {
1155         PetscInt colF = colsF[jF * NkFiber] / NkFiber;
1156 
1157         for (il = 0; il < dF; il++) {
1158           PetscReal val = 0.;
1159           for (jl = 0; jl < NkFiber; jl++) val += projFstar[il * NkFiber + jl] * PetscRealPart(valsF[jF * NkFiber + jl]);
1160           workF[il] = val;
1161         }
1162         if (kFiber < 0) {
1163           for (il = 0; il < dF; il++) workF2[il] = workF[il];
1164           ierr = PetscDTAltVStar(dim, PetscAbsInt(kFiber), 1, workF2, workF);CHKERRQ(ierr);
1165         }
1166         ierr = PetscDTAltVWedgeMatrix(dim, PetscAbsInt(kFiber), PetscAbsInt(kTrace), workF, wedgeMat);CHKERRQ(ierr);
1167         for (jT = 0; jT < nformsT; jT++, j++) {
1168           PetscInt colT = colsT[jT * NkTrace] / NkTrace;
1169           PetscInt col = colF * (nTrace / NkTrace) + colT;
1170           const PetscScalar *vals;
1171 
1172           for (il = 0; il < dT; il++) {
1173             PetscReal val = 0.;
1174             for (jl = 0; jl < NkTrace; jl++) val += projTstar[il * NkTrace + jl] * PetscRealPart(valsT[jT * NkTrace + jl]);
1175             workT[il] = val;
1176           }
1177           if (kTrace < 0) {
1178             for (il = 0; il < dT; il++) workT2[il] = workT[il];
1179             ierr = PetscDTAltVStar(dim, PetscAbsInt(kTrace), 1, workT2, workT);CHKERRQ(ierr);
1180           }
1181 
1182           for (il = 0; il < Nk; il++) {
1183             PetscReal val = 0.;
1184             for (jl = 0; jl < dT; jl++) val += sign * wedgeMat[il * dT + jl] * workT[jl];
1185             work[il] = val;
1186           }
1187           if (k < 0) {
1188             ierr = PetscDTAltVStar(dim, PetscAbsInt(k), -1, work, workstar);CHKERRQ(ierr);
1189 #if defined(PETSC_USE_COMPLEX)
1190             for (l = 0; l < Nk; l++) workS[l] = workstar[l];
1191             vals = &workS[0];
1192 #else
1193             vals = &workstar[0];
1194 #endif
1195           } else {
1196 #if defined(PETSC_USE_COMPLEX)
1197             for (l = 0; l < Nk; l++) workS[l] = work[l];
1198             vals = &workS[0];
1199 #else
1200             vals = &work[0];
1201 #endif
1202           }
1203           for (l = 0; l < Nk; l++) {
1204             ierr = MatSetValue(prod, i, col * Nk + l, vals[l], INSERT_VALUES);CHKERRQ(ierr);
1205           } /* Nk */
1206         } /* jT */
1207       } /* jF */
1208       ierr = MatRestoreRow(trace, iT, &ncolsT, &colsT, &valsT);CHKERRQ(ierr);
1209     } /* iT */
1210     ierr = MatRestoreRow(fiber, iF, &ncolsF, &colsF, &valsF);CHKERRQ(ierr);
1211   } /* iF */
1212   ierr = PetscFree(wedgeMat);CHKERRQ(ierr);
1213   ierr = PetscFree4(projT, projF, projTstar, projFstar);CHKERRQ(ierr);
1214   ierr = PetscFree2(workT2, workF2);CHKERRQ(ierr);
1215   ierr = PetscFree5(workT, workF, work, workstar, workS);CHKERRQ(ierr);
1216   ierr = MatAssemblyBegin(prod, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1217   ierr = MatAssemblyEnd(prod, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1218   *product = prod;
1219   PetscFunctionReturn(0);
1220 }
1221 
1222 /* Union of quadrature points, with an attempt to identify commont points in the two sets */
1223 static PetscErrorCode PetscQuadraturePointsMerge(PetscQuadrature quadA, PetscQuadrature quadB, PetscQuadrature *quadJoint, PetscInt *aToJoint[], PetscInt *bToJoint[])
1224 {
1225   PetscInt         dimA, dimB;
1226   PetscInt         nA, nB, nJoint, i, j, d;
1227   const PetscReal *pointsA;
1228   const PetscReal *pointsB;
1229   PetscReal       *pointsJoint;
1230   PetscInt        *aToJ, *bToJ;
1231   PetscQuadrature  qJ;
1232   PetscErrorCode   ierr;
1233 
1234   PetscFunctionBegin;
1235   ierr = PetscQuadratureGetData(quadA, &dimA, NULL, &nA, &pointsA, NULL);CHKERRQ(ierr);
1236   ierr = PetscQuadratureGetData(quadB, &dimB, NULL, &nB, &pointsB, NULL);CHKERRQ(ierr);
1237   PetscCheckFalse(dimA != dimB,PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Quadrature points must be in the same dimension");
1238   nJoint = nA;
1239   ierr = PetscMalloc1(nA, &aToJ);CHKERRQ(ierr);
1240   for (i = 0; i < nA; i++) aToJ[i] = i;
1241   ierr = PetscMalloc1(nB, &bToJ);CHKERRQ(ierr);
1242   for (i = 0; i < nB; i++) {
1243     for (j = 0; j < nA; j++) {
1244       bToJ[i] = -1;
1245       for (d = 0; d < dimA; d++) if (PetscAbsReal(pointsB[i * dimA + d] - pointsA[j * dimA + d]) > PETSC_SMALL) break;
1246       if (d == dimA) {
1247         bToJ[i] = j;
1248         break;
1249       }
1250     }
1251     if (bToJ[i] == -1) {
1252       bToJ[i] = nJoint++;
1253     }
1254   }
1255   *aToJoint = aToJ;
1256   *bToJoint = bToJ;
1257   ierr = PetscMalloc1(nJoint * dimA, &pointsJoint);CHKERRQ(ierr);
1258   ierr = PetscArraycpy(pointsJoint, pointsA, nA * dimA);CHKERRQ(ierr);
1259   for (i = 0; i < nB; i++) {
1260     if (bToJ[i] >= nA) {
1261       for (d = 0; d < dimA; d++) pointsJoint[bToJ[i] * dimA + d] = pointsB[i * dimA + d];
1262     }
1263   }
1264   ierr = PetscQuadratureCreate(PETSC_COMM_SELF, &qJ);CHKERRQ(ierr);
1265   ierr = PetscQuadratureSetData(qJ, dimA, 0, nJoint, pointsJoint, NULL);CHKERRQ(ierr);
1266   *quadJoint = qJ;
1267   PetscFunctionReturn(0);
1268 }
1269 
1270 /* Matrices matA and matB are both quadrature -> dof matrices: produce a matrix that is joint quadrature -> union of
1271  * dofs, where the joint quadrature was produced by PetscQuadraturePointsMerge */
1272 static PetscErrorCode MatricesMerge(Mat matA, Mat matB, PetscInt dim, PetscInt k, PetscInt numMerged, const PetscInt aToMerged[], const PetscInt bToMerged[], Mat *matMerged)
1273 {
1274   PetscInt m, n, mA, nA, mB, nB, Nk, i, j, l;
1275   Mat      M;
1276   PetscInt *nnz;
1277   PetscInt maxnnz;
1278   PetscInt *work;
1279   PetscErrorCode ierr;
1280 
1281   PetscFunctionBegin;
1282   ierr = PetscDTBinomialInt(dim, PetscAbsInt(k), &Nk);CHKERRQ(ierr);
1283   ierr = MatGetSize(matA, &mA, &nA);CHKERRQ(ierr);
1284   PetscCheckFalse(nA % Nk,PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "matA column space not a multiple of k-form size");
1285   ierr = MatGetSize(matB, &mB, &nB);CHKERRQ(ierr);
1286   PetscCheckFalse(nB % Nk,PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "matB column space not a multiple of k-form size");
1287   m = mA + mB;
1288   n = numMerged * Nk;
1289   ierr = PetscMalloc1(m, &nnz);CHKERRQ(ierr);
1290   maxnnz = 0;
1291   for (i = 0; i < mA; i++) {
1292     ierr = MatGetRow(matA, i, &(nnz[i]), NULL, NULL);CHKERRQ(ierr);
1293     PetscCheckFalse(nnz[i] % Nk,PETSC_COMM_SELF, PETSC_ERR_PLIB, "nonzeros in matA are not in k-form size blocks");
1294     maxnnz = PetscMax(maxnnz, nnz[i]);
1295   }
1296   for (i = 0; i < mB; i++) {
1297     ierr = MatGetRow(matB, i, &(nnz[i+mA]), NULL, NULL);CHKERRQ(ierr);
1298     PetscCheckFalse(nnz[i+mA] % Nk,PETSC_COMM_SELF, PETSC_ERR_PLIB, "nonzeros in matB are not in k-form size blocks");
1299     maxnnz = PetscMax(maxnnz, nnz[i+mA]);
1300   }
1301   ierr = MatCreateSeqAIJ(PETSC_COMM_SELF, m, n, 0, nnz, &M);CHKERRQ(ierr);
1302   ierr = PetscFree(nnz);CHKERRQ(ierr);
1303   /* reasoning about which points each dof needs depends on having zeros computed at points preserved */
1304   ierr = MatSetOption(M, MAT_IGNORE_ZERO_ENTRIES, PETSC_FALSE);CHKERRQ(ierr);
1305   ierr = PetscMalloc1(maxnnz, &work);CHKERRQ(ierr);
1306   for (i = 0; i < mA; i++) {
1307     const PetscInt *cols;
1308     const PetscScalar *vals;
1309     PetscInt nCols;
1310     ierr = MatGetRow(matA, i, &nCols, &cols, &vals);CHKERRQ(ierr);
1311     for (j = 0; j < nCols / Nk; j++) {
1312       PetscInt newCol = aToMerged[cols[j * Nk] / Nk];
1313       for (l = 0; l < Nk; l++) work[j * Nk + l] = newCol * Nk + l;
1314     }
1315     ierr = MatSetValuesBlocked(M, 1, &i, nCols, work, vals, INSERT_VALUES);CHKERRQ(ierr);
1316     ierr = MatRestoreRow(matA, i, &nCols, &cols, &vals);CHKERRQ(ierr);
1317   }
1318   for (i = 0; i < mB; i++) {
1319     const PetscInt *cols;
1320     const PetscScalar *vals;
1321 
1322     PetscInt row = i + mA;
1323     PetscInt nCols;
1324     ierr = MatGetRow(matB, i, &nCols, &cols, &vals);CHKERRQ(ierr);
1325     for (j = 0; j < nCols / Nk; j++) {
1326       PetscInt newCol = bToMerged[cols[j * Nk] / Nk];
1327       for (l = 0; l < Nk; l++) work[j * Nk + l] = newCol * Nk + l;
1328     }
1329     ierr = MatSetValuesBlocked(M, 1, &row, nCols, work, vals, INSERT_VALUES);CHKERRQ(ierr);
1330     ierr = MatRestoreRow(matB, i, &nCols, &cols, &vals);CHKERRQ(ierr);
1331   }
1332   ierr = PetscFree(work);CHKERRQ(ierr);
1333   ierr = MatAssemblyBegin(M, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1334   ierr = MatAssemblyEnd(M, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1335   *matMerged = M;
1336   PetscFunctionReturn(0);
1337 }
1338 
1339 /* Take a dual space and product a segment space that has all the same specifications (trimmed, continuous, order,
1340  * node set), except for the form degree.  For computing boundary dofs and for making tensor product spaces */
1341 static PetscErrorCode PetscDualSpaceCreateFacetSubspace_Lagrange(PetscDualSpace sp, DM K, PetscInt f, PetscInt k, PetscInt Ncopies, PetscBool interiorOnly, PetscDualSpace *bdsp)
1342 {
1343   PetscInt           Nknew, Ncnew;
1344   PetscInt           dim, pointDim = -1;
1345   PetscInt           depth;
1346   DM                 dm;
1347   PetscDualSpace_Lag *newlag;
1348   PetscErrorCode     ierr;
1349 
1350   PetscFunctionBegin;
1351   ierr = PetscDualSpaceGetDM(sp,&dm);CHKERRQ(ierr);
1352   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
1353   ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr);
1354   ierr = PetscDualSpaceDuplicate(sp,bdsp);CHKERRQ(ierr);
1355   ierr = PetscDualSpaceSetFormDegree(*bdsp,k);CHKERRQ(ierr);
1356   if (!K) {
1357     PetscBool isSimplex;
1358 
1359     if (depth == dim) {
1360       DMPolytopeType ct;
1361 
1362       pointDim = dim - 1;
1363       ierr = DMPlexGetCellType(dm, f, &ct);CHKERRQ(ierr);
1364       ierr = DMPlexCreateReferenceCell(PETSC_COMM_SELF, ct, &K);CHKERRQ(ierr);
1365     } else if (depth == 1) {
1366       pointDim = 0;
1367       ierr = DMPlexCreateReferenceCell(PETSC_COMM_SELF, DM_POLYTOPE_POINT, &K);CHKERRQ(ierr);
1368     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unsupported interpolation state of reference element");
1369   } else {
1370     ierr = PetscObjectReference((PetscObject)K);CHKERRQ(ierr);
1371     ierr = DMGetDimension(K, &pointDim);CHKERRQ(ierr);
1372   }
1373   ierr = PetscDualSpaceSetDM(*bdsp, K);CHKERRQ(ierr);
1374   ierr = DMDestroy(&K);CHKERRQ(ierr);
1375   ierr = PetscDTBinomialInt(pointDim, PetscAbsInt(k), &Nknew);CHKERRQ(ierr);
1376   Ncnew = Nknew * Ncopies;
1377   ierr = PetscDualSpaceSetNumComponents(*bdsp, Ncnew);CHKERRQ(ierr);
1378   newlag = (PetscDualSpace_Lag *) (*bdsp)->data;
1379   newlag->interiorOnly = interiorOnly;
1380   ierr = PetscDualSpaceSetUp(*bdsp);CHKERRQ(ierr);
1381   PetscFunctionReturn(0);
1382 }
1383 
1384 /* Construct simplex nodes from a nodefamily, add Nk dof vectors of length Nk at each node.
1385  * Return the (quadrature, matrix) form of the dofs and the nodeIndices form as well.
1386  *
1387  * Sometimes we want a set of nodes to be contained in the interior of the element,
1388  * even when the node scheme puts nodes on the boundaries.  numNodeSkip tells
1389  * the routine how many "layers" of nodes need to be skipped.
1390  * */
1391 static PetscErrorCode PetscDualSpaceLagrangeCreateSimplexNodeMat(Petsc1DNodeFamily nodeFamily, PetscInt dim, PetscInt sum, PetscInt Nk, PetscInt numNodeSkip, PetscQuadrature *iNodes, Mat *iMat, PetscLagNodeIndices *nodeIndices)
1392 {
1393   PetscReal *extraNodeCoords, *nodeCoords;
1394   PetscInt nNodes, nExtraNodes;
1395   PetscInt i, j, k, extraSum = sum + numNodeSkip * (1 + dim);
1396   PetscQuadrature intNodes;
1397   Mat intMat;
1398   PetscLagNodeIndices ni;
1399   PetscErrorCode ierr;
1400 
1401   PetscFunctionBegin;
1402   ierr = PetscDTBinomialInt(dim + sum, dim, &nNodes);CHKERRQ(ierr);
1403   ierr = PetscDTBinomialInt(dim + extraSum, dim, &nExtraNodes);CHKERRQ(ierr);
1404 
1405   ierr = PetscMalloc1(dim * nExtraNodes, &extraNodeCoords);CHKERRQ(ierr);
1406   ierr = PetscNew(&ni);CHKERRQ(ierr);
1407   ni->nodeIdxDim = dim + 1;
1408   ni->nodeVecDim = Nk;
1409   ni->nNodes = nNodes * Nk;
1410   ni->refct = 1;
1411   ierr = PetscMalloc1(nNodes * Nk * (dim + 1), &(ni->nodeIdx));CHKERRQ(ierr);
1412   ierr = PetscMalloc1(nNodes * Nk * Nk, &(ni->nodeVec));CHKERRQ(ierr);
1413   for (i = 0; i < nNodes; i++) for (j = 0; j < Nk; j++) for (k = 0; k < Nk; k++) ni->nodeVec[(i * Nk + j) * Nk + k] = (j == k) ? 1. : 0.;
1414   ierr = Petsc1DNodeFamilyComputeSimplexNodes(nodeFamily, dim, extraSum, extraNodeCoords);CHKERRQ(ierr);
1415   if (numNodeSkip) {
1416     PetscInt k;
1417     PetscInt *tup;
1418 
1419     ierr = PetscMalloc1(dim * nNodes, &nodeCoords);CHKERRQ(ierr);
1420     ierr = PetscMalloc1(dim + 1, &tup);CHKERRQ(ierr);
1421     for (k = 0; k < nNodes; k++) {
1422       PetscInt j, c;
1423       PetscInt index;
1424 
1425       ierr = PetscDTIndexToBary(dim + 1, sum, k, tup);CHKERRQ(ierr);
1426       for (j = 0; j < dim + 1; j++) tup[j] += numNodeSkip;
1427       for (c = 0; c < Nk; c++) {
1428         for (j = 0; j < dim + 1; j++) {
1429           ni->nodeIdx[(k * Nk + c) * (dim + 1) + j] = tup[j] + 1;
1430         }
1431       }
1432       ierr = PetscDTBaryToIndex(dim + 1, extraSum, tup, &index);CHKERRQ(ierr);
1433       for (j = 0; j < dim; j++) nodeCoords[k * dim + j] = extraNodeCoords[index * dim + j];
1434     }
1435     ierr = PetscFree(tup);CHKERRQ(ierr);
1436     ierr = PetscFree(extraNodeCoords);CHKERRQ(ierr);
1437   } else {
1438     PetscInt k;
1439     PetscInt *tup;
1440 
1441     nodeCoords = extraNodeCoords;
1442     ierr = PetscMalloc1(dim + 1, &tup);CHKERRQ(ierr);
1443     for (k = 0; k < nNodes; k++) {
1444       PetscInt j, c;
1445 
1446       ierr = PetscDTIndexToBary(dim + 1, sum, k, tup);CHKERRQ(ierr);
1447       for (c = 0; c < Nk; c++) {
1448         for (j = 0; j < dim + 1; j++) {
1449           /* barycentric indices can have zeros, but we don't want to push forward zeros because it makes it harder to
1450            * determine which nodes correspond to which under symmetries, so we increase by 1.  This is fine
1451            * because the nodeIdx coordinates don't have any meaning other than helping to identify symmetries */
1452           ni->nodeIdx[(k * Nk + c) * (dim + 1) + j] = tup[j] + 1;
1453         }
1454       }
1455     }
1456     ierr = PetscFree(tup);CHKERRQ(ierr);
1457   }
1458   ierr = PetscQuadratureCreate(PETSC_COMM_SELF, &intNodes);CHKERRQ(ierr);
1459   ierr = PetscQuadratureSetData(intNodes, dim, 0, nNodes, nodeCoords, NULL);CHKERRQ(ierr);
1460   ierr = MatCreateSeqAIJ(PETSC_COMM_SELF, nNodes * Nk, nNodes * Nk, Nk, NULL, &intMat);CHKERRQ(ierr);
1461   ierr = MatSetOption(intMat,MAT_IGNORE_ZERO_ENTRIES,PETSC_FALSE);CHKERRQ(ierr);
1462   for (j = 0; j < nNodes * Nk; j++) {
1463     PetscInt rem = j % Nk;
1464     PetscInt a, aprev = j - rem;
1465     PetscInt anext = aprev + Nk;
1466 
1467     for (a = aprev; a < anext; a++) {
1468       ierr = MatSetValue(intMat, j, a, (a == j) ? 1. : 0., INSERT_VALUES);CHKERRQ(ierr);
1469     }
1470   }
1471   ierr = MatAssemblyBegin(intMat, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1472   ierr = MatAssemblyEnd(intMat, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1473   *iNodes = intNodes;
1474   *iMat = intMat;
1475   *nodeIndices = ni;
1476   PetscFunctionReturn(0);
1477 }
1478 
1479 /* once the nodeIndices have been created for the interior of the reference cell, and for all of the boundary cells,
1480  * push forward the boundary dofs and concatenate them into the full node indices for the dual space */
1481 static PetscErrorCode PetscDualSpaceLagrangeCreateAllNodeIdx(PetscDualSpace sp)
1482 {
1483   DM             dm;
1484   PetscInt       dim, nDofs;
1485   PetscSection   section;
1486   PetscInt       pStart, pEnd, p;
1487   PetscInt       formDegree, Nk;
1488   PetscInt       nodeIdxDim, spintdim;
1489   PetscDualSpace_Lag *lag;
1490   PetscLagNodeIndices ni, verti;
1491   PetscErrorCode ierr;
1492 
1493   PetscFunctionBegin;
1494   lag = (PetscDualSpace_Lag *) sp->data;
1495   verti = lag->vertIndices;
1496   ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr);
1497   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1498   ierr = PetscDualSpaceGetFormDegree(sp, &formDegree);CHKERRQ(ierr);
1499   ierr = PetscDTBinomialInt(dim, PetscAbsInt(formDegree), &Nk);CHKERRQ(ierr);
1500   ierr = PetscDualSpaceGetSection(sp, &section);CHKERRQ(ierr);
1501   ierr = PetscSectionGetStorageSize(section, &nDofs);CHKERRQ(ierr);
1502   ierr = PetscNew(&ni);CHKERRQ(ierr);
1503   ni->nodeIdxDim = nodeIdxDim = verti->nodeIdxDim;
1504   ni->nodeVecDim = Nk;
1505   ni->nNodes = nDofs;
1506   ni->refct = 1;
1507   ierr = PetscMalloc1(nodeIdxDim * nDofs, &(ni->nodeIdx));CHKERRQ(ierr);
1508   ierr = PetscMalloc1(Nk * nDofs, &(ni->nodeVec));CHKERRQ(ierr);
1509   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
1510   ierr = PetscSectionGetDof(section, 0, &spintdim);CHKERRQ(ierr);
1511   if (spintdim) {
1512     ierr = PetscArraycpy(ni->nodeIdx, lag->intNodeIndices->nodeIdx, spintdim * nodeIdxDim);CHKERRQ(ierr);
1513     ierr = PetscArraycpy(ni->nodeVec, lag->intNodeIndices->nodeVec, spintdim * Nk);CHKERRQ(ierr);
1514   }
1515   for (p = pStart + 1; p < pEnd; p++) {
1516     PetscDualSpace psp = sp->pointSpaces[p];
1517     PetscDualSpace_Lag *plag;
1518     PetscInt dof, off;
1519 
1520     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
1521     if (!dof) continue;
1522     plag = (PetscDualSpace_Lag *) psp->data;
1523     ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
1524     ierr = PetscLagNodeIndicesPushForward(dm, verti, p, plag->vertIndices, plag->intNodeIndices, 0, formDegree, &(ni->nodeIdx[off * nodeIdxDim]), &(ni->nodeVec[off * Nk]));CHKERRQ(ierr);
1525   }
1526   lag->allNodeIndices = ni;
1527   PetscFunctionReturn(0);
1528 }
1529 
1530 /* once the (quadrature, Matrix) forms of the dofs have been created for the interior of the
1531  * reference cell and for the boundary cells, jk
1532  * push forward the boundary data and concatenate them into the full (quadrature, matrix) data
1533  * for the dual space */
1534 static PetscErrorCode PetscDualSpaceCreateAllDataFromInteriorData(PetscDualSpace sp)
1535 {
1536   DM               dm;
1537   PetscSection     section;
1538   PetscInt         pStart, pEnd, p, k, Nk, dim, Nc;
1539   PetscInt         nNodes;
1540   PetscInt         countNodes;
1541   Mat              allMat;
1542   PetscQuadrature  allNodes;
1543   PetscInt         nDofs;
1544   PetscInt         maxNzforms, j;
1545   PetscScalar      *work;
1546   PetscReal        *L, *J, *Jinv, *v0, *pv0;
1547   PetscInt         *iwork;
1548   PetscReal        *nodes;
1549   PetscErrorCode   ierr;
1550 
1551   PetscFunctionBegin;
1552   ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr);
1553   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1554   ierr = PetscDualSpaceGetSection(sp, &section);CHKERRQ(ierr);
1555   ierr = PetscSectionGetStorageSize(section, &nDofs);CHKERRQ(ierr);
1556   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
1557   ierr = PetscDualSpaceGetFormDegree(sp, &k);CHKERRQ(ierr);
1558   ierr = PetscDualSpaceGetNumComponents(sp, &Nc);CHKERRQ(ierr);
1559   ierr = PetscDTBinomialInt(dim, PetscAbsInt(k), &Nk);CHKERRQ(ierr);
1560   for (p = pStart, nNodes = 0, maxNzforms = 0; p < pEnd; p++) {
1561     PetscDualSpace  psp;
1562     DM              pdm;
1563     PetscInt        pdim, pNk;
1564     PetscQuadrature intNodes;
1565     Mat intMat;
1566 
1567     ierr = PetscDualSpaceGetPointSubspace(sp, p, &psp);CHKERRQ(ierr);
1568     if (!psp) continue;
1569     ierr = PetscDualSpaceGetDM(psp, &pdm);CHKERRQ(ierr);
1570     ierr = DMGetDimension(pdm, &pdim);CHKERRQ(ierr);
1571     if (pdim < PetscAbsInt(k)) continue;
1572     ierr = PetscDTBinomialInt(pdim, PetscAbsInt(k), &pNk);CHKERRQ(ierr);
1573     ierr = PetscDualSpaceGetInteriorData(psp, &intNodes, &intMat);CHKERRQ(ierr);
1574     if (intNodes) {
1575       PetscInt nNodesp;
1576 
1577       ierr = PetscQuadratureGetData(intNodes, NULL, NULL, &nNodesp, NULL, NULL);CHKERRQ(ierr);
1578       nNodes += nNodesp;
1579     }
1580     if (intMat) {
1581       PetscInt maxNzsp;
1582       PetscInt maxNzformsp;
1583 
1584       ierr = MatSeqAIJGetMaxRowNonzeros(intMat, &maxNzsp);CHKERRQ(ierr);
1585       PetscCheckFalse(maxNzsp % pNk,PETSC_COMM_SELF, PETSC_ERR_PLIB, "interior matrix is not laid out as blocks of k-forms");
1586       maxNzformsp = maxNzsp / pNk;
1587       maxNzforms = PetscMax(maxNzforms, maxNzformsp);
1588     }
1589   }
1590   ierr = MatCreateSeqAIJ(PETSC_COMM_SELF, nDofs, nNodes * Nc, maxNzforms * Nk, NULL, &allMat);CHKERRQ(ierr);
1591   ierr = MatSetOption(allMat,MAT_IGNORE_ZERO_ENTRIES,PETSC_FALSE);CHKERRQ(ierr);
1592   ierr = PetscMalloc7(dim, &v0, dim, &pv0, dim * dim, &J, dim * dim, &Jinv, Nk * Nk, &L, maxNzforms * Nk, &work, maxNzforms * Nk, &iwork);CHKERRQ(ierr);
1593   for (j = 0; j < dim; j++) pv0[j] = -1.;
1594   ierr = PetscMalloc1(dim * nNodes, &nodes);CHKERRQ(ierr);
1595   for (p = pStart, countNodes = 0; p < pEnd; p++) {
1596     PetscDualSpace  psp;
1597     PetscQuadrature intNodes;
1598     DM pdm;
1599     PetscInt pdim, pNk;
1600     PetscInt countNodesIn = countNodes;
1601     PetscReal detJ;
1602     Mat intMat;
1603 
1604     ierr = PetscDualSpaceGetPointSubspace(sp, p, &psp);CHKERRQ(ierr);
1605     if (!psp) continue;
1606     ierr = PetscDualSpaceGetDM(psp, &pdm);CHKERRQ(ierr);
1607     ierr = DMGetDimension(pdm, &pdim);CHKERRQ(ierr);
1608     if (pdim < PetscAbsInt(k)) continue;
1609     ierr = PetscDualSpaceGetInteriorData(psp, &intNodes, &intMat);CHKERRQ(ierr);
1610     if (intNodes == NULL && intMat == NULL) continue;
1611     ierr = PetscDTBinomialInt(pdim, PetscAbsInt(k), &pNk);CHKERRQ(ierr);
1612     if (p) {
1613       ierr = DMPlexComputeCellGeometryAffineFEM(dm, p, v0, J, Jinv, &detJ);CHKERRQ(ierr);
1614     } else { /* identity */
1615       PetscInt i,j;
1616 
1617       for (i = 0; i < dim; i++) for (j = 0; j < dim; j++) J[i * dim + j] = Jinv[i * dim + j] = 0.;
1618       for (i = 0; i < dim; i++) J[i * dim + i] = Jinv[i * dim + i] = 1.;
1619       for (i = 0; i < dim; i++) v0[i] = -1.;
1620     }
1621     if (pdim != dim) { /* compactify Jacobian */
1622       PetscInt i, j;
1623 
1624       for (i = 0; i < dim; i++) for (j = 0; j < pdim; j++) J[i * pdim + j] = J[i * dim + j];
1625     }
1626     ierr = PetscDTAltVPullbackMatrix(pdim, dim, J, k, L);CHKERRQ(ierr);
1627     if (intNodes) { /* push forward quadrature locations by the affine transformation */
1628       PetscInt nNodesp;
1629       const PetscReal *nodesp;
1630       PetscInt j;
1631 
1632       ierr = PetscQuadratureGetData(intNodes, NULL, NULL, &nNodesp, &nodesp, NULL);CHKERRQ(ierr);
1633       for (j = 0; j < nNodesp; j++, countNodes++) {
1634         PetscInt d, e;
1635 
1636         for (d = 0; d < dim; d++) {
1637           nodes[countNodes * dim + d] = v0[d];
1638           for (e = 0; e < pdim; e++) {
1639             nodes[countNodes * dim + d] += J[d * pdim + e] * (nodesp[j * pdim + e] - pv0[e]);
1640           }
1641         }
1642       }
1643     }
1644     if (intMat) {
1645       PetscInt nrows;
1646       PetscInt off;
1647 
1648       ierr = PetscSectionGetDof(section, p, &nrows);CHKERRQ(ierr);
1649       ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
1650       for (j = 0; j < nrows; j++) {
1651         PetscInt ncols;
1652         const PetscInt *cols;
1653         const PetscScalar *vals;
1654         PetscInt l, d, e;
1655         PetscInt row = j + off;
1656 
1657         ierr = MatGetRow(intMat, j, &ncols, &cols, &vals);CHKERRQ(ierr);
1658         PetscCheckFalse(ncols % pNk,PETSC_COMM_SELF, PETSC_ERR_PLIB, "interior matrix is not laid out as blocks of k-forms");
1659         for (l = 0; l < ncols / pNk; l++) {
1660           PetscInt blockcol;
1661 
1662           for (d = 0; d < pNk; d++) {
1663             PetscCheckFalse((cols[l * pNk + d] % pNk) != d,PETSC_COMM_SELF, PETSC_ERR_PLIB, "interior matrix is not laid out as blocks of k-forms");
1664           }
1665           blockcol = cols[l * pNk] / pNk;
1666           for (d = 0; d < Nk; d++) {
1667             iwork[l * Nk + d] = (blockcol + countNodesIn) * Nk + d;
1668           }
1669           for (d = 0; d < Nk; d++) work[l * Nk + d] = 0.;
1670           for (d = 0; d < Nk; d++) {
1671             for (e = 0; e < pNk; e++) {
1672               /* "push forward" dof by pulling back a k-form to be evaluated on the point: multiply on the right by L */
1673               work[l * Nk + d] += vals[l * pNk + e] * L[e * Nk + d];
1674             }
1675           }
1676         }
1677         ierr = MatSetValues(allMat, 1, &row, (ncols / pNk) * Nk, iwork, work, INSERT_VALUES);CHKERRQ(ierr);
1678         ierr = MatRestoreRow(intMat, j, &ncols, &cols, &vals);CHKERRQ(ierr);
1679       }
1680     }
1681   }
1682   ierr = MatAssemblyBegin(allMat, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1683   ierr = MatAssemblyEnd(allMat, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1684   ierr = PetscQuadratureCreate(PETSC_COMM_SELF, &allNodes);CHKERRQ(ierr);
1685   ierr = PetscQuadratureSetData(allNodes, dim, 0, nNodes, nodes, NULL);CHKERRQ(ierr);
1686   ierr = PetscFree7(v0, pv0, J, Jinv, L, work, iwork);CHKERRQ(ierr);
1687   ierr = MatDestroy(&(sp->allMat));CHKERRQ(ierr);
1688   sp->allMat = allMat;
1689   ierr = PetscQuadratureDestroy(&(sp->allNodes));CHKERRQ(ierr);
1690   sp->allNodes = allNodes;
1691   PetscFunctionReturn(0);
1692 }
1693 
1694 /* rather than trying to get all data from the functionals, we create
1695  * the functionals from rows of the quadrature -> dof matrix.
1696  *
1697  * Ideally most of the uses of PetscDualSpace in PetscFE will switch
1698  * to using intMat and allMat, so that the individual functionals
1699  * don't need to be constructed at all */
1700 static PetscErrorCode PetscDualSpaceComputeFunctionalsFromAllData(PetscDualSpace sp)
1701 {
1702   PetscQuadrature allNodes;
1703   Mat             allMat;
1704   PetscInt        nDofs;
1705   PetscInt        dim, k, Nk, Nc, f;
1706   DM              dm;
1707   PetscInt        nNodes, spdim;
1708   const PetscReal *nodes = NULL;
1709   PetscSection    section;
1710   PetscBool       useMoments;
1711   PetscErrorCode  ierr;
1712 
1713   PetscFunctionBegin;
1714   ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr);
1715   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1716   ierr = PetscDualSpaceGetNumComponents(sp, &Nc);CHKERRQ(ierr);
1717   ierr = PetscDualSpaceGetFormDegree(sp, &k);CHKERRQ(ierr);
1718   ierr = PetscDTBinomialInt(dim, PetscAbsInt(k), &Nk);CHKERRQ(ierr);
1719   ierr = PetscDualSpaceGetAllData(sp, &allNodes, &allMat);CHKERRQ(ierr);
1720   nNodes = 0;
1721   if (allNodes) {
1722     ierr = PetscQuadratureGetData(allNodes, NULL, NULL, &nNodes, &nodes, NULL);CHKERRQ(ierr);
1723   }
1724   ierr = MatGetSize(allMat, &nDofs, NULL);CHKERRQ(ierr);
1725   ierr = PetscDualSpaceGetSection(sp, &section);CHKERRQ(ierr);
1726   ierr = PetscSectionGetStorageSize(section, &spdim);CHKERRQ(ierr);
1727   PetscCheckFalse(spdim != nDofs,PETSC_COMM_SELF, PETSC_ERR_PLIB, "incompatible all matrix size");
1728   ierr = PetscMalloc1(nDofs, &(sp->functional));CHKERRQ(ierr);
1729   ierr = PetscDualSpaceLagrangeGetUseMoments(sp, &useMoments);CHKERRQ(ierr);
1730   if (useMoments) {
1731     Mat              allMat;
1732     PetscInt         momentOrder, i;
1733     PetscBool        tensor;
1734     const PetscReal *weights;
1735     PetscScalar     *array;
1736 
1737     PetscCheckFalse(nDofs != 1,PETSC_COMM_SELF, PETSC_ERR_SUP, "We do not yet support moments beyond P0, nDofs == %D", nDofs);
1738     ierr = PetscDualSpaceLagrangeGetMomentOrder(sp, &momentOrder);CHKERRQ(ierr);
1739     ierr = PetscDualSpaceLagrangeGetTensor(sp, &tensor);CHKERRQ(ierr);
1740     if (!tensor) {ierr = PetscDTStroudConicalQuadrature(dim, Nc, PetscMax(momentOrder + 1,1), -1.0, 1.0, &(sp->functional[0]));CHKERRQ(ierr);}
1741     else         {ierr = PetscDTGaussTensorQuadrature(dim, Nc, PetscMax(momentOrder + 1,1), -1.0, 1.0, &(sp->functional[0]));CHKERRQ(ierr);}
1742     /* Need to replace allNodes and allMat */
1743     ierr = PetscObjectReference((PetscObject) sp->functional[0]);CHKERRQ(ierr);
1744     ierr = PetscQuadratureDestroy(&(sp->allNodes));CHKERRQ(ierr);
1745     sp->allNodes = sp->functional[0];
1746     ierr = PetscQuadratureGetData(sp->allNodes, NULL, NULL, &nNodes, NULL, &weights);CHKERRQ(ierr);
1747     ierr = MatCreateSeqDense(PETSC_COMM_SELF, nDofs, nNodes * Nc, NULL, &allMat);CHKERRQ(ierr);
1748     ierr = MatDenseGetArrayWrite(allMat, &array);CHKERRQ(ierr);
1749     for (i = 0; i < nNodes * Nc; ++i) array[i] = weights[i];
1750     ierr = MatDenseRestoreArrayWrite(allMat, &array);CHKERRQ(ierr);
1751     ierr = MatAssemblyBegin(allMat, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1752     ierr = MatAssemblyEnd(allMat, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1753     ierr = MatDestroy(&(sp->allMat));CHKERRQ(ierr);
1754     sp->allMat = allMat;
1755     PetscFunctionReturn(0);
1756   }
1757   for (f = 0; f < nDofs; f++) {
1758     PetscInt ncols, c;
1759     const PetscInt *cols;
1760     const PetscScalar *vals;
1761     PetscReal *nodesf;
1762     PetscReal *weightsf;
1763     PetscInt nNodesf;
1764     PetscInt countNodes;
1765 
1766     ierr = MatGetRow(allMat, f, &ncols, &cols, &vals);CHKERRQ(ierr);
1767     PetscCheckFalse(ncols % Nk,PETSC_COMM_SELF, PETSC_ERR_PLIB, "all matrix is not laid out as blocks of k-forms");
1768     for (c = 1, nNodesf = 1; c < ncols; c++) {
1769       if ((cols[c] / Nc) != (cols[c-1] / Nc)) nNodesf++;
1770     }
1771     ierr = PetscMalloc1(dim * nNodesf, &nodesf);CHKERRQ(ierr);
1772     ierr = PetscMalloc1(Nc * nNodesf, &weightsf);CHKERRQ(ierr);
1773     for (c = 0, countNodes = 0; c < ncols; c++) {
1774       if (!c || ((cols[c] / Nc) != (cols[c-1] / Nc))) {
1775         PetscInt d;
1776 
1777         for (d = 0; d < Nc; d++) {
1778           weightsf[countNodes * Nc + d] = 0.;
1779         }
1780         for (d = 0; d < dim; d++) {
1781           nodesf[countNodes * dim + d] = nodes[(cols[c] / Nc) * dim + d];
1782         }
1783         countNodes++;
1784       }
1785       weightsf[(countNodes - 1) * Nc + (cols[c] % Nc)] = PetscRealPart(vals[c]);
1786     }
1787     ierr = PetscQuadratureCreate(PETSC_COMM_SELF, &(sp->functional[f]));CHKERRQ(ierr);
1788     ierr = PetscQuadratureSetData(sp->functional[f], dim, Nc, nNodesf, nodesf, weightsf);CHKERRQ(ierr);
1789     ierr = MatRestoreRow(allMat, f, &ncols, &cols, &vals);CHKERRQ(ierr);
1790   }
1791   PetscFunctionReturn(0);
1792 }
1793 
1794 /* take a matrix meant for k-forms and expand it to one for Ncopies */
1795 static PetscErrorCode PetscDualSpaceLagrangeMatrixCreateCopies(Mat A, PetscInt Nk, PetscInt Ncopies, Mat *Abs)
1796 {
1797   PetscInt       m, n, i, j, k;
1798   PetscInt       maxnnz, *nnz, *iwork;
1799   Mat            Ac;
1800   PetscErrorCode ierr;
1801 
1802   PetscFunctionBegin;
1803   ierr = MatGetSize(A, &m, &n);CHKERRQ(ierr);
1804   PetscCheckFalse(n % Nk,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of columns in A %D is not a multiple of Nk %D", n, Nk);
1805   ierr = PetscMalloc1(m * Ncopies, &nnz);CHKERRQ(ierr);
1806   for (i = 0, maxnnz = 0; i < m; i++) {
1807     PetscInt innz;
1808     ierr = MatGetRow(A, i, &innz, NULL, NULL);CHKERRQ(ierr);
1809     PetscCheckFalse(innz % Nk,PETSC_COMM_SELF, PETSC_ERR_PLIB, "A row %D nnzs is not a multiple of Nk %D", innz, Nk);
1810     for (j = 0; j < Ncopies; j++) nnz[i * Ncopies + j] = innz;
1811     maxnnz = PetscMax(maxnnz, innz);
1812   }
1813   ierr = MatCreateSeqAIJ(PETSC_COMM_SELF, m * Ncopies, n * Ncopies, 0, nnz, &Ac);CHKERRQ(ierr);
1814   ierr = MatSetOption(Ac, MAT_IGNORE_ZERO_ENTRIES, PETSC_FALSE);CHKERRQ(ierr);
1815   ierr = PetscFree(nnz);CHKERRQ(ierr);
1816   ierr = PetscMalloc1(maxnnz, &iwork);CHKERRQ(ierr);
1817   for (i = 0; i < m; i++) {
1818     PetscInt innz;
1819     const PetscInt    *cols;
1820     const PetscScalar *vals;
1821 
1822     ierr = MatGetRow(A, i, &innz, &cols, &vals);CHKERRQ(ierr);
1823     for (j = 0; j < innz; j++) iwork[j] = (cols[j] / Nk) * (Nk * Ncopies) + (cols[j] % Nk);
1824     for (j = 0; j < Ncopies; j++) {
1825       PetscInt row = i * Ncopies + j;
1826 
1827       ierr = MatSetValues(Ac, 1, &row, innz, iwork, vals, INSERT_VALUES);CHKERRQ(ierr);
1828       for (k = 0; k < innz; k++) iwork[k] += Nk;
1829     }
1830     ierr = MatRestoreRow(A, i, &innz, &cols, &vals);CHKERRQ(ierr);
1831   }
1832   ierr = PetscFree(iwork);CHKERRQ(ierr);
1833   ierr = MatAssemblyBegin(Ac, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1834   ierr = MatAssemblyEnd(Ac, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1835   *Abs = Ac;
1836   PetscFunctionReturn(0);
1837 }
1838 
1839 /* check if a cell is a tensor product of the segment with a facet,
1840  * specifically checking if f and f2 can be the "endpoints" (like the triangles
1841  * at either end of a wedge) */
1842 static PetscErrorCode DMPlexPointIsTensor_Internal_Given(DM dm, PetscInt p, PetscInt f, PetscInt f2, PetscBool *isTensor)
1843 {
1844   PetscInt        coneSize, c;
1845   const PetscInt *cone;
1846   const PetscInt *fCone;
1847   const PetscInt *f2Cone;
1848   PetscInt        fs[2];
1849   PetscInt        meetSize, nmeet;
1850   const PetscInt *meet;
1851   PetscErrorCode  ierr;
1852 
1853   PetscFunctionBegin;
1854   fs[0] = f;
1855   fs[1] = f2;
1856   ierr = DMPlexGetMeet(dm, 2, fs, &meetSize, &meet);CHKERRQ(ierr);
1857   nmeet = meetSize;
1858   ierr = DMPlexRestoreMeet(dm, 2, fs, &meetSize, &meet);CHKERRQ(ierr);
1859   /* two points that have a non-empty meet cannot be at opposite ends of a cell */
1860   if (nmeet) {
1861     *isTensor = PETSC_FALSE;
1862     PetscFunctionReturn(0);
1863   }
1864   ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
1865   ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
1866   ierr = DMPlexGetCone(dm, f, &fCone);CHKERRQ(ierr);
1867   ierr = DMPlexGetCone(dm, f2, &f2Cone);CHKERRQ(ierr);
1868   for (c = 0; c < coneSize; c++) {
1869     PetscInt e, ef;
1870     PetscInt d = -1, d2 = -1;
1871     PetscInt dcount, d2count;
1872     PetscInt t = cone[c];
1873     PetscInt tConeSize;
1874     PetscBool tIsTensor;
1875     const PetscInt *tCone;
1876 
1877     if (t == f || t == f2) continue;
1878     /* for every other facet in the cone, check that is has
1879      * one ridge in common with each end */
1880     ierr = DMPlexGetConeSize(dm, t, &tConeSize);CHKERRQ(ierr);
1881     ierr = DMPlexGetCone(dm, t, &tCone);CHKERRQ(ierr);
1882 
1883     dcount = 0;
1884     d2count = 0;
1885     for (e = 0; e < tConeSize; e++) {
1886       PetscInt q = tCone[e];
1887       for (ef = 0; ef < coneSize - 2; ef++) {
1888         if (fCone[ef] == q) {
1889           if (dcount) {
1890             *isTensor = PETSC_FALSE;
1891             PetscFunctionReturn(0);
1892           }
1893           d = q;
1894           dcount++;
1895         } else if (f2Cone[ef] == q) {
1896           if (d2count) {
1897             *isTensor = PETSC_FALSE;
1898             PetscFunctionReturn(0);
1899           }
1900           d2 = q;
1901           d2count++;
1902         }
1903       }
1904     }
1905     /* if the whole cell is a tensor with the segment, then this
1906      * facet should be a tensor with the segment */
1907     ierr = DMPlexPointIsTensor_Internal_Given(dm, t, d, d2, &tIsTensor);CHKERRQ(ierr);
1908     if (!tIsTensor) {
1909       *isTensor = PETSC_FALSE;
1910       PetscFunctionReturn(0);
1911     }
1912   }
1913   *isTensor = PETSC_TRUE;
1914   PetscFunctionReturn(0);
1915 }
1916 
1917 /* determine if a cell is a tensor with a segment by looping over pairs of facets to find a pair
1918  * that could be the opposite ends */
1919 static PetscErrorCode DMPlexPointIsTensor_Internal(DM dm, PetscInt p, PetscBool *isTensor, PetscInt *endA, PetscInt *endB)
1920 {
1921   PetscInt        coneSize, c, c2;
1922   const PetscInt *cone;
1923   PetscErrorCode  ierr;
1924 
1925   PetscFunctionBegin;
1926   ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
1927   if (!coneSize) {
1928     if (isTensor) *isTensor = PETSC_FALSE;
1929     if (endA) *endA = -1;
1930     if (endB) *endB = -1;
1931   }
1932   ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
1933   for (c = 0; c < coneSize; c++) {
1934     PetscInt f = cone[c];
1935     PetscInt fConeSize;
1936 
1937     ierr = DMPlexGetConeSize(dm, f, &fConeSize);CHKERRQ(ierr);
1938     if (fConeSize != coneSize - 2) continue;
1939 
1940     for (c2 = c + 1; c2 < coneSize; c2++) {
1941       PetscInt  f2 = cone[c2];
1942       PetscBool isTensorff2;
1943       PetscInt f2ConeSize;
1944 
1945       ierr = DMPlexGetConeSize(dm, f2, &f2ConeSize);CHKERRQ(ierr);
1946       if (f2ConeSize != coneSize - 2) continue;
1947 
1948       ierr = DMPlexPointIsTensor_Internal_Given(dm, p, f, f2, &isTensorff2);CHKERRQ(ierr);
1949       if (isTensorff2) {
1950         if (isTensor) *isTensor = PETSC_TRUE;
1951         if (endA) *endA = f;
1952         if (endB) *endB = f2;
1953         PetscFunctionReturn(0);
1954       }
1955     }
1956   }
1957   if (isTensor) *isTensor = PETSC_FALSE;
1958   if (endA) *endA = -1;
1959   if (endB) *endB = -1;
1960   PetscFunctionReturn(0);
1961 }
1962 
1963 /* determine if a cell is a tensor with a segment by looping over pairs of facets to find a pair
1964  * that could be the opposite ends */
1965 static PetscErrorCode DMPlexPointIsTensor(DM dm, PetscInt p, PetscBool *isTensor, PetscInt *endA, PetscInt *endB)
1966 {
1967   DMPlexInterpolatedFlag interpolated;
1968   PetscErrorCode ierr;
1969 
1970   PetscFunctionBegin;
1971   ierr = DMPlexIsInterpolated(dm, &interpolated);CHKERRQ(ierr);
1972   PetscCheckFalse(interpolated != DMPLEX_INTERPOLATED_FULL,PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Only for interpolated DMPlex's");
1973   ierr = DMPlexPointIsTensor_Internal(dm, p, isTensor, endA, endB);CHKERRQ(ierr);
1974   PetscFunctionReturn(0);
1975 }
1976 
1977 /* Let k = formDegree and k' = -sign(k) * dim + k.  Transform a symmetric frame for k-forms on the biunit simplex into
1978  * a symmetric frame for k'-forms on the biunit simplex.
1979  *
1980  * A frame is "symmetric" if the pullback of every symmetry of the biunit simplex is a permutation of the frame.
1981  *
1982  * forms in the symmetric frame are used as dofs in the untrimmed simplex spaces.  This way, symmetries of the
1983  * reference cell result in permutations of dofs grouped by node.
1984  *
1985  * Use T to transform dof matrices for k'-forms into dof matrices for k-forms as a block diagonal transformation on
1986  * the right.
1987  */
1988 static PetscErrorCode BiunitSimplexSymmetricFormTransformation(PetscInt dim, PetscInt formDegree, PetscReal T[])
1989 {
1990   PetscInt       k = formDegree;
1991   PetscInt       kd = k < 0 ? dim + k : k - dim;
1992   PetscInt       Nk;
1993   PetscReal      *biToEq, *eqToBi, *biToEqStar, *eqToBiStar;
1994   PetscInt       fact;
1995   PetscErrorCode ierr;
1996 
1997   PetscFunctionBegin;
1998   ierr = PetscDTBinomialInt(dim, PetscAbsInt(k), &Nk);CHKERRQ(ierr);
1999   ierr = PetscCalloc4(dim * dim, &biToEq, dim * dim, &eqToBi, Nk * Nk, &biToEqStar, Nk * Nk, &eqToBiStar);CHKERRQ(ierr);
2000   /* fill in biToEq: Jacobian of the transformation from the biunit simplex to the equilateral simplex */
2001   fact = 0;
2002   for (PetscInt i = 0; i < dim; i++) {
2003     biToEq[i * dim + i] = PetscSqrtReal(((PetscReal)i + 2.) / (2.*((PetscReal)i+1.)));
2004     fact += 4*(i+1);
2005     for (PetscInt j = i+1; j < dim; j++) {
2006       biToEq[i * dim + j] = PetscSqrtReal(1./(PetscReal)fact);
2007     }
2008   }
2009   /* fill in eqToBi: Jacobian of the transformation from the equilateral simplex to the biunit simplex */
2010   fact = 0;
2011   for (PetscInt j = 0; j < dim; j++) {
2012     eqToBi[j * dim + j] = PetscSqrtReal(2.*((PetscReal)j+1.)/((PetscReal)j+2));
2013     fact += j+1;
2014     for (PetscInt i = 0; i < j; i++) {
2015       eqToBi[i * dim + j] = -PetscSqrtReal(1./(PetscReal)fact);
2016     }
2017   }
2018   ierr = PetscDTAltVPullbackMatrix(dim, dim, biToEq, kd, biToEqStar);CHKERRQ(ierr);
2019   ierr = PetscDTAltVPullbackMatrix(dim, dim, eqToBi, k, eqToBiStar);CHKERRQ(ierr);
2020   /* product of pullbacks simulates the following steps
2021    *
2022    * 1. start with frame W = [w_1, w_2, ..., w_m] of k forms that is symmetric on the biunit simplex:
2023           if J is the Jacobian of a symmetry of the biunit simplex, then J_k* W = [J_k*w_1, ..., J_k*w_m]
2024           is a permutation of W.
2025           Even though a k' form --- a (dim - k) form represented by its Hodge star --- has the same geometric
2026           content as a k form, W is not a symmetric frame of k' forms on the biunit simplex.  That's because,
2027           for general Jacobian J, J_k* != J_k'*.
2028    * 2. pullback W to the equilateral triangle using the k pullback, W_eq = eqToBi_k* W.  All symmetries of the
2029           equilateral simplex have orthonormal Jacobians.  For an orthonormal Jacobian O, J_k* = J_k'*, so W_eq is
2030           also a symmetric frame for k' forms on the equilateral simplex.
2031      3. pullback W_eq back to the biunit simplex using the k' pulback, V = biToEq_k'* W_eq = biToEq_k'* eqToBi_k* W.
2032           V is a symmetric frame for k' forms on the biunit simplex.
2033    */
2034   for (PetscInt i = 0; i < Nk; i++) {
2035     for (PetscInt j = 0; j < Nk; j++) {
2036       PetscReal val = 0.;
2037       for (PetscInt k = 0; k < Nk; k++) val += biToEqStar[i * Nk + k] * eqToBiStar[k * Nk + j];
2038       T[i * Nk + j] = val;
2039     }
2040   }
2041   ierr = PetscFree4(biToEq, eqToBi, biToEqStar, eqToBiStar);CHKERRQ(ierr);
2042   PetscFunctionReturn(0);
2043 }
2044 
2045 /* permute a quadrature -> dof matrix so that its rows are in revlex order by nodeIdx */
2046 static PetscErrorCode MatPermuteByNodeIdx(Mat A, PetscLagNodeIndices ni, Mat *Aperm)
2047 {
2048   PetscInt       m, n, i, j;
2049   PetscInt       nodeIdxDim = ni->nodeIdxDim;
2050   PetscInt       nodeVecDim = ni->nodeVecDim;
2051   PetscInt       *perm;
2052   IS             permIS;
2053   IS             id;
2054   PetscInt       *nIdxPerm;
2055   PetscReal      *nVecPerm;
2056   PetscErrorCode ierr;
2057 
2058   PetscFunctionBegin;
2059   ierr = PetscLagNodeIndicesGetPermutation(ni, &perm);CHKERRQ(ierr);
2060   ierr = MatGetSize(A, &m, &n);CHKERRQ(ierr);
2061   ierr = PetscMalloc1(nodeIdxDim * m, &nIdxPerm);CHKERRQ(ierr);
2062   ierr = PetscMalloc1(nodeVecDim * m, &nVecPerm);CHKERRQ(ierr);
2063   for (i = 0; i < m; i++) for (j = 0; j < nodeIdxDim; j++) nIdxPerm[i * nodeIdxDim + j] = ni->nodeIdx[perm[i] * nodeIdxDim + j];
2064   for (i = 0; i < m; i++) for (j = 0; j < nodeVecDim; j++) nVecPerm[i * nodeVecDim + j] = ni->nodeVec[perm[i] * nodeVecDim + j];
2065   ierr = ISCreateGeneral(PETSC_COMM_SELF, m, perm, PETSC_USE_POINTER, &permIS);CHKERRQ(ierr);
2066   ierr = ISSetPermutation(permIS);CHKERRQ(ierr);
2067   ierr = ISCreateStride(PETSC_COMM_SELF, n, 0, 1, &id);CHKERRQ(ierr);
2068   ierr = ISSetPermutation(id);CHKERRQ(ierr);
2069   ierr = MatPermute(A, permIS, id, Aperm);CHKERRQ(ierr);
2070   ierr = ISDestroy(&permIS);CHKERRQ(ierr);
2071   ierr = ISDestroy(&id);CHKERRQ(ierr);
2072   for (i = 0; i < m; i++) perm[i] = i;
2073   ierr = PetscFree(ni->nodeIdx);CHKERRQ(ierr);
2074   ierr = PetscFree(ni->nodeVec);CHKERRQ(ierr);
2075   ni->nodeIdx = nIdxPerm;
2076   ni->nodeVec = nVecPerm;
2077   PetscFunctionReturn(0);
2078 }
2079 
2080 static PetscErrorCode PetscDualSpaceSetUp_Lagrange(PetscDualSpace sp)
2081 {
2082   PetscDualSpace_Lag *lag   = (PetscDualSpace_Lag *) sp->data;
2083   DM                  dm    = sp->dm;
2084   DM                  dmint = NULL;
2085   PetscInt            order;
2086   PetscInt            Nc    = sp->Nc;
2087   MPI_Comm            comm;
2088   PetscBool           continuous;
2089   PetscSection        section;
2090   PetscInt            depth, dim, pStart, pEnd, cStart, cEnd, p, *pStratStart, *pStratEnd, d;
2091   PetscInt            formDegree, Nk, Ncopies;
2092   PetscInt            tensorf = -1, tensorf2 = -1;
2093   PetscBool           tensorCell, tensorSpace;
2094   PetscBool           uniform, trimmed;
2095   Petsc1DNodeFamily   nodeFamily;
2096   PetscInt            numNodeSkip;
2097   DMPlexInterpolatedFlag interpolated;
2098   PetscBool           isbdm;
2099   PetscErrorCode      ierr;
2100 
2101   PetscFunctionBegin;
2102   /* step 1: sanitize input */
2103   ierr = PetscObjectGetComm((PetscObject) sp, &comm);CHKERRQ(ierr);
2104   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
2105   ierr = PetscObjectTypeCompare((PetscObject)sp, PETSCDUALSPACEBDM, &isbdm);CHKERRQ(ierr);
2106   if (isbdm) {
2107     sp->k = -(dim-1); /* form degree of H-div */
2108     ierr = PetscObjectChangeTypeName((PetscObject)sp, PETSCDUALSPACELAGRANGE);CHKERRQ(ierr);
2109   }
2110   ierr = PetscDualSpaceGetFormDegree(sp, &formDegree);CHKERRQ(ierr);
2111   PetscCheckFalse(PetscAbsInt(formDegree) > dim,comm, PETSC_ERR_ARG_OUTOFRANGE, "Form degree must be bounded by dimension");
2112   ierr = PetscDTBinomialInt(dim,PetscAbsInt(formDegree),&Nk);CHKERRQ(ierr);
2113   if (sp->Nc <= 0 && lag->numCopies > 0) sp->Nc = Nk * lag->numCopies;
2114   Nc = sp->Nc;
2115   PetscCheckFalse(Nc % Nk,comm, PETSC_ERR_ARG_INCOMP, "Number of components is not a multiple of form degree size");
2116   if (lag->numCopies <= 0) lag->numCopies = Nc / Nk;
2117   Ncopies = lag->numCopies;
2118   PetscCheckFalse(Nc / Nk != Ncopies,comm, PETSC_ERR_ARG_INCOMP, "Number of copies * (dim choose k) != Nc");
2119   if (!dim) sp->order = 0;
2120   order = sp->order;
2121   uniform = sp->uniform;
2122   PetscCheckFalse(!uniform,PETSC_COMM_SELF, PETSC_ERR_SUP, "Variable order not supported yet");
2123   if (lag->trimmed && !formDegree) lag->trimmed = PETSC_FALSE; /* trimmed spaces are the same as full spaces for 0-forms */
2124   if (lag->nodeType == PETSCDTNODES_DEFAULT) {
2125     lag->nodeType = PETSCDTNODES_GAUSSJACOBI;
2126     lag->nodeExponent = 0.;
2127     /* trimmed spaces don't include corner vertices, so don't use end nodes by default */
2128     lag->endNodes = lag->trimmed ? PETSC_FALSE : PETSC_TRUE;
2129   }
2130   /* If a trimmed space and the user did choose nodes with endpoints, skip them by default */
2131   if (lag->numNodeSkip < 0) lag->numNodeSkip = (lag->trimmed && lag->endNodes) ? 1 : 0;
2132   numNodeSkip = lag->numNodeSkip;
2133   PetscCheckFalse(lag->trimmed && !order,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot have zeroth order trimmed elements");
2134   if (lag->trimmed && PetscAbsInt(formDegree) == dim) { /* convert trimmed n-forms to untrimmed of one polynomial order less */
2135     sp->order--;
2136     order--;
2137     lag->trimmed = PETSC_FALSE;
2138   }
2139   trimmed = lag->trimmed;
2140   if (!order || PetscAbsInt(formDegree) == dim) lag->continuous = PETSC_FALSE;
2141   continuous = lag->continuous;
2142   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
2143   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
2144   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
2145   PetscCheckFalse(pStart != 0 || cStart != 0,PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Expect DM with chart starting at zero and cells first");
2146   PetscCheckFalse(cEnd != 1,PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Use PETSCDUALSPACEREFINED for multi-cell reference meshes");
2147   ierr = DMPlexIsInterpolated(dm, &interpolated);CHKERRQ(ierr);
2148   if (interpolated != DMPLEX_INTERPOLATED_FULL) {
2149     ierr = DMPlexInterpolate(dm, &dmint);CHKERRQ(ierr);
2150   } else {
2151     ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);
2152     dmint = dm;
2153   }
2154   tensorCell = PETSC_FALSE;
2155   if (dim > 1) {
2156     ierr = DMPlexPointIsTensor(dmint, 0, &tensorCell, &tensorf, &tensorf2);CHKERRQ(ierr);
2157   }
2158   lag->tensorCell = tensorCell;
2159   if (dim < 2 || !lag->tensorCell) lag->tensorSpace = PETSC_FALSE;
2160   tensorSpace = lag->tensorSpace;
2161   if (!lag->nodeFamily) {
2162     ierr = Petsc1DNodeFamilyCreate(lag->nodeType, lag->nodeExponent, lag->endNodes, &lag->nodeFamily);CHKERRQ(ierr);
2163   }
2164   nodeFamily = lag->nodeFamily;
2165   PetscCheckFalse(interpolated != DMPLEX_INTERPOLATED_FULL && continuous && (PetscAbsInt(formDegree) > 0 || order > 1),PETSC_COMM_SELF,PETSC_ERR_PLIB,"Reference element won't support all boundary nodes");
2166 
2167   /* step 2: construct the boundary spaces */
2168   ierr = PetscMalloc2(depth+1,&pStratStart,depth+1,&pStratEnd);CHKERRQ(ierr);
2169   ierr = PetscCalloc1(pEnd,&(sp->pointSpaces));CHKERRQ(ierr);
2170   for (d = 0; d <= depth; ++d) {ierr = DMPlexGetDepthStratum(dm, d, &pStratStart[d], &pStratEnd[d]);CHKERRQ(ierr);}
2171   ierr = PetscDualSpaceSectionCreate_Internal(sp, &section);CHKERRQ(ierr);
2172   sp->pointSection = section;
2173   if (continuous && !(lag->interiorOnly)) {
2174     PetscInt h;
2175 
2176     for (p = pStratStart[depth - 1]; p < pStratEnd[depth - 1]; p++) { /* calculate the facet dual spaces */
2177       PetscReal v0[3];
2178       DMPolytopeType ptype;
2179       PetscReal J[9], detJ;
2180       PetscInt  q;
2181 
2182       ierr = DMPlexComputeCellGeometryAffineFEM(dm, p, v0, J, NULL, &detJ);CHKERRQ(ierr);
2183       ierr = DMPlexGetCellType(dm, p, &ptype);CHKERRQ(ierr);
2184 
2185       /* compare to previous facets: if computed, reference that dualspace */
2186       for (q = pStratStart[depth - 1]; q < p; q++) {
2187         DMPolytopeType qtype;
2188 
2189         ierr = DMPlexGetCellType(dm, q, &qtype);CHKERRQ(ierr);
2190         if (qtype == ptype) break;
2191       }
2192       if (q < p) { /* this facet has the same dual space as that one */
2193         ierr = PetscObjectReference((PetscObject)sp->pointSpaces[q]);CHKERRQ(ierr);
2194         sp->pointSpaces[p] = sp->pointSpaces[q];
2195         continue;
2196       }
2197       /* if not, recursively compute this dual space */
2198       ierr = PetscDualSpaceCreateFacetSubspace_Lagrange(sp,NULL,p,formDegree,Ncopies,PETSC_FALSE,&sp->pointSpaces[p]);CHKERRQ(ierr);
2199     }
2200     for (h = 2; h <= depth; h++) { /* get the higher subspaces from the facet subspaces */
2201       PetscInt hd = depth - h;
2202       PetscInt hdim = dim - h;
2203 
2204       if (hdim < PetscAbsInt(formDegree)) break;
2205       for (p = pStratStart[hd]; p < pStratEnd[hd]; p++) {
2206         PetscInt suppSize, s;
2207         const PetscInt *supp;
2208 
2209         ierr = DMPlexGetSupportSize(dm, p, &suppSize);CHKERRQ(ierr);
2210         ierr = DMPlexGetSupport(dm, p, &supp);CHKERRQ(ierr);
2211         for (s = 0; s < suppSize; s++) {
2212           DM             qdm;
2213           PetscDualSpace qsp, psp;
2214           PetscInt c, coneSize, q;
2215           const PetscInt *cone;
2216           const PetscInt *refCone;
2217 
2218           q = supp[0];
2219           qsp = sp->pointSpaces[q];
2220           ierr = DMPlexGetConeSize(dm, q, &coneSize);CHKERRQ(ierr);
2221           ierr = DMPlexGetCone(dm, q, &cone);CHKERRQ(ierr);
2222           for (c = 0; c < coneSize; c++) if (cone[c] == p) break;
2223           PetscCheckFalse(c == coneSize,PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "cone/support mismatch");
2224           ierr = PetscDualSpaceGetDM(qsp, &qdm);CHKERRQ(ierr);
2225           ierr = DMPlexGetCone(qdm, 0, &refCone);CHKERRQ(ierr);
2226           /* get the equivalent dual space from the support dual space */
2227           ierr = PetscDualSpaceGetPointSubspace(qsp, refCone[c], &psp);CHKERRQ(ierr);
2228           if (!s) {
2229             ierr = PetscObjectReference((PetscObject)psp);CHKERRQ(ierr);
2230             sp->pointSpaces[p] = psp;
2231           }
2232         }
2233       }
2234     }
2235     for (p = 1; p < pEnd; p++) {
2236       PetscInt pspdim;
2237       if (!sp->pointSpaces[p]) continue;
2238       ierr = PetscDualSpaceGetInteriorDimension(sp->pointSpaces[p], &pspdim);CHKERRQ(ierr);
2239       ierr = PetscSectionSetDof(section, p, pspdim);CHKERRQ(ierr);
2240     }
2241   }
2242 
2243   if (Ncopies > 1) {
2244     Mat intMatScalar, allMatScalar;
2245     PetscDualSpace scalarsp;
2246     PetscDualSpace_Lag *scalarlag;
2247 
2248     ierr = PetscDualSpaceDuplicate(sp, &scalarsp);CHKERRQ(ierr);
2249     /* Setting the number of components to Nk is a space with 1 copy of each k-form */
2250     ierr = PetscDualSpaceSetNumComponents(scalarsp, Nk);CHKERRQ(ierr);
2251     ierr = PetscDualSpaceSetUp(scalarsp);CHKERRQ(ierr);
2252     ierr = PetscDualSpaceGetInteriorData(scalarsp, &(sp->intNodes), &intMatScalar);CHKERRQ(ierr);
2253     ierr = PetscObjectReference((PetscObject)(sp->intNodes));CHKERRQ(ierr);
2254     if (intMatScalar) {ierr = PetscDualSpaceLagrangeMatrixCreateCopies(intMatScalar, Nk, Ncopies, &(sp->intMat));CHKERRQ(ierr);}
2255     ierr = PetscDualSpaceGetAllData(scalarsp, &(sp->allNodes), &allMatScalar);CHKERRQ(ierr);
2256     ierr = PetscObjectReference((PetscObject)(sp->allNodes));CHKERRQ(ierr);
2257     ierr = PetscDualSpaceLagrangeMatrixCreateCopies(allMatScalar, Nk, Ncopies, &(sp->allMat));CHKERRQ(ierr);
2258     sp->spdim = scalarsp->spdim * Ncopies;
2259     sp->spintdim = scalarsp->spintdim * Ncopies;
2260     scalarlag = (PetscDualSpace_Lag *) scalarsp->data;
2261     ierr = PetscLagNodeIndicesReference(scalarlag->vertIndices);CHKERRQ(ierr);
2262     lag->vertIndices = scalarlag->vertIndices;
2263     ierr = PetscLagNodeIndicesReference(scalarlag->intNodeIndices);CHKERRQ(ierr);
2264     lag->intNodeIndices = scalarlag->intNodeIndices;
2265     ierr = PetscLagNodeIndicesReference(scalarlag->allNodeIndices);CHKERRQ(ierr);
2266     lag->allNodeIndices = scalarlag->allNodeIndices;
2267     ierr = PetscDualSpaceDestroy(&scalarsp);CHKERRQ(ierr);
2268     ierr = PetscSectionSetDof(section, 0, sp->spintdim);CHKERRQ(ierr);
2269     ierr = PetscDualSpaceSectionSetUp_Internal(sp, section);CHKERRQ(ierr);
2270     ierr = PetscDualSpaceComputeFunctionalsFromAllData(sp);CHKERRQ(ierr);
2271     ierr = PetscFree2(pStratStart, pStratEnd);CHKERRQ(ierr);
2272     ierr = DMDestroy(&dmint);CHKERRQ(ierr);
2273     PetscFunctionReturn(0);
2274   }
2275 
2276   if (trimmed && !continuous) {
2277     /* the dofs of a trimmed space don't have a nice tensor/lattice structure:
2278      * just construct the continuous dual space and copy all of the data over,
2279      * allocating it all to the cell instead of splitting it up between the boundaries */
2280     PetscDualSpace  spcont;
2281     PetscInt        spdim, f;
2282     PetscQuadrature allNodes;
2283     PetscDualSpace_Lag *lagc;
2284     Mat             allMat;
2285 
2286     ierr = PetscDualSpaceDuplicate(sp, &spcont);CHKERRQ(ierr);
2287     ierr = PetscDualSpaceLagrangeSetContinuity(spcont, PETSC_TRUE);CHKERRQ(ierr);
2288     ierr = PetscDualSpaceSetUp(spcont);CHKERRQ(ierr);
2289     ierr = PetscDualSpaceGetDimension(spcont, &spdim);CHKERRQ(ierr);
2290     sp->spdim = sp->spintdim = spdim;
2291     ierr = PetscSectionSetDof(section, 0, spdim);CHKERRQ(ierr);
2292     ierr = PetscDualSpaceSectionSetUp_Internal(sp, section);CHKERRQ(ierr);
2293     ierr = PetscMalloc1(spdim, &(sp->functional));CHKERRQ(ierr);
2294     for (f = 0; f < spdim; f++) {
2295       PetscQuadrature fn;
2296 
2297       ierr = PetscDualSpaceGetFunctional(spcont, f, &fn);CHKERRQ(ierr);
2298       ierr = PetscObjectReference((PetscObject)fn);CHKERRQ(ierr);
2299       sp->functional[f] = fn;
2300     }
2301     ierr = PetscDualSpaceGetAllData(spcont, &allNodes, &allMat);CHKERRQ(ierr);
2302     ierr = PetscObjectReference((PetscObject) allNodes);CHKERRQ(ierr);
2303     ierr = PetscObjectReference((PetscObject) allNodes);CHKERRQ(ierr);
2304     sp->allNodes = sp->intNodes = allNodes;
2305     ierr = PetscObjectReference((PetscObject) allMat);CHKERRQ(ierr);
2306     ierr = PetscObjectReference((PetscObject) allMat);CHKERRQ(ierr);
2307     sp->allMat = sp->intMat = allMat;
2308     lagc = (PetscDualSpace_Lag *) spcont->data;
2309     ierr = PetscLagNodeIndicesReference(lagc->vertIndices);CHKERRQ(ierr);
2310     lag->vertIndices = lagc->vertIndices;
2311     ierr = PetscLagNodeIndicesReference(lagc->allNodeIndices);CHKERRQ(ierr);
2312     ierr = PetscLagNodeIndicesReference(lagc->allNodeIndices);CHKERRQ(ierr);
2313     lag->intNodeIndices = lagc->allNodeIndices;
2314     lag->allNodeIndices = lagc->allNodeIndices;
2315     ierr = PetscDualSpaceDestroy(&spcont);CHKERRQ(ierr);
2316     ierr = PetscFree2(pStratStart, pStratEnd);CHKERRQ(ierr);
2317     ierr = DMDestroy(&dmint);CHKERRQ(ierr);
2318     PetscFunctionReturn(0);
2319   }
2320 
2321   /* step 3: construct intNodes, and intMat, and combine it with boundray data to make allNodes and allMat */
2322   if (!tensorSpace) {
2323     if (!tensorCell) {ierr = PetscLagNodeIndicesCreateSimplexVertices(dm, &(lag->vertIndices));CHKERRQ(ierr);}
2324 
2325     if (trimmed) {
2326       /* there is one dof in the interior of the a trimmed element for each full polynomial of with degree at most
2327        * order + k - dim - 1 */
2328       if (order + PetscAbsInt(formDegree) > dim) {
2329         PetscInt sum = order + PetscAbsInt(formDegree) - dim - 1;
2330         PetscInt nDofs;
2331 
2332         ierr = PetscDualSpaceLagrangeCreateSimplexNodeMat(nodeFamily, dim, sum, Nk, numNodeSkip, &sp->intNodes, &sp->intMat, &(lag->intNodeIndices));CHKERRQ(ierr);
2333         ierr = MatGetSize(sp->intMat, &nDofs, NULL);CHKERRQ(ierr);
2334         ierr = PetscSectionSetDof(section, 0, nDofs);CHKERRQ(ierr);
2335       }
2336       ierr = PetscDualSpaceSectionSetUp_Internal(sp, section);CHKERRQ(ierr);
2337       ierr = PetscDualSpaceCreateAllDataFromInteriorData(sp);CHKERRQ(ierr);
2338       ierr = PetscDualSpaceLagrangeCreateAllNodeIdx(sp);CHKERRQ(ierr);
2339     } else {
2340       if (!continuous) {
2341         /* if discontinuous just construct one node for each set of dofs (a set of dofs is a basis for the k-form
2342          * space) */
2343         PetscInt sum = order;
2344         PetscInt nDofs;
2345 
2346         ierr = PetscDualSpaceLagrangeCreateSimplexNodeMat(nodeFamily, dim, sum, Nk, numNodeSkip, &sp->intNodes, &sp->intMat, &(lag->intNodeIndices));CHKERRQ(ierr);
2347         ierr = MatGetSize(sp->intMat, &nDofs, NULL);CHKERRQ(ierr);
2348         ierr = PetscSectionSetDof(section, 0, nDofs);CHKERRQ(ierr);
2349         ierr = PetscDualSpaceSectionSetUp_Internal(sp, section);CHKERRQ(ierr);
2350         ierr = PetscObjectReference((PetscObject)(sp->intNodes));CHKERRQ(ierr);
2351         sp->allNodes = sp->intNodes;
2352         ierr = PetscObjectReference((PetscObject)(sp->intMat));CHKERRQ(ierr);
2353         sp->allMat = sp->intMat;
2354         ierr = PetscLagNodeIndicesReference(lag->intNodeIndices);CHKERRQ(ierr);
2355         lag->allNodeIndices = lag->intNodeIndices;
2356       } else {
2357         /* there is one dof in the interior of the a full element for each trimmed polynomial of with degree at most
2358          * order + k - dim, but with complementary form degree */
2359         if (order + PetscAbsInt(formDegree) > dim) {
2360           PetscDualSpace trimmedsp;
2361           PetscDualSpace_Lag *trimmedlag;
2362           PetscQuadrature intNodes;
2363           PetscInt trFormDegree = formDegree >= 0 ? formDegree - dim : dim - PetscAbsInt(formDegree);
2364           PetscInt nDofs;
2365           Mat intMat;
2366 
2367           ierr = PetscDualSpaceDuplicate(sp, &trimmedsp);CHKERRQ(ierr);
2368           ierr = PetscDualSpaceLagrangeSetTrimmed(trimmedsp, PETSC_TRUE);CHKERRQ(ierr);
2369           ierr = PetscDualSpaceSetOrder(trimmedsp, order + PetscAbsInt(formDegree) - dim);CHKERRQ(ierr);
2370           ierr = PetscDualSpaceSetFormDegree(trimmedsp, trFormDegree);CHKERRQ(ierr);
2371           trimmedlag = (PetscDualSpace_Lag *) trimmedsp->data;
2372           trimmedlag->numNodeSkip = numNodeSkip + 1;
2373           ierr = PetscDualSpaceSetUp(trimmedsp);CHKERRQ(ierr);
2374           ierr = PetscDualSpaceGetAllData(trimmedsp, &intNodes, &intMat);CHKERRQ(ierr);
2375           ierr = PetscObjectReference((PetscObject)intNodes);CHKERRQ(ierr);
2376           sp->intNodes = intNodes;
2377           ierr = PetscLagNodeIndicesReference(trimmedlag->allNodeIndices);CHKERRQ(ierr);
2378           lag->intNodeIndices = trimmedlag->allNodeIndices;
2379           ierr = PetscObjectReference((PetscObject)intMat);CHKERRQ(ierr);
2380           if (PetscAbsInt(formDegree) > 0 && PetscAbsInt(formDegree) < dim) {
2381             PetscReal *T;
2382             PetscScalar *work;
2383             PetscInt nCols, nRows;
2384             Mat intMatT;
2385 
2386             ierr = MatDuplicate(intMat, MAT_COPY_VALUES, &intMatT);CHKERRQ(ierr);
2387             ierr = MatGetSize(intMat, &nRows, &nCols);CHKERRQ(ierr);
2388             ierr = PetscMalloc2(Nk * Nk, &T, nCols, &work);CHKERRQ(ierr);
2389             ierr = BiunitSimplexSymmetricFormTransformation(dim, formDegree, T);CHKERRQ(ierr);
2390             for (PetscInt row = 0; row < nRows; row++) {
2391               PetscInt nrCols;
2392               const PetscInt *rCols;
2393               const PetscScalar *rVals;
2394 
2395               ierr = MatGetRow(intMat, row, &nrCols, &rCols, &rVals);CHKERRQ(ierr);
2396               PetscCheckFalse(nrCols % Nk,PETSC_COMM_SELF, PETSC_ERR_PLIB, "nonzeros in intMat matrix are not in k-form size blocks");
2397               for (PetscInt b = 0; b < nrCols; b += Nk) {
2398                 const PetscScalar *v = &rVals[b];
2399                 PetscScalar *w = &work[b];
2400                 for (PetscInt j = 0; j < Nk; j++) {
2401                   w[j] = 0.;
2402                   for (PetscInt i = 0; i < Nk; i++) {
2403                     w[j] += v[i] * T[i * Nk + j];
2404                   }
2405                 }
2406               }
2407               ierr = MatSetValuesBlocked(intMatT, 1, &row, nrCols, rCols, work, INSERT_VALUES);CHKERRQ(ierr);
2408               ierr = MatRestoreRow(intMat, row, &nrCols, &rCols, &rVals);CHKERRQ(ierr);
2409             }
2410             ierr = MatAssemblyBegin(intMatT, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2411             ierr = MatAssemblyEnd(intMatT, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2412             ierr = MatDestroy(&intMat);CHKERRQ(ierr);
2413             intMat = intMatT;
2414             ierr = PetscLagNodeIndicesDestroy(&(lag->intNodeIndices));CHKERRQ(ierr);
2415             ierr = PetscLagNodeIndicesDuplicate(trimmedlag->allNodeIndices, &(lag->intNodeIndices));CHKERRQ(ierr);
2416             {
2417               PetscInt nNodes = lag->intNodeIndices->nNodes;
2418               PetscReal *newNodeVec = lag->intNodeIndices->nodeVec;
2419               const PetscReal *oldNodeVec = trimmedlag->allNodeIndices->nodeVec;
2420 
2421               for (PetscInt n = 0; n < nNodes; n++) {
2422                 PetscReal *w = &newNodeVec[n * Nk];
2423                 const PetscReal *v = &oldNodeVec[n * Nk];
2424 
2425                 for (PetscInt j = 0; j < Nk; j++) {
2426                   w[j] = 0.;
2427                   for (PetscInt i = 0; i < Nk; i++) {
2428                     w[j] += v[i] * T[i * Nk + j];
2429                   }
2430                 }
2431               }
2432             }
2433             ierr = PetscFree2(T, work);CHKERRQ(ierr);
2434           }
2435           sp->intMat = intMat;
2436           ierr = MatGetSize(sp->intMat, &nDofs, NULL);CHKERRQ(ierr);
2437           ierr = PetscDualSpaceDestroy(&trimmedsp);CHKERRQ(ierr);
2438           ierr = PetscSectionSetDof(section, 0, nDofs);CHKERRQ(ierr);
2439         }
2440         ierr = PetscDualSpaceSectionSetUp_Internal(sp, section);CHKERRQ(ierr);
2441         ierr = PetscDualSpaceCreateAllDataFromInteriorData(sp);CHKERRQ(ierr);
2442         ierr = PetscDualSpaceLagrangeCreateAllNodeIdx(sp);CHKERRQ(ierr);
2443       }
2444     }
2445   } else {
2446     PetscQuadrature intNodesTrace = NULL;
2447     PetscQuadrature intNodesFiber = NULL;
2448     PetscQuadrature intNodes = NULL;
2449     PetscLagNodeIndices intNodeIndices = NULL;
2450     Mat             intMat = NULL;
2451 
2452     if (PetscAbsInt(formDegree) < dim) { /* get the trace k-forms on the first facet, and the 0-forms on the edge,
2453                                             and wedge them together to create some of the k-form dofs */
2454       PetscDualSpace  trace, fiber;
2455       PetscDualSpace_Lag *tracel, *fiberl;
2456       Mat             intMatTrace, intMatFiber;
2457 
2458       if (sp->pointSpaces[tensorf]) {
2459         ierr = PetscObjectReference((PetscObject)(sp->pointSpaces[tensorf]));CHKERRQ(ierr);
2460         trace = sp->pointSpaces[tensorf];
2461       } else {
2462         ierr = PetscDualSpaceCreateFacetSubspace_Lagrange(sp,NULL,tensorf,formDegree,Ncopies,PETSC_TRUE,&trace);CHKERRQ(ierr);
2463       }
2464       ierr = PetscDualSpaceCreateEdgeSubspace_Lagrange(sp,order,0,1,PETSC_TRUE,&fiber);CHKERRQ(ierr);
2465       tracel = (PetscDualSpace_Lag *) trace->data;
2466       fiberl = (PetscDualSpace_Lag *) fiber->data;
2467       ierr = PetscLagNodeIndicesCreateTensorVertices(dm, tracel->vertIndices, &(lag->vertIndices));CHKERRQ(ierr);
2468       ierr = PetscDualSpaceGetInteriorData(trace, &intNodesTrace, &intMatTrace);CHKERRQ(ierr);
2469       ierr = PetscDualSpaceGetInteriorData(fiber, &intNodesFiber, &intMatFiber);CHKERRQ(ierr);
2470       if (intNodesTrace && intNodesFiber) {
2471         ierr = PetscQuadratureCreateTensor(intNodesTrace, intNodesFiber, &intNodes);CHKERRQ(ierr);
2472         ierr = MatTensorAltV(intMatTrace, intMatFiber, dim-1, formDegree, 1, 0, &intMat);CHKERRQ(ierr);
2473         ierr = PetscLagNodeIndicesTensor(tracel->intNodeIndices, dim - 1, formDegree, fiberl->intNodeIndices, 1, 0, &intNodeIndices);CHKERRQ(ierr);
2474       }
2475       ierr = PetscObjectReference((PetscObject) intNodesTrace);CHKERRQ(ierr);
2476       ierr = PetscObjectReference((PetscObject) intNodesFiber);CHKERRQ(ierr);
2477       ierr = PetscDualSpaceDestroy(&fiber);CHKERRQ(ierr);
2478       ierr = PetscDualSpaceDestroy(&trace);CHKERRQ(ierr);
2479     }
2480     if (PetscAbsInt(formDegree) > 0) { /* get the trace (k-1)-forms on the first facet, and the 1-forms on the edge,
2481                                           and wedge them together to create the remaining k-form dofs */
2482       PetscDualSpace  trace, fiber;
2483       PetscDualSpace_Lag *tracel, *fiberl;
2484       PetscQuadrature intNodesTrace2, intNodesFiber2, intNodes2;
2485       PetscLagNodeIndices intNodeIndices2;
2486       Mat             intMatTrace, intMatFiber, intMat2;
2487       PetscInt        traceDegree = formDegree > 0 ? formDegree - 1 : formDegree + 1;
2488       PetscInt        fiberDegree = formDegree > 0 ? 1 : -1;
2489 
2490       ierr = PetscDualSpaceCreateFacetSubspace_Lagrange(sp,NULL,tensorf,traceDegree,Ncopies,PETSC_TRUE,&trace);CHKERRQ(ierr);
2491       ierr = PetscDualSpaceCreateEdgeSubspace_Lagrange(sp,order,fiberDegree,1,PETSC_TRUE,&fiber);CHKERRQ(ierr);
2492       tracel = (PetscDualSpace_Lag *) trace->data;
2493       fiberl = (PetscDualSpace_Lag *) fiber->data;
2494       if (!lag->vertIndices) {
2495         ierr = PetscLagNodeIndicesCreateTensorVertices(dm, tracel->vertIndices, &(lag->vertIndices));CHKERRQ(ierr);
2496       }
2497       ierr = PetscDualSpaceGetInteriorData(trace, &intNodesTrace2, &intMatTrace);CHKERRQ(ierr);
2498       ierr = PetscDualSpaceGetInteriorData(fiber, &intNodesFiber2, &intMatFiber);CHKERRQ(ierr);
2499       if (intNodesTrace2 && intNodesFiber2) {
2500         ierr = PetscQuadratureCreateTensor(intNodesTrace2, intNodesFiber2, &intNodes2);CHKERRQ(ierr);
2501         ierr = MatTensorAltV(intMatTrace, intMatFiber, dim-1, traceDegree, 1, fiberDegree, &intMat2);CHKERRQ(ierr);
2502         ierr = PetscLagNodeIndicesTensor(tracel->intNodeIndices, dim - 1, traceDegree, fiberl->intNodeIndices, 1, fiberDegree, &intNodeIndices2);CHKERRQ(ierr);
2503         if (!intMat) {
2504           intMat = intMat2;
2505           intNodes = intNodes2;
2506           intNodeIndices = intNodeIndices2;
2507         } else {
2508           /* merge the matrices, quadrature points, and nodes */
2509           PetscInt         nM;
2510           PetscInt         nDof, nDof2;
2511           PetscInt        *toMerged = NULL, *toMerged2 = NULL;
2512           PetscQuadrature  merged = NULL;
2513           PetscLagNodeIndices intNodeIndicesMerged = NULL;
2514           Mat              matMerged = NULL;
2515 
2516           ierr = MatGetSize(intMat, &nDof, NULL);CHKERRQ(ierr);
2517           ierr = MatGetSize(intMat2, &nDof2, NULL);CHKERRQ(ierr);
2518           ierr = PetscQuadraturePointsMerge(intNodes, intNodes2, &merged, &toMerged, &toMerged2);CHKERRQ(ierr);
2519           ierr = PetscQuadratureGetData(merged, NULL, NULL, &nM, NULL, NULL);CHKERRQ(ierr);
2520           ierr = MatricesMerge(intMat, intMat2, dim, formDegree, nM, toMerged, toMerged2, &matMerged);CHKERRQ(ierr);
2521           ierr = PetscLagNodeIndicesMerge(intNodeIndices, intNodeIndices2, &intNodeIndicesMerged);CHKERRQ(ierr);
2522           ierr = PetscFree(toMerged);CHKERRQ(ierr);
2523           ierr = PetscFree(toMerged2);CHKERRQ(ierr);
2524           ierr = MatDestroy(&intMat);CHKERRQ(ierr);
2525           ierr = MatDestroy(&intMat2);CHKERRQ(ierr);
2526           ierr = PetscQuadratureDestroy(&intNodes);CHKERRQ(ierr);
2527           ierr = PetscQuadratureDestroy(&intNodes2);CHKERRQ(ierr);
2528           ierr = PetscLagNodeIndicesDestroy(&intNodeIndices);CHKERRQ(ierr);
2529           ierr = PetscLagNodeIndicesDestroy(&intNodeIndices2);CHKERRQ(ierr);
2530           intNodes = merged;
2531           intMat = matMerged;
2532           intNodeIndices = intNodeIndicesMerged;
2533           if (!trimmed) {
2534             /* I think users expect that, when a node has a full basis for the k-forms,
2535              * they should be consecutive dofs.  That isn't the case for trimmed spaces,
2536              * but is for some of the nodes in untrimmed spaces, so in that case we
2537              * sort them to group them by node */
2538             Mat intMatPerm;
2539 
2540             ierr = MatPermuteByNodeIdx(intMat, intNodeIndices, &intMatPerm);CHKERRQ(ierr);
2541             ierr = MatDestroy(&intMat);CHKERRQ(ierr);
2542             intMat = intMatPerm;
2543           }
2544         }
2545       }
2546       ierr = PetscDualSpaceDestroy(&fiber);CHKERRQ(ierr);
2547       ierr = PetscDualSpaceDestroy(&trace);CHKERRQ(ierr);
2548     }
2549     ierr = PetscQuadratureDestroy(&intNodesTrace);CHKERRQ(ierr);
2550     ierr = PetscQuadratureDestroy(&intNodesFiber);CHKERRQ(ierr);
2551     sp->intNodes = intNodes;
2552     sp->intMat = intMat;
2553     lag->intNodeIndices = intNodeIndices;
2554     {
2555       PetscInt nDofs = 0;
2556 
2557       if (intMat) {
2558         ierr = MatGetSize(intMat, &nDofs, NULL);CHKERRQ(ierr);
2559       }
2560       ierr = PetscSectionSetDof(section, 0, nDofs);CHKERRQ(ierr);
2561     }
2562     ierr = PetscDualSpaceSectionSetUp_Internal(sp, section);CHKERRQ(ierr);
2563     if (continuous) {
2564       ierr = PetscDualSpaceCreateAllDataFromInteriorData(sp);CHKERRQ(ierr);
2565       ierr = PetscDualSpaceLagrangeCreateAllNodeIdx(sp);CHKERRQ(ierr);
2566     } else {
2567       ierr = PetscObjectReference((PetscObject) intNodes);CHKERRQ(ierr);
2568       sp->allNodes = intNodes;
2569       ierr = PetscObjectReference((PetscObject) intMat);CHKERRQ(ierr);
2570       sp->allMat = intMat;
2571       ierr = PetscLagNodeIndicesReference(intNodeIndices);CHKERRQ(ierr);
2572       lag->allNodeIndices = intNodeIndices;
2573     }
2574   }
2575   ierr = PetscSectionGetStorageSize(section, &sp->spdim);CHKERRQ(ierr);
2576   ierr = PetscSectionGetConstrainedStorageSize(section, &sp->spintdim);CHKERRQ(ierr);
2577   ierr = PetscDualSpaceComputeFunctionalsFromAllData(sp);CHKERRQ(ierr);
2578   ierr = PetscFree2(pStratStart, pStratEnd);CHKERRQ(ierr);
2579   ierr = DMDestroy(&dmint);CHKERRQ(ierr);
2580   PetscFunctionReturn(0);
2581 }
2582 
2583 /* Create a matrix that represents the transformation that DMPlexVecGetClosure() would need
2584  * to get the representation of the dofs for a mesh point if the mesh point had this orientation
2585  * relative to the cell */
2586 PetscErrorCode PetscDualSpaceCreateInteriorSymmetryMatrix_Lagrange(PetscDualSpace sp, PetscInt ornt, Mat *symMat)
2587 {
2588   PetscDualSpace_Lag *lag;
2589   DM dm;
2590   PetscLagNodeIndices vertIndices, intNodeIndices;
2591   PetscLagNodeIndices ni;
2592   PetscInt nodeIdxDim, nodeVecDim, nNodes;
2593   PetscInt formDegree;
2594   PetscInt *perm, *permOrnt;
2595   PetscInt *nnz;
2596   PetscInt n;
2597   PetscInt maxGroupSize;
2598   PetscScalar *V, *W, *work;
2599   Mat A;
2600   PetscErrorCode ierr;
2601 
2602   PetscFunctionBegin;
2603   if (!sp->spintdim) {
2604     *symMat = NULL;
2605     PetscFunctionReturn(0);
2606   }
2607   lag = (PetscDualSpace_Lag *) sp->data;
2608   vertIndices = lag->vertIndices;
2609   intNodeIndices = lag->intNodeIndices;
2610   ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr);
2611   ierr = PetscDualSpaceGetFormDegree(sp, &formDegree);CHKERRQ(ierr);
2612   ierr = PetscNew(&ni);CHKERRQ(ierr);
2613   ni->refct = 1;
2614   ni->nodeIdxDim = nodeIdxDim = intNodeIndices->nodeIdxDim;
2615   ni->nodeVecDim = nodeVecDim = intNodeIndices->nodeVecDim;
2616   ni->nNodes = nNodes = intNodeIndices->nNodes;
2617   ierr = PetscMalloc1(nNodes * nodeIdxDim, &(ni->nodeIdx));CHKERRQ(ierr);
2618   ierr = PetscMalloc1(nNodes * nodeVecDim, &(ni->nodeVec));CHKERRQ(ierr);
2619   /* push forward the dofs by the symmetry of the reference element induced by ornt */
2620   ierr = PetscLagNodeIndicesPushForward(dm, vertIndices, 0, vertIndices, intNodeIndices, ornt, formDegree, ni->nodeIdx, ni->nodeVec);CHKERRQ(ierr);
2621   /* get the revlex order for both the original and transformed dofs */
2622   ierr = PetscLagNodeIndicesGetPermutation(intNodeIndices, &perm);CHKERRQ(ierr);
2623   ierr = PetscLagNodeIndicesGetPermutation(ni, &permOrnt);CHKERRQ(ierr);
2624   ierr = PetscMalloc1(nNodes, &nnz);CHKERRQ(ierr);
2625   for (n = 0, maxGroupSize = 0; n < nNodes;) { /* incremented in the loop */
2626     PetscInt *nind = &(ni->nodeIdx[permOrnt[n] * nodeIdxDim]);
2627     PetscInt m, nEnd;
2628     PetscInt groupSize;
2629     /* for each group of dofs that have the same nodeIdx coordinate */
2630     for (nEnd = n + 1; nEnd < nNodes; nEnd++) {
2631       PetscInt *mind = &(ni->nodeIdx[permOrnt[nEnd] * nodeIdxDim]);
2632       PetscInt d;
2633 
2634       /* compare the oriented permutation indices */
2635       for (d = 0; d < nodeIdxDim; d++) if (mind[d] != nind[d]) break;
2636       if (d < nodeIdxDim) break;
2637     }
2638     /* permOrnt[[n, nEnd)] is a group of dofs that, under the symmetry are at the same location */
2639 
2640     /* the symmetry had better map the group of dofs with the same permuted nodeIdx
2641      * to a group of dofs with the same size, otherwise we messed up */
2642     if (PetscDefined(USE_DEBUG)) {
2643       PetscInt m;
2644       PetscInt *nind = &(intNodeIndices->nodeIdx[perm[n] * nodeIdxDim]);
2645 
2646       for (m = n + 1; m < nEnd; m++) {
2647         PetscInt *mind = &(intNodeIndices->nodeIdx[perm[m] * nodeIdxDim]);
2648         PetscInt d;
2649 
2650         /* compare the oriented permutation indices */
2651         for (d = 0; d < nodeIdxDim; d++) if (mind[d] != nind[d]) break;
2652         if (d < nodeIdxDim) break;
2653       }
2654       PetscCheckFalse(m < nEnd,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Dofs with same index after symmetry not same block size");
2655     }
2656     groupSize = nEnd - n;
2657     /* each pushforward dof vector will be expressed in a basis of the unpermuted dofs */
2658     for (m = n; m < nEnd; m++) nnz[permOrnt[m]] = groupSize;
2659 
2660     maxGroupSize = PetscMax(maxGroupSize, nEnd - n);
2661     n = nEnd;
2662   }
2663   PetscCheckFalse(maxGroupSize > nodeVecDim,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Dofs are not in blocks that can be solved");
2664   ierr = MatCreateSeqAIJ(PETSC_COMM_SELF, nNodes, nNodes, 0, nnz, &A);CHKERRQ(ierr);
2665   ierr = PetscFree(nnz);CHKERRQ(ierr);
2666   ierr = PetscMalloc3(maxGroupSize * nodeVecDim, &V, maxGroupSize * nodeVecDim, &W, nodeVecDim * 2, &work);CHKERRQ(ierr);
2667   for (n = 0; n < nNodes;) { /* incremented in the loop */
2668     PetscInt *nind = &(ni->nodeIdx[permOrnt[n] * nodeIdxDim]);
2669     PetscInt nEnd;
2670     PetscInt m;
2671     PetscInt groupSize;
2672     for (nEnd = n + 1; nEnd < nNodes; nEnd++) {
2673       PetscInt *mind = &(ni->nodeIdx[permOrnt[nEnd] * nodeIdxDim]);
2674       PetscInt d;
2675 
2676       /* compare the oriented permutation indices */
2677       for (d = 0; d < nodeIdxDim; d++) if (mind[d] != nind[d]) break;
2678       if (d < nodeIdxDim) break;
2679     }
2680     groupSize = nEnd - n;
2681     /* get all of the vectors from the original and all of the pushforward vectors */
2682     for (m = n; m < nEnd; m++) {
2683       PetscInt d;
2684 
2685       for (d = 0; d < nodeVecDim; d++) {
2686         V[(m - n) * nodeVecDim + d] = intNodeIndices->nodeVec[perm[m] * nodeVecDim + d];
2687         W[(m - n) * nodeVecDim + d] = ni->nodeVec[permOrnt[m] * nodeVecDim + d];
2688       }
2689     }
2690     /* now we have to solve for W in terms of V: the systems isn't always square, but the span
2691      * of V and W should always be the same, so the solution of the normal equations works */
2692     {
2693       char transpose = 'N';
2694       PetscBLASInt bm = nodeVecDim;
2695       PetscBLASInt bn = groupSize;
2696       PetscBLASInt bnrhs = groupSize;
2697       PetscBLASInt blda = bm;
2698       PetscBLASInt bldb = bm;
2699       PetscBLASInt blwork = 2 * nodeVecDim;
2700       PetscBLASInt info;
2701 
2702       PetscStackCallBLAS("LAPACKgels",LAPACKgels_(&transpose,&bm,&bn,&bnrhs,V,&blda,W,&bldb,work,&blwork, &info));
2703       PetscCheckFalse(info != 0,PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to GELS");
2704       /* repack */
2705       {
2706         PetscInt i, j;
2707 
2708         for (i = 0; i < groupSize; i++) {
2709           for (j = 0; j < groupSize; j++) {
2710             /* notice the different leading dimension */
2711             V[i * groupSize + j] = W[i * nodeVecDim + j];
2712           }
2713         }
2714       }
2715       if (PetscDefined(USE_DEBUG)) {
2716         PetscReal res;
2717 
2718         /* check that the normal error is 0 */
2719         for (m = n; m < nEnd; m++) {
2720           PetscInt d;
2721 
2722           for (d = 0; d < nodeVecDim; d++) {
2723             W[(m - n) * nodeVecDim + d] = ni->nodeVec[permOrnt[m] * nodeVecDim + d];
2724           }
2725         }
2726         res = 0.;
2727         for (PetscInt i = 0; i < groupSize; i++) {
2728           for (PetscInt j = 0; j < nodeVecDim; j++) {
2729             for (PetscInt k = 0; k < groupSize; k++) {
2730               W[i * nodeVecDim + j] -= V[i * groupSize + k] * intNodeIndices->nodeVec[perm[n+k] * nodeVecDim + j];
2731             }
2732             res += PetscAbsScalar(W[i * nodeVecDim + j]);
2733           }
2734         }
2735         PetscCheckFalse(res > PETSC_SMALL,PETSC_COMM_SELF,PETSC_ERR_LIB,"Dof block did not solve");
2736       }
2737     }
2738     ierr = MatSetValues(A, groupSize, &permOrnt[n], groupSize, &perm[n], V, INSERT_VALUES);CHKERRQ(ierr);
2739     n = nEnd;
2740   }
2741   ierr = MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2742   ierr = MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2743   *symMat = A;
2744   ierr = PetscFree3(V,W,work);CHKERRQ(ierr);
2745   ierr = PetscLagNodeIndicesDestroy(&ni);CHKERRQ(ierr);
2746   PetscFunctionReturn(0);
2747 }
2748 
2749 #define BaryIndex(perEdge,a,b,c) (((b)*(2*perEdge+1-(b)))/2)+(c)
2750 
2751 #define CartIndex(perEdge,a,b) (perEdge*(a)+b)
2752 
2753 /* the existing interface for symmetries is insufficient for all cases:
2754  * - it should be sufficient for form degrees that are scalar (0 and n)
2755  * - it should be sufficient for hypercube dofs
2756  * - it isn't sufficient for simplex cells with non-scalar form degrees if
2757  *   there are any dofs in the interior
2758  *
2759  * We compute the general transformation matrices, and if they fit, we return them,
2760  * otherwise we error (but we should probably change the interface to allow for
2761  * these symmetries)
2762  */
2763 static PetscErrorCode PetscDualSpaceGetSymmetries_Lagrange(PetscDualSpace sp, const PetscInt ****perms, const PetscScalar ****flips)
2764 {
2765   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data;
2766   PetscInt           dim, order, Nc;
2767   PetscErrorCode     ierr;
2768 
2769   PetscFunctionBegin;
2770   ierr = PetscDualSpaceGetOrder(sp,&order);CHKERRQ(ierr);
2771   ierr = PetscDualSpaceGetNumComponents(sp,&Nc);CHKERRQ(ierr);
2772   ierr = DMGetDimension(sp->dm,&dim);CHKERRQ(ierr);
2773   if (!lag->symComputed) { /* store symmetries */
2774     PetscInt       pStart, pEnd, p;
2775     PetscInt       numPoints;
2776     PetscInt       numFaces;
2777     PetscInt       spintdim;
2778     PetscInt       ***symperms;
2779     PetscScalar    ***symflips;
2780 
2781     ierr = DMPlexGetChart(sp->dm, &pStart, &pEnd);CHKERRQ(ierr);
2782     numPoints = pEnd - pStart;
2783     {
2784       DMPolytopeType ct;
2785       /* The number of arrangements is no longer based on the number of faces */
2786       ierr = DMPlexGetCellType(sp->dm, 0, &ct);CHKERRQ(ierr);
2787       numFaces = DMPolytopeTypeGetNumArrangments(ct) / 2;
2788     }
2789     ierr = PetscCalloc1(numPoints,&symperms);CHKERRQ(ierr);
2790     ierr = PetscCalloc1(numPoints,&symflips);CHKERRQ(ierr);
2791     spintdim = sp->spintdim;
2792     /* The nodal symmetry behavior is not present when tensorSpace != tensorCell: someone might want this for the "S"
2793      * family of FEEC spaces.  Most used in particular are discontinuous polynomial L2 spaces in tensor cells, where
2794      * the symmetries are not necessary for FE assembly.  So for now we assume this is the case and don't return
2795      * symmetries if tensorSpace != tensorCell */
2796     if (spintdim && 0 < dim && dim < 3 && (lag->tensorSpace == lag->tensorCell)) { /* compute self symmetries */
2797       PetscInt **cellSymperms;
2798       PetscScalar **cellSymflips;
2799       PetscInt ornt;
2800       PetscInt nCopies = Nc / lag->intNodeIndices->nodeVecDim;
2801       PetscInt nNodes = lag->intNodeIndices->nNodes;
2802 
2803       lag->numSelfSym = 2 * numFaces;
2804       lag->selfSymOff = numFaces;
2805       ierr = PetscCalloc1(2*numFaces,&cellSymperms);CHKERRQ(ierr);
2806       ierr = PetscCalloc1(2*numFaces,&cellSymflips);CHKERRQ(ierr);
2807       /* we want to be able to index symmetries directly with the orientations, which range from [-numFaces,numFaces) */
2808       symperms[0] = &cellSymperms[numFaces];
2809       symflips[0] = &cellSymflips[numFaces];
2810       PetscCheckFalse(lag->intNodeIndices->nodeVecDim * nCopies != Nc,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Node indices incompatible with dofs");
2811       PetscCheckFalse(nNodes * nCopies != spintdim,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Node indices incompatible with dofs");
2812       for (ornt = -numFaces; ornt < numFaces; ornt++) { /* for every symmetry, compute the symmetry matrix, and extract rows to see if it fits in the perm + flip framework */
2813         Mat symMat;
2814         PetscInt *perm;
2815         PetscScalar *flips;
2816         PetscInt i;
2817 
2818         if (!ornt) continue;
2819         ierr = PetscMalloc1(spintdim, &perm);CHKERRQ(ierr);
2820         ierr = PetscCalloc1(spintdim, &flips);CHKERRQ(ierr);
2821         for (i = 0; i < spintdim; i++) perm[i] = -1;
2822         ierr = PetscDualSpaceCreateInteriorSymmetryMatrix_Lagrange(sp, ornt, &symMat);CHKERRQ(ierr);
2823         for (i = 0; i < nNodes; i++) {
2824           PetscInt ncols;
2825           PetscInt j, k;
2826           const PetscInt *cols;
2827           const PetscScalar *vals;
2828           PetscBool nz_seen = PETSC_FALSE;
2829 
2830           ierr = MatGetRow(symMat, i, &ncols, &cols, &vals);CHKERRQ(ierr);
2831           for (j = 0; j < ncols; j++) {
2832             if (PetscAbsScalar(vals[j]) > PETSC_SMALL) {
2833               PetscCheckFalse(nz_seen,PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "This dual space has symmetries that can't be described as a permutation + sign flips");
2834               nz_seen = PETSC_TRUE;
2835               PetscCheckFalse(PetscAbsReal(PetscAbsScalar(vals[j]) - PetscRealConstant(1.)) > PETSC_SMALL,PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "This dual space has symmetries that can't be described as a permutation + sign flips");
2836               PetscCheckFalse(PetscAbsReal(PetscImaginaryPart(vals[j])) > PETSC_SMALL,PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "This dual space has symmetries that can't be described as a permutation + sign flips");
2837               PetscCheckFalse(perm[cols[j] * nCopies] >= 0,PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "This dual space has symmetries that can't be described as a permutation + sign flips");
2838               for (k = 0; k < nCopies; k++) {
2839                 perm[cols[j] * nCopies + k] = i * nCopies + k;
2840               }
2841               if (PetscRealPart(vals[j]) < 0.) {
2842                 for (k = 0; k < nCopies; k++) {
2843                   flips[i * nCopies + k] = -1.;
2844                 }
2845               } else {
2846                 for (k = 0; k < nCopies; k++) {
2847                   flips[i * nCopies + k] = 1.;
2848                 }
2849               }
2850             }
2851           }
2852           ierr = MatRestoreRow(symMat, i, &ncols, &cols, &vals);CHKERRQ(ierr);
2853         }
2854         ierr = MatDestroy(&symMat);CHKERRQ(ierr);
2855         /* if there were no sign flips, keep NULL */
2856         for (i = 0; i < spintdim; i++) if (flips[i] != 1.) break;
2857         if (i == spintdim) {
2858           ierr = PetscFree(flips);CHKERRQ(ierr);
2859           flips = NULL;
2860         }
2861         /* if the permutation is identity, keep NULL */
2862         for (i = 0; i < spintdim; i++) if (perm[i] != i) break;
2863         if (i == spintdim) {
2864           ierr = PetscFree(perm);CHKERRQ(ierr);
2865           perm = NULL;
2866         }
2867         symperms[0][ornt] = perm;
2868         symflips[0][ornt] = flips;
2869       }
2870       /* if no orientations produced non-identity permutations, keep NULL */
2871       for (ornt = -numFaces; ornt < numFaces; ornt++) if (symperms[0][ornt]) break;
2872       if (ornt == numFaces) {
2873         ierr = PetscFree(cellSymperms);CHKERRQ(ierr);
2874         symperms[0] = NULL;
2875       }
2876       /* if no orientations produced sign flips, keep NULL */
2877       for (ornt = -numFaces; ornt < numFaces; ornt++) if (symflips[0][ornt]) break;
2878       if (ornt == numFaces) {
2879         ierr = PetscFree(cellSymflips);CHKERRQ(ierr);
2880         symflips[0] = NULL;
2881       }
2882     }
2883     { /* get the symmetries of closure points */
2884       PetscInt closureSize = 0;
2885       PetscInt *closure = NULL;
2886       PetscInt r;
2887 
2888       ierr = DMPlexGetTransitiveClosure(sp->dm,0,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
2889       for (r = 0; r < closureSize; r++) {
2890         PetscDualSpace psp;
2891         PetscInt point = closure[2 * r];
2892         PetscInt pspintdim;
2893         const PetscInt ***psymperms = NULL;
2894         const PetscScalar ***psymflips = NULL;
2895 
2896         if (!point) continue;
2897         ierr = PetscDualSpaceGetPointSubspace(sp, point, &psp);CHKERRQ(ierr);
2898         if (!psp) continue;
2899         ierr = PetscDualSpaceGetInteriorDimension(psp, &pspintdim);CHKERRQ(ierr);
2900         if (!pspintdim) continue;
2901         ierr = PetscDualSpaceGetSymmetries(psp,&psymperms,&psymflips);CHKERRQ(ierr);
2902         symperms[r] = (PetscInt **) (psymperms ? psymperms[0] : NULL);
2903         symflips[r] = (PetscScalar **) (psymflips ? psymflips[0] : NULL);
2904       }
2905       ierr = DMPlexRestoreTransitiveClosure(sp->dm,0,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
2906     }
2907     for (p = 0; p < pEnd; p++) if (symperms[p]) break;
2908     if (p == pEnd) {
2909       ierr = PetscFree(symperms);CHKERRQ(ierr);
2910       symperms = NULL;
2911     }
2912     for (p = 0; p < pEnd; p++) if (symflips[p]) break;
2913     if (p == pEnd) {
2914       ierr = PetscFree(symflips);CHKERRQ(ierr);
2915       symflips = NULL;
2916     }
2917     lag->symperms = symperms;
2918     lag->symflips = symflips;
2919     lag->symComputed = PETSC_TRUE;
2920   }
2921   if (perms) *perms = (const PetscInt ***) lag->symperms;
2922   if (flips) *flips = (const PetscScalar ***) lag->symflips;
2923   PetscFunctionReturn(0);
2924 }
2925 
2926 static PetscErrorCode PetscDualSpaceLagrangeGetContinuity_Lagrange(PetscDualSpace sp, PetscBool *continuous)
2927 {
2928   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data;
2929 
2930   PetscFunctionBegin;
2931   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
2932   PetscValidPointer(continuous, 2);
2933   *continuous = lag->continuous;
2934   PetscFunctionReturn(0);
2935 }
2936 
2937 static PetscErrorCode PetscDualSpaceLagrangeSetContinuity_Lagrange(PetscDualSpace sp, PetscBool continuous)
2938 {
2939   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data;
2940 
2941   PetscFunctionBegin;
2942   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
2943   lag->continuous = continuous;
2944   PetscFunctionReturn(0);
2945 }
2946 
2947 /*@
2948   PetscDualSpaceLagrangeGetContinuity - Retrieves the flag for element continuity
2949 
2950   Not Collective
2951 
2952   Input Parameter:
2953 . sp         - the PetscDualSpace
2954 
2955   Output Parameter:
2956 . continuous - flag for element continuity
2957 
2958   Level: intermediate
2959 
2960 .seealso: PetscDualSpaceLagrangeSetContinuity()
2961 @*/
2962 PetscErrorCode PetscDualSpaceLagrangeGetContinuity(PetscDualSpace sp, PetscBool *continuous)
2963 {
2964   PetscErrorCode ierr;
2965 
2966   PetscFunctionBegin;
2967   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
2968   PetscValidPointer(continuous, 2);
2969   ierr = PetscTryMethod(sp, "PetscDualSpaceLagrangeGetContinuity_C", (PetscDualSpace,PetscBool*),(sp,continuous));CHKERRQ(ierr);
2970   PetscFunctionReturn(0);
2971 }
2972 
2973 /*@
2974   PetscDualSpaceLagrangeSetContinuity - Indicate whether the element is continuous
2975 
2976   Logically Collective on sp
2977 
2978   Input Parameters:
2979 + sp         - the PetscDualSpace
2980 - continuous - flag for element continuity
2981 
2982   Options Database:
2983 . -petscdualspace_lagrange_continuity <bool>
2984 
2985   Level: intermediate
2986 
2987 .seealso: PetscDualSpaceLagrangeGetContinuity()
2988 @*/
2989 PetscErrorCode PetscDualSpaceLagrangeSetContinuity(PetscDualSpace sp, PetscBool continuous)
2990 {
2991   PetscErrorCode ierr;
2992 
2993   PetscFunctionBegin;
2994   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
2995   PetscValidLogicalCollectiveBool(sp, continuous, 2);
2996   ierr = PetscTryMethod(sp, "PetscDualSpaceLagrangeSetContinuity_C", (PetscDualSpace,PetscBool),(sp,continuous));CHKERRQ(ierr);
2997   PetscFunctionReturn(0);
2998 }
2999 
3000 static PetscErrorCode PetscDualSpaceLagrangeGetTensor_Lagrange(PetscDualSpace sp, PetscBool *tensor)
3001 {
3002   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data;
3003 
3004   PetscFunctionBegin;
3005   *tensor = lag->tensorSpace;
3006   PetscFunctionReturn(0);
3007 }
3008 
3009 static PetscErrorCode PetscDualSpaceLagrangeSetTensor_Lagrange(PetscDualSpace sp, PetscBool tensor)
3010 {
3011   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data;
3012 
3013   PetscFunctionBegin;
3014   lag->tensorSpace = tensor;
3015   PetscFunctionReturn(0);
3016 }
3017 
3018 static PetscErrorCode PetscDualSpaceLagrangeGetTrimmed_Lagrange(PetscDualSpace sp, PetscBool *trimmed)
3019 {
3020   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data;
3021 
3022   PetscFunctionBegin;
3023   *trimmed = lag->trimmed;
3024   PetscFunctionReturn(0);
3025 }
3026 
3027 static PetscErrorCode PetscDualSpaceLagrangeSetTrimmed_Lagrange(PetscDualSpace sp, PetscBool trimmed)
3028 {
3029   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data;
3030 
3031   PetscFunctionBegin;
3032   lag->trimmed = trimmed;
3033   PetscFunctionReturn(0);
3034 }
3035 
3036 static PetscErrorCode PetscDualSpaceLagrangeGetNodeType_Lagrange(PetscDualSpace sp, PetscDTNodeType *nodeType, PetscBool *boundary, PetscReal *exponent)
3037 {
3038   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data;
3039 
3040   PetscFunctionBegin;
3041   if (nodeType) *nodeType = lag->nodeType;
3042   if (boundary) *boundary = lag->endNodes;
3043   if (exponent) *exponent = lag->nodeExponent;
3044   PetscFunctionReturn(0);
3045 }
3046 
3047 static PetscErrorCode PetscDualSpaceLagrangeSetNodeType_Lagrange(PetscDualSpace sp, PetscDTNodeType nodeType, PetscBool boundary, PetscReal exponent)
3048 {
3049   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data;
3050 
3051   PetscFunctionBegin;
3052   PetscCheckFalse(nodeType == PETSCDTNODES_GAUSSJACOBI && exponent <= -1.,PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_OUTOFRANGE, "Exponent must be > -1");
3053   lag->nodeType = nodeType;
3054   lag->endNodes = boundary;
3055   lag->nodeExponent = exponent;
3056   PetscFunctionReturn(0);
3057 }
3058 
3059 static PetscErrorCode PetscDualSpaceLagrangeGetUseMoments_Lagrange(PetscDualSpace sp, PetscBool *useMoments)
3060 {
3061   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data;
3062 
3063   PetscFunctionBegin;
3064   *useMoments = lag->useMoments;
3065   PetscFunctionReturn(0);
3066 }
3067 
3068 static PetscErrorCode PetscDualSpaceLagrangeSetUseMoments_Lagrange(PetscDualSpace sp, PetscBool useMoments)
3069 {
3070   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data;
3071 
3072   PetscFunctionBegin;
3073   lag->useMoments = useMoments;
3074   PetscFunctionReturn(0);
3075 }
3076 
3077 static PetscErrorCode PetscDualSpaceLagrangeGetMomentOrder_Lagrange(PetscDualSpace sp, PetscInt *momentOrder)
3078 {
3079   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data;
3080 
3081   PetscFunctionBegin;
3082   *momentOrder = lag->momentOrder;
3083   PetscFunctionReturn(0);
3084 }
3085 
3086 static PetscErrorCode PetscDualSpaceLagrangeSetMomentOrder_Lagrange(PetscDualSpace sp, PetscInt momentOrder)
3087 {
3088   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data;
3089 
3090   PetscFunctionBegin;
3091   lag->momentOrder = momentOrder;
3092   PetscFunctionReturn(0);
3093 }
3094 
3095 /*@
3096   PetscDualSpaceLagrangeGetTensor - Get the tensor nature of the dual space
3097 
3098   Not collective
3099 
3100   Input Parameter:
3101 . sp - The PetscDualSpace
3102 
3103   Output Parameter:
3104 . tensor - Whether the dual space has tensor layout (vs. simplicial)
3105 
3106   Level: intermediate
3107 
3108 .seealso: PetscDualSpaceLagrangeSetTensor(), PetscDualSpaceCreate()
3109 @*/
3110 PetscErrorCode PetscDualSpaceLagrangeGetTensor(PetscDualSpace sp, PetscBool *tensor)
3111 {
3112   PetscErrorCode ierr;
3113 
3114   PetscFunctionBegin;
3115   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
3116   PetscValidPointer(tensor, 2);
3117   ierr = PetscTryMethod(sp,"PetscDualSpaceLagrangeGetTensor_C",(PetscDualSpace,PetscBool *),(sp,tensor));CHKERRQ(ierr);
3118   PetscFunctionReturn(0);
3119 }
3120 
3121 /*@
3122   PetscDualSpaceLagrangeSetTensor - Set the tensor nature of the dual space
3123 
3124   Not collective
3125 
3126   Input Parameters:
3127 + sp - The PetscDualSpace
3128 - tensor - Whether the dual space has tensor layout (vs. simplicial)
3129 
3130   Level: intermediate
3131 
3132 .seealso: PetscDualSpaceLagrangeGetTensor(), PetscDualSpaceCreate()
3133 @*/
3134 PetscErrorCode PetscDualSpaceLagrangeSetTensor(PetscDualSpace sp, PetscBool tensor)
3135 {
3136   PetscErrorCode ierr;
3137 
3138   PetscFunctionBegin;
3139   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
3140   ierr = PetscTryMethod(sp,"PetscDualSpaceLagrangeSetTensor_C",(PetscDualSpace,PetscBool),(sp,tensor));CHKERRQ(ierr);
3141   PetscFunctionReturn(0);
3142 }
3143 
3144 /*@
3145   PetscDualSpaceLagrangeGetTrimmed - Get the trimmed nature of the dual space
3146 
3147   Not collective
3148 
3149   Input Parameter:
3150 . sp - The PetscDualSpace
3151 
3152   Output Parameter:
3153 . trimmed - Whether the dual space represents to dual basis of a trimmed polynomial space (e.g. Raviart-Thomas and higher order / other form degree variants)
3154 
3155   Level: intermediate
3156 
3157 .seealso: PetscDualSpaceLagrangeSetTrimmed(), PetscDualSpaceCreate()
3158 @*/
3159 PetscErrorCode PetscDualSpaceLagrangeGetTrimmed(PetscDualSpace sp, PetscBool *trimmed)
3160 {
3161   PetscErrorCode ierr;
3162 
3163   PetscFunctionBegin;
3164   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
3165   PetscValidPointer(trimmed, 2);
3166   ierr = PetscTryMethod(sp,"PetscDualSpaceLagrangeGetTrimmed_C",(PetscDualSpace,PetscBool *),(sp,trimmed));CHKERRQ(ierr);
3167   PetscFunctionReturn(0);
3168 }
3169 
3170 /*@
3171   PetscDualSpaceLagrangeSetTrimmed - Set the trimmed nature of the dual space
3172 
3173   Not collective
3174 
3175   Input Parameters:
3176 + sp - The PetscDualSpace
3177 - trimmed - Whether the dual space represents to dual basis of a trimmed polynomial space (e.g. Raviart-Thomas and higher order / other form degree variants)
3178 
3179   Level: intermediate
3180 
3181 .seealso: PetscDualSpaceLagrangeGetTrimmed(), PetscDualSpaceCreate()
3182 @*/
3183 PetscErrorCode PetscDualSpaceLagrangeSetTrimmed(PetscDualSpace sp, PetscBool trimmed)
3184 {
3185   PetscErrorCode ierr;
3186 
3187   PetscFunctionBegin;
3188   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
3189   ierr = PetscTryMethod(sp,"PetscDualSpaceLagrangeSetTrimmed_C",(PetscDualSpace,PetscBool),(sp,trimmed));CHKERRQ(ierr);
3190   PetscFunctionReturn(0);
3191 }
3192 
3193 /*@
3194   PetscDualSpaceLagrangeGetNodeType - Get a description of how nodes are laid out for Lagrange polynomials in this
3195   dual space
3196 
3197   Not collective
3198 
3199   Input Parameter:
3200 . sp - The PetscDualSpace
3201 
3202   Output Parameters:
3203 + nodeType - The type of nodes
3204 . boundary - Whether the node type is one that includes endpoints (if nodeType is PETSCDTNODES_GAUSSJACOBI, nodes that
3205              include the boundary are Gauss-Lobatto-Jacobi nodes)
3206 - exponent - If nodeType is PETSCDTNODES_GAUSJACOBI, indicates the exponent used for both ends of the 1D Jacobi weight function
3207              '0' is Gauss-Legendre, '-0.5' is Gauss-Chebyshev of the first type, '0.5' is Gauss-Chebyshev of the second type
3208 
3209   Level: advanced
3210 
3211 .seealso: PetscDTNodeType, PetscDualSpaceLagrangeSetNodeType()
3212 @*/
3213 PetscErrorCode PetscDualSpaceLagrangeGetNodeType(PetscDualSpace sp, PetscDTNodeType *nodeType, PetscBool *boundary, PetscReal *exponent)
3214 {
3215   PetscErrorCode ierr;
3216 
3217   PetscFunctionBegin;
3218   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
3219   if (nodeType) PetscValidPointer(nodeType, 2);
3220   if (boundary) PetscValidPointer(boundary, 3);
3221   if (exponent) PetscValidPointer(exponent, 4);
3222   ierr = PetscTryMethod(sp,"PetscDualSpaceLagrangeGetNodeType_C",(PetscDualSpace,PetscDTNodeType *,PetscBool *,PetscReal *),(sp,nodeType,boundary,exponent));CHKERRQ(ierr);
3223   PetscFunctionReturn(0);
3224 }
3225 
3226 /*@
3227   PetscDualSpaceLagrangeSetNodeType - Set a description of how nodes are laid out for Lagrange polynomials in this
3228   dual space
3229 
3230   Logically collective
3231 
3232   Input Parameters:
3233 + sp - The PetscDualSpace
3234 . nodeType - The type of nodes
3235 . boundary - Whether the node type is one that includes endpoints (if nodeType is PETSCDTNODES_GAUSSJACOBI, nodes that
3236              include the boundary are Gauss-Lobatto-Jacobi nodes)
3237 - exponent - If nodeType is PETSCDTNODES_GAUSJACOBI, indicates the exponent used for both ends of the 1D Jacobi weight function
3238              '0' is Gauss-Legendre, '-0.5' is Gauss-Chebyshev of the first type, '0.5' is Gauss-Chebyshev of the second type
3239 
3240   Level: advanced
3241 
3242 .seealso: PetscDTNodeType, PetscDualSpaceLagrangeGetNodeType()
3243 @*/
3244 PetscErrorCode PetscDualSpaceLagrangeSetNodeType(PetscDualSpace sp, PetscDTNodeType nodeType, PetscBool boundary, PetscReal exponent)
3245 {
3246   PetscErrorCode ierr;
3247 
3248   PetscFunctionBegin;
3249   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
3250   ierr = PetscTryMethod(sp,"PetscDualSpaceLagrangeSetNodeType_C",(PetscDualSpace,PetscDTNodeType,PetscBool,PetscReal),(sp,nodeType,boundary,exponent));CHKERRQ(ierr);
3251   PetscFunctionReturn(0);
3252 }
3253 
3254 /*@
3255   PetscDualSpaceLagrangeGetUseMoments - Get the flag for using moment functionals
3256 
3257   Not collective
3258 
3259   Input Parameter:
3260 . sp - The PetscDualSpace
3261 
3262   Output Parameter:
3263 . useMoments - Moment flag
3264 
3265   Level: advanced
3266 
3267 .seealso: PetscDualSpaceLagrangeSetUseMoments()
3268 @*/
3269 PetscErrorCode PetscDualSpaceLagrangeGetUseMoments(PetscDualSpace sp, PetscBool *useMoments)
3270 {
3271   PetscErrorCode ierr;
3272 
3273   PetscFunctionBegin;
3274   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
3275   PetscValidBoolPointer(useMoments, 2);
3276   ierr = PetscUseMethod(sp,"PetscDualSpaceLagrangeGetUseMoments_C",(PetscDualSpace,PetscBool *),(sp,useMoments));CHKERRQ(ierr);
3277   PetscFunctionReturn(0);
3278 }
3279 
3280 /*@
3281   PetscDualSpaceLagrangeSetUseMoments - Set the flag for moment functionals
3282 
3283   Logically collective
3284 
3285   Input Parameters:
3286 + sp - The PetscDualSpace
3287 - useMoments - The flag for moment functionals
3288 
3289   Level: advanced
3290 
3291 .seealso: PetscDualSpaceLagrangeGetUseMoments()
3292 @*/
3293 PetscErrorCode PetscDualSpaceLagrangeSetUseMoments(PetscDualSpace sp, PetscBool useMoments)
3294 {
3295   PetscErrorCode ierr;
3296 
3297   PetscFunctionBegin;
3298   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
3299   ierr = PetscTryMethod(sp,"PetscDualSpaceLagrangeSetUseMoments_C",(PetscDualSpace,PetscBool),(sp,useMoments));CHKERRQ(ierr);
3300   PetscFunctionReturn(0);
3301 }
3302 
3303 /*@
3304   PetscDualSpaceLagrangeGetMomentOrder - Get the order for moment integration
3305 
3306   Not collective
3307 
3308   Input Parameter:
3309 . sp - The PetscDualSpace
3310 
3311   Output Parameter:
3312 . order - Moment integration order
3313 
3314   Level: advanced
3315 
3316 .seealso: PetscDualSpaceLagrangeSetMomentOrder()
3317 @*/
3318 PetscErrorCode PetscDualSpaceLagrangeGetMomentOrder(PetscDualSpace sp, PetscInt *order)
3319 {
3320   PetscErrorCode ierr;
3321 
3322   PetscFunctionBegin;
3323   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
3324   PetscValidIntPointer(order, 2);
3325   ierr = PetscUseMethod(sp,"PetscDualSpaceLagrangeGetMomentOrder_C",(PetscDualSpace,PetscInt *),(sp,order));CHKERRQ(ierr);
3326   PetscFunctionReturn(0);
3327 }
3328 
3329 /*@
3330   PetscDualSpaceLagrangeSetMomentOrder - Set the order for moment integration
3331 
3332   Logically collective
3333 
3334   Input Parameters:
3335 + sp - The PetscDualSpace
3336 - order - The order for moment integration
3337 
3338   Level: advanced
3339 
3340 .seealso: PetscDualSpaceLagrangeGetMomentOrder()
3341 @*/
3342 PetscErrorCode PetscDualSpaceLagrangeSetMomentOrder(PetscDualSpace sp, PetscInt order)
3343 {
3344   PetscErrorCode ierr;
3345 
3346   PetscFunctionBegin;
3347   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
3348   ierr = PetscTryMethod(sp,"PetscDualSpaceLagrangeSetMomentOrder_C",(PetscDualSpace,PetscInt),(sp,order));CHKERRQ(ierr);
3349   PetscFunctionReturn(0);
3350 }
3351 
3352 static PetscErrorCode PetscDualSpaceInitialize_Lagrange(PetscDualSpace sp)
3353 {
3354   PetscFunctionBegin;
3355   sp->ops->destroy              = PetscDualSpaceDestroy_Lagrange;
3356   sp->ops->view                 = PetscDualSpaceView_Lagrange;
3357   sp->ops->setfromoptions       = PetscDualSpaceSetFromOptions_Lagrange;
3358   sp->ops->duplicate            = PetscDualSpaceDuplicate_Lagrange;
3359   sp->ops->setup                = PetscDualSpaceSetUp_Lagrange;
3360   sp->ops->createheightsubspace = NULL;
3361   sp->ops->createpointsubspace  = NULL;
3362   sp->ops->getsymmetries        = PetscDualSpaceGetSymmetries_Lagrange;
3363   sp->ops->apply                = PetscDualSpaceApplyDefault;
3364   sp->ops->applyall             = PetscDualSpaceApplyAllDefault;
3365   sp->ops->applyint             = PetscDualSpaceApplyInteriorDefault;
3366   sp->ops->createalldata        = PetscDualSpaceCreateAllDataDefault;
3367   sp->ops->createintdata        = PetscDualSpaceCreateInteriorDataDefault;
3368   PetscFunctionReturn(0);
3369 }
3370 
3371 /*MC
3372   PETSCDUALSPACELAGRANGE = "lagrange" - A PetscDualSpace object that encapsulates a dual space of pointwise evaluation functionals
3373 
3374   Level: intermediate
3375 
3376 .seealso: PetscDualSpaceType, PetscDualSpaceCreate(), PetscDualSpaceSetType()
3377 M*/
3378 PETSC_EXTERN PetscErrorCode PetscDualSpaceCreate_Lagrange(PetscDualSpace sp)
3379 {
3380   PetscDualSpace_Lag *lag;
3381   PetscErrorCode      ierr;
3382 
3383   PetscFunctionBegin;
3384   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
3385   ierr     = PetscNewLog(sp,&lag);CHKERRQ(ierr);
3386   sp->data = lag;
3387 
3388   lag->tensorCell  = PETSC_FALSE;
3389   lag->tensorSpace = PETSC_FALSE;
3390   lag->continuous  = PETSC_TRUE;
3391   lag->numCopies   = PETSC_DEFAULT;
3392   lag->numNodeSkip = PETSC_DEFAULT;
3393   lag->nodeType    = PETSCDTNODES_DEFAULT;
3394   lag->useMoments  = PETSC_FALSE;
3395   lag->momentOrder = 0;
3396 
3397   ierr = PetscDualSpaceInitialize_Lagrange(sp);CHKERRQ(ierr);
3398   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetContinuity_C", PetscDualSpaceLagrangeGetContinuity_Lagrange);CHKERRQ(ierr);
3399   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetContinuity_C", PetscDualSpaceLagrangeSetContinuity_Lagrange);CHKERRQ(ierr);
3400   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetTensor_C", PetscDualSpaceLagrangeGetTensor_Lagrange);CHKERRQ(ierr);
3401   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetTensor_C", PetscDualSpaceLagrangeSetTensor_Lagrange);CHKERRQ(ierr);
3402   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetTrimmed_C", PetscDualSpaceLagrangeGetTrimmed_Lagrange);CHKERRQ(ierr);
3403   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetTrimmed_C", PetscDualSpaceLagrangeSetTrimmed_Lagrange);CHKERRQ(ierr);
3404   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetNodeType_C", PetscDualSpaceLagrangeGetNodeType_Lagrange);CHKERRQ(ierr);
3405   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetNodeType_C", PetscDualSpaceLagrangeSetNodeType_Lagrange);CHKERRQ(ierr);
3406   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetUseMoments_C", PetscDualSpaceLagrangeGetUseMoments_Lagrange);CHKERRQ(ierr);
3407   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetUseMoments_C", PetscDualSpaceLagrangeSetUseMoments_Lagrange);CHKERRQ(ierr);
3408   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetMomentOrder_C", PetscDualSpaceLagrangeGetMomentOrder_Lagrange);CHKERRQ(ierr);
3409   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetMomentOrder_C", PetscDualSpaceLagrangeSetMomentOrder_Lagrange);CHKERRQ(ierr);
3410   PetscFunctionReturn(0);
3411 }
3412