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