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