xref: /petsc/include/petscfe.h (revision b44575274871f84dc4018b7bb4f8b6493c17eb57)
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 
114129dba9SToby Isaac typedef struct _n_PetscFEGeom {
124129dba9SToby Isaac   const PetscReal *xi;
1341418987SMatthew G. Knepley   PetscReal *v;           /* v[Nc*Np*dE]:           The first point in each each in real coordinates */
1441418987SMatthew 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) */
1541418987SMatthew 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) */
1641418987SMatthew G. Knepley   PetscReal *detJ;        /* detJ[Nc*Np]:           The determinant of J, and if it is non-square its the volume change */
1741418987SMatthew G. Knepley   PetscReal *n;           /* n[Nc*Np*dE]:           For faces, the normal to the face in real coordinates */
1841418987SMatthew G. Knepley   PetscInt  (*face)[2];   /* face[Nc][s]:           For faces, the local face number (cone index) for this face in each supporting cell s */
1941418987SMatthew G. Knepley   PetscReal *suppJ[2];    /* sJ[s][Nc*Np*dE*dE]:    For faces, the Jacobian for each supporting cell s */
2041418987SMatthew G. Knepley   PetscReal *suppInvJ[2]; /* sInvJ[s][Nc*Np*dE*dE]: For faces, the inverse Jacobian for each supporting cell s */
2141418987SMatthew G. Knepley   PetscReal *suppDetJ[2]; /* sInvJ[s][Nc*Np*dE*dE]: For faces, the Jacobian determinant for each supporting cell s */
2241418987SMatthew G. Knepley   PetscInt  dim;          /* Topological dimension */
2341418987SMatthew G. Knepley   PetscInt  dimEmbed;     /* Real coordinate dimension */
2441418987SMatthew G. Knepley   PetscInt  numCells;     /* Number of mesh points represented in the arrays */
2541418987SMatthew G. Knepley   PetscInt  numPoints;    /* Number of evaluation points represented in the arrays */
2641418987SMatthew G. Knepley   PetscBool isAffine;     /* Flag for affine transforms */
274129dba9SToby Isaac } PetscFEGeom;
284129dba9SToby Isaac 
29dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEInitializePackage(void);
30dbe77d9eSMatthew G. Knepley 
31dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscClassId PETSCSPACE_CLASSID;
32dbe77d9eSMatthew G. Knepley 
33dbe77d9eSMatthew G. Knepley /*J
34dbe77d9eSMatthew G. Knepley   PetscSpaceType - String with the name of a PETSc linear space
35dbe77d9eSMatthew G. Knepley 
36dbe77d9eSMatthew G. Knepley   Level: beginner
37dbe77d9eSMatthew G. Knepley 
38dbe77d9eSMatthew G. Knepley .seealso: PetscSpaceSetType(), PetscSpace
39dbe77d9eSMatthew G. Knepley J*/
40dbe77d9eSMatthew G. Knepley typedef const char *PetscSpaceType;
41dbe77d9eSMatthew G. Knepley #define PETSCSPACEPOLYNOMIAL "poly"
4236e5648fSToby Isaac #define PETSCSPACETENSOR     "tensor"
439c3cf19fSMatthew G. Knepley #define PETSCSPACEPOINT      "point"
442f5fb066SToby Isaac #define PETSCSPACESUBSPACE   "subspace"
45dbe77d9eSMatthew G. Knepley 
46dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscFunctionList PetscSpaceList;
47dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceCreate(MPI_Comm, PetscSpace *);
48568c31deSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceDestroy(PetscSpace *);
49dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceSetType(PetscSpace, PetscSpaceType);
50dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceGetType(PetscSpace, PetscSpaceType *);
519a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceSetUp(PetscSpace);
52dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceSetFromOptions(PetscSpace);
53fe2efc57SMark PETSC_EXTERN PetscErrorCode PetscSpaceViewFromOptions(PetscSpace,PetscObject,const char[]);
548aec7d55SBarry Smith 
55fd41f350SJed Brown PETSC_EXTERN PetscErrorCode PetscSpaceView(PetscSpace,PetscViewer);
56dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceRegister(const char [], PetscErrorCode (*)(PetscSpace));
57dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceRegisterDestroy(void);
58dbe77d9eSMatthew G. Knepley 
599a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceGetDimension(PetscSpace, PetscInt *);
609c3cf19fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceSetNumComponents(PetscSpace, PetscInt);
619c3cf19fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceGetNumComponents(PetscSpace, PetscInt *);
62157782e2SToby Isaac PETSC_EXTERN PetscErrorCode PetscSpaceSetNumVariables(PetscSpace, PetscInt);
63157782e2SToby Isaac PETSC_EXTERN PetscErrorCode PetscSpaceGetNumVariables(PetscSpace, PetscInt *);
64d39dd5f5SToby Isaac PETSC_EXTERN PetscErrorCode PetscSpaceSetDegree(PetscSpace, PetscInt, PetscInt);
65cdcc6362SToby Isaac PETSC_EXTERN PetscErrorCode PetscSpaceGetDegree(PetscSpace, PetscInt *, PetscInt *);
662bdb15eaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceEvaluate(PetscSpace, PetscInt, const PetscReal[], PetscReal[], PetscReal[], PetscReal[]);
670c16f5adSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceGetHeightSubspace(PetscSpace, PetscInt, PetscSpace *);
689a1a3eb8SMatthew G. Knepley 
692bdb15eaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpacePolynomialSetSymmetric(PetscSpace, PetscBool);
702bdb15eaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpacePolynomialGetSymmetric(PetscSpace, PetscBool *);
71c43d9e4aSToby Isaac PETSC_EXTERN PetscErrorCode PetscSpacePolynomialSetTensor(PetscSpace, PetscBool);
72c43d9e4aSToby Isaac PETSC_EXTERN PetscErrorCode PetscSpacePolynomialGetTensor(PetscSpace, PetscBool *);
732bdb15eaSMatthew G. Knepley 
7436e5648fSToby Isaac PETSC_EXTERN PetscErrorCode PetscSpaceTensorSetNumSubspaces(PetscSpace, PetscInt);
7536e5648fSToby Isaac PETSC_EXTERN PetscErrorCode PetscSpaceTensorGetNumSubspaces(PetscSpace, PetscInt *);
7636e5648fSToby Isaac PETSC_EXTERN PetscErrorCode PetscSpaceTensorSetSubspace(PetscSpace, PetscInt, PetscSpace);
7736e5648fSToby Isaac PETSC_EXTERN PetscErrorCode PetscSpaceTensorGetSubspace(PetscSpace, PetscInt, PetscSpace *);
7836e5648fSToby Isaac 
798049c7f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpacePointGetPoints(PetscSpace, PetscQuadrature *);
808049c7f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpacePointSetPoints(PetscSpace, PetscQuadrature);
819a1a3eb8SMatthew G. Knepley 
822f5fb066SToby Isaac PETSC_EXTERN PetscErrorCode PetscSpaceCreateSubspace(PetscSpace, PetscDualSpace, PetscReal *, PetscReal *, PetscReal *, PetscReal *, PetscCopyMode, PetscSpace *);
832f5fb066SToby Isaac 
84dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscClassId PETSCDUALSPACE_CLASSID;
85dbe77d9eSMatthew G. Knepley 
86dbe77d9eSMatthew G. Knepley /*J
87dbe77d9eSMatthew G. Knepley   PetscDualSpaceType - String with the name of a PETSc dual space
88dbe77d9eSMatthew G. Knepley 
89dbe77d9eSMatthew G. Knepley   Level: beginner
90dbe77d9eSMatthew G. Knepley 
91dbe77d9eSMatthew G. Knepley .seealso: PetscDualSpaceSetType(), PetscDualSpace
92dbe77d9eSMatthew G. Knepley J*/
93dbe77d9eSMatthew G. Knepley typedef const char *PetscDualSpaceType;
94dbe77d9eSMatthew G. Knepley #define PETSCDUALSPACELAGRANGE "lagrange"
95c2765ee2SMatthew G. Knepley #define PETSCDUALSPACESIMPLE   "simple"
96dbe77d9eSMatthew G. Knepley 
97dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscFunctionList PetscDualSpaceList;
98dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceCreate(MPI_Comm, PetscDualSpace *);
99568c31deSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceDestroy(PetscDualSpace *);
100907df5ccSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceDuplicate(PetscDualSpace, PetscDualSpace *);
101dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceSetType(PetscDualSpace, PetscDualSpaceType);
102dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceGetType(PetscDualSpace, PetscDualSpaceType *);
103*b4457527SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceGetUniform(PetscDualSpace, PetscBool *);
1048d2f55e7SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceGetNumDof(PetscDualSpace, const PetscInt **);
105*b4457527SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceGetSection(PetscDualSpace, PetscSection *);
1069a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceSetUp(PetscDualSpace);
107dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceSetFromOptions(PetscDualSpace);
108fe2efc57SMark PETSC_EXTERN PetscErrorCode PetscDualSpaceViewFromOptions(PetscDualSpace,PetscObject,const char[]);
1098aec7d55SBarry Smith 
110fd41f350SJed Brown PETSC_EXTERN PetscErrorCode PetscDualSpaceView(PetscDualSpace,PetscViewer);
111dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceRegister(const char [], PetscErrorCode (*)(PetscDualSpace));
112dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceRegisterDestroy(void);
113dbe77d9eSMatthew G. Knepley 
114dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceGetDimension(PetscDualSpace, PetscInt *);
115*b4457527SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceGetInteriorDimension(PetscDualSpace, PetscInt *);
1169c3cf19fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceSetNumComponents(PetscDualSpace, PetscInt);
1179c3cf19fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceGetNumComponents(PetscDualSpace, PetscInt *);
1189a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceSetOrder(PetscDualSpace, PetscInt);
1199a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceGetOrder(PetscDualSpace, PetscInt *);
1209a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceSetDM(PetscDualSpace, DM);
1219a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceGetDM(PetscDualSpace, DM *);
122ebac44aeSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceGetFunctional(PetscDualSpace, PetscInt, PetscQuadrature *);
1230ddb9b0bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceCreateReferenceCell(PetscDualSpace, PetscInt, PetscBool, DM *);
124fffc2ff8SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceGetSymmetries(PetscDualSpace, const PetscInt ****, const PetscScalar ****);
125dbe77d9eSMatthew G. Knepley 
1264129dba9SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEGeomCreate(PetscQuadrature,PetscInt,PetscInt,PetscBool,PetscFEGeom**);
1274129dba9SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEGeomGetChunk(PetscFEGeom*,PetscInt,PetscInt,PetscFEGeom**);
1284129dba9SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEGeomRestoreChunk(PetscFEGeom*,PetscInt,PetscInt,PetscFEGeom**);
1294129dba9SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEGeomComplete(PetscFEGeom*);
1304129dba9SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEGeomDestroy(PetscFEGeom**);
1314129dba9SToby Isaac 
132c330f8ffSToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceApply(PetscDualSpace, PetscInt, PetscReal, PetscFEGeom *, PetscInt, PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void *, PetscScalar *);
133c330f8ffSToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceApplyDefault(PetscDualSpace, PetscInt, PetscReal, PetscFEGeom *, PetscInt, PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void *, PetscScalar *);
134ebac44aeSMatthew G. Knepley 
135*b4457527SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceGetAllData(PetscDualSpace, PetscQuadrature *, Mat *);
136*b4457527SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceCreateAllDataDefault(PetscDualSpace, PetscQuadrature *, Mat *);
137*b4457527SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceGetInteriorData(PetscDualSpace, PetscQuadrature *, Mat *);
138*b4457527SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceCreateInteriorDataDefault(PetscDualSpace, PetscQuadrature *, Mat *);
139976f670dSToby Isaac 
140edfea91aSToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceApplyAll(PetscDualSpace, const PetscScalar *, PetscScalar *);
141edfea91aSToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceApplyAllDefault(PetscDualSpace, const PetscScalar *, PetscScalar *);
142*b4457527SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceApplyInterior(PetscDualSpace, const PetscScalar *, PetscScalar *);
143*b4457527SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceApplyInteriorDefault(PetscDualSpace, const PetscScalar *, PetscScalar *);
1443f4bc389SToby Isaac 
145*b4457527SToby Isaac 
146*b4457527SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceGetFormDegree(PetscDualSpace, PetscInt *);
147*b4457527SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceSetFormDegree(PetscDualSpace, PetscInt);
1484bee2e38SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceGetDeRahm(PetscDualSpace, PetscInt *);
1494bee2e38SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceTransform(PetscDualSpace, PetscDualSpaceTransformType, PetscBool, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]);
150*b4457527SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpacePullback(PetscDualSpace, PetscFEGeom *, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]);
151*b4457527SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpacePushforward(PetscDualSpace, PetscFEGeom *, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]);
152*b4457527SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpacePushforwardGradient(PetscDualSpace, PetscFEGeom *, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]);
1534bee2e38SMatthew G. Knepley 
15419ded50bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeGetContinuity(PetscDualSpace, PetscBool *);
15519ded50bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeSetContinuity(PetscDualSpace, PetscBool);
15689cad2ffSToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeGetTensor(PetscDualSpace, PetscBool *);
15789cad2ffSToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeSetTensor(PetscDualSpace, PetscBool);
15819ded50bSMatthew G. Knepley 
1596b594bd0SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceGetHeightSubspace(PetscDualSpace,PetscInt,PetscDualSpace *);
1602f5fb066SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceGetPointSubspace(PetscDualSpace,PetscInt,PetscDualSpace *);
161c2765ee2SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceSimpleSetDimension(PetscDualSpace, PetscInt);
162c2765ee2SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceSimpleSetFunctional(PetscDualSpace, PetscInt, PetscQuadrature);
1636b594bd0SToby Isaac 
164dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscClassId PETSCFE_CLASSID;
165dbe77d9eSMatthew G. Knepley 
1660483ade4SMatthew G. Knepley /*J
1670483ade4SMatthew G. Knepley   PetscFEType - String with the name of a PETSc finite element space
1680483ade4SMatthew G. Knepley 
1690483ade4SMatthew G. Knepley   Level: beginner
1700483ade4SMatthew G. Knepley 
1710483ade4SMatthew G. Knepley   Note: Currently, the classes are concerned with the implementation of element integration
1720483ade4SMatthew G. Knepley 
1730483ade4SMatthew G. Knepley .seealso: PetscFESetType(), PetscFE
1740483ade4SMatthew G. Knepley J*/
1750483ade4SMatthew G. Knepley typedef const char *PetscFEType;
1760483ade4SMatthew G. Knepley #define PETSCFEBASIC     "basic"
1770483ade4SMatthew G. Knepley #define PETSCFEOPENCL    "opencl"
178aaf1837cSMatthew G. Knepley #define PETSCFECOMPOSITE "composite"
1790483ade4SMatthew G. Knepley 
1800483ade4SMatthew G. Knepley PETSC_EXTERN PetscFunctionList PetscFEList;
181568c31deSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreate(MPI_Comm, PetscFE *);
182568c31deSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEDestroy(PetscFE *);
1830483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetType(PetscFE, PetscFEType);
1840483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetType(PetscFE, PetscFEType *);
1850483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetUp(PetscFE);
1860483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetFromOptions(PetscFE);
187fe2efc57SMark PETSC_EXTERN PetscErrorCode PetscFEViewFromOptions(PetscFE,PetscObject,const char[]);
1883f6b16c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetName(PetscFE, const char []);
1898aec7d55SBarry Smith 
190fd41f350SJed Brown PETSC_EXTERN PetscErrorCode PetscFEView(PetscFE,PetscViewer);
1910483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFERegister(const char [], PetscErrorCode (*)(PetscFE));
1920483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFERegisterDestroy(void);
193e8d98e54SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateDefault(MPI_Comm, PetscInt, PetscInt, PetscBool, const char [], PetscInt, PetscFE *);
194e703855dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateLagrange(MPI_Comm, PetscInt, PetscInt, PetscBool, PetscInt, PetscInt, PetscFE *);
1959a1a3eb8SMatthew G. Knepley 
1969a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetDimension(PetscFE, PetscInt *);
1970ddb9b0bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetSpatialDimension(PetscFE, PetscInt *);
1989a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetNumComponents(PetscFE, PetscInt);
1999a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetNumComponents(PetscFE, PetscInt *);
2009a769b6bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetTileSizes(PetscFE, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
2019a769b6bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetTileSizes(PetscFE, PetscInt, PetscInt, PetscInt, PetscInt);
2029a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetBasisSpace(PetscFE, PetscSpace);
2039a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetBasisSpace(PetscFE, PetscSpace *);
2049a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetDualSpace(PetscFE, PetscDualSpace);
2059a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetDualSpace(PetscFE, PetscDualSpace *);
206bfa639d9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetQuadrature(PetscFE, PetscQuadrature);
207bfa639d9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetQuadrature(PetscFE, PetscQuadrature *);
2081bafc85dSSander Arens PETSC_EXTERN PetscErrorCode PetscFESetFaceQuadrature(PetscFE, PetscQuadrature);
2091bafc85dSSander Arens PETSC_EXTERN PetscErrorCode PetscFEGetFaceQuadrature(PetscFE, PetscQuadrature *);
2105dc5c000SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECopyQuadrature(PetscFE, PetscFE);
2110ddb9b0bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetNumDof(PetscFE, const PetscInt **);
212ef0bb6c7SMatthew G. Knepley 
213ef0bb6c7SMatthew G. Knepley /* TODO: Need a function to reuse the memory when retabulating the same FE at different points */
214ef0bb6c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetCellTabulation(PetscFE, PetscTabulation *);
215ef0bb6c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetFaceTabulation(PetscFE, PetscTabulation *);
216ef0bb6c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetFaceCentroidTabulation(PetscFE, PetscTabulation *);
217ef0bb6c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateTabulation(PetscFE, PetscInt, PetscInt, const PetscReal[], PetscInt, PetscTabulation *);
218ef0bb6c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEComputeTabulation(PetscFE, PetscInt, const PetscReal[], PetscInt, PetscTabulation);
219ef0bb6c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscTabulationDestroy(PetscTabulation *);
220ef0bb6c7SMatthew G. Knepley 
221bceba477SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFERefine(PetscFE, PetscFE *);
222228cfb18SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetHeightSubspace(PetscFE, PetscInt, PetscFE *);
223a0845e3aSMatthew G. Knepley 
224c9ba7969SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateCellGeometry(PetscFE, PetscQuadrature, PetscFEGeom *);
225c9ba7969SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEDestroyCellGeometry(PetscFE, PetscFEGeom *);
226*b4457527SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEPushforward(PetscFE, PetscFEGeom *, PetscInt, const PetscScalar[], PetscScalar[]);
227*b4457527SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEPushforwardGradient(PetscFE, PetscFEGeom *, PetscInt, const PetscScalar[], PetscScalar[]);
2284bee2e38SMatthew G. Knepley 
2294bee2e38SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrate(PetscDS, PetscInt, PetscInt, PetscFEGeom *, const PetscScalar[], PetscDS, const PetscScalar[], PetscScalar[]);
2304bee2e38SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateBd(PetscDS, PetscInt,
23164c72086SMatthew G. Knepley                                                void (*)(PetscInt, PetscInt, PetscInt,
23264c72086SMatthew G. Knepley                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
23364c72086SMatthew G. Knepley                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
23464c72086SMatthew G. Knepley                                                         PetscReal, const PetscReal[], const PetscReal[], PetscInt, const
23564c72086SMatthew G. Knepley                                                         PetscScalar[], PetscScalar[]),
23664c72086SMatthew G. Knepley                                                PetscInt, PetscFEGeom *, const PetscScalar[], PetscDS, const PetscScalar[], PetscScalar[]);
2374bee2e38SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateResidual(PetscDS, PetscInt, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscScalar[]);
2384bee2e38SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateBdResidual(PetscDS, PetscInt, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscScalar[]);
2394bee2e38SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateJacobian(PetscDS, PetscFEJacobianType, PetscInt, PetscInt, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscReal, PetscScalar[]);
2404bee2e38SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateBdJacobian(PetscDS, PetscInt, PetscInt, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscReal, PetscScalar[]);
241a0845e3aSMatthew G. Knepley 
24289710940SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECompositeGetMapping(PetscFE, PetscInt *, const PetscReal *[], const PetscReal *[], const PetscReal *[]);
24389710940SMatthew G. Knepley 
2442f5fb066SToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreateHeightTrace(PetscFE, PetscInt, PetscFE *);
2452f5fb066SToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreatePointTrace(PetscFE, PetscInt, PetscFE *);
2462f5fb066SToby Isaac 
247855cd083SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEOpenCLSetRealType(PetscFE, PetscDataType);
248855cd083SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEOpenCLGetRealType(PetscFE, PetscDataType *);
249855cd083SMatthew G. Knepley 
250dbe77d9eSMatthew G. Knepley #endif
251