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