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