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