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