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