120cf1dd8SToby Isaac /* Basis Jet Tabulation 220cf1dd8SToby Isaac 320cf1dd8SToby Isaac We would like to tabulate the nodal basis functions and derivatives at a set of points, usually quadrature points. We 420cf1dd8SToby Isaac follow here the derviation in http://www.math.ttu.edu/~kirby/papers/fiat-toms-2004.pdf. The nodal basis $\psi_i$ can 520cf1dd8SToby Isaac be expressed in terms of a prime basis $\phi_i$ which can be stably evaluated. In PETSc, we will use the Legendre basis 620cf1dd8SToby Isaac as a prime basis. 720cf1dd8SToby Isaac 820cf1dd8SToby Isaac \psi_i = \sum_k \alpha_{ki} \phi_k 920cf1dd8SToby Isaac 1020cf1dd8SToby Isaac Our nodal basis is defined in terms of the dual basis $n_j$ 1120cf1dd8SToby Isaac 1220cf1dd8SToby Isaac n_j \cdot \psi_i = \delta_{ji} 1320cf1dd8SToby Isaac 1420cf1dd8SToby Isaac and we may act on the first equation to obtain 1520cf1dd8SToby Isaac 1620cf1dd8SToby Isaac n_j \cdot \psi_i = \sum_k \alpha_{ki} n_j \cdot \phi_k 1720cf1dd8SToby Isaac \delta_{ji} = \sum_k \alpha_{ki} V_{jk} 1820cf1dd8SToby Isaac I = V \alpha 1920cf1dd8SToby Isaac 2020cf1dd8SToby Isaac so the coefficients of the nodal basis in the prime basis are 2120cf1dd8SToby Isaac 2220cf1dd8SToby Isaac \alpha = V^{-1} 2320cf1dd8SToby Isaac 2420cf1dd8SToby Isaac We will define the dual basis vectors $n_j$ using a quadrature rule. 2520cf1dd8SToby Isaac 2620cf1dd8SToby Isaac Right now, we will just use the polynomial spaces P^k. I know some elements use the space of symmetric polynomials 2720cf1dd8SToby Isaac (I think Nedelec), but we will neglect this for now. Constraints in the space, e.g. Arnold-Winther elements, can 2820cf1dd8SToby Isaac be implemented exactly as in FIAT using functionals $L_j$. 2920cf1dd8SToby Isaac 3020cf1dd8SToby Isaac I will have to count the degrees correctly for the Legendre product when we are on simplices. 3120cf1dd8SToby Isaac 3220cf1dd8SToby Isaac We will have three objects: 3320cf1dd8SToby Isaac - Space, P: this just need point evaluation I think 3420cf1dd8SToby 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 3520cf1dd8SToby Isaac - FEM: This keeps {P, P', Q} 3620cf1dd8SToby Isaac */ 3720cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/ 3820cf1dd8SToby Isaac #include <petscdmplex.h> 3920cf1dd8SToby Isaac 4020cf1dd8SToby Isaac PetscBool FEcite = PETSC_FALSE; 4120cf1dd8SToby Isaac const char FECitation[] = "@article{kirby2004,\n" 4220cf1dd8SToby Isaac " title = {Algorithm 839: FIAT, a New Paradigm for Computing Finite Element Basis Functions},\n" 4320cf1dd8SToby Isaac " journal = {ACM Transactions on Mathematical Software},\n" 4420cf1dd8SToby Isaac " author = {Robert C. Kirby},\n" 4520cf1dd8SToby Isaac " volume = {30},\n" 4620cf1dd8SToby Isaac " number = {4},\n" 4720cf1dd8SToby Isaac " pages = {502--516},\n" 4820cf1dd8SToby Isaac " doi = {10.1145/1039813.1039820},\n" 4920cf1dd8SToby Isaac " year = {2004}\n}\n"; 5020cf1dd8SToby Isaac 5120cf1dd8SToby Isaac PetscClassId PETSCFE_CLASSID = 0; 5220cf1dd8SToby Isaac 5320cf1dd8SToby Isaac PetscFunctionList PetscFEList = NULL; 5420cf1dd8SToby Isaac PetscBool PetscFERegisterAllCalled = PETSC_FALSE; 5520cf1dd8SToby Isaac 5620cf1dd8SToby Isaac /*@C 5720cf1dd8SToby Isaac PetscFERegister - Adds a new PetscFE implementation 5820cf1dd8SToby Isaac 5920cf1dd8SToby Isaac Not Collective 6020cf1dd8SToby Isaac 6120cf1dd8SToby Isaac Input Parameters: 6220cf1dd8SToby Isaac + name - The name of a new user-defined creation routine 6320cf1dd8SToby Isaac - create_func - The creation routine itself 6420cf1dd8SToby Isaac 6520cf1dd8SToby Isaac Notes: 6620cf1dd8SToby Isaac PetscFERegister() may be called multiple times to add several user-defined PetscFEs 6720cf1dd8SToby Isaac 6820cf1dd8SToby Isaac Sample usage: 6920cf1dd8SToby Isaac .vb 7020cf1dd8SToby Isaac PetscFERegister("my_fe", MyPetscFECreate); 7120cf1dd8SToby Isaac .ve 7220cf1dd8SToby Isaac 7320cf1dd8SToby Isaac Then, your PetscFE type can be chosen with the procedural interface via 7420cf1dd8SToby Isaac .vb 7520cf1dd8SToby Isaac PetscFECreate(MPI_Comm, PetscFE *); 7620cf1dd8SToby Isaac PetscFESetType(PetscFE, "my_fe"); 7720cf1dd8SToby Isaac .ve 7820cf1dd8SToby Isaac or at runtime via the option 7920cf1dd8SToby Isaac .vb 8020cf1dd8SToby Isaac -petscfe_type my_fe 8120cf1dd8SToby Isaac .ve 8220cf1dd8SToby Isaac 8320cf1dd8SToby Isaac Level: advanced 8420cf1dd8SToby Isaac 8520cf1dd8SToby Isaac .keywords: PetscFE, register 8620cf1dd8SToby Isaac .seealso: PetscFERegisterAll(), PetscFERegisterDestroy() 8720cf1dd8SToby Isaac 8820cf1dd8SToby Isaac @*/ 8920cf1dd8SToby Isaac PetscErrorCode PetscFERegister(const char sname[], PetscErrorCode (*function)(PetscFE)) 9020cf1dd8SToby Isaac { 9120cf1dd8SToby Isaac PetscErrorCode ierr; 9220cf1dd8SToby Isaac 9320cf1dd8SToby Isaac PetscFunctionBegin; 9420cf1dd8SToby Isaac ierr = PetscFunctionListAdd(&PetscFEList, sname, function);CHKERRQ(ierr); 9520cf1dd8SToby Isaac PetscFunctionReturn(0); 9620cf1dd8SToby Isaac } 9720cf1dd8SToby Isaac 9820cf1dd8SToby Isaac /*@C 9920cf1dd8SToby Isaac PetscFESetType - Builds a particular PetscFE 10020cf1dd8SToby Isaac 10120cf1dd8SToby Isaac Collective on PetscFE 10220cf1dd8SToby Isaac 10320cf1dd8SToby Isaac Input Parameters: 10420cf1dd8SToby Isaac + fem - The PetscFE object 10520cf1dd8SToby Isaac - name - The kind of FEM space 10620cf1dd8SToby Isaac 10720cf1dd8SToby Isaac Options Database Key: 10820cf1dd8SToby Isaac . -petscfe_type <type> - Sets the PetscFE type; use -help for a list of available types 10920cf1dd8SToby Isaac 11020cf1dd8SToby Isaac Level: intermediate 11120cf1dd8SToby Isaac 11220cf1dd8SToby Isaac .keywords: PetscFE, set, type 11320cf1dd8SToby Isaac .seealso: PetscFEGetType(), PetscFECreate() 11420cf1dd8SToby Isaac @*/ 11520cf1dd8SToby Isaac PetscErrorCode PetscFESetType(PetscFE fem, PetscFEType name) 11620cf1dd8SToby Isaac { 11720cf1dd8SToby Isaac PetscErrorCode (*r)(PetscFE); 11820cf1dd8SToby Isaac PetscBool match; 11920cf1dd8SToby Isaac PetscErrorCode ierr; 12020cf1dd8SToby Isaac 12120cf1dd8SToby Isaac PetscFunctionBegin; 12220cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 12320cf1dd8SToby Isaac ierr = PetscObjectTypeCompare((PetscObject) fem, name, &match);CHKERRQ(ierr); 12420cf1dd8SToby Isaac if (match) PetscFunctionReturn(0); 12520cf1dd8SToby Isaac 12620cf1dd8SToby Isaac if (!PetscFERegisterAllCalled) {ierr = PetscFERegisterAll();CHKERRQ(ierr);} 12720cf1dd8SToby Isaac ierr = PetscFunctionListFind(PetscFEList, name, &r);CHKERRQ(ierr); 12820cf1dd8SToby Isaac if (!r) SETERRQ1(PetscObjectComm((PetscObject) fem), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscFE type: %s", name); 12920cf1dd8SToby Isaac 13020cf1dd8SToby Isaac if (fem->ops->destroy) { 13120cf1dd8SToby Isaac ierr = (*fem->ops->destroy)(fem);CHKERRQ(ierr); 13220cf1dd8SToby Isaac fem->ops->destroy = NULL; 13320cf1dd8SToby Isaac } 13420cf1dd8SToby Isaac ierr = (*r)(fem);CHKERRQ(ierr); 13520cf1dd8SToby Isaac ierr = PetscObjectChangeTypeName((PetscObject) fem, name);CHKERRQ(ierr); 13620cf1dd8SToby Isaac PetscFunctionReturn(0); 13720cf1dd8SToby Isaac } 13820cf1dd8SToby Isaac 13920cf1dd8SToby Isaac /*@C 14020cf1dd8SToby Isaac PetscFEGetType - Gets the PetscFE type name (as a string) from the object. 14120cf1dd8SToby Isaac 14220cf1dd8SToby Isaac Not Collective 14320cf1dd8SToby Isaac 14420cf1dd8SToby Isaac Input Parameter: 14520cf1dd8SToby Isaac . fem - The PetscFE 14620cf1dd8SToby Isaac 14720cf1dd8SToby Isaac Output Parameter: 14820cf1dd8SToby Isaac . name - The PetscFE type name 14920cf1dd8SToby Isaac 15020cf1dd8SToby Isaac Level: intermediate 15120cf1dd8SToby Isaac 15220cf1dd8SToby Isaac .keywords: PetscFE, get, type, name 15320cf1dd8SToby Isaac .seealso: PetscFESetType(), PetscFECreate() 15420cf1dd8SToby Isaac @*/ 15520cf1dd8SToby Isaac PetscErrorCode PetscFEGetType(PetscFE fem, PetscFEType *name) 15620cf1dd8SToby Isaac { 15720cf1dd8SToby Isaac PetscErrorCode ierr; 15820cf1dd8SToby Isaac 15920cf1dd8SToby Isaac PetscFunctionBegin; 16020cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 16120cf1dd8SToby Isaac PetscValidPointer(name, 2); 16220cf1dd8SToby Isaac if (!PetscFERegisterAllCalled) { 16320cf1dd8SToby Isaac ierr = PetscFERegisterAll();CHKERRQ(ierr); 16420cf1dd8SToby Isaac } 16520cf1dd8SToby Isaac *name = ((PetscObject) fem)->type_name; 16620cf1dd8SToby Isaac PetscFunctionReturn(0); 16720cf1dd8SToby Isaac } 16820cf1dd8SToby Isaac 16920cf1dd8SToby Isaac /*@C 17020cf1dd8SToby Isaac PetscFEView - Views a PetscFE 17120cf1dd8SToby Isaac 17220cf1dd8SToby Isaac Collective on PetscFE 17320cf1dd8SToby Isaac 17420cf1dd8SToby Isaac Input Parameter: 17520cf1dd8SToby Isaac + fem - the PetscFE object to view 176d9bac1caSLisandro Dalcin - viewer - the viewer 17720cf1dd8SToby Isaac 17820cf1dd8SToby Isaac Level: developer 17920cf1dd8SToby Isaac 18020cf1dd8SToby Isaac .seealso PetscFEDestroy() 18120cf1dd8SToby Isaac @*/ 182d9bac1caSLisandro Dalcin PetscErrorCode PetscFEView(PetscFE fem, PetscViewer viewer) 18320cf1dd8SToby Isaac { 184d9bac1caSLisandro Dalcin PetscBool iascii; 18520cf1dd8SToby Isaac PetscErrorCode ierr; 18620cf1dd8SToby Isaac 18720cf1dd8SToby Isaac PetscFunctionBegin; 18820cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 189d9bac1caSLisandro Dalcin if (viewer) PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 190d9bac1caSLisandro Dalcin if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) fem), &viewer);CHKERRQ(ierr);} 191d9bac1caSLisandro Dalcin ierr = PetscObjectPrintClassNamePrefixType((PetscObject)fem, viewer);CHKERRQ(ierr); 192d9bac1caSLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 193d9bac1caSLisandro Dalcin if (fem->ops->view) {ierr = (*fem->ops->view)(fem, viewer);CHKERRQ(ierr);} 19420cf1dd8SToby Isaac PetscFunctionReturn(0); 19520cf1dd8SToby Isaac } 19620cf1dd8SToby Isaac 19720cf1dd8SToby Isaac /*@ 19820cf1dd8SToby Isaac PetscFESetFromOptions - sets parameters in a PetscFE from the options database 19920cf1dd8SToby Isaac 20020cf1dd8SToby Isaac Collective on PetscFE 20120cf1dd8SToby Isaac 20220cf1dd8SToby Isaac Input Parameter: 20320cf1dd8SToby Isaac . fem - the PetscFE object to set options for 20420cf1dd8SToby Isaac 20520cf1dd8SToby Isaac Options Database: 206*a2b725a8SWilliam Gropp + -petscfe_num_blocks - the number of cell blocks to integrate concurrently 207*a2b725a8SWilliam Gropp - -petscfe_num_batches - the number of cell batches to integrate serially 20820cf1dd8SToby Isaac 20920cf1dd8SToby Isaac Level: developer 21020cf1dd8SToby Isaac 21120cf1dd8SToby Isaac .seealso PetscFEView() 21220cf1dd8SToby Isaac @*/ 21320cf1dd8SToby Isaac PetscErrorCode PetscFESetFromOptions(PetscFE fem) 21420cf1dd8SToby Isaac { 21520cf1dd8SToby Isaac const char *defaultType; 21620cf1dd8SToby Isaac char name[256]; 21720cf1dd8SToby Isaac PetscBool flg; 21820cf1dd8SToby Isaac PetscErrorCode ierr; 21920cf1dd8SToby Isaac 22020cf1dd8SToby Isaac PetscFunctionBegin; 22120cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 22220cf1dd8SToby Isaac if (!((PetscObject) fem)->type_name) { 22320cf1dd8SToby Isaac defaultType = PETSCFEBASIC; 22420cf1dd8SToby Isaac } else { 22520cf1dd8SToby Isaac defaultType = ((PetscObject) fem)->type_name; 22620cf1dd8SToby Isaac } 22720cf1dd8SToby Isaac if (!PetscFERegisterAllCalled) {ierr = PetscFERegisterAll();CHKERRQ(ierr);} 22820cf1dd8SToby Isaac 22920cf1dd8SToby Isaac ierr = PetscObjectOptionsBegin((PetscObject) fem);CHKERRQ(ierr); 23020cf1dd8SToby Isaac ierr = PetscOptionsFList("-petscfe_type", "Finite element space", "PetscFESetType", PetscFEList, defaultType, name, 256, &flg);CHKERRQ(ierr); 23120cf1dd8SToby Isaac if (flg) { 23220cf1dd8SToby Isaac ierr = PetscFESetType(fem, name);CHKERRQ(ierr); 23320cf1dd8SToby Isaac } else if (!((PetscObject) fem)->type_name) { 23420cf1dd8SToby Isaac ierr = PetscFESetType(fem, defaultType);CHKERRQ(ierr); 23520cf1dd8SToby Isaac } 23620cf1dd8SToby Isaac ierr = PetscOptionsInt("-petscfe_num_blocks", "The number of cell blocks to integrate concurrently", "PetscSpaceSetTileSizes", fem->numBlocks, &fem->numBlocks, NULL);CHKERRQ(ierr); 23720cf1dd8SToby Isaac ierr = PetscOptionsInt("-petscfe_num_batches", "The number of cell batches to integrate serially", "PetscSpaceSetTileSizes", fem->numBatches, &fem->numBatches, NULL);CHKERRQ(ierr); 23820cf1dd8SToby Isaac if (fem->ops->setfromoptions) { 23920cf1dd8SToby Isaac ierr = (*fem->ops->setfromoptions)(PetscOptionsObject,fem);CHKERRQ(ierr); 24020cf1dd8SToby Isaac } 24120cf1dd8SToby Isaac /* process any options handlers added with PetscObjectAddOptionsHandler() */ 24220cf1dd8SToby Isaac ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) fem);CHKERRQ(ierr); 24320cf1dd8SToby Isaac ierr = PetscOptionsEnd();CHKERRQ(ierr); 24420cf1dd8SToby Isaac ierr = PetscFEViewFromOptions(fem, NULL, "-petscfe_view");CHKERRQ(ierr); 24520cf1dd8SToby Isaac PetscFunctionReturn(0); 24620cf1dd8SToby Isaac } 24720cf1dd8SToby Isaac 24820cf1dd8SToby Isaac /*@C 24920cf1dd8SToby Isaac PetscFESetUp - Construct data structures for the PetscFE 25020cf1dd8SToby Isaac 25120cf1dd8SToby Isaac Collective on PetscFE 25220cf1dd8SToby Isaac 25320cf1dd8SToby Isaac Input Parameter: 25420cf1dd8SToby Isaac . fem - the PetscFE object to setup 25520cf1dd8SToby Isaac 25620cf1dd8SToby Isaac Level: developer 25720cf1dd8SToby Isaac 25820cf1dd8SToby Isaac .seealso PetscFEView(), PetscFEDestroy() 25920cf1dd8SToby Isaac @*/ 26020cf1dd8SToby Isaac PetscErrorCode PetscFESetUp(PetscFE fem) 26120cf1dd8SToby Isaac { 26220cf1dd8SToby Isaac PetscErrorCode ierr; 26320cf1dd8SToby Isaac 26420cf1dd8SToby Isaac PetscFunctionBegin; 26520cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 26620cf1dd8SToby Isaac if (fem->setupcalled) PetscFunctionReturn(0); 26720cf1dd8SToby Isaac fem->setupcalled = PETSC_TRUE; 26820cf1dd8SToby Isaac if (fem->ops->setup) {ierr = (*fem->ops->setup)(fem);CHKERRQ(ierr);} 26920cf1dd8SToby Isaac PetscFunctionReturn(0); 27020cf1dd8SToby Isaac } 27120cf1dd8SToby Isaac 27220cf1dd8SToby Isaac /*@ 27320cf1dd8SToby Isaac PetscFEDestroy - Destroys a PetscFE object 27420cf1dd8SToby Isaac 27520cf1dd8SToby Isaac Collective on PetscFE 27620cf1dd8SToby Isaac 27720cf1dd8SToby Isaac Input Parameter: 27820cf1dd8SToby Isaac . fem - the PetscFE object to destroy 27920cf1dd8SToby Isaac 28020cf1dd8SToby Isaac Level: developer 28120cf1dd8SToby Isaac 28220cf1dd8SToby Isaac .seealso PetscFEView() 28320cf1dd8SToby Isaac @*/ 28420cf1dd8SToby Isaac PetscErrorCode PetscFEDestroy(PetscFE *fem) 28520cf1dd8SToby Isaac { 28620cf1dd8SToby Isaac PetscErrorCode ierr; 28720cf1dd8SToby Isaac 28820cf1dd8SToby Isaac PetscFunctionBegin; 28920cf1dd8SToby Isaac if (!*fem) PetscFunctionReturn(0); 29020cf1dd8SToby Isaac PetscValidHeaderSpecific((*fem), PETSCFE_CLASSID, 1); 29120cf1dd8SToby Isaac 29220cf1dd8SToby Isaac if (--((PetscObject)(*fem))->refct > 0) {*fem = 0; PetscFunctionReturn(0);} 29320cf1dd8SToby Isaac ((PetscObject) (*fem))->refct = 0; 29420cf1dd8SToby Isaac 29520cf1dd8SToby Isaac if ((*fem)->subspaces) { 29620cf1dd8SToby Isaac PetscInt dim, d; 29720cf1dd8SToby Isaac 29820cf1dd8SToby Isaac ierr = PetscDualSpaceGetDimension((*fem)->dualSpace, &dim);CHKERRQ(ierr); 29920cf1dd8SToby Isaac for (d = 0; d < dim; ++d) {ierr = PetscFEDestroy(&(*fem)->subspaces[d]);CHKERRQ(ierr);} 30020cf1dd8SToby Isaac } 30120cf1dd8SToby Isaac ierr = PetscFree((*fem)->subspaces);CHKERRQ(ierr); 30220cf1dd8SToby Isaac ierr = PetscFree((*fem)->invV);CHKERRQ(ierr); 30320cf1dd8SToby Isaac ierr = PetscFERestoreTabulation((*fem), 0, NULL, &(*fem)->B, &(*fem)->D, NULL /*&(*fem)->H*/);CHKERRQ(ierr); 30420cf1dd8SToby Isaac ierr = PetscFERestoreTabulation((*fem), 0, NULL, &(*fem)->Bf, &(*fem)->Df, NULL /*&(*fem)->Hf*/);CHKERRQ(ierr); 30520cf1dd8SToby Isaac ierr = PetscFERestoreTabulation((*fem), 0, NULL, &(*fem)->F, NULL, NULL);CHKERRQ(ierr); 30620cf1dd8SToby Isaac ierr = PetscSpaceDestroy(&(*fem)->basisSpace);CHKERRQ(ierr); 30720cf1dd8SToby Isaac ierr = PetscDualSpaceDestroy(&(*fem)->dualSpace);CHKERRQ(ierr); 30820cf1dd8SToby Isaac ierr = PetscQuadratureDestroy(&(*fem)->quadrature);CHKERRQ(ierr); 30920cf1dd8SToby Isaac ierr = PetscQuadratureDestroy(&(*fem)->faceQuadrature);CHKERRQ(ierr); 31020cf1dd8SToby Isaac 31120cf1dd8SToby Isaac if ((*fem)->ops->destroy) {ierr = (*(*fem)->ops->destroy)(*fem);CHKERRQ(ierr);} 31220cf1dd8SToby Isaac ierr = PetscHeaderDestroy(fem);CHKERRQ(ierr); 31320cf1dd8SToby Isaac PetscFunctionReturn(0); 31420cf1dd8SToby Isaac } 31520cf1dd8SToby Isaac 31620cf1dd8SToby Isaac /*@ 31720cf1dd8SToby Isaac PetscFECreate - Creates an empty PetscFE object. The type can then be set with PetscFESetType(). 31820cf1dd8SToby Isaac 31920cf1dd8SToby Isaac Collective on MPI_Comm 32020cf1dd8SToby Isaac 32120cf1dd8SToby Isaac Input Parameter: 32220cf1dd8SToby Isaac . comm - The communicator for the PetscFE object 32320cf1dd8SToby Isaac 32420cf1dd8SToby Isaac Output Parameter: 32520cf1dd8SToby Isaac . fem - The PetscFE object 32620cf1dd8SToby Isaac 32720cf1dd8SToby Isaac Level: beginner 32820cf1dd8SToby Isaac 32920cf1dd8SToby Isaac .seealso: PetscFESetType(), PETSCFEGALERKIN 33020cf1dd8SToby Isaac @*/ 33120cf1dd8SToby Isaac PetscErrorCode PetscFECreate(MPI_Comm comm, PetscFE *fem) 33220cf1dd8SToby Isaac { 33320cf1dd8SToby Isaac PetscFE f; 33420cf1dd8SToby Isaac PetscErrorCode ierr; 33520cf1dd8SToby Isaac 33620cf1dd8SToby Isaac PetscFunctionBegin; 33720cf1dd8SToby Isaac PetscValidPointer(fem, 2); 33820cf1dd8SToby Isaac ierr = PetscCitationsRegister(FECitation,&FEcite);CHKERRQ(ierr); 33920cf1dd8SToby Isaac *fem = NULL; 34020cf1dd8SToby Isaac ierr = PetscFEInitializePackage();CHKERRQ(ierr); 34120cf1dd8SToby Isaac 34220cf1dd8SToby Isaac ierr = PetscHeaderCreate(f, PETSCFE_CLASSID, "PetscFE", "Finite Element", "PetscFE", comm, PetscFEDestroy, PetscFEView);CHKERRQ(ierr); 34320cf1dd8SToby Isaac 34420cf1dd8SToby Isaac f->basisSpace = NULL; 34520cf1dd8SToby Isaac f->dualSpace = NULL; 34620cf1dd8SToby Isaac f->numComponents = 1; 34720cf1dd8SToby Isaac f->subspaces = NULL; 34820cf1dd8SToby Isaac f->invV = NULL; 34920cf1dd8SToby Isaac f->B = NULL; 35020cf1dd8SToby Isaac f->D = NULL; 35120cf1dd8SToby Isaac f->H = NULL; 35220cf1dd8SToby Isaac f->Bf = NULL; 35320cf1dd8SToby Isaac f->Df = NULL; 35420cf1dd8SToby Isaac f->Hf = NULL; 35520cf1dd8SToby Isaac ierr = PetscMemzero(&f->quadrature, sizeof(PetscQuadrature));CHKERRQ(ierr); 35620cf1dd8SToby Isaac ierr = PetscMemzero(&f->faceQuadrature, sizeof(PetscQuadrature));CHKERRQ(ierr); 35720cf1dd8SToby Isaac f->blockSize = 0; 35820cf1dd8SToby Isaac f->numBlocks = 1; 35920cf1dd8SToby Isaac f->batchSize = 0; 36020cf1dd8SToby Isaac f->numBatches = 1; 36120cf1dd8SToby Isaac 36220cf1dd8SToby Isaac *fem = f; 36320cf1dd8SToby Isaac PetscFunctionReturn(0); 36420cf1dd8SToby Isaac } 36520cf1dd8SToby Isaac 36620cf1dd8SToby Isaac /*@ 36720cf1dd8SToby Isaac PetscFEGetSpatialDimension - Returns the spatial dimension of the element 36820cf1dd8SToby Isaac 36920cf1dd8SToby Isaac Not collective 37020cf1dd8SToby Isaac 37120cf1dd8SToby Isaac Input Parameter: 37220cf1dd8SToby Isaac . fem - The PetscFE object 37320cf1dd8SToby Isaac 37420cf1dd8SToby Isaac Output Parameter: 37520cf1dd8SToby Isaac . dim - The spatial dimension 37620cf1dd8SToby Isaac 37720cf1dd8SToby Isaac Level: intermediate 37820cf1dd8SToby Isaac 37920cf1dd8SToby Isaac .seealso: PetscFECreate() 38020cf1dd8SToby Isaac @*/ 38120cf1dd8SToby Isaac PetscErrorCode PetscFEGetSpatialDimension(PetscFE fem, PetscInt *dim) 38220cf1dd8SToby Isaac { 38320cf1dd8SToby Isaac DM dm; 38420cf1dd8SToby Isaac PetscErrorCode ierr; 38520cf1dd8SToby Isaac 38620cf1dd8SToby Isaac PetscFunctionBegin; 38720cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 38820cf1dd8SToby Isaac PetscValidPointer(dim, 2); 38920cf1dd8SToby Isaac ierr = PetscDualSpaceGetDM(fem->dualSpace, &dm);CHKERRQ(ierr); 39020cf1dd8SToby Isaac ierr = DMGetDimension(dm, dim);CHKERRQ(ierr); 39120cf1dd8SToby Isaac PetscFunctionReturn(0); 39220cf1dd8SToby Isaac } 39320cf1dd8SToby Isaac 39420cf1dd8SToby Isaac /*@ 39520cf1dd8SToby Isaac PetscFESetNumComponents - Sets the number of components in the element 39620cf1dd8SToby Isaac 39720cf1dd8SToby Isaac Not collective 39820cf1dd8SToby Isaac 39920cf1dd8SToby Isaac Input Parameters: 40020cf1dd8SToby Isaac + fem - The PetscFE object 40120cf1dd8SToby Isaac - comp - The number of field components 40220cf1dd8SToby Isaac 40320cf1dd8SToby Isaac Level: intermediate 40420cf1dd8SToby Isaac 40520cf1dd8SToby Isaac .seealso: PetscFECreate() 40620cf1dd8SToby Isaac @*/ 40720cf1dd8SToby Isaac PetscErrorCode PetscFESetNumComponents(PetscFE fem, PetscInt comp) 40820cf1dd8SToby Isaac { 40920cf1dd8SToby Isaac PetscFunctionBegin; 41020cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 41120cf1dd8SToby Isaac fem->numComponents = comp; 41220cf1dd8SToby Isaac PetscFunctionReturn(0); 41320cf1dd8SToby Isaac } 41420cf1dd8SToby Isaac 41520cf1dd8SToby Isaac /*@ 41620cf1dd8SToby Isaac PetscFEGetNumComponents - Returns the number of components in the element 41720cf1dd8SToby Isaac 41820cf1dd8SToby Isaac Not collective 41920cf1dd8SToby Isaac 42020cf1dd8SToby Isaac Input Parameter: 42120cf1dd8SToby Isaac . fem - The PetscFE object 42220cf1dd8SToby Isaac 42320cf1dd8SToby Isaac Output Parameter: 42420cf1dd8SToby Isaac . comp - The number of field components 42520cf1dd8SToby Isaac 42620cf1dd8SToby Isaac Level: intermediate 42720cf1dd8SToby Isaac 42820cf1dd8SToby Isaac .seealso: PetscFECreate() 42920cf1dd8SToby Isaac @*/ 43020cf1dd8SToby Isaac PetscErrorCode PetscFEGetNumComponents(PetscFE fem, PetscInt *comp) 43120cf1dd8SToby Isaac { 43220cf1dd8SToby Isaac PetscFunctionBegin; 43320cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 43420cf1dd8SToby Isaac PetscValidPointer(comp, 2); 43520cf1dd8SToby Isaac *comp = fem->numComponents; 43620cf1dd8SToby Isaac PetscFunctionReturn(0); 43720cf1dd8SToby Isaac } 43820cf1dd8SToby Isaac 43920cf1dd8SToby Isaac /*@ 44020cf1dd8SToby Isaac PetscFESetTileSizes - Sets the tile sizes for evaluation 44120cf1dd8SToby Isaac 44220cf1dd8SToby Isaac Not collective 44320cf1dd8SToby Isaac 44420cf1dd8SToby Isaac Input Parameters: 44520cf1dd8SToby Isaac + fem - The PetscFE object 44620cf1dd8SToby Isaac . blockSize - The number of elements in a block 44720cf1dd8SToby Isaac . numBlocks - The number of blocks in a batch 44820cf1dd8SToby Isaac . batchSize - The number of elements in a batch 44920cf1dd8SToby Isaac - numBatches - The number of batches in a chunk 45020cf1dd8SToby Isaac 45120cf1dd8SToby Isaac Level: intermediate 45220cf1dd8SToby Isaac 45320cf1dd8SToby Isaac .seealso: PetscFECreate() 45420cf1dd8SToby Isaac @*/ 45520cf1dd8SToby Isaac PetscErrorCode PetscFESetTileSizes(PetscFE fem, PetscInt blockSize, PetscInt numBlocks, PetscInt batchSize, PetscInt numBatches) 45620cf1dd8SToby Isaac { 45720cf1dd8SToby Isaac PetscFunctionBegin; 45820cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 45920cf1dd8SToby Isaac fem->blockSize = blockSize; 46020cf1dd8SToby Isaac fem->numBlocks = numBlocks; 46120cf1dd8SToby Isaac fem->batchSize = batchSize; 46220cf1dd8SToby Isaac fem->numBatches = numBatches; 46320cf1dd8SToby Isaac PetscFunctionReturn(0); 46420cf1dd8SToby Isaac } 46520cf1dd8SToby Isaac 46620cf1dd8SToby Isaac /*@ 46720cf1dd8SToby Isaac PetscFEGetTileSizes - Returns the tile sizes for evaluation 46820cf1dd8SToby Isaac 46920cf1dd8SToby Isaac Not collective 47020cf1dd8SToby Isaac 47120cf1dd8SToby Isaac Input Parameter: 47220cf1dd8SToby Isaac . fem - The PetscFE object 47320cf1dd8SToby Isaac 47420cf1dd8SToby Isaac Output Parameters: 47520cf1dd8SToby Isaac + blockSize - The number of elements in a block 47620cf1dd8SToby Isaac . numBlocks - The number of blocks in a batch 47720cf1dd8SToby Isaac . batchSize - The number of elements in a batch 47820cf1dd8SToby Isaac - numBatches - The number of batches in a chunk 47920cf1dd8SToby Isaac 48020cf1dd8SToby Isaac Level: intermediate 48120cf1dd8SToby Isaac 48220cf1dd8SToby Isaac .seealso: PetscFECreate() 48320cf1dd8SToby Isaac @*/ 48420cf1dd8SToby Isaac PetscErrorCode PetscFEGetTileSizes(PetscFE fem, PetscInt *blockSize, PetscInt *numBlocks, PetscInt *batchSize, PetscInt *numBatches) 48520cf1dd8SToby Isaac { 48620cf1dd8SToby Isaac PetscFunctionBegin; 48720cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 48820cf1dd8SToby Isaac if (blockSize) PetscValidPointer(blockSize, 2); 48920cf1dd8SToby Isaac if (numBlocks) PetscValidPointer(numBlocks, 3); 49020cf1dd8SToby Isaac if (batchSize) PetscValidPointer(batchSize, 4); 49120cf1dd8SToby Isaac if (numBatches) PetscValidPointer(numBatches, 5); 49220cf1dd8SToby Isaac if (blockSize) *blockSize = fem->blockSize; 49320cf1dd8SToby Isaac if (numBlocks) *numBlocks = fem->numBlocks; 49420cf1dd8SToby Isaac if (batchSize) *batchSize = fem->batchSize; 49520cf1dd8SToby Isaac if (numBatches) *numBatches = fem->numBatches; 49620cf1dd8SToby Isaac PetscFunctionReturn(0); 49720cf1dd8SToby Isaac } 49820cf1dd8SToby Isaac 49920cf1dd8SToby Isaac /*@ 50020cf1dd8SToby Isaac PetscFEGetBasisSpace - Returns the PetscSpace used for approximation of the solution 50120cf1dd8SToby Isaac 50220cf1dd8SToby Isaac Not collective 50320cf1dd8SToby Isaac 50420cf1dd8SToby Isaac Input Parameter: 50520cf1dd8SToby Isaac . fem - The PetscFE object 50620cf1dd8SToby Isaac 50720cf1dd8SToby Isaac Output Parameter: 50820cf1dd8SToby Isaac . sp - The PetscSpace object 50920cf1dd8SToby Isaac 51020cf1dd8SToby Isaac Level: intermediate 51120cf1dd8SToby Isaac 51220cf1dd8SToby Isaac .seealso: PetscFECreate() 51320cf1dd8SToby Isaac @*/ 51420cf1dd8SToby Isaac PetscErrorCode PetscFEGetBasisSpace(PetscFE fem, PetscSpace *sp) 51520cf1dd8SToby Isaac { 51620cf1dd8SToby Isaac PetscFunctionBegin; 51720cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 51820cf1dd8SToby Isaac PetscValidPointer(sp, 2); 51920cf1dd8SToby Isaac *sp = fem->basisSpace; 52020cf1dd8SToby Isaac PetscFunctionReturn(0); 52120cf1dd8SToby Isaac } 52220cf1dd8SToby Isaac 52320cf1dd8SToby Isaac /*@ 52420cf1dd8SToby Isaac PetscFESetBasisSpace - Sets the PetscSpace used for approximation of the solution 52520cf1dd8SToby Isaac 52620cf1dd8SToby Isaac Not collective 52720cf1dd8SToby Isaac 52820cf1dd8SToby Isaac Input Parameters: 52920cf1dd8SToby Isaac + fem - The PetscFE object 53020cf1dd8SToby Isaac - sp - The PetscSpace object 53120cf1dd8SToby Isaac 53220cf1dd8SToby Isaac Level: intermediate 53320cf1dd8SToby Isaac 53420cf1dd8SToby Isaac .seealso: PetscFECreate() 53520cf1dd8SToby Isaac @*/ 53620cf1dd8SToby Isaac PetscErrorCode PetscFESetBasisSpace(PetscFE fem, PetscSpace sp) 53720cf1dd8SToby Isaac { 53820cf1dd8SToby Isaac PetscErrorCode ierr; 53920cf1dd8SToby Isaac 54020cf1dd8SToby Isaac PetscFunctionBegin; 54120cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 54220cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 2); 54320cf1dd8SToby Isaac ierr = PetscSpaceDestroy(&fem->basisSpace);CHKERRQ(ierr); 54420cf1dd8SToby Isaac fem->basisSpace = sp; 54520cf1dd8SToby Isaac ierr = PetscObjectReference((PetscObject) fem->basisSpace);CHKERRQ(ierr); 54620cf1dd8SToby Isaac PetscFunctionReturn(0); 54720cf1dd8SToby Isaac } 54820cf1dd8SToby Isaac 54920cf1dd8SToby Isaac /*@ 55020cf1dd8SToby Isaac PetscFEGetDualSpace - Returns the PetscDualSpace used to define the inner product 55120cf1dd8SToby Isaac 55220cf1dd8SToby Isaac Not collective 55320cf1dd8SToby Isaac 55420cf1dd8SToby Isaac Input Parameter: 55520cf1dd8SToby Isaac . fem - The PetscFE object 55620cf1dd8SToby Isaac 55720cf1dd8SToby Isaac Output Parameter: 55820cf1dd8SToby Isaac . sp - The PetscDualSpace object 55920cf1dd8SToby Isaac 56020cf1dd8SToby Isaac Level: intermediate 56120cf1dd8SToby Isaac 56220cf1dd8SToby Isaac .seealso: PetscFECreate() 56320cf1dd8SToby Isaac @*/ 56420cf1dd8SToby Isaac PetscErrorCode PetscFEGetDualSpace(PetscFE fem, PetscDualSpace *sp) 56520cf1dd8SToby Isaac { 56620cf1dd8SToby Isaac PetscFunctionBegin; 56720cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 56820cf1dd8SToby Isaac PetscValidPointer(sp, 2); 56920cf1dd8SToby Isaac *sp = fem->dualSpace; 57020cf1dd8SToby Isaac PetscFunctionReturn(0); 57120cf1dd8SToby Isaac } 57220cf1dd8SToby Isaac 57320cf1dd8SToby Isaac /*@ 57420cf1dd8SToby Isaac PetscFESetDualSpace - Sets the PetscDualSpace used to define the inner product 57520cf1dd8SToby Isaac 57620cf1dd8SToby Isaac Not collective 57720cf1dd8SToby Isaac 57820cf1dd8SToby Isaac Input Parameters: 57920cf1dd8SToby Isaac + fem - The PetscFE object 58020cf1dd8SToby Isaac - sp - The PetscDualSpace object 58120cf1dd8SToby Isaac 58220cf1dd8SToby Isaac Level: intermediate 58320cf1dd8SToby Isaac 58420cf1dd8SToby Isaac .seealso: PetscFECreate() 58520cf1dd8SToby Isaac @*/ 58620cf1dd8SToby Isaac PetscErrorCode PetscFESetDualSpace(PetscFE fem, PetscDualSpace sp) 58720cf1dd8SToby Isaac { 58820cf1dd8SToby Isaac PetscErrorCode ierr; 58920cf1dd8SToby Isaac 59020cf1dd8SToby Isaac PetscFunctionBegin; 59120cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 59220cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 2); 59320cf1dd8SToby Isaac ierr = PetscDualSpaceDestroy(&fem->dualSpace);CHKERRQ(ierr); 59420cf1dd8SToby Isaac fem->dualSpace = sp; 59520cf1dd8SToby Isaac ierr = PetscObjectReference((PetscObject) fem->dualSpace);CHKERRQ(ierr); 59620cf1dd8SToby Isaac PetscFunctionReturn(0); 59720cf1dd8SToby Isaac } 59820cf1dd8SToby Isaac 59920cf1dd8SToby Isaac /*@ 60020cf1dd8SToby Isaac PetscFEGetQuadrature - Returns the PetscQuadrature used to calculate inner products 60120cf1dd8SToby Isaac 60220cf1dd8SToby Isaac Not collective 60320cf1dd8SToby Isaac 60420cf1dd8SToby Isaac Input Parameter: 60520cf1dd8SToby Isaac . fem - The PetscFE object 60620cf1dd8SToby Isaac 60720cf1dd8SToby Isaac Output Parameter: 60820cf1dd8SToby Isaac . q - The PetscQuadrature object 60920cf1dd8SToby Isaac 61020cf1dd8SToby Isaac Level: intermediate 61120cf1dd8SToby Isaac 61220cf1dd8SToby Isaac .seealso: PetscFECreate() 61320cf1dd8SToby Isaac @*/ 61420cf1dd8SToby Isaac PetscErrorCode PetscFEGetQuadrature(PetscFE fem, PetscQuadrature *q) 61520cf1dd8SToby Isaac { 61620cf1dd8SToby Isaac PetscFunctionBegin; 61720cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 61820cf1dd8SToby Isaac PetscValidPointer(q, 2); 61920cf1dd8SToby Isaac *q = fem->quadrature; 62020cf1dd8SToby Isaac PetscFunctionReturn(0); 62120cf1dd8SToby Isaac } 62220cf1dd8SToby Isaac 62320cf1dd8SToby Isaac /*@ 62420cf1dd8SToby Isaac PetscFESetQuadrature - Sets the PetscQuadrature used to calculate inner products 62520cf1dd8SToby Isaac 62620cf1dd8SToby Isaac Not collective 62720cf1dd8SToby Isaac 62820cf1dd8SToby Isaac Input Parameters: 62920cf1dd8SToby Isaac + fem - The PetscFE object 63020cf1dd8SToby Isaac - q - The PetscQuadrature object 63120cf1dd8SToby Isaac 63220cf1dd8SToby Isaac Level: intermediate 63320cf1dd8SToby Isaac 63420cf1dd8SToby Isaac .seealso: PetscFECreate() 63520cf1dd8SToby Isaac @*/ 63620cf1dd8SToby Isaac PetscErrorCode PetscFESetQuadrature(PetscFE fem, PetscQuadrature q) 63720cf1dd8SToby Isaac { 63820cf1dd8SToby Isaac PetscInt Nc, qNc; 63920cf1dd8SToby Isaac PetscErrorCode ierr; 64020cf1dd8SToby Isaac 64120cf1dd8SToby Isaac PetscFunctionBegin; 64220cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 64320cf1dd8SToby Isaac ierr = PetscFEGetNumComponents(fem, &Nc);CHKERRQ(ierr); 64420cf1dd8SToby Isaac ierr = PetscQuadratureGetNumComponents(q, &qNc);CHKERRQ(ierr); 64520cf1dd8SToby 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); 64620cf1dd8SToby Isaac ierr = PetscFERestoreTabulation(fem, 0, NULL, &fem->B, &fem->D, NULL /*&(*fem)->H*/);CHKERRQ(ierr); 64720cf1dd8SToby Isaac ierr = PetscQuadratureDestroy(&fem->quadrature);CHKERRQ(ierr); 64820cf1dd8SToby Isaac fem->quadrature = q; 64920cf1dd8SToby Isaac ierr = PetscObjectReference((PetscObject) q);CHKERRQ(ierr); 65020cf1dd8SToby Isaac PetscFunctionReturn(0); 65120cf1dd8SToby Isaac } 65220cf1dd8SToby Isaac 65320cf1dd8SToby Isaac /*@ 65420cf1dd8SToby Isaac PetscFEGetFaceQuadrature - Returns the PetscQuadrature used to calculate inner products on faces 65520cf1dd8SToby Isaac 65620cf1dd8SToby Isaac Not collective 65720cf1dd8SToby Isaac 65820cf1dd8SToby Isaac Input Parameter: 65920cf1dd8SToby Isaac . fem - The PetscFE object 66020cf1dd8SToby Isaac 66120cf1dd8SToby Isaac Output Parameter: 66220cf1dd8SToby Isaac . q - The PetscQuadrature object 66320cf1dd8SToby Isaac 66420cf1dd8SToby Isaac Level: intermediate 66520cf1dd8SToby Isaac 66620cf1dd8SToby Isaac .seealso: PetscFECreate() 66720cf1dd8SToby Isaac @*/ 66820cf1dd8SToby Isaac PetscErrorCode PetscFEGetFaceQuadrature(PetscFE fem, PetscQuadrature *q) 66920cf1dd8SToby Isaac { 67020cf1dd8SToby Isaac PetscFunctionBegin; 67120cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 67220cf1dd8SToby Isaac PetscValidPointer(q, 2); 67320cf1dd8SToby Isaac *q = fem->faceQuadrature; 67420cf1dd8SToby Isaac PetscFunctionReturn(0); 67520cf1dd8SToby Isaac } 67620cf1dd8SToby Isaac 67720cf1dd8SToby Isaac /*@ 67820cf1dd8SToby Isaac PetscFESetFaceQuadrature - Sets the PetscQuadrature used to calculate inner products on faces 67920cf1dd8SToby Isaac 68020cf1dd8SToby Isaac Not collective 68120cf1dd8SToby Isaac 68220cf1dd8SToby Isaac Input Parameters: 68320cf1dd8SToby Isaac + fem - The PetscFE object 68420cf1dd8SToby Isaac - q - The PetscQuadrature object 68520cf1dd8SToby Isaac 68620cf1dd8SToby Isaac Level: intermediate 68720cf1dd8SToby Isaac 68820cf1dd8SToby Isaac .seealso: PetscFECreate() 68920cf1dd8SToby Isaac @*/ 69020cf1dd8SToby Isaac PetscErrorCode PetscFESetFaceQuadrature(PetscFE fem, PetscQuadrature q) 69120cf1dd8SToby Isaac { 69220cf1dd8SToby Isaac PetscErrorCode ierr; 69320cf1dd8SToby Isaac 69420cf1dd8SToby Isaac PetscFunctionBegin; 69520cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 69620cf1dd8SToby Isaac ierr = PetscFERestoreTabulation(fem, 0, NULL, &fem->Bf, &fem->Df, NULL /*&(*fem)->Hf*/);CHKERRQ(ierr); 69720cf1dd8SToby Isaac ierr = PetscQuadratureDestroy(&fem->faceQuadrature);CHKERRQ(ierr); 69820cf1dd8SToby Isaac fem->faceQuadrature = q; 69920cf1dd8SToby Isaac ierr = PetscObjectReference((PetscObject) q);CHKERRQ(ierr); 70020cf1dd8SToby Isaac PetscFunctionReturn(0); 70120cf1dd8SToby Isaac } 70220cf1dd8SToby Isaac 70320cf1dd8SToby Isaac /*@C 70420cf1dd8SToby Isaac PetscFEGetNumDof - Returns the number of dofs (dual basis vectors) associated to mesh points on the reference cell of a given dimension 70520cf1dd8SToby Isaac 70620cf1dd8SToby Isaac Not collective 70720cf1dd8SToby Isaac 70820cf1dd8SToby Isaac Input Parameter: 70920cf1dd8SToby Isaac . fem - The PetscFE object 71020cf1dd8SToby Isaac 71120cf1dd8SToby Isaac Output Parameter: 71220cf1dd8SToby Isaac . numDof - Array with the number of dofs per dimension 71320cf1dd8SToby Isaac 71420cf1dd8SToby Isaac Level: intermediate 71520cf1dd8SToby Isaac 71620cf1dd8SToby Isaac .seealso: PetscFECreate() 71720cf1dd8SToby Isaac @*/ 71820cf1dd8SToby Isaac PetscErrorCode PetscFEGetNumDof(PetscFE fem, const PetscInt **numDof) 71920cf1dd8SToby Isaac { 72020cf1dd8SToby Isaac PetscErrorCode ierr; 72120cf1dd8SToby Isaac 72220cf1dd8SToby Isaac PetscFunctionBegin; 72320cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 72420cf1dd8SToby Isaac PetscValidPointer(numDof, 2); 72520cf1dd8SToby Isaac ierr = PetscDualSpaceGetNumDof(fem->dualSpace, numDof);CHKERRQ(ierr); 72620cf1dd8SToby Isaac PetscFunctionReturn(0); 72720cf1dd8SToby Isaac } 72820cf1dd8SToby Isaac 72920cf1dd8SToby Isaac /*@C 73020cf1dd8SToby Isaac PetscFEGetDefaultTabulation - Returns the tabulation of the basis functions at the quadrature points 73120cf1dd8SToby Isaac 73220cf1dd8SToby Isaac Not collective 73320cf1dd8SToby Isaac 73420cf1dd8SToby Isaac Input Parameter: 73520cf1dd8SToby Isaac . fem - The PetscFE object 73620cf1dd8SToby Isaac 73720cf1dd8SToby Isaac Output Parameters: 73820cf1dd8SToby Isaac + B - The basis function values at quadrature points 73920cf1dd8SToby Isaac . D - The basis function derivatives at quadrature points 74020cf1dd8SToby Isaac - H - The basis function second derivatives at quadrature points 74120cf1dd8SToby Isaac 74220cf1dd8SToby Isaac Note: 74320cf1dd8SToby Isaac $ B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c 74420cf1dd8SToby 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 74520cf1dd8SToby 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 74620cf1dd8SToby Isaac 74720cf1dd8SToby Isaac Level: intermediate 74820cf1dd8SToby Isaac 74920cf1dd8SToby Isaac .seealso: PetscFEGetTabulation(), PetscFERestoreTabulation() 75020cf1dd8SToby Isaac @*/ 75120cf1dd8SToby Isaac PetscErrorCode PetscFEGetDefaultTabulation(PetscFE fem, PetscReal **B, PetscReal **D, PetscReal **H) 75220cf1dd8SToby Isaac { 75320cf1dd8SToby Isaac PetscInt npoints; 75420cf1dd8SToby Isaac const PetscReal *points; 75520cf1dd8SToby Isaac PetscErrorCode ierr; 75620cf1dd8SToby Isaac 75720cf1dd8SToby Isaac PetscFunctionBegin; 75820cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 75920cf1dd8SToby Isaac if (B) PetscValidPointer(B, 2); 76020cf1dd8SToby Isaac if (D) PetscValidPointer(D, 3); 76120cf1dd8SToby Isaac if (H) PetscValidPointer(H, 4); 76220cf1dd8SToby Isaac ierr = PetscQuadratureGetData(fem->quadrature, NULL, NULL, &npoints, &points, NULL);CHKERRQ(ierr); 76320cf1dd8SToby Isaac if (!fem->B) {ierr = PetscFEGetTabulation(fem, npoints, points, &fem->B, &fem->D, NULL/*&fem->H*/);CHKERRQ(ierr);} 76420cf1dd8SToby Isaac if (B) *B = fem->B; 76520cf1dd8SToby Isaac if (D) *D = fem->D; 76620cf1dd8SToby Isaac if (H) *H = fem->H; 76720cf1dd8SToby Isaac PetscFunctionReturn(0); 76820cf1dd8SToby Isaac } 76920cf1dd8SToby Isaac 77020cf1dd8SToby Isaac PetscErrorCode PetscFEGetFaceTabulation(PetscFE fem, PetscReal **Bf, PetscReal **Df, PetscReal **Hf) 77120cf1dd8SToby Isaac { 77220cf1dd8SToby Isaac PetscErrorCode ierr; 77320cf1dd8SToby Isaac 77420cf1dd8SToby Isaac PetscFunctionBegin; 77520cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 77620cf1dd8SToby Isaac if (Bf) PetscValidPointer(Bf, 2); 77720cf1dd8SToby Isaac if (Df) PetscValidPointer(Df, 3); 77820cf1dd8SToby Isaac if (Hf) PetscValidPointer(Hf, 4); 77920cf1dd8SToby Isaac if (!fem->Bf) { 78020cf1dd8SToby Isaac const PetscReal xi0[3] = {-1., -1., -1.}; 78120cf1dd8SToby Isaac PetscReal v0[3], J[9], detJ; 78220cf1dd8SToby Isaac PetscQuadrature fq; 78320cf1dd8SToby Isaac PetscDualSpace sp; 78420cf1dd8SToby Isaac DM dm; 78520cf1dd8SToby Isaac const PetscInt *faces; 78620cf1dd8SToby Isaac PetscInt dim, numFaces, f, npoints, q; 78720cf1dd8SToby Isaac const PetscReal *points; 78820cf1dd8SToby Isaac PetscReal *facePoints; 78920cf1dd8SToby Isaac 79020cf1dd8SToby Isaac ierr = PetscFEGetDualSpace(fem, &sp);CHKERRQ(ierr); 79120cf1dd8SToby Isaac ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr); 79220cf1dd8SToby Isaac ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 79320cf1dd8SToby Isaac ierr = DMPlexGetConeSize(dm, 0, &numFaces);CHKERRQ(ierr); 79420cf1dd8SToby Isaac ierr = DMPlexGetCone(dm, 0, &faces);CHKERRQ(ierr); 79520cf1dd8SToby Isaac ierr = PetscFEGetFaceQuadrature(fem, &fq);CHKERRQ(ierr); 79620cf1dd8SToby Isaac if (fq) { 79720cf1dd8SToby Isaac ierr = PetscQuadratureGetData(fq, NULL, NULL, &npoints, &points, NULL);CHKERRQ(ierr); 79820cf1dd8SToby Isaac ierr = PetscMalloc1(numFaces*npoints*dim, &facePoints);CHKERRQ(ierr); 79920cf1dd8SToby Isaac for (f = 0; f < numFaces; ++f) { 80020cf1dd8SToby Isaac ierr = DMPlexComputeCellGeometryFEM(dm, faces[f], NULL, v0, J, NULL, &detJ);CHKERRQ(ierr); 80120cf1dd8SToby Isaac for (q = 0; q < npoints; ++q) CoordinatesRefToReal(dim, dim-1, xi0, v0, J, &points[q*(dim-1)], &facePoints[(f*npoints+q)*dim]); 80220cf1dd8SToby Isaac } 80320cf1dd8SToby Isaac ierr = PetscFEGetTabulation(fem, numFaces*npoints, facePoints, &fem->Bf, &fem->Df, NULL/*&fem->Hf*/);CHKERRQ(ierr); 80420cf1dd8SToby Isaac ierr = PetscFree(facePoints);CHKERRQ(ierr); 80520cf1dd8SToby Isaac } 80620cf1dd8SToby Isaac } 80720cf1dd8SToby Isaac if (Bf) *Bf = fem->Bf; 80820cf1dd8SToby Isaac if (Df) *Df = fem->Df; 80920cf1dd8SToby Isaac if (Hf) *Hf = fem->Hf; 81020cf1dd8SToby Isaac PetscFunctionReturn(0); 81120cf1dd8SToby Isaac } 81220cf1dd8SToby Isaac 81320cf1dd8SToby Isaac PetscErrorCode PetscFEGetFaceCentroidTabulation(PetscFE fem, PetscReal **F) 81420cf1dd8SToby Isaac { 81520cf1dd8SToby Isaac PetscErrorCode ierr; 81620cf1dd8SToby Isaac 81720cf1dd8SToby Isaac PetscFunctionBegin; 81820cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 81920cf1dd8SToby Isaac PetscValidPointer(F, 2); 82020cf1dd8SToby Isaac if (!fem->F) { 82120cf1dd8SToby Isaac PetscDualSpace sp; 82220cf1dd8SToby Isaac DM dm; 82320cf1dd8SToby Isaac const PetscInt *cone; 82420cf1dd8SToby Isaac PetscReal *centroids; 82520cf1dd8SToby Isaac PetscInt dim, numFaces, f; 82620cf1dd8SToby Isaac 82720cf1dd8SToby Isaac ierr = PetscFEGetDualSpace(fem, &sp);CHKERRQ(ierr); 82820cf1dd8SToby Isaac ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr); 82920cf1dd8SToby Isaac ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 83020cf1dd8SToby Isaac ierr = DMPlexGetConeSize(dm, 0, &numFaces);CHKERRQ(ierr); 83120cf1dd8SToby Isaac ierr = DMPlexGetCone(dm, 0, &cone);CHKERRQ(ierr); 83220cf1dd8SToby Isaac ierr = PetscMalloc1(numFaces*dim, ¢roids);CHKERRQ(ierr); 83320cf1dd8SToby Isaac for (f = 0; f < numFaces; ++f) {ierr = DMPlexComputeCellGeometryFVM(dm, cone[f], NULL, ¢roids[f*dim], NULL);CHKERRQ(ierr);} 83420cf1dd8SToby Isaac ierr = PetscFEGetTabulation(fem, numFaces, centroids, &fem->F, NULL, NULL);CHKERRQ(ierr); 83520cf1dd8SToby Isaac ierr = PetscFree(centroids);CHKERRQ(ierr); 83620cf1dd8SToby Isaac } 83720cf1dd8SToby Isaac *F = fem->F; 83820cf1dd8SToby Isaac PetscFunctionReturn(0); 83920cf1dd8SToby Isaac } 84020cf1dd8SToby Isaac 84120cf1dd8SToby Isaac /*@C 84220cf1dd8SToby Isaac PetscFEGetTabulation - Tabulates the basis functions, and perhaps derivatives, at the points provided. 84320cf1dd8SToby Isaac 84420cf1dd8SToby Isaac Not collective 84520cf1dd8SToby Isaac 84620cf1dd8SToby Isaac Input Parameters: 84720cf1dd8SToby Isaac + fem - The PetscFE object 84820cf1dd8SToby Isaac . npoints - The number of tabulation points 84920cf1dd8SToby Isaac - points - The tabulation point coordinates 85020cf1dd8SToby Isaac 85120cf1dd8SToby Isaac Output Parameters: 85220cf1dd8SToby Isaac + B - The basis function values at tabulation points 85320cf1dd8SToby Isaac . D - The basis function derivatives at tabulation points 85420cf1dd8SToby Isaac - H - The basis function second derivatives at tabulation points 85520cf1dd8SToby Isaac 85620cf1dd8SToby Isaac Note: 85720cf1dd8SToby Isaac $ B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c 85820cf1dd8SToby 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 85920cf1dd8SToby 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 86020cf1dd8SToby Isaac 86120cf1dd8SToby Isaac Level: intermediate 86220cf1dd8SToby Isaac 86320cf1dd8SToby Isaac .seealso: PetscFERestoreTabulation(), PetscFEGetDefaultTabulation() 86420cf1dd8SToby Isaac @*/ 86520cf1dd8SToby Isaac PetscErrorCode PetscFEGetTabulation(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscReal **B, PetscReal **D, PetscReal **H) 86620cf1dd8SToby Isaac { 86720cf1dd8SToby Isaac DM dm; 86820cf1dd8SToby Isaac PetscInt pdim; /* Dimension of FE space P */ 86920cf1dd8SToby Isaac PetscInt dim; /* Spatial dimension */ 87020cf1dd8SToby Isaac PetscInt comp; /* Field components */ 87120cf1dd8SToby Isaac PetscErrorCode ierr; 87220cf1dd8SToby Isaac 87320cf1dd8SToby Isaac PetscFunctionBegin; 87420cf1dd8SToby Isaac if (!npoints) { 87520cf1dd8SToby Isaac if (B) *B = NULL; 87620cf1dd8SToby Isaac if (D) *D = NULL; 87720cf1dd8SToby Isaac if (H) *H = NULL; 87820cf1dd8SToby Isaac PetscFunctionReturn(0); 87920cf1dd8SToby Isaac } 88020cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 88120cf1dd8SToby Isaac PetscValidPointer(points, 3); 88220cf1dd8SToby Isaac if (B) PetscValidPointer(B, 4); 88320cf1dd8SToby Isaac if (D) PetscValidPointer(D, 5); 88420cf1dd8SToby Isaac if (H) PetscValidPointer(H, 6); 88520cf1dd8SToby Isaac ierr = PetscDualSpaceGetDM(fem->dualSpace, &dm);CHKERRQ(ierr); 88620cf1dd8SToby Isaac ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 88720cf1dd8SToby Isaac ierr = PetscDualSpaceGetDimension(fem->dualSpace, &pdim);CHKERRQ(ierr); 88820cf1dd8SToby Isaac ierr = PetscFEGetNumComponents(fem, &comp);CHKERRQ(ierr); 88920cf1dd8SToby Isaac if (B) {ierr = DMGetWorkArray(dm, npoints*pdim*comp, MPIU_REAL, B);CHKERRQ(ierr);} 89020cf1dd8SToby Isaac if (!dim) { 89120cf1dd8SToby Isaac if (D) *D = NULL; 89220cf1dd8SToby Isaac if (H) *H = NULL; 89320cf1dd8SToby Isaac } else { 89420cf1dd8SToby Isaac if (D) {ierr = DMGetWorkArray(dm, npoints*pdim*comp*dim, MPIU_REAL, D);CHKERRQ(ierr);} 89520cf1dd8SToby Isaac if (H) {ierr = DMGetWorkArray(dm, npoints*pdim*comp*dim*dim, MPIU_REAL, H);CHKERRQ(ierr);} 89620cf1dd8SToby Isaac } 89720cf1dd8SToby Isaac ierr = (*fem->ops->gettabulation)(fem, npoints, points, B ? *B : NULL, D ? *D : NULL, H ? *H : NULL);CHKERRQ(ierr); 89820cf1dd8SToby Isaac PetscFunctionReturn(0); 89920cf1dd8SToby Isaac } 90020cf1dd8SToby Isaac 90120cf1dd8SToby Isaac PetscErrorCode PetscFERestoreTabulation(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscReal **B, PetscReal **D, PetscReal **H) 90220cf1dd8SToby Isaac { 90320cf1dd8SToby Isaac DM dm; 90420cf1dd8SToby Isaac PetscErrorCode ierr; 90520cf1dd8SToby Isaac 90620cf1dd8SToby Isaac PetscFunctionBegin; 90720cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 90820cf1dd8SToby Isaac ierr = PetscDualSpaceGetDM(fem->dualSpace, &dm);CHKERRQ(ierr); 90920cf1dd8SToby Isaac if (B && *B) {ierr = DMRestoreWorkArray(dm, 0, MPIU_REAL, B);CHKERRQ(ierr);} 91020cf1dd8SToby Isaac if (D && *D) {ierr = DMRestoreWorkArray(dm, 0, MPIU_REAL, D);CHKERRQ(ierr);} 91120cf1dd8SToby Isaac if (H && *H) {ierr = DMRestoreWorkArray(dm, 0, MPIU_REAL, H);CHKERRQ(ierr);} 91220cf1dd8SToby Isaac PetscFunctionReturn(0); 91320cf1dd8SToby Isaac } 91420cf1dd8SToby Isaac 91520cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreatePointTrace(PetscFE fe, PetscInt refPoint, PetscFE *trFE) 91620cf1dd8SToby Isaac { 91720cf1dd8SToby Isaac PetscSpace bsp, bsubsp; 91820cf1dd8SToby Isaac PetscDualSpace dsp, dsubsp; 91920cf1dd8SToby Isaac PetscInt dim, depth, numComp, i, j, coneSize, order; 92020cf1dd8SToby Isaac PetscFEType type; 92120cf1dd8SToby Isaac DM dm; 92220cf1dd8SToby Isaac DMLabel label; 92320cf1dd8SToby Isaac PetscReal *xi, *v, *J, detJ; 924db11e2ebSMatthew G. Knepley const char *name; 92520cf1dd8SToby Isaac PetscQuadrature origin, fullQuad, subQuad; 92620cf1dd8SToby Isaac PetscErrorCode ierr; 92720cf1dd8SToby Isaac 92820cf1dd8SToby Isaac PetscFunctionBegin; 92920cf1dd8SToby Isaac PetscValidHeaderSpecific(fe,PETSCFE_CLASSID,1); 93020cf1dd8SToby Isaac PetscValidPointer(trFE,3); 93120cf1dd8SToby Isaac ierr = PetscFEGetBasisSpace(fe,&bsp);CHKERRQ(ierr); 93220cf1dd8SToby Isaac ierr = PetscFEGetDualSpace(fe,&dsp);CHKERRQ(ierr); 93320cf1dd8SToby Isaac ierr = PetscDualSpaceGetDM(dsp,&dm);CHKERRQ(ierr); 93420cf1dd8SToby Isaac ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 93520cf1dd8SToby Isaac ierr = DMPlexGetDepthLabel(dm,&label);CHKERRQ(ierr); 93620cf1dd8SToby Isaac ierr = DMLabelGetValue(label,refPoint,&depth);CHKERRQ(ierr); 93720cf1dd8SToby Isaac ierr = PetscCalloc1(depth,&xi);CHKERRQ(ierr); 93820cf1dd8SToby Isaac ierr = PetscMalloc1(dim,&v);CHKERRQ(ierr); 93920cf1dd8SToby Isaac ierr = PetscMalloc1(dim*dim,&J);CHKERRQ(ierr); 94020cf1dd8SToby Isaac for (i = 0; i < depth; i++) xi[i] = 0.; 94120cf1dd8SToby Isaac ierr = PetscQuadratureCreate(PETSC_COMM_SELF,&origin);CHKERRQ(ierr); 94220cf1dd8SToby Isaac ierr = PetscQuadratureSetData(origin,depth,0,1,xi,NULL);CHKERRQ(ierr); 94320cf1dd8SToby Isaac ierr = DMPlexComputeCellGeometryFEM(dm,refPoint,origin,v,J,NULL,&detJ);CHKERRQ(ierr); 94420cf1dd8SToby Isaac /* CellGeometryFEM computes the expanded Jacobian, we want the true jacobian */ 94520cf1dd8SToby Isaac for (i = 1; i < dim; i++) { 94620cf1dd8SToby Isaac for (j = 0; j < depth; j++) { 94720cf1dd8SToby Isaac J[i * depth + j] = J[i * dim + j]; 94820cf1dd8SToby Isaac } 94920cf1dd8SToby Isaac } 95020cf1dd8SToby Isaac ierr = PetscQuadratureDestroy(&origin);CHKERRQ(ierr); 95120cf1dd8SToby Isaac ierr = PetscDualSpaceGetPointSubspace(dsp,refPoint,&dsubsp);CHKERRQ(ierr); 95220cf1dd8SToby Isaac ierr = PetscSpaceCreateSubspace(bsp,dsubsp,v,J,NULL,NULL,PETSC_OWN_POINTER,&bsubsp);CHKERRQ(ierr); 95320cf1dd8SToby Isaac ierr = PetscSpaceSetUp(bsubsp);CHKERRQ(ierr); 95420cf1dd8SToby Isaac ierr = PetscFECreate(PetscObjectComm((PetscObject)fe),trFE);CHKERRQ(ierr); 95520cf1dd8SToby Isaac ierr = PetscFEGetType(fe,&type);CHKERRQ(ierr); 95620cf1dd8SToby Isaac ierr = PetscFESetType(*trFE,type);CHKERRQ(ierr); 95720cf1dd8SToby Isaac ierr = PetscFEGetNumComponents(fe,&numComp);CHKERRQ(ierr); 95820cf1dd8SToby Isaac ierr = PetscFESetNumComponents(*trFE,numComp);CHKERRQ(ierr); 95920cf1dd8SToby Isaac ierr = PetscFESetBasisSpace(*trFE,bsubsp);CHKERRQ(ierr); 96020cf1dd8SToby Isaac ierr = PetscFESetDualSpace(*trFE,dsubsp);CHKERRQ(ierr); 961db11e2ebSMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) fe, &name);CHKERRQ(ierr); 962db11e2ebSMatthew G. Knepley if (name) {ierr = PetscFESetName(*trFE, name);CHKERRQ(ierr);} 96320cf1dd8SToby Isaac ierr = PetscFEGetQuadrature(fe,&fullQuad);CHKERRQ(ierr); 96420cf1dd8SToby Isaac ierr = PetscQuadratureGetOrder(fullQuad,&order);CHKERRQ(ierr); 96520cf1dd8SToby Isaac ierr = DMPlexGetConeSize(dm,refPoint,&coneSize);CHKERRQ(ierr); 96620cf1dd8SToby Isaac if (coneSize == 2 * depth) { 96720cf1dd8SToby Isaac ierr = PetscDTGaussTensorQuadrature(depth,1,(order + 1)/2,-1.,1.,&subQuad);CHKERRQ(ierr); 96820cf1dd8SToby Isaac } else { 96920cf1dd8SToby Isaac ierr = PetscDTGaussJacobiQuadrature(depth,1,(order + 1)/2,-1.,1.,&subQuad);CHKERRQ(ierr); 97020cf1dd8SToby Isaac } 97120cf1dd8SToby Isaac ierr = PetscFESetQuadrature(*trFE,subQuad);CHKERRQ(ierr); 97220cf1dd8SToby Isaac ierr = PetscFESetUp(*trFE);CHKERRQ(ierr); 97320cf1dd8SToby Isaac ierr = PetscQuadratureDestroy(&subQuad);CHKERRQ(ierr); 97420cf1dd8SToby Isaac ierr = PetscSpaceDestroy(&bsubsp);CHKERRQ(ierr); 97520cf1dd8SToby Isaac PetscFunctionReturn(0); 97620cf1dd8SToby Isaac } 97720cf1dd8SToby Isaac 97820cf1dd8SToby Isaac PetscErrorCode PetscFECreateHeightTrace(PetscFE fe, PetscInt height, PetscFE *trFE) 97920cf1dd8SToby Isaac { 98020cf1dd8SToby Isaac PetscInt hStart, hEnd; 98120cf1dd8SToby Isaac PetscDualSpace dsp; 98220cf1dd8SToby Isaac DM dm; 98320cf1dd8SToby Isaac PetscErrorCode ierr; 98420cf1dd8SToby Isaac 98520cf1dd8SToby Isaac PetscFunctionBegin; 98620cf1dd8SToby Isaac PetscValidHeaderSpecific(fe,PETSCFE_CLASSID,1); 98720cf1dd8SToby Isaac PetscValidPointer(trFE,3); 98820cf1dd8SToby Isaac *trFE = NULL; 98920cf1dd8SToby Isaac ierr = PetscFEGetDualSpace(fe,&dsp);CHKERRQ(ierr); 99020cf1dd8SToby Isaac ierr = PetscDualSpaceGetDM(dsp,&dm);CHKERRQ(ierr); 99120cf1dd8SToby Isaac ierr = DMPlexGetHeightStratum(dm,height,&hStart,&hEnd);CHKERRQ(ierr); 99220cf1dd8SToby Isaac if (hEnd <= hStart) PetscFunctionReturn(0); 99320cf1dd8SToby Isaac ierr = PetscFECreatePointTrace(fe,hStart,trFE);CHKERRQ(ierr); 99420cf1dd8SToby Isaac PetscFunctionReturn(0); 99520cf1dd8SToby Isaac } 99620cf1dd8SToby Isaac 99720cf1dd8SToby Isaac 99820cf1dd8SToby Isaac /*@ 99920cf1dd8SToby Isaac PetscFEGetDimension - Get the dimension of the finite element space on a cell 100020cf1dd8SToby Isaac 100120cf1dd8SToby Isaac Not collective 100220cf1dd8SToby Isaac 100320cf1dd8SToby Isaac Input Parameter: 100420cf1dd8SToby Isaac . fe - The PetscFE 100520cf1dd8SToby Isaac 100620cf1dd8SToby Isaac Output Parameter: 100720cf1dd8SToby Isaac . dim - The dimension 100820cf1dd8SToby Isaac 100920cf1dd8SToby Isaac Level: intermediate 101020cf1dd8SToby Isaac 101120cf1dd8SToby Isaac .seealso: PetscFECreate(), PetscSpaceGetDimension(), PetscDualSpaceGetDimension() 101220cf1dd8SToby Isaac @*/ 101320cf1dd8SToby Isaac PetscErrorCode PetscFEGetDimension(PetscFE fem, PetscInt *dim) 101420cf1dd8SToby Isaac { 101520cf1dd8SToby Isaac PetscErrorCode ierr; 101620cf1dd8SToby Isaac 101720cf1dd8SToby Isaac PetscFunctionBegin; 101820cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 101920cf1dd8SToby Isaac PetscValidPointer(dim, 2); 102020cf1dd8SToby Isaac if (fem->ops->getdimension) {ierr = (*fem->ops->getdimension)(fem, dim);CHKERRQ(ierr);} 102120cf1dd8SToby Isaac PetscFunctionReturn(0); 102220cf1dd8SToby Isaac } 102320cf1dd8SToby Isaac 102420cf1dd8SToby Isaac /* 102520cf1dd8SToby Isaac Purpose: Compute element vector for chunk of elements 102620cf1dd8SToby Isaac 102720cf1dd8SToby Isaac Input: 102820cf1dd8SToby Isaac Sizes: 102920cf1dd8SToby Isaac Ne: number of elements 103020cf1dd8SToby Isaac Nf: number of fields 103120cf1dd8SToby Isaac PetscFE 103220cf1dd8SToby Isaac dim: spatial dimension 103320cf1dd8SToby Isaac Nb: number of basis functions 103420cf1dd8SToby Isaac Nc: number of field components 103520cf1dd8SToby Isaac PetscQuadrature 103620cf1dd8SToby Isaac Nq: number of quadrature points 103720cf1dd8SToby Isaac 103820cf1dd8SToby Isaac Geometry: 103920cf1dd8SToby Isaac PetscFEGeom[Ne] possibly *Nq 104020cf1dd8SToby Isaac PetscReal v0s[dim] 104120cf1dd8SToby Isaac PetscReal n[dim] 104220cf1dd8SToby Isaac PetscReal jacobians[dim*dim] 104320cf1dd8SToby Isaac PetscReal jacobianInverses[dim*dim] 104420cf1dd8SToby Isaac PetscReal jacobianDeterminants 104520cf1dd8SToby Isaac FEM: 104620cf1dd8SToby Isaac PetscFE 104720cf1dd8SToby Isaac PetscQuadrature 104820cf1dd8SToby Isaac PetscReal quadPoints[Nq*dim] 104920cf1dd8SToby Isaac PetscReal quadWeights[Nq] 105020cf1dd8SToby Isaac PetscReal basis[Nq*Nb*Nc] 105120cf1dd8SToby Isaac PetscReal basisDer[Nq*Nb*Nc*dim] 105220cf1dd8SToby Isaac PetscScalar coefficients[Ne*Nb*Nc] 105320cf1dd8SToby Isaac PetscScalar elemVec[Ne*Nb*Nc] 105420cf1dd8SToby Isaac 105520cf1dd8SToby Isaac Problem: 105620cf1dd8SToby Isaac PetscInt f: the active field 105720cf1dd8SToby Isaac f0, f1 105820cf1dd8SToby Isaac 105920cf1dd8SToby Isaac Work Space: 106020cf1dd8SToby Isaac PetscFE 106120cf1dd8SToby Isaac PetscScalar f0[Nq*dim]; 106220cf1dd8SToby Isaac PetscScalar f1[Nq*dim*dim]; 106320cf1dd8SToby Isaac PetscScalar u[Nc]; 106420cf1dd8SToby Isaac PetscScalar gradU[Nc*dim]; 106520cf1dd8SToby Isaac PetscReal x[dim]; 106620cf1dd8SToby Isaac PetscScalar realSpaceDer[dim]; 106720cf1dd8SToby Isaac 106820cf1dd8SToby Isaac Purpose: Compute element vector for N_cb batches of elements 106920cf1dd8SToby Isaac 107020cf1dd8SToby Isaac Input: 107120cf1dd8SToby Isaac Sizes: 107220cf1dd8SToby Isaac N_cb: Number of serial cell batches 107320cf1dd8SToby Isaac 107420cf1dd8SToby Isaac Geometry: 107520cf1dd8SToby Isaac PetscReal v0s[Ne*dim] 107620cf1dd8SToby Isaac PetscReal jacobians[Ne*dim*dim] possibly *Nq 107720cf1dd8SToby Isaac PetscReal jacobianInverses[Ne*dim*dim] possibly *Nq 107820cf1dd8SToby Isaac PetscReal jacobianDeterminants[Ne] possibly *Nq 107920cf1dd8SToby Isaac FEM: 108020cf1dd8SToby Isaac static PetscReal quadPoints[Nq*dim] 108120cf1dd8SToby Isaac static PetscReal quadWeights[Nq] 108220cf1dd8SToby Isaac static PetscReal basis[Nq*Nb*Nc] 108320cf1dd8SToby Isaac static PetscReal basisDer[Nq*Nb*Nc*dim] 108420cf1dd8SToby Isaac PetscScalar coefficients[Ne*Nb*Nc] 108520cf1dd8SToby Isaac PetscScalar elemVec[Ne*Nb*Nc] 108620cf1dd8SToby Isaac 108720cf1dd8SToby Isaac ex62.c: 108820cf1dd8SToby Isaac PetscErrorCode PetscFEIntegrateResidualBatch(PetscInt Ne, PetscInt numFields, PetscInt field, PetscQuadrature quad[], const PetscScalar coefficients[], 108920cf1dd8SToby Isaac const PetscReal v0s[], const PetscReal jacobians[], const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], 109020cf1dd8SToby Isaac void (*f0_func)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar f0[]), 109120cf1dd8SToby Isaac void (*f1_func)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar f1[]), PetscScalar elemVec[]) 109220cf1dd8SToby Isaac 109320cf1dd8SToby Isaac ex52.c: 109420cf1dd8SToby 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) 109520cf1dd8SToby 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) 109620cf1dd8SToby Isaac 109720cf1dd8SToby Isaac ex52_integrateElement.cu 109820cf1dd8SToby Isaac __global__ void integrateElementQuadrature(int N_cb, realType *coefficients, realType *jacobianInverses, realType *jacobianDeterminants, realType *elemVec) 109920cf1dd8SToby Isaac 110020cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode IntegrateElementBatchGPU(PetscInt spatial_dim, PetscInt Ne, PetscInt Ncb, PetscInt Nbc, PetscInt Nbl, const PetscScalar coefficients[], 110120cf1dd8SToby Isaac const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], PetscScalar elemVec[], 110220cf1dd8SToby Isaac PetscLogEvent event, PetscInt debug, PetscInt pde_op) 110320cf1dd8SToby Isaac 110420cf1dd8SToby Isaac ex52_integrateElementOpenCL.c: 110520cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode IntegrateElementBatchGPU(PetscInt spatial_dim, PetscInt Ne, PetscInt Ncb, PetscInt Nbc, PetscInt N_bl, const PetscScalar coefficients[], 110620cf1dd8SToby Isaac const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], PetscScalar elemVec[], 110720cf1dd8SToby Isaac PetscLogEvent event, PetscInt debug, PetscInt pde_op) 110820cf1dd8SToby Isaac 110920cf1dd8SToby Isaac __kernel void integrateElementQuadrature(int N_cb, __global float *coefficients, __global float *jacobianInverses, __global float *jacobianDeterminants, __global float *elemVec) 111020cf1dd8SToby Isaac */ 111120cf1dd8SToby Isaac 111220cf1dd8SToby Isaac /*@C 111320cf1dd8SToby Isaac PetscFEIntegrate - Produce the integral for the given field for a chunk of elements by quadrature integration 111420cf1dd8SToby Isaac 111520cf1dd8SToby Isaac Not collective 111620cf1dd8SToby Isaac 111720cf1dd8SToby Isaac Input Parameters: 111820cf1dd8SToby Isaac + fem - The PetscFE object for the field being integrated 111920cf1dd8SToby Isaac . prob - The PetscDS specifying the discretizations and continuum functions 112020cf1dd8SToby Isaac . field - The field being integrated 112120cf1dd8SToby Isaac . Ne - The number of elements in the chunk 112220cf1dd8SToby Isaac . cgeom - The cell geometry for each cell in the chunk 112320cf1dd8SToby Isaac . coefficients - The array of FEM basis coefficients for the elements 112420cf1dd8SToby Isaac . probAux - The PetscDS specifying the auxiliary discretizations 112520cf1dd8SToby Isaac - coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 112620cf1dd8SToby Isaac 112720cf1dd8SToby Isaac Output Parameter 112820cf1dd8SToby Isaac . integral - the integral for this field 112920cf1dd8SToby Isaac 113020cf1dd8SToby Isaac Level: developer 113120cf1dd8SToby Isaac 113220cf1dd8SToby Isaac .seealso: PetscFEIntegrateResidual() 113320cf1dd8SToby Isaac @*/ 113420cf1dd8SToby Isaac PetscErrorCode PetscFEIntegrate(PetscFE fem, PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom, 113520cf1dd8SToby Isaac const PetscScalar coefficients[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscScalar integral[]) 113620cf1dd8SToby Isaac { 113720cf1dd8SToby Isaac PetscErrorCode ierr; 113820cf1dd8SToby Isaac 113920cf1dd8SToby Isaac PetscFunctionBegin; 114020cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 114120cf1dd8SToby Isaac PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 2); 114220cf1dd8SToby Isaac if (fem->ops->integrate) {ierr = (*fem->ops->integrate)(fem, prob, field, Ne, cgeom, coefficients, probAux, coefficientsAux, integral);CHKERRQ(ierr);} 114320cf1dd8SToby Isaac PetscFunctionReturn(0); 114420cf1dd8SToby Isaac } 114520cf1dd8SToby Isaac 114620cf1dd8SToby Isaac /*@C 1147afe6d6adSToby Isaac PetscFEIntegrateBd - Produce the integral for the given field for a chunk of elements by quadrature integration 1148afe6d6adSToby Isaac 1149afe6d6adSToby Isaac Not collective 1150afe6d6adSToby Isaac 1151afe6d6adSToby Isaac Input Parameters: 1152afe6d6adSToby Isaac + fem - The PetscFE object for the field being integrated 1153afe6d6adSToby Isaac . prob - The PetscDS specifying the discretizations and continuum functions 1154afe6d6adSToby Isaac . field - The field being integrated 1155afe6d6adSToby Isaac . obj_func - The function to be integrated 1156afe6d6adSToby Isaac . Ne - The number of elements in the chunk 1157afe6d6adSToby Isaac . fgeom - The face geometry for each face in the chunk 1158afe6d6adSToby Isaac . coefficients - The array of FEM basis coefficients for the elements 1159afe6d6adSToby Isaac . probAux - The PetscDS specifying the auxiliary discretizations 1160afe6d6adSToby Isaac - coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 1161afe6d6adSToby Isaac 1162afe6d6adSToby Isaac Output Parameter 1163afe6d6adSToby Isaac . integral - the integral for this field 1164afe6d6adSToby Isaac 1165afe6d6adSToby Isaac Level: developer 1166afe6d6adSToby Isaac 1167afe6d6adSToby Isaac .seealso: PetscFEIntegrateResidual() 1168afe6d6adSToby Isaac @*/ 1169afe6d6adSToby Isaac PetscErrorCode PetscFEIntegrateBd(PetscFE fem, PetscDS prob, PetscInt field, 1170afe6d6adSToby Isaac void (*obj_func)(PetscInt, PetscInt, PetscInt, 1171afe6d6adSToby Isaac const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 1172afe6d6adSToby Isaac const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 1173afe6d6adSToby Isaac PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 1174afe6d6adSToby Isaac PetscInt Ne, PetscFEGeom *geom, const PetscScalar coefficients[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscScalar integral[]) 1175afe6d6adSToby Isaac { 1176afe6d6adSToby Isaac PetscErrorCode ierr; 1177afe6d6adSToby Isaac 1178afe6d6adSToby Isaac PetscFunctionBegin; 1179afe6d6adSToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 1180afe6d6adSToby Isaac PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 2); 1181afe6d6adSToby Isaac if (fem->ops->integratebd) {ierr = (*fem->ops->integratebd)(fem, prob, field, obj_func, Ne, geom, coefficients, probAux, coefficientsAux, integral);CHKERRQ(ierr);} 1182afe6d6adSToby Isaac PetscFunctionReturn(0); 1183afe6d6adSToby Isaac } 1184afe6d6adSToby Isaac 1185afe6d6adSToby Isaac /*@C 118620cf1dd8SToby Isaac PetscFEIntegrateResidual - Produce the element residual vector for a chunk of elements by quadrature integration 118720cf1dd8SToby Isaac 118820cf1dd8SToby Isaac Not collective 118920cf1dd8SToby Isaac 119020cf1dd8SToby Isaac Input Parameters: 119120cf1dd8SToby Isaac + fem - The PetscFE object for the field being integrated 119220cf1dd8SToby Isaac . prob - The PetscDS specifying the discretizations and continuum functions 119320cf1dd8SToby Isaac . field - The field being integrated 119420cf1dd8SToby Isaac . Ne - The number of elements in the chunk 119520cf1dd8SToby Isaac . cgeom - The cell geometry for each cell in the chunk 119620cf1dd8SToby Isaac . coefficients - The array of FEM basis coefficients for the elements 119720cf1dd8SToby Isaac . coefficients_t - The array of FEM basis time derivative coefficients for the elements 119820cf1dd8SToby Isaac . probAux - The PetscDS specifying the auxiliary discretizations 119920cf1dd8SToby Isaac . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 120020cf1dd8SToby Isaac - t - The time 120120cf1dd8SToby Isaac 120220cf1dd8SToby Isaac Output Parameter 120320cf1dd8SToby Isaac . elemVec - the element residual vectors from each element 120420cf1dd8SToby Isaac 120520cf1dd8SToby Isaac Note: 120620cf1dd8SToby Isaac $ Loop over batch of elements (e): 120720cf1dd8SToby Isaac $ Loop over quadrature points (q): 120820cf1dd8SToby Isaac $ Make u_q and gradU_q (loops over fields,Nb,Ncomp) and x_q 120920cf1dd8SToby Isaac $ Call f_0 and f_1 121020cf1dd8SToby Isaac $ Loop over element vector entries (f,fc --> i): 121120cf1dd8SToby Isaac $ elemVec[i] += \psi^{fc}_f(q) f0_{fc}(u, \nabla u) + \nabla\psi^{fc}_f(q) \cdot f1_{fc,df}(u, \nabla u) 121220cf1dd8SToby Isaac 121320cf1dd8SToby Isaac Level: developer 121420cf1dd8SToby Isaac 121520cf1dd8SToby Isaac .seealso: PetscFEIntegrateResidual() 121620cf1dd8SToby Isaac @*/ 121720cf1dd8SToby Isaac PetscErrorCode PetscFEIntegrateResidual(PetscFE fem, PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom, 121820cf1dd8SToby Isaac const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[]) 121920cf1dd8SToby Isaac { 122020cf1dd8SToby Isaac PetscErrorCode ierr; 122120cf1dd8SToby Isaac 122220cf1dd8SToby Isaac PetscFunctionBegin; 122320cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 122420cf1dd8SToby Isaac PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 2); 122520cf1dd8SToby Isaac if (fem->ops->integrateresidual) {ierr = (*fem->ops->integrateresidual)(fem, prob, field, Ne, cgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec);CHKERRQ(ierr);} 122620cf1dd8SToby Isaac PetscFunctionReturn(0); 122720cf1dd8SToby Isaac } 122820cf1dd8SToby Isaac 122920cf1dd8SToby Isaac /*@C 123020cf1dd8SToby Isaac PetscFEIntegrateBdResidual - Produce the element residual vector for a chunk of elements by quadrature integration over a boundary 123120cf1dd8SToby Isaac 123220cf1dd8SToby Isaac Not collective 123320cf1dd8SToby Isaac 123420cf1dd8SToby Isaac Input Parameters: 123520cf1dd8SToby Isaac + fem - The PetscFE object for the field being integrated 123620cf1dd8SToby Isaac . prob - The PetscDS specifying the discretizations and continuum functions 123720cf1dd8SToby Isaac . field - The field being integrated 123820cf1dd8SToby Isaac . Ne - The number of elements in the chunk 123920cf1dd8SToby Isaac . fgeom - The face geometry for each cell in the chunk 124020cf1dd8SToby Isaac . coefficients - The array of FEM basis coefficients for the elements 124120cf1dd8SToby Isaac . coefficients_t - The array of FEM basis time derivative coefficients for the elements 124220cf1dd8SToby Isaac . probAux - The PetscDS specifying the auxiliary discretizations 124320cf1dd8SToby Isaac . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 124420cf1dd8SToby Isaac - t - The time 124520cf1dd8SToby Isaac 124620cf1dd8SToby Isaac Output Parameter 124720cf1dd8SToby Isaac . elemVec - the element residual vectors from each element 124820cf1dd8SToby Isaac 124920cf1dd8SToby Isaac Level: developer 125020cf1dd8SToby Isaac 125120cf1dd8SToby Isaac .seealso: PetscFEIntegrateResidual() 125220cf1dd8SToby Isaac @*/ 125320cf1dd8SToby Isaac PetscErrorCode PetscFEIntegrateBdResidual(PetscFE fem, PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *fgeom, 125420cf1dd8SToby Isaac const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[]) 125520cf1dd8SToby Isaac { 125620cf1dd8SToby Isaac PetscErrorCode ierr; 125720cf1dd8SToby Isaac 125820cf1dd8SToby Isaac PetscFunctionBegin; 125920cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 126020cf1dd8SToby Isaac if (fem->ops->integratebdresidual) {ierr = (*fem->ops->integratebdresidual)(fem, prob, field, Ne, fgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec);CHKERRQ(ierr);} 126120cf1dd8SToby Isaac PetscFunctionReturn(0); 126220cf1dd8SToby Isaac } 126320cf1dd8SToby Isaac 126420cf1dd8SToby Isaac /*@C 126520cf1dd8SToby Isaac PetscFEIntegrateJacobian - Produce the element Jacobian for a chunk of elements by quadrature integration 126620cf1dd8SToby Isaac 126720cf1dd8SToby Isaac Not collective 126820cf1dd8SToby Isaac 126920cf1dd8SToby Isaac Input Parameters: 127020cf1dd8SToby Isaac + fem - The PetscFE object for the field being integrated 127120cf1dd8SToby Isaac . prob - The PetscDS specifying the discretizations and continuum functions 127220cf1dd8SToby Isaac . jtype - The type of matrix pointwise functions that should be used 127320cf1dd8SToby Isaac . fieldI - The test field being integrated 127420cf1dd8SToby Isaac . fieldJ - The basis field being integrated 127520cf1dd8SToby Isaac . Ne - The number of elements in the chunk 127620cf1dd8SToby Isaac . cgeom - The cell geometry for each cell in the chunk 127720cf1dd8SToby Isaac . coefficients - The array of FEM basis coefficients for the elements for the Jacobian evaluation point 127820cf1dd8SToby Isaac . coefficients_t - The array of FEM basis time derivative coefficients for the elements 127920cf1dd8SToby Isaac . probAux - The PetscDS specifying the auxiliary discretizations 128020cf1dd8SToby Isaac . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 128120cf1dd8SToby Isaac . t - The time 128220cf1dd8SToby Isaac - u_tShift - A multiplier for the dF/du_t term (as opposed to the dF/du term) 128320cf1dd8SToby Isaac 128420cf1dd8SToby Isaac Output Parameter 128520cf1dd8SToby Isaac . elemMat - the element matrices for the Jacobian from each element 128620cf1dd8SToby Isaac 128720cf1dd8SToby Isaac Note: 128820cf1dd8SToby Isaac $ Loop over batch of elements (e): 128920cf1dd8SToby Isaac $ Loop over element matrix entries (f,fc,g,gc --> i,j): 129020cf1dd8SToby Isaac $ Loop over quadrature points (q): 129120cf1dd8SToby Isaac $ Make u_q and gradU_q (loops over fields,Nb,Ncomp) 129220cf1dd8SToby Isaac $ elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q) 129320cf1dd8SToby Isaac $ + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q) 129420cf1dd8SToby Isaac $ + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q) 129520cf1dd8SToby Isaac $ + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q) 129620cf1dd8SToby Isaac Level: developer 129720cf1dd8SToby Isaac 129820cf1dd8SToby Isaac .seealso: PetscFEIntegrateResidual() 129920cf1dd8SToby Isaac @*/ 130020cf1dd8SToby Isaac PetscErrorCode PetscFEIntegrateJacobian(PetscFE fem, PetscDS prob, PetscFEJacobianType jtype, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *cgeom, 130120cf1dd8SToby Isaac const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) 130220cf1dd8SToby Isaac { 130320cf1dd8SToby Isaac PetscErrorCode ierr; 130420cf1dd8SToby Isaac 130520cf1dd8SToby Isaac PetscFunctionBegin; 130620cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 130720cf1dd8SToby 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);} 130820cf1dd8SToby Isaac PetscFunctionReturn(0); 130920cf1dd8SToby Isaac } 131020cf1dd8SToby Isaac 131120cf1dd8SToby Isaac /*@C 131220cf1dd8SToby Isaac PetscFEIntegrateBdJacobian - Produce the boundary element Jacobian for a chunk of elements by quadrature integration 131320cf1dd8SToby Isaac 131420cf1dd8SToby Isaac Not collective 131520cf1dd8SToby Isaac 131620cf1dd8SToby Isaac Input Parameters: 131720cf1dd8SToby Isaac + fem = The PetscFE object for the field being integrated 131820cf1dd8SToby Isaac . prob - The PetscDS specifying the discretizations and continuum functions 131920cf1dd8SToby Isaac . fieldI - The test field being integrated 132020cf1dd8SToby Isaac . fieldJ - The basis field being integrated 132120cf1dd8SToby Isaac . Ne - The number of elements in the chunk 132220cf1dd8SToby Isaac . fgeom - The face geometry for each cell in the chunk 132320cf1dd8SToby Isaac . coefficients - The array of FEM basis coefficients for the elements for the Jacobian evaluation point 132420cf1dd8SToby Isaac . coefficients_t - The array of FEM basis time derivative coefficients for the elements 132520cf1dd8SToby Isaac . probAux - The PetscDS specifying the auxiliary discretizations 132620cf1dd8SToby Isaac . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 132720cf1dd8SToby Isaac . t - The time 132820cf1dd8SToby Isaac - u_tShift - A multiplier for the dF/du_t term (as opposed to the dF/du term) 132920cf1dd8SToby Isaac 133020cf1dd8SToby Isaac Output Parameter 133120cf1dd8SToby Isaac . elemMat - the element matrices for the Jacobian from each element 133220cf1dd8SToby Isaac 133320cf1dd8SToby Isaac Note: 133420cf1dd8SToby Isaac $ Loop over batch of elements (e): 133520cf1dd8SToby Isaac $ Loop over element matrix entries (f,fc,g,gc --> i,j): 133620cf1dd8SToby Isaac $ Loop over quadrature points (q): 133720cf1dd8SToby Isaac $ Make u_q and gradU_q (loops over fields,Nb,Ncomp) 133820cf1dd8SToby Isaac $ elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q) 133920cf1dd8SToby Isaac $ + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q) 134020cf1dd8SToby Isaac $ + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q) 134120cf1dd8SToby Isaac $ + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q) 134220cf1dd8SToby Isaac Level: developer 134320cf1dd8SToby Isaac 134420cf1dd8SToby Isaac .seealso: PetscFEIntegrateJacobian(), PetscFEIntegrateResidual() 134520cf1dd8SToby Isaac @*/ 134620cf1dd8SToby Isaac PetscErrorCode PetscFEIntegrateBdJacobian(PetscFE fem, PetscDS prob, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *fgeom, 134720cf1dd8SToby Isaac const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) 134820cf1dd8SToby Isaac { 134920cf1dd8SToby Isaac PetscErrorCode ierr; 135020cf1dd8SToby Isaac 135120cf1dd8SToby Isaac PetscFunctionBegin; 135220cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 135320cf1dd8SToby 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);} 135420cf1dd8SToby Isaac PetscFunctionReturn(0); 135520cf1dd8SToby Isaac } 135620cf1dd8SToby Isaac 135720cf1dd8SToby Isaac PetscErrorCode PetscFEGetHeightSubspace(PetscFE fe, PetscInt height, PetscFE *subfe) 135820cf1dd8SToby Isaac { 135920cf1dd8SToby Isaac PetscSpace P, subP; 136020cf1dd8SToby Isaac PetscDualSpace Q, subQ; 136120cf1dd8SToby Isaac PetscQuadrature subq; 136220cf1dd8SToby Isaac PetscFEType fetype; 136320cf1dd8SToby Isaac PetscInt dim, Nc; 136420cf1dd8SToby Isaac PetscErrorCode ierr; 136520cf1dd8SToby Isaac 136620cf1dd8SToby Isaac PetscFunctionBegin; 136720cf1dd8SToby Isaac PetscValidHeaderSpecific(fe, PETSCFE_CLASSID, 1); 136820cf1dd8SToby Isaac PetscValidPointer(subfe, 3); 136920cf1dd8SToby Isaac if (height == 0) { 137020cf1dd8SToby Isaac *subfe = fe; 137120cf1dd8SToby Isaac PetscFunctionReturn(0); 137220cf1dd8SToby Isaac } 137320cf1dd8SToby Isaac ierr = PetscFEGetBasisSpace(fe, &P);CHKERRQ(ierr); 137420cf1dd8SToby Isaac ierr = PetscFEGetDualSpace(fe, &Q);CHKERRQ(ierr); 137520cf1dd8SToby Isaac ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 137620cf1dd8SToby Isaac ierr = PetscFEGetFaceQuadrature(fe, &subq);CHKERRQ(ierr); 137720cf1dd8SToby Isaac ierr = PetscDualSpaceGetDimension(Q, &dim);CHKERRQ(ierr); 137820cf1dd8SToby 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);} 137920cf1dd8SToby Isaac if (!fe->subspaces) {ierr = PetscCalloc1(dim, &fe->subspaces);CHKERRQ(ierr);} 138020cf1dd8SToby Isaac if (height <= dim) { 138120cf1dd8SToby Isaac if (!fe->subspaces[height-1]) { 138220cf1dd8SToby Isaac PetscFE sub; 13833f6b16c7SMatthew G. Knepley const char *name; 138420cf1dd8SToby Isaac 138520cf1dd8SToby Isaac ierr = PetscSpaceGetHeightSubspace(P, height, &subP);CHKERRQ(ierr); 138620cf1dd8SToby Isaac ierr = PetscDualSpaceGetHeightSubspace(Q, height, &subQ);CHKERRQ(ierr); 138720cf1dd8SToby Isaac ierr = PetscFECreate(PetscObjectComm((PetscObject) fe), &sub);CHKERRQ(ierr); 13883f6b16c7SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) fe, &name);CHKERRQ(ierr); 13893f6b16c7SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) sub, name);CHKERRQ(ierr); 139020cf1dd8SToby Isaac ierr = PetscFEGetType(fe, &fetype);CHKERRQ(ierr); 139120cf1dd8SToby Isaac ierr = PetscFESetType(sub, fetype);CHKERRQ(ierr); 139220cf1dd8SToby Isaac ierr = PetscFESetBasisSpace(sub, subP);CHKERRQ(ierr); 139320cf1dd8SToby Isaac ierr = PetscFESetDualSpace(sub, subQ);CHKERRQ(ierr); 139420cf1dd8SToby Isaac ierr = PetscFESetNumComponents(sub, Nc);CHKERRQ(ierr); 139520cf1dd8SToby Isaac ierr = PetscFESetUp(sub);CHKERRQ(ierr); 139620cf1dd8SToby Isaac ierr = PetscFESetQuadrature(sub, subq);CHKERRQ(ierr); 139720cf1dd8SToby Isaac fe->subspaces[height-1] = sub; 139820cf1dd8SToby Isaac } 139920cf1dd8SToby Isaac *subfe = fe->subspaces[height-1]; 140020cf1dd8SToby Isaac } else { 140120cf1dd8SToby Isaac *subfe = NULL; 140220cf1dd8SToby Isaac } 140320cf1dd8SToby Isaac PetscFunctionReturn(0); 140420cf1dd8SToby Isaac } 140520cf1dd8SToby Isaac 140620cf1dd8SToby Isaac /*@ 140720cf1dd8SToby Isaac PetscFERefine - Create a "refined" PetscFE object that refines the reference cell into smaller copies. This is typically used 140820cf1dd8SToby Isaac to precondition a higher order method with a lower order method on a refined mesh having the same number of dofs (but more 140920cf1dd8SToby Isaac sparsity). It is also used to create an interpolation between regularly refined meshes. 141020cf1dd8SToby Isaac 141120cf1dd8SToby Isaac Collective on PetscFE 141220cf1dd8SToby Isaac 141320cf1dd8SToby Isaac Input Parameter: 141420cf1dd8SToby Isaac . fe - The initial PetscFE 141520cf1dd8SToby Isaac 141620cf1dd8SToby Isaac Output Parameter: 141720cf1dd8SToby Isaac . feRef - The refined PetscFE 141820cf1dd8SToby Isaac 141920cf1dd8SToby Isaac Level: developer 142020cf1dd8SToby Isaac 142120cf1dd8SToby Isaac .seealso: PetscFEType, PetscFECreate(), PetscFESetType() 142220cf1dd8SToby Isaac @*/ 142320cf1dd8SToby Isaac PetscErrorCode PetscFERefine(PetscFE fe, PetscFE *feRef) 142420cf1dd8SToby Isaac { 142520cf1dd8SToby Isaac PetscSpace P, Pref; 142620cf1dd8SToby Isaac PetscDualSpace Q, Qref; 142720cf1dd8SToby Isaac DM K, Kref; 142820cf1dd8SToby Isaac PetscQuadrature q, qref; 142920cf1dd8SToby Isaac const PetscReal *v0, *jac; 143020cf1dd8SToby Isaac PetscInt numComp, numSubelements; 143120cf1dd8SToby Isaac PetscErrorCode ierr; 143220cf1dd8SToby Isaac 143320cf1dd8SToby Isaac PetscFunctionBegin; 143420cf1dd8SToby Isaac ierr = PetscFEGetBasisSpace(fe, &P);CHKERRQ(ierr); 143520cf1dd8SToby Isaac ierr = PetscFEGetDualSpace(fe, &Q);CHKERRQ(ierr); 143620cf1dd8SToby Isaac ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr); 143720cf1dd8SToby Isaac ierr = PetscDualSpaceGetDM(Q, &K);CHKERRQ(ierr); 143820cf1dd8SToby Isaac /* Create space */ 143920cf1dd8SToby Isaac ierr = PetscObjectReference((PetscObject) P);CHKERRQ(ierr); 144020cf1dd8SToby Isaac Pref = P; 144120cf1dd8SToby Isaac /* Create dual space */ 144220cf1dd8SToby Isaac ierr = PetscDualSpaceDuplicate(Q, &Qref);CHKERRQ(ierr); 144320cf1dd8SToby Isaac ierr = DMRefine(K, PetscObjectComm((PetscObject) fe), &Kref);CHKERRQ(ierr); 144420cf1dd8SToby Isaac ierr = PetscDualSpaceSetDM(Qref, Kref);CHKERRQ(ierr); 144520cf1dd8SToby Isaac ierr = DMDestroy(&Kref);CHKERRQ(ierr); 144620cf1dd8SToby Isaac ierr = PetscDualSpaceSetUp(Qref);CHKERRQ(ierr); 144720cf1dd8SToby Isaac /* Create element */ 144820cf1dd8SToby Isaac ierr = PetscFECreate(PetscObjectComm((PetscObject) fe), feRef);CHKERRQ(ierr); 144920cf1dd8SToby Isaac ierr = PetscFESetType(*feRef, PETSCFECOMPOSITE);CHKERRQ(ierr); 145020cf1dd8SToby Isaac ierr = PetscFESetBasisSpace(*feRef, Pref);CHKERRQ(ierr); 145120cf1dd8SToby Isaac ierr = PetscFESetDualSpace(*feRef, Qref);CHKERRQ(ierr); 145220cf1dd8SToby Isaac ierr = PetscFEGetNumComponents(fe, &numComp);CHKERRQ(ierr); 145320cf1dd8SToby Isaac ierr = PetscFESetNumComponents(*feRef, numComp);CHKERRQ(ierr); 145420cf1dd8SToby Isaac ierr = PetscFESetUp(*feRef);CHKERRQ(ierr); 145520cf1dd8SToby Isaac ierr = PetscSpaceDestroy(&Pref);CHKERRQ(ierr); 145620cf1dd8SToby Isaac ierr = PetscDualSpaceDestroy(&Qref);CHKERRQ(ierr); 145720cf1dd8SToby Isaac /* Create quadrature */ 145820cf1dd8SToby Isaac ierr = PetscFECompositeGetMapping(*feRef, &numSubelements, &v0, &jac, NULL);CHKERRQ(ierr); 145920cf1dd8SToby Isaac ierr = PetscQuadratureExpandComposite(q, numSubelements, v0, jac, &qref);CHKERRQ(ierr); 146020cf1dd8SToby Isaac ierr = PetscFESetQuadrature(*feRef, qref);CHKERRQ(ierr); 146120cf1dd8SToby Isaac ierr = PetscQuadratureDestroy(&qref);CHKERRQ(ierr); 146220cf1dd8SToby Isaac PetscFunctionReturn(0); 146320cf1dd8SToby Isaac } 146420cf1dd8SToby Isaac 146520cf1dd8SToby Isaac /*@C 146620cf1dd8SToby Isaac PetscFECreateDefault - Create a PetscFE for basic FEM computation 146720cf1dd8SToby Isaac 146820cf1dd8SToby Isaac Collective on DM 146920cf1dd8SToby Isaac 147020cf1dd8SToby Isaac Input Parameters: 14717be5e748SToby Isaac + comm - The MPI comm 147220cf1dd8SToby Isaac . dim - The spatial dimension 147320cf1dd8SToby Isaac . Nc - The number of components 147420cf1dd8SToby Isaac . isSimplex - Flag for simplex reference cell, otherwise its a tensor product 147520cf1dd8SToby Isaac . prefix - The options prefix, or NULL 147620cf1dd8SToby Isaac - qorder - The quadrature order 147720cf1dd8SToby Isaac 147820cf1dd8SToby Isaac Output Parameter: 147920cf1dd8SToby Isaac . fem - The PetscFE object 148020cf1dd8SToby Isaac 148120cf1dd8SToby Isaac Level: beginner 148220cf1dd8SToby Isaac 148320cf1dd8SToby Isaac .keywords: PetscFE, finite element 148420cf1dd8SToby Isaac .seealso: PetscFECreate(), PetscSpaceCreate(), PetscDualSpaceCreate() 148520cf1dd8SToby Isaac @*/ 14867be5e748SToby Isaac PetscErrorCode PetscFECreateDefault(MPI_Comm comm, PetscInt dim, PetscInt Nc, PetscBool isSimplex, const char prefix[], PetscInt qorder, PetscFE *fem) 148720cf1dd8SToby Isaac { 148820cf1dd8SToby Isaac PetscQuadrature q, fq; 148920cf1dd8SToby Isaac DM K; 149020cf1dd8SToby Isaac PetscSpace P; 149120cf1dd8SToby Isaac PetscDualSpace Q; 149220cf1dd8SToby Isaac PetscInt order, quadPointsPerEdge; 149320cf1dd8SToby Isaac PetscBool tensor = isSimplex ? PETSC_FALSE : PETSC_TRUE; 149420cf1dd8SToby Isaac PetscErrorCode ierr; 149520cf1dd8SToby Isaac 149620cf1dd8SToby Isaac PetscFunctionBegin; 149720cf1dd8SToby Isaac /* Create space */ 14987be5e748SToby Isaac ierr = PetscSpaceCreate(comm, &P);CHKERRQ(ierr); 149920cf1dd8SToby Isaac ierr = PetscObjectSetOptionsPrefix((PetscObject) P, prefix);CHKERRQ(ierr); 150020cf1dd8SToby Isaac ierr = PetscSpacePolynomialSetTensor(P, tensor);CHKERRQ(ierr); 150120cf1dd8SToby Isaac ierr = PetscSpaceSetNumComponents(P, Nc);CHKERRQ(ierr); 150220cf1dd8SToby Isaac ierr = PetscSpaceSetNumVariables(P, dim);CHKERRQ(ierr); 1503028afddaSToby Isaac ierr = PetscSpaceSetFromOptions(P);CHKERRQ(ierr); 150420cf1dd8SToby Isaac ierr = PetscSpaceSetUp(P);CHKERRQ(ierr); 150520cf1dd8SToby Isaac ierr = PetscSpaceGetDegree(P, &order, NULL);CHKERRQ(ierr); 150620cf1dd8SToby Isaac ierr = PetscSpacePolynomialGetTensor(P, &tensor);CHKERRQ(ierr); 150720cf1dd8SToby Isaac /* Create dual space */ 15087be5e748SToby Isaac ierr = PetscDualSpaceCreate(comm, &Q);CHKERRQ(ierr); 150920cf1dd8SToby Isaac ierr = PetscDualSpaceSetType(Q,PETSCDUALSPACELAGRANGE);CHKERRQ(ierr); 151020cf1dd8SToby Isaac ierr = PetscObjectSetOptionsPrefix((PetscObject) Q, prefix);CHKERRQ(ierr); 151120cf1dd8SToby Isaac ierr = PetscDualSpaceCreateReferenceCell(Q, dim, isSimplex, &K);CHKERRQ(ierr); 151220cf1dd8SToby Isaac ierr = PetscDualSpaceSetDM(Q, K);CHKERRQ(ierr); 151320cf1dd8SToby Isaac ierr = DMDestroy(&K);CHKERRQ(ierr); 151420cf1dd8SToby Isaac ierr = PetscDualSpaceSetNumComponents(Q, Nc);CHKERRQ(ierr); 151520cf1dd8SToby Isaac ierr = PetscDualSpaceSetOrder(Q, order);CHKERRQ(ierr); 151620cf1dd8SToby Isaac ierr = PetscDualSpaceLagrangeSetTensor(Q, tensor);CHKERRQ(ierr); 151720cf1dd8SToby Isaac ierr = PetscDualSpaceSetFromOptions(Q);CHKERRQ(ierr); 151820cf1dd8SToby Isaac ierr = PetscDualSpaceSetUp(Q);CHKERRQ(ierr); 151920cf1dd8SToby Isaac /* Create element */ 15207be5e748SToby Isaac ierr = PetscFECreate(comm, fem);CHKERRQ(ierr); 152120cf1dd8SToby Isaac ierr = PetscObjectSetOptionsPrefix((PetscObject) *fem, prefix);CHKERRQ(ierr); 152220cf1dd8SToby Isaac ierr = PetscFESetFromOptions(*fem);CHKERRQ(ierr); 152320cf1dd8SToby Isaac ierr = PetscFESetBasisSpace(*fem, P);CHKERRQ(ierr); 152420cf1dd8SToby Isaac ierr = PetscFESetDualSpace(*fem, Q);CHKERRQ(ierr); 152520cf1dd8SToby Isaac ierr = PetscFESetNumComponents(*fem, Nc);CHKERRQ(ierr); 152620cf1dd8SToby Isaac ierr = PetscFESetUp(*fem);CHKERRQ(ierr); 152720cf1dd8SToby Isaac ierr = PetscSpaceDestroy(&P);CHKERRQ(ierr); 152820cf1dd8SToby Isaac ierr = PetscDualSpaceDestroy(&Q);CHKERRQ(ierr); 152920cf1dd8SToby Isaac /* Create quadrature (with specified order if given) */ 153020cf1dd8SToby Isaac qorder = qorder >= 0 ? qorder : order; 153120cf1dd8SToby Isaac ierr = PetscObjectOptionsBegin((PetscObject)*fem);CHKERRQ(ierr); 153220cf1dd8SToby Isaac ierr = PetscOptionsInt("-petscfe_default_quadrature_order","Quadrature order is one less than quadture points per edge","PetscFECreateDefault",qorder,&qorder,NULL);CHKERRQ(ierr); 153320cf1dd8SToby Isaac ierr = PetscOptionsEnd();CHKERRQ(ierr); 153420cf1dd8SToby Isaac quadPointsPerEdge = PetscMax(qorder + 1,1); 153520cf1dd8SToby Isaac if (isSimplex) { 153620cf1dd8SToby Isaac ierr = PetscDTGaussJacobiQuadrature(dim, 1, quadPointsPerEdge, -1.0, 1.0, &q);CHKERRQ(ierr); 153720cf1dd8SToby Isaac ierr = PetscDTGaussJacobiQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0, &fq);CHKERRQ(ierr); 15384ccfa306SStefano Zampini } else { 153920cf1dd8SToby Isaac ierr = PetscDTGaussTensorQuadrature(dim, 1, quadPointsPerEdge, -1.0, 1.0, &q);CHKERRQ(ierr); 154020cf1dd8SToby Isaac ierr = PetscDTGaussTensorQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0, &fq);CHKERRQ(ierr); 154120cf1dd8SToby Isaac } 154220cf1dd8SToby Isaac ierr = PetscFESetQuadrature(*fem, q);CHKERRQ(ierr); 154320cf1dd8SToby Isaac ierr = PetscFESetFaceQuadrature(*fem, fq);CHKERRQ(ierr); 154420cf1dd8SToby Isaac ierr = PetscQuadratureDestroy(&q);CHKERRQ(ierr); 154520cf1dd8SToby Isaac ierr = PetscQuadratureDestroy(&fq);CHKERRQ(ierr); 154620cf1dd8SToby Isaac PetscFunctionReturn(0); 154720cf1dd8SToby Isaac } 15483f6b16c7SMatthew G. Knepley 15493f6b16c7SMatthew G. Knepley /*@C 15503f6b16c7SMatthew G. Knepley PetscFESetName - Names the FE and its subobjects 15513f6b16c7SMatthew G. Knepley 15523f6b16c7SMatthew G. Knepley Not collective 15533f6b16c7SMatthew G. Knepley 15543f6b16c7SMatthew G. Knepley Input Parameters: 15553f6b16c7SMatthew G. Knepley + fe - The PetscFE 15563f6b16c7SMatthew G. Knepley - name - The name 15573f6b16c7SMatthew G. Knepley 15583f6b16c7SMatthew G. Knepley Level: beginner 15593f6b16c7SMatthew G. Knepley 15603f6b16c7SMatthew G. Knepley .keywords: PetscFE, finite element 15613f6b16c7SMatthew G. Knepley .seealso: PetscFECreate(), PetscSpaceCreate(), PetscDualSpaceCreate() 15623f6b16c7SMatthew G. Knepley @*/ 15633f6b16c7SMatthew G. Knepley PetscErrorCode PetscFESetName(PetscFE fe, const char name[]) 15643f6b16c7SMatthew G. Knepley { 15653f6b16c7SMatthew G. Knepley PetscSpace P; 15663f6b16c7SMatthew G. Knepley PetscDualSpace Q; 15673f6b16c7SMatthew G. Knepley PetscErrorCode ierr; 15683f6b16c7SMatthew G. Knepley 15693f6b16c7SMatthew G. Knepley PetscFunctionBegin; 15703f6b16c7SMatthew G. Knepley ierr = PetscFEGetBasisSpace(fe, &P);CHKERRQ(ierr); 15713f6b16c7SMatthew G. Knepley ierr = PetscFEGetDualSpace(fe, &Q);CHKERRQ(ierr); 15723f6b16c7SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) fe, name);CHKERRQ(ierr); 15733f6b16c7SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) P, name);CHKERRQ(ierr); 15743f6b16c7SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) Q, name);CHKERRQ(ierr); 15753f6b16c7SMatthew G. Knepley PetscFunctionReturn(0); 15763f6b16c7SMatthew G. Knepley } 1577