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 14*ac9d17c7SMatthew G. Knepley /*E 15*ac9d17c7SMatthew G. Knepley PetscFEGeomMode - Describes the type of geometry being encoded. 16*ac9d17c7SMatthew G. Knepley 17*ac9d17c7SMatthew G. Knepley Values: 18*ac9d17c7SMatthew G. Knepley + `PETSC_FEGEOM_BASIC` - These are normal dim-cells, with dim == dE, and only bulk data is stored. 19*ac9d17c7SMatthew G. Knepley . `PETSC_FEGEOM_EMBEDDED` - These are dim-cells embedded in a higher dimension, as an embedded manifold, where dim < dE and only bulk data is stored. 20*ac9d17c7SMatthew G. Knepley . `PETSC_FEGEOM_BOUNDARY` - These are dim-cells on the boundary of a dE-mesh, so that dim < dE, and both bulk and s = 1 face data are stored. 21*ac9d17c7SMatthew G. Knepley - `PETSC_FEGEOM_COHESIVE` - These are dim-cells in the interior of a dE-mesh, so that dim < dE, and both bulk and s = 2 face data are stored. 22*ac9d17c7SMatthew G. Knepley 23*ac9d17c7SMatthew G. Knepley Level: beginner 24*ac9d17c7SMatthew G. Knepley 25*ac9d17c7SMatthew G. Knepley Note: 26*ac9d17c7SMatthew G. Knepley .vb 27*ac9d17c7SMatthew G. Knepley dim - The topological dimension and reference coordinate dimension 28*ac9d17c7SMatthew G. Knepley dE - The real coordinate dimension 29*ac9d17c7SMatthew G. Knepley s - The number of supporting cells for a face 30*ac9d17c7SMatthew G. Knepley .ve 31*ac9d17c7SMatthew G. Knepley 32*ac9d17c7SMatthew G. Knepley .seealso: [](ch_dmbase), `PetscFEGeom`, `DM`, `DMPLEX`, `PetscFEGeomCreate()` 33*ac9d17c7SMatthew G. Knepley E*/ 34*ac9d17c7SMatthew G. Knepley typedef enum { 35*ac9d17c7SMatthew G. Knepley PETSC_FEGEOM_BASIC, 36*ac9d17c7SMatthew G. Knepley PETSC_FEGEOM_EMBEDDED, 37*ac9d17c7SMatthew G. Knepley PETSC_FEGEOM_BOUNDARY, 38*ac9d17c7SMatthew G. Knepley PETSC_FEGEOM_COHESIVE 39*ac9d17c7SMatthew G. Knepley } PetscFEGeomMode; 40*ac9d17c7SMatthew G. Knepley 41dce8aebaSBarry Smith /*MC 42dce8aebaSBarry Smith PetscFEGeom - Structure for geometric information for `PetscFE` 43dce8aebaSBarry Smith 44dce8aebaSBarry Smith Level: intermediate 45dce8aebaSBarry Smith 465d83a8b1SBarry Smith Note: 475d83a8b1SBarry Smith This is a struct, not a `PetscObject` 485d83a8b1SBarry Smith 4916a05f60SBarry Smith .seealso: `PetscFE`, `PetscFEGeomCreate()`, `PetscFEGeomDestroy()`, `PetscFEGeomGetChunk()`, `PetscFEGeomRestoreChunk()`, `PetscFEGeomGetPoint()`, `PetscFEGeomGetCellPoint()`, 50b24fb147SBarry Smith `PetscFEGeomComplete()`, `PetscSpace`, `PetscDualSpace` 51dce8aebaSBarry Smith M*/ 524129dba9SToby Isaac typedef struct _n_PetscFEGeom { 53*ac9d17c7SMatthew G. Knepley // We can represent several different types of geometry, which we call modes: 54*ac9d17c7SMatthew G. Knepley // basic: dim == dE, only bulk data 55*ac9d17c7SMatthew G. Knepley // These are normal dim-cells 56*ac9d17c7SMatthew G. Knepley // embedded: dim < dE, only bulk data 57*ac9d17c7SMatthew G. Knepley // These are dim-cells embedded in a higher dimension, as an embedded manifold 58*ac9d17c7SMatthew G. Knepley // boundary: dim < dE, bulk and face data 59*ac9d17c7SMatthew G. Knepley // These are dim-cells on the boundary of a dE-mesh 60*ac9d17c7SMatthew G. Knepley // cohesive: dim < dE, bulk and face data 61*ac9d17c7SMatthew G. Knepley // These are dim-cells in the interior of a dE-mesh 62*ac9d17c7SMatthew G. Knepley // affine: 63*ac9d17c7SMatthew G. Knepley // For all modes, the transforms between real and reference are affine 64*ac9d17c7SMatthew G. Knepley PetscFEGeomMode mode; // The type of geometric data stored 65*ac9d17c7SMatthew G. Knepley PetscBool isAffine; // Flag for affine transforms 66*ac9d17c7SMatthew G. Knepley // Sizes 67*ac9d17c7SMatthew G. Knepley PetscInt dim; // dim: topological dimension and reference coordinate dimension 68*ac9d17c7SMatthew G. Knepley PetscInt dimEmbed; // dE: real coordinate dimension 69*ac9d17c7SMatthew G. Knepley PetscInt numCells; // Nc: Number of mesh points represented in the arrays (points are assumed to be the same DMPolytopeType) 70*ac9d17c7SMatthew G. Knepley PetscInt numPoints; // Np: Number of evaluation points represented in the arrays 71*ac9d17c7SMatthew G. Knepley // Bulk data 72*ac9d17c7SMatthew G. Knepley const PetscReal *xi; // xi[dim] The first point in each cell in reference coordinates 73*ac9d17c7SMatthew G. Knepley PetscReal *v; // v[Nc*Np*dE]: The first point in each cell in real coordinates 74*ac9d17c7SMatthew 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) 75*ac9d17c7SMatthew 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) 76*ac9d17c7SMatthew G. Knepley PetscReal *detJ; // detJ[Nc*Np]: The determinant of J, and if J is non-square it is the volume change 77*ac9d17c7SMatthew G. Knepley // Face data 78*ac9d17c7SMatthew G. Knepley PetscReal *n; // n[Nc*Np*dE]: For faces, the normal to the face in real coordinates, outward for the first supporting cell 79*ac9d17c7SMatthew 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 80*ac9d17c7SMatthew G. Knepley PetscReal *suppJ[2]; // sJ[s][Nc*Np*dE*dE]: For faces, the Jacobian for each supporting cell 81*ac9d17c7SMatthew G. Knepley PetscReal *suppInvJ[2]; // sInvJ[s][Nc*Np*dE*dE]: For faces, the inverse Jacobian for each supporting cell 82*ac9d17c7SMatthew G. Knepley PetscReal *suppDetJ[2]; // sdetJ[s][Nc*Np]: For faces, the Jacobian determinant for each supporting cell 834129dba9SToby Isaac } PetscFEGeom; 844129dba9SToby Isaac 85dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEInitializePackage(void); 86dbe77d9eSMatthew G. Knepley 87*ac9d17c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGeomCreate(PetscQuadrature, PetscInt, PetscInt, PetscFEGeomMode, PetscFEGeom **); 884129dba9SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEGeomGetChunk(PetscFEGeom *, PetscInt, PetscInt, PetscFEGeom **); 894129dba9SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEGeomRestoreChunk(PetscFEGeom *, PetscInt, PetscInt, PetscFEGeom **); 906587ee25SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGeomGetPoint(PetscFEGeom *, PetscInt, PetscInt, const PetscReal[], PetscFEGeom *); 916587ee25SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGeomGetCellPoint(PetscFEGeom *, PetscInt, PetscInt, PetscFEGeom *); 924129dba9SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEGeomComplete(PetscFEGeom *); 934129dba9SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEGeomDestroy(PetscFEGeom **); 944129dba9SToby Isaac 95c330f8ffSToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceApply(PetscDualSpace, PetscInt, PetscReal, PetscFEGeom *, PetscInt, PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void *, PetscScalar *); 96c330f8ffSToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceApplyDefault(PetscDualSpace, PetscInt, PetscReal, PetscFEGeom *, PetscInt, PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void *, PetscScalar *); 97ebac44aeSMatthew G. Knepley 984bee2e38SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceTransform(PetscDualSpace, PetscDualSpaceTransformType, PetscBool, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]); 99f9244615SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceTransformGradient(PetscDualSpace, PetscDualSpaceTransformType, PetscBool, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]); 100f9244615SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceTransformHessian(PetscDualSpace, PetscDualSpaceTransformType, PetscBool, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]); 1012edcad52SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpacePullback(PetscDualSpace, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]); 1022edcad52SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpacePushforward(PetscDualSpace, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]); 1032edcad52SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpacePushforwardGradient(PetscDualSpace, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]); 104f9244615SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpacePushforwardHessian(PetscDualSpace, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]); 1054bee2e38SMatthew G. Knepley 106dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscClassId PETSCFE_CLASSID; 107dbe77d9eSMatthew G. Knepley 1080483ade4SMatthew G. Knepley /*J 1090483ade4SMatthew G. Knepley PetscFEType - String with the name of a PETSc finite element space 1100483ade4SMatthew G. Knepley 1110483ade4SMatthew G. Knepley Level: beginner 1120483ade4SMatthew G. Knepley 11316a05f60SBarry Smith Note: 11416a05f60SBarry Smith Currently, the classes are concerned with the implementation of element integration 1150483ade4SMatthew G. Knepley 116db781477SPatrick Sanan .seealso: `PetscFESetType()`, `PetscFE` 1170483ade4SMatthew G. Knepley J*/ 1180483ade4SMatthew G. Knepley typedef const char *PetscFEType; 1190483ade4SMatthew G. Knepley #define PETSCFEBASIC "basic" 1200483ade4SMatthew G. Knepley #define PETSCFEOPENCL "opencl" 121aaf1837cSMatthew G. Knepley #define PETSCFECOMPOSITE "composite" 1222dce792eSToby Isaac #define PETSCFEVECTOR "vector" 1230483ade4SMatthew G. Knepley 1240483ade4SMatthew G. Knepley PETSC_EXTERN PetscFunctionList PetscFEList; 125568c31deSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreate(MPI_Comm, PetscFE *); 126568c31deSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEDestroy(PetscFE *); 1270483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetType(PetscFE, PetscFEType); 1280483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetType(PetscFE, PetscFEType *); 1290483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetUp(PetscFE); 1300483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetFromOptions(PetscFE); 131fe2efc57SMark PETSC_EXTERN PetscErrorCode PetscFEViewFromOptions(PetscFE, PetscObject, const char[]); 1323f6b16c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetName(PetscFE, const char[]); 1332dce792eSToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreateVector(PetscFE, PetscInt, PetscBool, PetscBool, PetscFE *); 1348aec7d55SBarry Smith 135fd41f350SJed Brown PETSC_EXTERN PetscErrorCode PetscFEView(PetscFE, PetscViewer); 1360483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFERegister(const char[], PetscErrorCode (*)(PetscFE)); 1370483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFERegisterDestroy(void); 138e8d98e54SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateDefault(MPI_Comm, PetscInt, PetscInt, PetscBool, const char[], PetscInt, PetscFE *); 1392df84da0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateByCell(MPI_Comm, PetscInt, PetscInt, DMPolytopeType, const char[], PetscInt, PetscFE *); 140e703855dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateLagrange(MPI_Comm, PetscInt, PetscInt, PetscBool, PetscInt, PetscInt, PetscFE *); 1412df84da0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateLagrangeByCell(MPI_Comm, PetscInt, PetscInt, DMPolytopeType, PetscInt, PetscInt, PetscFE *); 1427c48043bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateFromSpaces(PetscSpace, PetscDualSpace, PetscQuadrature, PetscQuadrature, PetscFE *); 143bb4b53efSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFELimitDegree(PetscFE, PetscInt, PetscInt, PetscFE *); 1449a1a3eb8SMatthew G. Knepley 1459a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetDimension(PetscFE, PetscInt *); 1460ddb9b0bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetSpatialDimension(PetscFE, PetscInt *); 1479a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetNumComponents(PetscFE, PetscInt); 1489a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetNumComponents(PetscFE, PetscInt *); 1499a769b6bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetTileSizes(PetscFE, PetscInt *, PetscInt *, PetscInt *, PetscInt *); 1509a769b6bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetTileSizes(PetscFE, PetscInt, PetscInt, PetscInt, PetscInt); 1519a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetBasisSpace(PetscFE, PetscSpace); 1529a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetBasisSpace(PetscFE, PetscSpace *); 1539a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetDualSpace(PetscFE, PetscDualSpace); 1549a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetDualSpace(PetscFE, PetscDualSpace *); 155bfa639d9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetQuadrature(PetscFE, PetscQuadrature); 156bfa639d9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetQuadrature(PetscFE, PetscQuadrature *); 1571bafc85dSSander Arens PETSC_EXTERN PetscErrorCode PetscFESetFaceQuadrature(PetscFE, PetscQuadrature); 1581bafc85dSSander Arens PETSC_EXTERN PetscErrorCode PetscFEGetFaceQuadrature(PetscFE, PetscQuadrature *); 1595dc5c000SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECopyQuadrature(PetscFE, PetscFE); 1600ddb9b0bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetNumDof(PetscFE, const PetscInt **); 161ef0bb6c7SMatthew G. Knepley 162ef0bb6c7SMatthew G. Knepley /* TODO: Need a function to reuse the memory when retabulating the same FE at different points */ 163f9244615SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetCellTabulation(PetscFE, PetscInt, PetscTabulation *); 164f9244615SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetFaceTabulation(PetscFE, PetscInt, PetscTabulation *); 165ef0bb6c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetFaceCentroidTabulation(PetscFE, PetscTabulation *); 166ef0bb6c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateTabulation(PetscFE, PetscInt, PetscInt, const PetscReal[], PetscInt, PetscTabulation *); 167ef0bb6c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEComputeTabulation(PetscFE, PetscInt, const PetscReal[], PetscInt, PetscTabulation); 168ef0bb6c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscTabulationDestroy(PetscTabulation *); 169ef0bb6c7SMatthew G. Knepley 170bceba477SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFERefine(PetscFE, PetscFE *); 171228cfb18SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetHeightSubspace(PetscFE, PetscInt, PetscFE *); 172a0845e3aSMatthew G. Knepley 173c9ba7969SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateCellGeometry(PetscFE, PetscQuadrature, PetscFEGeom *); 174c9ba7969SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEDestroyCellGeometry(PetscFE, PetscFEGeom *); 1752edcad52SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEPushforward(PetscFE, PetscFEGeom *, PetscInt, PetscScalar[]); 1762edcad52SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEPushforwardGradient(PetscFE, PetscFEGeom *, PetscInt, PetscScalar[]); 177f9244615SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEPushforwardHessian(PetscFE, PetscFEGeom *, PetscInt, PetscScalar[]); 1784bee2e38SMatthew G. Knepley 1794bee2e38SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrate(PetscDS, PetscInt, PetscInt, PetscFEGeom *, const PetscScalar[], PetscDS, const PetscScalar[], PetscScalar[]); 1809371c9d4SSatish 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[]); 18106ad1575SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateResidual(PetscDS, PetscFormKey, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscScalar[]); 18206ad1575SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateBdResidual(PetscDS, PetscWeakForm, PetscFormKey, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscScalar[]); 18307218a29SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateHybridResidual(PetscDS, PetscDS, PetscFormKey, PetscInt, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscScalar[]); 18406ad1575SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateJacobian(PetscDS, PetscFEJacobianType, PetscFormKey, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscReal, PetscScalar[]); 185e3d591f2SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateBdJacobian(PetscDS, PetscWeakForm, PetscFEJacobianType, PetscFormKey, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscReal, PetscScalar[]); 18607218a29SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateHybridJacobian(PetscDS, PetscDS, PetscFEJacobianType, PetscFormKey, PetscInt, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscReal, PetscScalar[]); 187a0845e3aSMatthew G. Knepley 18889710940SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECompositeGetMapping(PetscFE, PetscInt *, const PetscReal *[], const PetscReal *[], const PetscReal *[]); 18989710940SMatthew G. Knepley 1902f5fb066SToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreateHeightTrace(PetscFE, PetscInt, PetscFE *); 1912f5fb066SToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreatePointTrace(PetscFE, PetscInt, PetscFE *); 1922f5fb066SToby Isaac 193855cd083SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEOpenCLSetRealType(PetscFE, PetscDataType); 194855cd083SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEOpenCLGetRealType(PetscFE, PetscDataType *); 195d2b2dc1eSMatthew G. Knepley 196d2b2dc1eSMatthew G. Knepley #ifdef PETSC_HAVE_LIBCEED 197d2b2dc1eSMatthew G. Knepley 198d2b2dc1eSMatthew G. Knepley #ifndef PLEXFE_QFUNCTION 199d2b2dc1eSMatthew G. Knepley #define PLEXFE_QFUNCTION(fname, f0_name, f1_name) \ 200d2b2dc1eSMatthew G. Knepley CEED_QFUNCTION(PlexQFunction##fname)(void *ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) \ 201d2b2dc1eSMatthew G. Knepley { \ 202d2b2dc1eSMatthew G. Knepley const CeedScalar *u = in[0], *du = in[1], *qdata = in[2]; \ 203d2b2dc1eSMatthew G. Knepley CeedScalar *v = out[0], *dv = out[1]; \ 204d2b2dc1eSMatthew G. Knepley const PetscInt Nc = 1; \ 205d2b2dc1eSMatthew G. Knepley const PetscInt cdim = 2; \ 206d2b2dc1eSMatthew G. Knepley \ 207d2b2dc1eSMatthew G. Knepley CeedPragmaSIMD for (CeedInt i = 0; i < Q; ++i) \ 208d2b2dc1eSMatthew G. Knepley { \ 209d2b2dc1eSMatthew G. Knepley const PetscInt uOff[2] = {0, Nc}; \ 210d2b2dc1eSMatthew G. Knepley const PetscInt uOff_x[2] = {0, Nc * cdim}; \ 211d2b2dc1eSMatthew G. Knepley const CeedScalar x[2] = {qdata[i + Q * 1], qdata[i + Q * 2]}; \ 212d2b2dc1eSMatthew G. Knepley const CeedScalar invJ[2][2] = { \ 213d2b2dc1eSMatthew G. Knepley {qdata[i + Q * 3], qdata[i + Q * 5]}, \ 214d2b2dc1eSMatthew G. Knepley {qdata[i + Q * 4], qdata[i + Q * 6]} \ 215d2b2dc1eSMatthew G. Knepley }; \ 216d2b2dc1eSMatthew 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]}; \ 217d2b2dc1eSMatthew G. Knepley PetscScalar f0[Nc]; \ 218d2b2dc1eSMatthew G. Knepley PetscScalar f1[Nc * cdim]; \ 219d2b2dc1eSMatthew G. Knepley \ 220d2b2dc1eSMatthew G. Knepley for (PetscInt k = 0; k < Nc; ++k) f0[k] = 0; \ 221d2b2dc1eSMatthew G. Knepley for (PetscInt k = 0; k < Nc * cdim; ++k) f1[k] = 0; \ 222d2b2dc1eSMatthew 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); \ 223d2b2dc1eSMatthew 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); \ 224d2b2dc1eSMatthew G. Knepley \ 225d2b2dc1eSMatthew G. Knepley dv[i + Q * 0] = qdata[i + Q * 0] * (invJ[0][0] * f1[0] + invJ[0][1] * f1[1]); \ 226d2b2dc1eSMatthew G. Knepley dv[i + Q * 1] = qdata[i + Q * 0] * (invJ[1][0] * f1[0] + invJ[1][1] * f1[1]); \ 227d2b2dc1eSMatthew G. Knepley v[i] = qdata[i + Q * 0] * f0[0]; \ 228d2b2dc1eSMatthew G. Knepley } \ 229d2b2dc1eSMatthew G. Knepley return CEED_ERROR_SUCCESS; \ 230d2b2dc1eSMatthew G. Knepley } 231d2b2dc1eSMatthew G. Knepley #endif 232d2b2dc1eSMatthew G. Knepley 233d2b2dc1eSMatthew G. Knepley #else 234d2b2dc1eSMatthew G. Knepley 235d2b2dc1eSMatthew G. Knepley #ifndef PLEXFE_QFUNCTION 236d2b2dc1eSMatthew G. Knepley #define PLEXFE_QFUNCTION(fname, f0_name, f1_name) 237d2b2dc1eSMatthew G. Knepley #endif 238d2b2dc1eSMatthew G. Knepley 239d2b2dc1eSMatthew G. Knepley #endif 240