xref: /petsc/include/petscfe.h (revision 9d47de495d3c23378050c1b4a410c12a375cb6c6)
1dbe77d9eSMatthew G. Knepley /*
2dbe77d9eSMatthew G. Knepley       Objects which encapsulate finite element spaces and operations
3dbe77d9eSMatthew G. Knepley */
4a4963045SJacob Faibussowitsch #pragma once
5dbe77d9eSMatthew G. Knepley #include <petscdm.h>
60ddb9b0bSMatthew G. Knepley #include <petscdt.h>
7dbe77d9eSMatthew G. Knepley #include <petscfetypes.h>
82764a2aaSMatthew G. Knepley #include <petscdstypes.h>
9660d4ad9SBarry Smith #include <petscspace.h>
10660d4ad9SBarry Smith #include <petscdualspace.h>
11dbe77d9eSMatthew G. Knepley 
12ce78bad3SBarry Smith /* MANSEC = DM */
13ac09b921SBarry Smith /* SUBMANSEC = FE */
14ac09b921SBarry Smith 
15dce8aebaSBarry Smith /*MC
16dce8aebaSBarry Smith     PetscFEGeom - Structure for geometric information for `PetscFE`
17dce8aebaSBarry Smith 
18dce8aebaSBarry Smith     Level: intermediate
19dce8aebaSBarry Smith 
205d83a8b1SBarry Smith     Note:
215d83a8b1SBarry Smith     This is a struct, not a `PetscObject`
225d83a8b1SBarry Smith 
2316a05f60SBarry Smith .seealso: `PetscFE`, `PetscFEGeomCreate()`, `PetscFEGeomDestroy()`, `PetscFEGeomGetChunk()`, `PetscFEGeomRestoreChunk()`, `PetscFEGeomGetPoint()`, `PetscFEGeomGetCellPoint()`,
24b24fb147SBarry Smith           `PetscFEGeomComplete()`, `PetscSpace`, `PetscDualSpace`
25dce8aebaSBarry Smith M*/
26ce78bad3SBarry Smith typedef struct {
27ac9d17c7SMatthew G. Knepley   // We can represent several different types of geometry, which we call modes:
28ac9d17c7SMatthew G. Knepley   //   basic:    dim == dE, only bulk data
29ac9d17c7SMatthew G. Knepley   //     These are normal dim-cells
30ac9d17c7SMatthew G. Knepley   //   embedded: dim < dE, only bulk data
31ac9d17c7SMatthew G. Knepley   //     These are dim-cells embedded in a higher dimension, as an embedded manifold
32ac9d17c7SMatthew G. Knepley   //   boundary: dim < dE, bulk and face data
33ac9d17c7SMatthew G. Knepley   //     These are dim-cells on the boundary of a dE-mesh
34ac9d17c7SMatthew G. Knepley   //   cohesive: dim < dE, bulk and face data
35ac9d17c7SMatthew G. Knepley   //     These are dim-cells in the interior of a dE-mesh
36ac9d17c7SMatthew G. Knepley   //   affine:
37ac9d17c7SMatthew G. Knepley   //     For all modes, the transforms between real and reference are affine
38ac9d17c7SMatthew G. Knepley   PetscFEGeomMode mode;     // The type of geometric data stored
39ac9d17c7SMatthew G. Knepley   PetscBool       isAffine; // Flag for affine transforms
40ac9d17c7SMatthew G. Knepley   // Sizes
41ac9d17c7SMatthew G. Knepley   PetscInt dim;       // dim: topological dimension and reference coordinate dimension
42ac9d17c7SMatthew G. Knepley   PetscInt dimEmbed;  // dE:  real coordinate dimension
43ac9d17c7SMatthew G. Knepley   PetscInt numCells;  // Nc:  Number of mesh points represented in the arrays (points are assumed to be the same DMPolytopeType)
44ac9d17c7SMatthew G. Knepley   PetscInt numPoints; // Np:  Number of evaluation points represented in the arrays
45ac9d17c7SMatthew G. Knepley   // Bulk data
46ac9d17c7SMatthew G. Knepley   const PetscReal *xi;   // xi[dim]                The first point in each cell in reference coordinates
47ac9d17c7SMatthew G. Knepley   PetscReal       *v;    // v[Nc*Np*dE]:           The first point in each cell in real coordinates
48ac9d17c7SMatthew G. Knepley   PetscReal       *J;    // J[Nc*Np*dE*dE]:        The Jacobian of the map from reference to real coordinates (if nonsquare it is completed with orthogonal columns)
49ac9d17c7SMatthew G. Knepley   PetscReal       *invJ; // invJ[Nc*Np*dE*dE]:     The inverse of the Jacobian of the map from reference to real coordinates (if nonsquare it is completed with orthogonal columns)
50ac9d17c7SMatthew G. Knepley   PetscReal       *detJ; // detJ[Nc*Np]:           The determinant of J, and if J is non-square it is the volume change
51ac9d17c7SMatthew G. Knepley   // Face data
52ac9d17c7SMatthew G. Knepley   PetscReal *n;           // n[Nc*Np*dE]:           For faces, the normal to the face in real coordinates, outward for the first supporting cell
53ac9d17c7SMatthew G. Knepley   PetscInt (*face)[4];    // face[Nc][s*2]:         For faces, the local face number (cone index) and orientation for this face in each supporting cell
54ac9d17c7SMatthew G. Knepley   PetscReal *suppJ[2];    // sJ[s][Nc*Np*dE*dE]:    For faces, the Jacobian for each supporting cell
55ac9d17c7SMatthew G. Knepley   PetscReal *suppInvJ[2]; // sInvJ[s][Nc*Np*dE*dE]: For faces, the inverse Jacobian for each supporting cell
56ac9d17c7SMatthew G. Knepley   PetscReal *suppDetJ[2]; // sdetJ[s][Nc*Np]:       For faces, the Jacobian determinant for each supporting cell
574129dba9SToby Isaac } PetscFEGeom;
584129dba9SToby Isaac 
59dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEInitializePackage(void);
60dbe77d9eSMatthew G. Knepley 
61ac9d17c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGeomCreate(PetscQuadrature, PetscInt, PetscInt, PetscFEGeomMode, PetscFEGeom **);
624129dba9SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEGeomGetChunk(PetscFEGeom *, PetscInt, PetscInt, PetscFEGeom **);
634129dba9SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEGeomRestoreChunk(PetscFEGeom *, PetscInt, PetscInt, PetscFEGeom **);
646587ee25SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGeomGetPoint(PetscFEGeom *, PetscInt, PetscInt, const PetscReal[], PetscFEGeom *);
656587ee25SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGeomGetCellPoint(PetscFEGeom *, PetscInt, PetscInt, PetscFEGeom *);
664129dba9SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEGeomComplete(PetscFEGeom *);
674129dba9SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEGeomDestroy(PetscFEGeom **);
684129dba9SToby Isaac 
69c330f8ffSToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceApply(PetscDualSpace, PetscInt, PetscReal, PetscFEGeom *, PetscInt, PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void *, PetscScalar *);
70c330f8ffSToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceApplyDefault(PetscDualSpace, PetscInt, PetscReal, PetscFEGeom *, PetscInt, PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void *, PetscScalar *);
71ebac44aeSMatthew G. Knepley 
724bee2e38SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceTransform(PetscDualSpace, PetscDualSpaceTransformType, PetscBool, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]);
73f9244615SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceTransformGradient(PetscDualSpace, PetscDualSpaceTransformType, PetscBool, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]);
74f9244615SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceTransformHessian(PetscDualSpace, PetscDualSpaceTransformType, PetscBool, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]);
752edcad52SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpacePullback(PetscDualSpace, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]);
762edcad52SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpacePushforward(PetscDualSpace, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]);
772edcad52SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpacePushforwardGradient(PetscDualSpace, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]);
78f9244615SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpacePushforwardHessian(PetscDualSpace, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]);
794bee2e38SMatthew G. Knepley 
80dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscClassId PETSCFE_CLASSID;
81dbe77d9eSMatthew G. Knepley 
820483ade4SMatthew G. Knepley /*J
830483ade4SMatthew G. Knepley   PetscFEType - String with the name of a PETSc finite element space
840483ade4SMatthew G. Knepley 
850483ade4SMatthew G. Knepley   Level: beginner
860483ade4SMatthew G. Knepley 
8716a05f60SBarry Smith   Note:
8816a05f60SBarry Smith   Currently, the classes are concerned with the implementation of element integration
890483ade4SMatthew G. Knepley 
90db781477SPatrick Sanan .seealso: `PetscFESetType()`, `PetscFE`
910483ade4SMatthew G. Knepley J*/
920483ade4SMatthew G. Knepley typedef const char *PetscFEType;
930483ade4SMatthew G. Knepley #define PETSCFEBASIC     "basic"
940483ade4SMatthew G. Knepley #define PETSCFEOPENCL    "opencl"
95aaf1837cSMatthew G. Knepley #define PETSCFECOMPOSITE "composite"
962dce792eSToby Isaac #define PETSCFEVECTOR    "vector"
970483ade4SMatthew G. Knepley 
980483ade4SMatthew G. Knepley PETSC_EXTERN PetscFunctionList PetscFEList;
99568c31deSMatthew G. Knepley PETSC_EXTERN PetscErrorCode    PetscFECreate(MPI_Comm, PetscFE *);
100568c31deSMatthew G. Knepley PETSC_EXTERN PetscErrorCode    PetscFEDestroy(PetscFE *);
1010483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode    PetscFESetType(PetscFE, PetscFEType);
1020483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode    PetscFEGetType(PetscFE, PetscFEType *);
1030483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode    PetscFESetUp(PetscFE);
1040483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode    PetscFESetFromOptions(PetscFE);
105fe2efc57SMark PETSC_EXTERN PetscErrorCode    PetscFEViewFromOptions(PetscFE, PetscObject, const char[]);
1063f6b16c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode    PetscFESetName(PetscFE, const char[]);
1072dce792eSToby Isaac PETSC_EXTERN PetscErrorCode    PetscFECreateVector(PetscFE, PetscInt, PetscBool, PetscBool, PetscFE *);
1088aec7d55SBarry Smith 
109fd41f350SJed Brown PETSC_EXTERN PetscErrorCode PetscFEView(PetscFE, PetscViewer);
1100483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFERegister(const char[], PetscErrorCode (*)(PetscFE));
1110483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFERegisterDestroy(void);
112e8d98e54SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateDefault(MPI_Comm, PetscInt, PetscInt, PetscBool, const char[], PetscInt, PetscFE *);
1132df84da0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateByCell(MPI_Comm, PetscInt, PetscInt, DMPolytopeType, const char[], PetscInt, PetscFE *);
114e703855dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateLagrange(MPI_Comm, PetscInt, PetscInt, PetscBool, PetscInt, PetscInt, PetscFE *);
1152df84da0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateLagrangeByCell(MPI_Comm, PetscInt, PetscInt, DMPolytopeType, PetscInt, PetscInt, PetscFE *);
1167c48043bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateFromSpaces(PetscSpace, PetscDualSpace, PetscQuadrature, PetscQuadrature, PetscFE *);
117bb4b53efSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFELimitDegree(PetscFE, PetscInt, PetscInt, PetscFE *);
1184c712d99Sksagiyam PETSC_EXTERN PetscErrorCode PetscFECreateBrokenElement(PetscFE, PetscFE *);
1199a1a3eb8SMatthew G. Knepley 
1209a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetDimension(PetscFE, PetscInt *);
1210ddb9b0bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetSpatialDimension(PetscFE, PetscInt *);
1229a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetNumComponents(PetscFE, PetscInt);
1239a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetNumComponents(PetscFE, PetscInt *);
1249a769b6bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetTileSizes(PetscFE, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
1259a769b6bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetTileSizes(PetscFE, PetscInt, PetscInt, PetscInt, PetscInt);
1269a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetBasisSpace(PetscFE, PetscSpace);
1279a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetBasisSpace(PetscFE, PetscSpace *);
1289a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetDualSpace(PetscFE, PetscDualSpace);
1299a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetDualSpace(PetscFE, PetscDualSpace *);
130bfa639d9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetQuadrature(PetscFE, PetscQuadrature);
131bfa639d9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetQuadrature(PetscFE, PetscQuadrature *);
1321bafc85dSSander Arens PETSC_EXTERN PetscErrorCode PetscFESetFaceQuadrature(PetscFE, PetscQuadrature);
1331bafc85dSSander Arens PETSC_EXTERN PetscErrorCode PetscFEGetFaceQuadrature(PetscFE, PetscQuadrature *);
134989fa639SBrad Aagaard PETSC_EXTERN PetscErrorCode PetscFEExpandFaceQuadrature(PetscFE, PetscQuadrature, PetscQuadrature *);
1355dc5c000SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECopyQuadrature(PetscFE, PetscFE);
1360ddb9b0bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetNumDof(PetscFE, const PetscInt **);
137ef0bb6c7SMatthew G. Knepley 
138ef0bb6c7SMatthew G. Knepley /* TODO: Need a function to reuse the memory when retabulating the same FE at different points */
139f9244615SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetCellTabulation(PetscFE, PetscInt, PetscTabulation *);
140f9244615SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetFaceTabulation(PetscFE, PetscInt, PetscTabulation *);
141ef0bb6c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetFaceCentroidTabulation(PetscFE, PetscTabulation *);
142ef0bb6c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateTabulation(PetscFE, PetscInt, PetscInt, const PetscReal[], PetscInt, PetscTabulation *);
143ef0bb6c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEComputeTabulation(PetscFE, PetscInt, const PetscReal[], PetscInt, PetscTabulation);
144ef0bb6c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscTabulationDestroy(PetscTabulation *);
145ef0bb6c7SMatthew G. Knepley 
146bceba477SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFERefine(PetscFE, PetscFE *);
147228cfb18SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetHeightSubspace(PetscFE, PetscInt, PetscFE *);
148a0845e3aSMatthew G. Knepley 
149c9ba7969SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateCellGeometry(PetscFE, PetscQuadrature, PetscFEGeom *);
150c9ba7969SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEDestroyCellGeometry(PetscFE, PetscFEGeom *);
1512edcad52SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEPushforward(PetscFE, PetscFEGeom *, PetscInt, PetscScalar[]);
1522edcad52SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEPushforwardGradient(PetscFE, PetscFEGeom *, PetscInt, PetscScalar[]);
153f9244615SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEPushforwardHessian(PetscFE, PetscFEGeom *, PetscInt, PetscScalar[]);
1544bee2e38SMatthew G. Knepley 
1554bee2e38SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrate(PetscDS, PetscInt, PetscInt, PetscFEGeom *, const PetscScalar[], PetscDS, const PetscScalar[], PetscScalar[]);
1569371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PetscFEIntegrateBd(PetscDS, PetscInt, void (*)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt, PetscFEGeom *, const PetscScalar[], PetscDS, const PetscScalar[], PetscScalar[]);
15706ad1575SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateResidual(PetscDS, PetscFormKey, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscScalar[]);
15806ad1575SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateBdResidual(PetscDS, PetscWeakForm, PetscFormKey, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscScalar[]);
159989fa639SBrad Aagaard PETSC_EXTERN PetscErrorCode PetscFEIntegrateHybridResidual(PetscDS, PetscDS, PetscFormKey, PetscInt, PetscInt, PetscFEGeom *, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscScalar[]);
1604561e6c9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateJacobian(PetscDS, PetscDS, PetscFEJacobianType, PetscFormKey, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscReal, PetscScalar[]);
161e3d591f2SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateBdJacobian(PetscDS, PetscWeakForm, PetscFEJacobianType, PetscFormKey, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscReal, PetscScalar[]);
162989fa639SBrad Aagaard PETSC_EXTERN PetscErrorCode PetscFEIntegrateHybridJacobian(PetscDS, PetscDS, PetscFEJacobianType, PetscFormKey, PetscInt, PetscInt, PetscFEGeom *, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscReal, PetscScalar[]);
163a0845e3aSMatthew G. Knepley 
16489710940SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECompositeGetMapping(PetscFE, PetscInt *, const PetscReal *[], const PetscReal *[], const PetscReal *[]);
16589710940SMatthew G. Knepley 
1662f5fb066SToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreateHeightTrace(PetscFE, PetscInt, PetscFE *);
1672f5fb066SToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreatePointTrace(PetscFE, PetscInt, PetscFE *);
1682f5fb066SToby Isaac 
169855cd083SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEOpenCLSetRealType(PetscFE, PetscDataType);
170855cd083SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEOpenCLGetRealType(PetscFE, PetscDataType *);
171d2b2dc1eSMatthew G. Knepley 
172d2b2dc1eSMatthew G. Knepley #ifdef PETSC_HAVE_LIBCEED
173d2b2dc1eSMatthew G. Knepley 
174*beceaeb6SBarry Smith   #if !defined(PLEXFE_QFUNCTION)
175d2b2dc1eSMatthew G. Knepley     #define PLEXFE_QFUNCTION(fname, f0_name, f1_name) \
1762a8381b2SBarry Smith       CEED_QFUNCTION(PlexQFunction##fname)(PetscCtx ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) \
177d2b2dc1eSMatthew G. Knepley       { \
178d2b2dc1eSMatthew G. Knepley         const CeedScalar *u = in[0], *du = in[1], *qdata = in[2]; \
179d2b2dc1eSMatthew G. Knepley         CeedScalar       *v = out[0], *dv = out[1]; \
180d2b2dc1eSMatthew G. Knepley         const PetscInt    Nc   = 1; \
181d2b2dc1eSMatthew G. Knepley         const PetscInt    cdim = 2; \
182d2b2dc1eSMatthew G. Knepley \
183d2b2dc1eSMatthew G. Knepley         CeedPragmaSIMD for (CeedInt i = 0; i < Q; ++i) \
184d2b2dc1eSMatthew G. Knepley         { \
185d2b2dc1eSMatthew G. Knepley           const PetscInt   uOff[2]    = {0, Nc}; \
186d2b2dc1eSMatthew G. Knepley           const PetscInt   uOff_x[2]  = {0, Nc * cdim}; \
187d2b2dc1eSMatthew G. Knepley           const CeedScalar x[2]       = {qdata[i + Q * 1], qdata[i + Q * 2]}; \
188d2b2dc1eSMatthew G. Knepley           const CeedScalar invJ[2][2] = { \
189d2b2dc1eSMatthew G. Knepley             {qdata[i + Q * 3], qdata[i + Q * 5]}, \
190d2b2dc1eSMatthew G. Knepley             {qdata[i + Q * 4], qdata[i + Q * 6]} \
191d2b2dc1eSMatthew G. Knepley           }; \
192d2b2dc1eSMatthew G. Knepley           const CeedScalar u_x[2] = {invJ[0][0] * du[i + Q * 0] + invJ[1][0] * du[i + Q * 1], invJ[0][1] * du[i + Q * 0] + invJ[1][1] * du[i + Q * 1]}; \
193d2b2dc1eSMatthew G. Knepley           PetscScalar      f0[Nc]; \
194d2b2dc1eSMatthew G. Knepley           PetscScalar      f1[Nc * cdim]; \
195d2b2dc1eSMatthew G. Knepley \
196d2b2dc1eSMatthew G. Knepley           for (PetscInt k = 0; k < Nc; ++k) f0[k] = 0; \
197d2b2dc1eSMatthew G. Knepley           for (PetscInt k = 0; k < Nc * cdim; ++k) f1[k] = 0; \
198d2b2dc1eSMatthew G. Knepley           f0_name(2, 1, 0, uOff, uOff_x, u, NULL, u_x, NULL, NULL, NULL, NULL, NULL, 0.0, x, 0, NULL, f0); \
199d2b2dc1eSMatthew G. Knepley           f1_name(2, 1, 0, uOff, uOff_x, u, NULL, u_x, NULL, NULL, NULL, NULL, NULL, 0.0, x, 0, NULL, f1); \
200d2b2dc1eSMatthew G. Knepley \
201d2b2dc1eSMatthew G. Knepley           dv[i + Q * 0] = qdata[i + Q * 0] * (invJ[0][0] * f1[0] + invJ[0][1] * f1[1]); \
202d2b2dc1eSMatthew G. Knepley           dv[i + Q * 1] = qdata[i + Q * 0] * (invJ[1][0] * f1[0] + invJ[1][1] * f1[1]); \
203d2b2dc1eSMatthew G. Knepley           v[i]          = qdata[i + Q * 0] * f0[0]; \
204d2b2dc1eSMatthew G. Knepley         } \
205d2b2dc1eSMatthew G. Knepley         return CEED_ERROR_SUCCESS; \
206d2b2dc1eSMatthew G. Knepley       }
207d2b2dc1eSMatthew G. Knepley   #endif
208d2b2dc1eSMatthew G. Knepley 
209d2b2dc1eSMatthew G. Knepley #else
210d2b2dc1eSMatthew G. Knepley 
211*beceaeb6SBarry Smith   #if !defined(PLEXFE_QFUNCTION)
212d2b2dc1eSMatthew G. Knepley     #define PLEXFE_QFUNCTION(fname, f0_name, f1_name)
213d2b2dc1eSMatthew G. Knepley   #endif
214d2b2dc1eSMatthew G. Knepley 
215d2b2dc1eSMatthew G. Knepley #endif
216