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 PetscQuadrature origin, fullQuad, subQuad; 925 PetscErrorCode ierr; 926 927 PetscFunctionBegin; 928 PetscValidHeaderSpecific(fe,PETSCFE_CLASSID,1); 929 PetscValidPointer(trFE,3); 930 ierr = PetscFEGetBasisSpace(fe,&bsp);CHKERRQ(ierr); 931 ierr = PetscFEGetDualSpace(fe,&dsp);CHKERRQ(ierr); 932 ierr = PetscDualSpaceGetDM(dsp,&dm);CHKERRQ(ierr); 933 ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 934 ierr = DMPlexGetDepthLabel(dm,&label);CHKERRQ(ierr); 935 ierr = DMLabelGetValue(label,refPoint,&depth);CHKERRQ(ierr); 936 ierr = PetscCalloc1(depth,&xi);CHKERRQ(ierr); 937 ierr = PetscMalloc1(dim,&v);CHKERRQ(ierr); 938 ierr = PetscMalloc1(dim*dim,&J);CHKERRQ(ierr); 939 for (i = 0; i < depth; i++) xi[i] = 0.; 940 ierr = PetscQuadratureCreate(PETSC_COMM_SELF,&origin);CHKERRQ(ierr); 941 ierr = PetscQuadratureSetData(origin,depth,0,1,xi,NULL);CHKERRQ(ierr); 942 ierr = DMPlexComputeCellGeometryFEM(dm,refPoint,origin,v,J,NULL,&detJ);CHKERRQ(ierr); 943 /* CellGeometryFEM computes the expanded Jacobian, we want the true jacobian */ 944 for (i = 1; i < dim; i++) { 945 for (j = 0; j < depth; j++) { 946 J[i * depth + j] = J[i * dim + j]; 947 } 948 } 949 ierr = PetscQuadratureDestroy(&origin);CHKERRQ(ierr); 950 ierr = PetscDualSpaceGetPointSubspace(dsp,refPoint,&dsubsp);CHKERRQ(ierr); 951 ierr = PetscSpaceCreateSubspace(bsp,dsubsp,v,J,NULL,NULL,PETSC_OWN_POINTER,&bsubsp);CHKERRQ(ierr); 952 ierr = PetscSpaceSetUp(bsubsp);CHKERRQ(ierr); 953 ierr = PetscFECreate(PetscObjectComm((PetscObject)fe),trFE);CHKERRQ(ierr); 954 ierr = PetscFEGetType(fe,&type);CHKERRQ(ierr); 955 ierr = PetscFESetType(*trFE,type);CHKERRQ(ierr); 956 ierr = PetscFEGetNumComponents(fe,&numComp);CHKERRQ(ierr); 957 ierr = PetscFESetNumComponents(*trFE,numComp);CHKERRQ(ierr); 958 ierr = PetscFESetBasisSpace(*trFE,bsubsp);CHKERRQ(ierr); 959 ierr = PetscFESetDualSpace(*trFE,dsubsp);CHKERRQ(ierr); 960 ierr = PetscFEGetQuadrature(fe,&fullQuad);CHKERRQ(ierr); 961 ierr = PetscQuadratureGetOrder(fullQuad,&order);CHKERRQ(ierr); 962 ierr = DMPlexGetConeSize(dm,refPoint,&coneSize);CHKERRQ(ierr); 963 if (coneSize == 2 * depth) { 964 ierr = PetscDTGaussTensorQuadrature(depth,1,(order + 1)/2,-1.,1.,&subQuad);CHKERRQ(ierr); 965 } else { 966 ierr = PetscDTGaussJacobiQuadrature(depth,1,(order + 1)/2,-1.,1.,&subQuad);CHKERRQ(ierr); 967 } 968 ierr = PetscFESetQuadrature(*trFE,subQuad);CHKERRQ(ierr); 969 ierr = PetscFESetUp(*trFE);CHKERRQ(ierr); 970 ierr = PetscQuadratureDestroy(&subQuad);CHKERRQ(ierr); 971 ierr = PetscSpaceDestroy(&bsubsp);CHKERRQ(ierr); 972 PetscFunctionReturn(0); 973 } 974 975 PetscErrorCode PetscFECreateHeightTrace(PetscFE fe, PetscInt height, PetscFE *trFE) 976 { 977 PetscInt hStart, hEnd; 978 PetscDualSpace dsp; 979 DM dm; 980 PetscErrorCode ierr; 981 982 PetscFunctionBegin; 983 PetscValidHeaderSpecific(fe,PETSCFE_CLASSID,1); 984 PetscValidPointer(trFE,3); 985 *trFE = NULL; 986 ierr = PetscFEGetDualSpace(fe,&dsp);CHKERRQ(ierr); 987 ierr = PetscDualSpaceGetDM(dsp,&dm);CHKERRQ(ierr); 988 ierr = DMPlexGetHeightStratum(dm,height,&hStart,&hEnd);CHKERRQ(ierr); 989 if (hEnd <= hStart) PetscFunctionReturn(0); 990 ierr = PetscFECreatePointTrace(fe,hStart,trFE);CHKERRQ(ierr); 991 PetscFunctionReturn(0); 992 } 993 994 995 /*@ 996 PetscFEGetDimension - Get the dimension of the finite element space on a cell 997 998 Not collective 999 1000 Input Parameter: 1001 . fe - The PetscFE 1002 1003 Output Parameter: 1004 . dim - The dimension 1005 1006 Level: intermediate 1007 1008 .seealso: PetscFECreate(), PetscSpaceGetDimension(), PetscDualSpaceGetDimension() 1009 @*/ 1010 PetscErrorCode PetscFEGetDimension(PetscFE fem, PetscInt *dim) 1011 { 1012 PetscErrorCode ierr; 1013 1014 PetscFunctionBegin; 1015 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 1016 PetscValidPointer(dim, 2); 1017 if (fem->ops->getdimension) {ierr = (*fem->ops->getdimension)(fem, dim);CHKERRQ(ierr);} 1018 PetscFunctionReturn(0); 1019 } 1020 1021 /* 1022 Purpose: Compute element vector for chunk of elements 1023 1024 Input: 1025 Sizes: 1026 Ne: number of elements 1027 Nf: number of fields 1028 PetscFE 1029 dim: spatial dimension 1030 Nb: number of basis functions 1031 Nc: number of field components 1032 PetscQuadrature 1033 Nq: number of quadrature points 1034 1035 Geometry: 1036 PetscFEGeom[Ne] possibly *Nq 1037 PetscReal v0s[dim] 1038 PetscReal n[dim] 1039 PetscReal jacobians[dim*dim] 1040 PetscReal jacobianInverses[dim*dim] 1041 PetscReal jacobianDeterminants 1042 FEM: 1043 PetscFE 1044 PetscQuadrature 1045 PetscReal quadPoints[Nq*dim] 1046 PetscReal quadWeights[Nq] 1047 PetscReal basis[Nq*Nb*Nc] 1048 PetscReal basisDer[Nq*Nb*Nc*dim] 1049 PetscScalar coefficients[Ne*Nb*Nc] 1050 PetscScalar elemVec[Ne*Nb*Nc] 1051 1052 Problem: 1053 PetscInt f: the active field 1054 f0, f1 1055 1056 Work Space: 1057 PetscFE 1058 PetscScalar f0[Nq*dim]; 1059 PetscScalar f1[Nq*dim*dim]; 1060 PetscScalar u[Nc]; 1061 PetscScalar gradU[Nc*dim]; 1062 PetscReal x[dim]; 1063 PetscScalar realSpaceDer[dim]; 1064 1065 Purpose: Compute element vector for N_cb batches of elements 1066 1067 Input: 1068 Sizes: 1069 N_cb: Number of serial cell batches 1070 1071 Geometry: 1072 PetscReal v0s[Ne*dim] 1073 PetscReal jacobians[Ne*dim*dim] possibly *Nq 1074 PetscReal jacobianInverses[Ne*dim*dim] possibly *Nq 1075 PetscReal jacobianDeterminants[Ne] possibly *Nq 1076 FEM: 1077 static PetscReal quadPoints[Nq*dim] 1078 static PetscReal quadWeights[Nq] 1079 static PetscReal basis[Nq*Nb*Nc] 1080 static PetscReal basisDer[Nq*Nb*Nc*dim] 1081 PetscScalar coefficients[Ne*Nb*Nc] 1082 PetscScalar elemVec[Ne*Nb*Nc] 1083 1084 ex62.c: 1085 PetscErrorCode PetscFEIntegrateResidualBatch(PetscInt Ne, PetscInt numFields, PetscInt field, PetscQuadrature quad[], const PetscScalar coefficients[], 1086 const PetscReal v0s[], const PetscReal jacobians[], const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], 1087 void (*f0_func)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar f0[]), 1088 void (*f1_func)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar f1[]), PetscScalar elemVec[]) 1089 1090 ex52.c: 1091 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) 1092 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) 1093 1094 ex52_integrateElement.cu 1095 __global__ void integrateElementQuadrature(int N_cb, realType *coefficients, realType *jacobianInverses, realType *jacobianDeterminants, realType *elemVec) 1096 1097 PETSC_EXTERN PetscErrorCode IntegrateElementBatchGPU(PetscInt spatial_dim, PetscInt Ne, PetscInt Ncb, PetscInt Nbc, PetscInt Nbl, const PetscScalar coefficients[], 1098 const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], PetscScalar elemVec[], 1099 PetscLogEvent event, PetscInt debug, PetscInt pde_op) 1100 1101 ex52_integrateElementOpenCL.c: 1102 PETSC_EXTERN PetscErrorCode IntegrateElementBatchGPU(PetscInt spatial_dim, PetscInt Ne, PetscInt Ncb, PetscInt Nbc, PetscInt N_bl, const PetscScalar coefficients[], 1103 const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], PetscScalar elemVec[], 1104 PetscLogEvent event, PetscInt debug, PetscInt pde_op) 1105 1106 __kernel void integrateElementQuadrature(int N_cb, __global float *coefficients, __global float *jacobianInverses, __global float *jacobianDeterminants, __global float *elemVec) 1107 */ 1108 1109 /*@C 1110 PetscFEIntegrate - Produce the integral for the given field for a chunk of elements by quadrature integration 1111 1112 Not collective 1113 1114 Input Parameters: 1115 + fem - The PetscFE object for the field being integrated 1116 . prob - The PetscDS specifying the discretizations and continuum functions 1117 . field - The field being integrated 1118 . Ne - The number of elements in the chunk 1119 . cgeom - The cell geometry for each cell in the chunk 1120 . coefficients - The array of FEM basis coefficients for the elements 1121 . probAux - The PetscDS specifying the auxiliary discretizations 1122 - coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 1123 1124 Output Parameter 1125 . integral - the integral for this field 1126 1127 Level: developer 1128 1129 .seealso: PetscFEIntegrateResidual() 1130 @*/ 1131 PetscErrorCode PetscFEIntegrate(PetscFE fem, PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom, 1132 const PetscScalar coefficients[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscScalar integral[]) 1133 { 1134 PetscErrorCode ierr; 1135 1136 PetscFunctionBegin; 1137 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 1138 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 2); 1139 if (fem->ops->integrate) {ierr = (*fem->ops->integrate)(fem, prob, field, Ne, cgeom, coefficients, probAux, coefficientsAux, integral);CHKERRQ(ierr);} 1140 PetscFunctionReturn(0); 1141 } 1142 1143 /*@C 1144 PetscFEIntegrateBd - Produce the integral for the given field for a chunk of elements by quadrature integration 1145 1146 Not collective 1147 1148 Input Parameters: 1149 + fem - The PetscFE object for the field being integrated 1150 . prob - The PetscDS specifying the discretizations and continuum functions 1151 . field - The field being integrated 1152 . obj_func - The function to be integrated 1153 . Ne - The number of elements in the chunk 1154 . fgeom - The face geometry for each face in the chunk 1155 . coefficients - The array of FEM basis coefficients for the elements 1156 . probAux - The PetscDS specifying the auxiliary discretizations 1157 - coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 1158 1159 Output Parameter 1160 . integral - the integral for this field 1161 1162 Level: developer 1163 1164 .seealso: PetscFEIntegrateResidual() 1165 @*/ 1166 PetscErrorCode PetscFEIntegrateBd(PetscFE fem, PetscDS prob, PetscInt field, 1167 void (*obj_func)(PetscInt, PetscInt, PetscInt, 1168 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 1169 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 1170 PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 1171 PetscInt Ne, PetscFEGeom *geom, const PetscScalar coefficients[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscScalar integral[]) 1172 { 1173 PetscErrorCode ierr; 1174 1175 PetscFunctionBegin; 1176 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 1177 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 2); 1178 if (fem->ops->integratebd) {ierr = (*fem->ops->integratebd)(fem, prob, field, obj_func, Ne, geom, coefficients, probAux, coefficientsAux, integral);CHKERRQ(ierr);} 1179 PetscFunctionReturn(0); 1180 } 1181 1182 /*@C 1183 PetscFEIntegrateResidual - Produce the element residual vector for a chunk of elements by quadrature integration 1184 1185 Not collective 1186 1187 Input Parameters: 1188 + fem - The PetscFE object for the field being integrated 1189 . prob - The PetscDS specifying the discretizations and continuum functions 1190 . field - The field being integrated 1191 . Ne - The number of elements in the chunk 1192 . cgeom - The cell geometry for each cell in the chunk 1193 . coefficients - The array of FEM basis coefficients for the elements 1194 . coefficients_t - The array of FEM basis time derivative coefficients for the elements 1195 . probAux - The PetscDS specifying the auxiliary discretizations 1196 . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 1197 - t - The time 1198 1199 Output Parameter 1200 . elemVec - the element residual vectors from each element 1201 1202 Note: 1203 $ Loop over batch of elements (e): 1204 $ Loop over quadrature points (q): 1205 $ Make u_q and gradU_q (loops over fields,Nb,Ncomp) and x_q 1206 $ Call f_0 and f_1 1207 $ Loop over element vector entries (f,fc --> i): 1208 $ elemVec[i] += \psi^{fc}_f(q) f0_{fc}(u, \nabla u) + \nabla\psi^{fc}_f(q) \cdot f1_{fc,df}(u, \nabla u) 1209 1210 Level: developer 1211 1212 .seealso: PetscFEIntegrateResidual() 1213 @*/ 1214 PetscErrorCode PetscFEIntegrateResidual(PetscFE fem, PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom, 1215 const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[]) 1216 { 1217 PetscErrorCode ierr; 1218 1219 PetscFunctionBegin; 1220 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 1221 PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 2); 1222 if (fem->ops->integrateresidual) {ierr = (*fem->ops->integrateresidual)(fem, prob, field, Ne, cgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec);CHKERRQ(ierr);} 1223 PetscFunctionReturn(0); 1224 } 1225 1226 /*@C 1227 PetscFEIntegrateBdResidual - Produce the element residual vector for a chunk of elements by quadrature integration over a boundary 1228 1229 Not collective 1230 1231 Input Parameters: 1232 + fem - The PetscFE object for the field being integrated 1233 . prob - The PetscDS specifying the discretizations and continuum functions 1234 . field - The field being integrated 1235 . Ne - The number of elements in the chunk 1236 . fgeom - The face geometry for each cell in the chunk 1237 . coefficients - The array of FEM basis coefficients for the elements 1238 . coefficients_t - The array of FEM basis time derivative coefficients for the elements 1239 . probAux - The PetscDS specifying the auxiliary discretizations 1240 . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 1241 - t - The time 1242 1243 Output Parameter 1244 . elemVec - the element residual vectors from each element 1245 1246 Level: developer 1247 1248 .seealso: PetscFEIntegrateResidual() 1249 @*/ 1250 PetscErrorCode PetscFEIntegrateBdResidual(PetscFE fem, PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *fgeom, 1251 const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[]) 1252 { 1253 PetscErrorCode ierr; 1254 1255 PetscFunctionBegin; 1256 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 1257 if (fem->ops->integratebdresidual) {ierr = (*fem->ops->integratebdresidual)(fem, prob, field, Ne, fgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec);CHKERRQ(ierr);} 1258 PetscFunctionReturn(0); 1259 } 1260 1261 /*@C 1262 PetscFEIntegrateJacobian - Produce the element Jacobian for a chunk of elements by quadrature integration 1263 1264 Not collective 1265 1266 Input Parameters: 1267 + fem - The PetscFE object for the field being integrated 1268 . prob - The PetscDS specifying the discretizations and continuum functions 1269 . jtype - The type of matrix pointwise functions that should be used 1270 . fieldI - The test field being integrated 1271 . fieldJ - The basis field being integrated 1272 . Ne - The number of elements in the chunk 1273 . cgeom - The cell geometry for each cell in the chunk 1274 . coefficients - The array of FEM basis coefficients for the elements for the Jacobian evaluation point 1275 . coefficients_t - The array of FEM basis time derivative coefficients for the elements 1276 . probAux - The PetscDS specifying the auxiliary discretizations 1277 . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 1278 . t - The time 1279 - u_tShift - A multiplier for the dF/du_t term (as opposed to the dF/du term) 1280 1281 Output Parameter 1282 . elemMat - the element matrices for the Jacobian from each element 1283 1284 Note: 1285 $ Loop over batch of elements (e): 1286 $ Loop over element matrix entries (f,fc,g,gc --> i,j): 1287 $ Loop over quadrature points (q): 1288 $ Make u_q and gradU_q (loops over fields,Nb,Ncomp) 1289 $ elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q) 1290 $ + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q) 1291 $ + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q) 1292 $ + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q) 1293 Level: developer 1294 1295 .seealso: PetscFEIntegrateResidual() 1296 @*/ 1297 PetscErrorCode PetscFEIntegrateJacobian(PetscFE fem, PetscDS prob, PetscFEJacobianType jtype, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *cgeom, 1298 const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) 1299 { 1300 PetscErrorCode ierr; 1301 1302 PetscFunctionBegin; 1303 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 1304 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);} 1305 PetscFunctionReturn(0); 1306 } 1307 1308 /*@C 1309 PetscFEIntegrateBdJacobian - Produce the boundary element Jacobian for a chunk of elements by quadrature integration 1310 1311 Not collective 1312 1313 Input Parameters: 1314 + fem = The PetscFE object for the field being integrated 1315 . prob - The PetscDS specifying the discretizations and continuum functions 1316 . fieldI - The test field being integrated 1317 . fieldJ - The basis field being integrated 1318 . Ne - The number of elements in the chunk 1319 . fgeom - The face geometry for each cell in the chunk 1320 . coefficients - The array of FEM basis coefficients for the elements for the Jacobian evaluation point 1321 . coefficients_t - The array of FEM basis time derivative coefficients for the elements 1322 . probAux - The PetscDS specifying the auxiliary discretizations 1323 . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 1324 . t - The time 1325 - u_tShift - A multiplier for the dF/du_t term (as opposed to the dF/du term) 1326 1327 Output Parameter 1328 . elemMat - the element matrices for the Jacobian from each element 1329 1330 Note: 1331 $ Loop over batch of elements (e): 1332 $ Loop over element matrix entries (f,fc,g,gc --> i,j): 1333 $ Loop over quadrature points (q): 1334 $ Make u_q and gradU_q (loops over fields,Nb,Ncomp) 1335 $ elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q) 1336 $ + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q) 1337 $ + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q) 1338 $ + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q) 1339 Level: developer 1340 1341 .seealso: PetscFEIntegrateJacobian(), PetscFEIntegrateResidual() 1342 @*/ 1343 PetscErrorCode PetscFEIntegrateBdJacobian(PetscFE fem, PetscDS prob, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *fgeom, 1344 const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) 1345 { 1346 PetscErrorCode ierr; 1347 1348 PetscFunctionBegin; 1349 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 1350 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);} 1351 PetscFunctionReturn(0); 1352 } 1353 1354 PetscErrorCode PetscFEGetHeightSubspace(PetscFE fe, PetscInt height, PetscFE *subfe) 1355 { 1356 PetscSpace P, subP; 1357 PetscDualSpace Q, subQ; 1358 PetscQuadrature subq; 1359 PetscFEType fetype; 1360 PetscInt dim, Nc; 1361 PetscErrorCode ierr; 1362 1363 PetscFunctionBegin; 1364 PetscValidHeaderSpecific(fe, PETSCFE_CLASSID, 1); 1365 PetscValidPointer(subfe, 3); 1366 if (height == 0) { 1367 *subfe = fe; 1368 PetscFunctionReturn(0); 1369 } 1370 ierr = PetscFEGetBasisSpace(fe, &P);CHKERRQ(ierr); 1371 ierr = PetscFEGetDualSpace(fe, &Q);CHKERRQ(ierr); 1372 ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 1373 ierr = PetscFEGetFaceQuadrature(fe, &subq);CHKERRQ(ierr); 1374 ierr = PetscDualSpaceGetDimension(Q, &dim);CHKERRQ(ierr); 1375 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);} 1376 if (!fe->subspaces) {ierr = PetscCalloc1(dim, &fe->subspaces);CHKERRQ(ierr);} 1377 if (height <= dim) { 1378 if (!fe->subspaces[height-1]) { 1379 PetscFE sub; 1380 1381 ierr = PetscSpaceGetHeightSubspace(P, height, &subP);CHKERRQ(ierr); 1382 ierr = PetscDualSpaceGetHeightSubspace(Q, height, &subQ);CHKERRQ(ierr); 1383 ierr = PetscFECreate(PetscObjectComm((PetscObject) fe), &sub);CHKERRQ(ierr); 1384 ierr = PetscFEGetType(fe, &fetype);CHKERRQ(ierr); 1385 ierr = PetscFESetType(sub, fetype);CHKERRQ(ierr); 1386 ierr = PetscFESetBasisSpace(sub, subP);CHKERRQ(ierr); 1387 ierr = PetscFESetDualSpace(sub, subQ);CHKERRQ(ierr); 1388 ierr = PetscFESetNumComponents(sub, Nc);CHKERRQ(ierr); 1389 ierr = PetscFESetUp(sub);CHKERRQ(ierr); 1390 ierr = PetscFESetQuadrature(sub, subq);CHKERRQ(ierr); 1391 fe->subspaces[height-1] = sub; 1392 } 1393 *subfe = fe->subspaces[height-1]; 1394 } else { 1395 *subfe = NULL; 1396 } 1397 PetscFunctionReturn(0); 1398 } 1399 1400 /*@ 1401 PetscFERefine - Create a "refined" PetscFE object that refines the reference cell into smaller copies. This is typically used 1402 to precondition a higher order method with a lower order method on a refined mesh having the same number of dofs (but more 1403 sparsity). It is also used to create an interpolation between regularly refined meshes. 1404 1405 Collective on PetscFE 1406 1407 Input Parameter: 1408 . fe - The initial PetscFE 1409 1410 Output Parameter: 1411 . feRef - The refined PetscFE 1412 1413 Level: developer 1414 1415 .seealso: PetscFEType, PetscFECreate(), PetscFESetType() 1416 @*/ 1417 PetscErrorCode PetscFERefine(PetscFE fe, PetscFE *feRef) 1418 { 1419 PetscSpace P, Pref; 1420 PetscDualSpace Q, Qref; 1421 DM K, Kref; 1422 PetscQuadrature q, qref; 1423 const PetscReal *v0, *jac; 1424 PetscInt numComp, numSubelements; 1425 PetscErrorCode ierr; 1426 1427 PetscFunctionBegin; 1428 ierr = PetscFEGetBasisSpace(fe, &P);CHKERRQ(ierr); 1429 ierr = PetscFEGetDualSpace(fe, &Q);CHKERRQ(ierr); 1430 ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr); 1431 ierr = PetscDualSpaceGetDM(Q, &K);CHKERRQ(ierr); 1432 /* Create space */ 1433 ierr = PetscObjectReference((PetscObject) P);CHKERRQ(ierr); 1434 Pref = P; 1435 /* Create dual space */ 1436 ierr = PetscDualSpaceDuplicate(Q, &Qref);CHKERRQ(ierr); 1437 ierr = DMRefine(K, PetscObjectComm((PetscObject) fe), &Kref);CHKERRQ(ierr); 1438 ierr = PetscDualSpaceSetDM(Qref, Kref);CHKERRQ(ierr); 1439 ierr = DMDestroy(&Kref);CHKERRQ(ierr); 1440 ierr = PetscDualSpaceSetUp(Qref);CHKERRQ(ierr); 1441 /* Create element */ 1442 ierr = PetscFECreate(PetscObjectComm((PetscObject) fe), feRef);CHKERRQ(ierr); 1443 ierr = PetscFESetType(*feRef, PETSCFECOMPOSITE);CHKERRQ(ierr); 1444 ierr = PetscFESetBasisSpace(*feRef, Pref);CHKERRQ(ierr); 1445 ierr = PetscFESetDualSpace(*feRef, Qref);CHKERRQ(ierr); 1446 ierr = PetscFEGetNumComponents(fe, &numComp);CHKERRQ(ierr); 1447 ierr = PetscFESetNumComponents(*feRef, numComp);CHKERRQ(ierr); 1448 ierr = PetscFESetUp(*feRef);CHKERRQ(ierr); 1449 ierr = PetscSpaceDestroy(&Pref);CHKERRQ(ierr); 1450 ierr = PetscDualSpaceDestroy(&Qref);CHKERRQ(ierr); 1451 /* Create quadrature */ 1452 ierr = PetscFECompositeGetMapping(*feRef, &numSubelements, &v0, &jac, NULL);CHKERRQ(ierr); 1453 ierr = PetscQuadratureExpandComposite(q, numSubelements, v0, jac, &qref);CHKERRQ(ierr); 1454 ierr = PetscFESetQuadrature(*feRef, qref);CHKERRQ(ierr); 1455 ierr = PetscQuadratureDestroy(&qref);CHKERRQ(ierr); 1456 PetscFunctionReturn(0); 1457 } 1458 1459 /*@C 1460 PetscFECreateDefault - Create a PetscFE for basic FEM computation 1461 1462 Collective on DM 1463 1464 Input Parameters: 1465 + comm - The MPI comm 1466 . dim - The spatial dimension 1467 . Nc - The number of components 1468 . isSimplex - Flag for simplex reference cell, otherwise its a tensor product 1469 . prefix - The options prefix, or NULL 1470 - qorder - The quadrature order 1471 1472 Output Parameter: 1473 . fem - The PetscFE object 1474 1475 Level: beginner 1476 1477 .keywords: PetscFE, finite element 1478 .seealso: PetscFECreate(), PetscSpaceCreate(), PetscDualSpaceCreate() 1479 @*/ 1480 PetscErrorCode PetscFECreateDefault(MPI_Comm comm, PetscInt dim, PetscInt Nc, PetscBool isSimplex, const char prefix[], PetscInt qorder, PetscFE *fem) 1481 { 1482 PetscQuadrature q, fq; 1483 DM K; 1484 PetscSpace P; 1485 PetscDualSpace Q; 1486 PetscInt order, quadPointsPerEdge; 1487 PetscBool tensor = isSimplex ? PETSC_FALSE : PETSC_TRUE; 1488 PetscErrorCode ierr; 1489 1490 PetscFunctionBegin; 1491 /* Create space */ 1492 ierr = PetscSpaceCreate(comm, &P);CHKERRQ(ierr); 1493 ierr = PetscObjectSetOptionsPrefix((PetscObject) P, prefix);CHKERRQ(ierr); 1494 ierr = PetscSpacePolynomialSetTensor(P, tensor);CHKERRQ(ierr); 1495 ierr = PetscSpaceSetNumComponents(P, Nc);CHKERRQ(ierr); 1496 ierr = PetscSpaceSetNumVariables(P, dim);CHKERRQ(ierr); 1497 ierr = PetscSpaceSetFromOptions(P);CHKERRQ(ierr); 1498 ierr = PetscSpaceSetUp(P);CHKERRQ(ierr); 1499 ierr = PetscSpaceGetDegree(P, &order, NULL);CHKERRQ(ierr); 1500 ierr = PetscSpacePolynomialGetTensor(P, &tensor);CHKERRQ(ierr); 1501 /* Create dual space */ 1502 ierr = PetscDualSpaceCreate(comm, &Q);CHKERRQ(ierr); 1503 ierr = PetscDualSpaceSetType(Q,PETSCDUALSPACELAGRANGE);CHKERRQ(ierr); 1504 ierr = PetscObjectSetOptionsPrefix((PetscObject) Q, prefix);CHKERRQ(ierr); 1505 ierr = PetscDualSpaceCreateReferenceCell(Q, dim, isSimplex, &K);CHKERRQ(ierr); 1506 ierr = PetscDualSpaceSetDM(Q, K);CHKERRQ(ierr); 1507 ierr = DMDestroy(&K);CHKERRQ(ierr); 1508 ierr = PetscDualSpaceSetNumComponents(Q, Nc);CHKERRQ(ierr); 1509 ierr = PetscDualSpaceSetOrder(Q, order);CHKERRQ(ierr); 1510 ierr = PetscDualSpaceLagrangeSetTensor(Q, tensor);CHKERRQ(ierr); 1511 ierr = PetscDualSpaceSetFromOptions(Q);CHKERRQ(ierr); 1512 ierr = PetscDualSpaceSetUp(Q);CHKERRQ(ierr); 1513 /* Create element */ 1514 ierr = PetscFECreate(comm, fem);CHKERRQ(ierr); 1515 ierr = PetscObjectSetOptionsPrefix((PetscObject) *fem, prefix);CHKERRQ(ierr); 1516 ierr = PetscFESetFromOptions(*fem);CHKERRQ(ierr); 1517 ierr = PetscFESetBasisSpace(*fem, P);CHKERRQ(ierr); 1518 ierr = PetscFESetDualSpace(*fem, Q);CHKERRQ(ierr); 1519 ierr = PetscFESetNumComponents(*fem, Nc);CHKERRQ(ierr); 1520 ierr = PetscFESetUp(*fem);CHKERRQ(ierr); 1521 ierr = PetscSpaceDestroy(&P);CHKERRQ(ierr); 1522 ierr = PetscDualSpaceDestroy(&Q);CHKERRQ(ierr); 1523 /* Create quadrature (with specified order if given) */ 1524 qorder = qorder >= 0 ? qorder : order; 1525 ierr = PetscObjectOptionsBegin((PetscObject)*fem);CHKERRQ(ierr); 1526 ierr = PetscOptionsInt("-petscfe_default_quadrature_order","Quadrature order is one less than quadture points per edge","PetscFECreateDefault",qorder,&qorder,NULL);CHKERRQ(ierr); 1527 ierr = PetscOptionsEnd();CHKERRQ(ierr); 1528 quadPointsPerEdge = PetscMax(qorder + 1,1); 1529 if (isSimplex) { 1530 ierr = PetscDTGaussJacobiQuadrature(dim, 1, quadPointsPerEdge, -1.0, 1.0, &q);CHKERRQ(ierr); 1531 ierr = PetscDTGaussJacobiQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0, &fq);CHKERRQ(ierr); 1532 } else { 1533 ierr = PetscDTGaussTensorQuadrature(dim, 1, quadPointsPerEdge, -1.0, 1.0, &q);CHKERRQ(ierr); 1534 ierr = PetscDTGaussTensorQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0, &fq);CHKERRQ(ierr); 1535 } 1536 ierr = PetscFESetQuadrature(*fem, q);CHKERRQ(ierr); 1537 ierr = PetscFESetFaceQuadrature(*fem, fq);CHKERRQ(ierr); 1538 ierr = PetscQuadratureDestroy(&q);CHKERRQ(ierr); 1539 ierr = PetscQuadratureDestroy(&fq);CHKERRQ(ierr); 1540 PetscFunctionReturn(0); 1541 } 1542