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 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 /*@C 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 /*@C 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 /*@C 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 167 - name - command line option name 168 169 Level: intermediate 170 171 .seealso: `PetscFE`, `PetscFEView()`, `PetscObjectViewFromOptions()`, `PetscFECreate()` 172 @*/ 173 PetscErrorCode PetscFEViewFromOptions(PetscFE A, 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 /*@C 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 iascii; 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, &iascii)); 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: `PetscFEV`, `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 /*@C 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 *fem = NULL; 352 PetscCall(PetscFEInitializePackage()); 353 354 PetscCall(PetscHeaderCreate(f, PETSCFE_CLASSID, "PetscFE", "Finite Element", "PetscFE", comm, PetscFEDestroy, PetscFEView)); 355 356 f->basisSpace = NULL; 357 f->dualSpace = NULL; 358 f->numComponents = 1; 359 f->subspaces = NULL; 360 f->invV = NULL; 361 f->T = NULL; 362 f->Tf = NULL; 363 f->Tc = NULL; 364 PetscCall(PetscArrayzero(&f->quadrature, 1)); 365 PetscCall(PetscArrayzero(&f->faceQuadrature, 1)); 366 f->blockSize = 0; 367 f->numBlocks = 1; 368 f->batchSize = 0; 369 f->numBatches = 1; 370 371 *fem = f; 372 PetscFunctionReturn(PETSC_SUCCESS); 373 } 374 375 /*@ 376 PetscFEGetSpatialDimension - Returns the spatial dimension of the element 377 378 Not Collective 379 380 Input Parameter: 381 . fem - The `PetscFE` object 382 383 Output Parameter: 384 . dim - The spatial dimension 385 386 Level: intermediate 387 388 .seealso: `PetscFE`, `PetscFECreate()` 389 @*/ 390 PetscErrorCode PetscFEGetSpatialDimension(PetscFE fem, PetscInt *dim) 391 { 392 DM dm; 393 394 PetscFunctionBegin; 395 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 396 PetscAssertPointer(dim, 2); 397 PetscCall(PetscDualSpaceGetDM(fem->dualSpace, &dm)); 398 PetscCall(DMGetDimension(dm, dim)); 399 PetscFunctionReturn(PETSC_SUCCESS); 400 } 401 402 /*@ 403 PetscFESetNumComponents - Sets the number of field components in the element 404 405 Not Collective 406 407 Input Parameters: 408 + fem - The `PetscFE` object 409 - comp - The number of field components 410 411 Level: intermediate 412 413 .seealso: `PetscFE`, `PetscFECreate()`, `PetscFEGetSpatialDimension()`, `PetscFEGetNumComponents()` 414 @*/ 415 PetscErrorCode PetscFESetNumComponents(PetscFE fem, PetscInt comp) 416 { 417 PetscFunctionBegin; 418 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 419 fem->numComponents = comp; 420 PetscFunctionReturn(PETSC_SUCCESS); 421 } 422 423 /*@ 424 PetscFEGetNumComponents - Returns the number of components in the element 425 426 Not Collective 427 428 Input Parameter: 429 . fem - The `PetscFE` object 430 431 Output Parameter: 432 . comp - The number of field components 433 434 Level: intermediate 435 436 .seealso: `PetscFE`, `PetscFECreate()`, `PetscFEGetSpatialDimension()` 437 @*/ 438 PetscErrorCode PetscFEGetNumComponents(PetscFE fem, PetscInt *comp) 439 { 440 PetscFunctionBegin; 441 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 442 PetscAssertPointer(comp, 2); 443 *comp = fem->numComponents; 444 PetscFunctionReturn(PETSC_SUCCESS); 445 } 446 447 /*@ 448 PetscFESetTileSizes - Sets the tile sizes for evaluation 449 450 Not Collective 451 452 Input Parameters: 453 + fem - The `PetscFE` object 454 . blockSize - The number of elements in a block 455 . numBlocks - The number of blocks in a batch 456 . batchSize - The number of elements in a batch 457 - numBatches - The number of batches in a chunk 458 459 Level: intermediate 460 461 .seealso: `PetscFE`, `PetscFECreate()`, `PetscFEGetTileSizes()` 462 @*/ 463 PetscErrorCode PetscFESetTileSizes(PetscFE fem, PetscInt blockSize, PetscInt numBlocks, PetscInt batchSize, PetscInt numBatches) 464 { 465 PetscFunctionBegin; 466 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 467 fem->blockSize = blockSize; 468 fem->numBlocks = numBlocks; 469 fem->batchSize = batchSize; 470 fem->numBatches = numBatches; 471 PetscFunctionReturn(PETSC_SUCCESS); 472 } 473 474 /*@ 475 PetscFEGetTileSizes - Returns the tile sizes for evaluation 476 477 Not Collective 478 479 Input Parameter: 480 . fem - The `PetscFE` object 481 482 Output Parameters: 483 + blockSize - The number of elements in a block 484 . numBlocks - The number of blocks in a batch 485 . batchSize - The number of elements in a batch 486 - numBatches - The number of batches in a chunk 487 488 Level: intermediate 489 490 .seealso: `PetscFE`, `PetscFECreate()`, `PetscFESetTileSizes()` 491 @*/ 492 PetscErrorCode PetscFEGetTileSizes(PetscFE fem, PetscInt *blockSize, PetscInt *numBlocks, PetscInt *batchSize, PetscInt *numBatches) 493 { 494 PetscFunctionBegin; 495 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 496 if (blockSize) PetscAssertPointer(blockSize, 2); 497 if (numBlocks) PetscAssertPointer(numBlocks, 3); 498 if (batchSize) PetscAssertPointer(batchSize, 4); 499 if (numBatches) PetscAssertPointer(numBatches, 5); 500 if (blockSize) *blockSize = fem->blockSize; 501 if (numBlocks) *numBlocks = fem->numBlocks; 502 if (batchSize) *batchSize = fem->batchSize; 503 if (numBatches) *numBatches = fem->numBatches; 504 PetscFunctionReturn(PETSC_SUCCESS); 505 } 506 507 /*@ 508 PetscFEGetBasisSpace - Returns the `PetscSpace` used for the approximation of the solution for the `PetscFE` 509 510 Not Collective 511 512 Input Parameter: 513 . fem - The `PetscFE` object 514 515 Output Parameter: 516 . sp - The `PetscSpace` object 517 518 Level: intermediate 519 520 .seealso: `PetscFE`, `PetscSpace`, `PetscFECreate()` 521 @*/ 522 PetscErrorCode PetscFEGetBasisSpace(PetscFE fem, PetscSpace *sp) 523 { 524 PetscFunctionBegin; 525 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 526 PetscAssertPointer(sp, 2); 527 *sp = fem->basisSpace; 528 PetscFunctionReturn(PETSC_SUCCESS); 529 } 530 531 /*@ 532 PetscFESetBasisSpace - Sets the `PetscSpace` used for the approximation of the solution 533 534 Not Collective 535 536 Input Parameters: 537 + fem - The `PetscFE` object 538 - sp - The `PetscSpace` object 539 540 Level: intermediate 541 542 Developer Notes: 543 There is `PetscFESetBasisSpace()` but the `PetscFESetDualSpace()`, likely the Basis is unneeded in the function name 544 545 .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscFECreate()`, `PetscFESetDualSpace()` 546 @*/ 547 PetscErrorCode PetscFESetBasisSpace(PetscFE fem, PetscSpace sp) 548 { 549 PetscFunctionBegin; 550 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 551 PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 2); 552 PetscCall(PetscSpaceDestroy(&fem->basisSpace)); 553 fem->basisSpace = sp; 554 PetscCall(PetscObjectReference((PetscObject)fem->basisSpace)); 555 PetscFunctionReturn(PETSC_SUCCESS); 556 } 557 558 /*@ 559 PetscFEGetDualSpace - Returns the `PetscDualSpace` used to define the inner product for a `PetscFE` 560 561 Not Collective 562 563 Input Parameter: 564 . fem - The `PetscFE` object 565 566 Output Parameter: 567 . sp - The `PetscDualSpace` object 568 569 Level: intermediate 570 571 .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscFECreate()` 572 @*/ 573 PetscErrorCode PetscFEGetDualSpace(PetscFE fem, PetscDualSpace *sp) 574 { 575 PetscFunctionBegin; 576 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 577 PetscAssertPointer(sp, 2); 578 *sp = fem->dualSpace; 579 PetscFunctionReturn(PETSC_SUCCESS); 580 } 581 582 /*@ 583 PetscFESetDualSpace - Sets the `PetscDualSpace` used to define the inner product 584 585 Not Collective 586 587 Input Parameters: 588 + fem - The `PetscFE` object 589 - sp - The `PetscDualSpace` object 590 591 Level: intermediate 592 593 .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscFECreate()`, `PetscFESetBasisSpace()` 594 @*/ 595 PetscErrorCode PetscFESetDualSpace(PetscFE fem, PetscDualSpace sp) 596 { 597 PetscFunctionBegin; 598 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 599 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 2); 600 PetscCall(PetscDualSpaceDestroy(&fem->dualSpace)); 601 fem->dualSpace = sp; 602 PetscCall(PetscObjectReference((PetscObject)fem->dualSpace)); 603 PetscFunctionReturn(PETSC_SUCCESS); 604 } 605 606 /*@ 607 PetscFEGetQuadrature - Returns the `PetscQuadrature` used to calculate inner products 608 609 Not Collective 610 611 Input Parameter: 612 . fem - The `PetscFE` object 613 614 Output Parameter: 615 . q - The `PetscQuadrature` object 616 617 Level: intermediate 618 619 .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscQuadrature`, `PetscFECreate()` 620 @*/ 621 PetscErrorCode PetscFEGetQuadrature(PetscFE fem, PetscQuadrature *q) 622 { 623 PetscFunctionBegin; 624 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 625 PetscAssertPointer(q, 2); 626 *q = fem->quadrature; 627 PetscFunctionReturn(PETSC_SUCCESS); 628 } 629 630 /*@ 631 PetscFESetQuadrature - Sets the `PetscQuadrature` used to calculate inner products 632 633 Not Collective 634 635 Input Parameters: 636 + fem - The `PetscFE` object 637 - q - The `PetscQuadrature` object 638 639 Level: intermediate 640 641 .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscQuadrature`, `PetscFECreate()`, `PetscFEGetFaceQuadrature()` 642 @*/ 643 PetscErrorCode PetscFESetQuadrature(PetscFE fem, PetscQuadrature q) 644 { 645 PetscInt Nc, qNc; 646 647 PetscFunctionBegin; 648 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 649 if (q == fem->quadrature) PetscFunctionReturn(PETSC_SUCCESS); 650 PetscCall(PetscFEGetNumComponents(fem, &Nc)); 651 PetscCall(PetscQuadratureGetNumComponents(q, &qNc)); 652 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); 653 PetscCall(PetscTabulationDestroy(&fem->T)); 654 PetscCall(PetscTabulationDestroy(&fem->Tc)); 655 PetscCall(PetscObjectReference((PetscObject)q)); 656 PetscCall(PetscQuadratureDestroy(&fem->quadrature)); 657 fem->quadrature = q; 658 PetscFunctionReturn(PETSC_SUCCESS); 659 } 660 661 /*@ 662 PetscFEGetFaceQuadrature - Returns the `PetscQuadrature` used to calculate inner products on faces 663 664 Not Collective 665 666 Input Parameter: 667 . fem - The `PetscFE` object 668 669 Output Parameter: 670 . q - The `PetscQuadrature` object 671 672 Level: intermediate 673 674 Developer Notes: 675 There is a special face quadrature but not edge, likely this API would benefit from a refactorization 676 677 .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscQuadrature`, `PetscFECreate()`, `PetscFESetQuadrature()`, `PetscFESetFaceQuadrature()` 678 @*/ 679 PetscErrorCode PetscFEGetFaceQuadrature(PetscFE fem, PetscQuadrature *q) 680 { 681 PetscFunctionBegin; 682 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 683 PetscAssertPointer(q, 2); 684 *q = fem->faceQuadrature; 685 PetscFunctionReturn(PETSC_SUCCESS); 686 } 687 688 /*@ 689 PetscFESetFaceQuadrature - Sets the `PetscQuadrature` used to calculate inner products on faces 690 691 Not Collective 692 693 Input Parameters: 694 + fem - The `PetscFE` object 695 - q - The `PetscQuadrature` object 696 697 Level: intermediate 698 699 .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscQuadrature`, `PetscFECreate()`, `PetscFESetQuadrature()` 700 @*/ 701 PetscErrorCode PetscFESetFaceQuadrature(PetscFE fem, PetscQuadrature q) 702 { 703 PetscInt Nc, qNc; 704 705 PetscFunctionBegin; 706 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 707 if (q == fem->faceQuadrature) PetscFunctionReturn(PETSC_SUCCESS); 708 PetscCall(PetscFEGetNumComponents(fem, &Nc)); 709 PetscCall(PetscQuadratureGetNumComponents(q, &qNc)); 710 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); 711 PetscCall(PetscTabulationDestroy(&fem->Tf)); 712 PetscCall(PetscObjectReference((PetscObject)q)); 713 PetscCall(PetscQuadratureDestroy(&fem->faceQuadrature)); 714 fem->faceQuadrature = q; 715 PetscFunctionReturn(PETSC_SUCCESS); 716 } 717 718 /*@ 719 PetscFECopyQuadrature - Copy both volumetric and surface quadrature to a new `PetscFE` 720 721 Not Collective 722 723 Input Parameters: 724 + sfe - The `PetscFE` source for the quadratures 725 - tfe - The `PetscFE` target for the quadratures 726 727 Level: intermediate 728 729 .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscQuadrature`, `PetscFECreate()`, `PetscFESetQuadrature()`, `PetscFESetFaceQuadrature()` 730 @*/ 731 PetscErrorCode PetscFECopyQuadrature(PetscFE sfe, PetscFE tfe) 732 { 733 PetscQuadrature q; 734 735 PetscFunctionBegin; 736 PetscValidHeaderSpecific(sfe, PETSCFE_CLASSID, 1); 737 PetscValidHeaderSpecific(tfe, PETSCFE_CLASSID, 2); 738 PetscCall(PetscFEGetQuadrature(sfe, &q)); 739 PetscCall(PetscFESetQuadrature(tfe, q)); 740 PetscCall(PetscFEGetFaceQuadrature(sfe, &q)); 741 PetscCall(PetscFESetFaceQuadrature(tfe, q)); 742 PetscFunctionReturn(PETSC_SUCCESS); 743 } 744 745 /*@C 746 PetscFEGetNumDof - Returns the number of dofs (dual basis vectors) associated to mesh points on the reference cell of a given dimension 747 748 Not Collective 749 750 Input Parameter: 751 . fem - The `PetscFE` object 752 753 Output Parameter: 754 . numDof - Array with the number of dofs per dimension 755 756 Level: intermediate 757 758 .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscFECreate()` 759 @*/ 760 PetscErrorCode PetscFEGetNumDof(PetscFE fem, const PetscInt **numDof) 761 { 762 PetscFunctionBegin; 763 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 764 PetscAssertPointer(numDof, 2); 765 PetscCall(PetscDualSpaceGetNumDof(fem->dualSpace, numDof)); 766 PetscFunctionReturn(PETSC_SUCCESS); 767 } 768 769 /*@C 770 PetscFEGetCellTabulation - Returns the tabulation of the basis functions at the quadrature points on the reference cell 771 772 Not Collective 773 774 Input Parameters: 775 + fem - The `PetscFE` object 776 - k - The highest derivative we need to tabulate, very often 1 777 778 Output Parameter: 779 . T - The basis function values and derivatives at quadrature points 780 781 Level: intermediate 782 783 Note: 784 .vb 785 T->T[0] = B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c 786 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 787 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 788 .ve 789 790 .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscTabulation`, `PetscFECreateTabulation()`, `PetscTabulationDestroy()` 791 @*/ 792 PetscErrorCode PetscFEGetCellTabulation(PetscFE fem, PetscInt k, PetscTabulation *T) 793 { 794 PetscInt npoints; 795 const PetscReal *points; 796 797 PetscFunctionBegin; 798 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 799 PetscAssertPointer(T, 3); 800 PetscCall(PetscQuadratureGetData(fem->quadrature, NULL, NULL, &npoints, &points, NULL)); 801 if (!fem->T) PetscCall(PetscFECreateTabulation(fem, 1, npoints, points, k, &fem->T)); 802 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); 803 *T = fem->T; 804 PetscFunctionReturn(PETSC_SUCCESS); 805 } 806 807 /*@C 808 PetscFEGetFaceTabulation - Returns the tabulation of the basis functions at the face quadrature points for each face of the reference cell 809 810 Not Collective 811 812 Input Parameters: 813 + fem - The `PetscFE` object 814 - k - The highest derivative we need to tabulate, very often 1 815 816 Output Parameter: 817 . Tf - The basis function values and derivatives at face quadrature points 818 819 Level: intermediate 820 821 Note: 822 .vb 823 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 824 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 825 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 826 .ve 827 828 .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscTabulation`, `PetscFEGetCellTabulation()`, `PetscFECreateTabulation()`, `PetscTabulationDestroy()` 829 @*/ 830 PetscErrorCode PetscFEGetFaceTabulation(PetscFE fem, PetscInt k, PetscTabulation *Tf) 831 { 832 PetscFunctionBegin; 833 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 834 PetscAssertPointer(Tf, 3); 835 if (!fem->Tf) { 836 const PetscReal xi0[3] = {-1., -1., -1.}; 837 PetscReal v0[3], J[9], detJ; 838 PetscQuadrature fq; 839 PetscDualSpace sp; 840 DM dm; 841 const PetscInt *faces; 842 PetscInt dim, numFaces, f, npoints, q; 843 const PetscReal *points; 844 PetscReal *facePoints; 845 846 PetscCall(PetscFEGetDualSpace(fem, &sp)); 847 PetscCall(PetscDualSpaceGetDM(sp, &dm)); 848 PetscCall(DMGetDimension(dm, &dim)); 849 PetscCall(DMPlexGetConeSize(dm, 0, &numFaces)); 850 PetscCall(DMPlexGetCone(dm, 0, &faces)); 851 PetscCall(PetscFEGetFaceQuadrature(fem, &fq)); 852 if (fq) { 853 PetscCall(PetscQuadratureGetData(fq, NULL, NULL, &npoints, &points, NULL)); 854 PetscCall(PetscMalloc1(numFaces * npoints * dim, &facePoints)); 855 for (f = 0; f < numFaces; ++f) { 856 PetscCall(DMPlexComputeCellGeometryFEM(dm, faces[f], NULL, v0, J, NULL, &detJ)); 857 for (q = 0; q < npoints; ++q) CoordinatesRefToReal(dim, dim - 1, xi0, v0, J, &points[q * (dim - 1)], &facePoints[(f * npoints + q) * dim]); 858 } 859 PetscCall(PetscFECreateTabulation(fem, numFaces, npoints, facePoints, k, &fem->Tf)); 860 PetscCall(PetscFree(facePoints)); 861 } 862 } 863 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); 864 *Tf = fem->Tf; 865 PetscFunctionReturn(PETSC_SUCCESS); 866 } 867 868 /*@C 869 PetscFEGetFaceCentroidTabulation - Returns the tabulation of the basis functions at the face centroid points 870 871 Not Collective 872 873 Input Parameter: 874 . fem - The `PetscFE` object 875 876 Output Parameter: 877 . Tc - The basis function values at face centroid points 878 879 Level: intermediate 880 881 Note: 882 .vb 883 T->T[0] = Bf[(f*pdim + i)*Nc + c] is the value at point f for basis function i and component c 884 .ve 885 886 .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscTabulation`, `PetscFEGetFaceTabulation()`, `PetscFEGetCellTabulation()`, `PetscFECreateTabulation()`, `PetscTabulationDestroy()` 887 @*/ 888 PetscErrorCode PetscFEGetFaceCentroidTabulation(PetscFE fem, PetscTabulation *Tc) 889 { 890 PetscFunctionBegin; 891 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 892 PetscAssertPointer(Tc, 2); 893 if (!fem->Tc) { 894 PetscDualSpace sp; 895 DM dm; 896 const PetscInt *cone; 897 PetscReal *centroids; 898 PetscInt dim, numFaces, f; 899 900 PetscCall(PetscFEGetDualSpace(fem, &sp)); 901 PetscCall(PetscDualSpaceGetDM(sp, &dm)); 902 PetscCall(DMGetDimension(dm, &dim)); 903 PetscCall(DMPlexGetConeSize(dm, 0, &numFaces)); 904 PetscCall(DMPlexGetCone(dm, 0, &cone)); 905 PetscCall(PetscMalloc1(numFaces * dim, ¢roids)); 906 for (f = 0; f < numFaces; ++f) PetscCall(DMPlexComputeCellGeometryFVM(dm, cone[f], NULL, ¢roids[f * dim], NULL)); 907 PetscCall(PetscFECreateTabulation(fem, 1, numFaces, centroids, 0, &fem->Tc)); 908 PetscCall(PetscFree(centroids)); 909 } 910 *Tc = fem->Tc; 911 PetscFunctionReturn(PETSC_SUCCESS); 912 } 913 914 /*@C 915 PetscFECreateTabulation - Tabulates the basis functions, and perhaps derivatives, at the points provided. 916 917 Not Collective 918 919 Input Parameters: 920 + fem - The `PetscFE` object 921 . nrepl - The number of replicas 922 . npoints - The number of tabulation points in a replica 923 . points - The tabulation point coordinates 924 - K - The number of derivatives calculated 925 926 Output Parameter: 927 . T - The basis function values and derivatives at tabulation points 928 929 Level: intermediate 930 931 Note: 932 .vb 933 T->T[0] = B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c 934 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 935 T->T[2] = H[(((p*pdim + i)*Nc + c)*dim + d)*dim + e] is the value at point p for basis 936 T->function i, component c, in directions d and e 937 .ve 938 939 .seealso: `PetscTabulation`, `PetscFEGetCellTabulation()`, `PetscTabulationDestroy()` 940 @*/ 941 PetscErrorCode PetscFECreateTabulation(PetscFE fem, PetscInt nrepl, PetscInt npoints, const PetscReal points[], PetscInt K, PetscTabulation *T) 942 { 943 DM dm; 944 PetscDualSpace Q; 945 PetscInt Nb; /* Dimension of FE space P */ 946 PetscInt Nc; /* Field components */ 947 PetscInt cdim; /* Reference coordinate dimension */ 948 PetscInt k; 949 950 PetscFunctionBegin; 951 if (!npoints || !fem->dualSpace || K < 0) { 952 *T = NULL; 953 PetscFunctionReturn(PETSC_SUCCESS); 954 } 955 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 956 PetscAssertPointer(points, 4); 957 PetscAssertPointer(T, 6); 958 PetscCall(PetscFEGetDualSpace(fem, &Q)); 959 PetscCall(PetscDualSpaceGetDM(Q, &dm)); 960 PetscCall(DMGetDimension(dm, &cdim)); 961 PetscCall(PetscDualSpaceGetDimension(Q, &Nb)); 962 PetscCall(PetscFEGetNumComponents(fem, &Nc)); 963 PetscCall(PetscMalloc1(1, T)); 964 (*T)->K = !cdim ? 0 : K; 965 (*T)->Nr = nrepl; 966 (*T)->Np = npoints; 967 (*T)->Nb = Nb; 968 (*T)->Nc = Nc; 969 (*T)->cdim = cdim; 970 PetscCall(PetscMalloc1((*T)->K + 1, &(*T)->T)); 971 for (k = 0; k <= (*T)->K; ++k) PetscCall(PetscMalloc1(nrepl * npoints * Nb * Nc * PetscPowInt(cdim, k), &(*T)->T[k])); 972 PetscUseTypeMethod(fem, createtabulation, nrepl * npoints, points, K, *T); 973 PetscFunctionReturn(PETSC_SUCCESS); 974 } 975 976 /*@C 977 PetscFEComputeTabulation - Tabulates the basis functions, and perhaps derivatives, at the points provided. 978 979 Not Collective 980 981 Input Parameters: 982 + fem - The `PetscFE` object 983 . npoints - The number of tabulation points 984 . points - The tabulation point coordinates 985 . K - The number of derivatives calculated 986 - T - An existing tabulation object with enough allocated space 987 988 Output Parameter: 989 . T - The basis function values and derivatives at tabulation points 990 991 Level: intermediate 992 993 Note: 994 .vb 995 T->T[0] = B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c 996 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 997 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 998 .ve 999 1000 .seealso: `PetscTabulation`, `PetscFEGetCellTabulation()`, `PetscTabulationDestroy()` 1001 @*/ 1002 PetscErrorCode PetscFEComputeTabulation(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscInt K, PetscTabulation T) 1003 { 1004 PetscFunctionBeginHot; 1005 if (!npoints || !fem->dualSpace || K < 0) PetscFunctionReturn(PETSC_SUCCESS); 1006 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 1007 PetscAssertPointer(points, 3); 1008 PetscAssertPointer(T, 5); 1009 if (PetscDefined(USE_DEBUG)) { 1010 DM dm; 1011 PetscDualSpace Q; 1012 PetscInt Nb; /* Dimension of FE space P */ 1013 PetscInt Nc; /* Field components */ 1014 PetscInt cdim; /* Reference coordinate dimension */ 1015 1016 PetscCall(PetscFEGetDualSpace(fem, &Q)); 1017 PetscCall(PetscDualSpaceGetDM(Q, &dm)); 1018 PetscCall(DMGetDimension(dm, &cdim)); 1019 PetscCall(PetscDualSpaceGetDimension(Q, &Nb)); 1020 PetscCall(PetscFEGetNumComponents(fem, &Nc)); 1021 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); 1022 PetscCheck(T->Nb == Nb, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Tabulation Nb %" PetscInt_FMT " must match requested Nb %" PetscInt_FMT, T->Nb, Nb); 1023 PetscCheck(T->Nc == Nc, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Tabulation Nc %" PetscInt_FMT " must match requested Nc %" PetscInt_FMT, T->Nc, Nc); 1024 PetscCheck(T->cdim == cdim, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Tabulation cdim %" PetscInt_FMT " must match requested cdim %" PetscInt_FMT, T->cdim, cdim); 1025 } 1026 T->Nr = 1; 1027 T->Np = npoints; 1028 PetscUseTypeMethod(fem, createtabulation, npoints, points, K, T); 1029 PetscFunctionReturn(PETSC_SUCCESS); 1030 } 1031 1032 /*@C 1033 PetscTabulationDestroy - Frees memory from the associated tabulation. 1034 1035 Not Collective 1036 1037 Input Parameter: 1038 . T - The tabulation 1039 1040 Level: intermediate 1041 1042 .seealso: `PetscTabulation`, `PetscFECreateTabulation()`, `PetscFEGetCellTabulation()` 1043 @*/ 1044 PetscErrorCode PetscTabulationDestroy(PetscTabulation *T) 1045 { 1046 PetscInt k; 1047 1048 PetscFunctionBegin; 1049 PetscAssertPointer(T, 1); 1050 if (!T || !(*T)) PetscFunctionReturn(PETSC_SUCCESS); 1051 for (k = 0; k <= (*T)->K; ++k) PetscCall(PetscFree((*T)->T[k])); 1052 PetscCall(PetscFree((*T)->T)); 1053 PetscCall(PetscFree(*T)); 1054 *T = NULL; 1055 PetscFunctionReturn(PETSC_SUCCESS); 1056 } 1057 1058 PETSC_EXTERN PetscErrorCode PetscFECreatePointTrace(PetscFE fe, PetscInt refPoint, PetscFE *trFE) 1059 { 1060 PetscSpace bsp, bsubsp; 1061 PetscDualSpace dsp, dsubsp; 1062 PetscInt dim, depth, numComp, i, j, coneSize, order; 1063 PetscFEType type; 1064 DM dm; 1065 DMLabel label; 1066 PetscReal *xi, *v, *J, detJ; 1067 const char *name; 1068 PetscQuadrature origin, fullQuad, subQuad; 1069 1070 PetscFunctionBegin; 1071 PetscValidHeaderSpecific(fe, PETSCFE_CLASSID, 1); 1072 PetscAssertPointer(trFE, 3); 1073 PetscCall(PetscFEGetBasisSpace(fe, &bsp)); 1074 PetscCall(PetscFEGetDualSpace(fe, &dsp)); 1075 PetscCall(PetscDualSpaceGetDM(dsp, &dm)); 1076 PetscCall(DMGetDimension(dm, &dim)); 1077 PetscCall(DMPlexGetDepthLabel(dm, &label)); 1078 PetscCall(DMLabelGetValue(label, refPoint, &depth)); 1079 PetscCall(PetscCalloc1(depth, &xi)); 1080 PetscCall(PetscMalloc1(dim, &v)); 1081 PetscCall(PetscMalloc1(dim * dim, &J)); 1082 for (i = 0; i < depth; i++) xi[i] = 0.; 1083 PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, &origin)); 1084 PetscCall(PetscQuadratureSetData(origin, depth, 0, 1, xi, NULL)); 1085 PetscCall(DMPlexComputeCellGeometryFEM(dm, refPoint, origin, v, J, NULL, &detJ)); 1086 /* CellGeometryFEM computes the expanded Jacobian, we want the true jacobian */ 1087 for (i = 1; i < dim; i++) { 1088 for (j = 0; j < depth; j++) J[i * depth + j] = J[i * dim + j]; 1089 } 1090 PetscCall(PetscQuadratureDestroy(&origin)); 1091 PetscCall(PetscDualSpaceGetPointSubspace(dsp, refPoint, &dsubsp)); 1092 PetscCall(PetscSpaceCreateSubspace(bsp, dsubsp, v, J, NULL, NULL, PETSC_OWN_POINTER, &bsubsp)); 1093 PetscCall(PetscSpaceSetUp(bsubsp)); 1094 PetscCall(PetscFECreate(PetscObjectComm((PetscObject)fe), trFE)); 1095 PetscCall(PetscFEGetType(fe, &type)); 1096 PetscCall(PetscFESetType(*trFE, type)); 1097 PetscCall(PetscFEGetNumComponents(fe, &numComp)); 1098 PetscCall(PetscFESetNumComponents(*trFE, numComp)); 1099 PetscCall(PetscFESetBasisSpace(*trFE, bsubsp)); 1100 PetscCall(PetscFESetDualSpace(*trFE, dsubsp)); 1101 PetscCall(PetscObjectGetName((PetscObject)fe, &name)); 1102 if (name) PetscCall(PetscFESetName(*trFE, name)); 1103 PetscCall(PetscFEGetQuadrature(fe, &fullQuad)); 1104 PetscCall(PetscQuadratureGetOrder(fullQuad, &order)); 1105 PetscCall(DMPlexGetConeSize(dm, refPoint, &coneSize)); 1106 if (coneSize == 2 * depth) PetscCall(PetscDTGaussTensorQuadrature(depth, 1, (order + 2) / 2, -1., 1., &subQuad)); 1107 else PetscCall(PetscDTSimplexQuadrature(depth, order, PETSCDTSIMPLEXQUAD_DEFAULT, &subQuad)); 1108 PetscCall(PetscFESetQuadrature(*trFE, subQuad)); 1109 PetscCall(PetscFESetUp(*trFE)); 1110 PetscCall(PetscQuadratureDestroy(&subQuad)); 1111 PetscCall(PetscSpaceDestroy(&bsubsp)); 1112 PetscFunctionReturn(PETSC_SUCCESS); 1113 } 1114 1115 PetscErrorCode PetscFECreateHeightTrace(PetscFE fe, PetscInt height, PetscFE *trFE) 1116 { 1117 PetscInt hStart, hEnd; 1118 PetscDualSpace dsp; 1119 DM dm; 1120 1121 PetscFunctionBegin; 1122 PetscValidHeaderSpecific(fe, PETSCFE_CLASSID, 1); 1123 PetscAssertPointer(trFE, 3); 1124 *trFE = NULL; 1125 PetscCall(PetscFEGetDualSpace(fe, &dsp)); 1126 PetscCall(PetscDualSpaceGetDM(dsp, &dm)); 1127 PetscCall(DMPlexGetHeightStratum(dm, height, &hStart, &hEnd)); 1128 if (hEnd <= hStart) PetscFunctionReturn(PETSC_SUCCESS); 1129 PetscCall(PetscFECreatePointTrace(fe, hStart, trFE)); 1130 PetscFunctionReturn(PETSC_SUCCESS); 1131 } 1132 1133 /*@ 1134 PetscFEGetDimension - Get the dimension of the finite element space on a cell 1135 1136 Not Collective 1137 1138 Input Parameter: 1139 . fem - The `PetscFE` 1140 1141 Output Parameter: 1142 . dim - The dimension 1143 1144 Level: intermediate 1145 1146 .seealso: `PetscFE`, `PetscFECreate()`, `PetscSpaceGetDimension()`, `PetscDualSpaceGetDimension()` 1147 @*/ 1148 PetscErrorCode PetscFEGetDimension(PetscFE fem, PetscInt *dim) 1149 { 1150 PetscFunctionBegin; 1151 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 1152 PetscAssertPointer(dim, 2); 1153 PetscTryTypeMethod(fem, getdimension, dim); 1154 PetscFunctionReturn(PETSC_SUCCESS); 1155 } 1156 1157 /*@C 1158 PetscFEPushforward - Map the reference element function to real space 1159 1160 Input Parameters: 1161 + fe - The `PetscFE` 1162 . fegeom - The cell geometry 1163 . Nv - The number of function values 1164 - vals - The function values 1165 1166 Output Parameter: 1167 . vals - The transformed function values 1168 1169 Level: advanced 1170 1171 Notes: 1172 This just forwards the call onto `PetscDualSpacePushforward()`. 1173 1174 It only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension. 1175 1176 .seealso: `PetscFE`, `PetscFEGeom`, `PetscDualSpace`, `PetscDualSpacePushforward()` 1177 @*/ 1178 PetscErrorCode PetscFEPushforward(PetscFE fe, PetscFEGeom *fegeom, PetscInt Nv, PetscScalar vals[]) 1179 { 1180 PetscFunctionBeginHot; 1181 PetscCall(PetscDualSpacePushforward(fe->dualSpace, fegeom, Nv, fe->numComponents, vals)); 1182 PetscFunctionReturn(PETSC_SUCCESS); 1183 } 1184 1185 /*@C 1186 PetscFEPushforwardGradient - Map the reference element function gradient to real space 1187 1188 Input Parameters: 1189 + fe - The `PetscFE` 1190 . fegeom - The cell geometry 1191 . Nv - The number of function gradient values 1192 - vals - The function gradient values 1193 1194 Output Parameter: 1195 . vals - The transformed function gradient values 1196 1197 Level: advanced 1198 1199 Notes: 1200 This just forwards the call onto `PetscDualSpacePushforwardGradient()`. 1201 1202 It only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension. 1203 1204 .seealso: `PetscFE`, `PetscFEGeom`, `PetscDualSpace`, `PetscFEPushforward()`, `PetscDualSpacePushforwardGradient()`, `PetscDualSpacePushforward()` 1205 @*/ 1206 PetscErrorCode PetscFEPushforwardGradient(PetscFE fe, PetscFEGeom *fegeom, PetscInt Nv, PetscScalar vals[]) 1207 { 1208 PetscFunctionBeginHot; 1209 PetscCall(PetscDualSpacePushforwardGradient(fe->dualSpace, fegeom, Nv, fe->numComponents, vals)); 1210 PetscFunctionReturn(PETSC_SUCCESS); 1211 } 1212 1213 /*@C 1214 PetscFEPushforwardHessian - Map the reference element function Hessian to real space 1215 1216 Input Parameters: 1217 + fe - The `PetscFE` 1218 . fegeom - The cell geometry 1219 . Nv - The number of function Hessian values 1220 - vals - The function Hessian values 1221 1222 Output Parameter: 1223 . vals - The transformed function Hessian values 1224 1225 Level: advanced 1226 1227 Notes: 1228 This just forwards the call onto `PetscDualSpacePushforwardHessian()`. 1229 1230 It only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension. 1231 1232 Developer Notes: 1233 It is unclear why all these one line convenience routines are desirable 1234 1235 .seealso: `PetscFE`, `PetscFEGeom`, `PetscDualSpace`, `PetscFEPushforward()`, `PetscDualSpacePushforwardHessian()`, `PetscDualSpacePushforward()` 1236 @*/ 1237 PetscErrorCode PetscFEPushforwardHessian(PetscFE fe, PetscFEGeom *fegeom, PetscInt Nv, PetscScalar vals[]) 1238 { 1239 PetscFunctionBeginHot; 1240 PetscCall(PetscDualSpacePushforwardHessian(fe->dualSpace, fegeom, Nv, fe->numComponents, vals)); 1241 PetscFunctionReturn(PETSC_SUCCESS); 1242 } 1243 1244 /* 1245 Purpose: Compute element vector for chunk of elements 1246 1247 Input: 1248 Sizes: 1249 Ne: number of elements 1250 Nf: number of fields 1251 PetscFE 1252 dim: spatial dimension 1253 Nb: number of basis functions 1254 Nc: number of field components 1255 PetscQuadrature 1256 Nq: number of quadrature points 1257 1258 Geometry: 1259 PetscFEGeom[Ne] possibly *Nq 1260 PetscReal v0s[dim] 1261 PetscReal n[dim] 1262 PetscReal jacobians[dim*dim] 1263 PetscReal jacobianInverses[dim*dim] 1264 PetscReal jacobianDeterminants 1265 FEM: 1266 PetscFE 1267 PetscQuadrature 1268 PetscReal quadPoints[Nq*dim] 1269 PetscReal quadWeights[Nq] 1270 PetscReal basis[Nq*Nb*Nc] 1271 PetscReal basisDer[Nq*Nb*Nc*dim] 1272 PetscScalar coefficients[Ne*Nb*Nc] 1273 PetscScalar elemVec[Ne*Nb*Nc] 1274 1275 Problem: 1276 PetscInt f: the active field 1277 f0, f1 1278 1279 Work Space: 1280 PetscFE 1281 PetscScalar f0[Nq*dim]; 1282 PetscScalar f1[Nq*dim*dim]; 1283 PetscScalar u[Nc]; 1284 PetscScalar gradU[Nc*dim]; 1285 PetscReal x[dim]; 1286 PetscScalar realSpaceDer[dim]; 1287 1288 Purpose: Compute element vector for N_cb batches of elements 1289 1290 Input: 1291 Sizes: 1292 N_cb: Number of serial cell batches 1293 1294 Geometry: 1295 PetscReal v0s[Ne*dim] 1296 PetscReal jacobians[Ne*dim*dim] possibly *Nq 1297 PetscReal jacobianInverses[Ne*dim*dim] possibly *Nq 1298 PetscReal jacobianDeterminants[Ne] possibly *Nq 1299 FEM: 1300 static PetscReal quadPoints[Nq*dim] 1301 static PetscReal quadWeights[Nq] 1302 static PetscReal basis[Nq*Nb*Nc] 1303 static PetscReal basisDer[Nq*Nb*Nc*dim] 1304 PetscScalar coefficients[Ne*Nb*Nc] 1305 PetscScalar elemVec[Ne*Nb*Nc] 1306 1307 ex62.c: 1308 PetscErrorCode PetscFEIntegrateResidualBatch(PetscInt Ne, PetscInt numFields, PetscInt field, PetscQuadrature quad[], const PetscScalar coefficients[], 1309 const PetscReal v0s[], const PetscReal jacobians[], const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], 1310 void (*f0_func)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar f0[]), 1311 void (*f1_func)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar f1[]), PetscScalar elemVec[]) 1312 1313 ex52.c: 1314 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) 1315 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) 1316 1317 ex52_integrateElement.cu 1318 __global__ void integrateElementQuadrature(int N_cb, realType *coefficients, realType *jacobianInverses, realType *jacobianDeterminants, realType *elemVec) 1319 1320 PETSC_EXTERN PetscErrorCode IntegrateElementBatchGPU(PetscInt spatial_dim, PetscInt Ne, PetscInt Ncb, PetscInt Nbc, PetscInt Nbl, const PetscScalar coefficients[], 1321 const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], PetscScalar elemVec[], 1322 PetscLogEvent event, PetscInt debug, PetscInt pde_op) 1323 1324 ex52_integrateElementOpenCL.c: 1325 PETSC_EXTERN PetscErrorCode IntegrateElementBatchGPU(PetscInt spatial_dim, PetscInt Ne, PetscInt Ncb, PetscInt Nbc, PetscInt N_bl, const PetscScalar coefficients[], 1326 const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], PetscScalar elemVec[], 1327 PetscLogEvent event, PetscInt debug, PetscInt pde_op) 1328 1329 __kernel void integrateElementQuadrature(int N_cb, __global float *coefficients, __global float *jacobianInverses, __global float *jacobianDeterminants, __global float *elemVec) 1330 */ 1331 1332 /*@C 1333 PetscFEIntegrate - Produce the integral for the given field for a chunk of elements by quadrature integration 1334 1335 Not Collective 1336 1337 Input Parameters: 1338 + prob - The `PetscDS` specifying the discretizations and continuum functions 1339 . field - The field being integrated 1340 . Ne - The number of elements in the chunk 1341 . cgeom - The cell geometry for each cell in the chunk 1342 . coefficients - The array of FEM basis coefficients for the elements 1343 . probAux - The `PetscDS` specifying the auxiliary discretizations 1344 - coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 1345 1346 Output Parameter: 1347 . integral - the integral for this field 1348 1349 Level: intermediate 1350 1351 Developer Notes: 1352 The function name begins with `PetscFE` and yet the first argument is `PetscDS` and it has no `PetscFE` arguments. 1353 1354 .seealso: `PetscFE`, `PetscDS`, `PetscFEIntegrateResidual()`, `PetscFEIntegrateBd()` 1355 @*/ 1356 PetscErrorCode PetscFEIntegrate(PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom, const PetscScalar coefficients[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscScalar integral[]) 1357 { 1358 PetscFE fe; 1359 1360 PetscFunctionBegin; 1361 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1362 PetscCall(PetscDSGetDiscretization(prob, field, (PetscObject *)&fe)); 1363 if (fe->ops->integrate) PetscCall((*fe->ops->integrate)(prob, field, Ne, cgeom, coefficients, probAux, coefficientsAux, integral)); 1364 PetscFunctionReturn(PETSC_SUCCESS); 1365 } 1366 1367 /*@C 1368 PetscFEIntegrateBd - Produce the integral for the given field for a chunk of elements by quadrature integration 1369 1370 Not Collective 1371 1372 Input Parameters: 1373 + prob - The `PetscDS` specifying the discretizations and continuum functions 1374 . field - The field being integrated 1375 . obj_func - The function to be integrated 1376 . Ne - The number of elements in the chunk 1377 . geom - The face geometry for each face in the chunk 1378 . coefficients - The array of FEM basis coefficients for the elements 1379 . probAux - The `PetscDS` specifying the auxiliary discretizations 1380 - coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 1381 1382 Output Parameter: 1383 . integral - the integral for this field 1384 1385 Level: intermediate 1386 1387 Developer Notes: 1388 The function name begins with `PetscFE` and yet the first argument is `PetscDS` and it has no `PetscFE` arguments. 1389 1390 .seealso: `PetscFE`, `PetscDS`, `PetscFEIntegrateResidual()`, `PetscFEIntegrate()` 1391 @*/ 1392 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[]) 1393 { 1394 PetscFE fe; 1395 1396 PetscFunctionBegin; 1397 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1398 PetscCall(PetscDSGetDiscretization(prob, field, (PetscObject *)&fe)); 1399 if (fe->ops->integratebd) PetscCall((*fe->ops->integratebd)(prob, field, obj_func, Ne, geom, coefficients, probAux, coefficientsAux, integral)); 1400 PetscFunctionReturn(PETSC_SUCCESS); 1401 } 1402 1403 /*@C 1404 PetscFEIntegrateResidual - Produce the element residual vector for a chunk of elements by quadrature integration 1405 1406 Not Collective 1407 1408 Input Parameters: 1409 + ds - The `PetscDS` specifying the discretizations and continuum functions 1410 . key - The (label+value, field) being integrated 1411 . Ne - The number of elements in the chunk 1412 . cgeom - The cell geometry for each cell in the chunk 1413 . coefficients - The array of FEM basis coefficients for the elements 1414 . coefficients_t - The array of FEM basis time derivative coefficients for the elements 1415 . probAux - The `PetscDS` specifying the auxiliary discretizations 1416 . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 1417 - t - The time 1418 1419 Output Parameter: 1420 . elemVec - the element residual vectors from each element 1421 1422 Level: intermediate 1423 1424 Note: 1425 .vb 1426 Loop over batch of elements (e): 1427 Loop over quadrature points (q): 1428 Make u_q and gradU_q (loops over fields,Nb,Ncomp) and x_q 1429 Call f_0 and f_1 1430 Loop over element vector entries (f,fc --> i): 1431 elemVec[i] += \psi^{fc}_f(q) f0_{fc}(u, \nabla u) + \nabla\psi^{fc}_f(q) \cdot f1_{fc,df}(u, \nabla u) 1432 .ve 1433 1434 .seealso: `PetscFEIntegrateBdResidual()` 1435 @*/ 1436 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[]) 1437 { 1438 PetscFE fe; 1439 1440 PetscFunctionBeginHot; 1441 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1442 PetscCall(PetscDSGetDiscretization(ds, key.field, (PetscObject *)&fe)); 1443 if (fe->ops->integrateresidual) PetscCall((*fe->ops->integrateresidual)(ds, key, Ne, cgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec)); 1444 PetscFunctionReturn(PETSC_SUCCESS); 1445 } 1446 1447 /*@C 1448 PetscFEIntegrateBdResidual - Produce the element residual vector for a chunk of elements by quadrature integration over a boundary 1449 1450 Not Collective 1451 1452 Input Parameters: 1453 + ds - The `PetscDS` specifying the discretizations and continuum functions 1454 . wf - The PetscWeakForm object holding the pointwise functions 1455 . key - The (label+value, field) being integrated 1456 . Ne - The number of elements in the chunk 1457 . fgeom - The face geometry for each cell in the chunk 1458 . coefficients - The array of FEM basis coefficients for the elements 1459 . coefficients_t - The array of FEM basis time derivative coefficients for the elements 1460 . probAux - The `PetscDS` specifying the auxiliary discretizations 1461 . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 1462 - t - The time 1463 1464 Output Parameter: 1465 . elemVec - the element residual vectors from each element 1466 1467 Level: intermediate 1468 1469 .seealso: `PetscFEIntegrateResidual()` 1470 @*/ 1471 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[]) 1472 { 1473 PetscFE fe; 1474 1475 PetscFunctionBegin; 1476 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1477 PetscCall(PetscDSGetDiscretization(ds, key.field, (PetscObject *)&fe)); 1478 if (fe->ops->integratebdresidual) PetscCall((*fe->ops->integratebdresidual)(ds, wf, key, Ne, fgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec)); 1479 PetscFunctionReturn(PETSC_SUCCESS); 1480 } 1481 1482 /*@C 1483 PetscFEIntegrateHybridResidual - Produce the element residual vector for a chunk of hybrid element faces by quadrature integration 1484 1485 Not Collective 1486 1487 Input Parameters: 1488 + ds - The `PetscDS` specifying the discretizations and continuum functions 1489 . dsIn - The `PetscDS` specifying the discretizations and continuum functions for input 1490 . key - The (label+value, field) being integrated 1491 . s - The side of the cell being integrated, 0 for negative and 1 for positive 1492 . Ne - The number of elements in the chunk 1493 . fgeom - The face geometry for each cell in the chunk 1494 . coefficients - The array of FEM basis coefficients for the elements 1495 . coefficients_t - The array of FEM basis time derivative coefficients for the elements 1496 . probAux - The `PetscDS` specifying the auxiliary discretizations 1497 . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 1498 - t - The time 1499 1500 Output Parameter: 1501 . elemVec - the element residual vectors from each element 1502 1503 Level: developer 1504 1505 .seealso: `PetscFEIntegrateResidual()` 1506 @*/ 1507 PetscErrorCode PetscFEIntegrateHybridResidual(PetscDS ds, PetscDS dsIn, PetscFormKey key, PetscInt s, PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[]) 1508 { 1509 PetscFE fe; 1510 1511 PetscFunctionBegin; 1512 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1513 PetscValidHeaderSpecific(dsIn, PETSCDS_CLASSID, 2); 1514 PetscCall(PetscDSGetDiscretization(ds, key.field, (PetscObject *)&fe)); 1515 if (fe->ops->integratehybridresidual) PetscCall((*fe->ops->integratehybridresidual)(ds, dsIn, key, s, Ne, fgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec)); 1516 PetscFunctionReturn(PETSC_SUCCESS); 1517 } 1518 1519 /*@C 1520 PetscFEIntegrateJacobian - Produce the element Jacobian for a chunk of elements by quadrature integration 1521 1522 Not Collective 1523 1524 Input Parameters: 1525 + ds - The `PetscDS` specifying the discretizations and continuum functions 1526 . jtype - The type of matrix pointwise functions that should be used 1527 . key - The (label+value, fieldI*Nf + fieldJ) being integrated 1528 . Ne - The number of elements in the chunk 1529 . cgeom - The cell geometry for each cell in the chunk 1530 . coefficients - The array of FEM basis coefficients for the elements for the Jacobian evaluation point 1531 . coefficients_t - The array of FEM basis time derivative coefficients for the elements 1532 . probAux - The `PetscDS` specifying the auxiliary discretizations 1533 . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 1534 . t - The time 1535 - u_tshift - A multiplier for the dF/du_t term (as opposed to the dF/du term) 1536 1537 Output Parameter: 1538 . elemMat - the element matrices for the Jacobian from each element 1539 1540 Level: intermediate 1541 1542 Note: 1543 .vb 1544 Loop over batch of elements (e): 1545 Loop over element matrix entries (f,fc,g,gc --> i,j): 1546 Loop over quadrature points (q): 1547 Make u_q and gradU_q (loops over fields,Nb,Ncomp) 1548 elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q) 1549 + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q) 1550 + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q) 1551 + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q) 1552 .ve 1553 1554 .seealso: `PetscFEIntegrateResidual()` 1555 @*/ 1556 PetscErrorCode PetscFEIntegrateJacobian(PetscDS ds, PetscFEJacobianType jtype, PetscFormKey key, PetscInt Ne, PetscFEGeom *cgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) 1557 { 1558 PetscFE fe; 1559 PetscInt Nf; 1560 1561 PetscFunctionBegin; 1562 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1563 PetscCall(PetscDSGetNumFields(ds, &Nf)); 1564 PetscCall(PetscDSGetDiscretization(ds, key.field / Nf, (PetscObject *)&fe)); 1565 if (fe->ops->integratejacobian) PetscCall((*fe->ops->integratejacobian)(ds, jtype, key, Ne, cgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, u_tshift, elemMat)); 1566 PetscFunctionReturn(PETSC_SUCCESS); 1567 } 1568 1569 /*@C 1570 PetscFEIntegrateBdJacobian - Produce the boundary element Jacobian for a chunk of elements by quadrature integration 1571 1572 Not Collective 1573 1574 Input Parameters: 1575 + ds - The `PetscDS` specifying the discretizations and continuum functions 1576 . wf - The PetscWeakForm holding the pointwise functions 1577 . key - The (label+value, fieldI*Nf + fieldJ) being integrated 1578 . Ne - The number of elements in the chunk 1579 . fgeom - The face geometry for each cell in the chunk 1580 . coefficients - The array of FEM basis coefficients for the elements for the Jacobian evaluation point 1581 . coefficients_t - The array of FEM basis time derivative coefficients for the elements 1582 . probAux - The `PetscDS` specifying the auxiliary discretizations 1583 . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 1584 . t - The time 1585 - u_tshift - A multiplier for the dF/du_t term (as opposed to the dF/du term) 1586 1587 Output Parameter: 1588 . elemMat - the element matrices for the Jacobian from each element 1589 1590 Level: intermediate 1591 1592 Note: 1593 .vb 1594 Loop over batch of elements (e): 1595 Loop over element matrix entries (f,fc,g,gc --> i,j): 1596 Loop over quadrature points (q): 1597 Make u_q and gradU_q (loops over fields,Nb,Ncomp) 1598 elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q) 1599 + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q) 1600 + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q) 1601 + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q) 1602 .ve 1603 1604 .seealso: `PetscFEIntegrateJacobian()`, `PetscFEIntegrateResidual()` 1605 @*/ 1606 PetscErrorCode PetscFEIntegrateBdJacobian(PetscDS ds, PetscWeakForm wf, 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[]) 1607 { 1608 PetscFE fe; 1609 PetscInt Nf; 1610 1611 PetscFunctionBegin; 1612 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1613 PetscCall(PetscDSGetNumFields(ds, &Nf)); 1614 PetscCall(PetscDSGetDiscretization(ds, key.field / Nf, (PetscObject *)&fe)); 1615 if (fe->ops->integratebdjacobian) PetscCall((*fe->ops->integratebdjacobian)(ds, wf, key, Ne, fgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, u_tshift, elemMat)); 1616 PetscFunctionReturn(PETSC_SUCCESS); 1617 } 1618 1619 /*@C 1620 PetscFEIntegrateHybridJacobian - Produce the boundary element Jacobian for a chunk of hybrid elements by quadrature integration 1621 1622 Not Collective 1623 1624 Input Parameters: 1625 + ds - The `PetscDS` specifying the discretizations and continuum functions for the output 1626 . dsIn - The `PetscDS` specifying the discretizations and continuum functions for the input 1627 . jtype - The type of matrix pointwise functions that should be used 1628 . key - The (label+value, fieldI*Nf + fieldJ) being integrated 1629 . s - The side of the cell being integrated, 0 for negative and 1 for positive 1630 . Ne - The number of elements in the chunk 1631 . fgeom - The face geometry for each cell in the chunk 1632 . coefficients - The array of FEM basis coefficients for the elements for the Jacobian evaluation point 1633 . coefficients_t - The array of FEM basis time derivative coefficients for the elements 1634 . probAux - The `PetscDS` specifying the auxiliary discretizations 1635 . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 1636 . t - The time 1637 - u_tshift - A multiplier for the dF/du_t term (as opposed to the dF/du term) 1638 1639 Output Parameter: 1640 . elemMat - the element matrices for the Jacobian from each element 1641 1642 Level: developer 1643 1644 Note: 1645 .vb 1646 Loop over batch of elements (e): 1647 Loop over element matrix entries (f,fc,g,gc --> i,j): 1648 Loop over quadrature points (q): 1649 Make u_q and gradU_q (loops over fields,Nb,Ncomp) 1650 elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q) 1651 + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q) 1652 + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q) 1653 + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q) 1654 .ve 1655 1656 .seealso: `PetscFEIntegrateJacobian()`, `PetscFEIntegrateResidual()` 1657 @*/ 1658 PetscErrorCode PetscFEIntegrateHybridJacobian(PetscDS ds, PetscDS dsIn, PetscFEJacobianType jtype, PetscFormKey key, PetscInt s, PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) 1659 { 1660 PetscFE fe; 1661 PetscInt Nf; 1662 1663 PetscFunctionBegin; 1664 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1665 PetscCall(PetscDSGetNumFields(ds, &Nf)); 1666 PetscCall(PetscDSGetDiscretization(ds, key.field / Nf, (PetscObject *)&fe)); 1667 if (fe->ops->integratehybridjacobian) PetscCall((*fe->ops->integratehybridjacobian)(ds, dsIn, jtype, key, s, Ne, fgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, u_tshift, elemMat)); 1668 PetscFunctionReturn(PETSC_SUCCESS); 1669 } 1670 1671 /*@ 1672 PetscFEGetHeightSubspace - Get the subspace of this space for a mesh point of a given height 1673 1674 Input Parameters: 1675 + fe - The finite element space 1676 - height - The height of the `DMPLEX` point 1677 1678 Output Parameter: 1679 . subfe - The subspace of this `PetscFE` space 1680 1681 Level: advanced 1682 1683 Note: 1684 For example, if we want the subspace of this space for a face, we would choose height = 1. 1685 1686 .seealso: `PetscFECreateDefault()` 1687 @*/ 1688 PetscErrorCode PetscFEGetHeightSubspace(PetscFE fe, PetscInt height, PetscFE *subfe) 1689 { 1690 PetscSpace P, subP; 1691 PetscDualSpace Q, subQ; 1692 PetscQuadrature subq; 1693 PetscFEType fetype; 1694 PetscInt dim, Nc; 1695 1696 PetscFunctionBegin; 1697 PetscValidHeaderSpecific(fe, PETSCFE_CLASSID, 1); 1698 PetscAssertPointer(subfe, 3); 1699 if (height == 0) { 1700 *subfe = fe; 1701 PetscFunctionReturn(PETSC_SUCCESS); 1702 } 1703 PetscCall(PetscFEGetBasisSpace(fe, &P)); 1704 PetscCall(PetscFEGetDualSpace(fe, &Q)); 1705 PetscCall(PetscFEGetNumComponents(fe, &Nc)); 1706 PetscCall(PetscFEGetFaceQuadrature(fe, &subq)); 1707 PetscCall(PetscDualSpaceGetDimension(Q, &dim)); 1708 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); 1709 if (!fe->subspaces) PetscCall(PetscCalloc1(dim, &fe->subspaces)); 1710 if (height <= dim) { 1711 if (!fe->subspaces[height - 1]) { 1712 PetscFE sub = NULL; 1713 const char *name; 1714 1715 PetscCall(PetscSpaceGetHeightSubspace(P, height, &subP)); 1716 PetscCall(PetscDualSpaceGetHeightSubspace(Q, height, &subQ)); 1717 if (subQ) { 1718 PetscCall(PetscFECreate(PetscObjectComm((PetscObject)fe), &sub)); 1719 PetscCall(PetscObjectGetName((PetscObject)fe, &name)); 1720 PetscCall(PetscObjectSetName((PetscObject)sub, name)); 1721 PetscCall(PetscFEGetType(fe, &fetype)); 1722 PetscCall(PetscFESetType(sub, fetype)); 1723 PetscCall(PetscFESetBasisSpace(sub, subP)); 1724 PetscCall(PetscFESetDualSpace(sub, subQ)); 1725 PetscCall(PetscFESetNumComponents(sub, Nc)); 1726 PetscCall(PetscFESetUp(sub)); 1727 PetscCall(PetscFESetQuadrature(sub, subq)); 1728 } 1729 fe->subspaces[height - 1] = sub; 1730 } 1731 *subfe = fe->subspaces[height - 1]; 1732 } else { 1733 *subfe = NULL; 1734 } 1735 PetscFunctionReturn(PETSC_SUCCESS); 1736 } 1737 1738 /*@ 1739 PetscFERefine - Create a "refined" `PetscFE` object that refines the reference cell into 1740 smaller copies. 1741 1742 Collective 1743 1744 Input Parameter: 1745 . fe - The initial `PetscFE` 1746 1747 Output Parameter: 1748 . feRef - The refined `PetscFE` 1749 1750 Level: advanced 1751 1752 Notes: 1753 This is typically used to generate a preconditioner for a higher order method from a lower order method on a 1754 refined mesh having the same number of dofs (but more sparsity). It is also used to create an 1755 interpolation between regularly refined meshes. 1756 1757 .seealso: `PetscFEType`, `PetscFECreate()`, `PetscFESetType()` 1758 @*/ 1759 PetscErrorCode PetscFERefine(PetscFE fe, PetscFE *feRef) 1760 { 1761 PetscSpace P, Pref; 1762 PetscDualSpace Q, Qref; 1763 DM K, Kref; 1764 PetscQuadrature q, qref; 1765 const PetscReal *v0, *jac; 1766 PetscInt numComp, numSubelements; 1767 PetscInt cStart, cEnd, c; 1768 PetscDualSpace *cellSpaces; 1769 1770 PetscFunctionBegin; 1771 PetscCall(PetscFEGetBasisSpace(fe, &P)); 1772 PetscCall(PetscFEGetDualSpace(fe, &Q)); 1773 PetscCall(PetscFEGetQuadrature(fe, &q)); 1774 PetscCall(PetscDualSpaceGetDM(Q, &K)); 1775 /* Create space */ 1776 PetscCall(PetscObjectReference((PetscObject)P)); 1777 Pref = P; 1778 /* Create dual space */ 1779 PetscCall(PetscDualSpaceDuplicate(Q, &Qref)); 1780 PetscCall(PetscDualSpaceSetType(Qref, PETSCDUALSPACEREFINED)); 1781 PetscCall(DMRefine(K, PetscObjectComm((PetscObject)fe), &Kref)); 1782 PetscCall(PetscDualSpaceSetDM(Qref, Kref)); 1783 PetscCall(DMPlexGetHeightStratum(Kref, 0, &cStart, &cEnd)); 1784 PetscCall(PetscMalloc1(cEnd - cStart, &cellSpaces)); 1785 /* TODO: fix for non-uniform refinement */ 1786 for (c = 0; c < cEnd - cStart; c++) cellSpaces[c] = Q; 1787 PetscCall(PetscDualSpaceRefinedSetCellSpaces(Qref, cellSpaces)); 1788 PetscCall(PetscFree(cellSpaces)); 1789 PetscCall(DMDestroy(&Kref)); 1790 PetscCall(PetscDualSpaceSetUp(Qref)); 1791 /* Create element */ 1792 PetscCall(PetscFECreate(PetscObjectComm((PetscObject)fe), feRef)); 1793 PetscCall(PetscFESetType(*feRef, PETSCFECOMPOSITE)); 1794 PetscCall(PetscFESetBasisSpace(*feRef, Pref)); 1795 PetscCall(PetscFESetDualSpace(*feRef, Qref)); 1796 PetscCall(PetscFEGetNumComponents(fe, &numComp)); 1797 PetscCall(PetscFESetNumComponents(*feRef, numComp)); 1798 PetscCall(PetscFESetUp(*feRef)); 1799 PetscCall(PetscSpaceDestroy(&Pref)); 1800 PetscCall(PetscDualSpaceDestroy(&Qref)); 1801 /* Create quadrature */ 1802 PetscCall(PetscFECompositeGetMapping(*feRef, &numSubelements, &v0, &jac, NULL)); 1803 PetscCall(PetscQuadratureExpandComposite(q, numSubelements, v0, jac, &qref)); 1804 PetscCall(PetscFESetQuadrature(*feRef, qref)); 1805 PetscCall(PetscQuadratureDestroy(&qref)); 1806 PetscFunctionReturn(PETSC_SUCCESS); 1807 } 1808 1809 static PetscErrorCode PetscFESetDefaultName_Private(PetscFE fe) 1810 { 1811 PetscSpace P; 1812 PetscDualSpace Q; 1813 DM K; 1814 DMPolytopeType ct; 1815 PetscInt degree; 1816 char name[64]; 1817 1818 PetscFunctionBegin; 1819 PetscCall(PetscFEGetBasisSpace(fe, &P)); 1820 PetscCall(PetscSpaceGetDegree(P, °ree, NULL)); 1821 PetscCall(PetscFEGetDualSpace(fe, &Q)); 1822 PetscCall(PetscDualSpaceGetDM(Q, &K)); 1823 PetscCall(DMPlexGetCellType(K, 0, &ct)); 1824 switch (ct) { 1825 case DM_POLYTOPE_SEGMENT: 1826 case DM_POLYTOPE_POINT_PRISM_TENSOR: 1827 case DM_POLYTOPE_QUADRILATERAL: 1828 case DM_POLYTOPE_SEG_PRISM_TENSOR: 1829 case DM_POLYTOPE_HEXAHEDRON: 1830 case DM_POLYTOPE_QUAD_PRISM_TENSOR: 1831 PetscCall(PetscSNPrintf(name, sizeof(name), "Q%" PetscInt_FMT, degree)); 1832 break; 1833 case DM_POLYTOPE_TRIANGLE: 1834 case DM_POLYTOPE_TETRAHEDRON: 1835 PetscCall(PetscSNPrintf(name, sizeof(name), "P%" PetscInt_FMT, degree)); 1836 break; 1837 case DM_POLYTOPE_TRI_PRISM: 1838 case DM_POLYTOPE_TRI_PRISM_TENSOR: 1839 PetscCall(PetscSNPrintf(name, sizeof(name), "P%" PetscInt_FMT "xQ%" PetscInt_FMT, degree, degree)); 1840 break; 1841 default: 1842 PetscCall(PetscSNPrintf(name, sizeof(name), "FE")); 1843 } 1844 PetscCall(PetscFESetName(fe, name)); 1845 PetscFunctionReturn(PETSC_SUCCESS); 1846 } 1847 1848 /*@ 1849 PetscFECreateFromSpaces - Create a `PetscFE` from the basis and dual spaces 1850 1851 Collective 1852 1853 Input Parameters: 1854 + P - The basis space 1855 . Q - The dual space 1856 . q - The cell quadrature 1857 - fq - The face quadrature 1858 1859 Output Parameter: 1860 . fem - The `PetscFE` object 1861 1862 Level: beginner 1863 1864 Note: 1865 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. 1866 1867 .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscQuadrature`, 1868 `PetscFECreateLagrangeByCell()`, `PetscFECreateDefault()`, `PetscFECreateByCell()`, `PetscFECreate()`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()` 1869 @*/ 1870 PetscErrorCode PetscFECreateFromSpaces(PetscSpace P, PetscDualSpace Q, PetscQuadrature q, PetscQuadrature fq, PetscFE *fem) 1871 { 1872 PetscInt Nc; 1873 const char *prefix; 1874 1875 PetscFunctionBegin; 1876 PetscCall(PetscFECreate(PetscObjectComm((PetscObject)P), fem)); 1877 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)P, &prefix)); 1878 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*fem, prefix)); 1879 PetscCall(PetscFESetType(*fem, PETSCFEBASIC)); 1880 PetscCall(PetscFESetBasisSpace(*fem, P)); 1881 PetscCall(PetscFESetDualSpace(*fem, Q)); 1882 PetscCall(PetscSpaceGetNumComponents(P, &Nc)); 1883 PetscCall(PetscFESetNumComponents(*fem, Nc)); 1884 PetscCall(PetscFESetUp(*fem)); 1885 PetscCall(PetscSpaceDestroy(&P)); 1886 PetscCall(PetscDualSpaceDestroy(&Q)); 1887 PetscCall(PetscFESetQuadrature(*fem, q)); 1888 PetscCall(PetscFESetFaceQuadrature(*fem, fq)); 1889 PetscCall(PetscQuadratureDestroy(&q)); 1890 PetscCall(PetscQuadratureDestroy(&fq)); 1891 PetscCall(PetscFESetDefaultName_Private(*fem)); 1892 PetscFunctionReturn(PETSC_SUCCESS); 1893 } 1894 1895 static PetscErrorCode PetscFECreate_Internal(MPI_Comm comm, PetscInt dim, PetscInt Nc, DMPolytopeType ct, const char prefix[], PetscInt degree, PetscInt qorder, PetscBool setFromOptions, PetscFE *fem) 1896 { 1897 DM K; 1898 PetscSpace P; 1899 PetscDualSpace Q; 1900 PetscQuadrature q, fq; 1901 PetscBool tensor; 1902 1903 PetscFunctionBegin; 1904 if (prefix) PetscAssertPointer(prefix, 5); 1905 PetscAssertPointer(fem, 9); 1906 switch (ct) { 1907 case DM_POLYTOPE_SEGMENT: 1908 case DM_POLYTOPE_POINT_PRISM_TENSOR: 1909 case DM_POLYTOPE_QUADRILATERAL: 1910 case DM_POLYTOPE_SEG_PRISM_TENSOR: 1911 case DM_POLYTOPE_HEXAHEDRON: 1912 case DM_POLYTOPE_QUAD_PRISM_TENSOR: 1913 tensor = PETSC_TRUE; 1914 break; 1915 default: 1916 tensor = PETSC_FALSE; 1917 } 1918 /* Create space */ 1919 PetscCall(PetscSpaceCreate(comm, &P)); 1920 PetscCall(PetscSpaceSetType(P, PETSCSPACEPOLYNOMIAL)); 1921 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)P, prefix)); 1922 PetscCall(PetscSpacePolynomialSetTensor(P, tensor)); 1923 PetscCall(PetscSpaceSetNumComponents(P, Nc)); 1924 PetscCall(PetscSpaceSetNumVariables(P, dim)); 1925 if (degree >= 0) { 1926 PetscCall(PetscSpaceSetDegree(P, degree, PETSC_DETERMINE)); 1927 if (ct == DM_POLYTOPE_TRI_PRISM || ct == DM_POLYTOPE_TRI_PRISM_TENSOR) { 1928 PetscSpace Pend, Pside; 1929 1930 PetscCall(PetscSpaceCreate(comm, &Pend)); 1931 PetscCall(PetscSpaceSetType(Pend, PETSCSPACEPOLYNOMIAL)); 1932 PetscCall(PetscSpacePolynomialSetTensor(Pend, PETSC_FALSE)); 1933 PetscCall(PetscSpaceSetNumComponents(Pend, Nc)); 1934 PetscCall(PetscSpaceSetNumVariables(Pend, dim - 1)); 1935 PetscCall(PetscSpaceSetDegree(Pend, degree, PETSC_DETERMINE)); 1936 PetscCall(PetscSpaceCreate(comm, &Pside)); 1937 PetscCall(PetscSpaceSetType(Pside, PETSCSPACEPOLYNOMIAL)); 1938 PetscCall(PetscSpacePolynomialSetTensor(Pside, PETSC_FALSE)); 1939 PetscCall(PetscSpaceSetNumComponents(Pside, 1)); 1940 PetscCall(PetscSpaceSetNumVariables(Pside, 1)); 1941 PetscCall(PetscSpaceSetDegree(Pside, degree, PETSC_DETERMINE)); 1942 PetscCall(PetscSpaceSetType(P, PETSCSPACETENSOR)); 1943 PetscCall(PetscSpaceTensorSetNumSubspaces(P, 2)); 1944 PetscCall(PetscSpaceTensorSetSubspace(P, 0, Pend)); 1945 PetscCall(PetscSpaceTensorSetSubspace(P, 1, Pside)); 1946 PetscCall(PetscSpaceDestroy(&Pend)); 1947 PetscCall(PetscSpaceDestroy(&Pside)); 1948 } 1949 } 1950 if (setFromOptions) PetscCall(PetscSpaceSetFromOptions(P)); 1951 PetscCall(PetscSpaceSetUp(P)); 1952 PetscCall(PetscSpaceGetDegree(P, °ree, NULL)); 1953 PetscCall(PetscSpacePolynomialGetTensor(P, &tensor)); 1954 PetscCall(PetscSpaceGetNumComponents(P, &Nc)); 1955 /* Create dual space */ 1956 PetscCall(PetscDualSpaceCreate(comm, &Q)); 1957 PetscCall(PetscDualSpaceSetType(Q, PETSCDUALSPACELAGRANGE)); 1958 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)Q, prefix)); 1959 PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, ct, &K)); 1960 PetscCall(PetscDualSpaceSetDM(Q, K)); 1961 PetscCall(DMDestroy(&K)); 1962 PetscCall(PetscDualSpaceSetNumComponents(Q, Nc)); 1963 PetscCall(PetscDualSpaceSetOrder(Q, degree)); 1964 /* TODO For some reason, we need a tensor dualspace with wedges */ 1965 PetscCall(PetscDualSpaceLagrangeSetTensor(Q, (tensor || (ct == DM_POLYTOPE_TRI_PRISM)) ? PETSC_TRUE : PETSC_FALSE)); 1966 if (setFromOptions) PetscCall(PetscDualSpaceSetFromOptions(Q)); 1967 PetscCall(PetscDualSpaceSetUp(Q)); 1968 /* Create quadrature */ 1969 qorder = qorder >= 0 ? qorder : degree; 1970 if (setFromOptions) { 1971 PetscObjectOptionsBegin((PetscObject)P); 1972 PetscCall(PetscOptionsBoundedInt("-petscfe_default_quadrature_order", "Quadrature order is one less than quadrature points per edge", "PetscFECreateDefault", qorder, &qorder, NULL, 0)); 1973 PetscOptionsEnd(); 1974 } 1975 PetscCall(PetscDTCreateDefaultQuadrature(ct, qorder, &q, &fq)); 1976 /* Create finite element */ 1977 PetscCall(PetscFECreateFromSpaces(P, Q, q, fq, fem)); 1978 if (setFromOptions) PetscCall(PetscFESetFromOptions(*fem)); 1979 PetscFunctionReturn(PETSC_SUCCESS); 1980 } 1981 1982 /*@C 1983 PetscFECreateDefault - Create a `PetscFE` for basic FEM computation 1984 1985 Collective 1986 1987 Input Parameters: 1988 + comm - The MPI comm 1989 . dim - The spatial dimension 1990 . Nc - The number of components 1991 . isSimplex - Flag for simplex reference cell, otherwise its a tensor product 1992 . prefix - The options prefix, or `NULL` 1993 - qorder - The quadrature order or `PETSC_DETERMINE` to use `PetscSpace` polynomial degree 1994 1995 Output Parameter: 1996 . fem - The `PetscFE` object 1997 1998 Level: beginner 1999 2000 Note: 2001 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. 2002 2003 .seealso: `PetscFECreateLagrange()`, `PetscFECreateByCell()`, `PetscSpaceSetFromOptions()`, `PetscDualSpaceSetFromOptions()`, `PetscFESetFromOptions()`, `PetscFECreate()`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()` 2004 @*/ 2005 PetscErrorCode PetscFECreateDefault(MPI_Comm comm, PetscInt dim, PetscInt Nc, PetscBool isSimplex, const char prefix[], PetscInt qorder, PetscFE *fem) 2006 { 2007 PetscFunctionBegin; 2008 PetscCall(PetscFECreate_Internal(comm, dim, Nc, DMPolytopeTypeSimpleShape(dim, isSimplex), prefix, PETSC_DECIDE, qorder, PETSC_TRUE, fem)); 2009 PetscFunctionReturn(PETSC_SUCCESS); 2010 } 2011 2012 /*@C 2013 PetscFECreateByCell - Create a `PetscFE` for basic FEM computation 2014 2015 Collective 2016 2017 Input Parameters: 2018 + comm - The MPI comm 2019 . dim - The spatial dimension 2020 . Nc - The number of components 2021 . ct - The celltype of the reference cell 2022 . prefix - The options prefix, or `NULL` 2023 - qorder - The quadrature order or `PETSC_DETERMINE` to use `PetscSpace` polynomial degree 2024 2025 Output Parameter: 2026 . fem - The `PetscFE` object 2027 2028 Level: beginner 2029 2030 Note: 2031 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. 2032 2033 .seealso: `PetscFECreateDefault()`, `PetscFECreateLagrange()`, `PetscSpaceSetFromOptions()`, `PetscDualSpaceSetFromOptions()`, `PetscFESetFromOptions()`, `PetscFECreate()`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()` 2034 @*/ 2035 PetscErrorCode PetscFECreateByCell(MPI_Comm comm, PetscInt dim, PetscInt Nc, DMPolytopeType ct, const char prefix[], PetscInt qorder, PetscFE *fem) 2036 { 2037 PetscFunctionBegin; 2038 PetscCall(PetscFECreate_Internal(comm, dim, Nc, ct, prefix, PETSC_DECIDE, qorder, PETSC_TRUE, fem)); 2039 PetscFunctionReturn(PETSC_SUCCESS); 2040 } 2041 2042 /*@ 2043 PetscFECreateLagrange - Create a `PetscFE` for the basic Lagrange space of degree k 2044 2045 Collective 2046 2047 Input Parameters: 2048 + comm - The MPI comm 2049 . dim - The spatial dimension 2050 . Nc - The number of components 2051 . isSimplex - Flag for simplex reference cell, otherwise its a tensor product 2052 . k - The degree k of the space 2053 - qorder - The quadrature order or `PETSC_DETERMINE` to use `PetscSpace` polynomial degree 2054 2055 Output Parameter: 2056 . fem - The `PetscFE` object 2057 2058 Level: beginner 2059 2060 Note: 2061 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. 2062 2063 .seealso: `PetscFECreateLagrangeByCell()`, `PetscFECreateDefault()`, `PetscFECreateByCell()`, `PetscFECreate()`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()` 2064 @*/ 2065 PetscErrorCode PetscFECreateLagrange(MPI_Comm comm, PetscInt dim, PetscInt Nc, PetscBool isSimplex, PetscInt k, PetscInt qorder, PetscFE *fem) 2066 { 2067 PetscFunctionBegin; 2068 PetscCall(PetscFECreate_Internal(comm, dim, Nc, DMPolytopeTypeSimpleShape(dim, isSimplex), NULL, k, qorder, PETSC_FALSE, fem)); 2069 PetscFunctionReturn(PETSC_SUCCESS); 2070 } 2071 2072 /*@ 2073 PetscFECreateLagrangeByCell - Create a `PetscFE` for the basic Lagrange space of degree k 2074 2075 Collective 2076 2077 Input Parameters: 2078 + comm - The MPI comm 2079 . dim - The spatial dimension 2080 . Nc - The number of components 2081 . ct - The celltype of the reference cell 2082 . k - The degree k of the space 2083 - qorder - The quadrature order or `PETSC_DETERMINE` to use `PetscSpace` polynomial degree 2084 2085 Output Parameter: 2086 . fem - The `PetscFE` object 2087 2088 Level: beginner 2089 2090 Note: 2091 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. 2092 2093 .seealso: `PetscFECreateLagrange()`, `PetscFECreateDefault()`, `PetscFECreateByCell()`, `PetscFECreate()`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()` 2094 @*/ 2095 PetscErrorCode PetscFECreateLagrangeByCell(MPI_Comm comm, PetscInt dim, PetscInt Nc, DMPolytopeType ct, PetscInt k, PetscInt qorder, PetscFE *fem) 2096 { 2097 PetscFunctionBegin; 2098 PetscCall(PetscFECreate_Internal(comm, dim, Nc, ct, NULL, k, qorder, PETSC_FALSE, fem)); 2099 PetscFunctionReturn(PETSC_SUCCESS); 2100 } 2101 2102 /*@C 2103 PetscFESetName - Names the `PetscFE` and its subobjects 2104 2105 Not Collective 2106 2107 Input Parameters: 2108 + fe - The `PetscFE` 2109 - name - The name 2110 2111 Level: intermediate 2112 2113 .seealso: `PetscFECreate()`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()` 2114 @*/ 2115 PetscErrorCode PetscFESetName(PetscFE fe, const char name[]) 2116 { 2117 PetscSpace P; 2118 PetscDualSpace Q; 2119 2120 PetscFunctionBegin; 2121 PetscCall(PetscFEGetBasisSpace(fe, &P)); 2122 PetscCall(PetscFEGetDualSpace(fe, &Q)); 2123 PetscCall(PetscObjectSetName((PetscObject)fe, name)); 2124 PetscCall(PetscObjectSetName((PetscObject)P, name)); 2125 PetscCall(PetscObjectSetName((PetscObject)Q, name)); 2126 PetscFunctionReturn(PETSC_SUCCESS); 2127 } 2128 2129 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[]) 2130 { 2131 PetscInt dOffset = 0, fOffset = 0, f, g; 2132 2133 for (f = 0; f < Nf; ++f) { 2134 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); 2135 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); 2136 PetscFE fe; 2137 const PetscInt k = ds->jetDegree[f]; 2138 const PetscInt cdim = T[f]->cdim; 2139 const PetscInt dE = fegeom->dimEmbed; 2140 const PetscInt Nq = T[f]->Np; 2141 const PetscInt Nbf = T[f]->Nb; 2142 const PetscInt Ncf = T[f]->Nc; 2143 const PetscReal *Bq = &T[f]->T[0][(r * Nq + q) * Nbf * Ncf]; 2144 const PetscReal *Dq = &T[f]->T[1][(r * Nq + q) * Nbf * Ncf * cdim]; 2145 const PetscReal *Hq = k > 1 ? &T[f]->T[2][(r * Nq + q) * Nbf * Ncf * cdim * cdim] : NULL; 2146 PetscInt hOffset = 0, b, c, d; 2147 2148 PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fe)); 2149 for (c = 0; c < Ncf; ++c) u[fOffset + c] = 0.0; 2150 for (d = 0; d < dE * Ncf; ++d) u_x[fOffset * dE + d] = 0.0; 2151 for (b = 0; b < Nbf; ++b) { 2152 for (c = 0; c < Ncf; ++c) { 2153 const PetscInt cidx = b * Ncf + c; 2154 2155 u[fOffset + c] += Bq[cidx] * coefficients[dOffset + b]; 2156 for (d = 0; d < cdim; ++d) u_x[(fOffset + c) * dE + d] += Dq[cidx * cdim + d] * coefficients[dOffset + b]; 2157 } 2158 } 2159 if (k > 1) { 2160 for (g = 0; g < Nf; ++g) hOffset += T[g]->Nc * dE; 2161 for (d = 0; d < dE * dE * Ncf; ++d) u_x[hOffset + fOffset * dE * dE + d] = 0.0; 2162 for (b = 0; b < Nbf; ++b) { 2163 for (c = 0; c < Ncf; ++c) { 2164 const PetscInt cidx = b * Ncf + c; 2165 2166 for (d = 0; d < cdim * cdim; ++d) u_x[hOffset + (fOffset + c) * dE * dE + d] += Hq[cidx * cdim * cdim + d] * coefficients[dOffset + b]; 2167 } 2168 } 2169 PetscCall(PetscFEPushforwardHessian(fe, fegeom, 1, &u_x[hOffset + fOffset * dE * dE])); 2170 } 2171 PetscCall(PetscFEPushforward(fe, fegeom, 1, &u[fOffset])); 2172 PetscCall(PetscFEPushforwardGradient(fe, fegeom, 1, &u_x[fOffset * dE])); 2173 if (u_t) { 2174 for (c = 0; c < Ncf; ++c) u_t[fOffset + c] = 0.0; 2175 for (b = 0; b < Nbf; ++b) { 2176 for (c = 0; c < Ncf; ++c) { 2177 const PetscInt cidx = b * Ncf + c; 2178 2179 u_t[fOffset + c] += Bq[cidx] * coefficients_t[dOffset + b]; 2180 } 2181 } 2182 PetscCall(PetscFEPushforward(fe, fegeom, 1, &u_t[fOffset])); 2183 } 2184 fOffset += Ncf; 2185 dOffset += Nbf; 2186 } 2187 return PETSC_SUCCESS; 2188 } 2189 2190 PetscErrorCode PetscFEEvaluateFieldJets_Hybrid_Internal(PetscDS ds, PetscInt Nf, PetscInt rc, PetscInt qc, PetscTabulation Tab[], const PetscInt rf[], const PetscInt qf[], PetscTabulation Tabf[], PetscFEGeom *fegeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscScalar u[], PetscScalar u_x[], PetscScalar u_t[]) 2191 { 2192 PetscInt dOffset = 0, fOffset = 0, f, g; 2193 2194 /* f is the field number in the DS, g is the field number in u[] */ 2195 for (f = 0, g = 0; f < Nf; ++f) { 2196 PetscBool isCohesive; 2197 PetscInt Ns, s; 2198 2199 if (!Tab[f]) continue; 2200 PetscCall(PetscDSGetCohesive(ds, f, &isCohesive)); 2201 Ns = isCohesive ? 1 : 2; 2202 { 2203 PetscTabulation T = isCohesive ? Tab[f] : Tabf[f]; 2204 PetscFE fe = (PetscFE)ds->disc[f]; 2205 const PetscInt dEt = T->cdim; 2206 const PetscInt dE = fegeom->dimEmbed; 2207 const PetscInt Nq = T->Np; 2208 const PetscInt Nbf = T->Nb; 2209 const PetscInt Ncf = T->Nc; 2210 2211 for (s = 0; s < Ns; ++s, ++g) { 2212 const PetscInt r = isCohesive ? rc : rf[s]; 2213 const PetscInt q = isCohesive ? qc : qf[s]; 2214 const PetscReal *Bq = &T->T[0][(r * Nq + q) * Nbf * Ncf]; 2215 const PetscReal *Dq = &T->T[1][(r * Nq + q) * Nbf * Ncf * dEt]; 2216 PetscInt b, c, d; 2217 2218 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); 2219 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); 2220 for (c = 0; c < Ncf; ++c) u[fOffset + c] = 0.0; 2221 for (d = 0; d < dE * Ncf; ++d) u_x[fOffset * dE + d] = 0.0; 2222 for (b = 0; b < Nbf; ++b) { 2223 for (c = 0; c < Ncf; ++c) { 2224 const PetscInt cidx = b * Ncf + c; 2225 2226 u[fOffset + c] += Bq[cidx] * coefficients[dOffset + b]; 2227 for (d = 0; d < dEt; ++d) u_x[(fOffset + c) * dE + d] += Dq[cidx * dEt + d] * coefficients[dOffset + b]; 2228 } 2229 } 2230 PetscCall(PetscFEPushforward(fe, fegeom, 1, &u[fOffset])); 2231 PetscCall(PetscFEPushforwardGradient(fe, fegeom, 1, &u_x[fOffset * dE])); 2232 if (u_t) { 2233 for (c = 0; c < Ncf; ++c) u_t[fOffset + c] = 0.0; 2234 for (b = 0; b < Nbf; ++b) { 2235 for (c = 0; c < Ncf; ++c) { 2236 const PetscInt cidx = b * Ncf + c; 2237 2238 u_t[fOffset + c] += Bq[cidx] * coefficients_t[dOffset + b]; 2239 } 2240 } 2241 PetscCall(PetscFEPushforward(fe, fegeom, 1, &u_t[fOffset])); 2242 } 2243 fOffset += Ncf; 2244 dOffset += Nbf; 2245 } 2246 } 2247 } 2248 return PETSC_SUCCESS; 2249 } 2250 2251 PetscErrorCode PetscFEEvaluateFaceFields_Internal(PetscDS prob, PetscInt field, PetscInt faceLoc, const PetscScalar coefficients[], PetscScalar u[]) 2252 { 2253 PetscFE fe; 2254 PetscTabulation Tc; 2255 PetscInt b, c; 2256 2257 if (!prob) return PETSC_SUCCESS; 2258 PetscCall(PetscDSGetDiscretization(prob, field, (PetscObject *)&fe)); 2259 PetscCall(PetscFEGetFaceCentroidTabulation(fe, &Tc)); 2260 { 2261 const PetscReal *faceBasis = Tc->T[0]; 2262 const PetscInt Nb = Tc->Nb; 2263 const PetscInt Nc = Tc->Nc; 2264 2265 for (c = 0; c < Nc; ++c) u[c] = 0.0; 2266 for (b = 0; b < Nb; ++b) { 2267 for (c = 0; c < Nc; ++c) u[c] += coefficients[b] * faceBasis[(faceLoc * Nb + b) * Nc + c]; 2268 } 2269 } 2270 return PETSC_SUCCESS; 2271 } 2272 2273 PetscErrorCode PetscFEUpdateElementVec_Internal(PetscFE fe, PetscTabulation T, PetscInt r, PetscScalar tmpBasis[], PetscScalar tmpBasisDer[], PetscInt e, PetscFEGeom *fegeom, PetscScalar f0[], PetscScalar f1[], PetscScalar elemVec[]) 2274 { 2275 PetscFEGeom pgeom; 2276 const PetscInt dEt = T->cdim; 2277 const PetscInt dE = fegeom->dimEmbed; 2278 const PetscInt Nq = T->Np; 2279 const PetscInt Nb = T->Nb; 2280 const PetscInt Nc = T->Nc; 2281 const PetscReal *basis = &T->T[0][r * Nq * Nb * Nc]; 2282 const PetscReal *basisDer = &T->T[1][r * Nq * Nb * Nc * dEt]; 2283 PetscInt q, b, c, d; 2284 2285 for (q = 0; q < Nq; ++q) { 2286 for (b = 0; b < Nb; ++b) { 2287 for (c = 0; c < Nc; ++c) { 2288 const PetscInt bcidx = b * Nc + c; 2289 2290 tmpBasis[bcidx] = basis[q * Nb * Nc + bcidx]; 2291 for (d = 0; d < dEt; ++d) tmpBasisDer[bcidx * dE + d] = basisDer[q * Nb * Nc * dEt + bcidx * dEt + d]; 2292 for (d = dEt; d < dE; ++d) tmpBasisDer[bcidx * dE + d] = 0.0; 2293 } 2294 } 2295 PetscCall(PetscFEGeomGetCellPoint(fegeom, e, q, &pgeom)); 2296 PetscCall(PetscFEPushforward(fe, &pgeom, Nb, tmpBasis)); 2297 PetscCall(PetscFEPushforwardGradient(fe, &pgeom, Nb, tmpBasisDer)); 2298 for (b = 0; b < Nb; ++b) { 2299 for (c = 0; c < Nc; ++c) { 2300 const PetscInt bcidx = b * Nc + c; 2301 const PetscInt qcidx = q * Nc + c; 2302 2303 elemVec[b] += tmpBasis[bcidx] * f0[qcidx]; 2304 for (d = 0; d < dE; ++d) elemVec[b] += tmpBasisDer[bcidx * dE + d] * f1[qcidx * dE + d]; 2305 } 2306 } 2307 } 2308 return PETSC_SUCCESS; 2309 } 2310 2311 PetscErrorCode PetscFEUpdateElementVec_Hybrid_Internal(PetscFE fe, PetscTabulation T, PetscInt r, PetscInt s, PetscScalar tmpBasis[], PetscScalar tmpBasisDer[], PetscFEGeom *fegeom, PetscScalar f0[], PetscScalar f1[], PetscScalar elemVec[]) 2312 { 2313 const PetscInt dE = T->cdim; 2314 const PetscInt Nq = T->Np; 2315 const PetscInt Nb = T->Nb; 2316 const PetscInt Nc = T->Nc; 2317 const PetscReal *basis = &T->T[0][r * Nq * Nb * Nc]; 2318 const PetscReal *basisDer = &T->T[1][r * Nq * Nb * Nc * dE]; 2319 PetscInt q, b, c, d; 2320 2321 for (q = 0; q < Nq; ++q) { 2322 for (b = 0; b < Nb; ++b) { 2323 for (c = 0; c < Nc; ++c) { 2324 const PetscInt bcidx = b * Nc + c; 2325 2326 tmpBasis[bcidx] = basis[q * Nb * Nc + bcidx]; 2327 for (d = 0; d < dE; ++d) tmpBasisDer[bcidx * dE + d] = basisDer[q * Nb * Nc * dE + bcidx * dE + d]; 2328 } 2329 } 2330 PetscCall(PetscFEPushforward(fe, fegeom, Nb, tmpBasis)); 2331 // TODO This is currently broken since we do not pull the geometry down to the lower dimension 2332 // PetscCall(PetscFEPushforwardGradient(fe, fegeom, Nb, tmpBasisDer)); 2333 for (b = 0; b < Nb; ++b) { 2334 for (c = 0; c < Nc; ++c) { 2335 const PetscInt bcidx = b * Nc + c; 2336 const PetscInt qcidx = q * Nc + c; 2337 2338 elemVec[Nb * s + b] += tmpBasis[bcidx] * f0[qcidx]; 2339 for (d = 0; d < dE; ++d) elemVec[Nb * s + b] += tmpBasisDer[bcidx * dE + d] * f1[qcidx * dE + d]; 2340 } 2341 } 2342 } 2343 return PETSC_SUCCESS; 2344 } 2345 2346 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 eOffset, PetscInt totDim, PetscInt offsetI, PetscInt offsetJ, PetscScalar elemMat[]) 2347 { 2348 const PetscInt cdim = TI->cdim; 2349 const PetscInt dE = fegeom->dimEmbed; 2350 const PetscInt NqI = TI->Np; 2351 const PetscInt NbI = TI->Nb; 2352 const PetscInt NcI = TI->Nc; 2353 const PetscReal *basisI = &TI->T[0][(r * NqI + q) * NbI * NcI]; 2354 const PetscReal *basisDerI = &TI->T[1][(r * NqI + q) * NbI * NcI * cdim]; 2355 const PetscInt NqJ = TJ->Np; 2356 const PetscInt NbJ = TJ->Nb; 2357 const PetscInt NcJ = TJ->Nc; 2358 const PetscReal *basisJ = &TJ->T[0][(r * NqJ + q) * NbJ * NcJ]; 2359 const PetscReal *basisDerJ = &TJ->T[1][(r * NqJ + q) * NbJ * NcJ * cdim]; 2360 PetscInt f, fc, g, gc, df, dg; 2361 2362 for (f = 0; f < NbI; ++f) { 2363 for (fc = 0; fc < NcI; ++fc) { 2364 const PetscInt fidx = f * NcI + fc; /* Test function basis index */ 2365 2366 tmpBasisI[fidx] = basisI[fidx]; 2367 for (df = 0; df < cdim; ++df) tmpBasisDerI[fidx * dE + df] = basisDerI[fidx * cdim + df]; 2368 } 2369 } 2370 PetscCall(PetscFEPushforward(feI, fegeom, NbI, tmpBasisI)); 2371 PetscCall(PetscFEPushforwardGradient(feI, fegeom, NbI, tmpBasisDerI)); 2372 for (g = 0; g < NbJ; ++g) { 2373 for (gc = 0; gc < NcJ; ++gc) { 2374 const PetscInt gidx = g * NcJ + gc; /* Trial function basis index */ 2375 2376 tmpBasisJ[gidx] = basisJ[gidx]; 2377 for (dg = 0; dg < cdim; ++dg) tmpBasisDerJ[gidx * dE + dg] = basisDerJ[gidx * cdim + dg]; 2378 } 2379 } 2380 PetscCall(PetscFEPushforward(feJ, fegeom, NbJ, tmpBasisJ)); 2381 PetscCall(PetscFEPushforwardGradient(feJ, fegeom, NbJ, tmpBasisDerJ)); 2382 for (f = 0; f < NbI; ++f) { 2383 for (fc = 0; fc < NcI; ++fc) { 2384 const PetscInt fidx = f * NcI + fc; /* Test function basis index */ 2385 const PetscInt i = offsetI + f; /* Element matrix row */ 2386 for (g = 0; g < NbJ; ++g) { 2387 for (gc = 0; gc < NcJ; ++gc) { 2388 const PetscInt gidx = g * NcJ + gc; /* Trial function basis index */ 2389 const PetscInt j = offsetJ + g; /* Element matrix column */ 2390 const PetscInt fOff = eOffset + i * totDim + j; 2391 2392 elemMat[fOff] += tmpBasisI[fidx] * g0[fc * NcJ + gc] * tmpBasisJ[gidx]; 2393 for (df = 0; df < dE; ++df) { 2394 elemMat[fOff] += tmpBasisI[fidx] * g1[(fc * NcJ + gc) * dE + df] * tmpBasisDerJ[gidx * dE + df]; 2395 elemMat[fOff] += tmpBasisDerI[fidx * dE + df] * g2[(fc * NcJ + gc) * dE + df] * tmpBasisJ[gidx]; 2396 for (dg = 0; dg < dE; ++dg) elemMat[fOff] += tmpBasisDerI[fidx * dE + df] * g3[((fc * NcJ + gc) * dE + df) * dE + dg] * tmpBasisDerJ[gidx * dE + dg]; 2397 } 2398 } 2399 } 2400 } 2401 } 2402 return PETSC_SUCCESS; 2403 } 2404 2405 PetscErrorCode PetscFEUpdateElementMat_Hybrid_Internal(PetscFE feI, PetscBool isHybridI, PetscFE feJ, PetscBool isHybridJ, PetscInt r, PetscInt s, 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[]) 2406 { 2407 const PetscInt dE = TI->cdim; 2408 const PetscInt NqI = TI->Np; 2409 const PetscInt NbI = TI->Nb; 2410 const PetscInt NcI = TI->Nc; 2411 const PetscReal *basisI = &TI->T[0][(r * NqI + q) * NbI * NcI]; 2412 const PetscReal *basisDerI = &TI->T[1][(r * NqI + q) * NbI * NcI * dE]; 2413 const PetscInt NqJ = TJ->Np; 2414 const PetscInt NbJ = TJ->Nb; 2415 const PetscInt NcJ = TJ->Nc; 2416 const PetscReal *basisJ = &TJ->T[0][(r * NqJ + q) * NbJ * NcJ]; 2417 const PetscReal *basisDerJ = &TJ->T[1][(r * NqJ + q) * NbJ * NcJ * dE]; 2418 const PetscInt so = isHybridI ? 0 : s; 2419 const PetscInt to = isHybridJ ? 0 : s; 2420 PetscInt f, fc, g, gc, df, dg; 2421 2422 for (f = 0; f < NbI; ++f) { 2423 for (fc = 0; fc < NcI; ++fc) { 2424 const PetscInt fidx = f * NcI + fc; /* Test function basis index */ 2425 2426 tmpBasisI[fidx] = basisI[fidx]; 2427 for (df = 0; df < dE; ++df) tmpBasisDerI[fidx * dE + df] = basisDerI[fidx * dE + df]; 2428 } 2429 } 2430 PetscCall(PetscFEPushforward(feI, fegeom, NbI, tmpBasisI)); 2431 PetscCall(PetscFEPushforwardGradient(feI, fegeom, NbI, tmpBasisDerI)); 2432 for (g = 0; g < NbJ; ++g) { 2433 for (gc = 0; gc < NcJ; ++gc) { 2434 const PetscInt gidx = g * NcJ + gc; /* Trial function basis index */ 2435 2436 tmpBasisJ[gidx] = basisJ[gidx]; 2437 for (dg = 0; dg < dE; ++dg) tmpBasisDerJ[gidx * dE + dg] = basisDerJ[gidx * dE + dg]; 2438 } 2439 } 2440 PetscCall(PetscFEPushforward(feJ, fegeom, NbJ, tmpBasisJ)); 2441 // TODO This is currently broken since we do not pull the geometry down to the lower dimension 2442 // PetscCall(PetscFEPushforwardGradient(feJ, fegeom, NbJ, tmpBasisDerJ)); 2443 for (f = 0; f < NbI; ++f) { 2444 for (fc = 0; fc < NcI; ++fc) { 2445 const PetscInt fidx = f * NcI + fc; /* Test function basis index */ 2446 const PetscInt i = offsetI + NbI * so + f; /* Element matrix row */ 2447 for (g = 0; g < NbJ; ++g) { 2448 for (gc = 0; gc < NcJ; ++gc) { 2449 const PetscInt gidx = g * NcJ + gc; /* Trial function basis index */ 2450 const PetscInt j = offsetJ + NbJ * to + g; /* Element matrix column */ 2451 const PetscInt fOff = eOffset + i * totDim + j; 2452 2453 elemMat[fOff] += tmpBasisI[fidx] * g0[fc * NcJ + gc] * tmpBasisJ[gidx]; 2454 for (df = 0; df < dE; ++df) { 2455 elemMat[fOff] += tmpBasisI[fidx] * g1[(fc * NcJ + gc) * dE + df] * tmpBasisDerJ[gidx * dE + df]; 2456 elemMat[fOff] += tmpBasisDerI[fidx * dE + df] * g2[(fc * NcJ + gc) * dE + df] * tmpBasisJ[gidx]; 2457 for (dg = 0; dg < dE; ++dg) elemMat[fOff] += tmpBasisDerI[fidx * dE + df] * g3[((fc * NcJ + gc) * dE + df) * dE + dg] * tmpBasisDerJ[gidx * dE + dg]; 2458 } 2459 } 2460 } 2461 } 2462 } 2463 return PETSC_SUCCESS; 2464 } 2465 2466 PetscErrorCode PetscFECreateCellGeometry(PetscFE fe, PetscQuadrature quad, PetscFEGeom *cgeom) 2467 { 2468 PetscDualSpace dsp; 2469 DM dm; 2470 PetscQuadrature quadDef; 2471 PetscInt dim, cdim, Nq; 2472 2473 PetscFunctionBegin; 2474 PetscCall(PetscFEGetDualSpace(fe, &dsp)); 2475 PetscCall(PetscDualSpaceGetDM(dsp, &dm)); 2476 PetscCall(DMGetDimension(dm, &dim)); 2477 PetscCall(DMGetCoordinateDim(dm, &cdim)); 2478 PetscCall(PetscFEGetQuadrature(fe, &quadDef)); 2479 quad = quad ? quad : quadDef; 2480 PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL)); 2481 PetscCall(PetscMalloc1(Nq * cdim, &cgeom->v)); 2482 PetscCall(PetscMalloc1(Nq * cdim * cdim, &cgeom->J)); 2483 PetscCall(PetscMalloc1(Nq * cdim * cdim, &cgeom->invJ)); 2484 PetscCall(PetscMalloc1(Nq, &cgeom->detJ)); 2485 cgeom->dim = dim; 2486 cgeom->dimEmbed = cdim; 2487 cgeom->numCells = 1; 2488 cgeom->numPoints = Nq; 2489 PetscCall(DMPlexComputeCellGeometryFEM(dm, 0, quad, cgeom->v, cgeom->J, cgeom->invJ, cgeom->detJ)); 2490 PetscFunctionReturn(PETSC_SUCCESS); 2491 } 2492 2493 PetscErrorCode PetscFEDestroyCellGeometry(PetscFE fe, PetscFEGeom *cgeom) 2494 { 2495 PetscFunctionBegin; 2496 PetscCall(PetscFree(cgeom->v)); 2497 PetscCall(PetscFree(cgeom->J)); 2498 PetscCall(PetscFree(cgeom->invJ)); 2499 PetscCall(PetscFree(cgeom->detJ)); 2500 PetscFunctionReturn(PETSC_SUCCESS); 2501 } 2502