xref: /petsc/src/dm/dt/fe/interface/fe.c (revision 5a856986583887c326abe5dfd149e8184a29cd80)
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 .seealso: PetscFERegisterAll(), PetscFERegisterDestroy()
8620cf1dd8SToby Isaac 
8720cf1dd8SToby Isaac @*/
8820cf1dd8SToby Isaac PetscErrorCode PetscFERegister(const char sname[], PetscErrorCode (*function)(PetscFE))
8920cf1dd8SToby Isaac {
9020cf1dd8SToby Isaac   PetscErrorCode ierr;
9120cf1dd8SToby Isaac 
9220cf1dd8SToby Isaac   PetscFunctionBegin;
9320cf1dd8SToby Isaac   ierr = PetscFunctionListAdd(&PetscFEList, sname, function);CHKERRQ(ierr);
9420cf1dd8SToby Isaac   PetscFunctionReturn(0);
9520cf1dd8SToby Isaac }
9620cf1dd8SToby Isaac 
9720cf1dd8SToby Isaac /*@C
9820cf1dd8SToby Isaac   PetscFESetType - Builds a particular PetscFE
9920cf1dd8SToby Isaac 
100d083f849SBarry Smith   Collective on fem
10120cf1dd8SToby Isaac 
10220cf1dd8SToby Isaac   Input Parameters:
10320cf1dd8SToby Isaac + fem  - The PetscFE object
10420cf1dd8SToby Isaac - name - The kind of FEM space
10520cf1dd8SToby Isaac 
10620cf1dd8SToby Isaac   Options Database Key:
10720cf1dd8SToby Isaac . -petscfe_type <type> - Sets the PetscFE type; use -help for a list of available types
10820cf1dd8SToby Isaac 
10920cf1dd8SToby Isaac   Level: intermediate
11020cf1dd8SToby Isaac 
11120cf1dd8SToby Isaac .seealso: PetscFEGetType(), PetscFECreate()
11220cf1dd8SToby Isaac @*/
11320cf1dd8SToby Isaac PetscErrorCode PetscFESetType(PetscFE fem, PetscFEType name)
11420cf1dd8SToby Isaac {
11520cf1dd8SToby Isaac   PetscErrorCode (*r)(PetscFE);
11620cf1dd8SToby Isaac   PetscBool      match;
11720cf1dd8SToby Isaac   PetscErrorCode ierr;
11820cf1dd8SToby Isaac 
11920cf1dd8SToby Isaac   PetscFunctionBegin;
12020cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
12120cf1dd8SToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) fem, name, &match);CHKERRQ(ierr);
12220cf1dd8SToby Isaac   if (match) PetscFunctionReturn(0);
12320cf1dd8SToby Isaac 
12420cf1dd8SToby Isaac   if (!PetscFERegisterAllCalled) {ierr = PetscFERegisterAll();CHKERRQ(ierr);}
12520cf1dd8SToby Isaac   ierr = PetscFunctionListFind(PetscFEList, name, &r);CHKERRQ(ierr);
12620cf1dd8SToby Isaac   if (!r) SETERRQ1(PetscObjectComm((PetscObject) fem), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscFE type: %s", name);
12720cf1dd8SToby Isaac 
12820cf1dd8SToby Isaac   if (fem->ops->destroy) {
12920cf1dd8SToby Isaac     ierr              = (*fem->ops->destroy)(fem);CHKERRQ(ierr);
13020cf1dd8SToby Isaac     fem->ops->destroy = NULL;
13120cf1dd8SToby Isaac   }
13220cf1dd8SToby Isaac   ierr = (*r)(fem);CHKERRQ(ierr);
13320cf1dd8SToby Isaac   ierr = PetscObjectChangeTypeName((PetscObject) fem, name);CHKERRQ(ierr);
13420cf1dd8SToby Isaac   PetscFunctionReturn(0);
13520cf1dd8SToby Isaac }
13620cf1dd8SToby Isaac 
13720cf1dd8SToby Isaac /*@C
13820cf1dd8SToby Isaac   PetscFEGetType - Gets the PetscFE type name (as a string) from the object.
13920cf1dd8SToby Isaac 
14020cf1dd8SToby Isaac   Not Collective
14120cf1dd8SToby Isaac 
14220cf1dd8SToby Isaac   Input Parameter:
14320cf1dd8SToby Isaac . fem  - The PetscFE
14420cf1dd8SToby Isaac 
14520cf1dd8SToby Isaac   Output Parameter:
14620cf1dd8SToby Isaac . name - The PetscFE type name
14720cf1dd8SToby Isaac 
14820cf1dd8SToby Isaac   Level: intermediate
14920cf1dd8SToby Isaac 
15020cf1dd8SToby Isaac .seealso: PetscFESetType(), PetscFECreate()
15120cf1dd8SToby Isaac @*/
15220cf1dd8SToby Isaac PetscErrorCode PetscFEGetType(PetscFE fem, PetscFEType *name)
15320cf1dd8SToby Isaac {
15420cf1dd8SToby Isaac   PetscErrorCode ierr;
15520cf1dd8SToby Isaac 
15620cf1dd8SToby Isaac   PetscFunctionBegin;
15720cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
15820cf1dd8SToby Isaac   PetscValidPointer(name, 2);
15920cf1dd8SToby Isaac   if (!PetscFERegisterAllCalled) {
16020cf1dd8SToby Isaac     ierr = PetscFERegisterAll();CHKERRQ(ierr);
16120cf1dd8SToby Isaac   }
16220cf1dd8SToby Isaac   *name = ((PetscObject) fem)->type_name;
16320cf1dd8SToby Isaac   PetscFunctionReturn(0);
16420cf1dd8SToby Isaac }
16520cf1dd8SToby Isaac 
16620cf1dd8SToby Isaac /*@C
16720cf1dd8SToby Isaac   PetscFEView - Views a PetscFE
16820cf1dd8SToby Isaac 
169d083f849SBarry Smith   Collective on fem
17020cf1dd8SToby Isaac 
17120cf1dd8SToby Isaac   Input Parameter:
17220cf1dd8SToby Isaac + fem - the PetscFE object to view
173d9bac1caSLisandro Dalcin - viewer   - the viewer
17420cf1dd8SToby Isaac 
17520cf1dd8SToby Isaac   Level: developer
17620cf1dd8SToby Isaac 
17720cf1dd8SToby Isaac .seealso PetscFEDestroy()
17820cf1dd8SToby Isaac @*/
179d9bac1caSLisandro Dalcin PetscErrorCode PetscFEView(PetscFE fem, PetscViewer viewer)
18020cf1dd8SToby Isaac {
181d9bac1caSLisandro Dalcin   PetscBool      iascii;
18220cf1dd8SToby Isaac   PetscErrorCode ierr;
18320cf1dd8SToby Isaac 
18420cf1dd8SToby Isaac   PetscFunctionBegin;
18520cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
186d9bac1caSLisandro Dalcin   if (viewer) PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
187d9bac1caSLisandro Dalcin   if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) fem), &viewer);CHKERRQ(ierr);}
188d9bac1caSLisandro Dalcin   ierr = PetscObjectPrintClassNamePrefixType((PetscObject)fem, viewer);CHKERRQ(ierr);
189d9bac1caSLisandro Dalcin   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
190d9bac1caSLisandro Dalcin   if (fem->ops->view) {ierr = (*fem->ops->view)(fem, viewer);CHKERRQ(ierr);}
19120cf1dd8SToby Isaac   PetscFunctionReturn(0);
19220cf1dd8SToby Isaac }
19320cf1dd8SToby Isaac 
19420cf1dd8SToby Isaac /*@
19520cf1dd8SToby Isaac   PetscFESetFromOptions - sets parameters in a PetscFE from the options database
19620cf1dd8SToby Isaac 
197d083f849SBarry Smith   Collective on fem
19820cf1dd8SToby Isaac 
19920cf1dd8SToby Isaac   Input Parameter:
20020cf1dd8SToby Isaac . fem - the PetscFE object to set options for
20120cf1dd8SToby Isaac 
20220cf1dd8SToby Isaac   Options Database:
20320cf1dd8SToby Isaac . -petscfe_num_blocks  the number of cell blocks to integrate concurrently
20420cf1dd8SToby Isaac . -petscfe_num_batches the number of cell batches to integrate serially
20520cf1dd8SToby Isaac 
20620cf1dd8SToby Isaac   Level: developer
20720cf1dd8SToby Isaac 
20820cf1dd8SToby Isaac .seealso PetscFEView()
20920cf1dd8SToby Isaac @*/
21020cf1dd8SToby Isaac PetscErrorCode PetscFESetFromOptions(PetscFE fem)
21120cf1dd8SToby Isaac {
21220cf1dd8SToby Isaac   const char    *defaultType;
21320cf1dd8SToby Isaac   char           name[256];
21420cf1dd8SToby Isaac   PetscBool      flg;
21520cf1dd8SToby Isaac   PetscErrorCode ierr;
21620cf1dd8SToby Isaac 
21720cf1dd8SToby Isaac   PetscFunctionBegin;
21820cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
21920cf1dd8SToby Isaac   if (!((PetscObject) fem)->type_name) {
22020cf1dd8SToby Isaac     defaultType = PETSCFEBASIC;
22120cf1dd8SToby Isaac   } else {
22220cf1dd8SToby Isaac     defaultType = ((PetscObject) fem)->type_name;
22320cf1dd8SToby Isaac   }
22420cf1dd8SToby Isaac   if (!PetscFERegisterAllCalled) {ierr = PetscFERegisterAll();CHKERRQ(ierr);}
22520cf1dd8SToby Isaac 
22620cf1dd8SToby Isaac   ierr = PetscObjectOptionsBegin((PetscObject) fem);CHKERRQ(ierr);
22720cf1dd8SToby Isaac   ierr = PetscOptionsFList("-petscfe_type", "Finite element space", "PetscFESetType", PetscFEList, defaultType, name, 256, &flg);CHKERRQ(ierr);
22820cf1dd8SToby Isaac   if (flg) {
22920cf1dd8SToby Isaac     ierr = PetscFESetType(fem, name);CHKERRQ(ierr);
23020cf1dd8SToby Isaac   } else if (!((PetscObject) fem)->type_name) {
23120cf1dd8SToby Isaac     ierr = PetscFESetType(fem, defaultType);CHKERRQ(ierr);
23220cf1dd8SToby Isaac   }
233*5a856986SBarry Smith   ierr = PetscOptionsBoundedInt("-petscfe_num_blocks", "The number of cell blocks to integrate concurrently", "PetscSpaceSetTileSizes", fem->numBlocks, &fem->numBlocks, NULL,1);CHKERRQ(ierr);
234*5a856986SBarry Smith   ierr = PetscOptionsBoundedInt("-petscfe_num_batches", "The number of cell batches to integrate serially", "PetscSpaceSetTileSizes", fem->numBatches, &fem->numBatches, NULL,1);CHKERRQ(ierr);
23520cf1dd8SToby Isaac   if (fem->ops->setfromoptions) {
23620cf1dd8SToby Isaac     ierr = (*fem->ops->setfromoptions)(PetscOptionsObject,fem);CHKERRQ(ierr);
23720cf1dd8SToby Isaac   }
23820cf1dd8SToby Isaac   /* process any options handlers added with PetscObjectAddOptionsHandler() */
23920cf1dd8SToby Isaac   ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) fem);CHKERRQ(ierr);
24020cf1dd8SToby Isaac   ierr = PetscOptionsEnd();CHKERRQ(ierr);
24120cf1dd8SToby Isaac   ierr = PetscFEViewFromOptions(fem, NULL, "-petscfe_view");CHKERRQ(ierr);
24220cf1dd8SToby Isaac   PetscFunctionReturn(0);
24320cf1dd8SToby Isaac }
24420cf1dd8SToby Isaac 
24520cf1dd8SToby Isaac /*@C
24620cf1dd8SToby Isaac   PetscFESetUp - Construct data structures for the PetscFE
24720cf1dd8SToby Isaac 
248d083f849SBarry Smith   Collective on fem
24920cf1dd8SToby Isaac 
25020cf1dd8SToby Isaac   Input Parameter:
25120cf1dd8SToby Isaac . fem - the PetscFE object to setup
25220cf1dd8SToby Isaac 
25320cf1dd8SToby Isaac   Level: developer
25420cf1dd8SToby Isaac 
25520cf1dd8SToby Isaac .seealso PetscFEView(), PetscFEDestroy()
25620cf1dd8SToby Isaac @*/
25720cf1dd8SToby Isaac PetscErrorCode PetscFESetUp(PetscFE fem)
25820cf1dd8SToby Isaac {
25920cf1dd8SToby Isaac   PetscErrorCode ierr;
26020cf1dd8SToby Isaac 
26120cf1dd8SToby Isaac   PetscFunctionBegin;
26220cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
26320cf1dd8SToby Isaac   if (fem->setupcalled) PetscFunctionReturn(0);
26420cf1dd8SToby Isaac   fem->setupcalled = PETSC_TRUE;
26520cf1dd8SToby Isaac   if (fem->ops->setup) {ierr = (*fem->ops->setup)(fem);CHKERRQ(ierr);}
26620cf1dd8SToby Isaac   PetscFunctionReturn(0);
26720cf1dd8SToby Isaac }
26820cf1dd8SToby Isaac 
26920cf1dd8SToby Isaac /*@
27020cf1dd8SToby Isaac   PetscFEDestroy - Destroys a PetscFE object
27120cf1dd8SToby Isaac 
272d083f849SBarry Smith   Collective on fem
27320cf1dd8SToby Isaac 
27420cf1dd8SToby Isaac   Input Parameter:
27520cf1dd8SToby Isaac . fem - the PetscFE object to destroy
27620cf1dd8SToby Isaac 
27720cf1dd8SToby Isaac   Level: developer
27820cf1dd8SToby Isaac 
27920cf1dd8SToby Isaac .seealso PetscFEView()
28020cf1dd8SToby Isaac @*/
28120cf1dd8SToby Isaac PetscErrorCode PetscFEDestroy(PetscFE *fem)
28220cf1dd8SToby Isaac {
28320cf1dd8SToby Isaac   PetscErrorCode ierr;
28420cf1dd8SToby Isaac 
28520cf1dd8SToby Isaac   PetscFunctionBegin;
28620cf1dd8SToby Isaac   if (!*fem) PetscFunctionReturn(0);
28720cf1dd8SToby Isaac   PetscValidHeaderSpecific((*fem), PETSCFE_CLASSID, 1);
28820cf1dd8SToby Isaac 
28920cf1dd8SToby Isaac   if (--((PetscObject)(*fem))->refct > 0) {*fem = 0; PetscFunctionReturn(0);}
29020cf1dd8SToby Isaac   ((PetscObject) (*fem))->refct = 0;
29120cf1dd8SToby Isaac 
29220cf1dd8SToby Isaac   if ((*fem)->subspaces) {
29320cf1dd8SToby Isaac     PetscInt dim, d;
29420cf1dd8SToby Isaac 
29520cf1dd8SToby Isaac     ierr = PetscDualSpaceGetDimension((*fem)->dualSpace, &dim);CHKERRQ(ierr);
29620cf1dd8SToby Isaac     for (d = 0; d < dim; ++d) {ierr = PetscFEDestroy(&(*fem)->subspaces[d]);CHKERRQ(ierr);}
29720cf1dd8SToby Isaac   }
29820cf1dd8SToby Isaac   ierr = PetscFree((*fem)->subspaces);CHKERRQ(ierr);
29920cf1dd8SToby Isaac   ierr = PetscFree((*fem)->invV);CHKERRQ(ierr);
30020cf1dd8SToby Isaac   ierr = PetscFERestoreTabulation((*fem), 0, NULL, &(*fem)->B, &(*fem)->D, NULL /*&(*fem)->H*/);CHKERRQ(ierr);
30120cf1dd8SToby Isaac   ierr = PetscFERestoreTabulation((*fem), 0, NULL, &(*fem)->Bf, &(*fem)->Df, NULL /*&(*fem)->Hf*/);CHKERRQ(ierr);
30220cf1dd8SToby Isaac   ierr = PetscFERestoreTabulation((*fem), 0, NULL, &(*fem)->F, NULL, NULL);CHKERRQ(ierr);
30320cf1dd8SToby Isaac   ierr = PetscSpaceDestroy(&(*fem)->basisSpace);CHKERRQ(ierr);
30420cf1dd8SToby Isaac   ierr = PetscDualSpaceDestroy(&(*fem)->dualSpace);CHKERRQ(ierr);
30520cf1dd8SToby Isaac   ierr = PetscQuadratureDestroy(&(*fem)->quadrature);CHKERRQ(ierr);
30620cf1dd8SToby Isaac   ierr = PetscQuadratureDestroy(&(*fem)->faceQuadrature);CHKERRQ(ierr);
30720cf1dd8SToby Isaac 
30820cf1dd8SToby Isaac   if ((*fem)->ops->destroy) {ierr = (*(*fem)->ops->destroy)(*fem);CHKERRQ(ierr);}
30920cf1dd8SToby Isaac   ierr = PetscHeaderDestroy(fem);CHKERRQ(ierr);
31020cf1dd8SToby Isaac   PetscFunctionReturn(0);
31120cf1dd8SToby Isaac }
31220cf1dd8SToby Isaac 
31320cf1dd8SToby Isaac /*@
31420cf1dd8SToby Isaac   PetscFECreate - Creates an empty PetscFE object. The type can then be set with PetscFESetType().
31520cf1dd8SToby Isaac 
316d083f849SBarry Smith   Collective
31720cf1dd8SToby Isaac 
31820cf1dd8SToby Isaac   Input Parameter:
31920cf1dd8SToby Isaac . comm - The communicator for the PetscFE object
32020cf1dd8SToby Isaac 
32120cf1dd8SToby Isaac   Output Parameter:
32220cf1dd8SToby Isaac . fem - The PetscFE object
32320cf1dd8SToby Isaac 
32420cf1dd8SToby Isaac   Level: beginner
32520cf1dd8SToby Isaac 
32620cf1dd8SToby Isaac .seealso: PetscFESetType(), PETSCFEGALERKIN
32720cf1dd8SToby Isaac @*/
32820cf1dd8SToby Isaac PetscErrorCode PetscFECreate(MPI_Comm comm, PetscFE *fem)
32920cf1dd8SToby Isaac {
33020cf1dd8SToby Isaac   PetscFE        f;
33120cf1dd8SToby Isaac   PetscErrorCode ierr;
33220cf1dd8SToby Isaac 
33320cf1dd8SToby Isaac   PetscFunctionBegin;
33420cf1dd8SToby Isaac   PetscValidPointer(fem, 2);
33520cf1dd8SToby Isaac   ierr = PetscCitationsRegister(FECitation,&FEcite);CHKERRQ(ierr);
33620cf1dd8SToby Isaac   *fem = NULL;
33720cf1dd8SToby Isaac   ierr = PetscFEInitializePackage();CHKERRQ(ierr);
33820cf1dd8SToby Isaac 
33920cf1dd8SToby Isaac   ierr = PetscHeaderCreate(f, PETSCFE_CLASSID, "PetscFE", "Finite Element", "PetscFE", comm, PetscFEDestroy, PetscFEView);CHKERRQ(ierr);
34020cf1dd8SToby Isaac 
34120cf1dd8SToby Isaac   f->basisSpace    = NULL;
34220cf1dd8SToby Isaac   f->dualSpace     = NULL;
34320cf1dd8SToby Isaac   f->numComponents = 1;
34420cf1dd8SToby Isaac   f->subspaces     = NULL;
34520cf1dd8SToby Isaac   f->invV          = NULL;
34620cf1dd8SToby Isaac   f->B             = NULL;
34720cf1dd8SToby Isaac   f->D             = NULL;
34820cf1dd8SToby Isaac   f->H             = NULL;
34920cf1dd8SToby Isaac   f->Bf            = NULL;
35020cf1dd8SToby Isaac   f->Df            = NULL;
35120cf1dd8SToby Isaac   f->Hf            = NULL;
352580bdb30SBarry Smith   ierr = PetscArrayzero(&f->quadrature, 1);CHKERRQ(ierr);
353580bdb30SBarry Smith   ierr = PetscArrayzero(&f->faceQuadrature, 1);CHKERRQ(ierr);
35420cf1dd8SToby Isaac   f->blockSize     = 0;
35520cf1dd8SToby Isaac   f->numBlocks     = 1;
35620cf1dd8SToby Isaac   f->batchSize     = 0;
35720cf1dd8SToby Isaac   f->numBatches    = 1;
35820cf1dd8SToby Isaac 
35920cf1dd8SToby Isaac   *fem = f;
36020cf1dd8SToby Isaac   PetscFunctionReturn(0);
36120cf1dd8SToby Isaac }
36220cf1dd8SToby Isaac 
36320cf1dd8SToby Isaac /*@
36420cf1dd8SToby Isaac   PetscFEGetSpatialDimension - Returns the spatial dimension of the element
36520cf1dd8SToby Isaac 
36620cf1dd8SToby Isaac   Not collective
36720cf1dd8SToby Isaac 
36820cf1dd8SToby Isaac   Input Parameter:
36920cf1dd8SToby Isaac . fem - The PetscFE object
37020cf1dd8SToby Isaac 
37120cf1dd8SToby Isaac   Output Parameter:
37220cf1dd8SToby Isaac . dim - The spatial dimension
37320cf1dd8SToby Isaac 
37420cf1dd8SToby Isaac   Level: intermediate
37520cf1dd8SToby Isaac 
37620cf1dd8SToby Isaac .seealso: PetscFECreate()
37720cf1dd8SToby Isaac @*/
37820cf1dd8SToby Isaac PetscErrorCode PetscFEGetSpatialDimension(PetscFE fem, PetscInt *dim)
37920cf1dd8SToby Isaac {
38020cf1dd8SToby Isaac   DM             dm;
38120cf1dd8SToby Isaac   PetscErrorCode ierr;
38220cf1dd8SToby Isaac 
38320cf1dd8SToby Isaac   PetscFunctionBegin;
38420cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
38520cf1dd8SToby Isaac   PetscValidPointer(dim, 2);
38620cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDM(fem->dualSpace, &dm);CHKERRQ(ierr);
38720cf1dd8SToby Isaac   ierr = DMGetDimension(dm, dim);CHKERRQ(ierr);
38820cf1dd8SToby Isaac   PetscFunctionReturn(0);
38920cf1dd8SToby Isaac }
39020cf1dd8SToby Isaac 
39120cf1dd8SToby Isaac /*@
39220cf1dd8SToby Isaac   PetscFESetNumComponents - Sets the number of components in the element
39320cf1dd8SToby Isaac 
39420cf1dd8SToby Isaac   Not collective
39520cf1dd8SToby Isaac 
39620cf1dd8SToby Isaac   Input Parameters:
39720cf1dd8SToby Isaac + fem - The PetscFE object
39820cf1dd8SToby Isaac - comp - The number of field components
39920cf1dd8SToby Isaac 
40020cf1dd8SToby Isaac   Level: intermediate
40120cf1dd8SToby Isaac 
40220cf1dd8SToby Isaac .seealso: PetscFECreate()
40320cf1dd8SToby Isaac @*/
40420cf1dd8SToby Isaac PetscErrorCode PetscFESetNumComponents(PetscFE fem, PetscInt comp)
40520cf1dd8SToby Isaac {
40620cf1dd8SToby Isaac   PetscFunctionBegin;
40720cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
40820cf1dd8SToby Isaac   fem->numComponents = comp;
40920cf1dd8SToby Isaac   PetscFunctionReturn(0);
41020cf1dd8SToby Isaac }
41120cf1dd8SToby Isaac 
41220cf1dd8SToby Isaac /*@
41320cf1dd8SToby Isaac   PetscFEGetNumComponents - Returns the number of components in the element
41420cf1dd8SToby Isaac 
41520cf1dd8SToby Isaac   Not collective
41620cf1dd8SToby Isaac 
41720cf1dd8SToby Isaac   Input Parameter:
41820cf1dd8SToby Isaac . fem - The PetscFE object
41920cf1dd8SToby Isaac 
42020cf1dd8SToby Isaac   Output Parameter:
42120cf1dd8SToby Isaac . comp - The number of field components
42220cf1dd8SToby Isaac 
42320cf1dd8SToby Isaac   Level: intermediate
42420cf1dd8SToby Isaac 
42520cf1dd8SToby Isaac .seealso: PetscFECreate()
42620cf1dd8SToby Isaac @*/
42720cf1dd8SToby Isaac PetscErrorCode PetscFEGetNumComponents(PetscFE fem, PetscInt *comp)
42820cf1dd8SToby Isaac {
42920cf1dd8SToby Isaac   PetscFunctionBegin;
43020cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
43120cf1dd8SToby Isaac   PetscValidPointer(comp, 2);
43220cf1dd8SToby Isaac   *comp = fem->numComponents;
43320cf1dd8SToby Isaac   PetscFunctionReturn(0);
43420cf1dd8SToby Isaac }
43520cf1dd8SToby Isaac 
43620cf1dd8SToby Isaac /*@
43720cf1dd8SToby Isaac   PetscFESetTileSizes - Sets the tile sizes for evaluation
43820cf1dd8SToby Isaac 
43920cf1dd8SToby Isaac   Not collective
44020cf1dd8SToby Isaac 
44120cf1dd8SToby Isaac   Input Parameters:
44220cf1dd8SToby Isaac + fem - The PetscFE object
44320cf1dd8SToby Isaac . blockSize - The number of elements in a block
44420cf1dd8SToby Isaac . numBlocks - The number of blocks in a batch
44520cf1dd8SToby Isaac . batchSize - The number of elements in a batch
44620cf1dd8SToby Isaac - numBatches - The number of batches in a chunk
44720cf1dd8SToby Isaac 
44820cf1dd8SToby Isaac   Level: intermediate
44920cf1dd8SToby Isaac 
45020cf1dd8SToby Isaac .seealso: PetscFECreate()
45120cf1dd8SToby Isaac @*/
45220cf1dd8SToby Isaac PetscErrorCode PetscFESetTileSizes(PetscFE fem, PetscInt blockSize, PetscInt numBlocks, PetscInt batchSize, PetscInt numBatches)
45320cf1dd8SToby Isaac {
45420cf1dd8SToby Isaac   PetscFunctionBegin;
45520cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
45620cf1dd8SToby Isaac   fem->blockSize  = blockSize;
45720cf1dd8SToby Isaac   fem->numBlocks  = numBlocks;
45820cf1dd8SToby Isaac   fem->batchSize  = batchSize;
45920cf1dd8SToby Isaac   fem->numBatches = numBatches;
46020cf1dd8SToby Isaac   PetscFunctionReturn(0);
46120cf1dd8SToby Isaac }
46220cf1dd8SToby Isaac 
46320cf1dd8SToby Isaac /*@
46420cf1dd8SToby Isaac   PetscFEGetTileSizes - Returns the tile sizes for evaluation
46520cf1dd8SToby Isaac 
46620cf1dd8SToby Isaac   Not collective
46720cf1dd8SToby Isaac 
46820cf1dd8SToby Isaac   Input Parameter:
46920cf1dd8SToby Isaac . fem - The PetscFE object
47020cf1dd8SToby Isaac 
47120cf1dd8SToby Isaac   Output Parameters:
47220cf1dd8SToby Isaac + blockSize - The number of elements in a block
47320cf1dd8SToby Isaac . numBlocks - The number of blocks in a batch
47420cf1dd8SToby Isaac . batchSize - The number of elements in a batch
47520cf1dd8SToby Isaac - numBatches - The number of batches in a chunk
47620cf1dd8SToby Isaac 
47720cf1dd8SToby Isaac   Level: intermediate
47820cf1dd8SToby Isaac 
47920cf1dd8SToby Isaac .seealso: PetscFECreate()
48020cf1dd8SToby Isaac @*/
48120cf1dd8SToby Isaac PetscErrorCode PetscFEGetTileSizes(PetscFE fem, PetscInt *blockSize, PetscInt *numBlocks, PetscInt *batchSize, PetscInt *numBatches)
48220cf1dd8SToby Isaac {
48320cf1dd8SToby Isaac   PetscFunctionBegin;
48420cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
48520cf1dd8SToby Isaac   if (blockSize)  PetscValidPointer(blockSize,  2);
48620cf1dd8SToby Isaac   if (numBlocks)  PetscValidPointer(numBlocks,  3);
48720cf1dd8SToby Isaac   if (batchSize)  PetscValidPointer(batchSize,  4);
48820cf1dd8SToby Isaac   if (numBatches) PetscValidPointer(numBatches, 5);
48920cf1dd8SToby Isaac   if (blockSize)  *blockSize  = fem->blockSize;
49020cf1dd8SToby Isaac   if (numBlocks)  *numBlocks  = fem->numBlocks;
49120cf1dd8SToby Isaac   if (batchSize)  *batchSize  = fem->batchSize;
49220cf1dd8SToby Isaac   if (numBatches) *numBatches = fem->numBatches;
49320cf1dd8SToby Isaac   PetscFunctionReturn(0);
49420cf1dd8SToby Isaac }
49520cf1dd8SToby Isaac 
49620cf1dd8SToby Isaac /*@
49720cf1dd8SToby Isaac   PetscFEGetBasisSpace - Returns the PetscSpace used for approximation of the solution
49820cf1dd8SToby Isaac 
49920cf1dd8SToby Isaac   Not collective
50020cf1dd8SToby Isaac 
50120cf1dd8SToby Isaac   Input Parameter:
50220cf1dd8SToby Isaac . fem - The PetscFE object
50320cf1dd8SToby Isaac 
50420cf1dd8SToby Isaac   Output Parameter:
50520cf1dd8SToby Isaac . sp - The PetscSpace object
50620cf1dd8SToby Isaac 
50720cf1dd8SToby Isaac   Level: intermediate
50820cf1dd8SToby Isaac 
50920cf1dd8SToby Isaac .seealso: PetscFECreate()
51020cf1dd8SToby Isaac @*/
51120cf1dd8SToby Isaac PetscErrorCode PetscFEGetBasisSpace(PetscFE fem, PetscSpace *sp)
51220cf1dd8SToby Isaac {
51320cf1dd8SToby Isaac   PetscFunctionBegin;
51420cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
51520cf1dd8SToby Isaac   PetscValidPointer(sp, 2);
51620cf1dd8SToby Isaac   *sp = fem->basisSpace;
51720cf1dd8SToby Isaac   PetscFunctionReturn(0);
51820cf1dd8SToby Isaac }
51920cf1dd8SToby Isaac 
52020cf1dd8SToby Isaac /*@
52120cf1dd8SToby Isaac   PetscFESetBasisSpace - Sets the PetscSpace used for approximation of the solution
52220cf1dd8SToby Isaac 
52320cf1dd8SToby Isaac   Not collective
52420cf1dd8SToby Isaac 
52520cf1dd8SToby Isaac   Input Parameters:
52620cf1dd8SToby Isaac + fem - The PetscFE object
52720cf1dd8SToby Isaac - sp - The PetscSpace object
52820cf1dd8SToby Isaac 
52920cf1dd8SToby Isaac   Level: intermediate
53020cf1dd8SToby Isaac 
53120cf1dd8SToby Isaac .seealso: PetscFECreate()
53220cf1dd8SToby Isaac @*/
53320cf1dd8SToby Isaac PetscErrorCode PetscFESetBasisSpace(PetscFE fem, PetscSpace sp)
53420cf1dd8SToby Isaac {
53520cf1dd8SToby Isaac   PetscErrorCode ierr;
53620cf1dd8SToby Isaac 
53720cf1dd8SToby Isaac   PetscFunctionBegin;
53820cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
53920cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 2);
54020cf1dd8SToby Isaac   ierr = PetscSpaceDestroy(&fem->basisSpace);CHKERRQ(ierr);
54120cf1dd8SToby Isaac   fem->basisSpace = sp;
54220cf1dd8SToby Isaac   ierr = PetscObjectReference((PetscObject) fem->basisSpace);CHKERRQ(ierr);
54320cf1dd8SToby Isaac   PetscFunctionReturn(0);
54420cf1dd8SToby Isaac }
54520cf1dd8SToby Isaac 
54620cf1dd8SToby Isaac /*@
54720cf1dd8SToby Isaac   PetscFEGetDualSpace - Returns the PetscDualSpace used to define the inner product
54820cf1dd8SToby Isaac 
54920cf1dd8SToby Isaac   Not collective
55020cf1dd8SToby Isaac 
55120cf1dd8SToby Isaac   Input Parameter:
55220cf1dd8SToby Isaac . fem - The PetscFE object
55320cf1dd8SToby Isaac 
55420cf1dd8SToby Isaac   Output Parameter:
55520cf1dd8SToby Isaac . sp - The PetscDualSpace object
55620cf1dd8SToby Isaac 
55720cf1dd8SToby Isaac   Level: intermediate
55820cf1dd8SToby Isaac 
55920cf1dd8SToby Isaac .seealso: PetscFECreate()
56020cf1dd8SToby Isaac @*/
56120cf1dd8SToby Isaac PetscErrorCode PetscFEGetDualSpace(PetscFE fem, PetscDualSpace *sp)
56220cf1dd8SToby Isaac {
56320cf1dd8SToby Isaac   PetscFunctionBegin;
56420cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
56520cf1dd8SToby Isaac   PetscValidPointer(sp, 2);
56620cf1dd8SToby Isaac   *sp = fem->dualSpace;
56720cf1dd8SToby Isaac   PetscFunctionReturn(0);
56820cf1dd8SToby Isaac }
56920cf1dd8SToby Isaac 
57020cf1dd8SToby Isaac /*@
57120cf1dd8SToby Isaac   PetscFESetDualSpace - Sets the PetscDualSpace used to define the inner product
57220cf1dd8SToby Isaac 
57320cf1dd8SToby Isaac   Not collective
57420cf1dd8SToby Isaac 
57520cf1dd8SToby Isaac   Input Parameters:
57620cf1dd8SToby Isaac + fem - The PetscFE object
57720cf1dd8SToby Isaac - sp - The PetscDualSpace object
57820cf1dd8SToby Isaac 
57920cf1dd8SToby Isaac   Level: intermediate
58020cf1dd8SToby Isaac 
58120cf1dd8SToby Isaac .seealso: PetscFECreate()
58220cf1dd8SToby Isaac @*/
58320cf1dd8SToby Isaac PetscErrorCode PetscFESetDualSpace(PetscFE fem, PetscDualSpace sp)
58420cf1dd8SToby Isaac {
58520cf1dd8SToby Isaac   PetscErrorCode ierr;
58620cf1dd8SToby Isaac 
58720cf1dd8SToby Isaac   PetscFunctionBegin;
58820cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
58920cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 2);
59020cf1dd8SToby Isaac   ierr = PetscDualSpaceDestroy(&fem->dualSpace);CHKERRQ(ierr);
59120cf1dd8SToby Isaac   fem->dualSpace = sp;
59220cf1dd8SToby Isaac   ierr = PetscObjectReference((PetscObject) fem->dualSpace);CHKERRQ(ierr);
59320cf1dd8SToby Isaac   PetscFunctionReturn(0);
59420cf1dd8SToby Isaac }
59520cf1dd8SToby Isaac 
59620cf1dd8SToby Isaac /*@
59720cf1dd8SToby Isaac   PetscFEGetQuadrature - Returns the PetscQuadrature used to calculate inner products
59820cf1dd8SToby Isaac 
59920cf1dd8SToby Isaac   Not collective
60020cf1dd8SToby Isaac 
60120cf1dd8SToby Isaac   Input Parameter:
60220cf1dd8SToby Isaac . fem - The PetscFE object
60320cf1dd8SToby Isaac 
60420cf1dd8SToby Isaac   Output Parameter:
60520cf1dd8SToby Isaac . q - The PetscQuadrature object
60620cf1dd8SToby Isaac 
60720cf1dd8SToby Isaac   Level: intermediate
60820cf1dd8SToby Isaac 
60920cf1dd8SToby Isaac .seealso: PetscFECreate()
61020cf1dd8SToby Isaac @*/
61120cf1dd8SToby Isaac PetscErrorCode PetscFEGetQuadrature(PetscFE fem, PetscQuadrature *q)
61220cf1dd8SToby Isaac {
61320cf1dd8SToby Isaac   PetscFunctionBegin;
61420cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
61520cf1dd8SToby Isaac   PetscValidPointer(q, 2);
61620cf1dd8SToby Isaac   *q = fem->quadrature;
61720cf1dd8SToby Isaac   PetscFunctionReturn(0);
61820cf1dd8SToby Isaac }
61920cf1dd8SToby Isaac 
62020cf1dd8SToby Isaac /*@
62120cf1dd8SToby Isaac   PetscFESetQuadrature - Sets the PetscQuadrature used to calculate inner products
62220cf1dd8SToby Isaac 
62320cf1dd8SToby Isaac   Not collective
62420cf1dd8SToby Isaac 
62520cf1dd8SToby Isaac   Input Parameters:
62620cf1dd8SToby Isaac + fem - The PetscFE object
62720cf1dd8SToby Isaac - q - The PetscQuadrature object
62820cf1dd8SToby Isaac 
62920cf1dd8SToby Isaac   Level: intermediate
63020cf1dd8SToby Isaac 
63120cf1dd8SToby Isaac .seealso: PetscFECreate()
63220cf1dd8SToby Isaac @*/
63320cf1dd8SToby Isaac PetscErrorCode PetscFESetQuadrature(PetscFE fem, PetscQuadrature q)
63420cf1dd8SToby Isaac {
63520cf1dd8SToby Isaac   PetscInt       Nc, qNc;
63620cf1dd8SToby Isaac   PetscErrorCode ierr;
63720cf1dd8SToby Isaac 
63820cf1dd8SToby Isaac   PetscFunctionBegin;
63920cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
64020cf1dd8SToby Isaac   ierr = PetscFEGetNumComponents(fem, &Nc);CHKERRQ(ierr);
64120cf1dd8SToby Isaac   ierr = PetscQuadratureGetNumComponents(q, &qNc);CHKERRQ(ierr);
64220cf1dd8SToby 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);
64320cf1dd8SToby Isaac   ierr = PetscFERestoreTabulation(fem, 0, NULL, &fem->B, &fem->D, NULL /*&(*fem)->H*/);CHKERRQ(ierr);
64420cf1dd8SToby Isaac   ierr = PetscQuadratureDestroy(&fem->quadrature);CHKERRQ(ierr);
64520cf1dd8SToby Isaac   fem->quadrature = q;
64620cf1dd8SToby Isaac   ierr = PetscObjectReference((PetscObject) q);CHKERRQ(ierr);
64720cf1dd8SToby Isaac   PetscFunctionReturn(0);
64820cf1dd8SToby Isaac }
64920cf1dd8SToby Isaac 
65020cf1dd8SToby Isaac /*@
65120cf1dd8SToby Isaac   PetscFEGetFaceQuadrature - Returns the PetscQuadrature used to calculate inner products on faces
65220cf1dd8SToby Isaac 
65320cf1dd8SToby Isaac   Not collective
65420cf1dd8SToby Isaac 
65520cf1dd8SToby Isaac   Input Parameter:
65620cf1dd8SToby Isaac . fem - The PetscFE object
65720cf1dd8SToby Isaac 
65820cf1dd8SToby Isaac   Output Parameter:
65920cf1dd8SToby Isaac . q - The PetscQuadrature object
66020cf1dd8SToby Isaac 
66120cf1dd8SToby Isaac   Level: intermediate
66220cf1dd8SToby Isaac 
66320cf1dd8SToby Isaac .seealso: PetscFECreate()
66420cf1dd8SToby Isaac @*/
66520cf1dd8SToby Isaac PetscErrorCode PetscFEGetFaceQuadrature(PetscFE fem, PetscQuadrature *q)
66620cf1dd8SToby Isaac {
66720cf1dd8SToby Isaac   PetscFunctionBegin;
66820cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
66920cf1dd8SToby Isaac   PetscValidPointer(q, 2);
67020cf1dd8SToby Isaac   *q = fem->faceQuadrature;
67120cf1dd8SToby Isaac   PetscFunctionReturn(0);
67220cf1dd8SToby Isaac }
67320cf1dd8SToby Isaac 
67420cf1dd8SToby Isaac /*@
67520cf1dd8SToby Isaac   PetscFESetFaceQuadrature - Sets the PetscQuadrature used to calculate inner products on faces
67620cf1dd8SToby Isaac 
67720cf1dd8SToby Isaac   Not collective
67820cf1dd8SToby Isaac 
67920cf1dd8SToby Isaac   Input Parameters:
68020cf1dd8SToby Isaac + fem - The PetscFE object
68120cf1dd8SToby Isaac - q - The PetscQuadrature object
68220cf1dd8SToby Isaac 
68320cf1dd8SToby Isaac   Level: intermediate
68420cf1dd8SToby Isaac 
68520cf1dd8SToby Isaac .seealso: PetscFECreate()
68620cf1dd8SToby Isaac @*/
68720cf1dd8SToby Isaac PetscErrorCode PetscFESetFaceQuadrature(PetscFE fem, PetscQuadrature q)
68820cf1dd8SToby Isaac {
68920cf1dd8SToby Isaac   PetscErrorCode ierr;
69020cf1dd8SToby Isaac 
69120cf1dd8SToby Isaac   PetscFunctionBegin;
69220cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
69320cf1dd8SToby Isaac   ierr = PetscFERestoreTabulation(fem, 0, NULL, &fem->Bf, &fem->Df, NULL /*&(*fem)->Hf*/);CHKERRQ(ierr);
69420cf1dd8SToby Isaac   ierr = PetscQuadratureDestroy(&fem->faceQuadrature);CHKERRQ(ierr);
69520cf1dd8SToby Isaac   fem->faceQuadrature = q;
69620cf1dd8SToby Isaac   ierr = PetscObjectReference((PetscObject) q);CHKERRQ(ierr);
69720cf1dd8SToby Isaac   PetscFunctionReturn(0);
69820cf1dd8SToby Isaac }
69920cf1dd8SToby Isaac 
70020cf1dd8SToby Isaac /*@C
70120cf1dd8SToby Isaac   PetscFEGetNumDof - Returns the number of dofs (dual basis vectors) associated to mesh points on the reference cell of a given dimension
70220cf1dd8SToby Isaac 
70320cf1dd8SToby Isaac   Not collective
70420cf1dd8SToby Isaac 
70520cf1dd8SToby Isaac   Input Parameter:
70620cf1dd8SToby Isaac . fem - The PetscFE object
70720cf1dd8SToby Isaac 
70820cf1dd8SToby Isaac   Output Parameter:
70920cf1dd8SToby Isaac . numDof - Array with the number of dofs per dimension
71020cf1dd8SToby Isaac 
71120cf1dd8SToby Isaac   Level: intermediate
71220cf1dd8SToby Isaac 
71320cf1dd8SToby Isaac .seealso: PetscFECreate()
71420cf1dd8SToby Isaac @*/
71520cf1dd8SToby Isaac PetscErrorCode PetscFEGetNumDof(PetscFE fem, const PetscInt **numDof)
71620cf1dd8SToby Isaac {
71720cf1dd8SToby Isaac   PetscErrorCode ierr;
71820cf1dd8SToby Isaac 
71920cf1dd8SToby Isaac   PetscFunctionBegin;
72020cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
72120cf1dd8SToby Isaac   PetscValidPointer(numDof, 2);
72220cf1dd8SToby Isaac   ierr = PetscDualSpaceGetNumDof(fem->dualSpace, numDof);CHKERRQ(ierr);
72320cf1dd8SToby Isaac   PetscFunctionReturn(0);
72420cf1dd8SToby Isaac }
72520cf1dd8SToby Isaac 
72620cf1dd8SToby Isaac /*@C
72720cf1dd8SToby Isaac   PetscFEGetDefaultTabulation - Returns the tabulation of the basis functions at the quadrature points
72820cf1dd8SToby Isaac 
72920cf1dd8SToby Isaac   Not collective
73020cf1dd8SToby Isaac 
73120cf1dd8SToby Isaac   Input Parameter:
73220cf1dd8SToby Isaac . fem - The PetscFE object
73320cf1dd8SToby Isaac 
73420cf1dd8SToby Isaac   Output Parameters:
73520cf1dd8SToby Isaac + B - The basis function values at quadrature points
73620cf1dd8SToby Isaac . D - The basis function derivatives at quadrature points
73720cf1dd8SToby Isaac - H - The basis function second derivatives at quadrature points
73820cf1dd8SToby Isaac 
73920cf1dd8SToby Isaac   Note:
74020cf1dd8SToby Isaac $ B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c
74120cf1dd8SToby 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
74220cf1dd8SToby 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
74320cf1dd8SToby Isaac 
74420cf1dd8SToby Isaac   Level: intermediate
74520cf1dd8SToby Isaac 
74620cf1dd8SToby Isaac .seealso: PetscFEGetTabulation(), PetscFERestoreTabulation()
74720cf1dd8SToby Isaac @*/
74820cf1dd8SToby Isaac PetscErrorCode PetscFEGetDefaultTabulation(PetscFE fem, PetscReal **B, PetscReal **D, PetscReal **H)
74920cf1dd8SToby Isaac {
75020cf1dd8SToby Isaac   PetscInt         npoints;
75120cf1dd8SToby Isaac   const PetscReal *points;
75220cf1dd8SToby Isaac   PetscErrorCode   ierr;
75320cf1dd8SToby Isaac 
75420cf1dd8SToby Isaac   PetscFunctionBegin;
75520cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
75620cf1dd8SToby Isaac   if (B) PetscValidPointer(B, 2);
75720cf1dd8SToby Isaac   if (D) PetscValidPointer(D, 3);
75820cf1dd8SToby Isaac   if (H) PetscValidPointer(H, 4);
75920cf1dd8SToby Isaac   ierr = PetscQuadratureGetData(fem->quadrature, NULL, NULL, &npoints, &points, NULL);CHKERRQ(ierr);
76020cf1dd8SToby Isaac   if (!fem->B) {ierr = PetscFEGetTabulation(fem, npoints, points, &fem->B, &fem->D, NULL/*&fem->H*/);CHKERRQ(ierr);}
76120cf1dd8SToby Isaac   if (B) *B = fem->B;
76220cf1dd8SToby Isaac   if (D) *D = fem->D;
76320cf1dd8SToby Isaac   if (H) *H = fem->H;
76420cf1dd8SToby Isaac   PetscFunctionReturn(0);
76520cf1dd8SToby Isaac }
76620cf1dd8SToby Isaac 
76720cf1dd8SToby Isaac PetscErrorCode PetscFEGetFaceTabulation(PetscFE fem, PetscReal **Bf, PetscReal **Df, PetscReal **Hf)
76820cf1dd8SToby Isaac {
76920cf1dd8SToby Isaac   PetscErrorCode   ierr;
77020cf1dd8SToby Isaac 
77120cf1dd8SToby Isaac   PetscFunctionBegin;
77220cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
77320cf1dd8SToby Isaac   if (Bf) PetscValidPointer(Bf, 2);
77420cf1dd8SToby Isaac   if (Df) PetscValidPointer(Df, 3);
77520cf1dd8SToby Isaac   if (Hf) PetscValidPointer(Hf, 4);
77620cf1dd8SToby Isaac   if (!fem->Bf) {
77720cf1dd8SToby Isaac     const PetscReal  xi0[3] = {-1., -1., -1.};
77820cf1dd8SToby Isaac     PetscReal        v0[3], J[9], detJ;
77920cf1dd8SToby Isaac     PetscQuadrature  fq;
78020cf1dd8SToby Isaac     PetscDualSpace   sp;
78120cf1dd8SToby Isaac     DM               dm;
78220cf1dd8SToby Isaac     const PetscInt  *faces;
78320cf1dd8SToby Isaac     PetscInt         dim, numFaces, f, npoints, q;
78420cf1dd8SToby Isaac     const PetscReal *points;
78520cf1dd8SToby Isaac     PetscReal       *facePoints;
78620cf1dd8SToby Isaac 
78720cf1dd8SToby Isaac     ierr = PetscFEGetDualSpace(fem, &sp);CHKERRQ(ierr);
78820cf1dd8SToby Isaac     ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr);
78920cf1dd8SToby Isaac     ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
79020cf1dd8SToby Isaac     ierr = DMPlexGetConeSize(dm, 0, &numFaces);CHKERRQ(ierr);
79120cf1dd8SToby Isaac     ierr = DMPlexGetCone(dm, 0, &faces);CHKERRQ(ierr);
79220cf1dd8SToby Isaac     ierr = PetscFEGetFaceQuadrature(fem, &fq);CHKERRQ(ierr);
79320cf1dd8SToby Isaac     if (fq) {
79420cf1dd8SToby Isaac       ierr = PetscQuadratureGetData(fq, NULL, NULL, &npoints, &points, NULL);CHKERRQ(ierr);
79520cf1dd8SToby Isaac       ierr = PetscMalloc1(numFaces*npoints*dim, &facePoints);CHKERRQ(ierr);
79620cf1dd8SToby Isaac       for (f = 0; f < numFaces; ++f) {
79720cf1dd8SToby Isaac         ierr = DMPlexComputeCellGeometryFEM(dm, faces[f], NULL, v0, J, NULL, &detJ);CHKERRQ(ierr);
79820cf1dd8SToby Isaac         for (q = 0; q < npoints; ++q) CoordinatesRefToReal(dim, dim-1, xi0, v0, J, &points[q*(dim-1)], &facePoints[(f*npoints+q)*dim]);
79920cf1dd8SToby Isaac       }
80020cf1dd8SToby Isaac       ierr = PetscFEGetTabulation(fem, numFaces*npoints, facePoints, &fem->Bf, &fem->Df, NULL/*&fem->Hf*/);CHKERRQ(ierr);
80120cf1dd8SToby Isaac       ierr = PetscFree(facePoints);CHKERRQ(ierr);
80220cf1dd8SToby Isaac     }
80320cf1dd8SToby Isaac   }
80420cf1dd8SToby Isaac   if (Bf) *Bf = fem->Bf;
80520cf1dd8SToby Isaac   if (Df) *Df = fem->Df;
80620cf1dd8SToby Isaac   if (Hf) *Hf = fem->Hf;
80720cf1dd8SToby Isaac   PetscFunctionReturn(0);
80820cf1dd8SToby Isaac }
80920cf1dd8SToby Isaac 
81020cf1dd8SToby Isaac PetscErrorCode PetscFEGetFaceCentroidTabulation(PetscFE fem, PetscReal **F)
81120cf1dd8SToby Isaac {
81220cf1dd8SToby Isaac   PetscErrorCode   ierr;
81320cf1dd8SToby Isaac 
81420cf1dd8SToby Isaac   PetscFunctionBegin;
81520cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
81620cf1dd8SToby Isaac   PetscValidPointer(F, 2);
81720cf1dd8SToby Isaac   if (!fem->F) {
81820cf1dd8SToby Isaac     PetscDualSpace  sp;
81920cf1dd8SToby Isaac     DM              dm;
82020cf1dd8SToby Isaac     const PetscInt *cone;
82120cf1dd8SToby Isaac     PetscReal      *centroids;
82220cf1dd8SToby Isaac     PetscInt        dim, numFaces, f;
82320cf1dd8SToby Isaac 
82420cf1dd8SToby Isaac     ierr = PetscFEGetDualSpace(fem, &sp);CHKERRQ(ierr);
82520cf1dd8SToby Isaac     ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr);
82620cf1dd8SToby Isaac     ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
82720cf1dd8SToby Isaac     ierr = DMPlexGetConeSize(dm, 0, &numFaces);CHKERRQ(ierr);
82820cf1dd8SToby Isaac     ierr = DMPlexGetCone(dm, 0, &cone);CHKERRQ(ierr);
82920cf1dd8SToby Isaac     ierr = PetscMalloc1(numFaces*dim, &centroids);CHKERRQ(ierr);
83020cf1dd8SToby Isaac     for (f = 0; f < numFaces; ++f) {ierr = DMPlexComputeCellGeometryFVM(dm, cone[f], NULL, &centroids[f*dim], NULL);CHKERRQ(ierr);}
83120cf1dd8SToby Isaac     ierr = PetscFEGetTabulation(fem, numFaces, centroids, &fem->F, NULL, NULL);CHKERRQ(ierr);
83220cf1dd8SToby Isaac     ierr = PetscFree(centroids);CHKERRQ(ierr);
83320cf1dd8SToby Isaac   }
83420cf1dd8SToby Isaac   *F = fem->F;
83520cf1dd8SToby Isaac   PetscFunctionReturn(0);
83620cf1dd8SToby Isaac }
83720cf1dd8SToby Isaac 
83820cf1dd8SToby Isaac /*@C
83920cf1dd8SToby Isaac   PetscFEGetTabulation - Tabulates the basis functions, and perhaps derivatives, at the points provided.
84020cf1dd8SToby Isaac 
84120cf1dd8SToby Isaac   Not collective
84220cf1dd8SToby Isaac 
84320cf1dd8SToby Isaac   Input Parameters:
84420cf1dd8SToby Isaac + fem     - The PetscFE object
84520cf1dd8SToby Isaac . npoints - The number of tabulation points
84620cf1dd8SToby Isaac - points  - The tabulation point coordinates
84720cf1dd8SToby Isaac 
84820cf1dd8SToby Isaac   Output Parameters:
84920cf1dd8SToby Isaac + B - The basis function values at tabulation points
85020cf1dd8SToby Isaac . D - The basis function derivatives at tabulation points
85120cf1dd8SToby Isaac - H - The basis function second derivatives at tabulation points
85220cf1dd8SToby Isaac 
85320cf1dd8SToby Isaac   Note:
85420cf1dd8SToby Isaac $ B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c
85520cf1dd8SToby 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
85620cf1dd8SToby 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
85720cf1dd8SToby Isaac 
85820cf1dd8SToby Isaac   Level: intermediate
85920cf1dd8SToby Isaac 
86020cf1dd8SToby Isaac .seealso: PetscFERestoreTabulation(), PetscFEGetDefaultTabulation()
86120cf1dd8SToby Isaac @*/
86220cf1dd8SToby Isaac PetscErrorCode PetscFEGetTabulation(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscReal **B, PetscReal **D, PetscReal **H)
86320cf1dd8SToby Isaac {
86420cf1dd8SToby Isaac   DM               dm;
86520cf1dd8SToby Isaac   PetscInt         pdim; /* Dimension of FE space P */
86620cf1dd8SToby Isaac   PetscInt         dim;  /* Spatial dimension */
86720cf1dd8SToby Isaac   PetscInt         comp; /* Field components */
86820cf1dd8SToby Isaac   PetscErrorCode   ierr;
86920cf1dd8SToby Isaac 
87020cf1dd8SToby Isaac   PetscFunctionBegin;
87120cf1dd8SToby Isaac   if (!npoints) {
87220cf1dd8SToby Isaac     if (B) *B = NULL;
87320cf1dd8SToby Isaac     if (D) *D = NULL;
87420cf1dd8SToby Isaac     if (H) *H = NULL;
87520cf1dd8SToby Isaac     PetscFunctionReturn(0);
87620cf1dd8SToby Isaac   }
87720cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
87820cf1dd8SToby Isaac   PetscValidPointer(points, 3);
87920cf1dd8SToby Isaac   if (B) PetscValidPointer(B, 4);
88020cf1dd8SToby Isaac   if (D) PetscValidPointer(D, 5);
88120cf1dd8SToby Isaac   if (H) PetscValidPointer(H, 6);
88220cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDM(fem->dualSpace, &dm);CHKERRQ(ierr);
88320cf1dd8SToby Isaac   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
88420cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDimension(fem->dualSpace, &pdim);CHKERRQ(ierr);
88520cf1dd8SToby Isaac   ierr = PetscFEGetNumComponents(fem, &comp);CHKERRQ(ierr);
88620cf1dd8SToby Isaac   if (B) {ierr = DMGetWorkArray(dm, npoints*pdim*comp, MPIU_REAL, B);CHKERRQ(ierr);}
88720cf1dd8SToby Isaac   if (!dim) {
88820cf1dd8SToby Isaac     if (D) *D = NULL;
88920cf1dd8SToby Isaac     if (H) *H = NULL;
89020cf1dd8SToby Isaac   } else {
89120cf1dd8SToby Isaac     if (D) {ierr = DMGetWorkArray(dm, npoints*pdim*comp*dim, MPIU_REAL, D);CHKERRQ(ierr);}
89220cf1dd8SToby Isaac     if (H) {ierr = DMGetWorkArray(dm, npoints*pdim*comp*dim*dim, MPIU_REAL, H);CHKERRQ(ierr);}
89320cf1dd8SToby Isaac   }
89420cf1dd8SToby Isaac   ierr = (*fem->ops->gettabulation)(fem, npoints, points, B ? *B : NULL, D ? *D : NULL, H ? *H : NULL);CHKERRQ(ierr);
89520cf1dd8SToby Isaac   PetscFunctionReturn(0);
89620cf1dd8SToby Isaac }
89720cf1dd8SToby Isaac 
89820cf1dd8SToby Isaac PetscErrorCode PetscFERestoreTabulation(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscReal **B, PetscReal **D, PetscReal **H)
89920cf1dd8SToby Isaac {
90020cf1dd8SToby Isaac   DM             dm;
90120cf1dd8SToby Isaac   PetscErrorCode ierr;
90220cf1dd8SToby Isaac 
90320cf1dd8SToby Isaac   PetscFunctionBegin;
90420cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
90520cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDM(fem->dualSpace, &dm);CHKERRQ(ierr);
90620cf1dd8SToby Isaac   if (B && *B) {ierr = DMRestoreWorkArray(dm, 0, MPIU_REAL, B);CHKERRQ(ierr);}
90720cf1dd8SToby Isaac   if (D && *D) {ierr = DMRestoreWorkArray(dm, 0, MPIU_REAL, D);CHKERRQ(ierr);}
90820cf1dd8SToby Isaac   if (H && *H) {ierr = DMRestoreWorkArray(dm, 0, MPIU_REAL, H);CHKERRQ(ierr);}
90920cf1dd8SToby Isaac   PetscFunctionReturn(0);
91020cf1dd8SToby Isaac }
91120cf1dd8SToby Isaac 
91220cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreatePointTrace(PetscFE fe, PetscInt refPoint, PetscFE *trFE)
91320cf1dd8SToby Isaac {
91420cf1dd8SToby Isaac   PetscSpace     bsp, bsubsp;
91520cf1dd8SToby Isaac   PetscDualSpace dsp, dsubsp;
91620cf1dd8SToby Isaac   PetscInt       dim, depth, numComp, i, j, coneSize, order;
91720cf1dd8SToby Isaac   PetscFEType    type;
91820cf1dd8SToby Isaac   DM             dm;
91920cf1dd8SToby Isaac   DMLabel        label;
92020cf1dd8SToby Isaac   PetscReal      *xi, *v, *J, detJ;
921db11e2ebSMatthew G. Knepley   const char     *name;
92220cf1dd8SToby Isaac   PetscQuadrature origin, fullQuad, subQuad;
92320cf1dd8SToby Isaac   PetscErrorCode ierr;
92420cf1dd8SToby Isaac 
92520cf1dd8SToby Isaac   PetscFunctionBegin;
92620cf1dd8SToby Isaac   PetscValidHeaderSpecific(fe,PETSCFE_CLASSID,1);
92720cf1dd8SToby Isaac   PetscValidPointer(trFE,3);
92820cf1dd8SToby Isaac   ierr = PetscFEGetBasisSpace(fe,&bsp);CHKERRQ(ierr);
92920cf1dd8SToby Isaac   ierr = PetscFEGetDualSpace(fe,&dsp);CHKERRQ(ierr);
93020cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDM(dsp,&dm);CHKERRQ(ierr);
93120cf1dd8SToby Isaac   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
93220cf1dd8SToby Isaac   ierr = DMPlexGetDepthLabel(dm,&label);CHKERRQ(ierr);
93320cf1dd8SToby Isaac   ierr = DMLabelGetValue(label,refPoint,&depth);CHKERRQ(ierr);
93420cf1dd8SToby Isaac   ierr = PetscCalloc1(depth,&xi);CHKERRQ(ierr);
93520cf1dd8SToby Isaac   ierr = PetscMalloc1(dim,&v);CHKERRQ(ierr);
93620cf1dd8SToby Isaac   ierr = PetscMalloc1(dim*dim,&J);CHKERRQ(ierr);
93720cf1dd8SToby Isaac   for (i = 0; i < depth; i++) xi[i] = 0.;
93820cf1dd8SToby Isaac   ierr = PetscQuadratureCreate(PETSC_COMM_SELF,&origin);CHKERRQ(ierr);
93920cf1dd8SToby Isaac   ierr = PetscQuadratureSetData(origin,depth,0,1,xi,NULL);CHKERRQ(ierr);
94020cf1dd8SToby Isaac   ierr = DMPlexComputeCellGeometryFEM(dm,refPoint,origin,v,J,NULL,&detJ);CHKERRQ(ierr);
94120cf1dd8SToby Isaac   /* CellGeometryFEM computes the expanded Jacobian, we want the true jacobian */
94220cf1dd8SToby Isaac   for (i = 1; i < dim; i++) {
94320cf1dd8SToby Isaac     for (j = 0; j < depth; j++) {
94420cf1dd8SToby Isaac       J[i * depth + j] = J[i * dim + j];
94520cf1dd8SToby Isaac     }
94620cf1dd8SToby Isaac   }
94720cf1dd8SToby Isaac   ierr = PetscQuadratureDestroy(&origin);CHKERRQ(ierr);
94820cf1dd8SToby Isaac   ierr = PetscDualSpaceGetPointSubspace(dsp,refPoint,&dsubsp);CHKERRQ(ierr);
94920cf1dd8SToby Isaac   ierr = PetscSpaceCreateSubspace(bsp,dsubsp,v,J,NULL,NULL,PETSC_OWN_POINTER,&bsubsp);CHKERRQ(ierr);
95020cf1dd8SToby Isaac   ierr = PetscSpaceSetUp(bsubsp);CHKERRQ(ierr);
95120cf1dd8SToby Isaac   ierr = PetscFECreate(PetscObjectComm((PetscObject)fe),trFE);CHKERRQ(ierr);
95220cf1dd8SToby Isaac   ierr = PetscFEGetType(fe,&type);CHKERRQ(ierr);
95320cf1dd8SToby Isaac   ierr = PetscFESetType(*trFE,type);CHKERRQ(ierr);
95420cf1dd8SToby Isaac   ierr = PetscFEGetNumComponents(fe,&numComp);CHKERRQ(ierr);
95520cf1dd8SToby Isaac   ierr = PetscFESetNumComponents(*trFE,numComp);CHKERRQ(ierr);
95620cf1dd8SToby Isaac   ierr = PetscFESetBasisSpace(*trFE,bsubsp);CHKERRQ(ierr);
95720cf1dd8SToby Isaac   ierr = PetscFESetDualSpace(*trFE,dsubsp);CHKERRQ(ierr);
958db11e2ebSMatthew G. Knepley   ierr = PetscObjectGetName((PetscObject) fe, &name);CHKERRQ(ierr);
959db11e2ebSMatthew G. Knepley   if (name) {ierr = PetscFESetName(*trFE, name);CHKERRQ(ierr);}
96020cf1dd8SToby Isaac   ierr = PetscFEGetQuadrature(fe,&fullQuad);CHKERRQ(ierr);
96120cf1dd8SToby Isaac   ierr = PetscQuadratureGetOrder(fullQuad,&order);CHKERRQ(ierr);
96220cf1dd8SToby Isaac   ierr = DMPlexGetConeSize(dm,refPoint,&coneSize);CHKERRQ(ierr);
96320cf1dd8SToby Isaac   if (coneSize == 2 * depth) {
96420cf1dd8SToby Isaac     ierr = PetscDTGaussTensorQuadrature(depth,1,(order + 1)/2,-1.,1.,&subQuad);CHKERRQ(ierr);
96520cf1dd8SToby Isaac   } else {
96620cf1dd8SToby Isaac     ierr = PetscDTGaussJacobiQuadrature(depth,1,(order + 1)/2,-1.,1.,&subQuad);CHKERRQ(ierr);
96720cf1dd8SToby Isaac   }
96820cf1dd8SToby Isaac   ierr = PetscFESetQuadrature(*trFE,subQuad);CHKERRQ(ierr);
96920cf1dd8SToby Isaac   ierr = PetscFESetUp(*trFE);CHKERRQ(ierr);
97020cf1dd8SToby Isaac   ierr = PetscQuadratureDestroy(&subQuad);CHKERRQ(ierr);
97120cf1dd8SToby Isaac   ierr = PetscSpaceDestroy(&bsubsp);CHKERRQ(ierr);
97220cf1dd8SToby Isaac   PetscFunctionReturn(0);
97320cf1dd8SToby Isaac }
97420cf1dd8SToby Isaac 
97520cf1dd8SToby Isaac PetscErrorCode PetscFECreateHeightTrace(PetscFE fe, PetscInt height, PetscFE *trFE)
97620cf1dd8SToby Isaac {
97720cf1dd8SToby Isaac   PetscInt       hStart, hEnd;
97820cf1dd8SToby Isaac   PetscDualSpace dsp;
97920cf1dd8SToby Isaac   DM             dm;
98020cf1dd8SToby Isaac   PetscErrorCode ierr;
98120cf1dd8SToby Isaac 
98220cf1dd8SToby Isaac   PetscFunctionBegin;
98320cf1dd8SToby Isaac   PetscValidHeaderSpecific(fe,PETSCFE_CLASSID,1);
98420cf1dd8SToby Isaac   PetscValidPointer(trFE,3);
98520cf1dd8SToby Isaac   *trFE = NULL;
98620cf1dd8SToby Isaac   ierr = PetscFEGetDualSpace(fe,&dsp);CHKERRQ(ierr);
98720cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDM(dsp,&dm);CHKERRQ(ierr);
98820cf1dd8SToby Isaac   ierr = DMPlexGetHeightStratum(dm,height,&hStart,&hEnd);CHKERRQ(ierr);
98920cf1dd8SToby Isaac   if (hEnd <= hStart) PetscFunctionReturn(0);
99020cf1dd8SToby Isaac   ierr = PetscFECreatePointTrace(fe,hStart,trFE);CHKERRQ(ierr);
99120cf1dd8SToby Isaac   PetscFunctionReturn(0);
99220cf1dd8SToby Isaac }
99320cf1dd8SToby Isaac 
99420cf1dd8SToby Isaac 
99520cf1dd8SToby Isaac /*@
99620cf1dd8SToby Isaac   PetscFEGetDimension - Get the dimension of the finite element space on a cell
99720cf1dd8SToby Isaac 
99820cf1dd8SToby Isaac   Not collective
99920cf1dd8SToby Isaac 
100020cf1dd8SToby Isaac   Input Parameter:
100120cf1dd8SToby Isaac . fe - The PetscFE
100220cf1dd8SToby Isaac 
100320cf1dd8SToby Isaac   Output Parameter:
100420cf1dd8SToby Isaac . dim - The dimension
100520cf1dd8SToby Isaac 
100620cf1dd8SToby Isaac   Level: intermediate
100720cf1dd8SToby Isaac 
100820cf1dd8SToby Isaac .seealso: PetscFECreate(), PetscSpaceGetDimension(), PetscDualSpaceGetDimension()
100920cf1dd8SToby Isaac @*/
101020cf1dd8SToby Isaac PetscErrorCode PetscFEGetDimension(PetscFE fem, PetscInt *dim)
101120cf1dd8SToby Isaac {
101220cf1dd8SToby Isaac   PetscErrorCode ierr;
101320cf1dd8SToby Isaac 
101420cf1dd8SToby Isaac   PetscFunctionBegin;
101520cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
101620cf1dd8SToby Isaac   PetscValidPointer(dim, 2);
101720cf1dd8SToby Isaac   if (fem->ops->getdimension) {ierr = (*fem->ops->getdimension)(fem, dim);CHKERRQ(ierr);}
101820cf1dd8SToby Isaac   PetscFunctionReturn(0);
101920cf1dd8SToby Isaac }
102020cf1dd8SToby Isaac 
10214bee2e38SMatthew G. Knepley /*@C
10224bee2e38SMatthew G. Knepley   PetscFEPushforward - Map the reference element function to real space
10234bee2e38SMatthew G. Knepley 
10244bee2e38SMatthew G. Knepley   Input Parameters:
10254bee2e38SMatthew G. Knepley + fe     - The PetscFE
10264bee2e38SMatthew G. Knepley . fegeom - The cell geometry
10274bee2e38SMatthew G. Knepley . Nv     - The number of function values
10284bee2e38SMatthew G. Knepley - vals   - The function values
10294bee2e38SMatthew G. Knepley 
10304bee2e38SMatthew G. Knepley   Output Parameter:
10314bee2e38SMatthew G. Knepley . vals   - The transformed function values
10324bee2e38SMatthew G. Knepley 
10334bee2e38SMatthew G. Knepley   Level: advanced
10344bee2e38SMatthew G. Knepley 
10354bee2e38SMatthew G. Knepley   Note: This just forwards the call onto PetscDualSpacePushforward().
10364bee2e38SMatthew G. Knepley 
10374bee2e38SMatthew G. Knepley .seealso: PetscDualSpacePushforward()
10384bee2e38SMatthew G. Knepley @*/
10394bee2e38SMatthew G. Knepley PetscErrorCode PetscFEPushforward(PetscFE fe, PetscFEGeom *fegeom, PetscInt Nv, PetscScalar vals[])
10404bee2e38SMatthew G. Knepley {
10414bee2e38SMatthew G. Knepley   PetscErrorCode ierr;
10424bee2e38SMatthew G. Knepley 
10432ae266adSMatthew G. Knepley   PetscFunctionBeginHot;
10442ae266adSMatthew G. Knepley   ierr = PetscDualSpacePushforward(fe->dualSpace, fegeom, Nv, fe->numComponents, vals);CHKERRQ(ierr);
10454bee2e38SMatthew G. Knepley   PetscFunctionReturn(0);
10464bee2e38SMatthew G. Knepley }
10474bee2e38SMatthew G. Knepley 
10484bee2e38SMatthew G. Knepley /*@C
10494bee2e38SMatthew G. Knepley   PetscFEPushforwardGradient - Map the reference element function gradient to real space
10504bee2e38SMatthew G. Knepley 
10514bee2e38SMatthew G. Knepley   Input Parameters:
10524bee2e38SMatthew G. Knepley + fe     - The PetscFE
10534bee2e38SMatthew G. Knepley . fegeom - The cell geometry
10544bee2e38SMatthew G. Knepley . Nv     - The number of function gradient values
10554bee2e38SMatthew G. Knepley - vals   - The function gradient values
10564bee2e38SMatthew G. Knepley 
10574bee2e38SMatthew G. Knepley   Output Parameter:
10584bee2e38SMatthew G. Knepley . vals   - The transformed function gradient values
10594bee2e38SMatthew G. Knepley 
10604bee2e38SMatthew G. Knepley   Level: advanced
10614bee2e38SMatthew G. Knepley 
10624bee2e38SMatthew G. Knepley   Note: This just forwards the call onto PetscDualSpacePushforwardGradient().
10634bee2e38SMatthew G. Knepley 
10644bee2e38SMatthew G. Knepley .seealso: PetscFEPushforward(), PetscDualSpacePushforwardGradient(), PetscDualSpacePushforward()
10654bee2e38SMatthew G. Knepley @*/
10664bee2e38SMatthew G. Knepley PetscErrorCode PetscFEPushforwardGradient(PetscFE fe, PetscFEGeom *fegeom, PetscInt Nv, PetscScalar vals[])
10674bee2e38SMatthew G. Knepley {
10684bee2e38SMatthew G. Knepley   PetscErrorCode ierr;
10694bee2e38SMatthew G. Knepley 
10702ae266adSMatthew G. Knepley   PetscFunctionBeginHot;
10712ae266adSMatthew G. Knepley   ierr = PetscDualSpacePushforwardGradient(fe->dualSpace, fegeom, Nv, fe->numComponents, vals);CHKERRQ(ierr);
10724bee2e38SMatthew G. Knepley   PetscFunctionReturn(0);
10734bee2e38SMatthew G. Knepley }
10744bee2e38SMatthew G. Knepley 
107520cf1dd8SToby Isaac /*
107620cf1dd8SToby Isaac Purpose: Compute element vector for chunk of elements
107720cf1dd8SToby Isaac 
107820cf1dd8SToby Isaac Input:
107920cf1dd8SToby Isaac   Sizes:
108020cf1dd8SToby Isaac      Ne:  number of elements
108120cf1dd8SToby Isaac      Nf:  number of fields
108220cf1dd8SToby Isaac      PetscFE
108320cf1dd8SToby Isaac        dim: spatial dimension
108420cf1dd8SToby Isaac        Nb:  number of basis functions
108520cf1dd8SToby Isaac        Nc:  number of field components
108620cf1dd8SToby Isaac        PetscQuadrature
108720cf1dd8SToby Isaac          Nq:  number of quadrature points
108820cf1dd8SToby Isaac 
108920cf1dd8SToby Isaac   Geometry:
109020cf1dd8SToby Isaac      PetscFEGeom[Ne] possibly *Nq
109120cf1dd8SToby Isaac        PetscReal v0s[dim]
109220cf1dd8SToby Isaac        PetscReal n[dim]
109320cf1dd8SToby Isaac        PetscReal jacobians[dim*dim]
109420cf1dd8SToby Isaac        PetscReal jacobianInverses[dim*dim]
109520cf1dd8SToby Isaac        PetscReal jacobianDeterminants
109620cf1dd8SToby Isaac   FEM:
109720cf1dd8SToby Isaac      PetscFE
109820cf1dd8SToby Isaac        PetscQuadrature
109920cf1dd8SToby Isaac          PetscReal   quadPoints[Nq*dim]
110020cf1dd8SToby Isaac          PetscReal   quadWeights[Nq]
110120cf1dd8SToby Isaac        PetscReal   basis[Nq*Nb*Nc]
110220cf1dd8SToby Isaac        PetscReal   basisDer[Nq*Nb*Nc*dim]
110320cf1dd8SToby Isaac      PetscScalar coefficients[Ne*Nb*Nc]
110420cf1dd8SToby Isaac      PetscScalar elemVec[Ne*Nb*Nc]
110520cf1dd8SToby Isaac 
110620cf1dd8SToby Isaac   Problem:
110720cf1dd8SToby Isaac      PetscInt f: the active field
110820cf1dd8SToby Isaac      f0, f1
110920cf1dd8SToby Isaac 
111020cf1dd8SToby Isaac   Work Space:
111120cf1dd8SToby Isaac      PetscFE
111220cf1dd8SToby Isaac        PetscScalar f0[Nq*dim];
111320cf1dd8SToby Isaac        PetscScalar f1[Nq*dim*dim];
111420cf1dd8SToby Isaac        PetscScalar u[Nc];
111520cf1dd8SToby Isaac        PetscScalar gradU[Nc*dim];
111620cf1dd8SToby Isaac        PetscReal   x[dim];
111720cf1dd8SToby Isaac        PetscScalar realSpaceDer[dim];
111820cf1dd8SToby Isaac 
111920cf1dd8SToby Isaac Purpose: Compute element vector for N_cb batches of elements
112020cf1dd8SToby Isaac 
112120cf1dd8SToby Isaac Input:
112220cf1dd8SToby Isaac   Sizes:
112320cf1dd8SToby Isaac      N_cb: Number of serial cell batches
112420cf1dd8SToby Isaac 
112520cf1dd8SToby Isaac   Geometry:
112620cf1dd8SToby Isaac      PetscReal v0s[Ne*dim]
112720cf1dd8SToby Isaac      PetscReal jacobians[Ne*dim*dim]        possibly *Nq
112820cf1dd8SToby Isaac      PetscReal jacobianInverses[Ne*dim*dim] possibly *Nq
112920cf1dd8SToby Isaac      PetscReal jacobianDeterminants[Ne]     possibly *Nq
113020cf1dd8SToby Isaac   FEM:
113120cf1dd8SToby Isaac      static PetscReal   quadPoints[Nq*dim]
113220cf1dd8SToby Isaac      static PetscReal   quadWeights[Nq]
113320cf1dd8SToby Isaac      static PetscReal   basis[Nq*Nb*Nc]
113420cf1dd8SToby Isaac      static PetscReal   basisDer[Nq*Nb*Nc*dim]
113520cf1dd8SToby Isaac      PetscScalar coefficients[Ne*Nb*Nc]
113620cf1dd8SToby Isaac      PetscScalar elemVec[Ne*Nb*Nc]
113720cf1dd8SToby Isaac 
113820cf1dd8SToby Isaac ex62.c:
113920cf1dd8SToby Isaac   PetscErrorCode PetscFEIntegrateResidualBatch(PetscInt Ne, PetscInt numFields, PetscInt field, PetscQuadrature quad[], const PetscScalar coefficients[],
114020cf1dd8SToby Isaac                                                const PetscReal v0s[], const PetscReal jacobians[], const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[],
114120cf1dd8SToby Isaac                                                void (*f0_func)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar f0[]),
114220cf1dd8SToby Isaac                                                void (*f1_func)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar f1[]), PetscScalar elemVec[])
114320cf1dd8SToby Isaac 
114420cf1dd8SToby Isaac ex52.c:
114520cf1dd8SToby 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)
114620cf1dd8SToby 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)
114720cf1dd8SToby Isaac 
114820cf1dd8SToby Isaac ex52_integrateElement.cu
114920cf1dd8SToby Isaac __global__ void integrateElementQuadrature(int N_cb, realType *coefficients, realType *jacobianInverses, realType *jacobianDeterminants, realType *elemVec)
115020cf1dd8SToby Isaac 
115120cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode IntegrateElementBatchGPU(PetscInt spatial_dim, PetscInt Ne, PetscInt Ncb, PetscInt Nbc, PetscInt Nbl, const PetscScalar coefficients[],
115220cf1dd8SToby Isaac                                                      const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], PetscScalar elemVec[],
115320cf1dd8SToby Isaac                                                      PetscLogEvent event, PetscInt debug, PetscInt pde_op)
115420cf1dd8SToby Isaac 
115520cf1dd8SToby Isaac ex52_integrateElementOpenCL.c:
115620cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode IntegrateElementBatchGPU(PetscInt spatial_dim, PetscInt Ne, PetscInt Ncb, PetscInt Nbc, PetscInt N_bl, const PetscScalar coefficients[],
115720cf1dd8SToby Isaac                                                      const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], PetscScalar elemVec[],
115820cf1dd8SToby Isaac                                                      PetscLogEvent event, PetscInt debug, PetscInt pde_op)
115920cf1dd8SToby Isaac 
116020cf1dd8SToby Isaac __kernel void integrateElementQuadrature(int N_cb, __global float *coefficients, __global float *jacobianInverses, __global float *jacobianDeterminants, __global float *elemVec)
116120cf1dd8SToby Isaac */
116220cf1dd8SToby Isaac 
116320cf1dd8SToby Isaac /*@C
116420cf1dd8SToby Isaac   PetscFEIntegrate - Produce the integral for the given field for a chunk of elements by quadrature integration
116520cf1dd8SToby Isaac 
116620cf1dd8SToby Isaac   Not collective
116720cf1dd8SToby Isaac 
116820cf1dd8SToby Isaac   Input Parameters:
116920cf1dd8SToby Isaac + fem          - The PetscFE object for the field being integrated
117020cf1dd8SToby Isaac . prob         - The PetscDS specifying the discretizations and continuum functions
117120cf1dd8SToby Isaac . field        - The field being integrated
117220cf1dd8SToby Isaac . Ne           - The number of elements in the chunk
117320cf1dd8SToby Isaac . cgeom        - The cell geometry for each cell in the chunk
117420cf1dd8SToby Isaac . coefficients - The array of FEM basis coefficients for the elements
117520cf1dd8SToby Isaac . probAux      - The PetscDS specifying the auxiliary discretizations
117620cf1dd8SToby Isaac - coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
117720cf1dd8SToby Isaac 
117820cf1dd8SToby Isaac   Output Parameter
117920cf1dd8SToby Isaac . integral     - the integral for this field
118020cf1dd8SToby Isaac 
118120cf1dd8SToby Isaac   Level: developer
118220cf1dd8SToby Isaac 
118320cf1dd8SToby Isaac .seealso: PetscFEIntegrateResidual()
118420cf1dd8SToby Isaac @*/
11854bee2e38SMatthew G. Knepley PetscErrorCode PetscFEIntegrate(PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom,
118620cf1dd8SToby Isaac                                 const PetscScalar coefficients[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscScalar integral[])
118720cf1dd8SToby Isaac {
11884bee2e38SMatthew G. Knepley   PetscFE        fe;
118920cf1dd8SToby Isaac   PetscErrorCode ierr;
119020cf1dd8SToby Isaac 
119120cf1dd8SToby Isaac   PetscFunctionBegin;
11924bee2e38SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
11934bee2e38SMatthew G. Knepley   ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);CHKERRQ(ierr);
11944bee2e38SMatthew G. Knepley   if (fe->ops->integrate) {ierr = (*fe->ops->integrate)(prob, field, Ne, cgeom, coefficients, probAux, coefficientsAux, integral);CHKERRQ(ierr);}
119520cf1dd8SToby Isaac   PetscFunctionReturn(0);
119620cf1dd8SToby Isaac }
119720cf1dd8SToby Isaac 
119820cf1dd8SToby Isaac /*@C
1199afe6d6adSToby Isaac   PetscFEIntegrateBd - Produce the integral for the given field for a chunk of elements by quadrature integration
1200afe6d6adSToby Isaac 
1201afe6d6adSToby Isaac   Not collective
1202afe6d6adSToby Isaac 
1203afe6d6adSToby Isaac   Input Parameters:
1204afe6d6adSToby Isaac + fem          - The PetscFE object for the field being integrated
1205afe6d6adSToby Isaac . prob         - The PetscDS specifying the discretizations and continuum functions
1206afe6d6adSToby Isaac . field        - The field being integrated
1207afe6d6adSToby Isaac . obj_func     - The function to be integrated
1208afe6d6adSToby Isaac . Ne           - The number of elements in the chunk
1209afe6d6adSToby Isaac . fgeom        - The face geometry for each face in the chunk
1210afe6d6adSToby Isaac . coefficients - The array of FEM basis coefficients for the elements
1211afe6d6adSToby Isaac . probAux      - The PetscDS specifying the auxiliary discretizations
1212afe6d6adSToby Isaac - coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1213afe6d6adSToby Isaac 
1214afe6d6adSToby Isaac   Output Parameter
1215afe6d6adSToby Isaac . integral     - the integral for this field
1216afe6d6adSToby Isaac 
1217afe6d6adSToby Isaac   Level: developer
1218afe6d6adSToby Isaac 
1219afe6d6adSToby Isaac .seealso: PetscFEIntegrateResidual()
1220afe6d6adSToby Isaac @*/
12214bee2e38SMatthew G. Knepley PetscErrorCode PetscFEIntegrateBd(PetscDS prob, PetscInt field,
1222afe6d6adSToby Isaac                                   void (*obj_func)(PetscInt, PetscInt, PetscInt,
1223afe6d6adSToby Isaac                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1224afe6d6adSToby Isaac                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1225afe6d6adSToby Isaac                                                    PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1226afe6d6adSToby Isaac                                   PetscInt Ne, PetscFEGeom *geom, const PetscScalar coefficients[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscScalar integral[])
1227afe6d6adSToby Isaac {
12284bee2e38SMatthew G. Knepley   PetscFE        fe;
1229afe6d6adSToby Isaac   PetscErrorCode ierr;
1230afe6d6adSToby Isaac 
1231afe6d6adSToby Isaac   PetscFunctionBegin;
12324bee2e38SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
12334bee2e38SMatthew G. Knepley   ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);CHKERRQ(ierr);
12344bee2e38SMatthew G. Knepley   if (fe->ops->integratebd) {ierr = (*fe->ops->integratebd)(prob, field, obj_func, Ne, geom, coefficients, probAux, coefficientsAux, integral);CHKERRQ(ierr);}
1235afe6d6adSToby Isaac   PetscFunctionReturn(0);
1236afe6d6adSToby Isaac }
1237afe6d6adSToby Isaac 
1238afe6d6adSToby Isaac /*@C
123920cf1dd8SToby Isaac   PetscFEIntegrateResidual - Produce the element residual vector for a chunk of elements by quadrature integration
124020cf1dd8SToby Isaac 
124120cf1dd8SToby Isaac   Not collective
124220cf1dd8SToby Isaac 
124320cf1dd8SToby Isaac   Input Parameters:
124420cf1dd8SToby Isaac + fem          - The PetscFE object for the field being integrated
124520cf1dd8SToby Isaac . prob         - The PetscDS specifying the discretizations and continuum functions
124620cf1dd8SToby Isaac . field        - The field being integrated
124720cf1dd8SToby Isaac . Ne           - The number of elements in the chunk
124820cf1dd8SToby Isaac . cgeom        - The cell geometry for each cell in the chunk
124920cf1dd8SToby Isaac . coefficients - The array of FEM basis coefficients for the elements
125020cf1dd8SToby Isaac . coefficients_t - The array of FEM basis time derivative coefficients for the elements
125120cf1dd8SToby Isaac . probAux      - The PetscDS specifying the auxiliary discretizations
125220cf1dd8SToby Isaac . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
125320cf1dd8SToby Isaac - t            - The time
125420cf1dd8SToby Isaac 
125520cf1dd8SToby Isaac   Output Parameter
125620cf1dd8SToby Isaac . elemVec      - the element residual vectors from each element
125720cf1dd8SToby Isaac 
125820cf1dd8SToby Isaac   Note:
125920cf1dd8SToby Isaac $ Loop over batch of elements (e):
126020cf1dd8SToby Isaac $   Loop over quadrature points (q):
126120cf1dd8SToby Isaac $     Make u_q and gradU_q (loops over fields,Nb,Ncomp) and x_q
126220cf1dd8SToby Isaac $     Call f_0 and f_1
126320cf1dd8SToby Isaac $   Loop over element vector entries (f,fc --> i):
126420cf1dd8SToby Isaac $     elemVec[i] += \psi^{fc}_f(q) f0_{fc}(u, \nabla u) + \nabla\psi^{fc}_f(q) \cdot f1_{fc,df}(u, \nabla u)
126520cf1dd8SToby Isaac 
126620cf1dd8SToby Isaac   Level: developer
126720cf1dd8SToby Isaac 
126820cf1dd8SToby Isaac .seealso: PetscFEIntegrateResidual()
126920cf1dd8SToby Isaac @*/
12704bee2e38SMatthew G. Knepley PetscErrorCode PetscFEIntegrateResidual(PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom,
127120cf1dd8SToby Isaac                                         const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
127220cf1dd8SToby Isaac {
12734bee2e38SMatthew G. Knepley   PetscFE        fe;
127420cf1dd8SToby Isaac   PetscErrorCode ierr;
127520cf1dd8SToby Isaac 
127620cf1dd8SToby Isaac   PetscFunctionBegin;
12774bee2e38SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
12784bee2e38SMatthew G. Knepley   ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);CHKERRQ(ierr);
12794bee2e38SMatthew G. Knepley   if (fe->ops->integrateresidual) {ierr = (*fe->ops->integrateresidual)(prob, field, Ne, cgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec);CHKERRQ(ierr);}
128020cf1dd8SToby Isaac   PetscFunctionReturn(0);
128120cf1dd8SToby Isaac }
128220cf1dd8SToby Isaac 
128320cf1dd8SToby Isaac /*@C
128420cf1dd8SToby Isaac   PetscFEIntegrateBdResidual - Produce the element residual vector for a chunk of elements by quadrature integration over a boundary
128520cf1dd8SToby Isaac 
128620cf1dd8SToby Isaac   Not collective
128720cf1dd8SToby Isaac 
128820cf1dd8SToby Isaac   Input Parameters:
128920cf1dd8SToby Isaac + fem          - The PetscFE object for the field being integrated
129020cf1dd8SToby Isaac . prob         - The PetscDS specifying the discretizations and continuum functions
129120cf1dd8SToby Isaac . field        - The field being integrated
129220cf1dd8SToby Isaac . Ne           - The number of elements in the chunk
129320cf1dd8SToby Isaac . fgeom        - The face geometry for each cell in the chunk
129420cf1dd8SToby Isaac . coefficients - The array of FEM basis coefficients for the elements
129520cf1dd8SToby Isaac . coefficients_t - The array of FEM basis time derivative coefficients for the elements
129620cf1dd8SToby Isaac . probAux      - The PetscDS specifying the auxiliary discretizations
129720cf1dd8SToby Isaac . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
129820cf1dd8SToby Isaac - t            - The time
129920cf1dd8SToby Isaac 
130020cf1dd8SToby Isaac   Output Parameter
130120cf1dd8SToby Isaac . elemVec      - the element residual vectors from each element
130220cf1dd8SToby Isaac 
130320cf1dd8SToby Isaac   Level: developer
130420cf1dd8SToby Isaac 
130520cf1dd8SToby Isaac .seealso: PetscFEIntegrateResidual()
130620cf1dd8SToby Isaac @*/
13074bee2e38SMatthew G. Knepley PetscErrorCode PetscFEIntegrateBdResidual(PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *fgeom,
130820cf1dd8SToby Isaac                                           const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
130920cf1dd8SToby Isaac {
13104bee2e38SMatthew G. Knepley   PetscFE        fe;
131120cf1dd8SToby Isaac   PetscErrorCode ierr;
131220cf1dd8SToby Isaac 
131320cf1dd8SToby Isaac   PetscFunctionBegin;
13144bee2e38SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
13154bee2e38SMatthew G. Knepley   ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);CHKERRQ(ierr);
13164bee2e38SMatthew G. Knepley   if (fe->ops->integratebdresidual) {ierr = (*fe->ops->integratebdresidual)(prob, field, Ne, fgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec);CHKERRQ(ierr);}
131720cf1dd8SToby Isaac   PetscFunctionReturn(0);
131820cf1dd8SToby Isaac }
131920cf1dd8SToby Isaac 
132020cf1dd8SToby Isaac /*@C
132120cf1dd8SToby Isaac   PetscFEIntegrateJacobian - Produce the element Jacobian for a chunk of elements by quadrature integration
132220cf1dd8SToby Isaac 
132320cf1dd8SToby Isaac   Not collective
132420cf1dd8SToby Isaac 
132520cf1dd8SToby Isaac   Input Parameters:
132620cf1dd8SToby Isaac + fem          - The PetscFE object for the field being integrated
132720cf1dd8SToby Isaac . prob         - The PetscDS specifying the discretizations and continuum functions
132820cf1dd8SToby Isaac . jtype        - The type of matrix pointwise functions that should be used
132920cf1dd8SToby Isaac . fieldI       - The test field being integrated
133020cf1dd8SToby Isaac . fieldJ       - The basis field being integrated
133120cf1dd8SToby Isaac . Ne           - The number of elements in the chunk
133220cf1dd8SToby Isaac . cgeom        - The cell geometry for each cell in the chunk
133320cf1dd8SToby Isaac . coefficients - The array of FEM basis coefficients for the elements for the Jacobian evaluation point
133420cf1dd8SToby Isaac . coefficients_t - The array of FEM basis time derivative coefficients for the elements
133520cf1dd8SToby Isaac . probAux      - The PetscDS specifying the auxiliary discretizations
133620cf1dd8SToby Isaac . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
133720cf1dd8SToby Isaac . t            - The time
133820cf1dd8SToby Isaac - u_tShift     - A multiplier for the dF/du_t term (as opposed to the dF/du term)
133920cf1dd8SToby Isaac 
134020cf1dd8SToby Isaac   Output Parameter
134120cf1dd8SToby Isaac . elemMat      - the element matrices for the Jacobian from each element
134220cf1dd8SToby Isaac 
134320cf1dd8SToby Isaac   Note:
134420cf1dd8SToby Isaac $ Loop over batch of elements (e):
134520cf1dd8SToby Isaac $   Loop over element matrix entries (f,fc,g,gc --> i,j):
134620cf1dd8SToby Isaac $     Loop over quadrature points (q):
134720cf1dd8SToby Isaac $       Make u_q and gradU_q (loops over fields,Nb,Ncomp)
134820cf1dd8SToby Isaac $         elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q)
134920cf1dd8SToby Isaac $                      + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
135020cf1dd8SToby Isaac $                      + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q)
135120cf1dd8SToby Isaac $                      + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
135220cf1dd8SToby Isaac   Level: developer
135320cf1dd8SToby Isaac 
135420cf1dd8SToby Isaac .seealso: PetscFEIntegrateResidual()
135520cf1dd8SToby Isaac @*/
13564bee2e38SMatthew G. Knepley PetscErrorCode PetscFEIntegrateJacobian(PetscDS prob, PetscFEJacobianType jtype, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *cgeom,
135720cf1dd8SToby Isaac                                         const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
135820cf1dd8SToby Isaac {
13594bee2e38SMatthew G. Knepley   PetscFE        fe;
136020cf1dd8SToby Isaac   PetscErrorCode ierr;
136120cf1dd8SToby Isaac 
136220cf1dd8SToby Isaac   PetscFunctionBegin;
13634bee2e38SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
13644bee2e38SMatthew G. Knepley   ierr = PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);CHKERRQ(ierr);
13654bee2e38SMatthew G. Knepley   if (fe->ops->integratejacobian) {ierr = (*fe->ops->integratejacobian)(prob, jtype, fieldI, fieldJ, Ne, cgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, u_tshift, elemMat);CHKERRQ(ierr);}
136620cf1dd8SToby Isaac   PetscFunctionReturn(0);
136720cf1dd8SToby Isaac }
136820cf1dd8SToby Isaac 
136920cf1dd8SToby Isaac /*@C
137020cf1dd8SToby Isaac   PetscFEIntegrateBdJacobian - Produce the boundary element Jacobian for a chunk of elements by quadrature integration
137120cf1dd8SToby Isaac 
137220cf1dd8SToby Isaac   Not collective
137320cf1dd8SToby Isaac 
137420cf1dd8SToby Isaac   Input Parameters:
137520cf1dd8SToby Isaac . prob         - The PetscDS specifying the discretizations and continuum functions
137620cf1dd8SToby Isaac . fieldI       - The test field being integrated
137720cf1dd8SToby Isaac . fieldJ       - The basis field being integrated
137820cf1dd8SToby Isaac . Ne           - The number of elements in the chunk
137920cf1dd8SToby Isaac . fgeom        - The face geometry for each cell in the chunk
138020cf1dd8SToby Isaac . coefficients - The array of FEM basis coefficients for the elements for the Jacobian evaluation point
138120cf1dd8SToby Isaac . coefficients_t - The array of FEM basis time derivative coefficients for the elements
138220cf1dd8SToby Isaac . probAux      - The PetscDS specifying the auxiliary discretizations
138320cf1dd8SToby Isaac . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
138420cf1dd8SToby Isaac . t            - The time
138520cf1dd8SToby Isaac - u_tShift     - A multiplier for the dF/du_t term (as opposed to the dF/du term)
138620cf1dd8SToby Isaac 
138720cf1dd8SToby Isaac   Output Parameter
138820cf1dd8SToby Isaac . elemMat              - the element matrices for the Jacobian from each element
138920cf1dd8SToby Isaac 
139020cf1dd8SToby Isaac   Note:
139120cf1dd8SToby Isaac $ Loop over batch of elements (e):
139220cf1dd8SToby Isaac $   Loop over element matrix entries (f,fc,g,gc --> i,j):
139320cf1dd8SToby Isaac $     Loop over quadrature points (q):
139420cf1dd8SToby Isaac $       Make u_q and gradU_q (loops over fields,Nb,Ncomp)
139520cf1dd8SToby Isaac $         elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q)
139620cf1dd8SToby Isaac $                      + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
139720cf1dd8SToby Isaac $                      + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q)
139820cf1dd8SToby Isaac $                      + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
139920cf1dd8SToby Isaac   Level: developer
140020cf1dd8SToby Isaac 
140120cf1dd8SToby Isaac .seealso: PetscFEIntegrateJacobian(), PetscFEIntegrateResidual()
140220cf1dd8SToby Isaac @*/
14034bee2e38SMatthew G. Knepley PetscErrorCode PetscFEIntegrateBdJacobian(PetscDS prob, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *fgeom,
140420cf1dd8SToby Isaac                                           const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
140520cf1dd8SToby Isaac {
14064bee2e38SMatthew G. Knepley   PetscFE        fe;
140720cf1dd8SToby Isaac   PetscErrorCode ierr;
140820cf1dd8SToby Isaac 
140920cf1dd8SToby Isaac   PetscFunctionBegin;
14104bee2e38SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
14114bee2e38SMatthew G. Knepley   ierr = PetscDSGetDiscretization(prob, fieldI, (PetscObject *) &fe);CHKERRQ(ierr);
14124bee2e38SMatthew G. Knepley   if (fe->ops->integratebdjacobian) {ierr = (*fe->ops->integratebdjacobian)(prob, fieldI, fieldJ, Ne, fgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, u_tshift, elemMat);CHKERRQ(ierr);}
141320cf1dd8SToby Isaac   PetscFunctionReturn(0);
141420cf1dd8SToby Isaac }
141520cf1dd8SToby Isaac 
141620cf1dd8SToby Isaac PetscErrorCode PetscFEGetHeightSubspace(PetscFE fe, PetscInt height, PetscFE *subfe)
141720cf1dd8SToby Isaac {
141820cf1dd8SToby Isaac   PetscSpace      P, subP;
141920cf1dd8SToby Isaac   PetscDualSpace  Q, subQ;
142020cf1dd8SToby Isaac   PetscQuadrature subq;
142120cf1dd8SToby Isaac   PetscFEType     fetype;
142220cf1dd8SToby Isaac   PetscInt        dim, Nc;
142320cf1dd8SToby Isaac   PetscErrorCode  ierr;
142420cf1dd8SToby Isaac 
142520cf1dd8SToby Isaac   PetscFunctionBegin;
142620cf1dd8SToby Isaac   PetscValidHeaderSpecific(fe, PETSCFE_CLASSID, 1);
142720cf1dd8SToby Isaac   PetscValidPointer(subfe, 3);
142820cf1dd8SToby Isaac   if (height == 0) {
142920cf1dd8SToby Isaac     *subfe = fe;
143020cf1dd8SToby Isaac     PetscFunctionReturn(0);
143120cf1dd8SToby Isaac   }
143220cf1dd8SToby Isaac   ierr = PetscFEGetBasisSpace(fe, &P);CHKERRQ(ierr);
143320cf1dd8SToby Isaac   ierr = PetscFEGetDualSpace(fe, &Q);CHKERRQ(ierr);
143420cf1dd8SToby Isaac   ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
143520cf1dd8SToby Isaac   ierr = PetscFEGetFaceQuadrature(fe, &subq);CHKERRQ(ierr);
143620cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDimension(Q, &dim);CHKERRQ(ierr);
143720cf1dd8SToby 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);}
143820cf1dd8SToby Isaac   if (!fe->subspaces) {ierr = PetscCalloc1(dim, &fe->subspaces);CHKERRQ(ierr);}
143920cf1dd8SToby Isaac   if (height <= dim) {
144020cf1dd8SToby Isaac     if (!fe->subspaces[height-1]) {
144120cf1dd8SToby Isaac       PetscFE     sub;
14423f6b16c7SMatthew G. Knepley       const char *name;
144320cf1dd8SToby Isaac 
144420cf1dd8SToby Isaac       ierr = PetscSpaceGetHeightSubspace(P, height, &subP);CHKERRQ(ierr);
144520cf1dd8SToby Isaac       ierr = PetscDualSpaceGetHeightSubspace(Q, height, &subQ);CHKERRQ(ierr);
144620cf1dd8SToby Isaac       ierr = PetscFECreate(PetscObjectComm((PetscObject) fe), &sub);CHKERRQ(ierr);
14473f6b16c7SMatthew G. Knepley       ierr = PetscObjectGetName((PetscObject) fe,  &name);CHKERRQ(ierr);
14483f6b16c7SMatthew G. Knepley       ierr = PetscObjectSetName((PetscObject) sub,  name);CHKERRQ(ierr);
144920cf1dd8SToby Isaac       ierr = PetscFEGetType(fe, &fetype);CHKERRQ(ierr);
145020cf1dd8SToby Isaac       ierr = PetscFESetType(sub, fetype);CHKERRQ(ierr);
145120cf1dd8SToby Isaac       ierr = PetscFESetBasisSpace(sub, subP);CHKERRQ(ierr);
145220cf1dd8SToby Isaac       ierr = PetscFESetDualSpace(sub, subQ);CHKERRQ(ierr);
145320cf1dd8SToby Isaac       ierr = PetscFESetNumComponents(sub, Nc);CHKERRQ(ierr);
145420cf1dd8SToby Isaac       ierr = PetscFESetUp(sub);CHKERRQ(ierr);
145520cf1dd8SToby Isaac       ierr = PetscFESetQuadrature(sub, subq);CHKERRQ(ierr);
145620cf1dd8SToby Isaac       fe->subspaces[height-1] = sub;
145720cf1dd8SToby Isaac     }
145820cf1dd8SToby Isaac     *subfe = fe->subspaces[height-1];
145920cf1dd8SToby Isaac   } else {
146020cf1dd8SToby Isaac     *subfe = NULL;
146120cf1dd8SToby Isaac   }
146220cf1dd8SToby Isaac   PetscFunctionReturn(0);
146320cf1dd8SToby Isaac }
146420cf1dd8SToby Isaac 
146520cf1dd8SToby Isaac /*@
146620cf1dd8SToby Isaac   PetscFERefine - Create a "refined" PetscFE object that refines the reference cell into smaller copies. This is typically used
146720cf1dd8SToby Isaac   to precondition a higher order method with a lower order method on a refined mesh having the same number of dofs (but more
146820cf1dd8SToby Isaac   sparsity). It is also used to create an interpolation between regularly refined meshes.
146920cf1dd8SToby Isaac 
1470d083f849SBarry Smith   Collective on fem
147120cf1dd8SToby Isaac 
147220cf1dd8SToby Isaac   Input Parameter:
147320cf1dd8SToby Isaac . fe - The initial PetscFE
147420cf1dd8SToby Isaac 
147520cf1dd8SToby Isaac   Output Parameter:
147620cf1dd8SToby Isaac . feRef - The refined PetscFE
147720cf1dd8SToby Isaac 
147820cf1dd8SToby Isaac   Level: developer
147920cf1dd8SToby Isaac 
148020cf1dd8SToby Isaac .seealso: PetscFEType, PetscFECreate(), PetscFESetType()
148120cf1dd8SToby Isaac @*/
148220cf1dd8SToby Isaac PetscErrorCode PetscFERefine(PetscFE fe, PetscFE *feRef)
148320cf1dd8SToby Isaac {
148420cf1dd8SToby Isaac   PetscSpace       P, Pref;
148520cf1dd8SToby Isaac   PetscDualSpace   Q, Qref;
148620cf1dd8SToby Isaac   DM               K, Kref;
148720cf1dd8SToby Isaac   PetscQuadrature  q, qref;
148820cf1dd8SToby Isaac   const PetscReal *v0, *jac;
148920cf1dd8SToby Isaac   PetscInt         numComp, numSubelements;
149020cf1dd8SToby Isaac   PetscErrorCode   ierr;
149120cf1dd8SToby Isaac 
149220cf1dd8SToby Isaac   PetscFunctionBegin;
149320cf1dd8SToby Isaac   ierr = PetscFEGetBasisSpace(fe, &P);CHKERRQ(ierr);
149420cf1dd8SToby Isaac   ierr = PetscFEGetDualSpace(fe, &Q);CHKERRQ(ierr);
149520cf1dd8SToby Isaac   ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
149620cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDM(Q, &K);CHKERRQ(ierr);
149720cf1dd8SToby Isaac   /* Create space */
149820cf1dd8SToby Isaac   ierr = PetscObjectReference((PetscObject) P);CHKERRQ(ierr);
149920cf1dd8SToby Isaac   Pref = P;
150020cf1dd8SToby Isaac   /* Create dual space */
150120cf1dd8SToby Isaac   ierr = PetscDualSpaceDuplicate(Q, &Qref);CHKERRQ(ierr);
150220cf1dd8SToby Isaac   ierr = DMRefine(K, PetscObjectComm((PetscObject) fe), &Kref);CHKERRQ(ierr);
150320cf1dd8SToby Isaac   ierr = PetscDualSpaceSetDM(Qref, Kref);CHKERRQ(ierr);
150420cf1dd8SToby Isaac   ierr = DMDestroy(&Kref);CHKERRQ(ierr);
150520cf1dd8SToby Isaac   ierr = PetscDualSpaceSetUp(Qref);CHKERRQ(ierr);
150620cf1dd8SToby Isaac   /* Create element */
150720cf1dd8SToby Isaac   ierr = PetscFECreate(PetscObjectComm((PetscObject) fe), feRef);CHKERRQ(ierr);
150820cf1dd8SToby Isaac   ierr = PetscFESetType(*feRef, PETSCFECOMPOSITE);CHKERRQ(ierr);
150920cf1dd8SToby Isaac   ierr = PetscFESetBasisSpace(*feRef, Pref);CHKERRQ(ierr);
151020cf1dd8SToby Isaac   ierr = PetscFESetDualSpace(*feRef, Qref);CHKERRQ(ierr);
151120cf1dd8SToby Isaac   ierr = PetscFEGetNumComponents(fe,    &numComp);CHKERRQ(ierr);
151220cf1dd8SToby Isaac   ierr = PetscFESetNumComponents(*feRef, numComp);CHKERRQ(ierr);
151320cf1dd8SToby Isaac   ierr = PetscFESetUp(*feRef);CHKERRQ(ierr);
151420cf1dd8SToby Isaac   ierr = PetscSpaceDestroy(&Pref);CHKERRQ(ierr);
151520cf1dd8SToby Isaac   ierr = PetscDualSpaceDestroy(&Qref);CHKERRQ(ierr);
151620cf1dd8SToby Isaac   /* Create quadrature */
151720cf1dd8SToby Isaac   ierr = PetscFECompositeGetMapping(*feRef, &numSubelements, &v0, &jac, NULL);CHKERRQ(ierr);
151820cf1dd8SToby Isaac   ierr = PetscQuadratureExpandComposite(q, numSubelements, v0, jac, &qref);CHKERRQ(ierr);
151920cf1dd8SToby Isaac   ierr = PetscFESetQuadrature(*feRef, qref);CHKERRQ(ierr);
152020cf1dd8SToby Isaac   ierr = PetscQuadratureDestroy(&qref);CHKERRQ(ierr);
152120cf1dd8SToby Isaac   PetscFunctionReturn(0);
152220cf1dd8SToby Isaac }
152320cf1dd8SToby Isaac 
152420cf1dd8SToby Isaac /*@C
152520cf1dd8SToby Isaac   PetscFECreateDefault - Create a PetscFE for basic FEM computation
152620cf1dd8SToby Isaac 
1527d083f849SBarry Smith   Collective
152820cf1dd8SToby Isaac 
152920cf1dd8SToby Isaac   Input Parameters:
15307be5e748SToby Isaac + comm      - The MPI comm
153120cf1dd8SToby Isaac . dim       - The spatial dimension
153220cf1dd8SToby Isaac . Nc        - The number of components
153320cf1dd8SToby Isaac . isSimplex - Flag for simplex reference cell, otherwise its a tensor product
153420cf1dd8SToby Isaac . prefix    - The options prefix, or NULL
153520cf1dd8SToby Isaac - qorder    - The quadrature order
153620cf1dd8SToby Isaac 
153720cf1dd8SToby Isaac   Output Parameter:
153820cf1dd8SToby Isaac . fem - The PetscFE object
153920cf1dd8SToby Isaac 
154020cf1dd8SToby Isaac   Level: beginner
154120cf1dd8SToby Isaac 
154220cf1dd8SToby Isaac .seealso: PetscFECreate(), PetscSpaceCreate(), PetscDualSpaceCreate()
154320cf1dd8SToby Isaac @*/
15447be5e748SToby Isaac PetscErrorCode PetscFECreateDefault(MPI_Comm comm, PetscInt dim, PetscInt Nc, PetscBool isSimplex, const char prefix[], PetscInt qorder, PetscFE *fem)
154520cf1dd8SToby Isaac {
154620cf1dd8SToby Isaac   PetscQuadrature q, fq;
154720cf1dd8SToby Isaac   DM              K;
154820cf1dd8SToby Isaac   PetscSpace      P;
154920cf1dd8SToby Isaac   PetscDualSpace  Q;
155020cf1dd8SToby Isaac   PetscInt        order, quadPointsPerEdge;
155120cf1dd8SToby Isaac   PetscBool       tensor = isSimplex ? PETSC_FALSE : PETSC_TRUE;
155220cf1dd8SToby Isaac   PetscErrorCode  ierr;
155320cf1dd8SToby Isaac 
155420cf1dd8SToby Isaac   PetscFunctionBegin;
155520cf1dd8SToby Isaac   /* Create space */
15567be5e748SToby Isaac   ierr = PetscSpaceCreate(comm, &P);CHKERRQ(ierr);
155720cf1dd8SToby Isaac   ierr = PetscObjectSetOptionsPrefix((PetscObject) P, prefix);CHKERRQ(ierr);
155820cf1dd8SToby Isaac   ierr = PetscSpacePolynomialSetTensor(P, tensor);CHKERRQ(ierr);
155920cf1dd8SToby Isaac   ierr = PetscSpaceSetNumComponents(P, Nc);CHKERRQ(ierr);
156020cf1dd8SToby Isaac   ierr = PetscSpaceSetNumVariables(P, dim);CHKERRQ(ierr);
1561028afddaSToby Isaac   ierr = PetscSpaceSetFromOptions(P);CHKERRQ(ierr);
156220cf1dd8SToby Isaac   ierr = PetscSpaceSetUp(P);CHKERRQ(ierr);
156320cf1dd8SToby Isaac   ierr = PetscSpaceGetDegree(P, &order, NULL);CHKERRQ(ierr);
156420cf1dd8SToby Isaac   ierr = PetscSpacePolynomialGetTensor(P, &tensor);CHKERRQ(ierr);
156520cf1dd8SToby Isaac   /* Create dual space */
15667be5e748SToby Isaac   ierr = PetscDualSpaceCreate(comm, &Q);CHKERRQ(ierr);
156720cf1dd8SToby Isaac   ierr = PetscDualSpaceSetType(Q,PETSCDUALSPACELAGRANGE);CHKERRQ(ierr);
156820cf1dd8SToby Isaac   ierr = PetscObjectSetOptionsPrefix((PetscObject) Q, prefix);CHKERRQ(ierr);
156920cf1dd8SToby Isaac   ierr = PetscDualSpaceCreateReferenceCell(Q, dim, isSimplex, &K);CHKERRQ(ierr);
157020cf1dd8SToby Isaac   ierr = PetscDualSpaceSetDM(Q, K);CHKERRQ(ierr);
157120cf1dd8SToby Isaac   ierr = DMDestroy(&K);CHKERRQ(ierr);
157220cf1dd8SToby Isaac   ierr = PetscDualSpaceSetNumComponents(Q, Nc);CHKERRQ(ierr);
157320cf1dd8SToby Isaac   ierr = PetscDualSpaceSetOrder(Q, order);CHKERRQ(ierr);
157420cf1dd8SToby Isaac   ierr = PetscDualSpaceLagrangeSetTensor(Q, tensor);CHKERRQ(ierr);
157520cf1dd8SToby Isaac   ierr = PetscDualSpaceSetFromOptions(Q);CHKERRQ(ierr);
157620cf1dd8SToby Isaac   ierr = PetscDualSpaceSetUp(Q);CHKERRQ(ierr);
157720cf1dd8SToby Isaac   /* Create element */
15787be5e748SToby Isaac   ierr = PetscFECreate(comm, fem);CHKERRQ(ierr);
157920cf1dd8SToby Isaac   ierr = PetscObjectSetOptionsPrefix((PetscObject) *fem, prefix);CHKERRQ(ierr);
158020cf1dd8SToby Isaac   ierr = PetscFESetBasisSpace(*fem, P);CHKERRQ(ierr);
158120cf1dd8SToby Isaac   ierr = PetscFESetDualSpace(*fem, Q);CHKERRQ(ierr);
158220cf1dd8SToby Isaac   ierr = PetscFESetNumComponents(*fem, Nc);CHKERRQ(ierr);
158391e89cf0SMatthew G. Knepley   ierr = PetscFESetFromOptions(*fem);CHKERRQ(ierr);
158420cf1dd8SToby Isaac   ierr = PetscFESetUp(*fem);CHKERRQ(ierr);
158520cf1dd8SToby Isaac   ierr = PetscSpaceDestroy(&P);CHKERRQ(ierr);
158620cf1dd8SToby Isaac   ierr = PetscDualSpaceDestroy(&Q);CHKERRQ(ierr);
158720cf1dd8SToby Isaac   /* Create quadrature (with specified order if given) */
158820cf1dd8SToby Isaac   qorder = qorder >= 0 ? qorder : order;
158920cf1dd8SToby Isaac   ierr = PetscObjectOptionsBegin((PetscObject)*fem);CHKERRQ(ierr);
1590*5a856986SBarry Smith   ierr = PetscOptionsBoundedInt("-petscfe_default_quadrature_order","Quadrature order is one less than quadrature points per edge","PetscFECreateDefault",qorder,&qorder,NULL,0);CHKERRQ(ierr);
159120cf1dd8SToby Isaac   ierr = PetscOptionsEnd();CHKERRQ(ierr);
159220cf1dd8SToby Isaac   quadPointsPerEdge = PetscMax(qorder + 1,1);
159320cf1dd8SToby Isaac   if (isSimplex) {
159420cf1dd8SToby Isaac     ierr = PetscDTGaussJacobiQuadrature(dim,   1, quadPointsPerEdge, -1.0, 1.0, &q);CHKERRQ(ierr);
159520cf1dd8SToby Isaac     ierr = PetscDTGaussJacobiQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0, &fq);CHKERRQ(ierr);
15964ccfa306SStefano Zampini   } else {
159720cf1dd8SToby Isaac     ierr = PetscDTGaussTensorQuadrature(dim,   1, quadPointsPerEdge, -1.0, 1.0, &q);CHKERRQ(ierr);
159820cf1dd8SToby Isaac     ierr = PetscDTGaussTensorQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0, &fq);CHKERRQ(ierr);
159920cf1dd8SToby Isaac   }
160020cf1dd8SToby Isaac   ierr = PetscFESetQuadrature(*fem, q);CHKERRQ(ierr);
160120cf1dd8SToby Isaac   ierr = PetscFESetFaceQuadrature(*fem, fq);CHKERRQ(ierr);
160220cf1dd8SToby Isaac   ierr = PetscQuadratureDestroy(&q);CHKERRQ(ierr);
160320cf1dd8SToby Isaac   ierr = PetscQuadratureDestroy(&fq);CHKERRQ(ierr);
160420cf1dd8SToby Isaac   PetscFunctionReturn(0);
160520cf1dd8SToby Isaac }
16063f6b16c7SMatthew G. Knepley 
16073f6b16c7SMatthew G. Knepley /*@C
16083f6b16c7SMatthew G. Knepley   PetscFESetName - Names the FE and its subobjects
16093f6b16c7SMatthew G. Knepley 
16103f6b16c7SMatthew G. Knepley   Not collective
16113f6b16c7SMatthew G. Knepley 
16123f6b16c7SMatthew G. Knepley   Input Parameters:
16133f6b16c7SMatthew G. Knepley + fe   - The PetscFE
16143f6b16c7SMatthew G. Knepley - name - The name
16153f6b16c7SMatthew G. Knepley 
16163f6b16c7SMatthew G. Knepley   Level: beginner
16173f6b16c7SMatthew G. Knepley 
16183f6b16c7SMatthew G. Knepley .seealso: PetscFECreate(), PetscSpaceCreate(), PetscDualSpaceCreate()
16193f6b16c7SMatthew G. Knepley @*/
16203f6b16c7SMatthew G. Knepley PetscErrorCode PetscFESetName(PetscFE fe, const char name[])
16213f6b16c7SMatthew G. Knepley {
16223f6b16c7SMatthew G. Knepley   PetscSpace     P;
16233f6b16c7SMatthew G. Knepley   PetscDualSpace Q;
16243f6b16c7SMatthew G. Knepley   PetscErrorCode ierr;
16253f6b16c7SMatthew G. Knepley 
16263f6b16c7SMatthew G. Knepley   PetscFunctionBegin;
16273f6b16c7SMatthew G. Knepley   ierr = PetscFEGetBasisSpace(fe, &P);CHKERRQ(ierr);
16283f6b16c7SMatthew G. Knepley   ierr = PetscFEGetDualSpace(fe, &Q);CHKERRQ(ierr);
16293f6b16c7SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) fe, name);CHKERRQ(ierr);
16303f6b16c7SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) P,  name);CHKERRQ(ierr);
16313f6b16c7SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) Q,  name);CHKERRQ(ierr);
16323f6b16c7SMatthew G. Knepley   PetscFunctionReturn(0);
16333f6b16c7SMatthew G. Knepley }
1634a8f1f9e5SMatthew G. Knepley 
1635a8f1f9e5SMatthew G. Knepley PetscErrorCode PetscFEEvaluateFieldJets_Internal(PetscDS ds, PetscInt dim, PetscInt Nf, const PetscInt Nb[], const PetscInt Nc[], PetscInt q, PetscReal *basisField[], PetscReal *basisFieldDer[], PetscFEGeom *fegeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscScalar u[], PetscScalar u_x[], PetscScalar u_t[])
1636a8f1f9e5SMatthew G. Knepley {
1637a8f1f9e5SMatthew G. Knepley   PetscInt       dOffset = 0, fOffset = 0, f;
1638a8f1f9e5SMatthew G. Knepley   PetscErrorCode ierr;
1639a8f1f9e5SMatthew G. Knepley 
1640a8f1f9e5SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
1641a8f1f9e5SMatthew G. Knepley     PetscFE          fe;
1642a8f1f9e5SMatthew G. Knepley     const PetscInt   Nbf = Nb[f], Ncf = Nc[f];
1643a8f1f9e5SMatthew G. Knepley     const PetscReal *Bq = &basisField[f][q*Nbf*Ncf];
1644a8f1f9e5SMatthew G. Knepley     const PetscReal *Dq = &basisFieldDer[f][q*Nbf*Ncf*dim];
1645a8f1f9e5SMatthew G. Knepley     PetscInt         b, c, d;
1646a8f1f9e5SMatthew G. Knepley 
1647a8f1f9e5SMatthew G. Knepley     ierr = PetscDSGetDiscretization(ds, f, (PetscObject *) &fe);CHKERRQ(ierr);
1648a8f1f9e5SMatthew G. Knepley     for (c = 0; c < Ncf; ++c)     u[fOffset+c] = 0.0;
1649a8f1f9e5SMatthew G. Knepley     for (d = 0; d < dim*Ncf; ++d) u_x[fOffset*dim+d] = 0.0;
1650a8f1f9e5SMatthew G. Knepley     for (b = 0; b < Nbf; ++b) {
1651a8f1f9e5SMatthew G. Knepley       for (c = 0; c < Ncf; ++c) {
1652a8f1f9e5SMatthew G. Knepley         const PetscInt cidx = b*Ncf+c;
1653a8f1f9e5SMatthew G. Knepley 
1654a8f1f9e5SMatthew G. Knepley         u[fOffset+c] += Bq[cidx]*coefficients[dOffset+b];
1655a8f1f9e5SMatthew G. Knepley         for (d = 0; d < dim; ++d) u_x[(fOffset+c)*dim+d] += Dq[cidx*dim+d]*coefficients[dOffset+b];
1656a8f1f9e5SMatthew G. Knepley       }
1657a8f1f9e5SMatthew G. Knepley     }
1658a8f1f9e5SMatthew G. Knepley     ierr = PetscFEPushforward(fe, fegeom, 1, &u[fOffset]);CHKERRQ(ierr);
1659a8f1f9e5SMatthew G. Knepley     ierr = PetscFEPushforwardGradient(fe, fegeom, 1, &u_x[fOffset*dim]);CHKERRQ(ierr);
1660a8f1f9e5SMatthew G. Knepley     if (u_t) {
1661a8f1f9e5SMatthew G. Knepley       for (c = 0; c < Ncf; ++c) u_t[fOffset+c] = 0.0;
1662a8f1f9e5SMatthew G. Knepley       for (b = 0; b < Nbf; ++b) {
1663a8f1f9e5SMatthew G. Knepley         for (c = 0; c < Ncf; ++c) {
1664a8f1f9e5SMatthew G. Knepley           const PetscInt cidx = b*Ncf+c;
1665a8f1f9e5SMatthew G. Knepley 
1666a8f1f9e5SMatthew G. Knepley           u_t[fOffset+c] += Bq[cidx]*coefficients_t[dOffset+b];
1667a8f1f9e5SMatthew G. Knepley         }
1668a8f1f9e5SMatthew G. Knepley       }
1669a8f1f9e5SMatthew G. Knepley       ierr = PetscFEPushforward(fe, fegeom, 1, &u_t[fOffset]);CHKERRQ(ierr);
1670a8f1f9e5SMatthew G. Knepley     }
1671a8f1f9e5SMatthew G. Knepley     fOffset += Ncf;
1672a8f1f9e5SMatthew G. Knepley     dOffset += Nbf;
1673a8f1f9e5SMatthew G. Knepley   }
1674a8f1f9e5SMatthew G. Knepley   return 0;
1675a8f1f9e5SMatthew G. Knepley }
1676a8f1f9e5SMatthew G. Knepley 
1677a8f1f9e5SMatthew G. Knepley PetscErrorCode PetscFEEvaluateFaceFields_Internal(PetscDS prob, PetscInt field, PetscInt faceLoc, const PetscScalar coefficients[], PetscScalar u[])
1678a8f1f9e5SMatthew G. Knepley {
1679a8f1f9e5SMatthew G. Knepley   PetscFE        fe;
1680a8f1f9e5SMatthew G. Knepley   PetscReal     *faceBasis;
1681a8f1f9e5SMatthew G. Knepley   PetscInt       Nb, Nc, b, c;
1682a8f1f9e5SMatthew G. Knepley   PetscErrorCode ierr;
1683a8f1f9e5SMatthew G. Knepley 
1684a8f1f9e5SMatthew G. Knepley   if (!prob) return 0;
1685a8f1f9e5SMatthew G. Knepley   ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fe);CHKERRQ(ierr);
1686a8f1f9e5SMatthew G. Knepley   ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
1687a8f1f9e5SMatthew G. Knepley   ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
1688a8f1f9e5SMatthew G. Knepley   ierr = PetscFEGetFaceCentroidTabulation(fe, &faceBasis);CHKERRQ(ierr);
1689a8f1f9e5SMatthew G. Knepley   for (c = 0; c < Nc; ++c) {u[c] = 0.0;}
1690a8f1f9e5SMatthew G. Knepley   for (b = 0; b < Nb; ++b) {
1691a8f1f9e5SMatthew G. Knepley     for (c = 0; c < Nc; ++c) {
1692a8f1f9e5SMatthew G. Knepley       const PetscInt cidx = b*Nc+c;
1693a8f1f9e5SMatthew G. Knepley 
1694a8f1f9e5SMatthew G. Knepley       u[c] += coefficients[cidx]*faceBasis[faceLoc*Nb*Nc+cidx];
1695a8f1f9e5SMatthew G. Knepley     }
1696a8f1f9e5SMatthew G. Knepley   }
1697a8f1f9e5SMatthew G. Knepley   return 0;
1698a8f1f9e5SMatthew G. Knepley }
1699a8f1f9e5SMatthew G. Knepley 
17006142fa51SMatthew G. Knepley PetscErrorCode PetscFEUpdateElementVec_Internal(PetscFE fe, PetscInt dim, PetscInt Nq, PetscInt Nb, PetscInt Nc, PetscReal basis[], PetscReal basisDer[], PetscScalar tmpBasis[], PetscScalar tmpBasisDer[], PetscFEGeom *fegeom, PetscScalar f0[], PetscScalar f1[], PetscScalar elemVec[])
1701a8f1f9e5SMatthew G. Knepley {
1702a8f1f9e5SMatthew G. Knepley   PetscInt       q, b, c, d;
1703a8f1f9e5SMatthew G. Knepley   PetscErrorCode ierr;
1704a8f1f9e5SMatthew G. Knepley 
1705a8f1f9e5SMatthew G. Knepley   for (b = 0; b < Nb; ++b) elemVec[b] = 0.0;
1706a8f1f9e5SMatthew G. Knepley   for (q = 0; q < Nq; ++q) {
1707a8f1f9e5SMatthew G. Knepley     for (b = 0; b < Nb; ++b) {
1708a8f1f9e5SMatthew G. Knepley       for (c = 0; c < Nc; ++c) {
1709a8f1f9e5SMatthew G. Knepley         const PetscInt bcidx = b*Nc+c;
1710a8f1f9e5SMatthew G. Knepley 
1711a8f1f9e5SMatthew G. Knepley         tmpBasis[bcidx] = basis[q*Nb*Nc+bcidx];
1712a8f1f9e5SMatthew G. Knepley         for (d = 0; d < dim; ++d) tmpBasisDer[bcidx*dim+d] = basisDer[q*Nb*Nc*dim+bcidx*dim+d];
1713a8f1f9e5SMatthew G. Knepley       }
1714a8f1f9e5SMatthew G. Knepley     }
1715a8f1f9e5SMatthew G. Knepley     ierr = PetscFEPushforward(fe, fegeom, Nb, tmpBasis);CHKERRQ(ierr);
1716a8f1f9e5SMatthew G. Knepley     ierr = PetscFEPushforwardGradient(fe, fegeom, Nb, tmpBasisDer);CHKERRQ(ierr);
1717a8f1f9e5SMatthew G. Knepley     for (b = 0; b < Nb; ++b) {
1718a8f1f9e5SMatthew G. Knepley       for (c = 0; c < Nc; ++c) {
1719a8f1f9e5SMatthew G. Knepley         const PetscInt bcidx = b*Nc+c;
1720a8f1f9e5SMatthew G. Knepley         const PetscInt qcidx = q*Nc+c;
1721a8f1f9e5SMatthew G. Knepley 
1722a8f1f9e5SMatthew G. Knepley         elemVec[b] += tmpBasis[bcidx]*f0[qcidx];
1723a8f1f9e5SMatthew G. Knepley         for (d = 0; d < dim; ++d) elemVec[b] += tmpBasisDer[bcidx*dim+d]*f1[qcidx*dim+d];
1724a8f1f9e5SMatthew G. Knepley       }
1725a8f1f9e5SMatthew G. Knepley     }
1726a8f1f9e5SMatthew G. Knepley   }
1727a8f1f9e5SMatthew G. Knepley   return(0);
1728a8f1f9e5SMatthew G. Knepley }
1729a8f1f9e5SMatthew G. Knepley 
17306142fa51SMatthew G. Knepley PetscErrorCode PetscFEUpdateElementMat_Internal(PetscFE feI, PetscFE feJ, PetscInt dim, PetscInt NbI, PetscInt NcI, const PetscReal basisI[], const PetscReal basisDerI[], PetscScalar tmpBasisI[], PetscScalar tmpBasisDerI[], PetscInt NbJ, PetscInt NcJ, const PetscReal basisJ[], const PetscReal basisDerJ[], PetscScalar tmpBasisJ[], PetscScalar tmpBasisDerJ[], PetscFEGeom *fegeom, const PetscScalar g0[], const PetscScalar g1[], const PetscScalar g2[], const PetscScalar g3[], PetscInt eOffset, PetscInt totDim, PetscInt offsetI, PetscInt offsetJ, PetscScalar elemMat[])
1731a8f1f9e5SMatthew G. Knepley {
1732a8f1f9e5SMatthew G. Knepley   PetscInt       f, fc, g, gc, df, dg;
1733a8f1f9e5SMatthew G. Knepley   PetscErrorCode ierr;
1734a8f1f9e5SMatthew G. Knepley 
1735a8f1f9e5SMatthew G. Knepley   for (f = 0; f < NbI; ++f) {
1736a8f1f9e5SMatthew G. Knepley     for (fc = 0; fc < NcI; ++fc) {
1737a8f1f9e5SMatthew G. Knepley       const PetscInt fidx = f*NcI+fc; /* Test function basis index */
1738a8f1f9e5SMatthew G. Knepley 
1739a8f1f9e5SMatthew G. Knepley       tmpBasisI[fidx] = basisI[fidx];
1740a8f1f9e5SMatthew G. Knepley       for (df = 0; df < dim; ++df) tmpBasisDerI[fidx*dim+df] = basisDerI[fidx*dim+df];
1741a8f1f9e5SMatthew G. Knepley     }
1742a8f1f9e5SMatthew G. Knepley   }
1743a8f1f9e5SMatthew G. Knepley   ierr = PetscFEPushforward(feI, fegeom, NbI, tmpBasisI);CHKERRQ(ierr);
1744a8f1f9e5SMatthew G. Knepley   ierr = PetscFEPushforwardGradient(feI, fegeom, NbI, tmpBasisDerI);CHKERRQ(ierr);
1745a8f1f9e5SMatthew G. Knepley   for (g = 0; g < NbJ; ++g) {
1746a8f1f9e5SMatthew G. Knepley     for (gc = 0; gc < NcJ; ++gc) {
1747a8f1f9e5SMatthew G. Knepley       const PetscInt gidx = g*NcJ+gc; /* Trial function basis index */
1748a8f1f9e5SMatthew G. Knepley 
1749a8f1f9e5SMatthew G. Knepley       tmpBasisJ[gidx] = basisJ[gidx];
1750a8f1f9e5SMatthew G. Knepley       for (dg = 0; dg < dim; ++dg) tmpBasisDerJ[gidx*dim+dg] = basisDerJ[gidx*dim+dg];
1751a8f1f9e5SMatthew G. Knepley     }
1752a8f1f9e5SMatthew G. Knepley   }
1753a8f1f9e5SMatthew G. Knepley   ierr = PetscFEPushforward(feJ, fegeom, NbJ, tmpBasisJ);CHKERRQ(ierr);
1754a8f1f9e5SMatthew G. Knepley   ierr = PetscFEPushforwardGradient(feJ, fegeom, NbJ, tmpBasisDerJ);CHKERRQ(ierr);
1755a8f1f9e5SMatthew G. Knepley   for (f = 0; f < NbI; ++f) {
1756a8f1f9e5SMatthew G. Knepley     for (fc = 0; fc < NcI; ++fc) {
1757a8f1f9e5SMatthew G. Knepley       const PetscInt fidx = f*NcI+fc; /* Test function basis index */
1758a8f1f9e5SMatthew G. Knepley       const PetscInt i    = offsetI+f; /* Element matrix row */
1759a8f1f9e5SMatthew G. Knepley       for (g = 0; g < NbJ; ++g) {
1760a8f1f9e5SMatthew G. Knepley         for (gc = 0; gc < NcJ; ++gc) {
1761a8f1f9e5SMatthew G. Knepley           const PetscInt gidx = g*NcJ+gc; /* Trial function basis index */
1762a8f1f9e5SMatthew G. Knepley           const PetscInt j    = offsetJ+g; /* Element matrix column */
1763a8f1f9e5SMatthew G. Knepley           const PetscInt fOff = eOffset+i*totDim+j;
1764a8f1f9e5SMatthew G. Knepley 
1765a8f1f9e5SMatthew G. Knepley           elemMat[fOff] += tmpBasisI[fidx]*g0[fc*NcJ+gc]*tmpBasisJ[gidx];
1766a8f1f9e5SMatthew G. Knepley           for (df = 0; df < dim; ++df) {
1767a8f1f9e5SMatthew G. Knepley             elemMat[fOff] += tmpBasisI[fidx]*g1[(fc*NcJ+gc)*dim+df]*tmpBasisDerJ[gidx*dim+df];
1768a8f1f9e5SMatthew G. Knepley             elemMat[fOff] += tmpBasisDerI[fidx*dim+df]*g2[(fc*NcJ+gc)*dim+df]*tmpBasisJ[gidx];
1769a8f1f9e5SMatthew G. Knepley             for (dg = 0; dg < dim; ++dg) {
1770a8f1f9e5SMatthew G. Knepley               elemMat[fOff] += tmpBasisDerI[fidx*dim+df]*g3[((fc*NcJ+gc)*dim+df)*dim+dg]*tmpBasisDerJ[gidx*dim+dg];
1771a8f1f9e5SMatthew G. Knepley             }
1772a8f1f9e5SMatthew G. Knepley           }
1773a8f1f9e5SMatthew G. Knepley         }
1774a8f1f9e5SMatthew G. Knepley       }
1775a8f1f9e5SMatthew G. Knepley     }
1776a8f1f9e5SMatthew G. Knepley   }
1777a8f1f9e5SMatthew G. Knepley   return(0);
1778a8f1f9e5SMatthew G. Knepley }
1779c9ba7969SMatthew G. Knepley 
1780c9ba7969SMatthew G. Knepley PetscErrorCode PetscFECreateCellGeometry(PetscFE fe, PetscQuadrature quad, PetscFEGeom *cgeom)
1781c9ba7969SMatthew G. Knepley {
1782c9ba7969SMatthew G. Knepley   PetscDualSpace  dsp;
1783c9ba7969SMatthew G. Knepley   DM              dm;
1784c9ba7969SMatthew G. Knepley   PetscQuadrature quadDef;
1785c9ba7969SMatthew G. Knepley   PetscInt        dim, cdim, Nq;
1786c9ba7969SMatthew G. Knepley   PetscErrorCode  ierr;
1787c9ba7969SMatthew G. Knepley 
1788c9ba7969SMatthew G. Knepley   PetscFunctionBegin;
1789c9ba7969SMatthew G. Knepley   ierr = PetscFEGetDualSpace(fe, &dsp);CHKERRQ(ierr);
1790c9ba7969SMatthew G. Knepley   ierr = PetscDualSpaceGetDM(dsp, &dm);CHKERRQ(ierr);
1791c9ba7969SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1792c9ba7969SMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &cdim);CHKERRQ(ierr);
1793c9ba7969SMatthew G. Knepley   ierr = PetscFEGetQuadrature(fe, &quadDef);CHKERRQ(ierr);
1794c9ba7969SMatthew G. Knepley   quad = quad ? quad : quadDef;
1795c9ba7969SMatthew G. Knepley   ierr = PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);
1796c9ba7969SMatthew G. Knepley   ierr = PetscMalloc1(Nq*cdim,      &cgeom->v);CHKERRQ(ierr);
1797c9ba7969SMatthew G. Knepley   ierr = PetscMalloc1(Nq*cdim*cdim, &cgeom->J);CHKERRQ(ierr);
1798c9ba7969SMatthew G. Knepley   ierr = PetscMalloc1(Nq*cdim*cdim, &cgeom->invJ);CHKERRQ(ierr);
1799c9ba7969SMatthew G. Knepley   ierr = PetscMalloc1(Nq,           &cgeom->detJ);CHKERRQ(ierr);
1800c9ba7969SMatthew G. Knepley   cgeom->dim       = dim;
1801c9ba7969SMatthew G. Knepley   cgeom->dimEmbed  = cdim;
1802c9ba7969SMatthew G. Knepley   cgeom->numCells  = 1;
1803c9ba7969SMatthew G. Knepley   cgeom->numPoints = Nq;
1804c9ba7969SMatthew G. Knepley   ierr = DMPlexComputeCellGeometryFEM(dm, 0, quad, cgeom->v, cgeom->J, cgeom->invJ, cgeom->detJ);CHKERRQ(ierr);
1805c9ba7969SMatthew G. Knepley   PetscFunctionReturn(0);
1806c9ba7969SMatthew G. Knepley }
1807c9ba7969SMatthew G. Knepley 
1808c9ba7969SMatthew G. Knepley PetscErrorCode PetscFEDestroyCellGeometry(PetscFE fe, PetscFEGeom *cgeom)
1809c9ba7969SMatthew G. Knepley {
1810c9ba7969SMatthew G. Knepley   PetscErrorCode ierr;
1811c9ba7969SMatthew G. Knepley 
1812c9ba7969SMatthew G. Knepley   PetscFunctionBegin;
1813c9ba7969SMatthew G. Knepley   ierr = PetscFree(cgeom->v);CHKERRQ(ierr);
1814c9ba7969SMatthew G. Knepley   ierr = PetscFree(cgeom->J);CHKERRQ(ierr);
1815c9ba7969SMatthew G. Knepley   ierr = PetscFree(cgeom->invJ);CHKERRQ(ierr);
1816c9ba7969SMatthew G. Knepley   ierr = PetscFree(cgeom->detJ);CHKERRQ(ierr);
1817c9ba7969SMatthew G. Knepley   PetscFunctionReturn(0);
1818c9ba7969SMatthew G. Knepley }
1819