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 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 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 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 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 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 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 if (q == fem->faceQuadrature) PetscFunctionReturn(PETSC_SUCCESS); 708 PetscCall(PetscFEGetNumComponents(fem, &Nc)); 709 PetscCall(PetscQuadratureGetNumComponents(q, &qNc)); 710 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); 711 PetscCall(PetscTabulationDestroy(&fem->Tf)); 712 PetscCall(PetscObjectReference((PetscObject)q)); 713 PetscCall(PetscQuadratureDestroy(&fem->faceQuadrature)); 714 fem->faceQuadrature = q; 715 PetscFunctionReturn(PETSC_SUCCESS); 716 } 717 718 /*@ 719 PetscFECopyQuadrature - Copy both volumetric and surface quadrature to a new `PetscFE` 720 721 Not Collective 722 723 Input Parameters: 724 + sfe - The `PetscFE` source for the quadratures 725 - tfe - The `PetscFE` target for the quadratures 726 727 Level: intermediate 728 729 .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscQuadrature`, `PetscFECreate()`, `PetscFESetQuadrature()`, `PetscFESetFaceQuadrature()` 730 @*/ 731 PetscErrorCode PetscFECopyQuadrature(PetscFE sfe, PetscFE tfe) 732 { 733 PetscQuadrature q; 734 735 PetscFunctionBegin; 736 PetscValidHeaderSpecific(sfe, PETSCFE_CLASSID, 1); 737 PetscValidHeaderSpecific(tfe, PETSCFE_CLASSID, 2); 738 PetscCall(PetscFEGetQuadrature(sfe, &q)); 739 PetscCall(PetscFESetQuadrature(tfe, q)); 740 PetscCall(PetscFEGetFaceQuadrature(sfe, &q)); 741 PetscCall(PetscFESetFaceQuadrature(tfe, q)); 742 PetscFunctionReturn(PETSC_SUCCESS); 743 } 744 745 /*@C 746 PetscFEGetNumDof - Returns the number of dofs (dual basis vectors) associated to mesh points on the reference cell of a given dimension 747 748 Not Collective 749 750 Input Parameter: 751 . fem - The `PetscFE` object 752 753 Output Parameter: 754 . numDof - Array with the number of dofs per dimension 755 756 Level: intermediate 757 758 .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscFECreate()` 759 @*/ 760 PetscErrorCode PetscFEGetNumDof(PetscFE fem, const PetscInt **numDof) 761 { 762 PetscFunctionBegin; 763 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 764 PetscValidPointer(numDof, 2); 765 PetscCall(PetscDualSpaceGetNumDof(fem->dualSpace, numDof)); 766 PetscFunctionReturn(PETSC_SUCCESS); 767 } 768 769 /*@C 770 PetscFEGetCellTabulation - Returns the tabulation of the basis functions at the quadrature points on the reference cell 771 772 Not Collective 773 774 Input Parameters: 775 + fem - The `PetscFE` object 776 - k - The highest derivative we need to tabulate, very often 1 777 778 Output Parameter: 779 . T - The basis function values and derivatives at quadrature points 780 781 Level: intermediate 782 783 Note: 784 .vb 785 T->T[0] = B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c 786 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 787 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 788 .ve 789 790 .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscTabulation`, `PetscFECreateTabulation()`, `PetscTabulationDestroy()` 791 @*/ 792 PetscErrorCode PetscFEGetCellTabulation(PetscFE fem, PetscInt k, PetscTabulation *T) 793 { 794 PetscInt npoints; 795 const PetscReal *points; 796 797 PetscFunctionBegin; 798 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 799 PetscValidPointer(T, 3); 800 PetscCall(PetscQuadratureGetData(fem->quadrature, NULL, NULL, &npoints, &points, NULL)); 801 if (!fem->T) PetscCall(PetscFECreateTabulation(fem, 1, npoints, points, k, &fem->T)); 802 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); 803 *T = fem->T; 804 PetscFunctionReturn(PETSC_SUCCESS); 805 } 806 807 /*@C 808 PetscFEGetFaceTabulation - Returns the tabulation of the basis functions at the face quadrature points for each face of the reference cell 809 810 Not Collective 811 812 Input Parameters: 813 + fem - The `PetscFE` object 814 - k - The highest derivative we need to tabulate, very often 1 815 816 Output Parameters: 817 . Tf - The basis function values and derivatives at face quadrature points 818 819 Level: intermediate 820 821 Note: 822 .vb 823 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 824 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 825 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 826 .ve 827 828 .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscTabulation`, `PetscFEGetCellTabulation()`, `PetscFECreateTabulation()`, `PetscTabulationDestroy()` 829 @*/ 830 PetscErrorCode PetscFEGetFaceTabulation(PetscFE fem, PetscInt k, PetscTabulation *Tf) 831 { 832 PetscFunctionBegin; 833 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 834 PetscValidPointer(Tf, 3); 835 if (!fem->Tf) { 836 const PetscReal xi0[3] = {-1., -1., -1.}; 837 PetscReal v0[3], J[9], detJ; 838 PetscQuadrature fq; 839 PetscDualSpace sp; 840 DM dm; 841 const PetscInt *faces; 842 PetscInt dim, numFaces, f, npoints, q; 843 const PetscReal *points; 844 PetscReal *facePoints; 845 846 PetscCall(PetscFEGetDualSpace(fem, &sp)); 847 PetscCall(PetscDualSpaceGetDM(sp, &dm)); 848 PetscCall(DMGetDimension(dm, &dim)); 849 PetscCall(DMPlexGetConeSize(dm, 0, &numFaces)); 850 PetscCall(DMPlexGetCone(dm, 0, &faces)); 851 PetscCall(PetscFEGetFaceQuadrature(fem, &fq)); 852 if (fq) { 853 PetscCall(PetscQuadratureGetData(fq, NULL, NULL, &npoints, &points, NULL)); 854 PetscCall(PetscMalloc1(numFaces * npoints * dim, &facePoints)); 855 for (f = 0; f < numFaces; ++f) { 856 PetscCall(DMPlexComputeCellGeometryFEM(dm, faces[f], NULL, v0, J, NULL, &detJ)); 857 for (q = 0; q < npoints; ++q) CoordinatesRefToReal(dim, dim - 1, xi0, v0, J, &points[q * (dim - 1)], &facePoints[(f * npoints + q) * dim]); 858 } 859 PetscCall(PetscFECreateTabulation(fem, numFaces, npoints, facePoints, k, &fem->Tf)); 860 PetscCall(PetscFree(facePoints)); 861 } 862 } 863 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); 864 *Tf = fem->Tf; 865 PetscFunctionReturn(PETSC_SUCCESS); 866 } 867 868 /*@C 869 PetscFEGetFaceCentroidTabulation - Returns the tabulation of the basis functions at the face centroid points 870 871 Not Collective 872 873 Input Parameter: 874 . fem - The `PetscFE` object 875 876 Output Parameters: 877 . Tc - The basis function values at face centroid points 878 879 Level: intermediate 880 881 Note: 882 .vb 883 T->T[0] = Bf[(f*pdim + i)*Nc + c] is the value at point f for basis function i and component c 884 .ve 885 886 .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscTabulation`, `PetscFEGetFaceTabulation()`, `PetscFEGetCellTabulation()`, `PetscFECreateTabulation()`, `PetscTabulationDestroy()` 887 @*/ 888 PetscErrorCode PetscFEGetFaceCentroidTabulation(PetscFE fem, PetscTabulation *Tc) 889 { 890 PetscFunctionBegin; 891 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 892 PetscValidPointer(Tc, 2); 893 if (!fem->Tc) { 894 PetscDualSpace sp; 895 DM dm; 896 const PetscInt *cone; 897 PetscReal *centroids; 898 PetscInt dim, numFaces, f; 899 900 PetscCall(PetscFEGetDualSpace(fem, &sp)); 901 PetscCall(PetscDualSpaceGetDM(sp, &dm)); 902 PetscCall(DMGetDimension(dm, &dim)); 903 PetscCall(DMPlexGetConeSize(dm, 0, &numFaces)); 904 PetscCall(DMPlexGetCone(dm, 0, &cone)); 905 PetscCall(PetscMalloc1(numFaces * dim, ¢roids)); 906 for (f = 0; f < numFaces; ++f) PetscCall(DMPlexComputeCellGeometryFVM(dm, cone[f], NULL, ¢roids[f * dim], NULL)); 907 PetscCall(PetscFECreateTabulation(fem, 1, numFaces, centroids, 0, &fem->Tc)); 908 PetscCall(PetscFree(centroids)); 909 } 910 *Tc = fem->Tc; 911 PetscFunctionReturn(PETSC_SUCCESS); 912 } 913 914 /*@C 915 PetscFECreateTabulation - Tabulates the basis functions, and perhaps derivatives, at the points provided. 916 917 Not Collective 918 919 Input Parameters: 920 + fem - The `PetscFE` object 921 . nrepl - The number of replicas 922 . npoints - The number of tabulation points in a replica 923 . points - The tabulation point coordinates 924 - K - The number of derivatives calculated 925 926 Output Parameter: 927 . T - The basis function values and derivatives at tabulation points 928 929 Level: intermediate 930 931 Note: 932 .vb 933 T->T[0] = B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c 934 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 935 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 936 937 .seealso: `PetscTabulation`, `PetscFEGetCellTabulation()`, `PetscTabulationDestroy()` 938 @*/ 939 PetscErrorCode PetscFECreateTabulation(PetscFE fem, PetscInt nrepl, PetscInt npoints, const PetscReal points[], PetscInt K, PetscTabulation *T) 940 { 941 DM dm; 942 PetscDualSpace Q; 943 PetscInt Nb; /* Dimension of FE space P */ 944 PetscInt Nc; /* Field components */ 945 PetscInt cdim; /* Reference coordinate dimension */ 946 PetscInt k; 947 948 PetscFunctionBegin; 949 if (!npoints || !fem->dualSpace || K < 0) { 950 *T = NULL; 951 PetscFunctionReturn(PETSC_SUCCESS); 952 } 953 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 954 PetscValidRealPointer(points, 4); 955 PetscValidPointer(T, 6); 956 PetscCall(PetscFEGetDualSpace(fem, &Q)); 957 PetscCall(PetscDualSpaceGetDM(Q, &dm)); 958 PetscCall(DMGetDimension(dm, &cdim)); 959 PetscCall(PetscDualSpaceGetDimension(Q, &Nb)); 960 PetscCall(PetscFEGetNumComponents(fem, &Nc)); 961 PetscCall(PetscMalloc1(1, T)); 962 (*T)->K = !cdim ? 0 : K; 963 (*T)->Nr = nrepl; 964 (*T)->Np = npoints; 965 (*T)->Nb = Nb; 966 (*T)->Nc = Nc; 967 (*T)->cdim = cdim; 968 PetscCall(PetscMalloc1((*T)->K + 1, &(*T)->T)); 969 for (k = 0; k <= (*T)->K; ++k) PetscCall(PetscMalloc1(nrepl * npoints * Nb * Nc * PetscPowInt(cdim, k), &(*T)->T[k])); 970 PetscUseTypeMethod(fem, createtabulation, nrepl * npoints, points, K, *T); 971 PetscFunctionReturn(PETSC_SUCCESS); 972 } 973 974 /*@C 975 PetscFEComputeTabulation - Tabulates the basis functions, and perhaps derivatives, at the points provided. 976 977 Not Collective 978 979 Input Parameters: 980 + fem - The `PetscFE` object 981 . npoints - The number of tabulation points 982 . points - The tabulation point coordinates 983 . K - The number of derivatives calculated 984 - T - An existing tabulation object with enough allocated space 985 986 Output Parameter: 987 . T - The basis function values and derivatives at tabulation points 988 989 Level: intermediate 990 991 Note: 992 .vb 993 T->T[0] = B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c 994 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 995 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 996 .ve 997 998 .seealso: `PetscTabulation`, `PetscFEGetCellTabulation()`, `PetscTabulationDestroy()` 999 @*/ 1000 PetscErrorCode PetscFEComputeTabulation(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscInt K, PetscTabulation T) 1001 { 1002 PetscFunctionBeginHot; 1003 if (!npoints || !fem->dualSpace || K < 0) PetscFunctionReturn(PETSC_SUCCESS); 1004 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 1005 PetscValidRealPointer(points, 3); 1006 PetscValidPointer(T, 5); 1007 if (PetscDefined(USE_DEBUG)) { 1008 DM dm; 1009 PetscDualSpace Q; 1010 PetscInt Nb; /* Dimension of FE space P */ 1011 PetscInt Nc; /* Field components */ 1012 PetscInt cdim; /* Reference coordinate dimension */ 1013 1014 PetscCall(PetscFEGetDualSpace(fem, &Q)); 1015 PetscCall(PetscDualSpaceGetDM(Q, &dm)); 1016 PetscCall(DMGetDimension(dm, &cdim)); 1017 PetscCall(PetscDualSpaceGetDimension(Q, &Nb)); 1018 PetscCall(PetscFEGetNumComponents(fem, &Nc)); 1019 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); 1020 PetscCheck(T->Nb == Nb, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Tabulation Nb %" PetscInt_FMT " must match requested Nb %" PetscInt_FMT, T->Nb, Nb); 1021 PetscCheck(T->Nc == Nc, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Tabulation Nc %" PetscInt_FMT " must match requested Nc %" PetscInt_FMT, T->Nc, Nc); 1022 PetscCheck(T->cdim == cdim, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Tabulation cdim %" PetscInt_FMT " must match requested cdim %" PetscInt_FMT, T->cdim, cdim); 1023 } 1024 T->Nr = 1; 1025 T->Np = npoints; 1026 PetscUseTypeMethod(fem, createtabulation, npoints, points, K, T); 1027 PetscFunctionReturn(PETSC_SUCCESS); 1028 } 1029 1030 /*@C 1031 PetscTabulationDestroy - Frees memory from the associated tabulation. 1032 1033 Not Collective 1034 1035 Input Parameter: 1036 . T - The tabulation 1037 1038 Level: intermediate 1039 1040 .seealso: `PetscTabulation`, `PetscFECreateTabulation()`, `PetscFEGetCellTabulation()` 1041 @*/ 1042 PetscErrorCode PetscTabulationDestroy(PetscTabulation *T) 1043 { 1044 PetscInt k; 1045 1046 PetscFunctionBegin; 1047 PetscValidPointer(T, 1); 1048 if (!T || !(*T)) PetscFunctionReturn(PETSC_SUCCESS); 1049 for (k = 0; k <= (*T)->K; ++k) PetscCall(PetscFree((*T)->T[k])); 1050 PetscCall(PetscFree((*T)->T)); 1051 PetscCall(PetscFree(*T)); 1052 *T = NULL; 1053 PetscFunctionReturn(PETSC_SUCCESS); 1054 } 1055 1056 PETSC_EXTERN PetscErrorCode PetscFECreatePointTrace(PetscFE fe, PetscInt refPoint, PetscFE *trFE) 1057 { 1058 PetscSpace bsp, bsubsp; 1059 PetscDualSpace dsp, dsubsp; 1060 PetscInt dim, depth, numComp, i, j, coneSize, order; 1061 PetscFEType type; 1062 DM dm; 1063 DMLabel label; 1064 PetscReal *xi, *v, *J, detJ; 1065 const char *name; 1066 PetscQuadrature origin, fullQuad, subQuad; 1067 1068 PetscFunctionBegin; 1069 PetscValidHeaderSpecific(fe, PETSCFE_CLASSID, 1); 1070 PetscValidPointer(trFE, 3); 1071 PetscCall(PetscFEGetBasisSpace(fe, &bsp)); 1072 PetscCall(PetscFEGetDualSpace(fe, &dsp)); 1073 PetscCall(PetscDualSpaceGetDM(dsp, &dm)); 1074 PetscCall(DMGetDimension(dm, &dim)); 1075 PetscCall(DMPlexGetDepthLabel(dm, &label)); 1076 PetscCall(DMLabelGetValue(label, refPoint, &depth)); 1077 PetscCall(PetscCalloc1(depth, &xi)); 1078 PetscCall(PetscMalloc1(dim, &v)); 1079 PetscCall(PetscMalloc1(dim * dim, &J)); 1080 for (i = 0; i < depth; i++) xi[i] = 0.; 1081 PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, &origin)); 1082 PetscCall(PetscQuadratureSetData(origin, depth, 0, 1, xi, NULL)); 1083 PetscCall(DMPlexComputeCellGeometryFEM(dm, refPoint, origin, v, J, NULL, &detJ)); 1084 /* CellGeometryFEM computes the expanded Jacobian, we want the true jacobian */ 1085 for (i = 1; i < dim; i++) { 1086 for (j = 0; j < depth; j++) J[i * depth + j] = J[i * dim + j]; 1087 } 1088 PetscCall(PetscQuadratureDestroy(&origin)); 1089 PetscCall(PetscDualSpaceGetPointSubspace(dsp, refPoint, &dsubsp)); 1090 PetscCall(PetscSpaceCreateSubspace(bsp, dsubsp, v, J, NULL, NULL, PETSC_OWN_POINTER, &bsubsp)); 1091 PetscCall(PetscSpaceSetUp(bsubsp)); 1092 PetscCall(PetscFECreate(PetscObjectComm((PetscObject)fe), trFE)); 1093 PetscCall(PetscFEGetType(fe, &type)); 1094 PetscCall(PetscFESetType(*trFE, type)); 1095 PetscCall(PetscFEGetNumComponents(fe, &numComp)); 1096 PetscCall(PetscFESetNumComponents(*trFE, numComp)); 1097 PetscCall(PetscFESetBasisSpace(*trFE, bsubsp)); 1098 PetscCall(PetscFESetDualSpace(*trFE, dsubsp)); 1099 PetscCall(PetscObjectGetName((PetscObject)fe, &name)); 1100 if (name) PetscCall(PetscFESetName(*trFE, name)); 1101 PetscCall(PetscFEGetQuadrature(fe, &fullQuad)); 1102 PetscCall(PetscQuadratureGetOrder(fullQuad, &order)); 1103 PetscCall(DMPlexGetConeSize(dm, refPoint, &coneSize)); 1104 if (coneSize == 2 * depth) PetscCall(PetscDTGaussTensorQuadrature(depth, 1, (order + 2) / 2, -1., 1., &subQuad)); 1105 else PetscCall(PetscDTSimplexQuadrature(depth, order, PETSCDTSIMPLEXQUAD_DEFAULT, &subQuad)); 1106 PetscCall(PetscFESetQuadrature(*trFE, subQuad)); 1107 PetscCall(PetscFESetUp(*trFE)); 1108 PetscCall(PetscQuadratureDestroy(&subQuad)); 1109 PetscCall(PetscSpaceDestroy(&bsubsp)); 1110 PetscFunctionReturn(PETSC_SUCCESS); 1111 } 1112 1113 PetscErrorCode PetscFECreateHeightTrace(PetscFE fe, PetscInt height, PetscFE *trFE) 1114 { 1115 PetscInt hStart, hEnd; 1116 PetscDualSpace dsp; 1117 DM dm; 1118 1119 PetscFunctionBegin; 1120 PetscValidHeaderSpecific(fe, PETSCFE_CLASSID, 1); 1121 PetscValidPointer(trFE, 3); 1122 *trFE = NULL; 1123 PetscCall(PetscFEGetDualSpace(fe, &dsp)); 1124 PetscCall(PetscDualSpaceGetDM(dsp, &dm)); 1125 PetscCall(DMPlexGetHeightStratum(dm, height, &hStart, &hEnd)); 1126 if (hEnd <= hStart) PetscFunctionReturn(PETSC_SUCCESS); 1127 PetscCall(PetscFECreatePointTrace(fe, hStart, trFE)); 1128 PetscFunctionReturn(PETSC_SUCCESS); 1129 } 1130 1131 /*@ 1132 PetscFEGetDimension - Get the dimension of the finite element space on a cell 1133 1134 Not Collective 1135 1136 Input Parameter: 1137 . fe - The `PetscFE` 1138 1139 Output Parameter: 1140 . dim - The dimension 1141 1142 Level: intermediate 1143 1144 .seealso: `PetscFE`, `PetscFECreate()`, `PetscSpaceGetDimension()`, `PetscDualSpaceGetDimension()` 1145 @*/ 1146 PetscErrorCode PetscFEGetDimension(PetscFE fem, PetscInt *dim) 1147 { 1148 PetscFunctionBegin; 1149 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 1150 PetscValidIntPointer(dim, 2); 1151 PetscTryTypeMethod(fem, getdimension, dim); 1152 PetscFunctionReturn(PETSC_SUCCESS); 1153 } 1154 1155 /*@C 1156 PetscFEPushforward - Map the reference element function to real space 1157 1158 Input Parameters: 1159 + fe - The `PetscFE` 1160 . fegeom - The cell geometry 1161 . Nv - The number of function values 1162 - vals - The function values 1163 1164 Output Parameter: 1165 . vals - The transformed function values 1166 1167 Level: advanced 1168 1169 Notes: 1170 This just forwards the call onto `PetscDualSpacePushforward()`. 1171 1172 It only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension. 1173 1174 .seealso: `PetscFE`, `PetscFEGeom`, `PetscDualSpace`, `PetscDualSpacePushforward()` 1175 @*/ 1176 PetscErrorCode PetscFEPushforward(PetscFE fe, PetscFEGeom *fegeom, PetscInt Nv, PetscScalar vals[]) 1177 { 1178 PetscFunctionBeginHot; 1179 PetscCall(PetscDualSpacePushforward(fe->dualSpace, fegeom, Nv, fe->numComponents, vals)); 1180 PetscFunctionReturn(PETSC_SUCCESS); 1181 } 1182 1183 /*@C 1184 PetscFEPushforwardGradient - Map the reference element function gradient to real space 1185 1186 Input Parameters: 1187 + fe - The `PetscFE` 1188 . fegeom - The cell geometry 1189 . Nv - The number of function gradient values 1190 - vals - The function gradient values 1191 1192 Output Parameter: 1193 . vals - The transformed function gradient values 1194 1195 Level: advanced 1196 1197 Notes: 1198 This just forwards the call onto `PetscDualSpacePushforwardGradient()`. 1199 1200 It only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension. 1201 1202 .seealso: `PetscFE`, `PetscFEGeom`, `PetscDualSpace`, `PetscFEPushforward()`, `PetscDualSpacePushforwardGradient()`, `PetscDualSpacePushforward()` 1203 @*/ 1204 PetscErrorCode PetscFEPushforwardGradient(PetscFE fe, PetscFEGeom *fegeom, PetscInt Nv, PetscScalar vals[]) 1205 { 1206 PetscFunctionBeginHot; 1207 PetscCall(PetscDualSpacePushforwardGradient(fe->dualSpace, fegeom, Nv, fe->numComponents, vals)); 1208 PetscFunctionReturn(PETSC_SUCCESS); 1209 } 1210 1211 /*@C 1212 PetscFEPushforwardHessian - Map the reference element function Hessian to real space 1213 1214 Input Parameters: 1215 + fe - The `PetscFE` 1216 . fegeom - The cell geometry 1217 . Nv - The number of function Hessian values 1218 - vals - The function Hessian values 1219 1220 Output Parameter: 1221 . vals - The transformed function Hessian values 1222 1223 Level: advanced 1224 1225 Notes: 1226 This just forwards the call onto `PetscDualSpacePushforwardHessian()`. 1227 1228 It only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension. 1229 1230 Developer Note: 1231 It is unclear why all these one line convenience routines are desirable 1232 1233 .seealso: `PetscFE`, `PetscFEGeom`, `PetscDualSpace`, `PetscFEPushforward()`, `PetscDualSpacePushforwardHessian()`, `PetscDualSpacePushforward()` 1234 @*/ 1235 PetscErrorCode PetscFEPushforwardHessian(PetscFE fe, PetscFEGeom *fegeom, PetscInt Nv, PetscScalar vals[]) 1236 { 1237 PetscFunctionBeginHot; 1238 PetscCall(PetscDualSpacePushforwardHessian(fe->dualSpace, fegeom, Nv, fe->numComponents, vals)); 1239 PetscFunctionReturn(PETSC_SUCCESS); 1240 } 1241 1242 /* 1243 Purpose: Compute element vector for chunk of elements 1244 1245 Input: 1246 Sizes: 1247 Ne: number of elements 1248 Nf: number of fields 1249 PetscFE 1250 dim: spatial dimension 1251 Nb: number of basis functions 1252 Nc: number of field components 1253 PetscQuadrature 1254 Nq: number of quadrature points 1255 1256 Geometry: 1257 PetscFEGeom[Ne] possibly *Nq 1258 PetscReal v0s[dim] 1259 PetscReal n[dim] 1260 PetscReal jacobians[dim*dim] 1261 PetscReal jacobianInverses[dim*dim] 1262 PetscReal jacobianDeterminants 1263 FEM: 1264 PetscFE 1265 PetscQuadrature 1266 PetscReal quadPoints[Nq*dim] 1267 PetscReal quadWeights[Nq] 1268 PetscReal basis[Nq*Nb*Nc] 1269 PetscReal basisDer[Nq*Nb*Nc*dim] 1270 PetscScalar coefficients[Ne*Nb*Nc] 1271 PetscScalar elemVec[Ne*Nb*Nc] 1272 1273 Problem: 1274 PetscInt f: the active field 1275 f0, f1 1276 1277 Work Space: 1278 PetscFE 1279 PetscScalar f0[Nq*dim]; 1280 PetscScalar f1[Nq*dim*dim]; 1281 PetscScalar u[Nc]; 1282 PetscScalar gradU[Nc*dim]; 1283 PetscReal x[dim]; 1284 PetscScalar realSpaceDer[dim]; 1285 1286 Purpose: Compute element vector for N_cb batches of elements 1287 1288 Input: 1289 Sizes: 1290 N_cb: Number of serial cell batches 1291 1292 Geometry: 1293 PetscReal v0s[Ne*dim] 1294 PetscReal jacobians[Ne*dim*dim] possibly *Nq 1295 PetscReal jacobianInverses[Ne*dim*dim] possibly *Nq 1296 PetscReal jacobianDeterminants[Ne] possibly *Nq 1297 FEM: 1298 static PetscReal quadPoints[Nq*dim] 1299 static PetscReal quadWeights[Nq] 1300 static PetscReal basis[Nq*Nb*Nc] 1301 static PetscReal basisDer[Nq*Nb*Nc*dim] 1302 PetscScalar coefficients[Ne*Nb*Nc] 1303 PetscScalar elemVec[Ne*Nb*Nc] 1304 1305 ex62.c: 1306 PetscErrorCode PetscFEIntegrateResidualBatch(PetscInt Ne, PetscInt numFields, PetscInt field, PetscQuadrature quad[], const PetscScalar coefficients[], 1307 const PetscReal v0s[], const PetscReal jacobians[], const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], 1308 void (*f0_func)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar f0[]), 1309 void (*f1_func)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar f1[]), PetscScalar elemVec[]) 1310 1311 ex52.c: 1312 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) 1313 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) 1314 1315 ex52_integrateElement.cu 1316 __global__ void integrateElementQuadrature(int N_cb, realType *coefficients, realType *jacobianInverses, realType *jacobianDeterminants, realType *elemVec) 1317 1318 PETSC_EXTERN PetscErrorCode IntegrateElementBatchGPU(PetscInt spatial_dim, PetscInt Ne, PetscInt Ncb, PetscInt Nbc, PetscInt Nbl, const PetscScalar coefficients[], 1319 const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], PetscScalar elemVec[], 1320 PetscLogEvent event, PetscInt debug, PetscInt pde_op) 1321 1322 ex52_integrateElementOpenCL.c: 1323 PETSC_EXTERN PetscErrorCode IntegrateElementBatchGPU(PetscInt spatial_dim, PetscInt Ne, PetscInt Ncb, PetscInt Nbc, PetscInt N_bl, const PetscScalar coefficients[], 1324 const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], PetscScalar elemVec[], 1325 PetscLogEvent event, PetscInt debug, PetscInt pde_op) 1326 1327 __kernel void integrateElementQuadrature(int N_cb, __global float *coefficients, __global float *jacobianInverses, __global float *jacobianDeterminants, __global float *elemVec) 1328 */ 1329 1330 /*@C 1331 PetscFEIntegrate - Produce the integral for the given field for a chunk of elements by quadrature integration 1332 1333 Not Collective 1334 1335 Input Parameters: 1336 + prob - The `PetscDS` specifying the discretizations and continuum functions 1337 . field - The field being integrated 1338 . Ne - The number of elements in the chunk 1339 . cgeom - The cell geometry for each cell in the chunk 1340 . coefficients - The array of FEM basis coefficients for the elements 1341 . probAux - The `PetscDS` specifying the auxiliary discretizations 1342 - coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 1343 1344 Output Parameter: 1345 . integral - the integral for this field 1346 1347 Level: intermediate 1348 1349 Developer Note: 1350 The function name begins with `PetscFE` and yet the first argument is `PetscDS` and it has no `PetscFE` arguments. 1351 1352 .seealso: `PetscFE`, `PetscDS`, `PetscFEIntegrateResidual()`, `PetscFEIntegrateBd()` 1353 @*/ 1354 PetscErrorCode PetscFEIntegrate(PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom, const PetscScalar coefficients[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscScalar integral[]) 1355 { 1356 PetscFE fe; 1357 1358 PetscFunctionBegin; 1359 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1360 PetscCall(PetscDSGetDiscretization(prob, field, (PetscObject *)&fe)); 1361 if (fe->ops->integrate) PetscCall((*fe->ops->integrate)(prob, field, Ne, cgeom, coefficients, probAux, coefficientsAux, integral)); 1362 PetscFunctionReturn(PETSC_SUCCESS); 1363 } 1364 1365 /*@C 1366 PetscFEIntegrateBd - Produce the integral for the given field for a chunk of elements by quadrature integration 1367 1368 Not Collective 1369 1370 Input Parameters: 1371 + prob - The `PetscDS` specifying the discretizations and continuum functions 1372 . field - The field being integrated 1373 . obj_func - The function to be integrated 1374 . Ne - The number of elements in the chunk 1375 . fgeom - The face geometry for each face in the chunk 1376 . coefficients - The array of FEM basis coefficients for the elements 1377 . probAux - The `PetscDS` specifying the auxiliary discretizations 1378 - coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 1379 1380 Output Parameter: 1381 . integral - the integral for this field 1382 1383 Level: intermediate 1384 1385 Developer Note: 1386 The function name begins with `PetscFE` and yet the first argument is `PetscDS` and it has no `PetscFE` arguments. 1387 1388 .seealso: `PetscFE`, `PetscDS`, `PetscFEIntegrateResidual()`, `PetscFEIntegrate()` 1389 @*/ 1390 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[]) 1391 { 1392 PetscFE fe; 1393 1394 PetscFunctionBegin; 1395 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1396 PetscCall(PetscDSGetDiscretization(prob, field, (PetscObject *)&fe)); 1397 if (fe->ops->integratebd) PetscCall((*fe->ops->integratebd)(prob, field, obj_func, Ne, geom, coefficients, probAux, coefficientsAux, integral)); 1398 PetscFunctionReturn(PETSC_SUCCESS); 1399 } 1400 1401 /*@C 1402 PetscFEIntegrateResidual - Produce the element residual vector for a chunk of elements by quadrature integration 1403 1404 Not Collective 1405 1406 Input Parameters: 1407 + ds - The `PetscDS` specifying the discretizations and continuum functions 1408 . key - The (label+value, field) being integrated 1409 . Ne - The number of elements in the chunk 1410 . cgeom - The cell geometry for each cell in the chunk 1411 . coefficients - The array of FEM basis coefficients for the elements 1412 . coefficients_t - The array of FEM basis time derivative coefficients for the elements 1413 . probAux - The `PetscDS` specifying the auxiliary discretizations 1414 . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 1415 - t - The time 1416 1417 Output Parameter: 1418 . elemVec - the element residual vectors from each element 1419 1420 Level: intermediate 1421 1422 Note: 1423 .vb 1424 Loop over batch of elements (e): 1425 Loop over quadrature points (q): 1426 Make u_q and gradU_q (loops over fields,Nb,Ncomp) and x_q 1427 Call f_0 and f_1 1428 Loop over element vector entries (f,fc --> i): 1429 elemVec[i] += \psi^{fc}_f(q) f0_{fc}(u, \nabla u) + \nabla\psi^{fc}_f(q) \cdot f1_{fc,df}(u, \nabla u) 1430 .ve 1431 1432 .seealso: `PetscFEIntegrateResidual()` 1433 @*/ 1434 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[]) 1435 { 1436 PetscFE fe; 1437 1438 PetscFunctionBeginHot; 1439 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1440 PetscCall(PetscDSGetDiscretization(ds, key.field, (PetscObject *)&fe)); 1441 if (fe->ops->integrateresidual) PetscCall((*fe->ops->integrateresidual)(ds, key, Ne, cgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec)); 1442 PetscFunctionReturn(PETSC_SUCCESS); 1443 } 1444 1445 /*@C 1446 PetscFEIntegrateBdResidual - Produce the element residual vector for a chunk of elements by quadrature integration over a boundary 1447 1448 Not Collective 1449 1450 Input Parameters: 1451 + ds - The `PetscDS` specifying the discretizations and continuum functions 1452 . wf - The PetscWeakForm object holding the pointwise functions 1453 . key - The (label+value, field) being integrated 1454 . Ne - The number of elements in the chunk 1455 . fgeom - The face geometry for each cell in the chunk 1456 . coefficients - The array of FEM basis coefficients for the elements 1457 . coefficients_t - The array of FEM basis time derivative coefficients for the elements 1458 . probAux - The `PetscDS` specifying the auxiliary discretizations 1459 . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 1460 - t - The time 1461 1462 Output Parameter: 1463 . elemVec - the element residual vectors from each element 1464 1465 Level: intermediate 1466 1467 .seealso: `PetscFEIntegrateResidual()` 1468 @*/ 1469 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[]) 1470 { 1471 PetscFE fe; 1472 1473 PetscFunctionBegin; 1474 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1475 PetscCall(PetscDSGetDiscretization(ds, key.field, (PetscObject *)&fe)); 1476 if (fe->ops->integratebdresidual) PetscCall((*fe->ops->integratebdresidual)(ds, wf, key, Ne, fgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec)); 1477 PetscFunctionReturn(PETSC_SUCCESS); 1478 } 1479 1480 /*@C 1481 PetscFEIntegrateHybridResidual - Produce the element residual vector for a chunk of hybrid element faces by quadrature integration 1482 1483 Not Collective 1484 1485 Input Parameters: 1486 + prob - The `PetscDS` specifying the discretizations and continuum functions 1487 . key - The (label+value, field) being integrated 1488 . s - The side of the cell being integrated, 0 for negative and 1 for positive 1489 . Ne - The number of elements in the chunk 1490 . fgeom - The face geometry for each cell in the chunk 1491 . coefficients - The array of FEM basis coefficients for the elements 1492 . coefficients_t - The array of FEM basis time derivative coefficients for the elements 1493 . probAux - The `PetscDS` specifying the auxiliary discretizations 1494 . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 1495 - t - The time 1496 1497 Output Parameter 1498 . elemVec - the element residual vectors from each element 1499 1500 Level: developer 1501 1502 .seealso: `PetscFEIntegrateResidual()` 1503 @*/ 1504 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[]) 1505 { 1506 PetscFE fe; 1507 1508 PetscFunctionBegin; 1509 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1); 1510 PetscCall(PetscDSGetDiscretization(prob, key.field, (PetscObject *)&fe)); 1511 if (fe->ops->integratehybridresidual) PetscCall((*fe->ops->integratehybridresidual)(prob, key, s, Ne, fgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec)); 1512 PetscFunctionReturn(PETSC_SUCCESS); 1513 } 1514 1515 /*@C 1516 PetscFEIntegrateJacobian - Produce the element Jacobian for a chunk of elements by quadrature integration 1517 1518 Not Collective 1519 1520 Input Parameters: 1521 + ds - The `PetscDS` specifying the discretizations and continuum functions 1522 . jtype - The type of matrix pointwise functions that should be used 1523 . key - The (label+value, fieldI*Nf + fieldJ) being integrated 1524 . s - The side of the cell being integrated, 0 for negative and 1 for positive 1525 . Ne - The number of elements in the chunk 1526 . cgeom - The cell geometry for each cell in the chunk 1527 . coefficients - The array of FEM basis coefficients for the elements for the Jacobian evaluation point 1528 . coefficients_t - The array of FEM basis time derivative coefficients for the elements 1529 . probAux - The `PetscDS` specifying the auxiliary discretizations 1530 . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 1531 . t - The time 1532 - u_tShift - A multiplier for the dF/du_t term (as opposed to the dF/du term) 1533 1534 Output Parameter: 1535 . elemMat - the element matrices for the Jacobian from each element 1536 1537 Level: intermediate 1538 1539 Note: 1540 .vb 1541 Loop over batch of elements (e): 1542 Loop over element matrix entries (f,fc,g,gc --> i,j): 1543 Loop over quadrature points (q): 1544 Make u_q and gradU_q (loops over fields,Nb,Ncomp) 1545 elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q) 1546 + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q) 1547 + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q) 1548 + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q) 1549 .ve 1550 1551 .seealso: `PetscFEIntegrateResidual()` 1552 @*/ 1553 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[]) 1554 { 1555 PetscFE fe; 1556 PetscInt Nf; 1557 1558 PetscFunctionBegin; 1559 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1560 PetscCall(PetscDSGetNumFields(ds, &Nf)); 1561 PetscCall(PetscDSGetDiscretization(ds, key.field / Nf, (PetscObject *)&fe)); 1562 if (fe->ops->integratejacobian) PetscCall((*fe->ops->integratejacobian)(ds, jtype, key, Ne, cgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, u_tshift, elemMat)); 1563 PetscFunctionReturn(PETSC_SUCCESS); 1564 } 1565 1566 /*@C 1567 PetscFEIntegrateBdJacobian - Produce the boundary element Jacobian for a chunk of elements by quadrature integration 1568 1569 Not Collective 1570 1571 Input Parameters: 1572 + ds - The `PetscDS` specifying the discretizations and continuum functions 1573 . wf - The PetscWeakForm holding the pointwise functions 1574 . key - The (label+value, fieldI*Nf + fieldJ) being integrated 1575 . Ne - The number of elements in the chunk 1576 . fgeom - The face geometry for each cell in the chunk 1577 . coefficients - The array of FEM basis coefficients for the elements for the Jacobian evaluation point 1578 . coefficients_t - The array of FEM basis time derivative coefficients for the elements 1579 . probAux - The `PetscDS` specifying the auxiliary discretizations 1580 . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 1581 . t - The time 1582 - u_tShift - A multiplier for the dF/du_t term (as opposed to the dF/du term) 1583 1584 Output Parameter: 1585 . elemMat - the element matrices for the Jacobian from each element 1586 1587 Level: intermediate 1588 1589 Note: 1590 .vb 1591 Loop over batch of elements (e): 1592 Loop over element matrix entries (f,fc,g,gc --> i,j): 1593 Loop over quadrature points (q): 1594 Make u_q and gradU_q (loops over fields,Nb,Ncomp) 1595 elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q) 1596 + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q) 1597 + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q) 1598 + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q) 1599 .ve 1600 1601 .seealso: `PetscFEIntegrateJacobian()`, `PetscFEIntegrateResidual()` 1602 @*/ 1603 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[]) 1604 { 1605 PetscFE fe; 1606 PetscInt Nf; 1607 1608 PetscFunctionBegin; 1609 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1610 PetscCall(PetscDSGetNumFields(ds, &Nf)); 1611 PetscCall(PetscDSGetDiscretization(ds, key.field / Nf, (PetscObject *)&fe)); 1612 if (fe->ops->integratebdjacobian) PetscCall((*fe->ops->integratebdjacobian)(ds, wf, key, Ne, fgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, u_tshift, elemMat)); 1613 PetscFunctionReturn(PETSC_SUCCESS); 1614 } 1615 1616 /*@C 1617 PetscFEIntegrateHybridJacobian - Produce the boundary element Jacobian for a chunk of hybrid elements by quadrature integration 1618 1619 Not Collective 1620 1621 Input Parameters: 1622 + ds - The `PetscDS` specifying the discretizations and continuum functions 1623 . jtype - The type of matrix pointwise functions that should be used 1624 . key - The (label+value, fieldI*Nf + fieldJ) being integrated 1625 . s - The side of the cell being integrated, 0 for negative and 1 for positive 1626 . Ne - The number of elements in the chunk 1627 . fgeom - The face geometry for each cell in the chunk 1628 . coefficients - The array of FEM basis coefficients for the elements for the Jacobian evaluation point 1629 . coefficients_t - The array of FEM basis time derivative coefficients for the elements 1630 . probAux - The `PetscDS` specifying the auxiliary discretizations 1631 . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 1632 . t - The time 1633 - u_tShift - A multiplier for the dF/du_t term (as opposed to the dF/du term) 1634 1635 Output Parameter 1636 . elemMat - the element matrices for the Jacobian from each element 1637 1638 Level: developer 1639 1640 Note: 1641 .vb 1642 Loop over batch of elements (e): 1643 Loop over element matrix entries (f,fc,g,gc --> i,j): 1644 Loop over quadrature points (q): 1645 Make u_q and gradU_q (loops over fields,Nb,Ncomp) 1646 elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q) 1647 + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q) 1648 + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q) 1649 + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q) 1650 .ve 1651 1652 .seealso: `PetscFEIntegrateJacobian()`, `PetscFEIntegrateResidual()` 1653 @*/ 1654 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[]) 1655 { 1656 PetscFE fe; 1657 PetscInt Nf; 1658 1659 PetscFunctionBegin; 1660 PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1); 1661 PetscCall(PetscDSGetNumFields(ds, &Nf)); 1662 PetscCall(PetscDSGetDiscretization(ds, key.field / Nf, (PetscObject *)&fe)); 1663 if (fe->ops->integratehybridjacobian) PetscCall((*fe->ops->integratehybridjacobian)(ds, jtype, key, s, Ne, fgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, u_tshift, elemMat)); 1664 PetscFunctionReturn(PETSC_SUCCESS); 1665 } 1666 1667 /*@ 1668 PetscFEGetHeightSubspace - Get the subspace of this space for a mesh point of a given height 1669 1670 Input Parameters: 1671 + fe - The finite element space 1672 - height - The height of the `DMPLEX` point 1673 1674 Output Parameter: 1675 . subfe - The subspace of this `PetscFE` space 1676 1677 Level: advanced 1678 1679 Note: 1680 For example, if we want the subspace of this space for a face, we would choose height = 1. 1681 1682 .seealso: `PetscFECreateDefault()` 1683 @*/ 1684 PetscErrorCode PetscFEGetHeightSubspace(PetscFE fe, PetscInt height, PetscFE *subfe) 1685 { 1686 PetscSpace P, subP; 1687 PetscDualSpace Q, subQ; 1688 PetscQuadrature subq; 1689 PetscFEType fetype; 1690 PetscInt dim, Nc; 1691 1692 PetscFunctionBegin; 1693 PetscValidHeaderSpecific(fe, PETSCFE_CLASSID, 1); 1694 PetscValidPointer(subfe, 3); 1695 if (height == 0) { 1696 *subfe = fe; 1697 PetscFunctionReturn(PETSC_SUCCESS); 1698 } 1699 PetscCall(PetscFEGetBasisSpace(fe, &P)); 1700 PetscCall(PetscFEGetDualSpace(fe, &Q)); 1701 PetscCall(PetscFEGetNumComponents(fe, &Nc)); 1702 PetscCall(PetscFEGetFaceQuadrature(fe, &subq)); 1703 PetscCall(PetscDualSpaceGetDimension(Q, &dim)); 1704 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); 1705 if (!fe->subspaces) PetscCall(PetscCalloc1(dim, &fe->subspaces)); 1706 if (height <= dim) { 1707 if (!fe->subspaces[height - 1]) { 1708 PetscFE sub = NULL; 1709 const char *name; 1710 1711 PetscCall(PetscSpaceGetHeightSubspace(P, height, &subP)); 1712 PetscCall(PetscDualSpaceGetHeightSubspace(Q, height, &subQ)); 1713 if (subQ) { 1714 PetscCall(PetscFECreate(PetscObjectComm((PetscObject)fe), &sub)); 1715 PetscCall(PetscObjectGetName((PetscObject)fe, &name)); 1716 PetscCall(PetscObjectSetName((PetscObject)sub, name)); 1717 PetscCall(PetscFEGetType(fe, &fetype)); 1718 PetscCall(PetscFESetType(sub, fetype)); 1719 PetscCall(PetscFESetBasisSpace(sub, subP)); 1720 PetscCall(PetscFESetDualSpace(sub, subQ)); 1721 PetscCall(PetscFESetNumComponents(sub, Nc)); 1722 PetscCall(PetscFESetUp(sub)); 1723 PetscCall(PetscFESetQuadrature(sub, subq)); 1724 } 1725 fe->subspaces[height - 1] = sub; 1726 } 1727 *subfe = fe->subspaces[height - 1]; 1728 } else { 1729 *subfe = NULL; 1730 } 1731 PetscFunctionReturn(PETSC_SUCCESS); 1732 } 1733 1734 /*@ 1735 PetscFERefine - Create a "refined" `PetscFE` object that refines the reference cell into smaller copies. This is typically used 1736 to precondition a higher order method with a lower order method on a refined mesh having the same number of dofs (but more 1737 sparsity). It is also used to create an interpolation between regularly refined meshes. 1738 1739 Collective 1740 1741 Input Parameter: 1742 . fe - The initial `PetscFE` 1743 1744 Output Parameter: 1745 . feRef - The refined `PetscFE` 1746 1747 Level: advanced 1748 1749 .seealso: `PetscFEType`, `PetscFECreate()`, `PetscFESetType()` 1750 @*/ 1751 PetscErrorCode PetscFERefine(PetscFE fe, PetscFE *feRef) 1752 { 1753 PetscSpace P, Pref; 1754 PetscDualSpace Q, Qref; 1755 DM K, Kref; 1756 PetscQuadrature q, qref; 1757 const PetscReal *v0, *jac; 1758 PetscInt numComp, numSubelements; 1759 PetscInt cStart, cEnd, c; 1760 PetscDualSpace *cellSpaces; 1761 1762 PetscFunctionBegin; 1763 PetscCall(PetscFEGetBasisSpace(fe, &P)); 1764 PetscCall(PetscFEGetDualSpace(fe, &Q)); 1765 PetscCall(PetscFEGetQuadrature(fe, &q)); 1766 PetscCall(PetscDualSpaceGetDM(Q, &K)); 1767 /* Create space */ 1768 PetscCall(PetscObjectReference((PetscObject)P)); 1769 Pref = P; 1770 /* Create dual space */ 1771 PetscCall(PetscDualSpaceDuplicate(Q, &Qref)); 1772 PetscCall(PetscDualSpaceSetType(Qref, PETSCDUALSPACEREFINED)); 1773 PetscCall(DMRefine(K, PetscObjectComm((PetscObject)fe), &Kref)); 1774 PetscCall(PetscDualSpaceSetDM(Qref, Kref)); 1775 PetscCall(DMPlexGetHeightStratum(Kref, 0, &cStart, &cEnd)); 1776 PetscCall(PetscMalloc1(cEnd - cStart, &cellSpaces)); 1777 /* TODO: fix for non-uniform refinement */ 1778 for (c = 0; c < cEnd - cStart; c++) cellSpaces[c] = Q; 1779 PetscCall(PetscDualSpaceRefinedSetCellSpaces(Qref, cellSpaces)); 1780 PetscCall(PetscFree(cellSpaces)); 1781 PetscCall(DMDestroy(&Kref)); 1782 PetscCall(PetscDualSpaceSetUp(Qref)); 1783 /* Create element */ 1784 PetscCall(PetscFECreate(PetscObjectComm((PetscObject)fe), feRef)); 1785 PetscCall(PetscFESetType(*feRef, PETSCFECOMPOSITE)); 1786 PetscCall(PetscFESetBasisSpace(*feRef, Pref)); 1787 PetscCall(PetscFESetDualSpace(*feRef, Qref)); 1788 PetscCall(PetscFEGetNumComponents(fe, &numComp)); 1789 PetscCall(PetscFESetNumComponents(*feRef, numComp)); 1790 PetscCall(PetscFESetUp(*feRef)); 1791 PetscCall(PetscSpaceDestroy(&Pref)); 1792 PetscCall(PetscDualSpaceDestroy(&Qref)); 1793 /* Create quadrature */ 1794 PetscCall(PetscFECompositeGetMapping(*feRef, &numSubelements, &v0, &jac, NULL)); 1795 PetscCall(PetscQuadratureExpandComposite(q, numSubelements, v0, jac, &qref)); 1796 PetscCall(PetscFESetQuadrature(*feRef, qref)); 1797 PetscCall(PetscQuadratureDestroy(&qref)); 1798 PetscFunctionReturn(PETSC_SUCCESS); 1799 } 1800 1801 static PetscErrorCode PetscFESetDefaultName_Private(PetscFE fe) 1802 { 1803 PetscSpace P; 1804 PetscDualSpace Q; 1805 DM K; 1806 DMPolytopeType ct; 1807 PetscInt degree; 1808 char name[64]; 1809 1810 PetscFunctionBegin; 1811 PetscCall(PetscFEGetBasisSpace(fe, &P)); 1812 PetscCall(PetscSpaceGetDegree(P, °ree, NULL)); 1813 PetscCall(PetscFEGetDualSpace(fe, &Q)); 1814 PetscCall(PetscDualSpaceGetDM(Q, &K)); 1815 PetscCall(DMPlexGetCellType(K, 0, &ct)); 1816 switch (ct) { 1817 case DM_POLYTOPE_SEGMENT: 1818 case DM_POLYTOPE_POINT_PRISM_TENSOR: 1819 case DM_POLYTOPE_QUADRILATERAL: 1820 case DM_POLYTOPE_SEG_PRISM_TENSOR: 1821 case DM_POLYTOPE_HEXAHEDRON: 1822 case DM_POLYTOPE_QUAD_PRISM_TENSOR: 1823 PetscCall(PetscSNPrintf(name, sizeof(name), "Q%" PetscInt_FMT, degree)); 1824 break; 1825 case DM_POLYTOPE_TRIANGLE: 1826 case DM_POLYTOPE_TETRAHEDRON: 1827 PetscCall(PetscSNPrintf(name, sizeof(name), "P%" PetscInt_FMT, degree)); 1828 break; 1829 case DM_POLYTOPE_TRI_PRISM: 1830 case DM_POLYTOPE_TRI_PRISM_TENSOR: 1831 PetscCall(PetscSNPrintf(name, sizeof(name), "P%" PetscInt_FMT "xQ%" PetscInt_FMT, degree, degree)); 1832 break; 1833 default: 1834 PetscCall(PetscSNPrintf(name, sizeof(name), "FE")); 1835 } 1836 PetscCall(PetscFESetName(fe, name)); 1837 PetscFunctionReturn(PETSC_SUCCESS); 1838 } 1839 1840 static PetscErrorCode PetscFECreateDefaultQuadrature_Private(PetscInt dim, DMPolytopeType ct, PetscInt qorder, PetscQuadrature *q, PetscQuadrature *fq) 1841 { 1842 const PetscInt quadPointsPerEdge = PetscMax(qorder + 1, 1); 1843 1844 PetscFunctionBegin; 1845 switch (ct) { 1846 case DM_POLYTOPE_SEGMENT: 1847 case DM_POLYTOPE_POINT_PRISM_TENSOR: 1848 case DM_POLYTOPE_QUADRILATERAL: 1849 case DM_POLYTOPE_SEG_PRISM_TENSOR: 1850 case DM_POLYTOPE_HEXAHEDRON: 1851 case DM_POLYTOPE_QUAD_PRISM_TENSOR: 1852 PetscCall(PetscDTGaussTensorQuadrature(dim, 1, quadPointsPerEdge, -1.0, 1.0, q)); 1853 PetscCall(PetscDTGaussTensorQuadrature(dim - 1, 1, quadPointsPerEdge, -1.0, 1.0, fq)); 1854 break; 1855 case DM_POLYTOPE_TRIANGLE: 1856 case DM_POLYTOPE_TETRAHEDRON: 1857 PetscCall(PetscDTSimplexQuadrature(dim, 2 * qorder, PETSCDTSIMPLEXQUAD_DEFAULT, q)); 1858 PetscCall(PetscDTSimplexQuadrature(dim - 1, 2 * qorder, PETSCDTSIMPLEXQUAD_DEFAULT, fq)); 1859 break; 1860 case DM_POLYTOPE_TRI_PRISM: 1861 case DM_POLYTOPE_TRI_PRISM_TENSOR: { 1862 PetscQuadrature q1, q2; 1863 1864 // TODO: this should be able to use symmetric rules, but doing so causes tests to fail 1865 PetscCall(PetscDTSimplexQuadrature(2, 2 * qorder, PETSCDTSIMPLEXQUAD_CONIC, &q1)); 1866 PetscCall(PetscDTGaussTensorQuadrature(1, 1, quadPointsPerEdge, -1.0, 1.0, &q2)); 1867 PetscCall(PetscDTTensorQuadratureCreate(q1, q2, q)); 1868 PetscCall(PetscQuadratureDestroy(&q2)); 1869 *fq = q1; 1870 /* TODO Need separate quadratures for each face */ 1871 } break; 1872 default: 1873 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No quadrature for celltype %s", DMPolytopeTypes[PetscMin(ct, DM_POLYTOPE_UNKNOWN)]); 1874 } 1875 PetscFunctionReturn(PETSC_SUCCESS); 1876 } 1877 1878 /*@ 1879 PetscFECreateFromSpaces - Create a `PetscFE` from the basis and dual spaces 1880 1881 Collective 1882 1883 Input Parameters: 1884 + P - The basis space 1885 . Q - The dual space 1886 . q - The cell quadrature 1887 - fq - The face quadrature 1888 1889 Output Parameter: 1890 . fem - The `PetscFE` object 1891 1892 Level: beginner 1893 1894 Note: 1895 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. 1896 1897 .seealso: `PetscFE`, `PetscSpace`, `PetscDualSpace`, `PetscQuadrature`, 1898 `PetscFECreateLagrangeByCell()`, `PetscFECreateDefault()`, `PetscFECreateByCell()`, `PetscFECreate()`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()` 1899 @*/ 1900 PetscErrorCode PetscFECreateFromSpaces(PetscSpace P, PetscDualSpace Q, PetscQuadrature q, PetscQuadrature fq, PetscFE *fem) 1901 { 1902 PetscInt Nc; 1903 const char *prefix; 1904 1905 PetscFunctionBegin; 1906 PetscCall(PetscFECreate(PetscObjectComm((PetscObject)P), fem)); 1907 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)P, &prefix)); 1908 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*fem, prefix)); 1909 PetscCall(PetscFESetType(*fem, PETSCFEBASIC)); 1910 PetscCall(PetscFESetBasisSpace(*fem, P)); 1911 PetscCall(PetscFESetDualSpace(*fem, Q)); 1912 PetscCall(PetscSpaceGetNumComponents(P, &Nc)); 1913 PetscCall(PetscFESetNumComponents(*fem, Nc)); 1914 PetscCall(PetscFESetUp(*fem)); 1915 PetscCall(PetscSpaceDestroy(&P)); 1916 PetscCall(PetscDualSpaceDestroy(&Q)); 1917 PetscCall(PetscFESetQuadrature(*fem, q)); 1918 PetscCall(PetscFESetFaceQuadrature(*fem, fq)); 1919 PetscCall(PetscQuadratureDestroy(&q)); 1920 PetscCall(PetscQuadratureDestroy(&fq)); 1921 PetscCall(PetscFESetDefaultName_Private(*fem)); 1922 PetscFunctionReturn(PETSC_SUCCESS); 1923 } 1924 1925 static PetscErrorCode PetscFECreate_Internal(MPI_Comm comm, PetscInt dim, PetscInt Nc, DMPolytopeType ct, const char prefix[], PetscInt degree, PetscInt qorder, PetscBool setFromOptions, PetscFE *fem) 1926 { 1927 DM K; 1928 PetscSpace P; 1929 PetscDualSpace Q; 1930 PetscQuadrature q, fq; 1931 PetscBool tensor; 1932 1933 PetscFunctionBegin; 1934 if (prefix) PetscValidCharPointer(prefix, 5); 1935 PetscValidPointer(fem, 9); 1936 switch (ct) { 1937 case DM_POLYTOPE_SEGMENT: 1938 case DM_POLYTOPE_POINT_PRISM_TENSOR: 1939 case DM_POLYTOPE_QUADRILATERAL: 1940 case DM_POLYTOPE_SEG_PRISM_TENSOR: 1941 case DM_POLYTOPE_HEXAHEDRON: 1942 case DM_POLYTOPE_QUAD_PRISM_TENSOR: 1943 tensor = PETSC_TRUE; 1944 break; 1945 default: 1946 tensor = PETSC_FALSE; 1947 } 1948 /* Create space */ 1949 PetscCall(PetscSpaceCreate(comm, &P)); 1950 PetscCall(PetscSpaceSetType(P, PETSCSPACEPOLYNOMIAL)); 1951 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)P, prefix)); 1952 PetscCall(PetscSpacePolynomialSetTensor(P, tensor)); 1953 PetscCall(PetscSpaceSetNumComponents(P, Nc)); 1954 PetscCall(PetscSpaceSetNumVariables(P, dim)); 1955 if (degree >= 0) { 1956 PetscCall(PetscSpaceSetDegree(P, degree, PETSC_DETERMINE)); 1957 if (ct == DM_POLYTOPE_TRI_PRISM || ct == DM_POLYTOPE_TRI_PRISM_TENSOR) { 1958 PetscSpace Pend, Pside; 1959 1960 PetscCall(PetscSpaceCreate(comm, &Pend)); 1961 PetscCall(PetscSpaceSetType(Pend, PETSCSPACEPOLYNOMIAL)); 1962 PetscCall(PetscSpacePolynomialSetTensor(Pend, PETSC_FALSE)); 1963 PetscCall(PetscSpaceSetNumComponents(Pend, Nc)); 1964 PetscCall(PetscSpaceSetNumVariables(Pend, dim - 1)); 1965 PetscCall(PetscSpaceSetDegree(Pend, degree, PETSC_DETERMINE)); 1966 PetscCall(PetscSpaceCreate(comm, &Pside)); 1967 PetscCall(PetscSpaceSetType(Pside, PETSCSPACEPOLYNOMIAL)); 1968 PetscCall(PetscSpacePolynomialSetTensor(Pside, PETSC_FALSE)); 1969 PetscCall(PetscSpaceSetNumComponents(Pside, 1)); 1970 PetscCall(PetscSpaceSetNumVariables(Pside, 1)); 1971 PetscCall(PetscSpaceSetDegree(Pside, degree, PETSC_DETERMINE)); 1972 PetscCall(PetscSpaceSetType(P, PETSCSPACETENSOR)); 1973 PetscCall(PetscSpaceTensorSetNumSubspaces(P, 2)); 1974 PetscCall(PetscSpaceTensorSetSubspace(P, 0, Pend)); 1975 PetscCall(PetscSpaceTensorSetSubspace(P, 1, Pside)); 1976 PetscCall(PetscSpaceDestroy(&Pend)); 1977 PetscCall(PetscSpaceDestroy(&Pside)); 1978 } 1979 } 1980 if (setFromOptions) PetscCall(PetscSpaceSetFromOptions(P)); 1981 PetscCall(PetscSpaceSetUp(P)); 1982 PetscCall(PetscSpaceGetDegree(P, °ree, NULL)); 1983 PetscCall(PetscSpacePolynomialGetTensor(P, &tensor)); 1984 PetscCall(PetscSpaceGetNumComponents(P, &Nc)); 1985 /* Create dual space */ 1986 PetscCall(PetscDualSpaceCreate(comm, &Q)); 1987 PetscCall(PetscDualSpaceSetType(Q, PETSCDUALSPACELAGRANGE)); 1988 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)Q, prefix)); 1989 PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, ct, &K)); 1990 PetscCall(PetscDualSpaceSetDM(Q, K)); 1991 PetscCall(DMDestroy(&K)); 1992 PetscCall(PetscDualSpaceSetNumComponents(Q, Nc)); 1993 PetscCall(PetscDualSpaceSetOrder(Q, degree)); 1994 /* TODO For some reason, we need a tensor dualspace with wedges */ 1995 PetscCall(PetscDualSpaceLagrangeSetTensor(Q, (tensor || (ct == DM_POLYTOPE_TRI_PRISM)) ? PETSC_TRUE : PETSC_FALSE)); 1996 if (setFromOptions) PetscCall(PetscDualSpaceSetFromOptions(Q)); 1997 PetscCall(PetscDualSpaceSetUp(Q)); 1998 /* Create quadrature */ 1999 qorder = qorder >= 0 ? qorder : degree; 2000 if (setFromOptions) { 2001 PetscObjectOptionsBegin((PetscObject)P); 2002 PetscCall(PetscOptionsBoundedInt("-petscfe_default_quadrature_order", "Quadrature order is one less than quadrature points per edge", "PetscFECreateDefault", qorder, &qorder, NULL, 0)); 2003 PetscOptionsEnd(); 2004 } 2005 PetscCall(PetscFECreateDefaultQuadrature_Private(dim, ct, qorder, &q, &fq)); 2006 /* Create finite element */ 2007 PetscCall(PetscFECreateFromSpaces(P, Q, q, fq, fem)); 2008 if (setFromOptions) PetscCall(PetscFESetFromOptions(*fem)); 2009 PetscFunctionReturn(PETSC_SUCCESS); 2010 } 2011 2012 /*@C 2013 PetscFECreateDefault - Create a `PetscFE` for basic FEM computation 2014 2015 Collective 2016 2017 Input Parameters: 2018 + comm - The MPI comm 2019 . dim - The spatial dimension 2020 . Nc - The number of components 2021 . isSimplex - Flag for simplex reference cell, otherwise its a tensor product 2022 . prefix - The options prefix, or `NULL` 2023 - qorder - The quadrature order or `PETSC_DETERMINE` to use `PetscSpace` polynomial degree 2024 2025 Output Parameter: 2026 . fem - The `PetscFE` object 2027 2028 Level: beginner 2029 2030 Note: 2031 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. 2032 2033 .seealso: `PetscFECreateLagrange()`, `PetscFECreateByCell()`, `PetscSpaceSetFromOptions()`, `PetscDualSpaceSetFromOptions()`, `PetscFESetFromOptions()`, `PetscFECreate()`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()` 2034 @*/ 2035 PetscErrorCode PetscFECreateDefault(MPI_Comm comm, PetscInt dim, PetscInt Nc, PetscBool isSimplex, const char prefix[], PetscInt qorder, PetscFE *fem) 2036 { 2037 PetscFunctionBegin; 2038 PetscCall(PetscFECreate_Internal(comm, dim, Nc, DMPolytopeTypeSimpleShape(dim, isSimplex), prefix, PETSC_DECIDE, qorder, PETSC_TRUE, fem)); 2039 PetscFunctionReturn(PETSC_SUCCESS); 2040 } 2041 2042 /*@C 2043 PetscFECreateByCell - Create a `PetscFE` for basic FEM computation 2044 2045 Collective 2046 2047 Input Parameters: 2048 + comm - The MPI comm 2049 . dim - The spatial dimension 2050 . Nc - The number of components 2051 . ct - The celltype of the reference cell 2052 . prefix - The options prefix, or `NULL` 2053 - qorder - The quadrature order or `PETSC_DETERMINE` to use `PetscSpace` polynomial degree 2054 2055 Output Parameter: 2056 . fem - The `PetscFE` object 2057 2058 Level: beginner 2059 2060 Note: 2061 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. 2062 2063 .seealso: `PetscFECreateDefault()`, `PetscFECreateLagrange()`, `PetscSpaceSetFromOptions()`, `PetscDualSpaceSetFromOptions()`, `PetscFESetFromOptions()`, `PetscFECreate()`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()` 2064 @*/ 2065 PetscErrorCode PetscFECreateByCell(MPI_Comm comm, PetscInt dim, PetscInt Nc, DMPolytopeType ct, const char prefix[], PetscInt qorder, PetscFE *fem) 2066 { 2067 PetscFunctionBegin; 2068 PetscCall(PetscFECreate_Internal(comm, dim, Nc, ct, prefix, PETSC_DECIDE, qorder, PETSC_TRUE, fem)); 2069 PetscFunctionReturn(PETSC_SUCCESS); 2070 } 2071 2072 /*@ 2073 PetscFECreateLagrange - Create a `PetscFE` for the basic Lagrange space of degree k 2074 2075 Collective 2076 2077 Input Parameters: 2078 + comm - The MPI comm 2079 . dim - The spatial dimension 2080 . Nc - The number of components 2081 . isSimplex - Flag for simplex reference cell, otherwise its a tensor product 2082 . k - The degree k of the space 2083 - qorder - The quadrature order or `PETSC_DETERMINE` to use `PetscSpace` polynomial degree 2084 2085 Output Parameter: 2086 . fem - The `PetscFE` object 2087 2088 Level: beginner 2089 2090 Note: 2091 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. 2092 2093 .seealso: `PetscFECreateLagrangeByCell()`, `PetscFECreateDefault()`, `PetscFECreateByCell()`, `PetscFECreate()`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()` 2094 @*/ 2095 PetscErrorCode PetscFECreateLagrange(MPI_Comm comm, PetscInt dim, PetscInt Nc, PetscBool isSimplex, PetscInt k, PetscInt qorder, PetscFE *fem) 2096 { 2097 PetscFunctionBegin; 2098 PetscCall(PetscFECreate_Internal(comm, dim, Nc, DMPolytopeTypeSimpleShape(dim, isSimplex), NULL, k, qorder, PETSC_FALSE, fem)); 2099 PetscFunctionReturn(PETSC_SUCCESS); 2100 } 2101 2102 /*@ 2103 PetscFECreateLagrangeByCell - Create a `PetscFE` for the basic Lagrange space of degree k 2104 2105 Collective 2106 2107 Input Parameters: 2108 + comm - The MPI comm 2109 . dim - The spatial dimension 2110 . Nc - The number of components 2111 . ct - The celltype of the reference cell 2112 . k - The degree k of the space 2113 - qorder - The quadrature order or `PETSC_DETERMINE` to use `PetscSpace` polynomial degree 2114 2115 Output Parameter: 2116 . fem - The `PetscFE` object 2117 2118 Level: beginner 2119 2120 Note: 2121 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. 2122 2123 .seealso: `PetscFECreateLagrange()`, `PetscFECreateDefault()`, `PetscFECreateByCell()`, `PetscFECreate()`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()` 2124 @*/ 2125 PetscErrorCode PetscFECreateLagrangeByCell(MPI_Comm comm, PetscInt dim, PetscInt Nc, DMPolytopeType ct, PetscInt k, PetscInt qorder, PetscFE *fem) 2126 { 2127 PetscFunctionBegin; 2128 PetscCall(PetscFECreate_Internal(comm, dim, Nc, ct, NULL, k, qorder, PETSC_FALSE, fem)); 2129 PetscFunctionReturn(PETSC_SUCCESS); 2130 } 2131 2132 /*@C 2133 PetscFESetName - Names the `PetscFE` and its subobjects 2134 2135 Not Collective 2136 2137 Input Parameters: 2138 + fe - The `PetscFE` 2139 - name - The name 2140 2141 Level: intermediate 2142 2143 .seealso: `PetscFECreate()`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()` 2144 @*/ 2145 PetscErrorCode PetscFESetName(PetscFE fe, const char name[]) 2146 { 2147 PetscSpace P; 2148 PetscDualSpace Q; 2149 2150 PetscFunctionBegin; 2151 PetscCall(PetscFEGetBasisSpace(fe, &P)); 2152 PetscCall(PetscFEGetDualSpace(fe, &Q)); 2153 PetscCall(PetscObjectSetName((PetscObject)fe, name)); 2154 PetscCall(PetscObjectSetName((PetscObject)P, name)); 2155 PetscCall(PetscObjectSetName((PetscObject)Q, name)); 2156 PetscFunctionReturn(PETSC_SUCCESS); 2157 } 2158 2159 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[]) 2160 { 2161 PetscInt dOffset = 0, fOffset = 0, f, g; 2162 2163 for (f = 0; f < Nf; ++f) { 2164 PetscCheck(r < T[f]->Nr, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Replica number %" PetscInt_FMT " should be in [0, %" PetscInt_FMT ")", r, T[f]->Nr); 2165 PetscCheck(q < T[f]->Np, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point number %" PetscInt_FMT " should be in [0, %" PetscInt_FMT ")", q, T[f]->Np); 2166 PetscFE fe; 2167 const PetscInt k = ds->jetDegree[f]; 2168 const PetscInt cdim = T[f]->cdim; 2169 const PetscInt Nq = T[f]->Np; 2170 const PetscInt Nbf = T[f]->Nb; 2171 const PetscInt Ncf = T[f]->Nc; 2172 const PetscReal *Bq = &T[f]->T[0][(r * Nq + q) * Nbf * Ncf]; 2173 const PetscReal *Dq = &T[f]->T[1][(r * Nq + q) * Nbf * Ncf * cdim]; 2174 const PetscReal *Hq = k > 1 ? &T[f]->T[2][(r * Nq + q) * Nbf * Ncf * cdim * cdim] : NULL; 2175 PetscInt hOffset = 0, b, c, d; 2176 2177 PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *)&fe)); 2178 for (c = 0; c < Ncf; ++c) u[fOffset + c] = 0.0; 2179 for (d = 0; d < cdim * Ncf; ++d) u_x[fOffset * cdim + d] = 0.0; 2180 for (b = 0; b < Nbf; ++b) { 2181 for (c = 0; c < Ncf; ++c) { 2182 const PetscInt cidx = b * Ncf + c; 2183 2184 u[fOffset + c] += Bq[cidx] * coefficients[dOffset + b]; 2185 for (d = 0; d < cdim; ++d) u_x[(fOffset + c) * cdim + d] += Dq[cidx * cdim + d] * coefficients[dOffset + b]; 2186 } 2187 } 2188 if (k > 1) { 2189 for (g = 0; g < Nf; ++g) hOffset += T[g]->Nc * cdim; 2190 for (d = 0; d < cdim * cdim * Ncf; ++d) u_x[hOffset + fOffset * cdim * cdim + d] = 0.0; 2191 for (b = 0; b < Nbf; ++b) { 2192 for (c = 0; c < Ncf; ++c) { 2193 const PetscInt cidx = b * Ncf + c; 2194 2195 for (d = 0; d < cdim * cdim; ++d) u_x[hOffset + (fOffset + c) * cdim * cdim + d] += Hq[cidx * cdim * cdim + d] * coefficients[dOffset + b]; 2196 } 2197 } 2198 PetscCall(PetscFEPushforwardHessian(fe, fegeom, 1, &u_x[hOffset + fOffset * cdim * cdim])); 2199 } 2200 PetscCall(PetscFEPushforward(fe, fegeom, 1, &u[fOffset])); 2201 PetscCall(PetscFEPushforwardGradient(fe, fegeom, 1, &u_x[fOffset * cdim])); 2202 if (u_t) { 2203 for (c = 0; c < Ncf; ++c) u_t[fOffset + c] = 0.0; 2204 for (b = 0; b < Nbf; ++b) { 2205 for (c = 0; c < Ncf; ++c) { 2206 const PetscInt cidx = b * Ncf + c; 2207 2208 u_t[fOffset + c] += Bq[cidx] * coefficients_t[dOffset + b]; 2209 } 2210 } 2211 PetscCall(PetscFEPushforward(fe, fegeom, 1, &u_t[fOffset])); 2212 } 2213 fOffset += Ncf; 2214 dOffset += Nbf; 2215 } 2216 return PETSC_SUCCESS; 2217 } 2218 2219 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[]) 2220 { 2221 PetscInt dOffset = 0, fOffset = 0, f, g; 2222 2223 /* f is the field number in the DS, g is the field number in u[] */ 2224 for (f = 0, g = 0; f < Nf; ++f) { 2225 PetscFE fe = (PetscFE)ds->disc[f]; 2226 const PetscInt dEt = T[f]->cdim; 2227 const PetscInt dE = fegeom->dimEmbed; 2228 const PetscInt Nq = T[f]->Np; 2229 const PetscInt Nbf = T[f]->Nb; 2230 const PetscInt Ncf = T[f]->Nc; 2231 const PetscReal *Bq = &T[f]->T[0][(r * Nq + q) * Nbf * Ncf]; 2232 const PetscReal *Dq = &T[f]->T[1][(r * Nq + q) * Nbf * Ncf * dEt]; 2233 PetscBool isCohesive; 2234 PetscInt Ns, s; 2235 2236 if (!T[f]) continue; 2237 PetscCall(PetscDSGetCohesive(ds, f, &isCohesive)); 2238 Ns = isCohesive ? 1 : 2; 2239 for (s = 0; s < Ns; ++s, ++g) { 2240 PetscInt b, c, d; 2241 2242 for (c = 0; c < Ncf; ++c) u[fOffset + c] = 0.0; 2243 for (d = 0; d < dE * Ncf; ++d) u_x[fOffset * dE + d] = 0.0; 2244 for (b = 0; b < Nbf; ++b) { 2245 for (c = 0; c < Ncf; ++c) { 2246 const PetscInt cidx = b * Ncf + c; 2247 2248 u[fOffset + c] += Bq[cidx] * coefficients[dOffset + b]; 2249 for (d = 0; d < dEt; ++d) u_x[(fOffset + c) * dE + d] += Dq[cidx * dEt + d] * coefficients[dOffset + b]; 2250 } 2251 } 2252 PetscCall(PetscFEPushforward(fe, fegeom, 1, &u[fOffset])); 2253 PetscCall(PetscFEPushforwardGradient(fe, fegeom, 1, &u_x[fOffset * dE])); 2254 if (u_t) { 2255 for (c = 0; c < Ncf; ++c) u_t[fOffset + c] = 0.0; 2256 for (b = 0; b < Nbf; ++b) { 2257 for (c = 0; c < Ncf; ++c) { 2258 const PetscInt cidx = b * Ncf + c; 2259 2260 u_t[fOffset + c] += Bq[cidx] * coefficients_t[dOffset + b]; 2261 } 2262 } 2263 PetscCall(PetscFEPushforward(fe, fegeom, 1, &u_t[fOffset])); 2264 } 2265 fOffset += Ncf; 2266 dOffset += Nbf; 2267 } 2268 } 2269 return PETSC_SUCCESS; 2270 } 2271 2272 PetscErrorCode PetscFEEvaluateFaceFields_Internal(PetscDS prob, PetscInt field, PetscInt faceLoc, const PetscScalar coefficients[], PetscScalar u[]) 2273 { 2274 PetscFE fe; 2275 PetscTabulation Tc; 2276 PetscInt b, c; 2277 2278 if (!prob) return PETSC_SUCCESS; 2279 PetscCall(PetscDSGetDiscretization(prob, field, (PetscObject *)&fe)); 2280 PetscCall(PetscFEGetFaceCentroidTabulation(fe, &Tc)); 2281 { 2282 const PetscReal *faceBasis = Tc->T[0]; 2283 const PetscInt Nb = Tc->Nb; 2284 const PetscInt Nc = Tc->Nc; 2285 2286 for (c = 0; c < Nc; ++c) u[c] = 0.0; 2287 for (b = 0; b < Nb; ++b) { 2288 for (c = 0; c < Nc; ++c) u[c] += coefficients[b] * faceBasis[(faceLoc * Nb + b) * Nc + c]; 2289 } 2290 } 2291 return PETSC_SUCCESS; 2292 } 2293 2294 PetscErrorCode PetscFEUpdateElementVec_Internal(PetscFE fe, PetscTabulation T, PetscInt r, PetscScalar tmpBasis[], PetscScalar tmpBasisDer[], PetscInt e, PetscFEGeom *fegeom, PetscScalar f0[], PetscScalar f1[], PetscScalar elemVec[]) 2295 { 2296 PetscFEGeom pgeom; 2297 const PetscInt dEt = T->cdim; 2298 const PetscInt dE = fegeom->dimEmbed; 2299 const PetscInt Nq = T->Np; 2300 const PetscInt Nb = T->Nb; 2301 const PetscInt Nc = T->Nc; 2302 const PetscReal *basis = &T->T[0][r * Nq * Nb * Nc]; 2303 const PetscReal *basisDer = &T->T[1][r * Nq * Nb * Nc * dEt]; 2304 PetscInt q, b, c, d; 2305 2306 for (q = 0; q < Nq; ++q) { 2307 for (b = 0; b < Nb; ++b) { 2308 for (c = 0; c < Nc; ++c) { 2309 const PetscInt bcidx = b * Nc + c; 2310 2311 tmpBasis[bcidx] = basis[q * Nb * Nc + bcidx]; 2312 for (d = 0; d < dEt; ++d) tmpBasisDer[bcidx * dE + d] = basisDer[q * Nb * Nc * dEt + bcidx * dEt + d]; 2313 for (d = dEt; d < dE; ++d) tmpBasisDer[bcidx * dE + d] = 0.0; 2314 } 2315 } 2316 PetscCall(PetscFEGeomGetCellPoint(fegeom, e, q, &pgeom)); 2317 PetscCall(PetscFEPushforward(fe, &pgeom, Nb, tmpBasis)); 2318 PetscCall(PetscFEPushforwardGradient(fe, &pgeom, Nb, tmpBasisDer)); 2319 for (b = 0; b < Nb; ++b) { 2320 for (c = 0; c < Nc; ++c) { 2321 const PetscInt bcidx = b * Nc + c; 2322 const PetscInt qcidx = q * Nc + c; 2323 2324 elemVec[b] += tmpBasis[bcidx] * f0[qcidx]; 2325 for (d = 0; d < dE; ++d) elemVec[b] += tmpBasisDer[bcidx * dE + d] * f1[qcidx * dE + d]; 2326 } 2327 } 2328 } 2329 return PETSC_SUCCESS; 2330 } 2331 2332 PetscErrorCode PetscFEUpdateElementVec_Hybrid_Internal(PetscFE fe, PetscTabulation T, PetscInt r, PetscInt s, PetscScalar tmpBasis[], PetscScalar tmpBasisDer[], PetscFEGeom *fegeom, PetscScalar f0[], PetscScalar f1[], PetscScalar elemVec[]) 2333 { 2334 const PetscInt dE = T->cdim; 2335 const PetscInt Nq = T->Np; 2336 const PetscInt Nb = T->Nb; 2337 const PetscInt Nc = T->Nc; 2338 const PetscReal *basis = &T->T[0][r * Nq * Nb * Nc]; 2339 const PetscReal *basisDer = &T->T[1][r * Nq * Nb * Nc * dE]; 2340 PetscInt q, b, c, d; 2341 2342 for (q = 0; q < Nq; ++q) { 2343 for (b = 0; b < Nb; ++b) { 2344 for (c = 0; c < Nc; ++c) { 2345 const PetscInt bcidx = b * Nc + c; 2346 2347 tmpBasis[bcidx] = basis[q * Nb * Nc + bcidx]; 2348 for (d = 0; d < dE; ++d) tmpBasisDer[bcidx * dE + d] = basisDer[q * Nb * Nc * dE + bcidx * dE + d]; 2349 } 2350 } 2351 PetscCall(PetscFEPushforward(fe, fegeom, Nb, tmpBasis)); 2352 PetscCall(PetscFEPushforwardGradient(fe, fegeom, Nb, tmpBasisDer)); 2353 for (b = 0; b < Nb; ++b) { 2354 for (c = 0; c < Nc; ++c) { 2355 const PetscInt bcidx = b * Nc + c; 2356 const PetscInt qcidx = q * Nc + c; 2357 2358 elemVec[Nb * s + b] += tmpBasis[bcidx] * f0[qcidx]; 2359 for (d = 0; d < dE; ++d) elemVec[Nb * s + b] += tmpBasisDer[bcidx * dE + d] * f1[qcidx * dE + d]; 2360 } 2361 } 2362 } 2363 return PETSC_SUCCESS; 2364 } 2365 2366 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[]) 2367 { 2368 const PetscInt dE = TI->cdim; 2369 const PetscInt NqI = TI->Np; 2370 const PetscInt NbI = TI->Nb; 2371 const PetscInt NcI = TI->Nc; 2372 const PetscReal *basisI = &TI->T[0][(r * NqI + q) * NbI * NcI]; 2373 const PetscReal *basisDerI = &TI->T[1][(r * NqI + q) * NbI * NcI * dE]; 2374 const PetscInt NqJ = TJ->Np; 2375 const PetscInt NbJ = TJ->Nb; 2376 const PetscInt NcJ = TJ->Nc; 2377 const PetscReal *basisJ = &TJ->T[0][(r * NqJ + q) * NbJ * NcJ]; 2378 const PetscReal *basisDerJ = &TJ->T[1][(r * NqJ + q) * NbJ * NcJ * dE]; 2379 PetscInt f, fc, g, gc, df, dg; 2380 2381 for (f = 0; f < NbI; ++f) { 2382 for (fc = 0; fc < NcI; ++fc) { 2383 const PetscInt fidx = f * NcI + fc; /* Test function basis index */ 2384 2385 tmpBasisI[fidx] = basisI[fidx]; 2386 for (df = 0; df < dE; ++df) tmpBasisDerI[fidx * dE + df] = basisDerI[fidx * dE + df]; 2387 } 2388 } 2389 PetscCall(PetscFEPushforward(feI, fegeom, NbI, tmpBasisI)); 2390 PetscCall(PetscFEPushforwardGradient(feI, fegeom, NbI, tmpBasisDerI)); 2391 for (g = 0; g < NbJ; ++g) { 2392 for (gc = 0; gc < NcJ; ++gc) { 2393 const PetscInt gidx = g * NcJ + gc; /* Trial function basis index */ 2394 2395 tmpBasisJ[gidx] = basisJ[gidx]; 2396 for (dg = 0; dg < dE; ++dg) tmpBasisDerJ[gidx * dE + dg] = basisDerJ[gidx * dE + dg]; 2397 } 2398 } 2399 PetscCall(PetscFEPushforward(feJ, fegeom, NbJ, tmpBasisJ)); 2400 PetscCall(PetscFEPushforwardGradient(feJ, fegeom, NbJ, tmpBasisDerJ)); 2401 for (f = 0; f < NbI; ++f) { 2402 for (fc = 0; fc < NcI; ++fc) { 2403 const PetscInt fidx = f * NcI + fc; /* Test function basis index */ 2404 const PetscInt i = offsetI + f; /* Element matrix row */ 2405 for (g = 0; g < NbJ; ++g) { 2406 for (gc = 0; gc < NcJ; ++gc) { 2407 const PetscInt gidx = g * NcJ + gc; /* Trial function basis index */ 2408 const PetscInt j = offsetJ + g; /* Element matrix column */ 2409 const PetscInt fOff = eOffset + i * totDim + j; 2410 2411 elemMat[fOff] += tmpBasisI[fidx] * g0[fc * NcJ + gc] * tmpBasisJ[gidx]; 2412 for (df = 0; df < dE; ++df) { 2413 elemMat[fOff] += tmpBasisI[fidx] * g1[(fc * NcJ + gc) * dE + df] * tmpBasisDerJ[gidx * dE + df]; 2414 elemMat[fOff] += tmpBasisDerI[fidx * dE + df] * g2[(fc * NcJ + gc) * dE + df] * tmpBasisJ[gidx]; 2415 for (dg = 0; dg < dE; ++dg) elemMat[fOff] += tmpBasisDerI[fidx * dE + df] * g3[((fc * NcJ + gc) * dE + df) * dE + dg] * tmpBasisDerJ[gidx * dE + dg]; 2416 } 2417 } 2418 } 2419 } 2420 } 2421 return PETSC_SUCCESS; 2422 } 2423 2424 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[]) 2425 { 2426 const PetscInt dE = TI->cdim; 2427 const PetscInt NqI = TI->Np; 2428 const PetscInt NbI = TI->Nb; 2429 const PetscInt NcI = TI->Nc; 2430 const PetscReal *basisI = &TI->T[0][(r * NqI + q) * NbI * NcI]; 2431 const PetscReal *basisDerI = &TI->T[1][(r * NqI + q) * NbI * NcI * dE]; 2432 const PetscInt NqJ = TJ->Np; 2433 const PetscInt NbJ = TJ->Nb; 2434 const PetscInt NcJ = TJ->Nc; 2435 const PetscReal *basisJ = &TJ->T[0][(r * NqJ + q) * NbJ * NcJ]; 2436 const PetscReal *basisDerJ = &TJ->T[1][(r * NqJ + q) * NbJ * NcJ * dE]; 2437 const PetscInt so = isHybridI ? 0 : s; 2438 const PetscInt to = isHybridJ ? 0 : s; 2439 PetscInt f, fc, g, gc, df, dg; 2440 2441 for (f = 0; f < NbI; ++f) { 2442 for (fc = 0; fc < NcI; ++fc) { 2443 const PetscInt fidx = f * NcI + fc; /* Test function basis index */ 2444 2445 tmpBasisI[fidx] = basisI[fidx]; 2446 for (df = 0; df < dE; ++df) tmpBasisDerI[fidx * dE + df] = basisDerI[fidx * dE + df]; 2447 } 2448 } 2449 PetscCall(PetscFEPushforward(feI, fegeom, NbI, tmpBasisI)); 2450 PetscCall(PetscFEPushforwardGradient(feI, fegeom, NbI, tmpBasisDerI)); 2451 for (g = 0; g < NbJ; ++g) { 2452 for (gc = 0; gc < NcJ; ++gc) { 2453 const PetscInt gidx = g * NcJ + gc; /* Trial function basis index */ 2454 2455 tmpBasisJ[gidx] = basisJ[gidx]; 2456 for (dg = 0; dg < dE; ++dg) tmpBasisDerJ[gidx * dE + dg] = basisDerJ[gidx * dE + dg]; 2457 } 2458 } 2459 PetscCall(PetscFEPushforward(feJ, fegeom, NbJ, tmpBasisJ)); 2460 PetscCall(PetscFEPushforwardGradient(feJ, fegeom, NbJ, tmpBasisDerJ)); 2461 for (f = 0; f < NbI; ++f) { 2462 for (fc = 0; fc < NcI; ++fc) { 2463 const PetscInt fidx = f * NcI + fc; /* Test function basis index */ 2464 const PetscInt i = offsetI + NbI * so + f; /* Element matrix row */ 2465 for (g = 0; g < NbJ; ++g) { 2466 for (gc = 0; gc < NcJ; ++gc) { 2467 const PetscInt gidx = g * NcJ + gc; /* Trial function basis index */ 2468 const PetscInt j = offsetJ + NbJ * to + g; /* Element matrix column */ 2469 const PetscInt fOff = eOffset + i * totDim + j; 2470 2471 elemMat[fOff] += tmpBasisI[fidx] * g0[fc * NcJ + gc] * tmpBasisJ[gidx]; 2472 for (df = 0; df < dE; ++df) { 2473 elemMat[fOff] += tmpBasisI[fidx] * g1[(fc * NcJ + gc) * dE + df] * tmpBasisDerJ[gidx * dE + df]; 2474 elemMat[fOff] += tmpBasisDerI[fidx * dE + df] * g2[(fc * NcJ + gc) * dE + df] * tmpBasisJ[gidx]; 2475 for (dg = 0; dg < dE; ++dg) elemMat[fOff] += tmpBasisDerI[fidx * dE + df] * g3[((fc * NcJ + gc) * dE + df) * dE + dg] * tmpBasisDerJ[gidx * dE + dg]; 2476 } 2477 } 2478 } 2479 } 2480 } 2481 return PETSC_SUCCESS; 2482 } 2483 2484 PetscErrorCode PetscFECreateCellGeometry(PetscFE fe, PetscQuadrature quad, PetscFEGeom *cgeom) 2485 { 2486 PetscDualSpace dsp; 2487 DM dm; 2488 PetscQuadrature quadDef; 2489 PetscInt dim, cdim, Nq; 2490 2491 PetscFunctionBegin; 2492 PetscCall(PetscFEGetDualSpace(fe, &dsp)); 2493 PetscCall(PetscDualSpaceGetDM(dsp, &dm)); 2494 PetscCall(DMGetDimension(dm, &dim)); 2495 PetscCall(DMGetCoordinateDim(dm, &cdim)); 2496 PetscCall(PetscFEGetQuadrature(fe, &quadDef)); 2497 quad = quad ? quad : quadDef; 2498 PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL)); 2499 PetscCall(PetscMalloc1(Nq * cdim, &cgeom->v)); 2500 PetscCall(PetscMalloc1(Nq * cdim * cdim, &cgeom->J)); 2501 PetscCall(PetscMalloc1(Nq * cdim * cdim, &cgeom->invJ)); 2502 PetscCall(PetscMalloc1(Nq, &cgeom->detJ)); 2503 cgeom->dim = dim; 2504 cgeom->dimEmbed = cdim; 2505 cgeom->numCells = 1; 2506 cgeom->numPoints = Nq; 2507 PetscCall(DMPlexComputeCellGeometryFEM(dm, 0, quad, cgeom->v, cgeom->J, cgeom->invJ, cgeom->detJ)); 2508 PetscFunctionReturn(PETSC_SUCCESS); 2509 } 2510 2511 PetscErrorCode PetscFEDestroyCellGeometry(PetscFE fe, PetscFEGeom *cgeom) 2512 { 2513 PetscFunctionBegin; 2514 PetscCall(PetscFree(cgeom->v)); 2515 PetscCall(PetscFree(cgeom->J)); 2516 PetscCall(PetscFree(cgeom->invJ)); 2517 PetscCall(PetscFree(cgeom->detJ)); 2518 PetscFunctionReturn(PETSC_SUCCESS); 2519 } 2520