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 73f27d899SToby Isaac struct _n_Petsc1DNodeFamily 83f27d899SToby Isaac { 93f27d899SToby Isaac PetscInt refct; 103f27d899SToby Isaac PetscDTNodeType nodeFamily; 113f27d899SToby Isaac PetscReal gaussJacobiExp; 123f27d899SToby Isaac PetscInt nComputed; 133f27d899SToby Isaac PetscReal **nodesets; 143f27d899SToby Isaac PetscBool endpoints; 153f27d899SToby Isaac }; 163f27d899SToby Isaac 1777f1a120SToby Isaac /* users set node families for PETSCDUALSPACELAGRANGE with just the inputs to this function, but internally we create 1877f1a120SToby Isaac * an object that can cache the computations across multiple dual spaces */ 193f27d899SToby Isaac static PetscErrorCode Petsc1DNodeFamilyCreate(PetscDTNodeType family, PetscReal gaussJacobiExp, PetscBool endpoints, Petsc1DNodeFamily *nf) 203f27d899SToby Isaac { 213f27d899SToby Isaac Petsc1DNodeFamily f; 223f27d899SToby Isaac 233f27d899SToby Isaac PetscFunctionBegin; 249566063dSJacob Faibussowitsch PetscCall(PetscNew(&f)); 253f27d899SToby Isaac switch (family) { 263f27d899SToby Isaac case PETSCDTNODES_GAUSSJACOBI: 273f27d899SToby Isaac case PETSCDTNODES_EQUISPACED: 283f27d899SToby Isaac f->nodeFamily = family; 293f27d899SToby Isaac break; 303f27d899SToby Isaac default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unknown 1D node family"); 313f27d899SToby Isaac } 323f27d899SToby Isaac f->endpoints = endpoints; 333f27d899SToby Isaac f->gaussJacobiExp = 0.; 343f27d899SToby Isaac if (family == PETSCDTNODES_GAUSSJACOBI) { 352c71b3e2SJacob Faibussowitsch PetscCheckFalse(gaussJacobiExp <= -1.,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Gauss-Jacobi exponent must be > -1."); 363f27d899SToby Isaac f->gaussJacobiExp = gaussJacobiExp; 373f27d899SToby Isaac } 383f27d899SToby Isaac f->refct = 1; 393f27d899SToby Isaac *nf = f; 403f27d899SToby Isaac PetscFunctionReturn(0); 413f27d899SToby Isaac } 423f27d899SToby Isaac 433f27d899SToby Isaac static PetscErrorCode Petsc1DNodeFamilyReference(Petsc1DNodeFamily nf) 443f27d899SToby Isaac { 453f27d899SToby Isaac PetscFunctionBegin; 463f27d899SToby Isaac if (nf) nf->refct++; 473f27d899SToby Isaac PetscFunctionReturn(0); 483f27d899SToby Isaac } 493f27d899SToby Isaac 50bdb10af2SPierre Jolivet static PetscErrorCode Petsc1DNodeFamilyDestroy(Petsc1DNodeFamily *nf) 51bdb10af2SPierre Jolivet { 523f27d899SToby Isaac PetscInt i, nc; 533f27d899SToby Isaac 543f27d899SToby Isaac PetscFunctionBegin; 553f27d899SToby Isaac if (!(*nf)) PetscFunctionReturn(0); 563f27d899SToby Isaac if (--(*nf)->refct > 0) { 573f27d899SToby Isaac *nf = NULL; 583f27d899SToby Isaac PetscFunctionReturn(0); 593f27d899SToby Isaac } 603f27d899SToby Isaac nc = (*nf)->nComputed; 613f27d899SToby Isaac for (i = 0; i < nc; i++) { 629566063dSJacob Faibussowitsch PetscCall(PetscFree((*nf)->nodesets[i])); 633f27d899SToby Isaac } 649566063dSJacob Faibussowitsch PetscCall(PetscFree((*nf)->nodesets)); 659566063dSJacob Faibussowitsch PetscCall(PetscFree(*nf)); 663f27d899SToby Isaac *nf = NULL; 673f27d899SToby Isaac PetscFunctionReturn(0); 683f27d899SToby Isaac } 693f27d899SToby Isaac 703f27d899SToby Isaac static PetscErrorCode Petsc1DNodeFamilyGetNodeSets(Petsc1DNodeFamily f, PetscInt degree, PetscReal ***nodesets) 713f27d899SToby Isaac { 723f27d899SToby Isaac PetscInt nc; 733f27d899SToby Isaac 743f27d899SToby Isaac PetscFunctionBegin; 753f27d899SToby Isaac nc = f->nComputed; 763f27d899SToby Isaac if (degree >= nc) { 773f27d899SToby Isaac PetscInt i, j; 783f27d899SToby Isaac PetscReal **new_nodesets; 793f27d899SToby Isaac PetscReal *w; 803f27d899SToby Isaac 819566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(degree + 1, &new_nodesets)); 829566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(new_nodesets, f->nodesets, nc)); 839566063dSJacob Faibussowitsch PetscCall(PetscFree(f->nodesets)); 843f27d899SToby Isaac f->nodesets = new_nodesets; 859566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(degree + 1, &w)); 863f27d899SToby Isaac for (i = nc; i < degree + 1; i++) { 879566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(i + 1, &(f->nodesets[i]))); 883f27d899SToby Isaac if (!i) { 893f27d899SToby Isaac f->nodesets[i][0] = 0.5; 903f27d899SToby Isaac } else { 913f27d899SToby Isaac switch (f->nodeFamily) { 923f27d899SToby Isaac case PETSCDTNODES_EQUISPACED: 933f27d899SToby Isaac if (f->endpoints) { 943f27d899SToby Isaac for (j = 0; j <= i; j++) f->nodesets[i][j] = (PetscReal) j / (PetscReal) i; 953f27d899SToby Isaac } else { 9677f1a120SToby Isaac /* these nodes are at the centroids of the small simplices created by the equispaced nodes that include 9777f1a120SToby Isaac * the endpoints */ 983f27d899SToby Isaac for (j = 0; j <= i; j++) f->nodesets[i][j] = ((PetscReal) j + 0.5) / ((PetscReal) i + 1.); 993f27d899SToby Isaac } 1003f27d899SToby Isaac break; 1013f27d899SToby Isaac case PETSCDTNODES_GAUSSJACOBI: 1023f27d899SToby Isaac if (f->endpoints) { 1039566063dSJacob Faibussowitsch PetscCall(PetscDTGaussLobattoJacobiQuadrature(i + 1, 0., 1., f->gaussJacobiExp, f->gaussJacobiExp, f->nodesets[i], w)); 1043f27d899SToby Isaac } else { 1059566063dSJacob Faibussowitsch PetscCall(PetscDTGaussJacobiQuadrature(i + 1, 0., 1., f->gaussJacobiExp, f->gaussJacobiExp, f->nodesets[i], w)); 1063f27d899SToby Isaac } 1073f27d899SToby Isaac break; 1083f27d899SToby Isaac default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unknown 1D node family"); 1093f27d899SToby Isaac } 1103f27d899SToby Isaac } 1113f27d899SToby Isaac } 1129566063dSJacob Faibussowitsch PetscCall(PetscFree(w)); 1133f27d899SToby Isaac f->nComputed = degree + 1; 1143f27d899SToby Isaac } 1153f27d899SToby Isaac *nodesets = f->nodesets; 1163f27d899SToby Isaac PetscFunctionReturn(0); 1173f27d899SToby Isaac } 1183f27d899SToby Isaac 11977f1a120SToby Isaac /* http://arxiv.org/abs/2002.09421 for details */ 1203f27d899SToby Isaac static PetscErrorCode PetscNodeRecursive_Internal(PetscInt dim, PetscInt degree, PetscReal **nodesets, PetscInt tup[], PetscReal node[]) 1213f27d899SToby Isaac { 1223f27d899SToby Isaac PetscReal w; 1233f27d899SToby Isaac PetscInt i, j; 1243f27d899SToby Isaac 1253f27d899SToby Isaac PetscFunctionBeginHot; 1263f27d899SToby Isaac w = 0.; 1273f27d899SToby Isaac if (dim == 1) { 1283f27d899SToby Isaac node[0] = nodesets[degree][tup[0]]; 1293f27d899SToby Isaac node[1] = nodesets[degree][tup[1]]; 1303f27d899SToby Isaac } else { 1313f27d899SToby Isaac for (i = 0; i < dim + 1; i++) node[i] = 0.; 1323f27d899SToby Isaac for (i = 0; i < dim + 1; i++) { 1333f27d899SToby Isaac PetscReal wi = nodesets[degree][degree-tup[i]]; 1343f27d899SToby Isaac 1353f27d899SToby Isaac for (j = 0; j < dim+1; j++) tup[dim+1+j] = tup[j+(j>=i)]; 1369566063dSJacob Faibussowitsch PetscCall(PetscNodeRecursive_Internal(dim-1,degree-tup[i],nodesets,&tup[dim+1],&node[dim+1])); 1373f27d899SToby Isaac for (j = 0; j < dim+1; j++) node[j+(j>=i)] += wi * node[dim+1+j]; 1383f27d899SToby Isaac w += wi; 1393f27d899SToby Isaac } 1403f27d899SToby Isaac for (i = 0; i < dim+1; i++) node[i] /= w; 1413f27d899SToby Isaac } 1423f27d899SToby Isaac PetscFunctionReturn(0); 1433f27d899SToby Isaac } 1443f27d899SToby Isaac 1453f27d899SToby Isaac /* compute simplex nodes for the biunit simplex from the 1D node family */ 1463f27d899SToby Isaac static PetscErrorCode Petsc1DNodeFamilyComputeSimplexNodes(Petsc1DNodeFamily f, PetscInt dim, PetscInt degree, PetscReal points[]) 1473f27d899SToby Isaac { 1483f27d899SToby Isaac PetscInt *tup; 1493f27d899SToby Isaac PetscInt k; 1503f27d899SToby Isaac PetscInt npoints; 1513f27d899SToby Isaac PetscReal **nodesets = NULL; 1523f27d899SToby Isaac PetscInt worksize; 1533f27d899SToby Isaac PetscReal *nodework; 1543f27d899SToby Isaac PetscInt *tupwork; 1553f27d899SToby Isaac 1563f27d899SToby Isaac PetscFunctionBegin; 1572c71b3e2SJacob Faibussowitsch PetscCheckFalse(dim < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must have non-negative dimension"); 1582c71b3e2SJacob Faibussowitsch PetscCheckFalse(degree < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must have non-negative degree"); 1593f27d899SToby Isaac if (!dim) PetscFunctionReturn(0); 1609566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(dim+2, &tup)); 1613f27d899SToby Isaac k = 0; 1629566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(degree + dim, dim, &npoints)); 1639566063dSJacob Faibussowitsch PetscCall(Petsc1DNodeFamilyGetNodeSets(f, degree, &nodesets)); 1643f27d899SToby Isaac worksize = ((dim + 2) * (dim + 3)) / 2; 1659566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(worksize, &nodework, worksize, &tupwork)); 16677f1a120SToby Isaac /* loop over the tuples of length dim with sum at most degree */ 1673f27d899SToby Isaac for (k = 0; k < npoints; k++) { 1683f27d899SToby Isaac PetscInt i; 1693f27d899SToby Isaac 17077f1a120SToby Isaac /* turn thm into tuples of length dim + 1 with sum equal to degree (barycentric indice) */ 1713f27d899SToby Isaac tup[0] = degree; 1723f27d899SToby Isaac for (i = 0; i < dim; i++) { 1733f27d899SToby Isaac tup[0] -= tup[i+1]; 1743f27d899SToby Isaac } 1753f27d899SToby Isaac switch(f->nodeFamily) { 1763f27d899SToby Isaac case PETSCDTNODES_EQUISPACED: 17777f1a120SToby Isaac /* compute equispaces nodes on the unit reference triangle */ 1783f27d899SToby Isaac if (f->endpoints) { 1793f27d899SToby Isaac for (i = 0; i < dim; i++) { 1803f27d899SToby Isaac points[dim*k + i] = (PetscReal) tup[i+1] / (PetscReal) degree; 1813f27d899SToby Isaac } 1823f27d899SToby Isaac } else { 1833f27d899SToby Isaac for (i = 0; i < dim; i++) { 18477f1a120SToby Isaac /* these nodes are at the centroids of the small simplices created by the equispaced nodes that include 18577f1a120SToby Isaac * the endpoints */ 1863f27d899SToby Isaac points[dim*k + i] = ((PetscReal) tup[i+1] + 1./(dim+1.)) / (PetscReal) (degree + 1.); 1873f27d899SToby Isaac } 1883f27d899SToby Isaac } 1893f27d899SToby Isaac break; 1903f27d899SToby Isaac default: 19177f1a120SToby Isaac /* compute equispaces nodes on the barycentric reference triangle (the trace on the first dim dimensions are the 19277f1a120SToby Isaac * unit reference triangle nodes */ 1933f27d899SToby Isaac for (i = 0; i < dim + 1; i++) tupwork[i] = tup[i]; 1949566063dSJacob Faibussowitsch PetscCall(PetscNodeRecursive_Internal(dim, degree, nodesets, tupwork, nodework)); 1953f27d899SToby Isaac for (i = 0; i < dim; i++) points[dim*k + i] = nodework[i + 1]; 1963f27d899SToby Isaac break; 1973f27d899SToby Isaac } 1989566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceLatticePointLexicographic_Internal(dim, degree, &tup[1])); 1993f27d899SToby Isaac } 2003f27d899SToby Isaac /* map from unit simplex to biunit simplex */ 2013f27d899SToby Isaac for (k = 0; k < npoints * dim; k++) points[k] = points[k] * 2. - 1.; 2029566063dSJacob Faibussowitsch PetscCall(PetscFree2(nodework, tupwork)); 2039566063dSJacob Faibussowitsch PetscCall(PetscFree(tup)); 2043f27d899SToby Isaac PetscFunctionReturn(0); 2053f27d899SToby Isaac } 2063f27d899SToby Isaac 20777f1a120SToby 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 20877f1a120SToby Isaac * on that mesh point, we have to be careful about getting/adding everything in the right place. 20977f1a120SToby Isaac * 21077f1a120SToby Isaac * With nodal dofs like PETSCDUALSPACELAGRANGE makes, the general approach to calculate the value of dofs associate 21177f1a120SToby Isaac * with a node A is 21277f1a120SToby Isaac * - transform the node locations x(A) by the map that takes the mesh point to its reorientation, x' = phi(x(A)) 21377f1a120SToby Isaac * - figure out which node was originally at the location of the transformed point, A' = idx(x') 21477f1a120SToby Isaac * - if the dofs are not scalars, figure out how to represent the transformed dofs in terms of the basis 21577f1a120SToby Isaac * of dofs at A' (using pushforward/pullback rules) 21677f1a120SToby Isaac * 21777f1a120SToby Isaac * The one sticky point with this approach is the "A' = idx(x')" step: trying to go from real valued coordinates 21877f1a120SToby Isaac * back to indices. I don't want to rely on floating point tolerances. Additionally, PETSCDUALSPACELAGRANGE may 21977f1a120SToby Isaac * eventually support quasi-Lagrangian dofs, which could involve quadrature at multiple points, so the location "x(A)" 22077f1a120SToby Isaac * would be ambiguous. 22177f1a120SToby Isaac * 22277f1a120SToby Isaac * So each dof gets an integer value coordinate (nodeIdx in the structure below). The choice of integer coordinates 22377f1a120SToby Isaac * is somewhat arbitrary, as long as all of the relevant symmetries of the mesh point correspond to *permutations* of 22477f1a120SToby Isaac * the integer coordinates, which do not depend on numerical precision. 22577f1a120SToby Isaac * 22677f1a120SToby Isaac * So 22777f1a120SToby Isaac * 22877f1a120SToby Isaac * - DMPlexGetTransitiveClosure_Internal() tells me how an orientation turns into a permutation of the vertices of a 22977f1a120SToby Isaac * mesh point 23077f1a120SToby Isaac * - The permutation of the vertices, and the nodeIdx values assigned to them, tells what permutation in index space 23177f1a120SToby Isaac * is associated with the orientation 23277f1a120SToby Isaac * - I uses that permutation to get xi' = phi(xi(A)), the integer coordinate of the transformed dof 23377f1a120SToby Isaac * - I can without numerical issues compute A' = idx(xi') 23477f1a120SToby Isaac * 23577f1a120SToby Isaac * Here are some examples of how the process works 23677f1a120SToby Isaac * 23777f1a120SToby Isaac * - With a triangle: 23877f1a120SToby Isaac * 23977f1a120SToby Isaac * The triangle has the following integer coordinates for vertices, taken from the barycentric triangle 24077f1a120SToby Isaac * 24177f1a120SToby Isaac * closure order 2 24277f1a120SToby Isaac * nodeIdx (0,0,1) 24377f1a120SToby Isaac * \ 24477f1a120SToby Isaac * + 24577f1a120SToby Isaac * |\ 24677f1a120SToby Isaac * | \ 24777f1a120SToby Isaac * | \ 24877f1a120SToby Isaac * | \ closure order 1 24977f1a120SToby Isaac * | \ / nodeIdx (0,1,0) 25077f1a120SToby Isaac * +-----+ 25177f1a120SToby Isaac * \ 25277f1a120SToby Isaac * closure order 0 25377f1a120SToby Isaac * nodeIdx (1,0,0) 25477f1a120SToby Isaac * 25577f1a120SToby Isaac * If I do DMPlexGetTransitiveClosure_Internal() with orientation 1, the vertices would appear 25677f1a120SToby Isaac * in the order (1, 2, 0) 25777f1a120SToby Isaac * 25877f1a120SToby 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 25977f1a120SToby Isaac * see 26077f1a120SToby Isaac * 26177f1a120SToby Isaac * orientation 0 | orientation 1 26277f1a120SToby Isaac * 26377f1a120SToby Isaac * [0] (1,0,0) [1] (0,1,0) 26477f1a120SToby Isaac * [1] (0,1,0) [2] (0,0,1) 26577f1a120SToby Isaac * [2] (0,0,1) [0] (1,0,0) 26677f1a120SToby Isaac * A B 26777f1a120SToby Isaac * 26877f1a120SToby Isaac * In other words, B is the result of a row permutation of A. But, there is also 26977f1a120SToby Isaac * a column permutation that accomplishes the same result, (2,0,1). 27077f1a120SToby Isaac * 27177f1a120SToby Isaac * So if a dof has nodeIdx coordinate (a,b,c), after the transformation its nodeIdx coordinate 27277f1a120SToby Isaac * is (c,a,b), and the transformed degree of freedom will be a linear combination of dofs 27377f1a120SToby Isaac * that originally had coordinate (c,a,b). 27477f1a120SToby Isaac * 27577f1a120SToby Isaac * - With a quadrilateral: 27677f1a120SToby Isaac * 27777f1a120SToby Isaac * The quadrilateral has the following integer coordinates for vertices, taken from concatenating barycentric 27877f1a120SToby Isaac * coordinates for two segments: 27977f1a120SToby Isaac * 28077f1a120SToby Isaac * closure order 3 closure order 2 28177f1a120SToby Isaac * nodeIdx (1,0,0,1) nodeIdx (0,1,0,1) 28277f1a120SToby Isaac * \ / 28377f1a120SToby Isaac * +----+ 28477f1a120SToby Isaac * | | 28577f1a120SToby Isaac * | | 28677f1a120SToby Isaac * +----+ 28777f1a120SToby Isaac * / \ 28877f1a120SToby Isaac * closure order 0 closure order 1 28977f1a120SToby Isaac * nodeIdx (1,0,1,0) nodeIdx (0,1,1,0) 29077f1a120SToby Isaac * 29177f1a120SToby Isaac * If I do DMPlexGetTransitiveClosure_Internal() with orientation 1, the vertices would appear 29277f1a120SToby Isaac * in the order (1, 2, 3, 0) 29377f1a120SToby Isaac * 29477f1a120SToby Isaac * If I list the nodeIdx of each vertex in closure order for orientation 0 (0, 1, 2, 3) and 29577f1a120SToby Isaac * orientation 1 (1, 2, 3, 0), I see 29677f1a120SToby Isaac * 29777f1a120SToby Isaac * orientation 0 | orientation 1 29877f1a120SToby Isaac * 29977f1a120SToby Isaac * [0] (1,0,1,0) [1] (0,1,1,0) 30077f1a120SToby Isaac * [1] (0,1,1,0) [2] (0,1,0,1) 30177f1a120SToby Isaac * [2] (0,1,0,1) [3] (1,0,0,1) 30277f1a120SToby Isaac * [3] (1,0,0,1) [0] (1,0,1,0) 30377f1a120SToby Isaac * A B 30477f1a120SToby Isaac * 30577f1a120SToby Isaac * The column permutation that accomplishes the same result is (3,2,0,1). 30677f1a120SToby Isaac * 30777f1a120SToby Isaac * So if a dof has nodeIdx coordinate (a,b,c,d), after the transformation its nodeIdx coordinate 30877f1a120SToby Isaac * is (d,c,a,b), and the transformed degree of freedom will be a linear combination of dofs 30977f1a120SToby Isaac * that originally had coordinate (d,c,a,b). 31077f1a120SToby Isaac * 31177f1a120SToby Isaac * Previously PETSCDUALSPACELAGRANGE had hardcoded symmetries for the triangle and quadrilateral, 31277f1a120SToby Isaac * but this approach will work for any polytope, such as the wedge (triangular prism). 31377f1a120SToby Isaac */ 3143f27d899SToby Isaac struct _n_PetscLagNodeIndices 3153f27d899SToby Isaac { 3163f27d899SToby Isaac PetscInt refct; 3173f27d899SToby Isaac PetscInt nodeIdxDim; 3183f27d899SToby Isaac PetscInt nodeVecDim; 3193f27d899SToby Isaac PetscInt nNodes; 3203f27d899SToby Isaac PetscInt *nodeIdx; /* for each node an index of size nodeIdxDim */ 3213f27d899SToby Isaac PetscReal *nodeVec; /* for each node a vector of size nodeVecDim */ 3223f27d899SToby Isaac PetscInt *perm; /* if these are vertices, perm takes DMPlex point index to closure order; 3233f27d899SToby Isaac if these are nodes, perm lists nodes in index revlex order */ 3243f27d899SToby Isaac }; 3253f27d899SToby Isaac 32677f1a120SToby Isaac /* this is just here so I can access the values in tests/ex1.c outside the library */ 3273f27d899SToby Isaac PetscErrorCode PetscLagNodeIndicesGetData_Internal(PetscLagNodeIndices ni, PetscInt *nodeIdxDim, PetscInt *nodeVecDim, PetscInt *nNodes, const PetscInt *nodeIdx[], const PetscReal *nodeVec[]) 3283f27d899SToby Isaac { 3293f27d899SToby Isaac PetscFunctionBegin; 3303f27d899SToby Isaac *nodeIdxDim = ni->nodeIdxDim; 3313f27d899SToby Isaac *nodeVecDim = ni->nodeVecDim; 3323f27d899SToby Isaac *nNodes = ni->nNodes; 3333f27d899SToby Isaac *nodeIdx = ni->nodeIdx; 3343f27d899SToby Isaac *nodeVec = ni->nodeVec; 3353f27d899SToby Isaac PetscFunctionReturn(0); 3363f27d899SToby Isaac } 3373f27d899SToby Isaac 3383f27d899SToby Isaac static PetscErrorCode PetscLagNodeIndicesReference(PetscLagNodeIndices ni) 3393f27d899SToby Isaac { 3403f27d899SToby Isaac PetscFunctionBegin; 3413f27d899SToby Isaac if (ni) ni->refct++; 3423f27d899SToby Isaac PetscFunctionReturn(0); 3433f27d899SToby Isaac } 3443f27d899SToby Isaac 3451f440fbeSToby Isaac static PetscErrorCode PetscLagNodeIndicesDuplicate(PetscLagNodeIndices ni, PetscLagNodeIndices *niNew) 3461f440fbeSToby Isaac { 3471f440fbeSToby Isaac PetscFunctionBegin; 3489566063dSJacob Faibussowitsch PetscCall(PetscNew(niNew)); 3491f440fbeSToby Isaac (*niNew)->refct = 1; 3501f440fbeSToby Isaac (*niNew)->nodeIdxDim = ni->nodeIdxDim; 3511f440fbeSToby Isaac (*niNew)->nodeVecDim = ni->nodeVecDim; 3521f440fbeSToby Isaac (*niNew)->nNodes = ni->nNodes; 3539566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ni->nNodes * ni->nodeIdxDim, &((*niNew)->nodeIdx))); 3549566063dSJacob Faibussowitsch PetscCall(PetscArraycpy((*niNew)->nodeIdx, ni->nodeIdx, ni->nNodes * ni->nodeIdxDim)); 3559566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ni->nNodes * ni->nodeVecDim, &((*niNew)->nodeVec))); 3569566063dSJacob Faibussowitsch PetscCall(PetscArraycpy((*niNew)->nodeVec, ni->nodeVec, ni->nNodes * ni->nodeVecDim)); 3571f440fbeSToby Isaac (*niNew)->perm = NULL; 3581f440fbeSToby Isaac PetscFunctionReturn(0); 3591f440fbeSToby Isaac } 3601f440fbeSToby Isaac 361bdb10af2SPierre Jolivet static PetscErrorCode PetscLagNodeIndicesDestroy(PetscLagNodeIndices *ni) 362bdb10af2SPierre Jolivet { 3633f27d899SToby Isaac PetscFunctionBegin; 3643f27d899SToby Isaac if (!(*ni)) PetscFunctionReturn(0); 3653f27d899SToby Isaac if (--(*ni)->refct > 0) { 3663f27d899SToby Isaac *ni = NULL; 3673f27d899SToby Isaac PetscFunctionReturn(0); 3683f27d899SToby Isaac } 3699566063dSJacob Faibussowitsch PetscCall(PetscFree((*ni)->nodeIdx)); 3709566063dSJacob Faibussowitsch PetscCall(PetscFree((*ni)->nodeVec)); 3719566063dSJacob Faibussowitsch PetscCall(PetscFree((*ni)->perm)); 3729566063dSJacob Faibussowitsch PetscCall(PetscFree(*ni)); 3733f27d899SToby Isaac *ni = NULL; 3743f27d899SToby Isaac PetscFunctionReturn(0); 3753f27d899SToby Isaac } 3763f27d899SToby Isaac 37777f1a120SToby Isaac /* The vertices are given nodeIdx coordinates (e.g. the corners of the barycentric triangle). Those coordinates are 37877f1a120SToby Isaac * in some other order, and to understand the effect of different symmetries, we need them to be in closure order. 37977f1a120SToby Isaac * 38077f1a120SToby Isaac * If sortIdx is PETSC_FALSE, the coordinates are already in revlex order, otherwise we must sort them 38177f1a120SToby Isaac * to that order before we do the real work of this function, which is 38277f1a120SToby Isaac * 38377f1a120SToby Isaac * - mark the vertices in closure order 38477f1a120SToby Isaac * - sort them in revlex order 38577f1a120SToby Isaac * - use the resulting permutation to list the vertex coordinates in closure order 38677f1a120SToby Isaac */ 3873f27d899SToby Isaac static PetscErrorCode PetscLagNodeIndicesComputeVertexOrder(DM dm, PetscLagNodeIndices ni, PetscBool sortIdx) 3883f27d899SToby Isaac { 3893f27d899SToby Isaac PetscInt v, w, vStart, vEnd, c, d; 3903f27d899SToby Isaac PetscInt nVerts; 3913f27d899SToby Isaac PetscInt closureSize = 0; 3923f27d899SToby Isaac PetscInt *closure = NULL; 3933f27d899SToby Isaac PetscInt *closureOrder; 3943f27d899SToby Isaac PetscInt *invClosureOrder; 3953f27d899SToby Isaac PetscInt *revlexOrder; 3963f27d899SToby Isaac PetscInt *newNodeIdx; 3973f27d899SToby Isaac PetscInt dim; 3983f27d899SToby Isaac Vec coordVec; 3993f27d899SToby Isaac const PetscScalar *coords; 4003f27d899SToby Isaac 4013f27d899SToby Isaac PetscFunctionBegin; 4029566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 4039566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 4043f27d899SToby Isaac nVerts = vEnd - vStart; 4059566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nVerts, &closureOrder)); 4069566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nVerts, &invClosureOrder)); 4079566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nVerts, &revlexOrder)); 40877f1a120SToby Isaac if (sortIdx) { /* bubble sort nodeIdx into revlex order */ 4093f27d899SToby Isaac PetscInt nodeIdxDim = ni->nodeIdxDim; 4103f27d899SToby Isaac PetscInt *idxOrder; 4113f27d899SToby Isaac 4129566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nVerts * nodeIdxDim, &newNodeIdx)); 4139566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nVerts, &idxOrder)); 4143f27d899SToby Isaac for (v = 0; v < nVerts; v++) idxOrder[v] = v; 4153f27d899SToby Isaac for (v = 0; v < nVerts; v++) { 4163f27d899SToby Isaac for (w = v + 1; w < nVerts; w++) { 4173f27d899SToby Isaac const PetscInt *iv = &(ni->nodeIdx[idxOrder[v] * nodeIdxDim]); 4183f27d899SToby Isaac const PetscInt *iw = &(ni->nodeIdx[idxOrder[w] * nodeIdxDim]); 4193f27d899SToby Isaac PetscInt diff = 0; 4203f27d899SToby Isaac 4213f27d899SToby Isaac for (d = nodeIdxDim - 1; d >= 0; d--) if ((diff = (iv[d] - iw[d]))) break; 4223f27d899SToby Isaac if (diff > 0) { 4233f27d899SToby Isaac PetscInt swap = idxOrder[v]; 4243f27d899SToby Isaac 4253f27d899SToby Isaac idxOrder[v] = idxOrder[w]; 4263f27d899SToby Isaac idxOrder[w] = swap; 4273f27d899SToby Isaac } 4283f27d899SToby Isaac } 4293f27d899SToby Isaac } 4303f27d899SToby Isaac for (v = 0; v < nVerts; v++) { 4313f27d899SToby Isaac for (d = 0; d < nodeIdxDim; d++) { 4323f27d899SToby Isaac newNodeIdx[v * ni->nodeIdxDim + d] = ni->nodeIdx[idxOrder[v] * nodeIdxDim + d]; 4333f27d899SToby Isaac } 4343f27d899SToby Isaac } 4359566063dSJacob Faibussowitsch PetscCall(PetscFree(ni->nodeIdx)); 4363f27d899SToby Isaac ni->nodeIdx = newNodeIdx; 4373f27d899SToby Isaac newNodeIdx = NULL; 4389566063dSJacob Faibussowitsch PetscCall(PetscFree(idxOrder)); 4393f27d899SToby Isaac } 4409566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, 0, PETSC_TRUE, &closureSize, &closure)); 4413f27d899SToby Isaac c = closureSize - nVerts; 4423f27d899SToby Isaac for (v = 0; v < nVerts; v++) closureOrder[v] = closure[2 * (c + v)] - vStart; 4433f27d899SToby Isaac for (v = 0; v < nVerts; v++) invClosureOrder[closureOrder[v]] = v; 4449566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, 0, PETSC_TRUE, &closureSize, &closure)); 4459566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordVec)); 4469566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(coordVec, &coords)); 4473f27d899SToby Isaac /* bubble sort closure vertices by coordinates in revlex order */ 4483f27d899SToby Isaac for (v = 0; v < nVerts; v++) revlexOrder[v] = v; 4493f27d899SToby Isaac for (v = 0; v < nVerts; v++) { 4503f27d899SToby Isaac for (w = v + 1; w < nVerts; w++) { 4513f27d899SToby Isaac const PetscScalar *cv = &coords[closureOrder[revlexOrder[v]] * dim]; 4523f27d899SToby Isaac const PetscScalar *cw = &coords[closureOrder[revlexOrder[w]] * dim]; 4533f27d899SToby Isaac PetscReal diff = 0; 4543f27d899SToby Isaac 4553f27d899SToby Isaac for (d = dim - 1; d >= 0; d--) if ((diff = PetscRealPart(cv[d] - cw[d])) != 0.) break; 4563f27d899SToby Isaac if (diff > 0.) { 4573f27d899SToby Isaac PetscInt swap = revlexOrder[v]; 4583f27d899SToby Isaac 4593f27d899SToby Isaac revlexOrder[v] = revlexOrder[w]; 4603f27d899SToby Isaac revlexOrder[w] = swap; 4613f27d899SToby Isaac } 4623f27d899SToby Isaac } 4633f27d899SToby Isaac } 4649566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(coordVec, &coords)); 4659566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ni->nodeIdxDim * nVerts, &newNodeIdx)); 4663f27d899SToby Isaac /* reorder nodeIdx to be in closure order */ 4673f27d899SToby Isaac for (v = 0; v < nVerts; v++) { 4683f27d899SToby Isaac for (d = 0; d < ni->nodeIdxDim; d++) { 4693f27d899SToby Isaac newNodeIdx[revlexOrder[v] * ni->nodeIdxDim + d] = ni->nodeIdx[v * ni->nodeIdxDim + d]; 4703f27d899SToby Isaac } 4713f27d899SToby Isaac } 4729566063dSJacob Faibussowitsch PetscCall(PetscFree(ni->nodeIdx)); 4733f27d899SToby Isaac ni->nodeIdx = newNodeIdx; 4743f27d899SToby Isaac ni->perm = invClosureOrder; 4759566063dSJacob Faibussowitsch PetscCall(PetscFree(revlexOrder)); 4769566063dSJacob Faibussowitsch PetscCall(PetscFree(closureOrder)); 4773f27d899SToby Isaac PetscFunctionReturn(0); 4783f27d899SToby Isaac } 4793f27d899SToby Isaac 48077f1a120SToby Isaac /* the coordinates of the simplex vertices are the corners of the barycentric simplex. 48177f1a120SToby Isaac * When we stack them on top of each other in revlex order, they look like the identity matrix */ 4823f27d899SToby Isaac static PetscErrorCode PetscLagNodeIndicesCreateSimplexVertices(DM dm, PetscLagNodeIndices *nodeIndices) 4833f27d899SToby Isaac { 4843f27d899SToby Isaac PetscLagNodeIndices ni; 4853f27d899SToby Isaac PetscInt dim, d; 4863f27d899SToby Isaac 4873f27d899SToby Isaac PetscFunctionBegin; 4889566063dSJacob Faibussowitsch PetscCall(PetscNew(&ni)); 4899566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 4903f27d899SToby Isaac ni->nodeIdxDim = dim + 1; 4913f27d899SToby Isaac ni->nodeVecDim = 0; 4923f27d899SToby Isaac ni->nNodes = dim + 1; 4933f27d899SToby Isaac ni->refct = 1; 4949566063dSJacob Faibussowitsch PetscCall(PetscCalloc1((dim + 1)*(dim + 1), &(ni->nodeIdx))); 4953f27d899SToby Isaac for (d = 0; d < dim + 1; d++) ni->nodeIdx[d*(dim + 2)] = 1; 4969566063dSJacob Faibussowitsch PetscCall(PetscLagNodeIndicesComputeVertexOrder(dm, ni, PETSC_FALSE)); 4973f27d899SToby Isaac *nodeIndices = ni; 4983f27d899SToby Isaac PetscFunctionReturn(0); 4993f27d899SToby Isaac } 5003f27d899SToby Isaac 50177f1a120SToby Isaac /* A polytope that is a tensor product of a facet and a segment. 50277f1a120SToby Isaac * We take whatever coordinate system was being used for the facet 5031f440fbeSToby Isaac * and we concatenate the barycentric coordinates for the vertices 50477f1a120SToby Isaac * at the end of the segment, (1,0) and (0,1), to get a coordinate 50577f1a120SToby Isaac * system for the tensor product element */ 5063f27d899SToby Isaac static PetscErrorCode PetscLagNodeIndicesCreateTensorVertices(DM dm, PetscLagNodeIndices facetni, PetscLagNodeIndices *nodeIndices) 5073f27d899SToby Isaac { 5083f27d899SToby Isaac PetscLagNodeIndices ni; 5093f27d899SToby Isaac PetscInt nodeIdxDim, subNodeIdxDim = facetni->nodeIdxDim; 5103f27d899SToby Isaac PetscInt nVerts, nSubVerts = facetni->nNodes; 5113f27d899SToby Isaac PetscInt dim, d, e, f, g; 5123f27d899SToby Isaac 5133f27d899SToby Isaac PetscFunctionBegin; 5149566063dSJacob Faibussowitsch PetscCall(PetscNew(&ni)); 5159566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 5163f27d899SToby Isaac ni->nodeIdxDim = nodeIdxDim = subNodeIdxDim + 2; 5173f27d899SToby Isaac ni->nodeVecDim = 0; 5183f27d899SToby Isaac ni->nNodes = nVerts = 2 * nSubVerts; 5193f27d899SToby Isaac ni->refct = 1; 5209566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(nodeIdxDim * nVerts, &(ni->nodeIdx))); 5213f27d899SToby Isaac for (f = 0, d = 0; d < 2; d++) { 5223f27d899SToby Isaac for (e = 0; e < nSubVerts; e++, f++) { 5233f27d899SToby Isaac for (g = 0; g < subNodeIdxDim; g++) { 5243f27d899SToby Isaac ni->nodeIdx[f * nodeIdxDim + g] = facetni->nodeIdx[e * subNodeIdxDim + g]; 5253f27d899SToby Isaac } 5263f27d899SToby Isaac ni->nodeIdx[f * nodeIdxDim + subNodeIdxDim] = (1 - d); 5273f27d899SToby Isaac ni->nodeIdx[f * nodeIdxDim + subNodeIdxDim + 1] = d; 5283f27d899SToby Isaac } 5293f27d899SToby Isaac } 5309566063dSJacob Faibussowitsch PetscCall(PetscLagNodeIndicesComputeVertexOrder(dm, ni, PETSC_TRUE)); 5313f27d899SToby Isaac *nodeIndices = ni; 5323f27d899SToby Isaac PetscFunctionReturn(0); 5333f27d899SToby Isaac } 5343f27d899SToby Isaac 53577f1a120SToby Isaac /* This helps us compute symmetries, and it also helps us compute coordinates for dofs that are being pushed 53677f1a120SToby Isaac * forward from a boundary mesh point. 53777f1a120SToby Isaac * 53877f1a120SToby Isaac * Input: 53977f1a120SToby Isaac * 54077f1a120SToby Isaac * dm - the target reference cell where we want new coordinates and dof directions to be valid 54177f1a120SToby Isaac * vert - the vertex coordinate system for the target reference cell 54277f1a120SToby Isaac * p - the point in the target reference cell that the dofs are coming from 54377f1a120SToby Isaac * vertp - the vertex coordinate system for p's reference cell 54477f1a120SToby Isaac * ornt - the resulting coordinates and dof vectors will be for p under this orientation 54577f1a120SToby Isaac * nodep - the node coordinates and dof vectors in p's reference cell 54677f1a120SToby Isaac * formDegree - the form degree that the dofs transform as 54777f1a120SToby Isaac * 54877f1a120SToby Isaac * Output: 54977f1a120SToby Isaac * 55077f1a120SToby Isaac * pfNodeIdx - the node coordinates for p's dofs, in the dm reference cell, from the ornt perspective 55177f1a120SToby Isaac * pfNodeVec - the node dof vectors for p's dofs, in the dm reference cell, from the ornt perspective 55277f1a120SToby Isaac */ 5533f27d899SToby Isaac static PetscErrorCode PetscLagNodeIndicesPushForward(DM dm, PetscLagNodeIndices vert, PetscInt p, PetscLagNodeIndices vertp, PetscLagNodeIndices nodep, PetscInt ornt, PetscInt formDegree, PetscInt pfNodeIdx[], PetscReal pfNodeVec[]) 5543f27d899SToby Isaac { 5553f27d899SToby Isaac PetscInt *closureVerts; 5563f27d899SToby Isaac PetscInt closureSize = 0; 5573f27d899SToby Isaac PetscInt *closure = NULL; 5583f27d899SToby Isaac PetscInt dim, pdim, c, i, j, k, n, v, vStart, vEnd; 5593f27d899SToby Isaac PetscInt nSubVert = vertp->nNodes; 5603f27d899SToby Isaac PetscInt nodeIdxDim = vert->nodeIdxDim; 5613f27d899SToby Isaac PetscInt subNodeIdxDim = vertp->nodeIdxDim; 5623f27d899SToby Isaac PetscInt nNodes = nodep->nNodes; 5633f27d899SToby Isaac const PetscInt *vertIdx = vert->nodeIdx; 5643f27d899SToby Isaac const PetscInt *subVertIdx = vertp->nodeIdx; 5653f27d899SToby Isaac const PetscInt *nodeIdx = nodep->nodeIdx; 5663f27d899SToby Isaac const PetscReal *nodeVec = nodep->nodeVec; 5673f27d899SToby Isaac PetscReal *J, *Jstar; 5683f27d899SToby Isaac PetscReal detJ; 5693f27d899SToby Isaac PetscInt depth, pdepth, Nk, pNk; 5703f27d899SToby Isaac Vec coordVec; 5713f27d899SToby Isaac PetscScalar *newCoords = NULL; 5723f27d899SToby Isaac const PetscScalar *oldCoords = NULL; 5733f27d899SToby Isaac 5743f27d899SToby Isaac PetscFunctionBegin; 5759566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 5769566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepth(dm, &depth)); 5779566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordVec)); 5789566063dSJacob Faibussowitsch PetscCall(DMPlexGetPointDepth(dm, p, &pdepth)); 5793f27d899SToby Isaac pdim = pdepth != depth ? pdepth != 0 ? pdepth : 0 : dim; 5809566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 5819566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(dm, nSubVert, MPIU_INT, &closureVerts)); 5829566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure_Internal(dm, p, ornt, PETSC_TRUE, &closureSize, &closure)); 5833f27d899SToby Isaac c = closureSize - nSubVert; 5843f27d899SToby Isaac /* we want which cell closure indices the closure of this point corresponds to */ 5853f27d899SToby Isaac for (v = 0; v < nSubVert; v++) closureVerts[v] = vert->perm[closure[2 * (c + v)] - vStart]; 5869566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, p, PETSC_TRUE, &closureSize, &closure)); 5873f27d899SToby Isaac /* push forward indices */ 5883f27d899SToby Isaac for (i = 0; i < nodeIdxDim; i++) { /* for every component of the target index space */ 5893f27d899SToby Isaac /* check if this is a component that all vertices around this point have in common */ 5903f27d899SToby Isaac for (j = 1; j < nSubVert; j++) { 5913f27d899SToby Isaac if (vertIdx[closureVerts[j] * nodeIdxDim + i] != vertIdx[closureVerts[0] * nodeIdxDim + i]) break; 5923f27d899SToby Isaac } 5933f27d899SToby Isaac if (j == nSubVert) { /* all vertices have this component in common, directly copy to output */ 5943f27d899SToby Isaac PetscInt val = vertIdx[closureVerts[0] * nodeIdxDim + i]; 5953f27d899SToby Isaac for (n = 0; n < nNodes; n++) pfNodeIdx[n * nodeIdxDim + i] = val; 5963f27d899SToby Isaac } else { 5973f27d899SToby Isaac PetscInt subi = -1; 5983f27d899SToby Isaac /* there must be a component in vertp that looks the same */ 5993f27d899SToby Isaac for (k = 0; k < subNodeIdxDim; k++) { 6003f27d899SToby Isaac for (j = 0; j < nSubVert; j++) { 6013f27d899SToby Isaac if (vertIdx[closureVerts[j] * nodeIdxDim + i] != subVertIdx[j * subNodeIdxDim + k]) break; 6023f27d899SToby Isaac } 6033f27d899SToby Isaac if (j == nSubVert) { 6043f27d899SToby Isaac subi = k; 6053f27d899SToby Isaac break; 6063f27d899SToby Isaac } 6073f27d899SToby Isaac } 6082c71b3e2SJacob Faibussowitsch PetscCheckFalse(subi < 0,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Did not find matching coordinate"); 60977f1a120SToby Isaac /* that component in the vertp system becomes component i in the vert system for each dof */ 6103f27d899SToby Isaac for (n = 0; n < nNodes; n++) pfNodeIdx[n * nodeIdxDim + i] = nodeIdx[n * subNodeIdxDim + subi]; 6113f27d899SToby Isaac } 6123f27d899SToby Isaac } 6133f27d899SToby Isaac /* push forward vectors */ 6149566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(dm, dim * dim, MPIU_REAL, &J)); 61577f1a120SToby Isaac if (ornt != 0) { /* temporarily change the coordinate vector so 61677f1a120SToby Isaac DMPlexComputeCellGeometryAffineFEM gives us the Jacobian we want */ 6173f27d899SToby Isaac PetscInt closureSize2 = 0; 6183f27d899SToby Isaac PetscInt *closure2 = NULL; 6193f27d899SToby Isaac 6209566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure_Internal(dm, p, 0, PETSC_TRUE, &closureSize2, &closure2)); 6219566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dim * nSubVert, &newCoords)); 6229566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(coordVec, &oldCoords)); 6233f27d899SToby Isaac for (v = 0; v < nSubVert; v++) { 6243f27d899SToby Isaac PetscInt d; 6253f27d899SToby Isaac for (d = 0; d < dim; d++) { 6263f27d899SToby Isaac newCoords[(closure2[2 * (c + v)] - vStart) * dim + d] = oldCoords[closureVerts[v] * dim + d]; 6273f27d899SToby Isaac } 6283f27d899SToby Isaac } 6299566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(coordVec, &oldCoords)); 6309566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, p, PETSC_TRUE, &closureSize2, &closure2)); 6319566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(coordVec, newCoords)); 6323f27d899SToby Isaac } 6339566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryAffineFEM(dm, p, NULL, J, NULL, &detJ)); 6343f27d899SToby Isaac if (ornt != 0) { 6359566063dSJacob Faibussowitsch PetscCall(VecResetArray(coordVec)); 6369566063dSJacob Faibussowitsch PetscCall(PetscFree(newCoords)); 6373f27d899SToby Isaac } 6389566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(dm, nSubVert, MPIU_INT, &closureVerts)); 6393f27d899SToby Isaac /* compactify */ 6403f27d899SToby Isaac for (i = 0; i < dim; i++) for (j = 0; j < pdim; j++) J[i * pdim + j] = J[i * dim + j]; 64177f1a120SToby Isaac /* We have the Jacobian mapping the point's reference cell to this reference cell: 64277f1a120SToby Isaac * pulling back a function to the point and applying the dof is what we want, 64377f1a120SToby Isaac * so we get the pullback matrix and multiply the dof by that matrix on the right */ 6449566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(formDegree), &Nk)); 6459566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(pdim, PetscAbsInt(formDegree), &pNk)); 6469566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(dm, pNk * Nk, MPIU_REAL, &Jstar)); 6479566063dSJacob Faibussowitsch PetscCall(PetscDTAltVPullbackMatrix(pdim, dim, J, formDegree, Jstar)); 6483f27d899SToby Isaac for (n = 0; n < nNodes; n++) { 6493f27d899SToby Isaac for (i = 0; i < Nk; i++) { 6503f27d899SToby Isaac PetscReal val = 0.; 6515efe5503SToby Isaac for (j = 0; j < pNk; j++) val += nodeVec[n * pNk + j] * Jstar[j * Nk + i]; 6523f27d899SToby Isaac pfNodeVec[n * Nk + i] = val; 6533f27d899SToby Isaac } 6543f27d899SToby Isaac } 6559566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(dm, pNk * Nk, MPIU_REAL, &Jstar)); 6569566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(dm, dim * dim, MPIU_REAL, &J)); 6573f27d899SToby Isaac PetscFunctionReturn(0); 6583f27d899SToby Isaac } 6593f27d899SToby Isaac 66077f1a120SToby Isaac /* given to sets of nodes, take the tensor product, where the product of the dof indices is concatenation and the 66177f1a120SToby Isaac * product of the dof vectors is the wedge product */ 6623f27d899SToby Isaac static PetscErrorCode PetscLagNodeIndicesTensor(PetscLagNodeIndices tracei, PetscInt dimT, PetscInt kT, PetscLagNodeIndices fiberi, PetscInt dimF, PetscInt kF, PetscLagNodeIndices *nodeIndices) 6633f27d899SToby Isaac { 6643f27d899SToby Isaac PetscInt dim = dimT + dimF; 6653f27d899SToby Isaac PetscInt nodeIdxDim, nNodes; 6663f27d899SToby Isaac PetscInt formDegree = kT + kF; 6673f27d899SToby Isaac PetscInt Nk, NkT, NkF; 6683f27d899SToby Isaac PetscInt MkT, MkF; 6693f27d899SToby Isaac PetscLagNodeIndices ni; 6703f27d899SToby Isaac PetscInt i, j, l; 6713f27d899SToby Isaac PetscReal *projF, *projT; 6723f27d899SToby Isaac PetscReal *projFstar, *projTstar; 6733f27d899SToby Isaac PetscReal *workF, *workF2, *workT, *workT2, *work, *work2; 6743f27d899SToby Isaac PetscReal *wedgeMat; 6753f27d899SToby Isaac PetscReal sign; 6763f27d899SToby Isaac 6773f27d899SToby Isaac PetscFunctionBegin; 6789566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(formDegree), &Nk)); 6799566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(dimT, PetscAbsInt(kT), &NkT)); 6809566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(dimF, PetscAbsInt(kF), &NkF)); 6819566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(kT), &MkT)); 6829566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(kF), &MkF)); 6839566063dSJacob Faibussowitsch PetscCall(PetscNew(&ni)); 6843f27d899SToby Isaac ni->nodeIdxDim = nodeIdxDim = tracei->nodeIdxDim + fiberi->nodeIdxDim; 6853f27d899SToby Isaac ni->nodeVecDim = Nk; 6863f27d899SToby Isaac ni->nNodes = nNodes = tracei->nNodes * fiberi->nNodes; 6873f27d899SToby Isaac ni->refct = 1; 6889566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nNodes * nodeIdxDim, &(ni->nodeIdx))); 6893f27d899SToby Isaac /* first concatenate the indices */ 6903f27d899SToby Isaac for (l = 0, j = 0; j < fiberi->nNodes; j++) { 6913f27d899SToby Isaac for (i = 0; i < tracei->nNodes; i++, l++) { 6923f27d899SToby Isaac PetscInt m, n = 0; 6933f27d899SToby Isaac 6943f27d899SToby Isaac for (m = 0; m < tracei->nodeIdxDim; m++) ni->nodeIdx[l * nodeIdxDim + n++] = tracei->nodeIdx[i * tracei->nodeIdxDim + m]; 6953f27d899SToby Isaac for (m = 0; m < fiberi->nodeIdxDim; m++) ni->nodeIdx[l * nodeIdxDim + n++] = fiberi->nodeIdx[j * fiberi->nodeIdxDim + m]; 6963f27d899SToby Isaac } 6973f27d899SToby Isaac } 6983f27d899SToby Isaac 6993f27d899SToby Isaac /* now wedge together the push-forward vectors */ 7009566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nNodes * Nk, &(ni->nodeVec))); 7019566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(dimT*dim, &projT, dimF*dim, &projF)); 7023f27d899SToby Isaac for (i = 0; i < dimT; i++) projT[i * (dim + 1)] = 1.; 7033f27d899SToby Isaac for (i = 0; i < dimF; i++) projF[i * (dim + dimT + 1) + dimT] = 1.; 7049566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(MkT*NkT, &projTstar, MkF*NkF, &projFstar)); 7059566063dSJacob Faibussowitsch PetscCall(PetscDTAltVPullbackMatrix(dim, dimT, projT, kT, projTstar)); 7069566063dSJacob Faibussowitsch PetscCall(PetscDTAltVPullbackMatrix(dim, dimF, projF, kF, projFstar)); 7079566063dSJacob Faibussowitsch PetscCall(PetscMalloc6(MkT, &workT, MkT, &workT2, MkF, &workF, MkF, &workF2, Nk, &work, Nk, &work2)); 7089566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Nk * MkT, &wedgeMat)); 7093f27d899SToby Isaac sign = (PetscAbsInt(kT * kF) & 1) ? -1. : 1.; 7103f27d899SToby Isaac for (l = 0, j = 0; j < fiberi->nNodes; j++) { 7113f27d899SToby Isaac PetscInt d, e; 7123f27d899SToby Isaac 7133f27d899SToby Isaac /* push forward fiber k-form */ 7143f27d899SToby Isaac for (d = 0; d < MkF; d++) { 7153f27d899SToby Isaac PetscReal val = 0.; 7163f27d899SToby Isaac for (e = 0; e < NkF; e++) val += projFstar[d * NkF + e] * fiberi->nodeVec[j * NkF + e]; 7173f27d899SToby Isaac workF[d] = val; 7183f27d899SToby Isaac } 7193f27d899SToby Isaac /* Hodge star to proper form if necessary */ 7203f27d899SToby Isaac if (kF < 0) { 7213f27d899SToby Isaac for (d = 0; d < MkF; d++) workF2[d] = workF[d]; 7229566063dSJacob Faibussowitsch PetscCall(PetscDTAltVStar(dim, PetscAbsInt(kF), 1, workF2, workF)); 7233f27d899SToby Isaac } 7243f27d899SToby Isaac /* Compute the matrix that wedges this form with one of the trace k-form */ 7259566063dSJacob Faibussowitsch PetscCall(PetscDTAltVWedgeMatrix(dim, PetscAbsInt(kF), PetscAbsInt(kT), workF, wedgeMat)); 7263f27d899SToby Isaac for (i = 0; i < tracei->nNodes; i++, l++) { 7273f27d899SToby Isaac /* push forward trace k-form */ 7283f27d899SToby Isaac for (d = 0; d < MkT; d++) { 7293f27d899SToby Isaac PetscReal val = 0.; 7303f27d899SToby Isaac for (e = 0; e < NkT; e++) val += projTstar[d * NkT + e] * tracei->nodeVec[i * NkT + e]; 7313f27d899SToby Isaac workT[d] = val; 7323f27d899SToby Isaac } 7333f27d899SToby Isaac /* Hodge star to proper form if necessary */ 7343f27d899SToby Isaac if (kT < 0) { 7353f27d899SToby Isaac for (d = 0; d < MkT; d++) workT2[d] = workT[d]; 7369566063dSJacob Faibussowitsch PetscCall(PetscDTAltVStar(dim, PetscAbsInt(kT), 1, workT2, workT)); 7373f27d899SToby Isaac } 7383f27d899SToby Isaac /* compute the wedge product of the push-forward trace form and firer forms */ 7393f27d899SToby Isaac for (d = 0; d < Nk; d++) { 7403f27d899SToby Isaac PetscReal val = 0.; 7413f27d899SToby Isaac for (e = 0; e < MkT; e++) val += wedgeMat[d * MkT + e] * workT[e]; 7423f27d899SToby Isaac work[d] = val; 7433f27d899SToby Isaac } 7443f27d899SToby Isaac /* inverse Hodge star from proper form if necessary */ 7453f27d899SToby Isaac if (formDegree < 0) { 7463f27d899SToby Isaac for (d = 0; d < Nk; d++) work2[d] = work[d]; 7479566063dSJacob Faibussowitsch PetscCall(PetscDTAltVStar(dim, PetscAbsInt(formDegree), -1, work2, work)); 7483f27d899SToby Isaac } 7493f27d899SToby Isaac /* insert into the array (adjusting for sign) */ 7503f27d899SToby Isaac for (d = 0; d < Nk; d++) ni->nodeVec[l * Nk + d] = sign * work[d]; 7513f27d899SToby Isaac } 7523f27d899SToby Isaac } 7539566063dSJacob Faibussowitsch PetscCall(PetscFree(wedgeMat)); 7549566063dSJacob Faibussowitsch PetscCall(PetscFree6(workT, workT2, workF, workF2, work, work2)); 7559566063dSJacob Faibussowitsch PetscCall(PetscFree2(projTstar, projFstar)); 7569566063dSJacob Faibussowitsch PetscCall(PetscFree2(projT, projF)); 7573f27d899SToby Isaac *nodeIndices = ni; 7583f27d899SToby Isaac PetscFunctionReturn(0); 7593f27d899SToby Isaac } 7603f27d899SToby Isaac 76177f1a120SToby Isaac /* simple union of two sets of nodes */ 7623f27d899SToby Isaac static PetscErrorCode PetscLagNodeIndicesMerge(PetscLagNodeIndices niA, PetscLagNodeIndices niB, PetscLagNodeIndices *nodeIndices) 7633f27d899SToby Isaac { 7643f27d899SToby Isaac PetscLagNodeIndices ni; 7653f27d899SToby Isaac PetscInt nodeIdxDim, nodeVecDim, nNodes; 7663f27d899SToby Isaac 7673f27d899SToby Isaac PetscFunctionBegin; 7689566063dSJacob Faibussowitsch PetscCall(PetscNew(&ni)); 7693f27d899SToby Isaac ni->nodeIdxDim = nodeIdxDim = niA->nodeIdxDim; 7702c71b3e2SJacob Faibussowitsch PetscCheckFalse(niB->nodeIdxDim != nodeIdxDim,PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Cannot merge PetscLagNodeIndices with different nodeIdxDim"); 7713f27d899SToby Isaac ni->nodeVecDim = nodeVecDim = niA->nodeVecDim; 7722c71b3e2SJacob Faibussowitsch PetscCheckFalse(niB->nodeVecDim != nodeVecDim,PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Cannot merge PetscLagNodeIndices with different nodeVecDim"); 7733f27d899SToby Isaac ni->nNodes = nNodes = niA->nNodes + niB->nNodes; 7743f27d899SToby Isaac ni->refct = 1; 7759566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nNodes * nodeIdxDim, &(ni->nodeIdx))); 7769566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nNodes * nodeVecDim, &(ni->nodeVec))); 7779566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(ni->nodeIdx, niA->nodeIdx, niA->nNodes * nodeIdxDim)); 7789566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(ni->nodeVec, niA->nodeVec, niA->nNodes * nodeVecDim)); 7799566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(&(ni->nodeIdx[niA->nNodes * nodeIdxDim]), niB->nodeIdx, niB->nNodes * nodeIdxDim)); 7809566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(&(ni->nodeVec[niA->nNodes * nodeVecDim]), niB->nodeVec, niB->nNodes * nodeVecDim)); 7813f27d899SToby Isaac *nodeIndices = ni; 7823f27d899SToby Isaac PetscFunctionReturn(0); 7833f27d899SToby Isaac } 7843f27d899SToby Isaac 7853f27d899SToby Isaac #define PETSCTUPINTCOMPREVLEX(N) \ 786d6a2e6abSJacob Faibussowitsch static int PetscConcat_(PetscTupIntCompRevlex_,N)(const void *a, const void *b) \ 7873f27d899SToby Isaac { \ 7883f27d899SToby Isaac const PetscInt *A = (const PetscInt *) a; \ 7893f27d899SToby Isaac const PetscInt *B = (const PetscInt *) b; \ 7903f27d899SToby Isaac int i; \ 7913f27d899SToby Isaac PetscInt diff = 0; \ 7923f27d899SToby Isaac for (i = 0; i < N; i++) { \ 7933f27d899SToby Isaac diff = A[N - i] - B[N - i]; \ 7943f27d899SToby Isaac if (diff) break; \ 7953f27d899SToby Isaac } \ 7963f27d899SToby Isaac return (diff <= 0) ? (diff < 0) ? -1 : 0 : 1; \ 7973f27d899SToby Isaac } 7983f27d899SToby Isaac 7993f27d899SToby Isaac PETSCTUPINTCOMPREVLEX(3) 8003f27d899SToby Isaac PETSCTUPINTCOMPREVLEX(4) 8013f27d899SToby Isaac PETSCTUPINTCOMPREVLEX(5) 8023f27d899SToby Isaac PETSCTUPINTCOMPREVLEX(6) 8033f27d899SToby Isaac PETSCTUPINTCOMPREVLEX(7) 8043f27d899SToby Isaac 8053f27d899SToby Isaac static int PetscTupIntCompRevlex_N(const void *a, const void *b) 8063f27d899SToby Isaac { 8073f27d899SToby Isaac const PetscInt *A = (const PetscInt *) a; 8083f27d899SToby Isaac const PetscInt *B = (const PetscInt *) b; 8093f27d899SToby Isaac int i; 8103f27d899SToby Isaac int N = A[0]; 8113f27d899SToby Isaac PetscInt diff = 0; 8123f27d899SToby Isaac for (i = 0; i < N; i++) { 8133f27d899SToby Isaac diff = A[N - i] - B[N - i]; 8143f27d899SToby Isaac if (diff) break; 8153f27d899SToby Isaac } 8163f27d899SToby Isaac return (diff <= 0) ? (diff < 0) ? -1 : 0 : 1; 8173f27d899SToby Isaac } 8183f27d899SToby Isaac 81977f1a120SToby Isaac /* The nodes are not necessarily in revlex order wrt nodeIdx: get the permutation 82077f1a120SToby Isaac * that puts them in that order */ 8213f27d899SToby Isaac static PetscErrorCode PetscLagNodeIndicesGetPermutation(PetscLagNodeIndices ni, PetscInt *perm[]) 8223f27d899SToby Isaac { 8233f27d899SToby Isaac PetscFunctionBegin; 8243f27d899SToby Isaac if (!(ni->perm)) { 8253f27d899SToby Isaac PetscInt *sorter; 8263f27d899SToby Isaac PetscInt m = ni->nNodes; 8273f27d899SToby Isaac PetscInt nodeIdxDim = ni->nodeIdxDim; 8283f27d899SToby Isaac PetscInt i, j, k, l; 8293f27d899SToby Isaac PetscInt *prm; 8303f27d899SToby Isaac int (*comp) (const void *, const void *); 8313f27d899SToby Isaac 8329566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((nodeIdxDim + 2) * m, &sorter)); 8333f27d899SToby Isaac for (k = 0, l = 0, i = 0; i < m; i++) { 8343f27d899SToby Isaac sorter[k++] = nodeIdxDim + 1; 8353f27d899SToby Isaac sorter[k++] = i; 8363f27d899SToby Isaac for (j = 0; j < nodeIdxDim; j++) sorter[k++] = ni->nodeIdx[l++]; 8373f27d899SToby Isaac } 8383f27d899SToby Isaac switch (nodeIdxDim) { 8393f27d899SToby Isaac case 2: 8403f27d899SToby Isaac comp = PetscTupIntCompRevlex_3; 8413f27d899SToby Isaac break; 8423f27d899SToby Isaac case 3: 8433f27d899SToby Isaac comp = PetscTupIntCompRevlex_4; 8443f27d899SToby Isaac break; 8453f27d899SToby Isaac case 4: 8463f27d899SToby Isaac comp = PetscTupIntCompRevlex_5; 8473f27d899SToby Isaac break; 8483f27d899SToby Isaac case 5: 8493f27d899SToby Isaac comp = PetscTupIntCompRevlex_6; 8503f27d899SToby Isaac break; 8513f27d899SToby Isaac case 6: 8523f27d899SToby Isaac comp = PetscTupIntCompRevlex_7; 8533f27d899SToby Isaac break; 8543f27d899SToby Isaac default: 8553f27d899SToby Isaac comp = PetscTupIntCompRevlex_N; 8563f27d899SToby Isaac break; 8573f27d899SToby Isaac } 8583f27d899SToby Isaac qsort(sorter, m, (nodeIdxDim + 2) * sizeof(PetscInt), comp); 8599566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m, &prm)); 8603f27d899SToby Isaac for (i = 0; i < m; i++) prm[i] = sorter[(nodeIdxDim + 2) * i + 1]; 8613f27d899SToby Isaac ni->perm = prm; 8629566063dSJacob Faibussowitsch PetscCall(PetscFree(sorter)); 8633f27d899SToby Isaac } 8643f27d899SToby Isaac *perm = ni->perm; 8653f27d899SToby Isaac PetscFunctionReturn(0); 8663f27d899SToby Isaac } 86720cf1dd8SToby Isaac 8686f905325SMatthew G. Knepley static PetscErrorCode PetscDualSpaceDestroy_Lagrange(PetscDualSpace sp) 86920cf1dd8SToby Isaac { 87020cf1dd8SToby Isaac PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 87120cf1dd8SToby Isaac 87220cf1dd8SToby Isaac PetscFunctionBegin; 8733f27d899SToby Isaac if (lag->symperms) { 8743f27d899SToby Isaac PetscInt **selfSyms = lag->symperms[0]; 8756f905325SMatthew G. Knepley 8766f905325SMatthew G. Knepley if (selfSyms) { 8776f905325SMatthew G. Knepley PetscInt i, **allocated = &selfSyms[-lag->selfSymOff]; 8786f905325SMatthew G. Knepley 8796f905325SMatthew G. Knepley for (i = 0; i < lag->numSelfSym; i++) { 8809566063dSJacob Faibussowitsch PetscCall(PetscFree(allocated[i])); 8816f905325SMatthew G. Knepley } 8829566063dSJacob Faibussowitsch PetscCall(PetscFree(allocated)); 8836f905325SMatthew G. Knepley } 8849566063dSJacob Faibussowitsch PetscCall(PetscFree(lag->symperms)); 8856f905325SMatthew G. Knepley } 8863f27d899SToby Isaac if (lag->symflips) { 8873f27d899SToby Isaac PetscScalar **selfSyms = lag->symflips[0]; 8883f27d899SToby Isaac 8893f27d899SToby Isaac if (selfSyms) { 8903f27d899SToby Isaac PetscInt i; 8913f27d899SToby Isaac PetscScalar **allocated = &selfSyms[-lag->selfSymOff]; 8923f27d899SToby Isaac 8933f27d899SToby Isaac for (i = 0; i < lag->numSelfSym; i++) { 8949566063dSJacob Faibussowitsch PetscCall(PetscFree(allocated[i])); 8956f905325SMatthew G. Knepley } 8969566063dSJacob Faibussowitsch PetscCall(PetscFree(allocated)); 8973f27d899SToby Isaac } 8989566063dSJacob Faibussowitsch PetscCall(PetscFree(lag->symflips)); 8993f27d899SToby Isaac } 9009566063dSJacob Faibussowitsch PetscCall(Petsc1DNodeFamilyDestroy(&(lag->nodeFamily))); 9019566063dSJacob Faibussowitsch PetscCall(PetscLagNodeIndicesDestroy(&(lag->vertIndices))); 9029566063dSJacob Faibussowitsch PetscCall(PetscLagNodeIndicesDestroy(&(lag->intNodeIndices))); 9039566063dSJacob Faibussowitsch PetscCall(PetscLagNodeIndicesDestroy(&(lag->allNodeIndices))); 9049566063dSJacob Faibussowitsch PetscCall(PetscFree(lag)); 9059566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetContinuity_C", NULL)); 9069566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetContinuity_C", NULL)); 9079566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetTensor_C", NULL)); 9089566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetTensor_C", NULL)); 9099566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetTrimmed_C", NULL)); 9109566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetTrimmed_C", NULL)); 9119566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetNodeType_C", NULL)); 9129566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetNodeType_C", NULL)); 9139566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetUseMoments_C", NULL)); 9149566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetUseMoments_C", NULL)); 9159566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetMomentOrder_C", NULL)); 9169566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetMomentOrder_C", NULL)); 91720cf1dd8SToby Isaac PetscFunctionReturn(0); 91820cf1dd8SToby Isaac } 91920cf1dd8SToby Isaac 9206f905325SMatthew G. Knepley static PetscErrorCode PetscDualSpaceLagrangeView_Ascii(PetscDualSpace sp, PetscViewer viewer) 92120cf1dd8SToby Isaac { 92220cf1dd8SToby Isaac PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 92320cf1dd8SToby Isaac 92420cf1dd8SToby Isaac PetscFunctionBegin; 9259566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%s %s%sLagrange dual space\n", lag->continuous ? "Continuous" : "Discontinuous", lag->tensorSpace ? "tensor " : "", lag->trimmed ? "trimmed " : "")); 92620cf1dd8SToby Isaac PetscFunctionReturn(0); 92720cf1dd8SToby Isaac } 92820cf1dd8SToby Isaac 9296f905325SMatthew G. Knepley static PetscErrorCode PetscDualSpaceView_Lagrange(PetscDualSpace sp, PetscViewer viewer) 93020cf1dd8SToby Isaac { 9316f905325SMatthew G. Knepley PetscBool iascii; 9326f905325SMatthew G. Knepley 93320cf1dd8SToby Isaac PetscFunctionBegin; 9346f905325SMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 9356f905325SMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 9369566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii)); 9379566063dSJacob Faibussowitsch if (iascii) PetscCall(PetscDualSpaceLagrangeView_Ascii(sp, viewer)); 93820cf1dd8SToby Isaac PetscFunctionReturn(0); 93920cf1dd8SToby Isaac } 94020cf1dd8SToby Isaac 9416f905325SMatthew G. Knepley static PetscErrorCode PetscDualSpaceSetFromOptions_Lagrange(PetscOptionItems *PetscOptionsObject,PetscDualSpace sp) 94220cf1dd8SToby Isaac { 9433f27d899SToby Isaac PetscBool continuous, tensor, trimmed, flg, flg2, flg3; 9443f27d899SToby Isaac PetscDTNodeType nodeType; 9453f27d899SToby Isaac PetscReal nodeExponent; 94666a6c23cSMatthew G. Knepley PetscInt momentOrder; 94766a6c23cSMatthew G. Knepley PetscBool nodeEndpoints, useMoments; 9486f905325SMatthew G. Knepley 9496f905325SMatthew G. Knepley PetscFunctionBegin; 9509566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceLagrangeGetContinuity(sp, &continuous)); 9519566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceLagrangeGetTensor(sp, &tensor)); 9529566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceLagrangeGetTrimmed(sp, &trimmed)); 9539566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceLagrangeGetNodeType(sp, &nodeType, &nodeEndpoints, &nodeExponent)); 9543f27d899SToby Isaac if (nodeType == PETSCDTNODES_DEFAULT) nodeType = PETSCDTNODES_GAUSSJACOBI; 9559566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceLagrangeGetUseMoments(sp, &useMoments)); 9569566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceLagrangeGetMomentOrder(sp, &momentOrder)); 9579566063dSJacob Faibussowitsch PetscCall(PetscOptionsHead(PetscOptionsObject,"PetscDualSpace Lagrange Options")); 9589566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-petscdualspace_lagrange_continuity", "Flag for continuous element", "PetscDualSpaceLagrangeSetContinuity", continuous, &continuous, &flg)); 9599566063dSJacob Faibussowitsch if (flg) PetscCall(PetscDualSpaceLagrangeSetContinuity(sp, continuous)); 9609566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-petscdualspace_lagrange_tensor", "Flag for tensor dual space", "PetscDualSpaceLagrangeSetTensor", tensor, &tensor, &flg)); 9619566063dSJacob Faibussowitsch if (flg) PetscCall(PetscDualSpaceLagrangeSetTensor(sp, tensor)); 9629566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-petscdualspace_lagrange_trimmed", "Flag for trimmed dual space", "PetscDualSpaceLagrangeSetTrimmed", trimmed, &trimmed, &flg)); 9639566063dSJacob Faibussowitsch if (flg) PetscCall(PetscDualSpaceLagrangeSetTrimmed(sp, trimmed)); 9649566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-petscdualspace_lagrange_node_type", "Lagrange node location type", "PetscDualSpaceLagrangeSetNodeType", PetscDTNodeTypes, (PetscEnum)nodeType, (PetscEnum *)&nodeType, &flg)); 9659566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-petscdualspace_lagrange_node_endpoints", "Flag for nodes that include endpoints", "PetscDualSpaceLagrangeSetNodeType", nodeEndpoints, &nodeEndpoints, &flg2)); 9663f27d899SToby Isaac flg3 = PETSC_FALSE; 9673f27d899SToby Isaac if (nodeType == PETSCDTNODES_GAUSSJACOBI) { 9689566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-petscdualspace_lagrange_node_exponent", "Gauss-Jacobi weight function exponent", "PetscDualSpaceLagrangeSetNodeType", nodeExponent, &nodeExponent, &flg3)); 9693f27d899SToby Isaac } 9709566063dSJacob Faibussowitsch if (flg || flg2 || flg3) PetscCall(PetscDualSpaceLagrangeSetNodeType(sp, nodeType, nodeEndpoints, nodeExponent)); 9719566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-petscdualspace_lagrange_use_moments", "Use moments (where appropriate) for functionals", "PetscDualSpaceLagrangeSetUseMoments", useMoments, &useMoments, &flg)); 9729566063dSJacob Faibussowitsch if (flg) PetscCall(PetscDualSpaceLagrangeSetUseMoments(sp, useMoments)); 9739566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-petscdualspace_lagrange_moment_order", "Quadrature order for moment functionals", "PetscDualSpaceLagrangeSetMomentOrder", momentOrder, &momentOrder, &flg)); 9749566063dSJacob Faibussowitsch if (flg) PetscCall(PetscDualSpaceLagrangeSetMomentOrder(sp, momentOrder)); 9759566063dSJacob Faibussowitsch PetscCall(PetscOptionsTail()); 9766f905325SMatthew G. Knepley PetscFunctionReturn(0); 9776f905325SMatthew G. Knepley } 9786f905325SMatthew G. Knepley 979b4457527SToby Isaac static PetscErrorCode PetscDualSpaceDuplicate_Lagrange(PetscDualSpace sp, PetscDualSpace spNew) 9806f905325SMatthew G. Knepley { 9813f27d899SToby Isaac PetscBool cont, tensor, trimmed, boundary; 9823f27d899SToby Isaac PetscDTNodeType nodeType; 9833f27d899SToby Isaac PetscReal exponent; 9843f27d899SToby Isaac PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 9856f905325SMatthew G. Knepley 9866f905325SMatthew G. Knepley PetscFunctionBegin; 9879566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceLagrangeGetContinuity(sp, &cont)); 9889566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceLagrangeSetContinuity(spNew, cont)); 9899566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceLagrangeGetTensor(sp, &tensor)); 9909566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceLagrangeSetTensor(spNew, tensor)); 9919566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceLagrangeGetTrimmed(sp, &trimmed)); 9929566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceLagrangeSetTrimmed(spNew, trimmed)); 9939566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceLagrangeGetNodeType(sp, &nodeType, &boundary, &exponent)); 9949566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceLagrangeSetNodeType(spNew, nodeType, boundary, exponent)); 9953f27d899SToby Isaac if (lag->nodeFamily) { 9963f27d899SToby Isaac PetscDualSpace_Lag *lagnew = (PetscDualSpace_Lag *) spNew->data; 9973f27d899SToby Isaac 9989566063dSJacob Faibussowitsch PetscCall(Petsc1DNodeFamilyReference(lag->nodeFamily)); 9993f27d899SToby Isaac lagnew->nodeFamily = lag->nodeFamily; 10003f27d899SToby Isaac } 10016f905325SMatthew G. Knepley PetscFunctionReturn(0); 10026f905325SMatthew G. Knepley } 10036f905325SMatthew G. Knepley 100477f1a120SToby Isaac /* for making tensor product spaces: take a dual space and product a segment space that has all the same 100577f1a120SToby Isaac * specifications (trimmed, continuous, order, node set), except for the form degree */ 10063f27d899SToby Isaac static PetscErrorCode PetscDualSpaceCreateEdgeSubspace_Lagrange(PetscDualSpace sp, PetscInt order, PetscInt k, PetscInt Nc, PetscBool interiorOnly, PetscDualSpace *bdsp) 10076f905325SMatthew G. Knepley { 10083f27d899SToby Isaac DM K; 10093f27d899SToby Isaac PetscDualSpace_Lag *newlag; 10106f905325SMatthew G. Knepley 10116f905325SMatthew G. Knepley PetscFunctionBegin; 10129566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceDuplicate(sp,bdsp)); 10139566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetFormDegree(*bdsp, k)); 10149566063dSJacob Faibussowitsch PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, DMPolytopeTypeSimpleShape(1, PETSC_TRUE), &K)); 10159566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetDM(*bdsp, K)); 10169566063dSJacob Faibussowitsch PetscCall(DMDestroy(&K)); 10179566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetOrder(*bdsp, order)); 10189566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetNumComponents(*bdsp, Nc)); 10193f27d899SToby Isaac newlag = (PetscDualSpace_Lag *) (*bdsp)->data; 10203f27d899SToby Isaac newlag->interiorOnly = interiorOnly; 10219566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetUp(*bdsp)); 10223f27d899SToby Isaac PetscFunctionReturn(0); 10236f905325SMatthew G. Knepley } 10243f27d899SToby Isaac 10253f27d899SToby Isaac /* just the points, weights aren't handled */ 10263f27d899SToby Isaac static PetscErrorCode PetscQuadratureCreateTensor(PetscQuadrature trace, PetscQuadrature fiber, PetscQuadrature *product) 10273f27d899SToby Isaac { 10283f27d899SToby Isaac PetscInt dimTrace, dimFiber; 10293f27d899SToby Isaac PetscInt numPointsTrace, numPointsFiber; 10303f27d899SToby Isaac PetscInt dim, numPoints; 10313f27d899SToby Isaac const PetscReal *pointsTrace; 10323f27d899SToby Isaac const PetscReal *pointsFiber; 10333f27d899SToby Isaac PetscReal *points; 10343f27d899SToby Isaac PetscInt i, j, k, p; 10353f27d899SToby Isaac 10363f27d899SToby Isaac PetscFunctionBegin; 10379566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(trace, &dimTrace, NULL, &numPointsTrace, &pointsTrace, NULL)); 10389566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(fiber, &dimFiber, NULL, &numPointsFiber, &pointsFiber, NULL)); 10393f27d899SToby Isaac dim = dimTrace + dimFiber; 10403f27d899SToby Isaac numPoints = numPointsFiber * numPointsTrace; 10419566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numPoints * dim, &points)); 10423f27d899SToby Isaac for (p = 0, j = 0; j < numPointsFiber; j++) { 10433f27d899SToby Isaac for (i = 0; i < numPointsTrace; i++, p++) { 10443f27d899SToby Isaac for (k = 0; k < dimTrace; k++) points[p * dim + k] = pointsTrace[i * dimTrace + k]; 10453f27d899SToby Isaac for (k = 0; k < dimFiber; k++) points[p * dim + dimTrace + k] = pointsFiber[j * dimFiber + k]; 10463f27d899SToby Isaac } 10473f27d899SToby Isaac } 10489566063dSJacob Faibussowitsch PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, product)); 10499566063dSJacob Faibussowitsch PetscCall(PetscQuadratureSetData(*product, dim, 0, numPoints, points, NULL)); 10503f27d899SToby Isaac PetscFunctionReturn(0); 10513f27d899SToby Isaac } 10523f27d899SToby Isaac 105377f1a120SToby Isaac /* Kronecker tensor product where matrix is considered a matrix of k-forms, so that 105477f1a120SToby Isaac * the entries in the product matrix are wedge products of the entries in the original matrices */ 10553f27d899SToby Isaac static PetscErrorCode MatTensorAltV(Mat trace, Mat fiber, PetscInt dimTrace, PetscInt kTrace, PetscInt dimFiber, PetscInt kFiber, Mat *product) 10563f27d899SToby Isaac { 10573f27d899SToby Isaac PetscInt mTrace, nTrace, mFiber, nFiber, m, n, k, i, j, l; 10583f27d899SToby Isaac PetscInt dim, NkTrace, NkFiber, Nk; 10593f27d899SToby Isaac PetscInt dT, dF; 10603f27d899SToby Isaac PetscInt *nnzTrace, *nnzFiber, *nnz; 10613f27d899SToby Isaac PetscInt iT, iF, jT, jF, il, jl; 10623f27d899SToby Isaac PetscReal *workT, *workT2, *workF, *workF2, *work, *workstar; 10633f27d899SToby Isaac PetscReal *projT, *projF; 10643f27d899SToby Isaac PetscReal *projTstar, *projFstar; 10653f27d899SToby Isaac PetscReal *wedgeMat; 10663f27d899SToby Isaac PetscReal sign; 10673f27d899SToby Isaac PetscScalar *workS; 10683f27d899SToby Isaac Mat prod; 10693f27d899SToby Isaac /* this produces dof groups that look like the identity */ 10703f27d899SToby Isaac 10713f27d899SToby Isaac PetscFunctionBegin; 10729566063dSJacob Faibussowitsch PetscCall(MatGetSize(trace, &mTrace, &nTrace)); 10739566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(dimTrace, PetscAbsInt(kTrace), &NkTrace)); 10742c71b3e2SJacob Faibussowitsch PetscCheckFalse(nTrace % NkTrace,PETSC_COMM_SELF, PETSC_ERR_PLIB, "point value space of trace matrix is not a multiple of k-form size"); 10759566063dSJacob Faibussowitsch PetscCall(MatGetSize(fiber, &mFiber, &nFiber)); 10769566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(dimFiber, PetscAbsInt(kFiber), &NkFiber)); 10772c71b3e2SJacob Faibussowitsch PetscCheckFalse(nFiber % NkFiber,PETSC_COMM_SELF, PETSC_ERR_PLIB, "point value space of fiber matrix is not a multiple of k-form size"); 10789566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(mTrace, &nnzTrace, mFiber, &nnzFiber)); 10793f27d899SToby Isaac for (i = 0; i < mTrace; i++) { 10809566063dSJacob Faibussowitsch PetscCall(MatGetRow(trace, i, &(nnzTrace[i]), NULL, NULL)); 10812c71b3e2SJacob Faibussowitsch PetscCheckFalse(nnzTrace[i] % NkTrace,PETSC_COMM_SELF, PETSC_ERR_PLIB, "nonzeros in trace matrix are not in k-form size blocks"); 10823f27d899SToby Isaac } 10833f27d899SToby Isaac for (i = 0; i < mFiber; i++) { 10849566063dSJacob Faibussowitsch PetscCall(MatGetRow(fiber, i, &(nnzFiber[i]), NULL, NULL)); 10852c71b3e2SJacob Faibussowitsch PetscCheckFalse(nnzFiber[i] % NkFiber,PETSC_COMM_SELF, PETSC_ERR_PLIB, "nonzeros in fiber matrix are not in k-form size blocks"); 10863f27d899SToby Isaac } 10873f27d899SToby Isaac dim = dimTrace + dimFiber; 10883f27d899SToby Isaac k = kFiber + kTrace; 10899566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(k), &Nk)); 10903f27d899SToby Isaac m = mTrace * mFiber; 10919566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m, &nnz)); 10923f27d899SToby Isaac for (l = 0, j = 0; j < mFiber; j++) for (i = 0; i < mTrace; i++, l++) nnz[l] = (nnzTrace[i] / NkTrace) * (nnzFiber[j] / NkFiber) * Nk; 10933f27d899SToby Isaac n = (nTrace / NkTrace) * (nFiber / NkFiber) * Nk; 10949566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, m, n, 0, nnz, &prod)); 10959566063dSJacob Faibussowitsch PetscCall(PetscFree(nnz)); 10969566063dSJacob Faibussowitsch PetscCall(PetscFree2(nnzTrace,nnzFiber)); 10973f27d899SToby Isaac /* reasoning about which points each dof needs depends on having zeros computed at points preserved */ 10989566063dSJacob Faibussowitsch PetscCall(MatSetOption(prod, MAT_IGNORE_ZERO_ENTRIES, PETSC_FALSE)); 10993f27d899SToby Isaac /* compute pullbacks */ 11009566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(kTrace), &dT)); 11019566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(kFiber), &dF)); 11029566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(dimTrace * dim, &projT, dimFiber * dim, &projF, dT * NkTrace, &projTstar, dF * NkFiber, &projFstar)); 11039566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(projT, dimTrace * dim)); 11043f27d899SToby Isaac for (i = 0; i < dimTrace; i++) projT[i * (dim + 1)] = 1.; 11059566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(projF, dimFiber * dim)); 11063f27d899SToby Isaac for (i = 0; i < dimFiber; i++) projF[i * (dim + 1) + dimTrace] = 1.; 11079566063dSJacob Faibussowitsch PetscCall(PetscDTAltVPullbackMatrix(dim, dimTrace, projT, kTrace, projTstar)); 11089566063dSJacob Faibussowitsch PetscCall(PetscDTAltVPullbackMatrix(dim, dimFiber, projF, kFiber, projFstar)); 11099566063dSJacob Faibussowitsch PetscCall(PetscMalloc5(dT, &workT, dF, &workF, Nk, &work, Nk, &workstar, Nk, &workS)); 11109566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(dT, &workT2, dF, &workF2)); 11119566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Nk * dT, &wedgeMat)); 11123f27d899SToby Isaac sign = (PetscAbsInt(kTrace * kFiber) & 1) ? -1. : 1.; 11133f27d899SToby Isaac for (i = 0, iF = 0; iF < mFiber; iF++) { 11143f27d899SToby Isaac PetscInt ncolsF, nformsF; 11153f27d899SToby Isaac const PetscInt *colsF; 11163f27d899SToby Isaac const PetscScalar *valsF; 11173f27d899SToby Isaac 11189566063dSJacob Faibussowitsch PetscCall(MatGetRow(fiber, iF, &ncolsF, &colsF, &valsF)); 11193f27d899SToby Isaac nformsF = ncolsF / NkFiber; 11203f27d899SToby Isaac for (iT = 0; iT < mTrace; iT++, i++) { 11213f27d899SToby Isaac PetscInt ncolsT, nformsT; 11223f27d899SToby Isaac const PetscInt *colsT; 11233f27d899SToby Isaac const PetscScalar *valsT; 11243f27d899SToby Isaac 11259566063dSJacob Faibussowitsch PetscCall(MatGetRow(trace, iT, &ncolsT, &colsT, &valsT)); 11263f27d899SToby Isaac nformsT = ncolsT / NkTrace; 11273f27d899SToby Isaac for (j = 0, jF = 0; jF < nformsF; jF++) { 11283f27d899SToby Isaac PetscInt colF = colsF[jF * NkFiber] / NkFiber; 11293f27d899SToby Isaac 11303f27d899SToby Isaac for (il = 0; il < dF; il++) { 11313f27d899SToby Isaac PetscReal val = 0.; 11323f27d899SToby Isaac for (jl = 0; jl < NkFiber; jl++) val += projFstar[il * NkFiber + jl] * PetscRealPart(valsF[jF * NkFiber + jl]); 11333f27d899SToby Isaac workF[il] = val; 11343f27d899SToby Isaac } 11353f27d899SToby Isaac if (kFiber < 0) { 11363f27d899SToby Isaac for (il = 0; il < dF; il++) workF2[il] = workF[il]; 11379566063dSJacob Faibussowitsch PetscCall(PetscDTAltVStar(dim, PetscAbsInt(kFiber), 1, workF2, workF)); 11383f27d899SToby Isaac } 11399566063dSJacob Faibussowitsch PetscCall(PetscDTAltVWedgeMatrix(dim, PetscAbsInt(kFiber), PetscAbsInt(kTrace), workF, wedgeMat)); 11403f27d899SToby Isaac for (jT = 0; jT < nformsT; jT++, j++) { 11413f27d899SToby Isaac PetscInt colT = colsT[jT * NkTrace] / NkTrace; 11423f27d899SToby Isaac PetscInt col = colF * (nTrace / NkTrace) + colT; 11433f27d899SToby Isaac const PetscScalar *vals; 11443f27d899SToby Isaac 11453f27d899SToby Isaac for (il = 0; il < dT; il++) { 11463f27d899SToby Isaac PetscReal val = 0.; 11473f27d899SToby Isaac for (jl = 0; jl < NkTrace; jl++) val += projTstar[il * NkTrace + jl] * PetscRealPart(valsT[jT * NkTrace + jl]); 11483f27d899SToby Isaac workT[il] = val; 11493f27d899SToby Isaac } 11503f27d899SToby Isaac if (kTrace < 0) { 11513f27d899SToby Isaac for (il = 0; il < dT; il++) workT2[il] = workT[il]; 11529566063dSJacob Faibussowitsch PetscCall(PetscDTAltVStar(dim, PetscAbsInt(kTrace), 1, workT2, workT)); 11533f27d899SToby Isaac } 11543f27d899SToby Isaac 11553f27d899SToby Isaac for (il = 0; il < Nk; il++) { 11563f27d899SToby Isaac PetscReal val = 0.; 11573f27d899SToby Isaac for (jl = 0; jl < dT; jl++) val += sign * wedgeMat[il * dT + jl] * workT[jl]; 11583f27d899SToby Isaac work[il] = val; 11593f27d899SToby Isaac } 11603f27d899SToby Isaac if (k < 0) { 11619566063dSJacob Faibussowitsch PetscCall(PetscDTAltVStar(dim, PetscAbsInt(k), -1, work, workstar)); 11623f27d899SToby Isaac #if defined(PETSC_USE_COMPLEX) 11633f27d899SToby Isaac for (l = 0; l < Nk; l++) workS[l] = workstar[l]; 11643f27d899SToby Isaac vals = &workS[0]; 11653f27d899SToby Isaac #else 11663f27d899SToby Isaac vals = &workstar[0]; 11673f27d899SToby Isaac #endif 11683f27d899SToby Isaac } else { 11693f27d899SToby Isaac #if defined(PETSC_USE_COMPLEX) 11703f27d899SToby Isaac for (l = 0; l < Nk; l++) workS[l] = work[l]; 11713f27d899SToby Isaac vals = &workS[0]; 11723f27d899SToby Isaac #else 11733f27d899SToby Isaac vals = &work[0]; 11743f27d899SToby Isaac #endif 11753f27d899SToby Isaac } 11763f27d899SToby Isaac for (l = 0; l < Nk; l++) { 11779566063dSJacob Faibussowitsch PetscCall(MatSetValue(prod, i, col * Nk + l, vals[l], INSERT_VALUES)); 11783f27d899SToby Isaac } /* Nk */ 11793f27d899SToby Isaac } /* jT */ 11803f27d899SToby Isaac } /* jF */ 11819566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(trace, iT, &ncolsT, &colsT, &valsT)); 11823f27d899SToby Isaac } /* iT */ 11839566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(fiber, iF, &ncolsF, &colsF, &valsF)); 11843f27d899SToby Isaac } /* iF */ 11859566063dSJacob Faibussowitsch PetscCall(PetscFree(wedgeMat)); 11869566063dSJacob Faibussowitsch PetscCall(PetscFree4(projT, projF, projTstar, projFstar)); 11879566063dSJacob Faibussowitsch PetscCall(PetscFree2(workT2, workF2)); 11889566063dSJacob Faibussowitsch PetscCall(PetscFree5(workT, workF, work, workstar, workS)); 11899566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(prod, MAT_FINAL_ASSEMBLY)); 11909566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(prod, MAT_FINAL_ASSEMBLY)); 11913f27d899SToby Isaac *product = prod; 11923f27d899SToby Isaac PetscFunctionReturn(0); 11933f27d899SToby Isaac } 11943f27d899SToby Isaac 119577f1a120SToby Isaac /* Union of quadrature points, with an attempt to identify commont points in the two sets */ 11963f27d899SToby Isaac static PetscErrorCode PetscQuadraturePointsMerge(PetscQuadrature quadA, PetscQuadrature quadB, PetscQuadrature *quadJoint, PetscInt *aToJoint[], PetscInt *bToJoint[]) 11973f27d899SToby Isaac { 11983f27d899SToby Isaac PetscInt dimA, dimB; 11993f27d899SToby Isaac PetscInt nA, nB, nJoint, i, j, d; 12003f27d899SToby Isaac const PetscReal *pointsA; 12013f27d899SToby Isaac const PetscReal *pointsB; 12023f27d899SToby Isaac PetscReal *pointsJoint; 12033f27d899SToby Isaac PetscInt *aToJ, *bToJ; 12043f27d899SToby Isaac PetscQuadrature qJ; 12053f27d899SToby Isaac 12063f27d899SToby Isaac PetscFunctionBegin; 12079566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quadA, &dimA, NULL, &nA, &pointsA, NULL)); 12089566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quadB, &dimB, NULL, &nB, &pointsB, NULL)); 12092c71b3e2SJacob Faibussowitsch PetscCheckFalse(dimA != dimB,PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Quadrature points must be in the same dimension"); 12103f27d899SToby Isaac nJoint = nA; 12119566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nA, &aToJ)); 12123f27d899SToby Isaac for (i = 0; i < nA; i++) aToJ[i] = i; 12139566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nB, &bToJ)); 12143f27d899SToby Isaac for (i = 0; i < nB; i++) { 12153f27d899SToby Isaac for (j = 0; j < nA; j++) { 12163f27d899SToby Isaac bToJ[i] = -1; 12176ff15688SToby Isaac for (d = 0; d < dimA; d++) if (PetscAbsReal(pointsB[i * dimA + d] - pointsA[j * dimA + d]) > PETSC_SMALL) break; 12183f27d899SToby Isaac if (d == dimA) { 12193f27d899SToby Isaac bToJ[i] = j; 12203f27d899SToby Isaac break; 12213f27d899SToby Isaac } 12223f27d899SToby Isaac } 12233f27d899SToby Isaac if (bToJ[i] == -1) { 12243f27d899SToby Isaac bToJ[i] = nJoint++; 12253f27d899SToby Isaac } 12263f27d899SToby Isaac } 12273f27d899SToby Isaac *aToJoint = aToJ; 12283f27d899SToby Isaac *bToJoint = bToJ; 12299566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nJoint * dimA, &pointsJoint)); 12309566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(pointsJoint, pointsA, nA * dimA)); 12313f27d899SToby Isaac for (i = 0; i < nB; i++) { 12323f27d899SToby Isaac if (bToJ[i] >= nA) { 12333f27d899SToby Isaac for (d = 0; d < dimA; d++) pointsJoint[bToJ[i] * dimA + d] = pointsB[i * dimA + d]; 12343f27d899SToby Isaac } 12353f27d899SToby Isaac } 12369566063dSJacob Faibussowitsch PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, &qJ)); 12379566063dSJacob Faibussowitsch PetscCall(PetscQuadratureSetData(qJ, dimA, 0, nJoint, pointsJoint, NULL)); 12383f27d899SToby Isaac *quadJoint = qJ; 12393f27d899SToby Isaac PetscFunctionReturn(0); 12403f27d899SToby Isaac } 12413f27d899SToby Isaac 124277f1a120SToby Isaac /* Matrices matA and matB are both quadrature -> dof matrices: produce a matrix that is joint quadrature -> union of 124377f1a120SToby Isaac * dofs, where the joint quadrature was produced by PetscQuadraturePointsMerge */ 12443f27d899SToby Isaac static PetscErrorCode MatricesMerge(Mat matA, Mat matB, PetscInt dim, PetscInt k, PetscInt numMerged, const PetscInt aToMerged[], const PetscInt bToMerged[], Mat *matMerged) 12453f27d899SToby Isaac { 12463f27d899SToby Isaac PetscInt m, n, mA, nA, mB, nB, Nk, i, j, l; 12473f27d899SToby Isaac Mat M; 12483f27d899SToby Isaac PetscInt *nnz; 12493f27d899SToby Isaac PetscInt maxnnz; 12503f27d899SToby Isaac PetscInt *work; 12513f27d899SToby Isaac 12523f27d899SToby Isaac PetscFunctionBegin; 12539566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(k), &Nk)); 12549566063dSJacob Faibussowitsch PetscCall(MatGetSize(matA, &mA, &nA)); 12552c71b3e2SJacob Faibussowitsch PetscCheckFalse(nA % Nk,PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "matA column space not a multiple of k-form size"); 12569566063dSJacob Faibussowitsch PetscCall(MatGetSize(matB, &mB, &nB)); 12572c71b3e2SJacob Faibussowitsch PetscCheckFalse(nB % Nk,PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "matB column space not a multiple of k-form size"); 12583f27d899SToby Isaac m = mA + mB; 12593f27d899SToby Isaac n = numMerged * Nk; 12609566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m, &nnz)); 12613f27d899SToby Isaac maxnnz = 0; 12623f27d899SToby Isaac for (i = 0; i < mA; i++) { 12639566063dSJacob Faibussowitsch PetscCall(MatGetRow(matA, i, &(nnz[i]), NULL, NULL)); 12642c71b3e2SJacob Faibussowitsch PetscCheckFalse(nnz[i] % Nk,PETSC_COMM_SELF, PETSC_ERR_PLIB, "nonzeros in matA are not in k-form size blocks"); 12653f27d899SToby Isaac maxnnz = PetscMax(maxnnz, nnz[i]); 12663f27d899SToby Isaac } 12673f27d899SToby Isaac for (i = 0; i < mB; i++) { 12689566063dSJacob Faibussowitsch PetscCall(MatGetRow(matB, i, &(nnz[i+mA]), NULL, NULL)); 12692c71b3e2SJacob Faibussowitsch PetscCheckFalse(nnz[i+mA] % Nk,PETSC_COMM_SELF, PETSC_ERR_PLIB, "nonzeros in matB are not in k-form size blocks"); 12703f27d899SToby Isaac maxnnz = PetscMax(maxnnz, nnz[i+mA]); 12713f27d899SToby Isaac } 12729566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, m, n, 0, nnz, &M)); 12739566063dSJacob Faibussowitsch PetscCall(PetscFree(nnz)); 12743f27d899SToby Isaac /* reasoning about which points each dof needs depends on having zeros computed at points preserved */ 12759566063dSJacob Faibussowitsch PetscCall(MatSetOption(M, MAT_IGNORE_ZERO_ENTRIES, PETSC_FALSE)); 12769566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(maxnnz, &work)); 12773f27d899SToby Isaac for (i = 0; i < mA; i++) { 12783f27d899SToby Isaac const PetscInt *cols; 12793f27d899SToby Isaac const PetscScalar *vals; 12803f27d899SToby Isaac PetscInt nCols; 12819566063dSJacob Faibussowitsch PetscCall(MatGetRow(matA, i, &nCols, &cols, &vals)); 12823f27d899SToby Isaac for (j = 0; j < nCols / Nk; j++) { 12833f27d899SToby Isaac PetscInt newCol = aToMerged[cols[j * Nk] / Nk]; 12843f27d899SToby Isaac for (l = 0; l < Nk; l++) work[j * Nk + l] = newCol * Nk + l; 12853f27d899SToby Isaac } 12869566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlocked(M, 1, &i, nCols, work, vals, INSERT_VALUES)); 12879566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(matA, i, &nCols, &cols, &vals)); 12883f27d899SToby Isaac } 12893f27d899SToby Isaac for (i = 0; i < mB; i++) { 12903f27d899SToby Isaac const PetscInt *cols; 12913f27d899SToby Isaac const PetscScalar *vals; 12923f27d899SToby Isaac 12933f27d899SToby Isaac PetscInt row = i + mA; 12943f27d899SToby Isaac PetscInt nCols; 12959566063dSJacob Faibussowitsch PetscCall(MatGetRow(matB, i, &nCols, &cols, &vals)); 12963f27d899SToby Isaac for (j = 0; j < nCols / Nk; j++) { 12973f27d899SToby Isaac PetscInt newCol = bToMerged[cols[j * Nk] / Nk]; 12983f27d899SToby Isaac for (l = 0; l < Nk; l++) work[j * Nk + l] = newCol * Nk + l; 12993f27d899SToby Isaac } 13009566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlocked(M, 1, &row, nCols, work, vals, INSERT_VALUES)); 13019566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(matB, i, &nCols, &cols, &vals)); 13023f27d899SToby Isaac } 13039566063dSJacob Faibussowitsch PetscCall(PetscFree(work)); 13049566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(M, MAT_FINAL_ASSEMBLY)); 13059566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(M, MAT_FINAL_ASSEMBLY)); 13063f27d899SToby Isaac *matMerged = M; 13073f27d899SToby Isaac PetscFunctionReturn(0); 13083f27d899SToby Isaac } 13093f27d899SToby Isaac 131077f1a120SToby Isaac /* Take a dual space and product a segment space that has all the same specifications (trimmed, continuous, order, 131177f1a120SToby Isaac * node set), except for the form degree. For computing boundary dofs and for making tensor product spaces */ 13123f27d899SToby Isaac static PetscErrorCode PetscDualSpaceCreateFacetSubspace_Lagrange(PetscDualSpace sp, DM K, PetscInt f, PetscInt k, PetscInt Ncopies, PetscBool interiorOnly, PetscDualSpace *bdsp) 13133f27d899SToby Isaac { 13143f27d899SToby Isaac PetscInt Nknew, Ncnew; 13153f27d899SToby Isaac PetscInt dim, pointDim = -1; 13163f27d899SToby Isaac PetscInt depth; 13173f27d899SToby Isaac DM dm; 13183f27d899SToby Isaac PetscDualSpace_Lag *newlag; 13193f27d899SToby Isaac 13203f27d899SToby Isaac PetscFunctionBegin; 13219566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDM(sp,&dm)); 13229566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm,&dim)); 13239566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepth(dm,&depth)); 13249566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceDuplicate(sp,bdsp)); 13259566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetFormDegree(*bdsp,k)); 13263f27d899SToby Isaac if (!K) { 13273f27d899SToby Isaac if (depth == dim) { 1328f783ec47SMatthew G. Knepley DMPolytopeType ct; 13293f27d899SToby Isaac 13306ff15688SToby Isaac pointDim = dim - 1; 13319566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, f, &ct)); 13329566063dSJacob Faibussowitsch PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, ct, &K)); 13333f27d899SToby Isaac } else if (depth == 1) { 13343f27d899SToby Isaac pointDim = 0; 13359566063dSJacob Faibussowitsch PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, DM_POLYTOPE_POINT, &K)); 13363f27d899SToby Isaac } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unsupported interpolation state of reference element"); 13373f27d899SToby Isaac } else { 13389566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)K)); 13399566063dSJacob Faibussowitsch PetscCall(DMGetDimension(K, &pointDim)); 13403f27d899SToby Isaac } 13419566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetDM(*bdsp, K)); 13429566063dSJacob Faibussowitsch PetscCall(DMDestroy(&K)); 13439566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(pointDim, PetscAbsInt(k), &Nknew)); 13443f27d899SToby Isaac Ncnew = Nknew * Ncopies; 13459566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetNumComponents(*bdsp, Ncnew)); 13463f27d899SToby Isaac newlag = (PetscDualSpace_Lag *) (*bdsp)->data; 13473f27d899SToby Isaac newlag->interiorOnly = interiorOnly; 13489566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetUp(*bdsp)); 13493f27d899SToby Isaac PetscFunctionReturn(0); 13503f27d899SToby Isaac } 13513f27d899SToby Isaac 135277f1a120SToby Isaac /* Construct simplex nodes from a nodefamily, add Nk dof vectors of length Nk at each node. 135377f1a120SToby Isaac * Return the (quadrature, matrix) form of the dofs and the nodeIndices form as well. 135477f1a120SToby Isaac * 135577f1a120SToby Isaac * Sometimes we want a set of nodes to be contained in the interior of the element, 135677f1a120SToby Isaac * even when the node scheme puts nodes on the boundaries. numNodeSkip tells 135777f1a120SToby Isaac * the routine how many "layers" of nodes need to be skipped. 135877f1a120SToby Isaac * */ 13593f27d899SToby Isaac static PetscErrorCode PetscDualSpaceLagrangeCreateSimplexNodeMat(Petsc1DNodeFamily nodeFamily, PetscInt dim, PetscInt sum, PetscInt Nk, PetscInt numNodeSkip, PetscQuadrature *iNodes, Mat *iMat, PetscLagNodeIndices *nodeIndices) 13603f27d899SToby Isaac { 13613f27d899SToby Isaac PetscReal *extraNodeCoords, *nodeCoords; 13623f27d899SToby Isaac PetscInt nNodes, nExtraNodes; 13633f27d899SToby Isaac PetscInt i, j, k, extraSum = sum + numNodeSkip * (1 + dim); 13643f27d899SToby Isaac PetscQuadrature intNodes; 13653f27d899SToby Isaac Mat intMat; 13663f27d899SToby Isaac PetscLagNodeIndices ni; 13673f27d899SToby Isaac 13683f27d899SToby Isaac PetscFunctionBegin; 13699566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(dim + sum, dim, &nNodes)); 13709566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(dim + extraSum, dim, &nExtraNodes)); 13713f27d899SToby Isaac 13729566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dim * nExtraNodes, &extraNodeCoords)); 13739566063dSJacob Faibussowitsch PetscCall(PetscNew(&ni)); 13743f27d899SToby Isaac ni->nodeIdxDim = dim + 1; 13753f27d899SToby Isaac ni->nodeVecDim = Nk; 13763f27d899SToby Isaac ni->nNodes = nNodes * Nk; 13773f27d899SToby Isaac ni->refct = 1; 13789566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nNodes * Nk * (dim + 1), &(ni->nodeIdx))); 13799566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nNodes * Nk * Nk, &(ni->nodeVec))); 13803f27d899SToby Isaac for (i = 0; i < nNodes; i++) for (j = 0; j < Nk; j++) for (k = 0; k < Nk; k++) ni->nodeVec[(i * Nk + j) * Nk + k] = (j == k) ? 1. : 0.; 13819566063dSJacob Faibussowitsch PetscCall(Petsc1DNodeFamilyComputeSimplexNodes(nodeFamily, dim, extraSum, extraNodeCoords)); 13823f27d899SToby Isaac if (numNodeSkip) { 13833f27d899SToby Isaac PetscInt k; 13843f27d899SToby Isaac PetscInt *tup; 13853f27d899SToby Isaac 13869566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dim * nNodes, &nodeCoords)); 13879566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dim + 1, &tup)); 13883f27d899SToby Isaac for (k = 0; k < nNodes; k++) { 13893f27d899SToby Isaac PetscInt j, c; 13903f27d899SToby Isaac PetscInt index; 13913f27d899SToby Isaac 13929566063dSJacob Faibussowitsch PetscCall(PetscDTIndexToBary(dim + 1, sum, k, tup)); 13933f27d899SToby Isaac for (j = 0; j < dim + 1; j++) tup[j] += numNodeSkip; 13943f27d899SToby Isaac for (c = 0; c < Nk; c++) { 13953f27d899SToby Isaac for (j = 0; j < dim + 1; j++) { 13963f27d899SToby Isaac ni->nodeIdx[(k * Nk + c) * (dim + 1) + j] = tup[j] + 1; 13973f27d899SToby Isaac } 13983f27d899SToby Isaac } 13999566063dSJacob Faibussowitsch PetscCall(PetscDTBaryToIndex(dim + 1, extraSum, tup, &index)); 14003f27d899SToby Isaac for (j = 0; j < dim; j++) nodeCoords[k * dim + j] = extraNodeCoords[index * dim + j]; 14013f27d899SToby Isaac } 14029566063dSJacob Faibussowitsch PetscCall(PetscFree(tup)); 14039566063dSJacob Faibussowitsch PetscCall(PetscFree(extraNodeCoords)); 14043f27d899SToby Isaac } else { 14053f27d899SToby Isaac PetscInt k; 14063f27d899SToby Isaac PetscInt *tup; 14073f27d899SToby Isaac 14083f27d899SToby Isaac nodeCoords = extraNodeCoords; 14099566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dim + 1, &tup)); 14103f27d899SToby Isaac for (k = 0; k < nNodes; k++) { 14113f27d899SToby Isaac PetscInt j, c; 14123f27d899SToby Isaac 14139566063dSJacob Faibussowitsch PetscCall(PetscDTIndexToBary(dim + 1, sum, k, tup)); 14143f27d899SToby Isaac for (c = 0; c < Nk; c++) { 14153f27d899SToby Isaac for (j = 0; j < dim + 1; j++) { 14163f27d899SToby Isaac /* barycentric indices can have zeros, but we don't want to push forward zeros because it makes it harder to 141777f1a120SToby Isaac * determine which nodes correspond to which under symmetries, so we increase by 1. This is fine 141877f1a120SToby Isaac * because the nodeIdx coordinates don't have any meaning other than helping to identify symmetries */ 14193f27d899SToby Isaac ni->nodeIdx[(k * Nk + c) * (dim + 1) + j] = tup[j] + 1; 14203f27d899SToby Isaac } 14213f27d899SToby Isaac } 14223f27d899SToby Isaac } 14239566063dSJacob Faibussowitsch PetscCall(PetscFree(tup)); 14243f27d899SToby Isaac } 14259566063dSJacob Faibussowitsch PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, &intNodes)); 14269566063dSJacob Faibussowitsch PetscCall(PetscQuadratureSetData(intNodes, dim, 0, nNodes, nodeCoords, NULL)); 14279566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, nNodes * Nk, nNodes * Nk, Nk, NULL, &intMat)); 14289566063dSJacob Faibussowitsch PetscCall(MatSetOption(intMat,MAT_IGNORE_ZERO_ENTRIES,PETSC_FALSE)); 14293f27d899SToby Isaac for (j = 0; j < nNodes * Nk; j++) { 14303f27d899SToby Isaac PetscInt rem = j % Nk; 14313f27d899SToby Isaac PetscInt a, aprev = j - rem; 14323f27d899SToby Isaac PetscInt anext = aprev + Nk; 14333f27d899SToby Isaac 14343f27d899SToby Isaac for (a = aprev; a < anext; a++) { 14359566063dSJacob Faibussowitsch PetscCall(MatSetValue(intMat, j, a, (a == j) ? 1. : 0., INSERT_VALUES)); 14363f27d899SToby Isaac } 14373f27d899SToby Isaac } 14389566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(intMat, MAT_FINAL_ASSEMBLY)); 14399566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(intMat, MAT_FINAL_ASSEMBLY)); 14403f27d899SToby Isaac *iNodes = intNodes; 14413f27d899SToby Isaac *iMat = intMat; 14423f27d899SToby Isaac *nodeIndices = ni; 14433f27d899SToby Isaac PetscFunctionReturn(0); 14443f27d899SToby Isaac } 14453f27d899SToby Isaac 144677f1a120SToby Isaac /* once the nodeIndices have been created for the interior of the reference cell, and for all of the boundary cells, 1447a5b23f4aSJose E. Roman * push forward the boundary dofs and concatenate them into the full node indices for the dual space */ 14483f27d899SToby Isaac static PetscErrorCode PetscDualSpaceLagrangeCreateAllNodeIdx(PetscDualSpace sp) 14493f27d899SToby Isaac { 14503f27d899SToby Isaac DM dm; 14513f27d899SToby Isaac PetscInt dim, nDofs; 14523f27d899SToby Isaac PetscSection section; 14533f27d899SToby Isaac PetscInt pStart, pEnd, p; 14543f27d899SToby Isaac PetscInt formDegree, Nk; 14553f27d899SToby Isaac PetscInt nodeIdxDim, spintdim; 14563f27d899SToby Isaac PetscDualSpace_Lag *lag; 14573f27d899SToby Isaac PetscLagNodeIndices ni, verti; 14583f27d899SToby Isaac 14593f27d899SToby Isaac PetscFunctionBegin; 14603f27d899SToby Isaac lag = (PetscDualSpace_Lag *) sp->data; 14613f27d899SToby Isaac verti = lag->vertIndices; 14629566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDM(sp, &dm)); 14639566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 14649566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetFormDegree(sp, &formDegree)); 14659566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(formDegree), &Nk)); 14669566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetSection(sp, §ion)); 14679566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(section, &nDofs)); 14689566063dSJacob Faibussowitsch PetscCall(PetscNew(&ni)); 14693f27d899SToby Isaac ni->nodeIdxDim = nodeIdxDim = verti->nodeIdxDim; 14703f27d899SToby Isaac ni->nodeVecDim = Nk; 14713f27d899SToby Isaac ni->nNodes = nDofs; 14723f27d899SToby Isaac ni->refct = 1; 14739566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nodeIdxDim * nDofs, &(ni->nodeIdx))); 14749566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Nk * nDofs, &(ni->nodeVec))); 14759566063dSJacob Faibussowitsch PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 14769566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, 0, &spintdim)); 14773f27d899SToby Isaac if (spintdim) { 14789566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(ni->nodeIdx, lag->intNodeIndices->nodeIdx, spintdim * nodeIdxDim)); 14799566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(ni->nodeVec, lag->intNodeIndices->nodeVec, spintdim * Nk)); 14803f27d899SToby Isaac } 14813f27d899SToby Isaac for (p = pStart + 1; p < pEnd; p++) { 14823f27d899SToby Isaac PetscDualSpace psp = sp->pointSpaces[p]; 14833f27d899SToby Isaac PetscDualSpace_Lag *plag; 14843f27d899SToby Isaac PetscInt dof, off; 14853f27d899SToby Isaac 14869566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, p, &dof)); 14873f27d899SToby Isaac if (!dof) continue; 14883f27d899SToby Isaac plag = (PetscDualSpace_Lag *) psp->data; 14899566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, p, &off)); 14909566063dSJacob Faibussowitsch PetscCall(PetscLagNodeIndicesPushForward(dm, verti, p, plag->vertIndices, plag->intNodeIndices, 0, formDegree, &(ni->nodeIdx[off * nodeIdxDim]), &(ni->nodeVec[off * Nk]))); 14913f27d899SToby Isaac } 14923f27d899SToby Isaac lag->allNodeIndices = ni; 14933f27d899SToby Isaac PetscFunctionReturn(0); 14943f27d899SToby Isaac } 14953f27d899SToby Isaac 149677f1a120SToby Isaac /* once the (quadrature, Matrix) forms of the dofs have been created for the interior of the 149777f1a120SToby Isaac * reference cell and for the boundary cells, jk 149877f1a120SToby Isaac * push forward the boundary data and concatenate them into the full (quadrature, matrix) data 149977f1a120SToby Isaac * for the dual space */ 15003f27d899SToby Isaac static PetscErrorCode PetscDualSpaceCreateAllDataFromInteriorData(PetscDualSpace sp) 15013f27d899SToby Isaac { 15023f27d899SToby Isaac DM dm; 15033f27d899SToby Isaac PetscSection section; 15043f27d899SToby Isaac PetscInt pStart, pEnd, p, k, Nk, dim, Nc; 15053f27d899SToby Isaac PetscInt nNodes; 15063f27d899SToby Isaac PetscInt countNodes; 15073f27d899SToby Isaac Mat allMat; 15083f27d899SToby Isaac PetscQuadrature allNodes; 15093f27d899SToby Isaac PetscInt nDofs; 15103f27d899SToby Isaac PetscInt maxNzforms, j; 15113f27d899SToby Isaac PetscScalar *work; 15123f27d899SToby Isaac PetscReal *L, *J, *Jinv, *v0, *pv0; 15133f27d899SToby Isaac PetscInt *iwork; 15143f27d899SToby Isaac PetscReal *nodes; 15153f27d899SToby Isaac 15163f27d899SToby Isaac PetscFunctionBegin; 15179566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDM(sp, &dm)); 15189566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 15199566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetSection(sp, §ion)); 15209566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(section, &nDofs)); 15219566063dSJacob Faibussowitsch PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 15229566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetFormDegree(sp, &k)); 15239566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetNumComponents(sp, &Nc)); 15249566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(k), &Nk)); 15253f27d899SToby Isaac for (p = pStart, nNodes = 0, maxNzforms = 0; p < pEnd; p++) { 15263f27d899SToby Isaac PetscDualSpace psp; 15273f27d899SToby Isaac DM pdm; 15283f27d899SToby Isaac PetscInt pdim, pNk; 15293f27d899SToby Isaac PetscQuadrature intNodes; 15303f27d899SToby Isaac Mat intMat; 15313f27d899SToby Isaac 15329566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetPointSubspace(sp, p, &psp)); 15333f27d899SToby Isaac if (!psp) continue; 15349566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDM(psp, &pdm)); 15359566063dSJacob Faibussowitsch PetscCall(DMGetDimension(pdm, &pdim)); 15363f27d899SToby Isaac if (pdim < PetscAbsInt(k)) continue; 15379566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(pdim, PetscAbsInt(k), &pNk)); 15389566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetInteriorData(psp, &intNodes, &intMat)); 15393f27d899SToby Isaac if (intNodes) { 15403f27d899SToby Isaac PetscInt nNodesp; 15413f27d899SToby Isaac 15429566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(intNodes, NULL, NULL, &nNodesp, NULL, NULL)); 15433f27d899SToby Isaac nNodes += nNodesp; 15443f27d899SToby Isaac } 15453f27d899SToby Isaac if (intMat) { 15463f27d899SToby Isaac PetscInt maxNzsp; 15473f27d899SToby Isaac PetscInt maxNzformsp; 15483f27d899SToby Isaac 15499566063dSJacob Faibussowitsch PetscCall(MatSeqAIJGetMaxRowNonzeros(intMat, &maxNzsp)); 15502c71b3e2SJacob Faibussowitsch PetscCheckFalse(maxNzsp % pNk,PETSC_COMM_SELF, PETSC_ERR_PLIB, "interior matrix is not laid out as blocks of k-forms"); 15513f27d899SToby Isaac maxNzformsp = maxNzsp / pNk; 15523f27d899SToby Isaac maxNzforms = PetscMax(maxNzforms, maxNzformsp); 15533f27d899SToby Isaac } 15543f27d899SToby Isaac } 15559566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, nDofs, nNodes * Nc, maxNzforms * Nk, NULL, &allMat)); 15569566063dSJacob Faibussowitsch PetscCall(MatSetOption(allMat,MAT_IGNORE_ZERO_ENTRIES,PETSC_FALSE)); 15579566063dSJacob Faibussowitsch PetscCall(PetscMalloc7(dim, &v0, dim, &pv0, dim * dim, &J, dim * dim, &Jinv, Nk * Nk, &L, maxNzforms * Nk, &work, maxNzforms * Nk, &iwork)); 15583f27d899SToby Isaac for (j = 0; j < dim; j++) pv0[j] = -1.; 15599566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dim * nNodes, &nodes)); 15603f27d899SToby Isaac for (p = pStart, countNodes = 0; p < pEnd; p++) { 15613f27d899SToby Isaac PetscDualSpace psp; 15623f27d899SToby Isaac PetscQuadrature intNodes; 15633f27d899SToby Isaac DM pdm; 15643f27d899SToby Isaac PetscInt pdim, pNk; 15653f27d899SToby Isaac PetscInt countNodesIn = countNodes; 15663f27d899SToby Isaac PetscReal detJ; 15673f27d899SToby Isaac Mat intMat; 15683f27d899SToby Isaac 15699566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetPointSubspace(sp, p, &psp)); 15703f27d899SToby Isaac if (!psp) continue; 15719566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDM(psp, &pdm)); 15729566063dSJacob Faibussowitsch PetscCall(DMGetDimension(pdm, &pdim)); 15733f27d899SToby Isaac if (pdim < PetscAbsInt(k)) continue; 15749566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetInteriorData(psp, &intNodes, &intMat)); 15753f27d899SToby Isaac if (intNodes == NULL && intMat == NULL) continue; 15769566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(pdim, PetscAbsInt(k), &pNk)); 15773f27d899SToby Isaac if (p) { 15789566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryAffineFEM(dm, p, v0, J, Jinv, &detJ)); 15793f27d899SToby Isaac } else { /* identity */ 15803f27d899SToby Isaac PetscInt i,j; 15813f27d899SToby Isaac 15823f27d899SToby Isaac for (i = 0; i < dim; i++) for (j = 0; j < dim; j++) J[i * dim + j] = Jinv[i * dim + j] = 0.; 15833f27d899SToby Isaac for (i = 0; i < dim; i++) J[i * dim + i] = Jinv[i * dim + i] = 1.; 15843f27d899SToby Isaac for (i = 0; i < dim; i++) v0[i] = -1.; 15853f27d899SToby Isaac } 15863f27d899SToby Isaac if (pdim != dim) { /* compactify Jacobian */ 15873f27d899SToby Isaac PetscInt i, j; 15883f27d899SToby Isaac 15893f27d899SToby Isaac for (i = 0; i < dim; i++) for (j = 0; j < pdim; j++) J[i * pdim + j] = J[i * dim + j]; 15903f27d899SToby Isaac } 15919566063dSJacob Faibussowitsch PetscCall(PetscDTAltVPullbackMatrix(pdim, dim, J, k, L)); 159277f1a120SToby Isaac if (intNodes) { /* push forward quadrature locations by the affine transformation */ 15933f27d899SToby Isaac PetscInt nNodesp; 15943f27d899SToby Isaac const PetscReal *nodesp; 15953f27d899SToby Isaac PetscInt j; 15963f27d899SToby Isaac 15979566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(intNodes, NULL, NULL, &nNodesp, &nodesp, NULL)); 15983f27d899SToby Isaac for (j = 0; j < nNodesp; j++, countNodes++) { 15993f27d899SToby Isaac PetscInt d, e; 16003f27d899SToby Isaac 16013f27d899SToby Isaac for (d = 0; d < dim; d++) { 16023f27d899SToby Isaac nodes[countNodes * dim + d] = v0[d]; 16033f27d899SToby Isaac for (e = 0; e < pdim; e++) { 16043f27d899SToby Isaac nodes[countNodes * dim + d] += J[d * pdim + e] * (nodesp[j * pdim + e] - pv0[e]); 16053f27d899SToby Isaac } 16063f27d899SToby Isaac } 16073f27d899SToby Isaac } 16083f27d899SToby Isaac } 16093f27d899SToby Isaac if (intMat) { 16103f27d899SToby Isaac PetscInt nrows; 16113f27d899SToby Isaac PetscInt off; 16123f27d899SToby Isaac 16139566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, p, &nrows)); 16149566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, p, &off)); 16153f27d899SToby Isaac for (j = 0; j < nrows; j++) { 16163f27d899SToby Isaac PetscInt ncols; 16173f27d899SToby Isaac const PetscInt *cols; 16183f27d899SToby Isaac const PetscScalar *vals; 16193f27d899SToby Isaac PetscInt l, d, e; 16203f27d899SToby Isaac PetscInt row = j + off; 16213f27d899SToby Isaac 16229566063dSJacob Faibussowitsch PetscCall(MatGetRow(intMat, j, &ncols, &cols, &vals)); 16232c71b3e2SJacob Faibussowitsch PetscCheckFalse(ncols % pNk,PETSC_COMM_SELF, PETSC_ERR_PLIB, "interior matrix is not laid out as blocks of k-forms"); 16243f27d899SToby Isaac for (l = 0; l < ncols / pNk; l++) { 16253f27d899SToby Isaac PetscInt blockcol; 16263f27d899SToby Isaac 16273f27d899SToby Isaac for (d = 0; d < pNk; d++) { 16282c71b3e2SJacob Faibussowitsch PetscCheckFalse((cols[l * pNk + d] % pNk) != d,PETSC_COMM_SELF, PETSC_ERR_PLIB, "interior matrix is not laid out as blocks of k-forms"); 16293f27d899SToby Isaac } 16303f27d899SToby Isaac blockcol = cols[l * pNk] / pNk; 16313f27d899SToby Isaac for (d = 0; d < Nk; d++) { 16323f27d899SToby Isaac iwork[l * Nk + d] = (blockcol + countNodesIn) * Nk + d; 16333f27d899SToby Isaac } 16343f27d899SToby Isaac for (d = 0; d < Nk; d++) work[l * Nk + d] = 0.; 16353f27d899SToby Isaac for (d = 0; d < Nk; d++) { 16363f27d899SToby Isaac for (e = 0; e < pNk; e++) { 16373f27d899SToby Isaac /* "push forward" dof by pulling back a k-form to be evaluated on the point: multiply on the right by L */ 16385efe5503SToby Isaac work[l * Nk + d] += vals[l * pNk + e] * L[e * Nk + d]; 16393f27d899SToby Isaac } 16403f27d899SToby Isaac } 16413f27d899SToby Isaac } 16429566063dSJacob Faibussowitsch PetscCall(MatSetValues(allMat, 1, &row, (ncols / pNk) * Nk, iwork, work, INSERT_VALUES)); 16439566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(intMat, j, &ncols, &cols, &vals)); 16443f27d899SToby Isaac } 16453f27d899SToby Isaac } 16463f27d899SToby Isaac } 16479566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(allMat, MAT_FINAL_ASSEMBLY)); 16489566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(allMat, MAT_FINAL_ASSEMBLY)); 16499566063dSJacob Faibussowitsch PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, &allNodes)); 16509566063dSJacob Faibussowitsch PetscCall(PetscQuadratureSetData(allNodes, dim, 0, nNodes, nodes, NULL)); 16519566063dSJacob Faibussowitsch PetscCall(PetscFree7(v0, pv0, J, Jinv, L, work, iwork)); 16529566063dSJacob Faibussowitsch PetscCall(MatDestroy(&(sp->allMat))); 16533f27d899SToby Isaac sp->allMat = allMat; 16549566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&(sp->allNodes))); 16553f27d899SToby Isaac sp->allNodes = allNodes; 16563f27d899SToby Isaac PetscFunctionReturn(0); 16573f27d899SToby Isaac } 16583f27d899SToby Isaac 165977f1a120SToby Isaac /* rather than trying to get all data from the functionals, we create 166077f1a120SToby Isaac * the functionals from rows of the quadrature -> dof matrix. 166177f1a120SToby Isaac * 166277f1a120SToby Isaac * Ideally most of the uses of PetscDualSpace in PetscFE will switch 166377f1a120SToby Isaac * to using intMat and allMat, so that the individual functionals 166477f1a120SToby Isaac * don't need to be constructed at all */ 16653f27d899SToby Isaac static PetscErrorCode PetscDualSpaceComputeFunctionalsFromAllData(PetscDualSpace sp) 16663f27d899SToby Isaac { 16673f27d899SToby Isaac PetscQuadrature allNodes; 16683f27d899SToby Isaac Mat allMat; 16693f27d899SToby Isaac PetscInt nDofs; 16703f27d899SToby Isaac PetscInt dim, k, Nk, Nc, f; 16713f27d899SToby Isaac DM dm; 16723f27d899SToby Isaac PetscInt nNodes, spdim; 16733f27d899SToby Isaac const PetscReal *nodes = NULL; 16743f27d899SToby Isaac PetscSection section; 167566a6c23cSMatthew G. Knepley PetscBool useMoments; 16763f27d899SToby Isaac 16773f27d899SToby Isaac PetscFunctionBegin; 16789566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDM(sp, &dm)); 16799566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 16809566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetNumComponents(sp, &Nc)); 16819566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetFormDegree(sp, &k)); 16829566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(k), &Nk)); 16839566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetAllData(sp, &allNodes, &allMat)); 16843f27d899SToby Isaac nNodes = 0; 16853f27d899SToby Isaac if (allNodes) { 16869566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(allNodes, NULL, NULL, &nNodes, &nodes, NULL)); 16873f27d899SToby Isaac } 16889566063dSJacob Faibussowitsch PetscCall(MatGetSize(allMat, &nDofs, NULL)); 16899566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetSection(sp, §ion)); 16909566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(section, &spdim)); 16912c71b3e2SJacob Faibussowitsch PetscCheckFalse(spdim != nDofs,PETSC_COMM_SELF, PETSC_ERR_PLIB, "incompatible all matrix size"); 16929566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nDofs, &(sp->functional))); 16939566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceLagrangeGetUseMoments(sp, &useMoments)); 169466a6c23cSMatthew G. Knepley if (useMoments) { 169566a6c23cSMatthew G. Knepley Mat allMat; 169666a6c23cSMatthew G. Knepley PetscInt momentOrder, i; 169766a6c23cSMatthew G. Knepley PetscBool tensor; 169866a6c23cSMatthew G. Knepley const PetscReal *weights; 169966a6c23cSMatthew G. Knepley PetscScalar *array; 170066a6c23cSMatthew G. Knepley 17012c71b3e2SJacob Faibussowitsch PetscCheckFalse(nDofs != 1,PETSC_COMM_SELF, PETSC_ERR_SUP, "We do not yet support moments beyond P0, nDofs == %D", nDofs); 17029566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceLagrangeGetMomentOrder(sp, &momentOrder)); 17039566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceLagrangeGetTensor(sp, &tensor)); 17049566063dSJacob Faibussowitsch if (!tensor) PetscCall(PetscDTStroudConicalQuadrature(dim, Nc, PetscMax(momentOrder + 1,1), -1.0, 1.0, &(sp->functional[0]))); 17059566063dSJacob Faibussowitsch else PetscCall(PetscDTGaussTensorQuadrature(dim, Nc, PetscMax(momentOrder + 1,1), -1.0, 1.0, &(sp->functional[0]))); 170666a6c23cSMatthew G. Knepley /* Need to replace allNodes and allMat */ 17079566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject) sp->functional[0])); 17089566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&(sp->allNodes))); 170966a6c23cSMatthew G. Knepley sp->allNodes = sp->functional[0]; 17109566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(sp->allNodes, NULL, NULL, &nNodes, NULL, &weights)); 17119566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, nDofs, nNodes * Nc, NULL, &allMat)); 17129566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite(allMat, &array)); 171366a6c23cSMatthew G. Knepley for (i = 0; i < nNodes * Nc; ++i) array[i] = weights[i]; 17149566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayWrite(allMat, &array)); 17159566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(allMat, MAT_FINAL_ASSEMBLY)); 17169566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(allMat, MAT_FINAL_ASSEMBLY)); 17179566063dSJacob Faibussowitsch PetscCall(MatDestroy(&(sp->allMat))); 171866a6c23cSMatthew G. Knepley sp->allMat = allMat; 171966a6c23cSMatthew G. Knepley PetscFunctionReturn(0); 172066a6c23cSMatthew G. Knepley } 17213f27d899SToby Isaac for (f = 0; f < nDofs; f++) { 17223f27d899SToby Isaac PetscInt ncols, c; 17233f27d899SToby Isaac const PetscInt *cols; 17243f27d899SToby Isaac const PetscScalar *vals; 17253f27d899SToby Isaac PetscReal *nodesf; 17263f27d899SToby Isaac PetscReal *weightsf; 17273f27d899SToby Isaac PetscInt nNodesf; 17283f27d899SToby Isaac PetscInt countNodes; 17293f27d899SToby Isaac 17309566063dSJacob Faibussowitsch PetscCall(MatGetRow(allMat, f, &ncols, &cols, &vals)); 17312c71b3e2SJacob Faibussowitsch PetscCheckFalse(ncols % Nk,PETSC_COMM_SELF, PETSC_ERR_PLIB, "all matrix is not laid out as blocks of k-forms"); 17323f27d899SToby Isaac for (c = 1, nNodesf = 1; c < ncols; c++) { 17333f27d899SToby Isaac if ((cols[c] / Nc) != (cols[c-1] / Nc)) nNodesf++; 17343f27d899SToby Isaac } 17359566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dim * nNodesf, &nodesf)); 17369566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Nc * nNodesf, &weightsf)); 17373f27d899SToby Isaac for (c = 0, countNodes = 0; c < ncols; c++) { 17383f27d899SToby Isaac if (!c || ((cols[c] / Nc) != (cols[c-1] / Nc))) { 17393f27d899SToby Isaac PetscInt d; 17403f27d899SToby Isaac 17413f27d899SToby Isaac for (d = 0; d < Nc; d++) { 17423f27d899SToby Isaac weightsf[countNodes * Nc + d] = 0.; 17433f27d899SToby Isaac } 17443f27d899SToby Isaac for (d = 0; d < dim; d++) { 17453f27d899SToby Isaac nodesf[countNodes * dim + d] = nodes[(cols[c] / Nc) * dim + d]; 17463f27d899SToby Isaac } 17473f27d899SToby Isaac countNodes++; 17483f27d899SToby Isaac } 17493f27d899SToby Isaac weightsf[(countNodes - 1) * Nc + (cols[c] % Nc)] = PetscRealPart(vals[c]); 17503f27d899SToby Isaac } 17519566063dSJacob Faibussowitsch PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, &(sp->functional[f]))); 17529566063dSJacob Faibussowitsch PetscCall(PetscQuadratureSetData(sp->functional[f], dim, Nc, nNodesf, nodesf, weightsf)); 17539566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(allMat, f, &ncols, &cols, &vals)); 17543f27d899SToby Isaac } 17553f27d899SToby Isaac PetscFunctionReturn(0); 17563f27d899SToby Isaac } 17573f27d899SToby Isaac 17583f27d899SToby Isaac /* take a matrix meant for k-forms and expand it to one for Ncopies */ 17593f27d899SToby Isaac static PetscErrorCode PetscDualSpaceLagrangeMatrixCreateCopies(Mat A, PetscInt Nk, PetscInt Ncopies, Mat *Abs) 17603f27d899SToby Isaac { 17613f27d899SToby Isaac PetscInt m, n, i, j, k; 17623f27d899SToby Isaac PetscInt maxnnz, *nnz, *iwork; 17633f27d899SToby Isaac Mat Ac; 17643f27d899SToby Isaac 17653f27d899SToby Isaac PetscFunctionBegin; 17669566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &m, &n)); 17672c71b3e2SJacob Faibussowitsch PetscCheckFalse(n % Nk,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of columns in A %D is not a multiple of Nk %D", n, Nk); 17689566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m * Ncopies, &nnz)); 17693f27d899SToby Isaac for (i = 0, maxnnz = 0; i < m; i++) { 17703f27d899SToby Isaac PetscInt innz; 17719566063dSJacob Faibussowitsch PetscCall(MatGetRow(A, i, &innz, NULL, NULL)); 17722c71b3e2SJacob Faibussowitsch PetscCheckFalse(innz % Nk,PETSC_COMM_SELF, PETSC_ERR_PLIB, "A row %D nnzs is not a multiple of Nk %D", innz, Nk); 17733f27d899SToby Isaac for (j = 0; j < Ncopies; j++) nnz[i * Ncopies + j] = innz; 17743f27d899SToby Isaac maxnnz = PetscMax(maxnnz, innz); 17753f27d899SToby Isaac } 17769566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, m * Ncopies, n * Ncopies, 0, nnz, &Ac)); 17779566063dSJacob Faibussowitsch PetscCall(MatSetOption(Ac, MAT_IGNORE_ZERO_ENTRIES, PETSC_FALSE)); 17789566063dSJacob Faibussowitsch PetscCall(PetscFree(nnz)); 17799566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(maxnnz, &iwork)); 17803f27d899SToby Isaac for (i = 0; i < m; i++) { 17813f27d899SToby Isaac PetscInt innz; 17823f27d899SToby Isaac const PetscInt *cols; 17833f27d899SToby Isaac const PetscScalar *vals; 17843f27d899SToby Isaac 17859566063dSJacob Faibussowitsch PetscCall(MatGetRow(A, i, &innz, &cols, &vals)); 17863f27d899SToby Isaac for (j = 0; j < innz; j++) iwork[j] = (cols[j] / Nk) * (Nk * Ncopies) + (cols[j] % Nk); 17873f27d899SToby Isaac for (j = 0; j < Ncopies; j++) { 17883f27d899SToby Isaac PetscInt row = i * Ncopies + j; 17893f27d899SToby Isaac 17909566063dSJacob Faibussowitsch PetscCall(MatSetValues(Ac, 1, &row, innz, iwork, vals, INSERT_VALUES)); 17913f27d899SToby Isaac for (k = 0; k < innz; k++) iwork[k] += Nk; 17923f27d899SToby Isaac } 17939566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A, i, &innz, &cols, &vals)); 17943f27d899SToby Isaac } 17959566063dSJacob Faibussowitsch PetscCall(PetscFree(iwork)); 17969566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Ac, MAT_FINAL_ASSEMBLY)); 17979566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Ac, MAT_FINAL_ASSEMBLY)); 17983f27d899SToby Isaac *Abs = Ac; 17993f27d899SToby Isaac PetscFunctionReturn(0); 18003f27d899SToby Isaac } 18013f27d899SToby Isaac 180277f1a120SToby Isaac /* check if a cell is a tensor product of the segment with a facet, 180377f1a120SToby Isaac * specifically checking if f and f2 can be the "endpoints" (like the triangles 180477f1a120SToby Isaac * at either end of a wedge) */ 18053f27d899SToby Isaac static PetscErrorCode DMPlexPointIsTensor_Internal_Given(DM dm, PetscInt p, PetscInt f, PetscInt f2, PetscBool *isTensor) 18063f27d899SToby Isaac { 18073f27d899SToby Isaac PetscInt coneSize, c; 18083f27d899SToby Isaac const PetscInt *cone; 18093f27d899SToby Isaac const PetscInt *fCone; 18103f27d899SToby Isaac const PetscInt *f2Cone; 18113f27d899SToby Isaac PetscInt fs[2]; 18123f27d899SToby Isaac PetscInt meetSize, nmeet; 18133f27d899SToby Isaac const PetscInt *meet; 18143f27d899SToby Isaac 18153f27d899SToby Isaac PetscFunctionBegin; 18163f27d899SToby Isaac fs[0] = f; 18173f27d899SToby Isaac fs[1] = f2; 18189566063dSJacob Faibussowitsch PetscCall(DMPlexGetMeet(dm, 2, fs, &meetSize, &meet)); 18193f27d899SToby Isaac nmeet = meetSize; 18209566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreMeet(dm, 2, fs, &meetSize, &meet)); 182177f1a120SToby Isaac /* two points that have a non-empty meet cannot be at opposite ends of a cell */ 18223f27d899SToby Isaac if (nmeet) { 18233f27d899SToby Isaac *isTensor = PETSC_FALSE; 18243f27d899SToby Isaac PetscFunctionReturn(0); 18253f27d899SToby Isaac } 18269566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, p, &coneSize)); 18279566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, p, &cone)); 18289566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, f, &fCone)); 18299566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, f2, &f2Cone)); 18303f27d899SToby Isaac for (c = 0; c < coneSize; c++) { 18313f27d899SToby Isaac PetscInt e, ef; 18323f27d899SToby Isaac PetscInt d = -1, d2 = -1; 18333f27d899SToby Isaac PetscInt dcount, d2count; 18343f27d899SToby Isaac PetscInt t = cone[c]; 18353f27d899SToby Isaac PetscInt tConeSize; 18363f27d899SToby Isaac PetscBool tIsTensor; 18373f27d899SToby Isaac const PetscInt *tCone; 18383f27d899SToby Isaac 18393f27d899SToby Isaac if (t == f || t == f2) continue; 184077f1a120SToby Isaac /* for every other facet in the cone, check that is has 184177f1a120SToby Isaac * one ridge in common with each end */ 18429566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, t, &tConeSize)); 18439566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, t, &tCone)); 18443f27d899SToby Isaac 18453f27d899SToby Isaac dcount = 0; 18463f27d899SToby Isaac d2count = 0; 18473f27d899SToby Isaac for (e = 0; e < tConeSize; e++) { 18483f27d899SToby Isaac PetscInt q = tCone[e]; 18493f27d899SToby Isaac for (ef = 0; ef < coneSize - 2; ef++) { 18503f27d899SToby Isaac if (fCone[ef] == q) { 18513f27d899SToby Isaac if (dcount) { 18523f27d899SToby Isaac *isTensor = PETSC_FALSE; 18533f27d899SToby Isaac PetscFunctionReturn(0); 18543f27d899SToby Isaac } 18553f27d899SToby Isaac d = q; 18563f27d899SToby Isaac dcount++; 18573f27d899SToby Isaac } else if (f2Cone[ef] == q) { 18583f27d899SToby Isaac if (d2count) { 18593f27d899SToby Isaac *isTensor = PETSC_FALSE; 18603f27d899SToby Isaac PetscFunctionReturn(0); 18613f27d899SToby Isaac } 18623f27d899SToby Isaac d2 = q; 18633f27d899SToby Isaac d2count++; 18643f27d899SToby Isaac } 18653f27d899SToby Isaac } 18663f27d899SToby Isaac } 186777f1a120SToby Isaac /* if the whole cell is a tensor with the segment, then this 186877f1a120SToby Isaac * facet should be a tensor with the segment */ 18699566063dSJacob Faibussowitsch PetscCall(DMPlexPointIsTensor_Internal_Given(dm, t, d, d2, &tIsTensor)); 18703f27d899SToby Isaac if (!tIsTensor) { 18713f27d899SToby Isaac *isTensor = PETSC_FALSE; 18723f27d899SToby Isaac PetscFunctionReturn(0); 18733f27d899SToby Isaac } 18743f27d899SToby Isaac } 18753f27d899SToby Isaac *isTensor = PETSC_TRUE; 18763f27d899SToby Isaac PetscFunctionReturn(0); 18773f27d899SToby Isaac } 18783f27d899SToby Isaac 187977f1a120SToby Isaac /* determine if a cell is a tensor with a segment by looping over pairs of facets to find a pair 188077f1a120SToby Isaac * that could be the opposite ends */ 18813f27d899SToby Isaac static PetscErrorCode DMPlexPointIsTensor_Internal(DM dm, PetscInt p, PetscBool *isTensor, PetscInt *endA, PetscInt *endB) 18823f27d899SToby Isaac { 18833f27d899SToby Isaac PetscInt coneSize, c, c2; 18843f27d899SToby Isaac const PetscInt *cone; 18853f27d899SToby Isaac 18863f27d899SToby Isaac PetscFunctionBegin; 18879566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, p, &coneSize)); 18883f27d899SToby Isaac if (!coneSize) { 18893f27d899SToby Isaac if (isTensor) *isTensor = PETSC_FALSE; 18903f27d899SToby Isaac if (endA) *endA = -1; 18913f27d899SToby Isaac if (endB) *endB = -1; 18923f27d899SToby Isaac } 18939566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, p, &cone)); 18943f27d899SToby Isaac for (c = 0; c < coneSize; c++) { 18953f27d899SToby Isaac PetscInt f = cone[c]; 18963f27d899SToby Isaac PetscInt fConeSize; 18973f27d899SToby Isaac 18989566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, f, &fConeSize)); 18993f27d899SToby Isaac if (fConeSize != coneSize - 2) continue; 19003f27d899SToby Isaac 19013f27d899SToby Isaac for (c2 = c + 1; c2 < coneSize; c2++) { 19023f27d899SToby Isaac PetscInt f2 = cone[c2]; 19033f27d899SToby Isaac PetscBool isTensorff2; 19043f27d899SToby Isaac PetscInt f2ConeSize; 19053f27d899SToby Isaac 19069566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, f2, &f2ConeSize)); 19073f27d899SToby Isaac if (f2ConeSize != coneSize - 2) continue; 19083f27d899SToby Isaac 19099566063dSJacob Faibussowitsch PetscCall(DMPlexPointIsTensor_Internal_Given(dm, p, f, f2, &isTensorff2)); 19103f27d899SToby Isaac if (isTensorff2) { 19113f27d899SToby Isaac if (isTensor) *isTensor = PETSC_TRUE; 19123f27d899SToby Isaac if (endA) *endA = f; 19133f27d899SToby Isaac if (endB) *endB = f2; 19143f27d899SToby Isaac PetscFunctionReturn(0); 19153f27d899SToby Isaac } 19163f27d899SToby Isaac } 19173f27d899SToby Isaac } 19183f27d899SToby Isaac if (isTensor) *isTensor = PETSC_FALSE; 19193f27d899SToby Isaac if (endA) *endA = -1; 19203f27d899SToby Isaac if (endB) *endB = -1; 19213f27d899SToby Isaac PetscFunctionReturn(0); 19223f27d899SToby Isaac } 19233f27d899SToby Isaac 192477f1a120SToby Isaac /* determine if a cell is a tensor with a segment by looping over pairs of facets to find a pair 192577f1a120SToby Isaac * that could be the opposite ends */ 19263f27d899SToby Isaac static PetscErrorCode DMPlexPointIsTensor(DM dm, PetscInt p, PetscBool *isTensor, PetscInt *endA, PetscInt *endB) 19273f27d899SToby Isaac { 19283f27d899SToby Isaac DMPlexInterpolatedFlag interpolated; 19293f27d899SToby Isaac 19303f27d899SToby Isaac PetscFunctionBegin; 19319566063dSJacob Faibussowitsch PetscCall(DMPlexIsInterpolated(dm, &interpolated)); 19322c71b3e2SJacob Faibussowitsch PetscCheckFalse(interpolated != DMPLEX_INTERPOLATED_FULL,PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Only for interpolated DMPlex's"); 19339566063dSJacob Faibussowitsch PetscCall(DMPlexPointIsTensor_Internal(dm, p, isTensor, endA, endB)); 19343f27d899SToby Isaac PetscFunctionReturn(0); 19353f27d899SToby Isaac } 19363f27d899SToby Isaac 19378f28b7bfSToby Isaac /* Let k = formDegree and k' = -sign(k) * dim + k. Transform a symmetric frame for k-forms on the biunit simplex into 19388f28b7bfSToby Isaac * a symmetric frame for k'-forms on the biunit simplex. 19391f440fbeSToby Isaac * 19408f28b7bfSToby Isaac * A frame is "symmetric" if the pullback of every symmetry of the biunit simplex is a permutation of the frame. 19411f440fbeSToby Isaac * 19428f28b7bfSToby Isaac * forms in the symmetric frame are used as dofs in the untrimmed simplex spaces. This way, symmetries of the 19438f28b7bfSToby Isaac * reference cell result in permutations of dofs grouped by node. 19441f440fbeSToby Isaac * 19458f28b7bfSToby Isaac * Use T to transform dof matrices for k'-forms into dof matrices for k-forms as a block diagonal transformation on 19468f28b7bfSToby Isaac * the right. 19471f440fbeSToby Isaac */ 19481f440fbeSToby Isaac static PetscErrorCode BiunitSimplexSymmetricFormTransformation(PetscInt dim, PetscInt formDegree, PetscReal T[]) 19491f440fbeSToby Isaac { 19501f440fbeSToby Isaac PetscInt k = formDegree; 19511f440fbeSToby Isaac PetscInt kd = k < 0 ? dim + k : k - dim; 19521f440fbeSToby Isaac PetscInt Nk; 19531f440fbeSToby Isaac PetscReal *biToEq, *eqToBi, *biToEqStar, *eqToBiStar; 19541f440fbeSToby Isaac PetscInt fact; 19551f440fbeSToby Isaac 19561f440fbeSToby Isaac PetscFunctionBegin; 19579566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(k), &Nk)); 19589566063dSJacob Faibussowitsch PetscCall(PetscCalloc4(dim * dim, &biToEq, dim * dim, &eqToBi, Nk * Nk, &biToEqStar, Nk * Nk, &eqToBiStar)); 19591f440fbeSToby Isaac /* fill in biToEq: Jacobian of the transformation from the biunit simplex to the equilateral simplex */ 19601f440fbeSToby Isaac fact = 0; 19611f440fbeSToby Isaac for (PetscInt i = 0; i < dim; i++) { 19621f440fbeSToby Isaac biToEq[i * dim + i] = PetscSqrtReal(((PetscReal)i + 2.) / (2.*((PetscReal)i+1.))); 19631f440fbeSToby Isaac fact += 4*(i+1); 19641f440fbeSToby Isaac for (PetscInt j = i+1; j < dim; j++) { 19651f440fbeSToby Isaac biToEq[i * dim + j] = PetscSqrtReal(1./(PetscReal)fact); 19661f440fbeSToby Isaac } 19671f440fbeSToby Isaac } 19688f28b7bfSToby Isaac /* fill in eqToBi: Jacobian of the transformation from the equilateral simplex to the biunit simplex */ 19691f440fbeSToby Isaac fact = 0; 19701f440fbeSToby Isaac for (PetscInt j = 0; j < dim; j++) { 19711f440fbeSToby Isaac eqToBi[j * dim + j] = PetscSqrtReal(2.*((PetscReal)j+1.)/((PetscReal)j+2)); 19721f440fbeSToby Isaac fact += j+1; 19731f440fbeSToby Isaac for (PetscInt i = 0; i < j; i++) { 19741f440fbeSToby Isaac eqToBi[i * dim + j] = -PetscSqrtReal(1./(PetscReal)fact); 19751f440fbeSToby Isaac } 19761f440fbeSToby Isaac } 19779566063dSJacob Faibussowitsch PetscCall(PetscDTAltVPullbackMatrix(dim, dim, biToEq, kd, biToEqStar)); 19789566063dSJacob Faibussowitsch PetscCall(PetscDTAltVPullbackMatrix(dim, dim, eqToBi, k, eqToBiStar)); 19798f28b7bfSToby Isaac /* product of pullbacks simulates the following steps 19808f28b7bfSToby Isaac * 19818f28b7bfSToby Isaac * 1. start with frame W = [w_1, w_2, ..., w_m] of k forms that is symmetric on the biunit simplex: 19828f28b7bfSToby 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] 19838f28b7bfSToby Isaac is a permutation of W. 19848f28b7bfSToby Isaac Even though a k' form --- a (dim - k) form represented by its Hodge star --- has the same geometric 19858f28b7bfSToby Isaac content as a k form, W is not a symmetric frame of k' forms on the biunit simplex. That's because, 19868f28b7bfSToby Isaac for general Jacobian J, J_k* != J_k'*. 19878f28b7bfSToby Isaac * 2. pullback W to the equilateral triangle using the k pullback, W_eq = eqToBi_k* W. All symmetries of the 19888f28b7bfSToby Isaac equilateral simplex have orthonormal Jacobians. For an orthonormal Jacobian O, J_k* = J_k'*, so W_eq is 19898f28b7bfSToby Isaac also a symmetric frame for k' forms on the equilateral simplex. 19908f28b7bfSToby Isaac 3. pullback W_eq back to the biunit simplex using the k' pulback, V = biToEq_k'* W_eq = biToEq_k'* eqToBi_k* W. 19918f28b7bfSToby Isaac V is a symmetric frame for k' forms on the biunit simplex. 19928f28b7bfSToby Isaac */ 19931f440fbeSToby Isaac for (PetscInt i = 0; i < Nk; i++) { 19941f440fbeSToby Isaac for (PetscInt j = 0; j < Nk; j++) { 19951f440fbeSToby Isaac PetscReal val = 0.; 19961f440fbeSToby Isaac for (PetscInt k = 0; k < Nk; k++) val += biToEqStar[i * Nk + k] * eqToBiStar[k * Nk + j]; 19971f440fbeSToby Isaac T[i * Nk + j] = val; 19981f440fbeSToby Isaac } 19991f440fbeSToby Isaac } 20009566063dSJacob Faibussowitsch PetscCall(PetscFree4(biToEq, eqToBi, biToEqStar, eqToBiStar)); 20011f440fbeSToby Isaac PetscFunctionReturn(0); 20021f440fbeSToby Isaac } 20031f440fbeSToby Isaac 200477f1a120SToby Isaac /* permute a quadrature -> dof matrix so that its rows are in revlex order by nodeIdx */ 20053f27d899SToby Isaac static PetscErrorCode MatPermuteByNodeIdx(Mat A, PetscLagNodeIndices ni, Mat *Aperm) 20063f27d899SToby Isaac { 20073f27d899SToby Isaac PetscInt m, n, i, j; 20083f27d899SToby Isaac PetscInt nodeIdxDim = ni->nodeIdxDim; 20093f27d899SToby Isaac PetscInt nodeVecDim = ni->nodeVecDim; 20103f27d899SToby Isaac PetscInt *perm; 20113f27d899SToby Isaac IS permIS; 20123f27d899SToby Isaac IS id; 20133f27d899SToby Isaac PetscInt *nIdxPerm; 20143f27d899SToby Isaac PetscReal *nVecPerm; 20153f27d899SToby Isaac 20163f27d899SToby Isaac PetscFunctionBegin; 20179566063dSJacob Faibussowitsch PetscCall(PetscLagNodeIndicesGetPermutation(ni, &perm)); 20189566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &m, &n)); 20199566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nodeIdxDim * m, &nIdxPerm)); 20209566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nodeVecDim * m, &nVecPerm)); 20213f27d899SToby Isaac for (i = 0; i < m; i++) for (j = 0; j < nodeIdxDim; j++) nIdxPerm[i * nodeIdxDim + j] = ni->nodeIdx[perm[i] * nodeIdxDim + j]; 20223f27d899SToby Isaac for (i = 0; i < m; i++) for (j = 0; j < nodeVecDim; j++) nVecPerm[i * nodeVecDim + j] = ni->nodeVec[perm[i] * nodeVecDim + j]; 20239566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, m, perm, PETSC_USE_POINTER, &permIS)); 20249566063dSJacob Faibussowitsch PetscCall(ISSetPermutation(permIS)); 20259566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, n, 0, 1, &id)); 20269566063dSJacob Faibussowitsch PetscCall(ISSetPermutation(id)); 20279566063dSJacob Faibussowitsch PetscCall(MatPermute(A, permIS, id, Aperm)); 20289566063dSJacob Faibussowitsch PetscCall(ISDestroy(&permIS)); 20299566063dSJacob Faibussowitsch PetscCall(ISDestroy(&id)); 20303f27d899SToby Isaac for (i = 0; i < m; i++) perm[i] = i; 20319566063dSJacob Faibussowitsch PetscCall(PetscFree(ni->nodeIdx)); 20329566063dSJacob Faibussowitsch PetscCall(PetscFree(ni->nodeVec)); 20333f27d899SToby Isaac ni->nodeIdx = nIdxPerm; 20343f27d899SToby Isaac ni->nodeVec = nVecPerm; 20356f905325SMatthew G. Knepley PetscFunctionReturn(0); 20366f905325SMatthew G. Knepley } 20376f905325SMatthew G. Knepley 20386f905325SMatthew G. Knepley static PetscErrorCode PetscDualSpaceSetUp_Lagrange(PetscDualSpace sp) 20396f905325SMatthew G. Knepley { 20406f905325SMatthew G. Knepley PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 20416f905325SMatthew G. Knepley DM dm = sp->dm; 20423f27d899SToby Isaac DM dmint = NULL; 20433f27d899SToby Isaac PetscInt order; 20446f905325SMatthew G. Knepley PetscInt Nc = sp->Nc; 20456f905325SMatthew G. Knepley MPI_Comm comm; 20466f905325SMatthew G. Knepley PetscBool continuous; 20473f27d899SToby Isaac PetscSection section; 20483f27d899SToby Isaac PetscInt depth, dim, pStart, pEnd, cStart, cEnd, p, *pStratStart, *pStratEnd, d; 20493f27d899SToby Isaac PetscInt formDegree, Nk, Ncopies; 20503f27d899SToby Isaac PetscInt tensorf = -1, tensorf2 = -1; 20513f27d899SToby Isaac PetscBool tensorCell, tensorSpace; 20523f27d899SToby Isaac PetscBool uniform, trimmed; 20533f27d899SToby Isaac Petsc1DNodeFamily nodeFamily; 20543f27d899SToby Isaac PetscInt numNodeSkip; 20553f27d899SToby Isaac DMPlexInterpolatedFlag interpolated; 20563f27d899SToby Isaac PetscBool isbdm; 20576f905325SMatthew G. Knepley 20586f905325SMatthew G. Knepley PetscFunctionBegin; 20593f27d899SToby Isaac /* step 1: sanitize input */ 20609566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject) sp, &comm)); 20619566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 20629566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)sp, PETSCDUALSPACEBDM, &isbdm)); 20633f27d899SToby Isaac if (isbdm) { 20643f27d899SToby Isaac sp->k = -(dim-1); /* form degree of H-div */ 20659566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)sp, PETSCDUALSPACELAGRANGE)); 20663f27d899SToby Isaac } 20679566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetFormDegree(sp, &formDegree)); 20682c71b3e2SJacob Faibussowitsch PetscCheckFalse(PetscAbsInt(formDegree) > dim,comm, PETSC_ERR_ARG_OUTOFRANGE, "Form degree must be bounded by dimension"); 20699566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(dim,PetscAbsInt(formDegree),&Nk)); 20703f27d899SToby Isaac if (sp->Nc <= 0 && lag->numCopies > 0) sp->Nc = Nk * lag->numCopies; 20713f27d899SToby Isaac Nc = sp->Nc; 20722c71b3e2SJacob Faibussowitsch PetscCheckFalse(Nc % Nk,comm, PETSC_ERR_ARG_INCOMP, "Number of components is not a multiple of form degree size"); 20733f27d899SToby Isaac if (lag->numCopies <= 0) lag->numCopies = Nc / Nk; 20743f27d899SToby Isaac Ncopies = lag->numCopies; 20752c71b3e2SJacob Faibussowitsch PetscCheckFalse(Nc / Nk != Ncopies,comm, PETSC_ERR_ARG_INCOMP, "Number of copies * (dim choose k) != Nc"); 20763f27d899SToby Isaac if (!dim) sp->order = 0; 20773f27d899SToby Isaac order = sp->order; 20783f27d899SToby Isaac uniform = sp->uniform; 207928b400f6SJacob Faibussowitsch PetscCheck(uniform,PETSC_COMM_SELF, PETSC_ERR_SUP, "Variable order not supported yet"); 20803f27d899SToby Isaac if (lag->trimmed && !formDegree) lag->trimmed = PETSC_FALSE; /* trimmed spaces are the same as full spaces for 0-forms */ 20813f27d899SToby Isaac if (lag->nodeType == PETSCDTNODES_DEFAULT) { 20823f27d899SToby Isaac lag->nodeType = PETSCDTNODES_GAUSSJACOBI; 20833f27d899SToby Isaac lag->nodeExponent = 0.; 20843f27d899SToby Isaac /* trimmed spaces don't include corner vertices, so don't use end nodes by default */ 20853f27d899SToby Isaac lag->endNodes = lag->trimmed ? PETSC_FALSE : PETSC_TRUE; 20863f27d899SToby Isaac } 20873f27d899SToby Isaac /* If a trimmed space and the user did choose nodes with endpoints, skip them by default */ 20883f27d899SToby Isaac if (lag->numNodeSkip < 0) lag->numNodeSkip = (lag->trimmed && lag->endNodes) ? 1 : 0; 20893f27d899SToby Isaac numNodeSkip = lag->numNodeSkip; 20902c71b3e2SJacob Faibussowitsch PetscCheckFalse(lag->trimmed && !order,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot have zeroth order trimmed elements"); 20913f27d899SToby Isaac if (lag->trimmed && PetscAbsInt(formDegree) == dim) { /* convert trimmed n-forms to untrimmed of one polynomial order less */ 20923f27d899SToby Isaac sp->order--; 20933f27d899SToby Isaac order--; 20943f27d899SToby Isaac lag->trimmed = PETSC_FALSE; 20953f27d899SToby Isaac } 20963f27d899SToby Isaac trimmed = lag->trimmed; 20973f27d899SToby Isaac if (!order || PetscAbsInt(formDegree) == dim) lag->continuous = PETSC_FALSE; 20983f27d899SToby Isaac continuous = lag->continuous; 20999566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepth(dm, &depth)); 21009566063dSJacob Faibussowitsch PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 21019566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 21022c71b3e2SJacob Faibussowitsch PetscCheckFalse(pStart != 0 || cStart != 0,PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Expect DM with chart starting at zero and cells first"); 21032c71b3e2SJacob Faibussowitsch PetscCheckFalse(cEnd != 1,PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Use PETSCDUALSPACEREFINED for multi-cell reference meshes"); 21049566063dSJacob Faibussowitsch PetscCall(DMPlexIsInterpolated(dm, &interpolated)); 21053f27d899SToby Isaac if (interpolated != DMPLEX_INTERPOLATED_FULL) { 21069566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(dm, &dmint)); 21073f27d899SToby Isaac } else { 21089566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)dm)); 21093f27d899SToby Isaac dmint = dm; 21103f27d899SToby Isaac } 21113f27d899SToby Isaac tensorCell = PETSC_FALSE; 21123f27d899SToby Isaac if (dim > 1) { 21139566063dSJacob Faibussowitsch PetscCall(DMPlexPointIsTensor(dmint, 0, &tensorCell, &tensorf, &tensorf2)); 21143f27d899SToby Isaac } 21153f27d899SToby Isaac lag->tensorCell = tensorCell; 21163f27d899SToby Isaac if (dim < 2 || !lag->tensorCell) lag->tensorSpace = PETSC_FALSE; 21176f905325SMatthew G. Knepley tensorSpace = lag->tensorSpace; 21183f27d899SToby Isaac if (!lag->nodeFamily) { 21199566063dSJacob Faibussowitsch PetscCall(Petsc1DNodeFamilyCreate(lag->nodeType, lag->nodeExponent, lag->endNodes, &lag->nodeFamily)); 21203f27d899SToby Isaac } 21213f27d899SToby Isaac nodeFamily = lag->nodeFamily; 21222c71b3e2SJacob Faibussowitsch PetscCheckFalse(interpolated != DMPLEX_INTERPOLATED_FULL && continuous && (PetscAbsInt(formDegree) > 0 || order > 1),PETSC_COMM_SELF,PETSC_ERR_PLIB,"Reference element won't support all boundary nodes"); 212320cf1dd8SToby Isaac 21243f27d899SToby Isaac /* step 2: construct the boundary spaces */ 21259566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(depth+1,&pStratStart,depth+1,&pStratEnd)); 21269566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(pEnd,&(sp->pointSpaces))); 21279566063dSJacob Faibussowitsch for (d = 0; d <= depth; ++d) PetscCall(DMPlexGetDepthStratum(dm, d, &pStratStart[d], &pStratEnd[d])); 21289566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSectionCreate_Internal(sp, §ion)); 21293f27d899SToby Isaac sp->pointSection = section; 21303f27d899SToby Isaac if (continuous && !(lag->interiorOnly)) { 21313f27d899SToby Isaac PetscInt h; 21326f905325SMatthew G. Knepley 21333f27d899SToby Isaac for (p = pStratStart[depth - 1]; p < pStratEnd[depth - 1]; p++) { /* calculate the facet dual spaces */ 21343f27d899SToby Isaac PetscReal v0[3]; 21353f27d899SToby Isaac DMPolytopeType ptype; 21363f27d899SToby Isaac PetscReal J[9], detJ; 21376f905325SMatthew G. Knepley PetscInt q; 21386f905325SMatthew G. Knepley 21399566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryAffineFEM(dm, p, v0, J, NULL, &detJ)); 21409566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, p, &ptype)); 21416f905325SMatthew G. Knepley 214277f1a120SToby Isaac /* compare to previous facets: if computed, reference that dualspace */ 21433f27d899SToby Isaac for (q = pStratStart[depth - 1]; q < p; q++) { 21443f27d899SToby Isaac DMPolytopeType qtype; 21456f905325SMatthew G. Knepley 21469566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(dm, q, &qtype)); 21473f27d899SToby Isaac if (qtype == ptype) break; 21486f905325SMatthew G. Knepley } 21493f27d899SToby Isaac if (q < p) { /* this facet has the same dual space as that one */ 21509566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)sp->pointSpaces[q])); 21513f27d899SToby Isaac sp->pointSpaces[p] = sp->pointSpaces[q]; 21523f27d899SToby Isaac continue; 21536f905325SMatthew G. Knepley } 21543f27d899SToby Isaac /* if not, recursively compute this dual space */ 21559566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceCreateFacetSubspace_Lagrange(sp,NULL,p,formDegree,Ncopies,PETSC_FALSE,&sp->pointSpaces[p])); 21566f905325SMatthew G. Knepley } 21573f27d899SToby Isaac for (h = 2; h <= depth; h++) { /* get the higher subspaces from the facet subspaces */ 21583f27d899SToby Isaac PetscInt hd = depth - h; 21593f27d899SToby Isaac PetscInt hdim = dim - h; 21606f905325SMatthew G. Knepley 21613f27d899SToby Isaac if (hdim < PetscAbsInt(formDegree)) break; 21623f27d899SToby Isaac for (p = pStratStart[hd]; p < pStratEnd[hd]; p++) { 21633f27d899SToby Isaac PetscInt suppSize, s; 21643f27d899SToby Isaac const PetscInt *supp; 21656f905325SMatthew G. Knepley 21669566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupportSize(dm, p, &suppSize)); 21679566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupport(dm, p, &supp)); 21683f27d899SToby Isaac for (s = 0; s < suppSize; s++) { 21693f27d899SToby Isaac DM qdm; 21703f27d899SToby Isaac PetscDualSpace qsp, psp; 21713f27d899SToby Isaac PetscInt c, coneSize, q; 21723f27d899SToby Isaac const PetscInt *cone; 21733f27d899SToby Isaac const PetscInt *refCone; 21746f905325SMatthew G. Knepley 21753f27d899SToby Isaac q = supp[0]; 21763f27d899SToby Isaac qsp = sp->pointSpaces[q]; 21779566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, q, &coneSize)); 21789566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, q, &cone)); 21793f27d899SToby Isaac for (c = 0; c < coneSize; c++) if (cone[c] == p) break; 21802c71b3e2SJacob Faibussowitsch PetscCheckFalse(c == coneSize,PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "cone/support mismatch"); 21819566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDM(qsp, &qdm)); 21829566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(qdm, 0, &refCone)); 21833f27d899SToby Isaac /* get the equivalent dual space from the support dual space */ 21849566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetPointSubspace(qsp, refCone[c], &psp)); 21853f27d899SToby Isaac if (!s) { 21869566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)psp)); 21873f27d899SToby Isaac sp->pointSpaces[p] = psp; 21883f27d899SToby Isaac } 21893f27d899SToby Isaac } 21903f27d899SToby Isaac } 21913f27d899SToby Isaac } 21923f27d899SToby Isaac for (p = 1; p < pEnd; p++) { 21933f27d899SToby Isaac PetscInt pspdim; 21943f27d899SToby Isaac if (!sp->pointSpaces[p]) continue; 21959566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetInteriorDimension(sp->pointSpaces[p], &pspdim)); 21969566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(section, p, pspdim)); 21973f27d899SToby Isaac } 21983f27d899SToby Isaac } 21996f905325SMatthew G. Knepley 22003f27d899SToby Isaac if (Ncopies > 1) { 22013f27d899SToby Isaac Mat intMatScalar, allMatScalar; 22023f27d899SToby Isaac PetscDualSpace scalarsp; 22033f27d899SToby Isaac PetscDualSpace_Lag *scalarlag; 22043f27d899SToby Isaac 22059566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceDuplicate(sp, &scalarsp)); 220677f1a120SToby Isaac /* Setting the number of components to Nk is a space with 1 copy of each k-form */ 22079566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetNumComponents(scalarsp, Nk)); 22089566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetUp(scalarsp)); 22099566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetInteriorData(scalarsp, &(sp->intNodes), &intMatScalar)); 22109566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)(sp->intNodes))); 22119566063dSJacob Faibussowitsch if (intMatScalar) PetscCall(PetscDualSpaceLagrangeMatrixCreateCopies(intMatScalar, Nk, Ncopies, &(sp->intMat))); 22129566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetAllData(scalarsp, &(sp->allNodes), &allMatScalar)); 22139566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)(sp->allNodes))); 22149566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceLagrangeMatrixCreateCopies(allMatScalar, Nk, Ncopies, &(sp->allMat))); 22153f27d899SToby Isaac sp->spdim = scalarsp->spdim * Ncopies; 22163f27d899SToby Isaac sp->spintdim = scalarsp->spintdim * Ncopies; 22173f27d899SToby Isaac scalarlag = (PetscDualSpace_Lag *) scalarsp->data; 22189566063dSJacob Faibussowitsch PetscCall(PetscLagNodeIndicesReference(scalarlag->vertIndices)); 22193f27d899SToby Isaac lag->vertIndices = scalarlag->vertIndices; 22209566063dSJacob Faibussowitsch PetscCall(PetscLagNodeIndicesReference(scalarlag->intNodeIndices)); 22213f27d899SToby Isaac lag->intNodeIndices = scalarlag->intNodeIndices; 22229566063dSJacob Faibussowitsch PetscCall(PetscLagNodeIndicesReference(scalarlag->allNodeIndices)); 22233f27d899SToby Isaac lag->allNodeIndices = scalarlag->allNodeIndices; 22249566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceDestroy(&scalarsp)); 22259566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(section, 0, sp->spintdim)); 22269566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSectionSetUp_Internal(sp, section)); 22279566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceComputeFunctionalsFromAllData(sp)); 22289566063dSJacob Faibussowitsch PetscCall(PetscFree2(pStratStart, pStratEnd)); 22299566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmint)); 223020cf1dd8SToby Isaac PetscFunctionReturn(0); 223120cf1dd8SToby Isaac } 223220cf1dd8SToby Isaac 22333f27d899SToby Isaac if (trimmed && !continuous) { 22343f27d899SToby Isaac /* the dofs of a trimmed space don't have a nice tensor/lattice structure: 22353f27d899SToby Isaac * just construct the continuous dual space and copy all of the data over, 22363f27d899SToby Isaac * allocating it all to the cell instead of splitting it up between the boundaries */ 22373f27d899SToby Isaac PetscDualSpace spcont; 22383f27d899SToby Isaac PetscInt spdim, f; 22393f27d899SToby Isaac PetscQuadrature allNodes; 22403f27d899SToby Isaac PetscDualSpace_Lag *lagc; 22413f27d899SToby Isaac Mat allMat; 22423f27d899SToby Isaac 22439566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceDuplicate(sp, &spcont)); 22449566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceLagrangeSetContinuity(spcont, PETSC_TRUE)); 22459566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetUp(spcont)); 22469566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDimension(spcont, &spdim)); 22473f27d899SToby Isaac sp->spdim = sp->spintdim = spdim; 22489566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(section, 0, spdim)); 22499566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSectionSetUp_Internal(sp, section)); 22509566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(spdim, &(sp->functional))); 22513f27d899SToby Isaac for (f = 0; f < spdim; f++) { 22523f27d899SToby Isaac PetscQuadrature fn; 22533f27d899SToby Isaac 22549566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetFunctional(spcont, f, &fn)); 22559566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)fn)); 22563f27d899SToby Isaac sp->functional[f] = fn; 22573f27d899SToby Isaac } 22589566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetAllData(spcont, &allNodes, &allMat)); 22599566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject) allNodes)); 22609566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject) allNodes)); 22613f27d899SToby Isaac sp->allNodes = sp->intNodes = allNodes; 22629566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject) allMat)); 22639566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject) allMat)); 22643f27d899SToby Isaac sp->allMat = sp->intMat = allMat; 22653f27d899SToby Isaac lagc = (PetscDualSpace_Lag *) spcont->data; 22669566063dSJacob Faibussowitsch PetscCall(PetscLagNodeIndicesReference(lagc->vertIndices)); 22673f27d899SToby Isaac lag->vertIndices = lagc->vertIndices; 22689566063dSJacob Faibussowitsch PetscCall(PetscLagNodeIndicesReference(lagc->allNodeIndices)); 22699566063dSJacob Faibussowitsch PetscCall(PetscLagNodeIndicesReference(lagc->allNodeIndices)); 22703f27d899SToby Isaac lag->intNodeIndices = lagc->allNodeIndices; 22713f27d899SToby Isaac lag->allNodeIndices = lagc->allNodeIndices; 22729566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceDestroy(&spcont)); 22739566063dSJacob Faibussowitsch PetscCall(PetscFree2(pStratStart, pStratEnd)); 22749566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmint)); 22753f27d899SToby Isaac PetscFunctionReturn(0); 22763f27d899SToby Isaac } 22773f27d899SToby Isaac 22783f27d899SToby Isaac /* step 3: construct intNodes, and intMat, and combine it with boundray data to make allNodes and allMat */ 22793f27d899SToby Isaac if (!tensorSpace) { 22809566063dSJacob Faibussowitsch if (!tensorCell) PetscCall(PetscLagNodeIndicesCreateSimplexVertices(dm, &(lag->vertIndices))); 22813f27d899SToby Isaac 22823f27d899SToby Isaac if (trimmed) { 228377f1a120SToby Isaac /* there is one dof in the interior of the a trimmed element for each full polynomial of with degree at most 228477f1a120SToby Isaac * order + k - dim - 1 */ 22853f27d899SToby Isaac if (order + PetscAbsInt(formDegree) > dim) { 22863f27d899SToby Isaac PetscInt sum = order + PetscAbsInt(formDegree) - dim - 1; 22873f27d899SToby Isaac PetscInt nDofs; 22883f27d899SToby Isaac 22899566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceLagrangeCreateSimplexNodeMat(nodeFamily, dim, sum, Nk, numNodeSkip, &sp->intNodes, &sp->intMat, &(lag->intNodeIndices))); 22909566063dSJacob Faibussowitsch PetscCall(MatGetSize(sp->intMat, &nDofs, NULL)); 22919566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(section, 0, nDofs)); 22923f27d899SToby Isaac } 22939566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSectionSetUp_Internal(sp, section)); 22949566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceCreateAllDataFromInteriorData(sp)); 22959566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceLagrangeCreateAllNodeIdx(sp)); 22963f27d899SToby Isaac } else { 22973f27d899SToby Isaac if (!continuous) { 229877f1a120SToby Isaac /* if discontinuous just construct one node for each set of dofs (a set of dofs is a basis for the k-form 229977f1a120SToby Isaac * space) */ 23003f27d899SToby Isaac PetscInt sum = order; 23013f27d899SToby Isaac PetscInt nDofs; 23023f27d899SToby Isaac 23039566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceLagrangeCreateSimplexNodeMat(nodeFamily, dim, sum, Nk, numNodeSkip, &sp->intNodes, &sp->intMat, &(lag->intNodeIndices))); 23049566063dSJacob Faibussowitsch PetscCall(MatGetSize(sp->intMat, &nDofs, NULL)); 23059566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(section, 0, nDofs)); 23069566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSectionSetUp_Internal(sp, section)); 23079566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)(sp->intNodes))); 23083f27d899SToby Isaac sp->allNodes = sp->intNodes; 23099566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)(sp->intMat))); 23103f27d899SToby Isaac sp->allMat = sp->intMat; 23119566063dSJacob Faibussowitsch PetscCall(PetscLagNodeIndicesReference(lag->intNodeIndices)); 23123f27d899SToby Isaac lag->allNodeIndices = lag->intNodeIndices; 23133f27d899SToby Isaac } else { 231477f1a120SToby Isaac /* there is one dof in the interior of the a full element for each trimmed polynomial of with degree at most 231577f1a120SToby Isaac * order + k - dim, but with complementary form degree */ 23163f27d899SToby Isaac if (order + PetscAbsInt(formDegree) > dim) { 23173f27d899SToby Isaac PetscDualSpace trimmedsp; 23183f27d899SToby Isaac PetscDualSpace_Lag *trimmedlag; 23193f27d899SToby Isaac PetscQuadrature intNodes; 23203f27d899SToby Isaac PetscInt trFormDegree = formDegree >= 0 ? formDegree - dim : dim - PetscAbsInt(formDegree); 23213f27d899SToby Isaac PetscInt nDofs; 23223f27d899SToby Isaac Mat intMat; 23233f27d899SToby Isaac 23249566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceDuplicate(sp, &trimmedsp)); 23259566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceLagrangeSetTrimmed(trimmedsp, PETSC_TRUE)); 23269566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetOrder(trimmedsp, order + PetscAbsInt(formDegree) - dim)); 23279566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetFormDegree(trimmedsp, trFormDegree)); 23283f27d899SToby Isaac trimmedlag = (PetscDualSpace_Lag *) trimmedsp->data; 23293f27d899SToby Isaac trimmedlag->numNodeSkip = numNodeSkip + 1; 23309566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetUp(trimmedsp)); 23319566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetAllData(trimmedsp, &intNodes, &intMat)); 23329566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)intNodes)); 23333f27d899SToby Isaac sp->intNodes = intNodes; 23349566063dSJacob Faibussowitsch PetscCall(PetscLagNodeIndicesReference(trimmedlag->allNodeIndices)); 23353f27d899SToby Isaac lag->intNodeIndices = trimmedlag->allNodeIndices; 23369566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)intMat)); 23371f440fbeSToby Isaac if (PetscAbsInt(formDegree) > 0 && PetscAbsInt(formDegree) < dim) { 23381f440fbeSToby Isaac PetscReal *T; 23391f440fbeSToby Isaac PetscScalar *work; 23401f440fbeSToby Isaac PetscInt nCols, nRows; 23411f440fbeSToby Isaac Mat intMatT; 23421f440fbeSToby Isaac 23439566063dSJacob Faibussowitsch PetscCall(MatDuplicate(intMat, MAT_COPY_VALUES, &intMatT)); 23449566063dSJacob Faibussowitsch PetscCall(MatGetSize(intMat, &nRows, &nCols)); 23459566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(Nk * Nk, &T, nCols, &work)); 23469566063dSJacob Faibussowitsch PetscCall(BiunitSimplexSymmetricFormTransformation(dim, formDegree, T)); 23471f440fbeSToby Isaac for (PetscInt row = 0; row < nRows; row++) { 23481f440fbeSToby Isaac PetscInt nrCols; 23491f440fbeSToby Isaac const PetscInt *rCols; 23501f440fbeSToby Isaac const PetscScalar *rVals; 23511f440fbeSToby Isaac 23529566063dSJacob Faibussowitsch PetscCall(MatGetRow(intMat, row, &nrCols, &rCols, &rVals)); 23532c71b3e2SJacob Faibussowitsch PetscCheckFalse(nrCols % Nk,PETSC_COMM_SELF, PETSC_ERR_PLIB, "nonzeros in intMat matrix are not in k-form size blocks"); 23541f440fbeSToby Isaac for (PetscInt b = 0; b < nrCols; b += Nk) { 23551f440fbeSToby Isaac const PetscScalar *v = &rVals[b]; 23561f440fbeSToby Isaac PetscScalar *w = &work[b]; 23571f440fbeSToby Isaac for (PetscInt j = 0; j < Nk; j++) { 23581f440fbeSToby Isaac w[j] = 0.; 23591f440fbeSToby Isaac for (PetscInt i = 0; i < Nk; i++) { 23601f440fbeSToby Isaac w[j] += v[i] * T[i * Nk + j]; 23611f440fbeSToby Isaac } 23621f440fbeSToby Isaac } 23631f440fbeSToby Isaac } 23649566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlocked(intMatT, 1, &row, nrCols, rCols, work, INSERT_VALUES)); 23659566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(intMat, row, &nrCols, &rCols, &rVals)); 23661f440fbeSToby Isaac } 23679566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(intMatT, MAT_FINAL_ASSEMBLY)); 23689566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(intMatT, MAT_FINAL_ASSEMBLY)); 23699566063dSJacob Faibussowitsch PetscCall(MatDestroy(&intMat)); 23701f440fbeSToby Isaac intMat = intMatT; 23719566063dSJacob Faibussowitsch PetscCall(PetscLagNodeIndicesDestroy(&(lag->intNodeIndices))); 23729566063dSJacob Faibussowitsch PetscCall(PetscLagNodeIndicesDuplicate(trimmedlag->allNodeIndices, &(lag->intNodeIndices))); 23731f440fbeSToby Isaac { 23741f440fbeSToby Isaac PetscInt nNodes = lag->intNodeIndices->nNodes; 23751f440fbeSToby Isaac PetscReal *newNodeVec = lag->intNodeIndices->nodeVec; 23761f440fbeSToby Isaac const PetscReal *oldNodeVec = trimmedlag->allNodeIndices->nodeVec; 23771f440fbeSToby Isaac 23781f440fbeSToby Isaac for (PetscInt n = 0; n < nNodes; n++) { 23791f440fbeSToby Isaac PetscReal *w = &newNodeVec[n * Nk]; 23801f440fbeSToby Isaac const PetscReal *v = &oldNodeVec[n * Nk]; 23811f440fbeSToby Isaac 23821f440fbeSToby Isaac for (PetscInt j = 0; j < Nk; j++) { 23831f440fbeSToby Isaac w[j] = 0.; 23841f440fbeSToby Isaac for (PetscInt i = 0; i < Nk; i++) { 23851f440fbeSToby Isaac w[j] += v[i] * T[i * Nk + j]; 23861f440fbeSToby Isaac } 23871f440fbeSToby Isaac } 23881f440fbeSToby Isaac } 23891f440fbeSToby Isaac } 23909566063dSJacob Faibussowitsch PetscCall(PetscFree2(T, work)); 23911f440fbeSToby Isaac } 23921f440fbeSToby Isaac sp->intMat = intMat; 23939566063dSJacob Faibussowitsch PetscCall(MatGetSize(sp->intMat, &nDofs, NULL)); 23949566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceDestroy(&trimmedsp)); 23959566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(section, 0, nDofs)); 23963f27d899SToby Isaac } 23979566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSectionSetUp_Internal(sp, section)); 23989566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceCreateAllDataFromInteriorData(sp)); 23999566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceLagrangeCreateAllNodeIdx(sp)); 24003f27d899SToby Isaac } 24013f27d899SToby Isaac } 24023f27d899SToby Isaac } else { 24033f27d899SToby Isaac PetscQuadrature intNodesTrace = NULL; 24043f27d899SToby Isaac PetscQuadrature intNodesFiber = NULL; 24053f27d899SToby Isaac PetscQuadrature intNodes = NULL; 24063f27d899SToby Isaac PetscLagNodeIndices intNodeIndices = NULL; 24073f27d899SToby Isaac Mat intMat = NULL; 24083f27d899SToby Isaac 240977f1a120SToby Isaac if (PetscAbsInt(formDegree) < dim) { /* get the trace k-forms on the first facet, and the 0-forms on the edge, 241077f1a120SToby Isaac and wedge them together to create some of the k-form dofs */ 24113f27d899SToby Isaac PetscDualSpace trace, fiber; 24123f27d899SToby Isaac PetscDualSpace_Lag *tracel, *fiberl; 24133f27d899SToby Isaac Mat intMatTrace, intMatFiber; 24143f27d899SToby Isaac 24153f27d899SToby Isaac if (sp->pointSpaces[tensorf]) { 24169566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)(sp->pointSpaces[tensorf]))); 24173f27d899SToby Isaac trace = sp->pointSpaces[tensorf]; 24183f27d899SToby Isaac } else { 24199566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceCreateFacetSubspace_Lagrange(sp,NULL,tensorf,formDegree,Ncopies,PETSC_TRUE,&trace)); 24203f27d899SToby Isaac } 24219566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceCreateEdgeSubspace_Lagrange(sp,order,0,1,PETSC_TRUE,&fiber)); 24223f27d899SToby Isaac tracel = (PetscDualSpace_Lag *) trace->data; 24233f27d899SToby Isaac fiberl = (PetscDualSpace_Lag *) fiber->data; 24249566063dSJacob Faibussowitsch PetscCall(PetscLagNodeIndicesCreateTensorVertices(dm, tracel->vertIndices, &(lag->vertIndices))); 24259566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetInteriorData(trace, &intNodesTrace, &intMatTrace)); 24269566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetInteriorData(fiber, &intNodesFiber, &intMatFiber)); 24273f27d899SToby Isaac if (intNodesTrace && intNodesFiber) { 24289566063dSJacob Faibussowitsch PetscCall(PetscQuadratureCreateTensor(intNodesTrace, intNodesFiber, &intNodes)); 24299566063dSJacob Faibussowitsch PetscCall(MatTensorAltV(intMatTrace, intMatFiber, dim-1, formDegree, 1, 0, &intMat)); 24309566063dSJacob Faibussowitsch PetscCall(PetscLagNodeIndicesTensor(tracel->intNodeIndices, dim - 1, formDegree, fiberl->intNodeIndices, 1, 0, &intNodeIndices)); 24313f27d899SToby Isaac } 24329566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject) intNodesTrace)); 24339566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject) intNodesFiber)); 24349566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceDestroy(&fiber)); 24359566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceDestroy(&trace)); 24363f27d899SToby Isaac } 243777f1a120SToby Isaac if (PetscAbsInt(formDegree) > 0) { /* get the trace (k-1)-forms on the first facet, and the 1-forms on the edge, 243877f1a120SToby Isaac and wedge them together to create the remaining k-form dofs */ 24393f27d899SToby Isaac PetscDualSpace trace, fiber; 24403f27d899SToby Isaac PetscDualSpace_Lag *tracel, *fiberl; 24413f27d899SToby Isaac PetscQuadrature intNodesTrace2, intNodesFiber2, intNodes2; 24423f27d899SToby Isaac PetscLagNodeIndices intNodeIndices2; 24433f27d899SToby Isaac Mat intMatTrace, intMatFiber, intMat2; 24443f27d899SToby Isaac PetscInt traceDegree = formDegree > 0 ? formDegree - 1 : formDegree + 1; 24453f27d899SToby Isaac PetscInt fiberDegree = formDegree > 0 ? 1 : -1; 24463f27d899SToby Isaac 24479566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceCreateFacetSubspace_Lagrange(sp,NULL,tensorf,traceDegree,Ncopies,PETSC_TRUE,&trace)); 24489566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceCreateEdgeSubspace_Lagrange(sp,order,fiberDegree,1,PETSC_TRUE,&fiber)); 24493f27d899SToby Isaac tracel = (PetscDualSpace_Lag *) trace->data; 24503f27d899SToby Isaac fiberl = (PetscDualSpace_Lag *) fiber->data; 24513f27d899SToby Isaac if (!lag->vertIndices) { 24529566063dSJacob Faibussowitsch PetscCall(PetscLagNodeIndicesCreateTensorVertices(dm, tracel->vertIndices, &(lag->vertIndices))); 24533f27d899SToby Isaac } 24549566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetInteriorData(trace, &intNodesTrace2, &intMatTrace)); 24559566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetInteriorData(fiber, &intNodesFiber2, &intMatFiber)); 24563f27d899SToby Isaac if (intNodesTrace2 && intNodesFiber2) { 24579566063dSJacob Faibussowitsch PetscCall(PetscQuadratureCreateTensor(intNodesTrace2, intNodesFiber2, &intNodes2)); 24589566063dSJacob Faibussowitsch PetscCall(MatTensorAltV(intMatTrace, intMatFiber, dim-1, traceDegree, 1, fiberDegree, &intMat2)); 24599566063dSJacob Faibussowitsch PetscCall(PetscLagNodeIndicesTensor(tracel->intNodeIndices, dim - 1, traceDegree, fiberl->intNodeIndices, 1, fiberDegree, &intNodeIndices2)); 24603f27d899SToby Isaac if (!intMat) { 24613f27d899SToby Isaac intMat = intMat2; 24623f27d899SToby Isaac intNodes = intNodes2; 24633f27d899SToby Isaac intNodeIndices = intNodeIndices2; 24643f27d899SToby Isaac } else { 246577f1a120SToby Isaac /* merge the matrices, quadrature points, and nodes */ 24663f27d899SToby Isaac PetscInt nM; 24673f27d899SToby Isaac PetscInt nDof, nDof2; 24686ff15688SToby Isaac PetscInt *toMerged = NULL, *toMerged2 = NULL; 24696ff15688SToby Isaac PetscQuadrature merged = NULL; 24703f27d899SToby Isaac PetscLagNodeIndices intNodeIndicesMerged = NULL; 24713f27d899SToby Isaac Mat matMerged = NULL; 24723f27d899SToby Isaac 24739566063dSJacob Faibussowitsch PetscCall(MatGetSize(intMat, &nDof, NULL)); 24749566063dSJacob Faibussowitsch PetscCall(MatGetSize(intMat2, &nDof2, NULL)); 24759566063dSJacob Faibussowitsch PetscCall(PetscQuadraturePointsMerge(intNodes, intNodes2, &merged, &toMerged, &toMerged2)); 24769566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(merged, NULL, NULL, &nM, NULL, NULL)); 24779566063dSJacob Faibussowitsch PetscCall(MatricesMerge(intMat, intMat2, dim, formDegree, nM, toMerged, toMerged2, &matMerged)); 24789566063dSJacob Faibussowitsch PetscCall(PetscLagNodeIndicesMerge(intNodeIndices, intNodeIndices2, &intNodeIndicesMerged)); 24799566063dSJacob Faibussowitsch PetscCall(PetscFree(toMerged)); 24809566063dSJacob Faibussowitsch PetscCall(PetscFree(toMerged2)); 24819566063dSJacob Faibussowitsch PetscCall(MatDestroy(&intMat)); 24829566063dSJacob Faibussowitsch PetscCall(MatDestroy(&intMat2)); 24839566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&intNodes)); 24849566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&intNodes2)); 24859566063dSJacob Faibussowitsch PetscCall(PetscLagNodeIndicesDestroy(&intNodeIndices)); 24869566063dSJacob Faibussowitsch PetscCall(PetscLagNodeIndicesDestroy(&intNodeIndices2)); 24873f27d899SToby Isaac intNodes = merged; 24883f27d899SToby Isaac intMat = matMerged; 24893f27d899SToby Isaac intNodeIndices = intNodeIndicesMerged; 24903f27d899SToby Isaac if (!trimmed) { 249177f1a120SToby Isaac /* I think users expect that, when a node has a full basis for the k-forms, 249277f1a120SToby Isaac * they should be consecutive dofs. That isn't the case for trimmed spaces, 249377f1a120SToby Isaac * but is for some of the nodes in untrimmed spaces, so in that case we 249477f1a120SToby Isaac * sort them to group them by node */ 24953f27d899SToby Isaac Mat intMatPerm; 24963f27d899SToby Isaac 24979566063dSJacob Faibussowitsch PetscCall(MatPermuteByNodeIdx(intMat, intNodeIndices, &intMatPerm)); 24989566063dSJacob Faibussowitsch PetscCall(MatDestroy(&intMat)); 24993f27d899SToby Isaac intMat = intMatPerm; 25003f27d899SToby Isaac } 25013f27d899SToby Isaac } 25023f27d899SToby Isaac } 25039566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceDestroy(&fiber)); 25049566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceDestroy(&trace)); 25053f27d899SToby Isaac } 25069566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&intNodesTrace)); 25079566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&intNodesFiber)); 25083f27d899SToby Isaac sp->intNodes = intNodes; 25093f27d899SToby Isaac sp->intMat = intMat; 25103f27d899SToby Isaac lag->intNodeIndices = intNodeIndices; 25116f905325SMatthew G. Knepley { 25123f27d899SToby Isaac PetscInt nDofs = 0; 25133f27d899SToby Isaac 25143f27d899SToby Isaac if (intMat) { 25159566063dSJacob Faibussowitsch PetscCall(MatGetSize(intMat, &nDofs, NULL)); 25163f27d899SToby Isaac } 25179566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(section, 0, nDofs)); 25183f27d899SToby Isaac } 25199566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSectionSetUp_Internal(sp, section)); 25203f27d899SToby Isaac if (continuous) { 25219566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceCreateAllDataFromInteriorData(sp)); 25229566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceLagrangeCreateAllNodeIdx(sp)); 25233f27d899SToby Isaac } else { 25249566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject) intNodes)); 25253f27d899SToby Isaac sp->allNodes = intNodes; 25269566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject) intMat)); 25273f27d899SToby Isaac sp->allMat = intMat; 25289566063dSJacob Faibussowitsch PetscCall(PetscLagNodeIndicesReference(intNodeIndices)); 25293f27d899SToby Isaac lag->allNodeIndices = intNodeIndices; 25303f27d899SToby Isaac } 25313f27d899SToby Isaac } 25329566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(section, &sp->spdim)); 25339566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstrainedStorageSize(section, &sp->spintdim)); 25349566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceComputeFunctionalsFromAllData(sp)); 25359566063dSJacob Faibussowitsch PetscCall(PetscFree2(pStratStart, pStratEnd)); 25369566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmint)); 25373f27d899SToby Isaac PetscFunctionReturn(0); 25383f27d899SToby Isaac } 25393f27d899SToby Isaac 254077f1a120SToby Isaac /* Create a matrix that represents the transformation that DMPlexVecGetClosure() would need 254177f1a120SToby Isaac * to get the representation of the dofs for a mesh point if the mesh point had this orientation 254277f1a120SToby Isaac * relative to the cell */ 25433f27d899SToby Isaac PetscErrorCode PetscDualSpaceCreateInteriorSymmetryMatrix_Lagrange(PetscDualSpace sp, PetscInt ornt, Mat *symMat) 25443f27d899SToby Isaac { 25453f27d899SToby Isaac PetscDualSpace_Lag *lag; 25463f27d899SToby Isaac DM dm; 25473f27d899SToby Isaac PetscLagNodeIndices vertIndices, intNodeIndices; 25483f27d899SToby Isaac PetscLagNodeIndices ni; 25493f27d899SToby Isaac PetscInt nodeIdxDim, nodeVecDim, nNodes; 25503f27d899SToby Isaac PetscInt formDegree; 25513f27d899SToby Isaac PetscInt *perm, *permOrnt; 25523f27d899SToby Isaac PetscInt *nnz; 25533f27d899SToby Isaac PetscInt n; 25543f27d899SToby Isaac PetscInt maxGroupSize; 25553f27d899SToby Isaac PetscScalar *V, *W, *work; 25563f27d899SToby Isaac Mat A; 25576f905325SMatthew G. Knepley 25586f905325SMatthew G. Knepley PetscFunctionBegin; 25593f27d899SToby Isaac if (!sp->spintdim) { 25603f27d899SToby Isaac *symMat = NULL; 25613f27d899SToby Isaac PetscFunctionReturn(0); 25626f905325SMatthew G. Knepley } 25633f27d899SToby Isaac lag = (PetscDualSpace_Lag *) sp->data; 25643f27d899SToby Isaac vertIndices = lag->vertIndices; 25653f27d899SToby Isaac intNodeIndices = lag->intNodeIndices; 25669566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDM(sp, &dm)); 25679566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetFormDegree(sp, &formDegree)); 25689566063dSJacob Faibussowitsch PetscCall(PetscNew(&ni)); 25693f27d899SToby Isaac ni->refct = 1; 25703f27d899SToby Isaac ni->nodeIdxDim = nodeIdxDim = intNodeIndices->nodeIdxDim; 25713f27d899SToby Isaac ni->nodeVecDim = nodeVecDim = intNodeIndices->nodeVecDim; 25723f27d899SToby Isaac ni->nNodes = nNodes = intNodeIndices->nNodes; 25739566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nNodes * nodeIdxDim, &(ni->nodeIdx))); 25749566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nNodes * nodeVecDim, &(ni->nodeVec))); 257577f1a120SToby Isaac /* push forward the dofs by the symmetry of the reference element induced by ornt */ 25769566063dSJacob Faibussowitsch PetscCall(PetscLagNodeIndicesPushForward(dm, vertIndices, 0, vertIndices, intNodeIndices, ornt, formDegree, ni->nodeIdx, ni->nodeVec)); 257777f1a120SToby Isaac /* get the revlex order for both the original and transformed dofs */ 25789566063dSJacob Faibussowitsch PetscCall(PetscLagNodeIndicesGetPermutation(intNodeIndices, &perm)); 25799566063dSJacob Faibussowitsch PetscCall(PetscLagNodeIndicesGetPermutation(ni, &permOrnt)); 25809566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nNodes, &nnz)); 25813f27d899SToby Isaac for (n = 0, maxGroupSize = 0; n < nNodes;) { /* incremented in the loop */ 25823f27d899SToby Isaac PetscInt *nind = &(ni->nodeIdx[permOrnt[n] * nodeIdxDim]); 25833f27d899SToby Isaac PetscInt m, nEnd; 25843f27d899SToby Isaac PetscInt groupSize; 258577f1a120SToby Isaac /* for each group of dofs that have the same nodeIdx coordinate */ 25863f27d899SToby Isaac for (nEnd = n + 1; nEnd < nNodes; nEnd++) { 25873f27d899SToby Isaac PetscInt *mind = &(ni->nodeIdx[permOrnt[nEnd] * nodeIdxDim]); 25883f27d899SToby Isaac PetscInt d; 25893f27d899SToby Isaac 25903f27d899SToby Isaac /* compare the oriented permutation indices */ 25913f27d899SToby Isaac for (d = 0; d < nodeIdxDim; d++) if (mind[d] != nind[d]) break; 25923f27d899SToby Isaac if (d < nodeIdxDim) break; 25933f27d899SToby Isaac } 259477f1a120SToby Isaac /* permOrnt[[n, nEnd)] is a group of dofs that, under the symmetry are at the same location */ 259576bd3646SJed Brown 259677f1a120SToby Isaac /* the symmetry had better map the group of dofs with the same permuted nodeIdx 259777f1a120SToby Isaac * to a group of dofs with the same size, otherwise we messed up */ 259876bd3646SJed Brown if (PetscDefined(USE_DEBUG)) { 25993f27d899SToby Isaac PetscInt m; 26003f27d899SToby Isaac PetscInt *nind = &(intNodeIndices->nodeIdx[perm[n] * nodeIdxDim]); 26013f27d899SToby Isaac 26023f27d899SToby Isaac for (m = n + 1; m < nEnd; m++) { 26033f27d899SToby Isaac PetscInt *mind = &(intNodeIndices->nodeIdx[perm[m] * nodeIdxDim]); 26043f27d899SToby Isaac PetscInt d; 26053f27d899SToby Isaac 26063f27d899SToby Isaac /* compare the oriented permutation indices */ 26073f27d899SToby Isaac for (d = 0; d < nodeIdxDim; d++) if (mind[d] != nind[d]) break; 26083f27d899SToby Isaac if (d < nodeIdxDim) break; 26093f27d899SToby Isaac } 26102c71b3e2SJacob Faibussowitsch PetscCheckFalse(m < nEnd,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Dofs with same index after symmetry not same block size"); 26113f27d899SToby Isaac } 26123f27d899SToby Isaac groupSize = nEnd - n; 261377f1a120SToby Isaac /* each pushforward dof vector will be expressed in a basis of the unpermuted dofs */ 26143f27d899SToby Isaac for (m = n; m < nEnd; m++) nnz[permOrnt[m]] = groupSize; 26153f27d899SToby Isaac 26163f27d899SToby Isaac maxGroupSize = PetscMax(maxGroupSize, nEnd - n); 26173f27d899SToby Isaac n = nEnd; 26183f27d899SToby Isaac } 26192c71b3e2SJacob Faibussowitsch PetscCheckFalse(maxGroupSize > nodeVecDim,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Dofs are not in blocks that can be solved"); 26209566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, nNodes, nNodes, 0, nnz, &A)); 26219566063dSJacob Faibussowitsch PetscCall(PetscFree(nnz)); 26229566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(maxGroupSize * nodeVecDim, &V, maxGroupSize * nodeVecDim, &W, nodeVecDim * 2, &work)); 26233f27d899SToby Isaac for (n = 0; n < nNodes;) { /* incremented in the loop */ 26243f27d899SToby Isaac PetscInt *nind = &(ni->nodeIdx[permOrnt[n] * nodeIdxDim]); 26253f27d899SToby Isaac PetscInt nEnd; 26263f27d899SToby Isaac PetscInt m; 26273f27d899SToby Isaac PetscInt groupSize; 26283f27d899SToby Isaac for (nEnd = n + 1; nEnd < nNodes; nEnd++) { 26293f27d899SToby Isaac PetscInt *mind = &(ni->nodeIdx[permOrnt[nEnd] * nodeIdxDim]); 26303f27d899SToby Isaac PetscInt d; 26313f27d899SToby Isaac 26323f27d899SToby Isaac /* compare the oriented permutation indices */ 26333f27d899SToby Isaac for (d = 0; d < nodeIdxDim; d++) if (mind[d] != nind[d]) break; 26343f27d899SToby Isaac if (d < nodeIdxDim) break; 26353f27d899SToby Isaac } 26363f27d899SToby Isaac groupSize = nEnd - n; 263777f1a120SToby Isaac /* get all of the vectors from the original and all of the pushforward vectors */ 26383f27d899SToby Isaac for (m = n; m < nEnd; m++) { 26393f27d899SToby Isaac PetscInt d; 26403f27d899SToby Isaac 26413f27d899SToby Isaac for (d = 0; d < nodeVecDim; d++) { 26423f27d899SToby Isaac V[(m - n) * nodeVecDim + d] = intNodeIndices->nodeVec[perm[m] * nodeVecDim + d]; 26433f27d899SToby Isaac W[(m - n) * nodeVecDim + d] = ni->nodeVec[permOrnt[m] * nodeVecDim + d]; 26443f27d899SToby Isaac } 26453f27d899SToby Isaac } 264677f1a120SToby Isaac /* now we have to solve for W in terms of V: the systems isn't always square, but the span 264777f1a120SToby Isaac * of V and W should always be the same, so the solution of the normal equations works */ 26483f27d899SToby Isaac { 26493f27d899SToby Isaac char transpose = 'N'; 26503f27d899SToby Isaac PetscBLASInt bm = nodeVecDim; 26513f27d899SToby Isaac PetscBLASInt bn = groupSize; 26523f27d899SToby Isaac PetscBLASInt bnrhs = groupSize; 26533f27d899SToby Isaac PetscBLASInt blda = bm; 26543f27d899SToby Isaac PetscBLASInt bldb = bm; 26553f27d899SToby Isaac PetscBLASInt blwork = 2 * nodeVecDim; 26563f27d899SToby Isaac PetscBLASInt info; 26573f27d899SToby Isaac 26583f27d899SToby Isaac PetscStackCallBLAS("LAPACKgels",LAPACKgels_(&transpose,&bm,&bn,&bnrhs,V,&blda,W,&bldb,work,&blwork, &info)); 26592c71b3e2SJacob Faibussowitsch PetscCheckFalse(info != 0,PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to GELS"); 26603f27d899SToby Isaac /* repack */ 26613f27d899SToby Isaac { 26623f27d899SToby Isaac PetscInt i, j; 26633f27d899SToby Isaac 26643f27d899SToby Isaac for (i = 0; i < groupSize; i++) { 26653f27d899SToby Isaac for (j = 0; j < groupSize; j++) { 266677f1a120SToby Isaac /* notice the different leading dimension */ 26673f27d899SToby Isaac V[i * groupSize + j] = W[i * nodeVecDim + j]; 26683f27d899SToby Isaac } 26693f27d899SToby Isaac } 26703f27d899SToby Isaac } 2671c5c386beSToby Isaac if (PetscDefined(USE_DEBUG)) { 2672c5c386beSToby Isaac PetscReal res; 2673c5c386beSToby Isaac 2674c5c386beSToby Isaac /* check that the normal error is 0 */ 2675c5c386beSToby Isaac for (m = n; m < nEnd; m++) { 2676c5c386beSToby Isaac PetscInt d; 2677c5c386beSToby Isaac 2678c5c386beSToby Isaac for (d = 0; d < nodeVecDim; d++) { 2679c5c386beSToby Isaac W[(m - n) * nodeVecDim + d] = ni->nodeVec[permOrnt[m] * nodeVecDim + d]; 2680c5c386beSToby Isaac } 2681c5c386beSToby Isaac } 2682c5c386beSToby Isaac res = 0.; 2683c5c386beSToby Isaac for (PetscInt i = 0; i < groupSize; i++) { 2684c5c386beSToby Isaac for (PetscInt j = 0; j < nodeVecDim; j++) { 2685c5c386beSToby Isaac for (PetscInt k = 0; k < groupSize; k++) { 2686c5c386beSToby Isaac W[i * nodeVecDim + j] -= V[i * groupSize + k] * intNodeIndices->nodeVec[perm[n+k] * nodeVecDim + j]; 2687c5c386beSToby Isaac } 2688c5c386beSToby Isaac res += PetscAbsScalar(W[i * nodeVecDim + j]); 2689c5c386beSToby Isaac } 2690c5c386beSToby Isaac } 26912c71b3e2SJacob Faibussowitsch PetscCheckFalse(res > PETSC_SMALL,PETSC_COMM_SELF,PETSC_ERR_LIB,"Dof block did not solve"); 2692c5c386beSToby Isaac } 26933f27d899SToby Isaac } 26949566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, groupSize, &permOrnt[n], groupSize, &perm[n], V, INSERT_VALUES)); 26953f27d899SToby Isaac n = nEnd; 26963f27d899SToby Isaac } 26979566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 26989566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 26993f27d899SToby Isaac *symMat = A; 27009566063dSJacob Faibussowitsch PetscCall(PetscFree3(V,W,work)); 27019566063dSJacob Faibussowitsch PetscCall(PetscLagNodeIndicesDestroy(&ni)); 27026f905325SMatthew G. Knepley PetscFunctionReturn(0); 27036f905325SMatthew G. Knepley } 270420cf1dd8SToby Isaac 270520cf1dd8SToby Isaac #define BaryIndex(perEdge,a,b,c) (((b)*(2*perEdge+1-(b)))/2)+(c) 270620cf1dd8SToby Isaac 270720cf1dd8SToby Isaac #define CartIndex(perEdge,a,b) (perEdge*(a)+b) 270820cf1dd8SToby Isaac 270977f1a120SToby Isaac /* the existing interface for symmetries is insufficient for all cases: 271077f1a120SToby Isaac * - it should be sufficient for form degrees that are scalar (0 and n) 271177f1a120SToby Isaac * - it should be sufficient for hypercube dofs 271277f1a120SToby Isaac * - it isn't sufficient for simplex cells with non-scalar form degrees if 271377f1a120SToby Isaac * there are any dofs in the interior 271477f1a120SToby Isaac * 271577f1a120SToby Isaac * We compute the general transformation matrices, and if they fit, we return them, 271677f1a120SToby Isaac * otherwise we error (but we should probably change the interface to allow for 271777f1a120SToby Isaac * these symmetries) 271877f1a120SToby Isaac */ 271920cf1dd8SToby Isaac static PetscErrorCode PetscDualSpaceGetSymmetries_Lagrange(PetscDualSpace sp, const PetscInt ****perms, const PetscScalar ****flips) 272020cf1dd8SToby Isaac { 272120cf1dd8SToby Isaac PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 27223f27d899SToby Isaac PetscInt dim, order, Nc; 272320cf1dd8SToby Isaac 272420cf1dd8SToby Isaac PetscFunctionBegin; 27259566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetOrder(sp,&order)); 27269566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetNumComponents(sp,&Nc)); 27279566063dSJacob Faibussowitsch PetscCall(DMGetDimension(sp->dm,&dim)); 27283f27d899SToby Isaac if (!lag->symComputed) { /* store symmetries */ 27293f27d899SToby Isaac PetscInt pStart, pEnd, p; 27303f27d899SToby Isaac PetscInt numPoints; 273120cf1dd8SToby Isaac PetscInt numFaces; 27323f27d899SToby Isaac PetscInt spintdim; 27333f27d899SToby Isaac PetscInt ***symperms; 27343f27d899SToby Isaac PetscScalar ***symflips; 273520cf1dd8SToby Isaac 27369566063dSJacob Faibussowitsch PetscCall(DMPlexGetChart(sp->dm, &pStart, &pEnd)); 27373f27d899SToby Isaac numPoints = pEnd - pStart; 2738b5a892a1SMatthew G. Knepley { 2739b5a892a1SMatthew G. Knepley DMPolytopeType ct; 2740b5a892a1SMatthew G. Knepley /* The number of arrangements is no longer based on the number of faces */ 27419566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(sp->dm, 0, &ct)); 2742b5a892a1SMatthew G. Knepley numFaces = DMPolytopeTypeGetNumArrangments(ct) / 2; 2743b5a892a1SMatthew G. Knepley } 27449566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(numPoints,&symperms)); 27459566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(numPoints,&symflips)); 27463f27d899SToby Isaac spintdim = sp->spintdim; 27473f27d899SToby Isaac /* The nodal symmetry behavior is not present when tensorSpace != tensorCell: someone might want this for the "S" 27483f27d899SToby Isaac * family of FEEC spaces. Most used in particular are discontinuous polynomial L2 spaces in tensor cells, where 27493f27d899SToby Isaac * the symmetries are not necessary for FE assembly. So for now we assume this is the case and don't return 27503f27d899SToby Isaac * symmetries if tensorSpace != tensorCell */ 27513f27d899SToby Isaac if (spintdim && 0 < dim && dim < 3 && (lag->tensorSpace == lag->tensorCell)) { /* compute self symmetries */ 27523f27d899SToby Isaac PetscInt **cellSymperms; 27533f27d899SToby Isaac PetscScalar **cellSymflips; 27543f27d899SToby Isaac PetscInt ornt; 27553f27d899SToby Isaac PetscInt nCopies = Nc / lag->intNodeIndices->nodeVecDim; 27563f27d899SToby Isaac PetscInt nNodes = lag->intNodeIndices->nNodes; 275720cf1dd8SToby Isaac 275820cf1dd8SToby Isaac lag->numSelfSym = 2 * numFaces; 275920cf1dd8SToby Isaac lag->selfSymOff = numFaces; 27609566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(2*numFaces,&cellSymperms)); 27619566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(2*numFaces,&cellSymflips)); 276220cf1dd8SToby Isaac /* we want to be able to index symmetries directly with the orientations, which range from [-numFaces,numFaces) */ 27633f27d899SToby Isaac symperms[0] = &cellSymperms[numFaces]; 27643f27d899SToby Isaac symflips[0] = &cellSymflips[numFaces]; 27652c71b3e2SJacob Faibussowitsch PetscCheckFalse(lag->intNodeIndices->nodeVecDim * nCopies != Nc,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Node indices incompatible with dofs"); 27662c71b3e2SJacob Faibussowitsch PetscCheckFalse(nNodes * nCopies != spintdim,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Node indices incompatible with dofs"); 27673f27d899SToby 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 */ 27683f27d899SToby Isaac Mat symMat; 27693f27d899SToby Isaac PetscInt *perm; 27703f27d899SToby Isaac PetscScalar *flips; 27713f27d899SToby Isaac PetscInt i; 277220cf1dd8SToby Isaac 27733f27d899SToby Isaac if (!ornt) continue; 27749566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(spintdim, &perm)); 27759566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(spintdim, &flips)); 27763f27d899SToby Isaac for (i = 0; i < spintdim; i++) perm[i] = -1; 27779566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceCreateInteriorSymmetryMatrix_Lagrange(sp, ornt, &symMat)); 27783f27d899SToby Isaac for (i = 0; i < nNodes; i++) { 27793f27d899SToby Isaac PetscInt ncols; 27803f27d899SToby Isaac PetscInt j, k; 27813f27d899SToby Isaac const PetscInt *cols; 27823f27d899SToby Isaac const PetscScalar *vals; 27833f27d899SToby Isaac PetscBool nz_seen = PETSC_FALSE; 278420cf1dd8SToby Isaac 27859566063dSJacob Faibussowitsch PetscCall(MatGetRow(symMat, i, &ncols, &cols, &vals)); 27863f27d899SToby Isaac for (j = 0; j < ncols; j++) { 27873f27d899SToby Isaac if (PetscAbsScalar(vals[j]) > PETSC_SMALL) { 278828b400f6SJacob 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"); 27893f27d899SToby Isaac nz_seen = PETSC_TRUE; 27902c71b3e2SJacob Faibussowitsch PetscCheckFalse(PetscAbsReal(PetscAbsScalar(vals[j]) - PetscRealConstant(1.)) > PETSC_SMALL,PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "This dual space has symmetries that can't be described as a permutation + sign flips"); 27912c71b3e2SJacob Faibussowitsch PetscCheckFalse(PetscAbsReal(PetscImaginaryPart(vals[j])) > PETSC_SMALL,PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "This dual space has symmetries that can't be described as a permutation + sign flips"); 27922c71b3e2SJacob Faibussowitsch PetscCheckFalse(perm[cols[j] * nCopies] >= 0,PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "This dual space has symmetries that can't be described as a permutation + sign flips"); 27933f27d899SToby Isaac for (k = 0; k < nCopies; k++) { 27943f27d899SToby Isaac perm[cols[j] * nCopies + k] = i * nCopies + k; 279520cf1dd8SToby Isaac } 27963f27d899SToby Isaac if (PetscRealPart(vals[j]) < 0.) { 27973f27d899SToby Isaac for (k = 0; k < nCopies; k++) { 27983f27d899SToby Isaac flips[i * nCopies + k] = -1.; 279920cf1dd8SToby Isaac } 280020cf1dd8SToby Isaac } else { 28013f27d899SToby Isaac for (k = 0; k < nCopies; k++) { 28023f27d899SToby Isaac flips[i * nCopies + k] = 1.; 28033f27d899SToby Isaac } 28043f27d899SToby Isaac } 28053f27d899SToby Isaac } 28063f27d899SToby Isaac } 28079566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(symMat, i, &ncols, &cols, &vals)); 28083f27d899SToby Isaac } 28099566063dSJacob Faibussowitsch PetscCall(MatDestroy(&symMat)); 28103f27d899SToby Isaac /* if there were no sign flips, keep NULL */ 28113f27d899SToby Isaac for (i = 0; i < spintdim; i++) if (flips[i] != 1.) break; 28123f27d899SToby Isaac if (i == spintdim) { 28139566063dSJacob Faibussowitsch PetscCall(PetscFree(flips)); 28143f27d899SToby Isaac flips = NULL; 28153f27d899SToby Isaac } 28163f27d899SToby Isaac /* if the permutation is identity, keep NULL */ 28173f27d899SToby Isaac for (i = 0; i < spintdim; i++) if (perm[i] != i) break; 28183f27d899SToby Isaac if (i == spintdim) { 28199566063dSJacob Faibussowitsch PetscCall(PetscFree(perm)); 28203f27d899SToby Isaac perm = NULL; 28213f27d899SToby Isaac } 28223f27d899SToby Isaac symperms[0][ornt] = perm; 28233f27d899SToby Isaac symflips[0][ornt] = flips; 28243f27d899SToby Isaac } 28253f27d899SToby Isaac /* if no orientations produced non-identity permutations, keep NULL */ 28263f27d899SToby Isaac for (ornt = -numFaces; ornt < numFaces; ornt++) if (symperms[0][ornt]) break; 28273f27d899SToby Isaac if (ornt == numFaces) { 28289566063dSJacob Faibussowitsch PetscCall(PetscFree(cellSymperms)); 28293f27d899SToby Isaac symperms[0] = NULL; 28303f27d899SToby Isaac } 28313f27d899SToby Isaac /* if no orientations produced sign flips, keep NULL */ 28323f27d899SToby Isaac for (ornt = -numFaces; ornt < numFaces; ornt++) if (symflips[0][ornt]) break; 28333f27d899SToby Isaac if (ornt == numFaces) { 28349566063dSJacob Faibussowitsch PetscCall(PetscFree(cellSymflips)); 28353f27d899SToby Isaac symflips[0] = NULL; 28363f27d899SToby Isaac } 28373f27d899SToby Isaac } 283877f1a120SToby Isaac { /* get the symmetries of closure points */ 28393f27d899SToby Isaac PetscInt closureSize = 0; 28403f27d899SToby Isaac PetscInt *closure = NULL; 28413f27d899SToby Isaac PetscInt r; 284220cf1dd8SToby Isaac 28439566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(sp->dm,0,PETSC_TRUE,&closureSize,&closure)); 28443f27d899SToby Isaac for (r = 0; r < closureSize; r++) { 28453f27d899SToby Isaac PetscDualSpace psp; 28463f27d899SToby Isaac PetscInt point = closure[2 * r]; 28473f27d899SToby Isaac PetscInt pspintdim; 28483f27d899SToby Isaac const PetscInt ***psymperms = NULL; 28493f27d899SToby Isaac const PetscScalar ***psymflips = NULL; 285020cf1dd8SToby Isaac 28513f27d899SToby Isaac if (!point) continue; 28529566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetPointSubspace(sp, point, &psp)); 28533f27d899SToby Isaac if (!psp) continue; 28549566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetInteriorDimension(psp, &pspintdim)); 28553f27d899SToby Isaac if (!pspintdim) continue; 28569566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetSymmetries(psp,&psymperms,&psymflips)); 28573f27d899SToby Isaac symperms[r] = (PetscInt **) (psymperms ? psymperms[0] : NULL); 28583f27d899SToby Isaac symflips[r] = (PetscScalar **) (psymflips ? psymflips[0] : NULL); 285920cf1dd8SToby Isaac } 28609566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(sp->dm,0,PETSC_TRUE,&closureSize,&closure)); 286120cf1dd8SToby Isaac } 28623f27d899SToby Isaac for (p = 0; p < pEnd; p++) if (symperms[p]) break; 28633f27d899SToby Isaac if (p == pEnd) { 28649566063dSJacob Faibussowitsch PetscCall(PetscFree(symperms)); 28653f27d899SToby Isaac symperms = NULL; 286620cf1dd8SToby Isaac } 28673f27d899SToby Isaac for (p = 0; p < pEnd; p++) if (symflips[p]) break; 28683f27d899SToby Isaac if (p == pEnd) { 28699566063dSJacob Faibussowitsch PetscCall(PetscFree(symflips)); 28703f27d899SToby Isaac symflips = NULL; 287120cf1dd8SToby Isaac } 28723f27d899SToby Isaac lag->symperms = symperms; 28733f27d899SToby Isaac lag->symflips = symflips; 28743f27d899SToby Isaac lag->symComputed = PETSC_TRUE; 287520cf1dd8SToby Isaac } 28763f27d899SToby Isaac if (perms) *perms = (const PetscInt ***) lag->symperms; 28773f27d899SToby Isaac if (flips) *flips = (const PetscScalar ***) lag->symflips; 287820cf1dd8SToby Isaac PetscFunctionReturn(0); 287920cf1dd8SToby Isaac } 288020cf1dd8SToby Isaac 288120cf1dd8SToby Isaac static PetscErrorCode PetscDualSpaceLagrangeGetContinuity_Lagrange(PetscDualSpace sp, PetscBool *continuous) 288220cf1dd8SToby Isaac { 288320cf1dd8SToby Isaac PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 288420cf1dd8SToby Isaac 288520cf1dd8SToby Isaac PetscFunctionBegin; 288620cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 2887dadcf809SJacob Faibussowitsch PetscValidBoolPointer(continuous, 2); 288820cf1dd8SToby Isaac *continuous = lag->continuous; 288920cf1dd8SToby Isaac PetscFunctionReturn(0); 289020cf1dd8SToby Isaac } 289120cf1dd8SToby Isaac 289220cf1dd8SToby Isaac static PetscErrorCode PetscDualSpaceLagrangeSetContinuity_Lagrange(PetscDualSpace sp, PetscBool continuous) 289320cf1dd8SToby Isaac { 289420cf1dd8SToby Isaac PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data; 289520cf1dd8SToby Isaac 289620cf1dd8SToby Isaac PetscFunctionBegin; 289720cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 289820cf1dd8SToby Isaac lag->continuous = continuous; 289920cf1dd8SToby Isaac PetscFunctionReturn(0); 290020cf1dd8SToby Isaac } 290120cf1dd8SToby Isaac 290220cf1dd8SToby Isaac /*@ 290320cf1dd8SToby Isaac PetscDualSpaceLagrangeGetContinuity - Retrieves the flag for element continuity 290420cf1dd8SToby Isaac 290520cf1dd8SToby Isaac Not Collective 290620cf1dd8SToby Isaac 290720cf1dd8SToby Isaac Input Parameter: 290820cf1dd8SToby Isaac . sp - the PetscDualSpace 290920cf1dd8SToby Isaac 291020cf1dd8SToby Isaac Output Parameter: 291120cf1dd8SToby Isaac . continuous - flag for element continuity 291220cf1dd8SToby Isaac 291320cf1dd8SToby Isaac Level: intermediate 291420cf1dd8SToby Isaac 291520cf1dd8SToby Isaac .seealso: PetscDualSpaceLagrangeSetContinuity() 291620cf1dd8SToby Isaac @*/ 291720cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceLagrangeGetContinuity(PetscDualSpace sp, PetscBool *continuous) 291820cf1dd8SToby Isaac { 291920cf1dd8SToby Isaac PetscFunctionBegin; 292020cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 2921dadcf809SJacob Faibussowitsch PetscValidBoolPointer(continuous, 2); 2922*cac4c232SBarry Smith PetscTryMethod(sp, "PetscDualSpaceLagrangeGetContinuity_C", (PetscDualSpace,PetscBool*),(sp,continuous)); 292320cf1dd8SToby Isaac PetscFunctionReturn(0); 292420cf1dd8SToby Isaac } 292520cf1dd8SToby Isaac 292620cf1dd8SToby Isaac /*@ 292720cf1dd8SToby Isaac PetscDualSpaceLagrangeSetContinuity - Indicate whether the element is continuous 292820cf1dd8SToby Isaac 2929d083f849SBarry Smith Logically Collective on sp 293020cf1dd8SToby Isaac 293120cf1dd8SToby Isaac Input Parameters: 293220cf1dd8SToby Isaac + sp - the PetscDualSpace 293320cf1dd8SToby Isaac - continuous - flag for element continuity 293420cf1dd8SToby Isaac 293520cf1dd8SToby Isaac Options Database: 2936147403d9SBarry Smith . -petscdualspace_lagrange_continuity <bool> - use a continuous element 293720cf1dd8SToby Isaac 293820cf1dd8SToby Isaac Level: intermediate 293920cf1dd8SToby Isaac 294020cf1dd8SToby Isaac .seealso: PetscDualSpaceLagrangeGetContinuity() 294120cf1dd8SToby Isaac @*/ 294220cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceLagrangeSetContinuity(PetscDualSpace sp, PetscBool continuous) 294320cf1dd8SToby Isaac { 294420cf1dd8SToby Isaac PetscFunctionBegin; 294520cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 294620cf1dd8SToby Isaac PetscValidLogicalCollectiveBool(sp, continuous, 2); 2947*cac4c232SBarry Smith PetscTryMethod(sp, "PetscDualSpaceLagrangeSetContinuity_C", (PetscDualSpace,PetscBool),(sp,continuous)); 294820cf1dd8SToby Isaac PetscFunctionReturn(0); 294920cf1dd8SToby Isaac } 295020cf1dd8SToby Isaac 29516f905325SMatthew G. Knepley static PetscErrorCode PetscDualSpaceLagrangeGetTensor_Lagrange(PetscDualSpace sp, PetscBool *tensor) 295220cf1dd8SToby Isaac { 295320cf1dd8SToby Isaac PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data; 29546f905325SMatthew G. Knepley 29556f905325SMatthew G. Knepley PetscFunctionBegin; 29566f905325SMatthew G. Knepley *tensor = lag->tensorSpace; 29576f905325SMatthew G. Knepley PetscFunctionReturn(0); 29586f905325SMatthew G. Knepley } 29596f905325SMatthew G. Knepley 29606f905325SMatthew G. Knepley static PetscErrorCode PetscDualSpaceLagrangeSetTensor_Lagrange(PetscDualSpace sp, PetscBool tensor) 29616f905325SMatthew G. Knepley { 29626f905325SMatthew G. Knepley PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data; 29636f905325SMatthew G. Knepley 29646f905325SMatthew G. Knepley PetscFunctionBegin; 29656f905325SMatthew G. Knepley lag->tensorSpace = tensor; 29666f905325SMatthew G. Knepley PetscFunctionReturn(0); 29676f905325SMatthew G. Knepley } 29686f905325SMatthew G. Knepley 29693f27d899SToby Isaac static PetscErrorCode PetscDualSpaceLagrangeGetTrimmed_Lagrange(PetscDualSpace sp, PetscBool *trimmed) 29703f27d899SToby Isaac { 29713f27d899SToby Isaac PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data; 29723f27d899SToby Isaac 29733f27d899SToby Isaac PetscFunctionBegin; 29743f27d899SToby Isaac *trimmed = lag->trimmed; 29753f27d899SToby Isaac PetscFunctionReturn(0); 29763f27d899SToby Isaac } 29773f27d899SToby Isaac 29783f27d899SToby Isaac static PetscErrorCode PetscDualSpaceLagrangeSetTrimmed_Lagrange(PetscDualSpace sp, PetscBool trimmed) 29793f27d899SToby Isaac { 29803f27d899SToby Isaac PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data; 29813f27d899SToby Isaac 29823f27d899SToby Isaac PetscFunctionBegin; 29833f27d899SToby Isaac lag->trimmed = trimmed; 29843f27d899SToby Isaac PetscFunctionReturn(0); 29853f27d899SToby Isaac } 29863f27d899SToby Isaac 29873f27d899SToby Isaac static PetscErrorCode PetscDualSpaceLagrangeGetNodeType_Lagrange(PetscDualSpace sp, PetscDTNodeType *nodeType, PetscBool *boundary, PetscReal *exponent) 29883f27d899SToby Isaac { 29893f27d899SToby Isaac PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data; 29903f27d899SToby Isaac 29913f27d899SToby Isaac PetscFunctionBegin; 29923f27d899SToby Isaac if (nodeType) *nodeType = lag->nodeType; 29933f27d899SToby Isaac if (boundary) *boundary = lag->endNodes; 29943f27d899SToby Isaac if (exponent) *exponent = lag->nodeExponent; 29953f27d899SToby Isaac PetscFunctionReturn(0); 29963f27d899SToby Isaac } 29973f27d899SToby Isaac 29983f27d899SToby Isaac static PetscErrorCode PetscDualSpaceLagrangeSetNodeType_Lagrange(PetscDualSpace sp, PetscDTNodeType nodeType, PetscBool boundary, PetscReal exponent) 29993f27d899SToby Isaac { 30003f27d899SToby Isaac PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data; 30013f27d899SToby Isaac 30023f27d899SToby Isaac PetscFunctionBegin; 30032c71b3e2SJacob Faibussowitsch PetscCheckFalse(nodeType == PETSCDTNODES_GAUSSJACOBI && exponent <= -1.,PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_OUTOFRANGE, "Exponent must be > -1"); 30043f27d899SToby Isaac lag->nodeType = nodeType; 30053f27d899SToby Isaac lag->endNodes = boundary; 30063f27d899SToby Isaac lag->nodeExponent = exponent; 30073f27d899SToby Isaac PetscFunctionReturn(0); 30083f27d899SToby Isaac } 30093f27d899SToby Isaac 301066a6c23cSMatthew G. Knepley static PetscErrorCode PetscDualSpaceLagrangeGetUseMoments_Lagrange(PetscDualSpace sp, PetscBool *useMoments) 301166a6c23cSMatthew G. Knepley { 301266a6c23cSMatthew G. Knepley PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data; 301366a6c23cSMatthew G. Knepley 301466a6c23cSMatthew G. Knepley PetscFunctionBegin; 301566a6c23cSMatthew G. Knepley *useMoments = lag->useMoments; 301666a6c23cSMatthew G. Knepley PetscFunctionReturn(0); 301766a6c23cSMatthew G. Knepley } 301866a6c23cSMatthew G. Knepley 301966a6c23cSMatthew G. Knepley static PetscErrorCode PetscDualSpaceLagrangeSetUseMoments_Lagrange(PetscDualSpace sp, PetscBool useMoments) 302066a6c23cSMatthew G. Knepley { 302166a6c23cSMatthew G. Knepley PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data; 302266a6c23cSMatthew G. Knepley 302366a6c23cSMatthew G. Knepley PetscFunctionBegin; 302466a6c23cSMatthew G. Knepley lag->useMoments = useMoments; 302566a6c23cSMatthew G. Knepley PetscFunctionReturn(0); 302666a6c23cSMatthew G. Knepley } 302766a6c23cSMatthew G. Knepley 302866a6c23cSMatthew G. Knepley static PetscErrorCode PetscDualSpaceLagrangeGetMomentOrder_Lagrange(PetscDualSpace sp, PetscInt *momentOrder) 302966a6c23cSMatthew G. Knepley { 303066a6c23cSMatthew G. Knepley PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data; 303166a6c23cSMatthew G. Knepley 303266a6c23cSMatthew G. Knepley PetscFunctionBegin; 303366a6c23cSMatthew G. Knepley *momentOrder = lag->momentOrder; 303466a6c23cSMatthew G. Knepley PetscFunctionReturn(0); 303566a6c23cSMatthew G. Knepley } 303666a6c23cSMatthew G. Knepley 303766a6c23cSMatthew G. Knepley static PetscErrorCode PetscDualSpaceLagrangeSetMomentOrder_Lagrange(PetscDualSpace sp, PetscInt momentOrder) 303866a6c23cSMatthew G. Knepley { 303966a6c23cSMatthew G. Knepley PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data; 304066a6c23cSMatthew G. Knepley 304166a6c23cSMatthew G. Knepley PetscFunctionBegin; 304266a6c23cSMatthew G. Knepley lag->momentOrder = momentOrder; 304366a6c23cSMatthew G. Knepley PetscFunctionReturn(0); 304466a6c23cSMatthew G. Knepley } 304566a6c23cSMatthew G. Knepley 30466f905325SMatthew G. Knepley /*@ 30476f905325SMatthew G. Knepley PetscDualSpaceLagrangeGetTensor - Get the tensor nature of the dual space 30486f905325SMatthew G. Knepley 30496f905325SMatthew G. Knepley Not collective 30506f905325SMatthew G. Knepley 30516f905325SMatthew G. Knepley Input Parameter: 30526f905325SMatthew G. Knepley . sp - The PetscDualSpace 30536f905325SMatthew G. Knepley 30546f905325SMatthew G. Knepley Output Parameter: 30556f905325SMatthew G. Knepley . tensor - Whether the dual space has tensor layout (vs. simplicial) 30566f905325SMatthew G. Knepley 30576f905325SMatthew G. Knepley Level: intermediate 30586f905325SMatthew G. Knepley 30596f905325SMatthew G. Knepley .seealso: PetscDualSpaceLagrangeSetTensor(), PetscDualSpaceCreate() 30606f905325SMatthew G. Knepley @*/ 30616f905325SMatthew G. Knepley PetscErrorCode PetscDualSpaceLagrangeGetTensor(PetscDualSpace sp, PetscBool *tensor) 30626f905325SMatthew G. Knepley { 306320cf1dd8SToby Isaac PetscFunctionBegin; 306420cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 3065dadcf809SJacob Faibussowitsch PetscValidBoolPointer(tensor, 2); 3066*cac4c232SBarry Smith PetscTryMethod(sp,"PetscDualSpaceLagrangeGetTensor_C",(PetscDualSpace,PetscBool *),(sp,tensor)); 306720cf1dd8SToby Isaac PetscFunctionReturn(0); 306820cf1dd8SToby Isaac } 306920cf1dd8SToby Isaac 30706f905325SMatthew G. Knepley /*@ 30716f905325SMatthew G. Knepley PetscDualSpaceLagrangeSetTensor - Set the tensor nature of the dual space 30726f905325SMatthew G. Knepley 30736f905325SMatthew G. Knepley Not collective 30746f905325SMatthew G. Knepley 30756f905325SMatthew G. Knepley Input Parameters: 30766f905325SMatthew G. Knepley + sp - The PetscDualSpace 30776f905325SMatthew G. Knepley - tensor - Whether the dual space has tensor layout (vs. simplicial) 30786f905325SMatthew G. Knepley 30796f905325SMatthew G. Knepley Level: intermediate 30806f905325SMatthew G. Knepley 30816f905325SMatthew G. Knepley .seealso: PetscDualSpaceLagrangeGetTensor(), PetscDualSpaceCreate() 30826f905325SMatthew G. Knepley @*/ 30836f905325SMatthew G. Knepley PetscErrorCode PetscDualSpaceLagrangeSetTensor(PetscDualSpace sp, PetscBool tensor) 30846f905325SMatthew G. Knepley { 30856f905325SMatthew G. Knepley PetscFunctionBegin; 30866f905325SMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 3087*cac4c232SBarry Smith PetscTryMethod(sp,"PetscDualSpaceLagrangeSetTensor_C",(PetscDualSpace,PetscBool),(sp,tensor)); 30886f905325SMatthew G. Knepley PetscFunctionReturn(0); 30896f905325SMatthew G. Knepley } 30906f905325SMatthew G. Knepley 30913f27d899SToby Isaac /*@ 30923f27d899SToby Isaac PetscDualSpaceLagrangeGetTrimmed - Get the trimmed nature of the dual space 30933f27d899SToby Isaac 30943f27d899SToby Isaac Not collective 30953f27d899SToby Isaac 30963f27d899SToby Isaac Input Parameter: 30973f27d899SToby Isaac . sp - The PetscDualSpace 30983f27d899SToby Isaac 30993f27d899SToby Isaac Output Parameter: 31003f27d899SToby 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) 31013f27d899SToby Isaac 31023f27d899SToby Isaac Level: intermediate 31033f27d899SToby Isaac 31043f27d899SToby Isaac .seealso: PetscDualSpaceLagrangeSetTrimmed(), PetscDualSpaceCreate() 31053f27d899SToby Isaac @*/ 31063f27d899SToby Isaac PetscErrorCode PetscDualSpaceLagrangeGetTrimmed(PetscDualSpace sp, PetscBool *trimmed) 31073f27d899SToby Isaac { 31083f27d899SToby Isaac PetscFunctionBegin; 31093f27d899SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 3110dadcf809SJacob Faibussowitsch PetscValidBoolPointer(trimmed, 2); 3111*cac4c232SBarry Smith PetscTryMethod(sp,"PetscDualSpaceLagrangeGetTrimmed_C",(PetscDualSpace,PetscBool *),(sp,trimmed)); 31123f27d899SToby Isaac PetscFunctionReturn(0); 31133f27d899SToby Isaac } 31143f27d899SToby Isaac 31153f27d899SToby Isaac /*@ 31163f27d899SToby Isaac PetscDualSpaceLagrangeSetTrimmed - Set the trimmed nature of the dual space 31173f27d899SToby Isaac 31183f27d899SToby Isaac Not collective 31193f27d899SToby Isaac 31203f27d899SToby Isaac Input Parameters: 31213f27d899SToby Isaac + sp - The PetscDualSpace 31223f27d899SToby 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) 31233f27d899SToby Isaac 31243f27d899SToby Isaac Level: intermediate 31253f27d899SToby Isaac 31263f27d899SToby Isaac .seealso: PetscDualSpaceLagrangeGetTrimmed(), PetscDualSpaceCreate() 31273f27d899SToby Isaac @*/ 31283f27d899SToby Isaac PetscErrorCode PetscDualSpaceLagrangeSetTrimmed(PetscDualSpace sp, PetscBool trimmed) 31293f27d899SToby Isaac { 31303f27d899SToby Isaac PetscFunctionBegin; 31313f27d899SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 3132*cac4c232SBarry Smith PetscTryMethod(sp,"PetscDualSpaceLagrangeSetTrimmed_C",(PetscDualSpace,PetscBool),(sp,trimmed)); 31333f27d899SToby Isaac PetscFunctionReturn(0); 31343f27d899SToby Isaac } 31353f27d899SToby Isaac 31363f27d899SToby Isaac /*@ 31373f27d899SToby Isaac PetscDualSpaceLagrangeGetNodeType - Get a description of how nodes are laid out for Lagrange polynomials in this 31383f27d899SToby Isaac dual space 31393f27d899SToby Isaac 31403f27d899SToby Isaac Not collective 31413f27d899SToby Isaac 31423f27d899SToby Isaac Input Parameter: 31433f27d899SToby Isaac . sp - The PetscDualSpace 31443f27d899SToby Isaac 31453f27d899SToby Isaac Output Parameters: 31463f27d899SToby Isaac + nodeType - The type of nodes 31473f27d899SToby Isaac . boundary - Whether the node type is one that includes endpoints (if nodeType is PETSCDTNODES_GAUSSJACOBI, nodes that 31483f27d899SToby Isaac include the boundary are Gauss-Lobatto-Jacobi nodes) 31493f27d899SToby Isaac - exponent - If nodeType is PETSCDTNODES_GAUSJACOBI, indicates the exponent used for both ends of the 1D Jacobi weight function 31503f27d899SToby Isaac '0' is Gauss-Legendre, '-0.5' is Gauss-Chebyshev of the first type, '0.5' is Gauss-Chebyshev of the second type 31513f27d899SToby Isaac 31523f27d899SToby Isaac Level: advanced 31533f27d899SToby Isaac 31543f27d899SToby Isaac .seealso: PetscDTNodeType, PetscDualSpaceLagrangeSetNodeType() 31553f27d899SToby Isaac @*/ 31563f27d899SToby Isaac PetscErrorCode PetscDualSpaceLagrangeGetNodeType(PetscDualSpace sp, PetscDTNodeType *nodeType, PetscBool *boundary, PetscReal *exponent) 31573f27d899SToby Isaac { 31583f27d899SToby Isaac PetscFunctionBegin; 31593f27d899SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 31603f27d899SToby Isaac if (nodeType) PetscValidPointer(nodeType, 2); 3161dadcf809SJacob Faibussowitsch if (boundary) PetscValidBoolPointer(boundary, 3); 3162dadcf809SJacob Faibussowitsch if (exponent) PetscValidRealPointer(exponent, 4); 3163*cac4c232SBarry Smith PetscTryMethod(sp,"PetscDualSpaceLagrangeGetNodeType_C",(PetscDualSpace,PetscDTNodeType *,PetscBool *,PetscReal *),(sp,nodeType,boundary,exponent)); 31643f27d899SToby Isaac PetscFunctionReturn(0); 31653f27d899SToby Isaac } 31663f27d899SToby Isaac 31673f27d899SToby Isaac /*@ 31683f27d899SToby Isaac PetscDualSpaceLagrangeSetNodeType - Set a description of how nodes are laid out for Lagrange polynomials in this 31693f27d899SToby Isaac dual space 31703f27d899SToby Isaac 31713f27d899SToby Isaac Logically collective 31723f27d899SToby Isaac 31733f27d899SToby Isaac Input Parameters: 31743f27d899SToby Isaac + sp - The PetscDualSpace 31753f27d899SToby Isaac . nodeType - The type of nodes 31763f27d899SToby Isaac . boundary - Whether the node type is one that includes endpoints (if nodeType is PETSCDTNODES_GAUSSJACOBI, nodes that 31773f27d899SToby Isaac include the boundary are Gauss-Lobatto-Jacobi nodes) 31783f27d899SToby Isaac - exponent - If nodeType is PETSCDTNODES_GAUSJACOBI, indicates the exponent used for both ends of the 1D Jacobi weight function 31793f27d899SToby Isaac '0' is Gauss-Legendre, '-0.5' is Gauss-Chebyshev of the first type, '0.5' is Gauss-Chebyshev of the second type 31803f27d899SToby Isaac 31813f27d899SToby Isaac Level: advanced 31823f27d899SToby Isaac 31833f27d899SToby Isaac .seealso: PetscDTNodeType, PetscDualSpaceLagrangeGetNodeType() 31843f27d899SToby Isaac @*/ 31853f27d899SToby Isaac PetscErrorCode PetscDualSpaceLagrangeSetNodeType(PetscDualSpace sp, PetscDTNodeType nodeType, PetscBool boundary, PetscReal exponent) 31863f27d899SToby Isaac { 31873f27d899SToby Isaac PetscFunctionBegin; 31883f27d899SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 3189*cac4c232SBarry Smith PetscTryMethod(sp,"PetscDualSpaceLagrangeSetNodeType_C",(PetscDualSpace,PetscDTNodeType,PetscBool,PetscReal),(sp,nodeType,boundary,exponent)); 31903f27d899SToby Isaac PetscFunctionReturn(0); 31913f27d899SToby Isaac } 31923f27d899SToby Isaac 319366a6c23cSMatthew G. Knepley /*@ 319466a6c23cSMatthew G. Knepley PetscDualSpaceLagrangeGetUseMoments - Get the flag for using moment functionals 319566a6c23cSMatthew G. Knepley 319666a6c23cSMatthew G. Knepley Not collective 319766a6c23cSMatthew G. Knepley 319866a6c23cSMatthew G. Knepley Input Parameter: 319966a6c23cSMatthew G. Knepley . sp - The PetscDualSpace 320066a6c23cSMatthew G. Knepley 320166a6c23cSMatthew G. Knepley Output Parameter: 320266a6c23cSMatthew G. Knepley . useMoments - Moment flag 320366a6c23cSMatthew G. Knepley 320466a6c23cSMatthew G. Knepley Level: advanced 320566a6c23cSMatthew G. Knepley 320666a6c23cSMatthew G. Knepley .seealso: PetscDualSpaceLagrangeSetUseMoments() 320766a6c23cSMatthew G. Knepley @*/ 320866a6c23cSMatthew G. Knepley PetscErrorCode PetscDualSpaceLagrangeGetUseMoments(PetscDualSpace sp, PetscBool *useMoments) 320966a6c23cSMatthew G. Knepley { 321066a6c23cSMatthew G. Knepley PetscFunctionBegin; 321166a6c23cSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 321266a6c23cSMatthew G. Knepley PetscValidBoolPointer(useMoments, 2); 3213*cac4c232SBarry Smith PetscUseMethod(sp,"PetscDualSpaceLagrangeGetUseMoments_C",(PetscDualSpace,PetscBool *),(sp,useMoments)); 321466a6c23cSMatthew G. Knepley PetscFunctionReturn(0); 321566a6c23cSMatthew G. Knepley } 321666a6c23cSMatthew G. Knepley 321766a6c23cSMatthew G. Knepley /*@ 321866a6c23cSMatthew G. Knepley PetscDualSpaceLagrangeSetUseMoments - Set the flag for moment functionals 321966a6c23cSMatthew G. Knepley 322066a6c23cSMatthew G. Knepley Logically collective 322166a6c23cSMatthew G. Knepley 322266a6c23cSMatthew G. Knepley Input Parameters: 322366a6c23cSMatthew G. Knepley + sp - The PetscDualSpace 322466a6c23cSMatthew G. Knepley - useMoments - The flag for moment functionals 322566a6c23cSMatthew G. Knepley 322666a6c23cSMatthew G. Knepley Level: advanced 322766a6c23cSMatthew G. Knepley 322866a6c23cSMatthew G. Knepley .seealso: PetscDualSpaceLagrangeGetUseMoments() 322966a6c23cSMatthew G. Knepley @*/ 323066a6c23cSMatthew G. Knepley PetscErrorCode PetscDualSpaceLagrangeSetUseMoments(PetscDualSpace sp, PetscBool useMoments) 323166a6c23cSMatthew G. Knepley { 323266a6c23cSMatthew G. Knepley PetscFunctionBegin; 323366a6c23cSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 3234*cac4c232SBarry Smith PetscTryMethod(sp,"PetscDualSpaceLagrangeSetUseMoments_C",(PetscDualSpace,PetscBool),(sp,useMoments)); 323566a6c23cSMatthew G. Knepley PetscFunctionReturn(0); 323666a6c23cSMatthew G. Knepley } 323766a6c23cSMatthew G. Knepley 323866a6c23cSMatthew G. Knepley /*@ 323966a6c23cSMatthew G. Knepley PetscDualSpaceLagrangeGetMomentOrder - Get the order for moment integration 324066a6c23cSMatthew G. Knepley 324166a6c23cSMatthew G. Knepley Not collective 324266a6c23cSMatthew G. Knepley 324366a6c23cSMatthew G. Knepley Input Parameter: 324466a6c23cSMatthew G. Knepley . sp - The PetscDualSpace 324566a6c23cSMatthew G. Knepley 324666a6c23cSMatthew G. Knepley Output Parameter: 324766a6c23cSMatthew G. Knepley . order - Moment integration order 324866a6c23cSMatthew G. Knepley 324966a6c23cSMatthew G. Knepley Level: advanced 325066a6c23cSMatthew G. Knepley 325166a6c23cSMatthew G. Knepley .seealso: PetscDualSpaceLagrangeSetMomentOrder() 325266a6c23cSMatthew G. Knepley @*/ 325366a6c23cSMatthew G. Knepley PetscErrorCode PetscDualSpaceLagrangeGetMomentOrder(PetscDualSpace sp, PetscInt *order) 325466a6c23cSMatthew G. Knepley { 325566a6c23cSMatthew G. Knepley PetscFunctionBegin; 325666a6c23cSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 325766a6c23cSMatthew G. Knepley PetscValidIntPointer(order, 2); 3258*cac4c232SBarry Smith PetscUseMethod(sp,"PetscDualSpaceLagrangeGetMomentOrder_C",(PetscDualSpace,PetscInt *),(sp,order)); 325966a6c23cSMatthew G. Knepley PetscFunctionReturn(0); 326066a6c23cSMatthew G. Knepley } 326166a6c23cSMatthew G. Knepley 326266a6c23cSMatthew G. Knepley /*@ 326366a6c23cSMatthew G. Knepley PetscDualSpaceLagrangeSetMomentOrder - Set the order for moment integration 326466a6c23cSMatthew G. Knepley 326566a6c23cSMatthew G. Knepley Logically collective 326666a6c23cSMatthew G. Knepley 326766a6c23cSMatthew G. Knepley Input Parameters: 326866a6c23cSMatthew G. Knepley + sp - The PetscDualSpace 326966a6c23cSMatthew G. Knepley - order - The order for moment integration 327066a6c23cSMatthew G. Knepley 327166a6c23cSMatthew G. Knepley Level: advanced 327266a6c23cSMatthew G. Knepley 327366a6c23cSMatthew G. Knepley .seealso: PetscDualSpaceLagrangeGetMomentOrder() 327466a6c23cSMatthew G. Knepley @*/ 327566a6c23cSMatthew G. Knepley PetscErrorCode PetscDualSpaceLagrangeSetMomentOrder(PetscDualSpace sp, PetscInt order) 327666a6c23cSMatthew G. Knepley { 327766a6c23cSMatthew G. Knepley PetscFunctionBegin; 327866a6c23cSMatthew G. Knepley PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 3279*cac4c232SBarry Smith PetscTryMethod(sp,"PetscDualSpaceLagrangeSetMomentOrder_C",(PetscDualSpace,PetscInt),(sp,order)); 328066a6c23cSMatthew G. Knepley PetscFunctionReturn(0); 328166a6c23cSMatthew G. Knepley } 32823f27d899SToby Isaac 32836f905325SMatthew G. Knepley static PetscErrorCode PetscDualSpaceInitialize_Lagrange(PetscDualSpace sp) 328420cf1dd8SToby Isaac { 328520cf1dd8SToby Isaac PetscFunctionBegin; 328620cf1dd8SToby Isaac sp->ops->destroy = PetscDualSpaceDestroy_Lagrange; 32876f905325SMatthew G. Knepley sp->ops->view = PetscDualSpaceView_Lagrange; 32886f905325SMatthew G. Knepley sp->ops->setfromoptions = PetscDualSpaceSetFromOptions_Lagrange; 328920cf1dd8SToby Isaac sp->ops->duplicate = PetscDualSpaceDuplicate_Lagrange; 32906f905325SMatthew G. Knepley sp->ops->setup = PetscDualSpaceSetUp_Lagrange; 32913f27d899SToby Isaac sp->ops->createheightsubspace = NULL; 32923f27d899SToby Isaac sp->ops->createpointsubspace = NULL; 329320cf1dd8SToby Isaac sp->ops->getsymmetries = PetscDualSpaceGetSymmetries_Lagrange; 329420cf1dd8SToby Isaac sp->ops->apply = PetscDualSpaceApplyDefault; 329520cf1dd8SToby Isaac sp->ops->applyall = PetscDualSpaceApplyAllDefault; 3296b4457527SToby Isaac sp->ops->applyint = PetscDualSpaceApplyInteriorDefault; 32973f27d899SToby Isaac sp->ops->createalldata = PetscDualSpaceCreateAllDataDefault; 3298b4457527SToby Isaac sp->ops->createintdata = PetscDualSpaceCreateInteriorDataDefault; 329920cf1dd8SToby Isaac PetscFunctionReturn(0); 330020cf1dd8SToby Isaac } 330120cf1dd8SToby Isaac 330220cf1dd8SToby Isaac /*MC 330320cf1dd8SToby Isaac PETSCDUALSPACELAGRANGE = "lagrange" - A PetscDualSpace object that encapsulates a dual space of pointwise evaluation functionals 330420cf1dd8SToby Isaac 330520cf1dd8SToby Isaac Level: intermediate 330620cf1dd8SToby Isaac 330720cf1dd8SToby Isaac .seealso: PetscDualSpaceType, PetscDualSpaceCreate(), PetscDualSpaceSetType() 330820cf1dd8SToby Isaac M*/ 330920cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceCreate_Lagrange(PetscDualSpace sp) 331020cf1dd8SToby Isaac { 331120cf1dd8SToby Isaac PetscDualSpace_Lag *lag; 331220cf1dd8SToby Isaac 331320cf1dd8SToby Isaac PetscFunctionBegin; 331420cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 33159566063dSJacob Faibussowitsch PetscCall(PetscNewLog(sp,&lag)); 331620cf1dd8SToby Isaac sp->data = lag; 331720cf1dd8SToby Isaac 33183f27d899SToby Isaac lag->tensorCell = PETSC_FALSE; 331920cf1dd8SToby Isaac lag->tensorSpace = PETSC_FALSE; 332020cf1dd8SToby Isaac lag->continuous = PETSC_TRUE; 33213f27d899SToby Isaac lag->numCopies = PETSC_DEFAULT; 33223f27d899SToby Isaac lag->numNodeSkip = PETSC_DEFAULT; 33233f27d899SToby Isaac lag->nodeType = PETSCDTNODES_DEFAULT; 332466a6c23cSMatthew G. Knepley lag->useMoments = PETSC_FALSE; 332566a6c23cSMatthew G. Knepley lag->momentOrder = 0; 332620cf1dd8SToby Isaac 33279566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceInitialize_Lagrange(sp)); 33289566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetContinuity_C", PetscDualSpaceLagrangeGetContinuity_Lagrange)); 33299566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetContinuity_C", PetscDualSpaceLagrangeSetContinuity_Lagrange)); 33309566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetTensor_C", PetscDualSpaceLagrangeGetTensor_Lagrange)); 33319566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetTensor_C", PetscDualSpaceLagrangeSetTensor_Lagrange)); 33329566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetTrimmed_C", PetscDualSpaceLagrangeGetTrimmed_Lagrange)); 33339566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetTrimmed_C", PetscDualSpaceLagrangeSetTrimmed_Lagrange)); 33349566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetNodeType_C", PetscDualSpaceLagrangeGetNodeType_Lagrange)); 33359566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetNodeType_C", PetscDualSpaceLagrangeSetNodeType_Lagrange)); 33369566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetUseMoments_C", PetscDualSpaceLagrangeGetUseMoments_Lagrange)); 33379566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetUseMoments_C", PetscDualSpaceLagrangeSetUseMoments_Lagrange)); 33389566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetMomentOrder_C", PetscDualSpaceLagrangeGetMomentOrder_Lagrange)); 33399566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetMomentOrder_C", PetscDualSpaceLagrangeSetMomentOrder_Lagrange)); 334020cf1dd8SToby Isaac PetscFunctionReturn(0); 334120cf1dd8SToby Isaac } 3342