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