xref: /petsc/src/dm/dt/fe/interface/fe.c (revision dbbe0bcd3f3a8fbab5a45420dc06f8387e5764c6)
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 
53ead873ccSMatthew G. Knepley PetscLogEvent PETSCFE_SetUp;
54ead873ccSMatthew G. Knepley 
5520cf1dd8SToby Isaac PetscFunctionList PetscFEList              = NULL;
5620cf1dd8SToby Isaac PetscBool         PetscFERegisterAllCalled = PETSC_FALSE;
5720cf1dd8SToby Isaac 
5820cf1dd8SToby Isaac /*@C
5920cf1dd8SToby Isaac   PetscFERegister - Adds a new PetscFE implementation
6020cf1dd8SToby Isaac 
6120cf1dd8SToby Isaac   Not Collective
6220cf1dd8SToby Isaac 
6320cf1dd8SToby Isaac   Input Parameters:
6420cf1dd8SToby Isaac + name        - The name of a new user-defined creation routine
6520cf1dd8SToby Isaac - create_func - The creation routine itself
6620cf1dd8SToby Isaac 
6720cf1dd8SToby Isaac   Notes:
6820cf1dd8SToby Isaac   PetscFERegister() may be called multiple times to add several user-defined PetscFEs
6920cf1dd8SToby Isaac 
7020cf1dd8SToby Isaac   Sample usage:
7120cf1dd8SToby Isaac .vb
7220cf1dd8SToby Isaac     PetscFERegister("my_fe", MyPetscFECreate);
7320cf1dd8SToby Isaac .ve
7420cf1dd8SToby Isaac 
7520cf1dd8SToby Isaac   Then, your PetscFE type can be chosen with the procedural interface via
7620cf1dd8SToby Isaac .vb
7720cf1dd8SToby Isaac     PetscFECreate(MPI_Comm, PetscFE *);
7820cf1dd8SToby Isaac     PetscFESetType(PetscFE, "my_fe");
7920cf1dd8SToby Isaac .ve
8020cf1dd8SToby Isaac    or at runtime via the option
8120cf1dd8SToby Isaac .vb
8220cf1dd8SToby Isaac     -petscfe_type my_fe
8320cf1dd8SToby Isaac .ve
8420cf1dd8SToby Isaac 
8520cf1dd8SToby Isaac   Level: advanced
8620cf1dd8SToby Isaac 
87db781477SPatrick Sanan .seealso: `PetscFERegisterAll()`, `PetscFERegisterDestroy()`
8820cf1dd8SToby Isaac 
8920cf1dd8SToby Isaac @*/
9020cf1dd8SToby Isaac PetscErrorCode PetscFERegister(const char sname[], PetscErrorCode (*function)(PetscFE))
9120cf1dd8SToby Isaac {
9220cf1dd8SToby Isaac   PetscFunctionBegin;
939566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&PetscFEList, sname, function));
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 
111db781477SPatrick Sanan .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 
11820cf1dd8SToby Isaac   PetscFunctionBegin;
11920cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
1209566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject) fem, name, &match));
12120cf1dd8SToby Isaac   if (match) PetscFunctionReturn(0);
12220cf1dd8SToby Isaac 
1239566063dSJacob Faibussowitsch   if (!PetscFERegisterAllCalled) PetscCall(PetscFERegisterAll());
1249566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(PetscFEList, name, &r));
12528b400f6SJacob Faibussowitsch   PetscCheck(r,PetscObjectComm((PetscObject) fem), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscFE type: %s", name);
12620cf1dd8SToby Isaac 
127*dbbe0bcdSBarry Smith   PetscTryTypeMethod(fem,destroy);
12820cf1dd8SToby Isaac   fem->ops->destroy = NULL;
129*dbbe0bcdSBarry Smith 
1309566063dSJacob Faibussowitsch   PetscCall((*r)(fem));
1319566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject) fem, name));
13220cf1dd8SToby Isaac   PetscFunctionReturn(0);
13320cf1dd8SToby Isaac }
13420cf1dd8SToby Isaac 
13520cf1dd8SToby Isaac /*@C
13620cf1dd8SToby Isaac   PetscFEGetType - Gets the PetscFE type name (as a string) from the object.
13720cf1dd8SToby Isaac 
13820cf1dd8SToby Isaac   Not Collective
13920cf1dd8SToby Isaac 
14020cf1dd8SToby Isaac   Input Parameter:
14120cf1dd8SToby Isaac . fem  - The PetscFE
14220cf1dd8SToby Isaac 
14320cf1dd8SToby Isaac   Output Parameter:
14420cf1dd8SToby Isaac . name - The PetscFE type name
14520cf1dd8SToby Isaac 
14620cf1dd8SToby Isaac   Level: intermediate
14720cf1dd8SToby Isaac 
148db781477SPatrick Sanan .seealso: `PetscFESetType()`, `PetscFECreate()`
14920cf1dd8SToby Isaac @*/
15020cf1dd8SToby Isaac PetscErrorCode PetscFEGetType(PetscFE fem, PetscFEType *name)
15120cf1dd8SToby Isaac {
15220cf1dd8SToby Isaac   PetscFunctionBegin;
15320cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
15420cf1dd8SToby Isaac   PetscValidPointer(name, 2);
15520cf1dd8SToby Isaac   if (!PetscFERegisterAllCalled) {
1569566063dSJacob Faibussowitsch     PetscCall(PetscFERegisterAll());
15720cf1dd8SToby Isaac   }
15820cf1dd8SToby Isaac   *name = ((PetscObject) fem)->type_name;
15920cf1dd8SToby Isaac   PetscFunctionReturn(0);
16020cf1dd8SToby Isaac }
16120cf1dd8SToby Isaac 
16220cf1dd8SToby Isaac /*@C
163fe2efc57SMark    PetscFEViewFromOptions - View from Options
164fe2efc57SMark 
165fe2efc57SMark    Collective on PetscFE
166fe2efc57SMark 
167fe2efc57SMark    Input Parameters:
168fe2efc57SMark +  A - the PetscFE object
169fe2efc57SMark .  obj - Optional object
170fe2efc57SMark -  name - command line option
171fe2efc57SMark 
172fe2efc57SMark    Level: intermediate
173db781477SPatrick Sanan .seealso: `PetscFE()`, `PetscFEView()`, `PetscObjectViewFromOptions()`, `PetscFECreate()`
174fe2efc57SMark @*/
175fe2efc57SMark PetscErrorCode  PetscFEViewFromOptions(PetscFE A,PetscObject obj,const char name[])
176fe2efc57SMark {
177fe2efc57SMark   PetscFunctionBegin;
178fe2efc57SMark   PetscValidHeaderSpecific(A,PETSCFE_CLASSID,1);
1799566063dSJacob Faibussowitsch   PetscCall(PetscObjectViewFromOptions((PetscObject)A,obj,name));
180fe2efc57SMark   PetscFunctionReturn(0);
181fe2efc57SMark }
182fe2efc57SMark 
183fe2efc57SMark /*@C
18420cf1dd8SToby Isaac   PetscFEView - Views a PetscFE
18520cf1dd8SToby Isaac 
186d083f849SBarry Smith   Collective on fem
18720cf1dd8SToby Isaac 
188d8d19677SJose E. Roman   Input Parameters:
18920cf1dd8SToby Isaac + fem - the PetscFE object to view
190d9bac1caSLisandro Dalcin - viewer   - the viewer
19120cf1dd8SToby Isaac 
1922b99622eSMatthew G. Knepley   Level: beginner
19320cf1dd8SToby Isaac 
194db781477SPatrick Sanan .seealso `PetscFEDestroy()`
19520cf1dd8SToby Isaac @*/
196d9bac1caSLisandro Dalcin PetscErrorCode PetscFEView(PetscFE fem, PetscViewer viewer)
19720cf1dd8SToby Isaac {
198d9bac1caSLisandro Dalcin   PetscBool      iascii;
19920cf1dd8SToby Isaac 
20020cf1dd8SToby Isaac   PetscFunctionBegin;
20120cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
202d9bac1caSLisandro Dalcin   if (viewer) PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
2039566063dSJacob Faibussowitsch   if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) fem), &viewer));
2049566063dSJacob Faibussowitsch   PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)fem, viewer));
2059566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii));
206*dbbe0bcdSBarry Smith   PetscTryTypeMethod(fem,view , viewer);
20720cf1dd8SToby Isaac   PetscFunctionReturn(0);
20820cf1dd8SToby Isaac }
20920cf1dd8SToby Isaac 
21020cf1dd8SToby Isaac /*@
21120cf1dd8SToby Isaac   PetscFESetFromOptions - sets parameters in a PetscFE from the options database
21220cf1dd8SToby Isaac 
213d083f849SBarry Smith   Collective on fem
21420cf1dd8SToby Isaac 
21520cf1dd8SToby Isaac   Input Parameter:
21620cf1dd8SToby Isaac . fem - the PetscFE object to set options for
21720cf1dd8SToby Isaac 
21820cf1dd8SToby Isaac   Options Database:
219a2b725a8SWilliam Gropp + -petscfe_num_blocks  - the number of cell blocks to integrate concurrently
220a2b725a8SWilliam Gropp - -petscfe_num_batches - the number of cell batches to integrate serially
22120cf1dd8SToby Isaac 
2222b99622eSMatthew G. Knepley   Level: intermediate
22320cf1dd8SToby Isaac 
224db781477SPatrick Sanan .seealso `PetscFEView()`
22520cf1dd8SToby Isaac @*/
22620cf1dd8SToby Isaac PetscErrorCode PetscFESetFromOptions(PetscFE fem)
22720cf1dd8SToby Isaac {
22820cf1dd8SToby Isaac   const char    *defaultType;
22920cf1dd8SToby Isaac   char           name[256];
23020cf1dd8SToby Isaac   PetscBool      flg;
23120cf1dd8SToby Isaac 
23220cf1dd8SToby Isaac   PetscFunctionBegin;
23320cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
23420cf1dd8SToby Isaac   if (!((PetscObject) fem)->type_name) {
23520cf1dd8SToby Isaac     defaultType = PETSCFEBASIC;
23620cf1dd8SToby Isaac   } else {
23720cf1dd8SToby Isaac     defaultType = ((PetscObject) fem)->type_name;
23820cf1dd8SToby Isaac   }
2399566063dSJacob Faibussowitsch   if (!PetscFERegisterAllCalled) PetscCall(PetscFERegisterAll());
24020cf1dd8SToby Isaac 
241d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject) fem);
2429566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-petscfe_type", "Finite element space", "PetscFESetType", PetscFEList, defaultType, name, 256, &flg));
24320cf1dd8SToby Isaac   if (flg) {
2449566063dSJacob Faibussowitsch     PetscCall(PetscFESetType(fem, name));
24520cf1dd8SToby Isaac   } else if (!((PetscObject) fem)->type_name) {
2469566063dSJacob Faibussowitsch     PetscCall(PetscFESetType(fem, defaultType));
24720cf1dd8SToby Isaac   }
2489566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-petscfe_num_blocks", "The number of cell blocks to integrate concurrently", "PetscSpaceSetTileSizes", fem->numBlocks, &fem->numBlocks, NULL,1));
2499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-petscfe_num_batches", "The number of cell batches to integrate serially", "PetscSpaceSetTileSizes", fem->numBatches, &fem->numBatches, NULL,1));
250*dbbe0bcdSBarry Smith   PetscTryTypeMethod(fem,setfromoptions,PetscOptionsObject);
25120cf1dd8SToby Isaac   /* process any options handlers added with PetscObjectAddOptionsHandler() */
252*dbbe0bcdSBarry Smith   PetscCall(PetscObjectProcessOptionsHandlers((PetscObject) fem,PetscOptionsObject));
253d0609cedSBarry Smith   PetscOptionsEnd();
2549566063dSJacob Faibussowitsch   PetscCall(PetscFEViewFromOptions(fem, NULL, "-petscfe_view"));
25520cf1dd8SToby Isaac   PetscFunctionReturn(0);
25620cf1dd8SToby Isaac }
25720cf1dd8SToby Isaac 
25820cf1dd8SToby Isaac /*@C
25920cf1dd8SToby Isaac   PetscFESetUp - Construct data structures for the PetscFE
26020cf1dd8SToby Isaac 
261d083f849SBarry Smith   Collective on fem
26220cf1dd8SToby Isaac 
26320cf1dd8SToby Isaac   Input Parameter:
26420cf1dd8SToby Isaac . fem - the PetscFE object to setup
26520cf1dd8SToby Isaac 
2662b99622eSMatthew G. Knepley   Level: intermediate
26720cf1dd8SToby Isaac 
268db781477SPatrick Sanan .seealso `PetscFEView()`, `PetscFEDestroy()`
26920cf1dd8SToby Isaac @*/
27020cf1dd8SToby Isaac PetscErrorCode PetscFESetUp(PetscFE fem)
27120cf1dd8SToby Isaac {
27220cf1dd8SToby Isaac   PetscFunctionBegin;
27320cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
27420cf1dd8SToby Isaac   if (fem->setupcalled) PetscFunctionReturn(0);
2759566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(PETSCFE_SetUp, fem, 0, 0, 0));
27620cf1dd8SToby Isaac   fem->setupcalled = PETSC_TRUE;
277*dbbe0bcdSBarry Smith   PetscTryTypeMethod(fem,setup);
2789566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(PETSCFE_SetUp, fem, 0, 0, 0));
27920cf1dd8SToby Isaac   PetscFunctionReturn(0);
28020cf1dd8SToby Isaac }
28120cf1dd8SToby Isaac 
28220cf1dd8SToby Isaac /*@
28320cf1dd8SToby Isaac   PetscFEDestroy - Destroys a PetscFE object
28420cf1dd8SToby Isaac 
285d083f849SBarry Smith   Collective on fem
28620cf1dd8SToby Isaac 
28720cf1dd8SToby Isaac   Input Parameter:
28820cf1dd8SToby Isaac . fem - the PetscFE object to destroy
28920cf1dd8SToby Isaac 
2902b99622eSMatthew G. Knepley   Level: beginner
29120cf1dd8SToby Isaac 
292db781477SPatrick Sanan .seealso `PetscFEView()`
29320cf1dd8SToby Isaac @*/
29420cf1dd8SToby Isaac PetscErrorCode PetscFEDestroy(PetscFE *fem)
29520cf1dd8SToby Isaac {
29620cf1dd8SToby Isaac   PetscFunctionBegin;
29720cf1dd8SToby Isaac   if (!*fem) PetscFunctionReturn(0);
29820cf1dd8SToby Isaac   PetscValidHeaderSpecific((*fem), PETSCFE_CLASSID, 1);
29920cf1dd8SToby Isaac 
300ea78f98cSLisandro Dalcin   if (--((PetscObject)(*fem))->refct > 0) {*fem = NULL; PetscFunctionReturn(0);}
30120cf1dd8SToby Isaac   ((PetscObject) (*fem))->refct = 0;
30220cf1dd8SToby Isaac 
30320cf1dd8SToby Isaac   if ((*fem)->subspaces) {
30420cf1dd8SToby Isaac     PetscInt dim, d;
30520cf1dd8SToby Isaac 
3069566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDimension((*fem)->dualSpace, &dim));
3079566063dSJacob Faibussowitsch     for (d = 0; d < dim; ++d) PetscCall(PetscFEDestroy(&(*fem)->subspaces[d]));
30820cf1dd8SToby Isaac   }
3099566063dSJacob Faibussowitsch   PetscCall(PetscFree((*fem)->subspaces));
3109566063dSJacob Faibussowitsch   PetscCall(PetscFree((*fem)->invV));
3119566063dSJacob Faibussowitsch   PetscCall(PetscTabulationDestroy(&(*fem)->T));
3129566063dSJacob Faibussowitsch   PetscCall(PetscTabulationDestroy(&(*fem)->Tf));
3139566063dSJacob Faibussowitsch   PetscCall(PetscTabulationDestroy(&(*fem)->Tc));
3149566063dSJacob Faibussowitsch   PetscCall(PetscSpaceDestroy(&(*fem)->basisSpace));
3159566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceDestroy(&(*fem)->dualSpace));
3169566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureDestroy(&(*fem)->quadrature));
3179566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureDestroy(&(*fem)->faceQuadrature));
318f918ec44SMatthew G. Knepley #ifdef PETSC_HAVE_LIBCEED
3199566063dSJacob Faibussowitsch   PetscCallCEED(CeedBasisDestroy(&(*fem)->ceedBasis));
3209566063dSJacob Faibussowitsch   PetscCallCEED(CeedDestroy(&(*fem)->ceed));
321f918ec44SMatthew G. Knepley #endif
32220cf1dd8SToby Isaac 
323*dbbe0bcdSBarry Smith   PetscTryTypeMethod((*fem),destroy);
3249566063dSJacob Faibussowitsch   PetscCall(PetscHeaderDestroy(fem));
32520cf1dd8SToby Isaac   PetscFunctionReturn(0);
32620cf1dd8SToby Isaac }
32720cf1dd8SToby Isaac 
32820cf1dd8SToby Isaac /*@
32920cf1dd8SToby Isaac   PetscFECreate - Creates an empty PetscFE object. The type can then be set with PetscFESetType().
33020cf1dd8SToby Isaac 
331d083f849SBarry Smith   Collective
33220cf1dd8SToby Isaac 
33320cf1dd8SToby Isaac   Input Parameter:
33420cf1dd8SToby Isaac . comm - The communicator for the PetscFE object
33520cf1dd8SToby Isaac 
33620cf1dd8SToby Isaac   Output Parameter:
33720cf1dd8SToby Isaac . fem - The PetscFE object
33820cf1dd8SToby Isaac 
33920cf1dd8SToby Isaac   Level: beginner
34020cf1dd8SToby Isaac 
341db781477SPatrick Sanan .seealso: `PetscFESetType()`, `PETSCFEGALERKIN`
34220cf1dd8SToby Isaac @*/
34320cf1dd8SToby Isaac PetscErrorCode PetscFECreate(MPI_Comm comm, PetscFE *fem)
34420cf1dd8SToby Isaac {
34520cf1dd8SToby Isaac   PetscFE        f;
34620cf1dd8SToby Isaac 
34720cf1dd8SToby Isaac   PetscFunctionBegin;
34820cf1dd8SToby Isaac   PetscValidPointer(fem, 2);
3499566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(FECitation,&FEcite));
35020cf1dd8SToby Isaac   *fem = NULL;
3519566063dSJacob Faibussowitsch   PetscCall(PetscFEInitializePackage());
35220cf1dd8SToby Isaac 
3539566063dSJacob Faibussowitsch   PetscCall(PetscHeaderCreate(f, PETSCFE_CLASSID, "PetscFE", "Finite Element", "PetscFE", comm, PetscFEDestroy, PetscFEView));
35420cf1dd8SToby Isaac 
35520cf1dd8SToby Isaac   f->basisSpace    = NULL;
35620cf1dd8SToby Isaac   f->dualSpace     = NULL;
35720cf1dd8SToby Isaac   f->numComponents = 1;
35820cf1dd8SToby Isaac   f->subspaces     = NULL;
35920cf1dd8SToby Isaac   f->invV          = NULL;
360ef0bb6c7SMatthew G. Knepley   f->T             = NULL;
361ef0bb6c7SMatthew G. Knepley   f->Tf            = NULL;
362ef0bb6c7SMatthew G. Knepley   f->Tc            = NULL;
3639566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(&f->quadrature, 1));
3649566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(&f->faceQuadrature, 1));
36520cf1dd8SToby Isaac   f->blockSize     = 0;
36620cf1dd8SToby Isaac   f->numBlocks     = 1;
36720cf1dd8SToby Isaac   f->batchSize     = 0;
36820cf1dd8SToby Isaac   f->numBatches    = 1;
36920cf1dd8SToby Isaac 
37020cf1dd8SToby Isaac   *fem = f;
37120cf1dd8SToby Isaac   PetscFunctionReturn(0);
37220cf1dd8SToby Isaac }
37320cf1dd8SToby Isaac 
37420cf1dd8SToby Isaac /*@
37520cf1dd8SToby Isaac   PetscFEGetSpatialDimension - Returns the spatial dimension of the element
37620cf1dd8SToby Isaac 
37720cf1dd8SToby Isaac   Not collective
37820cf1dd8SToby Isaac 
37920cf1dd8SToby Isaac   Input Parameter:
38020cf1dd8SToby Isaac . fem - The PetscFE object
38120cf1dd8SToby Isaac 
38220cf1dd8SToby Isaac   Output Parameter:
38320cf1dd8SToby Isaac . dim - The spatial dimension
38420cf1dd8SToby Isaac 
38520cf1dd8SToby Isaac   Level: intermediate
38620cf1dd8SToby Isaac 
387db781477SPatrick Sanan .seealso: `PetscFECreate()`
38820cf1dd8SToby Isaac @*/
38920cf1dd8SToby Isaac PetscErrorCode PetscFEGetSpatialDimension(PetscFE fem, PetscInt *dim)
39020cf1dd8SToby Isaac {
39120cf1dd8SToby Isaac   DM dm;
39220cf1dd8SToby Isaac 
39320cf1dd8SToby Isaac   PetscFunctionBegin;
39420cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
395dadcf809SJacob Faibussowitsch   PetscValidIntPointer(dim, 2);
3969566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDM(fem->dualSpace, &dm));
3979566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, dim));
39820cf1dd8SToby Isaac   PetscFunctionReturn(0);
39920cf1dd8SToby Isaac }
40020cf1dd8SToby Isaac 
40120cf1dd8SToby Isaac /*@
40220cf1dd8SToby Isaac   PetscFESetNumComponents - Sets the number of components in the element
40320cf1dd8SToby Isaac 
40420cf1dd8SToby Isaac   Not collective
40520cf1dd8SToby Isaac 
40620cf1dd8SToby Isaac   Input Parameters:
40720cf1dd8SToby Isaac + fem - The PetscFE object
40820cf1dd8SToby Isaac - comp - The number of field components
40920cf1dd8SToby Isaac 
41020cf1dd8SToby Isaac   Level: intermediate
41120cf1dd8SToby Isaac 
412db781477SPatrick Sanan .seealso: `PetscFECreate()`
41320cf1dd8SToby Isaac @*/
41420cf1dd8SToby Isaac PetscErrorCode PetscFESetNumComponents(PetscFE fem, PetscInt comp)
41520cf1dd8SToby Isaac {
41620cf1dd8SToby Isaac   PetscFunctionBegin;
41720cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
41820cf1dd8SToby Isaac   fem->numComponents = comp;
41920cf1dd8SToby Isaac   PetscFunctionReturn(0);
42020cf1dd8SToby Isaac }
42120cf1dd8SToby Isaac 
42220cf1dd8SToby Isaac /*@
42320cf1dd8SToby Isaac   PetscFEGetNumComponents - Returns the number of components in the element
42420cf1dd8SToby Isaac 
42520cf1dd8SToby Isaac   Not collective
42620cf1dd8SToby Isaac 
42720cf1dd8SToby Isaac   Input Parameter:
42820cf1dd8SToby Isaac . fem - The PetscFE object
42920cf1dd8SToby Isaac 
43020cf1dd8SToby Isaac   Output Parameter:
43120cf1dd8SToby Isaac . comp - The number of field components
43220cf1dd8SToby Isaac 
43320cf1dd8SToby Isaac   Level: intermediate
43420cf1dd8SToby Isaac 
435db781477SPatrick Sanan .seealso: `PetscFECreate()`
43620cf1dd8SToby Isaac @*/
43720cf1dd8SToby Isaac PetscErrorCode PetscFEGetNumComponents(PetscFE fem, PetscInt *comp)
43820cf1dd8SToby Isaac {
43920cf1dd8SToby Isaac   PetscFunctionBegin;
44020cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
441dadcf809SJacob Faibussowitsch   PetscValidIntPointer(comp, 2);
44220cf1dd8SToby Isaac   *comp = fem->numComponents;
44320cf1dd8SToby Isaac   PetscFunctionReturn(0);
44420cf1dd8SToby Isaac }
44520cf1dd8SToby Isaac 
44620cf1dd8SToby Isaac /*@
44720cf1dd8SToby Isaac   PetscFESetTileSizes - Sets the tile sizes for evaluation
44820cf1dd8SToby Isaac 
44920cf1dd8SToby Isaac   Not collective
45020cf1dd8SToby Isaac 
45120cf1dd8SToby Isaac   Input Parameters:
45220cf1dd8SToby Isaac + fem - The PetscFE object
45320cf1dd8SToby Isaac . blockSize - The number of elements in a block
45420cf1dd8SToby Isaac . numBlocks - The number of blocks in a batch
45520cf1dd8SToby Isaac . batchSize - The number of elements in a batch
45620cf1dd8SToby Isaac - numBatches - The number of batches in a chunk
45720cf1dd8SToby Isaac 
45820cf1dd8SToby Isaac   Level: intermediate
45920cf1dd8SToby Isaac 
460db781477SPatrick Sanan .seealso: `PetscFECreate()`
46120cf1dd8SToby Isaac @*/
46220cf1dd8SToby Isaac PetscErrorCode PetscFESetTileSizes(PetscFE fem, PetscInt blockSize, PetscInt numBlocks, PetscInt batchSize, PetscInt numBatches)
46320cf1dd8SToby Isaac {
46420cf1dd8SToby Isaac   PetscFunctionBegin;
46520cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
46620cf1dd8SToby Isaac   fem->blockSize  = blockSize;
46720cf1dd8SToby Isaac   fem->numBlocks  = numBlocks;
46820cf1dd8SToby Isaac   fem->batchSize  = batchSize;
46920cf1dd8SToby Isaac   fem->numBatches = numBatches;
47020cf1dd8SToby Isaac   PetscFunctionReturn(0);
47120cf1dd8SToby Isaac }
47220cf1dd8SToby Isaac 
47320cf1dd8SToby Isaac /*@
47420cf1dd8SToby Isaac   PetscFEGetTileSizes - Returns the tile sizes for evaluation
47520cf1dd8SToby Isaac 
47620cf1dd8SToby Isaac   Not collective
47720cf1dd8SToby Isaac 
47820cf1dd8SToby Isaac   Input Parameter:
47920cf1dd8SToby Isaac . fem - The PetscFE object
48020cf1dd8SToby Isaac 
48120cf1dd8SToby Isaac   Output Parameters:
48220cf1dd8SToby Isaac + blockSize - The number of elements in a block
48320cf1dd8SToby Isaac . numBlocks - The number of blocks in a batch
48420cf1dd8SToby Isaac . batchSize - The number of elements in a batch
48520cf1dd8SToby Isaac - numBatches - The number of batches in a chunk
48620cf1dd8SToby Isaac 
48720cf1dd8SToby Isaac   Level: intermediate
48820cf1dd8SToby Isaac 
489db781477SPatrick Sanan .seealso: `PetscFECreate()`
49020cf1dd8SToby Isaac @*/
49120cf1dd8SToby Isaac PetscErrorCode PetscFEGetTileSizes(PetscFE fem, PetscInt *blockSize, PetscInt *numBlocks, PetscInt *batchSize, PetscInt *numBatches)
49220cf1dd8SToby Isaac {
49320cf1dd8SToby Isaac   PetscFunctionBegin;
49420cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
495dadcf809SJacob Faibussowitsch   if (blockSize)  PetscValidIntPointer(blockSize,  2);
496dadcf809SJacob Faibussowitsch   if (numBlocks)  PetscValidIntPointer(numBlocks,  3);
497dadcf809SJacob Faibussowitsch   if (batchSize)  PetscValidIntPointer(batchSize,  4);
498dadcf809SJacob Faibussowitsch   if (numBatches) PetscValidIntPointer(numBatches, 5);
49920cf1dd8SToby Isaac   if (blockSize)  *blockSize  = fem->blockSize;
50020cf1dd8SToby Isaac   if (numBlocks)  *numBlocks  = fem->numBlocks;
50120cf1dd8SToby Isaac   if (batchSize)  *batchSize  = fem->batchSize;
50220cf1dd8SToby Isaac   if (numBatches) *numBatches = fem->numBatches;
50320cf1dd8SToby Isaac   PetscFunctionReturn(0);
50420cf1dd8SToby Isaac }
50520cf1dd8SToby Isaac 
50620cf1dd8SToby Isaac /*@
50720cf1dd8SToby Isaac   PetscFEGetBasisSpace - Returns the PetscSpace used for approximation of the solution
50820cf1dd8SToby Isaac 
50920cf1dd8SToby Isaac   Not collective
51020cf1dd8SToby Isaac 
51120cf1dd8SToby Isaac   Input Parameter:
51220cf1dd8SToby Isaac . fem - The PetscFE object
51320cf1dd8SToby Isaac 
51420cf1dd8SToby Isaac   Output Parameter:
51520cf1dd8SToby Isaac . sp - The PetscSpace object
51620cf1dd8SToby Isaac 
51720cf1dd8SToby Isaac   Level: intermediate
51820cf1dd8SToby Isaac 
519db781477SPatrick Sanan .seealso: `PetscFECreate()`
52020cf1dd8SToby Isaac @*/
52120cf1dd8SToby Isaac PetscErrorCode PetscFEGetBasisSpace(PetscFE fem, PetscSpace *sp)
52220cf1dd8SToby Isaac {
52320cf1dd8SToby Isaac   PetscFunctionBegin;
52420cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
52520cf1dd8SToby Isaac   PetscValidPointer(sp, 2);
52620cf1dd8SToby Isaac   *sp = fem->basisSpace;
52720cf1dd8SToby Isaac   PetscFunctionReturn(0);
52820cf1dd8SToby Isaac }
52920cf1dd8SToby Isaac 
53020cf1dd8SToby Isaac /*@
53120cf1dd8SToby Isaac   PetscFESetBasisSpace - Sets the PetscSpace used for approximation of the solution
53220cf1dd8SToby Isaac 
53320cf1dd8SToby Isaac   Not collective
53420cf1dd8SToby Isaac 
53520cf1dd8SToby Isaac   Input Parameters:
53620cf1dd8SToby Isaac + fem - The PetscFE object
53720cf1dd8SToby Isaac - sp - The PetscSpace object
53820cf1dd8SToby Isaac 
53920cf1dd8SToby Isaac   Level: intermediate
54020cf1dd8SToby Isaac 
541db781477SPatrick Sanan .seealso: `PetscFECreate()`
54220cf1dd8SToby Isaac @*/
54320cf1dd8SToby Isaac PetscErrorCode PetscFESetBasisSpace(PetscFE fem, PetscSpace sp)
54420cf1dd8SToby Isaac {
54520cf1dd8SToby Isaac   PetscFunctionBegin;
54620cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
54720cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 2);
5489566063dSJacob Faibussowitsch   PetscCall(PetscSpaceDestroy(&fem->basisSpace));
54920cf1dd8SToby Isaac   fem->basisSpace = sp;
5509566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject) fem->basisSpace));
55120cf1dd8SToby Isaac   PetscFunctionReturn(0);
55220cf1dd8SToby Isaac }
55320cf1dd8SToby Isaac 
55420cf1dd8SToby Isaac /*@
55520cf1dd8SToby Isaac   PetscFEGetDualSpace - Returns the PetscDualSpace used to define the inner product
55620cf1dd8SToby Isaac 
55720cf1dd8SToby Isaac   Not collective
55820cf1dd8SToby Isaac 
55920cf1dd8SToby Isaac   Input Parameter:
56020cf1dd8SToby Isaac . fem - The PetscFE object
56120cf1dd8SToby Isaac 
56220cf1dd8SToby Isaac   Output Parameter:
56320cf1dd8SToby Isaac . sp - The PetscDualSpace object
56420cf1dd8SToby Isaac 
56520cf1dd8SToby Isaac   Level: intermediate
56620cf1dd8SToby Isaac 
567db781477SPatrick Sanan .seealso: `PetscFECreate()`
56820cf1dd8SToby Isaac @*/
56920cf1dd8SToby Isaac PetscErrorCode PetscFEGetDualSpace(PetscFE fem, PetscDualSpace *sp)
57020cf1dd8SToby Isaac {
57120cf1dd8SToby Isaac   PetscFunctionBegin;
57220cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
57320cf1dd8SToby Isaac   PetscValidPointer(sp, 2);
57420cf1dd8SToby Isaac   *sp = fem->dualSpace;
57520cf1dd8SToby Isaac   PetscFunctionReturn(0);
57620cf1dd8SToby Isaac }
57720cf1dd8SToby Isaac 
57820cf1dd8SToby Isaac /*@
57920cf1dd8SToby Isaac   PetscFESetDualSpace - Sets the PetscDualSpace used to define the inner product
58020cf1dd8SToby Isaac 
58120cf1dd8SToby Isaac   Not collective
58220cf1dd8SToby Isaac 
58320cf1dd8SToby Isaac   Input Parameters:
58420cf1dd8SToby Isaac + fem - The PetscFE object
58520cf1dd8SToby Isaac - sp - The PetscDualSpace object
58620cf1dd8SToby Isaac 
58720cf1dd8SToby Isaac   Level: intermediate
58820cf1dd8SToby Isaac 
589db781477SPatrick Sanan .seealso: `PetscFECreate()`
59020cf1dd8SToby Isaac @*/
59120cf1dd8SToby Isaac PetscErrorCode PetscFESetDualSpace(PetscFE fem, PetscDualSpace sp)
59220cf1dd8SToby Isaac {
59320cf1dd8SToby Isaac   PetscFunctionBegin;
59420cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
59520cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 2);
5969566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceDestroy(&fem->dualSpace));
59720cf1dd8SToby Isaac   fem->dualSpace = sp;
5989566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject) fem->dualSpace));
59920cf1dd8SToby Isaac   PetscFunctionReturn(0);
60020cf1dd8SToby Isaac }
60120cf1dd8SToby Isaac 
60220cf1dd8SToby Isaac /*@
60320cf1dd8SToby Isaac   PetscFEGetQuadrature - Returns the PetscQuadrature used to calculate inner products
60420cf1dd8SToby Isaac 
60520cf1dd8SToby Isaac   Not collective
60620cf1dd8SToby Isaac 
60720cf1dd8SToby Isaac   Input Parameter:
60820cf1dd8SToby Isaac . fem - The PetscFE object
60920cf1dd8SToby Isaac 
61020cf1dd8SToby Isaac   Output Parameter:
61120cf1dd8SToby Isaac . q - The PetscQuadrature object
61220cf1dd8SToby Isaac 
61320cf1dd8SToby Isaac   Level: intermediate
61420cf1dd8SToby Isaac 
615db781477SPatrick Sanan .seealso: `PetscFECreate()`
61620cf1dd8SToby Isaac @*/
61720cf1dd8SToby Isaac PetscErrorCode PetscFEGetQuadrature(PetscFE fem, PetscQuadrature *q)
61820cf1dd8SToby Isaac {
61920cf1dd8SToby Isaac   PetscFunctionBegin;
62020cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
62120cf1dd8SToby Isaac   PetscValidPointer(q, 2);
62220cf1dd8SToby Isaac   *q = fem->quadrature;
62320cf1dd8SToby Isaac   PetscFunctionReturn(0);
62420cf1dd8SToby Isaac }
62520cf1dd8SToby Isaac 
62620cf1dd8SToby Isaac /*@
62720cf1dd8SToby Isaac   PetscFESetQuadrature - Sets the PetscQuadrature used to calculate inner products
62820cf1dd8SToby Isaac 
62920cf1dd8SToby Isaac   Not collective
63020cf1dd8SToby Isaac 
63120cf1dd8SToby Isaac   Input Parameters:
63220cf1dd8SToby Isaac + fem - The PetscFE object
63320cf1dd8SToby Isaac - q - The PetscQuadrature object
63420cf1dd8SToby Isaac 
63520cf1dd8SToby Isaac   Level: intermediate
63620cf1dd8SToby Isaac 
637db781477SPatrick Sanan .seealso: `PetscFECreate()`
63820cf1dd8SToby Isaac @*/
63920cf1dd8SToby Isaac PetscErrorCode PetscFESetQuadrature(PetscFE fem, PetscQuadrature q)
64020cf1dd8SToby Isaac {
64120cf1dd8SToby Isaac   PetscInt       Nc, qNc;
64220cf1dd8SToby Isaac 
64320cf1dd8SToby Isaac   PetscFunctionBegin;
64420cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
645fd2fdbddSMatthew G. Knepley   if (q == fem->quadrature) PetscFunctionReturn(0);
6469566063dSJacob Faibussowitsch   PetscCall(PetscFEGetNumComponents(fem, &Nc));
6479566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetNumComponents(q, &qNc));
64863a3b9bcSJacob Faibussowitsch   PetscCheck(!(qNc != 1) || !(Nc != qNc),PetscObjectComm((PetscObject) fem), PETSC_ERR_ARG_SIZ, "FE components %" PetscInt_FMT " != Quadrature components %" PetscInt_FMT " and non-scalar quadrature", Nc, qNc);
6499566063dSJacob Faibussowitsch   PetscCall(PetscTabulationDestroy(&fem->T));
6509566063dSJacob Faibussowitsch   PetscCall(PetscTabulationDestroy(&fem->Tc));
6519566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject) q));
6529566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureDestroy(&fem->quadrature));
65320cf1dd8SToby Isaac   fem->quadrature = q;
65420cf1dd8SToby Isaac   PetscFunctionReturn(0);
65520cf1dd8SToby Isaac }
65620cf1dd8SToby Isaac 
65720cf1dd8SToby Isaac /*@
65820cf1dd8SToby Isaac   PetscFEGetFaceQuadrature - Returns the PetscQuadrature used to calculate inner products on faces
65920cf1dd8SToby Isaac 
66020cf1dd8SToby Isaac   Not collective
66120cf1dd8SToby Isaac 
66220cf1dd8SToby Isaac   Input Parameter:
66320cf1dd8SToby Isaac . fem - The PetscFE object
66420cf1dd8SToby Isaac 
66520cf1dd8SToby Isaac   Output Parameter:
66620cf1dd8SToby Isaac . q - The PetscQuadrature object
66720cf1dd8SToby Isaac 
66820cf1dd8SToby Isaac   Level: intermediate
66920cf1dd8SToby Isaac 
670db781477SPatrick Sanan .seealso: `PetscFECreate()`
67120cf1dd8SToby Isaac @*/
67220cf1dd8SToby Isaac PetscErrorCode PetscFEGetFaceQuadrature(PetscFE fem, PetscQuadrature *q)
67320cf1dd8SToby Isaac {
67420cf1dd8SToby Isaac   PetscFunctionBegin;
67520cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
67620cf1dd8SToby Isaac   PetscValidPointer(q, 2);
67720cf1dd8SToby Isaac   *q = fem->faceQuadrature;
67820cf1dd8SToby Isaac   PetscFunctionReturn(0);
67920cf1dd8SToby Isaac }
68020cf1dd8SToby Isaac 
68120cf1dd8SToby Isaac /*@
68220cf1dd8SToby Isaac   PetscFESetFaceQuadrature - Sets the PetscQuadrature used to calculate inner products on faces
68320cf1dd8SToby Isaac 
68420cf1dd8SToby Isaac   Not collective
68520cf1dd8SToby Isaac 
68620cf1dd8SToby Isaac   Input Parameters:
68720cf1dd8SToby Isaac + fem - The PetscFE object
68820cf1dd8SToby Isaac - q - The PetscQuadrature object
68920cf1dd8SToby Isaac 
69020cf1dd8SToby Isaac   Level: intermediate
69120cf1dd8SToby Isaac 
692db781477SPatrick Sanan .seealso: `PetscFECreate()`
69320cf1dd8SToby Isaac @*/
69420cf1dd8SToby Isaac PetscErrorCode PetscFESetFaceQuadrature(PetscFE fem, PetscQuadrature q)
69520cf1dd8SToby Isaac {
696ef0bb6c7SMatthew G. Knepley   PetscInt       Nc, qNc;
69720cf1dd8SToby Isaac 
69820cf1dd8SToby Isaac   PetscFunctionBegin;
69920cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
7009566063dSJacob Faibussowitsch   PetscCall(PetscFEGetNumComponents(fem, &Nc));
7019566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetNumComponents(q, &qNc));
70263a3b9bcSJacob Faibussowitsch   PetscCheck(!(qNc != 1) || !(Nc != qNc),PetscObjectComm((PetscObject) fem), PETSC_ERR_ARG_SIZ, "FE components %" PetscInt_FMT " != Quadrature components %" PetscInt_FMT " and non-scalar quadrature", Nc, qNc);
7039566063dSJacob Faibussowitsch   PetscCall(PetscTabulationDestroy(&fem->Tf));
7049566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureDestroy(&fem->faceQuadrature));
70520cf1dd8SToby Isaac   fem->faceQuadrature = q;
7069566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject) q));
70720cf1dd8SToby Isaac   PetscFunctionReturn(0);
70820cf1dd8SToby Isaac }
70920cf1dd8SToby Isaac 
7105dc5c000SMatthew G. Knepley /*@
7115dc5c000SMatthew G. Knepley   PetscFECopyQuadrature - Copy both volumetric and surface quadrature
7125dc5c000SMatthew G. Knepley 
7135dc5c000SMatthew G. Knepley   Not collective
7145dc5c000SMatthew G. Knepley 
7155dc5c000SMatthew G. Knepley   Input Parameters:
7165dc5c000SMatthew G. Knepley + sfe - The PetscFE source for the quadratures
7175dc5c000SMatthew G. Knepley - tfe - The PetscFE target for the quadratures
7185dc5c000SMatthew G. Knepley 
7195dc5c000SMatthew G. Knepley   Level: intermediate
7205dc5c000SMatthew G. Knepley 
721db781477SPatrick Sanan .seealso: `PetscFECreate()`, `PetscFESetQuadrature()`, `PetscFESetFaceQuadrature()`
7225dc5c000SMatthew G. Knepley @*/
7235dc5c000SMatthew G. Knepley PetscErrorCode PetscFECopyQuadrature(PetscFE sfe, PetscFE tfe)
7245dc5c000SMatthew G. Knepley {
7255dc5c000SMatthew G. Knepley   PetscQuadrature q;
7265dc5c000SMatthew G. Knepley 
7275dc5c000SMatthew G. Knepley   PetscFunctionBegin;
7285dc5c000SMatthew G. Knepley   PetscValidHeaderSpecific(sfe, PETSCFE_CLASSID, 1);
7295dc5c000SMatthew G. Knepley   PetscValidHeaderSpecific(tfe, PETSCFE_CLASSID, 2);
7309566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(sfe, &q));
7319566063dSJacob Faibussowitsch   PetscCall(PetscFESetQuadrature(tfe,  q));
7329566063dSJacob Faibussowitsch   PetscCall(PetscFEGetFaceQuadrature(sfe, &q));
7339566063dSJacob Faibussowitsch   PetscCall(PetscFESetFaceQuadrature(tfe,  q));
7345dc5c000SMatthew G. Knepley   PetscFunctionReturn(0);
7355dc5c000SMatthew G. Knepley }
7365dc5c000SMatthew G. Knepley 
73720cf1dd8SToby Isaac /*@C
73820cf1dd8SToby Isaac   PetscFEGetNumDof - Returns the number of dofs (dual basis vectors) associated to mesh points on the reference cell of a given dimension
73920cf1dd8SToby Isaac 
74020cf1dd8SToby Isaac   Not collective
74120cf1dd8SToby Isaac 
74220cf1dd8SToby Isaac   Input Parameter:
74320cf1dd8SToby Isaac . fem - The PetscFE object
74420cf1dd8SToby Isaac 
74520cf1dd8SToby Isaac   Output Parameter:
74620cf1dd8SToby Isaac . numDof - Array with the number of dofs per dimension
74720cf1dd8SToby Isaac 
74820cf1dd8SToby Isaac   Level: intermediate
74920cf1dd8SToby Isaac 
750db781477SPatrick Sanan .seealso: `PetscFECreate()`
75120cf1dd8SToby Isaac @*/
75220cf1dd8SToby Isaac PetscErrorCode PetscFEGetNumDof(PetscFE fem, const PetscInt **numDof)
75320cf1dd8SToby Isaac {
75420cf1dd8SToby Isaac   PetscFunctionBegin;
75520cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
75620cf1dd8SToby Isaac   PetscValidPointer(numDof, 2);
7579566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetNumDof(fem->dualSpace, numDof));
75820cf1dd8SToby Isaac   PetscFunctionReturn(0);
75920cf1dd8SToby Isaac }
76020cf1dd8SToby Isaac 
76120cf1dd8SToby Isaac /*@C
762ef0bb6c7SMatthew G. Knepley   PetscFEGetCellTabulation - Returns the tabulation of the basis functions at the quadrature points on the reference cell
76320cf1dd8SToby Isaac 
76420cf1dd8SToby Isaac   Not collective
76520cf1dd8SToby Isaac 
766d8d19677SJose E. Roman   Input Parameters:
767f9244615SMatthew G. Knepley + fem - The PetscFE object
768f9244615SMatthew G. Knepley - k   - The highest derivative we need to tabulate, very often 1
76920cf1dd8SToby Isaac 
770ef0bb6c7SMatthew G. Knepley   Output Parameter:
771ef0bb6c7SMatthew G. Knepley . T - The basis function values and derivatives at quadrature points
77220cf1dd8SToby Isaac 
77320cf1dd8SToby Isaac   Note:
774ef0bb6c7SMatthew G. Knepley $ T->T[0] = B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c
775ef0bb6c7SMatthew G. Knepley $ T->T[1] = D[((p*pdim + i)*Nc + c)*dim + d] is the derivative value at point p for basis function i, component c, in direction d
776ef0bb6c7SMatthew G. Knepley $ T->T[2] = 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
77720cf1dd8SToby Isaac 
77820cf1dd8SToby Isaac   Level: intermediate
77920cf1dd8SToby Isaac 
780db781477SPatrick Sanan .seealso: `PetscFECreateTabulation()`, `PetscTabulationDestroy()`
78120cf1dd8SToby Isaac @*/
782f9244615SMatthew G. Knepley PetscErrorCode PetscFEGetCellTabulation(PetscFE fem, PetscInt k, PetscTabulation *T)
78320cf1dd8SToby Isaac {
78420cf1dd8SToby Isaac   PetscInt         npoints;
78520cf1dd8SToby Isaac   const PetscReal *points;
78620cf1dd8SToby Isaac 
78720cf1dd8SToby Isaac   PetscFunctionBegin;
78820cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
789064a246eSJacob Faibussowitsch   PetscValidPointer(T, 3);
7909566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(fem->quadrature, NULL, NULL, &npoints, &points, NULL));
7919566063dSJacob Faibussowitsch   if (!fem->T) PetscCall(PetscFECreateTabulation(fem, 1, npoints, points, k, &fem->T));
7921dca8a05SBarry Smith   PetscCheck(!fem->T || k <= fem->T->K,PetscObjectComm((PetscObject) fem), PETSC_ERR_ARG_OUTOFRANGE, "Requested %" PetscInt_FMT " derivatives, but only tabulated %" PetscInt_FMT, k, fem->T->K);
793ef0bb6c7SMatthew G. Knepley   *T = fem->T;
79420cf1dd8SToby Isaac   PetscFunctionReturn(0);
79520cf1dd8SToby Isaac }
79620cf1dd8SToby Isaac 
7972b99622eSMatthew G. Knepley /*@C
798ef0bb6c7SMatthew G. Knepley   PetscFEGetFaceTabulation - Returns the tabulation of the basis functions at the face quadrature points for each face of the reference cell
7992b99622eSMatthew G. Knepley 
8002b99622eSMatthew G. Knepley   Not collective
8012b99622eSMatthew G. Knepley 
802d8d19677SJose E. Roman   Input Parameters:
803f9244615SMatthew G. Knepley + fem - The PetscFE object
804f9244615SMatthew G. Knepley - k   - The highest derivative we need to tabulate, very often 1
8052b99622eSMatthew G. Knepley 
8062b99622eSMatthew G. Knepley   Output Parameters:
807a5b23f4aSJose E. Roman . Tf - The basis function values and derivatives at face quadrature points
8082b99622eSMatthew G. Knepley 
8092b99622eSMatthew G. Knepley   Note:
810ef0bb6c7SMatthew G. Knepley $ T->T[0] = Bf[((f*Nq + q)*pdim + i)*Nc + c] is the value at point f,q for basis function i and component c
811ef0bb6c7SMatthew G. Knepley $ T->T[1] = Df[(((f*Nq + q)*pdim + i)*Nc + c)*dim + d] is the derivative value at point f,q for basis function i, component c, in direction d
812ef0bb6c7SMatthew G. Knepley $ T->T[2] = Hf[((((f*Nq + q)*pdim + i)*Nc + c)*dim + d)*dim + e] is the value at point f,q for basis function i, component c, in directions d and e
8132b99622eSMatthew G. Knepley 
8142b99622eSMatthew G. Knepley   Level: intermediate
8152b99622eSMatthew G. Knepley 
816db781477SPatrick Sanan .seealso: `PetscFEGetCellTabulation()`, `PetscFECreateTabulation()`, `PetscTabulationDestroy()`
8172b99622eSMatthew G. Knepley @*/
818f9244615SMatthew G. Knepley PetscErrorCode PetscFEGetFaceTabulation(PetscFE fem, PetscInt k, PetscTabulation *Tf)
81920cf1dd8SToby Isaac {
82020cf1dd8SToby Isaac   PetscFunctionBegin;
82120cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
822064a246eSJacob Faibussowitsch   PetscValidPointer(Tf, 3);
823ef0bb6c7SMatthew G. Knepley   if (!fem->Tf) {
82420cf1dd8SToby Isaac     const PetscReal  xi0[3] = {-1., -1., -1.};
82520cf1dd8SToby Isaac     PetscReal        v0[3], J[9], detJ;
82620cf1dd8SToby Isaac     PetscQuadrature  fq;
82720cf1dd8SToby Isaac     PetscDualSpace   sp;
82820cf1dd8SToby Isaac     DM               dm;
82920cf1dd8SToby Isaac     const PetscInt  *faces;
83020cf1dd8SToby Isaac     PetscInt         dim, numFaces, f, npoints, q;
83120cf1dd8SToby Isaac     const PetscReal *points;
83220cf1dd8SToby Isaac     PetscReal       *facePoints;
83320cf1dd8SToby Isaac 
8349566063dSJacob Faibussowitsch     PetscCall(PetscFEGetDualSpace(fem, &sp));
8359566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDM(sp, &dm));
8369566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(dm, &dim));
8379566063dSJacob Faibussowitsch     PetscCall(DMPlexGetConeSize(dm, 0, &numFaces));
8389566063dSJacob Faibussowitsch     PetscCall(DMPlexGetCone(dm, 0, &faces));
8399566063dSJacob Faibussowitsch     PetscCall(PetscFEGetFaceQuadrature(fem, &fq));
84020cf1dd8SToby Isaac     if (fq) {
8419566063dSJacob Faibussowitsch       PetscCall(PetscQuadratureGetData(fq, NULL, NULL, &npoints, &points, NULL));
8429566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(numFaces*npoints*dim, &facePoints));
84320cf1dd8SToby Isaac       for (f = 0; f < numFaces; ++f) {
8449566063dSJacob Faibussowitsch         PetscCall(DMPlexComputeCellGeometryFEM(dm, faces[f], NULL, v0, J, NULL, &detJ));
84520cf1dd8SToby Isaac         for (q = 0; q < npoints; ++q) CoordinatesRefToReal(dim, dim-1, xi0, v0, J, &points[q*(dim-1)], &facePoints[(f*npoints+q)*dim]);
84620cf1dd8SToby Isaac       }
8479566063dSJacob Faibussowitsch       PetscCall(PetscFECreateTabulation(fem, numFaces, npoints, facePoints, k, &fem->Tf));
8489566063dSJacob Faibussowitsch       PetscCall(PetscFree(facePoints));
84920cf1dd8SToby Isaac     }
85020cf1dd8SToby Isaac   }
8511dca8a05SBarry Smith   PetscCheck(!fem->Tf || k <= fem->Tf->K,PetscObjectComm((PetscObject) fem), PETSC_ERR_ARG_OUTOFRANGE, "Requested %" PetscInt_FMT " derivatives, but only tabulated %" PetscInt_FMT, k, fem->Tf->K);
852ef0bb6c7SMatthew G. Knepley   *Tf = fem->Tf;
85320cf1dd8SToby Isaac   PetscFunctionReturn(0);
85420cf1dd8SToby Isaac }
85520cf1dd8SToby Isaac 
8562b99622eSMatthew G. Knepley /*@C
857ef0bb6c7SMatthew G. Knepley   PetscFEGetFaceCentroidTabulation - Returns the tabulation of the basis functions at the face centroid points
8582b99622eSMatthew G. Knepley 
8592b99622eSMatthew G. Knepley   Not collective
8602b99622eSMatthew G. Knepley 
8612b99622eSMatthew G. Knepley   Input Parameter:
8622b99622eSMatthew G. Knepley . fem - The PetscFE object
8632b99622eSMatthew G. Knepley 
8642b99622eSMatthew G. Knepley   Output Parameters:
865ef0bb6c7SMatthew G. Knepley . Tc - The basis function values at face centroid points
8662b99622eSMatthew G. Knepley 
8672b99622eSMatthew G. Knepley   Note:
868ef0bb6c7SMatthew G. Knepley $ T->T[0] = Bf[(f*pdim + i)*Nc + c] is the value at point f for basis function i and component c
8692b99622eSMatthew G. Knepley 
8702b99622eSMatthew G. Knepley   Level: intermediate
8712b99622eSMatthew G. Knepley 
872db781477SPatrick Sanan .seealso: `PetscFEGetFaceTabulation()`, `PetscFEGetCellTabulation()`, `PetscFECreateTabulation()`, `PetscTabulationDestroy()`
8732b99622eSMatthew G. Knepley @*/
874ef0bb6c7SMatthew G. Knepley PetscErrorCode PetscFEGetFaceCentroidTabulation(PetscFE fem, PetscTabulation *Tc)
87520cf1dd8SToby Isaac {
87620cf1dd8SToby Isaac   PetscFunctionBegin;
87720cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
878ef0bb6c7SMatthew G. Knepley   PetscValidPointer(Tc, 2);
879ef0bb6c7SMatthew G. Knepley   if (!fem->Tc) {
88020cf1dd8SToby Isaac     PetscDualSpace  sp;
88120cf1dd8SToby Isaac     DM              dm;
88220cf1dd8SToby Isaac     const PetscInt *cone;
88320cf1dd8SToby Isaac     PetscReal      *centroids;
88420cf1dd8SToby Isaac     PetscInt        dim, numFaces, f;
88520cf1dd8SToby Isaac 
8869566063dSJacob Faibussowitsch     PetscCall(PetscFEGetDualSpace(fem, &sp));
8879566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDM(sp, &dm));
8889566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(dm, &dim));
8899566063dSJacob Faibussowitsch     PetscCall(DMPlexGetConeSize(dm, 0, &numFaces));
8909566063dSJacob Faibussowitsch     PetscCall(DMPlexGetCone(dm, 0, &cone));
8919566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(numFaces*dim, &centroids));
8929566063dSJacob Faibussowitsch     for (f = 0; f < numFaces; ++f) PetscCall(DMPlexComputeCellGeometryFVM(dm, cone[f], NULL, &centroids[f*dim], NULL));
8939566063dSJacob Faibussowitsch     PetscCall(PetscFECreateTabulation(fem, 1, numFaces, centroids, 0, &fem->Tc));
8949566063dSJacob Faibussowitsch     PetscCall(PetscFree(centroids));
89520cf1dd8SToby Isaac   }
896ef0bb6c7SMatthew G. Knepley   *Tc = fem->Tc;
89720cf1dd8SToby Isaac   PetscFunctionReturn(0);
89820cf1dd8SToby Isaac }
89920cf1dd8SToby Isaac 
90020cf1dd8SToby Isaac /*@C
901ef0bb6c7SMatthew G. Knepley   PetscFECreateTabulation - Tabulates the basis functions, and perhaps derivatives, at the points provided.
90220cf1dd8SToby Isaac 
90320cf1dd8SToby Isaac   Not collective
90420cf1dd8SToby Isaac 
90520cf1dd8SToby Isaac   Input Parameters:
90620cf1dd8SToby Isaac + fem     - The PetscFE object
907ef0bb6c7SMatthew G. Knepley . nrepl   - The number of replicas
908ef0bb6c7SMatthew G. Knepley . npoints - The number of tabulation points in a replica
909ef0bb6c7SMatthew G. Knepley . points  - The tabulation point coordinates
910ef0bb6c7SMatthew G. Knepley - K       - The number of derivatives calculated
91120cf1dd8SToby Isaac 
912ef0bb6c7SMatthew G. Knepley   Output Parameter:
913ef0bb6c7SMatthew G. Knepley . T - The basis function values and derivatives at tabulation points
91420cf1dd8SToby Isaac 
91520cf1dd8SToby Isaac   Note:
916ef0bb6c7SMatthew G. Knepley $ T->T[0] = B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c
917ef0bb6c7SMatthew G. Knepley $ T->T[1] = D[((p*pdim + i)*Nc + c)*dim + d] is the derivative value at point p for basis function i, component c, in direction d
918ef0bb6c7SMatthew G. Knepley $ T->T[2] = 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
91920cf1dd8SToby Isaac 
92020cf1dd8SToby Isaac   Level: intermediate
92120cf1dd8SToby Isaac 
922db781477SPatrick Sanan .seealso: `PetscFEGetCellTabulation()`, `PetscTabulationDestroy()`
92320cf1dd8SToby Isaac @*/
924ef0bb6c7SMatthew G. Knepley PetscErrorCode PetscFECreateTabulation(PetscFE fem, PetscInt nrepl, PetscInt npoints, const PetscReal points[], PetscInt K, PetscTabulation *T)
92520cf1dd8SToby Isaac {
92620cf1dd8SToby Isaac   DM               dm;
927ef0bb6c7SMatthew G. Knepley   PetscDualSpace   Q;
928ef0bb6c7SMatthew G. Knepley   PetscInt         Nb;   /* Dimension of FE space P */
929ef0bb6c7SMatthew G. Knepley   PetscInt         Nc;   /* Field components */
930ef0bb6c7SMatthew G. Knepley   PetscInt         cdim; /* Reference coordinate dimension */
931ef0bb6c7SMatthew G. Knepley   PetscInt         k;
93220cf1dd8SToby Isaac 
93320cf1dd8SToby Isaac   PetscFunctionBegin;
934ef0bb6c7SMatthew G. Knepley   if (!npoints || !fem->dualSpace || K < 0) {
935ef0bb6c7SMatthew G. Knepley     *T = NULL;
93620cf1dd8SToby Isaac     PetscFunctionReturn(0);
93720cf1dd8SToby Isaac   }
93820cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
939dadcf809SJacob Faibussowitsch   PetscValidRealPointer(points, 4);
94040a2aa30SMatthew G. Knepley   PetscValidPointer(T, 6);
9419566063dSJacob Faibussowitsch   PetscCall(PetscFEGetDualSpace(fem, &Q));
9429566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDM(Q, &dm));
9439566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &cdim));
9449566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDimension(Q, &Nb));
9459566063dSJacob Faibussowitsch   PetscCall(PetscFEGetNumComponents(fem, &Nc));
9469566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(1, T));
947ef0bb6c7SMatthew G. Knepley   (*T)->K    = !cdim ? 0 : K;
948ef0bb6c7SMatthew G. Knepley   (*T)->Nr   = nrepl;
949ef0bb6c7SMatthew G. Knepley   (*T)->Np   = npoints;
950ef0bb6c7SMatthew G. Knepley   (*T)->Nb   = Nb;
951ef0bb6c7SMatthew G. Knepley   (*T)->Nc   = Nc;
952ef0bb6c7SMatthew G. Knepley   (*T)->cdim = cdim;
9539566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1((*T)->K+1, &(*T)->T));
954ef0bb6c7SMatthew G. Knepley   for (k = 0; k <= (*T)->K; ++k) {
9559566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nrepl*npoints*Nb*Nc*PetscPowInt(cdim, k), &(*T)->T[k]));
95620cf1dd8SToby Isaac   }
957*dbbe0bcdSBarry Smith   PetscUseTypeMethod(fem,createtabulation , nrepl*npoints, points, K, *T);
95820cf1dd8SToby Isaac   PetscFunctionReturn(0);
95920cf1dd8SToby Isaac }
96020cf1dd8SToby Isaac 
9612b99622eSMatthew G. Knepley /*@C
962ef0bb6c7SMatthew G. Knepley   PetscFEComputeTabulation - Tabulates the basis functions, and perhaps derivatives, at the points provided.
9632b99622eSMatthew G. Knepley 
9642b99622eSMatthew G. Knepley   Not collective
9652b99622eSMatthew G. Knepley 
9662b99622eSMatthew G. Knepley   Input Parameters:
9672b99622eSMatthew G. Knepley + fem     - The PetscFE object
9682b99622eSMatthew G. Knepley . npoints - The number of tabulation points
9692b99622eSMatthew G. Knepley . points  - The tabulation point coordinates
970ef0bb6c7SMatthew G. Knepley . K       - The number of derivatives calculated
971ef0bb6c7SMatthew G. Knepley - T       - An existing tabulation object with enough allocated space
972ef0bb6c7SMatthew G. Knepley 
973ef0bb6c7SMatthew G. Knepley   Output Parameter:
974ef0bb6c7SMatthew G. Knepley . T - The basis function values and derivatives at tabulation points
9752b99622eSMatthew G. Knepley 
9762b99622eSMatthew G. Knepley   Note:
977ef0bb6c7SMatthew G. Knepley $ T->T[0] = B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c
978ef0bb6c7SMatthew G. Knepley $ T->T[1] = D[((p*pdim + i)*Nc + c)*dim + d] is the derivative value at point p for basis function i, component c, in direction d
979ef0bb6c7SMatthew G. Knepley $ T->T[2] = 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
9802b99622eSMatthew G. Knepley 
9812b99622eSMatthew G. Knepley   Level: intermediate
9822b99622eSMatthew G. Knepley 
983db781477SPatrick Sanan .seealso: `PetscFEGetCellTabulation()`, `PetscTabulationDestroy()`
9842b99622eSMatthew G. Knepley @*/
985ef0bb6c7SMatthew G. Knepley PetscErrorCode PetscFEComputeTabulation(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscInt K, PetscTabulation T)
986ef0bb6c7SMatthew G. Knepley {
987ef0bb6c7SMatthew G. Knepley   PetscFunctionBeginHot;
988ef0bb6c7SMatthew G. Knepley   if (!npoints || !fem->dualSpace || K < 0) PetscFunctionReturn(0);
989ef0bb6c7SMatthew G. Knepley   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
990dadcf809SJacob Faibussowitsch   PetscValidRealPointer(points, 3);
991ef0bb6c7SMatthew G. Knepley   PetscValidPointer(T, 5);
99276bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {
99320cf1dd8SToby Isaac     DM               dm;
994ef0bb6c7SMatthew G. Knepley     PetscDualSpace   Q;
995ef0bb6c7SMatthew G. Knepley     PetscInt         Nb;   /* Dimension of FE space P */
996ef0bb6c7SMatthew G. Knepley     PetscInt         Nc;   /* Field components */
997ef0bb6c7SMatthew G. Knepley     PetscInt         cdim; /* Reference coordinate dimension */
998ef0bb6c7SMatthew G. Knepley 
9999566063dSJacob Faibussowitsch     PetscCall(PetscFEGetDualSpace(fem, &Q));
10009566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDM(Q, &dm));
10019566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(dm, &cdim));
10029566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetDimension(Q, &Nb));
10039566063dSJacob Faibussowitsch     PetscCall(PetscFEGetNumComponents(fem, &Nc));
100463a3b9bcSJacob Faibussowitsch     PetscCheck(T->K    == (!cdim ? 0 : K),PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Tabulation K %" PetscInt_FMT " must match requested K %" PetscInt_FMT, T->K, !cdim ? 0 : K);
100563a3b9bcSJacob Faibussowitsch     PetscCheck(T->Nb   == Nb,PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Tabulation Nb %" PetscInt_FMT " must match requested Nb %" PetscInt_FMT, T->Nb, Nb);
100663a3b9bcSJacob Faibussowitsch     PetscCheck(T->Nc   == Nc,PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Tabulation Nc %" PetscInt_FMT " must match requested Nc %" PetscInt_FMT, T->Nc, Nc);
100763a3b9bcSJacob Faibussowitsch     PetscCheck(T->cdim == cdim,PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Tabulation cdim %" PetscInt_FMT " must match requested cdim %" PetscInt_FMT, T->cdim, cdim);
1008ef0bb6c7SMatthew G. Knepley   }
1009ef0bb6c7SMatthew G. Knepley   T->Nr = 1;
1010ef0bb6c7SMatthew G. Knepley   T->Np = npoints;
1011*dbbe0bcdSBarry Smith   PetscUseTypeMethod(fem,createtabulation , npoints, points, K, T);
1012ef0bb6c7SMatthew G. Knepley   PetscFunctionReturn(0);
1013ef0bb6c7SMatthew G. Knepley }
1014ef0bb6c7SMatthew G. Knepley 
1015ef0bb6c7SMatthew G. Knepley /*@C
1016ef0bb6c7SMatthew G. Knepley   PetscTabulationDestroy - Frees memory from the associated tabulation.
1017ef0bb6c7SMatthew G. Knepley 
1018ef0bb6c7SMatthew G. Knepley   Not collective
1019ef0bb6c7SMatthew G. Knepley 
1020ef0bb6c7SMatthew G. Knepley   Input Parameter:
1021ef0bb6c7SMatthew G. Knepley . T - The tabulation
1022ef0bb6c7SMatthew G. Knepley 
1023ef0bb6c7SMatthew G. Knepley   Level: intermediate
1024ef0bb6c7SMatthew G. Knepley 
1025db781477SPatrick Sanan .seealso: `PetscFECreateTabulation()`, `PetscFEGetCellTabulation()`
1026ef0bb6c7SMatthew G. Knepley @*/
1027ef0bb6c7SMatthew G. Knepley PetscErrorCode PetscTabulationDestroy(PetscTabulation *T)
1028ef0bb6c7SMatthew G. Knepley {
1029ef0bb6c7SMatthew G. Knepley   PetscInt       k;
103020cf1dd8SToby Isaac 
103120cf1dd8SToby Isaac   PetscFunctionBegin;
1032ef0bb6c7SMatthew G. Knepley   PetscValidPointer(T, 1);
1033ef0bb6c7SMatthew G. Knepley   if (!T || !(*T)) PetscFunctionReturn(0);
10349566063dSJacob Faibussowitsch   for (k = 0; k <= (*T)->K; ++k) PetscCall(PetscFree((*T)->T[k]));
10359566063dSJacob Faibussowitsch   PetscCall(PetscFree((*T)->T));
10369566063dSJacob Faibussowitsch   PetscCall(PetscFree(*T));
1037ef0bb6c7SMatthew G. Knepley   *T = NULL;
103820cf1dd8SToby Isaac   PetscFunctionReturn(0);
103920cf1dd8SToby Isaac }
104020cf1dd8SToby Isaac 
104120cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreatePointTrace(PetscFE fe, PetscInt refPoint, PetscFE *trFE)
104220cf1dd8SToby Isaac {
104320cf1dd8SToby Isaac   PetscSpace     bsp, bsubsp;
104420cf1dd8SToby Isaac   PetscDualSpace dsp, dsubsp;
104520cf1dd8SToby Isaac   PetscInt       dim, depth, numComp, i, j, coneSize, order;
104620cf1dd8SToby Isaac   PetscFEType    type;
104720cf1dd8SToby Isaac   DM             dm;
104820cf1dd8SToby Isaac   DMLabel        label;
104920cf1dd8SToby Isaac   PetscReal      *xi, *v, *J, detJ;
1050db11e2ebSMatthew G. Knepley   const char     *name;
105120cf1dd8SToby Isaac   PetscQuadrature origin, fullQuad, subQuad;
105220cf1dd8SToby Isaac 
105320cf1dd8SToby Isaac   PetscFunctionBegin;
105420cf1dd8SToby Isaac   PetscValidHeaderSpecific(fe,PETSCFE_CLASSID,1);
105520cf1dd8SToby Isaac   PetscValidPointer(trFE,3);
10569566063dSJacob Faibussowitsch   PetscCall(PetscFEGetBasisSpace(fe,&bsp));
10579566063dSJacob Faibussowitsch   PetscCall(PetscFEGetDualSpace(fe,&dsp));
10589566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDM(dsp,&dm));
10599566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm,&dim));
10609566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthLabel(dm,&label));
10619566063dSJacob Faibussowitsch   PetscCall(DMLabelGetValue(label,refPoint,&depth));
10629566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(depth,&xi));
10639566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(dim,&v));
10649566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(dim*dim,&J));
106520cf1dd8SToby Isaac   for (i = 0; i < depth; i++) xi[i] = 0.;
10669566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF,&origin));
10679566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureSetData(origin,depth,0,1,xi,NULL));
10689566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeCellGeometryFEM(dm,refPoint,origin,v,J,NULL,&detJ));
106920cf1dd8SToby Isaac   /* CellGeometryFEM computes the expanded Jacobian, we want the true jacobian */
107020cf1dd8SToby Isaac   for (i = 1; i < dim; i++) {
107120cf1dd8SToby Isaac     for (j = 0; j < depth; j++) {
107220cf1dd8SToby Isaac       J[i * depth + j] = J[i * dim + j];
107320cf1dd8SToby Isaac     }
107420cf1dd8SToby Isaac   }
10759566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureDestroy(&origin));
10769566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetPointSubspace(dsp,refPoint,&dsubsp));
10779566063dSJacob Faibussowitsch   PetscCall(PetscSpaceCreateSubspace(bsp,dsubsp,v,J,NULL,NULL,PETSC_OWN_POINTER,&bsubsp));
10789566063dSJacob Faibussowitsch   PetscCall(PetscSpaceSetUp(bsubsp));
10799566063dSJacob Faibussowitsch   PetscCall(PetscFECreate(PetscObjectComm((PetscObject)fe),trFE));
10809566063dSJacob Faibussowitsch   PetscCall(PetscFEGetType(fe,&type));
10819566063dSJacob Faibussowitsch   PetscCall(PetscFESetType(*trFE,type));
10829566063dSJacob Faibussowitsch   PetscCall(PetscFEGetNumComponents(fe,&numComp));
10839566063dSJacob Faibussowitsch   PetscCall(PetscFESetNumComponents(*trFE,numComp));
10849566063dSJacob Faibussowitsch   PetscCall(PetscFESetBasisSpace(*trFE,bsubsp));
10859566063dSJacob Faibussowitsch   PetscCall(PetscFESetDualSpace(*trFE,dsubsp));
10869566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetName((PetscObject) fe, &name));
10879566063dSJacob Faibussowitsch   if (name) PetscCall(PetscFESetName(*trFE, name));
10889566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(fe,&fullQuad));
10899566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetOrder(fullQuad,&order));
10909566063dSJacob Faibussowitsch   PetscCall(DMPlexGetConeSize(dm,refPoint,&coneSize));
10911baa6e33SBarry Smith   if (coneSize == 2 * depth) PetscCall(PetscDTGaussTensorQuadrature(depth,1,(order + 1)/2,-1.,1.,&subQuad));
10921baa6e33SBarry Smith   else PetscCall(PetscDTStroudConicalQuadrature(depth,1,(order + 1)/2,-1.,1.,&subQuad));
10939566063dSJacob Faibussowitsch   PetscCall(PetscFESetQuadrature(*trFE,subQuad));
10949566063dSJacob Faibussowitsch   PetscCall(PetscFESetUp(*trFE));
10959566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureDestroy(&subQuad));
10969566063dSJacob Faibussowitsch   PetscCall(PetscSpaceDestroy(&bsubsp));
109720cf1dd8SToby Isaac   PetscFunctionReturn(0);
109820cf1dd8SToby Isaac }
109920cf1dd8SToby Isaac 
110020cf1dd8SToby Isaac PetscErrorCode PetscFECreateHeightTrace(PetscFE fe, PetscInt height, PetscFE *trFE)
110120cf1dd8SToby Isaac {
110220cf1dd8SToby Isaac   PetscInt       hStart, hEnd;
110320cf1dd8SToby Isaac   PetscDualSpace dsp;
110420cf1dd8SToby Isaac   DM             dm;
110520cf1dd8SToby Isaac 
110620cf1dd8SToby Isaac   PetscFunctionBegin;
110720cf1dd8SToby Isaac   PetscValidHeaderSpecific(fe,PETSCFE_CLASSID,1);
110820cf1dd8SToby Isaac   PetscValidPointer(trFE,3);
110920cf1dd8SToby Isaac   *trFE = NULL;
11109566063dSJacob Faibussowitsch   PetscCall(PetscFEGetDualSpace(fe,&dsp));
11119566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDM(dsp,&dm));
11129566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm,height,&hStart,&hEnd));
111320cf1dd8SToby Isaac   if (hEnd <= hStart) PetscFunctionReturn(0);
11149566063dSJacob Faibussowitsch   PetscCall(PetscFECreatePointTrace(fe,hStart,trFE));
111520cf1dd8SToby Isaac   PetscFunctionReturn(0);
111620cf1dd8SToby Isaac }
111720cf1dd8SToby Isaac 
111820cf1dd8SToby Isaac /*@
111920cf1dd8SToby Isaac   PetscFEGetDimension - Get the dimension of the finite element space on a cell
112020cf1dd8SToby Isaac 
112120cf1dd8SToby Isaac   Not collective
112220cf1dd8SToby Isaac 
112320cf1dd8SToby Isaac   Input Parameter:
112420cf1dd8SToby Isaac . fe - The PetscFE
112520cf1dd8SToby Isaac 
112620cf1dd8SToby Isaac   Output Parameter:
112720cf1dd8SToby Isaac . dim - The dimension
112820cf1dd8SToby Isaac 
112920cf1dd8SToby Isaac   Level: intermediate
113020cf1dd8SToby Isaac 
1131db781477SPatrick Sanan .seealso: `PetscFECreate()`, `PetscSpaceGetDimension()`, `PetscDualSpaceGetDimension()`
113220cf1dd8SToby Isaac @*/
113320cf1dd8SToby Isaac PetscErrorCode PetscFEGetDimension(PetscFE fem, PetscInt *dim)
113420cf1dd8SToby Isaac {
113520cf1dd8SToby Isaac   PetscFunctionBegin;
113620cf1dd8SToby Isaac   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
1137dadcf809SJacob Faibussowitsch   PetscValidIntPointer(dim, 2);
1138*dbbe0bcdSBarry Smith   PetscTryTypeMethod(fem,getdimension , dim);
113920cf1dd8SToby Isaac   PetscFunctionReturn(0);
114020cf1dd8SToby Isaac }
114120cf1dd8SToby Isaac 
11424bee2e38SMatthew G. Knepley /*@C
11434bee2e38SMatthew G. Knepley   PetscFEPushforward - Map the reference element function to real space
11444bee2e38SMatthew G. Knepley 
11454bee2e38SMatthew G. Knepley   Input Parameters:
11464bee2e38SMatthew G. Knepley + fe     - The PetscFE
11474bee2e38SMatthew G. Knepley . fegeom - The cell geometry
11484bee2e38SMatthew G. Knepley . Nv     - The number of function values
11494bee2e38SMatthew G. Knepley - vals   - The function values
11504bee2e38SMatthew G. Knepley 
11514bee2e38SMatthew G. Knepley   Output Parameter:
11524bee2e38SMatthew G. Knepley . vals   - The transformed function values
11534bee2e38SMatthew G. Knepley 
11544bee2e38SMatthew G. Knepley   Level: advanced
11554bee2e38SMatthew G. Knepley 
11564bee2e38SMatthew G. Knepley   Note: This just forwards the call onto PetscDualSpacePushforward().
11574bee2e38SMatthew G. Knepley 
1158f9244615SMatthew G. Knepley   Note: This only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension.
11592edcad52SToby Isaac 
1160db781477SPatrick Sanan .seealso: `PetscDualSpacePushforward()`
11614bee2e38SMatthew G. Knepley @*/
11622edcad52SToby Isaac PetscErrorCode PetscFEPushforward(PetscFE fe, PetscFEGeom *fegeom, PetscInt Nv, PetscScalar vals[])
11634bee2e38SMatthew G. Knepley {
11642ae266adSMatthew G. Knepley   PetscFunctionBeginHot;
11659566063dSJacob Faibussowitsch   PetscCall(PetscDualSpacePushforward(fe->dualSpace, fegeom, Nv, fe->numComponents, vals));
11664bee2e38SMatthew G. Knepley   PetscFunctionReturn(0);
11674bee2e38SMatthew G. Knepley }
11684bee2e38SMatthew G. Knepley 
11694bee2e38SMatthew G. Knepley /*@C
11704bee2e38SMatthew G. Knepley   PetscFEPushforwardGradient - Map the reference element function gradient to real space
11714bee2e38SMatthew G. Knepley 
11724bee2e38SMatthew G. Knepley   Input Parameters:
11734bee2e38SMatthew G. Knepley + fe     - The PetscFE
11744bee2e38SMatthew G. Knepley . fegeom - The cell geometry
11754bee2e38SMatthew G. Knepley . Nv     - The number of function gradient values
11764bee2e38SMatthew G. Knepley - vals   - The function gradient values
11774bee2e38SMatthew G. Knepley 
11784bee2e38SMatthew G. Knepley   Output Parameter:
11794bee2e38SMatthew G. Knepley . vals   - The transformed function gradient values
11804bee2e38SMatthew G. Knepley 
11814bee2e38SMatthew G. Knepley   Level: advanced
11824bee2e38SMatthew G. Knepley 
11834bee2e38SMatthew G. Knepley   Note: This just forwards the call onto PetscDualSpacePushforwardGradient().
11844bee2e38SMatthew G. Knepley 
1185f9244615SMatthew G. Knepley   Note: This only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension.
11862edcad52SToby Isaac 
1187db781477SPatrick Sanan .seealso: `PetscFEPushforward()`, `PetscDualSpacePushforwardGradient()`, `PetscDualSpacePushforward()`
11884bee2e38SMatthew G. Knepley @*/
11892edcad52SToby Isaac PetscErrorCode PetscFEPushforwardGradient(PetscFE fe, PetscFEGeom *fegeom, PetscInt Nv, PetscScalar vals[])
11904bee2e38SMatthew G. Knepley {
11912ae266adSMatthew G. Knepley   PetscFunctionBeginHot;
11929566063dSJacob Faibussowitsch   PetscCall(PetscDualSpacePushforwardGradient(fe->dualSpace, fegeom, Nv, fe->numComponents, vals));
11934bee2e38SMatthew G. Knepley   PetscFunctionReturn(0);
11944bee2e38SMatthew G. Knepley }
11954bee2e38SMatthew G. Knepley 
1196f9244615SMatthew G. Knepley /*@C
1197f9244615SMatthew G. Knepley   PetscFEPushforwardHessian - Map the reference element function Hessian to real space
1198f9244615SMatthew G. Knepley 
1199f9244615SMatthew G. Knepley   Input Parameters:
1200f9244615SMatthew G. Knepley + fe     - The PetscFE
1201f9244615SMatthew G. Knepley . fegeom - The cell geometry
1202f9244615SMatthew G. Knepley . Nv     - The number of function Hessian values
1203f9244615SMatthew G. Knepley - vals   - The function Hessian values
1204f9244615SMatthew G. Knepley 
1205f9244615SMatthew G. Knepley   Output Parameter:
1206f9244615SMatthew G. Knepley . vals   - The transformed function Hessian values
1207f9244615SMatthew G. Knepley 
1208f9244615SMatthew G. Knepley   Level: advanced
1209f9244615SMatthew G. Knepley 
1210f9244615SMatthew G. Knepley   Note: This just forwards the call onto PetscDualSpacePushforwardHessian().
1211f9244615SMatthew G. Knepley 
1212f9244615SMatthew G. Knepley   Note: This only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension.
1213f9244615SMatthew G. Knepley 
1214db781477SPatrick Sanan .seealso: `PetscFEPushforward()`, `PetscDualSpacePushforwardHessian()`, `PetscDualSpacePushforward()`
1215f9244615SMatthew G. Knepley @*/
1216f9244615SMatthew G. Knepley PetscErrorCode PetscFEPushforwardHessian(PetscFE fe, PetscFEGeom *fegeom, PetscInt Nv, PetscScalar vals[])
1217f9244615SMatthew G. Knepley {
1218f9244615SMatthew G. Knepley   PetscFunctionBeginHot;
12199566063dSJacob Faibussowitsch   PetscCall(PetscDualSpacePushforwardHessian(fe->dualSpace, fegeom, Nv, fe->numComponents, vals));
1220f9244615SMatthew G. Knepley   PetscFunctionReturn(0);
1221f9244615SMatthew G. Knepley }
1222f9244615SMatthew G. Knepley 
122320cf1dd8SToby Isaac /*
122420cf1dd8SToby Isaac Purpose: Compute element vector for chunk of elements
122520cf1dd8SToby Isaac 
122620cf1dd8SToby Isaac Input:
122720cf1dd8SToby Isaac   Sizes:
122820cf1dd8SToby Isaac      Ne:  number of elements
122920cf1dd8SToby Isaac      Nf:  number of fields
123020cf1dd8SToby Isaac      PetscFE
123120cf1dd8SToby Isaac        dim: spatial dimension
123220cf1dd8SToby Isaac        Nb:  number of basis functions
123320cf1dd8SToby Isaac        Nc:  number of field components
123420cf1dd8SToby Isaac        PetscQuadrature
123520cf1dd8SToby Isaac          Nq:  number of quadrature points
123620cf1dd8SToby Isaac 
123720cf1dd8SToby Isaac   Geometry:
123820cf1dd8SToby Isaac      PetscFEGeom[Ne] possibly *Nq
123920cf1dd8SToby Isaac        PetscReal v0s[dim]
124020cf1dd8SToby Isaac        PetscReal n[dim]
124120cf1dd8SToby Isaac        PetscReal jacobians[dim*dim]
124220cf1dd8SToby Isaac        PetscReal jacobianInverses[dim*dim]
124320cf1dd8SToby Isaac        PetscReal jacobianDeterminants
124420cf1dd8SToby Isaac   FEM:
124520cf1dd8SToby Isaac      PetscFE
124620cf1dd8SToby Isaac        PetscQuadrature
124720cf1dd8SToby Isaac          PetscReal   quadPoints[Nq*dim]
124820cf1dd8SToby Isaac          PetscReal   quadWeights[Nq]
124920cf1dd8SToby Isaac        PetscReal   basis[Nq*Nb*Nc]
125020cf1dd8SToby Isaac        PetscReal   basisDer[Nq*Nb*Nc*dim]
125120cf1dd8SToby Isaac      PetscScalar coefficients[Ne*Nb*Nc]
125220cf1dd8SToby Isaac      PetscScalar elemVec[Ne*Nb*Nc]
125320cf1dd8SToby Isaac 
125420cf1dd8SToby Isaac   Problem:
125520cf1dd8SToby Isaac      PetscInt f: the active field
125620cf1dd8SToby Isaac      f0, f1
125720cf1dd8SToby Isaac 
125820cf1dd8SToby Isaac   Work Space:
125920cf1dd8SToby Isaac      PetscFE
126020cf1dd8SToby Isaac        PetscScalar f0[Nq*dim];
126120cf1dd8SToby Isaac        PetscScalar f1[Nq*dim*dim];
126220cf1dd8SToby Isaac        PetscScalar u[Nc];
126320cf1dd8SToby Isaac        PetscScalar gradU[Nc*dim];
126420cf1dd8SToby Isaac        PetscReal   x[dim];
126520cf1dd8SToby Isaac        PetscScalar realSpaceDer[dim];
126620cf1dd8SToby Isaac 
126720cf1dd8SToby Isaac Purpose: Compute element vector for N_cb batches of elements
126820cf1dd8SToby Isaac 
126920cf1dd8SToby Isaac Input:
127020cf1dd8SToby Isaac   Sizes:
127120cf1dd8SToby Isaac      N_cb: Number of serial cell batches
127220cf1dd8SToby Isaac 
127320cf1dd8SToby Isaac   Geometry:
127420cf1dd8SToby Isaac      PetscReal v0s[Ne*dim]
127520cf1dd8SToby Isaac      PetscReal jacobians[Ne*dim*dim]        possibly *Nq
127620cf1dd8SToby Isaac      PetscReal jacobianInverses[Ne*dim*dim] possibly *Nq
127720cf1dd8SToby Isaac      PetscReal jacobianDeterminants[Ne]     possibly *Nq
127820cf1dd8SToby Isaac   FEM:
127920cf1dd8SToby Isaac      static PetscReal   quadPoints[Nq*dim]
128020cf1dd8SToby Isaac      static PetscReal   quadWeights[Nq]
128120cf1dd8SToby Isaac      static PetscReal   basis[Nq*Nb*Nc]
128220cf1dd8SToby Isaac      static PetscReal   basisDer[Nq*Nb*Nc*dim]
128320cf1dd8SToby Isaac      PetscScalar coefficients[Ne*Nb*Nc]
128420cf1dd8SToby Isaac      PetscScalar elemVec[Ne*Nb*Nc]
128520cf1dd8SToby Isaac 
128620cf1dd8SToby Isaac ex62.c:
128720cf1dd8SToby Isaac   PetscErrorCode PetscFEIntegrateResidualBatch(PetscInt Ne, PetscInt numFields, PetscInt field, PetscQuadrature quad[], const PetscScalar coefficients[],
128820cf1dd8SToby Isaac                                                const PetscReal v0s[], const PetscReal jacobians[], const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[],
128920cf1dd8SToby Isaac                                                void (*f0_func)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar f0[]),
129020cf1dd8SToby Isaac                                                void (*f1_func)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar f1[]), PetscScalar elemVec[])
129120cf1dd8SToby Isaac 
129220cf1dd8SToby Isaac ex52.c:
129320cf1dd8SToby 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)
129420cf1dd8SToby 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)
129520cf1dd8SToby Isaac 
129620cf1dd8SToby Isaac ex52_integrateElement.cu
129720cf1dd8SToby Isaac __global__ void integrateElementQuadrature(int N_cb, realType *coefficients, realType *jacobianInverses, realType *jacobianDeterminants, realType *elemVec)
129820cf1dd8SToby Isaac 
129920cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode IntegrateElementBatchGPU(PetscInt spatial_dim, PetscInt Ne, PetscInt Ncb, PetscInt Nbc, PetscInt Nbl, const PetscScalar coefficients[],
130020cf1dd8SToby Isaac                                                      const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], PetscScalar elemVec[],
130120cf1dd8SToby Isaac                                                      PetscLogEvent event, PetscInt debug, PetscInt pde_op)
130220cf1dd8SToby Isaac 
130320cf1dd8SToby Isaac ex52_integrateElementOpenCL.c:
130420cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode IntegrateElementBatchGPU(PetscInt spatial_dim, PetscInt Ne, PetscInt Ncb, PetscInt Nbc, PetscInt N_bl, const PetscScalar coefficients[],
130520cf1dd8SToby Isaac                                                      const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], PetscScalar elemVec[],
130620cf1dd8SToby Isaac                                                      PetscLogEvent event, PetscInt debug, PetscInt pde_op)
130720cf1dd8SToby Isaac 
130820cf1dd8SToby Isaac __kernel void integrateElementQuadrature(int N_cb, __global float *coefficients, __global float *jacobianInverses, __global float *jacobianDeterminants, __global float *elemVec)
130920cf1dd8SToby Isaac */
131020cf1dd8SToby Isaac 
131120cf1dd8SToby Isaac /*@C
131220cf1dd8SToby Isaac   PetscFEIntegrate - Produce the integral for the given field for a chunk of elements by quadrature integration
131320cf1dd8SToby Isaac 
131420cf1dd8SToby Isaac   Not collective
131520cf1dd8SToby Isaac 
131620cf1dd8SToby Isaac   Input Parameters:
1317360cf244SMatthew G. Knepley + prob         - The PetscDS specifying the discretizations and continuum functions
131820cf1dd8SToby Isaac . field        - The field being integrated
131920cf1dd8SToby Isaac . Ne           - The number of elements in the chunk
132020cf1dd8SToby Isaac . cgeom        - The cell geometry for each cell in the chunk
132120cf1dd8SToby Isaac . coefficients - The array of FEM basis coefficients for the elements
132220cf1dd8SToby Isaac . probAux      - The PetscDS specifying the auxiliary discretizations
132320cf1dd8SToby Isaac - coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
132420cf1dd8SToby Isaac 
13257a7aea1fSJed Brown   Output Parameter:
132620cf1dd8SToby Isaac . integral     - the integral for this field
132720cf1dd8SToby Isaac 
13282b99622eSMatthew G. Knepley   Level: intermediate
132920cf1dd8SToby Isaac 
1330db781477SPatrick Sanan .seealso: `PetscFEIntegrateResidual()`
133120cf1dd8SToby Isaac @*/
13324bee2e38SMatthew G. Knepley PetscErrorCode PetscFEIntegrate(PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom,
133320cf1dd8SToby Isaac                                 const PetscScalar coefficients[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscScalar integral[])
133420cf1dd8SToby Isaac {
13354bee2e38SMatthew G. Knepley   PetscFE        fe;
133620cf1dd8SToby Isaac 
133720cf1dd8SToby Isaac   PetscFunctionBegin;
13384bee2e38SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
13399566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(prob, field, (PetscObject *) &fe));
13409566063dSJacob Faibussowitsch   if (fe->ops->integrate) PetscCall((*fe->ops->integrate)(prob, field, Ne, cgeom, coefficients, probAux, coefficientsAux, integral));
134120cf1dd8SToby Isaac   PetscFunctionReturn(0);
134220cf1dd8SToby Isaac }
134320cf1dd8SToby Isaac 
134420cf1dd8SToby Isaac /*@C
1345afe6d6adSToby Isaac   PetscFEIntegrateBd - Produce the integral for the given field for a chunk of elements by quadrature integration
1346afe6d6adSToby Isaac 
1347afe6d6adSToby Isaac   Not collective
1348afe6d6adSToby Isaac 
1349afe6d6adSToby Isaac   Input Parameters:
1350360cf244SMatthew G. Knepley + prob         - The PetscDS specifying the discretizations and continuum functions
1351afe6d6adSToby Isaac . field        - The field being integrated
1352afe6d6adSToby Isaac . obj_func     - The function to be integrated
1353afe6d6adSToby Isaac . Ne           - The number of elements in the chunk
1354afe6d6adSToby Isaac . fgeom        - The face geometry for each face in the chunk
1355afe6d6adSToby Isaac . coefficients - The array of FEM basis coefficients for the elements
1356afe6d6adSToby Isaac . probAux      - The PetscDS specifying the auxiliary discretizations
1357afe6d6adSToby Isaac - coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
1358afe6d6adSToby Isaac 
13597a7aea1fSJed Brown   Output Parameter:
1360afe6d6adSToby Isaac . integral     - the integral for this field
1361afe6d6adSToby Isaac 
13622b99622eSMatthew G. Knepley   Level: intermediate
1363afe6d6adSToby Isaac 
1364db781477SPatrick Sanan .seealso: `PetscFEIntegrateResidual()`
1365afe6d6adSToby Isaac @*/
13664bee2e38SMatthew G. Knepley PetscErrorCode PetscFEIntegrateBd(PetscDS prob, PetscInt field,
1367afe6d6adSToby Isaac                                   void (*obj_func)(PetscInt, PetscInt, PetscInt,
1368afe6d6adSToby Isaac                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1369afe6d6adSToby Isaac                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1370afe6d6adSToby Isaac                                                    PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1371afe6d6adSToby Isaac                                   PetscInt Ne, PetscFEGeom *geom, const PetscScalar coefficients[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscScalar integral[])
1372afe6d6adSToby Isaac {
13734bee2e38SMatthew G. Knepley   PetscFE        fe;
1374afe6d6adSToby Isaac 
1375afe6d6adSToby Isaac   PetscFunctionBegin;
13764bee2e38SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
13779566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(prob, field, (PetscObject *) &fe));
13789566063dSJacob Faibussowitsch   if (fe->ops->integratebd) PetscCall((*fe->ops->integratebd)(prob, field, obj_func, Ne, geom, coefficients, probAux, coefficientsAux, integral));
1379afe6d6adSToby Isaac   PetscFunctionReturn(0);
1380afe6d6adSToby Isaac }
1381afe6d6adSToby Isaac 
1382afe6d6adSToby Isaac /*@C
138320cf1dd8SToby Isaac   PetscFEIntegrateResidual - Produce the element residual vector for a chunk of elements by quadrature integration
138420cf1dd8SToby Isaac 
138520cf1dd8SToby Isaac   Not collective
138620cf1dd8SToby Isaac 
138720cf1dd8SToby Isaac   Input Parameters:
13886528b96dSMatthew G. Knepley + ds           - The PetscDS specifying the discretizations and continuum functions
13896528b96dSMatthew G. Knepley . key          - The (label+value, field) being integrated
139020cf1dd8SToby Isaac . Ne           - The number of elements in the chunk
139120cf1dd8SToby Isaac . cgeom        - The cell geometry for each cell in the chunk
139220cf1dd8SToby Isaac . coefficients - The array of FEM basis coefficients for the elements
139320cf1dd8SToby Isaac . coefficients_t - The array of FEM basis time derivative coefficients for the elements
139420cf1dd8SToby Isaac . probAux      - The PetscDS specifying the auxiliary discretizations
139520cf1dd8SToby Isaac . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
139620cf1dd8SToby Isaac - t            - The time
139720cf1dd8SToby Isaac 
13987a7aea1fSJed Brown   Output Parameter:
139920cf1dd8SToby Isaac . elemVec      - the element residual vectors from each element
140020cf1dd8SToby Isaac 
140120cf1dd8SToby Isaac   Note:
140220cf1dd8SToby Isaac $ Loop over batch of elements (e):
140320cf1dd8SToby Isaac $   Loop over quadrature points (q):
140420cf1dd8SToby Isaac $     Make u_q and gradU_q (loops over fields,Nb,Ncomp) and x_q
140520cf1dd8SToby Isaac $     Call f_0 and f_1
140620cf1dd8SToby Isaac $   Loop over element vector entries (f,fc --> i):
140720cf1dd8SToby Isaac $     elemVec[i] += \psi^{fc}_f(q) f0_{fc}(u, \nabla u) + \nabla\psi^{fc}_f(q) \cdot f1_{fc,df}(u, \nabla u)
140820cf1dd8SToby Isaac 
14092b99622eSMatthew G. Knepley   Level: intermediate
141020cf1dd8SToby Isaac 
1411db781477SPatrick Sanan .seealso: `PetscFEIntegrateResidual()`
141220cf1dd8SToby Isaac @*/
141306ad1575SMatthew G. Knepley PetscErrorCode PetscFEIntegrateResidual(PetscDS ds, PetscFormKey key, PetscInt Ne, PetscFEGeom *cgeom,
141420cf1dd8SToby Isaac                                         const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
141520cf1dd8SToby Isaac {
14164bee2e38SMatthew G. Knepley   PetscFE        fe;
141720cf1dd8SToby Isaac 
14186528b96dSMatthew G. Knepley   PetscFunctionBeginHot;
14196528b96dSMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
14209566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, key.field, (PetscObject *) &fe));
14219566063dSJacob Faibussowitsch   if (fe->ops->integrateresidual) PetscCall((*fe->ops->integrateresidual)(ds, key, Ne, cgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec));
142220cf1dd8SToby Isaac   PetscFunctionReturn(0);
142320cf1dd8SToby Isaac }
142420cf1dd8SToby Isaac 
142520cf1dd8SToby Isaac /*@C
142620cf1dd8SToby Isaac   PetscFEIntegrateBdResidual - Produce the element residual vector for a chunk of elements by quadrature integration over a boundary
142720cf1dd8SToby Isaac 
142820cf1dd8SToby Isaac   Not collective
142920cf1dd8SToby Isaac 
143020cf1dd8SToby Isaac   Input Parameters:
143106d8a0d3SMatthew G. Knepley + ds           - The PetscDS specifying the discretizations and continuum functions
143245480ffeSMatthew G. Knepley . wf           - The PetscWeakForm object holding the pointwise functions
143306d8a0d3SMatthew G. Knepley . key          - The (label+value, field) being integrated
143420cf1dd8SToby Isaac . Ne           - The number of elements in the chunk
143520cf1dd8SToby Isaac . fgeom        - The face geometry for each cell in the chunk
143620cf1dd8SToby Isaac . coefficients - The array of FEM basis coefficients for the elements
143720cf1dd8SToby Isaac . coefficients_t - The array of FEM basis time derivative coefficients for the elements
143820cf1dd8SToby Isaac . probAux      - The PetscDS specifying the auxiliary discretizations
143920cf1dd8SToby Isaac . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
144020cf1dd8SToby Isaac - t            - The time
144120cf1dd8SToby Isaac 
14427a7aea1fSJed Brown   Output Parameter:
144320cf1dd8SToby Isaac . elemVec      - the element residual vectors from each element
144420cf1dd8SToby Isaac 
14452b99622eSMatthew G. Knepley   Level: intermediate
144620cf1dd8SToby Isaac 
1447db781477SPatrick Sanan .seealso: `PetscFEIntegrateResidual()`
144820cf1dd8SToby Isaac @*/
144906ad1575SMatthew G. Knepley PetscErrorCode PetscFEIntegrateBdResidual(PetscDS ds, PetscWeakForm wf, PetscFormKey key, PetscInt Ne, PetscFEGeom *fgeom,
145020cf1dd8SToby Isaac                                           const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
145120cf1dd8SToby Isaac {
14524bee2e38SMatthew G. Knepley   PetscFE        fe;
145320cf1dd8SToby Isaac 
145420cf1dd8SToby Isaac   PetscFunctionBegin;
145506d8a0d3SMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
14569566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, key.field, (PetscObject *) &fe));
14579566063dSJacob Faibussowitsch   if (fe->ops->integratebdresidual) PetscCall((*fe->ops->integratebdresidual)(ds, wf, key, Ne, fgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec));
145820cf1dd8SToby Isaac   PetscFunctionReturn(0);
145920cf1dd8SToby Isaac }
146020cf1dd8SToby Isaac 
146120cf1dd8SToby Isaac /*@C
146227f02ce8SMatthew G. Knepley   PetscFEIntegrateHybridResidual - Produce the element residual vector for a chunk of hybrid element faces by quadrature integration
146327f02ce8SMatthew G. Knepley 
146427f02ce8SMatthew G. Knepley   Not collective
146527f02ce8SMatthew G. Knepley 
146627f02ce8SMatthew G. Knepley   Input Parameters:
146727f02ce8SMatthew G. Knepley + prob         - The PetscDS specifying the discretizations and continuum functions
14686528b96dSMatthew G. Knepley . key          - The (label+value, field) being integrated
1469c2b7495fSMatthew G. Knepley . s            - The side of the cell being integrated, 0 for negative and 1 for positive
147027f02ce8SMatthew G. Knepley . Ne           - The number of elements in the chunk
147127f02ce8SMatthew G. Knepley . fgeom        - The face geometry for each cell in the chunk
147227f02ce8SMatthew G. Knepley . coefficients - The array of FEM basis coefficients for the elements
147327f02ce8SMatthew G. Knepley . coefficients_t - The array of FEM basis time derivative coefficients for the elements
147427f02ce8SMatthew G. Knepley . probAux      - The PetscDS specifying the auxiliary discretizations
147527f02ce8SMatthew G. Knepley . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
147627f02ce8SMatthew G. Knepley - t            - The time
147727f02ce8SMatthew G. Knepley 
147827f02ce8SMatthew G. Knepley   Output Parameter
147927f02ce8SMatthew G. Knepley . elemVec      - the element residual vectors from each element
148027f02ce8SMatthew G. Knepley 
148127f02ce8SMatthew G. Knepley   Level: developer
148227f02ce8SMatthew G. Knepley 
1483db781477SPatrick Sanan .seealso: `PetscFEIntegrateResidual()`
148427f02ce8SMatthew G. Knepley @*/
1485c2b7495fSMatthew G. Knepley PetscErrorCode PetscFEIntegrateHybridResidual(PetscDS prob, PetscFormKey key, PetscInt s, PetscInt Ne, PetscFEGeom *fgeom,
148627f02ce8SMatthew G. Knepley                                               const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
148727f02ce8SMatthew G. Knepley {
148827f02ce8SMatthew G. Knepley   PetscFE        fe;
148927f02ce8SMatthew G. Knepley 
149027f02ce8SMatthew G. Knepley   PetscFunctionBegin;
149127f02ce8SMatthew G. Knepley   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
14929566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(prob, key.field, (PetscObject *) &fe));
14939566063dSJacob Faibussowitsch   if (fe->ops->integratehybridresidual) PetscCall((*fe->ops->integratehybridresidual)(prob, key, s, Ne, fgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec));
149427f02ce8SMatthew G. Knepley   PetscFunctionReturn(0);
149527f02ce8SMatthew G. Knepley }
149627f02ce8SMatthew G. Knepley 
149727f02ce8SMatthew G. Knepley /*@C
149820cf1dd8SToby Isaac   PetscFEIntegrateJacobian - Produce the element Jacobian for a chunk of elements by quadrature integration
149920cf1dd8SToby Isaac 
150020cf1dd8SToby Isaac   Not collective
150120cf1dd8SToby Isaac 
150220cf1dd8SToby Isaac   Input Parameters:
15036528b96dSMatthew G. Knepley + ds           - The PetscDS specifying the discretizations and continuum functions
150420cf1dd8SToby Isaac . jtype        - The type of matrix pointwise functions that should be used
15056528b96dSMatthew G. Knepley . key          - The (label+value, fieldI*Nf + fieldJ) being integrated
15065fedec97SMatthew G. Knepley . s            - The side of the cell being integrated, 0 for negative and 1 for positive
150720cf1dd8SToby Isaac . Ne           - The number of elements in the chunk
150820cf1dd8SToby Isaac . cgeom        - The cell geometry for each cell in the chunk
150920cf1dd8SToby Isaac . coefficients - The array of FEM basis coefficients for the elements for the Jacobian evaluation point
151020cf1dd8SToby Isaac . coefficients_t - The array of FEM basis time derivative coefficients for the elements
151120cf1dd8SToby Isaac . probAux      - The PetscDS specifying the auxiliary discretizations
151220cf1dd8SToby Isaac . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
151320cf1dd8SToby Isaac . t            - The time
151420cf1dd8SToby Isaac - u_tShift     - A multiplier for the dF/du_t term (as opposed to the dF/du term)
151520cf1dd8SToby Isaac 
15167a7aea1fSJed Brown   Output Parameter:
151720cf1dd8SToby Isaac . elemMat      - the element matrices for the Jacobian from each element
151820cf1dd8SToby Isaac 
151920cf1dd8SToby Isaac   Note:
152020cf1dd8SToby Isaac $ Loop over batch of elements (e):
152120cf1dd8SToby Isaac $   Loop over element matrix entries (f,fc,g,gc --> i,j):
152220cf1dd8SToby Isaac $     Loop over quadrature points (q):
152320cf1dd8SToby Isaac $       Make u_q and gradU_q (loops over fields,Nb,Ncomp)
152420cf1dd8SToby Isaac $         elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q)
152520cf1dd8SToby Isaac $                      + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
152620cf1dd8SToby Isaac $                      + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q)
152720cf1dd8SToby Isaac $                      + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
15282b99622eSMatthew G. Knepley   Level: intermediate
152920cf1dd8SToby Isaac 
1530db781477SPatrick Sanan .seealso: `PetscFEIntegrateResidual()`
153120cf1dd8SToby Isaac @*/
153206ad1575SMatthew G. Knepley PetscErrorCode PetscFEIntegrateJacobian(PetscDS ds, PetscFEJacobianType jtype, PetscFormKey key, PetscInt Ne, PetscFEGeom *cgeom,
153320cf1dd8SToby Isaac                                         const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
153420cf1dd8SToby Isaac {
15354bee2e38SMatthew G. Knepley   PetscFE        fe;
15366528b96dSMatthew G. Knepley   PetscInt       Nf;
153720cf1dd8SToby Isaac 
153820cf1dd8SToby Isaac   PetscFunctionBegin;
15396528b96dSMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
15409566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
15419566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, key.field / Nf, (PetscObject *) &fe));
15429566063dSJacob Faibussowitsch   if (fe->ops->integratejacobian) PetscCall((*fe->ops->integratejacobian)(ds, jtype, key, Ne, cgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, u_tshift, elemMat));
154320cf1dd8SToby Isaac   PetscFunctionReturn(0);
154420cf1dd8SToby Isaac }
154520cf1dd8SToby Isaac 
154620cf1dd8SToby Isaac /*@C
154720cf1dd8SToby Isaac   PetscFEIntegrateBdJacobian - Produce the boundary element Jacobian for a chunk of elements by quadrature integration
154820cf1dd8SToby Isaac 
154920cf1dd8SToby Isaac   Not collective
155020cf1dd8SToby Isaac 
155120cf1dd8SToby Isaac   Input Parameters:
155245480ffeSMatthew G. Knepley + ds           - The PetscDS specifying the discretizations and continuum functions
155345480ffeSMatthew G. Knepley . wf           - The PetscWeakForm holding the pointwise functions
155445480ffeSMatthew G. Knepley . key          - The (label+value, fieldI*Nf + fieldJ) being integrated
155520cf1dd8SToby Isaac . Ne           - The number of elements in the chunk
155620cf1dd8SToby Isaac . fgeom        - The face geometry for each cell in the chunk
155720cf1dd8SToby Isaac . coefficients - The array of FEM basis coefficients for the elements for the Jacobian evaluation point
155820cf1dd8SToby Isaac . coefficients_t - The array of FEM basis time derivative coefficients for the elements
155920cf1dd8SToby Isaac . probAux      - The PetscDS specifying the auxiliary discretizations
156020cf1dd8SToby Isaac . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
156120cf1dd8SToby Isaac . t            - The time
156220cf1dd8SToby Isaac - u_tShift     - A multiplier for the dF/du_t term (as opposed to the dF/du term)
156320cf1dd8SToby Isaac 
15647a7aea1fSJed Brown   Output Parameter:
156520cf1dd8SToby Isaac . elemMat              - the element matrices for the Jacobian from each element
156620cf1dd8SToby Isaac 
156720cf1dd8SToby Isaac   Note:
156820cf1dd8SToby Isaac $ Loop over batch of elements (e):
156920cf1dd8SToby Isaac $   Loop over element matrix entries (f,fc,g,gc --> i,j):
157020cf1dd8SToby Isaac $     Loop over quadrature points (q):
157120cf1dd8SToby Isaac $       Make u_q and gradU_q (loops over fields,Nb,Ncomp)
157220cf1dd8SToby Isaac $         elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q)
157320cf1dd8SToby Isaac $                      + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
157420cf1dd8SToby Isaac $                      + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q)
157520cf1dd8SToby Isaac $                      + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
15762b99622eSMatthew G. Knepley   Level: intermediate
157720cf1dd8SToby Isaac 
1578db781477SPatrick Sanan .seealso: `PetscFEIntegrateJacobian()`, `PetscFEIntegrateResidual()`
157920cf1dd8SToby Isaac @*/
158006ad1575SMatthew G. Knepley PetscErrorCode PetscFEIntegrateBdJacobian(PetscDS ds, PetscWeakForm wf, PetscFormKey key, PetscInt Ne, PetscFEGeom *fgeom,
158120cf1dd8SToby Isaac                                           const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
158220cf1dd8SToby Isaac {
15834bee2e38SMatthew G. Knepley   PetscFE        fe;
158445480ffeSMatthew G. Knepley   PetscInt       Nf;
158520cf1dd8SToby Isaac 
158620cf1dd8SToby Isaac   PetscFunctionBegin;
158745480ffeSMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
15889566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
15899566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, key.field / Nf, (PetscObject *) &fe));
15909566063dSJacob Faibussowitsch   if (fe->ops->integratebdjacobian) PetscCall((*fe->ops->integratebdjacobian)(ds, wf, key, Ne, fgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, u_tshift, elemMat));
159120cf1dd8SToby Isaac   PetscFunctionReturn(0);
159220cf1dd8SToby Isaac }
159320cf1dd8SToby Isaac 
159427f02ce8SMatthew G. Knepley /*@C
159527f02ce8SMatthew G. Knepley   PetscFEIntegrateHybridJacobian - Produce the boundary element Jacobian for a chunk of hybrid elements by quadrature integration
159627f02ce8SMatthew G. Knepley 
159727f02ce8SMatthew G. Knepley   Not collective
159827f02ce8SMatthew G. Knepley 
159927f02ce8SMatthew G. Knepley   Input Parameters:
160045480ffeSMatthew G. Knepley + ds           - The PetscDS specifying the discretizations and continuum functions
160127f02ce8SMatthew G. Knepley . jtype        - The type of matrix pointwise functions that should be used
160245480ffeSMatthew G. Knepley . key          - The (label+value, fieldI*Nf + fieldJ) being integrated
16035fedec97SMatthew G. Knepley . s            - The side of the cell being integrated, 0 for negative and 1 for positive
160427f02ce8SMatthew G. Knepley . Ne           - The number of elements in the chunk
160527f02ce8SMatthew G. Knepley . fgeom        - The face geometry for each cell in the chunk
160627f02ce8SMatthew G. Knepley . coefficients - The array of FEM basis coefficients for the elements for the Jacobian evaluation point
160727f02ce8SMatthew G. Knepley . coefficients_t - The array of FEM basis time derivative coefficients for the elements
160827f02ce8SMatthew G. Knepley . probAux      - The PetscDS specifying the auxiliary discretizations
160927f02ce8SMatthew G. Knepley . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
161027f02ce8SMatthew G. Knepley . t            - The time
161127f02ce8SMatthew G. Knepley - u_tShift     - A multiplier for the dF/du_t term (as opposed to the dF/du term)
161227f02ce8SMatthew G. Knepley 
161327f02ce8SMatthew G. Knepley   Output Parameter
161427f02ce8SMatthew G. Knepley . elemMat              - the element matrices for the Jacobian from each element
161527f02ce8SMatthew G. Knepley 
161627f02ce8SMatthew G. Knepley   Note:
161727f02ce8SMatthew G. Knepley $ Loop over batch of elements (e):
161827f02ce8SMatthew G. Knepley $   Loop over element matrix entries (f,fc,g,gc --> i,j):
161927f02ce8SMatthew G. Knepley $     Loop over quadrature points (q):
162027f02ce8SMatthew G. Knepley $       Make u_q and gradU_q (loops over fields,Nb,Ncomp)
162127f02ce8SMatthew G. Knepley $         elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q)
162227f02ce8SMatthew G. Knepley $                      + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
162327f02ce8SMatthew G. Knepley $                      + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q)
162427f02ce8SMatthew G. Knepley $                      + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
162527f02ce8SMatthew G. Knepley   Level: developer
162627f02ce8SMatthew G. Knepley 
1627db781477SPatrick Sanan .seealso: `PetscFEIntegrateJacobian()`, `PetscFEIntegrateResidual()`
162827f02ce8SMatthew G. Knepley @*/
16295fedec97SMatthew G. Knepley PetscErrorCode PetscFEIntegrateHybridJacobian(PetscDS ds, PetscFEJacobianType jtype, PetscFormKey key, PetscInt s, PetscInt Ne, PetscFEGeom *fgeom,
163027f02ce8SMatthew G. Knepley                                               const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
163127f02ce8SMatthew G. Knepley {
163227f02ce8SMatthew G. Knepley   PetscFE        fe;
163345480ffeSMatthew G. Knepley   PetscInt       Nf;
163427f02ce8SMatthew G. Knepley 
163527f02ce8SMatthew G. Knepley   PetscFunctionBegin;
163645480ffeSMatthew G. Knepley   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
16379566063dSJacob Faibussowitsch   PetscCall(PetscDSGetNumFields(ds, &Nf));
16389566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(ds, key.field / Nf, (PetscObject *) &fe));
16399566063dSJacob Faibussowitsch   if (fe->ops->integratehybridjacobian) PetscCall((*fe->ops->integratehybridjacobian)(ds, jtype, key, s, Ne, fgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, u_tshift, elemMat));
164027f02ce8SMatthew G. Knepley   PetscFunctionReturn(0);
164127f02ce8SMatthew G. Knepley }
164227f02ce8SMatthew G. Knepley 
16432b99622eSMatthew G. Knepley /*@
16442b99622eSMatthew G. Knepley   PetscFEGetHeightSubspace - Get the subspace of this space for a mesh point of a given height
16452b99622eSMatthew G. Knepley 
16462b99622eSMatthew G. Knepley   Input Parameters:
16472b99622eSMatthew G. Knepley + fe     - The finite element space
16482b99622eSMatthew G. Knepley - height - The height of the Plex point
16492b99622eSMatthew G. Knepley 
16502b99622eSMatthew G. Knepley   Output Parameter:
16512b99622eSMatthew G. Knepley . subfe  - The subspace of this FE space
16522b99622eSMatthew G. Knepley 
16532b99622eSMatthew G. Knepley   Note: For example, if we want the subspace of this space for a face, we would choose height = 1.
16542b99622eSMatthew G. Knepley 
16552b99622eSMatthew G. Knepley   Level: advanced
16562b99622eSMatthew G. Knepley 
1657db781477SPatrick Sanan .seealso: `PetscFECreateDefault()`
16582b99622eSMatthew G. Knepley @*/
165920cf1dd8SToby Isaac PetscErrorCode PetscFEGetHeightSubspace(PetscFE fe, PetscInt height, PetscFE *subfe)
166020cf1dd8SToby Isaac {
166120cf1dd8SToby Isaac   PetscSpace      P, subP;
166220cf1dd8SToby Isaac   PetscDualSpace  Q, subQ;
166320cf1dd8SToby Isaac   PetscQuadrature subq;
166420cf1dd8SToby Isaac   PetscFEType     fetype;
166520cf1dd8SToby Isaac   PetscInt        dim, Nc;
166620cf1dd8SToby Isaac 
166720cf1dd8SToby Isaac   PetscFunctionBegin;
166820cf1dd8SToby Isaac   PetscValidHeaderSpecific(fe, PETSCFE_CLASSID, 1);
166920cf1dd8SToby Isaac   PetscValidPointer(subfe, 3);
167020cf1dd8SToby Isaac   if (height == 0) {
167120cf1dd8SToby Isaac     *subfe = fe;
167220cf1dd8SToby Isaac     PetscFunctionReturn(0);
167320cf1dd8SToby Isaac   }
16749566063dSJacob Faibussowitsch   PetscCall(PetscFEGetBasisSpace(fe, &P));
16759566063dSJacob Faibussowitsch   PetscCall(PetscFEGetDualSpace(fe, &Q));
16769566063dSJacob Faibussowitsch   PetscCall(PetscFEGetNumComponents(fe, &Nc));
16779566063dSJacob Faibussowitsch   PetscCall(PetscFEGetFaceQuadrature(fe, &subq));
16789566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDimension(Q, &dim));
16791dca8a05SBarry Smith   PetscCheck(height <= dim && height >= 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Asked for space at height %" PetscInt_FMT " for dimension %" PetscInt_FMT " space", height, dim);
16809566063dSJacob Faibussowitsch   if (!fe->subspaces) PetscCall(PetscCalloc1(dim, &fe->subspaces));
168120cf1dd8SToby Isaac   if (height <= dim) {
168220cf1dd8SToby Isaac     if (!fe->subspaces[height-1]) {
1683665f567fSMatthew G. Knepley       PetscFE     sub = NULL;
16843f6b16c7SMatthew G. Knepley       const char *name;
168520cf1dd8SToby Isaac 
16869566063dSJacob Faibussowitsch       PetscCall(PetscSpaceGetHeightSubspace(P, height, &subP));
16879566063dSJacob Faibussowitsch       PetscCall(PetscDualSpaceGetHeightSubspace(Q, height, &subQ));
1688665f567fSMatthew G. Knepley       if (subQ) {
16899566063dSJacob Faibussowitsch         PetscCall(PetscFECreate(PetscObjectComm((PetscObject) fe), &sub));
16909566063dSJacob Faibussowitsch         PetscCall(PetscObjectGetName((PetscObject) fe,  &name));
16919566063dSJacob Faibussowitsch         PetscCall(PetscObjectSetName((PetscObject) sub,  name));
16929566063dSJacob Faibussowitsch         PetscCall(PetscFEGetType(fe, &fetype));
16939566063dSJacob Faibussowitsch         PetscCall(PetscFESetType(sub, fetype));
16949566063dSJacob Faibussowitsch         PetscCall(PetscFESetBasisSpace(sub, subP));
16959566063dSJacob Faibussowitsch         PetscCall(PetscFESetDualSpace(sub, subQ));
16969566063dSJacob Faibussowitsch         PetscCall(PetscFESetNumComponents(sub, Nc));
16979566063dSJacob Faibussowitsch         PetscCall(PetscFESetUp(sub));
16989566063dSJacob Faibussowitsch         PetscCall(PetscFESetQuadrature(sub, subq));
1699665f567fSMatthew G. Knepley       }
170020cf1dd8SToby Isaac       fe->subspaces[height-1] = sub;
170120cf1dd8SToby Isaac     }
170220cf1dd8SToby Isaac     *subfe = fe->subspaces[height-1];
170320cf1dd8SToby Isaac   } else {
170420cf1dd8SToby Isaac     *subfe = NULL;
170520cf1dd8SToby Isaac   }
170620cf1dd8SToby Isaac   PetscFunctionReturn(0);
170720cf1dd8SToby Isaac }
170820cf1dd8SToby Isaac 
170920cf1dd8SToby Isaac /*@
171020cf1dd8SToby Isaac   PetscFERefine - Create a "refined" PetscFE object that refines the reference cell into smaller copies. This is typically used
171120cf1dd8SToby Isaac   to precondition a higher order method with a lower order method on a refined mesh having the same number of dofs (but more
171220cf1dd8SToby Isaac   sparsity). It is also used to create an interpolation between regularly refined meshes.
171320cf1dd8SToby Isaac 
1714d083f849SBarry Smith   Collective on fem
171520cf1dd8SToby Isaac 
171620cf1dd8SToby Isaac   Input Parameter:
171720cf1dd8SToby Isaac . fe - The initial PetscFE
171820cf1dd8SToby Isaac 
171920cf1dd8SToby Isaac   Output Parameter:
172020cf1dd8SToby Isaac . feRef - The refined PetscFE
172120cf1dd8SToby Isaac 
17222b99622eSMatthew G. Knepley   Level: advanced
172320cf1dd8SToby Isaac 
1724db781477SPatrick Sanan .seealso: `PetscFEType`, `PetscFECreate()`, `PetscFESetType()`
172520cf1dd8SToby Isaac @*/
172620cf1dd8SToby Isaac PetscErrorCode PetscFERefine(PetscFE fe, PetscFE *feRef)
172720cf1dd8SToby Isaac {
172820cf1dd8SToby Isaac   PetscSpace       P, Pref;
172920cf1dd8SToby Isaac   PetscDualSpace   Q, Qref;
173020cf1dd8SToby Isaac   DM               K, Kref;
173120cf1dd8SToby Isaac   PetscQuadrature  q, qref;
173220cf1dd8SToby Isaac   const PetscReal *v0, *jac;
173320cf1dd8SToby Isaac   PetscInt         numComp, numSubelements;
17341ac17e89SToby Isaac   PetscInt         cStart, cEnd, c;
17351ac17e89SToby Isaac   PetscDualSpace  *cellSpaces;
173620cf1dd8SToby Isaac 
173720cf1dd8SToby Isaac   PetscFunctionBegin;
17389566063dSJacob Faibussowitsch   PetscCall(PetscFEGetBasisSpace(fe, &P));
17399566063dSJacob Faibussowitsch   PetscCall(PetscFEGetDualSpace(fe, &Q));
17409566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(fe, &q));
17419566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDM(Q, &K));
174220cf1dd8SToby Isaac   /* Create space */
17439566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject) P));
174420cf1dd8SToby Isaac   Pref = P;
174520cf1dd8SToby Isaac   /* Create dual space */
17469566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceDuplicate(Q, &Qref));
17479566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSetType(Qref, PETSCDUALSPACEREFINED));
17489566063dSJacob Faibussowitsch   PetscCall(DMRefine(K, PetscObjectComm((PetscObject) fe), &Kref));
17499566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSetDM(Qref, Kref));
17509566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(Kref, 0, &cStart, &cEnd));
17519566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(cEnd - cStart, &cellSpaces));
17521ac17e89SToby Isaac   /* TODO: fix for non-uniform refinement */
17531ac17e89SToby Isaac   for (c = 0; c < cEnd - cStart; c++) cellSpaces[c] = Q;
17549566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceRefinedSetCellSpaces(Qref, cellSpaces));
17559566063dSJacob Faibussowitsch   PetscCall(PetscFree(cellSpaces));
17569566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&Kref));
17579566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSetUp(Qref));
175820cf1dd8SToby Isaac   /* Create element */
17599566063dSJacob Faibussowitsch   PetscCall(PetscFECreate(PetscObjectComm((PetscObject) fe), feRef));
17609566063dSJacob Faibussowitsch   PetscCall(PetscFESetType(*feRef, PETSCFECOMPOSITE));
17619566063dSJacob Faibussowitsch   PetscCall(PetscFESetBasisSpace(*feRef, Pref));
17629566063dSJacob Faibussowitsch   PetscCall(PetscFESetDualSpace(*feRef, Qref));
17639566063dSJacob Faibussowitsch   PetscCall(PetscFEGetNumComponents(fe,    &numComp));
17649566063dSJacob Faibussowitsch   PetscCall(PetscFESetNumComponents(*feRef, numComp));
17659566063dSJacob Faibussowitsch   PetscCall(PetscFESetUp(*feRef));
17669566063dSJacob Faibussowitsch   PetscCall(PetscSpaceDestroy(&Pref));
17679566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceDestroy(&Qref));
176820cf1dd8SToby Isaac   /* Create quadrature */
17699566063dSJacob Faibussowitsch   PetscCall(PetscFECompositeGetMapping(*feRef, &numSubelements, &v0, &jac, NULL));
17709566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureExpandComposite(q, numSubelements, v0, jac, &qref));
17719566063dSJacob Faibussowitsch   PetscCall(PetscFESetQuadrature(*feRef, qref));
17729566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureDestroy(&qref));
177320cf1dd8SToby Isaac   PetscFunctionReturn(0);
177420cf1dd8SToby Isaac }
177520cf1dd8SToby Isaac 
17767c48043bSMatthew G. Knepley static PetscErrorCode PetscFESetDefaultName_Private(PetscFE fe)
17777c48043bSMatthew G. Knepley {
17787c48043bSMatthew G. Knepley   PetscSpace     P;
17797c48043bSMatthew G. Knepley   PetscDualSpace Q;
17807c48043bSMatthew G. Knepley   DM             K;
17817c48043bSMatthew G. Knepley   DMPolytopeType ct;
17827c48043bSMatthew G. Knepley   PetscInt       degree;
17837c48043bSMatthew G. Knepley   char           name[64];
17847c48043bSMatthew G. Knepley 
17857c48043bSMatthew G. Knepley   PetscFunctionBegin;
17867c48043bSMatthew G. Knepley   PetscCall(PetscFEGetBasisSpace(fe, &P));
17877c48043bSMatthew G. Knepley   PetscCall(PetscSpaceGetDegree(P, &degree, NULL));
17887c48043bSMatthew G. Knepley   PetscCall(PetscFEGetDualSpace(fe, &Q));
17897c48043bSMatthew G. Knepley   PetscCall(PetscDualSpaceGetDM(Q, &K));
17907c48043bSMatthew G. Knepley   PetscCall(DMPlexGetCellType(K, 0, &ct));
17917c48043bSMatthew G. Knepley   switch (ct) {
17927c48043bSMatthew G. Knepley     case DM_POLYTOPE_SEGMENT:
17937c48043bSMatthew G. Knepley     case DM_POLYTOPE_POINT_PRISM_TENSOR:
17947c48043bSMatthew G. Knepley     case DM_POLYTOPE_QUADRILATERAL:
17957c48043bSMatthew G. Knepley     case DM_POLYTOPE_SEG_PRISM_TENSOR:
17967c48043bSMatthew G. Knepley     case DM_POLYTOPE_HEXAHEDRON:
17977c48043bSMatthew G. Knepley     case DM_POLYTOPE_QUAD_PRISM_TENSOR:
17987c48043bSMatthew G. Knepley       PetscCall(PetscSNPrintf(name, sizeof(name), "Q%" PetscInt_FMT, degree));
17997c48043bSMatthew G. Knepley       break;
18007c48043bSMatthew G. Knepley     case DM_POLYTOPE_TRIANGLE:
18017c48043bSMatthew G. Knepley     case DM_POLYTOPE_TETRAHEDRON:
18027c48043bSMatthew G. Knepley       PetscCall(PetscSNPrintf(name, sizeof(name), "P%" PetscInt_FMT, degree));
18037c48043bSMatthew G. Knepley       break;
18047c48043bSMatthew G. Knepley     case DM_POLYTOPE_TRI_PRISM:
18057c48043bSMatthew G. Knepley     case DM_POLYTOPE_TRI_PRISM_TENSOR:
18067c48043bSMatthew G. Knepley       PetscCall(PetscSNPrintf(name, sizeof(name), "P%" PetscInt_FMT "xQ%" PetscInt_FMT, degree, degree));
18077c48043bSMatthew G. Knepley       break;
18087c48043bSMatthew G. Knepley     default:
18097c48043bSMatthew G. Knepley       PetscCall(PetscSNPrintf(name, sizeof(name), "FE"));
18107c48043bSMatthew G. Knepley   }
18117c48043bSMatthew G. Knepley   PetscCall(PetscFESetName(fe, name));
18127c48043bSMatthew G. Knepley   PetscFunctionReturn(0);
18137c48043bSMatthew G. Knepley }
18147c48043bSMatthew G. Knepley 
18157c48043bSMatthew G. Knepley static PetscErrorCode PetscFECreateDefaultQuadrature_Private(PetscInt dim, DMPolytopeType ct, PetscInt qorder, PetscQuadrature *q, PetscQuadrature *fq)
18167c48043bSMatthew G. Knepley {
18177c48043bSMatthew G. Knepley   const PetscInt quadPointsPerEdge = PetscMax(qorder + 1, 1);
18187c48043bSMatthew G. Knepley 
18197c48043bSMatthew G. Knepley   PetscFunctionBegin;
18207c48043bSMatthew G. Knepley   switch (ct) {
18217c48043bSMatthew G. Knepley     case DM_POLYTOPE_SEGMENT:
18227c48043bSMatthew G. Knepley     case DM_POLYTOPE_POINT_PRISM_TENSOR:
18237c48043bSMatthew G. Knepley     case DM_POLYTOPE_QUADRILATERAL:
18247c48043bSMatthew G. Knepley     case DM_POLYTOPE_SEG_PRISM_TENSOR:
18257c48043bSMatthew G. Knepley     case DM_POLYTOPE_HEXAHEDRON:
18267c48043bSMatthew G. Knepley     case DM_POLYTOPE_QUAD_PRISM_TENSOR:
18277c48043bSMatthew G. Knepley       PetscCall(PetscDTGaussTensorQuadrature(dim,   1, quadPointsPerEdge, -1.0, 1.0, q));
18287c48043bSMatthew G. Knepley       PetscCall(PetscDTGaussTensorQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0, fq));
18297c48043bSMatthew G. Knepley       break;
18307c48043bSMatthew G. Knepley     case DM_POLYTOPE_TRIANGLE:
18317c48043bSMatthew G. Knepley     case DM_POLYTOPE_TETRAHEDRON:
18327c48043bSMatthew G. Knepley       PetscCall(PetscDTStroudConicalQuadrature(dim,   1, quadPointsPerEdge, -1.0, 1.0, q));
18337c48043bSMatthew G. Knepley       PetscCall(PetscDTStroudConicalQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0, fq));
18347c48043bSMatthew G. Knepley       break;
18357c48043bSMatthew G. Knepley     case DM_POLYTOPE_TRI_PRISM:
18367c48043bSMatthew G. Knepley     case DM_POLYTOPE_TRI_PRISM_TENSOR:
18377c48043bSMatthew G. Knepley       {
18387c48043bSMatthew G. Knepley         PetscQuadrature q1, q2;
18397c48043bSMatthew G. Knepley 
18407c48043bSMatthew G. Knepley         PetscCall(PetscDTStroudConicalQuadrature(2, 1, quadPointsPerEdge, -1.0, 1.0, &q1));
18417c48043bSMatthew G. Knepley         PetscCall(PetscDTGaussTensorQuadrature(1, 1, quadPointsPerEdge, -1.0, 1.0, &q2));
18427c48043bSMatthew G. Knepley         PetscCall(PetscDTTensorQuadratureCreate(q1, q2, q));
18437c48043bSMatthew G. Knepley         PetscCall(PetscQuadratureDestroy(&q1));
18447c48043bSMatthew G. Knepley         PetscCall(PetscQuadratureDestroy(&q2));
18457c48043bSMatthew G. Knepley       }
18467c48043bSMatthew G. Knepley       PetscCall(PetscDTStroudConicalQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0, fq));
18477c48043bSMatthew G. Knepley       /* TODO Need separate quadratures for each face */
18487c48043bSMatthew G. Knepley       break;
18497c48043bSMatthew G. Knepley     default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No quadrature for celltype %s", DMPolytopeTypes[PetscMin(ct, DM_POLYTOPE_UNKNOWN)]);
18507c48043bSMatthew G. Knepley   }
18517c48043bSMatthew G. Knepley   PetscFunctionReturn(0);
18527c48043bSMatthew G. Knepley }
18537c48043bSMatthew G. Knepley 
18547c48043bSMatthew G. Knepley /*@
18557c48043bSMatthew G. Knepley   PetscFECreateFromSpaces - Create a PetscFE from the basis and dual spaces
18567c48043bSMatthew G. Knepley 
18577c48043bSMatthew G. Knepley   Collective
18587c48043bSMatthew G. Knepley 
18597c48043bSMatthew G. Knepley   Input Parameters:
18607c48043bSMatthew G. Knepley + P  - The basis space
18617c48043bSMatthew G. Knepley . Q  - The dual space
18627c48043bSMatthew G. Knepley . q  - The cell quadrature
18637c48043bSMatthew G. Knepley - fq - The face quadrature
18647c48043bSMatthew G. Knepley 
18657c48043bSMatthew G. Knepley   Output Parameter:
18667c48043bSMatthew G. Knepley . fem    - The PetscFE object
18677c48043bSMatthew G. Knepley 
18687c48043bSMatthew G. Knepley   Note:
18697c48043bSMatthew G. Knepley   The PetscFE takes ownership of these spaces by calling destroy on each. They should not be used after this call, and for borrowed references from `PetscFEGetSpace()` and the like, the caller must use `PetscObjectReference` before this call.
18707c48043bSMatthew G. Knepley 
18717c48043bSMatthew G. Knepley   Level: beginner
18727c48043bSMatthew G. Knepley 
18737c48043bSMatthew G. Knepley .seealso: `PetscFECreateLagrangeByCell()`, `PetscFECreateDefault()`, `PetscFECreateByCell()`, `PetscFECreate()`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()`
18747c48043bSMatthew G. Knepley @*/
18757c48043bSMatthew G. Knepley PetscErrorCode PetscFECreateFromSpaces(PetscSpace P, PetscDualSpace Q, PetscQuadrature q, PetscQuadrature fq, PetscFE *fem)
18767c48043bSMatthew G. Knepley {
18777c48043bSMatthew G. Knepley   PetscInt    Nc;
18787c48043bSMatthew G. Knepley   const char *prefix;
18797c48043bSMatthew G. Knepley 
18807c48043bSMatthew G. Knepley   PetscFunctionBegin;
18817c48043bSMatthew G. Knepley   PetscCall(PetscFECreate(PetscObjectComm((PetscObject) P), fem));
18827c48043bSMatthew G. Knepley   PetscCall(PetscObjectGetOptionsPrefix((PetscObject) P, &prefix));
18837c48043bSMatthew G. Knepley   PetscCall(PetscObjectSetOptionsPrefix((PetscObject) *fem, prefix));
18847c48043bSMatthew G. Knepley   PetscCall(PetscFESetType(*fem, PETSCFEBASIC));
18857c48043bSMatthew G. Knepley   PetscCall(PetscFESetBasisSpace(*fem, P));
18867c48043bSMatthew G. Knepley   PetscCall(PetscFESetDualSpace(*fem, Q));
18877c48043bSMatthew G. Knepley   PetscCall(PetscSpaceGetNumComponents(P, &Nc));
18887c48043bSMatthew G. Knepley   PetscCall(PetscFESetNumComponents(*fem, Nc));
18897c48043bSMatthew G. Knepley   PetscCall(PetscFESetUp(*fem));
18907c48043bSMatthew G. Knepley   PetscCall(PetscSpaceDestroy(&P));
18917c48043bSMatthew G. Knepley   PetscCall(PetscDualSpaceDestroy(&Q));
18927c48043bSMatthew G. Knepley   PetscCall(PetscFESetQuadrature(*fem, q));
18937c48043bSMatthew G. Knepley   PetscCall(PetscFESetFaceQuadrature(*fem, fq));
18947c48043bSMatthew G. Knepley   PetscCall(PetscQuadratureDestroy(&q));
18957c48043bSMatthew G. Knepley   PetscCall(PetscQuadratureDestroy(&fq));
18967c48043bSMatthew G. Knepley   PetscCall(PetscFESetDefaultName_Private(*fem));
18977c48043bSMatthew G. Knepley   PetscFunctionReturn(0);
18987c48043bSMatthew G. Knepley }
18997c48043bSMatthew G. Knepley 
19002df84da0SMatthew G. Knepley static PetscErrorCode PetscFECreate_Internal(MPI_Comm comm, PetscInt dim, PetscInt Nc, DMPolytopeType ct, const char prefix[], PetscInt degree, PetscInt qorder, PetscBool setFromOptions, PetscFE *fem)
19012df84da0SMatthew G. Knepley {
19022df84da0SMatthew G. Knepley   DM              K;
19032df84da0SMatthew G. Knepley   PetscSpace      P;
19042df84da0SMatthew G. Knepley   PetscDualSpace  Q;
19057c48043bSMatthew G. Knepley   PetscQuadrature q, fq;
19062df84da0SMatthew G. Knepley   PetscBool       tensor;
19072df84da0SMatthew G. Knepley 
19082df84da0SMatthew G. Knepley   PetscFunctionBegin;
19092df84da0SMatthew G. Knepley   if (prefix) PetscValidCharPointer(prefix, 5);
19102df84da0SMatthew G. Knepley   PetscValidPointer(fem, 9);
19112df84da0SMatthew G. Knepley   switch (ct) {
19122df84da0SMatthew G. Knepley     case DM_POLYTOPE_SEGMENT:
19132df84da0SMatthew G. Knepley     case DM_POLYTOPE_POINT_PRISM_TENSOR:
19142df84da0SMatthew G. Knepley     case DM_POLYTOPE_QUADRILATERAL:
19152df84da0SMatthew G. Knepley     case DM_POLYTOPE_SEG_PRISM_TENSOR:
19162df84da0SMatthew G. Knepley     case DM_POLYTOPE_HEXAHEDRON:
19172df84da0SMatthew G. Knepley     case DM_POLYTOPE_QUAD_PRISM_TENSOR:
19182df84da0SMatthew G. Knepley       tensor = PETSC_TRUE;
19192df84da0SMatthew G. Knepley       break;
19202df84da0SMatthew G. Knepley     default: tensor = PETSC_FALSE;
19212df84da0SMatthew G. Knepley   }
19222df84da0SMatthew G. Knepley   /* Create space */
19239566063dSJacob Faibussowitsch   PetscCall(PetscSpaceCreate(comm, &P));
19249566063dSJacob Faibussowitsch   PetscCall(PetscSpaceSetType(P, PETSCSPACEPOLYNOMIAL));
19259566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject) P, prefix));
19269566063dSJacob Faibussowitsch   PetscCall(PetscSpacePolynomialSetTensor(P, tensor));
19279566063dSJacob Faibussowitsch   PetscCall(PetscSpaceSetNumComponents(P, Nc));
19289566063dSJacob Faibussowitsch   PetscCall(PetscSpaceSetNumVariables(P, dim));
19292df84da0SMatthew G. Knepley   if (degree >= 0) {
19309566063dSJacob Faibussowitsch     PetscCall(PetscSpaceSetDegree(P, degree, PETSC_DETERMINE));
1931cfd33b42SLisandro Dalcin     if (ct == DM_POLYTOPE_TRI_PRISM || ct == DM_POLYTOPE_TRI_PRISM_TENSOR) {
19322df84da0SMatthew G. Knepley       PetscSpace Pend, Pside;
19332df84da0SMatthew G. Knepley 
19349566063dSJacob Faibussowitsch       PetscCall(PetscSpaceCreate(comm, &Pend));
19359566063dSJacob Faibussowitsch       PetscCall(PetscSpaceSetType(Pend, PETSCSPACEPOLYNOMIAL));
19369566063dSJacob Faibussowitsch       PetscCall(PetscSpacePolynomialSetTensor(Pend, PETSC_FALSE));
19379566063dSJacob Faibussowitsch       PetscCall(PetscSpaceSetNumComponents(Pend, Nc));
19389566063dSJacob Faibussowitsch       PetscCall(PetscSpaceSetNumVariables(Pend, dim-1));
19399566063dSJacob Faibussowitsch       PetscCall(PetscSpaceSetDegree(Pend, degree, PETSC_DETERMINE));
19409566063dSJacob Faibussowitsch       PetscCall(PetscSpaceCreate(comm, &Pside));
19419566063dSJacob Faibussowitsch       PetscCall(PetscSpaceSetType(Pside, PETSCSPACEPOLYNOMIAL));
19429566063dSJacob Faibussowitsch       PetscCall(PetscSpacePolynomialSetTensor(Pside, PETSC_FALSE));
19439566063dSJacob Faibussowitsch       PetscCall(PetscSpaceSetNumComponents(Pside, 1));
19449566063dSJacob Faibussowitsch       PetscCall(PetscSpaceSetNumVariables(Pside, 1));
19459566063dSJacob Faibussowitsch       PetscCall(PetscSpaceSetDegree(Pside, degree, PETSC_DETERMINE));
19469566063dSJacob Faibussowitsch       PetscCall(PetscSpaceSetType(P, PETSCSPACETENSOR));
19479566063dSJacob Faibussowitsch       PetscCall(PetscSpaceTensorSetNumSubspaces(P, 2));
19489566063dSJacob Faibussowitsch       PetscCall(PetscSpaceTensorSetSubspace(P, 0, Pend));
19499566063dSJacob Faibussowitsch       PetscCall(PetscSpaceTensorSetSubspace(P, 1, Pside));
19509566063dSJacob Faibussowitsch       PetscCall(PetscSpaceDestroy(&Pend));
19519566063dSJacob Faibussowitsch       PetscCall(PetscSpaceDestroy(&Pside));
19522df84da0SMatthew G. Knepley     }
19532df84da0SMatthew G. Knepley   }
19549566063dSJacob Faibussowitsch   if (setFromOptions) PetscCall(PetscSpaceSetFromOptions(P));
19559566063dSJacob Faibussowitsch   PetscCall(PetscSpaceSetUp(P));
19569566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetDegree(P, &degree, NULL));
19579566063dSJacob Faibussowitsch   PetscCall(PetscSpacePolynomialGetTensor(P, &tensor));
19589566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetNumComponents(P, &Nc));
19592df84da0SMatthew G. Knepley   /* Create dual space */
19609566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceCreate(comm, &Q));
19619566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSetType(Q,PETSCDUALSPACELAGRANGE));
19629566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject) Q, prefix));
19639566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, ct, &K));
19649566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSetDM(Q, K));
19659566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&K));
19669566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSetNumComponents(Q, Nc));
19679566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSetOrder(Q, degree));
19682df84da0SMatthew G. Knepley   /* TODO For some reason, we need a tensor dualspace with wedges */
19699566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceLagrangeSetTensor(Q, (tensor || (ct == DM_POLYTOPE_TRI_PRISM)) ? PETSC_TRUE : PETSC_FALSE));
19709566063dSJacob Faibussowitsch   if (setFromOptions) PetscCall(PetscDualSpaceSetFromOptions(Q));
19719566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSetUp(Q));
19727c48043bSMatthew G. Knepley   /* Create quadrature */
19732df84da0SMatthew G. Knepley   qorder = qorder >= 0 ? qorder : degree;
19742df84da0SMatthew G. Knepley   if (setFromOptions) {
19757c48043bSMatthew G. Knepley     PetscObjectOptionsBegin((PetscObject) P);
19769566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBoundedInt("-petscfe_default_quadrature_order", "Quadrature order is one less than quadrature points per edge", "PetscFECreateDefault", qorder, &qorder, NULL, 0));
1977d0609cedSBarry Smith     PetscOptionsEnd();
19782df84da0SMatthew G. Knepley   }
19797c48043bSMatthew G. Knepley   PetscCall(PetscFECreateDefaultQuadrature_Private(dim, ct, qorder, &q, &fq));
19807c48043bSMatthew G. Knepley   /* Create finite element */
19817c48043bSMatthew G. Knepley   PetscCall(PetscFECreateFromSpaces(P, Q, q, fq, fem));
19827c48043bSMatthew G. Knepley   if (setFromOptions) PetscCall(PetscFESetFromOptions(*fem));
19832df84da0SMatthew G. Knepley   PetscFunctionReturn(0);
19842df84da0SMatthew G. Knepley }
19852df84da0SMatthew G. Knepley 
198620cf1dd8SToby Isaac /*@C
198720cf1dd8SToby Isaac   PetscFECreateDefault - Create a PetscFE for basic FEM computation
198820cf1dd8SToby Isaac 
1989d083f849SBarry Smith   Collective
199020cf1dd8SToby Isaac 
199120cf1dd8SToby Isaac   Input Parameters:
19927be5e748SToby Isaac + comm      - The MPI comm
199320cf1dd8SToby Isaac . dim       - The spatial dimension
199420cf1dd8SToby Isaac . Nc        - The number of components
199520cf1dd8SToby Isaac . isSimplex - Flag for simplex reference cell, otherwise its a tensor product
199620cf1dd8SToby Isaac . prefix    - The options prefix, or NULL
1997727cddd5SJacob Faibussowitsch - qorder    - The quadrature order or PETSC_DETERMINE to use PetscSpace polynomial degree
199820cf1dd8SToby Isaac 
199920cf1dd8SToby Isaac   Output Parameter:
200020cf1dd8SToby Isaac . fem - The PetscFE object
200120cf1dd8SToby Isaac 
2002e703855dSMatthew G. Knepley   Note:
20038f2aacc6SMatthew G. Knepley   Each subobject is SetFromOption() during creation, so that the object may be customized from the command line, using the prefix specified above. See the links below for the particular options available.
2004e703855dSMatthew G. Knepley 
200520cf1dd8SToby Isaac   Level: beginner
200620cf1dd8SToby Isaac 
2007db781477SPatrick Sanan .seealso: `PetscFECreateLagrange()`, `PetscFECreateByCell()`, `PetscSpaceSetFromOptions()`, `PetscDualSpaceSetFromOptions()`, `PetscFESetFromOptions()`, `PetscFECreate()`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()`
200820cf1dd8SToby Isaac @*/
20097be5e748SToby Isaac PetscErrorCode PetscFECreateDefault(MPI_Comm comm, PetscInt dim, PetscInt Nc, PetscBool isSimplex, const char prefix[], PetscInt qorder, PetscFE *fem)
201020cf1dd8SToby Isaac {
201120cf1dd8SToby Isaac   PetscFunctionBegin;
20129566063dSJacob Faibussowitsch   PetscCall(PetscFECreate_Internal(comm, dim, Nc, DMPolytopeTypeSimpleShape(dim, isSimplex), prefix, PETSC_DECIDE, qorder, PETSC_TRUE, fem));
20132df84da0SMatthew G. Knepley   PetscFunctionReturn(0);
201420cf1dd8SToby Isaac }
20152df84da0SMatthew G. Knepley 
20162df84da0SMatthew G. Knepley /*@C
20172df84da0SMatthew G. Knepley   PetscFECreateByCell - Create a PetscFE for basic FEM computation
20182df84da0SMatthew G. Knepley 
20192df84da0SMatthew G. Knepley   Collective
20202df84da0SMatthew G. Knepley 
20212df84da0SMatthew G. Knepley   Input Parameters:
20222df84da0SMatthew G. Knepley + comm   - The MPI comm
20232df84da0SMatthew G. Knepley . dim    - The spatial dimension
20242df84da0SMatthew G. Knepley . Nc     - The number of components
20252df84da0SMatthew G. Knepley . ct     - The celltype of the reference cell
20262df84da0SMatthew G. Knepley . prefix - The options prefix, or NULL
20272df84da0SMatthew G. Knepley - qorder - The quadrature order or PETSC_DETERMINE to use PetscSpace polynomial degree
20282df84da0SMatthew G. Knepley 
20292df84da0SMatthew G. Knepley   Output Parameter:
20302df84da0SMatthew G. Knepley . fem - The PetscFE object
20312df84da0SMatthew G. Knepley 
20322df84da0SMatthew G. Knepley   Note:
20332df84da0SMatthew G. Knepley   Each subobject is SetFromOption() during creation, so that the object may be customized from the command line, using the prefix specified above. See the links below for the particular options available.
20342df84da0SMatthew G. Knepley 
20352df84da0SMatthew G. Knepley   Level: beginner
20362df84da0SMatthew G. Knepley 
2037db781477SPatrick Sanan .seealso: `PetscFECreateDefault()`, `PetscFECreateLagrange()`, `PetscSpaceSetFromOptions()`, `PetscDualSpaceSetFromOptions()`, `PetscFESetFromOptions()`, `PetscFECreate()`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()`
20382df84da0SMatthew G. Knepley @*/
20392df84da0SMatthew G. Knepley PetscErrorCode PetscFECreateByCell(MPI_Comm comm, PetscInt dim, PetscInt Nc, DMPolytopeType ct, const char prefix[], PetscInt qorder, PetscFE *fem)
20402df84da0SMatthew G. Knepley {
20412df84da0SMatthew G. Knepley   PetscFunctionBegin;
20429566063dSJacob Faibussowitsch   PetscCall(PetscFECreate_Internal(comm, dim, Nc, ct, prefix, PETSC_DECIDE, qorder, PETSC_TRUE, fem));
204320cf1dd8SToby Isaac   PetscFunctionReturn(0);
204420cf1dd8SToby Isaac }
20453f6b16c7SMatthew G. Knepley 
2046e703855dSMatthew G. Knepley /*@
2047e703855dSMatthew G. Knepley   PetscFECreateLagrange - Create a PetscFE for the basic Lagrange space of degree k
2048e703855dSMatthew G. Knepley 
2049e703855dSMatthew G. Knepley   Collective
2050e703855dSMatthew G. Knepley 
2051e703855dSMatthew G. Knepley   Input Parameters:
2052e703855dSMatthew G. Knepley + comm      - The MPI comm
2053e703855dSMatthew G. Knepley . dim       - The spatial dimension
2054e703855dSMatthew G. Knepley . Nc        - The number of components
2055e703855dSMatthew G. Knepley . isSimplex - Flag for simplex reference cell, otherwise its a tensor product
2056e703855dSMatthew G. Knepley . k         - The degree k of the space
2057e703855dSMatthew G. Knepley - qorder    - The quadrature order or PETSC_DETERMINE to use PetscSpace polynomial degree
2058e703855dSMatthew G. Knepley 
2059e703855dSMatthew G. Knepley   Output Parameter:
2060e703855dSMatthew G. Knepley . fem       - The PetscFE object
2061e703855dSMatthew G. Knepley 
2062e703855dSMatthew G. Knepley   Level: beginner
2063e703855dSMatthew G. Knepley 
2064e703855dSMatthew G. Knepley   Notes:
2065e703855dSMatthew G. Knepley   For simplices, this element is the space of maximum polynomial degree k, otherwise it is a tensor product of 1D polynomials, each with maximal degree k.
2066e703855dSMatthew G. Knepley 
2067db781477SPatrick Sanan .seealso: `PetscFECreateLagrangeByCell()`, `PetscFECreateDefault()`, `PetscFECreateByCell()`, `PetscFECreate()`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()`
2068e703855dSMatthew G. Knepley @*/
2069e703855dSMatthew G. Knepley PetscErrorCode PetscFECreateLagrange(MPI_Comm comm, PetscInt dim, PetscInt Nc, PetscBool isSimplex, PetscInt k, PetscInt qorder, PetscFE *fem)
2070e703855dSMatthew G. Knepley {
2071e703855dSMatthew G. Knepley   PetscFunctionBegin;
20729566063dSJacob Faibussowitsch   PetscCall(PetscFECreate_Internal(comm, dim, Nc, DMPolytopeTypeSimpleShape(dim, isSimplex), NULL, k, qorder, PETSC_FALSE, fem));
20732df84da0SMatthew G. Knepley   PetscFunctionReturn(0);
2074e703855dSMatthew G. Knepley }
20752df84da0SMatthew G. Knepley 
20762df84da0SMatthew G. Knepley /*@
20772df84da0SMatthew G. Knepley   PetscFECreateLagrangeByCell - Create a PetscFE for the basic Lagrange space of degree k
20782df84da0SMatthew G. Knepley 
20792df84da0SMatthew G. Knepley   Collective
20802df84da0SMatthew G. Knepley 
20812df84da0SMatthew G. Knepley   Input Parameters:
20822df84da0SMatthew G. Knepley + comm      - The MPI comm
20832df84da0SMatthew G. Knepley . dim       - The spatial dimension
20842df84da0SMatthew G. Knepley . Nc        - The number of components
20852df84da0SMatthew G. Knepley . ct        - The celltype of the reference cell
20862df84da0SMatthew G. Knepley . k         - The degree k of the space
20872df84da0SMatthew G. Knepley - qorder    - The quadrature order or PETSC_DETERMINE to use PetscSpace polynomial degree
20882df84da0SMatthew G. Knepley 
20892df84da0SMatthew G. Knepley   Output Parameter:
20902df84da0SMatthew G. Knepley . fem       - The PetscFE object
20912df84da0SMatthew G. Knepley 
20922df84da0SMatthew G. Knepley   Level: beginner
20932df84da0SMatthew G. Knepley 
20942df84da0SMatthew G. Knepley   Notes:
20952df84da0SMatthew G. Knepley   For simplices, this element is the space of maximum polynomial degree k, otherwise it is a tensor product of 1D polynomials, each with maximal degree k.
20962df84da0SMatthew G. Knepley 
2097db781477SPatrick Sanan .seealso: `PetscFECreateLagrange()`, `PetscFECreateDefault()`, `PetscFECreateByCell()`, `PetscFECreate()`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()`
20982df84da0SMatthew G. Knepley @*/
20992df84da0SMatthew G. Knepley PetscErrorCode PetscFECreateLagrangeByCell(MPI_Comm comm, PetscInt dim, PetscInt Nc, DMPolytopeType ct, PetscInt k, PetscInt qorder, PetscFE *fem)
21002df84da0SMatthew G. Knepley {
21012df84da0SMatthew G. Knepley   PetscFunctionBegin;
21029566063dSJacob Faibussowitsch   PetscCall(PetscFECreate_Internal(comm, dim, Nc, ct, NULL, k, qorder, PETSC_FALSE, fem));
2103e703855dSMatthew G. Knepley   PetscFunctionReturn(0);
2104e703855dSMatthew G. Knepley }
2105e703855dSMatthew G. Knepley 
21063f6b16c7SMatthew G. Knepley /*@C
21073f6b16c7SMatthew G. Knepley   PetscFESetName - Names the FE and its subobjects
21083f6b16c7SMatthew G. Knepley 
21093f6b16c7SMatthew G. Knepley   Not collective
21103f6b16c7SMatthew G. Knepley 
21113f6b16c7SMatthew G. Knepley   Input Parameters:
21123f6b16c7SMatthew G. Knepley + fe   - The PetscFE
21133f6b16c7SMatthew G. Knepley - name - The name
21143f6b16c7SMatthew G. Knepley 
21152b99622eSMatthew G. Knepley   Level: intermediate
21163f6b16c7SMatthew G. Knepley 
2117db781477SPatrick Sanan .seealso: `PetscFECreate()`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()`
21183f6b16c7SMatthew G. Knepley @*/
21193f6b16c7SMatthew G. Knepley PetscErrorCode PetscFESetName(PetscFE fe, const char name[])
21203f6b16c7SMatthew G. Knepley {
21213f6b16c7SMatthew G. Knepley   PetscSpace     P;
21223f6b16c7SMatthew G. Knepley   PetscDualSpace Q;
21233f6b16c7SMatthew G. Knepley 
21243f6b16c7SMatthew G. Knepley   PetscFunctionBegin;
21259566063dSJacob Faibussowitsch   PetscCall(PetscFEGetBasisSpace(fe, &P));
21269566063dSJacob Faibussowitsch   PetscCall(PetscFEGetDualSpace(fe, &Q));
21279566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) fe, name));
21289566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) P,  name));
21299566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) Q,  name));
21303f6b16c7SMatthew G. Knepley   PetscFunctionReturn(0);
21313f6b16c7SMatthew G. Knepley }
2132a8f1f9e5SMatthew G. Knepley 
2133ef0bb6c7SMatthew G. Knepley PetscErrorCode PetscFEEvaluateFieldJets_Internal(PetscDS ds, PetscInt Nf, PetscInt r, PetscInt q, PetscTabulation T[], PetscFEGeom *fegeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscScalar u[], PetscScalar u_x[], PetscScalar u_t[])
2134a8f1f9e5SMatthew G. Knepley {
2135f9244615SMatthew G. Knepley   PetscInt       dOffset = 0, fOffset = 0, f, g;
2136a8f1f9e5SMatthew G. Knepley 
2137a8f1f9e5SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
2138a8f1f9e5SMatthew G. Knepley     PetscFE          fe;
2139f9244615SMatthew G. Knepley     const PetscInt   k    = ds->jetDegree[f];
2140ef0bb6c7SMatthew G. Knepley     const PetscInt   cdim = T[f]->cdim;
2141ef0bb6c7SMatthew G. Knepley     const PetscInt   Nq   = T[f]->Np;
2142ef0bb6c7SMatthew G. Knepley     const PetscInt   Nbf  = T[f]->Nb;
2143ef0bb6c7SMatthew G. Knepley     const PetscInt   Ncf  = T[f]->Nc;
2144ef0bb6c7SMatthew G. Knepley     const PetscReal *Bq   = &T[f]->T[0][(r*Nq+q)*Nbf*Ncf];
2145ef0bb6c7SMatthew G. Knepley     const PetscReal *Dq   = &T[f]->T[1][(r*Nq+q)*Nbf*Ncf*cdim];
2146f9244615SMatthew G. Knepley     const PetscReal *Hq   = k > 1 ? &T[f]->T[2][(r*Nq+q)*Nbf*Ncf*cdim*cdim] : NULL;
2147f9244615SMatthew G. Knepley     PetscInt         hOffset = 0, b, c, d;
2148a8f1f9e5SMatthew G. Knepley 
21499566063dSJacob Faibussowitsch     PetscCall(PetscDSGetDiscretization(ds, f, (PetscObject *) &fe));
2150a8f1f9e5SMatthew G. Knepley     for (c = 0; c < Ncf; ++c) u[fOffset+c] = 0.0;
2151ef0bb6c7SMatthew G. Knepley     for (d = 0; d < cdim*Ncf; ++d) u_x[fOffset*cdim+d] = 0.0;
2152a8f1f9e5SMatthew G. Knepley     for (b = 0; b < Nbf; ++b) {
2153a8f1f9e5SMatthew G. Knepley       for (c = 0; c < Ncf; ++c) {
2154a8f1f9e5SMatthew G. Knepley         const PetscInt cidx = b*Ncf+c;
2155a8f1f9e5SMatthew G. Knepley 
2156a8f1f9e5SMatthew G. Knepley         u[fOffset+c] += Bq[cidx]*coefficients[dOffset+b];
2157ef0bb6c7SMatthew G. Knepley         for (d = 0; d < cdim; ++d) u_x[(fOffset+c)*cdim+d] += Dq[cidx*cdim+d]*coefficients[dOffset+b];
2158a8f1f9e5SMatthew G. Knepley       }
2159a8f1f9e5SMatthew G. Knepley     }
2160f9244615SMatthew G. Knepley     if (k > 1) {
2161f9244615SMatthew G. Knepley       for (g = 0; g < Nf; ++g) hOffset += T[g]->Nc*cdim;
2162f9244615SMatthew G. Knepley       for (d = 0; d < cdim*cdim*Ncf; ++d) u_x[hOffset+fOffset*cdim*cdim+d] = 0.0;
2163f9244615SMatthew G. Knepley       for (b = 0; b < Nbf; ++b) {
2164f9244615SMatthew G. Knepley         for (c = 0; c < Ncf; ++c) {
2165f9244615SMatthew G. Knepley           const PetscInt cidx = b*Ncf+c;
2166f9244615SMatthew G. Knepley 
2167f9244615SMatthew G. Knepley           for (d = 0; d < cdim*cdim; ++d) u_x[hOffset+(fOffset+c)*cdim*cdim+d] += Hq[cidx*cdim*cdim+d]*coefficients[dOffset+b];
2168f9244615SMatthew G. Knepley         }
2169f9244615SMatthew G. Knepley       }
21709566063dSJacob Faibussowitsch       PetscCall(PetscFEPushforwardHessian(fe, fegeom, 1, &u_x[hOffset+fOffset*cdim*cdim]));
2171f9244615SMatthew G. Knepley     }
21729566063dSJacob Faibussowitsch     PetscCall(PetscFEPushforward(fe, fegeom, 1, &u[fOffset]));
21739566063dSJacob Faibussowitsch     PetscCall(PetscFEPushforwardGradient(fe, fegeom, 1, &u_x[fOffset*cdim]));
2174a8f1f9e5SMatthew G. Knepley     if (u_t) {
2175a8f1f9e5SMatthew G. Knepley       for (c = 0; c < Ncf; ++c) u_t[fOffset+c] = 0.0;
2176a8f1f9e5SMatthew G. Knepley       for (b = 0; b < Nbf; ++b) {
2177a8f1f9e5SMatthew G. Knepley         for (c = 0; c < Ncf; ++c) {
2178a8f1f9e5SMatthew G. Knepley           const PetscInt cidx = b*Ncf+c;
2179a8f1f9e5SMatthew G. Knepley 
2180a8f1f9e5SMatthew G. Knepley           u_t[fOffset+c] += Bq[cidx]*coefficients_t[dOffset+b];
2181a8f1f9e5SMatthew G. Knepley         }
2182a8f1f9e5SMatthew G. Knepley       }
21839566063dSJacob Faibussowitsch       PetscCall(PetscFEPushforward(fe, fegeom, 1, &u_t[fOffset]));
2184a8f1f9e5SMatthew G. Knepley     }
2185a8f1f9e5SMatthew G. Knepley     fOffset += Ncf;
2186a8f1f9e5SMatthew G. Knepley     dOffset += Nbf;
2187a8f1f9e5SMatthew G. Knepley   }
2188a8f1f9e5SMatthew G. Knepley   return 0;
2189a8f1f9e5SMatthew G. Knepley }
2190a8f1f9e5SMatthew G. Knepley 
2191665f567fSMatthew G. Knepley PetscErrorCode PetscFEEvaluateFieldJets_Hybrid_Internal(PetscDS ds, PetscInt Nf, PetscInt r, PetscInt q, PetscTabulation T[], PetscFEGeom *fegeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscScalar u[], PetscScalar u_x[], PetscScalar u_t[])
219227f02ce8SMatthew G. Knepley {
21935fedec97SMatthew G. Knepley   PetscInt       dOffset = 0, fOffset = 0, f, g;
219427f02ce8SMatthew G. Knepley 
21955fedec97SMatthew G. Knepley   /* f is the field number in the DS, g is the field number in u[] */
21965fedec97SMatthew G. Knepley   for (f = 0, g = 0; f < Nf; ++f) {
21975fedec97SMatthew G. Knepley     PetscFE          fe   = (PetscFE) ds->disc[f];
21989ee2af8cSMatthew G. Knepley     const PetscInt   dEt  = T[f]->cdim;
21999ee2af8cSMatthew G. Knepley     const PetscInt   dE   = fegeom->dimEmbed;
2200665f567fSMatthew G. Knepley     const PetscInt   Nq   = T[f]->Np;
2201665f567fSMatthew G. Knepley     const PetscInt   Nbf  = T[f]->Nb;
2202665f567fSMatthew G. Knepley     const PetscInt   Ncf  = T[f]->Nc;
2203665f567fSMatthew G. Knepley     const PetscReal *Bq   = &T[f]->T[0][(r*Nq+q)*Nbf*Ncf];
22049ee2af8cSMatthew G. Knepley     const PetscReal *Dq   = &T[f]->T[1][(r*Nq+q)*Nbf*Ncf*dEt];
22055fedec97SMatthew G. Knepley     PetscBool        isCohesive;
22065fedec97SMatthew G. Knepley     PetscInt         Ns, s;
22075fedec97SMatthew G. Knepley 
22085fedec97SMatthew G. Knepley     if (!T[f]) continue;
22099566063dSJacob Faibussowitsch     PetscCall(PetscDSGetCohesive(ds, f, &isCohesive));
22105fedec97SMatthew G. Knepley     Ns   = isCohesive ? 1 : 2;
22115fedec97SMatthew G. Knepley     for (s = 0; s < Ns; ++s, ++g) {
221227f02ce8SMatthew G. Knepley       PetscInt b, c, d;
221327f02ce8SMatthew G. Knepley 
221427f02ce8SMatthew G. Knepley       for (c = 0; c < Ncf; ++c)    u[fOffset+c]      = 0.0;
22159ee2af8cSMatthew G. Knepley       for (d = 0; d < dE*Ncf; ++d) u_x[fOffset*dE+d] = 0.0;
221627f02ce8SMatthew G. Knepley       for (b = 0; b < Nbf; ++b) {
221727f02ce8SMatthew G. Knepley         for (c = 0; c < Ncf; ++c) {
221827f02ce8SMatthew G. Knepley           const PetscInt cidx = b*Ncf+c;
221927f02ce8SMatthew G. Knepley 
222027f02ce8SMatthew G. Knepley           u[fOffset+c] += Bq[cidx]*coefficients[dOffset+b];
22219ee2af8cSMatthew G. Knepley           for (d = 0; d < dEt; ++d) u_x[(fOffset+c)*dE+d] += Dq[cidx*dEt+d]*coefficients[dOffset+b];
222227f02ce8SMatthew G. Knepley         }
222327f02ce8SMatthew G. Knepley       }
22249566063dSJacob Faibussowitsch       PetscCall(PetscFEPushforward(fe, fegeom, 1, &u[fOffset]));
22259566063dSJacob Faibussowitsch       PetscCall(PetscFEPushforwardGradient(fe, fegeom, 1, &u_x[fOffset*dE]));
222627f02ce8SMatthew G. Knepley       if (u_t) {
222727f02ce8SMatthew G. Knepley         for (c = 0; c < Ncf; ++c) u_t[fOffset+c] = 0.0;
222827f02ce8SMatthew G. Knepley         for (b = 0; b < Nbf; ++b) {
222927f02ce8SMatthew G. Knepley           for (c = 0; c < Ncf; ++c) {
223027f02ce8SMatthew G. Knepley             const PetscInt cidx = b*Ncf+c;
223127f02ce8SMatthew G. Knepley 
223227f02ce8SMatthew G. Knepley             u_t[fOffset+c] += Bq[cidx]*coefficients_t[dOffset+b];
223327f02ce8SMatthew G. Knepley           }
223427f02ce8SMatthew G. Knepley         }
22359566063dSJacob Faibussowitsch         PetscCall(PetscFEPushforward(fe, fegeom, 1, &u_t[fOffset]));
223627f02ce8SMatthew G. Knepley       }
223727f02ce8SMatthew G. Knepley       fOffset += Ncf;
223827f02ce8SMatthew G. Knepley       dOffset += Nbf;
223927f02ce8SMatthew G. Knepley     }
2240665f567fSMatthew G. Knepley   }
224127f02ce8SMatthew G. Knepley   return 0;
224227f02ce8SMatthew G. Knepley }
224327f02ce8SMatthew G. Knepley 
2244a8f1f9e5SMatthew G. Knepley PetscErrorCode PetscFEEvaluateFaceFields_Internal(PetscDS prob, PetscInt field, PetscInt faceLoc, const PetscScalar coefficients[], PetscScalar u[])
2245a8f1f9e5SMatthew G. Knepley {
2246a8f1f9e5SMatthew G. Knepley   PetscFE         fe;
2247ef0bb6c7SMatthew G. Knepley   PetscTabulation Tc;
2248ef0bb6c7SMatthew G. Knepley   PetscInt        b, c;
2249a8f1f9e5SMatthew G. Knepley 
2250a8f1f9e5SMatthew G. Knepley   if (!prob) return 0;
22519566063dSJacob Faibussowitsch   PetscCall(PetscDSGetDiscretization(prob, field, (PetscObject *) &fe));
22529566063dSJacob Faibussowitsch   PetscCall(PetscFEGetFaceCentroidTabulation(fe, &Tc));
2253ef0bb6c7SMatthew G. Knepley   {
2254ef0bb6c7SMatthew G. Knepley     const PetscReal *faceBasis = Tc->T[0];
2255ef0bb6c7SMatthew G. Knepley     const PetscInt   Nb        = Tc->Nb;
2256ef0bb6c7SMatthew G. Knepley     const PetscInt   Nc        = Tc->Nc;
2257ef0bb6c7SMatthew G. Knepley 
2258a8f1f9e5SMatthew G. Knepley     for (c = 0; c < Nc; ++c) {u[c] = 0.0;}
2259a8f1f9e5SMatthew G. Knepley     for (b = 0; b < Nb; ++b) {
2260a8f1f9e5SMatthew G. Knepley       for (c = 0; c < Nc; ++c) {
2261813a933aSJed Brown         u[c] += coefficients[b] * faceBasis[(faceLoc*Nb + b)*Nc + c];
2262a8f1f9e5SMatthew G. Knepley       }
2263a8f1f9e5SMatthew G. Knepley     }
2264ef0bb6c7SMatthew G. Knepley   }
2265a8f1f9e5SMatthew G. Knepley   return 0;
2266a8f1f9e5SMatthew G. Knepley }
2267a8f1f9e5SMatthew G. Knepley 
22686587ee25SMatthew G. Knepley PetscErrorCode PetscFEUpdateElementVec_Internal(PetscFE fe, PetscTabulation T, PetscInt r, PetscScalar tmpBasis[], PetscScalar tmpBasisDer[], PetscInt e, PetscFEGeom *fegeom, PetscScalar f0[], PetscScalar f1[], PetscScalar elemVec[])
2269a8f1f9e5SMatthew G. Knepley {
22706587ee25SMatthew G. Knepley   PetscFEGeom      pgeom;
2271bc3a64adSMatthew G. Knepley   const PetscInt   dEt      = T->cdim;
2272bc3a64adSMatthew G. Knepley   const PetscInt   dE       = fegeom->dimEmbed;
2273ef0bb6c7SMatthew G. Knepley   const PetscInt   Nq       = T->Np;
2274ef0bb6c7SMatthew G. Knepley   const PetscInt   Nb       = T->Nb;
2275ef0bb6c7SMatthew G. Knepley   const PetscInt   Nc       = T->Nc;
2276ef0bb6c7SMatthew G. Knepley   const PetscReal *basis    = &T->T[0][r*Nq*Nb*Nc];
2277bc3a64adSMatthew G. Knepley   const PetscReal *basisDer = &T->T[1][r*Nq*Nb*Nc*dEt];
2278a8f1f9e5SMatthew G. Knepley   PetscInt         q, b, c, d;
2279a8f1f9e5SMatthew G. Knepley 
2280a8f1f9e5SMatthew G. Knepley   for (q = 0; q < Nq; ++q) {
2281a8f1f9e5SMatthew G. Knepley     for (b = 0; b < Nb; ++b) {
2282a8f1f9e5SMatthew G. Knepley       for (c = 0; c < Nc; ++c) {
2283a8f1f9e5SMatthew G. Knepley         const PetscInt bcidx = b*Nc+c;
2284a8f1f9e5SMatthew G. Knepley 
2285a8f1f9e5SMatthew G. Knepley         tmpBasis[bcidx] = basis[q*Nb*Nc+bcidx];
2286bc3a64adSMatthew G. Knepley         for (d = 0; d < dEt; ++d) tmpBasisDer[bcidx*dE+d] = basisDer[q*Nb*Nc*dEt+bcidx*dEt+d];
22879ee2af8cSMatthew G. Knepley         for (d = dEt; d < dE; ++d) tmpBasisDer[bcidx*dE+d] = 0.0;
2288a8f1f9e5SMatthew G. Knepley       }
2289a8f1f9e5SMatthew G. Knepley     }
22909566063dSJacob Faibussowitsch     PetscCall(PetscFEGeomGetCellPoint(fegeom, e, q, &pgeom));
22919566063dSJacob Faibussowitsch     PetscCall(PetscFEPushforward(fe, &pgeom, Nb, tmpBasis));
22929566063dSJacob Faibussowitsch     PetscCall(PetscFEPushforwardGradient(fe, &pgeom, Nb, tmpBasisDer));
2293a8f1f9e5SMatthew G. Knepley     for (b = 0; b < Nb; ++b) {
2294a8f1f9e5SMatthew G. Knepley       for (c = 0; c < Nc; ++c) {
2295a8f1f9e5SMatthew G. Knepley         const PetscInt bcidx = b*Nc+c;
2296a8f1f9e5SMatthew G. Knepley         const PetscInt qcidx = q*Nc+c;
2297a8f1f9e5SMatthew G. Knepley 
2298a8f1f9e5SMatthew G. Knepley         elemVec[b] += tmpBasis[bcidx]*f0[qcidx];
229927f02ce8SMatthew G. Knepley         for (d = 0; d < dE; ++d) elemVec[b] += tmpBasisDer[bcidx*dE+d]*f1[qcidx*dE+d];
230027f02ce8SMatthew G. Knepley       }
230127f02ce8SMatthew G. Knepley     }
230227f02ce8SMatthew G. Knepley   }
230327f02ce8SMatthew G. Knepley   return(0);
230427f02ce8SMatthew G. Knepley }
230527f02ce8SMatthew G. Knepley 
2306c2b7495fSMatthew G. Knepley PetscErrorCode PetscFEUpdateElementVec_Hybrid_Internal(PetscFE fe, PetscTabulation T, PetscInt r, PetscInt s, PetscScalar tmpBasis[], PetscScalar tmpBasisDer[], PetscFEGeom *fegeom, PetscScalar f0[], PetscScalar f1[], PetscScalar elemVec[])
230727f02ce8SMatthew G. Knepley {
230827f02ce8SMatthew G. Knepley   const PetscInt   dE       = T->cdim;
230927f02ce8SMatthew G. Knepley   const PetscInt   Nq       = T->Np;
231027f02ce8SMatthew G. Knepley   const PetscInt   Nb       = T->Nb;
231127f02ce8SMatthew G. Knepley   const PetscInt   Nc       = T->Nc;
231227f02ce8SMatthew G. Knepley   const PetscReal *basis    = &T->T[0][r*Nq*Nb*Nc];
231327f02ce8SMatthew G. Knepley   const PetscReal *basisDer = &T->T[1][r*Nq*Nb*Nc*dE];
2314c2b7495fSMatthew G. Knepley   PetscInt         q, b, c, d;
231527f02ce8SMatthew G. Knepley 
231627f02ce8SMatthew G. Knepley   for (q = 0; q < Nq; ++q) {
231727f02ce8SMatthew G. Knepley     for (b = 0; b < Nb; ++b) {
231827f02ce8SMatthew G. Knepley       for (c = 0; c < Nc; ++c) {
231927f02ce8SMatthew G. Knepley         const PetscInt bcidx = b*Nc+c;
232027f02ce8SMatthew G. Knepley 
232127f02ce8SMatthew G. Knepley         tmpBasis[bcidx] = basis[q*Nb*Nc+bcidx];
232227f02ce8SMatthew G. Knepley         for (d = 0; d < dE; ++d) tmpBasisDer[bcidx*dE+d] = basisDer[q*Nb*Nc*dE+bcidx*dE+d];
232327f02ce8SMatthew G. Knepley       }
232427f02ce8SMatthew G. Knepley     }
23259566063dSJacob Faibussowitsch     PetscCall(PetscFEPushforward(fe, fegeom, Nb, tmpBasis));
23269566063dSJacob Faibussowitsch     PetscCall(PetscFEPushforwardGradient(fe, fegeom, Nb, tmpBasisDer));
232727f02ce8SMatthew G. Knepley     for (b = 0; b < Nb; ++b) {
232827f02ce8SMatthew G. Knepley       for (c = 0; c < Nc; ++c) {
232927f02ce8SMatthew G. Knepley         const PetscInt bcidx = b*Nc+c;
2330c2b7495fSMatthew G. Knepley         const PetscInt qcidx = q*Nc+c;
233127f02ce8SMatthew G. Knepley 
233227f02ce8SMatthew G. Knepley         elemVec[Nb*s+b] += tmpBasis[bcidx]*f0[qcidx];
233327f02ce8SMatthew G. Knepley         for (d = 0; d < dE; ++d) elemVec[Nb*s+b] += tmpBasisDer[bcidx*dE+d]*f1[qcidx*dE+d];
233427f02ce8SMatthew G. Knepley       }
2335a8f1f9e5SMatthew G. Knepley     }
2336a8f1f9e5SMatthew G. Knepley   }
2337a8f1f9e5SMatthew G. Knepley   return(0);
2338a8f1f9e5SMatthew G. Knepley }
2339a8f1f9e5SMatthew G. Knepley 
2340ef0bb6c7SMatthew G. Knepley PetscErrorCode PetscFEUpdateElementMat_Internal(PetscFE feI, PetscFE feJ, PetscInt r, PetscInt q, PetscTabulation TI, PetscScalar tmpBasisI[], PetscScalar tmpBasisDerI[], PetscTabulation TJ, 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[])
2341a8f1f9e5SMatthew G. Knepley {
234227f02ce8SMatthew G. Knepley   const PetscInt   dE        = TI->cdim;
2343ef0bb6c7SMatthew G. Knepley   const PetscInt   NqI       = TI->Np;
2344ef0bb6c7SMatthew G. Knepley   const PetscInt   NbI       = TI->Nb;
2345ef0bb6c7SMatthew G. Knepley   const PetscInt   NcI       = TI->Nc;
2346ef0bb6c7SMatthew G. Knepley   const PetscReal *basisI    = &TI->T[0][(r*NqI+q)*NbI*NcI];
2347665f567fSMatthew G. Knepley   const PetscReal *basisDerI = &TI->T[1][(r*NqI+q)*NbI*NcI*dE];
2348ef0bb6c7SMatthew G. Knepley   const PetscInt   NqJ       = TJ->Np;
2349ef0bb6c7SMatthew G. Knepley   const PetscInt   NbJ       = TJ->Nb;
2350ef0bb6c7SMatthew G. Knepley   const PetscInt   NcJ       = TJ->Nc;
2351ef0bb6c7SMatthew G. Knepley   const PetscReal *basisJ    = &TJ->T[0][(r*NqJ+q)*NbJ*NcJ];
2352665f567fSMatthew G. Knepley   const PetscReal *basisDerJ = &TJ->T[1][(r*NqJ+q)*NbJ*NcJ*dE];
2353a8f1f9e5SMatthew G. Knepley   PetscInt         f, fc, g, gc, df, dg;
2354a8f1f9e5SMatthew G. Knepley 
2355a8f1f9e5SMatthew G. Knepley   for (f = 0; f < NbI; ++f) {
2356a8f1f9e5SMatthew G. Knepley     for (fc = 0; fc < NcI; ++fc) {
2357a8f1f9e5SMatthew G. Knepley       const PetscInt fidx = f*NcI+fc; /* Test function basis index */
2358a8f1f9e5SMatthew G. Knepley 
2359a8f1f9e5SMatthew G. Knepley       tmpBasisI[fidx] = basisI[fidx];
236027f02ce8SMatthew G. Knepley       for (df = 0; df < dE; ++df) tmpBasisDerI[fidx*dE+df] = basisDerI[fidx*dE+df];
2361a8f1f9e5SMatthew G. Knepley     }
2362a8f1f9e5SMatthew G. Knepley   }
23639566063dSJacob Faibussowitsch   PetscCall(PetscFEPushforward(feI, fegeom, NbI, tmpBasisI));
23649566063dSJacob Faibussowitsch   PetscCall(PetscFEPushforwardGradient(feI, fegeom, NbI, tmpBasisDerI));
2365a8f1f9e5SMatthew G. Knepley   for (g = 0; g < NbJ; ++g) {
2366a8f1f9e5SMatthew G. Knepley     for (gc = 0; gc < NcJ; ++gc) {
2367a8f1f9e5SMatthew G. Knepley       const PetscInt gidx = g*NcJ+gc; /* Trial function basis index */
2368a8f1f9e5SMatthew G. Knepley 
2369a8f1f9e5SMatthew G. Knepley       tmpBasisJ[gidx] = basisJ[gidx];
237027f02ce8SMatthew G. Knepley       for (dg = 0; dg < dE; ++dg) tmpBasisDerJ[gidx*dE+dg] = basisDerJ[gidx*dE+dg];
2371a8f1f9e5SMatthew G. Knepley     }
2372a8f1f9e5SMatthew G. Knepley   }
23739566063dSJacob Faibussowitsch   PetscCall(PetscFEPushforward(feJ, fegeom, NbJ, tmpBasisJ));
23749566063dSJacob Faibussowitsch   PetscCall(PetscFEPushforwardGradient(feJ, fegeom, NbJ, tmpBasisDerJ));
2375a8f1f9e5SMatthew G. Knepley   for (f = 0; f < NbI; ++f) {
2376a8f1f9e5SMatthew G. Knepley     for (fc = 0; fc < NcI; ++fc) {
2377a8f1f9e5SMatthew G. Knepley       const PetscInt fidx = f*NcI+fc; /* Test function basis index */
2378a8f1f9e5SMatthew G. Knepley       const PetscInt i    = offsetI+f; /* Element matrix row */
2379a8f1f9e5SMatthew G. Knepley       for (g = 0; g < NbJ; ++g) {
2380a8f1f9e5SMatthew G. Knepley         for (gc = 0; gc < NcJ; ++gc) {
2381a8f1f9e5SMatthew G. Knepley           const PetscInt gidx = g*NcJ+gc; /* Trial function basis index */
2382a8f1f9e5SMatthew G. Knepley           const PetscInt j    = offsetJ+g; /* Element matrix column */
2383a8f1f9e5SMatthew G. Knepley           const PetscInt fOff = eOffset+i*totDim+j;
2384a8f1f9e5SMatthew G. Knepley 
2385a8f1f9e5SMatthew G. Knepley           elemMat[fOff] += tmpBasisI[fidx]*g0[fc*NcJ+gc]*tmpBasisJ[gidx];
238627f02ce8SMatthew G. Knepley           for (df = 0; df < dE; ++df) {
238727f02ce8SMatthew G. Knepley             elemMat[fOff] += tmpBasisI[fidx]*g1[(fc*NcJ+gc)*dE+df]*tmpBasisDerJ[gidx*dE+df];
238827f02ce8SMatthew G. Knepley             elemMat[fOff] += tmpBasisDerI[fidx*dE+df]*g2[(fc*NcJ+gc)*dE+df]*tmpBasisJ[gidx];
238927f02ce8SMatthew G. Knepley             for (dg = 0; dg < dE; ++dg) {
239027f02ce8SMatthew G. Knepley               elemMat[fOff] += tmpBasisDerI[fidx*dE+df]*g3[((fc*NcJ+gc)*dE+df)*dE+dg]*tmpBasisDerJ[gidx*dE+dg];
239127f02ce8SMatthew G. Knepley             }
239227f02ce8SMatthew G. Knepley           }
239327f02ce8SMatthew G. Knepley         }
239427f02ce8SMatthew G. Knepley       }
239527f02ce8SMatthew G. Knepley     }
239627f02ce8SMatthew G. Knepley   }
239727f02ce8SMatthew G. Knepley   return(0);
239827f02ce8SMatthew G. Knepley }
239927f02ce8SMatthew G. Knepley 
24005fedec97SMatthew G. Knepley PetscErrorCode PetscFEUpdateElementMat_Hybrid_Internal(PetscFE feI, PetscBool isHybridI, PetscFE feJ, PetscBool isHybridJ, PetscInt r, PetscInt s, PetscInt q, PetscTabulation TI, PetscScalar tmpBasisI[], PetscScalar tmpBasisDerI[], PetscTabulation TJ, 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[])
240127f02ce8SMatthew G. Knepley {
2402665f567fSMatthew G. Knepley   const PetscInt   dE        = TI->cdim;
2403665f567fSMatthew G. Knepley   const PetscInt   NqI       = TI->Np;
2404665f567fSMatthew G. Knepley   const PetscInt   NbI       = TI->Nb;
2405665f567fSMatthew G. Knepley   const PetscInt   NcI       = TI->Nc;
2406665f567fSMatthew G. Knepley   const PetscReal *basisI    = &TI->T[0][(r*NqI+q)*NbI*NcI];
2407665f567fSMatthew G. Knepley   const PetscReal *basisDerI = &TI->T[1][(r*NqI+q)*NbI*NcI*dE];
2408665f567fSMatthew G. Knepley   const PetscInt   NqJ       = TJ->Np;
2409665f567fSMatthew G. Knepley   const PetscInt   NbJ       = TJ->Nb;
2410665f567fSMatthew G. Knepley   const PetscInt   NcJ       = TJ->Nc;
2411665f567fSMatthew G. Knepley   const PetscReal *basisJ    = &TJ->T[0][(r*NqJ+q)*NbJ*NcJ];
2412665f567fSMatthew G. Knepley   const PetscReal *basisDerJ = &TJ->T[1][(r*NqJ+q)*NbJ*NcJ*dE];
24135fedec97SMatthew G. Knepley   const PetscInt   so        = isHybridI ? 0 : s;
24145fedec97SMatthew G. Knepley   const PetscInt   to        = isHybridJ ? 0 : s;
24155fedec97SMatthew G. Knepley   PetscInt         f, fc, g, gc, df, dg;
241627f02ce8SMatthew G. Knepley 
241727f02ce8SMatthew G. Knepley   for (f = 0; f < NbI; ++f) {
241827f02ce8SMatthew G. Knepley     for (fc = 0; fc < NcI; ++fc) {
241927f02ce8SMatthew G. Knepley       const PetscInt fidx = f*NcI+fc; /* Test function basis index */
242027f02ce8SMatthew G. Knepley 
242127f02ce8SMatthew G. Knepley       tmpBasisI[fidx] = basisI[fidx];
2422665f567fSMatthew G. Knepley       for (df = 0; df < dE; ++df) tmpBasisDerI[fidx*dE+df] = basisDerI[fidx*dE+df];
242327f02ce8SMatthew G. Knepley     }
242427f02ce8SMatthew G. Knepley   }
24259566063dSJacob Faibussowitsch   PetscCall(PetscFEPushforward(feI, fegeom, NbI, tmpBasisI));
24269566063dSJacob Faibussowitsch   PetscCall(PetscFEPushforwardGradient(feI, fegeom, NbI, tmpBasisDerI));
242727f02ce8SMatthew G. Knepley   for (g = 0; g < NbJ; ++g) {
242827f02ce8SMatthew G. Knepley     for (gc = 0; gc < NcJ; ++gc) {
242927f02ce8SMatthew G. Knepley       const PetscInt gidx = g*NcJ+gc; /* Trial function basis index */
243027f02ce8SMatthew G. Knepley 
243127f02ce8SMatthew G. Knepley       tmpBasisJ[gidx] = basisJ[gidx];
2432665f567fSMatthew G. Knepley       for (dg = 0; dg < dE; ++dg) tmpBasisDerJ[gidx*dE+dg] = basisDerJ[gidx*dE+dg];
243327f02ce8SMatthew G. Knepley     }
243427f02ce8SMatthew G. Knepley   }
24359566063dSJacob Faibussowitsch   PetscCall(PetscFEPushforward(feJ, fegeom, NbJ, tmpBasisJ));
24369566063dSJacob Faibussowitsch   PetscCall(PetscFEPushforwardGradient(feJ, fegeom, NbJ, tmpBasisDerJ));
243727f02ce8SMatthew G. Knepley   for (f = 0; f < NbI; ++f) {
243827f02ce8SMatthew G. Knepley     for (fc = 0; fc < NcI; ++fc) {
243927f02ce8SMatthew G. Knepley       const PetscInt fidx = f*NcI+fc;         /* Test function basis index */
24405fedec97SMatthew G. Knepley       const PetscInt i    = offsetI+NbI*so+f; /* Element matrix row */
244127f02ce8SMatthew G. Knepley       for (g = 0; g < NbJ; ++g) {
244227f02ce8SMatthew G. Knepley         for (gc = 0; gc < NcJ; ++gc) {
244327f02ce8SMatthew G. Knepley           const PetscInt gidx = g*NcJ+gc;         /* Trial function basis index */
24445fedec97SMatthew G. Knepley           const PetscInt j    = offsetJ+NbJ*to+g; /* Element matrix column */
244527f02ce8SMatthew G. Knepley           const PetscInt fOff = eOffset+i*totDim+j;
244627f02ce8SMatthew G. Knepley 
24475fedec97SMatthew G. Knepley           elemMat[fOff] += tmpBasisI[fidx]*g0[fc*NcJ+gc]*tmpBasisJ[gidx];
244827f02ce8SMatthew G. Knepley           for (df = 0; df < dE; ++df) {
24495fedec97SMatthew G. Knepley             elemMat[fOff] += tmpBasisI[fidx]*g1[(fc*NcJ+gc)*dE+df]*tmpBasisDerJ[gidx*dE+df];
24505fedec97SMatthew G. Knepley             elemMat[fOff] += tmpBasisDerI[fidx*dE+df]*g2[(fc*NcJ+gc)*dE+df]*tmpBasisJ[gidx];
245127f02ce8SMatthew G. Knepley             for (dg = 0; dg < dE; ++dg) {
24525fedec97SMatthew G. Knepley               elemMat[fOff] += tmpBasisDerI[fidx*dE+df]*g3[((fc*NcJ+gc)*dE+df)*dE+dg]*tmpBasisDerJ[gidx*dE+dg];
2453a8f1f9e5SMatthew G. Knepley             }
2454a8f1f9e5SMatthew G. Knepley           }
2455a8f1f9e5SMatthew G. Knepley         }
2456a8f1f9e5SMatthew G. Knepley       }
2457a8f1f9e5SMatthew G. Knepley     }
2458a8f1f9e5SMatthew G. Knepley   }
2459a8f1f9e5SMatthew G. Knepley   return(0);
2460a8f1f9e5SMatthew G. Knepley }
2461c9ba7969SMatthew G. Knepley 
2462c9ba7969SMatthew G. Knepley PetscErrorCode PetscFECreateCellGeometry(PetscFE fe, PetscQuadrature quad, PetscFEGeom *cgeom)
2463c9ba7969SMatthew G. Knepley {
2464c9ba7969SMatthew G. Knepley   PetscDualSpace  dsp;
2465c9ba7969SMatthew G. Knepley   DM              dm;
2466c9ba7969SMatthew G. Knepley   PetscQuadrature quadDef;
2467c9ba7969SMatthew G. Knepley   PetscInt        dim, cdim, Nq;
2468c9ba7969SMatthew G. Knepley 
2469c9ba7969SMatthew G. Knepley   PetscFunctionBegin;
24709566063dSJacob Faibussowitsch   PetscCall(PetscFEGetDualSpace(fe, &dsp));
24719566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDM(dsp, &dm));
24729566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
24739566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &cdim));
24749566063dSJacob Faibussowitsch   PetscCall(PetscFEGetQuadrature(fe, &quadDef));
2475c9ba7969SMatthew G. Knepley   quad = quad ? quad : quadDef;
24769566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL));
24779566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Nq*cdim,      &cgeom->v));
24789566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Nq*cdim*cdim, &cgeom->J));
24799566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Nq*cdim*cdim, &cgeom->invJ));
24809566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Nq,           &cgeom->detJ));
2481c9ba7969SMatthew G. Knepley   cgeom->dim       = dim;
2482c9ba7969SMatthew G. Knepley   cgeom->dimEmbed  = cdim;
2483c9ba7969SMatthew G. Knepley   cgeom->numCells  = 1;
2484c9ba7969SMatthew G. Knepley   cgeom->numPoints = Nq;
24859566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeCellGeometryFEM(dm, 0, quad, cgeom->v, cgeom->J, cgeom->invJ, cgeom->detJ));
2486c9ba7969SMatthew G. Knepley   PetscFunctionReturn(0);
2487c9ba7969SMatthew G. Knepley }
2488c9ba7969SMatthew G. Knepley 
2489c9ba7969SMatthew G. Knepley PetscErrorCode PetscFEDestroyCellGeometry(PetscFE fe, PetscFEGeom *cgeom)
2490c9ba7969SMatthew G. Knepley {
2491c9ba7969SMatthew G. Knepley   PetscFunctionBegin;
24929566063dSJacob Faibussowitsch   PetscCall(PetscFree(cgeom->v));
24939566063dSJacob Faibussowitsch   PetscCall(PetscFree(cgeom->J));
24949566063dSJacob Faibussowitsch   PetscCall(PetscFree(cgeom->invJ));
24959566063dSJacob Faibussowitsch   PetscCall(PetscFree(cgeom->detJ));
2496c9ba7969SMatthew G. Knepley   PetscFunctionReturn(0);
2497c9ba7969SMatthew G. Knepley }
2498