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 17620cf1dd8SToby Isaac - v - the viewer 17720cf1dd8SToby Isaac 17820cf1dd8SToby Isaac Level: developer 17920cf1dd8SToby Isaac 18020cf1dd8SToby Isaac .seealso PetscFEDestroy() 18120cf1dd8SToby Isaac @*/ 18220cf1dd8SToby Isaac PetscErrorCode PetscFEView(PetscFE fem, PetscViewer v) 18320cf1dd8SToby Isaac { 18420cf1dd8SToby Isaac PetscErrorCode ierr; 18520cf1dd8SToby Isaac 18620cf1dd8SToby Isaac PetscFunctionBegin; 18720cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 18820cf1dd8SToby Isaac if (!v) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) fem), &v);CHKERRQ(ierr);} 18920cf1dd8SToby Isaac if (fem->ops->view) {ierr = (*fem->ops->view)(fem, v);CHKERRQ(ierr);} 19020cf1dd8SToby Isaac PetscFunctionReturn(0); 19120cf1dd8SToby Isaac } 19220cf1dd8SToby Isaac 19320cf1dd8SToby Isaac /*@ 19420cf1dd8SToby Isaac PetscFESetFromOptions - sets parameters in a PetscFE from the options database 19520cf1dd8SToby Isaac 19620cf1dd8SToby Isaac Collective on PetscFE 19720cf1dd8SToby Isaac 19820cf1dd8SToby Isaac Input Parameter: 19920cf1dd8SToby Isaac . fem - the PetscFE object to set options for 20020cf1dd8SToby Isaac 20120cf1dd8SToby Isaac Options Database: 20220cf1dd8SToby Isaac . -petscfe_num_blocks the number of cell blocks to integrate concurrently 20320cf1dd8SToby Isaac . -petscfe_num_batches the number of cell batches to integrate serially 20420cf1dd8SToby Isaac 20520cf1dd8SToby Isaac Level: developer 20620cf1dd8SToby Isaac 20720cf1dd8SToby Isaac .seealso PetscFEView() 20820cf1dd8SToby Isaac @*/ 20920cf1dd8SToby Isaac PetscErrorCode PetscFESetFromOptions(PetscFE fem) 21020cf1dd8SToby Isaac { 21120cf1dd8SToby Isaac const char *defaultType; 21220cf1dd8SToby Isaac char name[256]; 21320cf1dd8SToby Isaac PetscBool flg; 21420cf1dd8SToby Isaac PetscErrorCode ierr; 21520cf1dd8SToby Isaac 21620cf1dd8SToby Isaac PetscFunctionBegin; 21720cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 21820cf1dd8SToby Isaac if (!((PetscObject) fem)->type_name) { 21920cf1dd8SToby Isaac defaultType = PETSCFEBASIC; 22020cf1dd8SToby Isaac } else { 22120cf1dd8SToby Isaac defaultType = ((PetscObject) fem)->type_name; 22220cf1dd8SToby Isaac } 22320cf1dd8SToby Isaac if (!PetscFERegisterAllCalled) {ierr = PetscFERegisterAll();CHKERRQ(ierr);} 22420cf1dd8SToby Isaac 22520cf1dd8SToby Isaac ierr = PetscObjectOptionsBegin((PetscObject) fem);CHKERRQ(ierr); 22620cf1dd8SToby Isaac ierr = PetscOptionsFList("-petscfe_type", "Finite element space", "PetscFESetType", PetscFEList, defaultType, name, 256, &flg);CHKERRQ(ierr); 22720cf1dd8SToby Isaac if (flg) { 22820cf1dd8SToby Isaac ierr = PetscFESetType(fem, name);CHKERRQ(ierr); 22920cf1dd8SToby Isaac } else if (!((PetscObject) fem)->type_name) { 23020cf1dd8SToby Isaac ierr = PetscFESetType(fem, defaultType);CHKERRQ(ierr); 23120cf1dd8SToby Isaac } 23220cf1dd8SToby Isaac ierr = PetscOptionsInt("-petscfe_num_blocks", "The number of cell blocks to integrate concurrently", "PetscSpaceSetTileSizes", fem->numBlocks, &fem->numBlocks, NULL);CHKERRQ(ierr); 23320cf1dd8SToby Isaac ierr = PetscOptionsInt("-petscfe_num_batches", "The number of cell batches to integrate serially", "PetscSpaceSetTileSizes", fem->numBatches, &fem->numBatches, NULL);CHKERRQ(ierr); 23420cf1dd8SToby Isaac if (fem->ops->setfromoptions) { 23520cf1dd8SToby Isaac ierr = (*fem->ops->setfromoptions)(PetscOptionsObject,fem);CHKERRQ(ierr); 23620cf1dd8SToby Isaac } 23720cf1dd8SToby Isaac /* process any options handlers added with PetscObjectAddOptionsHandler() */ 23820cf1dd8SToby Isaac ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) fem);CHKERRQ(ierr); 23920cf1dd8SToby Isaac ierr = PetscOptionsEnd();CHKERRQ(ierr); 24020cf1dd8SToby Isaac ierr = PetscFEViewFromOptions(fem, NULL, "-petscfe_view");CHKERRQ(ierr); 24120cf1dd8SToby Isaac PetscFunctionReturn(0); 24220cf1dd8SToby Isaac } 24320cf1dd8SToby Isaac 24420cf1dd8SToby Isaac /*@C 24520cf1dd8SToby Isaac PetscFESetUp - Construct data structures for the PetscFE 24620cf1dd8SToby Isaac 24720cf1dd8SToby Isaac Collective on PetscFE 24820cf1dd8SToby Isaac 24920cf1dd8SToby Isaac Input Parameter: 25020cf1dd8SToby Isaac . fem - the PetscFE object to setup 25120cf1dd8SToby Isaac 25220cf1dd8SToby Isaac Level: developer 25320cf1dd8SToby Isaac 25420cf1dd8SToby Isaac .seealso PetscFEView(), PetscFEDestroy() 25520cf1dd8SToby Isaac @*/ 25620cf1dd8SToby Isaac PetscErrorCode PetscFESetUp(PetscFE fem) 25720cf1dd8SToby Isaac { 25820cf1dd8SToby Isaac PetscErrorCode ierr; 25920cf1dd8SToby Isaac 26020cf1dd8SToby Isaac PetscFunctionBegin; 26120cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 26220cf1dd8SToby Isaac if (fem->setupcalled) PetscFunctionReturn(0); 26320cf1dd8SToby Isaac fem->setupcalled = PETSC_TRUE; 26420cf1dd8SToby Isaac if (fem->ops->setup) {ierr = (*fem->ops->setup)(fem);CHKERRQ(ierr);} 26520cf1dd8SToby Isaac PetscFunctionReturn(0); 26620cf1dd8SToby Isaac } 26720cf1dd8SToby Isaac 26820cf1dd8SToby Isaac /*@ 26920cf1dd8SToby Isaac PetscFEDestroy - Destroys a PetscFE object 27020cf1dd8SToby Isaac 27120cf1dd8SToby Isaac Collective on PetscFE 27220cf1dd8SToby Isaac 27320cf1dd8SToby Isaac Input Parameter: 27420cf1dd8SToby Isaac . fem - the PetscFE object to destroy 27520cf1dd8SToby Isaac 27620cf1dd8SToby Isaac Level: developer 27720cf1dd8SToby Isaac 27820cf1dd8SToby Isaac .seealso PetscFEView() 27920cf1dd8SToby Isaac @*/ 28020cf1dd8SToby Isaac PetscErrorCode PetscFEDestroy(PetscFE *fem) 28120cf1dd8SToby Isaac { 28220cf1dd8SToby Isaac PetscErrorCode ierr; 28320cf1dd8SToby Isaac 28420cf1dd8SToby Isaac PetscFunctionBegin; 28520cf1dd8SToby Isaac if (!*fem) PetscFunctionReturn(0); 28620cf1dd8SToby Isaac PetscValidHeaderSpecific((*fem), PETSCFE_CLASSID, 1); 28720cf1dd8SToby Isaac 28820cf1dd8SToby Isaac if (--((PetscObject)(*fem))->refct > 0) {*fem = 0; PetscFunctionReturn(0);} 28920cf1dd8SToby Isaac ((PetscObject) (*fem))->refct = 0; 29020cf1dd8SToby Isaac 29120cf1dd8SToby Isaac if ((*fem)->subspaces) { 29220cf1dd8SToby Isaac PetscInt dim, d; 29320cf1dd8SToby Isaac 29420cf1dd8SToby Isaac ierr = PetscDualSpaceGetDimension((*fem)->dualSpace, &dim);CHKERRQ(ierr); 29520cf1dd8SToby Isaac for (d = 0; d < dim; ++d) {ierr = PetscFEDestroy(&(*fem)->subspaces[d]);CHKERRQ(ierr);} 29620cf1dd8SToby Isaac } 29720cf1dd8SToby Isaac ierr = PetscFree((*fem)->subspaces);CHKERRQ(ierr); 29820cf1dd8SToby Isaac ierr = PetscFree((*fem)->invV);CHKERRQ(ierr); 29920cf1dd8SToby Isaac ierr = PetscFERestoreTabulation((*fem), 0, NULL, &(*fem)->B, &(*fem)->D, NULL /*&(*fem)->H*/);CHKERRQ(ierr); 30020cf1dd8SToby Isaac ierr = PetscFERestoreTabulation((*fem), 0, NULL, &(*fem)->Bf, &(*fem)->Df, NULL /*&(*fem)->Hf*/);CHKERRQ(ierr); 30120cf1dd8SToby Isaac ierr = PetscFERestoreTabulation((*fem), 0, NULL, &(*fem)->F, NULL, NULL);CHKERRQ(ierr); 30220cf1dd8SToby Isaac ierr = PetscSpaceDestroy(&(*fem)->basisSpace);CHKERRQ(ierr); 30320cf1dd8SToby Isaac ierr = PetscDualSpaceDestroy(&(*fem)->dualSpace);CHKERRQ(ierr); 30420cf1dd8SToby Isaac ierr = PetscQuadratureDestroy(&(*fem)->quadrature);CHKERRQ(ierr); 30520cf1dd8SToby Isaac ierr = PetscQuadratureDestroy(&(*fem)->faceQuadrature);CHKERRQ(ierr); 30620cf1dd8SToby Isaac 30720cf1dd8SToby Isaac if ((*fem)->ops->destroy) {ierr = (*(*fem)->ops->destroy)(*fem);CHKERRQ(ierr);} 30820cf1dd8SToby Isaac ierr = PetscHeaderDestroy(fem);CHKERRQ(ierr); 30920cf1dd8SToby Isaac PetscFunctionReturn(0); 31020cf1dd8SToby Isaac } 31120cf1dd8SToby Isaac 31220cf1dd8SToby Isaac /*@ 31320cf1dd8SToby Isaac PetscFECreate - Creates an empty PetscFE object. The type can then be set with PetscFESetType(). 31420cf1dd8SToby Isaac 31520cf1dd8SToby Isaac Collective on MPI_Comm 31620cf1dd8SToby Isaac 31720cf1dd8SToby Isaac Input Parameter: 31820cf1dd8SToby Isaac . comm - The communicator for the PetscFE object 31920cf1dd8SToby Isaac 32020cf1dd8SToby Isaac Output Parameter: 32120cf1dd8SToby Isaac . fem - The PetscFE object 32220cf1dd8SToby Isaac 32320cf1dd8SToby Isaac Level: beginner 32420cf1dd8SToby Isaac 32520cf1dd8SToby Isaac .seealso: PetscFESetType(), PETSCFEGALERKIN 32620cf1dd8SToby Isaac @*/ 32720cf1dd8SToby Isaac PetscErrorCode PetscFECreate(MPI_Comm comm, PetscFE *fem) 32820cf1dd8SToby Isaac { 32920cf1dd8SToby Isaac PetscFE f; 33020cf1dd8SToby Isaac PetscErrorCode ierr; 33120cf1dd8SToby Isaac 33220cf1dd8SToby Isaac PetscFunctionBegin; 33320cf1dd8SToby Isaac PetscValidPointer(fem, 2); 33420cf1dd8SToby Isaac ierr = PetscCitationsRegister(FECitation,&FEcite);CHKERRQ(ierr); 33520cf1dd8SToby Isaac *fem = NULL; 33620cf1dd8SToby Isaac ierr = PetscFEInitializePackage();CHKERRQ(ierr); 33720cf1dd8SToby Isaac 33820cf1dd8SToby Isaac ierr = PetscHeaderCreate(f, PETSCFE_CLASSID, "PetscFE", "Finite Element", "PetscFE", comm, PetscFEDestroy, PetscFEView);CHKERRQ(ierr); 33920cf1dd8SToby Isaac 34020cf1dd8SToby Isaac f->basisSpace = NULL; 34120cf1dd8SToby Isaac f->dualSpace = NULL; 34220cf1dd8SToby Isaac f->numComponents = 1; 34320cf1dd8SToby Isaac f->subspaces = NULL; 34420cf1dd8SToby Isaac f->invV = NULL; 34520cf1dd8SToby Isaac f->B = NULL; 34620cf1dd8SToby Isaac f->D = NULL; 34720cf1dd8SToby Isaac f->H = NULL; 34820cf1dd8SToby Isaac f->Bf = NULL; 34920cf1dd8SToby Isaac f->Df = NULL; 35020cf1dd8SToby Isaac f->Hf = NULL; 35120cf1dd8SToby Isaac ierr = PetscMemzero(&f->quadrature, sizeof(PetscQuadrature));CHKERRQ(ierr); 35220cf1dd8SToby Isaac ierr = PetscMemzero(&f->faceQuadrature, sizeof(PetscQuadrature));CHKERRQ(ierr); 35320cf1dd8SToby Isaac f->blockSize = 0; 35420cf1dd8SToby Isaac f->numBlocks = 1; 35520cf1dd8SToby Isaac f->batchSize = 0; 35620cf1dd8SToby Isaac f->numBatches = 1; 35720cf1dd8SToby Isaac 35820cf1dd8SToby Isaac *fem = f; 35920cf1dd8SToby Isaac PetscFunctionReturn(0); 36020cf1dd8SToby Isaac } 36120cf1dd8SToby Isaac 36220cf1dd8SToby Isaac /*@ 36320cf1dd8SToby Isaac PetscFEGetSpatialDimension - Returns the spatial dimension of the element 36420cf1dd8SToby Isaac 36520cf1dd8SToby Isaac Not collective 36620cf1dd8SToby Isaac 36720cf1dd8SToby Isaac Input Parameter: 36820cf1dd8SToby Isaac . fem - The PetscFE object 36920cf1dd8SToby Isaac 37020cf1dd8SToby Isaac Output Parameter: 37120cf1dd8SToby Isaac . dim - The spatial dimension 37220cf1dd8SToby Isaac 37320cf1dd8SToby Isaac Level: intermediate 37420cf1dd8SToby Isaac 37520cf1dd8SToby Isaac .seealso: PetscFECreate() 37620cf1dd8SToby Isaac @*/ 37720cf1dd8SToby Isaac PetscErrorCode PetscFEGetSpatialDimension(PetscFE fem, PetscInt *dim) 37820cf1dd8SToby Isaac { 37920cf1dd8SToby Isaac DM dm; 38020cf1dd8SToby Isaac PetscErrorCode ierr; 38120cf1dd8SToby Isaac 38220cf1dd8SToby Isaac PetscFunctionBegin; 38320cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 38420cf1dd8SToby Isaac PetscValidPointer(dim, 2); 38520cf1dd8SToby Isaac ierr = PetscDualSpaceGetDM(fem->dualSpace, &dm);CHKERRQ(ierr); 38620cf1dd8SToby Isaac ierr = DMGetDimension(dm, dim);CHKERRQ(ierr); 38720cf1dd8SToby Isaac PetscFunctionReturn(0); 38820cf1dd8SToby Isaac } 38920cf1dd8SToby Isaac 39020cf1dd8SToby Isaac /*@ 39120cf1dd8SToby Isaac PetscFESetNumComponents - Sets the number of components in the element 39220cf1dd8SToby Isaac 39320cf1dd8SToby Isaac Not collective 39420cf1dd8SToby Isaac 39520cf1dd8SToby Isaac Input Parameters: 39620cf1dd8SToby Isaac + fem - The PetscFE object 39720cf1dd8SToby Isaac - comp - The number of field components 39820cf1dd8SToby Isaac 39920cf1dd8SToby Isaac Level: intermediate 40020cf1dd8SToby Isaac 40120cf1dd8SToby Isaac .seealso: PetscFECreate() 40220cf1dd8SToby Isaac @*/ 40320cf1dd8SToby Isaac PetscErrorCode PetscFESetNumComponents(PetscFE fem, PetscInt comp) 40420cf1dd8SToby Isaac { 40520cf1dd8SToby Isaac PetscFunctionBegin; 40620cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 40720cf1dd8SToby Isaac fem->numComponents = comp; 40820cf1dd8SToby Isaac PetscFunctionReturn(0); 40920cf1dd8SToby Isaac } 41020cf1dd8SToby Isaac 41120cf1dd8SToby Isaac /*@ 41220cf1dd8SToby Isaac PetscFEGetNumComponents - Returns the number of components in the element 41320cf1dd8SToby Isaac 41420cf1dd8SToby Isaac Not collective 41520cf1dd8SToby Isaac 41620cf1dd8SToby Isaac Input Parameter: 41720cf1dd8SToby Isaac . fem - The PetscFE object 41820cf1dd8SToby Isaac 41920cf1dd8SToby Isaac Output Parameter: 42020cf1dd8SToby Isaac . comp - The number of field components 42120cf1dd8SToby Isaac 42220cf1dd8SToby Isaac Level: intermediate 42320cf1dd8SToby Isaac 42420cf1dd8SToby Isaac .seealso: PetscFECreate() 42520cf1dd8SToby Isaac @*/ 42620cf1dd8SToby Isaac PetscErrorCode PetscFEGetNumComponents(PetscFE fem, PetscInt *comp) 42720cf1dd8SToby Isaac { 42820cf1dd8SToby Isaac PetscFunctionBegin; 42920cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 43020cf1dd8SToby Isaac PetscValidPointer(comp, 2); 43120cf1dd8SToby Isaac *comp = fem->numComponents; 43220cf1dd8SToby Isaac PetscFunctionReturn(0); 43320cf1dd8SToby Isaac } 43420cf1dd8SToby Isaac 43520cf1dd8SToby Isaac /*@ 43620cf1dd8SToby Isaac PetscFESetTileSizes - Sets the tile sizes for evaluation 43720cf1dd8SToby Isaac 43820cf1dd8SToby Isaac Not collective 43920cf1dd8SToby Isaac 44020cf1dd8SToby Isaac Input Parameters: 44120cf1dd8SToby Isaac + fem - The PetscFE object 44220cf1dd8SToby Isaac . blockSize - The number of elements in a block 44320cf1dd8SToby Isaac . numBlocks - The number of blocks in a batch 44420cf1dd8SToby Isaac . batchSize - The number of elements in a batch 44520cf1dd8SToby Isaac - numBatches - The number of batches in a chunk 44620cf1dd8SToby Isaac 44720cf1dd8SToby Isaac Level: intermediate 44820cf1dd8SToby Isaac 44920cf1dd8SToby Isaac .seealso: PetscFECreate() 45020cf1dd8SToby Isaac @*/ 45120cf1dd8SToby Isaac PetscErrorCode PetscFESetTileSizes(PetscFE fem, PetscInt blockSize, PetscInt numBlocks, PetscInt batchSize, PetscInt numBatches) 45220cf1dd8SToby Isaac { 45320cf1dd8SToby Isaac PetscFunctionBegin; 45420cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 45520cf1dd8SToby Isaac fem->blockSize = blockSize; 45620cf1dd8SToby Isaac fem->numBlocks = numBlocks; 45720cf1dd8SToby Isaac fem->batchSize = batchSize; 45820cf1dd8SToby Isaac fem->numBatches = numBatches; 45920cf1dd8SToby Isaac PetscFunctionReturn(0); 46020cf1dd8SToby Isaac } 46120cf1dd8SToby Isaac 46220cf1dd8SToby Isaac /*@ 46320cf1dd8SToby Isaac PetscFEGetTileSizes - Returns the tile sizes for evaluation 46420cf1dd8SToby Isaac 46520cf1dd8SToby Isaac Not collective 46620cf1dd8SToby Isaac 46720cf1dd8SToby Isaac Input Parameter: 46820cf1dd8SToby Isaac . fem - The PetscFE object 46920cf1dd8SToby Isaac 47020cf1dd8SToby Isaac Output Parameters: 47120cf1dd8SToby Isaac + blockSize - The number of elements in a block 47220cf1dd8SToby Isaac . numBlocks - The number of blocks in a batch 47320cf1dd8SToby Isaac . batchSize - The number of elements in a batch 47420cf1dd8SToby Isaac - numBatches - The number of batches in a chunk 47520cf1dd8SToby Isaac 47620cf1dd8SToby Isaac Level: intermediate 47720cf1dd8SToby Isaac 47820cf1dd8SToby Isaac .seealso: PetscFECreate() 47920cf1dd8SToby Isaac @*/ 48020cf1dd8SToby Isaac PetscErrorCode PetscFEGetTileSizes(PetscFE fem, PetscInt *blockSize, PetscInt *numBlocks, PetscInt *batchSize, PetscInt *numBatches) 48120cf1dd8SToby Isaac { 48220cf1dd8SToby Isaac PetscFunctionBegin; 48320cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 48420cf1dd8SToby Isaac if (blockSize) PetscValidPointer(blockSize, 2); 48520cf1dd8SToby Isaac if (numBlocks) PetscValidPointer(numBlocks, 3); 48620cf1dd8SToby Isaac if (batchSize) PetscValidPointer(batchSize, 4); 48720cf1dd8SToby Isaac if (numBatches) PetscValidPointer(numBatches, 5); 48820cf1dd8SToby Isaac if (blockSize) *blockSize = fem->blockSize; 48920cf1dd8SToby Isaac if (numBlocks) *numBlocks = fem->numBlocks; 49020cf1dd8SToby Isaac if (batchSize) *batchSize = fem->batchSize; 49120cf1dd8SToby Isaac if (numBatches) *numBatches = fem->numBatches; 49220cf1dd8SToby Isaac PetscFunctionReturn(0); 49320cf1dd8SToby Isaac } 49420cf1dd8SToby Isaac 49520cf1dd8SToby Isaac /*@ 49620cf1dd8SToby Isaac PetscFEGetBasisSpace - Returns the PetscSpace used for approximation of the solution 49720cf1dd8SToby Isaac 49820cf1dd8SToby Isaac Not collective 49920cf1dd8SToby Isaac 50020cf1dd8SToby Isaac Input Parameter: 50120cf1dd8SToby Isaac . fem - The PetscFE object 50220cf1dd8SToby Isaac 50320cf1dd8SToby Isaac Output Parameter: 50420cf1dd8SToby Isaac . sp - The PetscSpace object 50520cf1dd8SToby Isaac 50620cf1dd8SToby Isaac Level: intermediate 50720cf1dd8SToby Isaac 50820cf1dd8SToby Isaac .seealso: PetscFECreate() 50920cf1dd8SToby Isaac @*/ 51020cf1dd8SToby Isaac PetscErrorCode PetscFEGetBasisSpace(PetscFE fem, PetscSpace *sp) 51120cf1dd8SToby Isaac { 51220cf1dd8SToby Isaac PetscFunctionBegin; 51320cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 51420cf1dd8SToby Isaac PetscValidPointer(sp, 2); 51520cf1dd8SToby Isaac *sp = fem->basisSpace; 51620cf1dd8SToby Isaac PetscFunctionReturn(0); 51720cf1dd8SToby Isaac } 51820cf1dd8SToby Isaac 51920cf1dd8SToby Isaac /*@ 52020cf1dd8SToby Isaac PetscFESetBasisSpace - Sets the PetscSpace used for approximation of the solution 52120cf1dd8SToby Isaac 52220cf1dd8SToby Isaac Not collective 52320cf1dd8SToby Isaac 52420cf1dd8SToby Isaac Input Parameters: 52520cf1dd8SToby Isaac + fem - The PetscFE object 52620cf1dd8SToby Isaac - sp - The PetscSpace object 52720cf1dd8SToby Isaac 52820cf1dd8SToby Isaac Level: intermediate 52920cf1dd8SToby Isaac 53020cf1dd8SToby Isaac .seealso: PetscFECreate() 53120cf1dd8SToby Isaac @*/ 53220cf1dd8SToby Isaac PetscErrorCode PetscFESetBasisSpace(PetscFE fem, PetscSpace sp) 53320cf1dd8SToby Isaac { 53420cf1dd8SToby Isaac PetscErrorCode ierr; 53520cf1dd8SToby Isaac 53620cf1dd8SToby Isaac PetscFunctionBegin; 53720cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 53820cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 2); 53920cf1dd8SToby Isaac ierr = PetscSpaceDestroy(&fem->basisSpace);CHKERRQ(ierr); 54020cf1dd8SToby Isaac fem->basisSpace = sp; 54120cf1dd8SToby Isaac ierr = PetscObjectReference((PetscObject) fem->basisSpace);CHKERRQ(ierr); 54220cf1dd8SToby Isaac PetscFunctionReturn(0); 54320cf1dd8SToby Isaac } 54420cf1dd8SToby Isaac 54520cf1dd8SToby Isaac /*@ 54620cf1dd8SToby Isaac PetscFEGetDualSpace - Returns the PetscDualSpace used to define the inner product 54720cf1dd8SToby Isaac 54820cf1dd8SToby Isaac Not collective 54920cf1dd8SToby Isaac 55020cf1dd8SToby Isaac Input Parameter: 55120cf1dd8SToby Isaac . fem - The PetscFE object 55220cf1dd8SToby Isaac 55320cf1dd8SToby Isaac Output Parameter: 55420cf1dd8SToby Isaac . sp - The PetscDualSpace object 55520cf1dd8SToby Isaac 55620cf1dd8SToby Isaac Level: intermediate 55720cf1dd8SToby Isaac 55820cf1dd8SToby Isaac .seealso: PetscFECreate() 55920cf1dd8SToby Isaac @*/ 56020cf1dd8SToby Isaac PetscErrorCode PetscFEGetDualSpace(PetscFE fem, PetscDualSpace *sp) 56120cf1dd8SToby Isaac { 56220cf1dd8SToby Isaac PetscFunctionBegin; 56320cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 56420cf1dd8SToby Isaac PetscValidPointer(sp, 2); 56520cf1dd8SToby Isaac *sp = fem->dualSpace; 56620cf1dd8SToby Isaac PetscFunctionReturn(0); 56720cf1dd8SToby Isaac } 56820cf1dd8SToby Isaac 56920cf1dd8SToby Isaac /*@ 57020cf1dd8SToby Isaac PetscFESetDualSpace - Sets the PetscDualSpace used to define the inner product 57120cf1dd8SToby Isaac 57220cf1dd8SToby Isaac Not collective 57320cf1dd8SToby Isaac 57420cf1dd8SToby Isaac Input Parameters: 57520cf1dd8SToby Isaac + fem - The PetscFE object 57620cf1dd8SToby Isaac - sp - The PetscDualSpace object 57720cf1dd8SToby Isaac 57820cf1dd8SToby Isaac Level: intermediate 57920cf1dd8SToby Isaac 58020cf1dd8SToby Isaac .seealso: PetscFECreate() 58120cf1dd8SToby Isaac @*/ 58220cf1dd8SToby Isaac PetscErrorCode PetscFESetDualSpace(PetscFE fem, PetscDualSpace sp) 58320cf1dd8SToby Isaac { 58420cf1dd8SToby Isaac PetscErrorCode ierr; 58520cf1dd8SToby Isaac 58620cf1dd8SToby Isaac PetscFunctionBegin; 58720cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 58820cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 2); 58920cf1dd8SToby Isaac ierr = PetscDualSpaceDestroy(&fem->dualSpace);CHKERRQ(ierr); 59020cf1dd8SToby Isaac fem->dualSpace = sp; 59120cf1dd8SToby Isaac ierr = PetscObjectReference((PetscObject) fem->dualSpace);CHKERRQ(ierr); 59220cf1dd8SToby Isaac PetscFunctionReturn(0); 59320cf1dd8SToby Isaac } 59420cf1dd8SToby Isaac 59520cf1dd8SToby Isaac /*@ 59620cf1dd8SToby Isaac PetscFEGetQuadrature - Returns the PetscQuadrature used to calculate inner products 59720cf1dd8SToby Isaac 59820cf1dd8SToby Isaac Not collective 59920cf1dd8SToby Isaac 60020cf1dd8SToby Isaac Input Parameter: 60120cf1dd8SToby Isaac . fem - The PetscFE object 60220cf1dd8SToby Isaac 60320cf1dd8SToby Isaac Output Parameter: 60420cf1dd8SToby Isaac . q - The PetscQuadrature object 60520cf1dd8SToby Isaac 60620cf1dd8SToby Isaac Level: intermediate 60720cf1dd8SToby Isaac 60820cf1dd8SToby Isaac .seealso: PetscFECreate() 60920cf1dd8SToby Isaac @*/ 61020cf1dd8SToby Isaac PetscErrorCode PetscFEGetQuadrature(PetscFE fem, PetscQuadrature *q) 61120cf1dd8SToby Isaac { 61220cf1dd8SToby Isaac PetscFunctionBegin; 61320cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 61420cf1dd8SToby Isaac PetscValidPointer(q, 2); 61520cf1dd8SToby Isaac *q = fem->quadrature; 61620cf1dd8SToby Isaac PetscFunctionReturn(0); 61720cf1dd8SToby Isaac } 61820cf1dd8SToby Isaac 61920cf1dd8SToby Isaac /*@ 62020cf1dd8SToby Isaac PetscFESetQuadrature - Sets the PetscQuadrature used to calculate inner products 62120cf1dd8SToby Isaac 62220cf1dd8SToby Isaac Not collective 62320cf1dd8SToby Isaac 62420cf1dd8SToby Isaac Input Parameters: 62520cf1dd8SToby Isaac + fem - The PetscFE object 62620cf1dd8SToby Isaac - q - The PetscQuadrature object 62720cf1dd8SToby Isaac 62820cf1dd8SToby Isaac Level: intermediate 62920cf1dd8SToby Isaac 63020cf1dd8SToby Isaac .seealso: PetscFECreate() 63120cf1dd8SToby Isaac @*/ 63220cf1dd8SToby Isaac PetscErrorCode PetscFESetQuadrature(PetscFE fem, PetscQuadrature q) 63320cf1dd8SToby Isaac { 63420cf1dd8SToby Isaac PetscInt Nc, qNc; 63520cf1dd8SToby Isaac PetscErrorCode ierr; 63620cf1dd8SToby Isaac 63720cf1dd8SToby Isaac PetscFunctionBegin; 63820cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 63920cf1dd8SToby Isaac ierr = PetscFEGetNumComponents(fem, &Nc);CHKERRQ(ierr); 64020cf1dd8SToby Isaac ierr = PetscQuadratureGetNumComponents(q, &qNc);CHKERRQ(ierr); 64120cf1dd8SToby 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); 64220cf1dd8SToby Isaac ierr = PetscFERestoreTabulation(fem, 0, NULL, &fem->B, &fem->D, NULL /*&(*fem)->H*/);CHKERRQ(ierr); 64320cf1dd8SToby Isaac ierr = PetscQuadratureDestroy(&fem->quadrature);CHKERRQ(ierr); 64420cf1dd8SToby Isaac fem->quadrature = q; 64520cf1dd8SToby Isaac ierr = PetscObjectReference((PetscObject) q);CHKERRQ(ierr); 64620cf1dd8SToby Isaac PetscFunctionReturn(0); 64720cf1dd8SToby Isaac } 64820cf1dd8SToby Isaac 64920cf1dd8SToby Isaac /*@ 65020cf1dd8SToby Isaac PetscFEGetFaceQuadrature - Returns the PetscQuadrature used to calculate inner products on faces 65120cf1dd8SToby Isaac 65220cf1dd8SToby Isaac Not collective 65320cf1dd8SToby Isaac 65420cf1dd8SToby Isaac Input Parameter: 65520cf1dd8SToby Isaac . fem - The PetscFE object 65620cf1dd8SToby Isaac 65720cf1dd8SToby Isaac Output Parameter: 65820cf1dd8SToby Isaac . q - The PetscQuadrature object 65920cf1dd8SToby Isaac 66020cf1dd8SToby Isaac Level: intermediate 66120cf1dd8SToby Isaac 66220cf1dd8SToby Isaac .seealso: PetscFECreate() 66320cf1dd8SToby Isaac @*/ 66420cf1dd8SToby Isaac PetscErrorCode PetscFEGetFaceQuadrature(PetscFE fem, PetscQuadrature *q) 66520cf1dd8SToby Isaac { 66620cf1dd8SToby Isaac PetscFunctionBegin; 66720cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 66820cf1dd8SToby Isaac PetscValidPointer(q, 2); 66920cf1dd8SToby Isaac *q = fem->faceQuadrature; 67020cf1dd8SToby Isaac PetscFunctionReturn(0); 67120cf1dd8SToby Isaac } 67220cf1dd8SToby Isaac 67320cf1dd8SToby Isaac /*@ 67420cf1dd8SToby Isaac PetscFESetFaceQuadrature - Sets the PetscQuadrature used to calculate inner products on faces 67520cf1dd8SToby Isaac 67620cf1dd8SToby Isaac Not collective 67720cf1dd8SToby Isaac 67820cf1dd8SToby Isaac Input Parameters: 67920cf1dd8SToby Isaac + fem - The PetscFE object 68020cf1dd8SToby Isaac - q - The PetscQuadrature object 68120cf1dd8SToby Isaac 68220cf1dd8SToby Isaac Level: intermediate 68320cf1dd8SToby Isaac 68420cf1dd8SToby Isaac .seealso: PetscFECreate() 68520cf1dd8SToby Isaac @*/ 68620cf1dd8SToby Isaac PetscErrorCode PetscFESetFaceQuadrature(PetscFE fem, PetscQuadrature q) 68720cf1dd8SToby Isaac { 68820cf1dd8SToby Isaac PetscErrorCode ierr; 68920cf1dd8SToby Isaac 69020cf1dd8SToby Isaac PetscFunctionBegin; 69120cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 69220cf1dd8SToby Isaac ierr = PetscFERestoreTabulation(fem, 0, NULL, &fem->Bf, &fem->Df, NULL /*&(*fem)->Hf*/);CHKERRQ(ierr); 69320cf1dd8SToby Isaac ierr = PetscQuadratureDestroy(&fem->faceQuadrature);CHKERRQ(ierr); 69420cf1dd8SToby Isaac fem->faceQuadrature = q; 69520cf1dd8SToby Isaac ierr = PetscObjectReference((PetscObject) q);CHKERRQ(ierr); 69620cf1dd8SToby Isaac PetscFunctionReturn(0); 69720cf1dd8SToby Isaac } 69820cf1dd8SToby Isaac 69920cf1dd8SToby Isaac /*@C 70020cf1dd8SToby Isaac PetscFEGetNumDof - Returns the number of dofs (dual basis vectors) associated to mesh points on the reference cell of a given dimension 70120cf1dd8SToby Isaac 70220cf1dd8SToby Isaac Not collective 70320cf1dd8SToby Isaac 70420cf1dd8SToby Isaac Input Parameter: 70520cf1dd8SToby Isaac . fem - The PetscFE object 70620cf1dd8SToby Isaac 70720cf1dd8SToby Isaac Output Parameter: 70820cf1dd8SToby Isaac . numDof - Array with the number of dofs per dimension 70920cf1dd8SToby Isaac 71020cf1dd8SToby Isaac Level: intermediate 71120cf1dd8SToby Isaac 71220cf1dd8SToby Isaac .seealso: PetscFECreate() 71320cf1dd8SToby Isaac @*/ 71420cf1dd8SToby Isaac PetscErrorCode PetscFEGetNumDof(PetscFE fem, const PetscInt **numDof) 71520cf1dd8SToby Isaac { 71620cf1dd8SToby Isaac PetscErrorCode ierr; 71720cf1dd8SToby Isaac 71820cf1dd8SToby Isaac PetscFunctionBegin; 71920cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 72020cf1dd8SToby Isaac PetscValidPointer(numDof, 2); 72120cf1dd8SToby Isaac ierr = PetscDualSpaceGetNumDof(fem->dualSpace, numDof);CHKERRQ(ierr); 72220cf1dd8SToby Isaac PetscFunctionReturn(0); 72320cf1dd8SToby Isaac } 72420cf1dd8SToby Isaac 72520cf1dd8SToby Isaac /*@C 72620cf1dd8SToby Isaac PetscFEGetDefaultTabulation - Returns the tabulation of the basis functions at the quadrature points 72720cf1dd8SToby Isaac 72820cf1dd8SToby Isaac Not collective 72920cf1dd8SToby Isaac 73020cf1dd8SToby Isaac Input Parameter: 73120cf1dd8SToby Isaac . fem - The PetscFE object 73220cf1dd8SToby Isaac 73320cf1dd8SToby Isaac Output Parameters: 73420cf1dd8SToby Isaac + B - The basis function values at quadrature points 73520cf1dd8SToby Isaac . D - The basis function derivatives at quadrature points 73620cf1dd8SToby Isaac - H - The basis function second derivatives at quadrature points 73720cf1dd8SToby Isaac 73820cf1dd8SToby Isaac Note: 73920cf1dd8SToby Isaac $ B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c 74020cf1dd8SToby 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 74120cf1dd8SToby 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 74220cf1dd8SToby Isaac 74320cf1dd8SToby Isaac Level: intermediate 74420cf1dd8SToby Isaac 74520cf1dd8SToby Isaac .seealso: PetscFEGetTabulation(), PetscFERestoreTabulation() 74620cf1dd8SToby Isaac @*/ 74720cf1dd8SToby Isaac PetscErrorCode PetscFEGetDefaultTabulation(PetscFE fem, PetscReal **B, PetscReal **D, PetscReal **H) 74820cf1dd8SToby Isaac { 74920cf1dd8SToby Isaac PetscInt npoints; 75020cf1dd8SToby Isaac const PetscReal *points; 75120cf1dd8SToby Isaac PetscErrorCode ierr; 75220cf1dd8SToby Isaac 75320cf1dd8SToby Isaac PetscFunctionBegin; 75420cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 75520cf1dd8SToby Isaac if (B) PetscValidPointer(B, 2); 75620cf1dd8SToby Isaac if (D) PetscValidPointer(D, 3); 75720cf1dd8SToby Isaac if (H) PetscValidPointer(H, 4); 75820cf1dd8SToby Isaac ierr = PetscQuadratureGetData(fem->quadrature, NULL, NULL, &npoints, &points, NULL);CHKERRQ(ierr); 75920cf1dd8SToby Isaac if (!fem->B) {ierr = PetscFEGetTabulation(fem, npoints, points, &fem->B, &fem->D, NULL/*&fem->H*/);CHKERRQ(ierr);} 76020cf1dd8SToby Isaac if (B) *B = fem->B; 76120cf1dd8SToby Isaac if (D) *D = fem->D; 76220cf1dd8SToby Isaac if (H) *H = fem->H; 76320cf1dd8SToby Isaac PetscFunctionReturn(0); 76420cf1dd8SToby Isaac } 76520cf1dd8SToby Isaac 76620cf1dd8SToby Isaac PetscErrorCode PetscFEGetFaceTabulation(PetscFE fem, PetscReal **Bf, PetscReal **Df, PetscReal **Hf) 76720cf1dd8SToby Isaac { 76820cf1dd8SToby Isaac PetscErrorCode ierr; 76920cf1dd8SToby Isaac 77020cf1dd8SToby Isaac PetscFunctionBegin; 77120cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 77220cf1dd8SToby Isaac if (Bf) PetscValidPointer(Bf, 2); 77320cf1dd8SToby Isaac if (Df) PetscValidPointer(Df, 3); 77420cf1dd8SToby Isaac if (Hf) PetscValidPointer(Hf, 4); 77520cf1dd8SToby Isaac if (!fem->Bf) { 77620cf1dd8SToby Isaac const PetscReal xi0[3] = {-1., -1., -1.}; 77720cf1dd8SToby Isaac PetscReal v0[3], J[9], detJ; 77820cf1dd8SToby Isaac PetscQuadrature fq; 77920cf1dd8SToby Isaac PetscDualSpace sp; 78020cf1dd8SToby Isaac DM dm; 78120cf1dd8SToby Isaac const PetscInt *faces; 78220cf1dd8SToby Isaac PetscInt dim, numFaces, f, npoints, q; 78320cf1dd8SToby Isaac const PetscReal *points; 78420cf1dd8SToby Isaac PetscReal *facePoints; 78520cf1dd8SToby Isaac 78620cf1dd8SToby Isaac ierr = PetscFEGetDualSpace(fem, &sp);CHKERRQ(ierr); 78720cf1dd8SToby Isaac ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr); 78820cf1dd8SToby Isaac ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 78920cf1dd8SToby Isaac ierr = DMPlexGetConeSize(dm, 0, &numFaces);CHKERRQ(ierr); 79020cf1dd8SToby Isaac ierr = DMPlexGetCone(dm, 0, &faces);CHKERRQ(ierr); 79120cf1dd8SToby Isaac ierr = PetscFEGetFaceQuadrature(fem, &fq);CHKERRQ(ierr); 79220cf1dd8SToby Isaac if (fq) { 79320cf1dd8SToby Isaac ierr = PetscQuadratureGetData(fq, NULL, NULL, &npoints, &points, NULL);CHKERRQ(ierr); 79420cf1dd8SToby Isaac ierr = PetscMalloc1(numFaces*npoints*dim, &facePoints);CHKERRQ(ierr); 79520cf1dd8SToby Isaac for (f = 0; f < numFaces; ++f) { 79620cf1dd8SToby Isaac ierr = DMPlexComputeCellGeometryFEM(dm, faces[f], NULL, v0, J, NULL, &detJ);CHKERRQ(ierr); 79720cf1dd8SToby Isaac for (q = 0; q < npoints; ++q) CoordinatesRefToReal(dim, dim-1, xi0, v0, J, &points[q*(dim-1)], &facePoints[(f*npoints+q)*dim]); 79820cf1dd8SToby Isaac } 79920cf1dd8SToby Isaac ierr = PetscFEGetTabulation(fem, numFaces*npoints, facePoints, &fem->Bf, &fem->Df, NULL/*&fem->Hf*/);CHKERRQ(ierr); 80020cf1dd8SToby Isaac ierr = PetscFree(facePoints);CHKERRQ(ierr); 80120cf1dd8SToby Isaac } 80220cf1dd8SToby Isaac } 80320cf1dd8SToby Isaac if (Bf) *Bf = fem->Bf; 80420cf1dd8SToby Isaac if (Df) *Df = fem->Df; 80520cf1dd8SToby Isaac if (Hf) *Hf = fem->Hf; 80620cf1dd8SToby Isaac PetscFunctionReturn(0); 80720cf1dd8SToby Isaac } 80820cf1dd8SToby Isaac 80920cf1dd8SToby Isaac PetscErrorCode PetscFEGetFaceCentroidTabulation(PetscFE fem, PetscReal **F) 81020cf1dd8SToby Isaac { 81120cf1dd8SToby Isaac PetscErrorCode ierr; 81220cf1dd8SToby Isaac 81320cf1dd8SToby Isaac PetscFunctionBegin; 81420cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 81520cf1dd8SToby Isaac PetscValidPointer(F, 2); 81620cf1dd8SToby Isaac if (!fem->F) { 81720cf1dd8SToby Isaac PetscDualSpace sp; 81820cf1dd8SToby Isaac DM dm; 81920cf1dd8SToby Isaac const PetscInt *cone; 82020cf1dd8SToby Isaac PetscReal *centroids; 82120cf1dd8SToby Isaac PetscInt dim, numFaces, f; 82220cf1dd8SToby Isaac 82320cf1dd8SToby Isaac ierr = PetscFEGetDualSpace(fem, &sp);CHKERRQ(ierr); 82420cf1dd8SToby Isaac ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr); 82520cf1dd8SToby Isaac ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 82620cf1dd8SToby Isaac ierr = DMPlexGetConeSize(dm, 0, &numFaces);CHKERRQ(ierr); 82720cf1dd8SToby Isaac ierr = DMPlexGetCone(dm, 0, &cone);CHKERRQ(ierr); 82820cf1dd8SToby Isaac ierr = PetscMalloc1(numFaces*dim, ¢roids);CHKERRQ(ierr); 82920cf1dd8SToby Isaac for (f = 0; f < numFaces; ++f) {ierr = DMPlexComputeCellGeometryFVM(dm, cone[f], NULL, ¢roids[f*dim], NULL);CHKERRQ(ierr);} 83020cf1dd8SToby Isaac ierr = PetscFEGetTabulation(fem, numFaces, centroids, &fem->F, NULL, NULL);CHKERRQ(ierr); 83120cf1dd8SToby Isaac ierr = PetscFree(centroids);CHKERRQ(ierr); 83220cf1dd8SToby Isaac } 83320cf1dd8SToby Isaac *F = fem->F; 83420cf1dd8SToby Isaac PetscFunctionReturn(0); 83520cf1dd8SToby Isaac } 83620cf1dd8SToby Isaac 83720cf1dd8SToby Isaac /*@C 83820cf1dd8SToby Isaac PetscFEGetTabulation - Tabulates the basis functions, and perhaps derivatives, at the points provided. 83920cf1dd8SToby Isaac 84020cf1dd8SToby Isaac Not collective 84120cf1dd8SToby Isaac 84220cf1dd8SToby Isaac Input Parameters: 84320cf1dd8SToby Isaac + fem - The PetscFE object 84420cf1dd8SToby Isaac . npoints - The number of tabulation points 84520cf1dd8SToby Isaac - points - The tabulation point coordinates 84620cf1dd8SToby Isaac 84720cf1dd8SToby Isaac Output Parameters: 84820cf1dd8SToby Isaac + B - The basis function values at tabulation points 84920cf1dd8SToby Isaac . D - The basis function derivatives at tabulation points 85020cf1dd8SToby Isaac - H - The basis function second derivatives at tabulation points 85120cf1dd8SToby Isaac 85220cf1dd8SToby Isaac Note: 85320cf1dd8SToby Isaac $ B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c 85420cf1dd8SToby 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 85520cf1dd8SToby 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 85620cf1dd8SToby Isaac 85720cf1dd8SToby Isaac Level: intermediate 85820cf1dd8SToby Isaac 85920cf1dd8SToby Isaac .seealso: PetscFERestoreTabulation(), PetscFEGetDefaultTabulation() 86020cf1dd8SToby Isaac @*/ 86120cf1dd8SToby Isaac PetscErrorCode PetscFEGetTabulation(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscReal **B, PetscReal **D, PetscReal **H) 86220cf1dd8SToby Isaac { 86320cf1dd8SToby Isaac DM dm; 86420cf1dd8SToby Isaac PetscInt pdim; /* Dimension of FE space P */ 86520cf1dd8SToby Isaac PetscInt dim; /* Spatial dimension */ 86620cf1dd8SToby Isaac PetscInt comp; /* Field components */ 86720cf1dd8SToby Isaac PetscErrorCode ierr; 86820cf1dd8SToby Isaac 86920cf1dd8SToby Isaac PetscFunctionBegin; 87020cf1dd8SToby Isaac if (!npoints) { 87120cf1dd8SToby Isaac if (B) *B = NULL; 87220cf1dd8SToby Isaac if (D) *D = NULL; 87320cf1dd8SToby Isaac if (H) *H = NULL; 87420cf1dd8SToby Isaac PetscFunctionReturn(0); 87520cf1dd8SToby Isaac } 87620cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 87720cf1dd8SToby Isaac PetscValidPointer(points, 3); 87820cf1dd8SToby Isaac if (B) PetscValidPointer(B, 4); 87920cf1dd8SToby Isaac if (D) PetscValidPointer(D, 5); 88020cf1dd8SToby Isaac if (H) PetscValidPointer(H, 6); 88120cf1dd8SToby Isaac ierr = PetscDualSpaceGetDM(fem->dualSpace, &dm);CHKERRQ(ierr); 88220cf1dd8SToby Isaac ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 88320cf1dd8SToby Isaac ierr = PetscDualSpaceGetDimension(fem->dualSpace, &pdim);CHKERRQ(ierr); 88420cf1dd8SToby Isaac ierr = PetscFEGetNumComponents(fem, &comp);CHKERRQ(ierr); 88520cf1dd8SToby Isaac if (B) {ierr = DMGetWorkArray(dm, npoints*pdim*comp, MPIU_REAL, B);CHKERRQ(ierr);} 88620cf1dd8SToby Isaac if (!dim) { 88720cf1dd8SToby Isaac if (D) *D = NULL; 88820cf1dd8SToby Isaac if (H) *H = NULL; 88920cf1dd8SToby Isaac } else { 89020cf1dd8SToby Isaac if (D) {ierr = DMGetWorkArray(dm, npoints*pdim*comp*dim, MPIU_REAL, D);CHKERRQ(ierr);} 89120cf1dd8SToby Isaac if (H) {ierr = DMGetWorkArray(dm, npoints*pdim*comp*dim*dim, MPIU_REAL, H);CHKERRQ(ierr);} 89220cf1dd8SToby Isaac } 89320cf1dd8SToby Isaac ierr = (*fem->ops->gettabulation)(fem, npoints, points, B ? *B : NULL, D ? *D : NULL, H ? *H : NULL);CHKERRQ(ierr); 89420cf1dd8SToby Isaac PetscFunctionReturn(0); 89520cf1dd8SToby Isaac } 89620cf1dd8SToby Isaac 89720cf1dd8SToby Isaac PetscErrorCode PetscFERestoreTabulation(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscReal **B, PetscReal **D, PetscReal **H) 89820cf1dd8SToby Isaac { 89920cf1dd8SToby Isaac DM dm; 90020cf1dd8SToby Isaac PetscErrorCode ierr; 90120cf1dd8SToby Isaac 90220cf1dd8SToby Isaac PetscFunctionBegin; 90320cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 90420cf1dd8SToby Isaac ierr = PetscDualSpaceGetDM(fem->dualSpace, &dm);CHKERRQ(ierr); 90520cf1dd8SToby Isaac if (B && *B) {ierr = DMRestoreWorkArray(dm, 0, MPIU_REAL, B);CHKERRQ(ierr);} 90620cf1dd8SToby Isaac if (D && *D) {ierr = DMRestoreWorkArray(dm, 0, MPIU_REAL, D);CHKERRQ(ierr);} 90720cf1dd8SToby Isaac if (H && *H) {ierr = DMRestoreWorkArray(dm, 0, MPIU_REAL, H);CHKERRQ(ierr);} 90820cf1dd8SToby Isaac PetscFunctionReturn(0); 90920cf1dd8SToby Isaac } 91020cf1dd8SToby Isaac 91120cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreatePointTrace(PetscFE fe, PetscInt refPoint, PetscFE *trFE) 91220cf1dd8SToby Isaac { 91320cf1dd8SToby Isaac PetscSpace bsp, bsubsp; 91420cf1dd8SToby Isaac PetscDualSpace dsp, dsubsp; 91520cf1dd8SToby Isaac PetscInt dim, depth, numComp, i, j, coneSize, order; 91620cf1dd8SToby Isaac PetscFEType type; 91720cf1dd8SToby Isaac DM dm; 91820cf1dd8SToby Isaac DMLabel label; 91920cf1dd8SToby Isaac PetscReal *xi, *v, *J, detJ; 92020cf1dd8SToby Isaac PetscQuadrature origin, fullQuad, subQuad; 92120cf1dd8SToby Isaac PetscErrorCode ierr; 92220cf1dd8SToby Isaac 92320cf1dd8SToby Isaac PetscFunctionBegin; 92420cf1dd8SToby Isaac PetscValidHeaderSpecific(fe,PETSCFE_CLASSID,1); 92520cf1dd8SToby Isaac PetscValidPointer(trFE,3); 92620cf1dd8SToby Isaac ierr = PetscFEGetBasisSpace(fe,&bsp);CHKERRQ(ierr); 92720cf1dd8SToby Isaac ierr = PetscFEGetDualSpace(fe,&dsp);CHKERRQ(ierr); 92820cf1dd8SToby Isaac ierr = PetscDualSpaceGetDM(dsp,&dm);CHKERRQ(ierr); 92920cf1dd8SToby Isaac ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 93020cf1dd8SToby Isaac ierr = DMPlexGetDepthLabel(dm,&label);CHKERRQ(ierr); 93120cf1dd8SToby Isaac ierr = DMLabelGetValue(label,refPoint,&depth);CHKERRQ(ierr); 93220cf1dd8SToby Isaac ierr = PetscCalloc1(depth,&xi);CHKERRQ(ierr); 93320cf1dd8SToby Isaac ierr = PetscMalloc1(dim,&v);CHKERRQ(ierr); 93420cf1dd8SToby Isaac ierr = PetscMalloc1(dim*dim,&J);CHKERRQ(ierr); 93520cf1dd8SToby Isaac for (i = 0; i < depth; i++) xi[i] = 0.; 93620cf1dd8SToby Isaac ierr = PetscQuadratureCreate(PETSC_COMM_SELF,&origin);CHKERRQ(ierr); 93720cf1dd8SToby Isaac ierr = PetscQuadratureSetData(origin,depth,0,1,xi,NULL);CHKERRQ(ierr); 93820cf1dd8SToby Isaac ierr = DMPlexComputeCellGeometryFEM(dm,refPoint,origin,v,J,NULL,&detJ);CHKERRQ(ierr); 93920cf1dd8SToby Isaac /* CellGeometryFEM computes the expanded Jacobian, we want the true jacobian */ 94020cf1dd8SToby Isaac for (i = 1; i < dim; i++) { 94120cf1dd8SToby Isaac for (j = 0; j < depth; j++) { 94220cf1dd8SToby Isaac J[i * depth + j] = J[i * dim + j]; 94320cf1dd8SToby Isaac } 94420cf1dd8SToby Isaac } 94520cf1dd8SToby Isaac ierr = PetscQuadratureDestroy(&origin);CHKERRQ(ierr); 94620cf1dd8SToby Isaac ierr = PetscDualSpaceGetPointSubspace(dsp,refPoint,&dsubsp);CHKERRQ(ierr); 94720cf1dd8SToby Isaac ierr = PetscSpaceCreateSubspace(bsp,dsubsp,v,J,NULL,NULL,PETSC_OWN_POINTER,&bsubsp);CHKERRQ(ierr); 94820cf1dd8SToby Isaac ierr = PetscSpaceSetUp(bsubsp);CHKERRQ(ierr); 94920cf1dd8SToby Isaac ierr = PetscFECreate(PetscObjectComm((PetscObject)fe),trFE);CHKERRQ(ierr); 95020cf1dd8SToby Isaac ierr = PetscFEGetType(fe,&type);CHKERRQ(ierr); 95120cf1dd8SToby Isaac ierr = PetscFESetType(*trFE,type);CHKERRQ(ierr); 95220cf1dd8SToby Isaac ierr = PetscFEGetNumComponents(fe,&numComp);CHKERRQ(ierr); 95320cf1dd8SToby Isaac ierr = PetscFESetNumComponents(*trFE,numComp);CHKERRQ(ierr); 95420cf1dd8SToby Isaac ierr = PetscFESetBasisSpace(*trFE,bsubsp);CHKERRQ(ierr); 95520cf1dd8SToby Isaac ierr = PetscFESetDualSpace(*trFE,dsubsp);CHKERRQ(ierr); 95620cf1dd8SToby Isaac ierr = PetscFEGetQuadrature(fe,&fullQuad);CHKERRQ(ierr); 95720cf1dd8SToby Isaac ierr = PetscQuadratureGetOrder(fullQuad,&order);CHKERRQ(ierr); 95820cf1dd8SToby Isaac ierr = DMPlexGetConeSize(dm,refPoint,&coneSize);CHKERRQ(ierr); 95920cf1dd8SToby Isaac if (coneSize == 2 * depth) { 96020cf1dd8SToby Isaac ierr = PetscDTGaussTensorQuadrature(depth,1,(order + 1)/2,-1.,1.,&subQuad);CHKERRQ(ierr); 96120cf1dd8SToby Isaac } else { 96220cf1dd8SToby Isaac ierr = PetscDTGaussJacobiQuadrature(depth,1,(order + 1)/2,-1.,1.,&subQuad);CHKERRQ(ierr); 96320cf1dd8SToby Isaac } 96420cf1dd8SToby Isaac ierr = PetscFESetQuadrature(*trFE,subQuad);CHKERRQ(ierr); 96520cf1dd8SToby Isaac ierr = PetscFESetUp(*trFE);CHKERRQ(ierr); 96620cf1dd8SToby Isaac ierr = PetscQuadratureDestroy(&subQuad);CHKERRQ(ierr); 96720cf1dd8SToby Isaac ierr = PetscSpaceDestroy(&bsubsp);CHKERRQ(ierr); 96820cf1dd8SToby Isaac PetscFunctionReturn(0); 96920cf1dd8SToby Isaac } 97020cf1dd8SToby Isaac 97120cf1dd8SToby Isaac PetscErrorCode PetscFECreateHeightTrace(PetscFE fe, PetscInt height, PetscFE *trFE) 97220cf1dd8SToby Isaac { 97320cf1dd8SToby Isaac PetscInt hStart, hEnd; 97420cf1dd8SToby Isaac PetscDualSpace dsp; 97520cf1dd8SToby Isaac DM dm; 97620cf1dd8SToby Isaac PetscErrorCode ierr; 97720cf1dd8SToby Isaac 97820cf1dd8SToby Isaac PetscFunctionBegin; 97920cf1dd8SToby Isaac PetscValidHeaderSpecific(fe,PETSCFE_CLASSID,1); 98020cf1dd8SToby Isaac PetscValidPointer(trFE,3); 98120cf1dd8SToby Isaac *trFE = NULL; 98220cf1dd8SToby Isaac ierr = PetscFEGetDualSpace(fe,&dsp);CHKERRQ(ierr); 98320cf1dd8SToby Isaac ierr = PetscDualSpaceGetDM(dsp,&dm);CHKERRQ(ierr); 98420cf1dd8SToby Isaac ierr = DMPlexGetHeightStratum(dm,height,&hStart,&hEnd);CHKERRQ(ierr); 98520cf1dd8SToby Isaac if (hEnd <= hStart) PetscFunctionReturn(0); 98620cf1dd8SToby Isaac ierr = PetscFECreatePointTrace(fe,hStart,trFE);CHKERRQ(ierr); 98720cf1dd8SToby Isaac PetscFunctionReturn(0); 98820cf1dd8SToby Isaac } 98920cf1dd8SToby Isaac 99020cf1dd8SToby Isaac 99120cf1dd8SToby Isaac /*@ 99220cf1dd8SToby Isaac PetscFEGetDimension - Get the dimension of the finite element space on a cell 99320cf1dd8SToby Isaac 99420cf1dd8SToby Isaac Not collective 99520cf1dd8SToby Isaac 99620cf1dd8SToby Isaac Input Parameter: 99720cf1dd8SToby Isaac . fe - The PetscFE 99820cf1dd8SToby Isaac 99920cf1dd8SToby Isaac Output Parameter: 100020cf1dd8SToby Isaac . dim - The dimension 100120cf1dd8SToby Isaac 100220cf1dd8SToby Isaac Level: intermediate 100320cf1dd8SToby Isaac 100420cf1dd8SToby Isaac .seealso: PetscFECreate(), PetscSpaceGetDimension(), PetscDualSpaceGetDimension() 100520cf1dd8SToby Isaac @*/ 100620cf1dd8SToby Isaac PetscErrorCode PetscFEGetDimension(PetscFE fem, PetscInt *dim) 100720cf1dd8SToby Isaac { 100820cf1dd8SToby Isaac PetscErrorCode ierr; 100920cf1dd8SToby Isaac 101020cf1dd8SToby Isaac PetscFunctionBegin; 101120cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 101220cf1dd8SToby Isaac PetscValidPointer(dim, 2); 101320cf1dd8SToby Isaac if (fem->ops->getdimension) {ierr = (*fem->ops->getdimension)(fem, dim);CHKERRQ(ierr);} 101420cf1dd8SToby Isaac PetscFunctionReturn(0); 101520cf1dd8SToby Isaac } 101620cf1dd8SToby Isaac 101720cf1dd8SToby Isaac /* 101820cf1dd8SToby Isaac Purpose: Compute element vector for chunk of elements 101920cf1dd8SToby Isaac 102020cf1dd8SToby Isaac Input: 102120cf1dd8SToby Isaac Sizes: 102220cf1dd8SToby Isaac Ne: number of elements 102320cf1dd8SToby Isaac Nf: number of fields 102420cf1dd8SToby Isaac PetscFE 102520cf1dd8SToby Isaac dim: spatial dimension 102620cf1dd8SToby Isaac Nb: number of basis functions 102720cf1dd8SToby Isaac Nc: number of field components 102820cf1dd8SToby Isaac PetscQuadrature 102920cf1dd8SToby Isaac Nq: number of quadrature points 103020cf1dd8SToby Isaac 103120cf1dd8SToby Isaac Geometry: 103220cf1dd8SToby Isaac PetscFEGeom[Ne] possibly *Nq 103320cf1dd8SToby Isaac PetscReal v0s[dim] 103420cf1dd8SToby Isaac PetscReal n[dim] 103520cf1dd8SToby Isaac PetscReal jacobians[dim*dim] 103620cf1dd8SToby Isaac PetscReal jacobianInverses[dim*dim] 103720cf1dd8SToby Isaac PetscReal jacobianDeterminants 103820cf1dd8SToby Isaac FEM: 103920cf1dd8SToby Isaac PetscFE 104020cf1dd8SToby Isaac PetscQuadrature 104120cf1dd8SToby Isaac PetscReal quadPoints[Nq*dim] 104220cf1dd8SToby Isaac PetscReal quadWeights[Nq] 104320cf1dd8SToby Isaac PetscReal basis[Nq*Nb*Nc] 104420cf1dd8SToby Isaac PetscReal basisDer[Nq*Nb*Nc*dim] 104520cf1dd8SToby Isaac PetscScalar coefficients[Ne*Nb*Nc] 104620cf1dd8SToby Isaac PetscScalar elemVec[Ne*Nb*Nc] 104720cf1dd8SToby Isaac 104820cf1dd8SToby Isaac Problem: 104920cf1dd8SToby Isaac PetscInt f: the active field 105020cf1dd8SToby Isaac f0, f1 105120cf1dd8SToby Isaac 105220cf1dd8SToby Isaac Work Space: 105320cf1dd8SToby Isaac PetscFE 105420cf1dd8SToby Isaac PetscScalar f0[Nq*dim]; 105520cf1dd8SToby Isaac PetscScalar f1[Nq*dim*dim]; 105620cf1dd8SToby Isaac PetscScalar u[Nc]; 105720cf1dd8SToby Isaac PetscScalar gradU[Nc*dim]; 105820cf1dd8SToby Isaac PetscReal x[dim]; 105920cf1dd8SToby Isaac PetscScalar realSpaceDer[dim]; 106020cf1dd8SToby Isaac 106120cf1dd8SToby Isaac Purpose: Compute element vector for N_cb batches of elements 106220cf1dd8SToby Isaac 106320cf1dd8SToby Isaac Input: 106420cf1dd8SToby Isaac Sizes: 106520cf1dd8SToby Isaac N_cb: Number of serial cell batches 106620cf1dd8SToby Isaac 106720cf1dd8SToby Isaac Geometry: 106820cf1dd8SToby Isaac PetscReal v0s[Ne*dim] 106920cf1dd8SToby Isaac PetscReal jacobians[Ne*dim*dim] possibly *Nq 107020cf1dd8SToby Isaac PetscReal jacobianInverses[Ne*dim*dim] possibly *Nq 107120cf1dd8SToby Isaac PetscReal jacobianDeterminants[Ne] possibly *Nq 107220cf1dd8SToby Isaac FEM: 107320cf1dd8SToby Isaac static PetscReal quadPoints[Nq*dim] 107420cf1dd8SToby Isaac static PetscReal quadWeights[Nq] 107520cf1dd8SToby Isaac static PetscReal basis[Nq*Nb*Nc] 107620cf1dd8SToby Isaac static PetscReal basisDer[Nq*Nb*Nc*dim] 107720cf1dd8SToby Isaac PetscScalar coefficients[Ne*Nb*Nc] 107820cf1dd8SToby Isaac PetscScalar elemVec[Ne*Nb*Nc] 107920cf1dd8SToby Isaac 108020cf1dd8SToby Isaac ex62.c: 108120cf1dd8SToby Isaac PetscErrorCode PetscFEIntegrateResidualBatch(PetscInt Ne, PetscInt numFields, PetscInt field, PetscQuadrature quad[], const PetscScalar coefficients[], 108220cf1dd8SToby Isaac const PetscReal v0s[], const PetscReal jacobians[], const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], 108320cf1dd8SToby Isaac void (*f0_func)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar f0[]), 108420cf1dd8SToby Isaac void (*f1_func)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar f1[]), PetscScalar elemVec[]) 108520cf1dd8SToby Isaac 108620cf1dd8SToby Isaac ex52.c: 108720cf1dd8SToby 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) 108820cf1dd8SToby 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) 108920cf1dd8SToby Isaac 109020cf1dd8SToby Isaac ex52_integrateElement.cu 109120cf1dd8SToby Isaac __global__ void integrateElementQuadrature(int N_cb, realType *coefficients, realType *jacobianInverses, realType *jacobianDeterminants, realType *elemVec) 109220cf1dd8SToby Isaac 109320cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode IntegrateElementBatchGPU(PetscInt spatial_dim, PetscInt Ne, PetscInt Ncb, PetscInt Nbc, PetscInt Nbl, const PetscScalar coefficients[], 109420cf1dd8SToby Isaac const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], PetscScalar elemVec[], 109520cf1dd8SToby Isaac PetscLogEvent event, PetscInt debug, PetscInt pde_op) 109620cf1dd8SToby Isaac 109720cf1dd8SToby Isaac ex52_integrateElementOpenCL.c: 109820cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode IntegrateElementBatchGPU(PetscInt spatial_dim, PetscInt Ne, PetscInt Ncb, PetscInt Nbc, PetscInt N_bl, const PetscScalar coefficients[], 109920cf1dd8SToby Isaac const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], PetscScalar elemVec[], 110020cf1dd8SToby Isaac PetscLogEvent event, PetscInt debug, PetscInt pde_op) 110120cf1dd8SToby Isaac 110220cf1dd8SToby Isaac __kernel void integrateElementQuadrature(int N_cb, __global float *coefficients, __global float *jacobianInverses, __global float *jacobianDeterminants, __global float *elemVec) 110320cf1dd8SToby Isaac */ 110420cf1dd8SToby Isaac 110520cf1dd8SToby Isaac /*@C 110620cf1dd8SToby Isaac PetscFEIntegrate - Produce the integral for the given field for a chunk of elements by quadrature integration 110720cf1dd8SToby Isaac 110820cf1dd8SToby Isaac Not collective 110920cf1dd8SToby Isaac 111020cf1dd8SToby Isaac Input Parameters: 111120cf1dd8SToby Isaac + fem - The PetscFE object for the field being integrated 111220cf1dd8SToby Isaac . prob - The PetscDS specifying the discretizations and continuum functions 111320cf1dd8SToby Isaac . field - The field being integrated 111420cf1dd8SToby Isaac . Ne - The number of elements in the chunk 111520cf1dd8SToby Isaac . cgeom - The cell geometry for each cell in the chunk 111620cf1dd8SToby Isaac . coefficients - The array of FEM basis coefficients for the elements 111720cf1dd8SToby Isaac . probAux - The PetscDS specifying the auxiliary discretizations 111820cf1dd8SToby Isaac - coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 111920cf1dd8SToby Isaac 112020cf1dd8SToby Isaac Output Parameter 112120cf1dd8SToby Isaac . integral - the integral for this field 112220cf1dd8SToby Isaac 112320cf1dd8SToby Isaac Level: developer 112420cf1dd8SToby Isaac 112520cf1dd8SToby Isaac .seealso: PetscFEIntegrateResidual() 112620cf1dd8SToby Isaac @*/ 112720cf1dd8SToby Isaac PetscErrorCode PetscFEIntegrate(PetscFE fem, PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom, 112820cf1dd8SToby Isaac const PetscScalar coefficients[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscScalar integral[]) 112920cf1dd8SToby Isaac { 113020cf1dd8SToby Isaac PetscErrorCode ierr; 113120cf1dd8SToby Isaac 113220cf1dd8SToby Isaac PetscFunctionBegin; 113320cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 113420cf1dd8SToby Isaac PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 2); 113520cf1dd8SToby Isaac if (fem->ops->integrate) {ierr = (*fem->ops->integrate)(fem, prob, field, Ne, cgeom, coefficients, probAux, coefficientsAux, integral);CHKERRQ(ierr);} 113620cf1dd8SToby Isaac PetscFunctionReturn(0); 113720cf1dd8SToby Isaac } 113820cf1dd8SToby Isaac 113920cf1dd8SToby Isaac /*@C 114020cf1dd8SToby Isaac PetscFEIntegrateResidual - Produce the element residual vector for a chunk of elements by quadrature integration 114120cf1dd8SToby Isaac 114220cf1dd8SToby Isaac Not collective 114320cf1dd8SToby Isaac 114420cf1dd8SToby Isaac Input Parameters: 114520cf1dd8SToby Isaac + fem - The PetscFE object for the field being integrated 114620cf1dd8SToby Isaac . prob - The PetscDS specifying the discretizations and continuum functions 114720cf1dd8SToby Isaac . field - The field being integrated 114820cf1dd8SToby Isaac . Ne - The number of elements in the chunk 114920cf1dd8SToby Isaac . cgeom - The cell geometry for each cell in the chunk 115020cf1dd8SToby Isaac . coefficients - The array of FEM basis coefficients for the elements 115120cf1dd8SToby Isaac . coefficients_t - The array of FEM basis time derivative coefficients for the elements 115220cf1dd8SToby Isaac . probAux - The PetscDS specifying the auxiliary discretizations 115320cf1dd8SToby Isaac . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 115420cf1dd8SToby Isaac - t - The time 115520cf1dd8SToby Isaac 115620cf1dd8SToby Isaac Output Parameter 115720cf1dd8SToby Isaac . elemVec - the element residual vectors from each element 115820cf1dd8SToby Isaac 115920cf1dd8SToby Isaac Note: 116020cf1dd8SToby Isaac $ Loop over batch of elements (e): 116120cf1dd8SToby Isaac $ Loop over quadrature points (q): 116220cf1dd8SToby Isaac $ Make u_q and gradU_q (loops over fields,Nb,Ncomp) and x_q 116320cf1dd8SToby Isaac $ Call f_0 and f_1 116420cf1dd8SToby Isaac $ Loop over element vector entries (f,fc --> i): 116520cf1dd8SToby Isaac $ elemVec[i] += \psi^{fc}_f(q) f0_{fc}(u, \nabla u) + \nabla\psi^{fc}_f(q) \cdot f1_{fc,df}(u, \nabla u) 116620cf1dd8SToby Isaac 116720cf1dd8SToby Isaac Level: developer 116820cf1dd8SToby Isaac 116920cf1dd8SToby Isaac .seealso: PetscFEIntegrateResidual() 117020cf1dd8SToby Isaac @*/ 117120cf1dd8SToby Isaac PetscErrorCode PetscFEIntegrateResidual(PetscFE fem, PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom, 117220cf1dd8SToby Isaac const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[]) 117320cf1dd8SToby Isaac { 117420cf1dd8SToby Isaac PetscErrorCode ierr; 117520cf1dd8SToby Isaac 117620cf1dd8SToby Isaac PetscFunctionBegin; 117720cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 117820cf1dd8SToby Isaac PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 2); 117920cf1dd8SToby Isaac if (fem->ops->integrateresidual) {ierr = (*fem->ops->integrateresidual)(fem, prob, field, Ne, cgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec);CHKERRQ(ierr);} 118020cf1dd8SToby Isaac PetscFunctionReturn(0); 118120cf1dd8SToby Isaac } 118220cf1dd8SToby Isaac 118320cf1dd8SToby Isaac /*@C 118420cf1dd8SToby Isaac PetscFEIntegrateBdResidual - Produce the element residual vector for a chunk of elements by quadrature integration over a boundary 118520cf1dd8SToby Isaac 118620cf1dd8SToby Isaac Not collective 118720cf1dd8SToby Isaac 118820cf1dd8SToby Isaac Input Parameters: 118920cf1dd8SToby Isaac + fem - The PetscFE object for the field being integrated 119020cf1dd8SToby Isaac . prob - The PetscDS specifying the discretizations and continuum functions 119120cf1dd8SToby Isaac . field - The field being integrated 119220cf1dd8SToby Isaac . Ne - The number of elements in the chunk 119320cf1dd8SToby Isaac . fgeom - The face geometry for each cell in the chunk 119420cf1dd8SToby Isaac . coefficients - The array of FEM basis coefficients for the elements 119520cf1dd8SToby Isaac . coefficients_t - The array of FEM basis time derivative coefficients for the elements 119620cf1dd8SToby Isaac . probAux - The PetscDS specifying the auxiliary discretizations 119720cf1dd8SToby Isaac . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 119820cf1dd8SToby Isaac - t - The time 119920cf1dd8SToby Isaac 120020cf1dd8SToby Isaac Output Parameter 120120cf1dd8SToby Isaac . elemVec - the element residual vectors from each element 120220cf1dd8SToby Isaac 120320cf1dd8SToby Isaac Level: developer 120420cf1dd8SToby Isaac 120520cf1dd8SToby Isaac .seealso: PetscFEIntegrateResidual() 120620cf1dd8SToby Isaac @*/ 120720cf1dd8SToby Isaac PetscErrorCode PetscFEIntegrateBdResidual(PetscFE fem, PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *fgeom, 120820cf1dd8SToby Isaac const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[]) 120920cf1dd8SToby Isaac { 121020cf1dd8SToby Isaac PetscErrorCode ierr; 121120cf1dd8SToby Isaac 121220cf1dd8SToby Isaac PetscFunctionBegin; 121320cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 121420cf1dd8SToby Isaac if (fem->ops->integratebdresidual) {ierr = (*fem->ops->integratebdresidual)(fem, prob, field, Ne, fgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec);CHKERRQ(ierr);} 121520cf1dd8SToby Isaac PetscFunctionReturn(0); 121620cf1dd8SToby Isaac } 121720cf1dd8SToby Isaac 121820cf1dd8SToby Isaac /*@C 121920cf1dd8SToby Isaac PetscFEIntegrateJacobian - Produce the element Jacobian for a chunk of elements by quadrature integration 122020cf1dd8SToby Isaac 122120cf1dd8SToby Isaac Not collective 122220cf1dd8SToby Isaac 122320cf1dd8SToby Isaac Input Parameters: 122420cf1dd8SToby Isaac + fem - The PetscFE object for the field being integrated 122520cf1dd8SToby Isaac . prob - The PetscDS specifying the discretizations and continuum functions 122620cf1dd8SToby Isaac . jtype - The type of matrix pointwise functions that should be used 122720cf1dd8SToby Isaac . fieldI - The test field being integrated 122820cf1dd8SToby Isaac . fieldJ - The basis field being integrated 122920cf1dd8SToby Isaac . Ne - The number of elements in the chunk 123020cf1dd8SToby Isaac . cgeom - The cell geometry for each cell in the chunk 123120cf1dd8SToby Isaac . coefficients - The array of FEM basis coefficients for the elements for the Jacobian evaluation point 123220cf1dd8SToby Isaac . coefficients_t - The array of FEM basis time derivative coefficients for the elements 123320cf1dd8SToby Isaac . probAux - The PetscDS specifying the auxiliary discretizations 123420cf1dd8SToby Isaac . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements 123520cf1dd8SToby Isaac . t - The time 123620cf1dd8SToby Isaac - u_tShift - A multiplier for the dF/du_t term (as opposed to the dF/du term) 123720cf1dd8SToby Isaac 123820cf1dd8SToby Isaac Output Parameter 123920cf1dd8SToby Isaac . elemMat - the element matrices for the Jacobian from each element 124020cf1dd8SToby Isaac 124120cf1dd8SToby Isaac Note: 124220cf1dd8SToby Isaac $ Loop over batch of elements (e): 124320cf1dd8SToby Isaac $ Loop over element matrix entries (f,fc,g,gc --> i,j): 124420cf1dd8SToby Isaac $ Loop over quadrature points (q): 124520cf1dd8SToby Isaac $ Make u_q and gradU_q (loops over fields,Nb,Ncomp) 124620cf1dd8SToby Isaac $ elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q) 124720cf1dd8SToby Isaac $ + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q) 124820cf1dd8SToby Isaac $ + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q) 124920cf1dd8SToby Isaac $ + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q) 125020cf1dd8SToby Isaac Level: developer 125120cf1dd8SToby Isaac 125220cf1dd8SToby Isaac .seealso: PetscFEIntegrateResidual() 125320cf1dd8SToby Isaac @*/ 125420cf1dd8SToby Isaac PetscErrorCode PetscFEIntegrateJacobian(PetscFE fem, PetscDS prob, PetscFEJacobianType jtype, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *cgeom, 125520cf1dd8SToby Isaac const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) 125620cf1dd8SToby Isaac { 125720cf1dd8SToby Isaac PetscErrorCode ierr; 125820cf1dd8SToby Isaac 125920cf1dd8SToby Isaac PetscFunctionBegin; 126020cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 126120cf1dd8SToby 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);} 126220cf1dd8SToby Isaac PetscFunctionReturn(0); 126320cf1dd8SToby Isaac } 126420cf1dd8SToby Isaac 126520cf1dd8SToby Isaac /*@C 126620cf1dd8SToby Isaac PetscFEIntegrateBdJacobian - Produce the boundary element Jacobian for a chunk of elements by quadrature integration 126720cf1dd8SToby Isaac 126820cf1dd8SToby Isaac Not collective 126920cf1dd8SToby Isaac 127020cf1dd8SToby Isaac Input Parameters: 127120cf1dd8SToby Isaac + fem = The PetscFE object for the field being integrated 127220cf1dd8SToby Isaac . prob - The PetscDS specifying the discretizations and continuum functions 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 . fgeom - The face 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: PetscFEIntegrateJacobian(), PetscFEIntegrateResidual() 129920cf1dd8SToby Isaac @*/ 130020cf1dd8SToby Isaac PetscErrorCode PetscFEIntegrateBdJacobian(PetscFE fem, PetscDS prob, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *fgeom, 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->integratebdjacobian) {ierr = (*fem->ops->integratebdjacobian)(fem, prob, fieldI, fieldJ, Ne, fgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, u_tshift, elemMat);CHKERRQ(ierr);} 130820cf1dd8SToby Isaac PetscFunctionReturn(0); 130920cf1dd8SToby Isaac } 131020cf1dd8SToby Isaac 131120cf1dd8SToby Isaac PetscErrorCode PetscFEGetHeightSubspace(PetscFE fe, PetscInt height, PetscFE *subfe) 131220cf1dd8SToby Isaac { 131320cf1dd8SToby Isaac PetscSpace P, subP; 131420cf1dd8SToby Isaac PetscDualSpace Q, subQ; 131520cf1dd8SToby Isaac PetscQuadrature subq; 131620cf1dd8SToby Isaac PetscFEType fetype; 131720cf1dd8SToby Isaac PetscInt dim, Nc; 131820cf1dd8SToby Isaac PetscErrorCode ierr; 131920cf1dd8SToby Isaac 132020cf1dd8SToby Isaac PetscFunctionBegin; 132120cf1dd8SToby Isaac PetscValidHeaderSpecific(fe, PETSCFE_CLASSID, 1); 132220cf1dd8SToby Isaac PetscValidPointer(subfe, 3); 132320cf1dd8SToby Isaac if (height == 0) { 132420cf1dd8SToby Isaac *subfe = fe; 132520cf1dd8SToby Isaac PetscFunctionReturn(0); 132620cf1dd8SToby Isaac } 132720cf1dd8SToby Isaac ierr = PetscFEGetBasisSpace(fe, &P);CHKERRQ(ierr); 132820cf1dd8SToby Isaac ierr = PetscFEGetDualSpace(fe, &Q);CHKERRQ(ierr); 132920cf1dd8SToby Isaac ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 133020cf1dd8SToby Isaac ierr = PetscFEGetFaceQuadrature(fe, &subq);CHKERRQ(ierr); 133120cf1dd8SToby Isaac ierr = PetscDualSpaceGetDimension(Q, &dim);CHKERRQ(ierr); 133220cf1dd8SToby 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);} 133320cf1dd8SToby Isaac if (!fe->subspaces) {ierr = PetscCalloc1(dim, &fe->subspaces);CHKERRQ(ierr);} 133420cf1dd8SToby Isaac if (height <= dim) { 133520cf1dd8SToby Isaac if (!fe->subspaces[height-1]) { 133620cf1dd8SToby Isaac PetscFE sub; 133720cf1dd8SToby Isaac 133820cf1dd8SToby Isaac ierr = PetscSpaceGetHeightSubspace(P, height, &subP);CHKERRQ(ierr); 133920cf1dd8SToby Isaac ierr = PetscDualSpaceGetHeightSubspace(Q, height, &subQ);CHKERRQ(ierr); 134020cf1dd8SToby Isaac ierr = PetscFECreate(PetscObjectComm((PetscObject) fe), &sub);CHKERRQ(ierr); 134120cf1dd8SToby Isaac ierr = PetscFEGetType(fe, &fetype);CHKERRQ(ierr); 134220cf1dd8SToby Isaac ierr = PetscFESetType(sub, fetype);CHKERRQ(ierr); 134320cf1dd8SToby Isaac ierr = PetscFESetBasisSpace(sub, subP);CHKERRQ(ierr); 134420cf1dd8SToby Isaac ierr = PetscFESetDualSpace(sub, subQ);CHKERRQ(ierr); 134520cf1dd8SToby Isaac ierr = PetscFESetNumComponents(sub, Nc);CHKERRQ(ierr); 134620cf1dd8SToby Isaac ierr = PetscFESetUp(sub);CHKERRQ(ierr); 134720cf1dd8SToby Isaac ierr = PetscFESetQuadrature(sub, subq);CHKERRQ(ierr); 134820cf1dd8SToby Isaac fe->subspaces[height-1] = sub; 134920cf1dd8SToby Isaac } 135020cf1dd8SToby Isaac *subfe = fe->subspaces[height-1]; 135120cf1dd8SToby Isaac } else { 135220cf1dd8SToby Isaac *subfe = NULL; 135320cf1dd8SToby Isaac } 135420cf1dd8SToby Isaac PetscFunctionReturn(0); 135520cf1dd8SToby Isaac } 135620cf1dd8SToby Isaac 135720cf1dd8SToby Isaac /*@ 135820cf1dd8SToby Isaac PetscFERefine - Create a "refined" PetscFE object that refines the reference cell into smaller copies. This is typically used 135920cf1dd8SToby Isaac to precondition a higher order method with a lower order method on a refined mesh having the same number of dofs (but more 136020cf1dd8SToby Isaac sparsity). It is also used to create an interpolation between regularly refined meshes. 136120cf1dd8SToby Isaac 136220cf1dd8SToby Isaac Collective on PetscFE 136320cf1dd8SToby Isaac 136420cf1dd8SToby Isaac Input Parameter: 136520cf1dd8SToby Isaac . fe - The initial PetscFE 136620cf1dd8SToby Isaac 136720cf1dd8SToby Isaac Output Parameter: 136820cf1dd8SToby Isaac . feRef - The refined PetscFE 136920cf1dd8SToby Isaac 137020cf1dd8SToby Isaac Level: developer 137120cf1dd8SToby Isaac 137220cf1dd8SToby Isaac .seealso: PetscFEType, PetscFECreate(), PetscFESetType() 137320cf1dd8SToby Isaac @*/ 137420cf1dd8SToby Isaac PetscErrorCode PetscFERefine(PetscFE fe, PetscFE *feRef) 137520cf1dd8SToby Isaac { 137620cf1dd8SToby Isaac PetscSpace P, Pref; 137720cf1dd8SToby Isaac PetscDualSpace Q, Qref; 137820cf1dd8SToby Isaac DM K, Kref; 137920cf1dd8SToby Isaac PetscQuadrature q, qref; 138020cf1dd8SToby Isaac const PetscReal *v0, *jac; 138120cf1dd8SToby Isaac PetscInt numComp, numSubelements; 138220cf1dd8SToby Isaac PetscErrorCode ierr; 138320cf1dd8SToby Isaac 138420cf1dd8SToby Isaac PetscFunctionBegin; 138520cf1dd8SToby Isaac ierr = PetscFEGetBasisSpace(fe, &P);CHKERRQ(ierr); 138620cf1dd8SToby Isaac ierr = PetscFEGetDualSpace(fe, &Q);CHKERRQ(ierr); 138720cf1dd8SToby Isaac ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr); 138820cf1dd8SToby Isaac ierr = PetscDualSpaceGetDM(Q, &K);CHKERRQ(ierr); 138920cf1dd8SToby Isaac /* Create space */ 139020cf1dd8SToby Isaac ierr = PetscObjectReference((PetscObject) P);CHKERRQ(ierr); 139120cf1dd8SToby Isaac Pref = P; 139220cf1dd8SToby Isaac /* Create dual space */ 139320cf1dd8SToby Isaac ierr = PetscDualSpaceDuplicate(Q, &Qref);CHKERRQ(ierr); 139420cf1dd8SToby Isaac ierr = DMRefine(K, PetscObjectComm((PetscObject) fe), &Kref);CHKERRQ(ierr); 139520cf1dd8SToby Isaac ierr = PetscDualSpaceSetDM(Qref, Kref);CHKERRQ(ierr); 139620cf1dd8SToby Isaac ierr = DMDestroy(&Kref);CHKERRQ(ierr); 139720cf1dd8SToby Isaac ierr = PetscDualSpaceSetUp(Qref);CHKERRQ(ierr); 139820cf1dd8SToby Isaac /* Create element */ 139920cf1dd8SToby Isaac ierr = PetscFECreate(PetscObjectComm((PetscObject) fe), feRef);CHKERRQ(ierr); 140020cf1dd8SToby Isaac ierr = PetscFESetType(*feRef, PETSCFECOMPOSITE);CHKERRQ(ierr); 140120cf1dd8SToby Isaac ierr = PetscFESetBasisSpace(*feRef, Pref);CHKERRQ(ierr); 140220cf1dd8SToby Isaac ierr = PetscFESetDualSpace(*feRef, Qref);CHKERRQ(ierr); 140320cf1dd8SToby Isaac ierr = PetscFEGetNumComponents(fe, &numComp);CHKERRQ(ierr); 140420cf1dd8SToby Isaac ierr = PetscFESetNumComponents(*feRef, numComp);CHKERRQ(ierr); 140520cf1dd8SToby Isaac ierr = PetscFESetUp(*feRef);CHKERRQ(ierr); 140620cf1dd8SToby Isaac ierr = PetscSpaceDestroy(&Pref);CHKERRQ(ierr); 140720cf1dd8SToby Isaac ierr = PetscDualSpaceDestroy(&Qref);CHKERRQ(ierr); 140820cf1dd8SToby Isaac /* Create quadrature */ 140920cf1dd8SToby Isaac ierr = PetscFECompositeGetMapping(*feRef, &numSubelements, &v0, &jac, NULL);CHKERRQ(ierr); 141020cf1dd8SToby Isaac ierr = PetscQuadratureExpandComposite(q, numSubelements, v0, jac, &qref);CHKERRQ(ierr); 141120cf1dd8SToby Isaac ierr = PetscFESetQuadrature(*feRef, qref);CHKERRQ(ierr); 141220cf1dd8SToby Isaac ierr = PetscQuadratureDestroy(&qref);CHKERRQ(ierr); 141320cf1dd8SToby Isaac PetscFunctionReturn(0); 141420cf1dd8SToby Isaac } 141520cf1dd8SToby Isaac 141620cf1dd8SToby Isaac /*@C 141720cf1dd8SToby Isaac PetscFECreateDefault - Create a PetscFE for basic FEM computation 141820cf1dd8SToby Isaac 141920cf1dd8SToby Isaac Collective on DM 142020cf1dd8SToby Isaac 142120cf1dd8SToby Isaac Input Parameters: 14227be5e748SToby Isaac + comm - The MPI comm 142320cf1dd8SToby Isaac . dim - The spatial dimension 142420cf1dd8SToby Isaac . Nc - The number of components 142520cf1dd8SToby Isaac . isSimplex - Flag for simplex reference cell, otherwise its a tensor product 142620cf1dd8SToby Isaac . prefix - The options prefix, or NULL 142720cf1dd8SToby Isaac - qorder - The quadrature order 142820cf1dd8SToby Isaac 142920cf1dd8SToby Isaac Output Parameter: 143020cf1dd8SToby Isaac . fem - The PetscFE object 143120cf1dd8SToby Isaac 143220cf1dd8SToby Isaac Level: beginner 143320cf1dd8SToby Isaac 143420cf1dd8SToby Isaac .keywords: PetscFE, finite element 143520cf1dd8SToby Isaac .seealso: PetscFECreate(), PetscSpaceCreate(), PetscDualSpaceCreate() 143620cf1dd8SToby Isaac @*/ 14377be5e748SToby Isaac PetscErrorCode PetscFECreateDefault(MPI_Comm comm, PetscInt dim, PetscInt Nc, PetscBool isSimplex, const char prefix[], PetscInt qorder, PetscFE *fem) 143820cf1dd8SToby Isaac { 143920cf1dd8SToby Isaac PetscQuadrature q, fq; 144020cf1dd8SToby Isaac DM K; 144120cf1dd8SToby Isaac PetscSpace P; 144220cf1dd8SToby Isaac PetscDualSpace Q; 144320cf1dd8SToby Isaac PetscInt order, quadPointsPerEdge; 144420cf1dd8SToby Isaac PetscBool tensor = isSimplex ? PETSC_FALSE : PETSC_TRUE; 144520cf1dd8SToby Isaac PetscErrorCode ierr; 144620cf1dd8SToby Isaac 144720cf1dd8SToby Isaac PetscFunctionBegin; 144820cf1dd8SToby Isaac /* Create space */ 14497be5e748SToby Isaac ierr = PetscSpaceCreate(comm, &P);CHKERRQ(ierr); 145020cf1dd8SToby Isaac ierr = PetscObjectSetOptionsPrefix((PetscObject) P, prefix);CHKERRQ(ierr); 145120cf1dd8SToby Isaac ierr = PetscSpacePolynomialSetTensor(P, tensor);CHKERRQ(ierr); 145220cf1dd8SToby Isaac ierr = PetscSpaceSetNumComponents(P, Nc);CHKERRQ(ierr); 145320cf1dd8SToby Isaac ierr = PetscSpaceSetNumVariables(P, dim);CHKERRQ(ierr); 1454*028afddaSToby Isaac ierr = PetscSpaceSetFromOptions(P);CHKERRQ(ierr); 145520cf1dd8SToby Isaac ierr = PetscSpaceSetUp(P);CHKERRQ(ierr); 145620cf1dd8SToby Isaac ierr = PetscSpaceGetDegree(P, &order, NULL);CHKERRQ(ierr); 145720cf1dd8SToby Isaac ierr = PetscSpacePolynomialGetTensor(P, &tensor);CHKERRQ(ierr); 145820cf1dd8SToby Isaac /* Create dual space */ 14597be5e748SToby Isaac ierr = PetscDualSpaceCreate(comm, &Q);CHKERRQ(ierr); 146020cf1dd8SToby Isaac ierr = PetscDualSpaceSetType(Q,PETSCDUALSPACELAGRANGE);CHKERRQ(ierr); 146120cf1dd8SToby Isaac ierr = PetscObjectSetOptionsPrefix((PetscObject) Q, prefix);CHKERRQ(ierr); 146220cf1dd8SToby Isaac ierr = PetscDualSpaceCreateReferenceCell(Q, dim, isSimplex, &K);CHKERRQ(ierr); 146320cf1dd8SToby Isaac ierr = PetscDualSpaceSetDM(Q, K);CHKERRQ(ierr); 146420cf1dd8SToby Isaac ierr = DMDestroy(&K);CHKERRQ(ierr); 146520cf1dd8SToby Isaac ierr = PetscDualSpaceSetNumComponents(Q, Nc);CHKERRQ(ierr); 146620cf1dd8SToby Isaac ierr = PetscDualSpaceSetOrder(Q, order);CHKERRQ(ierr); 146720cf1dd8SToby Isaac ierr = PetscDualSpaceLagrangeSetTensor(Q, tensor);CHKERRQ(ierr); 146820cf1dd8SToby Isaac ierr = PetscDualSpaceSetFromOptions(Q);CHKERRQ(ierr); 146920cf1dd8SToby Isaac ierr = PetscDualSpaceSetUp(Q);CHKERRQ(ierr); 147020cf1dd8SToby Isaac /* Create element */ 14717be5e748SToby Isaac ierr = PetscFECreate(comm, fem);CHKERRQ(ierr); 147220cf1dd8SToby Isaac ierr = PetscObjectSetOptionsPrefix((PetscObject) *fem, prefix);CHKERRQ(ierr); 147320cf1dd8SToby Isaac ierr = PetscFESetFromOptions(*fem);CHKERRQ(ierr); 147420cf1dd8SToby Isaac ierr = PetscFESetBasisSpace(*fem, P);CHKERRQ(ierr); 147520cf1dd8SToby Isaac ierr = PetscFESetDualSpace(*fem, Q);CHKERRQ(ierr); 147620cf1dd8SToby Isaac ierr = PetscFESetNumComponents(*fem, Nc);CHKERRQ(ierr); 147720cf1dd8SToby Isaac ierr = PetscFESetUp(*fem);CHKERRQ(ierr); 147820cf1dd8SToby Isaac ierr = PetscSpaceDestroy(&P);CHKERRQ(ierr); 147920cf1dd8SToby Isaac ierr = PetscDualSpaceDestroy(&Q);CHKERRQ(ierr); 148020cf1dd8SToby Isaac /* Create quadrature (with specified order if given) */ 148120cf1dd8SToby Isaac qorder = qorder >= 0 ? qorder : order; 148220cf1dd8SToby Isaac ierr = PetscObjectOptionsBegin((PetscObject)*fem);CHKERRQ(ierr); 148320cf1dd8SToby Isaac ierr = PetscOptionsInt("-petscfe_default_quadrature_order","Quadrature order is one less than quadture points per edge","PetscFECreateDefault",qorder,&qorder,NULL);CHKERRQ(ierr); 148420cf1dd8SToby Isaac ierr = PetscOptionsEnd();CHKERRQ(ierr); 148520cf1dd8SToby Isaac quadPointsPerEdge = PetscMax(qorder + 1,1); 148620cf1dd8SToby Isaac if (isSimplex) { 148720cf1dd8SToby Isaac ierr = PetscDTGaussJacobiQuadrature(dim, 1, quadPointsPerEdge, -1.0, 1.0, &q);CHKERRQ(ierr); 148820cf1dd8SToby Isaac ierr = PetscDTGaussJacobiQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0, &fq);CHKERRQ(ierr); 148920cf1dd8SToby Isaac } 149020cf1dd8SToby Isaac else { 149120cf1dd8SToby Isaac ierr = PetscDTGaussTensorQuadrature(dim, 1, quadPointsPerEdge, -1.0, 1.0, &q);CHKERRQ(ierr); 149220cf1dd8SToby Isaac ierr = PetscDTGaussTensorQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0, &fq);CHKERRQ(ierr); 149320cf1dd8SToby Isaac } 149420cf1dd8SToby Isaac ierr = PetscFESetQuadrature(*fem, q);CHKERRQ(ierr); 149520cf1dd8SToby Isaac ierr = PetscFESetFaceQuadrature(*fem, fq);CHKERRQ(ierr); 149620cf1dd8SToby Isaac ierr = PetscQuadratureDestroy(&q);CHKERRQ(ierr); 149720cf1dd8SToby Isaac ierr = PetscQuadratureDestroy(&fq);CHKERRQ(ierr); 149820cf1dd8SToby Isaac PetscFunctionReturn(0); 149920cf1dd8SToby Isaac } 150020cf1dd8SToby Isaac 1501