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 .keywords: PetscFE, register 86 .seealso: PetscFERegisterAll(), PetscFERegisterDestroy() 87 88 @*/ 89 PetscErrorCode PetscFERegister(const char sname[], PetscErrorCode (*function)(PetscFE)) 90 { 91 PetscErrorCode ierr; 92 93 PetscFunctionBegin; 94 ierr = PetscFunctionListAdd(&PetscFEList, sname, function);CHKERRQ(ierr); 95 PetscFunctionReturn(0); 96 } 97 98 /*@C 99 PetscFESetType - Builds a particular PetscFE 100 101 Collective on PetscFE 102 103 Input Parameters: 104 + fem - The PetscFE object 105 - name - The kind of FEM space 106 107 Options Database Key: 108 . -petscfe_type <type> - Sets the PetscFE type; use -help for a list of available types 109 110 Level: intermediate 111 112 .keywords: PetscFE, set, type 113 .seealso: PetscFEGetType(), PetscFECreate() 114 @*/ 115 PetscErrorCode PetscFESetType(PetscFE fem, PetscFEType name) 116 { 117 PetscErrorCode (*r)(PetscFE); 118 PetscBool match; 119 PetscErrorCode ierr; 120 121 PetscFunctionBegin; 122 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 123 ierr = PetscObjectTypeCompare((PetscObject) fem, name, &match);CHKERRQ(ierr); 124 if (match) PetscFunctionReturn(0); 125 126 if (!PetscFERegisterAllCalled) {ierr = PetscFERegisterAll();CHKERRQ(ierr);} 127 ierr = PetscFunctionListFind(PetscFEList, name, &r);CHKERRQ(ierr); 128 if (!r) SETERRQ1(PetscObjectComm((PetscObject) fem), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscFE type: %s", name); 129 130 if (fem->ops->destroy) { 131 ierr = (*fem->ops->destroy)(fem);CHKERRQ(ierr); 132 fem->ops->destroy = NULL; 133 } 134 ierr = (*r)(fem);CHKERRQ(ierr); 135 ierr = PetscObjectChangeTypeName((PetscObject) fem, name);CHKERRQ(ierr); 136 PetscFunctionReturn(0); 137 } 138 139 /*@C 140 PetscFEGetType - Gets the PetscFE type name (as a string) from the object. 141 142 Not Collective 143 144 Input Parameter: 145 . fem - The PetscFE 146 147 Output Parameter: 148 . name - The PetscFE type name 149 150 Level: intermediate 151 152 .keywords: PetscFE, get, type, name 153 .seealso: PetscFESetType(), PetscFECreate() 154 @*/ 155 PetscErrorCode PetscFEGetType(PetscFE fem, PetscFEType *name) 156 { 157 PetscErrorCode ierr; 158 159 PetscFunctionBegin; 160 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 161 PetscValidPointer(name, 2); 162 if (!PetscFERegisterAllCalled) { 163 ierr = PetscFERegisterAll();CHKERRQ(ierr); 164 } 165 *name = ((PetscObject) fem)->type_name; 166 PetscFunctionReturn(0); 167 } 168 169 /*@C 170 PetscFEView - Views a PetscFE 171 172 Collective on PetscFE 173 174 Input Parameter: 175 + fem - the PetscFE object to view 176 - viewer - the viewer 177 178 Level: developer 179 180 .seealso PetscFEDestroy() 181 @*/ 182 PetscErrorCode PetscFEView(PetscFE fem, PetscViewer viewer) 183 { 184 PetscBool iascii; 185 PetscErrorCode ierr; 186 187 PetscFunctionBegin; 188 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 189 if (viewer) PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 190 if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) fem), &viewer);CHKERRQ(ierr);} 191 ierr = PetscObjectPrintClassNamePrefixType((PetscObject)fem, viewer);CHKERRQ(ierr); 192 ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 193 if (fem->ops->view) {ierr = (*fem->ops->view)(fem, viewer);CHKERRQ(ierr);} 194 PetscFunctionReturn(0); 195 } 196 197 /*@ 198 PetscFESetFromOptions - sets parameters in a PetscFE from the options database 199 200 Collective on PetscFE 201 202 Input Parameter: 203 . fem - the PetscFE object to set options for 204 205 Options Database: 206 . -petscfe_num_blocks the number of cell blocks to integrate concurrently 207 . -petscfe_num_batches the number of cell batches to integrate serially 208 209 Level: developer 210 211 .seealso PetscFEView() 212 @*/ 213 PetscErrorCode PetscFESetFromOptions(PetscFE fem) 214 { 215 const char *defaultType; 216 char name[256]; 217 PetscBool flg; 218 PetscErrorCode ierr; 219 220 PetscFunctionBegin; 221 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 222 if (!((PetscObject) fem)->type_name) { 223 defaultType = PETSCFEBASIC; 224 } else { 225 defaultType = ((PetscObject) fem)->type_name; 226 } 227 if (!PetscFERegisterAllCalled) {ierr = PetscFERegisterAll();CHKERRQ(ierr);} 228 229 ierr = PetscObjectOptionsBegin((PetscObject) fem);CHKERRQ(ierr); 230 ierr = PetscOptionsFList("-petscfe_type", "Finite element space", "PetscFESetType", PetscFEList, defaultType, name, 256, &flg);CHKERRQ(ierr); 231 if (flg) { 232 ierr = PetscFESetType(fem, name);CHKERRQ(ierr); 233 } else if (!((PetscObject) fem)->type_name) { 234 ierr = PetscFESetType(fem, defaultType);CHKERRQ(ierr); 235 } 236 ierr = PetscOptionsInt("-petscfe_num_blocks", "The number of cell blocks to integrate concurrently", "PetscSpaceSetTileSizes", fem->numBlocks, &fem->numBlocks, NULL);CHKERRQ(ierr); 237 ierr = PetscOptionsInt("-petscfe_num_batches", "The number of cell batches to integrate serially", "PetscSpaceSetTileSizes", fem->numBatches, &fem->numBatches, NULL);CHKERRQ(ierr); 238 if (fem->ops->setfromoptions) { 239 ierr = (*fem->ops->setfromoptions)(PetscOptionsObject,fem);CHKERRQ(ierr); 240 } 241 /* process any options handlers added with PetscObjectAddOptionsHandler() */ 242 ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) fem);CHKERRQ(ierr); 243 ierr = PetscOptionsEnd();CHKERRQ(ierr); 244 ierr = PetscFEViewFromOptions(fem, NULL, "-petscfe_view");CHKERRQ(ierr); 245 PetscFunctionReturn(0); 246 } 247 248 /*@C 249 PetscFESetUp - Construct data structures for the PetscFE 250 251 Collective on PetscFE 252 253 Input Parameter: 254 . fem - the PetscFE object to setup 255 256 Level: developer 257 258 .seealso PetscFEView(), PetscFEDestroy() 259 @*/ 260 PetscErrorCode PetscFESetUp(PetscFE fem) 261 { 262 PetscErrorCode ierr; 263 264 PetscFunctionBegin; 265 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 266 if (fem->setupcalled) PetscFunctionReturn(0); 267 fem->setupcalled = PETSC_TRUE; 268 if (fem->ops->setup) {ierr = (*fem->ops->setup)(fem);CHKERRQ(ierr);} 269 PetscFunctionReturn(0); 270 } 271 272 /*@ 273 PetscFEDestroy - Destroys a PetscFE object 274 275 Collective on PetscFE 276 277 Input Parameter: 278 . fem - the PetscFE object to destroy 279 280 Level: developer 281 282 .seealso PetscFEView() 283 @*/ 284 PetscErrorCode PetscFEDestroy(PetscFE *fem) 285 { 286 PetscErrorCode ierr; 287 288 PetscFunctionBegin; 289 if (!*fem) PetscFunctionReturn(0); 290 PetscValidHeaderSpecific((*fem), PETSCFE_CLASSID, 1); 291 292 if (--((PetscObject)(*fem))->refct > 0) {*fem = 0; PetscFunctionReturn(0);} 293 ((PetscObject) (*fem))->refct = 0; 294 295 if ((*fem)->subspaces) { 296 PetscInt dim, d; 297 298 ierr = PetscDualSpaceGetDimension((*fem)->dualSpace, &dim);CHKERRQ(ierr); 299 for (d = 0; d < dim; ++d) {ierr = PetscFEDestroy(&(*fem)->subspaces[d]);CHKERRQ(ierr);} 300 } 301 ierr = PetscFree((*fem)->subspaces);CHKERRQ(ierr); 302 ierr = PetscFree((*fem)->invV);CHKERRQ(ierr); 303 ierr = PetscFERestoreTabulation((*fem), 0, NULL, &(*fem)->B, &(*fem)->D, NULL /*&(*fem)->H*/);CHKERRQ(ierr); 304 ierr = PetscFERestoreTabulation((*fem), 0, NULL, &(*fem)->Bf, &(*fem)->Df, NULL /*&(*fem)->Hf*/);CHKERRQ(ierr); 305 ierr = PetscFERestoreTabulation((*fem), 0, NULL, &(*fem)->F, NULL, NULL);CHKERRQ(ierr); 306 ierr = PetscSpaceDestroy(&(*fem)->basisSpace);CHKERRQ(ierr); 307 ierr = PetscDualSpaceDestroy(&(*fem)->dualSpace);CHKERRQ(ierr); 308 ierr = PetscQuadratureDestroy(&(*fem)->quadrature);CHKERRQ(ierr); 309 ierr = PetscQuadratureDestroy(&(*fem)->faceQuadrature);CHKERRQ(ierr); 310 311 if ((*fem)->ops->destroy) {ierr = (*(*fem)->ops->destroy)(*fem);CHKERRQ(ierr);} 312 ierr = PetscHeaderDestroy(fem);CHKERRQ(ierr); 313 PetscFunctionReturn(0); 314 } 315 316 /*@ 317 PetscFECreate - Creates an empty PetscFE object. The type can then be set with PetscFESetType(). 318 319 Collective on MPI_Comm 320 321 Input Parameter: 322 . comm - The communicator for the PetscFE object 323 324 Output Parameter: 325 . fem - The PetscFE object 326 327 Level: beginner 328 329 .seealso: PetscFESetType(), PETSCFEGALERKIN 330 @*/ 331 PetscErrorCode PetscFECreate(MPI_Comm comm, PetscFE *fem) 332 { 333 PetscFE f; 334 PetscErrorCode ierr; 335 336 PetscFunctionBegin; 337 PetscValidPointer(fem, 2); 338 ierr = PetscCitationsRegister(FECitation,&FEcite);CHKERRQ(ierr); 339 *fem = NULL; 340 ierr = PetscFEInitializePackage();CHKERRQ(ierr); 341 342 ierr = PetscHeaderCreate(f, PETSCFE_CLASSID, "PetscFE", "Finite Element", "PetscFE", comm, PetscFEDestroy, PetscFEView);CHKERRQ(ierr); 343 344 f->basisSpace = NULL; 345 f->dualSpace = NULL; 346 f->numComponents = 1; 347 f->subspaces = NULL; 348 f->invV = NULL; 349 f->B = NULL; 350 f->D = NULL; 351 f->H = NULL; 352 f->Bf = NULL; 353 f->Df = NULL; 354 f->Hf = NULL; 355 ierr = PetscMemzero(&f->quadrature, sizeof(PetscQuadrature));CHKERRQ(ierr); 356 ierr = PetscMemzero(&f->faceQuadrature, sizeof(PetscQuadrature));CHKERRQ(ierr); 357 f->blockSize = 0; 358 f->numBlocks = 1; 359 f->batchSize = 0; 360 f->numBatches = 1; 361 362 *fem = f; 363 PetscFunctionReturn(0); 364 } 365 366 /*@ 367 PetscFEGetSpatialDimension - Returns the spatial dimension of the element 368 369 Not collective 370 371 Input Parameter: 372 . fem - The PetscFE object 373 374 Output Parameter: 375 . dim - The spatial dimension 376 377 Level: intermediate 378 379 .seealso: PetscFECreate() 380 @*/ 381 PetscErrorCode PetscFEGetSpatialDimension(PetscFE fem, PetscInt *dim) 382 { 383 DM dm; 384 PetscErrorCode ierr; 385 386 PetscFunctionBegin; 387 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 388 PetscValidPointer(dim, 2); 389 ierr = PetscDualSpaceGetDM(fem->dualSpace, &dm);CHKERRQ(ierr); 390 ierr = DMGetDimension(dm, dim);CHKERRQ(ierr); 391 PetscFunctionReturn(0); 392 } 393 394 /*@ 395 PetscFESetNumComponents - Sets the number of components in the element 396 397 Not collective 398 399 Input Parameters: 400 + fem - The PetscFE object 401 - comp - The number of field components 402 403 Level: intermediate 404 405 .seealso: PetscFECreate() 406 @*/ 407 PetscErrorCode PetscFESetNumComponents(PetscFE fem, PetscInt comp) 408 { 409 PetscFunctionBegin; 410 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 411 fem->numComponents = comp; 412 PetscFunctionReturn(0); 413 } 414 415 /*@ 416 PetscFEGetNumComponents - Returns the number of components in the element 417 418 Not collective 419 420 Input Parameter: 421 . fem - The PetscFE object 422 423 Output Parameter: 424 . comp - The number of field components 425 426 Level: intermediate 427 428 .seealso: PetscFECreate() 429 @*/ 430 PetscErrorCode PetscFEGetNumComponents(PetscFE fem, PetscInt *comp) 431 { 432 PetscFunctionBegin; 433 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 434 PetscValidPointer(comp, 2); 435 *comp = fem->numComponents; 436 PetscFunctionReturn(0); 437 } 438 439 /*@ 440 PetscFESetTileSizes - Sets the tile sizes for evaluation 441 442 Not collective 443 444 Input Parameters: 445 + fem - The PetscFE object 446 . blockSize - The number of elements in a block 447 . numBlocks - The number of blocks in a batch 448 . batchSize - The number of elements in a batch 449 - numBatches - The number of batches in a chunk 450 451 Level: intermediate 452 453 .seealso: PetscFECreate() 454 @*/ 455 PetscErrorCode PetscFESetTileSizes(PetscFE fem, PetscInt blockSize, PetscInt numBlocks, PetscInt batchSize, PetscInt numBatches) 456 { 457 PetscFunctionBegin; 458 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 459 fem->blockSize = blockSize; 460 fem->numBlocks = numBlocks; 461 fem->batchSize = batchSize; 462 fem->numBatches = numBatches; 463 PetscFunctionReturn(0); 464 } 465 466 /*@ 467 PetscFEGetTileSizes - Returns the tile sizes for evaluation 468 469 Not collective 470 471 Input Parameter: 472 . fem - The PetscFE object 473 474 Output Parameters: 475 + blockSize - The number of elements in a block 476 . numBlocks - The number of blocks in a batch 477 . batchSize - The number of elements in a batch 478 - numBatches - The number of batches in a chunk 479 480 Level: intermediate 481 482 .seealso: PetscFECreate() 483 @*/ 484 PetscErrorCode PetscFEGetTileSizes(PetscFE fem, PetscInt *blockSize, PetscInt *numBlocks, PetscInt *batchSize, PetscInt *numBatches) 485 { 486 PetscFunctionBegin; 487 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 488 if (blockSize) PetscValidPointer(blockSize, 2); 489 if (numBlocks) PetscValidPointer(numBlocks, 3); 490 if (batchSize) PetscValidPointer(batchSize, 4); 491 if (numBatches) PetscValidPointer(numBatches, 5); 492 if (blockSize) *blockSize = fem->blockSize; 493 if (numBlocks) *numBlocks = fem->numBlocks; 494 if (batchSize) *batchSize = fem->batchSize; 495 if (numBatches) *numBatches = fem->numBatches; 496 PetscFunctionReturn(0); 497 } 498 499 /*@ 500 PetscFEGetBasisSpace - Returns the PetscSpace used for approximation of the solution 501 502 Not collective 503 504 Input Parameter: 505 . fem - The PetscFE object 506 507 Output Parameter: 508 . sp - The PetscSpace object 509 510 Level: intermediate 511 512 .seealso: PetscFECreate() 513 @*/ 514 PetscErrorCode PetscFEGetBasisSpace(PetscFE fem, PetscSpace *sp) 515 { 516 PetscFunctionBegin; 517 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 518 PetscValidPointer(sp, 2); 519 *sp = fem->basisSpace; 520 PetscFunctionReturn(0); 521 } 522 523 /*@ 524 PetscFESetBasisSpace - Sets the PetscSpace used for approximation of the solution 525 526 Not collective 527 528 Input Parameters: 529 + fem - The PetscFE object 530 - sp - The PetscSpace object 531 532 Level: intermediate 533 534 .seealso: PetscFECreate() 535 @*/ 536 PetscErrorCode PetscFESetBasisSpace(PetscFE fem, PetscSpace sp) 537 { 538 PetscErrorCode ierr; 539 540 PetscFunctionBegin; 541 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 542 PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 2); 543 ierr = PetscSpaceDestroy(&fem->basisSpace);CHKERRQ(ierr); 544 fem->basisSpace = sp; 545 ierr = PetscObjectReference((PetscObject) fem->basisSpace);CHKERRQ(ierr); 546 PetscFunctionReturn(0); 547 } 548 549 /*@ 550 PetscFEGetDualSpace - Returns the PetscDualSpace used to define the inner product 551 552 Not collective 553 554 Input Parameter: 555 . fem - The PetscFE object 556 557 Output Parameter: 558 . sp - The PetscDualSpace object 559 560 Level: intermediate 561 562 .seealso: PetscFECreate() 563 @*/ 564 PetscErrorCode PetscFEGetDualSpace(PetscFE fem, PetscDualSpace *sp) 565 { 566 PetscFunctionBegin; 567 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 568 PetscValidPointer(sp, 2); 569 *sp = fem->dualSpace; 570 PetscFunctionReturn(0); 571 } 572 573 /*@ 574 PetscFESetDualSpace - Sets the PetscDualSpace used to define the inner product 575 576 Not collective 577 578 Input Parameters: 579 + fem - The PetscFE object 580 - sp - The PetscDualSpace object 581 582 Level: intermediate 583 584 .seealso: PetscFECreate() 585 @*/ 586 PetscErrorCode PetscFESetDualSpace(PetscFE fem, PetscDualSpace sp) 587 { 588 PetscErrorCode ierr; 589 590 PetscFunctionBegin; 591 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 592 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 2); 593 ierr = PetscDualSpaceDestroy(&fem->dualSpace);CHKERRQ(ierr); 594 fem->dualSpace = sp; 595 ierr = PetscObjectReference((PetscObject) fem->dualSpace);CHKERRQ(ierr); 596 PetscFunctionReturn(0); 597 } 598 599 /*@ 600 PetscFEGetQuadrature - Returns the PetscQuadrature used to calculate inner products 601 602 Not collective 603 604 Input Parameter: 605 . fem - The PetscFE object 606 607 Output Parameter: 608 . q - The PetscQuadrature object 609 610 Level: intermediate 611 612 .seealso: PetscFECreate() 613 @*/ 614 PetscErrorCode PetscFEGetQuadrature(PetscFE fem, PetscQuadrature *q) 615 { 616 PetscFunctionBegin; 617 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 618 PetscValidPointer(q, 2); 619 *q = fem->quadrature; 620 PetscFunctionReturn(0); 621 } 622 623 /*@ 624 PetscFESetQuadrature - Sets the PetscQuadrature used to calculate inner products 625 626 Not collective 627 628 Input Parameters: 629 + fem - The PetscFE object 630 - q - The PetscQuadrature object 631 632 Level: intermediate 633 634 .seealso: PetscFECreate() 635 @*/ 636 PetscErrorCode PetscFESetQuadrature(PetscFE fem, PetscQuadrature q) 637 { 638 PetscInt Nc, qNc; 639 PetscErrorCode ierr; 640 641 PetscFunctionBegin; 642 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 643 ierr = PetscFEGetNumComponents(fem, &Nc);CHKERRQ(ierr); 644 ierr = PetscQuadratureGetNumComponents(q, &qNc);CHKERRQ(ierr); 645 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); 646 ierr = PetscFERestoreTabulation(fem, 0, NULL, &fem->B, &fem->D, NULL /*&(*fem)->H*/);CHKERRQ(ierr); 647 ierr = PetscQuadratureDestroy(&fem->quadrature);CHKERRQ(ierr); 648 fem->quadrature = q; 649 ierr = PetscObjectReference((PetscObject) q);CHKERRQ(ierr); 650 PetscFunctionReturn(0); 651 } 652 653 /*@ 654 PetscFEGetFaceQuadrature - Returns the PetscQuadrature used to calculate inner products on faces 655 656 Not collective 657 658 Input Parameter: 659 . fem - The PetscFE object 660 661 Output Parameter: 662 . q - The PetscQuadrature object 663 664 Level: intermediate 665 666 .seealso: PetscFECreate() 667 @*/ 668 PetscErrorCode PetscFEGetFaceQuadrature(PetscFE fem, PetscQuadrature *q) 669 { 670 PetscFunctionBegin; 671 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 672 PetscValidPointer(q, 2); 673 *q = fem->faceQuadrature; 674 PetscFunctionReturn(0); 675 } 676 677 /*@ 678 PetscFESetFaceQuadrature - Sets the PetscQuadrature used to calculate inner products on faces 679 680 Not collective 681 682 Input Parameters: 683 + fem - The PetscFE object 684 - q - The PetscQuadrature object 685 686 Level: intermediate 687 688 .seealso: PetscFECreate() 689 @*/ 690 PetscErrorCode PetscFESetFaceQuadrature(PetscFE fem, PetscQuadrature q) 691 { 692 PetscErrorCode ierr; 693 694 PetscFunctionBegin; 695 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 696 ierr = PetscFERestoreTabulation(fem, 0, NULL, &fem->Bf, &fem->Df, NULL /*&(*fem)->Hf*/);CHKERRQ(ierr); 697 ierr = PetscQuadratureDestroy(&fem->faceQuadrature);CHKERRQ(ierr); 698 fem->faceQuadrature = q; 699 ierr = PetscObjectReference((PetscObject) q);CHKERRQ(ierr); 700 PetscFunctionReturn(0); 701 } 702 703 /*@C 704 PetscFEGetNumDof - Returns the number of dofs (dual basis vectors) associated to mesh points on the reference cell of a given dimension 705 706 Not collective 707 708 Input Parameter: 709 . fem - The PetscFE object 710 711 Output Parameter: 712 . numDof - Array with the number of dofs per dimension 713 714 Level: intermediate 715 716 .seealso: PetscFECreate() 717 @*/ 718 PetscErrorCode PetscFEGetNumDof(PetscFE fem, const PetscInt **numDof) 719 { 720 PetscErrorCode ierr; 721 722 PetscFunctionBegin; 723 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 724 PetscValidPointer(numDof, 2); 725 ierr = PetscDualSpaceGetNumDof(fem->dualSpace, numDof);CHKERRQ(ierr); 726 PetscFunctionReturn(0); 727 } 728 729 /*@C 730 PetscFEGetDefaultTabulation - Returns the tabulation of the basis functions at the quadrature points 731 732 Not collective 733 734 Input Parameter: 735 . fem - The PetscFE object 736 737 Output Parameters: 738 + B - The basis function values at quadrature points 739 . D - The basis function derivatives at quadrature points 740 - H - The basis function second derivatives at quadrature points 741 742 Note: 743 $ B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c 744 $ D[((p*pdim + i)*Nc + c)*dim + d] is the derivative value at point p for basis function i, component c, in direction d 745 $ 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 746 747 Level: intermediate 748 749 .seealso: PetscFEGetTabulation(), PetscFERestoreTabulation() 750 @*/ 751 PetscErrorCode PetscFEGetDefaultTabulation(PetscFE fem, PetscReal **B, PetscReal **D, PetscReal **H) 752 { 753 PetscInt npoints; 754 const PetscReal *points; 755 PetscErrorCode ierr; 756 757 PetscFunctionBegin; 758 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 759 if (B) PetscValidPointer(B, 2); 760 if (D) PetscValidPointer(D, 3); 761 if (H) PetscValidPointer(H, 4); 762 ierr = PetscQuadratureGetData(fem->quadrature, NULL, NULL, &npoints, &points, NULL);CHKERRQ(ierr); 763 if (!fem->B) {ierr = PetscFEGetTabulation(fem, npoints, points, &fem->B, &fem->D, NULL/*&fem->H*/);CHKERRQ(ierr);} 764 if (B) *B = fem->B; 765 if (D) *D = fem->D; 766 if (H) *H = fem->H; 767 PetscFunctionReturn(0); 768 } 769 770 PetscErrorCode PetscFEGetFaceTabulation(PetscFE fem, PetscReal **Bf, PetscReal **Df, PetscReal **Hf) 771 { 772 PetscErrorCode ierr; 773 774 PetscFunctionBegin; 775 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 776 if (Bf) PetscValidPointer(Bf, 2); 777 if (Df) PetscValidPointer(Df, 3); 778 if (Hf) PetscValidPointer(Hf, 4); 779 if (!fem->Bf) { 780 const PetscReal xi0[3] = {-1., -1., -1.}; 781 PetscReal v0[3], J[9], detJ; 782 PetscQuadrature fq; 783 PetscDualSpace sp; 784 DM dm; 785 const PetscInt *faces; 786 PetscInt dim, numFaces, f, npoints, q; 787 const PetscReal *points; 788 PetscReal *facePoints; 789 790 ierr = PetscFEGetDualSpace(fem, &sp);CHKERRQ(ierr); 791 ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr); 792 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 793 ierr = DMPlexGetConeSize(dm, 0, &numFaces);CHKERRQ(ierr); 794 ierr = DMPlexGetCone(dm, 0, &faces);CHKERRQ(ierr); 795 ierr = PetscFEGetFaceQuadrature(fem, &fq);CHKERRQ(ierr); 796 if (fq) { 797 ierr = PetscQuadratureGetData(fq, NULL, NULL, &npoints, &points, NULL);CHKERRQ(ierr); 798 ierr = PetscMalloc1(numFaces*npoints*dim, &facePoints);CHKERRQ(ierr); 799 for (f = 0; f < numFaces; ++f) { 800 ierr = DMPlexComputeCellGeometryFEM(dm, faces[f], NULL, v0, J, NULL, &detJ);CHKERRQ(ierr); 801 for (q = 0; q < npoints; ++q) CoordinatesRefToReal(dim, dim-1, xi0, v0, J, &points[q*(dim-1)], &facePoints[(f*npoints+q)*dim]); 802 } 803 ierr = PetscFEGetTabulation(fem, numFaces*npoints, facePoints, &fem->Bf, &fem->Df, NULL/*&fem->Hf*/);CHKERRQ(ierr); 804 ierr = PetscFree(facePoints);CHKERRQ(ierr); 805 } 806 } 807 if (Bf) *Bf = fem->Bf; 808 if (Df) *Df = fem->Df; 809 if (Hf) *Hf = fem->Hf; 810 PetscFunctionReturn(0); 811 } 812 813 PetscErrorCode PetscFEGetFaceCentroidTabulation(PetscFE fem, PetscReal **F) 814 { 815 PetscErrorCode ierr; 816 817 PetscFunctionBegin; 818 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 819 PetscValidPointer(F, 2); 820 if (!fem->F) { 821 PetscDualSpace sp; 822 DM dm; 823 const PetscInt *cone; 824 PetscReal *centroids; 825 PetscInt dim, numFaces, f; 826 827 ierr = PetscFEGetDualSpace(fem, &sp);CHKERRQ(ierr); 828 ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr); 829 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 830 ierr = DMPlexGetConeSize(dm, 0, &numFaces);CHKERRQ(ierr); 831 ierr = DMPlexGetCone(dm, 0, &cone);CHKERRQ(ierr); 832 ierr = PetscMalloc1(numFaces*dim, ¢roids);CHKERRQ(ierr); 833 for (f = 0; f < numFaces; ++f) {ierr = DMPlexComputeCellGeometryFVM(dm, cone[f], NULL, ¢roids[f*dim], NULL);CHKERRQ(ierr);} 834 ierr = PetscFEGetTabulation(fem, numFaces, centroids, &fem->F, NULL, NULL);CHKERRQ(ierr); 835 ierr = PetscFree(centroids);CHKERRQ(ierr); 836 } 837 *F = fem->F; 838 PetscFunctionReturn(0); 839 } 840 841 /*@C 842 PetscFEGetTabulation - Tabulates the basis functions, and perhaps derivatives, at the points provided. 843 844 Not collective 845 846 Input Parameters: 847 + fem - The PetscFE object 848 . npoints - The number of tabulation points 849 - points - The tabulation point coordinates 850 851 Output Parameters: 852 + B - The basis function values at tabulation points 853 . D - The basis function derivatives at tabulation points 854 - H - The basis function second derivatives at tabulation points 855 856 Note: 857 $ B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c 858 $ D[((p*pdim + i)*Nc + c)*dim + d] is the derivative value at point p for basis function i, component c, in direction d 859 $ 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 860 861 Level: intermediate 862 863 .seealso: PetscFERestoreTabulation(), PetscFEGetDefaultTabulation() 864 @*/ 865 PetscErrorCode PetscFEGetTabulation(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscReal **B, PetscReal **D, PetscReal **H) 866 { 867 DM dm; 868 PetscInt pdim; /* Dimension of FE space P */ 869 PetscInt dim; /* Spatial dimension */ 870 PetscInt comp; /* Field components */ 871 PetscErrorCode ierr; 872 873 PetscFunctionBegin; 874 if (!npoints) { 875 if (B) *B = NULL; 876 if (D) *D = NULL; 877 if (H) *H = NULL; 878 PetscFunctionReturn(0); 879 } 880 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 881 PetscValidPointer(points, 3); 882 if (B) PetscValidPointer(B, 4); 883 if (D) PetscValidPointer(D, 5); 884 if (H) PetscValidPointer(H, 6); 885 ierr = PetscDualSpaceGetDM(fem->dualSpace, &dm);CHKERRQ(ierr); 886 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 887 ierr = PetscDualSpaceGetDimension(fem->dualSpace, &pdim);CHKERRQ(ierr); 888 ierr = PetscFEGetNumComponents(fem, &comp);CHKERRQ(ierr); 889 if (B) {ierr = DMGetWorkArray(dm, npoints*pdim*comp, MPIU_REAL, B);CHKERRQ(ierr);} 890 if (!dim) { 891 if (D) *D = NULL; 892 if (H) *H = NULL; 893 } else { 894 if (D) {ierr = DMGetWorkArray(dm, npoints*pdim*comp*dim, MPIU_REAL, D);CHKERRQ(ierr);} 895 if (H) {ierr = DMGetWorkArray(dm, npoints*pdim*comp*dim*dim, MPIU_REAL, H);CHKERRQ(ierr);} 896 } 897 ierr = (*fem->ops->gettabulation)(fem, npoints, points, B ? *B : NULL, D ? *D : NULL, H ? *H : NULL);CHKERRQ(ierr); 898 PetscFunctionReturn(0); 899 } 900 901 PetscErrorCode PetscFERestoreTabulation(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscReal **B, PetscReal **D, PetscReal **H) 902 { 903 DM dm; 904 PetscErrorCode ierr; 905 906 PetscFunctionBegin; 907 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 908 ierr = PetscDualSpaceGetDM(fem->dualSpace, &dm);CHKERRQ(ierr); 909 if (B && *B) {ierr = DMRestoreWorkArray(dm, 0, MPIU_REAL, B);CHKERRQ(ierr);} 910 if (D && *D) {ierr = DMRestoreWorkArray(dm, 0, MPIU_REAL, D);CHKERRQ(ierr);} 911 if (H && *H) {ierr = DMRestoreWorkArray(dm, 0, MPIU_REAL, H);CHKERRQ(ierr);} 912 PetscFunctionReturn(0); 913 } 914 915 PETSC_EXTERN PetscErrorCode PetscFECreatePointTrace(PetscFE fe, PetscInt refPoint, PetscFE *trFE) 916 { 917 PetscSpace bsp, bsubsp; 918 PetscDualSpace dsp, dsubsp; 919 PetscInt dim, depth, numComp, i, j, coneSize, order; 920 PetscFEType type; 921 DM dm; 922 DMLabel label; 923 PetscReal *xi, *v, *J, detJ; 924 const char *name; 925 PetscQuadrature origin, fullQuad, subQuad; 926 PetscErrorCode ierr; 927 928 PetscFunctionBegin; 929 PetscValidHeaderSpecific(fe,PETSCFE_CLASSID,1); 930 PetscValidPointer(trFE,3); 931 ierr = PetscFEGetBasisSpace(fe,&bsp);CHKERRQ(ierr); 932 ierr = PetscFEGetDualSpace(fe,&dsp);CHKERRQ(ierr); 933 ierr = PetscDualSpaceGetDM(dsp,&dm);CHKERRQ(ierr); 934 ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 935 ierr = DMPlexGetDepthLabel(dm,&label);CHKERRQ(ierr); 936 ierr = DMLabelGetValue(label,refPoint,&depth);CHKERRQ(ierr); 937 ierr = PetscCalloc1(depth,&xi);CHKERRQ(ierr); 938 ierr = PetscMalloc1(dim,&v);CHKERRQ(ierr); 939 ierr = PetscMalloc1(dim*dim,&J);CHKERRQ(ierr); 940 for (i = 0; i < depth; i++) xi[i] = 0.; 941 ierr = PetscQuadratureCreate(PETSC_COMM_SELF,&origin);CHKERRQ(ierr); 942 ierr = PetscQuadratureSetData(origin,depth,0,1,xi,NULL);CHKERRQ(ierr); 943 ierr = DMPlexComputeCellGeometryFEM(dm,refPoint,origin,v,J,NULL,&detJ);CHKERRQ(ierr); 944 /* CellGeometryFEM computes the expanded Jacobian, we want the true jacobian */ 945 for (i = 1; i < dim; i++) { 946 for (j = 0; j < depth; j++) { 947 J[i * depth + j] = J[i * dim + j]; 948 } 949 } 950 ierr = PetscQuadratureDestroy(&origin);CHKERRQ(ierr); 951 ierr = PetscDualSpaceGetPointSubspace(dsp,refPoint,&dsubsp);CHKERRQ(ierr); 952 ierr = PetscSpaceCreateSubspace(bsp,dsubsp,v,J,NULL,NULL,PETSC_OWN_POINTER,&bsubsp);CHKERRQ(ierr); 953 ierr = PetscSpaceSetUp(bsubsp);CHKERRQ(ierr); 954 ierr = PetscFECreate(PetscObjectComm((PetscObject)fe),trFE);CHKERRQ(ierr); 955 ierr = PetscFEGetType(fe,&type);CHKERRQ(ierr); 956 ierr = PetscFESetType(*trFE,type);CHKERRQ(ierr); 957 ierr = PetscFEGetNumComponents(fe,&numComp);CHKERRQ(ierr); 958 ierr = PetscFESetNumComponents(*trFE,numComp);CHKERRQ(ierr); 959 ierr = PetscFESetBasisSpace(*trFE,bsubsp);CHKERRQ(ierr); 960 ierr = PetscFESetDualSpace(*trFE,dsubsp);CHKERRQ(ierr); 961 ierr = PetscObjectGetName((PetscObject) fe, &name);CHKERRQ(ierr); 962 if (name) {ierr = PetscFESetName(*trFE, name);CHKERRQ(ierr);} 963 ierr = PetscFEGetQuadrature(fe,&fullQuad);CHKERRQ(ierr); 964 ierr = PetscQuadratureGetOrder(fullQuad,&order);CHKERRQ(ierr); 965 ierr = DMPlexGetConeSize(dm,refPoint,&coneSize);CHKERRQ(ierr); 966 if (coneSize == 2 * depth) { 967 ierr = PetscDTGaussTensorQuadrature(depth,1,(order + 1)/2,-1.,1.,&subQuad);CHKERRQ(ierr); 968 } else { 969 ierr = PetscDTGaussJacobiQuadrature(depth,1,(order + 1)/2,-1.,1.,&subQuad);CHKERRQ(ierr); 970 } 971 ierr = PetscFESetQuadrature(*trFE,subQuad);CHKERRQ(ierr); 972 ierr = PetscFESetUp(*trFE);CHKERRQ(ierr); 973 ierr = PetscQuadratureDestroy(&subQuad);CHKERRQ(ierr); 974 ierr = PetscSpaceDestroy(&bsubsp);CHKERRQ(ierr); 975 PetscFunctionReturn(0); 976 } 977 978 PetscErrorCode PetscFECreateHeightTrace(PetscFE fe, PetscInt height, PetscFE *trFE) 979 { 980 PetscInt hStart, hEnd; 981 PetscDualSpace dsp; 982 DM dm; 983 PetscErrorCode ierr; 984 985 PetscFunctionBegin; 986 PetscValidHeaderSpecific(fe,PETSCFE_CLASSID,1); 987 PetscValidPointer(trFE,3); 988 *trFE = NULL; 989 ierr = PetscFEGetDualSpace(fe,&dsp);CHKERRQ(ierr); 990 ierr = PetscDualSpaceGetDM(dsp,&dm);CHKERRQ(ierr); 991 ierr = DMPlexGetHeightStratum(dm,height,&hStart,&hEnd);CHKERRQ(ierr); 992 if (hEnd <= hStart) PetscFunctionReturn(0); 993 ierr = PetscFECreatePointTrace(fe,hStart,trFE);CHKERRQ(ierr); 994 PetscFunctionReturn(0); 995 } 996 997 998 /*@ 999 PetscFEGetDimension - Get the dimension of the finite element space on a cell 1000 1001 Not collective 1002 1003 Input Parameter: 1004 . fe - The PetscFE 1005 1006 Output Parameter: 1007 . dim - The dimension 1008 1009 Level: intermediate 1010 1011 .seealso: PetscFECreate(), PetscSpaceGetDimension(), PetscDualSpaceGetDimension() 1012 @*/ 1013 PetscErrorCode PetscFEGetDimension(PetscFE fem, PetscInt *dim) 1014 { 1015 PetscErrorCode ierr; 1016 1017 PetscFunctionBegin; 1018 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 1019 PetscValidPointer(dim, 2); 1020 if (fem->ops->getdimension) {ierr = (*fem->ops->getdimension)(fem, dim);CHKERRQ(ierr);} 1021 PetscFunctionReturn(0); 1022 } 1023 1024 /* 1025 Purpose: Compute element vector for chunk of elements 1026 1027 Input: 1028 Sizes: 1029 Ne: number of elements 1030 Nf: number of fields 1031 PetscFE 1032 dim: spatial dimension 1033 Nb: number of basis functions 1034 Nc: number of field components 1035 PetscQuadrature 1036 Nq: number of quadrature points 1037 1038 Geometry: 1039 PetscFEGeom[Ne] possibly *Nq 1040 PetscReal v0s[dim] 1041 PetscReal n[dim] 1042 PetscReal jacobians[dim*dim] 1043 PetscReal jacobianInverses[dim*dim] 1044 PetscReal jacobianDeterminants 1045 FEM: 1046 PetscFE 1047 PetscQuadrature 1048 PetscReal quadPoints[Nq*dim] 1049 PetscReal quadWeights[Nq] 1050 PetscReal basis[Nq*Nb*Nc] 1051 PetscReal basisDer[Nq*Nb*Nc*dim] 1052 PetscScalar coefficients[Ne*Nb*Nc] 1053 PetscScalar elemVec[Ne*Nb*Nc] 1054 1055 Problem: 1056 PetscInt f: the active field 1057 f0, f1 1058 1059 Work Space: 1060 PetscFE 1061 PetscScalar f0[Nq*dim]; 1062 PetscScalar f1[Nq*dim*dim]; 1063 PetscScalar u[Nc]; 1064 PetscScalar gradU[Nc*dim]; 1065 PetscReal x[dim]; 1066 PetscScalar realSpaceDer[dim]; 1067 1068 Purpose: Compute element vector for N_cb batches of elements 1069 1070 Input: 1071 Sizes: 1072 N_cb: Number of serial cell batches 1073 1074 Geometry: 1075 PetscReal v0s[Ne*dim] 1076 PetscReal jacobians[Ne*dim*dim] possibly *Nq 1077 PetscReal jacobianInverses[Ne*dim*dim] possibly *Nq 1078 PetscReal jacobianDeterminants[Ne] possibly *Nq 1079 FEM: 1080 static PetscReal quadPoints[Nq*dim] 1081 static PetscReal quadWeights[Nq] 1082 static PetscReal basis[Nq*Nb*Nc] 1083 static PetscReal basisDer[Nq*Nb*Nc*dim] 1084 PetscScalar coefficients[Ne*Nb*Nc] 1085 PetscScalar elemVec[Ne*Nb*Nc] 1086 1087 ex62.c: 1088 PetscErrorCode PetscFEIntegrateResidualBatch(PetscInt Ne, PetscInt numFields, PetscInt field, PetscQuadrature quad[], const PetscScalar coefficients[], 1089 const PetscReal v0s[], const PetscReal jacobians[], const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], 1090 void (*f0_func)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar f0[]), 1091 void (*f1_func)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar f1[]), PetscScalar elemVec[]) 1092 1093 ex52.c: 1094 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) 1095 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) 1096 1097 ex52_integrateElement.cu 1098 __global__ void integrateElementQuadrature(int N_cb, realType *coefficients, realType *jacobianInverses, realType *jacobianDeterminants, realType *elemVec) 1099 1100 PETSC_EXTERN PetscErrorCode IntegrateElementBatchGPU(PetscInt spatial_dim, PetscInt Ne, PetscInt Ncb, PetscInt Nbc, PetscInt Nbl, const PetscScalar coefficients[], 1101 const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], PetscScalar elemVec[], 1102 PetscLogEvent event, PetscInt debug, PetscInt pde_op) 1103 1104 ex52_integrateElementOpenCL.c: 1105 PETSC_EXTERN PetscErrorCode IntegrateElementBatchGPU(PetscInt spatial_dim, PetscInt Ne, PetscInt Ncb, PetscInt Nbc, PetscInt N_bl, const PetscScalar coefficients[], 1106 const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], PetscScalar elemVec[], 1107 PetscLogEvent event, PetscInt debug, PetscInt pde_op) 1108 1109 __kernel void integrateElementQuadrature(int N_cb, __global float *coefficients, __global float *jacobianInverses, __global float *jacobianDeterminants, __global float *elemVec) 1110 */ 1111 1112 /*@C 1113 PetscFEIntegrate - Produce the integral for the given field for a chunk of elements by quadrature integration 1114 1115 Not collective 1116 1117 Input Parameters: 1118 + fem - The PetscFE object for the field being integrated 1119 . prob - The PetscDS specifying the discretizations and continuum functions 1120 . field - The field being integrated 1121 . Ne - The number of elements in the chunk 1122 . cgeom - The cell geometry for each cell in the chunk 1123 . coefficients - The array of FEM basis coefficients for the elements 1124 . probAux - The PetscDS specifying the auxiliary discretizations 1125 - coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 1126 1127 Output Parameter 1128 . integral - the integral for this field 1129 1130 Level: developer 1131 1132 .seealso: PetscFEIntegrateResidual() 1133 @*/ 1134 PetscErrorCode PetscFEIntegrate(PetscFE fem, PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom, 1135 const PetscScalar coefficients[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscScalar integral[]) 1136 { 1137 PetscErrorCode ierr; 1138 1139 PetscFunctionBegin; 1140 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 1141 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 2); 1142 if (fem->ops->integrate) {ierr = (*fem->ops->integrate)(fem, prob, field, Ne, cgeom, coefficients, probAux, coefficientsAux, integral);CHKERRQ(ierr);} 1143 PetscFunctionReturn(0); 1144 } 1145 1146 /*@C 1147 PetscFEIntegrateBd - Produce the integral for the given field for a chunk of elements by quadrature integration 1148 1149 Not collective 1150 1151 Input Parameters: 1152 + fem - The PetscFE object for the field being integrated 1153 . prob - The PetscDS specifying the discretizations and continuum functions 1154 . field - The field being integrated 1155 . obj_func - The function to be integrated 1156 . Ne - The number of elements in the chunk 1157 . fgeom - The face geometry for each face in the chunk 1158 . coefficients - The array of FEM basis coefficients for the elements 1159 . probAux - The PetscDS specifying the auxiliary discretizations 1160 - coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 1161 1162 Output Parameter 1163 . integral - the integral for this field 1164 1165 Level: developer 1166 1167 .seealso: PetscFEIntegrateResidual() 1168 @*/ 1169 PetscErrorCode PetscFEIntegrateBd(PetscFE fem, PetscDS prob, PetscInt field, 1170 void (*obj_func)(PetscInt, PetscInt, PetscInt, 1171 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 1172 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 1173 PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 1174 PetscInt Ne, PetscFEGeom *geom, const PetscScalar coefficients[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscScalar integral[]) 1175 { 1176 PetscErrorCode ierr; 1177 1178 PetscFunctionBegin; 1179 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 1180 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 2); 1181 if (fem->ops->integratebd) {ierr = (*fem->ops->integratebd)(fem, prob, field, obj_func, Ne, geom, coefficients, probAux, coefficientsAux, integral);CHKERRQ(ierr);} 1182 PetscFunctionReturn(0); 1183 } 1184 1185 /*@C 1186 PetscFEIntegrateResidual - Produce the element residual vector for a chunk of elements by quadrature integration 1187 1188 Not collective 1189 1190 Input Parameters: 1191 + fem - The PetscFE object for the field being integrated 1192 . prob - The PetscDS specifying the discretizations and continuum functions 1193 . field - The field being integrated 1194 . Ne - The number of elements in the chunk 1195 . cgeom - The cell geometry for each cell in the chunk 1196 . coefficients - The array of FEM basis coefficients for the elements 1197 . coefficients_t - The array of FEM basis time derivative coefficients for the elements 1198 . probAux - The PetscDS specifying the auxiliary discretizations 1199 . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 1200 - t - The time 1201 1202 Output Parameter 1203 . elemVec - the element residual vectors from each element 1204 1205 Note: 1206 $ Loop over batch of elements (e): 1207 $ Loop over quadrature points (q): 1208 $ Make u_q and gradU_q (loops over fields,Nb,Ncomp) and x_q 1209 $ Call f_0 and f_1 1210 $ Loop over element vector entries (f,fc --> i): 1211 $ elemVec[i] += \psi^{fc}_f(q) f0_{fc}(u, \nabla u) + \nabla\psi^{fc}_f(q) \cdot f1_{fc,df}(u, \nabla u) 1212 1213 Level: developer 1214 1215 .seealso: PetscFEIntegrateResidual() 1216 @*/ 1217 PetscErrorCode PetscFEIntegrateResidual(PetscFE fem, PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom, 1218 const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[]) 1219 { 1220 PetscErrorCode ierr; 1221 1222 PetscFunctionBegin; 1223 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 1224 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 2); 1225 if (fem->ops->integrateresidual) {ierr = (*fem->ops->integrateresidual)(fem, prob, field, Ne, cgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec);CHKERRQ(ierr);} 1226 PetscFunctionReturn(0); 1227 } 1228 1229 /*@C 1230 PetscFEIntegrateBdResidual - Produce the element residual vector for a chunk of elements by quadrature integration over a boundary 1231 1232 Not collective 1233 1234 Input Parameters: 1235 + fem - The PetscFE object for the field being integrated 1236 . prob - The PetscDS specifying the discretizations and continuum functions 1237 . field - The field being integrated 1238 . Ne - The number of elements in the chunk 1239 . fgeom - The face geometry for each cell in the chunk 1240 . coefficients - The array of FEM basis coefficients for the elements 1241 . coefficients_t - The array of FEM basis time derivative coefficients for the elements 1242 . probAux - The PetscDS specifying the auxiliary discretizations 1243 . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 1244 - t - The time 1245 1246 Output Parameter 1247 . elemVec - the element residual vectors from each element 1248 1249 Level: developer 1250 1251 .seealso: PetscFEIntegrateResidual() 1252 @*/ 1253 PetscErrorCode PetscFEIntegrateBdResidual(PetscFE fem, PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *fgeom, 1254 const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[]) 1255 { 1256 PetscErrorCode ierr; 1257 1258 PetscFunctionBegin; 1259 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 1260 if (fem->ops->integratebdresidual) {ierr = (*fem->ops->integratebdresidual)(fem, prob, field, Ne, fgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec);CHKERRQ(ierr);} 1261 PetscFunctionReturn(0); 1262 } 1263 1264 /*@C 1265 PetscFEIntegrateJacobian - Produce the element Jacobian for a chunk of elements by quadrature integration 1266 1267 Not collective 1268 1269 Input Parameters: 1270 + fem - The PetscFE object for the field being integrated 1271 . prob - The PetscDS specifying the discretizations and continuum functions 1272 . jtype - The type of matrix pointwise functions that should be used 1273 . fieldI - The test field being integrated 1274 . fieldJ - The basis field being integrated 1275 . Ne - The number of elements in the chunk 1276 . cgeom - The cell geometry for each cell in the chunk 1277 . coefficients - The array of FEM basis coefficients for the elements for the Jacobian evaluation point 1278 . coefficients_t - The array of FEM basis time derivative coefficients for the elements 1279 . probAux - The PetscDS specifying the auxiliary discretizations 1280 . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 1281 . t - The time 1282 - u_tShift - A multiplier for the dF/du_t term (as opposed to the dF/du term) 1283 1284 Output Parameter 1285 . elemMat - the element matrices for the Jacobian from each element 1286 1287 Note: 1288 $ Loop over batch of elements (e): 1289 $ Loop over element matrix entries (f,fc,g,gc --> i,j): 1290 $ Loop over quadrature points (q): 1291 $ Make u_q and gradU_q (loops over fields,Nb,Ncomp) 1292 $ elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q) 1293 $ + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q) 1294 $ + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q) 1295 $ + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q) 1296 Level: developer 1297 1298 .seealso: PetscFEIntegrateResidual() 1299 @*/ 1300 PetscErrorCode PetscFEIntegrateJacobian(PetscFE fem, PetscDS prob, PetscFEJacobianType jtype, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *cgeom, 1301 const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) 1302 { 1303 PetscErrorCode ierr; 1304 1305 PetscFunctionBegin; 1306 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 1307 if (fem->ops->integratejacobian) {ierr = (*fem->ops->integratejacobian)(fem, prob, jtype, fieldI, fieldJ, Ne, cgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, u_tshift, elemMat);CHKERRQ(ierr);} 1308 PetscFunctionReturn(0); 1309 } 1310 1311 /*@C 1312 PetscFEIntegrateBdJacobian - Produce the boundary element Jacobian for a chunk of elements by quadrature integration 1313 1314 Not collective 1315 1316 Input Parameters: 1317 + fem = The PetscFE object for the field being integrated 1318 . prob - The PetscDS specifying the discretizations and continuum functions 1319 . fieldI - The test field being integrated 1320 . fieldJ - The basis field being integrated 1321 . Ne - The number of elements in the chunk 1322 . fgeom - The face geometry for each cell in the chunk 1323 . coefficients - The array of FEM basis coefficients for the elements for the Jacobian evaluation point 1324 . coefficients_t - The array of FEM basis time derivative coefficients for the elements 1325 . probAux - The PetscDS specifying the auxiliary discretizations 1326 . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 1327 . t - The time 1328 - u_tShift - A multiplier for the dF/du_t term (as opposed to the dF/du term) 1329 1330 Output Parameter 1331 . elemMat - the element matrices for the Jacobian from each element 1332 1333 Note: 1334 $ Loop over batch of elements (e): 1335 $ Loop over element matrix entries (f,fc,g,gc --> i,j): 1336 $ Loop over quadrature points (q): 1337 $ Make u_q and gradU_q (loops over fields,Nb,Ncomp) 1338 $ elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q) 1339 $ + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q) 1340 $ + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q) 1341 $ + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q) 1342 Level: developer 1343 1344 .seealso: PetscFEIntegrateJacobian(), PetscFEIntegrateResidual() 1345 @*/ 1346 PetscErrorCode PetscFEIntegrateBdJacobian(PetscFE fem, PetscDS prob, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *fgeom, 1347 const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) 1348 { 1349 PetscErrorCode ierr; 1350 1351 PetscFunctionBegin; 1352 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 1353 if (fem->ops->integratebdjacobian) {ierr = (*fem->ops->integratebdjacobian)(fem, prob, fieldI, fieldJ, Ne, fgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, u_tshift, elemMat);CHKERRQ(ierr);} 1354 PetscFunctionReturn(0); 1355 } 1356 1357 PetscErrorCode PetscFEGetHeightSubspace(PetscFE fe, PetscInt height, PetscFE *subfe) 1358 { 1359 PetscSpace P, subP; 1360 PetscDualSpace Q, subQ; 1361 PetscQuadrature subq; 1362 PetscFEType fetype; 1363 PetscInt dim, Nc; 1364 PetscErrorCode ierr; 1365 1366 PetscFunctionBegin; 1367 PetscValidHeaderSpecific(fe, PETSCFE_CLASSID, 1); 1368 PetscValidPointer(subfe, 3); 1369 if (height == 0) { 1370 *subfe = fe; 1371 PetscFunctionReturn(0); 1372 } 1373 ierr = PetscFEGetBasisSpace(fe, &P);CHKERRQ(ierr); 1374 ierr = PetscFEGetDualSpace(fe, &Q);CHKERRQ(ierr); 1375 ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 1376 ierr = PetscFEGetFaceQuadrature(fe, &subq);CHKERRQ(ierr); 1377 ierr = PetscDualSpaceGetDimension(Q, &dim);CHKERRQ(ierr); 1378 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);} 1379 if (!fe->subspaces) {ierr = PetscCalloc1(dim, &fe->subspaces);CHKERRQ(ierr);} 1380 if (height <= dim) { 1381 if (!fe->subspaces[height-1]) { 1382 PetscFE sub; 1383 const char *name; 1384 1385 ierr = PetscSpaceGetHeightSubspace(P, height, &subP);CHKERRQ(ierr); 1386 ierr = PetscDualSpaceGetHeightSubspace(Q, height, &subQ);CHKERRQ(ierr); 1387 ierr = PetscFECreate(PetscObjectComm((PetscObject) fe), &sub);CHKERRQ(ierr); 1388 ierr = PetscObjectGetName((PetscObject) fe, &name);CHKERRQ(ierr); 1389 ierr = PetscObjectSetName((PetscObject) sub, name);CHKERRQ(ierr); 1390 ierr = PetscFEGetType(fe, &fetype);CHKERRQ(ierr); 1391 ierr = PetscFESetType(sub, fetype);CHKERRQ(ierr); 1392 ierr = PetscFESetBasisSpace(sub, subP);CHKERRQ(ierr); 1393 ierr = PetscFESetDualSpace(sub, subQ);CHKERRQ(ierr); 1394 ierr = PetscFESetNumComponents(sub, Nc);CHKERRQ(ierr); 1395 ierr = PetscFESetUp(sub);CHKERRQ(ierr); 1396 ierr = PetscFESetQuadrature(sub, subq);CHKERRQ(ierr); 1397 fe->subspaces[height-1] = sub; 1398 } 1399 *subfe = fe->subspaces[height-1]; 1400 } else { 1401 *subfe = NULL; 1402 } 1403 PetscFunctionReturn(0); 1404 } 1405 1406 /*@ 1407 PetscFERefine - Create a "refined" PetscFE object that refines the reference cell into smaller copies. This is typically used 1408 to precondition a higher order method with a lower order method on a refined mesh having the same number of dofs (but more 1409 sparsity). It is also used to create an interpolation between regularly refined meshes. 1410 1411 Collective on PetscFE 1412 1413 Input Parameter: 1414 . fe - The initial PetscFE 1415 1416 Output Parameter: 1417 . feRef - The refined PetscFE 1418 1419 Level: developer 1420 1421 .seealso: PetscFEType, PetscFECreate(), PetscFESetType() 1422 @*/ 1423 PetscErrorCode PetscFERefine(PetscFE fe, PetscFE *feRef) 1424 { 1425 PetscSpace P, Pref; 1426 PetscDualSpace Q, Qref; 1427 DM K, Kref; 1428 PetscQuadrature q, qref; 1429 const PetscReal *v0, *jac; 1430 PetscInt numComp, numSubelements; 1431 PetscErrorCode ierr; 1432 1433 PetscFunctionBegin; 1434 ierr = PetscFEGetBasisSpace(fe, &P);CHKERRQ(ierr); 1435 ierr = PetscFEGetDualSpace(fe, &Q);CHKERRQ(ierr); 1436 ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr); 1437 ierr = PetscDualSpaceGetDM(Q, &K);CHKERRQ(ierr); 1438 /* Create space */ 1439 ierr = PetscObjectReference((PetscObject) P);CHKERRQ(ierr); 1440 Pref = P; 1441 /* Create dual space */ 1442 ierr = PetscDualSpaceDuplicate(Q, &Qref);CHKERRQ(ierr); 1443 ierr = DMRefine(K, PetscObjectComm((PetscObject) fe), &Kref);CHKERRQ(ierr); 1444 ierr = PetscDualSpaceSetDM(Qref, Kref);CHKERRQ(ierr); 1445 ierr = DMDestroy(&Kref);CHKERRQ(ierr); 1446 ierr = PetscDualSpaceSetUp(Qref);CHKERRQ(ierr); 1447 /* Create element */ 1448 ierr = PetscFECreate(PetscObjectComm((PetscObject) fe), feRef);CHKERRQ(ierr); 1449 ierr = PetscFESetType(*feRef, PETSCFECOMPOSITE);CHKERRQ(ierr); 1450 ierr = PetscFESetBasisSpace(*feRef, Pref);CHKERRQ(ierr); 1451 ierr = PetscFESetDualSpace(*feRef, Qref);CHKERRQ(ierr); 1452 ierr = PetscFEGetNumComponents(fe, &numComp);CHKERRQ(ierr); 1453 ierr = PetscFESetNumComponents(*feRef, numComp);CHKERRQ(ierr); 1454 ierr = PetscFESetUp(*feRef);CHKERRQ(ierr); 1455 ierr = PetscSpaceDestroy(&Pref);CHKERRQ(ierr); 1456 ierr = PetscDualSpaceDestroy(&Qref);CHKERRQ(ierr); 1457 /* Create quadrature */ 1458 ierr = PetscFECompositeGetMapping(*feRef, &numSubelements, &v0, &jac, NULL);CHKERRQ(ierr); 1459 ierr = PetscQuadratureExpandComposite(q, numSubelements, v0, jac, &qref);CHKERRQ(ierr); 1460 ierr = PetscFESetQuadrature(*feRef, qref);CHKERRQ(ierr); 1461 ierr = PetscQuadratureDestroy(&qref);CHKERRQ(ierr); 1462 PetscFunctionReturn(0); 1463 } 1464 1465 /*@C 1466 PetscFECreateDefault - Create a PetscFE for basic FEM computation 1467 1468 Collective on DM 1469 1470 Input Parameters: 1471 + comm - The MPI comm 1472 . dim - The spatial dimension 1473 . Nc - The number of components 1474 . isSimplex - Flag for simplex reference cell, otherwise its a tensor product 1475 . prefix - The options prefix, or NULL 1476 - qorder - The quadrature order 1477 1478 Output Parameter: 1479 . fem - The PetscFE object 1480 1481 Level: beginner 1482 1483 .keywords: PetscFE, finite element 1484 .seealso: PetscFECreate(), PetscSpaceCreate(), PetscDualSpaceCreate() 1485 @*/ 1486 PetscErrorCode PetscFECreateDefault(MPI_Comm comm, PetscInt dim, PetscInt Nc, PetscBool isSimplex, const char prefix[], PetscInt qorder, PetscFE *fem) 1487 { 1488 PetscQuadrature q, fq; 1489 DM K; 1490 PetscSpace P; 1491 PetscDualSpace Q; 1492 PetscInt order, quadPointsPerEdge; 1493 PetscBool tensor = isSimplex ? PETSC_FALSE : PETSC_TRUE; 1494 PetscErrorCode ierr; 1495 1496 PetscFunctionBegin; 1497 /* Create space */ 1498 ierr = PetscSpaceCreate(comm, &P);CHKERRQ(ierr); 1499 ierr = PetscObjectSetOptionsPrefix((PetscObject) P, prefix);CHKERRQ(ierr); 1500 ierr = PetscSpacePolynomialSetTensor(P, tensor);CHKERRQ(ierr); 1501 ierr = PetscSpaceSetNumComponents(P, Nc);CHKERRQ(ierr); 1502 ierr = PetscSpaceSetNumVariables(P, dim);CHKERRQ(ierr); 1503 ierr = PetscSpaceSetFromOptions(P);CHKERRQ(ierr); 1504 ierr = PetscSpaceSetUp(P);CHKERRQ(ierr); 1505 ierr = PetscSpaceGetDegree(P, &order, NULL);CHKERRQ(ierr); 1506 ierr = PetscSpacePolynomialGetTensor(P, &tensor);CHKERRQ(ierr); 1507 /* Create dual space */ 1508 ierr = PetscDualSpaceCreate(comm, &Q);CHKERRQ(ierr); 1509 ierr = PetscDualSpaceSetType(Q,PETSCDUALSPACELAGRANGE);CHKERRQ(ierr); 1510 ierr = PetscObjectSetOptionsPrefix((PetscObject) Q, prefix);CHKERRQ(ierr); 1511 ierr = PetscDualSpaceCreateReferenceCell(Q, dim, isSimplex, &K);CHKERRQ(ierr); 1512 ierr = PetscDualSpaceSetDM(Q, K);CHKERRQ(ierr); 1513 ierr = DMDestroy(&K);CHKERRQ(ierr); 1514 ierr = PetscDualSpaceSetNumComponents(Q, Nc);CHKERRQ(ierr); 1515 ierr = PetscDualSpaceSetOrder(Q, order);CHKERRQ(ierr); 1516 ierr = PetscDualSpaceLagrangeSetTensor(Q, tensor);CHKERRQ(ierr); 1517 ierr = PetscDualSpaceSetFromOptions(Q);CHKERRQ(ierr); 1518 ierr = PetscDualSpaceSetUp(Q);CHKERRQ(ierr); 1519 /* Create element */ 1520 ierr = PetscFECreate(comm, fem);CHKERRQ(ierr); 1521 ierr = PetscObjectSetOptionsPrefix((PetscObject) *fem, prefix);CHKERRQ(ierr); 1522 ierr = PetscFESetFromOptions(*fem);CHKERRQ(ierr); 1523 ierr = PetscFESetBasisSpace(*fem, P);CHKERRQ(ierr); 1524 ierr = PetscFESetDualSpace(*fem, Q);CHKERRQ(ierr); 1525 ierr = PetscFESetNumComponents(*fem, Nc);CHKERRQ(ierr); 1526 ierr = PetscFESetUp(*fem);CHKERRQ(ierr); 1527 ierr = PetscSpaceDestroy(&P);CHKERRQ(ierr); 1528 ierr = PetscDualSpaceDestroy(&Q);CHKERRQ(ierr); 1529 /* Create quadrature (with specified order if given) */ 1530 qorder = qorder >= 0 ? qorder : order; 1531 ierr = PetscObjectOptionsBegin((PetscObject)*fem);CHKERRQ(ierr); 1532 ierr = PetscOptionsInt("-petscfe_default_quadrature_order","Quadrature order is one less than quadture points per edge","PetscFECreateDefault",qorder,&qorder,NULL);CHKERRQ(ierr); 1533 ierr = PetscOptionsEnd();CHKERRQ(ierr); 1534 quadPointsPerEdge = PetscMax(qorder + 1,1); 1535 if (isSimplex) { 1536 ierr = PetscDTGaussJacobiQuadrature(dim, 1, quadPointsPerEdge, -1.0, 1.0, &q);CHKERRQ(ierr); 1537 ierr = PetscDTGaussJacobiQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0, &fq);CHKERRQ(ierr); 1538 } else { 1539 ierr = PetscDTGaussTensorQuadrature(dim, 1, quadPointsPerEdge, -1.0, 1.0, &q);CHKERRQ(ierr); 1540 ierr = PetscDTGaussTensorQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0, &fq);CHKERRQ(ierr); 1541 } 1542 ierr = PetscFESetQuadrature(*fem, q);CHKERRQ(ierr); 1543 ierr = PetscFESetFaceQuadrature(*fem, fq);CHKERRQ(ierr); 1544 ierr = PetscQuadratureDestroy(&q);CHKERRQ(ierr); 1545 ierr = PetscQuadratureDestroy(&fq);CHKERRQ(ierr); 1546 PetscFunctionReturn(0); 1547 } 1548 1549 /*@C 1550 PetscFESetName - Names the FE and its subobjects 1551 1552 Not collective 1553 1554 Input Parameters: 1555 + fe - The PetscFE 1556 - name - The name 1557 1558 Level: beginner 1559 1560 .keywords: PetscFE, finite element 1561 .seealso: PetscFECreate(), PetscSpaceCreate(), PetscDualSpaceCreate() 1562 @*/ 1563 PetscErrorCode PetscFESetName(PetscFE fe, const char name[]) 1564 { 1565 PetscSpace P; 1566 PetscDualSpace Q; 1567 PetscErrorCode ierr; 1568 1569 PetscFunctionBegin; 1570 ierr = PetscFEGetBasisSpace(fe, &P);CHKERRQ(ierr); 1571 ierr = PetscFEGetDualSpace(fe, &Q);CHKERRQ(ierr); 1572 ierr = PetscObjectSetName((PetscObject) fe, name);CHKERRQ(ierr); 1573 ierr = PetscObjectSetName((PetscObject) P, name);CHKERRQ(ierr); 1574 ierr = PetscObjectSetName((PetscObject) Q, name);CHKERRQ(ierr); 1575 PetscFunctionReturn(0); 1576 } 1577