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