xref: /petsc/include/petscfe.h (revision fffc2ff84c6db5f18eba2c5d4a8d1460a3f2e795)
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;
17c0d900a5SMatthew G. Knepley   PetscReal n[3];
187a1a1bd4SToby Isaac   PetscInt  dim;
197a1a1bd4SToby Isaac   PetscInt  dimEmbed;
20c0d900a5SMatthew G. Knepley } PetscFECellGeom;
21c0d900a5SMatthew G. Knepley 
22dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEInitializePackage(void);
23dbe77d9eSMatthew G. Knepley 
24dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscClassId PETSCSPACE_CLASSID;
25dbe77d9eSMatthew G. Knepley 
26dbe77d9eSMatthew G. Knepley /*J
27dbe77d9eSMatthew G. Knepley   PetscSpaceType - String with the name of a PETSc linear space
28dbe77d9eSMatthew G. Knepley 
29dbe77d9eSMatthew G. Knepley   Level: beginner
30dbe77d9eSMatthew G. Knepley 
31dbe77d9eSMatthew G. Knepley .seealso: PetscSpaceSetType(), PetscSpace
32dbe77d9eSMatthew G. Knepley J*/
33dbe77d9eSMatthew G. Knepley typedef const char *PetscSpaceType;
34dbe77d9eSMatthew G. Knepley #define PETSCSPACEPOLYNOMIAL "poly"
352bdb15eaSMatthew G. Knepley #define PETSCSPACEDG         "dg"
36dbe77d9eSMatthew G. Knepley 
37dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscFunctionList PetscSpaceList;
38dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceCreate(MPI_Comm, PetscSpace *);
39568c31deSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceDestroy(PetscSpace *);
40dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceSetType(PetscSpace, PetscSpaceType);
41dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceGetType(PetscSpace, PetscSpaceType *);
429a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceSetUp(PetscSpace);
43dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceSetFromOptions(PetscSpace);
448aec7d55SBarry Smith PETSC_STATIC_INLINE PetscErrorCode PetscSpaceViewFromOptions(PetscSpace A,PetscObject B,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,B,name);}
458aec7d55SBarry Smith 
46fd41f350SJed Brown PETSC_EXTERN PetscErrorCode PetscSpaceView(PetscSpace,PetscViewer);
47dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceRegister(const char [], PetscErrorCode (*)(PetscSpace));
48dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceRegisterDestroy(void);
49dbe77d9eSMatthew G. Knepley 
509a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceGetDimension(PetscSpace, PetscInt *);
519a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceSetOrder(PetscSpace, PetscInt);
529a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceGetOrder(PetscSpace, PetscInt *);
532bdb15eaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceEvaluate(PetscSpace, PetscInt, const PetscReal[], PetscReal[], PetscReal[], PetscReal[]);
549a1a3eb8SMatthew G. Knepley 
559a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpacePolynomialSetNumVariables(PetscSpace, PetscInt);
569a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpacePolynomialGetNumVariables(PetscSpace, PetscInt *);
572bdb15eaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpacePolynomialSetSymmetric(PetscSpace, PetscBool);
582bdb15eaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpacePolynomialGetSymmetric(PetscSpace, PetscBool *);
59c43d9e4aSToby Isaac PETSC_EXTERN PetscErrorCode PetscSpacePolynomialSetTensor(PetscSpace, PetscBool);
60c43d9e4aSToby Isaac PETSC_EXTERN PetscErrorCode PetscSpacePolynomialGetTensor(PetscSpace, PetscBool *);
612bdb15eaSMatthew G. Knepley 
622bdb15eaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceDGSetQuadrature(PetscSpace, PetscQuadrature);
632bdb15eaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceDGGetQuadrature(PetscSpace, PetscQuadrature *);
649a1a3eb8SMatthew G. Knepley 
65dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscClassId PETSCDUALSPACE_CLASSID;
66dbe77d9eSMatthew G. Knepley 
67dbe77d9eSMatthew G. Knepley /*J
68dbe77d9eSMatthew G. Knepley   PetscDualSpaceType - String with the name of a PETSc dual space
69dbe77d9eSMatthew G. Knepley 
70dbe77d9eSMatthew G. Knepley   Level: beginner
71dbe77d9eSMatthew G. Knepley 
72dbe77d9eSMatthew G. Knepley .seealso: PetscDualSpaceSetType(), PetscDualSpace
73dbe77d9eSMatthew G. Knepley J*/
74dbe77d9eSMatthew G. Knepley typedef const char *PetscDualSpaceType;
75dbe77d9eSMatthew G. Knepley #define PETSCDUALSPACELAGRANGE "lagrange"
76c2765ee2SMatthew G. Knepley #define PETSCDUALSPACESIMPLE   "simple"
77dbe77d9eSMatthew G. Knepley 
78dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscFunctionList PetscDualSpaceList;
79dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceCreate(MPI_Comm, PetscDualSpace *);
80568c31deSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceDestroy(PetscDualSpace *);
81907df5ccSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceDuplicate(PetscDualSpace, PetscDualSpace *);
82dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceSetType(PetscDualSpace, PetscDualSpaceType);
83dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceGetType(PetscDualSpace, PetscDualSpaceType *);
848d2f55e7SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceGetNumDof(PetscDualSpace, const PetscInt **);
859a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceSetUp(PetscDualSpace);
86dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceSetFromOptions(PetscDualSpace);
878aec7d55SBarry Smith PETSC_STATIC_INLINE PetscErrorCode PetscDualSpaceViewFromOptions(PetscDualSpace A,PetscObject B,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,B,name);}
888aec7d55SBarry Smith 
89fd41f350SJed Brown PETSC_EXTERN PetscErrorCode PetscDualSpaceView(PetscDualSpace,PetscViewer);
90dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceRegister(const char [], PetscErrorCode (*)(PetscDualSpace));
91dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceRegisterDestroy(void);
92dbe77d9eSMatthew G. Knepley 
93dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceGetDimension(PetscDualSpace, PetscInt *);
949a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceSetOrder(PetscDualSpace, PetscInt);
959a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceGetOrder(PetscDualSpace, PetscInt *);
969a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceSetDM(PetscDualSpace, DM);
979a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceGetDM(PetscDualSpace, DM *);
98ebac44aeSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceGetFunctional(PetscDualSpace, PetscInt, PetscQuadrature *);
990ddb9b0bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceCreateReferenceCell(PetscDualSpace, PetscInt, PetscBool, DM *);
100*fffc2ff8SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceGetSymmetries(PetscDualSpace, const PetscInt ****, const PetscScalar ****);
101dbe77d9eSMatthew G. Knepley 
1020163fd50SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceApply(PetscDualSpace, PetscInt, PetscReal, PetscFECellGeom *, PetscInt, PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void *, PetscScalar *);
103ebac44aeSMatthew G. Knepley 
10419ded50bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeGetContinuity(PetscDualSpace, PetscBool *);
10519ded50bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeSetContinuity(PetscDualSpace, PetscBool);
10689cad2ffSToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeGetTensor(PetscDualSpace, PetscBool *);
10789cad2ffSToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeSetTensor(PetscDualSpace, PetscBool);
10819ded50bSMatthew G. Knepley 
1096b594bd0SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceGetHeightSubspace(PetscDualSpace,PetscInt,PetscDualSpace *);
110c2765ee2SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceSimpleSetDimension(PetscDualSpace, PetscInt);
111c2765ee2SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceSimpleSetFunctional(PetscDualSpace, PetscInt, PetscQuadrature);
1126b594bd0SToby Isaac 
113dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscClassId PETSCFE_CLASSID;
114dbe77d9eSMatthew G. Knepley 
1150483ade4SMatthew G. Knepley /*J
1160483ade4SMatthew G. Knepley   PetscFEType - String with the name of a PETSc finite element space
1170483ade4SMatthew G. Knepley 
1180483ade4SMatthew G. Knepley   Level: beginner
1190483ade4SMatthew G. Knepley 
1200483ade4SMatthew G. Knepley   Note: Currently, the classes are concerned with the implementation of element integration
1210483ade4SMatthew G. Knepley 
1220483ade4SMatthew G. Knepley .seealso: PetscFESetType(), PetscFE
1230483ade4SMatthew G. Knepley J*/
1240483ade4SMatthew G. Knepley typedef const char *PetscFEType;
1250483ade4SMatthew G. Knepley #define PETSCFEBASIC     "basic"
1263dff1f76SMatthew G. Knepley #define PETSCFENONAFFINE "nonaffine"
1270483ade4SMatthew G. Knepley #define PETSCFEOPENCL    "opencl"
128aaf1837cSMatthew G. Knepley #define PETSCFECOMPOSITE "composite"
1290483ade4SMatthew G. Knepley 
1300483ade4SMatthew G. Knepley PETSC_EXTERN PetscFunctionList PetscFEList;
131568c31deSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreate(MPI_Comm, PetscFE *);
132568c31deSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEDestroy(PetscFE *);
1330483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetType(PetscFE, PetscFEType);
1340483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetType(PetscFE, PetscFEType *);
1350483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetUp(PetscFE);
1360483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetFromOptions(PetscFE);
1378aec7d55SBarry Smith PETSC_STATIC_INLINE PetscErrorCode PetscFEViewFromOptions(PetscFE A,PetscObject B,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,B,name);}
1388aec7d55SBarry Smith 
139fd41f350SJed Brown PETSC_EXTERN PetscErrorCode PetscFEView(PetscFE,PetscViewer);
1400483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFERegister(const char [], PetscErrorCode (*)(PetscFE));
1410483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFERegisterDestroy(void);
142a8e86d1aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateDefault(DM, PetscInt, PetscInt, PetscBool, const char [], PetscInt, PetscFE *);
1439a1a3eb8SMatthew G. Knepley 
1449a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetDimension(PetscFE, PetscInt *);
1450ddb9b0bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetSpatialDimension(PetscFE, PetscInt *);
1469a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetNumComponents(PetscFE, PetscInt);
1479a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetNumComponents(PetscFE, PetscInt *);
1489a769b6bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetTileSizes(PetscFE, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
1499a769b6bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetTileSizes(PetscFE, PetscInt, PetscInt, PetscInt, PetscInt);
1509a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetBasisSpace(PetscFE, PetscSpace);
1519a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetBasisSpace(PetscFE, PetscSpace *);
1529a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetDualSpace(PetscFE, PetscDualSpace);
1539a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetDualSpace(PetscFE, PetscDualSpace *);
154bfa639d9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetQuadrature(PetscFE, PetscQuadrature);
155bfa639d9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetQuadrature(PetscFE, PetscQuadrature *);
1560ddb9b0bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetNumDof(PetscFE, const PetscInt **);
157a319912fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetDefaultTabulation(PetscFE, PetscReal **, PetscReal **, PetscReal **);
1580966a4bfSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetFaceTabulation(PetscFE, PetscReal **);
159a319912fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetTabulation(PetscFE, PetscInt, const PetscReal[], PetscReal **, PetscReal **, PetscReal **);
160a319912fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFERestoreTabulation(PetscFE, PetscInt, const PetscReal[], PetscReal **, PetscReal **, PetscReal **);
161bceba477SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFERefine(PetscFE, PetscFE *);
162a0845e3aSMatthew G. Knepley 
163bbce034cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrate(PetscFE, PetscDS, PetscInt, PetscInt, PetscFECellGeom *, const PetscScalar[], PetscDS, const PetscScalar[], PetscReal[]);
16411dd639bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateResidual(PetscFE, PetscDS, PetscInt, PetscInt, PetscFECellGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscScalar[]);
16511dd639bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateBdResidual(PetscFE, PetscDS, PetscInt, PetscInt, PetscFECellGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscScalar[]);
16611dd639bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateJacobian(PetscFE, PetscDS, PetscFEJacobianType, PetscInt, PetscInt, PetscInt, PetscFECellGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscReal, PetscScalar[]);
16711dd639bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateBdJacobian(PetscFE, PetscDS, PetscInt, PetscInt, PetscInt, PetscFECellGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscReal, PetscScalar[]);
168a0845e3aSMatthew G. Knepley 
16989710940SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECompositeGetMapping(PetscFE, PetscInt *, const PetscReal *[], const PetscReal *[], const PetscReal *[]);
17089710940SMatthew G. Knepley 
171855cd083SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEOpenCLSetRealType(PetscFE, PetscDataType);
172855cd083SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEOpenCLGetRealType(PetscFE, PetscDataType *);
173855cd083SMatthew G. Knepley 
174dbe77d9eSMatthew G. Knepley #endif
175