1dbe77d9eSMatthew G. Knepley /* 2dbe77d9eSMatthew G. Knepley Objects which encapsulate finite element spaces and operations 3dbe77d9eSMatthew G. Knepley */ 4a4963045SJacob Faibussowitsch #pragma once 5*bb4b53efSMatthew G. Knepley #include "petscmacros.h" 6dbe77d9eSMatthew G. Knepley #include <petscdm.h> 70ddb9b0bSMatthew G. Knepley #include <petscdt.h> 8dbe77d9eSMatthew G. Knepley #include <petscfetypes.h> 92764a2aaSMatthew G. Knepley #include <petscdstypes.h> 10660d4ad9SBarry Smith #include <petscspace.h> 11660d4ad9SBarry Smith #include <petscdualspace.h> 12dbe77d9eSMatthew G. Knepley 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*/ 264129dba9SToby Isaac typedef struct _n_PetscFEGeom { 274129dba9SToby Isaac const PetscReal *xi; 2841418987SMatthew G. Knepley PetscReal *v; /* v[Nc*Np*dE]: The first point in each each in real coordinates */ 2941418987SMatthew 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) */ 3041418987SMatthew 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) */ 3141418987SMatthew G. Knepley PetscReal *detJ; /* detJ[Nc*Np]: The determinant of J, and if it is non-square its the volume change */ 32f15274beSMatthew Knepley PetscReal *n; /* n[Nc*Np*dE]: For faces, the normal to the face in real coordinates, outward for the first supporting cell */ 330e18dc48SMatthew 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 s */ 3441418987SMatthew G. Knepley PetscReal *suppJ[2]; /* sJ[s][Nc*Np*dE*dE]: For faces, the Jacobian for each supporting cell s */ 3541418987SMatthew G. Knepley PetscReal *suppInvJ[2]; /* sInvJ[s][Nc*Np*dE*dE]: For faces, the inverse Jacobian for each supporting cell s */ 3641418987SMatthew G. Knepley PetscReal *suppDetJ[2]; /* sInvJ[s][Nc*Np*dE*dE]: For faces, the Jacobian determinant for each supporting cell s */ 3762bd480fSMatthew G. Knepley PetscInt dim; /* dim: Topological dimension */ 3862bd480fSMatthew G. Knepley PetscInt dimEmbed; /* dE: coordinate dimension */ 3962bd480fSMatthew G. Knepley PetscInt numCells; /* Nc: Number of mesh points represented in the arrays */ 4062bd480fSMatthew G. Knepley PetscInt numPoints; /* Np: Number of evaluation points represented in the arrays */ 4141418987SMatthew G. Knepley PetscBool isAffine; /* Flag for affine transforms */ 425fedec97SMatthew G. Knepley PetscBool isCohesive; /* Flag for a cohesive cell */ 434129dba9SToby Isaac } PetscFEGeom; 444129dba9SToby Isaac 45dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEInitializePackage(void); 46dbe77d9eSMatthew G. Knepley 474129dba9SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEGeomCreate(PetscQuadrature, PetscInt, PetscInt, PetscBool, PetscFEGeom **); 484129dba9SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEGeomGetChunk(PetscFEGeom *, PetscInt, PetscInt, PetscFEGeom **); 494129dba9SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEGeomRestoreChunk(PetscFEGeom *, PetscInt, PetscInt, PetscFEGeom **); 506587ee25SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGeomGetPoint(PetscFEGeom *, PetscInt, PetscInt, const PetscReal[], PetscFEGeom *); 516587ee25SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGeomGetCellPoint(PetscFEGeom *, PetscInt, PetscInt, PetscFEGeom *); 524129dba9SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEGeomComplete(PetscFEGeom *); 534129dba9SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEGeomDestroy(PetscFEGeom **); 544129dba9SToby Isaac 55c330f8ffSToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceApply(PetscDualSpace, PetscInt, PetscReal, PetscFEGeom *, PetscInt, PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void *, PetscScalar *); 56c330f8ffSToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceApplyDefault(PetscDualSpace, PetscInt, PetscReal, PetscFEGeom *, PetscInt, PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void *, PetscScalar *); 57ebac44aeSMatthew G. Knepley 584bee2e38SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceTransform(PetscDualSpace, PetscDualSpaceTransformType, PetscBool, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]); 59f9244615SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceTransformGradient(PetscDualSpace, PetscDualSpaceTransformType, PetscBool, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]); 60f9244615SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceTransformHessian(PetscDualSpace, PetscDualSpaceTransformType, PetscBool, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]); 612edcad52SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpacePullback(PetscDualSpace, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]); 622edcad52SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpacePushforward(PetscDualSpace, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]); 632edcad52SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpacePushforwardGradient(PetscDualSpace, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]); 64f9244615SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpacePushforwardHessian(PetscDualSpace, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]); 654bee2e38SMatthew G. Knepley 66dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscClassId PETSCFE_CLASSID; 67dbe77d9eSMatthew G. Knepley 680483ade4SMatthew G. Knepley /*J 690483ade4SMatthew G. Knepley PetscFEType - String with the name of a PETSc finite element space 700483ade4SMatthew G. Knepley 710483ade4SMatthew G. Knepley Level: beginner 720483ade4SMatthew G. Knepley 7316a05f60SBarry Smith Note: 7416a05f60SBarry Smith Currently, the classes are concerned with the implementation of element integration 750483ade4SMatthew G. Knepley 76db781477SPatrick Sanan .seealso: `PetscFESetType()`, `PetscFE` 770483ade4SMatthew G. Knepley J*/ 780483ade4SMatthew G. Knepley typedef const char *PetscFEType; 790483ade4SMatthew G. Knepley #define PETSCFEBASIC "basic" 800483ade4SMatthew G. Knepley #define PETSCFEOPENCL "opencl" 81aaf1837cSMatthew G. Knepley #define PETSCFECOMPOSITE "composite" 822dce792eSToby Isaac #define PETSCFEVECTOR "vector" 830483ade4SMatthew G. Knepley 840483ade4SMatthew G. Knepley PETSC_EXTERN PetscFunctionList PetscFEList; 85568c31deSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreate(MPI_Comm, PetscFE *); 86568c31deSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEDestroy(PetscFE *); 870483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetType(PetscFE, PetscFEType); 880483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetType(PetscFE, PetscFEType *); 890483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetUp(PetscFE); 900483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetFromOptions(PetscFE); 91fe2efc57SMark PETSC_EXTERN PetscErrorCode PetscFEViewFromOptions(PetscFE, PetscObject, const char[]); 923f6b16c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetName(PetscFE, const char[]); 932dce792eSToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreateVector(PetscFE, PetscInt, PetscBool, PetscBool, PetscFE *); 948aec7d55SBarry Smith 95fd41f350SJed Brown PETSC_EXTERN PetscErrorCode PetscFEView(PetscFE, PetscViewer); 960483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFERegister(const char[], PetscErrorCode (*)(PetscFE)); 970483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFERegisterDestroy(void); 98e8d98e54SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateDefault(MPI_Comm, PetscInt, PetscInt, PetscBool, const char[], PetscInt, PetscFE *); 992df84da0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateByCell(MPI_Comm, PetscInt, PetscInt, DMPolytopeType, const char[], PetscInt, PetscFE *); 100e703855dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateLagrange(MPI_Comm, PetscInt, PetscInt, PetscBool, PetscInt, PetscInt, PetscFE *); 1012df84da0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateLagrangeByCell(MPI_Comm, PetscInt, PetscInt, DMPolytopeType, PetscInt, PetscInt, PetscFE *); 1027c48043bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateFromSpaces(PetscSpace, PetscDualSpace, PetscQuadrature, PetscQuadrature, PetscFE *); 103*bb4b53efSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFELimitDegree(PetscFE, PetscInt, PetscInt, PetscFE *); 1049a1a3eb8SMatthew G. Knepley 1059a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetDimension(PetscFE, PetscInt *); 1060ddb9b0bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetSpatialDimension(PetscFE, PetscInt *); 1079a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetNumComponents(PetscFE, PetscInt); 1089a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetNumComponents(PetscFE, PetscInt *); 1099a769b6bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetTileSizes(PetscFE, PetscInt *, PetscInt *, PetscInt *, PetscInt *); 1109a769b6bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetTileSizes(PetscFE, PetscInt, PetscInt, PetscInt, PetscInt); 1119a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetBasisSpace(PetscFE, PetscSpace); 1129a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetBasisSpace(PetscFE, PetscSpace *); 1139a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetDualSpace(PetscFE, PetscDualSpace); 1149a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetDualSpace(PetscFE, PetscDualSpace *); 115bfa639d9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetQuadrature(PetscFE, PetscQuadrature); 116bfa639d9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetQuadrature(PetscFE, PetscQuadrature *); 1171bafc85dSSander Arens PETSC_EXTERN PetscErrorCode PetscFESetFaceQuadrature(PetscFE, PetscQuadrature); 1181bafc85dSSander Arens PETSC_EXTERN PetscErrorCode PetscFEGetFaceQuadrature(PetscFE, PetscQuadrature *); 1195dc5c000SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECopyQuadrature(PetscFE, PetscFE); 1200ddb9b0bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetNumDof(PetscFE, const PetscInt **); 121ef0bb6c7SMatthew G. Knepley 122ef0bb6c7SMatthew G. Knepley /* TODO: Need a function to reuse the memory when retabulating the same FE at different points */ 123f9244615SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetCellTabulation(PetscFE, PetscInt, PetscTabulation *); 124f9244615SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetFaceTabulation(PetscFE, PetscInt, PetscTabulation *); 125ef0bb6c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetFaceCentroidTabulation(PetscFE, PetscTabulation *); 126ef0bb6c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateTabulation(PetscFE, PetscInt, PetscInt, const PetscReal[], PetscInt, PetscTabulation *); 127ef0bb6c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEComputeTabulation(PetscFE, PetscInt, const PetscReal[], PetscInt, PetscTabulation); 128ef0bb6c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscTabulationDestroy(PetscTabulation *); 129ef0bb6c7SMatthew G. Knepley 130bceba477SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFERefine(PetscFE, PetscFE *); 131228cfb18SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetHeightSubspace(PetscFE, PetscInt, PetscFE *); 132a0845e3aSMatthew G. Knepley 133c9ba7969SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateCellGeometry(PetscFE, PetscQuadrature, PetscFEGeom *); 134c9ba7969SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEDestroyCellGeometry(PetscFE, PetscFEGeom *); 1352edcad52SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEPushforward(PetscFE, PetscFEGeom *, PetscInt, PetscScalar[]); 1362edcad52SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEPushforwardGradient(PetscFE, PetscFEGeom *, PetscInt, PetscScalar[]); 137f9244615SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEPushforwardHessian(PetscFE, PetscFEGeom *, PetscInt, PetscScalar[]); 1384bee2e38SMatthew G. Knepley 1394bee2e38SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrate(PetscDS, PetscInt, PetscInt, PetscFEGeom *, const PetscScalar[], PetscDS, const PetscScalar[], PetscScalar[]); 1409371c9d4SSatish 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[]); 14106ad1575SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateResidual(PetscDS, PetscFormKey, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscScalar[]); 14206ad1575SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateBdResidual(PetscDS, PetscWeakForm, PetscFormKey, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscScalar[]); 14307218a29SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateHybridResidual(PetscDS, PetscDS, PetscFormKey, PetscInt, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscScalar[]); 14406ad1575SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateJacobian(PetscDS, PetscFEJacobianType, PetscFormKey, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscReal, PetscScalar[]); 145e3d591f2SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateBdJacobian(PetscDS, PetscWeakForm, PetscFEJacobianType, PetscFormKey, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscReal, PetscScalar[]); 14607218a29SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateHybridJacobian(PetscDS, PetscDS, PetscFEJacobianType, PetscFormKey, PetscInt, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscReal, PetscScalar[]); 147a0845e3aSMatthew G. Knepley 14889710940SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECompositeGetMapping(PetscFE, PetscInt *, const PetscReal *[], const PetscReal *[], const PetscReal *[]); 14989710940SMatthew G. Knepley 1502f5fb066SToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreateHeightTrace(PetscFE, PetscInt, PetscFE *); 1512f5fb066SToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreatePointTrace(PetscFE, PetscInt, PetscFE *); 1522f5fb066SToby Isaac 153855cd083SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEOpenCLSetRealType(PetscFE, PetscDataType); 154855cd083SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEOpenCLGetRealType(PetscFE, PetscDataType *); 155d2b2dc1eSMatthew G. Knepley 156d2b2dc1eSMatthew G. Knepley #ifdef PETSC_HAVE_LIBCEED 157d2b2dc1eSMatthew G. Knepley 158d2b2dc1eSMatthew G. Knepley #ifndef PLEXFE_QFUNCTION 159d2b2dc1eSMatthew G. Knepley #define PLEXFE_QFUNCTION(fname, f0_name, f1_name) \ 160d2b2dc1eSMatthew G. Knepley CEED_QFUNCTION(PlexQFunction##fname)(void *ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) \ 161d2b2dc1eSMatthew G. Knepley { \ 162d2b2dc1eSMatthew G. Knepley const CeedScalar *u = in[0], *du = in[1], *qdata = in[2]; \ 163d2b2dc1eSMatthew G. Knepley CeedScalar *v = out[0], *dv = out[1]; \ 164d2b2dc1eSMatthew G. Knepley const PetscInt Nc = 1; \ 165d2b2dc1eSMatthew G. Knepley const PetscInt cdim = 2; \ 166d2b2dc1eSMatthew G. Knepley \ 167d2b2dc1eSMatthew G. Knepley CeedPragmaSIMD for (CeedInt i = 0; i < Q; ++i) \ 168d2b2dc1eSMatthew G. Knepley { \ 169d2b2dc1eSMatthew G. Knepley const PetscInt uOff[2] = {0, Nc}; \ 170d2b2dc1eSMatthew G. Knepley const PetscInt uOff_x[2] = {0, Nc * cdim}; \ 171d2b2dc1eSMatthew G. Knepley const CeedScalar x[2] = {qdata[i + Q * 1], qdata[i + Q * 2]}; \ 172d2b2dc1eSMatthew G. Knepley const CeedScalar invJ[2][2] = { \ 173d2b2dc1eSMatthew G. Knepley {qdata[i + Q * 3], qdata[i + Q * 5]}, \ 174d2b2dc1eSMatthew G. Knepley {qdata[i + Q * 4], qdata[i + Q * 6]} \ 175d2b2dc1eSMatthew G. Knepley }; \ 176d2b2dc1eSMatthew 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]}; \ 177d2b2dc1eSMatthew G. Knepley PetscScalar f0[Nc]; \ 178d2b2dc1eSMatthew G. Knepley PetscScalar f1[Nc * cdim]; \ 179d2b2dc1eSMatthew G. Knepley \ 180d2b2dc1eSMatthew G. Knepley for (PetscInt k = 0; k < Nc; ++k) f0[k] = 0; \ 181d2b2dc1eSMatthew G. Knepley for (PetscInt k = 0; k < Nc * cdim; ++k) f1[k] = 0; \ 182d2b2dc1eSMatthew 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); \ 183d2b2dc1eSMatthew 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); \ 184d2b2dc1eSMatthew G. Knepley \ 185d2b2dc1eSMatthew G. Knepley dv[i + Q * 0] = qdata[i + Q * 0] * (invJ[0][0] * f1[0] + invJ[0][1] * f1[1]); \ 186d2b2dc1eSMatthew G. Knepley dv[i + Q * 1] = qdata[i + Q * 0] * (invJ[1][0] * f1[0] + invJ[1][1] * f1[1]); \ 187d2b2dc1eSMatthew G. Knepley v[i] = qdata[i + Q * 0] * f0[0]; \ 188d2b2dc1eSMatthew G. Knepley } \ 189d2b2dc1eSMatthew G. Knepley return CEED_ERROR_SUCCESS; \ 190d2b2dc1eSMatthew G. Knepley } 191d2b2dc1eSMatthew G. Knepley #endif 192d2b2dc1eSMatthew G. Knepley 193d2b2dc1eSMatthew G. Knepley #else 194d2b2dc1eSMatthew G. Knepley 195d2b2dc1eSMatthew G. Knepley #ifndef PLEXFE_QFUNCTION 196d2b2dc1eSMatthew G. Knepley #define PLEXFE_QFUNCTION(fname, f0_name, f1_name) 197d2b2dc1eSMatthew G. Knepley #endif 198d2b2dc1eSMatthew G. Knepley 199d2b2dc1eSMatthew G. Knepley #endif 200