xref: /petsc/include/petscfe.h (revision ac09b9214d23ea9ad238aa607de9fa447fd4e91b)
1dbe77d9eSMatthew G. Knepley /*
2dbe77d9eSMatthew G. Knepley       Objects which encapsulate finite element spaces and operations
3dbe77d9eSMatthew G. Knepley */
426bd1501SBarry Smith #if !defined(PETSCFE_H)
526bd1501SBarry Smith #define PETSCFE_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>
10dbe77d9eSMatthew G. Knepley 
11*ac09b921SBarry Smith /* SUBMANSEC = FE */
12*ac09b921SBarry Smith 
134129dba9SToby Isaac typedef struct _n_PetscFEGeom {
144129dba9SToby Isaac   const PetscReal *xi;
1541418987SMatthew G. Knepley   PetscReal *v;           /* v[Nc*Np*dE]:           The first point in each each in real coordinates */
1641418987SMatthew 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) */
1741418987SMatthew 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) */
1841418987SMatthew G. Knepley   PetscReal *detJ;        /* detJ[Nc*Np]:           The determinant of J, and if it is non-square its the volume change */
19f15274beSMatthew Knepley   PetscReal *n;           /* n[Nc*Np*dE]:           For faces, the normal to the face in real coordinates, outward for the first supporting cell */
2041418987SMatthew G. Knepley   PetscInt  (*face)[2];   /* face[Nc][s]:           For faces, the local face number (cone index) for this face in each supporting cell s */
2141418987SMatthew G. Knepley   PetscReal *suppJ[2];    /* sJ[s][Nc*Np*dE*dE]:    For faces, the Jacobian for each supporting cell s */
2241418987SMatthew G. Knepley   PetscReal *suppInvJ[2]; /* sInvJ[s][Nc*Np*dE*dE]: For faces, the inverse Jacobian for each supporting cell s */
2341418987SMatthew G. Knepley   PetscReal *suppDetJ[2]; /* sInvJ[s][Nc*Np*dE*dE]: For faces, the Jacobian determinant for each supporting cell s */
2462bd480fSMatthew G. Knepley   PetscInt  dim;          /* dim: Topological dimension */
2562bd480fSMatthew G. Knepley   PetscInt  dimEmbed;     /* dE:  coordinate dimension */
2662bd480fSMatthew G. Knepley   PetscInt  numCells;     /* Nc:  Number of mesh points represented in the arrays */
2762bd480fSMatthew G. Knepley   PetscInt  numPoints;    /* Np:  Number of evaluation points represented in the arrays */
2841418987SMatthew G. Knepley   PetscBool isAffine;     /* Flag for affine transforms */
295fedec97SMatthew G. Knepley   PetscBool isCohesive;   /* Flag for a cohesive cell */
304129dba9SToby Isaac } PetscFEGeom;
314129dba9SToby Isaac 
32dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEInitializePackage(void);
33dbe77d9eSMatthew G. Knepley 
34dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscClassId PETSCSPACE_CLASSID;
35dbe77d9eSMatthew G. Knepley 
36dbe77d9eSMatthew G. Knepley /*J
37dbe77d9eSMatthew G. Knepley   PetscSpaceType - String with the name of a PETSc linear space
38dbe77d9eSMatthew G. Knepley 
39dbe77d9eSMatthew G. Knepley   Level: beginner
40dbe77d9eSMatthew G. Knepley 
41db781477SPatrick Sanan .seealso: `PetscSpaceSetType()`, `PetscSpace`
42dbe77d9eSMatthew G. Knepley J*/
43dbe77d9eSMatthew G. Knepley typedef const char* PetscSpaceType;
44dbe77d9eSMatthew G. Knepley #define PETSCSPACEPOLYNOMIAL "poly"
45130d5748SToby Isaac #define PETSCSPACEPTRIMMED   "ptrimmed"
4636e5648fSToby Isaac #define PETSCSPACETENSOR     "tensor"
47d092c84bSBrandon Whitchurch #define PETSCSPACESUM        "sum"
489c3cf19fSMatthew G. Knepley #define PETSCSPACEPOINT      "point"
492f5fb066SToby Isaac #define PETSCSPACESUBSPACE   "subspace"
5069cc43acSMatthew G. Knepley #define PETSCSPACEWXY        "wxy"
51dbe77d9eSMatthew G. Knepley 
52dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscFunctionList PetscSpaceList;
53dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceCreate(MPI_Comm, PetscSpace *);
54568c31deSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceDestroy(PetscSpace *);
55dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceSetType(PetscSpace, PetscSpaceType);
56dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceGetType(PetscSpace, PetscSpaceType *);
579a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceSetUp(PetscSpace);
58dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceSetFromOptions(PetscSpace);
59fe2efc57SMark PETSC_EXTERN PetscErrorCode PetscSpaceViewFromOptions(PetscSpace,PetscObject,const char[]);
608aec7d55SBarry Smith 
61fd41f350SJed Brown PETSC_EXTERN PetscErrorCode PetscSpaceView(PetscSpace,PetscViewer);
62dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceRegister(const char [], PetscErrorCode (*)(PetscSpace));
63dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceRegisterDestroy(void);
64dbe77d9eSMatthew G. Knepley 
659a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceGetDimension(PetscSpace, PetscInt *);
669c3cf19fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceSetNumComponents(PetscSpace, PetscInt);
679c3cf19fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceGetNumComponents(PetscSpace, PetscInt *);
68157782e2SToby Isaac PETSC_EXTERN PetscErrorCode PetscSpaceSetNumVariables(PetscSpace, PetscInt);
69157782e2SToby Isaac PETSC_EXTERN PetscErrorCode PetscSpaceGetNumVariables(PetscSpace, PetscInt *);
70d39dd5f5SToby Isaac PETSC_EXTERN PetscErrorCode PetscSpaceSetDegree(PetscSpace, PetscInt, PetscInt);
71cdcc6362SToby Isaac PETSC_EXTERN PetscErrorCode PetscSpaceGetDegree(PetscSpace, PetscInt *, PetscInt *);
722bdb15eaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceEvaluate(PetscSpace, PetscInt, const PetscReal[], PetscReal[], PetscReal[], PetscReal[]);
730c16f5adSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceGetHeightSubspace(PetscSpace, PetscInt, PetscSpace *);
749a1a3eb8SMatthew G. Knepley 
759fbee547SJacob Faibussowitsch static inline PETSC_DEPRECATED_FUNCTION("Property not used (since v3.17)") PetscErrorCode PetscSpacePolynomialSetSymmetric(PetscSpace sp, PetscBool s)
76f1436e55SToby Isaac {
772c71b3e2SJacob Faibussowitsch   PetscCheck(!s,PetscObjectComm((PetscObject)sp),PETSC_ERR_SUP,"PETSCSPACEPOLYNOMIAL does not support symmetric polynomials");
78f1436e55SToby Isaac   return 0;
79f1436e55SToby Isaac }
809fbee547SJacob Faibussowitsch static inline PETSC_DEPRECATED_FUNCTION("Property not used (since v3.17)") PetscErrorCode PetscSpacePolynomialGetSymmetric(PETSC_UNUSED PetscSpace sp, PetscBool *s) {*s = PETSC_FALSE; return 0;}
81c43d9e4aSToby Isaac PETSC_EXTERN PetscErrorCode PetscSpacePolynomialSetTensor(PetscSpace, PetscBool);
82c43d9e4aSToby Isaac PETSC_EXTERN PetscErrorCode PetscSpacePolynomialGetTensor(PetscSpace, PetscBool *);
832bdb15eaSMatthew G. Knepley 
84130d5748SToby Isaac PETSC_EXTERN PetscErrorCode PetscSpacePTrimmedSetFormDegree(PetscSpace, PetscInt);
85130d5748SToby Isaac PETSC_EXTERN PetscErrorCode PetscSpacePTrimmedGetFormDegree(PetscSpace, PetscInt *);
86130d5748SToby Isaac 
8736e5648fSToby Isaac PETSC_EXTERN PetscErrorCode PetscSpaceTensorSetNumSubspaces(PetscSpace, PetscInt);
8836e5648fSToby Isaac PETSC_EXTERN PetscErrorCode PetscSpaceTensorGetNumSubspaces(PetscSpace, PetscInt *);
8936e5648fSToby Isaac PETSC_EXTERN PetscErrorCode PetscSpaceTensorSetSubspace(PetscSpace, PetscInt, PetscSpace);
9036e5648fSToby Isaac PETSC_EXTERN PetscErrorCode PetscSpaceTensorGetSubspace(PetscSpace, PetscInt, PetscSpace *);
9136e5648fSToby Isaac 
92d092c84bSBrandon Whitchurch PETSC_EXTERN PetscErrorCode PetscSpaceSumSetNumSubspaces(PetscSpace, PetscInt);
93d092c84bSBrandon Whitchurch PETSC_EXTERN PetscErrorCode PetscSpaceSumGetNumSubspaces(PetscSpace, PetscInt *);
94d092c84bSBrandon Whitchurch PETSC_EXTERN PetscErrorCode PetscSpaceSumSetSubspace(PetscSpace, PetscInt, PetscSpace);
95d092c84bSBrandon Whitchurch PETSC_EXTERN PetscErrorCode PetscSpaceSumGetSubspace(PetscSpace, PetscInt, PetscSpace *);
96d092c84bSBrandon Whitchurch PETSC_EXTERN PetscErrorCode PetscSpaceSumSetConcatenate(PetscSpace, PetscBool);
97d092c84bSBrandon Whitchurch PETSC_EXTERN PetscErrorCode PetscSpaceSumGetConcatenate(PetscSpace, PetscBool *);
98d092c84bSBrandon Whitchurch PETSC_EXTERN PetscErrorCode PetscSpaceCreateSum(PetscInt numSubspaces,const PetscSpace subspaces[],PetscBool concatenate,PetscSpace *sumSpace);
99d092c84bSBrandon Whitchurch 
1008049c7f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpacePointGetPoints(PetscSpace, PetscQuadrature *);
1018049c7f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpacePointSetPoints(PetscSpace, PetscQuadrature);
1029a1a3eb8SMatthew G. Knepley 
1032f5fb066SToby Isaac PETSC_EXTERN PetscErrorCode PetscSpaceCreateSubspace(PetscSpace, PetscDualSpace, PetscReal *, PetscReal *, PetscReal *, PetscReal *, PetscCopyMode, PetscSpace *);
1042f5fb066SToby Isaac 
105dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscClassId PETSCDUALSPACE_CLASSID;
106dbe77d9eSMatthew G. Knepley 
107dbe77d9eSMatthew G. Knepley /*J
108dbe77d9eSMatthew G. Knepley   PetscDualSpaceType - String with the name of a PETSc dual space
109dbe77d9eSMatthew G. Knepley 
110dbe77d9eSMatthew G. Knepley   Level: beginner
111dbe77d9eSMatthew G. Knepley 
112db781477SPatrick Sanan .seealso: `PetscDualSpaceSetType()`, `PetscDualSpace`
113dbe77d9eSMatthew G. Knepley J*/
114dbe77d9eSMatthew G. Knepley typedef const char *PetscDualSpaceType;
115dbe77d9eSMatthew G. Knepley #define PETSCDUALSPACELAGRANGE "lagrange"
116c2765ee2SMatthew G. Knepley #define PETSCDUALSPACESIMPLE   "simple"
1171ac17e89SToby Isaac #define PETSCDUALSPACEREFINED  "refined"
11870715647SToby Isaac #define PETSCDUALSPACEBDM      "bdm"
11970715647SToby Isaac 
12070715647SToby Isaac /*MC
12170715647SToby Isaac   PETSCDUALSPACEBDM = "bdm" - A PetscDualSpace object that encapsulates a dual space for Brezzi-Douglas-Marini elements
12270715647SToby Isaac 
12370715647SToby Isaac   Level: intermediate
12470715647SToby Isaac 
12570715647SToby Isaac   Note: This type is a constructor alias of PETSCDUALSPACELAGRANGE.  During
12670715647SToby Isaac   PetscDualSpaceSetUp(), the correct value of PetscDualSpaceSetFormDegree() is
12770715647SToby Isaac   set for H-div conforming spaces. The type of the dual space is then changed to
12870715647SToby Isaac   to PETSCDUALSPACELAGRANGE.
12970715647SToby Isaac 
130db781477SPatrick Sanan .seealso: `PetscDualSpaceType`, `PetscDualSpaceCreate()`, `PetscDualSpaceSetType()`, `PETSCDUALSPACELAGRANGE`, `PetscDualSpaceSetFormDegree()`
13170715647SToby Isaac M*/
132dbe77d9eSMatthew G. Knepley 
133dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscFunctionList PetscDualSpaceList;
134dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceCreate(MPI_Comm, PetscDualSpace *);
135568c31deSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceDestroy(PetscDualSpace *);
136907df5ccSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceDuplicate(PetscDualSpace, PetscDualSpace *);
137dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceSetType(PetscDualSpace, PetscDualSpaceType);
138dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceGetType(PetscDualSpace, PetscDualSpaceType *);
139b4457527SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceGetUniform(PetscDualSpace, PetscBool *);
1408d2f55e7SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceGetNumDof(PetscDualSpace, const PetscInt **);
141b4457527SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceGetSection(PetscDualSpace, PetscSection *);
1429a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceSetUp(PetscDualSpace);
143dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceSetFromOptions(PetscDualSpace);
144fe2efc57SMark PETSC_EXTERN PetscErrorCode PetscDualSpaceViewFromOptions(PetscDualSpace,PetscObject,const char[]);
1458aec7d55SBarry Smith 
146fd41f350SJed Brown PETSC_EXTERN PetscErrorCode PetscDualSpaceView(PetscDualSpace,PetscViewer);
147dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceRegister(const char [], PetscErrorCode (*)(PetscDualSpace));
148dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceRegisterDestroy(void);
149dbe77d9eSMatthew G. Knepley 
150dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceGetDimension(PetscDualSpace, PetscInt *);
151b4457527SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceGetInteriorDimension(PetscDualSpace, PetscInt *);
1529c3cf19fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceSetNumComponents(PetscDualSpace, PetscInt);
1539c3cf19fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceGetNumComponents(PetscDualSpace, PetscInt *);
1549a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceSetOrder(PetscDualSpace, PetscInt);
1559a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceGetOrder(PetscDualSpace, PetscInt *);
1569a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceSetDM(PetscDualSpace, DM);
1579a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceGetDM(PetscDualSpace, DM *);
158ebac44aeSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceGetFunctional(PetscDualSpace, PetscInt, PetscQuadrature *);
159fffc2ff8SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceGetSymmetries(PetscDualSpace, const PetscInt ****, const PetscScalar ****);
160dbe77d9eSMatthew G. Knepley 
1614129dba9SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEGeomCreate(PetscQuadrature,PetscInt,PetscInt,PetscBool,PetscFEGeom**);
16262bd480fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGeomGetQuadrature(PetscFEGeom*,PetscQuadrature*);
16362bd480fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGeomSetQuadrature(PetscFEGeom*,PetscQuadrature);
1644129dba9SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEGeomGetChunk(PetscFEGeom*,PetscInt,PetscInt,PetscFEGeom**);
1654129dba9SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEGeomRestoreChunk(PetscFEGeom*,PetscInt,PetscInt,PetscFEGeom**);
1666587ee25SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGeomGetPoint(PetscFEGeom*,PetscInt,PetscInt,const PetscReal[],PetscFEGeom*);
1676587ee25SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGeomGetCellPoint(PetscFEGeom*,PetscInt,PetscInt,PetscFEGeom*);
1684129dba9SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEGeomComplete(PetscFEGeom*);
1694129dba9SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEGeomDestroy(PetscFEGeom**);
1704129dba9SToby Isaac 
171c330f8ffSToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceApply(PetscDualSpace, PetscInt, PetscReal, PetscFEGeom *, PetscInt, PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void *, PetscScalar *);
172c330f8ffSToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceApplyDefault(PetscDualSpace, PetscInt, PetscReal, PetscFEGeom *, PetscInt, PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void *, PetscScalar *);
173ebac44aeSMatthew G. Knepley 
174b4457527SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceGetAllData(PetscDualSpace, PetscQuadrature *, Mat *);
175b4457527SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceCreateAllDataDefault(PetscDualSpace, PetscQuadrature *, Mat *);
176b4457527SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceGetInteriorData(PetscDualSpace, PetscQuadrature *, Mat *);
177b4457527SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceCreateInteriorDataDefault(PetscDualSpace, PetscQuadrature *, Mat *);
1784f9ab2b4SJed Brown PETSC_EXTERN PetscErrorCode PetscDualSpaceEqual(PetscDualSpace, PetscDualSpace, PetscBool *);
179976f670dSToby Isaac 
180edfea91aSToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceApplyAll(PetscDualSpace, const PetscScalar *, PetscScalar *);
181edfea91aSToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceApplyAllDefault(PetscDualSpace, const PetscScalar *, PetscScalar *);
182b4457527SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceApplyInterior(PetscDualSpace, const PetscScalar *, PetscScalar *);
183b4457527SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceApplyInteriorDefault(PetscDualSpace, const PetscScalar *, PetscScalar *);
1843f4bc389SToby Isaac 
185b4457527SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceGetFormDegree(PetscDualSpace, PetscInt *);
186b4457527SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceSetFormDegree(PetscDualSpace, PetscInt);
1874bee2e38SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceGetDeRahm(PetscDualSpace, PetscInt *);
1884bee2e38SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceTransform(PetscDualSpace, PetscDualSpaceTransformType, PetscBool, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]);
189f9244615SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceTransformGradient(PetscDualSpace, PetscDualSpaceTransformType, PetscBool, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]);
190f9244615SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceTransformHessian(PetscDualSpace, PetscDualSpaceTransformType, PetscBool, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]);
1912edcad52SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpacePullback(PetscDualSpace, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]);
1922edcad52SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpacePushforward(PetscDualSpace, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]);
1932edcad52SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpacePushforwardGradient(PetscDualSpace, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]);
194f9244615SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpacePushforwardHessian(PetscDualSpace, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]);
1954bee2e38SMatthew G. Knepley 
19619ded50bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeGetContinuity(PetscDualSpace, PetscBool *);
19719ded50bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeSetContinuity(PetscDualSpace, PetscBool);
19889cad2ffSToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeGetTensor(PetscDualSpace, PetscBool *);
19989cad2ffSToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeSetTensor(PetscDualSpace, PetscBool);
2003f27d899SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeGetTrimmed(PetscDualSpace, PetscBool *);
2013f27d899SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeSetTrimmed(PetscDualSpace, PetscBool);
2023f27d899SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeGetNodeType(PetscDualSpace, PetscDTNodeType *, PetscBool *, PetscReal *);
2033f27d899SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeSetNodeType(PetscDualSpace, PetscDTNodeType, PetscBool, PetscReal);
20466a6c23cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeGetUseMoments(PetscDualSpace, PetscBool *);
20566a6c23cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeSetUseMoments(PetscDualSpace, PetscBool);
20666a6c23cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeGetMomentOrder(PetscDualSpace, PetscInt *);
20766a6c23cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeSetMomentOrder(PetscDualSpace, PetscInt);
20819ded50bSMatthew G. Knepley 
2096b594bd0SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceGetHeightSubspace(PetscDualSpace,PetscInt,PetscDualSpace *);
2102f5fb066SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceGetPointSubspace(PetscDualSpace,PetscInt,PetscDualSpace *);
211c2765ee2SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceSimpleSetDimension(PetscDualSpace, PetscInt);
212c2765ee2SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceSimpleSetFunctional(PetscDualSpace, PetscInt, PetscQuadrature);
2136b594bd0SToby Isaac 
2141ac17e89SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceRefinedSetCellSpaces(PetscDualSpace, const PetscDualSpace[]);
2151ac17e89SToby Isaac 
216dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscClassId PETSCFE_CLASSID;
217dbe77d9eSMatthew G. Knepley 
2180483ade4SMatthew G. Knepley /*J
2190483ade4SMatthew G. Knepley   PetscFEType - String with the name of a PETSc finite element space
2200483ade4SMatthew G. Knepley 
2210483ade4SMatthew G. Knepley   Level: beginner
2220483ade4SMatthew G. Knepley 
2230483ade4SMatthew G. Knepley   Note: Currently, the classes are concerned with the implementation of element integration
2240483ade4SMatthew G. Knepley 
225db781477SPatrick Sanan .seealso: `PetscFESetType()`, `PetscFE`
2260483ade4SMatthew G. Knepley J*/
2270483ade4SMatthew G. Knepley typedef const char *PetscFEType;
2280483ade4SMatthew G. Knepley #define PETSCFEBASIC     "basic"
2290483ade4SMatthew G. Knepley #define PETSCFEOPENCL    "opencl"
230aaf1837cSMatthew G. Knepley #define PETSCFECOMPOSITE "composite"
2310483ade4SMatthew G. Knepley 
2320483ade4SMatthew G. Knepley PETSC_EXTERN PetscFunctionList PetscFEList;
233568c31deSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreate(MPI_Comm, PetscFE *);
234568c31deSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEDestroy(PetscFE *);
2350483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetType(PetscFE, PetscFEType);
2360483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetType(PetscFE, PetscFEType *);
2370483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetUp(PetscFE);
2380483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetFromOptions(PetscFE);
239fe2efc57SMark PETSC_EXTERN PetscErrorCode PetscFEViewFromOptions(PetscFE,PetscObject,const char[]);
2403f6b16c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetName(PetscFE, const char []);
2418aec7d55SBarry Smith 
242fd41f350SJed Brown PETSC_EXTERN PetscErrorCode PetscFEView(PetscFE,PetscViewer);
2430483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFERegister(const char [], PetscErrorCode (*)(PetscFE));
2440483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFERegisterDestroy(void);
245e8d98e54SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateDefault(MPI_Comm, PetscInt, PetscInt, PetscBool, const char[], PetscInt, PetscFE *);
2462df84da0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateByCell(MPI_Comm, PetscInt, PetscInt, DMPolytopeType, const char[], PetscInt, PetscFE *);
247e703855dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateLagrange(MPI_Comm, PetscInt, PetscInt, PetscBool, PetscInt, PetscInt, PetscFE *);
2482df84da0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateLagrangeByCell(MPI_Comm, PetscInt, PetscInt, DMPolytopeType, PetscInt, PetscInt, PetscFE *);
2499a1a3eb8SMatthew G. Knepley 
2509a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetDimension(PetscFE, PetscInt *);
2510ddb9b0bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetSpatialDimension(PetscFE, PetscInt *);
2529a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetNumComponents(PetscFE, PetscInt);
2539a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetNumComponents(PetscFE, PetscInt *);
2549a769b6bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetTileSizes(PetscFE, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
2559a769b6bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetTileSizes(PetscFE, PetscInt, PetscInt, PetscInt, PetscInt);
2569a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetBasisSpace(PetscFE, PetscSpace);
2579a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetBasisSpace(PetscFE, PetscSpace *);
2589a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetDualSpace(PetscFE, PetscDualSpace);
2599a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetDualSpace(PetscFE, PetscDualSpace *);
260bfa639d9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetQuadrature(PetscFE, PetscQuadrature);
261bfa639d9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetQuadrature(PetscFE, PetscQuadrature *);
2621bafc85dSSander Arens PETSC_EXTERN PetscErrorCode PetscFESetFaceQuadrature(PetscFE, PetscQuadrature);
2631bafc85dSSander Arens PETSC_EXTERN PetscErrorCode PetscFEGetFaceQuadrature(PetscFE, PetscQuadrature *);
2645dc5c000SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECopyQuadrature(PetscFE, PetscFE);
2650ddb9b0bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetNumDof(PetscFE, const PetscInt **);
266ef0bb6c7SMatthew G. Knepley 
267ef0bb6c7SMatthew G. Knepley /* TODO: Need a function to reuse the memory when retabulating the same FE at different points */
268f9244615SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetCellTabulation(PetscFE, PetscInt, PetscTabulation *);
269f9244615SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetFaceTabulation(PetscFE, PetscInt, PetscTabulation *);
270ef0bb6c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetFaceCentroidTabulation(PetscFE, PetscTabulation *);
271ef0bb6c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateTabulation(PetscFE, PetscInt, PetscInt, const PetscReal[], PetscInt, PetscTabulation *);
272ef0bb6c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEComputeTabulation(PetscFE, PetscInt, const PetscReal[], PetscInt, PetscTabulation);
273ef0bb6c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscTabulationDestroy(PetscTabulation *);
274ef0bb6c7SMatthew G. Knepley 
275bceba477SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFERefine(PetscFE, PetscFE *);
276228cfb18SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetHeightSubspace(PetscFE, PetscInt, PetscFE *);
277a0845e3aSMatthew G. Knepley 
278c9ba7969SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateCellGeometry(PetscFE, PetscQuadrature, PetscFEGeom *);
279c9ba7969SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEDestroyCellGeometry(PetscFE, PetscFEGeom *);
2802edcad52SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEPushforward(PetscFE, PetscFEGeom *, PetscInt, PetscScalar[]);
2812edcad52SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEPushforwardGradient(PetscFE, PetscFEGeom *, PetscInt, PetscScalar[]);
282f9244615SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEPushforwardHessian(PetscFE, PetscFEGeom *, PetscInt, PetscScalar[]);
2834bee2e38SMatthew G. Knepley 
2844bee2e38SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrate(PetscDS, PetscInt, PetscInt, PetscFEGeom *, const PetscScalar[], PetscDS, const PetscScalar[], PetscScalar[]);
2854bee2e38SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateBd(PetscDS, PetscInt,
28664c72086SMatthew G. Knepley                                                void (*)(PetscInt, PetscInt, PetscInt,
28764c72086SMatthew G. Knepley                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
28864c72086SMatthew G. Knepley                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
28964c72086SMatthew G. Knepley                                                         PetscReal, const PetscReal[], const PetscReal[], PetscInt, const
29064c72086SMatthew G. Knepley                                                         PetscScalar[], PetscScalar[]),
29164c72086SMatthew G. Knepley                                                PetscInt, PetscFEGeom *, const PetscScalar[], PetscDS, const PetscScalar[], PetscScalar[]);
29206ad1575SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateResidual(PetscDS, PetscFormKey, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscScalar[]);
29306ad1575SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateBdResidual(PetscDS, PetscWeakForm, PetscFormKey, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscScalar[]);
294c2b7495fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateHybridResidual(PetscDS, PetscFormKey, PetscInt, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscScalar[]);
29506ad1575SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateJacobian(PetscDS, PetscFEJacobianType, PetscFormKey, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscReal, PetscScalar[]);
29606ad1575SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateBdJacobian(PetscDS, PetscWeakForm, PetscFormKey, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscReal, PetscScalar[]);
2975fedec97SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateHybridJacobian(PetscDS, PetscFEJacobianType, PetscFormKey, PetscInt, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscReal, PetscScalar[]);
298a0845e3aSMatthew G. Knepley 
29989710940SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECompositeGetMapping(PetscFE, PetscInt *, const PetscReal *[], const PetscReal *[], const PetscReal *[]);
30089710940SMatthew G. Knepley 
3012f5fb066SToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreateHeightTrace(PetscFE, PetscInt, PetscFE *);
3022f5fb066SToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreatePointTrace(PetscFE, PetscInt, PetscFE *);
3032f5fb066SToby Isaac 
304855cd083SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEOpenCLSetRealType(PetscFE, PetscDataType);
305855cd083SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEOpenCLGetRealType(PetscFE, PetscDataType *);
306855cd083SMatthew G. Knepley 
307dbe77d9eSMatthew G. Knepley #endif
308