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 12ac09b921SBarry Smith /* SUBMANSEC = FE */ 13ac09b921SBarry Smith 14dce8aebaSBarry Smith /*MC 15dce8aebaSBarry Smith PetscFEGeom - Structure for geometric information for `PetscFE` 16dce8aebaSBarry Smith 17dce8aebaSBarry Smith Level: intermediate 18dce8aebaSBarry Smith 19*5d83a8b1SBarry Smith Note: 20*5d83a8b1SBarry Smith This is a struct, not a `PetscObject` 21*5d83a8b1SBarry Smith 2216a05f60SBarry Smith .seealso: `PetscFE`, `PetscFEGeomCreate()`, `PetscFEGeomDestroy()`, `PetscFEGeomGetChunk()`, `PetscFEGeomRestoreChunk()`, `PetscFEGeomGetPoint()`, `PetscFEGeomGetCellPoint()`, 23b24fb147SBarry Smith `PetscFEGeomComplete()`, `PetscSpace`, `PetscDualSpace` 24dce8aebaSBarry Smith M*/ 254129dba9SToby Isaac typedef struct _n_PetscFEGeom { 264129dba9SToby Isaac const PetscReal *xi; 2741418987SMatthew G. Knepley PetscReal *v; /* v[Nc*Np*dE]: The first point in each each in real coordinates */ 2841418987SMatthew 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) */ 2941418987SMatthew 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) */ 3041418987SMatthew G. Knepley PetscReal *detJ; /* detJ[Nc*Np]: The determinant of J, and if it is non-square its the volume change */ 31f15274beSMatthew Knepley PetscReal *n; /* n[Nc*Np*dE]: For faces, the normal to the face in real coordinates, outward for the first supporting cell */ 320e18dc48SMatthew 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 */ 3341418987SMatthew G. Knepley PetscReal *suppJ[2]; /* sJ[s][Nc*Np*dE*dE]: For faces, the Jacobian for each supporting cell s */ 3441418987SMatthew G. Knepley PetscReal *suppInvJ[2]; /* sInvJ[s][Nc*Np*dE*dE]: For faces, the inverse Jacobian for each supporting cell s */ 3541418987SMatthew G. Knepley PetscReal *suppDetJ[2]; /* sInvJ[s][Nc*Np*dE*dE]: For faces, the Jacobian determinant for each supporting cell s */ 3662bd480fSMatthew G. Knepley PetscInt dim; /* dim: Topological dimension */ 3762bd480fSMatthew G. Knepley PetscInt dimEmbed; /* dE: coordinate dimension */ 3862bd480fSMatthew G. Knepley PetscInt numCells; /* Nc: Number of mesh points represented in the arrays */ 3962bd480fSMatthew G. Knepley PetscInt numPoints; /* Np: Number of evaluation points represented in the arrays */ 4041418987SMatthew G. Knepley PetscBool isAffine; /* Flag for affine transforms */ 415fedec97SMatthew G. Knepley PetscBool isCohesive; /* Flag for a cohesive cell */ 424129dba9SToby Isaac } PetscFEGeom; 434129dba9SToby Isaac 44dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEInitializePackage(void); 45dbe77d9eSMatthew G. Knepley 464129dba9SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEGeomCreate(PetscQuadrature, PetscInt, PetscInt, PetscBool, PetscFEGeom **); 474129dba9SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEGeomGetChunk(PetscFEGeom *, PetscInt, PetscInt, PetscFEGeom **); 484129dba9SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEGeomRestoreChunk(PetscFEGeom *, PetscInt, PetscInt, PetscFEGeom **); 496587ee25SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGeomGetPoint(PetscFEGeom *, PetscInt, PetscInt, const PetscReal[], PetscFEGeom *); 506587ee25SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGeomGetCellPoint(PetscFEGeom *, PetscInt, PetscInt, PetscFEGeom *); 514129dba9SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEGeomComplete(PetscFEGeom *); 524129dba9SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEGeomDestroy(PetscFEGeom **); 534129dba9SToby Isaac 54c330f8ffSToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceApply(PetscDualSpace, PetscInt, PetscReal, PetscFEGeom *, PetscInt, PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void *, PetscScalar *); 55c330f8ffSToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceApplyDefault(PetscDualSpace, PetscInt, PetscReal, PetscFEGeom *, PetscInt, PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void *, PetscScalar *); 56ebac44aeSMatthew G. Knepley 574bee2e38SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceTransform(PetscDualSpace, PetscDualSpaceTransformType, PetscBool, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]); 58f9244615SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceTransformGradient(PetscDualSpace, PetscDualSpaceTransformType, PetscBool, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]); 59f9244615SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceTransformHessian(PetscDualSpace, PetscDualSpaceTransformType, PetscBool, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]); 602edcad52SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpacePullback(PetscDualSpace, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]); 612edcad52SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpacePushforward(PetscDualSpace, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]); 622edcad52SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpacePushforwardGradient(PetscDualSpace, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]); 63f9244615SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpacePushforwardHessian(PetscDualSpace, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]); 644bee2e38SMatthew G. Knepley 65dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscClassId PETSCFE_CLASSID; 66dbe77d9eSMatthew G. Knepley 670483ade4SMatthew G. Knepley /*J 680483ade4SMatthew G. Knepley PetscFEType - String with the name of a PETSc finite element space 690483ade4SMatthew G. Knepley 700483ade4SMatthew G. Knepley Level: beginner 710483ade4SMatthew G. Knepley 7216a05f60SBarry Smith Note: 7316a05f60SBarry Smith Currently, the classes are concerned with the implementation of element integration 740483ade4SMatthew G. Knepley 75db781477SPatrick Sanan .seealso: `PetscFESetType()`, `PetscFE` 760483ade4SMatthew G. Knepley J*/ 770483ade4SMatthew G. Knepley typedef const char *PetscFEType; 780483ade4SMatthew G. Knepley #define PETSCFEBASIC "basic" 790483ade4SMatthew G. Knepley #define PETSCFEOPENCL "opencl" 80aaf1837cSMatthew G. Knepley #define PETSCFECOMPOSITE "composite" 812dce792eSToby Isaac #define PETSCFEVECTOR "vector" 820483ade4SMatthew G. Knepley 830483ade4SMatthew G. Knepley PETSC_EXTERN PetscFunctionList PetscFEList; 84568c31deSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreate(MPI_Comm, PetscFE *); 85568c31deSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEDestroy(PetscFE *); 860483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetType(PetscFE, PetscFEType); 870483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetType(PetscFE, PetscFEType *); 880483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetUp(PetscFE); 890483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetFromOptions(PetscFE); 90fe2efc57SMark PETSC_EXTERN PetscErrorCode PetscFEViewFromOptions(PetscFE, PetscObject, const char[]); 913f6b16c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetName(PetscFE, const char[]); 922dce792eSToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreateVector(PetscFE, PetscInt, PetscBool, PetscBool, PetscFE *); 938aec7d55SBarry Smith 94fd41f350SJed Brown PETSC_EXTERN PetscErrorCode PetscFEView(PetscFE, PetscViewer); 950483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFERegister(const char[], PetscErrorCode (*)(PetscFE)); 960483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFERegisterDestroy(void); 97e8d98e54SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateDefault(MPI_Comm, PetscInt, PetscInt, PetscBool, const char[], PetscInt, PetscFE *); 982df84da0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateByCell(MPI_Comm, PetscInt, PetscInt, DMPolytopeType, const char[], PetscInt, PetscFE *); 99e703855dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateLagrange(MPI_Comm, PetscInt, PetscInt, PetscBool, PetscInt, PetscInt, PetscFE *); 1002df84da0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateLagrangeByCell(MPI_Comm, PetscInt, PetscInt, DMPolytopeType, PetscInt, PetscInt, PetscFE *); 1017c48043bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateFromSpaces(PetscSpace, PetscDualSpace, PetscQuadrature, PetscQuadrature, PetscFE *); 1029a1a3eb8SMatthew G. Knepley 1039a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetDimension(PetscFE, PetscInt *); 1040ddb9b0bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetSpatialDimension(PetscFE, PetscInt *); 1059a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetNumComponents(PetscFE, PetscInt); 1069a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetNumComponents(PetscFE, PetscInt *); 1079a769b6bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetTileSizes(PetscFE, PetscInt *, PetscInt *, PetscInt *, PetscInt *); 1089a769b6bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetTileSizes(PetscFE, PetscInt, PetscInt, PetscInt, PetscInt); 1099a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetBasisSpace(PetscFE, PetscSpace); 1109a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetBasisSpace(PetscFE, PetscSpace *); 1119a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetDualSpace(PetscFE, PetscDualSpace); 1129a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetDualSpace(PetscFE, PetscDualSpace *); 113bfa639d9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetQuadrature(PetscFE, PetscQuadrature); 114bfa639d9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetQuadrature(PetscFE, PetscQuadrature *); 1151bafc85dSSander Arens PETSC_EXTERN PetscErrorCode PetscFESetFaceQuadrature(PetscFE, PetscQuadrature); 1161bafc85dSSander Arens PETSC_EXTERN PetscErrorCode PetscFEGetFaceQuadrature(PetscFE, PetscQuadrature *); 1175dc5c000SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECopyQuadrature(PetscFE, PetscFE); 1180ddb9b0bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetNumDof(PetscFE, const PetscInt **); 119ef0bb6c7SMatthew G. Knepley 120ef0bb6c7SMatthew G. Knepley /* TODO: Need a function to reuse the memory when retabulating the same FE at different points */ 121f9244615SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetCellTabulation(PetscFE, PetscInt, PetscTabulation *); 122f9244615SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetFaceTabulation(PetscFE, PetscInt, PetscTabulation *); 123ef0bb6c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetFaceCentroidTabulation(PetscFE, PetscTabulation *); 124ef0bb6c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateTabulation(PetscFE, PetscInt, PetscInt, const PetscReal[], PetscInt, PetscTabulation *); 125ef0bb6c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEComputeTabulation(PetscFE, PetscInt, const PetscReal[], PetscInt, PetscTabulation); 126ef0bb6c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscTabulationDestroy(PetscTabulation *); 127ef0bb6c7SMatthew G. Knepley 128bceba477SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFERefine(PetscFE, PetscFE *); 129228cfb18SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetHeightSubspace(PetscFE, PetscInt, PetscFE *); 130a0845e3aSMatthew G. Knepley 131c9ba7969SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateCellGeometry(PetscFE, PetscQuadrature, PetscFEGeom *); 132c9ba7969SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEDestroyCellGeometry(PetscFE, PetscFEGeom *); 1332edcad52SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEPushforward(PetscFE, PetscFEGeom *, PetscInt, PetscScalar[]); 1342edcad52SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEPushforwardGradient(PetscFE, PetscFEGeom *, PetscInt, PetscScalar[]); 135f9244615SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEPushforwardHessian(PetscFE, PetscFEGeom *, PetscInt, PetscScalar[]); 1364bee2e38SMatthew G. Knepley 1374bee2e38SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrate(PetscDS, PetscInt, PetscInt, PetscFEGeom *, const PetscScalar[], PetscDS, const PetscScalar[], PetscScalar[]); 1389371c9d4SSatish 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[]); 13906ad1575SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateResidual(PetscDS, PetscFormKey, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscScalar[]); 14006ad1575SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateBdResidual(PetscDS, PetscWeakForm, PetscFormKey, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscScalar[]); 14107218a29SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateHybridResidual(PetscDS, PetscDS, PetscFormKey, PetscInt, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscScalar[]); 14206ad1575SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateJacobian(PetscDS, PetscFEJacobianType, PetscFormKey, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscReal, PetscScalar[]); 143e3d591f2SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateBdJacobian(PetscDS, PetscWeakForm, PetscFEJacobianType, PetscFormKey, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscReal, PetscScalar[]); 14407218a29SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateHybridJacobian(PetscDS, PetscDS, PetscFEJacobianType, PetscFormKey, PetscInt, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscReal, PetscScalar[]); 145a0845e3aSMatthew G. Knepley 14689710940SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECompositeGetMapping(PetscFE, PetscInt *, const PetscReal *[], const PetscReal *[], const PetscReal *[]); 14789710940SMatthew G. Knepley 1482f5fb066SToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreateHeightTrace(PetscFE, PetscInt, PetscFE *); 1492f5fb066SToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreatePointTrace(PetscFE, PetscInt, PetscFE *); 1502f5fb066SToby Isaac 151855cd083SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEOpenCLSetRealType(PetscFE, PetscDataType); 152855cd083SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEOpenCLGetRealType(PetscFE, PetscDataType *); 153d2b2dc1eSMatthew G. Knepley 154d2b2dc1eSMatthew G. Knepley #ifdef PETSC_HAVE_LIBCEED 155d2b2dc1eSMatthew G. Knepley 156d2b2dc1eSMatthew G. Knepley #ifndef PLEXFE_QFUNCTION 157d2b2dc1eSMatthew G. Knepley #define PLEXFE_QFUNCTION(fname, f0_name, f1_name) \ 158d2b2dc1eSMatthew G. Knepley CEED_QFUNCTION(PlexQFunction##fname)(void *ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) \ 159d2b2dc1eSMatthew G. Knepley { \ 160d2b2dc1eSMatthew G. Knepley const CeedScalar *u = in[0], *du = in[1], *qdata = in[2]; \ 161d2b2dc1eSMatthew G. Knepley CeedScalar *v = out[0], *dv = out[1]; \ 162d2b2dc1eSMatthew G. Knepley const PetscInt Nc = 1; \ 163d2b2dc1eSMatthew G. Knepley const PetscInt cdim = 2; \ 164d2b2dc1eSMatthew G. Knepley \ 165d2b2dc1eSMatthew G. Knepley CeedPragmaSIMD for (CeedInt i = 0; i < Q; ++i) \ 166d2b2dc1eSMatthew G. Knepley { \ 167d2b2dc1eSMatthew G. Knepley const PetscInt uOff[2] = {0, Nc}; \ 168d2b2dc1eSMatthew G. Knepley const PetscInt uOff_x[2] = {0, Nc * cdim}; \ 169d2b2dc1eSMatthew G. Knepley const CeedScalar x[2] = {qdata[i + Q * 1], qdata[i + Q * 2]}; \ 170d2b2dc1eSMatthew G. Knepley const CeedScalar invJ[2][2] = { \ 171d2b2dc1eSMatthew G. Knepley {qdata[i + Q * 3], qdata[i + Q * 5]}, \ 172d2b2dc1eSMatthew G. Knepley {qdata[i + Q * 4], qdata[i + Q * 6]} \ 173d2b2dc1eSMatthew G. Knepley }; \ 174d2b2dc1eSMatthew 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]}; \ 175d2b2dc1eSMatthew G. Knepley PetscScalar f0[Nc]; \ 176d2b2dc1eSMatthew G. Knepley PetscScalar f1[Nc * cdim]; \ 177d2b2dc1eSMatthew G. Knepley \ 178d2b2dc1eSMatthew G. Knepley for (PetscInt k = 0; k < Nc; ++k) f0[k] = 0; \ 179d2b2dc1eSMatthew G. Knepley for (PetscInt k = 0; k < Nc * cdim; ++k) f1[k] = 0; \ 180d2b2dc1eSMatthew 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); \ 181d2b2dc1eSMatthew 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); \ 182d2b2dc1eSMatthew G. Knepley \ 183d2b2dc1eSMatthew G. Knepley dv[i + Q * 0] = qdata[i + Q * 0] * (invJ[0][0] * f1[0] + invJ[0][1] * f1[1]); \ 184d2b2dc1eSMatthew G. Knepley dv[i + Q * 1] = qdata[i + Q * 0] * (invJ[1][0] * f1[0] + invJ[1][1] * f1[1]); \ 185d2b2dc1eSMatthew G. Knepley v[i] = qdata[i + Q * 0] * f0[0]; \ 186d2b2dc1eSMatthew G. Knepley } \ 187d2b2dc1eSMatthew G. Knepley return CEED_ERROR_SUCCESS; \ 188d2b2dc1eSMatthew G. Knepley } 189d2b2dc1eSMatthew G. Knepley #endif 190d2b2dc1eSMatthew G. Knepley 191d2b2dc1eSMatthew G. Knepley #else 192d2b2dc1eSMatthew G. Knepley 193d2b2dc1eSMatthew G. Knepley #ifndef PLEXFE_QFUNCTION 194d2b2dc1eSMatthew G. Knepley #define PLEXFE_QFUNCTION(fname, f0_name, f1_name) 195d2b2dc1eSMatthew G. Knepley #endif 196d2b2dc1eSMatthew G. Knepley 197d2b2dc1eSMatthew G. Knepley #endif 198