1 /* Basis Jet Tabulation 2 3 We would like to tabulate the nodal basis functions and derivatives at a set of points, usually quadrature points. We 4 follow here the derviation in http://www.math.ttu.edu/~kirby/papers/fiat-toms-2004.pdf. The nodal basis $\psi_i$ can 5 be expressed in terms of a prime basis $\phi_i$ which can be stably evaluated. In PETSc, we will use the Legendre basis 6 as a prime basis. 7 8 \psi_i = \sum_k \alpha_{ki} \phi_k 9 10 Our nodal basis is defined in terms of the dual basis $n_j$ 11 12 n_j \cdot \psi_i = \delta_{ji} 13 14 and we may act on the first equation to obtain 15 16 n_j \cdot \psi_i = \sum_k \alpha_{ki} n_j \cdot \phi_k 17 \delta_{ji} = \sum_k \alpha_{ki} V_{jk} 18 I = V \alpha 19 20 so the coefficients of the nodal basis in the prime basis are 21 22 \alpha = V^{-1} 23 24 We will define the dual basis vectors $n_j$ using a quadrature rule. 25 26 Right now, we will just use the polynomial spaces P^k. I know some elements use the space of symmetric polynomials 27 (I think Nedelec), but we will neglect this for now. Constraints in the space, e.g. Arnold-Winther elements, can 28 be implemented exactly as in FIAT using functionals $L_j$. 29 30 I will have to count the degrees correctly for the Legendre product when we are on simplices. 31 32 We will have three objects: 33 - Space, P: this just need point evaluation I think 34 - Dual Space, P'+K: This looks like a set of functionals that can act on members of P, each n is defined by a Q 35 - FEM: This keeps {P, P', Q} 36 */ 37 #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/ 38 #include <petscdmplex.h> 39 40 PetscBool FEcite = PETSC_FALSE; 41 const char FECitation[] = "@article{kirby2004,\n" 42 " title = {Algorithm 839: FIAT, a New Paradigm for Computing Finite Element Basis Functions},\n" 43 " journal = {ACM Transactions on Mathematical Software},\n" 44 " author = {Robert C. Kirby},\n" 45 " volume = {30},\n" 46 " number = {4},\n" 47 " pages = {502--516},\n" 48 " doi = {10.1145/1039813.1039820},\n" 49 " year = {2004}\n}\n"; 50 51 PetscClassId PETSCFE_CLASSID = 0; 52 53 PetscLogEvent PETSCFE_SetUp; 54 55 PetscFunctionList PetscFEList = NULL; 56 PetscBool PetscFERegisterAllCalled = PETSC_FALSE; 57 58 /*@C 59 PetscFERegister - Adds a new `PetscFEType` 60 61 Not Collective, No Fortran Support 62 63 Input Parameters: 64 + sname - The name of a new user-defined creation routine 65 - function - The creation routine 66 67 Example Usage: 68 .vb 69 PetscFERegister("my_fe", MyPetscFECreate); 70 .ve 71 72 Then, your PetscFE type can be chosen with the procedural interface via 73 .vb 74 PetscFECreate(MPI_Comm, PetscFE *); 75 PetscFESetType(PetscFE, "my_fe"); 76 .ve 77 or at runtime via the option 78 .vb 79 -petscfe_type my_fe 80 .ve 81 82 Level: advanced 83 84 Note: 85 `PetscFERegister()` may be called multiple times to add several user-defined `PetscFE`s 86 87 .seealso: `PetscFE`, `PetscFEType`, `PetscFERegisterAll()`, `PetscFERegisterDestroy()` 88 @*/ 89 PetscErrorCode PetscFERegister(const char sname[], PetscErrorCode (*function)(PetscFE)) 90 { 91 PetscFunctionBegin; 92 PetscCall(PetscFunctionListAdd(&PetscFEList, sname, function)); 93 PetscFunctionReturn(PETSC_SUCCESS); 94 } 95 96 /*@ 97 PetscFESetType - Builds a particular `PetscFE` 98 99 Collective 100 101 Input Parameters: 102 + fem - The `PetscFE` object 103 - name - The kind of FEM space 104 105 Options Database Key: 106 . -petscfe_type <type> - Sets the `PetscFE` type; use -help for a list of available types 107 108 Level: intermediate 109 110 .seealso: `PetscFEType`, `PetscFE`, `PetscFEGetType()`, `PetscFECreate()` 111 @*/ 112 PetscErrorCode PetscFESetType(PetscFE fem, PetscFEType name) 113 { 114 PetscErrorCode (*r)(PetscFE); 115 PetscBool match; 116 117 PetscFunctionBegin; 118 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 119 PetscCall(PetscObjectTypeCompare((PetscObject)fem, name, &match)); 120 if (match) PetscFunctionReturn(PETSC_SUCCESS); 121 122 if (!PetscFERegisterAllCalled) PetscCall(PetscFERegisterAll()); 123 PetscCall(PetscFunctionListFind(PetscFEList, name, &r)); 124 PetscCheck(r, PetscObjectComm((PetscObject)fem), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscFE type: %s", name); 125 126 PetscTryTypeMethod(fem, destroy); 127 fem->ops->destroy = NULL; 128 129 PetscCall((*r)(fem)); 130 PetscCall(PetscObjectChangeTypeName((PetscObject)fem, name)); 131 PetscFunctionReturn(PETSC_SUCCESS); 132 } 133 134 /*@ 135 PetscFEGetType - Gets the `PetscFEType` (as a string) from the `PetscFE` object. 136 137 Not Collective 138 139 Input Parameter: 140 . fem - The `PetscFE` 141 142 Output Parameter: 143 . name - The `PetscFEType` name 144 145 Level: intermediate 146 147 .seealso: `PetscFEType`, `PetscFE`, `PetscFESetType()`, `PetscFECreate()` 148 @*/ 149 PetscErrorCode PetscFEGetType(PetscFE fem, PetscFEType *name) 150 { 151 PetscFunctionBegin; 152 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 153 PetscAssertPointer(name, 2); 154 if (!PetscFERegisterAllCalled) PetscCall(PetscFERegisterAll()); 155 *name = ((PetscObject)fem)->type_name; 156 PetscFunctionReturn(PETSC_SUCCESS); 157 } 158 159 /*@ 160 PetscFEViewFromOptions - View from a `PetscFE` based on values in the options database 161 162 Collective 163 164 Input Parameters: 165 + A - the `PetscFE` object 166 . obj - Optional object that provides the options prefix, pass `NULL` to use the options prefix of `A` 167 - name - command line option name 168 169 Level: intermediate 170 171 .seealso: `PetscFE`, `PetscFEView()`, `PetscObjectViewFromOptions()`, `PetscFECreate()` 172 @*/ 173 PetscErrorCode PetscFEViewFromOptions(PetscFE A, PeOp PetscObject obj, const char name[]) 174 { 175 PetscFunctionBegin; 176 PetscValidHeaderSpecific(A, PETSCFE_CLASSID, 1); 177 PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name)); 178 PetscFunctionReturn(PETSC_SUCCESS); 179 } 180 181 /*@ 182 PetscFEView - Views a `PetscFE` 183 184 Collective 185 186 Input Parameters: 187 + fem - the `PetscFE` object to view 188 - viewer - the viewer 189 190 Level: beginner 191 192 .seealso: `PetscFE`, `PetscViewer`, `PetscFEDestroy()`, `PetscFEViewFromOptions()` 193 @*/ 194 PetscErrorCode PetscFEView(PetscFE fem, PetscViewer viewer) 195 { 196 PetscBool isascii; 197 198 PetscFunctionBegin; 199 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 200 if (viewer) PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 201 if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)fem), &viewer)); 202 PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)fem, viewer)); 203 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 204 PetscTryTypeMethod(fem, view, viewer); 205 PetscFunctionReturn(PETSC_SUCCESS); 206 } 207 208 /*@ 209 PetscFESetFromOptions - sets parameters in a `PetscFE` from the options database 210 211 Collective 212 213 Input Parameter: 214 . fem - the `PetscFE` object to set options for 215 216 Options Database Keys: 217 + -petscfe_num_blocks - the number of cell blocks to integrate concurrently 218 - -petscfe_num_batches - the number of cell batches to integrate serially 219 220 Level: intermediate 221 222 .seealso: `PetscFE`, `PetscFEView()` 223 @*/ 224 PetscErrorCode PetscFESetFromOptions(PetscFE fem) 225 { 226 const char *defaultType; 227 char name[256]; 228 PetscBool flg; 229 230 PetscFunctionBegin; 231 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 232 if (!((PetscObject)fem)->type_name) { 233 defaultType = PETSCFEBASIC; 234 } else { 235 defaultType = ((PetscObject)fem)->type_name; 236 } 237 if (!PetscFERegisterAllCalled) PetscCall(PetscFERegisterAll()); 238 239 PetscObjectOptionsBegin((PetscObject)fem); 240 PetscCall(PetscOptionsFList("-petscfe_type", "Finite element space", "PetscFESetType", PetscFEList, defaultType, name, 256, &flg)); 241 if (flg) { 242 PetscCall(PetscFESetType(fem, name)); 243 } else if (!((PetscObject)fem)->type_name) { 244 PetscCall(PetscFESetType(fem, defaultType)); 245 } 246 PetscCall(PetscOptionsBoundedInt("-petscfe_num_blocks", "The number of cell blocks to integrate concurrently", "PetscSpaceSetTileSizes", fem->numBlocks, &fem->numBlocks, NULL, 1)); 247 PetscCall(PetscOptionsBoundedInt("-petscfe_num_batches", "The number of cell batches to integrate serially", "PetscSpaceSetTileSizes", fem->numBatches, &fem->numBatches, NULL, 1)); 248 PetscTryTypeMethod(fem, setfromoptions, PetscOptionsObject); 249 /* process any options handlers added with PetscObjectAddOptionsHandler() */ 250 PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)fem, PetscOptionsObject)); 251 PetscOptionsEnd(); 252 PetscCall(PetscFEViewFromOptions(fem, NULL, "-petscfe_view")); 253 PetscFunctionReturn(PETSC_SUCCESS); 254 } 255 256 /*@ 257 PetscFESetUp - Construct data structures for the `PetscFE` after the `PetscFEType` has been set 258 259 Collective 260 261 Input Parameter: 262 . fem - the `PetscFE` object to setup 263 264 Level: intermediate 265 266 .seealso: `PetscFE`, `PetscFEView()`, `PetscFEDestroy()` 267 @*/ 268 PetscErrorCode PetscFESetUp(PetscFE fem) 269 { 270 PetscFunctionBegin; 271 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 272 if (fem->setupcalled) PetscFunctionReturn(PETSC_SUCCESS); 273 PetscCall(PetscLogEventBegin(PETSCFE_SetUp, fem, 0, 0, 0)); 274 fem->setupcalled = PETSC_TRUE; 275 PetscTryTypeMethod(fem, setup); 276 PetscCall(PetscLogEventEnd(PETSCFE_SetUp, fem, 0, 0, 0)); 277 PetscFunctionReturn(PETSC_SUCCESS); 278 } 279 280 /*@ 281 PetscFEDestroy - Destroys a `PetscFE` object 282 283 Collective 284 285 Input Parameter: 286 . fem - the `PetscFE` object to destroy 287 288 Level: beginner 289 290 .seealso: `PetscFE`, `PetscFEView()` 291 @*/ 292 PetscErrorCode PetscFEDestroy(PetscFE *fem) 293 { 294 PetscFunctionBegin; 295 if (!*fem) PetscFunctionReturn(PETSC_SUCCESS); 296 PetscValidHeaderSpecific(*fem, PETSCFE_CLASSID, 1); 297 298 if (--((PetscObject)*fem)->refct > 0) { 299 *fem = NULL; 300 PetscFunctionReturn(PETSC_SUCCESS); 301 } 302 ((PetscObject)*fem)->refct = 0; 303 304 if ((*fem)->subspaces) { 305 PetscInt dim, d; 306 307 PetscCall(PetscDualSpaceGetDimension((*fem)->dualSpace, &dim)); 308 for (d = 0; d < dim; ++d) PetscCall(PetscFEDestroy(&(*fem)->subspaces[d])); 309 } 310 PetscCall(PetscFree((*fem)->subspaces)); 311 PetscCall(PetscFree((*fem)->invV)); 312 PetscCall(PetscTabulationDestroy(&(*fem)->T)); 313 PetscCall(PetscTabulationDestroy(&(*fem)->Tf)); 314 PetscCall(PetscTabulationDestroy(&(*fem)->Tc)); 315 PetscCall(PetscSpaceDestroy(&(*fem)->basisSpace)); 316 PetscCall(PetscDualSpaceDestroy(&(*fem)->dualSpace)); 317 PetscCall(PetscQuadratureDestroy(&(*fem)->quadrature)); 318 PetscCall(PetscQuadratureDestroy(&(*fem)->faceQuadrature)); 319 #ifdef PETSC_HAVE_LIBCEED 320 PetscCallCEED(CeedBasisDestroy(&(*fem)->ceedBasis)); 321 PetscCallCEED(CeedDestroy(&(*fem)->ceed)); 322 #endif 323 324 PetscTryTypeMethod(*fem, destroy); 325 PetscCall(PetscHeaderDestroy(fem)); 326 PetscFunctionReturn(PETSC_SUCCESS); 327 } 328 329 /*@ 330 PetscFECreate - Creates an empty `PetscFE` object. The type can then be set with `PetscFESetType()`. 331 332 Collective 333 334 Input Parameter: 335 . comm - The communicator for the `PetscFE` object 336 337 Output Parameter: 338 . fem - The `PetscFE` object 339 340 Level: beginner 341 342 .seealso: `PetscFE`, `PetscFEType`, `PetscFESetType()`, `PetscFECreateDefault()`, `PETSCFEGALERKIN` 343 @*/ 344 PetscErrorCode PetscFECreate(MPI_Comm comm, PetscFE *fem) 345 { 346 PetscFE f; 347 348 PetscFunctionBegin; 349 PetscAssertPointer(fem, 2); 350 PetscCall(PetscCitationsRegister(FECitation, &FEcite)); 351 PetscCall(PetscFEInitializePackage()); 352 353 PetscCall(PetscHeaderCreate(f, PETSCFE_CLASSID, "PetscFE", "Finite Element", "PetscFE", comm, PetscFEDestroy, PetscFEView)); 354 355 f->basisSpace = NULL; 356 f->dualSpace = NULL; 357 f->numComponents = 1; 358 f->subspaces = NULL; 359 f->invV = NULL; 360 f->T = NULL; 361 f->Tf = NULL; 362 f->Tc = NULL; 363 PetscCall(PetscArrayzero(&f->quadrature, 1)); 364 PetscCall(PetscArrayzero(&f->faceQuadrature, 1)); 365 f->blockSize = 0; 366 f->numBlocks = 1; 367 f->batchSize = 0; 368 f->numBatches = 1; 369 370 *fem = f; 371 PetscFunctionReturn(PETSC_SUCCESS); 372 } 373 374 /*@ 375 PetscFEGetSpatialDimension - Returns the spatial dimension of the element 376 377 Not Collective 378 379 Input Parameter: 380 . fem - The `PetscFE` object 381 382 Output Parameter: 383 . dim - The spatial dimension 384 385 Level: intermediate 386 387 .seealso: `PetscFE`, `PetscFECreate()` 388 @*/ 389 PetscErrorCode PetscFEGetSpatialDimension(PetscFE fem, PetscInt *dim) 390 { 391 DM dm; 392 393 PetscFunctionBegin; 394 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 395 PetscAssertPointer(dim, 2); 396 PetscCall(PetscDualSpaceGetDM(fem->dualSpace, &dm)); 397 PetscCall(DMGetDimension(dm, dim)); 398 PetscFunctionReturn(PETSC_SUCCESS); 399 } 400 401 /*@ 402 PetscFESetNumComponents - Sets the number of field components in the element 403 404 Not Collective 405 406 Input Parameters: 407 + fem - The `PetscFE` object 408 - comp - The number of field components 409 410 Level: intermediate 411 412 .seealso: `PetscFE`, `PetscFECreate()`, `PetscFEGetSpatialDimension()`, `PetscFEGetNumComponents()` 413 @*/ 414 PetscErrorCode PetscFESetNumComponents(PetscFE fem, PetscInt comp) 415 { 416 PetscFunctionBegin; 417 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 418 fem->numComponents = comp; 419 PetscFunctionReturn(PETSC_SUCCESS); 420 } 421 422 /*@ 423 PetscFEGetNumComponents - Returns the number of components in the element 424 425 Not Collective 426 427 Input Parameter: 428 . fem - The `PetscFE` object 429 430 Output Parameter: 431 . comp - The number of field components 432 433 Level: intermediate 434 435 .seealso: `PetscFE`, `PetscFECreate()`, `PetscFEGetSpatialDimension()` 436 @*/ 437 PetscErrorCode PetscFEGetNumComponents(PetscFE fem, PetscInt *comp) 438 { 439 PetscFunctionBegin; 440 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 441 PetscAssertPointer(comp, 2); 442 *comp = fem->numComponents; 443 PetscFunctionReturn(PETSC_SUCCESS); 444 } 445 446 /*@ 447 PetscFESetTileSizes - Sets the tile sizes for evaluation 448 449 Not Collective 450 451 Input Parameters: 452 + fem - The `PetscFE` object 453 . blockSize - The number of elements in a block 454 . numBlocks - The number of blocks in a batch 455 . batchSize - The number of elements in a batch 456 - numBatches - The number of batches in a chunk 457 458 Level: intermediate 459 460 .seealso: `PetscFE`, `PetscFECreate()`, `PetscFEGetTileSizes()` 461 @*/ 462 PetscErrorCode PetscFESetTileSizes(PetscFE fem, PetscInt blockSize, PetscInt numBlocks, PetscInt batchSize, PetscInt numBatches) 463 { 464 PetscFunctionBegin; 465 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 466 fem->blockSize = blockSize; 467 fem->numBlocks = numBlocks; 468 fem->batchSize = batchSize; 469 fem->numBatches = numBatches; 470 PetscFunctionReturn(PETSC_SUCCESS); 471 } 472 473 /*@ 474 PetscFEGetTileSizes - Returns the tile sizes for evaluation 475 476 Not Collective 477 478 Input Parameter: 479 . fem - The `PetscFE` object 480 481 Output Parameters: 482 + blockSize - The number of elements in a block, pass `NULL` if not needed 483 . numBlocks - The number of blocks in a batch, pass `NULL` if not needed 484 . batchSize - The number of elements in a batch, pass `NULL` if not needed 485 - numBatches - The number of batches in a chunk, pass `NULL` if not needed 486 487 Level: intermediate 488 489 .seealso: `PetscFE`, `PetscFECreate()`, `PetscFESetTileSizes()` 490 @*/ 491 PetscErrorCode PetscFEGetTileSizes(PetscFE fem, PeOp PetscInt *blockSize, PeOp PetscInt *numBlocks, PeOp PetscInt *batchSize, PeOp PetscInt *numBatches) 492 { 493 PetscFunctionBegin; 494 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 495 if (blockSize) PetscAssertPointer(blockSize, 2); 496 if (numBlocks) PetscAssertPointer(numBlocks, 3); 497 if (batchSize) PetscAssertPointer(batchSize, 4); 498 if (numBatches) PetscAssertPointer(numBatches, 5); 499 if (blockSize) *blockSize = fem->blockSize; 500 if (numBlocks) *numBlocks = fem->numBlocks; 501 if (batchSize) *batchSize = fem->batchSize; 502 if (numBatches) *numBatches = fem->numBatches; 503 PetscFunctionReturn(PETSC_SUCCESS); 504 } 505 506 /*@ 507 PetscFEGetBasisSpace - Returns the `PetscSpace` used for the approximation of the solution for the `PetscFE` 508 509 Not Collective 510 511 Input Parameter: 512 . fem - The `PetscFE` object 513 514 Output Parameter: 515 . sp - The `PetscSpace` object 516 517 Level: intermediate 518 519 .seealso: `PetscFE`, `PetscSpace`, `PetscFECreate()` 520 @*/ 521 PetscErrorCode PetscFEGetBasisSpace(PetscFE fem, PetscSpace *sp) 522 { 523 PetscFunctionBegin; 524 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 525 PetscAssertPointer(sp, 2); 526 *sp = fem->basisSpace; 527 PetscFunctionReturn(PETSC_SUCCESS); 528 } 529 530 /*@ 531 PetscFESetBasisSpace - Sets the `PetscSpace` used for the approximation of the solution 532 533 Not Collective 534 535 Input Parameters: 536 + fem - The `PetscFE` object 537 - sp - The `PetscSpace` object 538 539 Level: intermediate 540 541 Developer Notes: 542 There is `PetscFESetBasisSpace()` but the `PetscFESetDualSpace()`, likely the Basis is unneeded in the function name 543 544 .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscFECreate()`, `PetscFESetDualSpace()` 545 @*/ 546 PetscErrorCode PetscFESetBasisSpace(PetscFE fem, PetscSpace sp) 547 { 548 PetscFunctionBegin; 549 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 550 PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 2); 551 PetscCall(PetscSpaceDestroy(&fem->basisSpace)); 552 fem->basisSpace = sp; 553 PetscCall(PetscObjectReference((PetscObject)fem->basisSpace)); 554 PetscFunctionReturn(PETSC_SUCCESS); 555 } 556 557 /*@ 558 PetscFEGetDualSpace - Returns the `PetscDualSpace` used to define the inner product for a `PetscFE` 559 560 Not Collective 561 562 Input Parameter: 563 . fem - The `PetscFE` object 564 565 Output Parameter: 566 . sp - The `PetscDualSpace` object 567 568 Level: intermediate 569 570 .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscFECreate()` 571 @*/ 572 PetscErrorCode PetscFEGetDualSpace(PetscFE fem, PetscDualSpace *sp) 573 { 574 PetscFunctionBegin; 575 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 576 PetscAssertPointer(sp, 2); 577 *sp = fem->dualSpace; 578 PetscFunctionReturn(PETSC_SUCCESS); 579 } 580 581 /*@ 582 PetscFESetDualSpace - Sets the `PetscDualSpace` used to define the inner product 583 584 Not Collective 585 586 Input Parameters: 587 + fem - The `PetscFE` object 588 - sp - The `PetscDualSpace` object 589 590 Level: intermediate 591 592 .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscFECreate()`, `PetscFESetBasisSpace()` 593 @*/ 594 PetscErrorCode PetscFESetDualSpace(PetscFE fem, PetscDualSpace sp) 595 { 596 PetscFunctionBegin; 597 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 598 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 2); 599 PetscCall(PetscDualSpaceDestroy(&fem->dualSpace)); 600 fem->dualSpace = sp; 601 PetscCall(PetscObjectReference((PetscObject)fem->dualSpace)); 602 PetscFunctionReturn(PETSC_SUCCESS); 603 } 604 605 /*@ 606 PetscFEGetQuadrature - Returns the `PetscQuadrature` used to calculate inner products 607 608 Not Collective 609 610 Input Parameter: 611 . fem - The `PetscFE` object 612 613 Output Parameter: 614 . q - The `PetscQuadrature` object 615 616 Level: intermediate 617 618 .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscQuadrature`, `PetscFECreate()` 619 @*/ 620 PetscErrorCode PetscFEGetQuadrature(PetscFE fem, PetscQuadrature *q) 621 { 622 PetscFunctionBegin; 623 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 624 PetscAssertPointer(q, 2); 625 *q = fem->quadrature; 626 PetscFunctionReturn(PETSC_SUCCESS); 627 } 628 629 /*@ 630 PetscFESetQuadrature - Sets the `PetscQuadrature` used to calculate inner products 631 632 Not Collective 633 634 Input Parameters: 635 + fem - The `PetscFE` object 636 - q - The `PetscQuadrature` object 637 638 Level: intermediate 639 640 .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscQuadrature`, `PetscFECreate()`, `PetscFEGetFaceQuadrature()` 641 @*/ 642 PetscErrorCode PetscFESetQuadrature(PetscFE fem, PetscQuadrature q) 643 { 644 PetscInt Nc, qNc; 645 646 PetscFunctionBegin; 647 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 648 if (q == fem->quadrature) PetscFunctionReturn(PETSC_SUCCESS); 649 PetscCall(PetscFEGetNumComponents(fem, &Nc)); 650 PetscCall(PetscQuadratureGetNumComponents(q, &qNc)); 651 PetscCheck(!(qNc != 1) || !(Nc != qNc), PetscObjectComm((PetscObject)fem), PETSC_ERR_ARG_SIZ, "FE components %" PetscInt_FMT " != Quadrature components %" PetscInt_FMT " and non-scalar quadrature", Nc, qNc); 652 PetscCall(PetscTabulationDestroy(&fem->T)); 653 PetscCall(PetscTabulationDestroy(&fem->Tc)); 654 PetscCall(PetscObjectReference((PetscObject)q)); 655 PetscCall(PetscQuadratureDestroy(&fem->quadrature)); 656 fem->quadrature = q; 657 PetscFunctionReturn(PETSC_SUCCESS); 658 } 659 660 /*@ 661 PetscFEGetFaceQuadrature - Returns the `PetscQuadrature` used to calculate inner products on faces 662 663 Not Collective 664 665 Input Parameter: 666 . fem - The `PetscFE` object 667 668 Output Parameter: 669 . q - The `PetscQuadrature` object 670 671 Level: intermediate 672 673 Developer Notes: 674 There is a special face quadrature but not edge, likely this API would benefit from a refactorization 675 676 .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscQuadrature`, `PetscFECreate()`, `PetscFESetQuadrature()`, `PetscFESetFaceQuadrature()` 677 @*/ 678 PetscErrorCode PetscFEGetFaceQuadrature(PetscFE fem, PetscQuadrature *q) 679 { 680 PetscFunctionBegin; 681 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 682 PetscAssertPointer(q, 2); 683 *q = fem->faceQuadrature; 684 PetscFunctionReturn(PETSC_SUCCESS); 685 } 686 687 /*@ 688 PetscFESetFaceQuadrature - Sets the `PetscQuadrature` used to calculate inner products on faces 689 690 Not Collective 691 692 Input Parameters: 693 + fem - The `PetscFE` object 694 - q - The `PetscQuadrature` object 695 696 Level: intermediate 697 698 .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscQuadrature`, `PetscFECreate()`, `PetscFESetQuadrature()` 699 @*/ 700 PetscErrorCode PetscFESetFaceQuadrature(PetscFE fem, PetscQuadrature q) 701 { 702 PetscInt Nc, qNc; 703 704 PetscFunctionBegin; 705 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 706 if (q == fem->faceQuadrature) PetscFunctionReturn(PETSC_SUCCESS); 707 PetscCall(PetscFEGetNumComponents(fem, &Nc)); 708 PetscCall(PetscQuadratureGetNumComponents(q, &qNc)); 709 PetscCheck(!(qNc != 1) || !(Nc != qNc), PetscObjectComm((PetscObject)fem), PETSC_ERR_ARG_SIZ, "FE components %" PetscInt_FMT " != Quadrature components %" PetscInt_FMT " and non-scalar quadrature", Nc, qNc); 710 PetscCall(PetscTabulationDestroy(&fem->Tf)); 711 PetscCall(PetscObjectReference((PetscObject)q)); 712 PetscCall(PetscQuadratureDestroy(&fem->faceQuadrature)); 713 fem->faceQuadrature = q; 714 PetscFunctionReturn(PETSC_SUCCESS); 715 } 716 717 /*@ 718 PetscFECopyQuadrature - Copy both volumetric and surface quadrature to a new `PetscFE` 719 720 Not Collective 721 722 Input Parameters: 723 + sfe - The `PetscFE` source for the quadratures 724 - tfe - The `PetscFE` target for the quadratures 725 726 Level: intermediate 727 728 .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscQuadrature`, `PetscFECreate()`, `PetscFESetQuadrature()`, `PetscFESetFaceQuadrature()` 729 @*/ 730 PetscErrorCode PetscFECopyQuadrature(PetscFE sfe, PetscFE tfe) 731 { 732 PetscQuadrature q; 733 734 PetscFunctionBegin; 735 PetscValidHeaderSpecific(sfe, PETSCFE_CLASSID, 1); 736 PetscValidHeaderSpecific(tfe, PETSCFE_CLASSID, 2); 737 PetscCall(PetscFEGetQuadrature(sfe, &q)); 738 PetscCall(PetscFESetQuadrature(tfe, q)); 739 PetscCall(PetscFEGetFaceQuadrature(sfe, &q)); 740 PetscCall(PetscFESetFaceQuadrature(tfe, q)); 741 PetscFunctionReturn(PETSC_SUCCESS); 742 } 743 744 /*@C 745 PetscFEGetNumDof - Returns the number of dofs (dual basis vectors) associated to mesh points on the reference cell of a given dimension 746 747 Not Collective 748 749 Input Parameter: 750 . fem - The `PetscFE` object 751 752 Output Parameter: 753 . numDof - Array of length `dim` with the number of dofs in each dimension 754 755 Level: intermediate 756 757 .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscFECreate()` 758 @*/ 759 PetscErrorCode PetscFEGetNumDof(PetscFE fem, const PetscInt *numDof[]) 760 { 761 PetscFunctionBegin; 762 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 763 PetscAssertPointer(numDof, 2); 764 PetscCall(PetscDualSpaceGetNumDof(fem->dualSpace, numDof)); 765 PetscFunctionReturn(PETSC_SUCCESS); 766 } 767 768 /*@C 769 PetscFEGetCellTabulation - Returns the tabulation of the basis functions at the quadrature points on the reference cell 770 771 Not Collective 772 773 Input Parameters: 774 + fem - The `PetscFE` object 775 - k - The highest derivative we need to tabulate, very often 1 776 777 Output Parameter: 778 . T - The basis function values and derivatives at quadrature points 779 780 Level: intermediate 781 782 Note: 783 .vb 784 T->T[0] = B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c 785 T->T[1] = D[((p*pdim + i)*Nc + c)*dim + d] is the derivative value at point p for basis function i, component c, in direction d 786 T->T[2] = H[(((p*pdim + i)*Nc + c)*dim + d)*dim + e] is the value at point p for basis function i, component c, in directions d and e 787 .ve 788 789 .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscTabulation`, `PetscFECreateTabulation()`, `PetscTabulationDestroy()` 790 @*/ 791 PetscErrorCode PetscFEGetCellTabulation(PetscFE fem, PetscInt k, PetscTabulation *T) 792 { 793 PetscInt npoints; 794 const PetscReal *points; 795 796 PetscFunctionBegin; 797 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 798 PetscAssertPointer(T, 3); 799 PetscCall(PetscQuadratureGetData(fem->quadrature, NULL, NULL, &npoints, &points, NULL)); 800 if (!fem->T) PetscCall(PetscFECreateTabulation(fem, 1, npoints, points, k, &fem->T)); 801 PetscCheck(!fem->T || k <= fem->T->K || (!fem->T->cdim && !fem->T->K), PetscObjectComm((PetscObject)fem), PETSC_ERR_ARG_OUTOFRANGE, "Requested %" PetscInt_FMT " derivatives, but only tabulated %" PetscInt_FMT, k, fem->T->K); 802 *T = fem->T; 803 PetscFunctionReturn(PETSC_SUCCESS); 804 } 805 806 PetscErrorCode PetscFEExpandFaceQuadrature(PetscFE fe, PetscQuadrature fq, PetscQuadrature *efq) 807 { 808 DM dm; 809 PetscDualSpace sp; 810 const PetscInt *faces; 811 const PetscReal *points, *weights; 812 DMPolytopeType ct; 813 PetscReal *facePoints, *faceWeights; 814 PetscInt dim, cStart, Nf, Nc, Np, order; 815 816 PetscFunctionBegin; 817 PetscCall(PetscFEGetDualSpace(fe, &sp)); 818 PetscCall(PetscDualSpaceGetDM(sp, &dm)); 819 PetscCall(DMGetDimension(dm, &dim)); 820 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL)); 821 PetscCall(DMPlexGetConeSize(dm, cStart, &Nf)); 822 PetscCall(DMPlexGetCone(dm, cStart, &faces)); 823 PetscCall(PetscQuadratureGetData(fq, NULL, &Nc, &Np, &points, &weights)); 824 PetscCall(PetscMalloc1(Nf * Np * dim, &facePoints)); 825 PetscCall(PetscMalloc1(Nf * Np * Nc, &faceWeights)); 826 for (PetscInt f = 0; f < Nf; ++f) { 827 const PetscReal xi0[3] = {-1., -1., -1.}; 828 PetscReal v0[3], J[9], detJ; 829 830 PetscCall(DMPlexComputeCellGeometryFEM(dm, faces[f], NULL, v0, J, NULL, &detJ)); 831 for (PetscInt q = 0; q < Np; ++q) { 832 CoordinatesRefToReal(dim, dim - 1, xi0, v0, J, &points[q * (dim - 1)], &facePoints[(f * Np + q) * dim]); 833 for (PetscInt c = 0; c < Nc; ++c) faceWeights[(f * Np + q) * Nc + c] = weights[q * Nc + c]; 834 } 835 } 836 PetscCall(PetscQuadratureCreate(PetscObjectComm((PetscObject)fq), efq)); 837 PetscCall(PetscQuadratureGetCellType(fq, &ct)); 838 PetscCall(PetscQuadratureSetCellType(*efq, ct)); 839 PetscCall(PetscQuadratureGetOrder(fq, &order)); 840 PetscCall(PetscQuadratureSetOrder(*efq, order)); 841 PetscCall(PetscQuadratureSetData(*efq, dim, Nc, Nf * Np, facePoints, faceWeights)); 842 PetscFunctionReturn(PETSC_SUCCESS); 843 } 844 845 /*@C 846 PetscFEGetFaceTabulation - Returns the tabulation of the basis functions at the face quadrature points for each face of the reference cell 847 848 Not Collective 849 850 Input Parameters: 851 + fem - The `PetscFE` object 852 - k - The highest derivative we need to tabulate, very often 1 853 854 Output Parameter: 855 . Tf - The basis function values and derivatives at face quadrature points 856 857 Level: intermediate 858 859 Note: 860 .vb 861 T->T[0] = Bf[((f*Nq + q)*pdim + i)*Nc + c] is the value at point f,q for basis function i and component c 862 T->T[1] = Df[(((f*Nq + q)*pdim + i)*Nc + c)*dim + d] is the derivative value at point f,q for basis function i, component c, in direction d 863 T->T[2] = Hf[((((f*Nq + q)*pdim + i)*Nc + c)*dim + d)*dim + e] is the value at point f,q for basis function i, component c, in directions d and e 864 .ve 865 866 .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscTabulation`, `PetscFEGetCellTabulation()`, `PetscFECreateTabulation()`, `PetscTabulationDestroy()` 867 @*/ 868 PetscErrorCode PetscFEGetFaceTabulation(PetscFE fem, PetscInt k, PetscTabulation *Tf) 869 { 870 PetscFunctionBegin; 871 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 872 PetscAssertPointer(Tf, 3); 873 if (!fem->Tf) { 874 PetscQuadrature fq; 875 876 PetscCall(PetscFEGetFaceQuadrature(fem, &fq)); 877 if (fq) { 878 PetscQuadrature efq; 879 const PetscReal *facePoints; 880 PetscInt Np, eNp; 881 882 PetscCall(PetscFEExpandFaceQuadrature(fem, fq, &efq)); 883 PetscCall(PetscQuadratureGetData(fq, NULL, NULL, &Np, NULL, NULL)); 884 PetscCall(PetscQuadratureGetData(efq, NULL, NULL, &eNp, &facePoints, NULL)); 885 if (PetscDefined(USE_DEBUG)) { 886 PetscDualSpace sp; 887 DM dm; 888 PetscInt cStart, Nf; 889 890 PetscCall(PetscFEGetDualSpace(fem, &sp)); 891 PetscCall(PetscDualSpaceGetDM(sp, &dm)); 892 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL)); 893 PetscCall(DMPlexGetConeSize(dm, cStart, &Nf)); 894 PetscCheck(Nf == eNp / Np, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of faces %" PetscInt_FMT " != %" PetscInt_FMT " number of quadrature replicas", Nf, eNp / Np); 895 } 896 PetscCall(PetscFECreateTabulation(fem, eNp / Np, Np, facePoints, k, &fem->Tf)); 897 PetscCall(PetscQuadratureDestroy(&efq)); 898 } 899 } 900 PetscCheck(!fem->Tf || k <= fem->Tf->K, PetscObjectComm((PetscObject)fem), PETSC_ERR_ARG_OUTOFRANGE, "Requested %" PetscInt_FMT " derivatives, but only tabulated %" PetscInt_FMT, k, fem->Tf->K); 901 *Tf = fem->Tf; 902 PetscFunctionReturn(PETSC_SUCCESS); 903 } 904 905 /*@C 906 PetscFEGetFaceCentroidTabulation - Returns the tabulation of the basis functions at the face centroid points 907 908 Not Collective 909 910 Input Parameter: 911 . fem - The `PetscFE` object 912 913 Output Parameter: 914 . Tc - The basis function values at face centroid points 915 916 Level: intermediate 917 918 Note: 919 .vb 920 T->T[0] = Bf[(f*pdim + i)*Nc + c] is the value at point f for basis function i and component c 921 .ve 922 923 .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscTabulation`, `PetscFEGetFaceTabulation()`, `PetscFEGetCellTabulation()`, `PetscFECreateTabulation()`, `PetscTabulationDestroy()` 924 @*/ 925 PetscErrorCode PetscFEGetFaceCentroidTabulation(PetscFE fem, PetscTabulation *Tc) 926 { 927 PetscFunctionBegin; 928 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 929 PetscAssertPointer(Tc, 2); 930 if (!fem->Tc) { 931 PetscDualSpace sp; 932 DM dm; 933 const PetscInt *cone; 934 PetscReal *centroids; 935 PetscInt dim, numFaces, f; 936 937 PetscCall(PetscFEGetDualSpace(fem, &sp)); 938 PetscCall(PetscDualSpaceGetDM(sp, &dm)); 939 PetscCall(DMGetDimension(dm, &dim)); 940 PetscCall(DMPlexGetConeSize(dm, 0, &numFaces)); 941 PetscCall(DMPlexGetCone(dm, 0, &cone)); 942 PetscCall(PetscMalloc1(numFaces * dim, ¢roids)); 943 for (f = 0; f < numFaces; ++f) PetscCall(DMPlexComputeCellGeometryFVM(dm, cone[f], NULL, ¢roids[f * dim], NULL)); 944 PetscCall(PetscFECreateTabulation(fem, 1, numFaces, centroids, 0, &fem->Tc)); 945 PetscCall(PetscFree(centroids)); 946 } 947 *Tc = fem->Tc; 948 PetscFunctionReturn(PETSC_SUCCESS); 949 } 950 951 /*@C 952 PetscFECreateTabulation - Tabulates the basis functions, and perhaps derivatives, at the points provided. 953 954 Not Collective 955 956 Input Parameters: 957 + fem - The `PetscFE` object 958 . nrepl - The number of replicas 959 . npoints - The number of tabulation points in a replica 960 . points - The tabulation point coordinates 961 - K - The number of derivatives calculated 962 963 Output Parameter: 964 . T - The basis function values and derivatives at tabulation points 965 966 Level: intermediate 967 968 Note: 969 .vb 970 T->T[0] = B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c 971 T->T[1] = D[((p*pdim + i)*Nc + c)*dim + d] is the derivative value at point p for basis function i, component c, in direction d 972 T->T[2] = H[(((p*pdim + i)*Nc + c)*dim + d)*dim + e] is the value at point p for basis 973 T->function i, component c, in directions d and e 974 .ve 975 976 .seealso: `PetscTabulation`, `PetscFEGetCellTabulation()`, `PetscTabulationDestroy()` 977 @*/ 978 PetscErrorCode PetscFECreateTabulation(PetscFE fem, PetscInt nrepl, PetscInt npoints, const PetscReal points[], PetscInt K, PetscTabulation *T) 979 { 980 DM dm; 981 PetscDualSpace Q; 982 PetscInt Nb; /* Dimension of FE space P */ 983 PetscInt Nc; /* Field components */ 984 PetscInt cdim; /* Reference coordinate dimension */ 985 PetscInt k; 986 987 PetscFunctionBegin; 988 if (!npoints || !fem->dualSpace || K < 0) { 989 *T = NULL; 990 PetscFunctionReturn(PETSC_SUCCESS); 991 } 992 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 993 PetscAssertPointer(points, 4); 994 PetscAssertPointer(T, 6); 995 PetscCall(PetscFEGetDualSpace(fem, &Q)); 996 PetscCall(PetscDualSpaceGetDM(Q, &dm)); 997 PetscCall(DMGetDimension(dm, &cdim)); 998 PetscCall(PetscDualSpaceGetDimension(Q, &Nb)); 999 PetscCall(PetscFEGetNumComponents(fem, &Nc)); 1000 PetscCall(PetscMalloc1(1, T)); 1001 (*T)->K = !cdim ? 0 : K; 1002 (*T)->Nr = nrepl; 1003 (*T)->Np = npoints; 1004 (*T)->Nb = Nb; 1005 (*T)->Nc = Nc; 1006 (*T)->cdim = cdim; 1007 PetscCall(PetscMalloc1((*T)->K + 1, &(*T)->T)); 1008 for (k = 0; k <= (*T)->K; ++k) PetscCall(PetscCalloc1(nrepl * npoints * Nb * Nc * PetscPowInt(cdim, k), &(*T)->T[k])); 1009 PetscUseTypeMethod(fem, computetabulation, nrepl * npoints, points, K, *T); 1010 PetscFunctionReturn(PETSC_SUCCESS); 1011 } 1012 1013 /*@C 1014 PetscFEComputeTabulation - Tabulates the basis functions, and perhaps derivatives, at the points provided. 1015 1016 Not Collective 1017 1018 Input Parameters: 1019 + fem - The `PetscFE` object 1020 . npoints - The number of tabulation points 1021 . points - The tabulation point coordinates 1022 . K - The number of derivatives calculated 1023 - T - An existing tabulation object with enough allocated space 1024 1025 Output Parameter: 1026 . T - The basis function values and derivatives at tabulation points 1027 1028 Level: intermediate 1029 1030 Note: 1031 .vb 1032 T->T[0] = B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c 1033 T->T[1] = D[((p*pdim + i)*Nc + c)*dim + d] is the derivative value at point p for basis function i, component c, in direction d 1034 T->T[2] = H[(((p*pdim + i)*Nc + c)*dim + d)*dim + e] is the value at point p for basis function i, component c, in directions d and e 1035 .ve 1036 1037 .seealso: `PetscTabulation`, `PetscFEGetCellTabulation()`, `PetscTabulationDestroy()` 1038 @*/ 1039 PetscErrorCode PetscFEComputeTabulation(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscInt K, PetscTabulation T) 1040 { 1041 PetscFunctionBeginHot; 1042 if (!npoints || !fem->dualSpace || K < 0) PetscFunctionReturn(PETSC_SUCCESS); 1043 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 1044 PetscAssertPointer(points, 3); 1045 PetscAssertPointer(T, 5); 1046 if (PetscDefined(USE_DEBUG)) { 1047 DM dm; 1048 PetscDualSpace Q; 1049 PetscInt Nb; /* Dimension of FE space P */ 1050 PetscInt Nc; /* Field components */ 1051 PetscInt cdim; /* Reference coordinate dimension */ 1052 1053 PetscCall(PetscFEGetDualSpace(fem, &Q)); 1054 PetscCall(PetscDualSpaceGetDM(Q, &dm)); 1055 PetscCall(DMGetDimension(dm, &cdim)); 1056 PetscCall(PetscDualSpaceGetDimension(Q, &Nb)); 1057 PetscCall(PetscFEGetNumComponents(fem, &Nc)); 1058 PetscCheck(T->K == (!cdim ? 0 : K), PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Tabulation K %" PetscInt_FMT " must match requested K %" PetscInt_FMT, T->K, !cdim ? 0 : K); 1059 PetscCheck(T->Nb == Nb, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Tabulation Nb %" PetscInt_FMT " must match requested Nb %" PetscInt_FMT, T->Nb, Nb); 1060 PetscCheck(T->Nc == Nc, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Tabulation Nc %" PetscInt_FMT " must match requested Nc %" PetscInt_FMT, T->Nc, Nc); 1061 PetscCheck(T->cdim == cdim, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Tabulation cdim %" PetscInt_FMT " must match requested cdim %" PetscInt_FMT, T->cdim, cdim); 1062 } 1063 T->Nr = 1; 1064 T->Np = npoints; 1065 PetscUseTypeMethod(fem, computetabulation, npoints, points, K, T); 1066 PetscFunctionReturn(PETSC_SUCCESS); 1067 } 1068 1069 /*@ 1070 PetscTabulationDestroy - Frees memory from the associated tabulation. 1071 1072 Not Collective 1073 1074 Input Parameter: 1075 . T - The tabulation 1076 1077 Level: intermediate 1078 1079 .seealso: `PetscTabulation`, `PetscFECreateTabulation()`, `PetscFEGetCellTabulation()` 1080 @*/ 1081 PetscErrorCode PetscTabulationDestroy(PetscTabulation *T) 1082 { 1083 PetscInt k; 1084 1085 PetscFunctionBegin; 1086 PetscAssertPointer(T, 1); 1087 if (!T || !*T) PetscFunctionReturn(PETSC_SUCCESS); 1088 for (k = 0; k <= (*T)->K; ++k) PetscCall(PetscFree((*T)->T[k])); 1089 PetscCall(PetscFree((*T)->T)); 1090 PetscCall(PetscFree(*T)); 1091 *T = NULL; 1092 PetscFunctionReturn(PETSC_SUCCESS); 1093 } 1094 1095 static PetscErrorCode PetscFECreatePointTraceDefault_Internal(PetscFE fe, PetscInt refPoint, PetscFE *trFE) 1096 { 1097 PetscSpace bsp, bsubsp; 1098 PetscDualSpace dsp, dsubsp; 1099 PetscInt dim, depth, numComp, i, j, coneSize, order; 1100 DM dm; 1101 DMLabel label; 1102 PetscReal *xi, *v, *J, detJ; 1103 const char *name; 1104 PetscQuadrature origin, fullQuad, subQuad; 1105 1106 PetscFunctionBegin; 1107 PetscCall(PetscFEGetBasisSpace(fe, &bsp)); 1108 PetscCall(PetscFEGetDualSpace(fe, &dsp)); 1109 PetscCall(PetscDualSpaceGetDM(dsp, &dm)); 1110 PetscCall(DMGetDimension(dm, &dim)); 1111 PetscCall(DMPlexGetDepthLabel(dm, &label)); 1112 PetscCall(DMLabelGetValue(label, refPoint, &depth)); 1113 PetscCall(PetscCalloc1(depth, &xi)); 1114 PetscCall(PetscMalloc1(dim, &v)); 1115 PetscCall(PetscMalloc1(dim * dim, &J)); 1116 for (i = 0; i < depth; i++) xi[i] = 0.; 1117 PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, &origin)); 1118 PetscCall(PetscQuadratureSetData(origin, depth, 0, 1, xi, NULL)); 1119 PetscCall(DMPlexComputeCellGeometryFEM(dm, refPoint, origin, v, J, NULL, &detJ)); 1120 /* CellGeometryFEM computes the expanded Jacobian, we want the true jacobian */ 1121 for (i = 1; i < dim; i++) { 1122 for (j = 0; j < depth; j++) J[i * depth + j] = J[i * dim + j]; 1123 } 1124 PetscCall(PetscQuadratureDestroy(&origin)); 1125 PetscCall(PetscDualSpaceGetPointSubspace(dsp, refPoint, &dsubsp)); 1126 PetscCall(PetscSpaceCreateSubspace(bsp, dsubsp, v, J, NULL, NULL, PETSC_OWN_POINTER, &bsubsp)); 1127 PetscCall(PetscSpaceSetUp(bsubsp)); 1128 PetscCall(PetscFECreate(PetscObjectComm((PetscObject)fe), trFE)); 1129 PetscCall(PetscFESetType(*trFE, PETSCFEBASIC)); 1130 PetscCall(PetscFEGetNumComponents(fe, &numComp)); 1131 PetscCall(PetscFESetNumComponents(*trFE, numComp)); 1132 PetscCall(PetscFESetBasisSpace(*trFE, bsubsp)); 1133 PetscCall(PetscFESetDualSpace(*trFE, dsubsp)); 1134 PetscCall(PetscObjectGetName((PetscObject)fe, &name)); 1135 if (name) PetscCall(PetscFESetName(*trFE, name)); 1136 PetscCall(PetscFEGetQuadrature(fe, &fullQuad)); 1137 PetscCall(PetscQuadratureGetOrder(fullQuad, &order)); 1138 PetscCall(DMPlexGetConeSize(dm, refPoint, &coneSize)); 1139 if (coneSize == 2 * depth) PetscCall(PetscDTGaussTensorQuadrature(depth, 1, (order + 2) / 2, -1., 1., &subQuad)); 1140 else PetscCall(PetscDTSimplexQuadrature(depth, order, PETSCDTSIMPLEXQUAD_DEFAULT, &subQuad)); 1141 PetscCall(PetscFESetQuadrature(*trFE, subQuad)); 1142 PetscCall(PetscFESetUp(*trFE)); 1143 PetscCall(PetscQuadratureDestroy(&subQuad)); 1144 PetscCall(PetscSpaceDestroy(&bsubsp)); 1145 PetscFunctionReturn(PETSC_SUCCESS); 1146 } 1147 1148 PETSC_EXTERN PetscErrorCode PetscFECreatePointTrace(PetscFE fe, PetscInt refPoint, PetscFE *trFE) 1149 { 1150 PetscFunctionBegin; 1151 PetscValidHeaderSpecific(fe, PETSCFE_CLASSID, 1); 1152 PetscAssertPointer(trFE, 3); 1153 if (fe->ops->createpointtrace) PetscUseTypeMethod(fe, createpointtrace, refPoint, trFE); 1154 else PetscCall(PetscFECreatePointTraceDefault_Internal(fe, refPoint, trFE)); 1155 PetscFunctionReturn(PETSC_SUCCESS); 1156 } 1157 1158 PetscErrorCode PetscFECreateHeightTrace(PetscFE fe, PetscInt height, PetscFE *trFE) 1159 { 1160 PetscInt hStart, hEnd; 1161 PetscDualSpace dsp; 1162 DM dm; 1163 1164 PetscFunctionBegin; 1165 PetscValidHeaderSpecific(fe, PETSCFE_CLASSID, 1); 1166 PetscAssertPointer(trFE, 3); 1167 *trFE = NULL; 1168 PetscCall(PetscFEGetDualSpace(fe, &dsp)); 1169 PetscCall(PetscDualSpaceGetDM(dsp, &dm)); 1170 PetscCall(DMPlexGetHeightStratum(dm, height, &hStart, &hEnd)); 1171 if (hEnd <= hStart) PetscFunctionReturn(PETSC_SUCCESS); 1172 PetscCall(PetscFECreatePointTrace(fe, hStart, trFE)); 1173 PetscFunctionReturn(PETSC_SUCCESS); 1174 } 1175 1176 /*@ 1177 PetscFEGetDimension - Get the dimension of the finite element space on a cell 1178 1179 Not Collective 1180 1181 Input Parameter: 1182 . fem - The `PetscFE` 1183 1184 Output Parameter: 1185 . dim - The dimension 1186 1187 Level: intermediate 1188 1189 .seealso: `PetscFE`, `PetscFECreate()`, `PetscSpaceGetDimension()`, `PetscDualSpaceGetDimension()` 1190 @*/ 1191 PetscErrorCode PetscFEGetDimension(PetscFE fem, PetscInt *dim) 1192 { 1193 PetscFunctionBegin; 1194 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 1195 PetscAssertPointer(dim, 2); 1196 PetscTryTypeMethod(fem, getdimension, dim); 1197 PetscFunctionReturn(PETSC_SUCCESS); 1198 } 1199 1200 /*@ 1201 PetscFEPushforward - Map the reference element function to real space 1202 1203 Input Parameters: 1204 + fe - The `PetscFE` 1205 . fegeom - The cell geometry 1206 . Nv - The number of function values 1207 - vals - The function values 1208 1209 Output Parameter: 1210 . vals - The transformed function values 1211 1212 Level: advanced 1213 1214 Notes: 1215 This just forwards the call onto `PetscDualSpacePushforward()`. 1216 1217 It only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension. 1218 1219 .seealso: `PetscFE`, `PetscFEGeom`, `PetscDualSpace`, `PetscDualSpacePushforward()` 1220 @*/ 1221 PetscErrorCode PetscFEPushforward(PetscFE fe, PetscFEGeom *fegeom, PetscInt Nv, PetscScalar vals[]) 1222 { 1223 PetscFunctionBeginHot; 1224 PetscCall(PetscDualSpacePushforward(fe->dualSpace, fegeom, Nv, fe->numComponents, vals)); 1225 PetscFunctionReturn(PETSC_SUCCESS); 1226 } 1227 1228 /*@ 1229 PetscFEPushforwardGradient - Map the reference element function gradient to real space 1230 1231 Input Parameters: 1232 + fe - The `PetscFE` 1233 . fegeom - The cell geometry 1234 . Nv - The number of function gradient values 1235 - vals - The function gradient values 1236 1237 Output Parameter: 1238 . vals - The transformed function gradient values 1239 1240 Level: advanced 1241 1242 Notes: 1243 This just forwards the call onto `PetscDualSpacePushforwardGradient()`. 1244 1245 It only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension. 1246 1247 .seealso: `PetscFE`, `PetscFEGeom`, `PetscDualSpace`, `PetscFEPushforward()`, `PetscDualSpacePushforwardGradient()`, `PetscDualSpacePushforward()` 1248 @*/ 1249 PetscErrorCode PetscFEPushforwardGradient(PetscFE fe, PetscFEGeom *fegeom, PetscInt Nv, PetscScalar vals[]) 1250 { 1251 PetscFunctionBeginHot; 1252 PetscCall(PetscDualSpacePushforwardGradient(fe->dualSpace, fegeom, Nv, fe->numComponents, vals)); 1253 PetscFunctionReturn(PETSC_SUCCESS); 1254 } 1255 1256 /*@ 1257 PetscFEPushforwardHessian - Map the reference element function Hessian to real space 1258 1259 Input Parameters: 1260 + fe - The `PetscFE` 1261 . fegeom - The cell geometry 1262 . Nv - The number of function Hessian values 1263 - vals - The function Hessian values 1264 1265 Output Parameter: 1266 . vals - The transformed function Hessian values 1267 1268 Level: advanced 1269 1270 Notes: 1271 This just forwards the call onto `PetscDualSpacePushforwardHessian()`. 1272 1273 It only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension. 1274 1275 Developer Notes: 1276 It is unclear why all these one line convenience routines are desirable 1277 1278 .seealso: `PetscFE`, `PetscFEGeom`, `PetscDualSpace`, `PetscFEPushforward()`, `PetscDualSpacePushforwardHessian()`, `PetscDualSpacePushforward()` 1279 @*/ 1280 PetscErrorCode PetscFEPushforwardHessian(PetscFE fe, PetscFEGeom *fegeom, PetscInt Nv, PetscScalar vals[]) 1281 { 1282 PetscFunctionBeginHot; 1283 PetscCall(PetscDualSpacePushforwardHessian(fe->dualSpace, fegeom, Nv, fe->numComponents, vals)); 1284 PetscFunctionReturn(PETSC_SUCCESS); 1285 } 1286 1287 /* 1288 Purpose: Compute element vector for chunk of elements 1289 1290 Input: 1291 Sizes: 1292 Ne: number of elements 1293 Nf: number of fields 1294 PetscFE 1295 dim: spatial dimension 1296 Nb: number of basis functions 1297 Nc: number of field components 1298 PetscQuadrature 1299 Nq: number of quadrature points 1300 1301 Geometry: 1302 PetscFEGeom[Ne] possibly *Nq 1303 PetscReal v0s[dim] 1304 PetscReal n[dim] 1305 PetscReal jacobians[dim*dim] 1306 PetscReal jacobianInverses[dim*dim] 1307 PetscReal jacobianDeterminants 1308 FEM: 1309 PetscFE 1310 PetscQuadrature 1311 PetscReal quadPoints[Nq*dim] 1312 PetscReal quadWeights[Nq] 1313 PetscReal basis[Nq*Nb*Nc] 1314 PetscReal basisDer[Nq*Nb*Nc*dim] 1315 PetscScalar coefficients[Ne*Nb*Nc] 1316 PetscScalar elemVec[Ne*Nb*Nc] 1317 1318 Problem: 1319 PetscInt f: the active field 1320 f0, f1 1321 1322 Work Space: 1323 PetscFE 1324 PetscScalar f0[Nq*dim]; 1325 PetscScalar f1[Nq*dim*dim]; 1326 PetscScalar u[Nc]; 1327 PetscScalar gradU[Nc*dim]; 1328 PetscReal x[dim]; 1329 PetscScalar realSpaceDer[dim]; 1330 1331 Purpose: Compute element vector for N_cb batches of elements 1332 1333 Input: 1334 Sizes: 1335 N_cb: Number of serial cell batches 1336 1337 Geometry: 1338 PetscReal v0s[Ne*dim] 1339 PetscReal jacobians[Ne*dim*dim] possibly *Nq 1340 PetscReal jacobianInverses[Ne*dim*dim] possibly *Nq 1341 PetscReal jacobianDeterminants[Ne] possibly *Nq 1342 FEM: 1343 static PetscReal quadPoints[Nq*dim] 1344 static PetscReal quadWeights[Nq] 1345 static PetscReal basis[Nq*Nb*Nc] 1346 static PetscReal basisDer[Nq*Nb*Nc*dim] 1347 PetscScalar coefficients[Ne*Nb*Nc] 1348 PetscScalar elemVec[Ne*Nb*Nc] 1349 1350 ex62.c: 1351 PetscErrorCode PetscFEIntegrateResidualBatch(PetscInt Ne, PetscInt numFields, PetscInt field, PetscQuadrature quad[], const PetscScalar coefficients[], 1352 const PetscReal v0s[], const PetscReal jacobians[], const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], 1353 void (*f0_func)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar f0[]), 1354 void (*f1_func)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar f1[]), PetscScalar elemVec[]) 1355 1356 ex52.c: 1357 PetscErrorCode IntegrateLaplacianBatchCPU(PetscInt Ne, PetscInt Nb, const PetscScalar coefficients[], const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], PetscInt Nq, const PetscReal quadPoints[], const PetscReal quadWeights[], const PetscReal basisTabulation[], const PetscReal basisDerTabulation[], PetscScalar elemVec[], AppCtx *user) 1358 PetscErrorCode IntegrateElasticityBatchCPU(PetscInt Ne, PetscInt Nb, PetscInt Ncomp, const PetscScalar coefficients[], const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], PetscInt Nq, const PetscReal quadPoints[], const PetscReal quadWeights[], const PetscReal basisTabulation[], const PetscReal basisDerTabulation[], PetscScalar elemVec[], AppCtx *user) 1359 1360 ex52_integrateElement.cu 1361 __global__ void integrateElementQuadrature(int N_cb, realType *coefficients, realType *jacobianInverses, realType *jacobianDeterminants, realType *elemVec) 1362 1363 PETSC_EXTERN PetscErrorCode IntegrateElementBatchGPU(PetscInt spatial_dim, PetscInt Ne, PetscInt Ncb, PetscInt Nbc, PetscInt Nbl, const PetscScalar coefficients[], 1364 const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], PetscScalar elemVec[], 1365 PetscLogEvent event, PetscInt debug, PetscInt pde_op) 1366 1367 ex52_integrateElementOpenCL.c: 1368 PETSC_EXTERN PetscErrorCode IntegrateElementBatchGPU(PetscInt spatial_dim, PetscInt Ne, PetscInt Ncb, PetscInt Nbc, PetscInt N_bl, const PetscScalar coefficients[], 1369 const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], PetscScalar elemVec[], 1370 PetscLogEvent event, PetscInt debug, PetscInt pde_op) 1371 1372 __kernel void integrateElementQuadrature(int N_cb, __global float *coefficients, __global float *jacobianInverses, __global float *jacobianDeterminants, __global float *elemVec) 1373 */ 1374 1375 /*@ 1376 PetscFEIntegrate - Produce the integral for the given field for a chunk of elements by quadrature integration 1377 1378 Not Collective 1379 1380 Input Parameters: 1381 + prob - The `PetscDS` specifying the discretizations and continuum functions 1382 . field - The field being integrated 1383 . Ne - The number of elements in the chunk 1384 . cgeom - The cell geometry for each cell in the chunk 1385 . coefficients - The array of FEM basis coefficients for the elements 1386 . probAux - The `PetscDS` specifying the auxiliary discretizations 1387 - coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 1388 1389 Output Parameter: 1390 . integral - the integral for this field 1391 1392 Level: intermediate 1393 1394 Developer Notes: 1395 The function name begins with `PetscFE` and yet the first argument is `PetscDS` and it has no `PetscFE` arguments. 1396 1397 .seealso: `PetscFE`, `PetscDS`, `PetscFEIntegrateResidual()`, `PetscFEIntegrateBd()` 1398 @*/ 1399 PetscErrorCode PetscFEIntegrate(PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom, const PetscScalar coefficients[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscScalar integral[]) 1400 { 1401 PetscFE fe; 1402 1403 PetscFunctionBegin; 1404 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1405 PetscCall(PetscDSGetDiscretization(prob, field, (PetscObject *)&fe)); 1406 if (fe->ops->integrate) PetscCall((*fe->ops->integrate)(prob, field, Ne, cgeom, coefficients, probAux, coefficientsAux, integral)); 1407 PetscFunctionReturn(PETSC_SUCCESS); 1408 } 1409 1410 /*@C 1411 PetscFEIntegrateBd - Produce the integral for the given field for a chunk of elements by quadrature integration 1412 1413 Not Collective 1414 1415 Input Parameters: 1416 + prob - The `PetscDS` specifying the discretizations and continuum functions 1417 . field - The field being integrated 1418 . obj_func - The function to be integrated 1419 . Ne - The number of elements in the chunk 1420 . geom - The face geometry for each face in the chunk 1421 . coefficients - The array of FEM basis coefficients for the elements 1422 . probAux - The `PetscDS` specifying the auxiliary discretizations 1423 - coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 1424 1425 Output Parameter: 1426 . integral - the integral for this field 1427 1428 Level: intermediate 1429 1430 Developer Notes: 1431 The function name begins with `PetscFE` and yet the first argument is `PetscDS` and it has no `PetscFE` arguments. 1432 1433 .seealso: `PetscFE`, `PetscDS`, `PetscFEIntegrateResidual()`, `PetscFEIntegrate()` 1434 @*/ 1435 PetscErrorCode PetscFEIntegrateBd(PetscDS prob, PetscInt field, void (*obj_func)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt Ne, PetscFEGeom *geom, const PetscScalar coefficients[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscScalar integral[]) 1436 { 1437 PetscFE fe; 1438 1439 PetscFunctionBegin; 1440 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1441 PetscCall(PetscDSGetDiscretization(prob, field, (PetscObject *)&fe)); 1442 if (fe->ops->integratebd) PetscCall((*fe->ops->integratebd)(prob, field, obj_func, Ne, geom, coefficients, probAux, coefficientsAux, integral)); 1443 PetscFunctionReturn(PETSC_SUCCESS); 1444 } 1445 1446 /*@ 1447 PetscFEIntegrateResidual - Produce the element residual vector for a chunk of elements by quadrature integration 1448 1449 Not Collective 1450 1451 Input Parameters: 1452 + ds - The `PetscDS` specifying the discretizations and continuum functions 1453 . key - The (label+value, field) being integrated 1454 . Ne - The number of elements in the chunk 1455 . cgeom - The cell geometry for each cell in the chunk 1456 . coefficients - The array of FEM basis coefficients for the elements 1457 . coefficients_t - The array of FEM basis time derivative coefficients for the elements 1458 . probAux - The `PetscDS` specifying the auxiliary discretizations 1459 . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 1460 - t - The time 1461 1462 Output Parameter: 1463 . elemVec - the element residual vectors from each element 1464 1465 Level: intermediate 1466 1467 Note: 1468 .vb 1469 Loop over batch of elements (e): 1470 Loop over quadrature points (q): 1471 Make u_q and gradU_q (loops over fields,Nb,Ncomp) and x_q 1472 Call f_0 and f_1 1473 Loop over element vector entries (f,fc --> i): 1474 elemVec[i] += \psi^{fc}_f(q) f0_{fc}(u, \nabla u) + \nabla\psi^{fc}_f(q) \cdot f1_{fc,df}(u, \nabla u) 1475 .ve 1476 1477 .seealso: `PetscFEIntegrateBdResidual()` 1478 @*/ 1479 PetscErrorCode PetscFEIntegrateResidual(PetscDS ds, PetscFormKey key, PetscInt Ne, PetscFEGeom *cgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[]) 1480 { 1481 PetscFE fe; 1482 1483 PetscFunctionBeginHot; 1484 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1485 PetscCall(PetscDSGetDiscretization(ds, key.field, (PetscObject *)&fe)); 1486 if (fe->ops->integrateresidual) PetscCall((*fe->ops->integrateresidual)(ds, key, Ne, cgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec)); 1487 PetscFunctionReturn(PETSC_SUCCESS); 1488 } 1489 1490 /*@ 1491 PetscFEIntegrateBdResidual - Produce the element residual vector for a chunk of elements by quadrature integration over a boundary 1492 1493 Not Collective 1494 1495 Input Parameters: 1496 + ds - The `PetscDS` specifying the discretizations and continuum functions 1497 . wf - The PetscWeakForm object holding the pointwise functions 1498 . key - The (label+value, field) being integrated 1499 . Ne - The number of elements in the chunk 1500 . fgeom - The face geometry for each cell in the chunk 1501 . coefficients - The array of FEM basis coefficients for the elements 1502 . coefficients_t - The array of FEM basis time derivative coefficients for the elements 1503 . probAux - The `PetscDS` specifying the auxiliary discretizations 1504 . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 1505 - t - The time 1506 1507 Output Parameter: 1508 . elemVec - the element residual vectors from each element 1509 1510 Level: intermediate 1511 1512 .seealso: `PetscFEIntegrateResidual()` 1513 @*/ 1514 PetscErrorCode PetscFEIntegrateBdResidual(PetscDS ds, PetscWeakForm wf, PetscFormKey key, PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[]) 1515 { 1516 PetscFE fe; 1517 1518 PetscFunctionBegin; 1519 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1520 PetscCall(PetscDSGetDiscretization(ds, key.field, (PetscObject *)&fe)); 1521 if (fe->ops->integratebdresidual) PetscCall((*fe->ops->integratebdresidual)(ds, wf, key, Ne, fgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec)); 1522 PetscFunctionReturn(PETSC_SUCCESS); 1523 } 1524 1525 /*@ 1526 PetscFEIntegrateHybridResidual - Produce the element residual vector for a chunk of hybrid element faces by quadrature integration 1527 1528 Not Collective 1529 1530 Input Parameters: 1531 + ds - The `PetscDS` specifying the discretizations and continuum functions 1532 . dsIn - The `PetscDS` specifying the discretizations and continuum functions for input 1533 . key - The (label+value, field) being integrated 1534 . s - The side of the cell being integrated, 0 for negative and 1 for positive 1535 . Ne - The number of elements in the chunk 1536 . fgeom - The face geometry for each cell in the chunk 1537 . cgeom - The cell geometry for each neighbor cell in the chunk 1538 . coefficients - The array of FEM basis coefficients for the elements 1539 . coefficients_t - The array of FEM basis time derivative coefficients for the elements 1540 . probAux - The `PetscDS` specifying the auxiliary discretizations 1541 . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 1542 - t - The time 1543 1544 Output Parameter: 1545 . elemVec - the element residual vectors from each element 1546 1547 Level: developer 1548 1549 .seealso: `PetscFEIntegrateResidual()` 1550 @*/ 1551 PetscErrorCode PetscFEIntegrateHybridResidual(PetscDS ds, PetscDS dsIn, PetscFormKey key, PetscInt s, PetscInt Ne, PetscFEGeom *fgeom, PetscFEGeom *cgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[]) 1552 { 1553 PetscFE fe; 1554 1555 PetscFunctionBegin; 1556 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1557 PetscValidHeaderSpecific(dsIn, PETSCDS_CLASSID, 2); 1558 PetscCall(PetscDSGetDiscretization(ds, key.field, (PetscObject *)&fe)); 1559 if (fe->ops->integratehybridresidual) PetscCall((*fe->ops->integratehybridresidual)(ds, dsIn, key, s, Ne, fgeom, cgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec)); 1560 PetscFunctionReturn(PETSC_SUCCESS); 1561 } 1562 1563 /*@ 1564 PetscFEIntegrateJacobian - Produce the element Jacobian for a chunk of elements by quadrature integration 1565 1566 Not Collective 1567 1568 Input Parameters: 1569 + rds - The `PetscDS` specifying the row discretizations and continuum functions 1570 . cds - The `PetscDS` specifying the column discretizations 1571 . jtype - The type of matrix pointwise functions that should be used 1572 . key - The (label+value, fieldI*Nf + fieldJ) being integrated 1573 . Ne - The number of elements in the chunk 1574 . cgeom - The cell geometry for each cell in the chunk 1575 . coefficients - The array of FEM basis coefficients for the elements for the Jacobian evaluation point 1576 . coefficients_t - The array of FEM basis time derivative coefficients for the elements 1577 . dsAux - The `PetscDS` specifying the auxiliary discretizations 1578 . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 1579 . t - The time 1580 - u_tshift - A multiplier for the dF/du_t term (as opposed to the dF/du term) 1581 1582 Output Parameter: 1583 . elemMat - the element matrices for the Jacobian from each element 1584 1585 Level: intermediate 1586 1587 Note: 1588 .vb 1589 Loop over batch of elements (e): 1590 Loop over element matrix entries (f,fc,g,gc --> i,j): 1591 Loop over quadrature points (q): 1592 Make u_q and gradU_q (loops over fields,Nb,Ncomp) 1593 elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q) 1594 + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q) 1595 + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q) 1596 + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q) 1597 .ve 1598 1599 .seealso: `PetscFEIntegrateResidual()` 1600 @*/ 1601 PetscErrorCode PetscFEIntegrateJacobian(PetscDS rds, PetscDS cds, PetscFEJacobianType jtype, PetscFormKey key, PetscInt Ne, PetscFEGeom *cgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) 1602 { 1603 PetscFE fe; 1604 PetscInt Nf; 1605 1606 PetscFunctionBegin; 1607 PetscValidHeaderSpecific(rds, PETSCDS_CLASSID, 1); 1608 PetscValidHeaderSpecific(cds, PETSCDS_CLASSID, 2); 1609 PetscCall(PetscDSGetNumFields(rds, &Nf)); 1610 PetscCall(PetscDSGetDiscretization(rds, key.field / Nf, (PetscObject *)&fe)); 1611 if (fe->ops->integratejacobian) PetscCall((*fe->ops->integratejacobian)(rds, cds, jtype, key, Ne, cgeom, coefficients, coefficients_t, dsAux, coefficientsAux, t, u_tshift, elemMat)); 1612 PetscFunctionReturn(PETSC_SUCCESS); 1613 } 1614 1615 /*@ 1616 PetscFEIntegrateBdJacobian - Produce the boundary element Jacobian for a chunk of elements by quadrature integration 1617 1618 Not Collective 1619 1620 Input Parameters: 1621 + ds - The `PetscDS` specifying the discretizations and continuum functions 1622 . wf - The PetscWeakForm holding the pointwise functions 1623 . jtype - The type of matrix pointwise functions that should be used 1624 . key - The (label+value, fieldI*Nf + fieldJ) being integrated 1625 . Ne - The number of elements in the chunk 1626 . fgeom - The face geometry for each cell in the chunk 1627 . coefficients - The array of FEM basis coefficients for the elements for the Jacobian evaluation point 1628 . coefficients_t - The array of FEM basis time derivative coefficients for the elements 1629 . probAux - The `PetscDS` specifying the auxiliary discretizations 1630 . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 1631 . t - The time 1632 - u_tshift - A multiplier for the dF/du_t term (as opposed to the dF/du term) 1633 1634 Output Parameter: 1635 . elemMat - the element matrices for the Jacobian from each element 1636 1637 Level: intermediate 1638 1639 Note: 1640 .vb 1641 Loop over batch of elements (e): 1642 Loop over element matrix entries (f,fc,g,gc --> i,j): 1643 Loop over quadrature points (q): 1644 Make u_q and gradU_q (loops over fields,Nb,Ncomp) 1645 elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q) 1646 + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q) 1647 + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q) 1648 + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q) 1649 .ve 1650 1651 .seealso: `PetscFEIntegrateJacobian()`, `PetscFEIntegrateResidual()` 1652 @*/ 1653 PetscErrorCode PetscFEIntegrateBdJacobian(PetscDS ds, PetscWeakForm wf, PetscFEJacobianType jtype, PetscFormKey key, PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) 1654 { 1655 PetscFE fe; 1656 PetscInt Nf; 1657 1658 PetscFunctionBegin; 1659 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1660 PetscCall(PetscDSGetNumFields(ds, &Nf)); 1661 PetscCall(PetscDSGetDiscretization(ds, key.field / Nf, (PetscObject *)&fe)); 1662 if (fe->ops->integratebdjacobian) PetscCall((*fe->ops->integratebdjacobian)(ds, wf, jtype, key, Ne, fgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, u_tshift, elemMat)); 1663 PetscFunctionReturn(PETSC_SUCCESS); 1664 } 1665 1666 /*@ 1667 PetscFEIntegrateHybridJacobian - Produce the boundary element Jacobian for a chunk of hybrid elements by quadrature integration 1668 1669 Not Collective 1670 1671 Input Parameters: 1672 + ds - The `PetscDS` specifying the discretizations and continuum functions for the output 1673 . dsIn - The `PetscDS` specifying the discretizations and continuum functions for the input 1674 . jtype - The type of matrix pointwise functions that should be used 1675 . key - The (label+value, fieldI*Nf + fieldJ) being integrated 1676 . s - The side of the cell being integrated, 0 for negative and 1 for positive 1677 . Ne - The number of elements in the chunk 1678 . fgeom - The face geometry for each cell in the chunk 1679 . cgeom - The cell geometry for each neighbor cell in the chunk 1680 . coefficients - The array of FEM basis coefficients for the elements for the Jacobian evaluation point 1681 . coefficients_t - The array of FEM basis time derivative coefficients for the elements 1682 . probAux - The `PetscDS` specifying the auxiliary discretizations 1683 . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 1684 . t - The time 1685 - u_tshift - A multiplier for the dF/du_t term (as opposed to the dF/du term) 1686 1687 Output Parameter: 1688 . elemMat - the element matrices for the Jacobian from each element 1689 1690 Level: developer 1691 1692 Note: 1693 .vb 1694 Loop over batch of elements (e): 1695 Loop over element matrix entries (f,fc,g,gc --> i,j): 1696 Loop over quadrature points (q): 1697 Make u_q and gradU_q (loops over fields,Nb,Ncomp) 1698 elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q) 1699 + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q) 1700 + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q) 1701 + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q) 1702 .ve 1703 1704 .seealso: `PetscFEIntegrateJacobian()`, `PetscFEIntegrateResidual()` 1705 @*/ 1706 PetscErrorCode PetscFEIntegrateHybridJacobian(PetscDS ds, PetscDS dsIn, PetscFEJacobianType jtype, PetscFormKey key, PetscInt s, PetscInt Ne, PetscFEGeom *fgeom, PetscFEGeom *cgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) 1707 { 1708 PetscFE fe; 1709 PetscInt Nf; 1710 1711 PetscFunctionBegin; 1712 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1713 PetscCall(PetscDSGetNumFields(ds, &Nf)); 1714 PetscCall(PetscDSGetDiscretization(ds, key.field / Nf, (PetscObject *)&fe)); 1715 if (fe->ops->integratehybridjacobian) PetscCall((*fe->ops->integratehybridjacobian)(ds, dsIn, jtype, key, s, Ne, fgeom, cgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, u_tshift, elemMat)); 1716 PetscFunctionReturn(PETSC_SUCCESS); 1717 } 1718 1719 /*@ 1720 PetscFEGetHeightSubspace - Get the subspace of this space for a mesh point of a given height 1721 1722 Input Parameters: 1723 + fe - The finite element space 1724 - height - The height of the `DMPLEX` point 1725 1726 Output Parameter: 1727 . subfe - The subspace of this `PetscFE` space 1728 1729 Level: advanced 1730 1731 Note: 1732 For example, if we want the subspace of this space for a face, we would choose height = 1. 1733 1734 .seealso: `PetscFECreateDefault()` 1735 @*/ 1736 PetscErrorCode PetscFEGetHeightSubspace(PetscFE fe, PetscInt height, PetscFE *subfe) 1737 { 1738 PetscSpace P, subP; 1739 PetscDualSpace Q, subQ; 1740 PetscQuadrature subq; 1741 PetscInt dim, Nc; 1742 1743 PetscFunctionBegin; 1744 PetscValidHeaderSpecific(fe, PETSCFE_CLASSID, 1); 1745 PetscAssertPointer(subfe, 3); 1746 if (height == 0) { 1747 *subfe = fe; 1748 PetscFunctionReturn(PETSC_SUCCESS); 1749 } 1750 PetscCall(PetscFEGetBasisSpace(fe, &P)); 1751 PetscCall(PetscFEGetDualSpace(fe, &Q)); 1752 PetscCall(PetscFEGetNumComponents(fe, &Nc)); 1753 PetscCall(PetscFEGetFaceQuadrature(fe, &subq)); 1754 PetscCall(PetscDualSpaceGetDimension(Q, &dim)); 1755 PetscCheck(height <= dim && height >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Asked for space at height %" PetscInt_FMT " for dimension %" PetscInt_FMT " space", height, dim); 1756 if (!fe->subspaces) PetscCall(PetscCalloc1(dim, &fe->subspaces)); 1757 if (height <= dim) { 1758 if (!fe->subspaces[height - 1]) { 1759 PetscFE sub = NULL; 1760 const char *name; 1761 1762 PetscCall(PetscSpaceGetHeightSubspace(P, height, &subP)); 1763 PetscCall(PetscDualSpaceGetHeightSubspace(Q, height, &subQ)); 1764 if (subQ) { 1765 PetscCall(PetscObjectReference((PetscObject)subP)); 1766 PetscCall(PetscObjectReference((PetscObject)subQ)); 1767 PetscCall(PetscObjectReference((PetscObject)subq)); 1768 PetscCall(PetscFECreateFromSpaces(subP, subQ, subq, NULL, &sub)); 1769 } 1770 if (sub) { 1771 PetscCall(PetscObjectGetName((PetscObject)fe, &name)); 1772 if (name) PetscCall(PetscFESetName(sub, name)); 1773 } 1774 fe->subspaces[height - 1] = sub; 1775 } 1776 *subfe = fe->subspaces[height - 1]; 1777 } else { 1778 *subfe = NULL; 1779 } 1780 PetscFunctionReturn(PETSC_SUCCESS); 1781 } 1782 1783 /*@ 1784 PetscFERefine - Create a "refined" `PetscFE` object that refines the reference cell into 1785 smaller copies. 1786 1787 Collective 1788 1789 Input Parameter: 1790 . fe - The initial `PetscFE` 1791 1792 Output Parameter: 1793 . feRef - The refined `PetscFE` 1794 1795 Level: advanced 1796 1797 Notes: 1798 This is typically used to generate a preconditioner for a higher order method from a lower order method on a 1799 refined mesh having the same number of dofs (but more sparsity). It is also used to create an 1800 interpolation between regularly refined meshes. 1801 1802 .seealso: `PetscFEType`, `PetscFECreate()`, `PetscFESetType()` 1803 @*/ 1804 PetscErrorCode PetscFERefine(PetscFE fe, PetscFE *feRef) 1805 { 1806 PetscSpace P, Pref; 1807 PetscDualSpace Q, Qref; 1808 DM K, Kref; 1809 PetscQuadrature q, qref; 1810 const PetscReal *v0, *jac; 1811 PetscInt numComp, numSubelements; 1812 PetscInt cStart, cEnd, c; 1813 PetscDualSpace *cellSpaces; 1814 1815 PetscFunctionBegin; 1816 PetscCall(PetscFEGetBasisSpace(fe, &P)); 1817 PetscCall(PetscFEGetDualSpace(fe, &Q)); 1818 PetscCall(PetscFEGetQuadrature(fe, &q)); 1819 PetscCall(PetscDualSpaceGetDM(Q, &K)); 1820 /* Create space */ 1821 PetscCall(PetscObjectReference((PetscObject)P)); 1822 Pref = P; 1823 /* Create dual space */ 1824 PetscCall(PetscDualSpaceDuplicate(Q, &Qref)); 1825 PetscCall(PetscDualSpaceSetType(Qref, PETSCDUALSPACEREFINED)); 1826 PetscCall(DMRefine(K, PetscObjectComm((PetscObject)fe), &Kref)); 1827 PetscCall(DMGetCoordinatesLocalSetUp(Kref)); 1828 PetscCall(PetscDualSpaceSetDM(Qref, Kref)); 1829 PetscCall(DMPlexGetHeightStratum(Kref, 0, &cStart, &cEnd)); 1830 PetscCall(PetscMalloc1(cEnd - cStart, &cellSpaces)); 1831 /* TODO: fix for non-uniform refinement */ 1832 for (c = 0; c < cEnd - cStart; c++) cellSpaces[c] = Q; 1833 PetscCall(PetscDualSpaceRefinedSetCellSpaces(Qref, cellSpaces)); 1834 PetscCall(PetscFree(cellSpaces)); 1835 PetscCall(DMDestroy(&Kref)); 1836 PetscCall(PetscDualSpaceSetUp(Qref)); 1837 /* Create element */ 1838 PetscCall(PetscFECreate(PetscObjectComm((PetscObject)fe), feRef)); 1839 PetscCall(PetscFESetType(*feRef, PETSCFECOMPOSITE)); 1840 PetscCall(PetscFESetBasisSpace(*feRef, Pref)); 1841 PetscCall(PetscFESetDualSpace(*feRef, Qref)); 1842 PetscCall(PetscFEGetNumComponents(fe, &numComp)); 1843 PetscCall(PetscFESetNumComponents(*feRef, numComp)); 1844 PetscCall(PetscFESetUp(*feRef)); 1845 PetscCall(PetscSpaceDestroy(&Pref)); 1846 PetscCall(PetscDualSpaceDestroy(&Qref)); 1847 /* Create quadrature */ 1848 PetscCall(PetscFECompositeGetMapping(*feRef, &numSubelements, &v0, &jac, NULL)); 1849 PetscCall(PetscQuadratureExpandComposite(q, numSubelements, v0, jac, &qref)); 1850 PetscCall(PetscFESetQuadrature(*feRef, qref)); 1851 PetscCall(PetscQuadratureDestroy(&qref)); 1852 PetscFunctionReturn(PETSC_SUCCESS); 1853 } 1854 1855 static PetscErrorCode PetscFESetDefaultName_Private(PetscFE fe) 1856 { 1857 PetscSpace P; 1858 PetscDualSpace Q; 1859 DM K; 1860 DMPolytopeType ct; 1861 PetscInt degree; 1862 char name[64]; 1863 1864 PetscFunctionBegin; 1865 PetscCall(PetscFEGetBasisSpace(fe, &P)); 1866 PetscCall(PetscSpaceGetDegree(P, °ree, NULL)); 1867 PetscCall(PetscFEGetDualSpace(fe, &Q)); 1868 PetscCall(PetscDualSpaceGetDM(Q, &K)); 1869 PetscCall(DMPlexGetCellType(K, 0, &ct)); 1870 switch (ct) { 1871 case DM_POLYTOPE_SEGMENT: 1872 case DM_POLYTOPE_POINT_PRISM_TENSOR: 1873 case DM_POLYTOPE_QUADRILATERAL: 1874 case DM_POLYTOPE_SEG_PRISM_TENSOR: 1875 case DM_POLYTOPE_HEXAHEDRON: 1876 case DM_POLYTOPE_QUAD_PRISM_TENSOR: 1877 PetscCall(PetscSNPrintf(name, sizeof(name), "Q%" PetscInt_FMT, degree)); 1878 break; 1879 case DM_POLYTOPE_TRIANGLE: 1880 case DM_POLYTOPE_TETRAHEDRON: 1881 PetscCall(PetscSNPrintf(name, sizeof(name), "P%" PetscInt_FMT, degree)); 1882 break; 1883 case DM_POLYTOPE_TRI_PRISM: 1884 case DM_POLYTOPE_TRI_PRISM_TENSOR: 1885 PetscCall(PetscSNPrintf(name, sizeof(name), "P%" PetscInt_FMT "xQ%" PetscInt_FMT, degree, degree)); 1886 break; 1887 default: 1888 PetscCall(PetscSNPrintf(name, sizeof(name), "FE")); 1889 } 1890 PetscCall(PetscFESetName(fe, name)); 1891 PetscFunctionReturn(PETSC_SUCCESS); 1892 } 1893 1894 /*@ 1895 PetscFECreateFromSpaces - Create a `PetscFE` from the basis and dual spaces 1896 1897 Collective 1898 1899 Input Parameters: 1900 + P - The basis space 1901 . Q - The dual space 1902 . q - The cell quadrature 1903 - fq - The face quadrature 1904 1905 Output Parameter: 1906 . fem - The `PetscFE` object 1907 1908 Level: beginner 1909 1910 Note: 1911 The `PetscFE` takes ownership of these spaces by calling destroy on each. They should not be used after this call, and for borrowed references from `PetscFEGetSpace()` and the like, the caller must use `PetscObjectReference` before this call. 1912 1913 .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscQuadrature`, 1914 `PetscFECreateLagrangeByCell()`, `PetscFECreateDefault()`, `PetscFECreateByCell()`, `PetscFECreate()`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()` 1915 @*/ 1916 PetscErrorCode PetscFECreateFromSpaces(PetscSpace P, PetscDualSpace Q, PetscQuadrature q, PetscQuadrature fq, PetscFE *fem) 1917 { 1918 PetscInt Nc; 1919 PetscInt p_Ns = -1, p_Nc = -1, q_Ns = -1, q_Nc = -1; 1920 PetscBool p_is_uniform_sum = PETSC_FALSE, p_interleave_basis = PETSC_FALSE, p_interleave_components = PETSC_FALSE; 1921 PetscBool q_is_uniform_sum = PETSC_FALSE, q_interleave_basis = PETSC_FALSE, q_interleave_components = PETSC_FALSE; 1922 const char *prefix; 1923 1924 PetscFunctionBegin; 1925 PetscCall(PetscObjectTypeCompare((PetscObject)P, PETSCSPACESUM, &p_is_uniform_sum)); 1926 if (p_is_uniform_sum) { 1927 PetscSpace subsp_0 = NULL; 1928 PetscCall(PetscSpaceSumGetNumSubspaces(P, &p_Ns)); 1929 PetscCall(PetscSpaceGetNumComponents(P, &p_Nc)); 1930 PetscCall(PetscSpaceSumGetConcatenate(P, &p_is_uniform_sum)); 1931 PetscCall(PetscSpaceSumGetInterleave(P, &p_interleave_basis, &p_interleave_components)); 1932 for (PetscInt s = 0; s < p_Ns; s++) { 1933 PetscSpace subsp; 1934 1935 PetscCall(PetscSpaceSumGetSubspace(P, s, &subsp)); 1936 if (!s) { 1937 subsp_0 = subsp; 1938 } else if (subsp != subsp_0) { 1939 p_is_uniform_sum = PETSC_FALSE; 1940 } 1941 } 1942 } 1943 PetscCall(PetscObjectTypeCompare((PetscObject)Q, PETSCDUALSPACESUM, &q_is_uniform_sum)); 1944 if (q_is_uniform_sum) { 1945 PetscDualSpace subsp_0 = NULL; 1946 PetscCall(PetscDualSpaceSumGetNumSubspaces(Q, &q_Ns)); 1947 PetscCall(PetscDualSpaceGetNumComponents(Q, &q_Nc)); 1948 PetscCall(PetscDualSpaceSumGetConcatenate(Q, &q_is_uniform_sum)); 1949 PetscCall(PetscDualSpaceSumGetInterleave(Q, &q_interleave_basis, &q_interleave_components)); 1950 for (PetscInt s = 0; s < q_Ns; s++) { 1951 PetscDualSpace subsp; 1952 1953 PetscCall(PetscDualSpaceSumGetSubspace(Q, s, &subsp)); 1954 if (!s) { 1955 subsp_0 = subsp; 1956 } else if (subsp != subsp_0) { 1957 q_is_uniform_sum = PETSC_FALSE; 1958 } 1959 } 1960 } 1961 if (p_is_uniform_sum && q_is_uniform_sum && (p_interleave_basis == q_interleave_basis) && (p_interleave_components == q_interleave_components) && (p_Ns == q_Ns) && (p_Nc == q_Nc)) { 1962 PetscSpace scalar_space; 1963 PetscDualSpace scalar_dspace; 1964 PetscFE scalar_fe; 1965 1966 PetscCall(PetscSpaceSumGetSubspace(P, 0, &scalar_space)); 1967 PetscCall(PetscDualSpaceSumGetSubspace(Q, 0, &scalar_dspace)); 1968 PetscCall(PetscObjectReference((PetscObject)scalar_space)); 1969 PetscCall(PetscObjectReference((PetscObject)scalar_dspace)); 1970 PetscCall(PetscObjectReference((PetscObject)q)); 1971 PetscCall(PetscObjectReference((PetscObject)fq)); 1972 PetscCall(PetscFECreateFromSpaces(scalar_space, scalar_dspace, q, fq, &scalar_fe)); 1973 PetscCall(PetscFECreateVector(scalar_fe, p_Ns, p_interleave_basis, p_interleave_components, fem)); 1974 PetscCall(PetscFEDestroy(&scalar_fe)); 1975 } else { 1976 PetscCall(PetscFECreate(PetscObjectComm((PetscObject)P), fem)); 1977 PetscCall(PetscFESetType(*fem, PETSCFEBASIC)); 1978 } 1979 PetscCall(PetscSpaceGetNumComponents(P, &Nc)); 1980 PetscCall(PetscFESetNumComponents(*fem, Nc)); 1981 PetscCall(PetscFESetBasisSpace(*fem, P)); 1982 PetscCall(PetscFESetDualSpace(*fem, Q)); 1983 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)P, &prefix)); 1984 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*fem, prefix)); 1985 PetscCall(PetscFESetUp(*fem)); 1986 PetscCall(PetscSpaceDestroy(&P)); 1987 PetscCall(PetscDualSpaceDestroy(&Q)); 1988 PetscCall(PetscFESetQuadrature(*fem, q)); 1989 PetscCall(PetscFESetFaceQuadrature(*fem, fq)); 1990 PetscCall(PetscQuadratureDestroy(&q)); 1991 PetscCall(PetscQuadratureDestroy(&fq)); 1992 PetscCall(PetscFESetDefaultName_Private(*fem)); 1993 PetscFunctionReturn(PETSC_SUCCESS); 1994 } 1995 1996 static PetscErrorCode PetscFECreate_Internal(MPI_Comm comm, PetscInt dim, PetscInt Nc, DMPolytopeType ct, const char prefix[], PetscInt degree, PetscInt qorder, PetscBool setFromOptions, PetscFE *fem) 1997 { 1998 DM K; 1999 PetscSpace P; 2000 PetscDualSpace Q; 2001 PetscQuadrature q, fq; 2002 PetscBool tensor; 2003 2004 PetscFunctionBegin; 2005 if (prefix) PetscAssertPointer(prefix, 5); 2006 PetscAssertPointer(fem, 9); 2007 switch (ct) { 2008 case DM_POLYTOPE_SEGMENT: 2009 case DM_POLYTOPE_POINT_PRISM_TENSOR: 2010 case DM_POLYTOPE_QUADRILATERAL: 2011 case DM_POLYTOPE_SEG_PRISM_TENSOR: 2012 case DM_POLYTOPE_HEXAHEDRON: 2013 case DM_POLYTOPE_QUAD_PRISM_TENSOR: 2014 tensor = PETSC_TRUE; 2015 break; 2016 default: 2017 tensor = PETSC_FALSE; 2018 } 2019 /* Create space */ 2020 PetscCall(PetscSpaceCreate(comm, &P)); 2021 PetscCall(PetscSpaceSetType(P, PETSCSPACEPOLYNOMIAL)); 2022 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)P, prefix)); 2023 PetscCall(PetscSpacePolynomialSetTensor(P, tensor)); 2024 PetscCall(PetscSpaceSetNumComponents(P, Nc)); 2025 PetscCall(PetscSpaceSetNumVariables(P, dim)); 2026 if (degree >= 0) { 2027 PetscCall(PetscSpaceSetDegree(P, degree, PETSC_DETERMINE)); 2028 if (ct == DM_POLYTOPE_TRI_PRISM || ct == DM_POLYTOPE_TRI_PRISM_TENSOR) { 2029 PetscSpace Pend, Pside; 2030 2031 PetscCall(PetscSpaceSetNumComponents(P, 1)); 2032 PetscCall(PetscSpaceCreate(comm, &Pend)); 2033 PetscCall(PetscSpaceSetType(Pend, PETSCSPACEPOLYNOMIAL)); 2034 PetscCall(PetscSpacePolynomialSetTensor(Pend, PETSC_FALSE)); 2035 PetscCall(PetscSpaceSetNumComponents(Pend, 1)); 2036 PetscCall(PetscSpaceSetNumVariables(Pend, dim - 1)); 2037 PetscCall(PetscSpaceSetDegree(Pend, degree, PETSC_DETERMINE)); 2038 PetscCall(PetscSpaceCreate(comm, &Pside)); 2039 PetscCall(PetscSpaceSetType(Pside, PETSCSPACEPOLYNOMIAL)); 2040 PetscCall(PetscSpacePolynomialSetTensor(Pside, PETSC_FALSE)); 2041 PetscCall(PetscSpaceSetNumComponents(Pside, 1)); 2042 PetscCall(PetscSpaceSetNumVariables(Pside, 1)); 2043 PetscCall(PetscSpaceSetDegree(Pside, degree, PETSC_DETERMINE)); 2044 PetscCall(PetscSpaceSetType(P, PETSCSPACETENSOR)); 2045 PetscCall(PetscSpaceTensorSetNumSubspaces(P, 2)); 2046 PetscCall(PetscSpaceTensorSetSubspace(P, 0, Pend)); 2047 PetscCall(PetscSpaceTensorSetSubspace(P, 1, Pside)); 2048 PetscCall(PetscSpaceDestroy(&Pend)); 2049 PetscCall(PetscSpaceDestroy(&Pside)); 2050 2051 if (Nc > 1) { 2052 PetscSpace scalar_P = P; 2053 2054 PetscCall(PetscSpaceCreate(comm, &P)); 2055 PetscCall(PetscSpaceSetNumVariables(P, dim)); 2056 PetscCall(PetscSpaceSetNumComponents(P, Nc)); 2057 PetscCall(PetscSpaceSetType(P, PETSCSPACESUM)); 2058 PetscCall(PetscSpaceSumSetNumSubspaces(P, Nc)); 2059 PetscCall(PetscSpaceSumSetConcatenate(P, PETSC_TRUE)); 2060 PetscCall(PetscSpaceSumSetInterleave(P, PETSC_TRUE, PETSC_FALSE)); 2061 for (PetscInt i = 0; i < Nc; i++) PetscCall(PetscSpaceSumSetSubspace(P, i, scalar_P)); 2062 PetscCall(PetscSpaceDestroy(&scalar_P)); 2063 } 2064 } 2065 } 2066 if (setFromOptions) PetscCall(PetscSpaceSetFromOptions(P)); 2067 PetscCall(PetscSpaceSetUp(P)); 2068 PetscCall(PetscSpaceGetDegree(P, °ree, NULL)); 2069 PetscCall(PetscSpacePolynomialGetTensor(P, &tensor)); 2070 PetscCall(PetscSpaceGetNumComponents(P, &Nc)); 2071 /* Create dual space */ 2072 PetscCall(PetscDualSpaceCreate(comm, &Q)); 2073 PetscCall(PetscDualSpaceSetType(Q, PETSCDUALSPACELAGRANGE)); 2074 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)Q, prefix)); 2075 PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, ct, &K)); 2076 PetscCall(PetscDualSpaceSetDM(Q, K)); 2077 PetscCall(DMDestroy(&K)); 2078 PetscCall(PetscDualSpaceSetNumComponents(Q, Nc)); 2079 PetscCall(PetscDualSpaceSetOrder(Q, degree)); 2080 PetscCall(PetscDualSpaceLagrangeSetTensor(Q, (tensor || (ct == DM_POLYTOPE_TRI_PRISM)) ? PETSC_TRUE : PETSC_FALSE)); 2081 if (setFromOptions) PetscCall(PetscDualSpaceSetFromOptions(Q)); 2082 PetscCall(PetscDualSpaceSetUp(Q)); 2083 /* Create quadrature */ 2084 PetscDTSimplexQuadratureType qtype = PETSCDTSIMPLEXQUAD_DEFAULT; 2085 2086 qorder = qorder >= 0 ? qorder : degree; 2087 if (setFromOptions) { 2088 PetscObjectOptionsBegin((PetscObject)P); 2089 PetscCall(PetscOptionsBoundedInt("-petscfe_default_quadrature_order", "Quadrature order is one less than quadrature points per edge", "PetscFECreateDefault", qorder, &qorder, NULL, 0)); 2090 PetscCall(PetscOptionsEnum("-petscfe_default_quadrature_type", "Simplex quadrature type", "PetscDTSimplexQuadratureType", PetscDTSimplexQuadratureTypes, (PetscEnum)qtype, (PetscEnum *)&qtype, NULL)); 2091 PetscOptionsEnd(); 2092 } 2093 PetscCall(PetscDTCreateQuadratureByCell(ct, qorder, qtype, &q, &fq)); 2094 /* Create finite element */ 2095 PetscCall(PetscFECreateFromSpaces(P, Q, q, fq, fem)); 2096 if (setFromOptions) PetscCall(PetscFESetFromOptions(*fem)); 2097 PetscFunctionReturn(PETSC_SUCCESS); 2098 } 2099 2100 /*@ 2101 PetscFECreateDefault - Create a `PetscFE` for basic FEM computation 2102 2103 Collective 2104 2105 Input Parameters: 2106 + comm - The MPI comm 2107 . dim - The spatial dimension 2108 . Nc - The number of components 2109 . isSimplex - Flag for simplex reference cell, otherwise its a tensor product 2110 . prefix - The options prefix, or `NULL` 2111 - qorder - The quadrature order or `PETSC_DETERMINE` to use `PetscSpace` polynomial degree 2112 2113 Output Parameter: 2114 . fem - The `PetscFE` object 2115 2116 Level: beginner 2117 2118 Note: 2119 Each subobject is SetFromOption() during creation, so that the object may be customized from the command line, using the prefix specified above. See the links below for the particular options available. 2120 2121 .seealso: `PetscFECreateLagrange()`, `PetscFECreateByCell()`, `PetscSpaceSetFromOptions()`, `PetscDualSpaceSetFromOptions()`, `PetscFESetFromOptions()`, `PetscFECreate()`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()` 2122 @*/ 2123 PetscErrorCode PetscFECreateDefault(MPI_Comm comm, PetscInt dim, PetscInt Nc, PetscBool isSimplex, const char prefix[], PetscInt qorder, PetscFE *fem) 2124 { 2125 PetscFunctionBegin; 2126 PetscCall(PetscFECreate_Internal(comm, dim, Nc, DMPolytopeTypeSimpleShape(dim, isSimplex), prefix, PETSC_DECIDE, qorder, PETSC_TRUE, fem)); 2127 PetscFunctionReturn(PETSC_SUCCESS); 2128 } 2129 2130 /*@ 2131 PetscFECreateByCell - Create a `PetscFE` for basic FEM computation 2132 2133 Collective 2134 2135 Input Parameters: 2136 + comm - The MPI comm 2137 . dim - The spatial dimension 2138 . Nc - The number of components 2139 . ct - The celltype of the reference cell 2140 . prefix - The options prefix, or `NULL` 2141 - qorder - The quadrature order or `PETSC_DETERMINE` to use `PetscSpace` polynomial degree 2142 2143 Output Parameter: 2144 . fem - The `PetscFE` object 2145 2146 Level: beginner 2147 2148 Note: 2149 Each subobject is SetFromOption() during creation, so that the object may be customized from the command line, using the prefix specified above. See the links below for the particular options available. 2150 2151 .seealso: `PetscFECreateDefault()`, `PetscFECreateLagrange()`, `PetscSpaceSetFromOptions()`, `PetscDualSpaceSetFromOptions()`, `PetscFESetFromOptions()`, `PetscFECreate()`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()` 2152 @*/ 2153 PetscErrorCode PetscFECreateByCell(MPI_Comm comm, PetscInt dim, PetscInt Nc, DMPolytopeType ct, const char prefix[], PetscInt qorder, PetscFE *fem) 2154 { 2155 PetscFunctionBegin; 2156 PetscCall(PetscFECreate_Internal(comm, dim, Nc, ct, prefix, PETSC_DECIDE, qorder, PETSC_TRUE, fem)); 2157 PetscFunctionReturn(PETSC_SUCCESS); 2158 } 2159 2160 /*@ 2161 PetscFECreateLagrange - Create a `PetscFE` for the basic Lagrange space of degree k 2162 2163 Collective 2164 2165 Input Parameters: 2166 + comm - The MPI comm 2167 . dim - The spatial dimension 2168 . Nc - The number of components 2169 . isSimplex - Flag for simplex reference cell, otherwise its a tensor product 2170 . k - The degree k of the space 2171 - qorder - The quadrature order or `PETSC_DETERMINE` to use `PetscSpace` polynomial degree 2172 2173 Output Parameter: 2174 . fem - The `PetscFE` object 2175 2176 Level: beginner 2177 2178 Note: 2179 For simplices, this element is the space of maximum polynomial degree k, otherwise it is a tensor product of 1D polynomials, each with maximal degree k. 2180 2181 .seealso: `PetscFECreateLagrangeByCell()`, `PetscFECreateDefault()`, `PetscFECreateByCell()`, `PetscFECreate()`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()` 2182 @*/ 2183 PetscErrorCode PetscFECreateLagrange(MPI_Comm comm, PetscInt dim, PetscInt Nc, PetscBool isSimplex, PetscInt k, PetscInt qorder, PetscFE *fem) 2184 { 2185 PetscFunctionBegin; 2186 PetscCall(PetscFECreate_Internal(comm, dim, Nc, DMPolytopeTypeSimpleShape(dim, isSimplex), NULL, k, qorder, PETSC_FALSE, fem)); 2187 PetscFunctionReturn(PETSC_SUCCESS); 2188 } 2189 2190 /*@ 2191 PetscFECreateLagrangeByCell - Create a `PetscFE` for the basic Lagrange space of degree k 2192 2193 Collective 2194 2195 Input Parameters: 2196 + comm - The MPI comm 2197 . dim - The spatial dimension 2198 . Nc - The number of components 2199 . ct - The celltype of the reference cell 2200 . k - The degree k of the space 2201 - qorder - The quadrature order or `PETSC_DETERMINE` to use `PetscSpace` polynomial degree 2202 2203 Output Parameter: 2204 . fem - The `PetscFE` object 2205 2206 Level: beginner 2207 2208 Note: 2209 For simplices, this element is the space of maximum polynomial degree k, otherwise it is a tensor product of 1D polynomials, each with maximal degree k. 2210 2211 .seealso: `PetscFECreateLagrange()`, `PetscFECreateDefault()`, `PetscFECreateByCell()`, `PetscFECreate()`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()` 2212 @*/ 2213 PetscErrorCode PetscFECreateLagrangeByCell(MPI_Comm comm, PetscInt dim, PetscInt Nc, DMPolytopeType ct, PetscInt k, PetscInt qorder, PetscFE *fem) 2214 { 2215 PetscFunctionBegin; 2216 PetscCall(PetscFECreate_Internal(comm, dim, Nc, ct, NULL, k, qorder, PETSC_FALSE, fem)); 2217 PetscFunctionReturn(PETSC_SUCCESS); 2218 } 2219 2220 /*@ 2221 PetscFELimitDegree - Copy a `PetscFE` but limit the degree to be in the given range 2222 2223 Collective 2224 2225 Input Parameters: 2226 + fe - The `PetscFE` 2227 . minDegree - The minimum degree, or `PETSC_DETERMINE` for no limit 2228 - maxDegree - The maximum degree, or `PETSC_DETERMINE` for no limit 2229 2230 Output Parameter: 2231 . newfe - The `PetscFE` object 2232 2233 Level: advanced 2234 2235 Note: 2236 This currently only works for Lagrange elements. 2237 2238 .seealso: `PetscFECreateLagrange()`, `PetscFECreateDefault()`, `PetscFECreateByCell()`, `PetscFECreate()`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()` 2239 @*/ 2240 PetscErrorCode PetscFELimitDegree(PetscFE fe, PetscInt minDegree, PetscInt maxDegree, PetscFE *newfe) 2241 { 2242 PetscDualSpace Q; 2243 PetscBool islag, issum; 2244 PetscInt oldk = 0, k; 2245 2246 PetscFunctionBegin; 2247 PetscCall(PetscFEGetDualSpace(fe, &Q)); 2248 PetscCall(PetscObjectTypeCompare((PetscObject)Q, PETSCDUALSPACELAGRANGE, &islag)); 2249 PetscCall(PetscObjectTypeCompare((PetscObject)Q, PETSCDUALSPACESUM, &issum)); 2250 if (islag) { 2251 PetscCall(PetscDualSpaceGetOrder(Q, &oldk)); 2252 } else if (issum) { 2253 PetscDualSpace subQ; 2254 2255 PetscCall(PetscDualSpaceSumGetSubspace(Q, 0, &subQ)); 2256 PetscCall(PetscDualSpaceGetOrder(subQ, &oldk)); 2257 } else { 2258 PetscCall(PetscObjectReference((PetscObject)fe)); 2259 *newfe = fe; 2260 PetscFunctionReturn(PETSC_SUCCESS); 2261 } 2262 k = oldk; 2263 if (minDegree >= 0) k = PetscMax(k, minDegree); 2264 if (maxDegree >= 0) k = PetscMin(k, maxDegree); 2265 if (k != oldk) { 2266 DM K; 2267 PetscSpace P; 2268 PetscQuadrature q; 2269 DMPolytopeType ct; 2270 PetscInt dim, Nc; 2271 2272 PetscCall(PetscFEGetBasisSpace(fe, &P)); 2273 PetscCall(PetscSpaceGetNumVariables(P, &dim)); 2274 PetscCall(PetscSpaceGetNumComponents(P, &Nc)); 2275 PetscCall(PetscDualSpaceGetDM(Q, &K)); 2276 PetscCall(DMPlexGetCellType(K, 0, &ct)); 2277 PetscCall(PetscFECreateLagrangeByCell(PetscObjectComm((PetscObject)fe), dim, Nc, ct, k, PETSC_DETERMINE, newfe)); 2278 PetscCall(PetscFEGetQuadrature(fe, &q)); 2279 PetscCall(PetscFESetQuadrature(*newfe, q)); 2280 } else { 2281 PetscCall(PetscObjectReference((PetscObject)fe)); 2282 *newfe = fe; 2283 } 2284 PetscFunctionReturn(PETSC_SUCCESS); 2285 } 2286 2287 /*@ 2288 PetscFECreateBrokenElement - Create a discontinuous version of the input `PetscFE` 2289 2290 Collective 2291 2292 Input Parameters: 2293 . cgfe - The continuous `PetscFE` object 2294 2295 Output Parameter: 2296 . dgfe - The discontinuous `PetscFE` object 2297 2298 Level: advanced 2299 2300 Note: 2301 This only works for Lagrange elements. 2302 2303 .seealso: `PetscFECreate()`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()`, `PetscFECreateLagrange()`, `PetscFECreateLagrangeByCell()`, `PetscDualSpaceLagrangeSetContinuity()` 2304 @*/ 2305 PetscErrorCode PetscFECreateBrokenElement(PetscFE cgfe, PetscFE *dgfe) 2306 { 2307 PetscSpace P; 2308 PetscDualSpace Q, dgQ; 2309 PetscQuadrature q, fq; 2310 PetscBool is_lagrange, is_sum; 2311 2312 PetscFunctionBegin; 2313 PetscCall(PetscFEGetBasisSpace(cgfe, &P)); 2314 PetscCall(PetscObjectReference((PetscObject)P)); 2315 PetscCall(PetscFEGetDualSpace(cgfe, &Q)); 2316 PetscCall(PetscObjectTypeCompare((PetscObject)Q, PETSCDUALSPACELAGRANGE, &is_lagrange)); 2317 PetscCall(PetscObjectTypeCompare((PetscObject)Q, PETSCDUALSPACESUM, &is_sum)); 2318 PetscCheck(is_lagrange || is_sum, PETSC_COMM_SELF, PETSC_ERR_SUP, "Can only create broken elements of Lagrange elements"); 2319 PetscCall(PetscDualSpaceDuplicate(Q, &dgQ)); 2320 PetscCall(PetscDualSpaceLagrangeSetContinuity(dgQ, PETSC_FALSE)); 2321 PetscCall(PetscDualSpaceSetUp(dgQ)); 2322 PetscCall(PetscFEGetQuadrature(cgfe, &q)); 2323 PetscCall(PetscObjectReference((PetscObject)q)); 2324 PetscCall(PetscFEGetFaceQuadrature(cgfe, &fq)); 2325 PetscCall(PetscObjectReference((PetscObject)fq)); 2326 PetscCall(PetscFECreateFromSpaces(P, dgQ, q, fq, dgfe)); 2327 PetscFunctionReturn(PETSC_SUCCESS); 2328 } 2329 2330 /*@ 2331 PetscFESetName - Names the `PetscFE` and its subobjects 2332 2333 Not Collective 2334 2335 Input Parameters: 2336 + fe - The `PetscFE` 2337 - name - The name 2338 2339 Level: intermediate 2340 2341 .seealso: `PetscFECreate()`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()` 2342 @*/ 2343 PetscErrorCode PetscFESetName(PetscFE fe, const char name[]) 2344 { 2345 PetscSpace P; 2346 PetscDualSpace Q; 2347 2348 PetscFunctionBegin; 2349 PetscCall(PetscFEGetBasisSpace(fe, &P)); 2350 PetscCall(PetscFEGetDualSpace(fe, &Q)); 2351 PetscCall(PetscObjectSetName((PetscObject)fe, name)); 2352 PetscCall(PetscObjectSetName((PetscObject)P, name)); 2353 PetscCall(PetscObjectSetName((PetscObject)Q, name)); 2354 PetscFunctionReturn(PETSC_SUCCESS); 2355 } 2356 2357 PetscErrorCode PetscFEEvaluateFieldJets_Internal(PetscDS ds, PetscInt Nf, PetscInt r, PetscInt q, PetscTabulation T[], PetscFEGeom *fegeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscScalar u[], PetscScalar u_x[], PetscScalar u_t[]) 2358 { 2359 PetscInt dOffset = 0, fOffset = 0, f, g; 2360 2361 for (f = 0; f < Nf; ++f) { 2362 PetscCheck(r < T[f]->Nr, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Replica number %" PetscInt_FMT " should be in [0, %" PetscInt_FMT ")", r, T[f]->Nr); 2363 PetscCheck(q < T[f]->Np, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point number %" PetscInt_FMT " should be in [0, %" PetscInt_FMT ")", q, T[f]->Np); 2364 PetscFE fe; 2365 const PetscInt k = ds->jetDegree[f]; 2366 const PetscInt cdim = T[f]->cdim; 2367 const PetscInt dE = fegeom->dimEmbed; 2368 const PetscInt Nq = T[f]->Np; 2369 const PetscInt Nbf = T[f]->Nb; 2370 const PetscInt Ncf = T[f]->Nc; 2371 const PetscReal *Bq = &T[f]->T[0][(r * Nq + q) * Nbf * Ncf]; 2372 const PetscReal *Dq = &T[f]->T[1][(r * Nq + q) * Nbf * Ncf * cdim]; 2373 const PetscReal *Hq = k > 1 ? &T[f]->T[2][(r * Nq + q) * Nbf * Ncf * cdim * cdim] : NULL; 2374 PetscInt hOffset = 0, b, c, d; 2375 2376 PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fe)); 2377 for (c = 0; c < Ncf; ++c) u[fOffset + c] = 0.0; 2378 for (d = 0; d < dE * Ncf; ++d) u_x[fOffset * dE + d] = 0.0; 2379 for (b = 0; b < Nbf; ++b) { 2380 for (c = 0; c < Ncf; ++c) { 2381 const PetscInt cidx = b * Ncf + c; 2382 2383 u[fOffset + c] += Bq[cidx] * coefficients[dOffset + b]; 2384 for (d = 0; d < cdim; ++d) u_x[(fOffset + c) * dE + d] += Dq[cidx * cdim + d] * coefficients[dOffset + b]; 2385 } 2386 } 2387 if (k > 1) { 2388 for (g = 0; g < Nf; ++g) hOffset += T[g]->Nc * dE; 2389 for (d = 0; d < dE * dE * Ncf; ++d) u_x[hOffset + fOffset * dE * dE + d] = 0.0; 2390 for (b = 0; b < Nbf; ++b) { 2391 for (c = 0; c < Ncf; ++c) { 2392 const PetscInt cidx = b * Ncf + c; 2393 2394 for (d = 0; d < cdim * cdim; ++d) u_x[hOffset + (fOffset + c) * dE * dE + d] += Hq[cidx * cdim * cdim + d] * coefficients[dOffset + b]; 2395 } 2396 } 2397 PetscCall(PetscFEPushforwardHessian(fe, fegeom, 1, &u_x[hOffset + fOffset * dE * dE])); 2398 } 2399 PetscCall(PetscFEPushforward(fe, fegeom, 1, &u[fOffset])); 2400 PetscCall(PetscFEPushforwardGradient(fe, fegeom, 1, &u_x[fOffset * dE])); 2401 if (u_t) { 2402 for (c = 0; c < Ncf; ++c) u_t[fOffset + c] = 0.0; 2403 for (b = 0; b < Nbf; ++b) { 2404 for (c = 0; c < Ncf; ++c) { 2405 const PetscInt cidx = b * Ncf + c; 2406 2407 u_t[fOffset + c] += Bq[cidx] * coefficients_t[dOffset + b]; 2408 } 2409 } 2410 PetscCall(PetscFEPushforward(fe, fegeom, 1, &u_t[fOffset])); 2411 } 2412 fOffset += Ncf; 2413 dOffset += Nbf; 2414 } 2415 return PETSC_SUCCESS; 2416 } 2417 2418 PetscErrorCode PetscFEEvaluateFieldJets_Hybrid_Internal(PetscDS ds, PetscInt Nf, PetscInt rc, PetscInt qc, PetscTabulation Tab[], const PetscInt rf[], const PetscInt qf[], PetscTabulation Tabf[], PetscFEGeom *fegeom, PetscFEGeom *fegeomNbr, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscScalar u[], PetscScalar u_x[], PetscScalar u_t[]) 2419 { 2420 PetscInt dOffset = 0, fOffset = 0, f; 2421 2422 /* f is the field number in the DS */ 2423 for (f = 0; f < Nf; ++f) { 2424 PetscBool isCohesive; 2425 PetscInt Ns, s; 2426 2427 if (!Tab[f]) continue; 2428 PetscCall(PetscDSGetCohesive(ds, f, &isCohesive)); 2429 Ns = isCohesive ? 1 : 2; 2430 { 2431 PetscTabulation T = isCohesive ? Tab[f] : Tabf[f]; 2432 PetscFE fe = (PetscFE)ds->disc[f]; 2433 const PetscInt dEt = T->cdim; 2434 const PetscInt dE = fegeom->dimEmbed; 2435 const PetscInt Nq = T->Np; 2436 const PetscInt Nbf = T->Nb; 2437 const PetscInt Ncf = T->Nc; 2438 2439 for (s = 0; s < Ns; ++s) { 2440 const PetscInt r = isCohesive ? rc : rf[s]; 2441 const PetscInt q = isCohesive ? qc : qf[s]; 2442 const PetscReal *Bq = &T->T[0][(r * Nq + q) * Nbf * Ncf]; 2443 const PetscReal *Dq = &T->T[1][(r * Nq + q) * Nbf * Ncf * dEt]; 2444 PetscInt b, c, d; 2445 2446 PetscCheck(r < T->Nr, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field %" PetscInt_FMT " Side %" PetscInt_FMT " Replica number %" PetscInt_FMT " should be in [0, %" PetscInt_FMT ")", f, s, r, T->Nr); 2447 PetscCheck(q < T->Np, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field %" PetscInt_FMT " Side %" PetscInt_FMT " Point number %" PetscInt_FMT " should be in [0, %" PetscInt_FMT ")", f, s, q, T->Np); 2448 for (c = 0; c < Ncf; ++c) u[fOffset + c] = 0.0; 2449 for (d = 0; d < dE * Ncf; ++d) u_x[fOffset * dE + d] = 0.0; 2450 for (b = 0; b < Nbf; ++b) { 2451 for (c = 0; c < Ncf; ++c) { 2452 const PetscInt cidx = b * Ncf + c; 2453 2454 u[fOffset + c] += Bq[cidx] * coefficients[dOffset + b]; 2455 for (d = 0; d < dEt; ++d) u_x[(fOffset + c) * dE + d] += Dq[cidx * dEt + d] * coefficients[dOffset + b]; 2456 } 2457 } 2458 PetscCall(PetscFEPushforward(fe, isCohesive ? fegeom : &fegeomNbr[s], 1, &u[fOffset])); 2459 PetscCall(PetscFEPushforwardGradient(fe, isCohesive ? fegeom : &fegeomNbr[s], 1, &u_x[fOffset * dE])); 2460 if (u_t) { 2461 for (c = 0; c < Ncf; ++c) u_t[fOffset + c] = 0.0; 2462 for (b = 0; b < Nbf; ++b) { 2463 for (c = 0; c < Ncf; ++c) { 2464 const PetscInt cidx = b * Ncf + c; 2465 2466 u_t[fOffset + c] += Bq[cidx] * coefficients_t[dOffset + b]; 2467 } 2468 } 2469 PetscCall(PetscFEPushforward(fe, fegeom, 1, &u_t[fOffset])); 2470 } 2471 fOffset += Ncf; 2472 dOffset += Nbf; 2473 } 2474 } 2475 } 2476 return PETSC_SUCCESS; 2477 } 2478 2479 PetscErrorCode PetscFEEvaluateFaceFields_Internal(PetscDS prob, PetscInt field, PetscInt faceLoc, const PetscScalar coefficients[], PetscScalar u[]) 2480 { 2481 PetscFE fe; 2482 PetscTabulation Tc; 2483 PetscInt b, c; 2484 2485 if (!prob) return PETSC_SUCCESS; 2486 PetscCall(PetscDSGetDiscretization(prob, field, (PetscObject *)&fe)); 2487 PetscCall(PetscFEGetFaceCentroidTabulation(fe, &Tc)); 2488 { 2489 const PetscReal *faceBasis = Tc->T[0]; 2490 const PetscInt Nb = Tc->Nb; 2491 const PetscInt Nc = Tc->Nc; 2492 2493 for (c = 0; c < Nc; ++c) u[c] = 0.0; 2494 for (b = 0; b < Nb; ++b) { 2495 for (c = 0; c < Nc; ++c) u[c] += coefficients[b] * faceBasis[(faceLoc * Nb + b) * Nc + c]; 2496 } 2497 } 2498 return PETSC_SUCCESS; 2499 } 2500 2501 PetscErrorCode PetscFEUpdateElementVec_Internal(PetscFE fe, PetscTabulation T, PetscInt r, PetscScalar tmpBasis[], PetscScalar tmpBasisDer[], PetscInt e, PetscFEGeom *fegeom, PetscScalar f0[], PetscScalar f1[], PetscScalar elemVec[]) 2502 { 2503 PetscFEGeom pgeom; 2504 const PetscInt dEt = T->cdim; 2505 const PetscInt dE = fegeom->dimEmbed; 2506 const PetscInt Nq = T->Np; 2507 const PetscInt Nb = T->Nb; 2508 const PetscInt Nc = T->Nc; 2509 const PetscReal *basis = &T->T[0][r * Nq * Nb * Nc]; 2510 const PetscReal *basisDer = &T->T[1][r * Nq * Nb * Nc * dEt]; 2511 PetscInt q, b, c, d; 2512 2513 for (q = 0; q < Nq; ++q) { 2514 for (b = 0; b < Nb; ++b) { 2515 for (c = 0; c < Nc; ++c) { 2516 const PetscInt bcidx = b * Nc + c; 2517 2518 tmpBasis[bcidx] = basis[q * Nb * Nc + bcidx]; 2519 for (d = 0; d < dEt; ++d) tmpBasisDer[bcidx * dE + d] = basisDer[q * Nb * Nc * dEt + bcidx * dEt + d]; 2520 for (d = dEt; d < dE; ++d) tmpBasisDer[bcidx * dE + d] = 0.0; 2521 } 2522 } 2523 PetscCall(PetscFEGeomGetCellPoint(fegeom, e, q, &pgeom)); 2524 PetscCall(PetscFEPushforward(fe, &pgeom, Nb, tmpBasis)); 2525 PetscCall(PetscFEPushforwardGradient(fe, &pgeom, Nb, tmpBasisDer)); 2526 for (b = 0; b < Nb; ++b) { 2527 for (c = 0; c < Nc; ++c) { 2528 const PetscInt bcidx = b * Nc + c; 2529 const PetscInt qcidx = q * Nc + c; 2530 2531 elemVec[b] += tmpBasis[bcidx] * f0[qcidx]; 2532 for (d = 0; d < dE; ++d) elemVec[b] += tmpBasisDer[bcidx * dE + d] * f1[qcidx * dE + d]; 2533 } 2534 } 2535 } 2536 return PETSC_SUCCESS; 2537 } 2538 2539 PetscErrorCode PetscFEUpdateElementVec_Hybrid_Internal(PetscFE fe, PetscTabulation T, PetscInt r, PetscInt side, PetscScalar tmpBasis[], PetscScalar tmpBasisDer[], PetscFEGeom *fegeom, PetscScalar f0[], PetscScalar f1[], PetscScalar elemVec[]) 2540 { 2541 const PetscInt dE = T->cdim; 2542 const PetscInt Nq = T->Np; 2543 const PetscInt Nb = T->Nb; 2544 const PetscInt Nc = T->Nc; 2545 const PetscReal *basis = &T->T[0][r * Nq * Nb * Nc]; 2546 const PetscReal *basisDer = &T->T[1][r * Nq * Nb * Nc * dE]; 2547 2548 for (PetscInt q = 0; q < Nq; ++q) { 2549 for (PetscInt b = 0; b < Nb; ++b) { 2550 for (PetscInt c = 0; c < Nc; ++c) { 2551 const PetscInt bcidx = b * Nc + c; 2552 2553 tmpBasis[bcidx] = basis[q * Nb * Nc + bcidx]; 2554 for (PetscInt d = 0; d < dE; ++d) tmpBasisDer[bcidx * dE + d] = basisDer[q * Nb * Nc * dE + bcidx * dE + d]; 2555 } 2556 } 2557 PetscCall(PetscFEPushforward(fe, fegeom, Nb, tmpBasis)); 2558 // TODO This is currently broken since we do not pull the geometry down to the lower dimension 2559 // PetscCall(PetscFEPushforwardGradient(fe, fegeom, Nb, tmpBasisDer)); 2560 if (side == 2) { 2561 // Integrating over whole cohesive cell, so insert for both sides 2562 for (PetscInt s = 0; s < 2; ++s) { 2563 for (PetscInt b = 0; b < Nb; ++b) { 2564 for (PetscInt c = 0; c < Nc; ++c) { 2565 const PetscInt bcidx = b * Nc + c; 2566 const PetscInt qcidx = (q * 2 + s) * Nc + c; 2567 2568 elemVec[Nb * s + b] += tmpBasis[bcidx] * f0[qcidx]; 2569 for (PetscInt d = 0; d < dE; ++d) elemVec[Nb * s + b] += tmpBasisDer[bcidx * dE + d] * f1[qcidx * dE + d]; 2570 } 2571 } 2572 } 2573 } else { 2574 // Integrating over endcaps of cohesive cell, so insert for correct side 2575 for (PetscInt b = 0; b < Nb; ++b) { 2576 for (PetscInt c = 0; c < Nc; ++c) { 2577 const PetscInt bcidx = b * Nc + c; 2578 const PetscInt qcidx = q * Nc + c; 2579 2580 elemVec[Nb * side + b] += tmpBasis[bcidx] * f0[qcidx]; 2581 for (PetscInt d = 0; d < dE; ++d) elemVec[Nb * side + b] += tmpBasisDer[bcidx * dE + d] * f1[qcidx * dE + d]; 2582 } 2583 } 2584 } 2585 } 2586 return PETSC_SUCCESS; 2587 } 2588 2589 #define petsc_elemmat_kernel_g1(_NbI, _NcI, _NbJ, _NcJ, _dE) \ 2590 do { \ 2591 for (PetscInt fc = 0; fc < (_NcI); ++fc) { \ 2592 for (PetscInt gc = 0; gc < (_NcJ); ++gc) { \ 2593 const PetscScalar *G = g1 + (fc * (_NcJ) + gc) * _dE; \ 2594 for (PetscInt f = 0; f < (_NbI); ++f) { \ 2595 const PetscScalar tBIv = tmpBasisI[f * (_NcI) + fc]; \ 2596 for (PetscInt g = 0; g < (_NbJ); ++g) { \ 2597 const PetscScalar *tBDJ = tmpBasisDerJ + (g * (_NcJ) + gc) * (_dE); \ 2598 PetscScalar s = 0.0; \ 2599 for (PetscInt df = 0; df < _dE; ++df) s += G[df] * tBDJ[df]; \ 2600 elemMat[(offsetI + f) * totDim + (offsetJ + g)] += s * tBIv; \ 2601 } \ 2602 } \ 2603 } \ 2604 } \ 2605 } while (0) 2606 2607 #define petsc_elemmat_kernel_g2(_NbI, _NcI, _NbJ, _NcJ, _dE) \ 2608 do { \ 2609 for (PetscInt gc = 0; gc < (_NcJ); ++gc) { \ 2610 for (PetscInt fc = 0; fc < (_NcI); ++fc) { \ 2611 const PetscScalar *G = g2 + (fc * (_NcJ) + gc) * _dE; \ 2612 for (PetscInt g = 0; g < (_NbJ); ++g) { \ 2613 const PetscScalar tBJv = tmpBasisJ[g * (_NcJ) + gc]; \ 2614 for (PetscInt f = 0; f < (_NbI); ++f) { \ 2615 const PetscScalar *tBDI = tmpBasisDerI + (f * (_NcI) + fc) * (_dE); \ 2616 PetscScalar s = 0.0; \ 2617 for (PetscInt df = 0; df < _dE; ++df) s += tBDI[df] * G[df]; \ 2618 elemMat[(offsetI + f) * totDim + (offsetJ + g)] += s * tBJv; \ 2619 } \ 2620 } \ 2621 } \ 2622 } \ 2623 } while (0) 2624 2625 #define petsc_elemmat_kernel_g3(_NbI, _NcI, _NbJ, _NcJ, _dE) \ 2626 do { \ 2627 for (PetscInt fc = 0; fc < (_NcI); ++fc) { \ 2628 for (PetscInt gc = 0; gc < (_NcJ); ++gc) { \ 2629 const PetscScalar *G = g3 + (fc * (_NcJ) + gc) * (_dE) * (_dE); \ 2630 for (PetscInt f = 0; f < (_NbI); ++f) { \ 2631 const PetscScalar *tBDI = tmpBasisDerI + (f * (_NcI) + fc) * (_dE); \ 2632 for (PetscInt g = 0; g < (_NbJ); ++g) { \ 2633 PetscScalar s = 0.0; \ 2634 const PetscScalar *tBDJ = tmpBasisDerJ + (g * (_NcJ) + gc) * (_dE); \ 2635 for (PetscInt df = 0; df < (_dE); ++df) { \ 2636 for (PetscInt dg = 0; dg < (_dE); ++dg) s += tBDI[df] * G[df * (_dE) + dg] * tBDJ[dg]; \ 2637 } \ 2638 elemMat[(offsetI + f) * totDim + (offsetJ + g)] += s; \ 2639 } \ 2640 } \ 2641 } \ 2642 } \ 2643 } while (0) 2644 2645 PetscErrorCode PetscFEUpdateElementMat_Internal(PetscFE feI, PetscFE feJ, PetscInt r, PetscInt q, PetscTabulation TI, PetscScalar tmpBasisI[], PetscScalar tmpBasisDerI[], PetscTabulation TJ, PetscScalar tmpBasisJ[], PetscScalar tmpBasisDerJ[], PetscFEGeom *fegeom, const PetscScalar g0[], const PetscScalar g1[], const PetscScalar g2[], const PetscScalar g3[], PetscInt totDim, PetscInt offsetI, PetscInt offsetJ, PetscScalar elemMat[]) 2646 { 2647 const PetscInt cdim = TI->cdim; 2648 const PetscInt dE = fegeom->dimEmbed; 2649 const PetscInt NqI = TI->Np; 2650 const PetscInt NbI = TI->Nb; 2651 const PetscInt NcI = TI->Nc; 2652 const PetscReal *basisI = &TI->T[0][(r * NqI + q) * NbI * NcI]; 2653 const PetscReal *basisDerI = &TI->T[1][(r * NqI + q) * NbI * NcI * cdim]; 2654 const PetscInt NqJ = TJ->Np; 2655 const PetscInt NbJ = TJ->Nb; 2656 const PetscInt NcJ = TJ->Nc; 2657 const PetscReal *basisJ = &TJ->T[0][(r * NqJ + q) * NbJ * NcJ]; 2658 const PetscReal *basisDerJ = &TJ->T[1][(r * NqJ + q) * NbJ * NcJ * cdim]; 2659 2660 for (PetscInt f = 0; f < NbI; ++f) { 2661 for (PetscInt fc = 0; fc < NcI; ++fc) { 2662 const PetscInt fidx = f * NcI + fc; /* Test function basis index */ 2663 2664 tmpBasisI[fidx] = basisI[fidx]; 2665 for (PetscInt df = 0; df < cdim; ++df) tmpBasisDerI[fidx * dE + df] = basisDerI[fidx * cdim + df]; 2666 } 2667 } 2668 PetscCall(PetscFEPushforward(feI, fegeom, NbI, tmpBasisI)); 2669 PetscCall(PetscFEPushforwardGradient(feI, fegeom, NbI, tmpBasisDerI)); 2670 if (feI != feJ) { 2671 for (PetscInt g = 0; g < NbJ; ++g) { 2672 for (PetscInt gc = 0; gc < NcJ; ++gc) { 2673 const PetscInt gidx = g * NcJ + gc; /* Trial function basis index */ 2674 2675 tmpBasisJ[gidx] = basisJ[gidx]; 2676 for (PetscInt dg = 0; dg < cdim; ++dg) tmpBasisDerJ[gidx * dE + dg] = basisDerJ[gidx * cdim + dg]; 2677 } 2678 } 2679 PetscCall(PetscFEPushforward(feJ, fegeom, NbJ, tmpBasisJ)); 2680 PetscCall(PetscFEPushforwardGradient(feJ, fegeom, NbJ, tmpBasisDerJ)); 2681 } else { 2682 tmpBasisJ = tmpBasisI; 2683 tmpBasisDerJ = tmpBasisDerI; 2684 } 2685 if (PetscUnlikely(g0)) { 2686 for (PetscInt f = 0; f < NbI; ++f) { 2687 const PetscInt i = offsetI + f; /* Element matrix row */ 2688 2689 for (PetscInt fc = 0; fc < NcI; ++fc) { 2690 const PetscScalar bI = tmpBasisI[f * NcI + fc]; /* Test function basis value */ 2691 2692 for (PetscInt g = 0; g < NbJ; ++g) { 2693 const PetscInt j = offsetJ + g; /* Element matrix column */ 2694 const PetscInt fOff = i * totDim + j; 2695 2696 for (PetscInt gc = 0; gc < NcJ; ++gc) elemMat[fOff] += bI * g0[fc * NcJ + gc] * tmpBasisJ[g * NcJ + gc]; 2697 } 2698 } 2699 } 2700 } 2701 if (PetscUnlikely(g1)) { 2702 #if 1 2703 if (dE == 2) { 2704 petsc_elemmat_kernel_g1(NbI, NcI, NbJ, NcJ, 2); 2705 } else if (dE == 3) { 2706 petsc_elemmat_kernel_g1(NbI, NcI, NbJ, NcJ, 3); 2707 } else { 2708 petsc_elemmat_kernel_g1(NbI, NcI, NbJ, NcJ, dE); 2709 } 2710 #else 2711 for (PetscInt f = 0; f < NbI; ++f) { 2712 const PetscInt i = offsetI + f; /* Element matrix row */ 2713 2714 for (PetscInt fc = 0; fc < NcI; ++fc) { 2715 const PetscScalar bI = tmpBasisI[f * NcI + fc]; /* Test function basis value */ 2716 2717 for (PetscInt g = 0; g < NbJ; ++g) { 2718 const PetscInt j = offsetJ + g; /* Element matrix column */ 2719 const PetscInt fOff = i * totDim + j; 2720 2721 for (PetscInt gc = 0; gc < NcJ; ++gc) { 2722 const PetscInt gidx = g * NcJ + gc; /* Trial function basis index */ 2723 2724 for (PetscInt df = 0; df < dE; ++df) elemMat[fOff] += bI * g1[(fc * NcJ + gc) * dE + df] * tmpBasisDerJ[gidx * dE + df]; 2725 } 2726 } 2727 } 2728 } 2729 #endif 2730 } 2731 if (PetscUnlikely(g2)) { 2732 #if 1 2733 if (dE == 2) { 2734 petsc_elemmat_kernel_g2(NbI, NcI, NbJ, NcJ, 2); 2735 } else if (dE == 3) { 2736 petsc_elemmat_kernel_g2(NbI, NcI, NbJ, NcJ, 3); 2737 } else { 2738 petsc_elemmat_kernel_g2(NbI, NcI, NbJ, NcJ, dE); 2739 } 2740 #else 2741 for (PetscInt g = 0; g < NbJ; ++g) { 2742 const PetscInt j = offsetJ + g; /* Element matrix column */ 2743 2744 for (PetscInt gc = 0; gc < NcJ; ++gc) { 2745 const PetscScalar bJ = tmpBasisJ[g * NcJ + gc]; /* Trial function basis value */ 2746 2747 for (PetscInt f = 0; f < NbI; ++f) { 2748 const PetscInt i = offsetI + f; /* Element matrix row */ 2749 const PetscInt fOff = i * totDim + j; 2750 2751 for (PetscInt fc = 0; fc < NcI; ++fc) { 2752 const PetscInt fidx = f * NcI + fc; /* Test function basis index */ 2753 2754 for (PetscInt df = 0; df < dE; ++df) elemMat[fOff] += tmpBasisDerI[fidx * dE + df] * g2[(fc * NcJ + gc) * dE + df] * bJ; 2755 } 2756 } 2757 } 2758 } 2759 #endif 2760 } 2761 if (PetscUnlikely(g3)) { 2762 #if 1 2763 if (dE == 2) { 2764 petsc_elemmat_kernel_g3(NbI, NcI, NbJ, NcJ, 2); 2765 } else if (dE == 3) { 2766 petsc_elemmat_kernel_g3(NbI, NcI, NbJ, NcJ, 3); 2767 } else { 2768 petsc_elemmat_kernel_g3(NbI, NcI, NbJ, NcJ, dE); 2769 } 2770 #else 2771 for (PetscInt f = 0; f < NbI; ++f) { 2772 const PetscInt i = offsetI + f; /* Element matrix row */ 2773 2774 for (PetscInt fc = 0; fc < NcI; ++fc) { 2775 const PetscInt fidx = f * NcI + fc; /* Test function basis index */ 2776 2777 for (PetscInt g = 0; g < NbJ; ++g) { 2778 const PetscInt j = offsetJ + g; /* Element matrix column */ 2779 const PetscInt fOff = i * totDim + j; 2780 2781 for (PetscInt gc = 0; gc < NcJ; ++gc) { 2782 const PetscInt gidx = g * NcJ + gc; /* Trial function basis index */ 2783 2784 for (PetscInt df = 0; df < dE; ++df) { 2785 for (PetscInt dg = 0; dg < dE; ++dg) elemMat[fOff] += tmpBasisDerI[fidx * dE + df] * g3[((fc * NcJ + gc) * dE + df) * dE + dg] * tmpBasisDerJ[gidx * dE + dg]; 2786 } 2787 } 2788 } 2789 } 2790 } 2791 #endif 2792 } 2793 return PETSC_SUCCESS; 2794 } 2795 2796 #undef petsc_elemmat_kernel_g1 2797 #undef petsc_elemmat_kernel_g2 2798 #undef petsc_elemmat_kernel_g3 2799 2800 PetscErrorCode PetscFEUpdateElementMat_Hybrid_Internal(PetscFE feI, PetscBool isHybridI, PetscFE feJ, PetscBool isHybridJ, PetscInt r, PetscInt s, PetscInt t, PetscInt q, PetscTabulation TI, PetscScalar tmpBasisI[], PetscScalar tmpBasisDerI[], PetscTabulation TJ, PetscScalar tmpBasisJ[], PetscScalar tmpBasisDerJ[], PetscFEGeom *fegeom, const PetscScalar g0[], const PetscScalar g1[], const PetscScalar g2[], const PetscScalar g3[], PetscInt eOffset, PetscInt totDim, PetscInt offsetI, PetscInt offsetJ, PetscScalar elemMat[]) 2801 { 2802 const PetscInt dE = TI->cdim; 2803 const PetscInt NqI = TI->Np; 2804 const PetscInt NbI = TI->Nb; 2805 const PetscInt NcI = TI->Nc; 2806 const PetscReal *basisI = &TI->T[0][(r * NqI + q) * NbI * NcI]; 2807 const PetscReal *basisDerI = &TI->T[1][(r * NqI + q) * NbI * NcI * dE]; 2808 const PetscInt NqJ = TJ->Np; 2809 const PetscInt NbJ = TJ->Nb; 2810 const PetscInt NcJ = TJ->Nc; 2811 const PetscReal *basisJ = &TJ->T[0][(r * NqJ + q) * NbJ * NcJ]; 2812 const PetscReal *basisDerJ = &TJ->T[1][(r * NqJ + q) * NbJ * NcJ * dE]; 2813 const PetscInt so = isHybridI ? 0 : s; 2814 const PetscInt to = isHybridJ ? 0 : t; 2815 PetscInt f, fc, g, gc, df, dg; 2816 2817 for (f = 0; f < NbI; ++f) { 2818 for (fc = 0; fc < NcI; ++fc) { 2819 const PetscInt fidx = f * NcI + fc; /* Test function basis index */ 2820 2821 tmpBasisI[fidx] = basisI[fidx]; 2822 for (df = 0; df < dE; ++df) tmpBasisDerI[fidx * dE + df] = basisDerI[fidx * dE + df]; 2823 } 2824 } 2825 PetscCall(PetscFEPushforward(feI, fegeom, NbI, tmpBasisI)); 2826 PetscCall(PetscFEPushforwardGradient(feI, fegeom, NbI, tmpBasisDerI)); 2827 for (g = 0; g < NbJ; ++g) { 2828 for (gc = 0; gc < NcJ; ++gc) { 2829 const PetscInt gidx = g * NcJ + gc; /* Trial function basis index */ 2830 2831 tmpBasisJ[gidx] = basisJ[gidx]; 2832 for (dg = 0; dg < dE; ++dg) tmpBasisDerJ[gidx * dE + dg] = basisDerJ[gidx * dE + dg]; 2833 } 2834 } 2835 PetscCall(PetscFEPushforward(feJ, fegeom, NbJ, tmpBasisJ)); 2836 // TODO This is currently broken since we do not pull the geometry down to the lower dimension 2837 // PetscCall(PetscFEPushforwardGradient(feJ, fegeom, NbJ, tmpBasisDerJ)); 2838 for (f = 0; f < NbI; ++f) { 2839 for (fc = 0; fc < NcI; ++fc) { 2840 const PetscInt fidx = f * NcI + fc; /* Test function basis index */ 2841 const PetscInt i = offsetI + NbI * so + f; /* Element matrix row */ 2842 for (g = 0; g < NbJ; ++g) { 2843 for (gc = 0; gc < NcJ; ++gc) { 2844 const PetscInt gidx = g * NcJ + gc; /* Trial function basis index */ 2845 const PetscInt j = offsetJ + NbJ * to + g; /* Element matrix column */ 2846 const PetscInt fOff = eOffset + i * totDim + j; 2847 2848 elemMat[fOff] += tmpBasisI[fidx] * g0[fc * NcJ + gc] * tmpBasisJ[gidx]; 2849 for (df = 0; df < dE; ++df) { 2850 elemMat[fOff] += tmpBasisI[fidx] * g1[(fc * NcJ + gc) * dE + df] * tmpBasisDerJ[gidx * dE + df]; 2851 elemMat[fOff] += tmpBasisDerI[fidx * dE + df] * g2[(fc * NcJ + gc) * dE + df] * tmpBasisJ[gidx]; 2852 for (dg = 0; dg < dE; ++dg) elemMat[fOff] += tmpBasisDerI[fidx * dE + df] * g3[((fc * NcJ + gc) * dE + df) * dE + dg] * tmpBasisDerJ[gidx * dE + dg]; 2853 } 2854 } 2855 } 2856 } 2857 } 2858 return PETSC_SUCCESS; 2859 } 2860 2861 PetscErrorCode PetscFECreateCellGeometry(PetscFE fe, PetscQuadrature quad, PetscFEGeom *cgeom) 2862 { 2863 PetscDualSpace dsp; 2864 DM dm; 2865 PetscQuadrature quadDef; 2866 PetscInt dim, cdim, Nq; 2867 2868 PetscFunctionBegin; 2869 PetscCall(PetscFEGetDualSpace(fe, &dsp)); 2870 PetscCall(PetscDualSpaceGetDM(dsp, &dm)); 2871 PetscCall(DMGetDimension(dm, &dim)); 2872 PetscCall(DMGetCoordinateDim(dm, &cdim)); 2873 PetscCall(PetscFEGetQuadrature(fe, &quadDef)); 2874 quad = quad ? quad : quadDef; 2875 PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL)); 2876 PetscCall(PetscMalloc1(Nq * cdim, &cgeom->v)); 2877 PetscCall(PetscMalloc1(Nq * cdim * cdim, &cgeom->J)); 2878 PetscCall(PetscMalloc1(Nq * cdim * cdim, &cgeom->invJ)); 2879 PetscCall(PetscMalloc1(Nq, &cgeom->detJ)); 2880 cgeom->dim = dim; 2881 cgeom->dimEmbed = cdim; 2882 cgeom->numCells = 1; 2883 cgeom->numPoints = Nq; 2884 PetscCall(DMPlexComputeCellGeometryFEM(dm, 0, quad, cgeom->v, cgeom->J, cgeom->invJ, cgeom->detJ)); 2885 PetscFunctionReturn(PETSC_SUCCESS); 2886 } 2887 2888 PetscErrorCode PetscFEDestroyCellGeometry(PetscFE fe, PetscFEGeom *cgeom) 2889 { 2890 PetscFunctionBegin; 2891 PetscCall(PetscFree(cgeom->v)); 2892 PetscCall(PetscFree(cgeom->J)); 2893 PetscCall(PetscFree(cgeom->invJ)); 2894 PetscCall(PetscFree(cgeom->detJ)); 2895 PetscFunctionReturn(PETSC_SUCCESS); 2896 } 2897 2898 #if 0 2899 PetscErrorCode PetscFEUpdateElementMat_Internal_SparseIndices(PetscTabulation TI, PetscTabulation TJ, PetscInt dimEmbed, const PetscInt g0[], const PetscInt g1[], const PetscInt g2[], const PetscInt g3[], PetscInt totDim, PetscInt offsetI, PetscInt offsetJ, PetscInt *n_g0, PetscInt **g0_idxs_out, PetscInt *n_g1, PetscInt **g1_idxs_out, PetscInt *n_g2, PetscInt **g2_idxs_out, PetscInt *n_g3, PetscInt **g3_idxs_out) 2900 { 2901 const PetscInt dE = dimEmbed; 2902 const PetscInt NbI = TI->Nb; 2903 const PetscInt NcI = TI->Nc; 2904 const PetscInt NbJ = TJ->Nb; 2905 const PetscInt NcJ = TJ->Nc; 2906 PetscBool has_g0 = g0 ? PETSC_TRUE : PETSC_FALSE; 2907 PetscBool has_g1 = g1 ? PETSC_TRUE : PETSC_FALSE; 2908 PetscBool has_g2 = g2 ? PETSC_TRUE : PETSC_FALSE; 2909 PetscBool has_g3 = g3 ? PETSC_TRUE : PETSC_FALSE; 2910 PetscInt *g0_idxs = NULL, *g1_idxs = NULL, *g2_idxs = NULL, *g3_idxs = NULL; 2911 PetscInt g0_i, g1_i, g2_i, g3_i; 2912 2913 PetscFunctionBegin; 2914 g0_i = g1_i = g2_i = g3_i = 0; 2915 if (has_g0) 2916 for (PetscInt i = 0; i < NcI * NcJ; i++) 2917 if (g0[i]) g0_i += NbI * NbJ; 2918 if (has_g1) 2919 for (PetscInt i = 0; i < NcI * NcJ * dE; i++) 2920 if (g1[i]) g1_i += NbI * NbJ; 2921 if (has_g2) 2922 for (PetscInt i = 0; i < NcI * NcJ * dE; i++) 2923 if (g2[i]) g2_i += NbI * NbJ; 2924 if (has_g3) 2925 for (PetscInt i = 0; i < NcI * NcJ * dE * dE; i++) 2926 if (g3[i]) g3_i += NbI * NbJ; 2927 if (g0_i == NbI * NbJ * NcI * NcJ) g0_i = 0; 2928 if (g1_i == NbI * NbJ * NcI * NcJ * dE) g1_i = 0; 2929 if (g2_i == NbI * NbJ * NcI * NcJ * dE) g2_i = 0; 2930 if (g3_i == NbI * NbJ * NcI * NcJ * dE * dE) g3_i = 0; 2931 has_g0 = g0_i ? PETSC_TRUE : PETSC_FALSE; 2932 has_g1 = g1_i ? PETSC_TRUE : PETSC_FALSE; 2933 has_g2 = g2_i ? PETSC_TRUE : PETSC_FALSE; 2934 has_g3 = g3_i ? PETSC_TRUE : PETSC_FALSE; 2935 if (has_g0) PetscCall(PetscMalloc1(4 * g0_i, &g0_idxs)); 2936 if (has_g1) PetscCall(PetscMalloc1(4 * g1_i, &g1_idxs)); 2937 if (has_g2) PetscCall(PetscMalloc1(4 * g2_i, &g2_idxs)); 2938 if (has_g3) PetscCall(PetscMalloc1(4 * g3_i, &g3_idxs)); 2939 g0_i = g1_i = g2_i = g3_i = 0; 2940 2941 for (PetscInt f = 0; f < NbI; ++f) { 2942 const PetscInt i = offsetI + f; /* Element matrix row */ 2943 for (PetscInt fc = 0; fc < NcI; ++fc) { 2944 const PetscInt fidx = f * NcI + fc; /* Test function basis index */ 2945 2946 for (PetscInt g = 0; g < NbJ; ++g) { 2947 const PetscInt j = offsetJ + g; /* Element matrix column */ 2948 const PetscInt fOff = i * totDim + j; 2949 for (PetscInt gc = 0; gc < NcJ; ++gc) { 2950 const PetscInt gidx = g * NcJ + gc; /* Trial function basis index */ 2951 2952 if (has_g0) { 2953 if (g0[fc * NcJ + gc]) { 2954 g0_idxs[4 * g0_i + 0] = fidx; 2955 g0_idxs[4 * g0_i + 1] = fc * NcJ + gc; 2956 g0_idxs[4 * g0_i + 2] = gidx; 2957 g0_idxs[4 * g0_i + 3] = fOff; 2958 g0_i++; 2959 } 2960 } 2961 2962 for (PetscInt df = 0; df < dE; ++df) { 2963 if (has_g1) { 2964 if (g1[(fc * NcJ + gc) * dE + df]) { 2965 g1_idxs[4 * g1_i + 0] = fidx; 2966 g1_idxs[4 * g1_i + 1] = (fc * NcJ + gc) * dE + df; 2967 g1_idxs[4 * g1_i + 2] = gidx * dE + df; 2968 g1_idxs[4 * g1_i + 3] = fOff; 2969 g1_i++; 2970 } 2971 } 2972 if (has_g2) { 2973 if (g2[(fc * NcJ + gc) * dE + df]) { 2974 g2_idxs[4 * g2_i + 0] = fidx * dE + df; 2975 g2_idxs[4 * g2_i + 1] = (fc * NcJ + gc) * dE + df; 2976 g2_idxs[4 * g2_i + 2] = gidx; 2977 g2_idxs[4 * g2_i + 3] = fOff; 2978 g2_i++; 2979 } 2980 } 2981 if (has_g3) { 2982 for (PetscInt dg = 0; dg < dE; ++dg) { 2983 if (g3[((fc * NcJ + gc) * dE + df) * dE + dg]) { 2984 g3_idxs[4 * g3_i + 0] = fidx * dE + df; 2985 g3_idxs[4 * g3_i + 1] = ((fc * NcJ + gc) * dE + df) * dE + dg; 2986 g3_idxs[4 * g3_i + 2] = gidx * dE + dg; 2987 g3_idxs[4 * g3_i + 3] = fOff; 2988 g3_i++; 2989 } 2990 } 2991 } 2992 } 2993 } 2994 } 2995 } 2996 } 2997 *n_g0 = g0_i; 2998 *n_g1 = g1_i; 2999 *n_g2 = g2_i; 3000 *n_g3 = g3_i; 3001 3002 *g0_idxs_out = g0_idxs; 3003 *g1_idxs_out = g1_idxs; 3004 *g2_idxs_out = g2_idxs; 3005 *g3_idxs_out = g3_idxs; 3006 PetscFunctionReturn(PETSC_SUCCESS); 3007 } 3008 3009 //example HOW TO USE 3010 for (PetscInt i = 0; i < g0_sparse_n; i++) { 3011 PetscInt bM = g0_sparse_idxs[4 * i + 0]; 3012 PetscInt bN = g0_sparse_idxs[4 * i + 1]; 3013 PetscInt bK = g0_sparse_idxs[4 * i + 2]; 3014 PetscInt bO = g0_sparse_idxs[4 * i + 3]; 3015 elemMat[bO] += tmpBasisI[bM] * g0[bN] * tmpBasisJ[bK]; 3016 } 3017 #endif 3018