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