xref: /petsc/src/dm/dt/dualspace/impls/lagrange/dspacelagrange.c (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
120cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
220cf1dd8SToby Isaac #include <petscdmplex.h>
33f27d899SToby Isaac #include <petscblaslapack.h>
43f27d899SToby Isaac 
53f27d899SToby Isaac PetscErrorCode DMPlexGetTransitiveClosure_Internal(DM, PetscInt, PetscInt, PetscBool, PetscInt *, PetscInt *[]);
63f27d899SToby Isaac 
79371c9d4SSatish Balay struct _n_Petsc1DNodeFamily {
83f27d899SToby Isaac   PetscInt        refct;
93f27d899SToby Isaac   PetscDTNodeType nodeFamily;
103f27d899SToby Isaac   PetscReal       gaussJacobiExp;
113f27d899SToby Isaac   PetscInt        nComputed;
123f27d899SToby Isaac   PetscReal     **nodesets;
133f27d899SToby Isaac   PetscBool       endpoints;
143f27d899SToby Isaac };
153f27d899SToby Isaac 
1677f1a120SToby Isaac /* users set node families for PETSCDUALSPACELAGRANGE with just the inputs to this function, but internally we create
1777f1a120SToby Isaac  * an object that can cache the computations across multiple dual spaces */
189371c9d4SSatish Balay static PetscErrorCode Petsc1DNodeFamilyCreate(PetscDTNodeType family, PetscReal gaussJacobiExp, PetscBool endpoints, Petsc1DNodeFamily *nf) {
193f27d899SToby Isaac   Petsc1DNodeFamily f;
203f27d899SToby Isaac 
213f27d899SToby Isaac   PetscFunctionBegin;
229566063dSJacob Faibussowitsch   PetscCall(PetscNew(&f));
233f27d899SToby Isaac   switch (family) {
243f27d899SToby Isaac   case PETSCDTNODES_GAUSSJACOBI:
259371c9d4SSatish Balay   case PETSCDTNODES_EQUISPACED: f->nodeFamily = family; break;
263f27d899SToby Isaac   default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unknown 1D node family");
273f27d899SToby Isaac   }
283f27d899SToby Isaac   f->endpoints      = endpoints;
293f27d899SToby Isaac   f->gaussJacobiExp = 0.;
303f27d899SToby Isaac   if (family == PETSCDTNODES_GAUSSJACOBI) {
3108401ef6SPierre Jolivet     PetscCheck(gaussJacobiExp > -1., PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Gauss-Jacobi exponent must be > -1.");
323f27d899SToby Isaac     f->gaussJacobiExp = gaussJacobiExp;
333f27d899SToby Isaac   }
343f27d899SToby Isaac   f->refct = 1;
353f27d899SToby Isaac   *nf      = f;
363f27d899SToby Isaac   PetscFunctionReturn(0);
373f27d899SToby Isaac }
383f27d899SToby Isaac 
399371c9d4SSatish Balay static PetscErrorCode Petsc1DNodeFamilyReference(Petsc1DNodeFamily nf) {
403f27d899SToby Isaac   PetscFunctionBegin;
413f27d899SToby Isaac   if (nf) nf->refct++;
423f27d899SToby Isaac   PetscFunctionReturn(0);
433f27d899SToby Isaac }
443f27d899SToby Isaac 
459371c9d4SSatish Balay static PetscErrorCode Petsc1DNodeFamilyDestroy(Petsc1DNodeFamily *nf) {
463f27d899SToby Isaac   PetscInt i, nc;
473f27d899SToby Isaac 
483f27d899SToby Isaac   PetscFunctionBegin;
493f27d899SToby Isaac   if (!(*nf)) PetscFunctionReturn(0);
503f27d899SToby Isaac   if (--(*nf)->refct > 0) {
513f27d899SToby Isaac     *nf = NULL;
523f27d899SToby Isaac     PetscFunctionReturn(0);
533f27d899SToby Isaac   }
543f27d899SToby Isaac   nc = (*nf)->nComputed;
5548a46eb9SPierre Jolivet   for (i = 0; i < nc; i++) PetscCall(PetscFree((*nf)->nodesets[i]));
569566063dSJacob Faibussowitsch   PetscCall(PetscFree((*nf)->nodesets));
579566063dSJacob Faibussowitsch   PetscCall(PetscFree(*nf));
583f27d899SToby Isaac   *nf = NULL;
593f27d899SToby Isaac   PetscFunctionReturn(0);
603f27d899SToby Isaac }
613f27d899SToby Isaac 
629371c9d4SSatish Balay static PetscErrorCode Petsc1DNodeFamilyGetNodeSets(Petsc1DNodeFamily f, PetscInt degree, PetscReal ***nodesets) {
633f27d899SToby Isaac   PetscInt nc;
643f27d899SToby Isaac 
653f27d899SToby Isaac   PetscFunctionBegin;
663f27d899SToby Isaac   nc = f->nComputed;
673f27d899SToby Isaac   if (degree >= nc) {
683f27d899SToby Isaac     PetscInt    i, j;
693f27d899SToby Isaac     PetscReal **new_nodesets;
703f27d899SToby Isaac     PetscReal  *w;
713f27d899SToby Isaac 
729566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(degree + 1, &new_nodesets));
739566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(new_nodesets, f->nodesets, nc));
749566063dSJacob Faibussowitsch     PetscCall(PetscFree(f->nodesets));
753f27d899SToby Isaac     f->nodesets = new_nodesets;
769566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(degree + 1, &w));
773f27d899SToby Isaac     for (i = nc; i < degree + 1; i++) {
789566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(i + 1, &(f->nodesets[i])));
793f27d899SToby Isaac       if (!i) {
803f27d899SToby Isaac         f->nodesets[i][0] = 0.5;
813f27d899SToby Isaac       } else {
823f27d899SToby Isaac         switch (f->nodeFamily) {
833f27d899SToby Isaac         case PETSCDTNODES_EQUISPACED:
843f27d899SToby Isaac           if (f->endpoints) {
853f27d899SToby Isaac             for (j = 0; j <= i; j++) f->nodesets[i][j] = (PetscReal)j / (PetscReal)i;
863f27d899SToby Isaac           } else {
8777f1a120SToby Isaac             /* these nodes are at the centroids of the small simplices created by the equispaced nodes that include
8877f1a120SToby Isaac              * the endpoints */
893f27d899SToby Isaac             for (j = 0; j <= i; j++) f->nodesets[i][j] = ((PetscReal)j + 0.5) / ((PetscReal)i + 1.);
903f27d899SToby Isaac           }
913f27d899SToby Isaac           break;
923f27d899SToby Isaac         case PETSCDTNODES_GAUSSJACOBI:
933f27d899SToby Isaac           if (f->endpoints) {
949566063dSJacob Faibussowitsch             PetscCall(PetscDTGaussLobattoJacobiQuadrature(i + 1, 0., 1., f->gaussJacobiExp, f->gaussJacobiExp, f->nodesets[i], w));
953f27d899SToby Isaac           } else {
969566063dSJacob Faibussowitsch             PetscCall(PetscDTGaussJacobiQuadrature(i + 1, 0., 1., f->gaussJacobiExp, f->gaussJacobiExp, f->nodesets[i], w));
973f27d899SToby Isaac           }
983f27d899SToby Isaac           break;
993f27d899SToby Isaac         default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unknown 1D node family");
1003f27d899SToby Isaac         }
1013f27d899SToby Isaac       }
1023f27d899SToby Isaac     }
1039566063dSJacob Faibussowitsch     PetscCall(PetscFree(w));
1043f27d899SToby Isaac     f->nComputed = degree + 1;
1053f27d899SToby Isaac   }
1063f27d899SToby Isaac   *nodesets = f->nodesets;
1073f27d899SToby Isaac   PetscFunctionReturn(0);
1083f27d899SToby Isaac }
1093f27d899SToby Isaac 
11077f1a120SToby Isaac /* http://arxiv.org/abs/2002.09421 for details */
1119371c9d4SSatish Balay static PetscErrorCode PetscNodeRecursive_Internal(PetscInt dim, PetscInt degree, PetscReal **nodesets, PetscInt tup[], PetscReal node[]) {
1123f27d899SToby Isaac   PetscReal w;
1133f27d899SToby Isaac   PetscInt  i, j;
1143f27d899SToby Isaac 
1153f27d899SToby Isaac   PetscFunctionBeginHot;
1163f27d899SToby Isaac   w = 0.;
1173f27d899SToby Isaac   if (dim == 1) {
1183f27d899SToby Isaac     node[0] = nodesets[degree][tup[0]];
1193f27d899SToby Isaac     node[1] = nodesets[degree][tup[1]];
1203f27d899SToby Isaac   } else {
1213f27d899SToby Isaac     for (i = 0; i < dim + 1; i++) node[i] = 0.;
1223f27d899SToby Isaac     for (i = 0; i < dim + 1; i++) {
1233f27d899SToby Isaac       PetscReal wi = nodesets[degree][degree - tup[i]];
1243f27d899SToby Isaac 
1253f27d899SToby Isaac       for (j = 0; j < dim + 1; j++) tup[dim + 1 + j] = tup[j + (j >= i)];
1269566063dSJacob Faibussowitsch       PetscCall(PetscNodeRecursive_Internal(dim - 1, degree - tup[i], nodesets, &tup[dim + 1], &node[dim + 1]));
1273f27d899SToby Isaac       for (j = 0; j < dim + 1; j++) node[j + (j >= i)] += wi * node[dim + 1 + j];
1283f27d899SToby Isaac       w += wi;
1293f27d899SToby Isaac     }
1303f27d899SToby Isaac     for (i = 0; i < dim + 1; i++) node[i] /= w;
1313f27d899SToby Isaac   }
1323f27d899SToby Isaac   PetscFunctionReturn(0);
1333f27d899SToby Isaac }
1343f27d899SToby Isaac 
1353f27d899SToby Isaac /* compute simplex nodes for the biunit simplex from the 1D node family */
1369371c9d4SSatish Balay static PetscErrorCode Petsc1DNodeFamilyComputeSimplexNodes(Petsc1DNodeFamily f, PetscInt dim, PetscInt degree, PetscReal points[]) {
1373f27d899SToby Isaac   PetscInt   *tup;
1383f27d899SToby Isaac   PetscInt    k;
1393f27d899SToby Isaac   PetscInt    npoints;
1403f27d899SToby Isaac   PetscReal **nodesets = NULL;
1413f27d899SToby Isaac   PetscInt    worksize;
1423f27d899SToby Isaac   PetscReal  *nodework;
1433f27d899SToby Isaac   PetscInt   *tupwork;
1443f27d899SToby Isaac 
1453f27d899SToby Isaac   PetscFunctionBegin;
14608401ef6SPierre Jolivet   PetscCheck(dim >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must have non-negative dimension");
14708401ef6SPierre Jolivet   PetscCheck(degree >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must have non-negative degree");
1483f27d899SToby Isaac   if (!dim) PetscFunctionReturn(0);
1499566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(dim + 2, &tup));
1503f27d899SToby Isaac   k = 0;
1519566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(degree + dim, dim, &npoints));
1529566063dSJacob Faibussowitsch   PetscCall(Petsc1DNodeFamilyGetNodeSets(f, degree, &nodesets));
1533f27d899SToby Isaac   worksize = ((dim + 2) * (dim + 3)) / 2;
1549566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(worksize, &nodework, worksize, &tupwork));
15577f1a120SToby Isaac   /* loop over the tuples of length dim with sum at most degree */
1563f27d899SToby Isaac   for (k = 0; k < npoints; k++) {
1573f27d899SToby Isaac     PetscInt i;
1583f27d899SToby Isaac 
15977f1a120SToby Isaac     /* turn thm into tuples of length dim + 1 with sum equal to degree (barycentric indice) */
1603f27d899SToby Isaac     tup[0] = degree;
161ad540459SPierre Jolivet     for (i = 0; i < dim; i++) tup[0] -= tup[i + 1];
1623f27d899SToby Isaac     switch (f->nodeFamily) {
1633f27d899SToby Isaac     case PETSCDTNODES_EQUISPACED:
16477f1a120SToby Isaac       /* compute equispaces nodes on the unit reference triangle */
1653f27d899SToby Isaac       if (f->endpoints) {
166ad540459SPierre Jolivet         for (i = 0; i < dim; i++) points[dim * k + i] = (PetscReal)tup[i + 1] / (PetscReal)degree;
1673f27d899SToby Isaac       } else {
1683f27d899SToby Isaac         for (i = 0; i < dim; i++) {
16977f1a120SToby Isaac           /* these nodes are at the centroids of the small simplices created by the equispaced nodes that include
17077f1a120SToby Isaac            * the endpoints */
1713f27d899SToby Isaac           points[dim * k + i] = ((PetscReal)tup[i + 1] + 1. / (dim + 1.)) / (PetscReal)(degree + 1.);
1723f27d899SToby Isaac         }
1733f27d899SToby Isaac       }
1743f27d899SToby Isaac       break;
1753f27d899SToby Isaac     default:
17677f1a120SToby Isaac       /* compute equispaces nodes on the barycentric reference triangle (the trace on the first dim dimensions are the
17777f1a120SToby Isaac        * unit reference triangle nodes */
1783f27d899SToby Isaac       for (i = 0; i < dim + 1; i++) tupwork[i] = tup[i];
1799566063dSJacob Faibussowitsch       PetscCall(PetscNodeRecursive_Internal(dim, degree, nodesets, tupwork, nodework));
1803f27d899SToby Isaac       for (i = 0; i < dim; i++) points[dim * k + i] = nodework[i + 1];
1813f27d899SToby Isaac       break;
1823f27d899SToby Isaac     }
1839566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceLatticePointLexicographic_Internal(dim, degree, &tup[1]));
1843f27d899SToby Isaac   }
1853f27d899SToby Isaac   /* map from unit simplex to biunit simplex */
1863f27d899SToby Isaac   for (k = 0; k < npoints * dim; k++) points[k] = points[k] * 2. - 1.;
1879566063dSJacob Faibussowitsch   PetscCall(PetscFree2(nodework, tupwork));
1889566063dSJacob Faibussowitsch   PetscCall(PetscFree(tup));
1893f27d899SToby Isaac   PetscFunctionReturn(0);
1903f27d899SToby Isaac }
1913f27d899SToby Isaac 
19277f1a120SToby Isaac /* 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
19377f1a120SToby Isaac  * on that mesh point, we have to be careful about getting/adding everything in the right place.
19477f1a120SToby Isaac  *
19577f1a120SToby Isaac  * With nodal dofs like PETSCDUALSPACELAGRANGE makes, the general approach to calculate the value of dofs associate
19677f1a120SToby Isaac  * with a node A is
19777f1a120SToby Isaac  * - transform the node locations x(A) by the map that takes the mesh point to its reorientation, x' = phi(x(A))
19877f1a120SToby Isaac  * - figure out which node was originally at the location of the transformed point, A' = idx(x')
19977f1a120SToby Isaac  * - if the dofs are not scalars, figure out how to represent the transformed dofs in terms of the basis
20077f1a120SToby Isaac  *   of dofs at A' (using pushforward/pullback rules)
20177f1a120SToby Isaac  *
20277f1a120SToby Isaac  * The one sticky point with this approach is the "A' = idx(x')" step: trying to go from real valued coordinates
20377f1a120SToby Isaac  * back to indices.  I don't want to rely on floating point tolerances.  Additionally, PETSCDUALSPACELAGRANGE may
20477f1a120SToby Isaac  * eventually support quasi-Lagrangian dofs, which could involve quadrature at multiple points, so the location "x(A)"
20577f1a120SToby Isaac  * would be ambiguous.
20677f1a120SToby Isaac  *
20777f1a120SToby Isaac  * So each dof gets an integer value coordinate (nodeIdx in the structure below).  The choice of integer coordinates
20877f1a120SToby Isaac  * is somewhat arbitrary, as long as all of the relevant symmetries of the mesh point correspond to *permutations* of
20977f1a120SToby Isaac  * the integer coordinates, which do not depend on numerical precision.
21077f1a120SToby Isaac  *
21177f1a120SToby Isaac  * So
21277f1a120SToby Isaac  *
21377f1a120SToby Isaac  * - DMPlexGetTransitiveClosure_Internal() tells me how an orientation turns into a permutation of the vertices of a
21477f1a120SToby Isaac  *   mesh point
21577f1a120SToby Isaac  * - The permutation of the vertices, and the nodeIdx values assigned to them, tells what permutation in index space
21677f1a120SToby Isaac  *   is associated with the orientation
21777f1a120SToby Isaac  * - I uses that permutation to get xi' = phi(xi(A)), the integer coordinate of the transformed dof
21877f1a120SToby Isaac  * - I can without numerical issues compute A' = idx(xi')
21977f1a120SToby Isaac  *
22077f1a120SToby Isaac  * Here are some examples of how the process works
22177f1a120SToby Isaac  *
22277f1a120SToby Isaac  * - With a triangle:
22377f1a120SToby Isaac  *
22477f1a120SToby Isaac  *   The triangle has the following integer coordinates for vertices, taken from the barycentric triangle
22577f1a120SToby Isaac  *
22677f1a120SToby Isaac  *     closure order 2
22777f1a120SToby Isaac  *     nodeIdx (0,0,1)
22877f1a120SToby Isaac  *      \
22977f1a120SToby Isaac  *       +
23077f1a120SToby Isaac  *       |\
23177f1a120SToby Isaac  *       | \
23277f1a120SToby Isaac  *       |  \
23377f1a120SToby Isaac  *       |   \    closure order 1
23477f1a120SToby Isaac  *       |    \ / nodeIdx (0,1,0)
23577f1a120SToby Isaac  *       +-----+
23677f1a120SToby Isaac  *        \
23777f1a120SToby Isaac  *      closure order 0
23877f1a120SToby Isaac  *      nodeIdx (1,0,0)
23977f1a120SToby Isaac  *
24077f1a120SToby Isaac  *   If I do DMPlexGetTransitiveClosure_Internal() with orientation 1, the vertices would appear
24177f1a120SToby Isaac  *   in the order (1, 2, 0)
24277f1a120SToby Isaac  *
24377f1a120SToby Isaac  *   If I list the nodeIdx of each vertex in closure order for orientation 0 (0, 1, 2) and orientation 1 (1, 2, 0), I
24477f1a120SToby Isaac  *   see
24577f1a120SToby Isaac  *
24677f1a120SToby Isaac  *   orientation 0  | orientation 1
24777f1a120SToby Isaac  *
24877f1a120SToby Isaac  *   [0] (1,0,0)      [1] (0,1,0)
24977f1a120SToby Isaac  *   [1] (0,1,0)      [2] (0,0,1)
25077f1a120SToby Isaac  *   [2] (0,0,1)      [0] (1,0,0)
25177f1a120SToby Isaac  *          A                B
25277f1a120SToby Isaac  *
25377f1a120SToby Isaac  *   In other words, B is the result of a row permutation of A.  But, there is also
25477f1a120SToby Isaac  *   a column permutation that accomplishes the same result, (2,0,1).
25577f1a120SToby Isaac  *
25677f1a120SToby Isaac  *   So if a dof has nodeIdx coordinate (a,b,c), after the transformation its nodeIdx coordinate
25777f1a120SToby Isaac  *   is (c,a,b), and the transformed degree of freedom will be a linear combination of dofs
25877f1a120SToby Isaac  *   that originally had coordinate (c,a,b).
25977f1a120SToby Isaac  *
26077f1a120SToby Isaac  * - With a quadrilateral:
26177f1a120SToby Isaac  *
26277f1a120SToby Isaac  *   The quadrilateral has the following integer coordinates for vertices, taken from concatenating barycentric
26377f1a120SToby Isaac  *   coordinates for two segments:
26477f1a120SToby Isaac  *
26577f1a120SToby Isaac  *     closure order 3      closure order 2
26677f1a120SToby Isaac  *     nodeIdx (1,0,0,1)    nodeIdx (0,1,0,1)
26777f1a120SToby Isaac  *                   \      /
26877f1a120SToby Isaac  *                    +----+
26977f1a120SToby Isaac  *                    |    |
27077f1a120SToby Isaac  *                    |    |
27177f1a120SToby Isaac  *                    +----+
27277f1a120SToby Isaac  *                   /      \
27377f1a120SToby Isaac  *     closure order 0      closure order 1
27477f1a120SToby Isaac  *     nodeIdx (1,0,1,0)    nodeIdx (0,1,1,0)
27577f1a120SToby Isaac  *
27677f1a120SToby Isaac  *   If I do DMPlexGetTransitiveClosure_Internal() with orientation 1, the vertices would appear
27777f1a120SToby Isaac  *   in the order (1, 2, 3, 0)
27877f1a120SToby Isaac  *
27977f1a120SToby Isaac  *   If I list the nodeIdx of each vertex in closure order for orientation 0 (0, 1, 2, 3) and
28077f1a120SToby Isaac  *   orientation 1 (1, 2, 3, 0), I see
28177f1a120SToby Isaac  *
28277f1a120SToby Isaac  *   orientation 0  | orientation 1
28377f1a120SToby Isaac  *
28477f1a120SToby Isaac  *   [0] (1,0,1,0)    [1] (0,1,1,0)
28577f1a120SToby Isaac  *   [1] (0,1,1,0)    [2] (0,1,0,1)
28677f1a120SToby Isaac  *   [2] (0,1,0,1)    [3] (1,0,0,1)
28777f1a120SToby Isaac  *   [3] (1,0,0,1)    [0] (1,0,1,0)
28877f1a120SToby Isaac  *          A                B
28977f1a120SToby Isaac  *
29077f1a120SToby Isaac  *   The column permutation that accomplishes the same result is (3,2,0,1).
29177f1a120SToby Isaac  *
29277f1a120SToby Isaac  *   So if a dof has nodeIdx coordinate (a,b,c,d), after the transformation its nodeIdx coordinate
29377f1a120SToby Isaac  *   is (d,c,a,b), and the transformed degree of freedom will be a linear combination of dofs
29477f1a120SToby Isaac  *   that originally had coordinate (d,c,a,b).
29577f1a120SToby Isaac  *
29677f1a120SToby Isaac  * Previously PETSCDUALSPACELAGRANGE had hardcoded symmetries for the triangle and quadrilateral,
29777f1a120SToby Isaac  * but this approach will work for any polytope, such as the wedge (triangular prism).
29877f1a120SToby Isaac  */
2999371c9d4SSatish Balay struct _n_PetscLagNodeIndices {
3003f27d899SToby Isaac   PetscInt   refct;
3013f27d899SToby Isaac   PetscInt   nodeIdxDim;
3023f27d899SToby Isaac   PetscInt   nodeVecDim;
3033f27d899SToby Isaac   PetscInt   nNodes;
3043f27d899SToby Isaac   PetscInt  *nodeIdx; /* for each node an index of size nodeIdxDim */
3053f27d899SToby Isaac   PetscReal *nodeVec; /* for each node a vector of size nodeVecDim */
3063f27d899SToby Isaac   PetscInt  *perm;    /* if these are vertices, perm takes DMPlex point index to closure order;
3073f27d899SToby Isaac                               if these are nodes, perm lists nodes in index revlex order */
3083f27d899SToby Isaac };
3093f27d899SToby Isaac 
31077f1a120SToby Isaac /* this is just here so I can access the values in tests/ex1.c outside the library */
3119371c9d4SSatish Balay PetscErrorCode PetscLagNodeIndicesGetData_Internal(PetscLagNodeIndices ni, PetscInt *nodeIdxDim, PetscInt *nodeVecDim, PetscInt *nNodes, const PetscInt *nodeIdx[], const PetscReal *nodeVec[]) {
3123f27d899SToby Isaac   PetscFunctionBegin;
3133f27d899SToby Isaac   *nodeIdxDim = ni->nodeIdxDim;
3143f27d899SToby Isaac   *nodeVecDim = ni->nodeVecDim;
3153f27d899SToby Isaac   *nNodes     = ni->nNodes;
3163f27d899SToby Isaac   *nodeIdx    = ni->nodeIdx;
3173f27d899SToby Isaac   *nodeVec    = ni->nodeVec;
3183f27d899SToby Isaac   PetscFunctionReturn(0);
3193f27d899SToby Isaac }
3203f27d899SToby Isaac 
3219371c9d4SSatish Balay static PetscErrorCode PetscLagNodeIndicesReference(PetscLagNodeIndices ni) {
3223f27d899SToby Isaac   PetscFunctionBegin;
3233f27d899SToby Isaac   if (ni) ni->refct++;
3243f27d899SToby Isaac   PetscFunctionReturn(0);
3253f27d899SToby Isaac }
3263f27d899SToby Isaac 
3279371c9d4SSatish Balay static PetscErrorCode PetscLagNodeIndicesDuplicate(PetscLagNodeIndices ni, PetscLagNodeIndices *niNew) {
3281f440fbeSToby Isaac   PetscFunctionBegin;
3299566063dSJacob Faibussowitsch   PetscCall(PetscNew(niNew));
3301f440fbeSToby Isaac   (*niNew)->refct      = 1;
3311f440fbeSToby Isaac   (*niNew)->nodeIdxDim = ni->nodeIdxDim;
3321f440fbeSToby Isaac   (*niNew)->nodeVecDim = ni->nodeVecDim;
3331f440fbeSToby Isaac   (*niNew)->nNodes     = ni->nNodes;
3349566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ni->nNodes * ni->nodeIdxDim, &((*niNew)->nodeIdx)));
3359566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy((*niNew)->nodeIdx, ni->nodeIdx, ni->nNodes * ni->nodeIdxDim));
3369566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ni->nNodes * ni->nodeVecDim, &((*niNew)->nodeVec)));
3379566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy((*niNew)->nodeVec, ni->nodeVec, ni->nNodes * ni->nodeVecDim));
3381f440fbeSToby Isaac   (*niNew)->perm = NULL;
3391f440fbeSToby Isaac   PetscFunctionReturn(0);
3401f440fbeSToby Isaac }
3411f440fbeSToby Isaac 
3429371c9d4SSatish Balay static PetscErrorCode PetscLagNodeIndicesDestroy(PetscLagNodeIndices *ni) {
3433f27d899SToby Isaac   PetscFunctionBegin;
3443f27d899SToby Isaac   if (!(*ni)) PetscFunctionReturn(0);
3453f27d899SToby Isaac   if (--(*ni)->refct > 0) {
3463f27d899SToby Isaac     *ni = NULL;
3473f27d899SToby Isaac     PetscFunctionReturn(0);
3483f27d899SToby Isaac   }
3499566063dSJacob Faibussowitsch   PetscCall(PetscFree((*ni)->nodeIdx));
3509566063dSJacob Faibussowitsch   PetscCall(PetscFree((*ni)->nodeVec));
3519566063dSJacob Faibussowitsch   PetscCall(PetscFree((*ni)->perm));
3529566063dSJacob Faibussowitsch   PetscCall(PetscFree(*ni));
3533f27d899SToby Isaac   *ni = NULL;
3543f27d899SToby Isaac   PetscFunctionReturn(0);
3553f27d899SToby Isaac }
3563f27d899SToby Isaac 
35777f1a120SToby Isaac /* The vertices are given nodeIdx coordinates (e.g. the corners of the barycentric triangle).  Those coordinates are
35877f1a120SToby Isaac  * in some other order, and to understand the effect of different symmetries, we need them to be in closure order.
35977f1a120SToby Isaac  *
36077f1a120SToby Isaac  * If sortIdx is PETSC_FALSE, the coordinates are already in revlex order, otherwise we must sort them
36177f1a120SToby Isaac  * to that order before we do the real work of this function, which is
36277f1a120SToby Isaac  *
36377f1a120SToby Isaac  * - mark the vertices in closure order
36477f1a120SToby Isaac  * - sort them in revlex order
36577f1a120SToby Isaac  * - use the resulting permutation to list the vertex coordinates in closure order
36677f1a120SToby Isaac  */
3679371c9d4SSatish Balay static PetscErrorCode PetscLagNodeIndicesComputeVertexOrder(DM dm, PetscLagNodeIndices ni, PetscBool sortIdx) {
3683f27d899SToby Isaac   PetscInt           v, w, vStart, vEnd, c, d;
3693f27d899SToby Isaac   PetscInt           nVerts;
3703f27d899SToby Isaac   PetscInt           closureSize = 0;
3713f27d899SToby Isaac   PetscInt          *closure     = NULL;
3723f27d899SToby Isaac   PetscInt          *closureOrder;
3733f27d899SToby Isaac   PetscInt          *invClosureOrder;
3743f27d899SToby Isaac   PetscInt          *revlexOrder;
3753f27d899SToby Isaac   PetscInt          *newNodeIdx;
3763f27d899SToby Isaac   PetscInt           dim;
3773f27d899SToby Isaac   Vec                coordVec;
3783f27d899SToby Isaac   const PetscScalar *coords;
3793f27d899SToby Isaac 
3803f27d899SToby Isaac   PetscFunctionBegin;
3819566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
3829566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
3833f27d899SToby Isaac   nVerts = vEnd - vStart;
3849566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nVerts, &closureOrder));
3859566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nVerts, &invClosureOrder));
3869566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nVerts, &revlexOrder));
38777f1a120SToby Isaac   if (sortIdx) { /* bubble sort nodeIdx into revlex order */
3883f27d899SToby Isaac     PetscInt  nodeIdxDim = ni->nodeIdxDim;
3893f27d899SToby Isaac     PetscInt *idxOrder;
3903f27d899SToby Isaac 
3919566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nVerts * nodeIdxDim, &newNodeIdx));
3929566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nVerts, &idxOrder));
3933f27d899SToby Isaac     for (v = 0; v < nVerts; v++) idxOrder[v] = v;
3943f27d899SToby Isaac     for (v = 0; v < nVerts; v++) {
3953f27d899SToby Isaac       for (w = v + 1; w < nVerts; w++) {
3963f27d899SToby Isaac         const PetscInt *iv   = &(ni->nodeIdx[idxOrder[v] * nodeIdxDim]);
3973f27d899SToby Isaac         const PetscInt *iw   = &(ni->nodeIdx[idxOrder[w] * nodeIdxDim]);
3983f27d899SToby Isaac         PetscInt        diff = 0;
3993f27d899SToby Isaac 
4009371c9d4SSatish Balay         for (d = nodeIdxDim - 1; d >= 0; d--)
4019371c9d4SSatish Balay           if ((diff = (iv[d] - iw[d]))) break;
4023f27d899SToby Isaac         if (diff > 0) {
4033f27d899SToby Isaac           PetscInt swap = idxOrder[v];
4043f27d899SToby Isaac 
4053f27d899SToby Isaac           idxOrder[v] = idxOrder[w];
4063f27d899SToby Isaac           idxOrder[w] = swap;
4073f27d899SToby Isaac         }
4083f27d899SToby Isaac       }
4093f27d899SToby Isaac     }
4103f27d899SToby Isaac     for (v = 0; v < nVerts; v++) {
411ad540459SPierre Jolivet       for (d = 0; d < nodeIdxDim; d++) newNodeIdx[v * ni->nodeIdxDim + d] = ni->nodeIdx[idxOrder[v] * nodeIdxDim + d];
4123f27d899SToby Isaac     }
4139566063dSJacob Faibussowitsch     PetscCall(PetscFree(ni->nodeIdx));
4143f27d899SToby Isaac     ni->nodeIdx = newNodeIdx;
4153f27d899SToby Isaac     newNodeIdx  = NULL;
4169566063dSJacob Faibussowitsch     PetscCall(PetscFree(idxOrder));
4173f27d899SToby Isaac   }
4189566063dSJacob Faibussowitsch   PetscCall(DMPlexGetTransitiveClosure(dm, 0, PETSC_TRUE, &closureSize, &closure));
4193f27d899SToby Isaac   c = closureSize - nVerts;
4203f27d899SToby Isaac   for (v = 0; v < nVerts; v++) closureOrder[v] = closure[2 * (c + v)] - vStart;
4213f27d899SToby Isaac   for (v = 0; v < nVerts; v++) invClosureOrder[closureOrder[v]] = v;
4229566063dSJacob Faibussowitsch   PetscCall(DMPlexRestoreTransitiveClosure(dm, 0, PETSC_TRUE, &closureSize, &closure));
4239566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dm, &coordVec));
4249566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(coordVec, &coords));
4253f27d899SToby Isaac   /* bubble sort closure vertices by coordinates in revlex order */
4263f27d899SToby Isaac   for (v = 0; v < nVerts; v++) revlexOrder[v] = v;
4273f27d899SToby Isaac   for (v = 0; v < nVerts; v++) {
4283f27d899SToby Isaac     for (w = v + 1; w < nVerts; w++) {
4293f27d899SToby Isaac       const PetscScalar *cv   = &coords[closureOrder[revlexOrder[v]] * dim];
4303f27d899SToby Isaac       const PetscScalar *cw   = &coords[closureOrder[revlexOrder[w]] * dim];
4313f27d899SToby Isaac       PetscReal          diff = 0;
4323f27d899SToby Isaac 
4339371c9d4SSatish Balay       for (d = dim - 1; d >= 0; d--)
4349371c9d4SSatish Balay         if ((diff = PetscRealPart(cv[d] - cw[d])) != 0.) break;
4353f27d899SToby Isaac       if (diff > 0.) {
4363f27d899SToby Isaac         PetscInt swap = revlexOrder[v];
4373f27d899SToby Isaac 
4383f27d899SToby Isaac         revlexOrder[v] = revlexOrder[w];
4393f27d899SToby Isaac         revlexOrder[w] = swap;
4403f27d899SToby Isaac       }
4413f27d899SToby Isaac     }
4423f27d899SToby Isaac   }
4439566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(coordVec, &coords));
4449566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ni->nodeIdxDim * nVerts, &newNodeIdx));
4453f27d899SToby Isaac   /* reorder nodeIdx to be in closure order */
4463f27d899SToby Isaac   for (v = 0; v < nVerts; v++) {
447ad540459SPierre Jolivet     for (d = 0; d < ni->nodeIdxDim; d++) newNodeIdx[revlexOrder[v] * ni->nodeIdxDim + d] = ni->nodeIdx[v * ni->nodeIdxDim + d];
4483f27d899SToby Isaac   }
4499566063dSJacob Faibussowitsch   PetscCall(PetscFree(ni->nodeIdx));
4503f27d899SToby Isaac   ni->nodeIdx = newNodeIdx;
4513f27d899SToby Isaac   ni->perm    = invClosureOrder;
4529566063dSJacob Faibussowitsch   PetscCall(PetscFree(revlexOrder));
4539566063dSJacob Faibussowitsch   PetscCall(PetscFree(closureOrder));
4543f27d899SToby Isaac   PetscFunctionReturn(0);
4553f27d899SToby Isaac }
4563f27d899SToby Isaac 
45777f1a120SToby Isaac /* the coordinates of the simplex vertices are the corners of the barycentric simplex.
45877f1a120SToby Isaac  * When we stack them on top of each other in revlex order, they look like the identity matrix */
4599371c9d4SSatish Balay static PetscErrorCode PetscLagNodeIndicesCreateSimplexVertices(DM dm, PetscLagNodeIndices *nodeIndices) {
4603f27d899SToby Isaac   PetscLagNodeIndices ni;
4613f27d899SToby Isaac   PetscInt            dim, d;
4623f27d899SToby Isaac 
4633f27d899SToby Isaac   PetscFunctionBegin;
4649566063dSJacob Faibussowitsch   PetscCall(PetscNew(&ni));
4659566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
4663f27d899SToby Isaac   ni->nodeIdxDim = dim + 1;
4673f27d899SToby Isaac   ni->nodeVecDim = 0;
4683f27d899SToby Isaac   ni->nNodes     = dim + 1;
4693f27d899SToby Isaac   ni->refct      = 1;
4709566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1((dim + 1) * (dim + 1), &(ni->nodeIdx)));
4713f27d899SToby Isaac   for (d = 0; d < dim + 1; d++) ni->nodeIdx[d * (dim + 2)] = 1;
4729566063dSJacob Faibussowitsch   PetscCall(PetscLagNodeIndicesComputeVertexOrder(dm, ni, PETSC_FALSE));
4733f27d899SToby Isaac   *nodeIndices = ni;
4743f27d899SToby Isaac   PetscFunctionReturn(0);
4753f27d899SToby Isaac }
4763f27d899SToby Isaac 
47777f1a120SToby Isaac /* A polytope that is a tensor product of a facet and a segment.
47877f1a120SToby Isaac  * We take whatever coordinate system was being used for the facet
4791f440fbeSToby Isaac  * and we concatenate the barycentric coordinates for the vertices
48077f1a120SToby Isaac  * at the end of the segment, (1,0) and (0,1), to get a coordinate
48177f1a120SToby Isaac  * system for the tensor product element */
4829371c9d4SSatish Balay static PetscErrorCode PetscLagNodeIndicesCreateTensorVertices(DM dm, PetscLagNodeIndices facetni, PetscLagNodeIndices *nodeIndices) {
4833f27d899SToby Isaac   PetscLagNodeIndices ni;
4843f27d899SToby Isaac   PetscInt            nodeIdxDim, subNodeIdxDim = facetni->nodeIdxDim;
4853f27d899SToby Isaac   PetscInt            nVerts, nSubVerts         = facetni->nNodes;
4863f27d899SToby Isaac   PetscInt            dim, d, e, f, g;
4873f27d899SToby Isaac 
4883f27d899SToby Isaac   PetscFunctionBegin;
4899566063dSJacob Faibussowitsch   PetscCall(PetscNew(&ni));
4909566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
4913f27d899SToby Isaac   ni->nodeIdxDim = nodeIdxDim = subNodeIdxDim + 2;
4923f27d899SToby Isaac   ni->nodeVecDim              = 0;
4933f27d899SToby Isaac   ni->nNodes = nVerts = 2 * nSubVerts;
4943f27d899SToby Isaac   ni->refct           = 1;
4959566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(nodeIdxDim * nVerts, &(ni->nodeIdx)));
4963f27d899SToby Isaac   for (f = 0, d = 0; d < 2; d++) {
4973f27d899SToby Isaac     for (e = 0; e < nSubVerts; e++, f++) {
498ad540459SPierre Jolivet       for (g = 0; g < subNodeIdxDim; g++) ni->nodeIdx[f * nodeIdxDim + g] = facetni->nodeIdx[e * subNodeIdxDim + g];
4993f27d899SToby Isaac       ni->nodeIdx[f * nodeIdxDim + subNodeIdxDim]     = (1 - d);
5003f27d899SToby Isaac       ni->nodeIdx[f * nodeIdxDim + subNodeIdxDim + 1] = d;
5013f27d899SToby Isaac     }
5023f27d899SToby Isaac   }
5039566063dSJacob Faibussowitsch   PetscCall(PetscLagNodeIndicesComputeVertexOrder(dm, ni, PETSC_TRUE));
5043f27d899SToby Isaac   *nodeIndices = ni;
5053f27d899SToby Isaac   PetscFunctionReturn(0);
5063f27d899SToby Isaac }
5073f27d899SToby Isaac 
50877f1a120SToby Isaac /* This helps us compute symmetries, and it also helps us compute coordinates for dofs that are being pushed
50977f1a120SToby Isaac  * forward from a boundary mesh point.
51077f1a120SToby Isaac  *
51177f1a120SToby Isaac  * Input:
51277f1a120SToby Isaac  *
51377f1a120SToby Isaac  * dm - the target reference cell where we want new coordinates and dof directions to be valid
51477f1a120SToby Isaac  * vert - the vertex coordinate system for the target reference cell
51577f1a120SToby Isaac  * p - the point in the target reference cell that the dofs are coming from
51677f1a120SToby Isaac  * vertp - the vertex coordinate system for p's reference cell
51777f1a120SToby Isaac  * ornt - the resulting coordinates and dof vectors will be for p under this orientation
51877f1a120SToby Isaac  * nodep - the node coordinates and dof vectors in p's reference cell
51977f1a120SToby Isaac  * formDegree - the form degree that the dofs transform as
52077f1a120SToby Isaac  *
52177f1a120SToby Isaac  * Output:
52277f1a120SToby Isaac  *
52377f1a120SToby Isaac  * pfNodeIdx - the node coordinates for p's dofs, in the dm reference cell, from the ornt perspective
52477f1a120SToby Isaac  * pfNodeVec - the node dof vectors for p's dofs, in the dm reference cell, from the ornt perspective
52577f1a120SToby Isaac  */
5269371c9d4SSatish Balay static PetscErrorCode PetscLagNodeIndicesPushForward(DM dm, PetscLagNodeIndices vert, PetscInt p, PetscLagNodeIndices vertp, PetscLagNodeIndices nodep, PetscInt ornt, PetscInt formDegree, PetscInt pfNodeIdx[], PetscReal pfNodeVec[]) {
5273f27d899SToby Isaac   PetscInt          *closureVerts;
5283f27d899SToby Isaac   PetscInt           closureSize = 0;
5293f27d899SToby Isaac   PetscInt          *closure     = NULL;
5303f27d899SToby Isaac   PetscInt           dim, pdim, c, i, j, k, n, v, vStart, vEnd;
5313f27d899SToby Isaac   PetscInt           nSubVert      = vertp->nNodes;
5323f27d899SToby Isaac   PetscInt           nodeIdxDim    = vert->nodeIdxDim;
5333f27d899SToby Isaac   PetscInt           subNodeIdxDim = vertp->nodeIdxDim;
5343f27d899SToby Isaac   PetscInt           nNodes        = nodep->nNodes;
5353f27d899SToby Isaac   const PetscInt    *vertIdx       = vert->nodeIdx;
5363f27d899SToby Isaac   const PetscInt    *subVertIdx    = vertp->nodeIdx;
5373f27d899SToby Isaac   const PetscInt    *nodeIdx       = nodep->nodeIdx;
5383f27d899SToby Isaac   const PetscReal   *nodeVec       = nodep->nodeVec;
5393f27d899SToby Isaac   PetscReal         *J, *Jstar;
5403f27d899SToby Isaac   PetscReal          detJ;
5413f27d899SToby Isaac   PetscInt           depth, pdepth, Nk, pNk;
5423f27d899SToby Isaac   Vec                coordVec;
5433f27d899SToby Isaac   PetscScalar       *newCoords = NULL;
5443f27d899SToby Isaac   const PetscScalar *oldCoords = NULL;
5453f27d899SToby Isaac 
5463f27d899SToby Isaac   PetscFunctionBegin;
5479566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
5489566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(dm, &depth));
5499566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dm, &coordVec));
5509566063dSJacob Faibussowitsch   PetscCall(DMPlexGetPointDepth(dm, p, &pdepth));
5513f27d899SToby Isaac   pdim = pdepth != depth ? pdepth != 0 ? pdepth : 0 : dim;
5529566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
5539566063dSJacob Faibussowitsch   PetscCall(DMGetWorkArray(dm, nSubVert, MPIU_INT, &closureVerts));
5549566063dSJacob Faibussowitsch   PetscCall(DMPlexGetTransitiveClosure_Internal(dm, p, ornt, PETSC_TRUE, &closureSize, &closure));
5553f27d899SToby Isaac   c = closureSize - nSubVert;
5563f27d899SToby Isaac   /* we want which cell closure indices the closure of this point corresponds to */
5573f27d899SToby Isaac   for (v = 0; v < nSubVert; v++) closureVerts[v] = vert->perm[closure[2 * (c + v)] - vStart];
5589566063dSJacob Faibussowitsch   PetscCall(DMPlexRestoreTransitiveClosure(dm, p, PETSC_TRUE, &closureSize, &closure));
5593f27d899SToby Isaac   /* push forward indices */
5603f27d899SToby Isaac   for (i = 0; i < nodeIdxDim; i++) { /* for every component of the target index space */
5613f27d899SToby Isaac     /* check if this is a component that all vertices around this point have in common */
5623f27d899SToby Isaac     for (j = 1; j < nSubVert; j++) {
5633f27d899SToby Isaac       if (vertIdx[closureVerts[j] * nodeIdxDim + i] != vertIdx[closureVerts[0] * nodeIdxDim + i]) break;
5643f27d899SToby Isaac     }
5653f27d899SToby Isaac     if (j == nSubVert) { /* all vertices have this component in common, directly copy to output */
5663f27d899SToby Isaac       PetscInt val = vertIdx[closureVerts[0] * nodeIdxDim + i];
5673f27d899SToby Isaac       for (n = 0; n < nNodes; n++) pfNodeIdx[n * nodeIdxDim + i] = val;
5683f27d899SToby Isaac     } else {
5693f27d899SToby Isaac       PetscInt subi = -1;
5703f27d899SToby Isaac       /* there must be a component in vertp that looks the same */
5713f27d899SToby Isaac       for (k = 0; k < subNodeIdxDim; k++) {
5723f27d899SToby Isaac         for (j = 0; j < nSubVert; j++) {
5733f27d899SToby Isaac           if (vertIdx[closureVerts[j] * nodeIdxDim + i] != subVertIdx[j * subNodeIdxDim + k]) break;
5743f27d899SToby Isaac         }
5753f27d899SToby Isaac         if (j == nSubVert) {
5763f27d899SToby Isaac           subi = k;
5773f27d899SToby Isaac           break;
5783f27d899SToby Isaac         }
5793f27d899SToby Isaac       }
58008401ef6SPierre Jolivet       PetscCheck(subi >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Did not find matching coordinate");
58177f1a120SToby Isaac       /* that component in the vertp system becomes component i in the vert system for each dof */
5823f27d899SToby Isaac       for (n = 0; n < nNodes; n++) pfNodeIdx[n * nodeIdxDim + i] = nodeIdx[n * subNodeIdxDim + subi];
5833f27d899SToby Isaac     }
5843f27d899SToby Isaac   }
5853f27d899SToby Isaac   /* push forward vectors */
5869566063dSJacob Faibussowitsch   PetscCall(DMGetWorkArray(dm, dim * dim, MPIU_REAL, &J));
58777f1a120SToby Isaac   if (ornt != 0) { /* temporarily change the coordinate vector so
58877f1a120SToby Isaac                       DMPlexComputeCellGeometryAffineFEM gives us the Jacobian we want */
5893f27d899SToby Isaac     PetscInt  closureSize2 = 0;
5903f27d899SToby Isaac     PetscInt *closure2     = NULL;
5913f27d899SToby Isaac 
5929566063dSJacob Faibussowitsch     PetscCall(DMPlexGetTransitiveClosure_Internal(dm, p, 0, PETSC_TRUE, &closureSize2, &closure2));
5939566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dim * nSubVert, &newCoords));
5949566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(coordVec, &oldCoords));
5953f27d899SToby Isaac     for (v = 0; v < nSubVert; v++) {
5963f27d899SToby Isaac       PetscInt d;
597ad540459SPierre Jolivet       for (d = 0; d < dim; d++) newCoords[(closure2[2 * (c + v)] - vStart) * dim + d] = oldCoords[closureVerts[v] * dim + d];
5983f27d899SToby Isaac     }
5999566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(coordVec, &oldCoords));
6009566063dSJacob Faibussowitsch     PetscCall(DMPlexRestoreTransitiveClosure(dm, p, PETSC_TRUE, &closureSize2, &closure2));
6019566063dSJacob Faibussowitsch     PetscCall(VecPlaceArray(coordVec, newCoords));
6023f27d899SToby Isaac   }
6039566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeCellGeometryAffineFEM(dm, p, NULL, J, NULL, &detJ));
6043f27d899SToby Isaac   if (ornt != 0) {
6059566063dSJacob Faibussowitsch     PetscCall(VecResetArray(coordVec));
6069566063dSJacob Faibussowitsch     PetscCall(PetscFree(newCoords));
6073f27d899SToby Isaac   }
6089566063dSJacob Faibussowitsch   PetscCall(DMRestoreWorkArray(dm, nSubVert, MPIU_INT, &closureVerts));
6093f27d899SToby Isaac   /* compactify */
6109371c9d4SSatish Balay   for (i = 0; i < dim; i++)
6119371c9d4SSatish Balay     for (j = 0; j < pdim; j++) J[i * pdim + j] = J[i * dim + j];
61277f1a120SToby Isaac   /* We have the Jacobian mapping the point's reference cell to this reference cell:
61377f1a120SToby Isaac    * pulling back a function to the point and applying the dof is what we want,
61477f1a120SToby Isaac    * so we get the pullback matrix and multiply the dof by that matrix on the right */
6159566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(formDegree), &Nk));
6169566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(pdim, PetscAbsInt(formDegree), &pNk));
6179566063dSJacob Faibussowitsch   PetscCall(DMGetWorkArray(dm, pNk * Nk, MPIU_REAL, &Jstar));
6189566063dSJacob Faibussowitsch   PetscCall(PetscDTAltVPullbackMatrix(pdim, dim, J, formDegree, Jstar));
6193f27d899SToby Isaac   for (n = 0; n < nNodes; n++) {
6203f27d899SToby Isaac     for (i = 0; i < Nk; i++) {
6213f27d899SToby Isaac       PetscReal val = 0.;
6225efe5503SToby Isaac       for (j = 0; j < pNk; j++) val += nodeVec[n * pNk + j] * Jstar[j * Nk + i];
6233f27d899SToby Isaac       pfNodeVec[n * Nk + i] = val;
6243f27d899SToby Isaac     }
6253f27d899SToby Isaac   }
6269566063dSJacob Faibussowitsch   PetscCall(DMRestoreWorkArray(dm, pNk * Nk, MPIU_REAL, &Jstar));
6279566063dSJacob Faibussowitsch   PetscCall(DMRestoreWorkArray(dm, dim * dim, MPIU_REAL, &J));
6283f27d899SToby Isaac   PetscFunctionReturn(0);
6293f27d899SToby Isaac }
6303f27d899SToby Isaac 
63177f1a120SToby Isaac /* given to sets of nodes, take the tensor product, where the product of the dof indices is concatenation and the
63277f1a120SToby Isaac  * product of the dof vectors is the wedge product */
6339371c9d4SSatish Balay static PetscErrorCode PetscLagNodeIndicesTensor(PetscLagNodeIndices tracei, PetscInt dimT, PetscInt kT, PetscLagNodeIndices fiberi, PetscInt dimF, PetscInt kF, PetscLagNodeIndices *nodeIndices) {
6343f27d899SToby Isaac   PetscInt            dim = dimT + dimF;
6353f27d899SToby Isaac   PetscInt            nodeIdxDim, nNodes;
6363f27d899SToby Isaac   PetscInt            formDegree = kT + kF;
6373f27d899SToby Isaac   PetscInt            Nk, NkT, NkF;
6383f27d899SToby Isaac   PetscInt            MkT, MkF;
6393f27d899SToby Isaac   PetscLagNodeIndices ni;
6403f27d899SToby Isaac   PetscInt            i, j, l;
6413f27d899SToby Isaac   PetscReal          *projF, *projT;
6423f27d899SToby Isaac   PetscReal          *projFstar, *projTstar;
6433f27d899SToby Isaac   PetscReal          *workF, *workF2, *workT, *workT2, *work, *work2;
6443f27d899SToby Isaac   PetscReal          *wedgeMat;
6453f27d899SToby Isaac   PetscReal           sign;
6463f27d899SToby Isaac 
6473f27d899SToby Isaac   PetscFunctionBegin;
6489566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(formDegree), &Nk));
6499566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(dimT, PetscAbsInt(kT), &NkT));
6509566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(dimF, PetscAbsInt(kF), &NkF));
6519566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(kT), &MkT));
6529566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(kF), &MkF));
6539566063dSJacob Faibussowitsch   PetscCall(PetscNew(&ni));
6543f27d899SToby Isaac   ni->nodeIdxDim = nodeIdxDim = tracei->nodeIdxDim + fiberi->nodeIdxDim;
6553f27d899SToby Isaac   ni->nodeVecDim              = Nk;
6563f27d899SToby Isaac   ni->nNodes = nNodes = tracei->nNodes * fiberi->nNodes;
6573f27d899SToby Isaac   ni->refct           = 1;
6589566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nNodes * nodeIdxDim, &(ni->nodeIdx)));
6593f27d899SToby Isaac   /* first concatenate the indices */
6603f27d899SToby Isaac   for (l = 0, j = 0; j < fiberi->nNodes; j++) {
6613f27d899SToby Isaac     for (i = 0; i < tracei->nNodes; i++, l++) {
6623f27d899SToby Isaac       PetscInt m, n = 0;
6633f27d899SToby Isaac 
6643f27d899SToby Isaac       for (m = 0; m < tracei->nodeIdxDim; m++) ni->nodeIdx[l * nodeIdxDim + n++] = tracei->nodeIdx[i * tracei->nodeIdxDim + m];
6653f27d899SToby Isaac       for (m = 0; m < fiberi->nodeIdxDim; m++) ni->nodeIdx[l * nodeIdxDim + n++] = fiberi->nodeIdx[j * fiberi->nodeIdxDim + m];
6663f27d899SToby Isaac     }
6673f27d899SToby Isaac   }
6683f27d899SToby Isaac 
6693f27d899SToby Isaac   /* now wedge together the push-forward vectors */
6709566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nNodes * Nk, &(ni->nodeVec)));
6719566063dSJacob Faibussowitsch   PetscCall(PetscCalloc2(dimT * dim, &projT, dimF * dim, &projF));
6723f27d899SToby Isaac   for (i = 0; i < dimT; i++) projT[i * (dim + 1)] = 1.;
6733f27d899SToby Isaac   for (i = 0; i < dimF; i++) projF[i * (dim + dimT + 1) + dimT] = 1.;
6749566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(MkT * NkT, &projTstar, MkF * NkF, &projFstar));
6759566063dSJacob Faibussowitsch   PetscCall(PetscDTAltVPullbackMatrix(dim, dimT, projT, kT, projTstar));
6769566063dSJacob Faibussowitsch   PetscCall(PetscDTAltVPullbackMatrix(dim, dimF, projF, kF, projFstar));
6779566063dSJacob Faibussowitsch   PetscCall(PetscMalloc6(MkT, &workT, MkT, &workT2, MkF, &workF, MkF, &workF2, Nk, &work, Nk, &work2));
6789566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Nk * MkT, &wedgeMat));
6793f27d899SToby Isaac   sign = (PetscAbsInt(kT * kF) & 1) ? -1. : 1.;
6803f27d899SToby Isaac   for (l = 0, j = 0; j < fiberi->nNodes; j++) {
6813f27d899SToby Isaac     PetscInt d, e;
6823f27d899SToby Isaac 
6833f27d899SToby Isaac     /* push forward fiber k-form */
6843f27d899SToby Isaac     for (d = 0; d < MkF; d++) {
6853f27d899SToby Isaac       PetscReal val = 0.;
6863f27d899SToby Isaac       for (e = 0; e < NkF; e++) val += projFstar[d * NkF + e] * fiberi->nodeVec[j * NkF + e];
6873f27d899SToby Isaac       workF[d] = val;
6883f27d899SToby Isaac     }
6893f27d899SToby Isaac     /* Hodge star to proper form if necessary */
6903f27d899SToby Isaac     if (kF < 0) {
6913f27d899SToby Isaac       for (d = 0; d < MkF; d++) workF2[d] = workF[d];
6929566063dSJacob Faibussowitsch       PetscCall(PetscDTAltVStar(dim, PetscAbsInt(kF), 1, workF2, workF));
6933f27d899SToby Isaac     }
6943f27d899SToby Isaac     /* Compute the matrix that wedges this form with one of the trace k-form */
6959566063dSJacob Faibussowitsch     PetscCall(PetscDTAltVWedgeMatrix(dim, PetscAbsInt(kF), PetscAbsInt(kT), workF, wedgeMat));
6963f27d899SToby Isaac     for (i = 0; i < tracei->nNodes; i++, l++) {
6973f27d899SToby Isaac       /* push forward trace k-form */
6983f27d899SToby Isaac       for (d = 0; d < MkT; d++) {
6993f27d899SToby Isaac         PetscReal val = 0.;
7003f27d899SToby Isaac         for (e = 0; e < NkT; e++) val += projTstar[d * NkT + e] * tracei->nodeVec[i * NkT + e];
7013f27d899SToby Isaac         workT[d] = val;
7023f27d899SToby Isaac       }
7033f27d899SToby Isaac       /* Hodge star to proper form if necessary */
7043f27d899SToby Isaac       if (kT < 0) {
7053f27d899SToby Isaac         for (d = 0; d < MkT; d++) workT2[d] = workT[d];
7069566063dSJacob Faibussowitsch         PetscCall(PetscDTAltVStar(dim, PetscAbsInt(kT), 1, workT2, workT));
7073f27d899SToby Isaac       }
7083f27d899SToby Isaac       /* compute the wedge product of the push-forward trace form and firer forms */
7093f27d899SToby Isaac       for (d = 0; d < Nk; d++) {
7103f27d899SToby Isaac         PetscReal val = 0.;
7113f27d899SToby Isaac         for (e = 0; e < MkT; e++) val += wedgeMat[d * MkT + e] * workT[e];
7123f27d899SToby Isaac         work[d] = val;
7133f27d899SToby Isaac       }
7143f27d899SToby Isaac       /* inverse Hodge star from proper form if necessary */
7153f27d899SToby Isaac       if (formDegree < 0) {
7163f27d899SToby Isaac         for (d = 0; d < Nk; d++) work2[d] = work[d];
7179566063dSJacob Faibussowitsch         PetscCall(PetscDTAltVStar(dim, PetscAbsInt(formDegree), -1, work2, work));
7183f27d899SToby Isaac       }
7193f27d899SToby Isaac       /* insert into the array (adjusting for sign) */
7203f27d899SToby Isaac       for (d = 0; d < Nk; d++) ni->nodeVec[l * Nk + d] = sign * work[d];
7213f27d899SToby Isaac     }
7223f27d899SToby Isaac   }
7239566063dSJacob Faibussowitsch   PetscCall(PetscFree(wedgeMat));
7249566063dSJacob Faibussowitsch   PetscCall(PetscFree6(workT, workT2, workF, workF2, work, work2));
7259566063dSJacob Faibussowitsch   PetscCall(PetscFree2(projTstar, projFstar));
7269566063dSJacob Faibussowitsch   PetscCall(PetscFree2(projT, projF));
7273f27d899SToby Isaac   *nodeIndices = ni;
7283f27d899SToby Isaac   PetscFunctionReturn(0);
7293f27d899SToby Isaac }
7303f27d899SToby Isaac 
73177f1a120SToby Isaac /* simple union of two sets of nodes */
7329371c9d4SSatish Balay static PetscErrorCode PetscLagNodeIndicesMerge(PetscLagNodeIndices niA, PetscLagNodeIndices niB, PetscLagNodeIndices *nodeIndices) {
7333f27d899SToby Isaac   PetscLagNodeIndices ni;
7343f27d899SToby Isaac   PetscInt            nodeIdxDim, nodeVecDim, nNodes;
7353f27d899SToby Isaac 
7363f27d899SToby Isaac   PetscFunctionBegin;
7379566063dSJacob Faibussowitsch   PetscCall(PetscNew(&ni));
7383f27d899SToby Isaac   ni->nodeIdxDim = nodeIdxDim = niA->nodeIdxDim;
73908401ef6SPierre Jolivet   PetscCheck(niB->nodeIdxDim == nodeIdxDim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Cannot merge PetscLagNodeIndices with different nodeIdxDim");
7403f27d899SToby Isaac   ni->nodeVecDim = nodeVecDim = niA->nodeVecDim;
74108401ef6SPierre Jolivet   PetscCheck(niB->nodeVecDim == nodeVecDim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Cannot merge PetscLagNodeIndices with different nodeVecDim");
7423f27d899SToby Isaac   ni->nNodes = nNodes = niA->nNodes + niB->nNodes;
7433f27d899SToby Isaac   ni->refct           = 1;
7449566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nNodes * nodeIdxDim, &(ni->nodeIdx)));
7459566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nNodes * nodeVecDim, &(ni->nodeVec)));
7469566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(ni->nodeIdx, niA->nodeIdx, niA->nNodes * nodeIdxDim));
7479566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(ni->nodeVec, niA->nodeVec, niA->nNodes * nodeVecDim));
7489566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(&(ni->nodeIdx[niA->nNodes * nodeIdxDim]), niB->nodeIdx, niB->nNodes * nodeIdxDim));
7499566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(&(ni->nodeVec[niA->nNodes * nodeVecDim]), niB->nodeVec, niB->nNodes * nodeVecDim));
7503f27d899SToby Isaac   *nodeIndices = ni;
7513f27d899SToby Isaac   PetscFunctionReturn(0);
7523f27d899SToby Isaac }
7533f27d899SToby Isaac 
7543f27d899SToby Isaac #define PETSCTUPINTCOMPREVLEX(N) \
7559371c9d4SSatish Balay   static int PetscConcat_(PetscTupIntCompRevlex_, N)(const void *a, const void *b) { \
7563f27d899SToby Isaac     const PetscInt *A = (const PetscInt *)a; \
7573f27d899SToby Isaac     const PetscInt *B = (const PetscInt *)b; \
7583f27d899SToby Isaac     int             i; \
7593f27d899SToby Isaac     PetscInt        diff = 0; \
7603f27d899SToby Isaac     for (i = 0; i < N; i++) { \
7613f27d899SToby Isaac       diff = A[N - i] - B[N - i]; \
7623f27d899SToby Isaac       if (diff) break; \
7633f27d899SToby Isaac     } \
7643f27d899SToby Isaac     return (diff <= 0) ? (diff < 0) ? -1 : 0 : 1; \
7653f27d899SToby Isaac   }
7663f27d899SToby Isaac 
7673f27d899SToby Isaac PETSCTUPINTCOMPREVLEX(3)
7683f27d899SToby Isaac PETSCTUPINTCOMPREVLEX(4)
7693f27d899SToby Isaac PETSCTUPINTCOMPREVLEX(5)
7703f27d899SToby Isaac PETSCTUPINTCOMPREVLEX(6)
7713f27d899SToby Isaac PETSCTUPINTCOMPREVLEX(7)
7723f27d899SToby Isaac 
7739371c9d4SSatish Balay static int PetscTupIntCompRevlex_N(const void *a, const void *b) {
7743f27d899SToby Isaac   const PetscInt *A = (const PetscInt *)a;
7753f27d899SToby Isaac   const PetscInt *B = (const PetscInt *)b;
7763f27d899SToby Isaac   int             i;
7773f27d899SToby Isaac   int             N    = A[0];
7783f27d899SToby Isaac   PetscInt        diff = 0;
7793f27d899SToby Isaac   for (i = 0; i < N; i++) {
7803f27d899SToby Isaac     diff = A[N - i] - B[N - i];
7813f27d899SToby Isaac     if (diff) break;
7823f27d899SToby Isaac   }
7833f27d899SToby Isaac   return (diff <= 0) ? (diff < 0) ? -1 : 0 : 1;
7843f27d899SToby Isaac }
7853f27d899SToby Isaac 
78677f1a120SToby Isaac /* The nodes are not necessarily in revlex order wrt nodeIdx: get the permutation
78777f1a120SToby Isaac  * that puts them in that order */
7889371c9d4SSatish Balay static PetscErrorCode PetscLagNodeIndicesGetPermutation(PetscLagNodeIndices ni, PetscInt *perm[]) {
7893f27d899SToby Isaac   PetscFunctionBegin;
7903f27d899SToby Isaac   if (!(ni->perm)) {
7913f27d899SToby Isaac     PetscInt *sorter;
7923f27d899SToby Isaac     PetscInt  m          = ni->nNodes;
7933f27d899SToby Isaac     PetscInt  nodeIdxDim = ni->nodeIdxDim;
7943f27d899SToby Isaac     PetscInt  i, j, k, l;
7953f27d899SToby Isaac     PetscInt *prm;
7963f27d899SToby Isaac     int (*comp)(const void *, const void *);
7973f27d899SToby Isaac 
7989566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1((nodeIdxDim + 2) * m, &sorter));
7993f27d899SToby Isaac     for (k = 0, l = 0, i = 0; i < m; i++) {
8003f27d899SToby Isaac       sorter[k++] = nodeIdxDim + 1;
8013f27d899SToby Isaac       sorter[k++] = i;
8023f27d899SToby Isaac       for (j = 0; j < nodeIdxDim; j++) sorter[k++] = ni->nodeIdx[l++];
8033f27d899SToby Isaac     }
8043f27d899SToby Isaac     switch (nodeIdxDim) {
8059371c9d4SSatish Balay     case 2: comp = PetscTupIntCompRevlex_3; break;
8069371c9d4SSatish Balay     case 3: comp = PetscTupIntCompRevlex_4; break;
8079371c9d4SSatish Balay     case 4: comp = PetscTupIntCompRevlex_5; break;
8089371c9d4SSatish Balay     case 5: comp = PetscTupIntCompRevlex_6; break;
8099371c9d4SSatish Balay     case 6: comp = PetscTupIntCompRevlex_7; break;
8109371c9d4SSatish Balay     default: comp = PetscTupIntCompRevlex_N; break;
8113f27d899SToby Isaac     }
8123f27d899SToby Isaac     qsort(sorter, m, (nodeIdxDim + 2) * sizeof(PetscInt), comp);
8139566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(m, &prm));
8143f27d899SToby Isaac     for (i = 0; i < m; i++) prm[i] = sorter[(nodeIdxDim + 2) * i + 1];
8153f27d899SToby Isaac     ni->perm = prm;
8169566063dSJacob Faibussowitsch     PetscCall(PetscFree(sorter));
8173f27d899SToby Isaac   }
8183f27d899SToby Isaac   *perm = ni->perm;
8193f27d899SToby Isaac   PetscFunctionReturn(0);
8203f27d899SToby Isaac }
82120cf1dd8SToby Isaac 
8229371c9d4SSatish Balay static PetscErrorCode PetscDualSpaceDestroy_Lagrange(PetscDualSpace sp) {
82320cf1dd8SToby Isaac   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data;
82420cf1dd8SToby Isaac 
82520cf1dd8SToby Isaac   PetscFunctionBegin;
8263f27d899SToby Isaac   if (lag->symperms) {
8273f27d899SToby Isaac     PetscInt **selfSyms = lag->symperms[0];
8286f905325SMatthew G. Knepley 
8296f905325SMatthew G. Knepley     if (selfSyms) {
8306f905325SMatthew G. Knepley       PetscInt i, **allocated = &selfSyms[-lag->selfSymOff];
8316f905325SMatthew G. Knepley 
83248a46eb9SPierre Jolivet       for (i = 0; i < lag->numSelfSym; i++) PetscCall(PetscFree(allocated[i]));
8339566063dSJacob Faibussowitsch       PetscCall(PetscFree(allocated));
8346f905325SMatthew G. Knepley     }
8359566063dSJacob Faibussowitsch     PetscCall(PetscFree(lag->symperms));
8366f905325SMatthew G. Knepley   }
8373f27d899SToby Isaac   if (lag->symflips) {
8383f27d899SToby Isaac     PetscScalar **selfSyms = lag->symflips[0];
8393f27d899SToby Isaac 
8403f27d899SToby Isaac     if (selfSyms) {
8413f27d899SToby Isaac       PetscInt      i;
8423f27d899SToby Isaac       PetscScalar **allocated = &selfSyms[-lag->selfSymOff];
8433f27d899SToby Isaac 
84448a46eb9SPierre Jolivet       for (i = 0; i < lag->numSelfSym; i++) PetscCall(PetscFree(allocated[i]));
8459566063dSJacob Faibussowitsch       PetscCall(PetscFree(allocated));
8463f27d899SToby Isaac     }
8479566063dSJacob Faibussowitsch     PetscCall(PetscFree(lag->symflips));
8483f27d899SToby Isaac   }
8499566063dSJacob Faibussowitsch   PetscCall(Petsc1DNodeFamilyDestroy(&(lag->nodeFamily)));
8509566063dSJacob Faibussowitsch   PetscCall(PetscLagNodeIndicesDestroy(&(lag->vertIndices)));
8519566063dSJacob Faibussowitsch   PetscCall(PetscLagNodeIndicesDestroy(&(lag->intNodeIndices)));
8529566063dSJacob Faibussowitsch   PetscCall(PetscLagNodeIndicesDestroy(&(lag->allNodeIndices)));
8539566063dSJacob Faibussowitsch   PetscCall(PetscFree(lag));
8549566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetContinuity_C", NULL));
8559566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetContinuity_C", NULL));
8569566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetTensor_C", NULL));
8579566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetTensor_C", NULL));
8589566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetTrimmed_C", NULL));
8599566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetTrimmed_C", NULL));
8609566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetNodeType_C", NULL));
8619566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetNodeType_C", NULL));
8629566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetUseMoments_C", NULL));
8639566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetUseMoments_C", NULL));
8649566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetMomentOrder_C", NULL));
8659566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetMomentOrder_C", NULL));
86620cf1dd8SToby Isaac   PetscFunctionReturn(0);
86720cf1dd8SToby Isaac }
86820cf1dd8SToby Isaac 
8699371c9d4SSatish Balay static PetscErrorCode PetscDualSpaceLagrangeView_Ascii(PetscDualSpace sp, PetscViewer viewer) {
87020cf1dd8SToby Isaac   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data;
87120cf1dd8SToby Isaac 
87220cf1dd8SToby Isaac   PetscFunctionBegin;
8739566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "%s %s%sLagrange dual space\n", lag->continuous ? "Continuous" : "Discontinuous", lag->tensorSpace ? "tensor " : "", lag->trimmed ? "trimmed " : ""));
87420cf1dd8SToby Isaac   PetscFunctionReturn(0);
87520cf1dd8SToby Isaac }
87620cf1dd8SToby Isaac 
8779371c9d4SSatish Balay static PetscErrorCode PetscDualSpaceView_Lagrange(PetscDualSpace sp, PetscViewer viewer) {
8786f905325SMatthew G. Knepley   PetscBool iascii;
8796f905325SMatthew G. Knepley 
88020cf1dd8SToby Isaac   PetscFunctionBegin;
8816f905325SMatthew G. Knepley   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
8826f905325SMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
8839566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
8849566063dSJacob Faibussowitsch   if (iascii) PetscCall(PetscDualSpaceLagrangeView_Ascii(sp, viewer));
88520cf1dd8SToby Isaac   PetscFunctionReturn(0);
88620cf1dd8SToby Isaac }
88720cf1dd8SToby Isaac 
8889371c9d4SSatish Balay static PetscErrorCode PetscDualSpaceSetFromOptions_Lagrange(PetscDualSpace sp, PetscOptionItems *PetscOptionsObject) {
8893f27d899SToby Isaac   PetscBool       continuous, tensor, trimmed, flg, flg2, flg3;
8903f27d899SToby Isaac   PetscDTNodeType nodeType;
8913f27d899SToby Isaac   PetscReal       nodeExponent;
89266a6c23cSMatthew G. Knepley   PetscInt        momentOrder;
89366a6c23cSMatthew G. Knepley   PetscBool       nodeEndpoints, useMoments;
8946f905325SMatthew G. Knepley 
8956f905325SMatthew G. Knepley   PetscFunctionBegin;
8969566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceLagrangeGetContinuity(sp, &continuous));
8979566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceLagrangeGetTensor(sp, &tensor));
8989566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceLagrangeGetTrimmed(sp, &trimmed));
8999566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceLagrangeGetNodeType(sp, &nodeType, &nodeEndpoints, &nodeExponent));
9003f27d899SToby Isaac   if (nodeType == PETSCDTNODES_DEFAULT) nodeType = PETSCDTNODES_GAUSSJACOBI;
9019566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceLagrangeGetUseMoments(sp, &useMoments));
9029566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceLagrangeGetMomentOrder(sp, &momentOrder));
903d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "PetscDualSpace Lagrange Options");
9049566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-petscdualspace_lagrange_continuity", "Flag for continuous element", "PetscDualSpaceLagrangeSetContinuity", continuous, &continuous, &flg));
9059566063dSJacob Faibussowitsch   if (flg) PetscCall(PetscDualSpaceLagrangeSetContinuity(sp, continuous));
9069566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-petscdualspace_lagrange_tensor", "Flag for tensor dual space", "PetscDualSpaceLagrangeSetTensor", tensor, &tensor, &flg));
9079566063dSJacob Faibussowitsch   if (flg) PetscCall(PetscDualSpaceLagrangeSetTensor(sp, tensor));
9089566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-petscdualspace_lagrange_trimmed", "Flag for trimmed dual space", "PetscDualSpaceLagrangeSetTrimmed", trimmed, &trimmed, &flg));
9099566063dSJacob Faibussowitsch   if (flg) PetscCall(PetscDualSpaceLagrangeSetTrimmed(sp, trimmed));
9109566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-petscdualspace_lagrange_node_type", "Lagrange node location type", "PetscDualSpaceLagrangeSetNodeType", PetscDTNodeTypes, (PetscEnum)nodeType, (PetscEnum *)&nodeType, &flg));
9119566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-petscdualspace_lagrange_node_endpoints", "Flag for nodes that include endpoints", "PetscDualSpaceLagrangeSetNodeType", nodeEndpoints, &nodeEndpoints, &flg2));
9123f27d899SToby Isaac   flg3 = PETSC_FALSE;
91348a46eb9SPierre Jolivet   if (nodeType == PETSCDTNODES_GAUSSJACOBI) PetscCall(PetscOptionsReal("-petscdualspace_lagrange_node_exponent", "Gauss-Jacobi weight function exponent", "PetscDualSpaceLagrangeSetNodeType", nodeExponent, &nodeExponent, &flg3));
9149566063dSJacob Faibussowitsch   if (flg || flg2 || flg3) PetscCall(PetscDualSpaceLagrangeSetNodeType(sp, nodeType, nodeEndpoints, nodeExponent));
9159566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-petscdualspace_lagrange_use_moments", "Use moments (where appropriate) for functionals", "PetscDualSpaceLagrangeSetUseMoments", useMoments, &useMoments, &flg));
9169566063dSJacob Faibussowitsch   if (flg) PetscCall(PetscDualSpaceLagrangeSetUseMoments(sp, useMoments));
9179566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-petscdualspace_lagrange_moment_order", "Quadrature order for moment functionals", "PetscDualSpaceLagrangeSetMomentOrder", momentOrder, &momentOrder, &flg));
9189566063dSJacob Faibussowitsch   if (flg) PetscCall(PetscDualSpaceLagrangeSetMomentOrder(sp, momentOrder));
919d0609cedSBarry Smith   PetscOptionsHeadEnd();
9206f905325SMatthew G. Knepley   PetscFunctionReturn(0);
9216f905325SMatthew G. Knepley }
9226f905325SMatthew G. Knepley 
9239371c9d4SSatish Balay static PetscErrorCode PetscDualSpaceDuplicate_Lagrange(PetscDualSpace sp, PetscDualSpace spNew) {
9243f27d899SToby Isaac   PetscBool           cont, tensor, trimmed, boundary;
9253f27d899SToby Isaac   PetscDTNodeType     nodeType;
9263f27d899SToby Isaac   PetscReal           exponent;
9273f27d899SToby Isaac   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data;
9286f905325SMatthew G. Knepley 
9296f905325SMatthew G. Knepley   PetscFunctionBegin;
9309566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceLagrangeGetContinuity(sp, &cont));
9319566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceLagrangeSetContinuity(spNew, cont));
9329566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceLagrangeGetTensor(sp, &tensor));
9339566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceLagrangeSetTensor(spNew, tensor));
9349566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceLagrangeGetTrimmed(sp, &trimmed));
9359566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceLagrangeSetTrimmed(spNew, trimmed));
9369566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceLagrangeGetNodeType(sp, &nodeType, &boundary, &exponent));
9379566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceLagrangeSetNodeType(spNew, nodeType, boundary, exponent));
9383f27d899SToby Isaac   if (lag->nodeFamily) {
9393f27d899SToby Isaac     PetscDualSpace_Lag *lagnew = (PetscDualSpace_Lag *)spNew->data;
9403f27d899SToby Isaac 
9419566063dSJacob Faibussowitsch     PetscCall(Petsc1DNodeFamilyReference(lag->nodeFamily));
9423f27d899SToby Isaac     lagnew->nodeFamily = lag->nodeFamily;
9433f27d899SToby Isaac   }
9446f905325SMatthew G. Knepley   PetscFunctionReturn(0);
9456f905325SMatthew G. Knepley }
9466f905325SMatthew G. Knepley 
94777f1a120SToby Isaac /* for making tensor product spaces: take a dual space and product a segment space that has all the same
94877f1a120SToby Isaac  * specifications (trimmed, continuous, order, node set), except for the form degree */
9499371c9d4SSatish Balay static PetscErrorCode PetscDualSpaceCreateEdgeSubspace_Lagrange(PetscDualSpace sp, PetscInt order, PetscInt k, PetscInt Nc, PetscBool interiorOnly, PetscDualSpace *bdsp) {
9503f27d899SToby Isaac   DM                  K;
9513f27d899SToby Isaac   PetscDualSpace_Lag *newlag;
9526f905325SMatthew G. Knepley 
9536f905325SMatthew G. Knepley   PetscFunctionBegin;
9549566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceDuplicate(sp, bdsp));
9559566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSetFormDegree(*bdsp, k));
9569566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, DMPolytopeTypeSimpleShape(1, PETSC_TRUE), &K));
9579566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSetDM(*bdsp, K));
9589566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&K));
9599566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSetOrder(*bdsp, order));
9609566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSetNumComponents(*bdsp, Nc));
9613f27d899SToby Isaac   newlag               = (PetscDualSpace_Lag *)(*bdsp)->data;
9623f27d899SToby Isaac   newlag->interiorOnly = interiorOnly;
9639566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSetUp(*bdsp));
9643f27d899SToby Isaac   PetscFunctionReturn(0);
9656f905325SMatthew G. Knepley }
9663f27d899SToby Isaac 
9673f27d899SToby Isaac /* just the points, weights aren't handled */
9689371c9d4SSatish Balay static PetscErrorCode PetscQuadratureCreateTensor(PetscQuadrature trace, PetscQuadrature fiber, PetscQuadrature *product) {
9693f27d899SToby Isaac   PetscInt         dimTrace, dimFiber;
9703f27d899SToby Isaac   PetscInt         numPointsTrace, numPointsFiber;
9713f27d899SToby Isaac   PetscInt         dim, numPoints;
9723f27d899SToby Isaac   const PetscReal *pointsTrace;
9733f27d899SToby Isaac   const PetscReal *pointsFiber;
9743f27d899SToby Isaac   PetscReal       *points;
9753f27d899SToby Isaac   PetscInt         i, j, k, p;
9763f27d899SToby Isaac 
9773f27d899SToby Isaac   PetscFunctionBegin;
9789566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(trace, &dimTrace, NULL, &numPointsTrace, &pointsTrace, NULL));
9799566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(fiber, &dimFiber, NULL, &numPointsFiber, &pointsFiber, NULL));
9803f27d899SToby Isaac   dim       = dimTrace + dimFiber;
9813f27d899SToby Isaac   numPoints = numPointsFiber * numPointsTrace;
9829566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numPoints * dim, &points));
9833f27d899SToby Isaac   for (p = 0, j = 0; j < numPointsFiber; j++) {
9843f27d899SToby Isaac     for (i = 0; i < numPointsTrace; i++, p++) {
9853f27d899SToby Isaac       for (k = 0; k < dimTrace; k++) points[p * dim + k] = pointsTrace[i * dimTrace + k];
9863f27d899SToby Isaac       for (k = 0; k < dimFiber; k++) points[p * dim + dimTrace + k] = pointsFiber[j * dimFiber + k];
9873f27d899SToby Isaac     }
9883f27d899SToby Isaac   }
9899566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, product));
9909566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureSetData(*product, dim, 0, numPoints, points, NULL));
9913f27d899SToby Isaac   PetscFunctionReturn(0);
9923f27d899SToby Isaac }
9933f27d899SToby Isaac 
99477f1a120SToby Isaac /* Kronecker tensor product where matrix is considered a matrix of k-forms, so that
99577f1a120SToby Isaac  * the entries in the product matrix are wedge products of the entries in the original matrices */
9969371c9d4SSatish Balay static PetscErrorCode MatTensorAltV(Mat trace, Mat fiber, PetscInt dimTrace, PetscInt kTrace, PetscInt dimFiber, PetscInt kFiber, Mat *product) {
9973f27d899SToby Isaac   PetscInt     mTrace, nTrace, mFiber, nFiber, m, n, k, i, j, l;
9983f27d899SToby Isaac   PetscInt     dim, NkTrace, NkFiber, Nk;
9993f27d899SToby Isaac   PetscInt     dT, dF;
10003f27d899SToby Isaac   PetscInt    *nnzTrace, *nnzFiber, *nnz;
10013f27d899SToby Isaac   PetscInt     iT, iF, jT, jF, il, jl;
10023f27d899SToby Isaac   PetscReal   *workT, *workT2, *workF, *workF2, *work, *workstar;
10033f27d899SToby Isaac   PetscReal   *projT, *projF;
10043f27d899SToby Isaac   PetscReal   *projTstar, *projFstar;
10053f27d899SToby Isaac   PetscReal   *wedgeMat;
10063f27d899SToby Isaac   PetscReal    sign;
10073f27d899SToby Isaac   PetscScalar *workS;
10083f27d899SToby Isaac   Mat          prod;
10093f27d899SToby Isaac   /* this produces dof groups that look like the identity */
10103f27d899SToby Isaac 
10113f27d899SToby Isaac   PetscFunctionBegin;
10129566063dSJacob Faibussowitsch   PetscCall(MatGetSize(trace, &mTrace, &nTrace));
10139566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(dimTrace, PetscAbsInt(kTrace), &NkTrace));
101408401ef6SPierre Jolivet   PetscCheck(nTrace % NkTrace == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "point value space of trace matrix is not a multiple of k-form size");
10159566063dSJacob Faibussowitsch   PetscCall(MatGetSize(fiber, &mFiber, &nFiber));
10169566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(dimFiber, PetscAbsInt(kFiber), &NkFiber));
101708401ef6SPierre Jolivet   PetscCheck(nFiber % NkFiber == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "point value space of fiber matrix is not a multiple of k-form size");
10189566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(mTrace, &nnzTrace, mFiber, &nnzFiber));
10193f27d899SToby Isaac   for (i = 0; i < mTrace; i++) {
10209566063dSJacob Faibussowitsch     PetscCall(MatGetRow(trace, i, &(nnzTrace[i]), NULL, NULL));
102108401ef6SPierre Jolivet     PetscCheck(nnzTrace[i] % NkTrace == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "nonzeros in trace matrix are not in k-form size blocks");
10223f27d899SToby Isaac   }
10233f27d899SToby Isaac   for (i = 0; i < mFiber; i++) {
10249566063dSJacob Faibussowitsch     PetscCall(MatGetRow(fiber, i, &(nnzFiber[i]), NULL, NULL));
102508401ef6SPierre Jolivet     PetscCheck(nnzFiber[i] % NkFiber == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "nonzeros in fiber matrix are not in k-form size blocks");
10263f27d899SToby Isaac   }
10273f27d899SToby Isaac   dim = dimTrace + dimFiber;
10283f27d899SToby Isaac   k   = kFiber + kTrace;
10299566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(k), &Nk));
10303f27d899SToby Isaac   m = mTrace * mFiber;
10319566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(m, &nnz));
10329371c9d4SSatish Balay   for (l = 0, j = 0; j < mFiber; j++)
10339371c9d4SSatish Balay     for (i = 0; i < mTrace; i++, l++) nnz[l] = (nnzTrace[i] / NkTrace) * (nnzFiber[j] / NkFiber) * Nk;
10343f27d899SToby Isaac   n = (nTrace / NkTrace) * (nFiber / NkFiber) * Nk;
10359566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, m, n, 0, nnz, &prod));
10369566063dSJacob Faibussowitsch   PetscCall(PetscFree(nnz));
10379566063dSJacob Faibussowitsch   PetscCall(PetscFree2(nnzTrace, nnzFiber));
10383f27d899SToby Isaac   /* reasoning about which points each dof needs depends on having zeros computed at points preserved */
10399566063dSJacob Faibussowitsch   PetscCall(MatSetOption(prod, MAT_IGNORE_ZERO_ENTRIES, PETSC_FALSE));
10403f27d899SToby Isaac   /* compute pullbacks */
10419566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(kTrace), &dT));
10429566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(kFiber), &dF));
10439566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(dimTrace * dim, &projT, dimFiber * dim, &projF, dT * NkTrace, &projTstar, dF * NkFiber, &projFstar));
10449566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(projT, dimTrace * dim));
10453f27d899SToby Isaac   for (i = 0; i < dimTrace; i++) projT[i * (dim + 1)] = 1.;
10469566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(projF, dimFiber * dim));
10473f27d899SToby Isaac   for (i = 0; i < dimFiber; i++) projF[i * (dim + 1) + dimTrace] = 1.;
10489566063dSJacob Faibussowitsch   PetscCall(PetscDTAltVPullbackMatrix(dim, dimTrace, projT, kTrace, projTstar));
10499566063dSJacob Faibussowitsch   PetscCall(PetscDTAltVPullbackMatrix(dim, dimFiber, projF, kFiber, projFstar));
10509566063dSJacob Faibussowitsch   PetscCall(PetscMalloc5(dT, &workT, dF, &workF, Nk, &work, Nk, &workstar, Nk, &workS));
10519566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(dT, &workT2, dF, &workF2));
10529566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Nk * dT, &wedgeMat));
10533f27d899SToby Isaac   sign = (PetscAbsInt(kTrace * kFiber) & 1) ? -1. : 1.;
10543f27d899SToby Isaac   for (i = 0, iF = 0; iF < mFiber; iF++) {
10553f27d899SToby Isaac     PetscInt           ncolsF, nformsF;
10563f27d899SToby Isaac     const PetscInt    *colsF;
10573f27d899SToby Isaac     const PetscScalar *valsF;
10583f27d899SToby Isaac 
10599566063dSJacob Faibussowitsch     PetscCall(MatGetRow(fiber, iF, &ncolsF, &colsF, &valsF));
10603f27d899SToby Isaac     nformsF = ncolsF / NkFiber;
10613f27d899SToby Isaac     for (iT = 0; iT < mTrace; iT++, i++) {
10623f27d899SToby Isaac       PetscInt           ncolsT, nformsT;
10633f27d899SToby Isaac       const PetscInt    *colsT;
10643f27d899SToby Isaac       const PetscScalar *valsT;
10653f27d899SToby Isaac 
10669566063dSJacob Faibussowitsch       PetscCall(MatGetRow(trace, iT, &ncolsT, &colsT, &valsT));
10673f27d899SToby Isaac       nformsT = ncolsT / NkTrace;
10683f27d899SToby Isaac       for (j = 0, jF = 0; jF < nformsF; jF++) {
10693f27d899SToby Isaac         PetscInt colF = colsF[jF * NkFiber] / NkFiber;
10703f27d899SToby Isaac 
10713f27d899SToby Isaac         for (il = 0; il < dF; il++) {
10723f27d899SToby Isaac           PetscReal val = 0.;
10733f27d899SToby Isaac           for (jl = 0; jl < NkFiber; jl++) val += projFstar[il * NkFiber + jl] * PetscRealPart(valsF[jF * NkFiber + jl]);
10743f27d899SToby Isaac           workF[il] = val;
10753f27d899SToby Isaac         }
10763f27d899SToby Isaac         if (kFiber < 0) {
10773f27d899SToby Isaac           for (il = 0; il < dF; il++) workF2[il] = workF[il];
10789566063dSJacob Faibussowitsch           PetscCall(PetscDTAltVStar(dim, PetscAbsInt(kFiber), 1, workF2, workF));
10793f27d899SToby Isaac         }
10809566063dSJacob Faibussowitsch         PetscCall(PetscDTAltVWedgeMatrix(dim, PetscAbsInt(kFiber), PetscAbsInt(kTrace), workF, wedgeMat));
10813f27d899SToby Isaac         for (jT = 0; jT < nformsT; jT++, j++) {
10823f27d899SToby Isaac           PetscInt           colT = colsT[jT * NkTrace] / NkTrace;
10833f27d899SToby Isaac           PetscInt           col  = colF * (nTrace / NkTrace) + colT;
10843f27d899SToby Isaac           const PetscScalar *vals;
10853f27d899SToby Isaac 
10863f27d899SToby Isaac           for (il = 0; il < dT; il++) {
10873f27d899SToby Isaac             PetscReal val = 0.;
10883f27d899SToby Isaac             for (jl = 0; jl < NkTrace; jl++) val += projTstar[il * NkTrace + jl] * PetscRealPart(valsT[jT * NkTrace + jl]);
10893f27d899SToby Isaac             workT[il] = val;
10903f27d899SToby Isaac           }
10913f27d899SToby Isaac           if (kTrace < 0) {
10923f27d899SToby Isaac             for (il = 0; il < dT; il++) workT2[il] = workT[il];
10939566063dSJacob Faibussowitsch             PetscCall(PetscDTAltVStar(dim, PetscAbsInt(kTrace), 1, workT2, workT));
10943f27d899SToby Isaac           }
10953f27d899SToby Isaac 
10963f27d899SToby Isaac           for (il = 0; il < Nk; il++) {
10973f27d899SToby Isaac             PetscReal val = 0.;
10983f27d899SToby Isaac             for (jl = 0; jl < dT; jl++) val += sign * wedgeMat[il * dT + jl] * workT[jl];
10993f27d899SToby Isaac             work[il] = val;
11003f27d899SToby Isaac           }
11013f27d899SToby Isaac           if (k < 0) {
11029566063dSJacob Faibussowitsch             PetscCall(PetscDTAltVStar(dim, PetscAbsInt(k), -1, work, workstar));
11033f27d899SToby Isaac #if defined(PETSC_USE_COMPLEX)
11043f27d899SToby Isaac             for (l = 0; l < Nk; l++) workS[l] = workstar[l];
11053f27d899SToby Isaac             vals = &workS[0];
11063f27d899SToby Isaac #else
11073f27d899SToby Isaac             vals = &workstar[0];
11083f27d899SToby Isaac #endif
11093f27d899SToby Isaac           } else {
11103f27d899SToby Isaac #if defined(PETSC_USE_COMPLEX)
11113f27d899SToby Isaac             for (l = 0; l < Nk; l++) workS[l] = work[l];
11123f27d899SToby Isaac             vals = &workS[0];
11133f27d899SToby Isaac #else
11143f27d899SToby Isaac             vals = &work[0];
11153f27d899SToby Isaac #endif
11163f27d899SToby Isaac           }
111748a46eb9SPierre Jolivet           for (l = 0; l < Nk; l++) PetscCall(MatSetValue(prod, i, col * Nk + l, vals[l], INSERT_VALUES)); /* Nk */
11183f27d899SToby Isaac         }                                                                                                 /* jT */
11193f27d899SToby Isaac       }                                                                                                   /* jF */
11209566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(trace, iT, &ncolsT, &colsT, &valsT));
11213f27d899SToby Isaac     } /* iT */
11229566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(fiber, iF, &ncolsF, &colsF, &valsF));
11233f27d899SToby Isaac   } /* iF */
11249566063dSJacob Faibussowitsch   PetscCall(PetscFree(wedgeMat));
11259566063dSJacob Faibussowitsch   PetscCall(PetscFree4(projT, projF, projTstar, projFstar));
11269566063dSJacob Faibussowitsch   PetscCall(PetscFree2(workT2, workF2));
11279566063dSJacob Faibussowitsch   PetscCall(PetscFree5(workT, workF, work, workstar, workS));
11289566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(prod, MAT_FINAL_ASSEMBLY));
11299566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(prod, MAT_FINAL_ASSEMBLY));
11303f27d899SToby Isaac   *product = prod;
11313f27d899SToby Isaac   PetscFunctionReturn(0);
11323f27d899SToby Isaac }
11333f27d899SToby Isaac 
113477f1a120SToby Isaac /* Union of quadrature points, with an attempt to identify commont points in the two sets */
11359371c9d4SSatish Balay static PetscErrorCode PetscQuadraturePointsMerge(PetscQuadrature quadA, PetscQuadrature quadB, PetscQuadrature *quadJoint, PetscInt *aToJoint[], PetscInt *bToJoint[]) {
11363f27d899SToby Isaac   PetscInt         dimA, dimB;
11373f27d899SToby Isaac   PetscInt         nA, nB, nJoint, i, j, d;
11383f27d899SToby Isaac   const PetscReal *pointsA;
11393f27d899SToby Isaac   const PetscReal *pointsB;
11403f27d899SToby Isaac   PetscReal       *pointsJoint;
11413f27d899SToby Isaac   PetscInt        *aToJ, *bToJ;
11423f27d899SToby Isaac   PetscQuadrature  qJ;
11433f27d899SToby Isaac 
11443f27d899SToby Isaac   PetscFunctionBegin;
11459566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quadA, &dimA, NULL, &nA, &pointsA, NULL));
11469566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quadB, &dimB, NULL, &nB, &pointsB, NULL));
114708401ef6SPierre Jolivet   PetscCheck(dimA == dimB, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Quadrature points must be in the same dimension");
11483f27d899SToby Isaac   nJoint = nA;
11499566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nA, &aToJ));
11503f27d899SToby Isaac   for (i = 0; i < nA; i++) aToJ[i] = i;
11519566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nB, &bToJ));
11523f27d899SToby Isaac   for (i = 0; i < nB; i++) {
11533f27d899SToby Isaac     for (j = 0; j < nA; j++) {
11543f27d899SToby Isaac       bToJ[i] = -1;
11559371c9d4SSatish Balay       for (d = 0; d < dimA; d++)
11569371c9d4SSatish Balay         if (PetscAbsReal(pointsB[i * dimA + d] - pointsA[j * dimA + d]) > PETSC_SMALL) break;
11573f27d899SToby Isaac       if (d == dimA) {
11583f27d899SToby Isaac         bToJ[i] = j;
11593f27d899SToby Isaac         break;
11603f27d899SToby Isaac       }
11613f27d899SToby Isaac     }
1162ad540459SPierre Jolivet     if (bToJ[i] == -1) bToJ[i] = nJoint++;
11633f27d899SToby Isaac   }
11643f27d899SToby Isaac   *aToJoint = aToJ;
11653f27d899SToby Isaac   *bToJoint = bToJ;
11669566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nJoint * dimA, &pointsJoint));
11679566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(pointsJoint, pointsA, nA * dimA));
11683f27d899SToby Isaac   for (i = 0; i < nB; i++) {
11693f27d899SToby Isaac     if (bToJ[i] >= nA) {
11703f27d899SToby Isaac       for (d = 0; d < dimA; d++) pointsJoint[bToJ[i] * dimA + d] = pointsB[i * dimA + d];
11713f27d899SToby Isaac     }
11723f27d899SToby Isaac   }
11739566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, &qJ));
11749566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureSetData(qJ, dimA, 0, nJoint, pointsJoint, NULL));
11753f27d899SToby Isaac   *quadJoint = qJ;
11763f27d899SToby Isaac   PetscFunctionReturn(0);
11773f27d899SToby Isaac }
11783f27d899SToby Isaac 
117977f1a120SToby Isaac /* Matrices matA and matB are both quadrature -> dof matrices: produce a matrix that is joint quadrature -> union of
118077f1a120SToby Isaac  * dofs, where the joint quadrature was produced by PetscQuadraturePointsMerge */
11819371c9d4SSatish Balay static PetscErrorCode MatricesMerge(Mat matA, Mat matB, PetscInt dim, PetscInt k, PetscInt numMerged, const PetscInt aToMerged[], const PetscInt bToMerged[], Mat *matMerged) {
11823f27d899SToby Isaac   PetscInt  m, n, mA, nA, mB, nB, Nk, i, j, l;
11833f27d899SToby Isaac   Mat       M;
11843f27d899SToby Isaac   PetscInt *nnz;
11853f27d899SToby Isaac   PetscInt  maxnnz;
11863f27d899SToby Isaac   PetscInt *work;
11873f27d899SToby Isaac 
11883f27d899SToby Isaac   PetscFunctionBegin;
11899566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(k), &Nk));
11909566063dSJacob Faibussowitsch   PetscCall(MatGetSize(matA, &mA, &nA));
119108401ef6SPierre Jolivet   PetscCheck(nA % Nk == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "matA column space not a multiple of k-form size");
11929566063dSJacob Faibussowitsch   PetscCall(MatGetSize(matB, &mB, &nB));
119308401ef6SPierre Jolivet   PetscCheck(nB % Nk == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "matB column space not a multiple of k-form size");
11943f27d899SToby Isaac   m = mA + mB;
11953f27d899SToby Isaac   n = numMerged * Nk;
11969566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(m, &nnz));
11973f27d899SToby Isaac   maxnnz = 0;
11983f27d899SToby Isaac   for (i = 0; i < mA; i++) {
11999566063dSJacob Faibussowitsch     PetscCall(MatGetRow(matA, i, &(nnz[i]), NULL, NULL));
120008401ef6SPierre Jolivet     PetscCheck(nnz[i] % Nk == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "nonzeros in matA are not in k-form size blocks");
12013f27d899SToby Isaac     maxnnz = PetscMax(maxnnz, nnz[i]);
12023f27d899SToby Isaac   }
12033f27d899SToby Isaac   for (i = 0; i < mB; i++) {
12049566063dSJacob Faibussowitsch     PetscCall(MatGetRow(matB, i, &(nnz[i + mA]), NULL, NULL));
120508401ef6SPierre Jolivet     PetscCheck(nnz[i + mA] % Nk == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "nonzeros in matB are not in k-form size blocks");
12063f27d899SToby Isaac     maxnnz = PetscMax(maxnnz, nnz[i + mA]);
12073f27d899SToby Isaac   }
12089566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, m, n, 0, nnz, &M));
12099566063dSJacob Faibussowitsch   PetscCall(PetscFree(nnz));
12103f27d899SToby Isaac   /* reasoning about which points each dof needs depends on having zeros computed at points preserved */
12119566063dSJacob Faibussowitsch   PetscCall(MatSetOption(M, MAT_IGNORE_ZERO_ENTRIES, PETSC_FALSE));
12129566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(maxnnz, &work));
12133f27d899SToby Isaac   for (i = 0; i < mA; i++) {
12143f27d899SToby Isaac     const PetscInt    *cols;
12153f27d899SToby Isaac     const PetscScalar *vals;
12163f27d899SToby Isaac     PetscInt           nCols;
12179566063dSJacob Faibussowitsch     PetscCall(MatGetRow(matA, i, &nCols, &cols, &vals));
12183f27d899SToby Isaac     for (j = 0; j < nCols / Nk; j++) {
12193f27d899SToby Isaac       PetscInt newCol = aToMerged[cols[j * Nk] / Nk];
12203f27d899SToby Isaac       for (l = 0; l < Nk; l++) work[j * Nk + l] = newCol * Nk + l;
12213f27d899SToby Isaac     }
12229566063dSJacob Faibussowitsch     PetscCall(MatSetValuesBlocked(M, 1, &i, nCols, work, vals, INSERT_VALUES));
12239566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(matA, i, &nCols, &cols, &vals));
12243f27d899SToby Isaac   }
12253f27d899SToby Isaac   for (i = 0; i < mB; i++) {
12263f27d899SToby Isaac     const PetscInt    *cols;
12273f27d899SToby Isaac     const PetscScalar *vals;
12283f27d899SToby Isaac 
12293f27d899SToby Isaac     PetscInt row = i + mA;
12303f27d899SToby Isaac     PetscInt nCols;
12319566063dSJacob Faibussowitsch     PetscCall(MatGetRow(matB, i, &nCols, &cols, &vals));
12323f27d899SToby Isaac     for (j = 0; j < nCols / Nk; j++) {
12333f27d899SToby Isaac       PetscInt newCol = bToMerged[cols[j * Nk] / Nk];
12343f27d899SToby Isaac       for (l = 0; l < Nk; l++) work[j * Nk + l] = newCol * Nk + l;
12353f27d899SToby Isaac     }
12369566063dSJacob Faibussowitsch     PetscCall(MatSetValuesBlocked(M, 1, &row, nCols, work, vals, INSERT_VALUES));
12379566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(matB, i, &nCols, &cols, &vals));
12383f27d899SToby Isaac   }
12399566063dSJacob Faibussowitsch   PetscCall(PetscFree(work));
12409566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(M, MAT_FINAL_ASSEMBLY));
12419566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(M, MAT_FINAL_ASSEMBLY));
12423f27d899SToby Isaac   *matMerged = M;
12433f27d899SToby Isaac   PetscFunctionReturn(0);
12443f27d899SToby Isaac }
12453f27d899SToby Isaac 
124677f1a120SToby Isaac /* Take a dual space and product a segment space that has all the same specifications (trimmed, continuous, order,
124777f1a120SToby Isaac  * node set), except for the form degree.  For computing boundary dofs and for making tensor product spaces */
12489371c9d4SSatish Balay static PetscErrorCode PetscDualSpaceCreateFacetSubspace_Lagrange(PetscDualSpace sp, DM K, PetscInt f, PetscInt k, PetscInt Ncopies, PetscBool interiorOnly, PetscDualSpace *bdsp) {
12493f27d899SToby Isaac   PetscInt            Nknew, Ncnew;
12503f27d899SToby Isaac   PetscInt            dim, pointDim = -1;
12513f27d899SToby Isaac   PetscInt            depth;
12523f27d899SToby Isaac   DM                  dm;
12533f27d899SToby Isaac   PetscDualSpace_Lag *newlag;
12543f27d899SToby Isaac 
12553f27d899SToby Isaac   PetscFunctionBegin;
12569566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDM(sp, &dm));
12579566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
12589566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(dm, &depth));
12599566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceDuplicate(sp, bdsp));
12609566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSetFormDegree(*bdsp, k));
12613f27d899SToby Isaac   if (!K) {
12623f27d899SToby Isaac     if (depth == dim) {
1263f783ec47SMatthew G. Knepley       DMPolytopeType ct;
12643f27d899SToby Isaac 
12656ff15688SToby Isaac       pointDim = dim - 1;
12669566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCellType(dm, f, &ct));
12679566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, ct, &K));
12683f27d899SToby Isaac     } else if (depth == 1) {
12693f27d899SToby Isaac       pointDim = 0;
12709566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, DM_POLYTOPE_POINT, &K));
12713f27d899SToby Isaac     } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unsupported interpolation state of reference element");
12723f27d899SToby Isaac   } else {
12739566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)K));
12749566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(K, &pointDim));
12753f27d899SToby Isaac   }
12769566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSetDM(*bdsp, K));
12779566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&K));
12789566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(pointDim, PetscAbsInt(k), &Nknew));
12793f27d899SToby Isaac   Ncnew = Nknew * Ncopies;
12809566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSetNumComponents(*bdsp, Ncnew));
12813f27d899SToby Isaac   newlag               = (PetscDualSpace_Lag *)(*bdsp)->data;
12823f27d899SToby Isaac   newlag->interiorOnly = interiorOnly;
12839566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSetUp(*bdsp));
12843f27d899SToby Isaac   PetscFunctionReturn(0);
12853f27d899SToby Isaac }
12863f27d899SToby Isaac 
128777f1a120SToby Isaac /* Construct simplex nodes from a nodefamily, add Nk dof vectors of length Nk at each node.
128877f1a120SToby Isaac  * Return the (quadrature, matrix) form of the dofs and the nodeIndices form as well.
128977f1a120SToby Isaac  *
129077f1a120SToby Isaac  * Sometimes we want a set of nodes to be contained in the interior of the element,
129177f1a120SToby Isaac  * even when the node scheme puts nodes on the boundaries.  numNodeSkip tells
129277f1a120SToby Isaac  * the routine how many "layers" of nodes need to be skipped.
129377f1a120SToby Isaac  * */
12949371c9d4SSatish Balay static PetscErrorCode PetscDualSpaceLagrangeCreateSimplexNodeMat(Petsc1DNodeFamily nodeFamily, PetscInt dim, PetscInt sum, PetscInt Nk, PetscInt numNodeSkip, PetscQuadrature *iNodes, Mat *iMat, PetscLagNodeIndices *nodeIndices) {
12953f27d899SToby Isaac   PetscReal          *extraNodeCoords, *nodeCoords;
12963f27d899SToby Isaac   PetscInt            nNodes, nExtraNodes;
12973f27d899SToby Isaac   PetscInt            i, j, k, extraSum = sum + numNodeSkip * (1 + dim);
12983f27d899SToby Isaac   PetscQuadrature     intNodes;
12993f27d899SToby Isaac   Mat                 intMat;
13003f27d899SToby Isaac   PetscLagNodeIndices ni;
13013f27d899SToby Isaac 
13023f27d899SToby Isaac   PetscFunctionBegin;
13039566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(dim + sum, dim, &nNodes));
13049566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(dim + extraSum, dim, &nExtraNodes));
13053f27d899SToby Isaac 
13069566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(dim * nExtraNodes, &extraNodeCoords));
13079566063dSJacob Faibussowitsch   PetscCall(PetscNew(&ni));
13083f27d899SToby Isaac   ni->nodeIdxDim = dim + 1;
13093f27d899SToby Isaac   ni->nodeVecDim = Nk;
13103f27d899SToby Isaac   ni->nNodes     = nNodes * Nk;
13113f27d899SToby Isaac   ni->refct      = 1;
13129566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nNodes * Nk * (dim + 1), &(ni->nodeIdx)));
13139566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nNodes * Nk * Nk, &(ni->nodeVec)));
13149371c9d4SSatish Balay   for (i = 0; i < nNodes; i++)
13159371c9d4SSatish Balay     for (j = 0; j < Nk; j++)
13169371c9d4SSatish Balay       for (k = 0; k < Nk; k++) ni->nodeVec[(i * Nk + j) * Nk + k] = (j == k) ? 1. : 0.;
13179566063dSJacob Faibussowitsch   PetscCall(Petsc1DNodeFamilyComputeSimplexNodes(nodeFamily, dim, extraSum, extraNodeCoords));
13183f27d899SToby Isaac   if (numNodeSkip) {
13193f27d899SToby Isaac     PetscInt  k;
13203f27d899SToby Isaac     PetscInt *tup;
13213f27d899SToby Isaac 
13229566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dim * nNodes, &nodeCoords));
13239566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dim + 1, &tup));
13243f27d899SToby Isaac     for (k = 0; k < nNodes; k++) {
13253f27d899SToby Isaac       PetscInt j, c;
13263f27d899SToby Isaac       PetscInt index;
13273f27d899SToby Isaac 
13289566063dSJacob Faibussowitsch       PetscCall(PetscDTIndexToBary(dim + 1, sum, k, tup));
13293f27d899SToby Isaac       for (j = 0; j < dim + 1; j++) tup[j] += numNodeSkip;
13303f27d899SToby Isaac       for (c = 0; c < Nk; c++) {
1331ad540459SPierre Jolivet         for (j = 0; j < dim + 1; j++) ni->nodeIdx[(k * Nk + c) * (dim + 1) + j] = tup[j] + 1;
13323f27d899SToby Isaac       }
13339566063dSJacob Faibussowitsch       PetscCall(PetscDTBaryToIndex(dim + 1, extraSum, tup, &index));
13343f27d899SToby Isaac       for (j = 0; j < dim; j++) nodeCoords[k * dim + j] = extraNodeCoords[index * dim + j];
13353f27d899SToby Isaac     }
13369566063dSJacob Faibussowitsch     PetscCall(PetscFree(tup));
13379566063dSJacob Faibussowitsch     PetscCall(PetscFree(extraNodeCoords));
13383f27d899SToby Isaac   } else {
13393f27d899SToby Isaac     PetscInt  k;
13403f27d899SToby Isaac     PetscInt *tup;
13413f27d899SToby Isaac 
13423f27d899SToby Isaac     nodeCoords = extraNodeCoords;
13439566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dim + 1, &tup));
13443f27d899SToby Isaac     for (k = 0; k < nNodes; k++) {
13453f27d899SToby Isaac       PetscInt j, c;
13463f27d899SToby Isaac 
13479566063dSJacob Faibussowitsch       PetscCall(PetscDTIndexToBary(dim + 1, sum, k, tup));
13483f27d899SToby Isaac       for (c = 0; c < Nk; c++) {
13493f27d899SToby Isaac         for (j = 0; j < dim + 1; j++) {
13503f27d899SToby Isaac           /* barycentric indices can have zeros, but we don't want to push forward zeros because it makes it harder to
135177f1a120SToby Isaac            * determine which nodes correspond to which under symmetries, so we increase by 1.  This is fine
135277f1a120SToby Isaac            * because the nodeIdx coordinates don't have any meaning other than helping to identify symmetries */
13533f27d899SToby Isaac           ni->nodeIdx[(k * Nk + c) * (dim + 1) + j] = tup[j] + 1;
13543f27d899SToby Isaac         }
13553f27d899SToby Isaac       }
13563f27d899SToby Isaac     }
13579566063dSJacob Faibussowitsch     PetscCall(PetscFree(tup));
13583f27d899SToby Isaac   }
13599566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, &intNodes));
13609566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureSetData(intNodes, dim, 0, nNodes, nodeCoords, NULL));
13619566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, nNodes * Nk, nNodes * Nk, Nk, NULL, &intMat));
13629566063dSJacob Faibussowitsch   PetscCall(MatSetOption(intMat, MAT_IGNORE_ZERO_ENTRIES, PETSC_FALSE));
13633f27d899SToby Isaac   for (j = 0; j < nNodes * Nk; j++) {
13643f27d899SToby Isaac     PetscInt rem = j % Nk;
13653f27d899SToby Isaac     PetscInt a, aprev = j - rem;
13663f27d899SToby Isaac     PetscInt anext = aprev + Nk;
13673f27d899SToby Isaac 
136848a46eb9SPierre Jolivet     for (a = aprev; a < anext; a++) PetscCall(MatSetValue(intMat, j, a, (a == j) ? 1. : 0., INSERT_VALUES));
13693f27d899SToby Isaac   }
13709566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(intMat, MAT_FINAL_ASSEMBLY));
13719566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(intMat, MAT_FINAL_ASSEMBLY));
13723f27d899SToby Isaac   *iNodes      = intNodes;
13733f27d899SToby Isaac   *iMat        = intMat;
13743f27d899SToby Isaac   *nodeIndices = ni;
13753f27d899SToby Isaac   PetscFunctionReturn(0);
13763f27d899SToby Isaac }
13773f27d899SToby Isaac 
137877f1a120SToby Isaac /* once the nodeIndices have been created for the interior of the reference cell, and for all of the boundary cells,
1379a5b23f4aSJose E. Roman  * push forward the boundary dofs and concatenate them into the full node indices for the dual space */
13809371c9d4SSatish Balay static PetscErrorCode PetscDualSpaceLagrangeCreateAllNodeIdx(PetscDualSpace sp) {
13813f27d899SToby Isaac   DM                  dm;
13823f27d899SToby Isaac   PetscInt            dim, nDofs;
13833f27d899SToby Isaac   PetscSection        section;
13843f27d899SToby Isaac   PetscInt            pStart, pEnd, p;
13853f27d899SToby Isaac   PetscInt            formDegree, Nk;
13863f27d899SToby Isaac   PetscInt            nodeIdxDim, spintdim;
13873f27d899SToby Isaac   PetscDualSpace_Lag *lag;
13883f27d899SToby Isaac   PetscLagNodeIndices ni, verti;
13893f27d899SToby Isaac 
13903f27d899SToby Isaac   PetscFunctionBegin;
13913f27d899SToby Isaac   lag   = (PetscDualSpace_Lag *)sp->data;
13923f27d899SToby Isaac   verti = lag->vertIndices;
13939566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDM(sp, &dm));
13949566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
13959566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetFormDegree(sp, &formDegree));
13969566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(formDegree), &Nk));
13979566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetSection(sp, &section));
13989566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(section, &nDofs));
13999566063dSJacob Faibussowitsch   PetscCall(PetscNew(&ni));
14003f27d899SToby Isaac   ni->nodeIdxDim = nodeIdxDim = verti->nodeIdxDim;
14013f27d899SToby Isaac   ni->nodeVecDim              = Nk;
14023f27d899SToby Isaac   ni->nNodes                  = nDofs;
14033f27d899SToby Isaac   ni->refct                   = 1;
14049566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nodeIdxDim * nDofs, &(ni->nodeIdx)));
14059566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Nk * nDofs, &(ni->nodeVec)));
14069566063dSJacob Faibussowitsch   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
14079566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetDof(section, 0, &spintdim));
14083f27d899SToby Isaac   if (spintdim) {
14099566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(ni->nodeIdx, lag->intNodeIndices->nodeIdx, spintdim * nodeIdxDim));
14109566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(ni->nodeVec, lag->intNodeIndices->nodeVec, spintdim * Nk));
14113f27d899SToby Isaac   }
14123f27d899SToby Isaac   for (p = pStart + 1; p < pEnd; p++) {
14133f27d899SToby Isaac     PetscDualSpace      psp = sp->pointSpaces[p];
14143f27d899SToby Isaac     PetscDualSpace_Lag *plag;
14153f27d899SToby Isaac     PetscInt            dof, off;
14163f27d899SToby Isaac 
14179566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, p, &dof));
14183f27d899SToby Isaac     if (!dof) continue;
14193f27d899SToby Isaac     plag = (PetscDualSpace_Lag *)psp->data;
14209566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(section, p, &off));
14219566063dSJacob Faibussowitsch     PetscCall(PetscLagNodeIndicesPushForward(dm, verti, p, plag->vertIndices, plag->intNodeIndices, 0, formDegree, &(ni->nodeIdx[off * nodeIdxDim]), &(ni->nodeVec[off * Nk])));
14223f27d899SToby Isaac   }
14233f27d899SToby Isaac   lag->allNodeIndices = ni;
14243f27d899SToby Isaac   PetscFunctionReturn(0);
14253f27d899SToby Isaac }
14263f27d899SToby Isaac 
142777f1a120SToby Isaac /* once the (quadrature, Matrix) forms of the dofs have been created for the interior of the
142877f1a120SToby Isaac  * reference cell and for the boundary cells, jk
142977f1a120SToby Isaac  * push forward the boundary data and concatenate them into the full (quadrature, matrix) data
143077f1a120SToby Isaac  * for the dual space */
14319371c9d4SSatish Balay static PetscErrorCode PetscDualSpaceCreateAllDataFromInteriorData(PetscDualSpace sp) {
14323f27d899SToby Isaac   DM              dm;
14333f27d899SToby Isaac   PetscSection    section;
14343f27d899SToby Isaac   PetscInt        pStart, pEnd, p, k, Nk, dim, Nc;
14353f27d899SToby Isaac   PetscInt        nNodes;
14363f27d899SToby Isaac   PetscInt        countNodes;
14373f27d899SToby Isaac   Mat             allMat;
14383f27d899SToby Isaac   PetscQuadrature allNodes;
14393f27d899SToby Isaac   PetscInt        nDofs;
14403f27d899SToby Isaac   PetscInt        maxNzforms, j;
14413f27d899SToby Isaac   PetscScalar    *work;
14423f27d899SToby Isaac   PetscReal      *L, *J, *Jinv, *v0, *pv0;
14433f27d899SToby Isaac   PetscInt       *iwork;
14443f27d899SToby Isaac   PetscReal      *nodes;
14453f27d899SToby Isaac 
14463f27d899SToby Isaac   PetscFunctionBegin;
14479566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDM(sp, &dm));
14489566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
14499566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetSection(sp, &section));
14509566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(section, &nDofs));
14519566063dSJacob Faibussowitsch   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
14529566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetFormDegree(sp, &k));
14539566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetNumComponents(sp, &Nc));
14549566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(k), &Nk));
14553f27d899SToby Isaac   for (p = pStart, nNodes = 0, maxNzforms = 0; p < pEnd; p++) {
14563f27d899SToby Isaac     PetscDualSpace  psp;
14573f27d899SToby Isaac     DM              pdm;
14583f27d899SToby Isaac     PetscInt        pdim, pNk;
14593f27d899SToby Isaac     PetscQuadrature intNodes;
14603f27d899SToby Isaac     Mat             intMat;
14613f27d899SToby Isaac 
14629566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetPointSubspace(sp, p, &psp));
14633f27d899SToby Isaac     if (!psp) continue;
14649566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDM(psp, &pdm));
14659566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(pdm, &pdim));
14663f27d899SToby Isaac     if (pdim < PetscAbsInt(k)) continue;
14679566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(pdim, PetscAbsInt(k), &pNk));
14689566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetInteriorData(psp, &intNodes, &intMat));
14693f27d899SToby Isaac     if (intNodes) {
14703f27d899SToby Isaac       PetscInt nNodesp;
14713f27d899SToby Isaac 
14729566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureGetData(intNodes, NULL, NULL, &nNodesp, NULL, NULL));
14733f27d899SToby Isaac       nNodes += nNodesp;
14743f27d899SToby Isaac     }
14753f27d899SToby Isaac     if (intMat) {
14763f27d899SToby Isaac       PetscInt maxNzsp;
14773f27d899SToby Isaac       PetscInt maxNzformsp;
14783f27d899SToby Isaac 
14799566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJGetMaxRowNonzeros(intMat, &maxNzsp));
148008401ef6SPierre Jolivet       PetscCheck(maxNzsp % pNk == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "interior matrix is not laid out as blocks of k-forms");
14813f27d899SToby Isaac       maxNzformsp = maxNzsp / pNk;
14823f27d899SToby Isaac       maxNzforms  = PetscMax(maxNzforms, maxNzformsp);
14833f27d899SToby Isaac     }
14843f27d899SToby Isaac   }
14859566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, nDofs, nNodes * Nc, maxNzforms * Nk, NULL, &allMat));
14869566063dSJacob Faibussowitsch   PetscCall(MatSetOption(allMat, MAT_IGNORE_ZERO_ENTRIES, PETSC_FALSE));
14879566063dSJacob Faibussowitsch   PetscCall(PetscMalloc7(dim, &v0, dim, &pv0, dim * dim, &J, dim * dim, &Jinv, Nk * Nk, &L, maxNzforms * Nk, &work, maxNzforms * Nk, &iwork));
14883f27d899SToby Isaac   for (j = 0; j < dim; j++) pv0[j] = -1.;
14899566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(dim * nNodes, &nodes));
14903f27d899SToby Isaac   for (p = pStart, countNodes = 0; p < pEnd; p++) {
14913f27d899SToby Isaac     PetscDualSpace  psp;
14923f27d899SToby Isaac     PetscQuadrature intNodes;
14933f27d899SToby Isaac     DM              pdm;
14943f27d899SToby Isaac     PetscInt        pdim, pNk;
14953f27d899SToby Isaac     PetscInt        countNodesIn = countNodes;
14963f27d899SToby Isaac     PetscReal       detJ;
14973f27d899SToby Isaac     Mat             intMat;
14983f27d899SToby Isaac 
14999566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetPointSubspace(sp, p, &psp));
15003f27d899SToby Isaac     if (!psp) continue;
15019566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDM(psp, &pdm));
15029566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(pdm, &pdim));
15033f27d899SToby Isaac     if (pdim < PetscAbsInt(k)) continue;
15049566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetInteriorData(psp, &intNodes, &intMat));
15053f27d899SToby Isaac     if (intNodes == NULL && intMat == NULL) continue;
15069566063dSJacob Faibussowitsch     PetscCall(PetscDTBinomialInt(pdim, PetscAbsInt(k), &pNk));
15073f27d899SToby Isaac     if (p) {
15089566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryAffineFEM(dm, p, v0, J, Jinv, &detJ));
15093f27d899SToby Isaac     } else { /* identity */
15103f27d899SToby Isaac       PetscInt i, j;
15113f27d899SToby Isaac 
15129371c9d4SSatish Balay       for (i = 0; i < dim; i++)
15139371c9d4SSatish Balay         for (j = 0; j < dim; j++) J[i * dim + j] = Jinv[i * dim + j] = 0.;
15143f27d899SToby Isaac       for (i = 0; i < dim; i++) J[i * dim + i] = Jinv[i * dim + i] = 1.;
15153f27d899SToby Isaac       for (i = 0; i < dim; i++) v0[i] = -1.;
15163f27d899SToby Isaac     }
15173f27d899SToby Isaac     if (pdim != dim) { /* compactify Jacobian */
15183f27d899SToby Isaac       PetscInt i, j;
15193f27d899SToby Isaac 
15209371c9d4SSatish Balay       for (i = 0; i < dim; i++)
15219371c9d4SSatish Balay         for (j = 0; j < pdim; j++) J[i * pdim + j] = J[i * dim + j];
15223f27d899SToby Isaac     }
15239566063dSJacob Faibussowitsch     PetscCall(PetscDTAltVPullbackMatrix(pdim, dim, J, k, L));
152477f1a120SToby Isaac     if (intNodes) { /* push forward quadrature locations by the affine transformation */
15253f27d899SToby Isaac       PetscInt         nNodesp;
15263f27d899SToby Isaac       const PetscReal *nodesp;
15273f27d899SToby Isaac       PetscInt         j;
15283f27d899SToby Isaac 
15299566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureGetData(intNodes, NULL, NULL, &nNodesp, &nodesp, NULL));
15303f27d899SToby Isaac       for (j = 0; j < nNodesp; j++, countNodes++) {
15313f27d899SToby Isaac         PetscInt d, e;
15323f27d899SToby Isaac 
15333f27d899SToby Isaac         for (d = 0; d < dim; d++) {
15343f27d899SToby Isaac           nodes[countNodes * dim + d] = v0[d];
1535ad540459SPierre Jolivet           for (e = 0; e < pdim; e++) nodes[countNodes * dim + d] += J[d * pdim + e] * (nodesp[j * pdim + e] - pv0[e]);
15363f27d899SToby Isaac         }
15373f27d899SToby Isaac       }
15383f27d899SToby Isaac     }
15393f27d899SToby Isaac     if (intMat) {
15403f27d899SToby Isaac       PetscInt nrows;
15413f27d899SToby Isaac       PetscInt off;
15423f27d899SToby Isaac 
15439566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(section, p, &nrows));
15449566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(section, p, &off));
15453f27d899SToby Isaac       for (j = 0; j < nrows; j++) {
15463f27d899SToby Isaac         PetscInt           ncols;
15473f27d899SToby Isaac         const PetscInt    *cols;
15483f27d899SToby Isaac         const PetscScalar *vals;
15493f27d899SToby Isaac         PetscInt           l, d, e;
15503f27d899SToby Isaac         PetscInt           row = j + off;
15513f27d899SToby Isaac 
15529566063dSJacob Faibussowitsch         PetscCall(MatGetRow(intMat, j, &ncols, &cols, &vals));
155308401ef6SPierre Jolivet         PetscCheck(ncols % pNk == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "interior matrix is not laid out as blocks of k-forms");
15543f27d899SToby Isaac         for (l = 0; l < ncols / pNk; l++) {
15553f27d899SToby Isaac           PetscInt blockcol;
15563f27d899SToby Isaac 
1557ad540459SPierre Jolivet           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");
15583f27d899SToby Isaac           blockcol = cols[l * pNk] / pNk;
1559ad540459SPierre Jolivet           for (d = 0; d < Nk; d++) iwork[l * Nk + d] = (blockcol + countNodesIn) * Nk + d;
15603f27d899SToby Isaac           for (d = 0; d < Nk; d++) work[l * Nk + d] = 0.;
15613f27d899SToby Isaac           for (d = 0; d < Nk; d++) {
15623f27d899SToby Isaac             for (e = 0; e < pNk; e++) {
15633f27d899SToby Isaac               /* "push forward" dof by pulling back a k-form to be evaluated on the point: multiply on the right by L */
15645efe5503SToby Isaac               work[l * Nk + d] += vals[l * pNk + e] * L[e * Nk + d];
15653f27d899SToby Isaac             }
15663f27d899SToby Isaac           }
15673f27d899SToby Isaac         }
15689566063dSJacob Faibussowitsch         PetscCall(MatSetValues(allMat, 1, &row, (ncols / pNk) * Nk, iwork, work, INSERT_VALUES));
15699566063dSJacob Faibussowitsch         PetscCall(MatRestoreRow(intMat, j, &ncols, &cols, &vals));
15703f27d899SToby Isaac       }
15713f27d899SToby Isaac     }
15723f27d899SToby Isaac   }
15739566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(allMat, MAT_FINAL_ASSEMBLY));
15749566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(allMat, MAT_FINAL_ASSEMBLY));
15759566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, &allNodes));
15769566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureSetData(allNodes, dim, 0, nNodes, nodes, NULL));
15779566063dSJacob Faibussowitsch   PetscCall(PetscFree7(v0, pv0, J, Jinv, L, work, iwork));
15789566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&(sp->allMat)));
15793f27d899SToby Isaac   sp->allMat = allMat;
15809566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureDestroy(&(sp->allNodes)));
15813f27d899SToby Isaac   sp->allNodes = allNodes;
15823f27d899SToby Isaac   PetscFunctionReturn(0);
15833f27d899SToby Isaac }
15843f27d899SToby Isaac 
158577f1a120SToby Isaac /* rather than trying to get all data from the functionals, we create
158677f1a120SToby Isaac  * the functionals from rows of the quadrature -> dof matrix.
158777f1a120SToby Isaac  *
158877f1a120SToby Isaac  * Ideally most of the uses of PetscDualSpace in PetscFE will switch
158977f1a120SToby Isaac  * to using intMat and allMat, so that the individual functionals
159077f1a120SToby Isaac  * don't need to be constructed at all */
15919371c9d4SSatish Balay static PetscErrorCode PetscDualSpaceComputeFunctionalsFromAllData(PetscDualSpace sp) {
15923f27d899SToby Isaac   PetscQuadrature  allNodes;
15933f27d899SToby Isaac   Mat              allMat;
15943f27d899SToby Isaac   PetscInt         nDofs;
15953f27d899SToby Isaac   PetscInt         dim, k, Nk, Nc, f;
15963f27d899SToby Isaac   DM               dm;
15973f27d899SToby Isaac   PetscInt         nNodes, spdim;
15983f27d899SToby Isaac   const PetscReal *nodes = NULL;
15993f27d899SToby Isaac   PetscSection     section;
160066a6c23cSMatthew G. Knepley   PetscBool        useMoments;
16013f27d899SToby Isaac 
16023f27d899SToby Isaac   PetscFunctionBegin;
16039566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDM(sp, &dm));
16049566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
16059566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetNumComponents(sp, &Nc));
16069566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetFormDegree(sp, &k));
16079566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(k), &Nk));
16089566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetAllData(sp, &allNodes, &allMat));
16093f27d899SToby Isaac   nNodes = 0;
161048a46eb9SPierre Jolivet   if (allNodes) PetscCall(PetscQuadratureGetData(allNodes, NULL, NULL, &nNodes, &nodes, NULL));
16119566063dSJacob Faibussowitsch   PetscCall(MatGetSize(allMat, &nDofs, NULL));
16129566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetSection(sp, &section));
16139566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(section, &spdim));
161408401ef6SPierre Jolivet   PetscCheck(spdim == nDofs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "incompatible all matrix size");
16159566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nDofs, &(sp->functional)));
16169566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceLagrangeGetUseMoments(sp, &useMoments));
161766a6c23cSMatthew G. Knepley   if (useMoments) {
161866a6c23cSMatthew G. Knepley     Mat              allMat;
161966a6c23cSMatthew G. Knepley     PetscInt         momentOrder, i;
162066a6c23cSMatthew G. Knepley     PetscBool        tensor;
162166a6c23cSMatthew G. Knepley     const PetscReal *weights;
162266a6c23cSMatthew G. Knepley     PetscScalar     *array;
162366a6c23cSMatthew G. Knepley 
162463a3b9bcSJacob Faibussowitsch     PetscCheck(nDofs == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "We do not yet support moments beyond P0, nDofs == %" PetscInt_FMT, nDofs);
16259566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceLagrangeGetMomentOrder(sp, &momentOrder));
16269566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceLagrangeGetTensor(sp, &tensor));
16279566063dSJacob Faibussowitsch     if (!tensor) PetscCall(PetscDTStroudConicalQuadrature(dim, Nc, PetscMax(momentOrder + 1, 1), -1.0, 1.0, &(sp->functional[0])));
16289566063dSJacob Faibussowitsch     else PetscCall(PetscDTGaussTensorQuadrature(dim, Nc, PetscMax(momentOrder + 1, 1), -1.0, 1.0, &(sp->functional[0])));
162966a6c23cSMatthew G. Knepley     /* Need to replace allNodes and allMat */
16309566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)sp->functional[0]));
16319566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureDestroy(&(sp->allNodes)));
163266a6c23cSMatthew G. Knepley     sp->allNodes = sp->functional[0];
16339566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureGetData(sp->allNodes, NULL, NULL, &nNodes, NULL, &weights));
16349566063dSJacob Faibussowitsch     PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, nDofs, nNodes * Nc, NULL, &allMat));
16359566063dSJacob Faibussowitsch     PetscCall(MatDenseGetArrayWrite(allMat, &array));
163666a6c23cSMatthew G. Knepley     for (i = 0; i < nNodes * Nc; ++i) array[i] = weights[i];
16379566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreArrayWrite(allMat, &array));
16389566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(allMat, MAT_FINAL_ASSEMBLY));
16399566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(allMat, MAT_FINAL_ASSEMBLY));
16409566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&(sp->allMat)));
164166a6c23cSMatthew G. Knepley     sp->allMat = allMat;
164266a6c23cSMatthew G. Knepley     PetscFunctionReturn(0);
164366a6c23cSMatthew G. Knepley   }
16443f27d899SToby Isaac   for (f = 0; f < nDofs; f++) {
16453f27d899SToby Isaac     PetscInt           ncols, c;
16463f27d899SToby Isaac     const PetscInt    *cols;
16473f27d899SToby Isaac     const PetscScalar *vals;
16483f27d899SToby Isaac     PetscReal         *nodesf;
16493f27d899SToby Isaac     PetscReal         *weightsf;
16503f27d899SToby Isaac     PetscInt           nNodesf;
16513f27d899SToby Isaac     PetscInt           countNodes;
16523f27d899SToby Isaac 
16539566063dSJacob Faibussowitsch     PetscCall(MatGetRow(allMat, f, &ncols, &cols, &vals));
165408401ef6SPierre Jolivet     PetscCheck(ncols % Nk == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "all matrix is not laid out as blocks of k-forms");
16553f27d899SToby Isaac     for (c = 1, nNodesf = 1; c < ncols; c++) {
16563f27d899SToby Isaac       if ((cols[c] / Nc) != (cols[c - 1] / Nc)) nNodesf++;
16573f27d899SToby Isaac     }
16589566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(dim * nNodesf, &nodesf));
16599566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(Nc * nNodesf, &weightsf));
16603f27d899SToby Isaac     for (c = 0, countNodes = 0; c < ncols; c++) {
16613f27d899SToby Isaac       if (!c || ((cols[c] / Nc) != (cols[c - 1] / Nc))) {
16623f27d899SToby Isaac         PetscInt d;
16633f27d899SToby Isaac 
1664ad540459SPierre Jolivet         for (d = 0; d < Nc; d++) weightsf[countNodes * Nc + d] = 0.;
1665ad540459SPierre Jolivet         for (d = 0; d < dim; d++) nodesf[countNodes * dim + d] = nodes[(cols[c] / Nc) * dim + d];
16663f27d899SToby Isaac         countNodes++;
16673f27d899SToby Isaac       }
16683f27d899SToby Isaac       weightsf[(countNodes - 1) * Nc + (cols[c] % Nc)] = PetscRealPart(vals[c]);
16693f27d899SToby Isaac     }
16709566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, &(sp->functional[f])));
16719566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureSetData(sp->functional[f], dim, Nc, nNodesf, nodesf, weightsf));
16729566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(allMat, f, &ncols, &cols, &vals));
16733f27d899SToby Isaac   }
16743f27d899SToby Isaac   PetscFunctionReturn(0);
16753f27d899SToby Isaac }
16763f27d899SToby Isaac 
16773f27d899SToby Isaac /* take a matrix meant for k-forms and expand it to one for Ncopies */
16789371c9d4SSatish Balay static PetscErrorCode PetscDualSpaceLagrangeMatrixCreateCopies(Mat A, PetscInt Nk, PetscInt Ncopies, Mat *Abs) {
16793f27d899SToby Isaac   PetscInt m, n, i, j, k;
16803f27d899SToby Isaac   PetscInt maxnnz, *nnz, *iwork;
16813f27d899SToby Isaac   Mat      Ac;
16823f27d899SToby Isaac 
16833f27d899SToby Isaac   PetscFunctionBegin;
16849566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &m, &n));
168563a3b9bcSJacob Faibussowitsch   PetscCheck(n % Nk == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of columns in A %" PetscInt_FMT " is not a multiple of Nk %" PetscInt_FMT, n, Nk);
16869566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(m * Ncopies, &nnz));
16873f27d899SToby Isaac   for (i = 0, maxnnz = 0; i < m; i++) {
16883f27d899SToby Isaac     PetscInt innz;
16899566063dSJacob Faibussowitsch     PetscCall(MatGetRow(A, i, &innz, NULL, NULL));
169063a3b9bcSJacob Faibussowitsch     PetscCheck(innz % Nk == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "A row %" PetscInt_FMT " nnzs is not a multiple of Nk %" PetscInt_FMT, innz, Nk);
16913f27d899SToby Isaac     for (j = 0; j < Ncopies; j++) nnz[i * Ncopies + j] = innz;
16923f27d899SToby Isaac     maxnnz = PetscMax(maxnnz, innz);
16933f27d899SToby Isaac   }
16949566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, m * Ncopies, n * Ncopies, 0, nnz, &Ac));
16959566063dSJacob Faibussowitsch   PetscCall(MatSetOption(Ac, MAT_IGNORE_ZERO_ENTRIES, PETSC_FALSE));
16969566063dSJacob Faibussowitsch   PetscCall(PetscFree(nnz));
16979566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(maxnnz, &iwork));
16983f27d899SToby Isaac   for (i = 0; i < m; i++) {
16993f27d899SToby Isaac     PetscInt           innz;
17003f27d899SToby Isaac     const PetscInt    *cols;
17013f27d899SToby Isaac     const PetscScalar *vals;
17023f27d899SToby Isaac 
17039566063dSJacob Faibussowitsch     PetscCall(MatGetRow(A, i, &innz, &cols, &vals));
17043f27d899SToby Isaac     for (j = 0; j < innz; j++) iwork[j] = (cols[j] / Nk) * (Nk * Ncopies) + (cols[j] % Nk);
17053f27d899SToby Isaac     for (j = 0; j < Ncopies; j++) {
17063f27d899SToby Isaac       PetscInt row = i * Ncopies + j;
17073f27d899SToby Isaac 
17089566063dSJacob Faibussowitsch       PetscCall(MatSetValues(Ac, 1, &row, innz, iwork, vals, INSERT_VALUES));
17093f27d899SToby Isaac       for (k = 0; k < innz; k++) iwork[k] += Nk;
17103f27d899SToby Isaac     }
17119566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(A, i, &innz, &cols, &vals));
17123f27d899SToby Isaac   }
17139566063dSJacob Faibussowitsch   PetscCall(PetscFree(iwork));
17149566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Ac, MAT_FINAL_ASSEMBLY));
17159566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Ac, MAT_FINAL_ASSEMBLY));
17163f27d899SToby Isaac   *Abs = Ac;
17173f27d899SToby Isaac   PetscFunctionReturn(0);
17183f27d899SToby Isaac }
17193f27d899SToby Isaac 
172077f1a120SToby Isaac /* check if a cell is a tensor product of the segment with a facet,
172177f1a120SToby Isaac  * specifically checking if f and f2 can be the "endpoints" (like the triangles
172277f1a120SToby Isaac  * at either end of a wedge) */
17239371c9d4SSatish Balay static PetscErrorCode DMPlexPointIsTensor_Internal_Given(DM dm, PetscInt p, PetscInt f, PetscInt f2, PetscBool *isTensor) {
17243f27d899SToby Isaac   PetscInt        coneSize, c;
17253f27d899SToby Isaac   const PetscInt *cone;
17263f27d899SToby Isaac   const PetscInt *fCone;
17273f27d899SToby Isaac   const PetscInt *f2Cone;
17283f27d899SToby Isaac   PetscInt        fs[2];
17293f27d899SToby Isaac   PetscInt        meetSize, nmeet;
17303f27d899SToby Isaac   const PetscInt *meet;
17313f27d899SToby Isaac 
17323f27d899SToby Isaac   PetscFunctionBegin;
17333f27d899SToby Isaac   fs[0] = f;
17343f27d899SToby Isaac   fs[1] = f2;
17359566063dSJacob Faibussowitsch   PetscCall(DMPlexGetMeet(dm, 2, fs, &meetSize, &meet));
17363f27d899SToby Isaac   nmeet = meetSize;
17379566063dSJacob Faibussowitsch   PetscCall(DMPlexRestoreMeet(dm, 2, fs, &meetSize, &meet));
173877f1a120SToby Isaac   /* two points that have a non-empty meet cannot be at opposite ends of a cell */
17393f27d899SToby Isaac   if (nmeet) {
17403f27d899SToby Isaac     *isTensor = PETSC_FALSE;
17413f27d899SToby Isaac     PetscFunctionReturn(0);
17423f27d899SToby Isaac   }
17439566063dSJacob Faibussowitsch   PetscCall(DMPlexGetConeSize(dm, p, &coneSize));
17449566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCone(dm, p, &cone));
17459566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCone(dm, f, &fCone));
17469566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCone(dm, f2, &f2Cone));
17473f27d899SToby Isaac   for (c = 0; c < coneSize; c++) {
17483f27d899SToby Isaac     PetscInt        e, ef;
17493f27d899SToby Isaac     PetscInt        d = -1, d2 = -1;
17503f27d899SToby Isaac     PetscInt        dcount, d2count;
17513f27d899SToby Isaac     PetscInt        t = cone[c];
17523f27d899SToby Isaac     PetscInt        tConeSize;
17533f27d899SToby Isaac     PetscBool       tIsTensor;
17543f27d899SToby Isaac     const PetscInt *tCone;
17553f27d899SToby Isaac 
17563f27d899SToby Isaac     if (t == f || t == f2) continue;
175777f1a120SToby Isaac     /* for every other facet in the cone, check that is has
175877f1a120SToby Isaac      * one ridge in common with each end */
17599566063dSJacob Faibussowitsch     PetscCall(DMPlexGetConeSize(dm, t, &tConeSize));
17609566063dSJacob Faibussowitsch     PetscCall(DMPlexGetCone(dm, t, &tCone));
17613f27d899SToby Isaac 
17623f27d899SToby Isaac     dcount  = 0;
17633f27d899SToby Isaac     d2count = 0;
17643f27d899SToby Isaac     for (e = 0; e < tConeSize; e++) {
17653f27d899SToby Isaac       PetscInt q = tCone[e];
17663f27d899SToby Isaac       for (ef = 0; ef < coneSize - 2; ef++) {
17673f27d899SToby Isaac         if (fCone[ef] == q) {
17683f27d899SToby Isaac           if (dcount) {
17693f27d899SToby Isaac             *isTensor = PETSC_FALSE;
17703f27d899SToby Isaac             PetscFunctionReturn(0);
17713f27d899SToby Isaac           }
17723f27d899SToby Isaac           d = q;
17733f27d899SToby Isaac           dcount++;
17743f27d899SToby Isaac         } else if (f2Cone[ef] == q) {
17753f27d899SToby Isaac           if (d2count) {
17763f27d899SToby Isaac             *isTensor = PETSC_FALSE;
17773f27d899SToby Isaac             PetscFunctionReturn(0);
17783f27d899SToby Isaac           }
17793f27d899SToby Isaac           d2 = q;
17803f27d899SToby Isaac           d2count++;
17813f27d899SToby Isaac         }
17823f27d899SToby Isaac       }
17833f27d899SToby Isaac     }
178477f1a120SToby Isaac     /* if the whole cell is a tensor with the segment, then this
178577f1a120SToby Isaac      * facet should be a tensor with the segment */
17869566063dSJacob Faibussowitsch     PetscCall(DMPlexPointIsTensor_Internal_Given(dm, t, d, d2, &tIsTensor));
17873f27d899SToby Isaac     if (!tIsTensor) {
17883f27d899SToby Isaac       *isTensor = PETSC_FALSE;
17893f27d899SToby Isaac       PetscFunctionReturn(0);
17903f27d899SToby Isaac     }
17913f27d899SToby Isaac   }
17923f27d899SToby Isaac   *isTensor = PETSC_TRUE;
17933f27d899SToby Isaac   PetscFunctionReturn(0);
17943f27d899SToby Isaac }
17953f27d899SToby Isaac 
179677f1a120SToby Isaac /* determine if a cell is a tensor with a segment by looping over pairs of facets to find a pair
179777f1a120SToby Isaac  * that could be the opposite ends */
17989371c9d4SSatish Balay static PetscErrorCode DMPlexPointIsTensor_Internal(DM dm, PetscInt p, PetscBool *isTensor, PetscInt *endA, PetscInt *endB) {
17993f27d899SToby Isaac   PetscInt        coneSize, c, c2;
18003f27d899SToby Isaac   const PetscInt *cone;
18013f27d899SToby Isaac 
18023f27d899SToby Isaac   PetscFunctionBegin;
18039566063dSJacob Faibussowitsch   PetscCall(DMPlexGetConeSize(dm, p, &coneSize));
18043f27d899SToby Isaac   if (!coneSize) {
18053f27d899SToby Isaac     if (isTensor) *isTensor = PETSC_FALSE;
18063f27d899SToby Isaac     if (endA) *endA = -1;
18073f27d899SToby Isaac     if (endB) *endB = -1;
18083f27d899SToby Isaac   }
18099566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCone(dm, p, &cone));
18103f27d899SToby Isaac   for (c = 0; c < coneSize; c++) {
18113f27d899SToby Isaac     PetscInt f = cone[c];
18123f27d899SToby Isaac     PetscInt fConeSize;
18133f27d899SToby Isaac 
18149566063dSJacob Faibussowitsch     PetscCall(DMPlexGetConeSize(dm, f, &fConeSize));
18153f27d899SToby Isaac     if (fConeSize != coneSize - 2) continue;
18163f27d899SToby Isaac 
18173f27d899SToby Isaac     for (c2 = c + 1; c2 < coneSize; c2++) {
18183f27d899SToby Isaac       PetscInt  f2 = cone[c2];
18193f27d899SToby Isaac       PetscBool isTensorff2;
18203f27d899SToby Isaac       PetscInt  f2ConeSize;
18213f27d899SToby Isaac 
18229566063dSJacob Faibussowitsch       PetscCall(DMPlexGetConeSize(dm, f2, &f2ConeSize));
18233f27d899SToby Isaac       if (f2ConeSize != coneSize - 2) continue;
18243f27d899SToby Isaac 
18259566063dSJacob Faibussowitsch       PetscCall(DMPlexPointIsTensor_Internal_Given(dm, p, f, f2, &isTensorff2));
18263f27d899SToby Isaac       if (isTensorff2) {
18273f27d899SToby Isaac         if (isTensor) *isTensor = PETSC_TRUE;
18283f27d899SToby Isaac         if (endA) *endA = f;
18293f27d899SToby Isaac         if (endB) *endB = f2;
18303f27d899SToby Isaac         PetscFunctionReturn(0);
18313f27d899SToby Isaac       }
18323f27d899SToby Isaac     }
18333f27d899SToby Isaac   }
18343f27d899SToby Isaac   if (isTensor) *isTensor = PETSC_FALSE;
18353f27d899SToby Isaac   if (endA) *endA = -1;
18363f27d899SToby Isaac   if (endB) *endB = -1;
18373f27d899SToby Isaac   PetscFunctionReturn(0);
18383f27d899SToby Isaac }
18393f27d899SToby Isaac 
184077f1a120SToby Isaac /* determine if a cell is a tensor with a segment by looping over pairs of facets to find a pair
184177f1a120SToby Isaac  * that could be the opposite ends */
18429371c9d4SSatish Balay static PetscErrorCode DMPlexPointIsTensor(DM dm, PetscInt p, PetscBool *isTensor, PetscInt *endA, PetscInt *endB) {
18433f27d899SToby Isaac   DMPlexInterpolatedFlag interpolated;
18443f27d899SToby Isaac 
18453f27d899SToby Isaac   PetscFunctionBegin;
18469566063dSJacob Faibussowitsch   PetscCall(DMPlexIsInterpolated(dm, &interpolated));
184708401ef6SPierre Jolivet   PetscCheck(interpolated == DMPLEX_INTERPOLATED_FULL, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Only for interpolated DMPlex's");
18489566063dSJacob Faibussowitsch   PetscCall(DMPlexPointIsTensor_Internal(dm, p, isTensor, endA, endB));
18493f27d899SToby Isaac   PetscFunctionReturn(0);
18503f27d899SToby Isaac }
18513f27d899SToby Isaac 
18528f28b7bfSToby Isaac /* Let k = formDegree and k' = -sign(k) * dim + k.  Transform a symmetric frame for k-forms on the biunit simplex into
18538f28b7bfSToby Isaac  * a symmetric frame for k'-forms on the biunit simplex.
18541f440fbeSToby Isaac  *
18558f28b7bfSToby Isaac  * A frame is "symmetric" if the pullback of every symmetry of the biunit simplex is a permutation of the frame.
18561f440fbeSToby Isaac  *
18578f28b7bfSToby Isaac  * forms in the symmetric frame are used as dofs in the untrimmed simplex spaces.  This way, symmetries of the
18588f28b7bfSToby Isaac  * reference cell result in permutations of dofs grouped by node.
18591f440fbeSToby Isaac  *
18608f28b7bfSToby Isaac  * Use T to transform dof matrices for k'-forms into dof matrices for k-forms as a block diagonal transformation on
18618f28b7bfSToby Isaac  * the right.
18621f440fbeSToby Isaac  */
18639371c9d4SSatish Balay static PetscErrorCode BiunitSimplexSymmetricFormTransformation(PetscInt dim, PetscInt formDegree, PetscReal T[]) {
18641f440fbeSToby Isaac   PetscInt   k  = formDegree;
18651f440fbeSToby Isaac   PetscInt   kd = k < 0 ? dim + k : k - dim;
18661f440fbeSToby Isaac   PetscInt   Nk;
18671f440fbeSToby Isaac   PetscReal *biToEq, *eqToBi, *biToEqStar, *eqToBiStar;
18681f440fbeSToby Isaac   PetscInt   fact;
18691f440fbeSToby Isaac 
18701f440fbeSToby Isaac   PetscFunctionBegin;
18719566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(k), &Nk));
18729566063dSJacob Faibussowitsch   PetscCall(PetscCalloc4(dim * dim, &biToEq, dim * dim, &eqToBi, Nk * Nk, &biToEqStar, Nk * Nk, &eqToBiStar));
18731f440fbeSToby Isaac   /* fill in biToEq: Jacobian of the transformation from the biunit simplex to the equilateral simplex */
18741f440fbeSToby Isaac   fact = 0;
18751f440fbeSToby Isaac   for (PetscInt i = 0; i < dim; i++) {
18761f440fbeSToby Isaac     biToEq[i * dim + i] = PetscSqrtReal(((PetscReal)i + 2.) / (2. * ((PetscReal)i + 1.)));
18771f440fbeSToby Isaac     fact += 4 * (i + 1);
1878ad540459SPierre Jolivet     for (PetscInt j = i + 1; j < dim; j++) biToEq[i * dim + j] = PetscSqrtReal(1. / (PetscReal)fact);
18791f440fbeSToby Isaac   }
18808f28b7bfSToby Isaac   /* fill in eqToBi: Jacobian of the transformation from the equilateral simplex to the biunit simplex */
18811f440fbeSToby Isaac   fact = 0;
18821f440fbeSToby Isaac   for (PetscInt j = 0; j < dim; j++) {
18831f440fbeSToby Isaac     eqToBi[j * dim + j] = PetscSqrtReal(2. * ((PetscReal)j + 1.) / ((PetscReal)j + 2));
18841f440fbeSToby Isaac     fact += j + 1;
1885ad540459SPierre Jolivet     for (PetscInt i = 0; i < j; i++) eqToBi[i * dim + j] = -PetscSqrtReal(1. / (PetscReal)fact);
18861f440fbeSToby Isaac   }
18879566063dSJacob Faibussowitsch   PetscCall(PetscDTAltVPullbackMatrix(dim, dim, biToEq, kd, biToEqStar));
18889566063dSJacob Faibussowitsch   PetscCall(PetscDTAltVPullbackMatrix(dim, dim, eqToBi, k, eqToBiStar));
18898f28b7bfSToby Isaac   /* product of pullbacks simulates the following steps
18908f28b7bfSToby Isaac    *
18918f28b7bfSToby Isaac    * 1. start with frame W = [w_1, w_2, ..., w_m] of k forms that is symmetric on the biunit simplex:
18928f28b7bfSToby Isaac           if J is the Jacobian of a symmetry of the biunit simplex, then J_k* W = [J_k*w_1, ..., J_k*w_m]
18938f28b7bfSToby Isaac           is a permutation of W.
18948f28b7bfSToby Isaac           Even though a k' form --- a (dim - k) form represented by its Hodge star --- has the same geometric
18958f28b7bfSToby Isaac           content as a k form, W is not a symmetric frame of k' forms on the biunit simplex.  That's because,
18968f28b7bfSToby Isaac           for general Jacobian J, J_k* != J_k'*.
18978f28b7bfSToby Isaac    * 2. pullback W to the equilateral triangle using the k pullback, W_eq = eqToBi_k* W.  All symmetries of the
18988f28b7bfSToby Isaac           equilateral simplex have orthonormal Jacobians.  For an orthonormal Jacobian O, J_k* = J_k'*, so W_eq is
18998f28b7bfSToby Isaac           also a symmetric frame for k' forms on the equilateral simplex.
19008f28b7bfSToby Isaac      3. pullback W_eq back to the biunit simplex using the k' pulback, V = biToEq_k'* W_eq = biToEq_k'* eqToBi_k* W.
19018f28b7bfSToby Isaac           V is a symmetric frame for k' forms on the biunit simplex.
19028f28b7bfSToby Isaac    */
19031f440fbeSToby Isaac   for (PetscInt i = 0; i < Nk; i++) {
19041f440fbeSToby Isaac     for (PetscInt j = 0; j < Nk; j++) {
19051f440fbeSToby Isaac       PetscReal val = 0.;
19061f440fbeSToby Isaac       for (PetscInt k = 0; k < Nk; k++) val += biToEqStar[i * Nk + k] * eqToBiStar[k * Nk + j];
19071f440fbeSToby Isaac       T[i * Nk + j] = val;
19081f440fbeSToby Isaac     }
19091f440fbeSToby Isaac   }
19109566063dSJacob Faibussowitsch   PetscCall(PetscFree4(biToEq, eqToBi, biToEqStar, eqToBiStar));
19111f440fbeSToby Isaac   PetscFunctionReturn(0);
19121f440fbeSToby Isaac }
19131f440fbeSToby Isaac 
191477f1a120SToby Isaac /* permute a quadrature -> dof matrix so that its rows are in revlex order by nodeIdx */
19159371c9d4SSatish Balay static PetscErrorCode MatPermuteByNodeIdx(Mat A, PetscLagNodeIndices ni, Mat *Aperm) {
19163f27d899SToby Isaac   PetscInt   m, n, i, j;
19173f27d899SToby Isaac   PetscInt   nodeIdxDim = ni->nodeIdxDim;
19183f27d899SToby Isaac   PetscInt   nodeVecDim = ni->nodeVecDim;
19193f27d899SToby Isaac   PetscInt  *perm;
19203f27d899SToby Isaac   IS         permIS;
19213f27d899SToby Isaac   IS         id;
19223f27d899SToby Isaac   PetscInt  *nIdxPerm;
19233f27d899SToby Isaac   PetscReal *nVecPerm;
19243f27d899SToby Isaac 
19253f27d899SToby Isaac   PetscFunctionBegin;
19269566063dSJacob Faibussowitsch   PetscCall(PetscLagNodeIndicesGetPermutation(ni, &perm));
19279566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &m, &n));
19289566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nodeIdxDim * m, &nIdxPerm));
19299566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nodeVecDim * m, &nVecPerm));
19309371c9d4SSatish Balay   for (i = 0; i < m; i++)
19319371c9d4SSatish Balay     for (j = 0; j < nodeIdxDim; j++) nIdxPerm[i * nodeIdxDim + j] = ni->nodeIdx[perm[i] * nodeIdxDim + j];
19329371c9d4SSatish Balay   for (i = 0; i < m; i++)
19339371c9d4SSatish Balay     for (j = 0; j < nodeVecDim; j++) nVecPerm[i * nodeVecDim + j] = ni->nodeVec[perm[i] * nodeVecDim + j];
19349566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, m, perm, PETSC_USE_POINTER, &permIS));
19359566063dSJacob Faibussowitsch   PetscCall(ISSetPermutation(permIS));
19369566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF, n, 0, 1, &id));
19379566063dSJacob Faibussowitsch   PetscCall(ISSetPermutation(id));
19389566063dSJacob Faibussowitsch   PetscCall(MatPermute(A, permIS, id, Aperm));
19399566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&permIS));
19409566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&id));
19413f27d899SToby Isaac   for (i = 0; i < m; i++) perm[i] = i;
19429566063dSJacob Faibussowitsch   PetscCall(PetscFree(ni->nodeIdx));
19439566063dSJacob Faibussowitsch   PetscCall(PetscFree(ni->nodeVec));
19443f27d899SToby Isaac   ni->nodeIdx = nIdxPerm;
19453f27d899SToby Isaac   ni->nodeVec = nVecPerm;
19466f905325SMatthew G. Knepley   PetscFunctionReturn(0);
19476f905325SMatthew G. Knepley }
19486f905325SMatthew G. Knepley 
19499371c9d4SSatish Balay static PetscErrorCode PetscDualSpaceSetUp_Lagrange(PetscDualSpace sp) {
19506f905325SMatthew G. Knepley   PetscDualSpace_Lag    *lag   = (PetscDualSpace_Lag *)sp->data;
19516f905325SMatthew G. Knepley   DM                     dm    = sp->dm;
19523f27d899SToby Isaac   DM                     dmint = NULL;
19533f27d899SToby Isaac   PetscInt               order;
19546f905325SMatthew G. Knepley   PetscInt               Nc = sp->Nc;
19556f905325SMatthew G. Knepley   MPI_Comm               comm;
19566f905325SMatthew G. Knepley   PetscBool              continuous;
19573f27d899SToby Isaac   PetscSection           section;
19583f27d899SToby Isaac   PetscInt               depth, dim, pStart, pEnd, cStart, cEnd, p, *pStratStart, *pStratEnd, d;
19593f27d899SToby Isaac   PetscInt               formDegree, Nk, Ncopies;
19603f27d899SToby Isaac   PetscInt               tensorf = -1, tensorf2 = -1;
19613f27d899SToby Isaac   PetscBool              tensorCell, tensorSpace;
19623f27d899SToby Isaac   PetscBool              uniform, trimmed;
19633f27d899SToby Isaac   Petsc1DNodeFamily      nodeFamily;
19643f27d899SToby Isaac   PetscInt               numNodeSkip;
19653f27d899SToby Isaac   DMPlexInterpolatedFlag interpolated;
19663f27d899SToby Isaac   PetscBool              isbdm;
19676f905325SMatthew G. Knepley 
19686f905325SMatthew G. Knepley   PetscFunctionBegin;
19693f27d899SToby Isaac   /* step 1: sanitize input */
19709566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)sp, &comm));
19719566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
19729566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)sp, PETSCDUALSPACEBDM, &isbdm));
19733f27d899SToby Isaac   if (isbdm) {
19743f27d899SToby Isaac     sp->k = -(dim - 1); /* form degree of H-div */
19759566063dSJacob Faibussowitsch     PetscCall(PetscObjectChangeTypeName((PetscObject)sp, PETSCDUALSPACELAGRANGE));
19763f27d899SToby Isaac   }
19779566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetFormDegree(sp, &formDegree));
197808401ef6SPierre Jolivet   PetscCheck(PetscAbsInt(formDegree) <= dim, comm, PETSC_ERR_ARG_OUTOFRANGE, "Form degree must be bounded by dimension");
19799566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(formDegree), &Nk));
19803f27d899SToby Isaac   if (sp->Nc <= 0 && lag->numCopies > 0) sp->Nc = Nk * lag->numCopies;
19813f27d899SToby Isaac   Nc = sp->Nc;
198208401ef6SPierre Jolivet   PetscCheck(Nc % Nk == 0, comm, PETSC_ERR_ARG_INCOMP, "Number of components is not a multiple of form degree size");
19833f27d899SToby Isaac   if (lag->numCopies <= 0) lag->numCopies = Nc / Nk;
19843f27d899SToby Isaac   Ncopies = lag->numCopies;
19851dca8a05SBarry Smith   PetscCheck(Nc / Nk == Ncopies, comm, PETSC_ERR_ARG_INCOMP, "Number of copies * (dim choose k) != Nc");
19863f27d899SToby Isaac   if (!dim) sp->order = 0;
19873f27d899SToby Isaac   order   = sp->order;
19883f27d899SToby Isaac   uniform = sp->uniform;
198928b400f6SJacob Faibussowitsch   PetscCheck(uniform, PETSC_COMM_SELF, PETSC_ERR_SUP, "Variable order not supported yet");
19903f27d899SToby Isaac   if (lag->trimmed && !formDegree) lag->trimmed = PETSC_FALSE; /* trimmed spaces are the same as full spaces for 0-forms */
19913f27d899SToby Isaac   if (lag->nodeType == PETSCDTNODES_DEFAULT) {
19923f27d899SToby Isaac     lag->nodeType     = PETSCDTNODES_GAUSSJACOBI;
19933f27d899SToby Isaac     lag->nodeExponent = 0.;
19943f27d899SToby Isaac     /* trimmed spaces don't include corner vertices, so don't use end nodes by default */
19953f27d899SToby Isaac     lag->endNodes     = lag->trimmed ? PETSC_FALSE : PETSC_TRUE;
19963f27d899SToby Isaac   }
19973f27d899SToby Isaac   /* If a trimmed space and the user did choose nodes with endpoints, skip them by default */
19983f27d899SToby Isaac   if (lag->numNodeSkip < 0) lag->numNodeSkip = (lag->trimmed && lag->endNodes) ? 1 : 0;
19993f27d899SToby Isaac   numNodeSkip = lag->numNodeSkip;
200008401ef6SPierre Jolivet   PetscCheck(!lag->trimmed || order, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot have zeroth order trimmed elements");
20013f27d899SToby Isaac   if (lag->trimmed && PetscAbsInt(formDegree) == dim) { /* convert trimmed n-forms to untrimmed of one polynomial order less */
20023f27d899SToby Isaac     sp->order--;
20033f27d899SToby Isaac     order--;
20043f27d899SToby Isaac     lag->trimmed = PETSC_FALSE;
20053f27d899SToby Isaac   }
20063f27d899SToby Isaac   trimmed = lag->trimmed;
20073f27d899SToby Isaac   if (!order || PetscAbsInt(formDegree) == dim) lag->continuous = PETSC_FALSE;
20083f27d899SToby Isaac   continuous = lag->continuous;
20099566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(dm, &depth));
20109566063dSJacob Faibussowitsch   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
20119566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
20121dca8a05SBarry Smith   PetscCheck(pStart == 0 && cStart == 0, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Expect DM with chart starting at zero and cells first");
201308401ef6SPierre Jolivet   PetscCheck(cEnd == 1, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Use PETSCDUALSPACEREFINED for multi-cell reference meshes");
20149566063dSJacob Faibussowitsch   PetscCall(DMPlexIsInterpolated(dm, &interpolated));
20153f27d899SToby Isaac   if (interpolated != DMPLEX_INTERPOLATED_FULL) {
20169566063dSJacob Faibussowitsch     PetscCall(DMPlexInterpolate(dm, &dmint));
20173f27d899SToby Isaac   } else {
20189566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)dm));
20193f27d899SToby Isaac     dmint = dm;
20203f27d899SToby Isaac   }
20213f27d899SToby Isaac   tensorCell = PETSC_FALSE;
202248a46eb9SPierre Jolivet   if (dim > 1) PetscCall(DMPlexPointIsTensor(dmint, 0, &tensorCell, &tensorf, &tensorf2));
20233f27d899SToby Isaac   lag->tensorCell = tensorCell;
20243f27d899SToby Isaac   if (dim < 2 || !lag->tensorCell) lag->tensorSpace = PETSC_FALSE;
20256f905325SMatthew G. Knepley   tensorSpace = lag->tensorSpace;
202648a46eb9SPierre Jolivet   if (!lag->nodeFamily) PetscCall(Petsc1DNodeFamilyCreate(lag->nodeType, lag->nodeExponent, lag->endNodes, &lag->nodeFamily));
20273f27d899SToby Isaac   nodeFamily = lag->nodeFamily;
20281dca8a05SBarry Smith   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");
202920cf1dd8SToby Isaac 
20303f27d899SToby Isaac   /* step 2: construct the boundary spaces */
20319566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(depth + 1, &pStratStart, depth + 1, &pStratEnd));
20329566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(pEnd, &(sp->pointSpaces)));
20339566063dSJacob Faibussowitsch   for (d = 0; d <= depth; ++d) PetscCall(DMPlexGetDepthStratum(dm, d, &pStratStart[d], &pStratEnd[d]));
20349566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSectionCreate_Internal(sp, &section));
20353f27d899SToby Isaac   sp->pointSection = section;
20363f27d899SToby Isaac   if (continuous && !(lag->interiorOnly)) {
20373f27d899SToby Isaac     PetscInt h;
20386f905325SMatthew G. Knepley 
20393f27d899SToby Isaac     for (p = pStratStart[depth - 1]; p < pStratEnd[depth - 1]; p++) { /* calculate the facet dual spaces */
20403f27d899SToby Isaac       PetscReal      v0[3];
20413f27d899SToby Isaac       DMPolytopeType ptype;
20423f27d899SToby Isaac       PetscReal      J[9], detJ;
20436f905325SMatthew G. Knepley       PetscInt       q;
20446f905325SMatthew G. Knepley 
20459566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryAffineFEM(dm, p, v0, J, NULL, &detJ));
20469566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCellType(dm, p, &ptype));
20476f905325SMatthew G. Knepley 
204877f1a120SToby Isaac       /* compare to previous facets: if computed, reference that dualspace */
20493f27d899SToby Isaac       for (q = pStratStart[depth - 1]; q < p; q++) {
20503f27d899SToby Isaac         DMPolytopeType qtype;
20516f905325SMatthew G. Knepley 
20529566063dSJacob Faibussowitsch         PetscCall(DMPlexGetCellType(dm, q, &qtype));
20533f27d899SToby Isaac         if (qtype == ptype) break;
20546f905325SMatthew G. Knepley       }
20553f27d899SToby Isaac       if (q < p) { /* this facet has the same dual space as that one */
20569566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)sp->pointSpaces[q]));
20573f27d899SToby Isaac         sp->pointSpaces[p] = sp->pointSpaces[q];
20583f27d899SToby Isaac         continue;
20596f905325SMatthew G. Knepley       }
20603f27d899SToby Isaac       /* if not, recursively compute this dual space */
20619566063dSJacob Faibussowitsch       PetscCall(PetscDualSpaceCreateFacetSubspace_Lagrange(sp, NULL, p, formDegree, Ncopies, PETSC_FALSE, &sp->pointSpaces[p]));
20626f905325SMatthew G. Knepley     }
20633f27d899SToby Isaac     for (h = 2; h <= depth; h++) { /* get the higher subspaces from the facet subspaces */
20643f27d899SToby Isaac       PetscInt hd   = depth - h;
20653f27d899SToby Isaac       PetscInt hdim = dim - h;
20666f905325SMatthew G. Knepley 
20673f27d899SToby Isaac       if (hdim < PetscAbsInt(formDegree)) break;
20683f27d899SToby Isaac       for (p = pStratStart[hd]; p < pStratEnd[hd]; p++) {
20693f27d899SToby Isaac         PetscInt        suppSize, s;
20703f27d899SToby Isaac         const PetscInt *supp;
20716f905325SMatthew G. Knepley 
20729566063dSJacob Faibussowitsch         PetscCall(DMPlexGetSupportSize(dm, p, &suppSize));
20739566063dSJacob Faibussowitsch         PetscCall(DMPlexGetSupport(dm, p, &supp));
20743f27d899SToby Isaac         for (s = 0; s < suppSize; s++) {
20753f27d899SToby Isaac           DM              qdm;
20763f27d899SToby Isaac           PetscDualSpace  qsp, psp;
20773f27d899SToby Isaac           PetscInt        c, coneSize, q;
20783f27d899SToby Isaac           const PetscInt *cone;
20793f27d899SToby Isaac           const PetscInt *refCone;
20806f905325SMatthew G. Knepley 
20813f27d899SToby Isaac           q   = supp[0];
20823f27d899SToby Isaac           qsp = sp->pointSpaces[q];
20839566063dSJacob Faibussowitsch           PetscCall(DMPlexGetConeSize(dm, q, &coneSize));
20849566063dSJacob Faibussowitsch           PetscCall(DMPlexGetCone(dm, q, &cone));
20859371c9d4SSatish Balay           for (c = 0; c < coneSize; c++)
20869371c9d4SSatish Balay             if (cone[c] == p) break;
208708401ef6SPierre Jolivet           PetscCheck(c != coneSize, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "cone/support mismatch");
20889566063dSJacob Faibussowitsch           PetscCall(PetscDualSpaceGetDM(qsp, &qdm));
20899566063dSJacob Faibussowitsch           PetscCall(DMPlexGetCone(qdm, 0, &refCone));
20903f27d899SToby Isaac           /* get the equivalent dual space from the support dual space */
20919566063dSJacob Faibussowitsch           PetscCall(PetscDualSpaceGetPointSubspace(qsp, refCone[c], &psp));
20923f27d899SToby Isaac           if (!s) {
20939566063dSJacob Faibussowitsch             PetscCall(PetscObjectReference((PetscObject)psp));
20943f27d899SToby Isaac             sp->pointSpaces[p] = psp;
20953f27d899SToby Isaac           }
20963f27d899SToby Isaac         }
20973f27d899SToby Isaac       }
20983f27d899SToby Isaac     }
20993f27d899SToby Isaac     for (p = 1; p < pEnd; p++) {
21003f27d899SToby Isaac       PetscInt pspdim;
21013f27d899SToby Isaac       if (!sp->pointSpaces[p]) continue;
21029566063dSJacob Faibussowitsch       PetscCall(PetscDualSpaceGetInteriorDimension(sp->pointSpaces[p], &pspdim));
21039566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetDof(section, p, pspdim));
21043f27d899SToby Isaac     }
21053f27d899SToby Isaac   }
21066f905325SMatthew G. Knepley 
21073f27d899SToby Isaac   if (Ncopies > 1) {
21083f27d899SToby Isaac     Mat                 intMatScalar, allMatScalar;
21093f27d899SToby Isaac     PetscDualSpace      scalarsp;
21103f27d899SToby Isaac     PetscDualSpace_Lag *scalarlag;
21113f27d899SToby Isaac 
21129566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceDuplicate(sp, &scalarsp));
211377f1a120SToby Isaac     /* Setting the number of components to Nk is a space with 1 copy of each k-form */
21149566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceSetNumComponents(scalarsp, Nk));
21159566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceSetUp(scalarsp));
21169566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetInteriorData(scalarsp, &(sp->intNodes), &intMatScalar));
21179566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)(sp->intNodes)));
21189566063dSJacob Faibussowitsch     if (intMatScalar) PetscCall(PetscDualSpaceLagrangeMatrixCreateCopies(intMatScalar, Nk, Ncopies, &(sp->intMat)));
21199566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetAllData(scalarsp, &(sp->allNodes), &allMatScalar));
21209566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)(sp->allNodes)));
21219566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceLagrangeMatrixCreateCopies(allMatScalar, Nk, Ncopies, &(sp->allMat)));
21223f27d899SToby Isaac     sp->spdim    = scalarsp->spdim * Ncopies;
21233f27d899SToby Isaac     sp->spintdim = scalarsp->spintdim * Ncopies;
21243f27d899SToby Isaac     scalarlag    = (PetscDualSpace_Lag *)scalarsp->data;
21259566063dSJacob Faibussowitsch     PetscCall(PetscLagNodeIndicesReference(scalarlag->vertIndices));
21263f27d899SToby Isaac     lag->vertIndices = scalarlag->vertIndices;
21279566063dSJacob Faibussowitsch     PetscCall(PetscLagNodeIndicesReference(scalarlag->intNodeIndices));
21283f27d899SToby Isaac     lag->intNodeIndices = scalarlag->intNodeIndices;
21299566063dSJacob Faibussowitsch     PetscCall(PetscLagNodeIndicesReference(scalarlag->allNodeIndices));
21303f27d899SToby Isaac     lag->allNodeIndices = scalarlag->allNodeIndices;
21319566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceDestroy(&scalarsp));
21329566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetDof(section, 0, sp->spintdim));
21339566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceSectionSetUp_Internal(sp, section));
21349566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceComputeFunctionalsFromAllData(sp));
21359566063dSJacob Faibussowitsch     PetscCall(PetscFree2(pStratStart, pStratEnd));
21369566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&dmint));
213720cf1dd8SToby Isaac     PetscFunctionReturn(0);
213820cf1dd8SToby Isaac   }
213920cf1dd8SToby Isaac 
21403f27d899SToby Isaac   if (trimmed && !continuous) {
21413f27d899SToby Isaac     /* the dofs of a trimmed space don't have a nice tensor/lattice structure:
21423f27d899SToby Isaac      * just construct the continuous dual space and copy all of the data over,
21433f27d899SToby Isaac      * allocating it all to the cell instead of splitting it up between the boundaries */
21443f27d899SToby Isaac     PetscDualSpace      spcont;
21453f27d899SToby Isaac     PetscInt            spdim, f;
21463f27d899SToby Isaac     PetscQuadrature     allNodes;
21473f27d899SToby Isaac     PetscDualSpace_Lag *lagc;
21483f27d899SToby Isaac     Mat                 allMat;
21493f27d899SToby Isaac 
21509566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceDuplicate(sp, &spcont));
21519566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceLagrangeSetContinuity(spcont, PETSC_TRUE));
21529566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceSetUp(spcont));
21539566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDimension(spcont, &spdim));
21543f27d899SToby Isaac     sp->spdim = sp->spintdim = spdim;
21559566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetDof(section, 0, spdim));
21569566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceSectionSetUp_Internal(sp, section));
21579566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(spdim, &(sp->functional)));
21583f27d899SToby Isaac     for (f = 0; f < spdim; f++) {
21593f27d899SToby Isaac       PetscQuadrature fn;
21603f27d899SToby Isaac 
21619566063dSJacob Faibussowitsch       PetscCall(PetscDualSpaceGetFunctional(spcont, f, &fn));
21629566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)fn));
21633f27d899SToby Isaac       sp->functional[f] = fn;
21643f27d899SToby Isaac     }
21659566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetAllData(spcont, &allNodes, &allMat));
21669566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)allNodes));
21679566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)allNodes));
21683f27d899SToby Isaac     sp->allNodes = sp->intNodes = allNodes;
21699566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)allMat));
21709566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)allMat));
21713f27d899SToby Isaac     sp->allMat = sp->intMat = allMat;
21723f27d899SToby Isaac     lagc                    = (PetscDualSpace_Lag *)spcont->data;
21739566063dSJacob Faibussowitsch     PetscCall(PetscLagNodeIndicesReference(lagc->vertIndices));
21743f27d899SToby Isaac     lag->vertIndices = lagc->vertIndices;
21759566063dSJacob Faibussowitsch     PetscCall(PetscLagNodeIndicesReference(lagc->allNodeIndices));
21769566063dSJacob Faibussowitsch     PetscCall(PetscLagNodeIndicesReference(lagc->allNodeIndices));
21773f27d899SToby Isaac     lag->intNodeIndices = lagc->allNodeIndices;
21783f27d899SToby Isaac     lag->allNodeIndices = lagc->allNodeIndices;
21799566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceDestroy(&spcont));
21809566063dSJacob Faibussowitsch     PetscCall(PetscFree2(pStratStart, pStratEnd));
21819566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&dmint));
21823f27d899SToby Isaac     PetscFunctionReturn(0);
21833f27d899SToby Isaac   }
21843f27d899SToby Isaac 
21853f27d899SToby Isaac   /* step 3: construct intNodes, and intMat, and combine it with boundray data to make allNodes and allMat */
21863f27d899SToby Isaac   if (!tensorSpace) {
21879566063dSJacob Faibussowitsch     if (!tensorCell) PetscCall(PetscLagNodeIndicesCreateSimplexVertices(dm, &(lag->vertIndices)));
21883f27d899SToby Isaac 
21893f27d899SToby Isaac     if (trimmed) {
219077f1a120SToby Isaac       /* there is one dof in the interior of the a trimmed element for each full polynomial of with degree at most
219177f1a120SToby Isaac        * order + k - dim - 1 */
21923f27d899SToby Isaac       if (order + PetscAbsInt(formDegree) > dim) {
21933f27d899SToby Isaac         PetscInt sum = order + PetscAbsInt(formDegree) - dim - 1;
21943f27d899SToby Isaac         PetscInt nDofs;
21953f27d899SToby Isaac 
21969566063dSJacob Faibussowitsch         PetscCall(PetscDualSpaceLagrangeCreateSimplexNodeMat(nodeFamily, dim, sum, Nk, numNodeSkip, &sp->intNodes, &sp->intMat, &(lag->intNodeIndices)));
21979566063dSJacob Faibussowitsch         PetscCall(MatGetSize(sp->intMat, &nDofs, NULL));
21989566063dSJacob Faibussowitsch         PetscCall(PetscSectionSetDof(section, 0, nDofs));
21993f27d899SToby Isaac       }
22009566063dSJacob Faibussowitsch       PetscCall(PetscDualSpaceSectionSetUp_Internal(sp, section));
22019566063dSJacob Faibussowitsch       PetscCall(PetscDualSpaceCreateAllDataFromInteriorData(sp));
22029566063dSJacob Faibussowitsch       PetscCall(PetscDualSpaceLagrangeCreateAllNodeIdx(sp));
22033f27d899SToby Isaac     } else {
22043f27d899SToby Isaac       if (!continuous) {
220577f1a120SToby Isaac         /* if discontinuous just construct one node for each set of dofs (a set of dofs is a basis for the k-form
220677f1a120SToby Isaac          * space) */
22073f27d899SToby Isaac         PetscInt sum = order;
22083f27d899SToby Isaac         PetscInt nDofs;
22093f27d899SToby Isaac 
22109566063dSJacob Faibussowitsch         PetscCall(PetscDualSpaceLagrangeCreateSimplexNodeMat(nodeFamily, dim, sum, Nk, numNodeSkip, &sp->intNodes, &sp->intMat, &(lag->intNodeIndices)));
22119566063dSJacob Faibussowitsch         PetscCall(MatGetSize(sp->intMat, &nDofs, NULL));
22129566063dSJacob Faibussowitsch         PetscCall(PetscSectionSetDof(section, 0, nDofs));
22139566063dSJacob Faibussowitsch         PetscCall(PetscDualSpaceSectionSetUp_Internal(sp, section));
22149566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)(sp->intNodes)));
22153f27d899SToby Isaac         sp->allNodes = sp->intNodes;
22169566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)(sp->intMat)));
22173f27d899SToby Isaac         sp->allMat = sp->intMat;
22189566063dSJacob Faibussowitsch         PetscCall(PetscLagNodeIndicesReference(lag->intNodeIndices));
22193f27d899SToby Isaac         lag->allNodeIndices = lag->intNodeIndices;
22203f27d899SToby Isaac       } else {
222177f1a120SToby Isaac         /* there is one dof in the interior of the a full element for each trimmed polynomial of with degree at most
222277f1a120SToby Isaac          * order + k - dim, but with complementary form degree */
22233f27d899SToby Isaac         if (order + PetscAbsInt(formDegree) > dim) {
22243f27d899SToby Isaac           PetscDualSpace      trimmedsp;
22253f27d899SToby Isaac           PetscDualSpace_Lag *trimmedlag;
22263f27d899SToby Isaac           PetscQuadrature     intNodes;
22273f27d899SToby Isaac           PetscInt            trFormDegree = formDegree >= 0 ? formDegree - dim : dim - PetscAbsInt(formDegree);
22283f27d899SToby Isaac           PetscInt            nDofs;
22293f27d899SToby Isaac           Mat                 intMat;
22303f27d899SToby Isaac 
22319566063dSJacob Faibussowitsch           PetscCall(PetscDualSpaceDuplicate(sp, &trimmedsp));
22329566063dSJacob Faibussowitsch           PetscCall(PetscDualSpaceLagrangeSetTrimmed(trimmedsp, PETSC_TRUE));
22339566063dSJacob Faibussowitsch           PetscCall(PetscDualSpaceSetOrder(trimmedsp, order + PetscAbsInt(formDegree) - dim));
22349566063dSJacob Faibussowitsch           PetscCall(PetscDualSpaceSetFormDegree(trimmedsp, trFormDegree));
22353f27d899SToby Isaac           trimmedlag              = (PetscDualSpace_Lag *)trimmedsp->data;
22363f27d899SToby Isaac           trimmedlag->numNodeSkip = numNodeSkip + 1;
22379566063dSJacob Faibussowitsch           PetscCall(PetscDualSpaceSetUp(trimmedsp));
22389566063dSJacob Faibussowitsch           PetscCall(PetscDualSpaceGetAllData(trimmedsp, &intNodes, &intMat));
22399566063dSJacob Faibussowitsch           PetscCall(PetscObjectReference((PetscObject)intNodes));
22403f27d899SToby Isaac           sp->intNodes = intNodes;
22419566063dSJacob Faibussowitsch           PetscCall(PetscLagNodeIndicesReference(trimmedlag->allNodeIndices));
22423f27d899SToby Isaac           lag->intNodeIndices = trimmedlag->allNodeIndices;
22439566063dSJacob Faibussowitsch           PetscCall(PetscObjectReference((PetscObject)intMat));
22441f440fbeSToby Isaac           if (PetscAbsInt(formDegree) > 0 && PetscAbsInt(formDegree) < dim) {
22451f440fbeSToby Isaac             PetscReal   *T;
22461f440fbeSToby Isaac             PetscScalar *work;
22471f440fbeSToby Isaac             PetscInt     nCols, nRows;
22481f440fbeSToby Isaac             Mat          intMatT;
22491f440fbeSToby Isaac 
22509566063dSJacob Faibussowitsch             PetscCall(MatDuplicate(intMat, MAT_COPY_VALUES, &intMatT));
22519566063dSJacob Faibussowitsch             PetscCall(MatGetSize(intMat, &nRows, &nCols));
22529566063dSJacob Faibussowitsch             PetscCall(PetscMalloc2(Nk * Nk, &T, nCols, &work));
22539566063dSJacob Faibussowitsch             PetscCall(BiunitSimplexSymmetricFormTransformation(dim, formDegree, T));
22541f440fbeSToby Isaac             for (PetscInt row = 0; row < nRows; row++) {
22551f440fbeSToby Isaac               PetscInt           nrCols;
22561f440fbeSToby Isaac               const PetscInt    *rCols;
22571f440fbeSToby Isaac               const PetscScalar *rVals;
22581f440fbeSToby Isaac 
22599566063dSJacob Faibussowitsch               PetscCall(MatGetRow(intMat, row, &nrCols, &rCols, &rVals));
226008401ef6SPierre Jolivet               PetscCheck(nrCols % Nk == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "nonzeros in intMat matrix are not in k-form size blocks");
22611f440fbeSToby Isaac               for (PetscInt b = 0; b < nrCols; b += Nk) {
22621f440fbeSToby Isaac                 const PetscScalar *v = &rVals[b];
22631f440fbeSToby Isaac                 PetscScalar       *w = &work[b];
22641f440fbeSToby Isaac                 for (PetscInt j = 0; j < Nk; j++) {
22651f440fbeSToby Isaac                   w[j] = 0.;
2266ad540459SPierre Jolivet                   for (PetscInt i = 0; i < Nk; i++) w[j] += v[i] * T[i * Nk + j];
22671f440fbeSToby Isaac                 }
22681f440fbeSToby Isaac               }
22699566063dSJacob Faibussowitsch               PetscCall(MatSetValuesBlocked(intMatT, 1, &row, nrCols, rCols, work, INSERT_VALUES));
22709566063dSJacob Faibussowitsch               PetscCall(MatRestoreRow(intMat, row, &nrCols, &rCols, &rVals));
22711f440fbeSToby Isaac             }
22729566063dSJacob Faibussowitsch             PetscCall(MatAssemblyBegin(intMatT, MAT_FINAL_ASSEMBLY));
22739566063dSJacob Faibussowitsch             PetscCall(MatAssemblyEnd(intMatT, MAT_FINAL_ASSEMBLY));
22749566063dSJacob Faibussowitsch             PetscCall(MatDestroy(&intMat));
22751f440fbeSToby Isaac             intMat = intMatT;
22769566063dSJacob Faibussowitsch             PetscCall(PetscLagNodeIndicesDestroy(&(lag->intNodeIndices)));
22779566063dSJacob Faibussowitsch             PetscCall(PetscLagNodeIndicesDuplicate(trimmedlag->allNodeIndices, &(lag->intNodeIndices)));
22781f440fbeSToby Isaac             {
22791f440fbeSToby Isaac               PetscInt         nNodes     = lag->intNodeIndices->nNodes;
22801f440fbeSToby Isaac               PetscReal       *newNodeVec = lag->intNodeIndices->nodeVec;
22811f440fbeSToby Isaac               const PetscReal *oldNodeVec = trimmedlag->allNodeIndices->nodeVec;
22821f440fbeSToby Isaac 
22831f440fbeSToby Isaac               for (PetscInt n = 0; n < nNodes; n++) {
22841f440fbeSToby Isaac                 PetscReal       *w = &newNodeVec[n * Nk];
22851f440fbeSToby Isaac                 const PetscReal *v = &oldNodeVec[n * Nk];
22861f440fbeSToby Isaac 
22871f440fbeSToby Isaac                 for (PetscInt j = 0; j < Nk; j++) {
22881f440fbeSToby Isaac                   w[j] = 0.;
2289ad540459SPierre Jolivet                   for (PetscInt i = 0; i < Nk; i++) w[j] += v[i] * T[i * Nk + j];
22901f440fbeSToby Isaac                 }
22911f440fbeSToby Isaac               }
22921f440fbeSToby Isaac             }
22939566063dSJacob Faibussowitsch             PetscCall(PetscFree2(T, work));
22941f440fbeSToby Isaac           }
22951f440fbeSToby Isaac           sp->intMat = intMat;
22969566063dSJacob Faibussowitsch           PetscCall(MatGetSize(sp->intMat, &nDofs, NULL));
22979566063dSJacob Faibussowitsch           PetscCall(PetscDualSpaceDestroy(&trimmedsp));
22989566063dSJacob Faibussowitsch           PetscCall(PetscSectionSetDof(section, 0, nDofs));
22993f27d899SToby Isaac         }
23009566063dSJacob Faibussowitsch         PetscCall(PetscDualSpaceSectionSetUp_Internal(sp, section));
23019566063dSJacob Faibussowitsch         PetscCall(PetscDualSpaceCreateAllDataFromInteriorData(sp));
23029566063dSJacob Faibussowitsch         PetscCall(PetscDualSpaceLagrangeCreateAllNodeIdx(sp));
23033f27d899SToby Isaac       }
23043f27d899SToby Isaac     }
23053f27d899SToby Isaac   } else {
23063f27d899SToby Isaac     PetscQuadrature     intNodesTrace  = NULL;
23073f27d899SToby Isaac     PetscQuadrature     intNodesFiber  = NULL;
23083f27d899SToby Isaac     PetscQuadrature     intNodes       = NULL;
23093f27d899SToby Isaac     PetscLagNodeIndices intNodeIndices = NULL;
23103f27d899SToby Isaac     Mat                 intMat         = NULL;
23113f27d899SToby Isaac 
231277f1a120SToby Isaac     if (PetscAbsInt(formDegree) < dim) { /* get the trace k-forms on the first facet, and the 0-forms on the edge,
231377f1a120SToby Isaac                                             and wedge them together to create some of the k-form dofs */
23143f27d899SToby Isaac       PetscDualSpace      trace, fiber;
23153f27d899SToby Isaac       PetscDualSpace_Lag *tracel, *fiberl;
23163f27d899SToby Isaac       Mat                 intMatTrace, intMatFiber;
23173f27d899SToby Isaac 
23183f27d899SToby Isaac       if (sp->pointSpaces[tensorf]) {
23199566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)(sp->pointSpaces[tensorf])));
23203f27d899SToby Isaac         trace = sp->pointSpaces[tensorf];
23213f27d899SToby Isaac       } else {
23229566063dSJacob Faibussowitsch         PetscCall(PetscDualSpaceCreateFacetSubspace_Lagrange(sp, NULL, tensorf, formDegree, Ncopies, PETSC_TRUE, &trace));
23233f27d899SToby Isaac       }
23249566063dSJacob Faibussowitsch       PetscCall(PetscDualSpaceCreateEdgeSubspace_Lagrange(sp, order, 0, 1, PETSC_TRUE, &fiber));
23253f27d899SToby Isaac       tracel = (PetscDualSpace_Lag *)trace->data;
23263f27d899SToby Isaac       fiberl = (PetscDualSpace_Lag *)fiber->data;
23279566063dSJacob Faibussowitsch       PetscCall(PetscLagNodeIndicesCreateTensorVertices(dm, tracel->vertIndices, &(lag->vertIndices)));
23289566063dSJacob Faibussowitsch       PetscCall(PetscDualSpaceGetInteriorData(trace, &intNodesTrace, &intMatTrace));
23299566063dSJacob Faibussowitsch       PetscCall(PetscDualSpaceGetInteriorData(fiber, &intNodesFiber, &intMatFiber));
23303f27d899SToby Isaac       if (intNodesTrace && intNodesFiber) {
23319566063dSJacob Faibussowitsch         PetscCall(PetscQuadratureCreateTensor(intNodesTrace, intNodesFiber, &intNodes));
23329566063dSJacob Faibussowitsch         PetscCall(MatTensorAltV(intMatTrace, intMatFiber, dim - 1, formDegree, 1, 0, &intMat));
23339566063dSJacob Faibussowitsch         PetscCall(PetscLagNodeIndicesTensor(tracel->intNodeIndices, dim - 1, formDegree, fiberl->intNodeIndices, 1, 0, &intNodeIndices));
23343f27d899SToby Isaac       }
23359566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)intNodesTrace));
23369566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)intNodesFiber));
23379566063dSJacob Faibussowitsch       PetscCall(PetscDualSpaceDestroy(&fiber));
23389566063dSJacob Faibussowitsch       PetscCall(PetscDualSpaceDestroy(&trace));
23393f27d899SToby Isaac     }
234077f1a120SToby Isaac     if (PetscAbsInt(formDegree) > 0) { /* get the trace (k-1)-forms on the first facet, and the 1-forms on the edge,
234177f1a120SToby Isaac                                           and wedge them together to create the remaining k-form dofs */
23423f27d899SToby Isaac       PetscDualSpace      trace, fiber;
23433f27d899SToby Isaac       PetscDualSpace_Lag *tracel, *fiberl;
23443f27d899SToby Isaac       PetscQuadrature     intNodesTrace2, intNodesFiber2, intNodes2;
23453f27d899SToby Isaac       PetscLagNodeIndices intNodeIndices2;
23463f27d899SToby Isaac       Mat                 intMatTrace, intMatFiber, intMat2;
23473f27d899SToby Isaac       PetscInt            traceDegree = formDegree > 0 ? formDegree - 1 : formDegree + 1;
23483f27d899SToby Isaac       PetscInt            fiberDegree = formDegree > 0 ? 1 : -1;
23493f27d899SToby Isaac 
23509566063dSJacob Faibussowitsch       PetscCall(PetscDualSpaceCreateFacetSubspace_Lagrange(sp, NULL, tensorf, traceDegree, Ncopies, PETSC_TRUE, &trace));
23519566063dSJacob Faibussowitsch       PetscCall(PetscDualSpaceCreateEdgeSubspace_Lagrange(sp, order, fiberDegree, 1, PETSC_TRUE, &fiber));
23523f27d899SToby Isaac       tracel = (PetscDualSpace_Lag *)trace->data;
23533f27d899SToby Isaac       fiberl = (PetscDualSpace_Lag *)fiber->data;
235448a46eb9SPierre Jolivet       if (!lag->vertIndices) PetscCall(PetscLagNodeIndicesCreateTensorVertices(dm, tracel->vertIndices, &(lag->vertIndices)));
23559566063dSJacob Faibussowitsch       PetscCall(PetscDualSpaceGetInteriorData(trace, &intNodesTrace2, &intMatTrace));
23569566063dSJacob Faibussowitsch       PetscCall(PetscDualSpaceGetInteriorData(fiber, &intNodesFiber2, &intMatFiber));
23573f27d899SToby Isaac       if (intNodesTrace2 && intNodesFiber2) {
23589566063dSJacob Faibussowitsch         PetscCall(PetscQuadratureCreateTensor(intNodesTrace2, intNodesFiber2, &intNodes2));
23599566063dSJacob Faibussowitsch         PetscCall(MatTensorAltV(intMatTrace, intMatFiber, dim - 1, traceDegree, 1, fiberDegree, &intMat2));
23609566063dSJacob Faibussowitsch         PetscCall(PetscLagNodeIndicesTensor(tracel->intNodeIndices, dim - 1, traceDegree, fiberl->intNodeIndices, 1, fiberDegree, &intNodeIndices2));
23613f27d899SToby Isaac         if (!intMat) {
23623f27d899SToby Isaac           intMat         = intMat2;
23633f27d899SToby Isaac           intNodes       = intNodes2;
23643f27d899SToby Isaac           intNodeIndices = intNodeIndices2;
23653f27d899SToby Isaac         } else {
236677f1a120SToby Isaac           /* merge the matrices, quadrature points, and nodes */
23673f27d899SToby Isaac           PetscInt            nM;
23683f27d899SToby Isaac           PetscInt            nDof, nDof2;
23696ff15688SToby Isaac           PetscInt           *toMerged = NULL, *toMerged2 = NULL;
23706ff15688SToby Isaac           PetscQuadrature     merged               = NULL;
23713f27d899SToby Isaac           PetscLagNodeIndices intNodeIndicesMerged = NULL;
23723f27d899SToby Isaac           Mat                 matMerged            = NULL;
23733f27d899SToby Isaac 
23749566063dSJacob Faibussowitsch           PetscCall(MatGetSize(intMat, &nDof, NULL));
23759566063dSJacob Faibussowitsch           PetscCall(MatGetSize(intMat2, &nDof2, NULL));
23769566063dSJacob Faibussowitsch           PetscCall(PetscQuadraturePointsMerge(intNodes, intNodes2, &merged, &toMerged, &toMerged2));
23779566063dSJacob Faibussowitsch           PetscCall(PetscQuadratureGetData(merged, NULL, NULL, &nM, NULL, NULL));
23789566063dSJacob Faibussowitsch           PetscCall(MatricesMerge(intMat, intMat2, dim, formDegree, nM, toMerged, toMerged2, &matMerged));
23799566063dSJacob Faibussowitsch           PetscCall(PetscLagNodeIndicesMerge(intNodeIndices, intNodeIndices2, &intNodeIndicesMerged));
23809566063dSJacob Faibussowitsch           PetscCall(PetscFree(toMerged));
23819566063dSJacob Faibussowitsch           PetscCall(PetscFree(toMerged2));
23829566063dSJacob Faibussowitsch           PetscCall(MatDestroy(&intMat));
23839566063dSJacob Faibussowitsch           PetscCall(MatDestroy(&intMat2));
23849566063dSJacob Faibussowitsch           PetscCall(PetscQuadratureDestroy(&intNodes));
23859566063dSJacob Faibussowitsch           PetscCall(PetscQuadratureDestroy(&intNodes2));
23869566063dSJacob Faibussowitsch           PetscCall(PetscLagNodeIndicesDestroy(&intNodeIndices));
23879566063dSJacob Faibussowitsch           PetscCall(PetscLagNodeIndicesDestroy(&intNodeIndices2));
23883f27d899SToby Isaac           intNodes       = merged;
23893f27d899SToby Isaac           intMat         = matMerged;
23903f27d899SToby Isaac           intNodeIndices = intNodeIndicesMerged;
23913f27d899SToby Isaac           if (!trimmed) {
239277f1a120SToby Isaac             /* I think users expect that, when a node has a full basis for the k-forms,
239377f1a120SToby Isaac              * they should be consecutive dofs.  That isn't the case for trimmed spaces,
239477f1a120SToby Isaac              * but is for some of the nodes in untrimmed spaces, so in that case we
239577f1a120SToby Isaac              * sort them to group them by node */
23963f27d899SToby Isaac             Mat intMatPerm;
23973f27d899SToby Isaac 
23989566063dSJacob Faibussowitsch             PetscCall(MatPermuteByNodeIdx(intMat, intNodeIndices, &intMatPerm));
23999566063dSJacob Faibussowitsch             PetscCall(MatDestroy(&intMat));
24003f27d899SToby Isaac             intMat = intMatPerm;
24013f27d899SToby Isaac           }
24023f27d899SToby Isaac         }
24033f27d899SToby Isaac       }
24049566063dSJacob Faibussowitsch       PetscCall(PetscDualSpaceDestroy(&fiber));
24059566063dSJacob Faibussowitsch       PetscCall(PetscDualSpaceDestroy(&trace));
24063f27d899SToby Isaac     }
24079566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureDestroy(&intNodesTrace));
24089566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureDestroy(&intNodesFiber));
24093f27d899SToby Isaac     sp->intNodes        = intNodes;
24103f27d899SToby Isaac     sp->intMat          = intMat;
24113f27d899SToby Isaac     lag->intNodeIndices = intNodeIndices;
24126f905325SMatthew G. Knepley     {
24133f27d899SToby Isaac       PetscInt nDofs = 0;
24143f27d899SToby Isaac 
241548a46eb9SPierre Jolivet       if (intMat) PetscCall(MatGetSize(intMat, &nDofs, NULL));
24169566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetDof(section, 0, nDofs));
24173f27d899SToby Isaac     }
24189566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceSectionSetUp_Internal(sp, section));
24193f27d899SToby Isaac     if (continuous) {
24209566063dSJacob Faibussowitsch       PetscCall(PetscDualSpaceCreateAllDataFromInteriorData(sp));
24219566063dSJacob Faibussowitsch       PetscCall(PetscDualSpaceLagrangeCreateAllNodeIdx(sp));
24223f27d899SToby Isaac     } else {
24239566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)intNodes));
24243f27d899SToby Isaac       sp->allNodes = intNodes;
24259566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)intMat));
24263f27d899SToby Isaac       sp->allMat = intMat;
24279566063dSJacob Faibussowitsch       PetscCall(PetscLagNodeIndicesReference(intNodeIndices));
24283f27d899SToby Isaac       lag->allNodeIndices = intNodeIndices;
24293f27d899SToby Isaac     }
24303f27d899SToby Isaac   }
24319566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(section, &sp->spdim));
24329566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetConstrainedStorageSize(section, &sp->spintdim));
24339566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceComputeFunctionalsFromAllData(sp));
24349566063dSJacob Faibussowitsch   PetscCall(PetscFree2(pStratStart, pStratEnd));
24359566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmint));
24363f27d899SToby Isaac   PetscFunctionReturn(0);
24373f27d899SToby Isaac }
24383f27d899SToby Isaac 
243977f1a120SToby Isaac /* Create a matrix that represents the transformation that DMPlexVecGetClosure() would need
244077f1a120SToby Isaac  * to get the representation of the dofs for a mesh point if the mesh point had this orientation
244177f1a120SToby Isaac  * relative to the cell */
24429371c9d4SSatish Balay PetscErrorCode PetscDualSpaceCreateInteriorSymmetryMatrix_Lagrange(PetscDualSpace sp, PetscInt ornt, Mat *symMat) {
24433f27d899SToby Isaac   PetscDualSpace_Lag *lag;
24443f27d899SToby Isaac   DM                  dm;
24453f27d899SToby Isaac   PetscLagNodeIndices vertIndices, intNodeIndices;
24463f27d899SToby Isaac   PetscLagNodeIndices ni;
24473f27d899SToby Isaac   PetscInt            nodeIdxDim, nodeVecDim, nNodes;
24483f27d899SToby Isaac   PetscInt            formDegree;
24493f27d899SToby Isaac   PetscInt           *perm, *permOrnt;
24503f27d899SToby Isaac   PetscInt           *nnz;
24513f27d899SToby Isaac   PetscInt            n;
24523f27d899SToby Isaac   PetscInt            maxGroupSize;
24533f27d899SToby Isaac   PetscScalar        *V, *W, *work;
24543f27d899SToby Isaac   Mat                 A;
24556f905325SMatthew G. Knepley 
24566f905325SMatthew G. Knepley   PetscFunctionBegin;
24573f27d899SToby Isaac   if (!sp->spintdim) {
24583f27d899SToby Isaac     *symMat = NULL;
24593f27d899SToby Isaac     PetscFunctionReturn(0);
24606f905325SMatthew G. Knepley   }
24613f27d899SToby Isaac   lag            = (PetscDualSpace_Lag *)sp->data;
24623f27d899SToby Isaac   vertIndices    = lag->vertIndices;
24633f27d899SToby Isaac   intNodeIndices = lag->intNodeIndices;
24649566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDM(sp, &dm));
24659566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetFormDegree(sp, &formDegree));
24669566063dSJacob Faibussowitsch   PetscCall(PetscNew(&ni));
24673f27d899SToby Isaac   ni->refct      = 1;
24683f27d899SToby Isaac   ni->nodeIdxDim = nodeIdxDim = intNodeIndices->nodeIdxDim;
24693f27d899SToby Isaac   ni->nodeVecDim = nodeVecDim = intNodeIndices->nodeVecDim;
24703f27d899SToby Isaac   ni->nNodes = nNodes = intNodeIndices->nNodes;
24719566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nNodes * nodeIdxDim, &(ni->nodeIdx)));
24729566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nNodes * nodeVecDim, &(ni->nodeVec)));
247377f1a120SToby Isaac   /* push forward the dofs by the symmetry of the reference element induced by ornt */
24749566063dSJacob Faibussowitsch   PetscCall(PetscLagNodeIndicesPushForward(dm, vertIndices, 0, vertIndices, intNodeIndices, ornt, formDegree, ni->nodeIdx, ni->nodeVec));
247577f1a120SToby Isaac   /* get the revlex order for both the original and transformed dofs */
24769566063dSJacob Faibussowitsch   PetscCall(PetscLagNodeIndicesGetPermutation(intNodeIndices, &perm));
24779566063dSJacob Faibussowitsch   PetscCall(PetscLagNodeIndicesGetPermutation(ni, &permOrnt));
24789566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nNodes, &nnz));
24793f27d899SToby Isaac   for (n = 0, maxGroupSize = 0; n < nNodes;) { /* incremented in the loop */
24803f27d899SToby Isaac     PetscInt *nind = &(ni->nodeIdx[permOrnt[n] * nodeIdxDim]);
24813f27d899SToby Isaac     PetscInt  m, nEnd;
24823f27d899SToby Isaac     PetscInt  groupSize;
248377f1a120SToby Isaac     /* for each group of dofs that have the same nodeIdx coordinate */
24843f27d899SToby Isaac     for (nEnd = n + 1; nEnd < nNodes; nEnd++) {
24853f27d899SToby Isaac       PetscInt *mind = &(ni->nodeIdx[permOrnt[nEnd] * nodeIdxDim]);
24863f27d899SToby Isaac       PetscInt  d;
24873f27d899SToby Isaac 
24883f27d899SToby Isaac       /* compare the oriented permutation indices */
24899371c9d4SSatish Balay       for (d = 0; d < nodeIdxDim; d++)
24909371c9d4SSatish Balay         if (mind[d] != nind[d]) break;
24913f27d899SToby Isaac       if (d < nodeIdxDim) break;
24923f27d899SToby Isaac     }
249377f1a120SToby Isaac     /* permOrnt[[n, nEnd)] is a group of dofs that, under the symmetry are at the same location */
249476bd3646SJed Brown 
249577f1a120SToby Isaac     /* the symmetry had better map the group of dofs with the same permuted nodeIdx
249677f1a120SToby Isaac      * to a group of dofs with the same size, otherwise we messed up */
249776bd3646SJed Brown     if (PetscDefined(USE_DEBUG)) {
24983f27d899SToby Isaac       PetscInt  m;
24993f27d899SToby Isaac       PetscInt *nind = &(intNodeIndices->nodeIdx[perm[n] * nodeIdxDim]);
25003f27d899SToby Isaac 
25013f27d899SToby Isaac       for (m = n + 1; m < nEnd; m++) {
25023f27d899SToby Isaac         PetscInt *mind = &(intNodeIndices->nodeIdx[perm[m] * nodeIdxDim]);
25033f27d899SToby Isaac         PetscInt  d;
25043f27d899SToby Isaac 
25053f27d899SToby Isaac         /* compare the oriented permutation indices */
25069371c9d4SSatish Balay         for (d = 0; d < nodeIdxDim; d++)
25079371c9d4SSatish Balay           if (mind[d] != nind[d]) break;
25083f27d899SToby Isaac         if (d < nodeIdxDim) break;
25093f27d899SToby Isaac       }
251008401ef6SPierre Jolivet       PetscCheck(m >= nEnd, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Dofs with same index after symmetry not same block size");
25113f27d899SToby Isaac     }
25123f27d899SToby Isaac     groupSize = nEnd - n;
251377f1a120SToby Isaac     /* each pushforward dof vector will be expressed in a basis of the unpermuted dofs */
25143f27d899SToby Isaac     for (m = n; m < nEnd; m++) nnz[permOrnt[m]] = groupSize;
25153f27d899SToby Isaac 
25163f27d899SToby Isaac     maxGroupSize = PetscMax(maxGroupSize, nEnd - n);
25173f27d899SToby Isaac     n            = nEnd;
25183f27d899SToby Isaac   }
251908401ef6SPierre Jolivet   PetscCheck(maxGroupSize <= nodeVecDim, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Dofs are not in blocks that can be solved");
25209566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, nNodes, nNodes, 0, nnz, &A));
25219566063dSJacob Faibussowitsch   PetscCall(PetscFree(nnz));
25229566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(maxGroupSize * nodeVecDim, &V, maxGroupSize * nodeVecDim, &W, nodeVecDim * 2, &work));
25233f27d899SToby Isaac   for (n = 0; n < nNodes;) { /* incremented in the loop */
25243f27d899SToby Isaac     PetscInt *nind = &(ni->nodeIdx[permOrnt[n] * nodeIdxDim]);
25253f27d899SToby Isaac     PetscInt  nEnd;
25263f27d899SToby Isaac     PetscInt  m;
25273f27d899SToby Isaac     PetscInt  groupSize;
25283f27d899SToby Isaac     for (nEnd = n + 1; nEnd < nNodes; nEnd++) {
25293f27d899SToby Isaac       PetscInt *mind = &(ni->nodeIdx[permOrnt[nEnd] * nodeIdxDim]);
25303f27d899SToby Isaac       PetscInt  d;
25313f27d899SToby Isaac 
25323f27d899SToby Isaac       /* compare the oriented permutation indices */
25339371c9d4SSatish Balay       for (d = 0; d < nodeIdxDim; d++)
25349371c9d4SSatish Balay         if (mind[d] != nind[d]) break;
25353f27d899SToby Isaac       if (d < nodeIdxDim) break;
25363f27d899SToby Isaac     }
25373f27d899SToby Isaac     groupSize = nEnd - n;
253877f1a120SToby Isaac     /* get all of the vectors from the original and all of the pushforward vectors */
25393f27d899SToby Isaac     for (m = n; m < nEnd; m++) {
25403f27d899SToby Isaac       PetscInt d;
25413f27d899SToby Isaac 
25423f27d899SToby Isaac       for (d = 0; d < nodeVecDim; d++) {
25433f27d899SToby Isaac         V[(m - n) * nodeVecDim + d] = intNodeIndices->nodeVec[perm[m] * nodeVecDim + d];
25443f27d899SToby Isaac         W[(m - n) * nodeVecDim + d] = ni->nodeVec[permOrnt[m] * nodeVecDim + d];
25453f27d899SToby Isaac       }
25463f27d899SToby Isaac     }
254777f1a120SToby Isaac     /* now we have to solve for W in terms of V: the systems isn't always square, but the span
254877f1a120SToby Isaac      * of V and W should always be the same, so the solution of the normal equations works */
25493f27d899SToby Isaac     {
25503f27d899SToby Isaac       char         transpose = 'N';
25513f27d899SToby Isaac       PetscBLASInt bm        = nodeVecDim;
25523f27d899SToby Isaac       PetscBLASInt bn        = groupSize;
25533f27d899SToby Isaac       PetscBLASInt bnrhs     = groupSize;
25543f27d899SToby Isaac       PetscBLASInt blda      = bm;
25553f27d899SToby Isaac       PetscBLASInt bldb      = bm;
25563f27d899SToby Isaac       PetscBLASInt blwork    = 2 * nodeVecDim;
25573f27d899SToby Isaac       PetscBLASInt info;
25583f27d899SToby Isaac 
2559792fecdfSBarry Smith       PetscCallBLAS("LAPACKgels", LAPACKgels_(&transpose, &bm, &bn, &bnrhs, V, &blda, W, &bldb, work, &blwork, &info));
256008401ef6SPierre Jolivet       PetscCheck(info == 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Bad argument to GELS");
25613f27d899SToby Isaac       /* repack */
25623f27d899SToby Isaac       {
25633f27d899SToby Isaac         PetscInt i, j;
25643f27d899SToby Isaac 
25653f27d899SToby Isaac         for (i = 0; i < groupSize; i++) {
25663f27d899SToby Isaac           for (j = 0; j < groupSize; j++) {
256777f1a120SToby Isaac             /* notice the different leading dimension */
25683f27d899SToby Isaac             V[i * groupSize + j] = W[i * nodeVecDim + j];
25693f27d899SToby Isaac           }
25703f27d899SToby Isaac         }
25713f27d899SToby Isaac       }
2572c5c386beSToby Isaac       if (PetscDefined(USE_DEBUG)) {
2573c5c386beSToby Isaac         PetscReal res;
2574c5c386beSToby Isaac 
2575c5c386beSToby Isaac         /* check that the normal error is 0 */
2576c5c386beSToby Isaac         for (m = n; m < nEnd; m++) {
2577c5c386beSToby Isaac           PetscInt d;
2578c5c386beSToby Isaac 
2579ad540459SPierre Jolivet           for (d = 0; d < nodeVecDim; d++) W[(m - n) * nodeVecDim + d] = ni->nodeVec[permOrnt[m] * nodeVecDim + d];
2580c5c386beSToby Isaac         }
2581c5c386beSToby Isaac         res = 0.;
2582c5c386beSToby Isaac         for (PetscInt i = 0; i < groupSize; i++) {
2583c5c386beSToby Isaac           for (PetscInt j = 0; j < nodeVecDim; j++) {
2584ad540459SPierre Jolivet             for (PetscInt k = 0; k < groupSize; k++) W[i * nodeVecDim + j] -= V[i * groupSize + k] * intNodeIndices->nodeVec[perm[n + k] * nodeVecDim + j];
2585c5c386beSToby Isaac             res += PetscAbsScalar(W[i * nodeVecDim + j]);
2586c5c386beSToby Isaac           }
2587c5c386beSToby Isaac         }
258808401ef6SPierre Jolivet         PetscCheck(res <= PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_LIB, "Dof block did not solve");
2589c5c386beSToby Isaac       }
25903f27d899SToby Isaac     }
25919566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, groupSize, &permOrnt[n], groupSize, &perm[n], V, INSERT_VALUES));
25923f27d899SToby Isaac     n = nEnd;
25933f27d899SToby Isaac   }
25949566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
25959566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
25963f27d899SToby Isaac   *symMat = A;
25979566063dSJacob Faibussowitsch   PetscCall(PetscFree3(V, W, work));
25989566063dSJacob Faibussowitsch   PetscCall(PetscLagNodeIndicesDestroy(&ni));
25996f905325SMatthew G. Knepley   PetscFunctionReturn(0);
26006f905325SMatthew G. Knepley }
260120cf1dd8SToby Isaac 
260220cf1dd8SToby Isaac #define BaryIndex(perEdge, a, b, c) (((b) * (2 * perEdge + 1 - (b))) / 2) + (c)
260320cf1dd8SToby Isaac 
260420cf1dd8SToby Isaac #define CartIndex(perEdge, a, b) (perEdge * (a) + b)
260520cf1dd8SToby Isaac 
260677f1a120SToby Isaac /* the existing interface for symmetries is insufficient for all cases:
260777f1a120SToby Isaac  * - it should be sufficient for form degrees that are scalar (0 and n)
260877f1a120SToby Isaac  * - it should be sufficient for hypercube dofs
260977f1a120SToby Isaac  * - it isn't sufficient for simplex cells with non-scalar form degrees if
261077f1a120SToby Isaac  *   there are any dofs in the interior
261177f1a120SToby Isaac  *
261277f1a120SToby Isaac  * We compute the general transformation matrices, and if they fit, we return them,
261377f1a120SToby Isaac  * otherwise we error (but we should probably change the interface to allow for
261477f1a120SToby Isaac  * these symmetries)
261577f1a120SToby Isaac  */
26169371c9d4SSatish Balay static PetscErrorCode PetscDualSpaceGetSymmetries_Lagrange(PetscDualSpace sp, const PetscInt ****perms, const PetscScalar ****flips) {
261720cf1dd8SToby Isaac   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data;
26183f27d899SToby Isaac   PetscInt            dim, order, Nc;
261920cf1dd8SToby Isaac 
262020cf1dd8SToby Isaac   PetscFunctionBegin;
26219566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetOrder(sp, &order));
26229566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetNumComponents(sp, &Nc));
26239566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(sp->dm, &dim));
26243f27d899SToby Isaac   if (!lag->symComputed) { /* store symmetries */
26253f27d899SToby Isaac     PetscInt       pStart, pEnd, p;
26263f27d899SToby Isaac     PetscInt       numPoints;
262720cf1dd8SToby Isaac     PetscInt       numFaces;
26283f27d899SToby Isaac     PetscInt       spintdim;
26293f27d899SToby Isaac     PetscInt    ***symperms;
26303f27d899SToby Isaac     PetscScalar ***symflips;
263120cf1dd8SToby Isaac 
26329566063dSJacob Faibussowitsch     PetscCall(DMPlexGetChart(sp->dm, &pStart, &pEnd));
26333f27d899SToby Isaac     numPoints = pEnd - pStart;
2634b5a892a1SMatthew G. Knepley     {
2635b5a892a1SMatthew G. Knepley       DMPolytopeType ct;
2636b5a892a1SMatthew G. Knepley       /* The number of arrangements is no longer based on the number of faces */
26379566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCellType(sp->dm, 0, &ct));
2638b5a892a1SMatthew G. Knepley       numFaces = DMPolytopeTypeGetNumArrangments(ct) / 2;
2639b5a892a1SMatthew G. Knepley     }
26409566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(numPoints, &symperms));
26419566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(numPoints, &symflips));
26423f27d899SToby Isaac     spintdim = sp->spintdim;
26433f27d899SToby Isaac     /* The nodal symmetry behavior is not present when tensorSpace != tensorCell: someone might want this for the "S"
26443f27d899SToby Isaac      * family of FEEC spaces.  Most used in particular are discontinuous polynomial L2 spaces in tensor cells, where
26453f27d899SToby Isaac      * the symmetries are not necessary for FE assembly.  So for now we assume this is the case and don't return
26463f27d899SToby Isaac      * symmetries if tensorSpace != tensorCell */
26473f27d899SToby Isaac     if (spintdim && 0 < dim && dim < 3 && (lag->tensorSpace == lag->tensorCell)) { /* compute self symmetries */
26483f27d899SToby Isaac       PetscInt    **cellSymperms;
26493f27d899SToby Isaac       PetscScalar **cellSymflips;
26503f27d899SToby Isaac       PetscInt      ornt;
26513f27d899SToby Isaac       PetscInt      nCopies = Nc / lag->intNodeIndices->nodeVecDim;
26523f27d899SToby Isaac       PetscInt      nNodes  = lag->intNodeIndices->nNodes;
265320cf1dd8SToby Isaac 
265420cf1dd8SToby Isaac       lag->numSelfSym = 2 * numFaces;
265520cf1dd8SToby Isaac       lag->selfSymOff = numFaces;
26569566063dSJacob Faibussowitsch       PetscCall(PetscCalloc1(2 * numFaces, &cellSymperms));
26579566063dSJacob Faibussowitsch       PetscCall(PetscCalloc1(2 * numFaces, &cellSymflips));
265820cf1dd8SToby Isaac       /* we want to be able to index symmetries directly with the orientations, which range from [-numFaces,numFaces) */
26593f27d899SToby Isaac       symperms[0] = &cellSymperms[numFaces];
26603f27d899SToby Isaac       symflips[0] = &cellSymflips[numFaces];
26611dca8a05SBarry Smith       PetscCheck(lag->intNodeIndices->nodeVecDim * nCopies == Nc, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Node indices incompatible with dofs");
26621dca8a05SBarry Smith       PetscCheck(nNodes * nCopies == spintdim, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Node indices incompatible with dofs");
26633f27d899SToby Isaac       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 */
26643f27d899SToby Isaac         Mat          symMat;
26653f27d899SToby Isaac         PetscInt    *perm;
26663f27d899SToby Isaac         PetscScalar *flips;
26673f27d899SToby Isaac         PetscInt     i;
266820cf1dd8SToby Isaac 
26693f27d899SToby Isaac         if (!ornt) continue;
26709566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(spintdim, &perm));
26719566063dSJacob Faibussowitsch         PetscCall(PetscCalloc1(spintdim, &flips));
26723f27d899SToby Isaac         for (i = 0; i < spintdim; i++) perm[i] = -1;
26739566063dSJacob Faibussowitsch         PetscCall(PetscDualSpaceCreateInteriorSymmetryMatrix_Lagrange(sp, ornt, &symMat));
26743f27d899SToby Isaac         for (i = 0; i < nNodes; i++) {
26753f27d899SToby Isaac           PetscInt           ncols;
26763f27d899SToby Isaac           PetscInt           j, k;
26773f27d899SToby Isaac           const PetscInt    *cols;
26783f27d899SToby Isaac           const PetscScalar *vals;
26793f27d899SToby Isaac           PetscBool          nz_seen = PETSC_FALSE;
268020cf1dd8SToby Isaac 
26819566063dSJacob Faibussowitsch           PetscCall(MatGetRow(symMat, i, &ncols, &cols, &vals));
26823f27d899SToby Isaac           for (j = 0; j < ncols; j++) {
26833f27d899SToby Isaac             if (PetscAbsScalar(vals[j]) > PETSC_SMALL) {
268428b400f6SJacob Faibussowitsch               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");
26853f27d899SToby Isaac               nz_seen = PETSC_TRUE;
26861dca8a05SBarry Smith               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");
26871dca8a05SBarry Smith               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");
26881dca8a05SBarry Smith               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");
2689ad540459SPierre Jolivet               for (k = 0; k < nCopies; k++) perm[cols[j] * nCopies + k] = i * nCopies + k;
26903f27d899SToby Isaac               if (PetscRealPart(vals[j]) < 0.) {
2691ad540459SPierre Jolivet                 for (k = 0; k < nCopies; k++) flips[i * nCopies + k] = -1.;
269220cf1dd8SToby Isaac               } else {
2693ad540459SPierre Jolivet                 for (k = 0; k < nCopies; k++) flips[i * nCopies + k] = 1.;
26943f27d899SToby Isaac               }
26953f27d899SToby Isaac             }
26963f27d899SToby Isaac           }
26979566063dSJacob Faibussowitsch           PetscCall(MatRestoreRow(symMat, i, &ncols, &cols, &vals));
26983f27d899SToby Isaac         }
26999566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&symMat));
27003f27d899SToby Isaac         /* if there were no sign flips, keep NULL */
27019371c9d4SSatish Balay         for (i = 0; i < spintdim; i++)
27029371c9d4SSatish Balay           if (flips[i] != 1.) break;
27033f27d899SToby Isaac         if (i == spintdim) {
27049566063dSJacob Faibussowitsch           PetscCall(PetscFree(flips));
27053f27d899SToby Isaac           flips = NULL;
27063f27d899SToby Isaac         }
27073f27d899SToby Isaac         /* if the permutation is identity, keep NULL */
27089371c9d4SSatish Balay         for (i = 0; i < spintdim; i++)
27099371c9d4SSatish Balay           if (perm[i] != i) break;
27103f27d899SToby Isaac         if (i == spintdim) {
27119566063dSJacob Faibussowitsch           PetscCall(PetscFree(perm));
27123f27d899SToby Isaac           perm = NULL;
27133f27d899SToby Isaac         }
27143f27d899SToby Isaac         symperms[0][ornt] = perm;
27153f27d899SToby Isaac         symflips[0][ornt] = flips;
27163f27d899SToby Isaac       }
27173f27d899SToby Isaac       /* if no orientations produced non-identity permutations, keep NULL */
27189371c9d4SSatish Balay       for (ornt = -numFaces; ornt < numFaces; ornt++)
27199371c9d4SSatish Balay         if (symperms[0][ornt]) break;
27203f27d899SToby Isaac       if (ornt == numFaces) {
27219566063dSJacob Faibussowitsch         PetscCall(PetscFree(cellSymperms));
27223f27d899SToby Isaac         symperms[0] = NULL;
27233f27d899SToby Isaac       }
27243f27d899SToby Isaac       /* if no orientations produced sign flips, keep NULL */
27259371c9d4SSatish Balay       for (ornt = -numFaces; ornt < numFaces; ornt++)
27269371c9d4SSatish Balay         if (symflips[0][ornt]) break;
27273f27d899SToby Isaac       if (ornt == numFaces) {
27289566063dSJacob Faibussowitsch         PetscCall(PetscFree(cellSymflips));
27293f27d899SToby Isaac         symflips[0] = NULL;
27303f27d899SToby Isaac       }
27313f27d899SToby Isaac     }
273277f1a120SToby Isaac     { /* get the symmetries of closure points */
27333f27d899SToby Isaac       PetscInt  closureSize = 0;
27343f27d899SToby Isaac       PetscInt *closure     = NULL;
27353f27d899SToby Isaac       PetscInt  r;
273620cf1dd8SToby Isaac 
27379566063dSJacob Faibussowitsch       PetscCall(DMPlexGetTransitiveClosure(sp->dm, 0, PETSC_TRUE, &closureSize, &closure));
27383f27d899SToby Isaac       for (r = 0; r < closureSize; r++) {
27393f27d899SToby Isaac         PetscDualSpace       psp;
27403f27d899SToby Isaac         PetscInt             point = closure[2 * r];
27413f27d899SToby Isaac         PetscInt             pspintdim;
27423f27d899SToby Isaac         const PetscInt    ***psymperms = NULL;
27433f27d899SToby Isaac         const PetscScalar ***psymflips = NULL;
274420cf1dd8SToby Isaac 
27453f27d899SToby Isaac         if (!point) continue;
27469566063dSJacob Faibussowitsch         PetscCall(PetscDualSpaceGetPointSubspace(sp, point, &psp));
27473f27d899SToby Isaac         if (!psp) continue;
27489566063dSJacob Faibussowitsch         PetscCall(PetscDualSpaceGetInteriorDimension(psp, &pspintdim));
27493f27d899SToby Isaac         if (!pspintdim) continue;
27509566063dSJacob Faibussowitsch         PetscCall(PetscDualSpaceGetSymmetries(psp, &psymperms, &psymflips));
27513f27d899SToby Isaac         symperms[r] = (PetscInt **)(psymperms ? psymperms[0] : NULL);
27523f27d899SToby Isaac         symflips[r] = (PetscScalar **)(psymflips ? psymflips[0] : NULL);
275320cf1dd8SToby Isaac       }
27549566063dSJacob Faibussowitsch       PetscCall(DMPlexRestoreTransitiveClosure(sp->dm, 0, PETSC_TRUE, &closureSize, &closure));
275520cf1dd8SToby Isaac     }
27569371c9d4SSatish Balay     for (p = 0; p < pEnd; p++)
27579371c9d4SSatish Balay       if (symperms[p]) break;
27583f27d899SToby Isaac     if (p == pEnd) {
27599566063dSJacob Faibussowitsch       PetscCall(PetscFree(symperms));
27603f27d899SToby Isaac       symperms = NULL;
276120cf1dd8SToby Isaac     }
27629371c9d4SSatish Balay     for (p = 0; p < pEnd; p++)
27639371c9d4SSatish Balay       if (symflips[p]) break;
27643f27d899SToby Isaac     if (p == pEnd) {
27659566063dSJacob Faibussowitsch       PetscCall(PetscFree(symflips));
27663f27d899SToby Isaac       symflips = NULL;
276720cf1dd8SToby Isaac     }
27683f27d899SToby Isaac     lag->symperms    = symperms;
27693f27d899SToby Isaac     lag->symflips    = symflips;
27703f27d899SToby Isaac     lag->symComputed = PETSC_TRUE;
277120cf1dd8SToby Isaac   }
27723f27d899SToby Isaac   if (perms) *perms = (const PetscInt ***)lag->symperms;
27733f27d899SToby Isaac   if (flips) *flips = (const PetscScalar ***)lag->symflips;
277420cf1dd8SToby Isaac   PetscFunctionReturn(0);
277520cf1dd8SToby Isaac }
277620cf1dd8SToby Isaac 
27779371c9d4SSatish Balay static PetscErrorCode PetscDualSpaceLagrangeGetContinuity_Lagrange(PetscDualSpace sp, PetscBool *continuous) {
277820cf1dd8SToby Isaac   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data;
277920cf1dd8SToby Isaac 
278020cf1dd8SToby Isaac   PetscFunctionBegin;
278120cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
2782dadcf809SJacob Faibussowitsch   PetscValidBoolPointer(continuous, 2);
278320cf1dd8SToby Isaac   *continuous = lag->continuous;
278420cf1dd8SToby Isaac   PetscFunctionReturn(0);
278520cf1dd8SToby Isaac }
278620cf1dd8SToby Isaac 
27879371c9d4SSatish Balay static PetscErrorCode PetscDualSpaceLagrangeSetContinuity_Lagrange(PetscDualSpace sp, PetscBool continuous) {
278820cf1dd8SToby Isaac   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data;
278920cf1dd8SToby Isaac 
279020cf1dd8SToby Isaac   PetscFunctionBegin;
279120cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
279220cf1dd8SToby Isaac   lag->continuous = continuous;
279320cf1dd8SToby Isaac   PetscFunctionReturn(0);
279420cf1dd8SToby Isaac }
279520cf1dd8SToby Isaac 
279620cf1dd8SToby Isaac /*@
279720cf1dd8SToby Isaac   PetscDualSpaceLagrangeGetContinuity - Retrieves the flag for element continuity
279820cf1dd8SToby Isaac 
279920cf1dd8SToby Isaac   Not Collective
280020cf1dd8SToby Isaac 
280120cf1dd8SToby Isaac   Input Parameter:
280220cf1dd8SToby Isaac . sp         - the PetscDualSpace
280320cf1dd8SToby Isaac 
280420cf1dd8SToby Isaac   Output Parameter:
280520cf1dd8SToby Isaac . continuous - flag for element continuity
280620cf1dd8SToby Isaac 
280720cf1dd8SToby Isaac   Level: intermediate
280820cf1dd8SToby Isaac 
2809db781477SPatrick Sanan .seealso: `PetscDualSpaceLagrangeSetContinuity()`
281020cf1dd8SToby Isaac @*/
28119371c9d4SSatish Balay PetscErrorCode PetscDualSpaceLagrangeGetContinuity(PetscDualSpace sp, PetscBool *continuous) {
281220cf1dd8SToby Isaac   PetscFunctionBegin;
281320cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
2814dadcf809SJacob Faibussowitsch   PetscValidBoolPointer(continuous, 2);
2815cac4c232SBarry Smith   PetscTryMethod(sp, "PetscDualSpaceLagrangeGetContinuity_C", (PetscDualSpace, PetscBool *), (sp, continuous));
281620cf1dd8SToby Isaac   PetscFunctionReturn(0);
281720cf1dd8SToby Isaac }
281820cf1dd8SToby Isaac 
281920cf1dd8SToby Isaac /*@
282020cf1dd8SToby Isaac   PetscDualSpaceLagrangeSetContinuity - Indicate whether the element is continuous
282120cf1dd8SToby Isaac 
2822d083f849SBarry Smith   Logically Collective on sp
282320cf1dd8SToby Isaac 
282420cf1dd8SToby Isaac   Input Parameters:
282520cf1dd8SToby Isaac + sp         - the PetscDualSpace
282620cf1dd8SToby Isaac - continuous - flag for element continuity
282720cf1dd8SToby Isaac 
282820cf1dd8SToby Isaac   Options Database:
2829147403d9SBarry Smith . -petscdualspace_lagrange_continuity <bool> - use a continuous element
283020cf1dd8SToby Isaac 
283120cf1dd8SToby Isaac   Level: intermediate
283220cf1dd8SToby Isaac 
2833db781477SPatrick Sanan .seealso: `PetscDualSpaceLagrangeGetContinuity()`
283420cf1dd8SToby Isaac @*/
28359371c9d4SSatish Balay PetscErrorCode PetscDualSpaceLagrangeSetContinuity(PetscDualSpace sp, PetscBool continuous) {
283620cf1dd8SToby Isaac   PetscFunctionBegin;
283720cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
283820cf1dd8SToby Isaac   PetscValidLogicalCollectiveBool(sp, continuous, 2);
2839cac4c232SBarry Smith   PetscTryMethod(sp, "PetscDualSpaceLagrangeSetContinuity_C", (PetscDualSpace, PetscBool), (sp, continuous));
284020cf1dd8SToby Isaac   PetscFunctionReturn(0);
284120cf1dd8SToby Isaac }
284220cf1dd8SToby Isaac 
28439371c9d4SSatish Balay static PetscErrorCode PetscDualSpaceLagrangeGetTensor_Lagrange(PetscDualSpace sp, PetscBool *tensor) {
284420cf1dd8SToby Isaac   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data;
28456f905325SMatthew G. Knepley 
28466f905325SMatthew G. Knepley   PetscFunctionBegin;
28476f905325SMatthew G. Knepley   *tensor = lag->tensorSpace;
28486f905325SMatthew G. Knepley   PetscFunctionReturn(0);
28496f905325SMatthew G. Knepley }
28506f905325SMatthew G. Knepley 
28519371c9d4SSatish Balay static PetscErrorCode PetscDualSpaceLagrangeSetTensor_Lagrange(PetscDualSpace sp, PetscBool tensor) {
28526f905325SMatthew G. Knepley   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data;
28536f905325SMatthew G. Knepley 
28546f905325SMatthew G. Knepley   PetscFunctionBegin;
28556f905325SMatthew G. Knepley   lag->tensorSpace = tensor;
28566f905325SMatthew G. Knepley   PetscFunctionReturn(0);
28576f905325SMatthew G. Knepley }
28586f905325SMatthew G. Knepley 
28599371c9d4SSatish Balay static PetscErrorCode PetscDualSpaceLagrangeGetTrimmed_Lagrange(PetscDualSpace sp, PetscBool *trimmed) {
28603f27d899SToby Isaac   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data;
28613f27d899SToby Isaac 
28623f27d899SToby Isaac   PetscFunctionBegin;
28633f27d899SToby Isaac   *trimmed = lag->trimmed;
28643f27d899SToby Isaac   PetscFunctionReturn(0);
28653f27d899SToby Isaac }
28663f27d899SToby Isaac 
28679371c9d4SSatish Balay static PetscErrorCode PetscDualSpaceLagrangeSetTrimmed_Lagrange(PetscDualSpace sp, PetscBool trimmed) {
28683f27d899SToby Isaac   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data;
28693f27d899SToby Isaac 
28703f27d899SToby Isaac   PetscFunctionBegin;
28713f27d899SToby Isaac   lag->trimmed = trimmed;
28723f27d899SToby Isaac   PetscFunctionReturn(0);
28733f27d899SToby Isaac }
28743f27d899SToby Isaac 
28759371c9d4SSatish Balay static PetscErrorCode PetscDualSpaceLagrangeGetNodeType_Lagrange(PetscDualSpace sp, PetscDTNodeType *nodeType, PetscBool *boundary, PetscReal *exponent) {
28763f27d899SToby Isaac   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data;
28773f27d899SToby Isaac 
28783f27d899SToby Isaac   PetscFunctionBegin;
28793f27d899SToby Isaac   if (nodeType) *nodeType = lag->nodeType;
28803f27d899SToby Isaac   if (boundary) *boundary = lag->endNodes;
28813f27d899SToby Isaac   if (exponent) *exponent = lag->nodeExponent;
28823f27d899SToby Isaac   PetscFunctionReturn(0);
28833f27d899SToby Isaac }
28843f27d899SToby Isaac 
28859371c9d4SSatish Balay static PetscErrorCode PetscDualSpaceLagrangeSetNodeType_Lagrange(PetscDualSpace sp, PetscDTNodeType nodeType, PetscBool boundary, PetscReal exponent) {
28863f27d899SToby Isaac   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data;
28873f27d899SToby Isaac 
28883f27d899SToby Isaac   PetscFunctionBegin;
28891dca8a05SBarry Smith   PetscCheck(nodeType != PETSCDTNODES_GAUSSJACOBI || exponent > -1., PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_OUTOFRANGE, "Exponent must be > -1");
28903f27d899SToby Isaac   lag->nodeType     = nodeType;
28913f27d899SToby Isaac   lag->endNodes     = boundary;
28923f27d899SToby Isaac   lag->nodeExponent = exponent;
28933f27d899SToby Isaac   PetscFunctionReturn(0);
28943f27d899SToby Isaac }
28953f27d899SToby Isaac 
28969371c9d4SSatish Balay static PetscErrorCode PetscDualSpaceLagrangeGetUseMoments_Lagrange(PetscDualSpace sp, PetscBool *useMoments) {
289766a6c23cSMatthew G. Knepley   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data;
289866a6c23cSMatthew G. Knepley 
289966a6c23cSMatthew G. Knepley   PetscFunctionBegin;
290066a6c23cSMatthew G. Knepley   *useMoments = lag->useMoments;
290166a6c23cSMatthew G. Knepley   PetscFunctionReturn(0);
290266a6c23cSMatthew G. Knepley }
290366a6c23cSMatthew G. Knepley 
29049371c9d4SSatish Balay static PetscErrorCode PetscDualSpaceLagrangeSetUseMoments_Lagrange(PetscDualSpace sp, PetscBool useMoments) {
290566a6c23cSMatthew G. Knepley   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data;
290666a6c23cSMatthew G. Knepley 
290766a6c23cSMatthew G. Knepley   PetscFunctionBegin;
290866a6c23cSMatthew G. Knepley   lag->useMoments = useMoments;
290966a6c23cSMatthew G. Knepley   PetscFunctionReturn(0);
291066a6c23cSMatthew G. Knepley }
291166a6c23cSMatthew G. Knepley 
29129371c9d4SSatish Balay static PetscErrorCode PetscDualSpaceLagrangeGetMomentOrder_Lagrange(PetscDualSpace sp, PetscInt *momentOrder) {
291366a6c23cSMatthew G. Knepley   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data;
291466a6c23cSMatthew G. Knepley 
291566a6c23cSMatthew G. Knepley   PetscFunctionBegin;
291666a6c23cSMatthew G. Knepley   *momentOrder = lag->momentOrder;
291766a6c23cSMatthew G. Knepley   PetscFunctionReturn(0);
291866a6c23cSMatthew G. Knepley }
291966a6c23cSMatthew G. Knepley 
29209371c9d4SSatish Balay static PetscErrorCode PetscDualSpaceLagrangeSetMomentOrder_Lagrange(PetscDualSpace sp, PetscInt momentOrder) {
292166a6c23cSMatthew G. Knepley   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data;
292266a6c23cSMatthew G. Knepley 
292366a6c23cSMatthew G. Knepley   PetscFunctionBegin;
292466a6c23cSMatthew G. Knepley   lag->momentOrder = momentOrder;
292566a6c23cSMatthew G. Knepley   PetscFunctionReturn(0);
292666a6c23cSMatthew G. Knepley }
292766a6c23cSMatthew G. Knepley 
29286f905325SMatthew G. Knepley /*@
29296f905325SMatthew G. Knepley   PetscDualSpaceLagrangeGetTensor - Get the tensor nature of the dual space
29306f905325SMatthew G. Knepley 
29316f905325SMatthew G. Knepley   Not collective
29326f905325SMatthew G. Knepley 
29336f905325SMatthew G. Knepley   Input Parameter:
29346f905325SMatthew G. Knepley . sp - The PetscDualSpace
29356f905325SMatthew G. Knepley 
29366f905325SMatthew G. Knepley   Output Parameter:
29376f905325SMatthew G. Knepley . tensor - Whether the dual space has tensor layout (vs. simplicial)
29386f905325SMatthew G. Knepley 
29396f905325SMatthew G. Knepley   Level: intermediate
29406f905325SMatthew G. Knepley 
2941db781477SPatrick Sanan .seealso: `PetscDualSpaceLagrangeSetTensor()`, `PetscDualSpaceCreate()`
29426f905325SMatthew G. Knepley @*/
29439371c9d4SSatish Balay PetscErrorCode PetscDualSpaceLagrangeGetTensor(PetscDualSpace sp, PetscBool *tensor) {
294420cf1dd8SToby Isaac   PetscFunctionBegin;
294520cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
2946dadcf809SJacob Faibussowitsch   PetscValidBoolPointer(tensor, 2);
2947cac4c232SBarry Smith   PetscTryMethod(sp, "PetscDualSpaceLagrangeGetTensor_C", (PetscDualSpace, PetscBool *), (sp, tensor));
294820cf1dd8SToby Isaac   PetscFunctionReturn(0);
294920cf1dd8SToby Isaac }
295020cf1dd8SToby Isaac 
29516f905325SMatthew G. Knepley /*@
29526f905325SMatthew G. Knepley   PetscDualSpaceLagrangeSetTensor - Set the tensor nature of the dual space
29536f905325SMatthew G. Knepley 
29546f905325SMatthew G. Knepley   Not collective
29556f905325SMatthew G. Knepley 
29566f905325SMatthew G. Knepley   Input Parameters:
29576f905325SMatthew G. Knepley + sp - The PetscDualSpace
29586f905325SMatthew G. Knepley - tensor - Whether the dual space has tensor layout (vs. simplicial)
29596f905325SMatthew G. Knepley 
29606f905325SMatthew G. Knepley   Level: intermediate
29616f905325SMatthew G. Knepley 
2962db781477SPatrick Sanan .seealso: `PetscDualSpaceLagrangeGetTensor()`, `PetscDualSpaceCreate()`
29636f905325SMatthew G. Knepley @*/
29649371c9d4SSatish Balay PetscErrorCode PetscDualSpaceLagrangeSetTensor(PetscDualSpace sp, PetscBool tensor) {
29656f905325SMatthew G. Knepley   PetscFunctionBegin;
29666f905325SMatthew G. Knepley   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
2967cac4c232SBarry Smith   PetscTryMethod(sp, "PetscDualSpaceLagrangeSetTensor_C", (PetscDualSpace, PetscBool), (sp, tensor));
29686f905325SMatthew G. Knepley   PetscFunctionReturn(0);
29696f905325SMatthew G. Knepley }
29706f905325SMatthew G. Knepley 
29713f27d899SToby Isaac /*@
29723f27d899SToby Isaac   PetscDualSpaceLagrangeGetTrimmed - Get the trimmed nature of the dual space
29733f27d899SToby Isaac 
29743f27d899SToby Isaac   Not collective
29753f27d899SToby Isaac 
29763f27d899SToby Isaac   Input Parameter:
29773f27d899SToby Isaac . sp - The PetscDualSpace
29783f27d899SToby Isaac 
29793f27d899SToby Isaac   Output Parameter:
29803f27d899SToby Isaac . 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)
29813f27d899SToby Isaac 
29823f27d899SToby Isaac   Level: intermediate
29833f27d899SToby Isaac 
2984db781477SPatrick Sanan .seealso: `PetscDualSpaceLagrangeSetTrimmed()`, `PetscDualSpaceCreate()`
29853f27d899SToby Isaac @*/
29869371c9d4SSatish Balay PetscErrorCode PetscDualSpaceLagrangeGetTrimmed(PetscDualSpace sp, PetscBool *trimmed) {
29873f27d899SToby Isaac   PetscFunctionBegin;
29883f27d899SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
2989dadcf809SJacob Faibussowitsch   PetscValidBoolPointer(trimmed, 2);
2990cac4c232SBarry Smith   PetscTryMethod(sp, "PetscDualSpaceLagrangeGetTrimmed_C", (PetscDualSpace, PetscBool *), (sp, trimmed));
29913f27d899SToby Isaac   PetscFunctionReturn(0);
29923f27d899SToby Isaac }
29933f27d899SToby Isaac 
29943f27d899SToby Isaac /*@
29953f27d899SToby Isaac   PetscDualSpaceLagrangeSetTrimmed - Set the trimmed nature of the dual space
29963f27d899SToby Isaac 
29973f27d899SToby Isaac   Not collective
29983f27d899SToby Isaac 
29993f27d899SToby Isaac   Input Parameters:
30003f27d899SToby Isaac + sp - The PetscDualSpace
30013f27d899SToby Isaac - 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)
30023f27d899SToby Isaac 
30033f27d899SToby Isaac   Level: intermediate
30043f27d899SToby Isaac 
3005db781477SPatrick Sanan .seealso: `PetscDualSpaceLagrangeGetTrimmed()`, `PetscDualSpaceCreate()`
30063f27d899SToby Isaac @*/
30079371c9d4SSatish Balay PetscErrorCode PetscDualSpaceLagrangeSetTrimmed(PetscDualSpace sp, PetscBool trimmed) {
30083f27d899SToby Isaac   PetscFunctionBegin;
30093f27d899SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
3010cac4c232SBarry Smith   PetscTryMethod(sp, "PetscDualSpaceLagrangeSetTrimmed_C", (PetscDualSpace, PetscBool), (sp, trimmed));
30113f27d899SToby Isaac   PetscFunctionReturn(0);
30123f27d899SToby Isaac }
30133f27d899SToby Isaac 
30143f27d899SToby Isaac /*@
30153f27d899SToby Isaac   PetscDualSpaceLagrangeGetNodeType - Get a description of how nodes are laid out for Lagrange polynomials in this
30163f27d899SToby Isaac   dual space
30173f27d899SToby Isaac 
30183f27d899SToby Isaac   Not collective
30193f27d899SToby Isaac 
30203f27d899SToby Isaac   Input Parameter:
30213f27d899SToby Isaac . sp - The PetscDualSpace
30223f27d899SToby Isaac 
30233f27d899SToby Isaac   Output Parameters:
30243f27d899SToby Isaac + nodeType - The type of nodes
30253f27d899SToby Isaac . boundary - Whether the node type is one that includes endpoints (if nodeType is PETSCDTNODES_GAUSSJACOBI, nodes that
30263f27d899SToby Isaac              include the boundary are Gauss-Lobatto-Jacobi nodes)
30273f27d899SToby Isaac - exponent - If nodeType is PETSCDTNODES_GAUSJACOBI, indicates the exponent used for both ends of the 1D Jacobi weight function
30283f27d899SToby Isaac              '0' is Gauss-Legendre, '-0.5' is Gauss-Chebyshev of the first type, '0.5' is Gauss-Chebyshev of the second type
30293f27d899SToby Isaac 
30303f27d899SToby Isaac   Level: advanced
30313f27d899SToby Isaac 
3032db781477SPatrick Sanan .seealso: `PetscDTNodeType`, `PetscDualSpaceLagrangeSetNodeType()`
30333f27d899SToby Isaac @*/
30349371c9d4SSatish Balay PetscErrorCode PetscDualSpaceLagrangeGetNodeType(PetscDualSpace sp, PetscDTNodeType *nodeType, PetscBool *boundary, PetscReal *exponent) {
30353f27d899SToby Isaac   PetscFunctionBegin;
30363f27d899SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
30373f27d899SToby Isaac   if (nodeType) PetscValidPointer(nodeType, 2);
3038dadcf809SJacob Faibussowitsch   if (boundary) PetscValidBoolPointer(boundary, 3);
3039dadcf809SJacob Faibussowitsch   if (exponent) PetscValidRealPointer(exponent, 4);
3040cac4c232SBarry Smith   PetscTryMethod(sp, "PetscDualSpaceLagrangeGetNodeType_C", (PetscDualSpace, PetscDTNodeType *, PetscBool *, PetscReal *), (sp, nodeType, boundary, exponent));
30413f27d899SToby Isaac   PetscFunctionReturn(0);
30423f27d899SToby Isaac }
30433f27d899SToby Isaac 
30443f27d899SToby Isaac /*@
30453f27d899SToby Isaac   PetscDualSpaceLagrangeSetNodeType - Set a description of how nodes are laid out for Lagrange polynomials in this
30463f27d899SToby Isaac   dual space
30473f27d899SToby Isaac 
30483f27d899SToby Isaac   Logically collective
30493f27d899SToby Isaac 
30503f27d899SToby Isaac   Input Parameters:
30513f27d899SToby Isaac + sp - The PetscDualSpace
30523f27d899SToby Isaac . nodeType - The type of nodes
30533f27d899SToby Isaac . boundary - Whether the node type is one that includes endpoints (if nodeType is PETSCDTNODES_GAUSSJACOBI, nodes that
30543f27d899SToby Isaac              include the boundary are Gauss-Lobatto-Jacobi nodes)
30553f27d899SToby Isaac - exponent - If nodeType is PETSCDTNODES_GAUSJACOBI, indicates the exponent used for both ends of the 1D Jacobi weight function
30563f27d899SToby Isaac              '0' is Gauss-Legendre, '-0.5' is Gauss-Chebyshev of the first type, '0.5' is Gauss-Chebyshev of the second type
30573f27d899SToby Isaac 
30583f27d899SToby Isaac   Level: advanced
30593f27d899SToby Isaac 
3060db781477SPatrick Sanan .seealso: `PetscDTNodeType`, `PetscDualSpaceLagrangeGetNodeType()`
30613f27d899SToby Isaac @*/
30629371c9d4SSatish Balay PetscErrorCode PetscDualSpaceLagrangeSetNodeType(PetscDualSpace sp, PetscDTNodeType nodeType, PetscBool boundary, PetscReal exponent) {
30633f27d899SToby Isaac   PetscFunctionBegin;
30643f27d899SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
3065cac4c232SBarry Smith   PetscTryMethod(sp, "PetscDualSpaceLagrangeSetNodeType_C", (PetscDualSpace, PetscDTNodeType, PetscBool, PetscReal), (sp, nodeType, boundary, exponent));
30663f27d899SToby Isaac   PetscFunctionReturn(0);
30673f27d899SToby Isaac }
30683f27d899SToby Isaac 
306966a6c23cSMatthew G. Knepley /*@
307066a6c23cSMatthew G. Knepley   PetscDualSpaceLagrangeGetUseMoments - Get the flag for using moment functionals
307166a6c23cSMatthew G. Knepley 
307266a6c23cSMatthew G. Knepley   Not collective
307366a6c23cSMatthew G. Knepley 
307466a6c23cSMatthew G. Knepley   Input Parameter:
307566a6c23cSMatthew G. Knepley . sp - The PetscDualSpace
307666a6c23cSMatthew G. Knepley 
307766a6c23cSMatthew G. Knepley   Output Parameter:
307866a6c23cSMatthew G. Knepley . useMoments - Moment flag
307966a6c23cSMatthew G. Knepley 
308066a6c23cSMatthew G. Knepley   Level: advanced
308166a6c23cSMatthew G. Knepley 
3082db781477SPatrick Sanan .seealso: `PetscDualSpaceLagrangeSetUseMoments()`
308366a6c23cSMatthew G. Knepley @*/
30849371c9d4SSatish Balay PetscErrorCode PetscDualSpaceLagrangeGetUseMoments(PetscDualSpace sp, PetscBool *useMoments) {
308566a6c23cSMatthew G. Knepley   PetscFunctionBegin;
308666a6c23cSMatthew G. Knepley   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
308766a6c23cSMatthew G. Knepley   PetscValidBoolPointer(useMoments, 2);
3088cac4c232SBarry Smith   PetscUseMethod(sp, "PetscDualSpaceLagrangeGetUseMoments_C", (PetscDualSpace, PetscBool *), (sp, useMoments));
308966a6c23cSMatthew G. Knepley   PetscFunctionReturn(0);
309066a6c23cSMatthew G. Knepley }
309166a6c23cSMatthew G. Knepley 
309266a6c23cSMatthew G. Knepley /*@
309366a6c23cSMatthew G. Knepley   PetscDualSpaceLagrangeSetUseMoments - Set the flag for moment functionals
309466a6c23cSMatthew G. Knepley 
309566a6c23cSMatthew G. Knepley   Logically collective
309666a6c23cSMatthew G. Knepley 
309766a6c23cSMatthew G. Knepley   Input Parameters:
309866a6c23cSMatthew G. Knepley + sp - The PetscDualSpace
309966a6c23cSMatthew G. Knepley - useMoments - The flag for moment functionals
310066a6c23cSMatthew G. Knepley 
310166a6c23cSMatthew G. Knepley   Level: advanced
310266a6c23cSMatthew G. Knepley 
3103db781477SPatrick Sanan .seealso: `PetscDualSpaceLagrangeGetUseMoments()`
310466a6c23cSMatthew G. Knepley @*/
31059371c9d4SSatish Balay PetscErrorCode PetscDualSpaceLagrangeSetUseMoments(PetscDualSpace sp, PetscBool useMoments) {
310666a6c23cSMatthew G. Knepley   PetscFunctionBegin;
310766a6c23cSMatthew G. Knepley   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
3108cac4c232SBarry Smith   PetscTryMethod(sp, "PetscDualSpaceLagrangeSetUseMoments_C", (PetscDualSpace, PetscBool), (sp, useMoments));
310966a6c23cSMatthew G. Knepley   PetscFunctionReturn(0);
311066a6c23cSMatthew G. Knepley }
311166a6c23cSMatthew G. Knepley 
311266a6c23cSMatthew G. Knepley /*@
311366a6c23cSMatthew G. Knepley   PetscDualSpaceLagrangeGetMomentOrder - Get the order for moment integration
311466a6c23cSMatthew G. Knepley 
311566a6c23cSMatthew G. Knepley   Not collective
311666a6c23cSMatthew G. Knepley 
311766a6c23cSMatthew G. Knepley   Input Parameter:
311866a6c23cSMatthew G. Knepley . sp - The PetscDualSpace
311966a6c23cSMatthew G. Knepley 
312066a6c23cSMatthew G. Knepley   Output Parameter:
312166a6c23cSMatthew G. Knepley . order - Moment integration order
312266a6c23cSMatthew G. Knepley 
312366a6c23cSMatthew G. Knepley   Level: advanced
312466a6c23cSMatthew G. Knepley 
3125db781477SPatrick Sanan .seealso: `PetscDualSpaceLagrangeSetMomentOrder()`
312666a6c23cSMatthew G. Knepley @*/
31279371c9d4SSatish Balay PetscErrorCode PetscDualSpaceLagrangeGetMomentOrder(PetscDualSpace sp, PetscInt *order) {
312866a6c23cSMatthew G. Knepley   PetscFunctionBegin;
312966a6c23cSMatthew G. Knepley   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
313066a6c23cSMatthew G. Knepley   PetscValidIntPointer(order, 2);
3131cac4c232SBarry Smith   PetscUseMethod(sp, "PetscDualSpaceLagrangeGetMomentOrder_C", (PetscDualSpace, PetscInt *), (sp, order));
313266a6c23cSMatthew G. Knepley   PetscFunctionReturn(0);
313366a6c23cSMatthew G. Knepley }
313466a6c23cSMatthew G. Knepley 
313566a6c23cSMatthew G. Knepley /*@
313666a6c23cSMatthew G. Knepley   PetscDualSpaceLagrangeSetMomentOrder - Set the order for moment integration
313766a6c23cSMatthew G. Knepley 
313866a6c23cSMatthew G. Knepley   Logically collective
313966a6c23cSMatthew G. Knepley 
314066a6c23cSMatthew G. Knepley   Input Parameters:
314166a6c23cSMatthew G. Knepley + sp - The PetscDualSpace
314266a6c23cSMatthew G. Knepley - order - The order for moment integration
314366a6c23cSMatthew G. Knepley 
314466a6c23cSMatthew G. Knepley   Level: advanced
314566a6c23cSMatthew G. Knepley 
3146db781477SPatrick Sanan .seealso: `PetscDualSpaceLagrangeGetMomentOrder()`
314766a6c23cSMatthew G. Knepley @*/
31489371c9d4SSatish Balay PetscErrorCode PetscDualSpaceLagrangeSetMomentOrder(PetscDualSpace sp, PetscInt order) {
314966a6c23cSMatthew G. Knepley   PetscFunctionBegin;
315066a6c23cSMatthew G. Knepley   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
3151cac4c232SBarry Smith   PetscTryMethod(sp, "PetscDualSpaceLagrangeSetMomentOrder_C", (PetscDualSpace, PetscInt), (sp, order));
315266a6c23cSMatthew G. Knepley   PetscFunctionReturn(0);
315366a6c23cSMatthew G. Knepley }
31543f27d899SToby Isaac 
31559371c9d4SSatish Balay static PetscErrorCode PetscDualSpaceInitialize_Lagrange(PetscDualSpace sp) {
315620cf1dd8SToby Isaac   PetscFunctionBegin;
315720cf1dd8SToby Isaac   sp->ops->destroy              = PetscDualSpaceDestroy_Lagrange;
31586f905325SMatthew G. Knepley   sp->ops->view                 = PetscDualSpaceView_Lagrange;
31596f905325SMatthew G. Knepley   sp->ops->setfromoptions       = PetscDualSpaceSetFromOptions_Lagrange;
316020cf1dd8SToby Isaac   sp->ops->duplicate            = PetscDualSpaceDuplicate_Lagrange;
31616f905325SMatthew G. Knepley   sp->ops->setup                = PetscDualSpaceSetUp_Lagrange;
31623f27d899SToby Isaac   sp->ops->createheightsubspace = NULL;
31633f27d899SToby Isaac   sp->ops->createpointsubspace  = NULL;
316420cf1dd8SToby Isaac   sp->ops->getsymmetries        = PetscDualSpaceGetSymmetries_Lagrange;
316520cf1dd8SToby Isaac   sp->ops->apply                = PetscDualSpaceApplyDefault;
316620cf1dd8SToby Isaac   sp->ops->applyall             = PetscDualSpaceApplyAllDefault;
3167b4457527SToby Isaac   sp->ops->applyint             = PetscDualSpaceApplyInteriorDefault;
31683f27d899SToby Isaac   sp->ops->createalldata        = PetscDualSpaceCreateAllDataDefault;
3169b4457527SToby Isaac   sp->ops->createintdata        = PetscDualSpaceCreateInteriorDataDefault;
317020cf1dd8SToby Isaac   PetscFunctionReturn(0);
317120cf1dd8SToby Isaac }
317220cf1dd8SToby Isaac 
317320cf1dd8SToby Isaac /*MC
317420cf1dd8SToby Isaac   PETSCDUALSPACELAGRANGE = "lagrange" - A PetscDualSpace object that encapsulates a dual space of pointwise evaluation functionals
317520cf1dd8SToby Isaac 
317620cf1dd8SToby Isaac   Level: intermediate
317720cf1dd8SToby Isaac 
3178db781477SPatrick Sanan .seealso: `PetscDualSpaceType`, `PetscDualSpaceCreate()`, `PetscDualSpaceSetType()`
317920cf1dd8SToby Isaac M*/
31809371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PetscDualSpaceCreate_Lagrange(PetscDualSpace sp) {
318120cf1dd8SToby Isaac   PetscDualSpace_Lag *lag;
318220cf1dd8SToby Isaac 
318320cf1dd8SToby Isaac   PetscFunctionBegin;
318420cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
3185*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&lag));
318620cf1dd8SToby Isaac   sp->data = lag;
318720cf1dd8SToby Isaac 
31883f27d899SToby Isaac   lag->tensorCell  = PETSC_FALSE;
318920cf1dd8SToby Isaac   lag->tensorSpace = PETSC_FALSE;
319020cf1dd8SToby Isaac   lag->continuous  = PETSC_TRUE;
31913f27d899SToby Isaac   lag->numCopies   = PETSC_DEFAULT;
31923f27d899SToby Isaac   lag->numNodeSkip = PETSC_DEFAULT;
31933f27d899SToby Isaac   lag->nodeType    = PETSCDTNODES_DEFAULT;
319466a6c23cSMatthew G. Knepley   lag->useMoments  = PETSC_FALSE;
319566a6c23cSMatthew G. Knepley   lag->momentOrder = 0;
319620cf1dd8SToby Isaac 
31979566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceInitialize_Lagrange(sp));
31989566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetContinuity_C", PetscDualSpaceLagrangeGetContinuity_Lagrange));
31999566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetContinuity_C", PetscDualSpaceLagrangeSetContinuity_Lagrange));
32009566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetTensor_C", PetscDualSpaceLagrangeGetTensor_Lagrange));
32019566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetTensor_C", PetscDualSpaceLagrangeSetTensor_Lagrange));
32029566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetTrimmed_C", PetscDualSpaceLagrangeGetTrimmed_Lagrange));
32039566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetTrimmed_C", PetscDualSpaceLagrangeSetTrimmed_Lagrange));
32049566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetNodeType_C", PetscDualSpaceLagrangeGetNodeType_Lagrange));
32059566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetNodeType_C", PetscDualSpaceLagrangeSetNodeType_Lagrange));
32069566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetUseMoments_C", PetscDualSpaceLagrangeGetUseMoments_Lagrange));
32079566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetUseMoments_C", PetscDualSpaceLagrangeSetUseMoments_Lagrange));
32089566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetMomentOrder_C", PetscDualSpaceLagrangeGetMomentOrder_Lagrange));
32099566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetMomentOrder_C", PetscDualSpaceLagrangeSetMomentOrder_Lagrange));
321020cf1dd8SToby Isaac   PetscFunctionReturn(0);
321120cf1dd8SToby Isaac }
3212