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