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; 964 PetscDTNodeType nodeType; 965 PetscReal exponent; 966 PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data; 967 968 PetscFunctionBegin; 969 PetscCall(PetscDualSpaceLagrangeGetContinuity(sp, &cont)); 970 PetscCall(PetscDualSpaceLagrangeSetContinuity(spNew, cont)); 971 PetscCall(PetscDualSpaceLagrangeGetTensor(sp, &tensor)); 972 PetscCall(PetscDualSpaceLagrangeSetTensor(spNew, tensor)); 973 PetscCall(PetscDualSpaceLagrangeGetTrimmed(sp, &trimmed)); 974 PetscCall(PetscDualSpaceLagrangeSetTrimmed(spNew, trimmed)); 975 PetscCall(PetscDualSpaceLagrangeGetNodeType(sp, &nodeType, &boundary, &exponent)); 976 PetscCall(PetscDualSpaceLagrangeSetNodeType(spNew, nodeType, boundary, exponent)); 977 if (lag->nodeFamily) { 978 PetscDualSpace_Lag *lagnew = (PetscDualSpace_Lag *)spNew->data; 979 980 PetscCall(Petsc1DNodeFamilyReference(lag->nodeFamily)); 981 lagnew->nodeFamily = lag->nodeFamily; 982 } 983 PetscFunctionReturn(PETSC_SUCCESS); 984 } 985 986 /* for making tensor product spaces: take a dual space and product a segment space that has all the same 987 * specifications (trimmed, continuous, order, node set), except for the form degree */ 988 static PetscErrorCode PetscDualSpaceCreateEdgeSubspace_Lagrange(PetscDualSpace sp, PetscInt order, PetscInt k, PetscInt Nc, PetscBool interiorOnly, PetscDualSpace *bdsp) 989 { 990 DM K; 991 PetscDualSpace_Lag *newlag; 992 993 PetscFunctionBegin; 994 PetscCall(PetscDualSpaceDuplicate(sp, bdsp)); 995 PetscCall(PetscDualSpaceSetFormDegree(*bdsp, k)); 996 PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, DMPolytopeTypeSimpleShape(1, PETSC_TRUE), &K)); 997 PetscCall(PetscDualSpaceSetDM(*bdsp, K)); 998 PetscCall(DMDestroy(&K)); 999 PetscCall(PetscDualSpaceSetOrder(*bdsp, order)); 1000 PetscCall(PetscDualSpaceSetNumComponents(*bdsp, Nc)); 1001 newlag = (PetscDualSpace_Lag *)(*bdsp)->data; 1002 newlag->interiorOnly = interiorOnly; 1003 PetscCall(PetscDualSpaceSetUp(*bdsp)); 1004 PetscFunctionReturn(PETSC_SUCCESS); 1005 } 1006 1007 /* just the points, weights aren't handled */ 1008 static PetscErrorCode PetscQuadratureCreateTensor(PetscQuadrature trace, PetscQuadrature fiber, PetscQuadrature *product) 1009 { 1010 PetscInt dimTrace, dimFiber; 1011 PetscInt numPointsTrace, numPointsFiber; 1012 PetscInt dim, numPoints; 1013 const PetscReal *pointsTrace; 1014 const PetscReal *pointsFiber; 1015 PetscReal *points; 1016 PetscInt i, j, k, p; 1017 1018 PetscFunctionBegin; 1019 PetscCall(PetscQuadratureGetData(trace, &dimTrace, NULL, &numPointsTrace, &pointsTrace, NULL)); 1020 PetscCall(PetscQuadratureGetData(fiber, &dimFiber, NULL, &numPointsFiber, &pointsFiber, NULL)); 1021 dim = dimTrace + dimFiber; 1022 numPoints = numPointsFiber * numPointsTrace; 1023 PetscCall(PetscMalloc1(numPoints * dim, &points)); 1024 for (p = 0, j = 0; j < numPointsFiber; j++) { 1025 for (i = 0; i < numPointsTrace; i++, p++) { 1026 for (k = 0; k < dimTrace; k++) points[p * dim + k] = pointsTrace[i * dimTrace + k]; 1027 for (k = 0; k < dimFiber; k++) points[p * dim + dimTrace + k] = pointsFiber[j * dimFiber + k]; 1028 } 1029 } 1030 PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, product)); 1031 PetscCall(PetscQuadratureSetData(*product, dim, 0, numPoints, points, NULL)); 1032 PetscFunctionReturn(PETSC_SUCCESS); 1033 } 1034 1035 /* Kronecker tensor product where matrix is considered a matrix of k-forms, so that 1036 * the entries in the product matrix are wedge products of the entries in the original matrices */ 1037 static PetscErrorCode MatTensorAltV(Mat trace, Mat fiber, PetscInt dimTrace, PetscInt kTrace, PetscInt dimFiber, PetscInt kFiber, Mat *product) 1038 { 1039 PetscInt mTrace, nTrace, mFiber, nFiber, m, n, k, i, j, l; 1040 PetscInt dim, NkTrace, NkFiber, Nk; 1041 PetscInt dT, dF; 1042 PetscInt *nnzTrace, *nnzFiber, *nnz; 1043 PetscInt iT, iF, jT, jF, il, jl; 1044 PetscReal *workT, *workT2, *workF, *workF2, *work, *workstar; 1045 PetscReal *projT, *projF; 1046 PetscReal *projTstar, *projFstar; 1047 PetscReal *wedgeMat; 1048 PetscReal sign; 1049 PetscScalar *workS; 1050 Mat prod; 1051 /* this produces dof groups that look like the identity */ 1052 1053 PetscFunctionBegin; 1054 PetscCall(MatGetSize(trace, &mTrace, &nTrace)); 1055 PetscCall(PetscDTBinomialInt(dimTrace, PetscAbsInt(kTrace), &NkTrace)); 1056 PetscCheck(nTrace % NkTrace == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "point value space of trace matrix is not a multiple of k-form size"); 1057 PetscCall(MatGetSize(fiber, &mFiber, &nFiber)); 1058 PetscCall(PetscDTBinomialInt(dimFiber, PetscAbsInt(kFiber), &NkFiber)); 1059 PetscCheck(nFiber % NkFiber == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "point value space of fiber matrix is not a multiple of k-form size"); 1060 PetscCall(PetscMalloc2(mTrace, &nnzTrace, mFiber, &nnzFiber)); 1061 for (i = 0; i < mTrace; i++) { 1062 PetscCall(MatGetRow(trace, i, &nnzTrace[i], NULL, NULL)); 1063 PetscCheck(nnzTrace[i] % NkTrace == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "nonzeros in trace matrix are not in k-form size blocks"); 1064 } 1065 for (i = 0; i < mFiber; i++) { 1066 PetscCall(MatGetRow(fiber, i, &nnzFiber[i], NULL, NULL)); 1067 PetscCheck(nnzFiber[i] % NkFiber == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "nonzeros in fiber matrix are not in k-form size blocks"); 1068 } 1069 dim = dimTrace + dimFiber; 1070 k = kFiber + kTrace; 1071 PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(k), &Nk)); 1072 m = mTrace * mFiber; 1073 PetscCall(PetscMalloc1(m, &nnz)); 1074 for (l = 0, j = 0; j < mFiber; j++) 1075 for (i = 0; i < mTrace; i++, l++) nnz[l] = (nnzTrace[i] / NkTrace) * (nnzFiber[j] / NkFiber) * Nk; 1076 n = (nTrace / NkTrace) * (nFiber / NkFiber) * Nk; 1077 PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, m, n, 0, nnz, &prod)); 1078 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)prod, "altv_")); 1079 PetscCall(PetscFree(nnz)); 1080 PetscCall(PetscFree2(nnzTrace, nnzFiber)); 1081 /* reasoning about which points each dof needs depends on having zeros computed at points preserved */ 1082 PetscCall(MatSetOption(prod, MAT_IGNORE_ZERO_ENTRIES, PETSC_FALSE)); 1083 /* compute pullbacks */ 1084 PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(kTrace), &dT)); 1085 PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(kFiber), &dF)); 1086 PetscCall(PetscMalloc4(dimTrace * dim, &projT, dimFiber * dim, &projF, dT * NkTrace, &projTstar, dF * NkFiber, &projFstar)); 1087 PetscCall(PetscArrayzero(projT, dimTrace * dim)); 1088 for (i = 0; i < dimTrace; i++) projT[i * (dim + 1)] = 1.; 1089 PetscCall(PetscArrayzero(projF, dimFiber * dim)); 1090 for (i = 0; i < dimFiber; i++) projF[i * (dim + 1) + dimTrace] = 1.; 1091 PetscCall(PetscDTAltVPullbackMatrix(dim, dimTrace, projT, kTrace, projTstar)); 1092 PetscCall(PetscDTAltVPullbackMatrix(dim, dimFiber, projF, kFiber, projFstar)); 1093 PetscCall(PetscMalloc5(dT, &workT, dF, &workF, Nk, &work, Nk, &workstar, Nk, &workS)); 1094 PetscCall(PetscMalloc2(dT, &workT2, dF, &workF2)); 1095 PetscCall(PetscMalloc1(Nk * dT, &wedgeMat)); 1096 sign = (PetscAbsInt(kTrace * kFiber) & 1) ? -1. : 1.; 1097 for (i = 0, iF = 0; iF < mFiber; iF++) { 1098 PetscInt ncolsF, nformsF; 1099 const PetscInt *colsF; 1100 const PetscScalar *valsF; 1101 1102 PetscCall(MatGetRow(fiber, iF, &ncolsF, &colsF, &valsF)); 1103 nformsF = ncolsF / NkFiber; 1104 for (iT = 0; iT < mTrace; iT++, i++) { 1105 PetscInt ncolsT, nformsT; 1106 const PetscInt *colsT; 1107 const PetscScalar *valsT; 1108 1109 PetscCall(MatGetRow(trace, iT, &ncolsT, &colsT, &valsT)); 1110 nformsT = ncolsT / NkTrace; 1111 for (j = 0, jF = 0; jF < nformsF; jF++) { 1112 PetscInt colF = colsF[jF * NkFiber] / NkFiber; 1113 1114 for (il = 0; il < dF; il++) { 1115 PetscReal val = 0.; 1116 for (jl = 0; jl < NkFiber; jl++) val += projFstar[il * NkFiber + jl] * PetscRealPart(valsF[jF * NkFiber + jl]); 1117 workF[il] = val; 1118 } 1119 if (kFiber < 0) { 1120 for (il = 0; il < dF; il++) workF2[il] = workF[il]; 1121 PetscCall(PetscDTAltVStar(dim, PetscAbsInt(kFiber), 1, workF2, workF)); 1122 } 1123 PetscCall(PetscDTAltVWedgeMatrix(dim, PetscAbsInt(kFiber), PetscAbsInt(kTrace), workF, wedgeMat)); 1124 for (jT = 0; jT < nformsT; jT++, j++) { 1125 PetscInt colT = colsT[jT * NkTrace] / NkTrace; 1126 PetscInt col = colF * (nTrace / NkTrace) + colT; 1127 const PetscScalar *vals; 1128 1129 for (il = 0; il < dT; il++) { 1130 PetscReal val = 0.; 1131 for (jl = 0; jl < NkTrace; jl++) val += projTstar[il * NkTrace + jl] * PetscRealPart(valsT[jT * NkTrace + jl]); 1132 workT[il] = val; 1133 } 1134 if (kTrace < 0) { 1135 for (il = 0; il < dT; il++) workT2[il] = workT[il]; 1136 PetscCall(PetscDTAltVStar(dim, PetscAbsInt(kTrace), 1, workT2, workT)); 1137 } 1138 1139 for (il = 0; il < Nk; il++) { 1140 PetscReal val = 0.; 1141 for (jl = 0; jl < dT; jl++) val += sign * wedgeMat[il * dT + jl] * workT[jl]; 1142 work[il] = val; 1143 } 1144 if (k < 0) { 1145 PetscCall(PetscDTAltVStar(dim, PetscAbsInt(k), -1, work, workstar)); 1146 #if defined(PETSC_USE_COMPLEX) 1147 for (l = 0; l < Nk; l++) workS[l] = workstar[l]; 1148 vals = &workS[0]; 1149 #else 1150 vals = &workstar[0]; 1151 #endif 1152 } else { 1153 #if defined(PETSC_USE_COMPLEX) 1154 for (l = 0; l < Nk; l++) workS[l] = work[l]; 1155 vals = &workS[0]; 1156 #else 1157 vals = &work[0]; 1158 #endif 1159 } 1160 for (l = 0; l < Nk; l++) PetscCall(MatSetValue(prod, i, col * Nk + l, vals[l], INSERT_VALUES)); /* Nk */ 1161 } /* jT */ 1162 } /* jF */ 1163 PetscCall(MatRestoreRow(trace, iT, &ncolsT, &colsT, &valsT)); 1164 } /* iT */ 1165 PetscCall(MatRestoreRow(fiber, iF, &ncolsF, &colsF, &valsF)); 1166 } /* iF */ 1167 PetscCall(PetscFree(wedgeMat)); 1168 PetscCall(PetscFree4(projT, projF, projTstar, projFstar)); 1169 PetscCall(PetscFree2(workT2, workF2)); 1170 PetscCall(PetscFree5(workT, workF, work, workstar, workS)); 1171 PetscCall(MatAssemblyBegin(prod, MAT_FINAL_ASSEMBLY)); 1172 PetscCall(MatAssemblyEnd(prod, MAT_FINAL_ASSEMBLY)); 1173 *product = prod; 1174 PetscFunctionReturn(PETSC_SUCCESS); 1175 } 1176 1177 /* Union of quadrature points, with an attempt to identify common points in the two sets */ 1178 static PetscErrorCode PetscQuadraturePointsMerge(PetscQuadrature quadA, PetscQuadrature quadB, PetscQuadrature *quadJoint, PetscInt *aToJoint[], PetscInt *bToJoint[]) 1179 { 1180 PetscInt dimA, dimB; 1181 PetscInt nA, nB, nJoint, i, j, d; 1182 const PetscReal *pointsA; 1183 const PetscReal *pointsB; 1184 PetscReal *pointsJoint; 1185 PetscInt *aToJ, *bToJ; 1186 PetscQuadrature qJ; 1187 1188 PetscFunctionBegin; 1189 PetscCall(PetscQuadratureGetData(quadA, &dimA, NULL, &nA, &pointsA, NULL)); 1190 PetscCall(PetscQuadratureGetData(quadB, &dimB, NULL, &nB, &pointsB, NULL)); 1191 PetscCheck(dimA == dimB, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Quadrature points must be in the same dimension"); 1192 nJoint = nA; 1193 PetscCall(PetscMalloc1(nA, &aToJ)); 1194 for (i = 0; i < nA; i++) aToJ[i] = i; 1195 PetscCall(PetscMalloc1(nB, &bToJ)); 1196 for (i = 0; i < nB; i++) { 1197 for (j = 0; j < nA; j++) { 1198 bToJ[i] = -1; 1199 for (d = 0; d < dimA; d++) 1200 if (PetscAbsReal(pointsB[i * dimA + d] - pointsA[j * dimA + d]) > PETSC_SMALL) break; 1201 if (d == dimA) { 1202 bToJ[i] = j; 1203 break; 1204 } 1205 } 1206 if (bToJ[i] == -1) bToJ[i] = nJoint++; 1207 } 1208 *aToJoint = aToJ; 1209 *bToJoint = bToJ; 1210 PetscCall(PetscMalloc1(nJoint * dimA, &pointsJoint)); 1211 PetscCall(PetscArraycpy(pointsJoint, pointsA, nA * dimA)); 1212 for (i = 0; i < nB; i++) { 1213 if (bToJ[i] >= nA) { 1214 for (d = 0; d < dimA; d++) pointsJoint[bToJ[i] * dimA + d] = pointsB[i * dimA + d]; 1215 } 1216 } 1217 PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, &qJ)); 1218 PetscCall(PetscQuadratureSetData(qJ, dimA, 0, nJoint, pointsJoint, NULL)); 1219 *quadJoint = qJ; 1220 PetscFunctionReturn(PETSC_SUCCESS); 1221 } 1222 1223 /* Matrices matA and matB are both quadrature -> dof matrices: produce a matrix that is joint quadrature -> union of 1224 * dofs, where the joint quadrature was produced by PetscQuadraturePointsMerge */ 1225 static PetscErrorCode MatricesMerge(Mat matA, Mat matB, PetscInt dim, PetscInt k, PetscInt numMerged, const PetscInt aToMerged[], const PetscInt bToMerged[], Mat *matMerged) 1226 { 1227 PetscInt m, n, mA, nA, mB, nB, Nk, i, j, l; 1228 Mat M; 1229 PetscInt *nnz; 1230 PetscInt maxnnz; 1231 PetscInt *work; 1232 1233 PetscFunctionBegin; 1234 PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(k), &Nk)); 1235 PetscCall(MatGetSize(matA, &mA, &nA)); 1236 PetscCheck(nA % Nk == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "matA column space not a multiple of k-form size"); 1237 PetscCall(MatGetSize(matB, &mB, &nB)); 1238 PetscCheck(nB % Nk == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "matB column space not a multiple of k-form size"); 1239 m = mA + mB; 1240 n = numMerged * Nk; 1241 PetscCall(PetscMalloc1(m, &nnz)); 1242 maxnnz = 0; 1243 for (i = 0; i < mA; i++) { 1244 PetscCall(MatGetRow(matA, i, &nnz[i], NULL, NULL)); 1245 PetscCheck(nnz[i] % Nk == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "nonzeros in matA are not in k-form size blocks"); 1246 maxnnz = PetscMax(maxnnz, nnz[i]); 1247 } 1248 for (i = 0; i < mB; i++) { 1249 PetscCall(MatGetRow(matB, i, &nnz[i + mA], NULL, NULL)); 1250 PetscCheck(nnz[i + mA] % Nk == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "nonzeros in matB are not in k-form size blocks"); 1251 maxnnz = PetscMax(maxnnz, nnz[i + mA]); 1252 } 1253 PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, m, n, 0, nnz, &M)); 1254 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)M, "altv_")); 1255 PetscCall(PetscFree(nnz)); 1256 /* reasoning about which points each dof needs depends on having zeros computed at points preserved */ 1257 PetscCall(MatSetOption(M, MAT_IGNORE_ZERO_ENTRIES, PETSC_FALSE)); 1258 PetscCall(PetscMalloc1(maxnnz, &work)); 1259 for (i = 0; i < mA; i++) { 1260 const PetscInt *cols; 1261 const PetscScalar *vals; 1262 PetscInt nCols; 1263 PetscCall(MatGetRow(matA, i, &nCols, &cols, &vals)); 1264 for (j = 0; j < nCols / Nk; j++) { 1265 PetscInt newCol = aToMerged[cols[j * Nk] / Nk]; 1266 for (l = 0; l < Nk; l++) work[j * Nk + l] = newCol * Nk + l; 1267 } 1268 PetscCall(MatSetValuesBlocked(M, 1, &i, nCols, work, vals, INSERT_VALUES)); 1269 PetscCall(MatRestoreRow(matA, i, &nCols, &cols, &vals)); 1270 } 1271 for (i = 0; i < mB; i++) { 1272 const PetscInt *cols; 1273 const PetscScalar *vals; 1274 1275 PetscInt row = i + mA; 1276 PetscInt nCols; 1277 PetscCall(MatGetRow(matB, i, &nCols, &cols, &vals)); 1278 for (j = 0; j < nCols / Nk; j++) { 1279 PetscInt newCol = bToMerged[cols[j * Nk] / Nk]; 1280 for (l = 0; l < Nk; l++) work[j * Nk + l] = newCol * Nk + l; 1281 } 1282 PetscCall(MatSetValuesBlocked(M, 1, &row, nCols, work, vals, INSERT_VALUES)); 1283 PetscCall(MatRestoreRow(matB, i, &nCols, &cols, &vals)); 1284 } 1285 PetscCall(PetscFree(work)); 1286 PetscCall(MatAssemblyBegin(M, MAT_FINAL_ASSEMBLY)); 1287 PetscCall(MatAssemblyEnd(M, MAT_FINAL_ASSEMBLY)); 1288 *matMerged = M; 1289 PetscFunctionReturn(PETSC_SUCCESS); 1290 } 1291 1292 /* Take a dual space and product a segment space that has all the same specifications (trimmed, continuous, order, 1293 * node set), except for the form degree. For computing boundary dofs and for making tensor product spaces */ 1294 static PetscErrorCode PetscDualSpaceCreateFacetSubspace_Lagrange(PetscDualSpace sp, DM K, PetscInt f, PetscInt k, PetscInt Ncopies, PetscBool interiorOnly, PetscDualSpace *bdsp) 1295 { 1296 PetscInt Nknew, Ncnew; 1297 PetscInt dim, pointDim = -1; 1298 PetscInt depth; 1299 DM dm; 1300 PetscDualSpace_Lag *newlag; 1301 1302 PetscFunctionBegin; 1303 PetscCall(PetscDualSpaceGetDM(sp, &dm)); 1304 PetscCall(DMGetDimension(dm, &dim)); 1305 PetscCall(DMPlexGetDepth(dm, &depth)); 1306 PetscCall(PetscDualSpaceDuplicate(sp, bdsp)); 1307 PetscCall(PetscDualSpaceSetFormDegree(*bdsp, k)); 1308 if (!K) { 1309 if (depth == dim) { 1310 DMPolytopeType ct; 1311 1312 pointDim = dim - 1; 1313 PetscCall(DMPlexGetCellType(dm, f, &ct)); 1314 PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, ct, &K)); 1315 } else if (depth == 1) { 1316 pointDim = 0; 1317 PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, DM_POLYTOPE_POINT, &K)); 1318 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unsupported interpolation state of reference element"); 1319 } else { 1320 PetscCall(PetscObjectReference((PetscObject)K)); 1321 PetscCall(DMGetDimension(K, &pointDim)); 1322 } 1323 PetscCall(PetscDualSpaceSetDM(*bdsp, K)); 1324 PetscCall(DMDestroy(&K)); 1325 PetscCall(PetscDTBinomialInt(pointDim, PetscAbsInt(k), &Nknew)); 1326 Ncnew = Nknew * Ncopies; 1327 PetscCall(PetscDualSpaceSetNumComponents(*bdsp, Ncnew)); 1328 newlag = (PetscDualSpace_Lag *)(*bdsp)->data; 1329 newlag->interiorOnly = interiorOnly; 1330 PetscCall(PetscDualSpaceSetUp(*bdsp)); 1331 PetscFunctionReturn(PETSC_SUCCESS); 1332 } 1333 1334 /* Construct simplex nodes from a nodefamily, add Nk dof vectors of length Nk at each node. 1335 * Return the (quadrature, matrix) form of the dofs and the nodeIndices form as well. 1336 * 1337 * Sometimes we want a set of nodes to be contained in the interior of the element, 1338 * even when the node scheme puts nodes on the boundaries. numNodeSkip tells 1339 * the routine how many "layers" of nodes need to be skipped. 1340 * */ 1341 static PetscErrorCode PetscDualSpaceLagrangeCreateSimplexNodeMat(Petsc1DNodeFamily nodeFamily, PetscInt dim, PetscInt sum, PetscInt Nk, PetscInt numNodeSkip, PetscQuadrature *iNodes, Mat *iMat, PetscLagNodeIndices *nodeIndices) 1342 { 1343 PetscReal *extraNodeCoords, *nodeCoords; 1344 PetscInt nNodes, nExtraNodes; 1345 PetscInt i, j, k, extraSum = sum + numNodeSkip * (1 + dim); 1346 PetscQuadrature intNodes; 1347 Mat intMat; 1348 PetscLagNodeIndices ni; 1349 1350 PetscFunctionBegin; 1351 PetscCall(PetscDTBinomialInt(dim + sum, dim, &nNodes)); 1352 PetscCall(PetscDTBinomialInt(dim + extraSum, dim, &nExtraNodes)); 1353 1354 PetscCall(PetscMalloc1(dim * nExtraNodes, &extraNodeCoords)); 1355 PetscCall(PetscNew(&ni)); 1356 ni->nodeIdxDim = dim + 1; 1357 ni->nodeVecDim = Nk; 1358 ni->nNodes = nNodes * Nk; 1359 ni->refct = 1; 1360 PetscCall(PetscMalloc1(nNodes * Nk * (dim + 1), &ni->nodeIdx)); 1361 PetscCall(PetscMalloc1(nNodes * Nk * Nk, &ni->nodeVec)); 1362 for (i = 0; i < nNodes; i++) 1363 for (j = 0; j < Nk; j++) 1364 for (k = 0; k < Nk; k++) ni->nodeVec[(i * Nk + j) * Nk + k] = (j == k) ? 1. : 0.; 1365 PetscCall(Petsc1DNodeFamilyComputeSimplexNodes(nodeFamily, dim, extraSum, extraNodeCoords)); 1366 if (numNodeSkip) { 1367 PetscInt k; 1368 PetscInt *tup; 1369 1370 PetscCall(PetscMalloc1(dim * nNodes, &nodeCoords)); 1371 PetscCall(PetscMalloc1(dim + 1, &tup)); 1372 for (k = 0; k < nNodes; k++) { 1373 PetscInt j, c; 1374 PetscInt index; 1375 1376 PetscCall(PetscDTIndexToBary(dim + 1, sum, k, tup)); 1377 for (j = 0; j < dim + 1; j++) tup[j] += numNodeSkip; 1378 for (c = 0; c < Nk; c++) { 1379 for (j = 0; j < dim + 1; j++) ni->nodeIdx[(k * Nk + c) * (dim + 1) + j] = tup[j] + 1; 1380 } 1381 PetscCall(PetscDTBaryToIndex(dim + 1, extraSum, tup, &index)); 1382 for (j = 0; j < dim; j++) nodeCoords[k * dim + j] = extraNodeCoords[index * dim + j]; 1383 } 1384 PetscCall(PetscFree(tup)); 1385 PetscCall(PetscFree(extraNodeCoords)); 1386 } else { 1387 PetscInt k; 1388 PetscInt *tup; 1389 1390 nodeCoords = extraNodeCoords; 1391 PetscCall(PetscMalloc1(dim + 1, &tup)); 1392 for (k = 0; k < nNodes; k++) { 1393 PetscInt j, c; 1394 1395 PetscCall(PetscDTIndexToBary(dim + 1, sum, k, tup)); 1396 for (c = 0; c < Nk; c++) { 1397 for (j = 0; j < dim + 1; j++) { 1398 /* barycentric indices can have zeros, but we don't want to push forward zeros because it makes it harder to 1399 * determine which nodes correspond to which under symmetries, so we increase by 1. This is fine 1400 * because the nodeIdx coordinates don't have any meaning other than helping to identify symmetries */ 1401 ni->nodeIdx[(k * Nk + c) * (dim + 1) + j] = tup[j] + 1; 1402 } 1403 } 1404 } 1405 PetscCall(PetscFree(tup)); 1406 } 1407 PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, &intNodes)); 1408 PetscCall(PetscQuadratureSetData(intNodes, dim, 0, nNodes, nodeCoords, NULL)); 1409 PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, nNodes * Nk, nNodes * Nk, Nk, NULL, &intMat)); 1410 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)intMat, "lag_")); 1411 PetscCall(MatSetOption(intMat, MAT_IGNORE_ZERO_ENTRIES, PETSC_FALSE)); 1412 for (j = 0; j < nNodes * Nk; j++) { 1413 PetscInt rem = j % Nk; 1414 PetscInt a, aprev = j - rem; 1415 PetscInt anext = aprev + Nk; 1416 1417 for (a = aprev; a < anext; a++) PetscCall(MatSetValue(intMat, j, a, (a == j) ? 1. : 0., INSERT_VALUES)); 1418 } 1419 PetscCall(MatAssemblyBegin(intMat, MAT_FINAL_ASSEMBLY)); 1420 PetscCall(MatAssemblyEnd(intMat, MAT_FINAL_ASSEMBLY)); 1421 *iNodes = intNodes; 1422 *iMat = intMat; 1423 *nodeIndices = ni; 1424 PetscFunctionReturn(PETSC_SUCCESS); 1425 } 1426 1427 /* once the nodeIndices have been created for the interior of the reference cell, and for all of the boundary cells, 1428 * push forward the boundary dofs and concatenate them into the full node indices for the dual space */ 1429 static PetscErrorCode PetscDualSpaceLagrangeCreateAllNodeIdx(PetscDualSpace sp) 1430 { 1431 DM dm; 1432 PetscInt dim, nDofs; 1433 PetscSection section; 1434 PetscInt pStart, pEnd, p; 1435 PetscInt formDegree, Nk; 1436 PetscInt nodeIdxDim, spintdim; 1437 PetscDualSpace_Lag *lag; 1438 PetscLagNodeIndices ni, verti; 1439 1440 PetscFunctionBegin; 1441 lag = (PetscDualSpace_Lag *)sp->data; 1442 verti = lag->vertIndices; 1443 PetscCall(PetscDualSpaceGetDM(sp, &dm)); 1444 PetscCall(DMGetDimension(dm, &dim)); 1445 PetscCall(PetscDualSpaceGetFormDegree(sp, &formDegree)); 1446 PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(formDegree), &Nk)); 1447 PetscCall(PetscDualSpaceGetSection(sp, §ion)); 1448 PetscCall(PetscSectionGetStorageSize(section, &nDofs)); 1449 PetscCall(PetscNew(&ni)); 1450 ni->nodeIdxDim = nodeIdxDim = verti->nodeIdxDim; 1451 ni->nodeVecDim = Nk; 1452 ni->nNodes = nDofs; 1453 ni->refct = 1; 1454 PetscCall(PetscMalloc1(nodeIdxDim * nDofs, &ni->nodeIdx)); 1455 PetscCall(PetscMalloc1(Nk * nDofs, &ni->nodeVec)); 1456 PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 1457 PetscCall(PetscSectionGetDof(section, 0, &spintdim)); 1458 if (spintdim) { 1459 PetscCall(PetscArraycpy(ni->nodeIdx, lag->intNodeIndices->nodeIdx, spintdim * nodeIdxDim)); 1460 PetscCall(PetscArraycpy(ni->nodeVec, lag->intNodeIndices->nodeVec, spintdim * Nk)); 1461 } 1462 for (p = pStart + 1; p < pEnd; p++) { 1463 PetscDualSpace psp = sp->pointSpaces[p]; 1464 PetscDualSpace_Lag *plag; 1465 PetscInt dof, off; 1466 1467 PetscCall(PetscSectionGetDof(section, p, &dof)); 1468 if (!dof) continue; 1469 plag = (PetscDualSpace_Lag *)psp->data; 1470 PetscCall(PetscSectionGetOffset(section, p, &off)); 1471 PetscCall(PetscLagNodeIndicesPushForward(dm, verti, p, plag->vertIndices, plag->intNodeIndices, 0, formDegree, &ni->nodeIdx[off * nodeIdxDim], &ni->nodeVec[off * Nk])); 1472 } 1473 lag->allNodeIndices = ni; 1474 PetscFunctionReturn(PETSC_SUCCESS); 1475 } 1476 1477 /* once the (quadrature, Matrix) forms of the dofs have been created for the interior of the 1478 * reference cell and for the boundary cells, jk 1479 * push forward the boundary data and concatenate them into the full (quadrature, matrix) data 1480 * for the dual space */ 1481 static PetscErrorCode PetscDualSpaceCreateAllDataFromInteriorData(PetscDualSpace sp) 1482 { 1483 DM dm; 1484 PetscSection section; 1485 PetscInt pStart, pEnd, p, k, Nk, dim, Nc; 1486 PetscInt nNodes; 1487 PetscInt countNodes; 1488 Mat allMat; 1489 PetscQuadrature allNodes; 1490 PetscInt nDofs; 1491 PetscInt maxNzforms, j; 1492 PetscScalar *work; 1493 PetscReal *L, *J, *Jinv, *v0, *pv0; 1494 PetscInt *iwork; 1495 PetscReal *nodes; 1496 1497 PetscFunctionBegin; 1498 PetscCall(PetscDualSpaceGetDM(sp, &dm)); 1499 PetscCall(DMGetDimension(dm, &dim)); 1500 PetscCall(PetscDualSpaceGetSection(sp, §ion)); 1501 PetscCall(PetscSectionGetStorageSize(section, &nDofs)); 1502 PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 1503 PetscCall(PetscDualSpaceGetFormDegree(sp, &k)); 1504 PetscCall(PetscDualSpaceGetNumComponents(sp, &Nc)); 1505 PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(k), &Nk)); 1506 for (p = pStart, nNodes = 0, maxNzforms = 0; p < pEnd; p++) { 1507 PetscDualSpace psp; 1508 DM pdm; 1509 PetscInt pdim, pNk; 1510 PetscQuadrature intNodes; 1511 Mat intMat; 1512 1513 PetscCall(PetscDualSpaceGetPointSubspace(sp, p, &psp)); 1514 if (!psp) continue; 1515 PetscCall(PetscDualSpaceGetDM(psp, &pdm)); 1516 PetscCall(DMGetDimension(pdm, &pdim)); 1517 if (pdim < PetscAbsInt(k)) continue; 1518 PetscCall(PetscDTBinomialInt(pdim, PetscAbsInt(k), &pNk)); 1519 PetscCall(PetscDualSpaceGetInteriorData(psp, &intNodes, &intMat)); 1520 if (intNodes) { 1521 PetscInt nNodesp; 1522 1523 PetscCall(PetscQuadratureGetData(intNodes, NULL, NULL, &nNodesp, NULL, NULL)); 1524 nNodes += nNodesp; 1525 } 1526 if (intMat) { 1527 PetscInt maxNzsp; 1528 PetscInt maxNzformsp; 1529 1530 PetscCall(MatSeqAIJGetMaxRowNonzeros(intMat, &maxNzsp)); 1531 PetscCheck(maxNzsp % pNk == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "interior matrix is not laid out as blocks of k-forms"); 1532 maxNzformsp = maxNzsp / pNk; 1533 maxNzforms = PetscMax(maxNzforms, maxNzformsp); 1534 } 1535 } 1536 PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, nDofs, nNodes * Nc, maxNzforms * Nk, NULL, &allMat)); 1537 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)allMat, "ds_")); 1538 PetscCall(MatSetOption(allMat, MAT_IGNORE_ZERO_ENTRIES, PETSC_FALSE)); 1539 PetscCall(PetscMalloc7(dim, &v0, dim, &pv0, dim * dim, &J, dim * dim, &Jinv, Nk * Nk, &L, maxNzforms * Nk, &work, maxNzforms * Nk, &iwork)); 1540 for (j = 0; j < dim; j++) pv0[j] = -1.; 1541 PetscCall(PetscMalloc1(dim * nNodes, &nodes)); 1542 for (p = pStart, countNodes = 0; p < pEnd; p++) { 1543 PetscDualSpace psp; 1544 PetscQuadrature intNodes; 1545 DM pdm; 1546 PetscInt pdim, pNk; 1547 PetscInt countNodesIn = countNodes; 1548 PetscReal detJ; 1549 Mat intMat; 1550 1551 PetscCall(PetscDualSpaceGetPointSubspace(sp, p, &psp)); 1552 if (!psp) continue; 1553 PetscCall(PetscDualSpaceGetDM(psp, &pdm)); 1554 PetscCall(DMGetDimension(pdm, &pdim)); 1555 if (pdim < PetscAbsInt(k)) continue; 1556 PetscCall(PetscDualSpaceGetInteriorData(psp, &intNodes, &intMat)); 1557 if (intNodes == NULL && intMat == NULL) continue; 1558 PetscCall(PetscDTBinomialInt(pdim, PetscAbsInt(k), &pNk)); 1559 if (p) { 1560 PetscCall(DMPlexComputeCellGeometryAffineFEM(dm, p, v0, J, Jinv, &detJ)); 1561 } else { /* identity */ 1562 PetscInt i, j; 1563 1564 for (i = 0; i < dim; i++) 1565 for (j = 0; j < dim; j++) J[i * dim + j] = Jinv[i * dim + j] = 0.; 1566 for (i = 0; i < dim; i++) J[i * dim + i] = Jinv[i * dim + i] = 1.; 1567 for (i = 0; i < dim; i++) v0[i] = -1.; 1568 } 1569 if (pdim != dim) { /* compactify Jacobian */ 1570 PetscInt i, j; 1571 1572 for (i = 0; i < dim; i++) 1573 for (j = 0; j < pdim; j++) J[i * pdim + j] = J[i * dim + j]; 1574 } 1575 PetscCall(PetscDTAltVPullbackMatrix(pdim, dim, J, k, L)); 1576 if (intNodes) { /* push forward quadrature locations by the affine transformation */ 1577 PetscInt nNodesp; 1578 const PetscReal *nodesp; 1579 PetscInt j; 1580 1581 PetscCall(PetscQuadratureGetData(intNodes, NULL, NULL, &nNodesp, &nodesp, NULL)); 1582 for (j = 0; j < nNodesp; j++, countNodes++) { 1583 PetscInt d, e; 1584 1585 for (d = 0; d < dim; d++) { 1586 nodes[countNodes * dim + d] = v0[d]; 1587 for (e = 0; e < pdim; e++) nodes[countNodes * dim + d] += J[d * pdim + e] * (nodesp[j * pdim + e] - pv0[e]); 1588 } 1589 } 1590 } 1591 if (intMat) { 1592 PetscInt nrows; 1593 PetscInt off; 1594 1595 PetscCall(PetscSectionGetDof(section, p, &nrows)); 1596 PetscCall(PetscSectionGetOffset(section, p, &off)); 1597 for (j = 0; j < nrows; j++) { 1598 PetscInt ncols; 1599 const PetscInt *cols; 1600 const PetscScalar *vals; 1601 PetscInt l, d, e; 1602 PetscInt row = j + off; 1603 1604 PetscCall(MatGetRow(intMat, j, &ncols, &cols, &vals)); 1605 PetscCheck(ncols % pNk == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "interior matrix is not laid out as blocks of k-forms"); 1606 for (l = 0; l < ncols / pNk; l++) { 1607 PetscInt blockcol; 1608 1609 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"); 1610 blockcol = cols[l * pNk] / pNk; 1611 for (d = 0; d < Nk; d++) iwork[l * Nk + d] = (blockcol + countNodesIn) * Nk + d; 1612 for (d = 0; d < Nk; d++) work[l * Nk + d] = 0.; 1613 for (d = 0; d < Nk; d++) { 1614 for (e = 0; e < pNk; e++) { 1615 /* "push forward" dof by pulling back a k-form to be evaluated on the point: multiply on the right by L */ 1616 work[l * Nk + d] += vals[l * pNk + e] * L[e * Nk + d]; 1617 } 1618 } 1619 } 1620 PetscCall(MatSetValues(allMat, 1, &row, (ncols / pNk) * Nk, iwork, work, INSERT_VALUES)); 1621 PetscCall(MatRestoreRow(intMat, j, &ncols, &cols, &vals)); 1622 } 1623 } 1624 } 1625 PetscCall(MatAssemblyBegin(allMat, MAT_FINAL_ASSEMBLY)); 1626 PetscCall(MatAssemblyEnd(allMat, MAT_FINAL_ASSEMBLY)); 1627 PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, &allNodes)); 1628 PetscCall(PetscQuadratureSetData(allNodes, dim, 0, nNodes, nodes, NULL)); 1629 PetscCall(PetscFree7(v0, pv0, J, Jinv, L, work, iwork)); 1630 PetscCall(MatDestroy(&sp->allMat)); 1631 sp->allMat = allMat; 1632 PetscCall(PetscQuadratureDestroy(&sp->allNodes)); 1633 sp->allNodes = allNodes; 1634 PetscFunctionReturn(PETSC_SUCCESS); 1635 } 1636 1637 static PetscErrorCode PetscDualSpaceComputeFunctionalsFromAllData_Moments(PetscDualSpace sp) 1638 { 1639 Mat allMat; 1640 PetscInt momentOrder, i; 1641 PetscBool tensor = PETSC_FALSE; 1642 const PetscReal *weights; 1643 PetscScalar *array; 1644 PetscInt nDofs; 1645 PetscInt dim, Nc; 1646 DM dm; 1647 PetscQuadrature allNodes; 1648 PetscInt nNodes; 1649 1650 PetscFunctionBegin; 1651 PetscCall(PetscDualSpaceGetDM(sp, &dm)); 1652 PetscCall(DMGetDimension(dm, &dim)); 1653 PetscCall(PetscDualSpaceGetNumComponents(sp, &Nc)); 1654 PetscCall(PetscDualSpaceGetAllData(sp, &allNodes, &allMat)); 1655 PetscCall(MatGetSize(allMat, &nDofs, NULL)); 1656 PetscCheck(nDofs == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "We do not yet support moments beyond P0, nDofs == %" PetscInt_FMT, nDofs); 1657 PetscCall(PetscMalloc1(nDofs, &sp->functional)); 1658 PetscCall(PetscDualSpaceLagrangeGetMomentOrder(sp, &momentOrder)); 1659 PetscCall(PetscDualSpaceLagrangeGetTensor(sp, &tensor)); 1660 if (!tensor) PetscCall(PetscDTStroudConicalQuadrature(dim, Nc, PetscMax(momentOrder + 1, 1), -1.0, 1.0, &sp->functional[0])); 1661 else PetscCall(PetscDTGaussTensorQuadrature(dim, Nc, PetscMax(momentOrder + 1, 1), -1.0, 1.0, &sp->functional[0])); 1662 /* Need to replace allNodes and allMat */ 1663 PetscCall(PetscObjectReference((PetscObject)sp->functional[0])); 1664 PetscCall(PetscQuadratureDestroy(&sp->allNodes)); 1665 sp->allNodes = sp->functional[0]; 1666 PetscCall(PetscQuadratureGetData(sp->allNodes, NULL, NULL, &nNodes, NULL, &weights)); 1667 PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, nDofs, nNodes * Nc, NULL, &allMat)); 1668 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)allMat, "ds_")); 1669 PetscCall(MatDenseGetArrayWrite(allMat, &array)); 1670 for (i = 0; i < nNodes * Nc; ++i) array[i] = weights[i]; 1671 PetscCall(MatDenseRestoreArrayWrite(allMat, &array)); 1672 PetscCall(MatAssemblyBegin(allMat, MAT_FINAL_ASSEMBLY)); 1673 PetscCall(MatAssemblyEnd(allMat, MAT_FINAL_ASSEMBLY)); 1674 PetscCall(MatDestroy(&sp->allMat)); 1675 sp->allMat = allMat; 1676 PetscFunctionReturn(PETSC_SUCCESS); 1677 } 1678 1679 /* rather than trying to get all data from the functionals, we create 1680 * the functionals from rows of the quadrature -> dof matrix. 1681 * 1682 * Ideally most of the uses of PetscDualSpace in PetscFE will switch 1683 * to using intMat and allMat, so that the individual functionals 1684 * don't need to be constructed at all */ 1685 PETSC_INTERN PetscErrorCode PetscDualSpaceComputeFunctionalsFromAllData(PetscDualSpace sp) 1686 { 1687 PetscQuadrature allNodes; 1688 Mat allMat; 1689 PetscInt nDofs; 1690 PetscInt dim, Nc, f; 1691 DM dm; 1692 PetscInt nNodes, spdim; 1693 const PetscReal *nodes = NULL; 1694 PetscSection section; 1695 PetscBool useMoments; 1696 1697 PetscFunctionBegin; 1698 PetscCall(PetscDualSpaceGetDM(sp, &dm)); 1699 PetscCall(DMGetDimension(dm, &dim)); 1700 PetscCall(PetscDualSpaceGetNumComponents(sp, &Nc)); 1701 PetscCall(PetscDualSpaceGetAllData(sp, &allNodes, &allMat)); 1702 nNodes = 0; 1703 if (allNodes) PetscCall(PetscQuadratureGetData(allNodes, NULL, NULL, &nNodes, &nodes, NULL)); 1704 PetscCall(MatGetSize(allMat, &nDofs, NULL)); 1705 PetscCall(PetscDualSpaceGetSection(sp, §ion)); 1706 PetscCall(PetscSectionGetStorageSize(section, &spdim)); 1707 PetscCheck(spdim == nDofs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "incompatible all matrix size"); 1708 PetscCall(PetscMalloc1(nDofs, &sp->functional)); 1709 PetscCall(PetscDualSpaceLagrangeGetUseMoments(sp, &useMoments)); 1710 for (f = 0; f < nDofs; f++) { 1711 PetscInt ncols, c; 1712 const PetscInt *cols; 1713 const PetscScalar *vals; 1714 PetscReal *nodesf; 1715 PetscReal *weightsf; 1716 PetscInt nNodesf; 1717 PetscInt countNodes; 1718 1719 PetscCall(MatGetRow(allMat, f, &ncols, &cols, &vals)); 1720 for (c = 1, nNodesf = 1; c < ncols; c++) { 1721 if ((cols[c] / Nc) != (cols[c - 1] / Nc)) nNodesf++; 1722 } 1723 PetscCall(PetscMalloc1(dim * nNodesf, &nodesf)); 1724 PetscCall(PetscMalloc1(Nc * nNodesf, &weightsf)); 1725 for (c = 0, countNodes = 0; c < ncols; c++) { 1726 if (!c || ((cols[c] / Nc) != (cols[c - 1] / Nc))) { 1727 PetscInt d; 1728 1729 for (d = 0; d < Nc; d++) weightsf[countNodes * Nc + d] = 0.; 1730 for (d = 0; d < dim; d++) nodesf[countNodes * dim + d] = nodes[(cols[c] / Nc) * dim + d]; 1731 countNodes++; 1732 } 1733 weightsf[(countNodes - 1) * Nc + (cols[c] % Nc)] = PetscRealPart(vals[c]); 1734 } 1735 PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, &sp->functional[f])); 1736 PetscCall(PetscQuadratureSetData(sp->functional[f], dim, Nc, nNodesf, nodesf, weightsf)); 1737 PetscCall(MatRestoreRow(allMat, f, &ncols, &cols, &vals)); 1738 } 1739 PetscFunctionReturn(PETSC_SUCCESS); 1740 } 1741 1742 /* check if a cell is a tensor product of the segment with a facet, 1743 * specifically checking if f and f2 can be the "endpoints" (like the triangles 1744 * at either end of a wedge) */ 1745 static PetscErrorCode DMPlexPointIsTensor_Internal_Given(DM dm, PetscInt p, PetscInt f, PetscInt f2, PetscBool *isTensor) 1746 { 1747 PetscInt coneSize, c; 1748 const PetscInt *cone; 1749 const PetscInt *fCone; 1750 const PetscInt *f2Cone; 1751 PetscInt fs[2]; 1752 PetscInt meetSize, nmeet; 1753 const PetscInt *meet; 1754 1755 PetscFunctionBegin; 1756 fs[0] = f; 1757 fs[1] = f2; 1758 PetscCall(DMPlexGetMeet(dm, 2, fs, &meetSize, &meet)); 1759 nmeet = meetSize; 1760 PetscCall(DMPlexRestoreMeet(dm, 2, fs, &meetSize, &meet)); 1761 /* two points that have a non-empty meet cannot be at opposite ends of a cell */ 1762 if (nmeet) { 1763 *isTensor = PETSC_FALSE; 1764 PetscFunctionReturn(PETSC_SUCCESS); 1765 } 1766 PetscCall(DMPlexGetConeSize(dm, p, &coneSize)); 1767 PetscCall(DMPlexGetCone(dm, p, &cone)); 1768 PetscCall(DMPlexGetCone(dm, f, &fCone)); 1769 PetscCall(DMPlexGetCone(dm, f2, &f2Cone)); 1770 for (c = 0; c < coneSize; c++) { 1771 PetscInt e, ef; 1772 PetscInt d = -1, d2 = -1; 1773 PetscInt dcount, d2count; 1774 PetscInt t = cone[c]; 1775 PetscInt tConeSize; 1776 PetscBool tIsTensor; 1777 const PetscInt *tCone; 1778 1779 if (t == f || t == f2) continue; 1780 /* for every other facet in the cone, check that is has 1781 * one ridge in common with each end */ 1782 PetscCall(DMPlexGetConeSize(dm, t, &tConeSize)); 1783 PetscCall(DMPlexGetCone(dm, t, &tCone)); 1784 1785 dcount = 0; 1786 d2count = 0; 1787 for (e = 0; e < tConeSize; e++) { 1788 PetscInt q = tCone[e]; 1789 for (ef = 0; ef < coneSize - 2; ef++) { 1790 if (fCone[ef] == q) { 1791 if (dcount) { 1792 *isTensor = PETSC_FALSE; 1793 PetscFunctionReturn(PETSC_SUCCESS); 1794 } 1795 d = q; 1796 dcount++; 1797 } else if (f2Cone[ef] == q) { 1798 if (d2count) { 1799 *isTensor = PETSC_FALSE; 1800 PetscFunctionReturn(PETSC_SUCCESS); 1801 } 1802 d2 = q; 1803 d2count++; 1804 } 1805 } 1806 } 1807 /* if the whole cell is a tensor with the segment, then this 1808 * facet should be a tensor with the segment */ 1809 PetscCall(DMPlexPointIsTensor_Internal_Given(dm, t, d, d2, &tIsTensor)); 1810 if (!tIsTensor) { 1811 *isTensor = PETSC_FALSE; 1812 PetscFunctionReturn(PETSC_SUCCESS); 1813 } 1814 } 1815 *isTensor = PETSC_TRUE; 1816 PetscFunctionReturn(PETSC_SUCCESS); 1817 } 1818 1819 /* determine if a cell is a tensor with a segment by looping over pairs of facets to find a pair 1820 * that could be the opposite ends */ 1821 static PetscErrorCode DMPlexPointIsTensor_Internal(DM dm, PetscInt p, PetscBool *isTensor, PetscInt *endA, PetscInt *endB) 1822 { 1823 PetscInt coneSize, c, c2; 1824 const PetscInt *cone; 1825 1826 PetscFunctionBegin; 1827 PetscCall(DMPlexGetConeSize(dm, p, &coneSize)); 1828 if (!coneSize) { 1829 if (isTensor) *isTensor = PETSC_FALSE; 1830 if (endA) *endA = -1; 1831 if (endB) *endB = -1; 1832 } 1833 PetscCall(DMPlexGetCone(dm, p, &cone)); 1834 for (c = 0; c < coneSize; c++) { 1835 PetscInt f = cone[c]; 1836 PetscInt fConeSize; 1837 1838 PetscCall(DMPlexGetConeSize(dm, f, &fConeSize)); 1839 if (fConeSize != coneSize - 2) continue; 1840 1841 for (c2 = c + 1; c2 < coneSize; c2++) { 1842 PetscInt f2 = cone[c2]; 1843 PetscBool isTensorff2; 1844 PetscInt f2ConeSize; 1845 1846 PetscCall(DMPlexGetConeSize(dm, f2, &f2ConeSize)); 1847 if (f2ConeSize != coneSize - 2) continue; 1848 1849 PetscCall(DMPlexPointIsTensor_Internal_Given(dm, p, f, f2, &isTensorff2)); 1850 if (isTensorff2) { 1851 if (isTensor) *isTensor = PETSC_TRUE; 1852 if (endA) *endA = f; 1853 if (endB) *endB = f2; 1854 PetscFunctionReturn(PETSC_SUCCESS); 1855 } 1856 } 1857 } 1858 if (isTensor) *isTensor = PETSC_FALSE; 1859 if (endA) *endA = -1; 1860 if (endB) *endB = -1; 1861 PetscFunctionReturn(PETSC_SUCCESS); 1862 } 1863 1864 /* determine if a cell is a tensor with a segment by looping over pairs of facets to find a pair 1865 * that could be the opposite ends */ 1866 static PetscErrorCode DMPlexPointIsTensor(DM dm, PetscInt p, PetscBool *isTensor, PetscInt *endA, PetscInt *endB) 1867 { 1868 DMPlexInterpolatedFlag interpolated; 1869 1870 PetscFunctionBegin; 1871 PetscCall(DMPlexIsInterpolated(dm, &interpolated)); 1872 PetscCheck(interpolated == DMPLEX_INTERPOLATED_FULL, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Only for interpolated DMPlex's"); 1873 PetscCall(DMPlexPointIsTensor_Internal(dm, p, isTensor, endA, endB)); 1874 PetscFunctionReturn(PETSC_SUCCESS); 1875 } 1876 1877 /* Let k = formDegree and k' = -sign(k) * dim + k. Transform a symmetric frame for k-forms on the biunit simplex into 1878 * a symmetric frame for k'-forms on the biunit simplex. 1879 * 1880 * A frame is "symmetric" if the pullback of every symmetry of the biunit simplex is a permutation of the frame. 1881 * 1882 * forms in the symmetric frame are used as dofs in the untrimmed simplex spaces. This way, symmetries of the 1883 * reference cell result in permutations of dofs grouped by node. 1884 * 1885 * Use T to transform dof matrices for k'-forms into dof matrices for k-forms as a block diagonal transformation on 1886 * the right. 1887 */ 1888 static PetscErrorCode BiunitSimplexSymmetricFormTransformation(PetscInt dim, PetscInt formDegree, PetscReal T[]) 1889 { 1890 PetscInt k = formDegree; 1891 PetscInt kd = k < 0 ? dim + k : k - dim; 1892 PetscInt Nk; 1893 PetscReal *biToEq, *eqToBi, *biToEqStar, *eqToBiStar; 1894 PetscInt fact; 1895 1896 PetscFunctionBegin; 1897 PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(k), &Nk)); 1898 PetscCall(PetscCalloc4(dim * dim, &biToEq, dim * dim, &eqToBi, Nk * Nk, &biToEqStar, Nk * Nk, &eqToBiStar)); 1899 /* fill in biToEq: Jacobian of the transformation from the biunit simplex to the equilateral simplex */ 1900 fact = 0; 1901 for (PetscInt i = 0; i < dim; i++) { 1902 biToEq[i * dim + i] = PetscSqrtReal(((PetscReal)i + 2.) / (2. * ((PetscReal)i + 1.))); 1903 fact += 4 * (i + 1); 1904 for (PetscInt j = i + 1; j < dim; j++) biToEq[i * dim + j] = PetscSqrtReal(1. / (PetscReal)fact); 1905 } 1906 /* fill in eqToBi: Jacobian of the transformation from the equilateral simplex to the biunit simplex */ 1907 fact = 0; 1908 for (PetscInt j = 0; j < dim; j++) { 1909 eqToBi[j * dim + j] = PetscSqrtReal(2. * ((PetscReal)j + 1.) / ((PetscReal)j + 2)); 1910 fact += j + 1; 1911 for (PetscInt i = 0; i < j; i++) eqToBi[i * dim + j] = -PetscSqrtReal(1. / (PetscReal)fact); 1912 } 1913 PetscCall(PetscDTAltVPullbackMatrix(dim, dim, biToEq, kd, biToEqStar)); 1914 PetscCall(PetscDTAltVPullbackMatrix(dim, dim, eqToBi, k, eqToBiStar)); 1915 /* product of pullbacks simulates the following steps 1916 * 1917 * 1. start with frame W = [w_1, w_2, ..., w_m] of k forms that is symmetric on the biunit simplex: 1918 if J is the Jacobian of a symmetry of the biunit simplex, then J_k* W = [J_k*w_1, ..., J_k*w_m] 1919 is a permutation of W. 1920 Even though a k' form --- a (dim - k) form represented by its Hodge star --- has the same geometric 1921 content as a k form, W is not a symmetric frame of k' forms on the biunit simplex. That's because, 1922 for general Jacobian J, J_k* != J_k'*. 1923 * 2. pullback W to the equilateral triangle using the k pullback, W_eq = eqToBi_k* W. All symmetries of the 1924 equilateral simplex have orthonormal Jacobians. For an orthonormal Jacobian O, J_k* = J_k'*, so W_eq is 1925 also a symmetric frame for k' forms on the equilateral simplex. 1926 3. pullback W_eq back to the biunit simplex using the k' pulback, V = biToEq_k'* W_eq = biToEq_k'* eqToBi_k* W. 1927 V is a symmetric frame for k' forms on the biunit simplex. 1928 */ 1929 for (PetscInt i = 0; i < Nk; i++) { 1930 for (PetscInt j = 0; j < Nk; j++) { 1931 PetscReal val = 0.; 1932 for (PetscInt k = 0; k < Nk; k++) val += biToEqStar[i * Nk + k] * eqToBiStar[k * Nk + j]; 1933 T[i * Nk + j] = val; 1934 } 1935 } 1936 PetscCall(PetscFree4(biToEq, eqToBi, biToEqStar, eqToBiStar)); 1937 PetscFunctionReturn(PETSC_SUCCESS); 1938 } 1939 1940 /* permute a quadrature -> dof matrix so that its rows are in revlex order by nodeIdx */ 1941 static PetscErrorCode MatPermuteByNodeIdx(Mat A, PetscLagNodeIndices ni, Mat *Aperm) 1942 { 1943 PetscInt m, n, i, j; 1944 PetscInt nodeIdxDim = ni->nodeIdxDim; 1945 PetscInt nodeVecDim = ni->nodeVecDim; 1946 PetscInt *perm; 1947 IS permIS; 1948 IS id; 1949 PetscInt *nIdxPerm; 1950 PetscReal *nVecPerm; 1951 1952 PetscFunctionBegin; 1953 PetscCall(PetscLagNodeIndicesGetPermutation(ni, &perm)); 1954 PetscCall(MatGetSize(A, &m, &n)); 1955 PetscCall(PetscMalloc1(nodeIdxDim * m, &nIdxPerm)); 1956 PetscCall(PetscMalloc1(nodeVecDim * m, &nVecPerm)); 1957 for (i = 0; i < m; i++) 1958 for (j = 0; j < nodeIdxDim; j++) nIdxPerm[i * nodeIdxDim + j] = ni->nodeIdx[perm[i] * nodeIdxDim + j]; 1959 for (i = 0; i < m; i++) 1960 for (j = 0; j < nodeVecDim; j++) nVecPerm[i * nodeVecDim + j] = ni->nodeVec[perm[i] * nodeVecDim + j]; 1961 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, m, perm, PETSC_USE_POINTER, &permIS)); 1962 PetscCall(ISSetPermutation(permIS)); 1963 PetscCall(ISCreateStride(PETSC_COMM_SELF, n, 0, 1, &id)); 1964 PetscCall(ISSetPermutation(id)); 1965 PetscCall(MatPermute(A, permIS, id, Aperm)); 1966 PetscCall(ISDestroy(&permIS)); 1967 PetscCall(ISDestroy(&id)); 1968 for (i = 0; i < m; i++) perm[i] = i; 1969 PetscCall(PetscFree(ni->nodeIdx)); 1970 PetscCall(PetscFree(ni->nodeVec)); 1971 ni->nodeIdx = nIdxPerm; 1972 ni->nodeVec = nVecPerm; 1973 PetscFunctionReturn(PETSC_SUCCESS); 1974 } 1975 1976 static PetscErrorCode PetscDualSpaceSetUp_Lagrange(PetscDualSpace sp) 1977 { 1978 PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data; 1979 DM dm = sp->dm; 1980 DM dmint = NULL; 1981 PetscInt order; 1982 PetscInt Nc; 1983 MPI_Comm comm; 1984 PetscBool continuous; 1985 PetscSection section; 1986 PetscInt depth, dim, pStart, pEnd, cStart, cEnd, p, *pStratStart, *pStratEnd, d; 1987 PetscInt formDegree, Nk, Ncopies; 1988 PetscInt tensorf = -1, tensorf2 = -1; 1989 PetscBool tensorCell, tensorSpace; 1990 PetscBool uniform, trimmed; 1991 Petsc1DNodeFamily nodeFamily; 1992 PetscInt numNodeSkip; 1993 DMPlexInterpolatedFlag interpolated; 1994 PetscBool isbdm; 1995 1996 PetscFunctionBegin; 1997 /* step 1: sanitize input */ 1998 PetscCall(PetscObjectGetComm((PetscObject)sp, &comm)); 1999 PetscCall(DMGetDimension(dm, &dim)); 2000 PetscCall(PetscObjectTypeCompare((PetscObject)sp, PETSCDUALSPACEBDM, &isbdm)); 2001 if (isbdm) { 2002 sp->k = -(dim - 1); /* form degree of H-div */ 2003 PetscCall(PetscObjectChangeTypeName((PetscObject)sp, PETSCDUALSPACELAGRANGE)); 2004 } 2005 PetscCall(PetscDualSpaceGetFormDegree(sp, &formDegree)); 2006 PetscCheck(PetscAbsInt(formDegree) <= dim, comm, PETSC_ERR_ARG_OUTOFRANGE, "Form degree must be bounded by dimension"); 2007 PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(formDegree), &Nk)); 2008 if (sp->Nc <= 0 && lag->numCopies > 0) sp->Nc = Nk * lag->numCopies; 2009 Nc = sp->Nc; 2010 PetscCheck(Nc % Nk == 0, comm, PETSC_ERR_ARG_INCOMP, "Number of components is not a multiple of form degree size"); 2011 if (lag->numCopies <= 0) lag->numCopies = Nc / Nk; 2012 Ncopies = lag->numCopies; 2013 PetscCheck(Nc / Nk == Ncopies, comm, PETSC_ERR_ARG_INCOMP, "Number of copies * (dim choose k) != Nc"); 2014 if (!dim) sp->order = 0; 2015 order = sp->order; 2016 uniform = sp->uniform; 2017 PetscCheck(uniform, PETSC_COMM_SELF, PETSC_ERR_SUP, "Variable order not supported yet"); 2018 if (lag->trimmed && !formDegree) lag->trimmed = PETSC_FALSE; /* trimmed spaces are the same as full spaces for 0-forms */ 2019 if (lag->nodeType == PETSCDTNODES_DEFAULT) { 2020 lag->nodeType = PETSCDTNODES_GAUSSJACOBI; 2021 lag->nodeExponent = 0.; 2022 /* trimmed spaces don't include corner vertices, so don't use end nodes by default */ 2023 lag->endNodes = lag->trimmed ? PETSC_FALSE : PETSC_TRUE; 2024 } 2025 /* If a trimmed space and the user did choose nodes with endpoints, skip them by default */ 2026 if (lag->numNodeSkip < 0) lag->numNodeSkip = (lag->trimmed && lag->endNodes) ? 1 : 0; 2027 numNodeSkip = lag->numNodeSkip; 2028 PetscCheck(!lag->trimmed || order, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot have zeroth order trimmed elements"); 2029 if (lag->trimmed && PetscAbsInt(formDegree) == dim) { /* convert trimmed n-forms to untrimmed of one polynomial order less */ 2030 sp->order--; 2031 order--; 2032 lag->trimmed = PETSC_FALSE; 2033 } 2034 trimmed = lag->trimmed; 2035 if (!order || PetscAbsInt(formDegree) == dim) lag->continuous = PETSC_FALSE; 2036 continuous = lag->continuous; 2037 PetscCall(DMPlexGetDepth(dm, &depth)); 2038 PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 2039 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 2040 PetscCheck(pStart == 0 && cStart == 0, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Expect DM with chart starting at zero and cells first"); 2041 PetscCheck(cEnd == 1, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Use PETSCDUALSPACEREFINED for multi-cell reference meshes"); 2042 PetscCall(DMPlexIsInterpolated(dm, &interpolated)); 2043 if (interpolated != DMPLEX_INTERPOLATED_FULL) { 2044 PetscCall(DMPlexInterpolate(dm, &dmint)); 2045 } else { 2046 PetscCall(PetscObjectReference((PetscObject)dm)); 2047 dmint = dm; 2048 } 2049 tensorCell = PETSC_FALSE; 2050 if (dim > 1) PetscCall(DMPlexPointIsTensor(dmint, 0, &tensorCell, &tensorf, &tensorf2)); 2051 lag->tensorCell = tensorCell; 2052 if (dim < 2 || !lag->tensorCell) lag->tensorSpace = PETSC_FALSE; 2053 tensorSpace = lag->tensorSpace; 2054 if (!lag->nodeFamily) PetscCall(Petsc1DNodeFamilyCreate(lag->nodeType, lag->nodeExponent, lag->endNodes, &lag->nodeFamily)); 2055 nodeFamily = lag->nodeFamily; 2056 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"); 2057 2058 if (Ncopies > 1) { 2059 PetscDualSpace scalarsp; 2060 2061 PetscCall(PetscDualSpaceDuplicate(sp, &scalarsp)); 2062 /* Setting the number of components to Nk is a space with 1 copy of each k-form */ 2063 sp->setupcalled = PETSC_FALSE; 2064 PetscCall(PetscDualSpaceSetNumComponents(scalarsp, Nk)); 2065 PetscCall(PetscDualSpaceSetUp(scalarsp)); 2066 PetscCall(PetscDualSpaceSetType(sp, PETSCDUALSPACESUM)); 2067 PetscCall(PetscDualSpaceSumSetNumSubspaces(sp, Ncopies)); 2068 PetscCall(PetscDualSpaceSumSetConcatenate(sp, PETSC_TRUE)); 2069 PetscCall(PetscDualSpaceSumSetInterleave(sp, PETSC_TRUE, PETSC_FALSE)); 2070 for (PetscInt i = 0; i < Ncopies; i++) PetscCall(PetscDualSpaceSumSetSubspace(sp, i, scalarsp)); 2071 PetscCall(PetscDualSpaceSetUp(sp)); 2072 PetscCall(PetscDualSpaceDestroy(&scalarsp)); 2073 PetscCall(DMDestroy(&dmint)); 2074 PetscFunctionReturn(PETSC_SUCCESS); 2075 } 2076 2077 /* step 2: construct the boundary spaces */ 2078 PetscCall(PetscMalloc2(depth + 1, &pStratStart, depth + 1, &pStratEnd)); 2079 PetscCall(PetscCalloc1(pEnd, &sp->pointSpaces)); 2080 for (d = 0; d <= depth; ++d) PetscCall(DMPlexGetDepthStratum(dm, d, &pStratStart[d], &pStratEnd[d])); 2081 PetscCall(PetscDualSpaceSectionCreate_Internal(sp, §ion)); 2082 sp->pointSection = section; 2083 if (continuous && !lag->interiorOnly) { 2084 PetscInt h; 2085 2086 for (p = pStratStart[depth - 1]; p < pStratEnd[depth - 1]; p++) { /* calculate the facet dual spaces */ 2087 PetscReal v0[3]; 2088 DMPolytopeType ptype; 2089 PetscReal J[9], detJ; 2090 PetscInt q; 2091 2092 PetscCall(DMPlexComputeCellGeometryAffineFEM(dm, p, v0, J, NULL, &detJ)); 2093 PetscCall(DMPlexGetCellType(dm, p, &ptype)); 2094 2095 /* compare to previous facets: if computed, reference that dualspace */ 2096 for (q = pStratStart[depth - 1]; q < p; q++) { 2097 DMPolytopeType qtype; 2098 2099 PetscCall(DMPlexGetCellType(dm, q, &qtype)); 2100 if (qtype == ptype) break; 2101 } 2102 if (q < p) { /* this facet has the same dual space as that one */ 2103 PetscCall(PetscObjectReference((PetscObject)sp->pointSpaces[q])); 2104 sp->pointSpaces[p] = sp->pointSpaces[q]; 2105 continue; 2106 } 2107 /* if not, recursively compute this dual space */ 2108 PetscCall(PetscDualSpaceCreateFacetSubspace_Lagrange(sp, NULL, p, formDegree, Ncopies, PETSC_FALSE, &sp->pointSpaces[p])); 2109 } 2110 for (h = 2; h <= depth; h++) { /* get the higher subspaces from the facet subspaces */ 2111 PetscInt hd = depth - h; 2112 PetscInt hdim = dim - h; 2113 2114 if (hdim < PetscAbsInt(formDegree)) break; 2115 for (p = pStratStart[hd]; p < pStratEnd[hd]; p++) { 2116 PetscInt suppSize, s; 2117 const PetscInt *supp; 2118 2119 PetscCall(DMPlexGetSupportSize(dm, p, &suppSize)); 2120 PetscCall(DMPlexGetSupport(dm, p, &supp)); 2121 for (s = 0; s < suppSize; s++) { 2122 DM qdm; 2123 PetscDualSpace qsp, psp; 2124 PetscInt c, coneSize, q; 2125 const PetscInt *cone; 2126 const PetscInt *refCone; 2127 2128 q = supp[s]; 2129 qsp = sp->pointSpaces[q]; 2130 PetscCall(DMPlexGetConeSize(dm, q, &coneSize)); 2131 PetscCall(DMPlexGetCone(dm, q, &cone)); 2132 for (c = 0; c < coneSize; c++) 2133 if (cone[c] == p) break; 2134 PetscCheck(c != coneSize, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "cone/support mismatch"); 2135 PetscCall(PetscDualSpaceGetDM(qsp, &qdm)); 2136 PetscCall(DMPlexGetCone(qdm, 0, &refCone)); 2137 /* get the equivalent dual space from the support dual space */ 2138 PetscCall(PetscDualSpaceGetPointSubspace(qsp, refCone[c], &psp)); 2139 if (!s) { 2140 PetscCall(PetscObjectReference((PetscObject)psp)); 2141 sp->pointSpaces[p] = psp; 2142 } 2143 } 2144 } 2145 } 2146 for (p = 1; p < pEnd; p++) { 2147 PetscInt pspdim; 2148 if (!sp->pointSpaces[p]) continue; 2149 PetscCall(PetscDualSpaceGetInteriorDimension(sp->pointSpaces[p], &pspdim)); 2150 PetscCall(PetscSectionSetDof(section, p, pspdim)); 2151 } 2152 } 2153 2154 if (trimmed && !continuous) { 2155 /* the dofs of a trimmed space don't have a nice tensor/lattice structure: 2156 * just construct the continuous dual space and copy all of the data over, 2157 * allocating it all to the cell instead of splitting it up between the boundaries */ 2158 PetscDualSpace spcont; 2159 PetscInt spdim, f; 2160 PetscQuadrature allNodes; 2161 PetscDualSpace_Lag *lagc; 2162 Mat allMat; 2163 2164 PetscCall(PetscDualSpaceDuplicate(sp, &spcont)); 2165 PetscCall(PetscDualSpaceLagrangeSetContinuity(spcont, PETSC_TRUE)); 2166 PetscCall(PetscDualSpaceSetUp(spcont)); 2167 PetscCall(PetscDualSpaceGetDimension(spcont, &spdim)); 2168 sp->spdim = sp->spintdim = spdim; 2169 PetscCall(PetscSectionSetDof(section, 0, spdim)); 2170 PetscCall(PetscDualSpaceSectionSetUp_Internal(sp, section)); 2171 PetscCall(PetscMalloc1(spdim, &sp->functional)); 2172 for (f = 0; f < spdim; f++) { 2173 PetscQuadrature fn; 2174 2175 PetscCall(PetscDualSpaceGetFunctional(spcont, f, &fn)); 2176 PetscCall(PetscObjectReference((PetscObject)fn)); 2177 sp->functional[f] = fn; 2178 } 2179 PetscCall(PetscDualSpaceGetAllData(spcont, &allNodes, &allMat)); 2180 PetscCall(PetscObjectReference((PetscObject)allNodes)); 2181 PetscCall(PetscObjectReference((PetscObject)allNodes)); 2182 sp->allNodes = sp->intNodes = allNodes; 2183 PetscCall(PetscObjectReference((PetscObject)allMat)); 2184 PetscCall(PetscObjectReference((PetscObject)allMat)); 2185 sp->allMat = sp->intMat = allMat; 2186 lagc = (PetscDualSpace_Lag *)spcont->data; 2187 PetscCall(PetscLagNodeIndicesReference(lagc->vertIndices)); 2188 lag->vertIndices = lagc->vertIndices; 2189 PetscCall(PetscLagNodeIndicesReference(lagc->allNodeIndices)); 2190 PetscCall(PetscLagNodeIndicesReference(lagc->allNodeIndices)); 2191 lag->intNodeIndices = lagc->allNodeIndices; 2192 lag->allNodeIndices = lagc->allNodeIndices; 2193 PetscCall(PetscDualSpaceDestroy(&spcont)); 2194 PetscCall(PetscFree2(pStratStart, pStratEnd)); 2195 PetscCall(DMDestroy(&dmint)); 2196 PetscFunctionReturn(PETSC_SUCCESS); 2197 } 2198 2199 /* step 3: construct intNodes, and intMat, and combine it with boundray data to make allNodes and allMat */ 2200 if (!tensorSpace) { 2201 if (!tensorCell) PetscCall(PetscLagNodeIndicesCreateSimplexVertices(dm, &lag->vertIndices)); 2202 2203 if (trimmed) { 2204 /* there is one dof in the interior of the a trimmed element for each full polynomial of with degree at most 2205 * order + k - dim - 1 */ 2206 if (order + PetscAbsInt(formDegree) > dim) { 2207 PetscInt sum = order + PetscAbsInt(formDegree) - dim - 1; 2208 PetscInt nDofs; 2209 2210 PetscCall(PetscDualSpaceLagrangeCreateSimplexNodeMat(nodeFamily, dim, sum, Nk, numNodeSkip, &sp->intNodes, &sp->intMat, &lag->intNodeIndices)); 2211 PetscCall(MatGetSize(sp->intMat, &nDofs, NULL)); 2212 PetscCall(PetscSectionSetDof(section, 0, nDofs)); 2213 } 2214 PetscCall(PetscDualSpaceSectionSetUp_Internal(sp, section)); 2215 PetscCall(PetscDualSpaceCreateAllDataFromInteriorData(sp)); 2216 PetscCall(PetscDualSpaceLagrangeCreateAllNodeIdx(sp)); 2217 } else { 2218 if (!continuous) { 2219 /* if discontinuous just construct one node for each set of dofs (a set of dofs is a basis for the k-form 2220 * space) */ 2221 PetscInt sum = order; 2222 PetscInt nDofs; 2223 2224 PetscCall(PetscDualSpaceLagrangeCreateSimplexNodeMat(nodeFamily, dim, sum, Nk, numNodeSkip, &sp->intNodes, &sp->intMat, &lag->intNodeIndices)); 2225 PetscCall(MatGetSize(sp->intMat, &nDofs, NULL)); 2226 PetscCall(PetscSectionSetDof(section, 0, nDofs)); 2227 PetscCall(PetscDualSpaceSectionSetUp_Internal(sp, section)); 2228 PetscCall(PetscObjectReference((PetscObject)sp->intNodes)); 2229 sp->allNodes = sp->intNodes; 2230 PetscCall(PetscObjectReference((PetscObject)sp->intMat)); 2231 sp->allMat = sp->intMat; 2232 PetscCall(PetscLagNodeIndicesReference(lag->intNodeIndices)); 2233 lag->allNodeIndices = lag->intNodeIndices; 2234 } else { 2235 /* there is one dof in the interior of the a full element for each trimmed polynomial of with degree at most 2236 * order + k - dim, but with complementary form degree */ 2237 if (order + PetscAbsInt(formDegree) > dim) { 2238 PetscDualSpace trimmedsp; 2239 PetscDualSpace_Lag *trimmedlag; 2240 PetscQuadrature intNodes; 2241 PetscInt trFormDegree = formDegree >= 0 ? formDegree - dim : dim - PetscAbsInt(formDegree); 2242 PetscInt nDofs; 2243 Mat intMat; 2244 2245 PetscCall(PetscDualSpaceDuplicate(sp, &trimmedsp)); 2246 PetscCall(PetscDualSpaceLagrangeSetTrimmed(trimmedsp, PETSC_TRUE)); 2247 PetscCall(PetscDualSpaceSetOrder(trimmedsp, order + PetscAbsInt(formDegree) - dim)); 2248 PetscCall(PetscDualSpaceSetFormDegree(trimmedsp, trFormDegree)); 2249 trimmedlag = (PetscDualSpace_Lag *)trimmedsp->data; 2250 trimmedlag->numNodeSkip = numNodeSkip + 1; 2251 PetscCall(PetscDualSpaceSetUp(trimmedsp)); 2252 PetscCall(PetscDualSpaceGetAllData(trimmedsp, &intNodes, &intMat)); 2253 PetscCall(PetscObjectReference((PetscObject)intNodes)); 2254 sp->intNodes = intNodes; 2255 PetscCall(PetscLagNodeIndicesReference(trimmedlag->allNodeIndices)); 2256 lag->intNodeIndices = trimmedlag->allNodeIndices; 2257 PetscCall(PetscObjectReference((PetscObject)intMat)); 2258 if (PetscAbsInt(formDegree) > 0 && PetscAbsInt(formDegree) < dim) { 2259 PetscReal *T; 2260 PetscScalar *work; 2261 PetscInt nCols, nRows; 2262 Mat intMatT; 2263 2264 PetscCall(MatDuplicate(intMat, MAT_COPY_VALUES, &intMatT)); 2265 PetscCall(MatGetSize(intMat, &nRows, &nCols)); 2266 PetscCall(PetscMalloc2(Nk * Nk, &T, nCols, &work)); 2267 PetscCall(BiunitSimplexSymmetricFormTransformation(dim, formDegree, T)); 2268 for (PetscInt row = 0; row < nRows; row++) { 2269 PetscInt nrCols; 2270 const PetscInt *rCols; 2271 const PetscScalar *rVals; 2272 2273 PetscCall(MatGetRow(intMat, row, &nrCols, &rCols, &rVals)); 2274 PetscCheck(nrCols % Nk == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "nonzeros in intMat matrix are not in k-form size blocks"); 2275 for (PetscInt b = 0; b < nrCols; b += Nk) { 2276 const PetscScalar *v = &rVals[b]; 2277 PetscScalar *w = &work[b]; 2278 for (PetscInt j = 0; j < Nk; j++) { 2279 w[j] = 0.; 2280 for (PetscInt i = 0; i < Nk; i++) w[j] += v[i] * T[i * Nk + j]; 2281 } 2282 } 2283 PetscCall(MatSetValuesBlocked(intMatT, 1, &row, nrCols, rCols, work, INSERT_VALUES)); 2284 PetscCall(MatRestoreRow(intMat, row, &nrCols, &rCols, &rVals)); 2285 } 2286 PetscCall(MatAssemblyBegin(intMatT, MAT_FINAL_ASSEMBLY)); 2287 PetscCall(MatAssemblyEnd(intMatT, MAT_FINAL_ASSEMBLY)); 2288 PetscCall(MatDestroy(&intMat)); 2289 intMat = intMatT; 2290 PetscCall(PetscLagNodeIndicesDestroy(&lag->intNodeIndices)); 2291 PetscCall(PetscLagNodeIndicesDuplicate(trimmedlag->allNodeIndices, &lag->intNodeIndices)); 2292 { 2293 PetscInt nNodes = lag->intNodeIndices->nNodes; 2294 PetscReal *newNodeVec = lag->intNodeIndices->nodeVec; 2295 const PetscReal *oldNodeVec = trimmedlag->allNodeIndices->nodeVec; 2296 2297 for (PetscInt n = 0; n < nNodes; n++) { 2298 PetscReal *w = &newNodeVec[n * Nk]; 2299 const PetscReal *v = &oldNodeVec[n * Nk]; 2300 2301 for (PetscInt j = 0; j < Nk; j++) { 2302 w[j] = 0.; 2303 for (PetscInt i = 0; i < Nk; i++) w[j] += v[i] * T[i * Nk + j]; 2304 } 2305 } 2306 } 2307 PetscCall(PetscFree2(T, work)); 2308 } 2309 sp->intMat = intMat; 2310 PetscCall(MatGetSize(sp->intMat, &nDofs, NULL)); 2311 PetscCall(PetscDualSpaceDestroy(&trimmedsp)); 2312 PetscCall(PetscSectionSetDof(section, 0, nDofs)); 2313 } 2314 PetscCall(PetscDualSpaceSectionSetUp_Internal(sp, section)); 2315 PetscCall(PetscDualSpaceCreateAllDataFromInteriorData(sp)); 2316 PetscCall(PetscDualSpaceLagrangeCreateAllNodeIdx(sp)); 2317 } 2318 } 2319 } else { 2320 PetscQuadrature intNodesTrace = NULL; 2321 PetscQuadrature intNodesFiber = NULL; 2322 PetscQuadrature intNodes = NULL; 2323 PetscLagNodeIndices intNodeIndices = NULL; 2324 Mat intMat = NULL; 2325 2326 if (PetscAbsInt(formDegree) < dim) { /* get the trace k-forms on the first facet, and the 0-forms on the edge, 2327 and wedge them together to create some of the k-form dofs */ 2328 PetscDualSpace trace, fiber; 2329 PetscDualSpace_Lag *tracel, *fiberl; 2330 Mat intMatTrace, intMatFiber; 2331 2332 if (sp->pointSpaces[tensorf]) { 2333 PetscCall(PetscObjectReference((PetscObject)sp->pointSpaces[tensorf])); 2334 trace = sp->pointSpaces[tensorf]; 2335 } else { 2336 PetscCall(PetscDualSpaceCreateFacetSubspace_Lagrange(sp, NULL, tensorf, formDegree, Ncopies, PETSC_TRUE, &trace)); 2337 } 2338 PetscCall(PetscDualSpaceCreateEdgeSubspace_Lagrange(sp, order, 0, 1, PETSC_TRUE, &fiber)); 2339 tracel = (PetscDualSpace_Lag *)trace->data; 2340 fiberl = (PetscDualSpace_Lag *)fiber->data; 2341 PetscCall(PetscLagNodeIndicesCreateTensorVertices(dm, tracel->vertIndices, &lag->vertIndices)); 2342 PetscCall(PetscDualSpaceGetInteriorData(trace, &intNodesTrace, &intMatTrace)); 2343 PetscCall(PetscDualSpaceGetInteriorData(fiber, &intNodesFiber, &intMatFiber)); 2344 if (intNodesTrace && intNodesFiber) { 2345 PetscCall(PetscQuadratureCreateTensor(intNodesTrace, intNodesFiber, &intNodes)); 2346 PetscCall(MatTensorAltV(intMatTrace, intMatFiber, dim - 1, formDegree, 1, 0, &intMat)); 2347 PetscCall(PetscLagNodeIndicesTensor(tracel->intNodeIndices, dim - 1, formDegree, fiberl->intNodeIndices, 1, 0, &intNodeIndices)); 2348 } 2349 PetscCall(PetscObjectReference((PetscObject)intNodesTrace)); 2350 PetscCall(PetscObjectReference((PetscObject)intNodesFiber)); 2351 PetscCall(PetscDualSpaceDestroy(&fiber)); 2352 PetscCall(PetscDualSpaceDestroy(&trace)); 2353 } 2354 if (PetscAbsInt(formDegree) > 0) { /* get the trace (k-1)-forms on the first facet, and the 1-forms on the edge, 2355 and wedge them together to create the remaining k-form dofs */ 2356 PetscDualSpace trace, fiber; 2357 PetscDualSpace_Lag *tracel, *fiberl; 2358 PetscQuadrature intNodesTrace2, intNodesFiber2, intNodes2; 2359 PetscLagNodeIndices intNodeIndices2; 2360 Mat intMatTrace, intMatFiber, intMat2; 2361 PetscInt traceDegree = formDegree > 0 ? formDegree - 1 : formDegree + 1; 2362 PetscInt fiberDegree = formDegree > 0 ? 1 : -1; 2363 2364 PetscCall(PetscDualSpaceCreateFacetSubspace_Lagrange(sp, NULL, tensorf, traceDegree, Ncopies, PETSC_TRUE, &trace)); 2365 PetscCall(PetscDualSpaceCreateEdgeSubspace_Lagrange(sp, order, fiberDegree, 1, PETSC_TRUE, &fiber)); 2366 tracel = (PetscDualSpace_Lag *)trace->data; 2367 fiberl = (PetscDualSpace_Lag *)fiber->data; 2368 if (!lag->vertIndices) PetscCall(PetscLagNodeIndicesCreateTensorVertices(dm, tracel->vertIndices, &lag->vertIndices)); 2369 PetscCall(PetscDualSpaceGetInteriorData(trace, &intNodesTrace2, &intMatTrace)); 2370 PetscCall(PetscDualSpaceGetInteriorData(fiber, &intNodesFiber2, &intMatFiber)); 2371 if (intNodesTrace2 && intNodesFiber2) { 2372 PetscCall(PetscQuadratureCreateTensor(intNodesTrace2, intNodesFiber2, &intNodes2)); 2373 PetscCall(MatTensorAltV(intMatTrace, intMatFiber, dim - 1, traceDegree, 1, fiberDegree, &intMat2)); 2374 PetscCall(PetscLagNodeIndicesTensor(tracel->intNodeIndices, dim - 1, traceDegree, fiberl->intNodeIndices, 1, fiberDegree, &intNodeIndices2)); 2375 if (!intMat) { 2376 intMat = intMat2; 2377 intNodes = intNodes2; 2378 intNodeIndices = intNodeIndices2; 2379 } else { 2380 /* merge the matrices, quadrature points, and nodes */ 2381 PetscInt nM; 2382 PetscInt nDof, nDof2; 2383 PetscInt *toMerged = NULL, *toMerged2 = NULL; 2384 PetscQuadrature merged = NULL; 2385 PetscLagNodeIndices intNodeIndicesMerged = NULL; 2386 Mat matMerged = NULL; 2387 2388 PetscCall(MatGetSize(intMat, &nDof, NULL)); 2389 PetscCall(MatGetSize(intMat2, &nDof2, NULL)); 2390 PetscCall(PetscQuadraturePointsMerge(intNodes, intNodes2, &merged, &toMerged, &toMerged2)); 2391 PetscCall(PetscQuadratureGetData(merged, NULL, NULL, &nM, NULL, NULL)); 2392 PetscCall(MatricesMerge(intMat, intMat2, dim, formDegree, nM, toMerged, toMerged2, &matMerged)); 2393 PetscCall(PetscLagNodeIndicesMerge(intNodeIndices, intNodeIndices2, &intNodeIndicesMerged)); 2394 PetscCall(PetscFree(toMerged)); 2395 PetscCall(PetscFree(toMerged2)); 2396 PetscCall(MatDestroy(&intMat)); 2397 PetscCall(MatDestroy(&intMat2)); 2398 PetscCall(PetscQuadratureDestroy(&intNodes)); 2399 PetscCall(PetscQuadratureDestroy(&intNodes2)); 2400 PetscCall(PetscLagNodeIndicesDestroy(&intNodeIndices)); 2401 PetscCall(PetscLagNodeIndicesDestroy(&intNodeIndices2)); 2402 intNodes = merged; 2403 intMat = matMerged; 2404 intNodeIndices = intNodeIndicesMerged; 2405 if (!trimmed) { 2406 /* I think users expect that, when a node has a full basis for the k-forms, 2407 * they should be consecutive dofs. That isn't the case for trimmed spaces, 2408 * but is for some of the nodes in untrimmed spaces, so in that case we 2409 * sort them to group them by node */ 2410 Mat intMatPerm; 2411 2412 PetscCall(MatPermuteByNodeIdx(intMat, intNodeIndices, &intMatPerm)); 2413 PetscCall(MatDestroy(&intMat)); 2414 intMat = intMatPerm; 2415 } 2416 } 2417 } 2418 PetscCall(PetscDualSpaceDestroy(&fiber)); 2419 PetscCall(PetscDualSpaceDestroy(&trace)); 2420 } 2421 PetscCall(PetscQuadratureDestroy(&intNodesTrace)); 2422 PetscCall(PetscQuadratureDestroy(&intNodesFiber)); 2423 sp->intNodes = intNodes; 2424 sp->intMat = intMat; 2425 lag->intNodeIndices = intNodeIndices; 2426 { 2427 PetscInt nDofs = 0; 2428 2429 if (intMat) PetscCall(MatGetSize(intMat, &nDofs, NULL)); 2430 PetscCall(PetscSectionSetDof(section, 0, nDofs)); 2431 } 2432 PetscCall(PetscDualSpaceSectionSetUp_Internal(sp, section)); 2433 if (continuous) { 2434 PetscCall(PetscDualSpaceCreateAllDataFromInteriorData(sp)); 2435 PetscCall(PetscDualSpaceLagrangeCreateAllNodeIdx(sp)); 2436 } else { 2437 PetscCall(PetscObjectReference((PetscObject)intNodes)); 2438 sp->allNodes = intNodes; 2439 PetscCall(PetscObjectReference((PetscObject)intMat)); 2440 sp->allMat = intMat; 2441 PetscCall(PetscLagNodeIndicesReference(intNodeIndices)); 2442 lag->allNodeIndices = intNodeIndices; 2443 } 2444 } 2445 PetscCall(PetscSectionGetStorageSize(section, &sp->spdim)); 2446 PetscCall(PetscSectionGetConstrainedStorageSize(section, &sp->spintdim)); 2447 // TODO: fix this, computing functionals from moments should be no different for nodal vs modal 2448 if (lag->useMoments) { 2449 PetscCall(PetscDualSpaceComputeFunctionalsFromAllData_Moments(sp)); 2450 } else { 2451 PetscCall(PetscDualSpaceComputeFunctionalsFromAllData(sp)); 2452 } 2453 PetscCall(PetscFree2(pStratStart, pStratEnd)); 2454 PetscCall(DMDestroy(&dmint)); 2455 PetscFunctionReturn(PETSC_SUCCESS); 2456 } 2457 2458 /* Create a matrix that represents the transformation that DMPlexVecGetClosure() would need 2459 * to get the representation of the dofs for a mesh point if the mesh point had this orientation 2460 * relative to the cell */ 2461 PetscErrorCode PetscDualSpaceCreateInteriorSymmetryMatrix_Lagrange(PetscDualSpace sp, PetscInt ornt, Mat *symMat) 2462 { 2463 PetscDualSpace_Lag *lag; 2464 DM dm; 2465 PetscLagNodeIndices vertIndices, intNodeIndices; 2466 PetscLagNodeIndices ni; 2467 PetscInt nodeIdxDim, nodeVecDim, nNodes; 2468 PetscInt formDegree; 2469 PetscInt *perm, *permOrnt; 2470 PetscInt *nnz; 2471 PetscInt n; 2472 PetscInt maxGroupSize; 2473 PetscScalar *V, *W, *work; 2474 Mat A; 2475 2476 PetscFunctionBegin; 2477 if (!sp->spintdim) { 2478 *symMat = NULL; 2479 PetscFunctionReturn(PETSC_SUCCESS); 2480 } 2481 lag = (PetscDualSpace_Lag *)sp->data; 2482 vertIndices = lag->vertIndices; 2483 intNodeIndices = lag->intNodeIndices; 2484 PetscCall(PetscDualSpaceGetDM(sp, &dm)); 2485 PetscCall(PetscDualSpaceGetFormDegree(sp, &formDegree)); 2486 PetscCall(PetscNew(&ni)); 2487 ni->refct = 1; 2488 ni->nodeIdxDim = nodeIdxDim = intNodeIndices->nodeIdxDim; 2489 ni->nodeVecDim = nodeVecDim = intNodeIndices->nodeVecDim; 2490 ni->nNodes = nNodes = intNodeIndices->nNodes; 2491 PetscCall(PetscMalloc1(nNodes * nodeIdxDim, &ni->nodeIdx)); 2492 PetscCall(PetscMalloc1(nNodes * nodeVecDim, &ni->nodeVec)); 2493 /* push forward the dofs by the symmetry of the reference element induced by ornt */ 2494 PetscCall(PetscLagNodeIndicesPushForward(dm, vertIndices, 0, vertIndices, intNodeIndices, ornt, formDegree, ni->nodeIdx, ni->nodeVec)); 2495 /* get the revlex order for both the original and transformed dofs */ 2496 PetscCall(PetscLagNodeIndicesGetPermutation(intNodeIndices, &perm)); 2497 PetscCall(PetscLagNodeIndicesGetPermutation(ni, &permOrnt)); 2498 PetscCall(PetscMalloc1(nNodes, &nnz)); 2499 for (n = 0, maxGroupSize = 0; n < nNodes;) { /* incremented in the loop */ 2500 PetscInt *nind = &(ni->nodeIdx[permOrnt[n] * nodeIdxDim]); 2501 PetscInt m, nEnd; 2502 PetscInt groupSize; 2503 /* for each group of dofs that have the same nodeIdx coordinate */ 2504 for (nEnd = n + 1; nEnd < nNodes; nEnd++) { 2505 PetscInt *mind = &(ni->nodeIdx[permOrnt[nEnd] * nodeIdxDim]); 2506 PetscInt d; 2507 2508 /* compare the oriented permutation indices */ 2509 for (d = 0; d < nodeIdxDim; d++) 2510 if (mind[d] != nind[d]) break; 2511 if (d < nodeIdxDim) break; 2512 } 2513 /* permOrnt[[n, nEnd)] is a group of dofs that, under the symmetry are at the same location */ 2514 2515 /* the symmetry had better map the group of dofs with the same permuted nodeIdx 2516 * to a group of dofs with the same size, otherwise we messed up */ 2517 if (PetscDefined(USE_DEBUG)) { 2518 PetscInt m; 2519 PetscInt *nind = &(intNodeIndices->nodeIdx[perm[n] * nodeIdxDim]); 2520 2521 for (m = n + 1; m < nEnd; m++) { 2522 PetscInt *mind = &(intNodeIndices->nodeIdx[perm[m] * nodeIdxDim]); 2523 PetscInt d; 2524 2525 /* compare the oriented permutation indices */ 2526 for (d = 0; d < nodeIdxDim; d++) 2527 if (mind[d] != nind[d]) break; 2528 if (d < nodeIdxDim) break; 2529 } 2530 PetscCheck(m >= nEnd, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Dofs with same index after symmetry not same block size"); 2531 } 2532 groupSize = nEnd - n; 2533 /* each pushforward dof vector will be expressed in a basis of the unpermuted dofs */ 2534 for (m = n; m < nEnd; m++) nnz[permOrnt[m]] = groupSize; 2535 2536 maxGroupSize = PetscMax(maxGroupSize, nEnd - n); 2537 n = nEnd; 2538 } 2539 PetscCheck(maxGroupSize <= nodeVecDim, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Dofs are not in blocks that can be solved"); 2540 PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, nNodes, nNodes, 0, nnz, &A)); 2541 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)A, "lag_")); 2542 PetscCall(PetscFree(nnz)); 2543 PetscCall(PetscMalloc3(maxGroupSize * nodeVecDim, &V, maxGroupSize * nodeVecDim, &W, nodeVecDim * 2, &work)); 2544 for (n = 0; n < nNodes;) { /* incremented in the loop */ 2545 PetscInt *nind = &(ni->nodeIdx[permOrnt[n] * nodeIdxDim]); 2546 PetscInt nEnd; 2547 PetscInt m; 2548 PetscInt groupSize; 2549 for (nEnd = n + 1; nEnd < nNodes; nEnd++) { 2550 PetscInt *mind = &(ni->nodeIdx[permOrnt[nEnd] * nodeIdxDim]); 2551 PetscInt d; 2552 2553 /* compare the oriented permutation indices */ 2554 for (d = 0; d < nodeIdxDim; d++) 2555 if (mind[d] != nind[d]) break; 2556 if (d < nodeIdxDim) break; 2557 } 2558 groupSize = nEnd - n; 2559 /* get all of the vectors from the original and all of the pushforward vectors */ 2560 for (m = n; m < nEnd; m++) { 2561 PetscInt d; 2562 2563 for (d = 0; d < nodeVecDim; d++) { 2564 V[(m - n) * nodeVecDim + d] = intNodeIndices->nodeVec[perm[m] * nodeVecDim + d]; 2565 W[(m - n) * nodeVecDim + d] = ni->nodeVec[permOrnt[m] * nodeVecDim + d]; 2566 } 2567 } 2568 /* now we have to solve for W in terms of V: the systems isn't always square, but the span 2569 * of V and W should always be the same, so the solution of the normal equations works */ 2570 { 2571 char transpose = 'N'; 2572 PetscBLASInt bm, bn, bnrhs, blda, bldb, blwork, info; 2573 2574 PetscCall(PetscBLASIntCast(nodeVecDim, &bm)); 2575 PetscCall(PetscBLASIntCast(groupSize, &bn)); 2576 PetscCall(PetscBLASIntCast(groupSize, &bnrhs)); 2577 PetscCall(PetscBLASIntCast(bm, &blda)); 2578 PetscCall(PetscBLASIntCast(bm, &bldb)); 2579 PetscCall(PetscBLASIntCast(2 * nodeVecDim, &blwork)); 2580 PetscCallBLAS("LAPACKgels", LAPACKgels_(&transpose, &bm, &bn, &bnrhs, V, &blda, W, &bldb, work, &blwork, &info)); 2581 PetscCheck(info == 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Bad argument to GELS"); 2582 /* repack */ 2583 { 2584 PetscInt i, j; 2585 2586 for (i = 0; i < groupSize; i++) { 2587 for (j = 0; j < groupSize; j++) { 2588 /* notice the different leading dimension */ 2589 V[i * groupSize + j] = W[i * nodeVecDim + j]; 2590 } 2591 } 2592 } 2593 if (PetscDefined(USE_DEBUG)) { 2594 PetscReal res; 2595 2596 /* check that the normal error is 0 */ 2597 for (m = n; m < nEnd; m++) { 2598 PetscInt d; 2599 2600 for (d = 0; d < nodeVecDim; d++) W[(m - n) * nodeVecDim + d] = ni->nodeVec[permOrnt[m] * nodeVecDim + d]; 2601 } 2602 res = 0.; 2603 for (PetscInt i = 0; i < groupSize; i++) { 2604 for (PetscInt j = 0; j < nodeVecDim; j++) { 2605 for (PetscInt k = 0; k < groupSize; k++) W[i * nodeVecDim + j] -= V[i * groupSize + k] * intNodeIndices->nodeVec[perm[n + k] * nodeVecDim + j]; 2606 res += PetscAbsScalar(W[i * nodeVecDim + j]); 2607 } 2608 } 2609 PetscCheck(res <= PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_LIB, "Dof block did not solve"); 2610 } 2611 } 2612 PetscCall(MatSetValues(A, groupSize, &permOrnt[n], groupSize, &perm[n], V, INSERT_VALUES)); 2613 n = nEnd; 2614 } 2615 PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 2616 PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 2617 *symMat = A; 2618 PetscCall(PetscFree3(V, W, work)); 2619 PetscCall(PetscLagNodeIndicesDestroy(&ni)); 2620 PetscFunctionReturn(PETSC_SUCCESS); 2621 } 2622 2623 // get the symmetries of closure points 2624 PETSC_INTERN PetscErrorCode PetscDualSpaceGetBoundarySymmetries_Internal(PetscDualSpace sp, PetscInt ***symperms, PetscScalar ***symflips) 2625 { 2626 PetscInt closureSize = 0; 2627 PetscInt *closure = NULL; 2628 PetscInt r; 2629 2630 PetscFunctionBegin; 2631 PetscCall(DMPlexGetTransitiveClosure(sp->dm, 0, PETSC_TRUE, &closureSize, &closure)); 2632 for (r = 0; r < closureSize; r++) { 2633 PetscDualSpace psp; 2634 PetscInt point = closure[2 * r]; 2635 PetscInt pspintdim; 2636 const PetscInt ***psymperms = NULL; 2637 const PetscScalar ***psymflips = NULL; 2638 2639 if (!point) continue; 2640 PetscCall(PetscDualSpaceGetPointSubspace(sp, point, &psp)); 2641 if (!psp) continue; 2642 PetscCall(PetscDualSpaceGetInteriorDimension(psp, &pspintdim)); 2643 if (!pspintdim) continue; 2644 PetscCall(PetscDualSpaceGetSymmetries(psp, &psymperms, &psymflips)); 2645 symperms[r] = (PetscInt **)(psymperms ? psymperms[0] : NULL); 2646 symflips[r] = (PetscScalar **)(psymflips ? psymflips[0] : NULL); 2647 } 2648 PetscCall(DMPlexRestoreTransitiveClosure(sp->dm, 0, PETSC_TRUE, &closureSize, &closure)); 2649 PetscFunctionReturn(PETSC_SUCCESS); 2650 } 2651 2652 #define BaryIndex(perEdge, a, b, c) (((b) * (2 * perEdge + 1 - (b))) / 2) + (c) 2653 2654 #define CartIndex(perEdge, a, b) (perEdge * (a) + b) 2655 2656 /* the existing interface for symmetries is insufficient for all cases: 2657 * - it should be sufficient for form degrees that are scalar (0 and n) 2658 * - it should be sufficient for hypercube dofs 2659 * - it isn't sufficient for simplex cells with non-scalar form degrees if 2660 * there are any dofs in the interior 2661 * 2662 * We compute the general transformation matrices, and if they fit, we return them, 2663 * otherwise we error (but we should probably change the interface to allow for 2664 * these symmetries) 2665 */ 2666 static PetscErrorCode PetscDualSpaceGetSymmetries_Lagrange(PetscDualSpace sp, const PetscInt ****perms, const PetscScalar ****flips) 2667 { 2668 PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data; 2669 PetscInt dim, order, Nc; 2670 2671 PetscFunctionBegin; 2672 PetscCall(PetscDualSpaceGetOrder(sp, &order)); 2673 PetscCall(PetscDualSpaceGetNumComponents(sp, &Nc)); 2674 PetscCall(DMGetDimension(sp->dm, &dim)); 2675 if (!lag->symComputed) { /* store symmetries */ 2676 PetscInt pStart, pEnd, p; 2677 PetscInt numPoints; 2678 PetscInt numFaces; 2679 PetscInt spintdim; 2680 PetscInt ***symperms; 2681 PetscScalar ***symflips; 2682 2683 PetscCall(DMPlexGetChart(sp->dm, &pStart, &pEnd)); 2684 numPoints = pEnd - pStart; 2685 { 2686 DMPolytopeType ct; 2687 /* The number of arrangements is no longer based on the number of faces */ 2688 PetscCall(DMPlexGetCellType(sp->dm, 0, &ct)); 2689 numFaces = DMPolytopeTypeGetNumArrangements(ct) / 2; 2690 } 2691 PetscCall(PetscCalloc1(numPoints, &symperms)); 2692 PetscCall(PetscCalloc1(numPoints, &symflips)); 2693 spintdim = sp->spintdim; 2694 /* The nodal symmetry behavior is not present when tensorSpace != tensorCell: someone might want this for the "S" 2695 * family of FEEC spaces. Most used in particular are discontinuous polynomial L2 spaces in tensor cells, where 2696 * the symmetries are not necessary for FE assembly. So for now we assume this is the case and don't return 2697 * symmetries if tensorSpace != tensorCell */ 2698 if (spintdim && 0 < dim && dim < 3 && (lag->tensorSpace == lag->tensorCell)) { /* compute self symmetries */ 2699 PetscInt **cellSymperms; 2700 PetscScalar **cellSymflips; 2701 PetscInt ornt; 2702 PetscInt nCopies = Nc / lag->intNodeIndices->nodeVecDim; 2703 PetscInt nNodes = lag->intNodeIndices->nNodes; 2704 2705 lag->numSelfSym = 2 * numFaces; 2706 lag->selfSymOff = numFaces; 2707 PetscCall(PetscCalloc1(2 * numFaces, &cellSymperms)); 2708 PetscCall(PetscCalloc1(2 * numFaces, &cellSymflips)); 2709 /* we want to be able to index symmetries directly with the orientations, which range from [-numFaces,numFaces) */ 2710 symperms[0] = &cellSymperms[numFaces]; 2711 symflips[0] = &cellSymflips[numFaces]; 2712 PetscCheck(lag->intNodeIndices->nodeVecDim * nCopies == Nc, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Node indices incompatible with dofs"); 2713 PetscCheck(nNodes * nCopies == spintdim, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Node indices incompatible with dofs"); 2714 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 */ 2715 Mat symMat; 2716 PetscInt *perm; 2717 PetscScalar *flips; 2718 PetscInt i; 2719 2720 if (!ornt) continue; 2721 PetscCall(PetscMalloc1(spintdim, &perm)); 2722 PetscCall(PetscCalloc1(spintdim, &flips)); 2723 for (i = 0; i < spintdim; i++) perm[i] = -1; 2724 PetscCall(PetscDualSpaceCreateInteriorSymmetryMatrix_Lagrange(sp, ornt, &symMat)); 2725 for (i = 0; i < nNodes; i++) { 2726 PetscInt ncols; 2727 PetscInt j, k; 2728 const PetscInt *cols; 2729 const PetscScalar *vals; 2730 PetscBool nz_seen = PETSC_FALSE; 2731 2732 PetscCall(MatGetRow(symMat, i, &ncols, &cols, &vals)); 2733 for (j = 0; j < ncols; j++) { 2734 if (PetscAbsScalar(vals[j]) > PETSC_SMALL) { 2735 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"); 2736 nz_seen = PETSC_TRUE; 2737 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"); 2738 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"); 2739 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"); 2740 for (k = 0; k < nCopies; k++) perm[cols[j] * nCopies + k] = i * nCopies + k; 2741 if (PetscRealPart(vals[j]) < 0.) { 2742 for (k = 0; k < nCopies; k++) flips[i * nCopies + k] = -1.; 2743 } else { 2744 for (k = 0; k < nCopies; k++) flips[i * nCopies + k] = 1.; 2745 } 2746 } 2747 } 2748 PetscCall(MatRestoreRow(symMat, i, &ncols, &cols, &vals)); 2749 } 2750 PetscCall(MatDestroy(&symMat)); 2751 /* if there were no sign flips, keep NULL */ 2752 for (i = 0; i < spintdim; i++) 2753 if (flips[i] != 1.) break; 2754 if (i == spintdim) { 2755 PetscCall(PetscFree(flips)); 2756 flips = NULL; 2757 } 2758 /* if the permutation is identity, keep NULL */ 2759 for (i = 0; i < spintdim; i++) 2760 if (perm[i] != i) break; 2761 if (i == spintdim) { 2762 PetscCall(PetscFree(perm)); 2763 perm = NULL; 2764 } 2765 symperms[0][ornt] = perm; 2766 symflips[0][ornt] = flips; 2767 } 2768 /* if no orientations produced non-identity permutations, keep NULL */ 2769 for (ornt = -numFaces; ornt < numFaces; ornt++) 2770 if (symperms[0][ornt]) break; 2771 if (ornt == numFaces) { 2772 PetscCall(PetscFree(cellSymperms)); 2773 symperms[0] = NULL; 2774 } 2775 /* if no orientations produced sign flips, keep NULL */ 2776 for (ornt = -numFaces; ornt < numFaces; ornt++) 2777 if (symflips[0][ornt]) break; 2778 if (ornt == numFaces) { 2779 PetscCall(PetscFree(cellSymflips)); 2780 symflips[0] = NULL; 2781 } 2782 } 2783 PetscCall(PetscDualSpaceGetBoundarySymmetries_Internal(sp, symperms, symflips)); 2784 for (p = 0; p < pEnd; p++) 2785 if (symperms[p]) break; 2786 if (p == pEnd) { 2787 PetscCall(PetscFree(symperms)); 2788 symperms = NULL; 2789 } 2790 for (p = 0; p < pEnd; p++) 2791 if (symflips[p]) break; 2792 if (p == pEnd) { 2793 PetscCall(PetscFree(symflips)); 2794 symflips = NULL; 2795 } 2796 lag->symperms = symperms; 2797 lag->symflips = symflips; 2798 lag->symComputed = PETSC_TRUE; 2799 } 2800 if (perms) *perms = (const PetscInt ***)lag->symperms; 2801 if (flips) *flips = (const PetscScalar ***)lag->symflips; 2802 PetscFunctionReturn(PETSC_SUCCESS); 2803 } 2804 2805 static PetscErrorCode PetscDualSpaceLagrangeGetContinuity_Lagrange(PetscDualSpace sp, PetscBool *continuous) 2806 { 2807 PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data; 2808 2809 PetscFunctionBegin; 2810 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 2811 PetscAssertPointer(continuous, 2); 2812 *continuous = lag->continuous; 2813 PetscFunctionReturn(PETSC_SUCCESS); 2814 } 2815 2816 static PetscErrorCode PetscDualSpaceLagrangeSetContinuity_Lagrange(PetscDualSpace sp, PetscBool continuous) 2817 { 2818 PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data; 2819 2820 PetscFunctionBegin; 2821 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 2822 lag->continuous = continuous; 2823 PetscFunctionReturn(PETSC_SUCCESS); 2824 } 2825 2826 /*@ 2827 PetscDualSpaceLagrangeGetContinuity - Retrieves the flag for element continuity 2828 2829 Not Collective 2830 2831 Input Parameter: 2832 . sp - the `PetscDualSpace` 2833 2834 Output Parameter: 2835 . continuous - flag for element continuity 2836 2837 Level: intermediate 2838 2839 .seealso: `PETSCDUALSPACELAGRANGE`, `PetscDualSpace`, `PetscDualSpaceLagrangeSetContinuity()` 2840 @*/ 2841 PetscErrorCode PetscDualSpaceLagrangeGetContinuity(PetscDualSpace sp, PetscBool *continuous) 2842 { 2843 PetscFunctionBegin; 2844 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 2845 PetscAssertPointer(continuous, 2); 2846 PetscTryMethod(sp, "PetscDualSpaceLagrangeGetContinuity_C", (PetscDualSpace, PetscBool *), (sp, continuous)); 2847 PetscFunctionReturn(PETSC_SUCCESS); 2848 } 2849 2850 /*@ 2851 PetscDualSpaceLagrangeSetContinuity - Indicate whether the element is continuous 2852 2853 Logically Collective 2854 2855 Input Parameters: 2856 + sp - the `PetscDualSpace` 2857 - continuous - flag for element continuity 2858 2859 Options Database Key: 2860 . -petscdualspace_lagrange_continuity <bool> - use a continuous element 2861 2862 Level: intermediate 2863 2864 .seealso: `PETSCDUALSPACELAGRANGE`, `PetscDualSpace`, `PetscDualSpaceLagrangeGetContinuity()` 2865 @*/ 2866 PetscErrorCode PetscDualSpaceLagrangeSetContinuity(PetscDualSpace sp, PetscBool continuous) 2867 { 2868 PetscFunctionBegin; 2869 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 2870 PetscValidLogicalCollectiveBool(sp, continuous, 2); 2871 PetscTryMethod(sp, "PetscDualSpaceLagrangeSetContinuity_C", (PetscDualSpace, PetscBool), (sp, continuous)); 2872 PetscFunctionReturn(PETSC_SUCCESS); 2873 } 2874 2875 static PetscErrorCode PetscDualSpaceLagrangeGetTensor_Lagrange(PetscDualSpace sp, PetscBool *tensor) 2876 { 2877 PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data; 2878 2879 PetscFunctionBegin; 2880 *tensor = lag->tensorSpace; 2881 PetscFunctionReturn(PETSC_SUCCESS); 2882 } 2883 2884 static PetscErrorCode PetscDualSpaceLagrangeSetTensor_Lagrange(PetscDualSpace sp, PetscBool tensor) 2885 { 2886 PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data; 2887 2888 PetscFunctionBegin; 2889 lag->tensorSpace = tensor; 2890 PetscFunctionReturn(PETSC_SUCCESS); 2891 } 2892 2893 static PetscErrorCode PetscDualSpaceLagrangeGetTrimmed_Lagrange(PetscDualSpace sp, PetscBool *trimmed) 2894 { 2895 PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data; 2896 2897 PetscFunctionBegin; 2898 *trimmed = lag->trimmed; 2899 PetscFunctionReturn(PETSC_SUCCESS); 2900 } 2901 2902 static PetscErrorCode PetscDualSpaceLagrangeSetTrimmed_Lagrange(PetscDualSpace sp, PetscBool trimmed) 2903 { 2904 PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data; 2905 2906 PetscFunctionBegin; 2907 lag->trimmed = trimmed; 2908 PetscFunctionReturn(PETSC_SUCCESS); 2909 } 2910 2911 static PetscErrorCode PetscDualSpaceLagrangeGetNodeType_Lagrange(PetscDualSpace sp, PetscDTNodeType *nodeType, PetscBool *boundary, PetscReal *exponent) 2912 { 2913 PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data; 2914 2915 PetscFunctionBegin; 2916 if (nodeType) *nodeType = lag->nodeType; 2917 if (boundary) *boundary = lag->endNodes; 2918 if (exponent) *exponent = lag->nodeExponent; 2919 PetscFunctionReturn(PETSC_SUCCESS); 2920 } 2921 2922 static PetscErrorCode PetscDualSpaceLagrangeSetNodeType_Lagrange(PetscDualSpace sp, PetscDTNodeType nodeType, PetscBool boundary, PetscReal exponent) 2923 { 2924 PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data; 2925 2926 PetscFunctionBegin; 2927 PetscCheck(nodeType != PETSCDTNODES_GAUSSJACOBI || exponent > -1., PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_OUTOFRANGE, "Exponent must be > -1"); 2928 lag->nodeType = nodeType; 2929 lag->endNodes = boundary; 2930 lag->nodeExponent = exponent; 2931 PetscFunctionReturn(PETSC_SUCCESS); 2932 } 2933 2934 static PetscErrorCode PetscDualSpaceLagrangeGetUseMoments_Lagrange(PetscDualSpace sp, PetscBool *useMoments) 2935 { 2936 PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data; 2937 2938 PetscFunctionBegin; 2939 *useMoments = lag->useMoments; 2940 PetscFunctionReturn(PETSC_SUCCESS); 2941 } 2942 2943 static PetscErrorCode PetscDualSpaceLagrangeSetUseMoments_Lagrange(PetscDualSpace sp, PetscBool useMoments) 2944 { 2945 PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data; 2946 2947 PetscFunctionBegin; 2948 lag->useMoments = useMoments; 2949 PetscFunctionReturn(PETSC_SUCCESS); 2950 } 2951 2952 static PetscErrorCode PetscDualSpaceLagrangeGetMomentOrder_Lagrange(PetscDualSpace sp, PetscInt *momentOrder) 2953 { 2954 PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data; 2955 2956 PetscFunctionBegin; 2957 *momentOrder = lag->momentOrder; 2958 PetscFunctionReturn(PETSC_SUCCESS); 2959 } 2960 2961 static PetscErrorCode PetscDualSpaceLagrangeSetMomentOrder_Lagrange(PetscDualSpace sp, PetscInt momentOrder) 2962 { 2963 PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data; 2964 2965 PetscFunctionBegin; 2966 lag->momentOrder = momentOrder; 2967 PetscFunctionReturn(PETSC_SUCCESS); 2968 } 2969 2970 /*@ 2971 PetscDualSpaceLagrangeGetTensor - Get the tensor nature of the dual space 2972 2973 Not Collective 2974 2975 Input Parameter: 2976 . sp - The `PetscDualSpace` 2977 2978 Output Parameter: 2979 . tensor - Whether the dual space has tensor layout (vs. simplicial) 2980 2981 Level: intermediate 2982 2983 .seealso: `PETSCDUALSPACELAGRANGE`, `PetscDualSpace`, `PetscDualSpaceLagrangeSetTensor()`, `PetscDualSpaceCreate()` 2984 @*/ 2985 PetscErrorCode PetscDualSpaceLagrangeGetTensor(PetscDualSpace sp, PetscBool *tensor) 2986 { 2987 PetscFunctionBegin; 2988 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 2989 PetscAssertPointer(tensor, 2); 2990 PetscTryMethod(sp, "PetscDualSpaceLagrangeGetTensor_C", (PetscDualSpace, PetscBool *), (sp, tensor)); 2991 PetscFunctionReturn(PETSC_SUCCESS); 2992 } 2993 2994 /*@ 2995 PetscDualSpaceLagrangeSetTensor - Set the tensor nature of the dual space 2996 2997 Not Collective 2998 2999 Input Parameters: 3000 + sp - The `PetscDualSpace` 3001 - tensor - Whether the dual space has tensor layout (vs. simplicial) 3002 3003 Level: intermediate 3004 3005 .seealso: `PETSCDUALSPACELAGRANGE`, `PetscDualSpace`, `PetscDualSpaceLagrangeGetTensor()`, `PetscDualSpaceCreate()` 3006 @*/ 3007 PetscErrorCode PetscDualSpaceLagrangeSetTensor(PetscDualSpace sp, PetscBool tensor) 3008 { 3009 PetscFunctionBegin; 3010 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 3011 PetscTryMethod(sp, "PetscDualSpaceLagrangeSetTensor_C", (PetscDualSpace, PetscBool), (sp, tensor)); 3012 PetscFunctionReturn(PETSC_SUCCESS); 3013 } 3014 3015 /*@ 3016 PetscDualSpaceLagrangeGetTrimmed - Get the trimmed nature of the dual space 3017 3018 Not Collective 3019 3020 Input Parameter: 3021 . sp - The `PetscDualSpace` 3022 3023 Output Parameter: 3024 . 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) 3025 3026 Level: intermediate 3027 3028 .seealso: `PETSCDUALSPACELAGRANGE`, `PetscDualSpace`, `PetscDualSpaceLagrangeSetTrimmed()`, `PetscDualSpaceCreate()` 3029 @*/ 3030 PetscErrorCode PetscDualSpaceLagrangeGetTrimmed(PetscDualSpace sp, PetscBool *trimmed) 3031 { 3032 PetscFunctionBegin; 3033 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 3034 PetscAssertPointer(trimmed, 2); 3035 PetscTryMethod(sp, "PetscDualSpaceLagrangeGetTrimmed_C", (PetscDualSpace, PetscBool *), (sp, trimmed)); 3036 PetscFunctionReturn(PETSC_SUCCESS); 3037 } 3038 3039 /*@ 3040 PetscDualSpaceLagrangeSetTrimmed - Set the trimmed nature of the dual space 3041 3042 Not Collective 3043 3044 Input Parameters: 3045 + sp - The `PetscDualSpace` 3046 - 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) 3047 3048 Level: intermediate 3049 3050 .seealso: `PETSCDUALSPACELAGRANGE`, `PetscDualSpace`, `PetscDualSpaceLagrangeGetTrimmed()`, `PetscDualSpaceCreate()` 3051 @*/ 3052 PetscErrorCode PetscDualSpaceLagrangeSetTrimmed(PetscDualSpace sp, PetscBool trimmed) 3053 { 3054 PetscFunctionBegin; 3055 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 3056 PetscTryMethod(sp, "PetscDualSpaceLagrangeSetTrimmed_C", (PetscDualSpace, PetscBool), (sp, trimmed)); 3057 PetscFunctionReturn(PETSC_SUCCESS); 3058 } 3059 3060 /*@ 3061 PetscDualSpaceLagrangeGetNodeType - Get a description of how nodes are laid out for Lagrange polynomials in this 3062 dual space 3063 3064 Not Collective 3065 3066 Input Parameter: 3067 . sp - The `PetscDualSpace` 3068 3069 Output Parameters: 3070 + nodeType - The type of nodes 3071 . boundary - Whether the node type is one that includes endpoints (if nodeType is `PETSCDTNODES_GAUSSJACOBI`, nodes that 3072 include the boundary are Gauss-Lobatto-Jacobi nodes) 3073 - exponent - If nodeType is `PETSCDTNODES_GAUSJACOBI`, indicates the exponent used for both ends of the 1D Jacobi weight function 3074 '0' is Gauss-Legendre, '-0.5' is Gauss-Chebyshev of the first type, '0.5' is Gauss-Chebyshev of the second type 3075 3076 Level: advanced 3077 3078 .seealso: `PETSCDUALSPACELAGRANGE`, `PetscDualSpace`, `PetscDTNodeType`, `PetscDualSpaceLagrangeSetNodeType()` 3079 @*/ 3080 PetscErrorCode PetscDualSpaceLagrangeGetNodeType(PetscDualSpace sp, PetscDTNodeType *nodeType, PetscBool *boundary, PetscReal *exponent) 3081 { 3082 PetscFunctionBegin; 3083 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 3084 if (nodeType) PetscAssertPointer(nodeType, 2); 3085 if (boundary) PetscAssertPointer(boundary, 3); 3086 if (exponent) PetscAssertPointer(exponent, 4); 3087 PetscTryMethod(sp, "PetscDualSpaceLagrangeGetNodeType_C", (PetscDualSpace, PetscDTNodeType *, PetscBool *, PetscReal *), (sp, nodeType, boundary, exponent)); 3088 PetscFunctionReturn(PETSC_SUCCESS); 3089 } 3090 3091 /*@ 3092 PetscDualSpaceLagrangeSetNodeType - Set a description of how nodes are laid out for Lagrange polynomials in this 3093 dual space 3094 3095 Logically Collective 3096 3097 Input Parameters: 3098 + sp - The `PetscDualSpace` 3099 . nodeType - The type of nodes 3100 . boundary - Whether the node type is one that includes endpoints (if nodeType is `PETSCDTNODES_GAUSSJACOBI`, nodes that 3101 include the boundary are Gauss-Lobatto-Jacobi nodes) 3102 - exponent - If nodeType is `PETSCDTNODES_GAUSJACOBI`, indicates the exponent used for both ends of the 1D Jacobi weight function 3103 '0' is Gauss-Legendre, '-0.5' is Gauss-Chebyshev of the first type, '0.5' is Gauss-Chebyshev of the second type 3104 3105 Level: advanced 3106 3107 .seealso: `PETSCDUALSPACELAGRANGE`, `PetscDualSpace`, `PetscDTNodeType`, `PetscDualSpaceLagrangeGetNodeType()` 3108 @*/ 3109 PetscErrorCode PetscDualSpaceLagrangeSetNodeType(PetscDualSpace sp, PetscDTNodeType nodeType, PetscBool boundary, PetscReal exponent) 3110 { 3111 PetscFunctionBegin; 3112 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 3113 PetscTryMethod(sp, "PetscDualSpaceLagrangeSetNodeType_C", (PetscDualSpace, PetscDTNodeType, PetscBool, PetscReal), (sp, nodeType, boundary, exponent)); 3114 PetscFunctionReturn(PETSC_SUCCESS); 3115 } 3116 3117 /*@ 3118 PetscDualSpaceLagrangeGetUseMoments - Get the flag for using moment functionals 3119 3120 Not Collective 3121 3122 Input Parameter: 3123 . sp - The `PetscDualSpace` 3124 3125 Output Parameter: 3126 . useMoments - Moment flag 3127 3128 Level: advanced 3129 3130 .seealso: `PETSCDUALSPACELAGRANGE`, `PetscDualSpace`, `PetscDualSpaceLagrangeSetUseMoments()` 3131 @*/ 3132 PetscErrorCode PetscDualSpaceLagrangeGetUseMoments(PetscDualSpace sp, PetscBool *useMoments) 3133 { 3134 PetscFunctionBegin; 3135 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 3136 PetscAssertPointer(useMoments, 2); 3137 PetscUseMethod(sp, "PetscDualSpaceLagrangeGetUseMoments_C", (PetscDualSpace, PetscBool *), (sp, useMoments)); 3138 PetscFunctionReturn(PETSC_SUCCESS); 3139 } 3140 3141 /*@ 3142 PetscDualSpaceLagrangeSetUseMoments - Set the flag for moment functionals 3143 3144 Logically Collective 3145 3146 Input Parameters: 3147 + sp - The `PetscDualSpace` 3148 - useMoments - The flag for moment functionals 3149 3150 Level: advanced 3151 3152 .seealso: `PETSCDUALSPACELAGRANGE`, `PetscDualSpace`, `PetscDualSpaceLagrangeGetUseMoments()` 3153 @*/ 3154 PetscErrorCode PetscDualSpaceLagrangeSetUseMoments(PetscDualSpace sp, PetscBool useMoments) 3155 { 3156 PetscFunctionBegin; 3157 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 3158 PetscTryMethod(sp, "PetscDualSpaceLagrangeSetUseMoments_C", (PetscDualSpace, PetscBool), (sp, useMoments)); 3159 PetscFunctionReturn(PETSC_SUCCESS); 3160 } 3161 3162 /*@ 3163 PetscDualSpaceLagrangeGetMomentOrder - Get the order for moment integration 3164 3165 Not Collective 3166 3167 Input Parameter: 3168 . sp - The `PetscDualSpace` 3169 3170 Output Parameter: 3171 . order - Moment integration order 3172 3173 Level: advanced 3174 3175 .seealso: `PETSCDUALSPACELAGRANGE`, `PetscDualSpace`, `PetscDualSpaceLagrangeSetMomentOrder()` 3176 @*/ 3177 PetscErrorCode PetscDualSpaceLagrangeGetMomentOrder(PetscDualSpace sp, PetscInt *order) 3178 { 3179 PetscFunctionBegin; 3180 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 3181 PetscAssertPointer(order, 2); 3182 PetscUseMethod(sp, "PetscDualSpaceLagrangeGetMomentOrder_C", (PetscDualSpace, PetscInt *), (sp, order)); 3183 PetscFunctionReturn(PETSC_SUCCESS); 3184 } 3185 3186 /*@ 3187 PetscDualSpaceLagrangeSetMomentOrder - Set the order for moment integration 3188 3189 Logically Collective 3190 3191 Input Parameters: 3192 + sp - The `PetscDualSpace` 3193 - order - The order for moment integration 3194 3195 Level: advanced 3196 3197 .seealso: `PETSCDUALSPACELAGRANGE`, `PetscDualSpace`, `PetscDualSpaceLagrangeGetMomentOrder()` 3198 @*/ 3199 PetscErrorCode PetscDualSpaceLagrangeSetMomentOrder(PetscDualSpace sp, PetscInt order) 3200 { 3201 PetscFunctionBegin; 3202 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 3203 PetscTryMethod(sp, "PetscDualSpaceLagrangeSetMomentOrder_C", (PetscDualSpace, PetscInt), (sp, order)); 3204 PetscFunctionReturn(PETSC_SUCCESS); 3205 } 3206 3207 static PetscErrorCode PetscDualSpaceInitialize_Lagrange(PetscDualSpace sp) 3208 { 3209 PetscFunctionBegin; 3210 sp->ops->destroy = PetscDualSpaceDestroy_Lagrange; 3211 sp->ops->view = PetscDualSpaceView_Lagrange; 3212 sp->ops->setfromoptions = PetscDualSpaceSetFromOptions_Lagrange; 3213 sp->ops->duplicate = PetscDualSpaceDuplicate_Lagrange; 3214 sp->ops->setup = PetscDualSpaceSetUp_Lagrange; 3215 sp->ops->createheightsubspace = NULL; 3216 sp->ops->createpointsubspace = NULL; 3217 sp->ops->getsymmetries = PetscDualSpaceGetSymmetries_Lagrange; 3218 sp->ops->apply = PetscDualSpaceApplyDefault; 3219 sp->ops->applyall = PetscDualSpaceApplyAllDefault; 3220 sp->ops->applyint = PetscDualSpaceApplyInteriorDefault; 3221 sp->ops->createalldata = PetscDualSpaceCreateAllDataDefault; 3222 sp->ops->createintdata = PetscDualSpaceCreateInteriorDataDefault; 3223 PetscFunctionReturn(PETSC_SUCCESS); 3224 } 3225 3226 /*MC 3227 PETSCDUALSPACELAGRANGE = "lagrange" - A `PetscDualSpaceType` that encapsulates a dual space of pointwise evaluation functionals 3228 3229 Level: intermediate 3230 3231 Developer Note: 3232 This `PetscDualSpace` seems to manage directly trimmed and untrimmed polynomials as well as tensor and non-tensor polynomials while for `PetscSpace` there seems to 3233 be different `PetscSpaceType` for them. 3234 3235 .seealso: `PetscDualSpace`, `PetscDualSpaceType`, `PetscDualSpaceCreate()`, `PetscDualSpaceSetType()`, 3236 `PetscDualSpaceLagrangeSetMomentOrder()`, `PetscDualSpaceLagrangeGetMomentOrder()`, `PetscDualSpaceLagrangeSetUseMoments()`, `PetscDualSpaceLagrangeGetUseMoments()`, 3237 `PetscDualSpaceLagrangeSetNodeType, PetscDualSpaceLagrangeGetNodeType, PetscDualSpaceLagrangeGetContinuity, PetscDualSpaceLagrangeSetContinuity, 3238 `PetscDualSpaceLagrangeGetTensor()`, `PetscDualSpaceLagrangeSetTensor()`, `PetscDualSpaceLagrangeGetTrimmed()`, `PetscDualSpaceLagrangeSetTrimmed()` 3239 M*/ 3240 PETSC_EXTERN PetscErrorCode PetscDualSpaceCreate_Lagrange(PetscDualSpace sp) 3241 { 3242 PetscDualSpace_Lag *lag; 3243 3244 PetscFunctionBegin; 3245 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 3246 PetscCall(PetscNew(&lag)); 3247 sp->data = lag; 3248 3249 lag->tensorCell = PETSC_FALSE; 3250 lag->tensorSpace = PETSC_FALSE; 3251 lag->continuous = PETSC_TRUE; 3252 lag->numCopies = PETSC_DEFAULT; 3253 lag->numNodeSkip = PETSC_DEFAULT; 3254 lag->nodeType = PETSCDTNODES_DEFAULT; 3255 lag->useMoments = PETSC_FALSE; 3256 lag->momentOrder = 0; 3257 3258 PetscCall(PetscDualSpaceInitialize_Lagrange(sp)); 3259 PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetContinuity_C", PetscDualSpaceLagrangeGetContinuity_Lagrange)); 3260 PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetContinuity_C", PetscDualSpaceLagrangeSetContinuity_Lagrange)); 3261 PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetTensor_C", PetscDualSpaceLagrangeGetTensor_Lagrange)); 3262 PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetTensor_C", PetscDualSpaceLagrangeSetTensor_Lagrange)); 3263 PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetTrimmed_C", PetscDualSpaceLagrangeGetTrimmed_Lagrange)); 3264 PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetTrimmed_C", PetscDualSpaceLagrangeSetTrimmed_Lagrange)); 3265 PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetNodeType_C", PetscDualSpaceLagrangeGetNodeType_Lagrange)); 3266 PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetNodeType_C", PetscDualSpaceLagrangeSetNodeType_Lagrange)); 3267 PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetUseMoments_C", PetscDualSpaceLagrangeGetUseMoments_Lagrange)); 3268 PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetUseMoments_C", PetscDualSpaceLagrangeSetUseMoments_Lagrange)); 3269 PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetMomentOrder_C", PetscDualSpaceLagrangeGetMomentOrder_Lagrange)); 3270 PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetMomentOrder_C", PetscDualSpaceLagrangeSetMomentOrder_Lagrange)); 3271 PetscFunctionReturn(PETSC_SUCCESS); 3272 } 3273