xref: /petsc/include/petscfe.h (revision 4d0b9603208f03c394dfdee779423082bcb528af)
1dbe77d9eSMatthew G. Knepley /*
2dbe77d9eSMatthew G. Knepley       Objects which encapsulate finite element spaces and operations
3dbe77d9eSMatthew G. Knepley */
4dbe77d9eSMatthew G. Knepley #if !defined(__PETSCFE_H)
5dbe77d9eSMatthew G. Knepley #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 
11c0d900a5SMatthew G. Knepley /* Assuming dim <= 3 */
12c0d900a5SMatthew G. Knepley typedef struct {
13c0d900a5SMatthew G. Knepley   PetscReal v0[3];
14c0d900a5SMatthew G. Knepley   PetscReal J[9];
15c0d900a5SMatthew G. Knepley   PetscReal invJ[9];
16c0d900a5SMatthew G. Knepley   PetscReal detJ;
177a1a1bd4SToby Isaac   PetscInt  dim;
187a1a1bd4SToby Isaac   PetscInt  dimEmbed;
19c0d900a5SMatthew G. Knepley } PetscFECellGeom;
20c0d900a5SMatthew G. Knepley 
21*4d0b9603SSander Arens typedef struct {
22*4d0b9603SSander Arens   PetscReal v0[3];
23*4d0b9603SSander Arens   PetscReal J[9];       /* Size could be reduced to 6 */
24*4d0b9603SSander Arens   PetscReal invJ[2][9]; /* invJ in the adjacent cells */
25*4d0b9603SSander Arens   PetscReal detJ;       /* Not really the determinant of J, but \sqrt{\det{J^T J}} */
26*4d0b9603SSander Arens   PetscReal n[3];
27*4d0b9603SSander Arens   PetscInt  dim;
28*4d0b9603SSander Arens   PetscInt  dimEmbed;
29*4d0b9603SSander Arens   PetscInt  face[2];    /* The local face numbers in the adjacent cells */
30*4d0b9603SSander Arens } PetscFEFaceGeom;
31*4d0b9603SSander Arens 
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 
41dbe77d9eSMatthew G. Knepley .seealso: PetscSpaceSetType(), PetscSpace
42dbe77d9eSMatthew G. Knepley J*/
43dbe77d9eSMatthew G. Knepley typedef const char *PetscSpaceType;
44dbe77d9eSMatthew G. Knepley #define PETSCSPACEPOLYNOMIAL "poly"
452bdb15eaSMatthew G. Knepley #define PETSCSPACEDG         "dg"
46dbe77d9eSMatthew G. Knepley 
47dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscFunctionList PetscSpaceList;
48dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceCreate(MPI_Comm, PetscSpace *);
49568c31deSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceDestroy(PetscSpace *);
50dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceSetType(PetscSpace, PetscSpaceType);
51dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceGetType(PetscSpace, PetscSpaceType *);
529a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceSetUp(PetscSpace);
53dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceSetFromOptions(PetscSpace);
548aec7d55SBarry Smith PETSC_STATIC_INLINE PetscErrorCode PetscSpaceViewFromOptions(PetscSpace A,PetscObject B,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,B,name);}
558aec7d55SBarry Smith 
56fd41f350SJed Brown PETSC_EXTERN PetscErrorCode PetscSpaceView(PetscSpace,PetscViewer);
57dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceRegister(const char [], PetscErrorCode (*)(PetscSpace));
58dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceRegisterDestroy(void);
59dbe77d9eSMatthew G. Knepley 
609a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceGetDimension(PetscSpace, PetscInt *);
619a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceSetOrder(PetscSpace, PetscInt);
629a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceGetOrder(PetscSpace, PetscInt *);
632bdb15eaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceEvaluate(PetscSpace, PetscInt, const PetscReal[], PetscReal[], PetscReal[], PetscReal[]);
649a1a3eb8SMatthew G. Knepley 
659a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpacePolynomialSetNumVariables(PetscSpace, PetscInt);
669a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpacePolynomialGetNumVariables(PetscSpace, PetscInt *);
672bdb15eaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpacePolynomialSetSymmetric(PetscSpace, PetscBool);
682bdb15eaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpacePolynomialGetSymmetric(PetscSpace, PetscBool *);
69c43d9e4aSToby Isaac PETSC_EXTERN PetscErrorCode PetscSpacePolynomialSetTensor(PetscSpace, PetscBool);
70c43d9e4aSToby Isaac PETSC_EXTERN PetscErrorCode PetscSpacePolynomialGetTensor(PetscSpace, PetscBool *);
712bdb15eaSMatthew G. Knepley 
722bdb15eaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceDGSetQuadrature(PetscSpace, PetscQuadrature);
732bdb15eaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceDGGetQuadrature(PetscSpace, PetscQuadrature *);
749a1a3eb8SMatthew G. Knepley 
75dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscClassId PETSCDUALSPACE_CLASSID;
76dbe77d9eSMatthew G. Knepley 
77dbe77d9eSMatthew G. Knepley /*J
78dbe77d9eSMatthew G. Knepley   PetscDualSpaceType - String with the name of a PETSc dual space
79dbe77d9eSMatthew G. Knepley 
80dbe77d9eSMatthew G. Knepley   Level: beginner
81dbe77d9eSMatthew G. Knepley 
82dbe77d9eSMatthew G. Knepley .seealso: PetscDualSpaceSetType(), PetscDualSpace
83dbe77d9eSMatthew G. Knepley J*/
84dbe77d9eSMatthew G. Knepley typedef const char *PetscDualSpaceType;
85dbe77d9eSMatthew G. Knepley #define PETSCDUALSPACELAGRANGE "lagrange"
86c2765ee2SMatthew G. Knepley #define PETSCDUALSPACESIMPLE   "simple"
87dbe77d9eSMatthew G. Knepley 
88dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscFunctionList PetscDualSpaceList;
89dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceCreate(MPI_Comm, PetscDualSpace *);
90568c31deSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceDestroy(PetscDualSpace *);
91907df5ccSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceDuplicate(PetscDualSpace, PetscDualSpace *);
92dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceSetType(PetscDualSpace, PetscDualSpaceType);
93dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceGetType(PetscDualSpace, PetscDualSpaceType *);
948d2f55e7SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceGetNumDof(PetscDualSpace, const PetscInt **);
959a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceSetUp(PetscDualSpace);
96dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceSetFromOptions(PetscDualSpace);
978aec7d55SBarry Smith PETSC_STATIC_INLINE PetscErrorCode PetscDualSpaceViewFromOptions(PetscDualSpace A,PetscObject B,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,B,name);}
988aec7d55SBarry Smith 
99fd41f350SJed Brown PETSC_EXTERN PetscErrorCode PetscDualSpaceView(PetscDualSpace,PetscViewer);
100dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceRegister(const char [], PetscErrorCode (*)(PetscDualSpace));
101dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceRegisterDestroy(void);
102dbe77d9eSMatthew G. Knepley 
103dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceGetDimension(PetscDualSpace, PetscInt *);
1049a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceSetOrder(PetscDualSpace, PetscInt);
1059a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceGetOrder(PetscDualSpace, PetscInt *);
1069a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceSetDM(PetscDualSpace, DM);
1079a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceGetDM(PetscDualSpace, DM *);
108ebac44aeSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceGetFunctional(PetscDualSpace, PetscInt, PetscQuadrature *);
1090ddb9b0bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceCreateReferenceCell(PetscDualSpace, PetscInt, PetscBool, DM *);
110fffc2ff8SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceGetSymmetries(PetscDualSpace, const PetscInt ****, const PetscScalar ****);
111dbe77d9eSMatthew G. Knepley 
1120163fd50SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceApply(PetscDualSpace, PetscInt, PetscReal, PetscFECellGeom *, PetscInt, PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void *, PetscScalar *);
113ebac44aeSMatthew G. Knepley 
11419ded50bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeGetContinuity(PetscDualSpace, PetscBool *);
11519ded50bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeSetContinuity(PetscDualSpace, PetscBool);
11689cad2ffSToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeGetTensor(PetscDualSpace, PetscBool *);
11789cad2ffSToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeSetTensor(PetscDualSpace, PetscBool);
11819ded50bSMatthew G. Knepley 
1196b594bd0SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceGetHeightSubspace(PetscDualSpace,PetscInt,PetscDualSpace *);
120c2765ee2SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceSimpleSetDimension(PetscDualSpace, PetscInt);
121c2765ee2SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceSimpleSetFunctional(PetscDualSpace, PetscInt, PetscQuadrature);
1226b594bd0SToby Isaac 
123dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscClassId PETSCFE_CLASSID;
124dbe77d9eSMatthew G. Knepley 
1250483ade4SMatthew G. Knepley /*J
1260483ade4SMatthew G. Knepley   PetscFEType - String with the name of a PETSc finite element space
1270483ade4SMatthew G. Knepley 
1280483ade4SMatthew G. Knepley   Level: beginner
1290483ade4SMatthew G. Knepley 
1300483ade4SMatthew G. Knepley   Note: Currently, the classes are concerned with the implementation of element integration
1310483ade4SMatthew G. Knepley 
1320483ade4SMatthew G. Knepley .seealso: PetscFESetType(), PetscFE
1330483ade4SMatthew G. Knepley J*/
1340483ade4SMatthew G. Knepley typedef const char *PetscFEType;
1350483ade4SMatthew G. Knepley #define PETSCFEBASIC     "basic"
1363dff1f76SMatthew G. Knepley #define PETSCFENONAFFINE "nonaffine"
1370483ade4SMatthew G. Knepley #define PETSCFEOPENCL    "opencl"
138aaf1837cSMatthew G. Knepley #define PETSCFECOMPOSITE "composite"
1390483ade4SMatthew G. Knepley 
1400483ade4SMatthew G. Knepley PETSC_EXTERN PetscFunctionList PetscFEList;
141568c31deSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreate(MPI_Comm, PetscFE *);
142568c31deSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEDestroy(PetscFE *);
1430483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetType(PetscFE, PetscFEType);
1440483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetType(PetscFE, PetscFEType *);
1450483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetUp(PetscFE);
1460483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetFromOptions(PetscFE);
1478aec7d55SBarry Smith PETSC_STATIC_INLINE PetscErrorCode PetscFEViewFromOptions(PetscFE A,PetscObject B,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,B,name);}
1488aec7d55SBarry Smith 
149fd41f350SJed Brown PETSC_EXTERN PetscErrorCode PetscFEView(PetscFE,PetscViewer);
1500483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFERegister(const char [], PetscErrorCode (*)(PetscFE));
1510483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFERegisterDestroy(void);
152a8e86d1aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateDefault(DM, PetscInt, PetscInt, PetscBool, const char [], PetscInt, PetscFE *);
1539a1a3eb8SMatthew G. Knepley 
1549a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetDimension(PetscFE, PetscInt *);
1550ddb9b0bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetSpatialDimension(PetscFE, PetscInt *);
1569a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetNumComponents(PetscFE, PetscInt);
1579a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetNumComponents(PetscFE, PetscInt *);
1589a769b6bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetTileSizes(PetscFE, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
1599a769b6bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetTileSizes(PetscFE, PetscInt, PetscInt, PetscInt, PetscInt);
1609a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetBasisSpace(PetscFE, PetscSpace);
1619a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetBasisSpace(PetscFE, PetscSpace *);
1629a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetDualSpace(PetscFE, PetscDualSpace);
1639a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetDualSpace(PetscFE, PetscDualSpace *);
164bfa639d9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetQuadrature(PetscFE, PetscQuadrature);
165bfa639d9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetQuadrature(PetscFE, PetscQuadrature *);
1661bafc85dSSander Arens PETSC_EXTERN PetscErrorCode PetscFESetFaceQuadrature(PetscFE, PetscQuadrature);
1671bafc85dSSander Arens PETSC_EXTERN PetscErrorCode PetscFEGetFaceQuadrature(PetscFE, PetscQuadrature *);
1680ddb9b0bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetNumDof(PetscFE, const PetscInt **);
169a319912fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetDefaultTabulation(PetscFE, PetscReal **, PetscReal **, PetscReal **);
1701bafc85dSSander Arens PETSC_EXTERN PetscErrorCode PetscFEGetFaceTabulation(PetscFE, PetscReal **, PetscReal **, PetscReal **);
1711bafc85dSSander Arens PETSC_EXTERN PetscErrorCode PetscFEGetFaceCentroidTabulation(PetscFE, PetscReal **);
172a319912fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetTabulation(PetscFE, PetscInt, const PetscReal[], PetscReal **, PetscReal **, PetscReal **);
173a319912fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFERestoreTabulation(PetscFE, PetscInt, const PetscReal[], PetscReal **, PetscReal **, PetscReal **);
174bceba477SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFERefine(PetscFE, PetscFE *);
175a0845e3aSMatthew G. Knepley 
176bbce034cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrate(PetscFE, PetscDS, PetscInt, PetscInt, PetscFECellGeom *, const PetscScalar[], PetscDS, const PetscScalar[], PetscReal[]);
17711dd639bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateResidual(PetscFE, PetscDS, PetscInt, PetscInt, PetscFECellGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscScalar[]);
178*4d0b9603SSander Arens PETSC_EXTERN PetscErrorCode PetscFEIntegrateBdResidual(PetscFE, PetscDS, PetscInt, PetscInt, PetscFEFaceGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscScalar[]);
17911dd639bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateJacobian(PetscFE, PetscDS, PetscFEJacobianType, PetscInt, PetscInt, PetscInt, PetscFECellGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscReal, PetscScalar[]);
180*4d0b9603SSander Arens PETSC_EXTERN PetscErrorCode PetscFEIntegrateBdJacobian(PetscFE, PetscDS, PetscInt, PetscInt, PetscInt, PetscFEFaceGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscReal, PetscScalar[]);
181a0845e3aSMatthew G. Knepley 
18289710940SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECompositeGetMapping(PetscFE, PetscInt *, const PetscReal *[], const PetscReal *[], const PetscReal *[]);
18389710940SMatthew G. Knepley 
184855cd083SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEOpenCLSetRealType(PetscFE, PetscDataType);
185855cd083SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEOpenCLGetRealType(PetscFE, PetscDataType *);
186855cd083SMatthew G. Knepley 
187dbe77d9eSMatthew G. Knepley #endif
188