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