xref: /petsc/include/petscfe.h (revision 785e854f82a3c614b452fca2cf5ad4f2afe8bdde)
1 /*
2       Objects which encapsulate finite element spaces and operations
3 */
4 #if !defined(__PETSCFE_H)
5 #define __PETSCFE_H
6 #include <petscdm.h>
7 #include <petscdt.h>
8 #include <petscfetypes.h>
9 
10 PETSC_EXTERN PetscErrorCode PetscFEInitializePackage(void);
11 
12 PETSC_EXTERN PetscClassId PETSCSPACE_CLASSID;
13 
14 /*J
15   PetscSpaceType - String with the name of a PETSc linear space
16 
17   Level: beginner
18 
19 .seealso: PetscSpaceSetType(), PetscSpace
20 J*/
21 typedef const char *PetscSpaceType;
22 #define PETSCSPACEPOLYNOMIAL "poly"
23 #define PETSCSPACEDG         "dg"
24 
25 PETSC_EXTERN PetscFunctionList PetscSpaceList;
26 PETSC_EXTERN PetscBool         PetscSpaceRegisterAllCalled;
27 PETSC_EXTERN PetscErrorCode PetscSpaceCreate(MPI_Comm, PetscSpace *);
28 PETSC_EXTERN PetscErrorCode PetscSpaceDestroy(PetscSpace *);
29 PETSC_EXTERN PetscErrorCode PetscSpaceSetType(PetscSpace, PetscSpaceType);
30 PETSC_EXTERN PetscErrorCode PetscSpaceGetType(PetscSpace, PetscSpaceType *);
31 PETSC_EXTERN PetscErrorCode PetscSpaceSetUp(PetscSpace);
32 PETSC_EXTERN PetscErrorCode PetscSpaceSetFromOptions(PetscSpace);
33 PETSC_EXTERN PetscErrorCode PetscSpaceRegister(const char [], PetscErrorCode (*)(PetscSpace));
34 PETSC_EXTERN PetscErrorCode PetscSpaceRegisterAll(void);
35 PETSC_EXTERN PetscErrorCode PetscSpaceRegisterDestroy(void);
36 
37 PETSC_EXTERN PetscErrorCode PetscSpaceGetDimension(PetscSpace, PetscInt *);
38 PETSC_EXTERN PetscErrorCode PetscSpaceSetOrder(PetscSpace, PetscInt);
39 PETSC_EXTERN PetscErrorCode PetscSpaceGetOrder(PetscSpace, PetscInt *);
40 PETSC_EXTERN PetscErrorCode PetscSpaceEvaluate(PetscSpace, PetscInt, const PetscReal[], PetscReal[], PetscReal[], PetscReal[]);
41 
42 PETSC_EXTERN PetscErrorCode PetscSpacePolynomialSetNumVariables(PetscSpace, PetscInt);
43 PETSC_EXTERN PetscErrorCode PetscSpacePolynomialGetNumVariables(PetscSpace, PetscInt *);
44 PETSC_EXTERN PetscErrorCode PetscSpacePolynomialSetSymmetric(PetscSpace, PetscBool);
45 PETSC_EXTERN PetscErrorCode PetscSpacePolynomialGetSymmetric(PetscSpace, PetscBool *);
46 
47 PETSC_EXTERN PetscErrorCode PetscSpaceDGSetQuadrature(PetscSpace, PetscQuadrature);
48 PETSC_EXTERN PetscErrorCode PetscSpaceDGGetQuadrature(PetscSpace, PetscQuadrature *);
49 
50 PETSC_EXTERN PetscClassId PETSCDUALSPACE_CLASSID;
51 
52 /*J
53   PetscDualSpaceType - String with the name of a PETSc dual space
54 
55   Level: beginner
56 
57 .seealso: PetscDualSpaceSetType(), PetscDualSpace
58 J*/
59 typedef const char *PetscDualSpaceType;
60 #define PETSCDUALSPACELAGRANGE "lagrange"
61 
62 PETSC_EXTERN PetscFunctionList PetscDualSpaceList;
63 PETSC_EXTERN PetscBool         PetscDualSpaceRegisterAllCalled;
64 PETSC_EXTERN PetscErrorCode PetscDualSpaceCreate(MPI_Comm, PetscDualSpace *);
65 PETSC_EXTERN PetscErrorCode PetscDualSpaceDestroy(PetscDualSpace *);
66 PETSC_EXTERN PetscErrorCode PetscDualSpaceSetType(PetscDualSpace, PetscDualSpaceType);
67 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetType(PetscDualSpace, PetscDualSpaceType *);
68 PETSC_EXTERN PetscErrorCode PetscDualSpaceSetUp(PetscDualSpace);
69 PETSC_EXTERN PetscErrorCode PetscDualSpaceSetFromOptions(PetscDualSpace);
70 PETSC_EXTERN PetscErrorCode PetscDualSpaceRegister(const char [], PetscErrorCode (*)(PetscDualSpace));
71 PETSC_EXTERN PetscErrorCode PetscDualSpaceRegisterAll(void);
72 PETSC_EXTERN PetscErrorCode PetscDualSpaceRegisterDestroy(void);
73 
74 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetDimension(PetscDualSpace, PetscInt *);
75 PETSC_EXTERN PetscErrorCode PetscDualSpaceSetOrder(PetscDualSpace, PetscInt);
76 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetOrder(PetscDualSpace, PetscInt *);
77 PETSC_EXTERN PetscErrorCode PetscDualSpaceSetDM(PetscDualSpace, DM);
78 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetDM(PetscDualSpace, DM *);
79 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetFunctional(PetscDualSpace, PetscInt, PetscQuadrature *);
80 PETSC_EXTERN PetscErrorCode PetscDualSpaceCreateReferenceCell(PetscDualSpace, PetscInt, PetscBool, DM *);
81 
82 PETSC_EXTERN PetscErrorCode PetscDualSpaceApply(PetscDualSpace, PetscInt, PetscCellGeometry, PetscInt, void (*)(const PetscReal [], PetscScalar *), PetscScalar *);
83 
84 PETSC_EXTERN PetscClassId PETSCFE_CLASSID;
85 
86 /*J
87   PetscFEType - String with the name of a PETSc finite element space
88 
89   Level: beginner
90 
91   Note: Currently, the classes are concerned with the implementation of element integration
92 
93 .seealso: PetscFESetType(), PetscFE
94 J*/
95 typedef const char *PetscFEType;
96 #define PETSCFEBASIC  "basic"
97 #define PETSCFEOPENCL "opencl"
98 
99 PETSC_EXTERN PetscFunctionList PetscFEList;
100 PETSC_EXTERN PetscBool         PetscFERegisterAllCalled;
101 PETSC_EXTERN PetscErrorCode PetscFECreate(MPI_Comm, PetscFE *);
102 PETSC_EXTERN PetscErrorCode PetscFEDestroy(PetscFE *);
103 PETSC_EXTERN PetscErrorCode PetscFESetType(PetscFE, PetscFEType);
104 PETSC_EXTERN PetscErrorCode PetscFEGetType(PetscFE, PetscFEType *);
105 PETSC_EXTERN PetscErrorCode PetscFESetUp(PetscFE);
106 PETSC_EXTERN PetscErrorCode PetscFESetFromOptions(PetscFE);
107 PETSC_EXTERN PetscErrorCode PetscFERegister(const char [], PetscErrorCode (*)(PetscFE));
108 PETSC_EXTERN PetscErrorCode PetscFERegisterAll(void);
109 PETSC_EXTERN PetscErrorCode PetscFERegisterDestroy(void);
110 
111 PETSC_EXTERN PetscErrorCode PetscFEGetDimension(PetscFE, PetscInt *);
112 PETSC_EXTERN PetscErrorCode PetscFEGetSpatialDimension(PetscFE, PetscInt *);
113 PETSC_EXTERN PetscErrorCode PetscFESetNumComponents(PetscFE, PetscInt);
114 PETSC_EXTERN PetscErrorCode PetscFEGetNumComponents(PetscFE, PetscInt *);
115 PETSC_EXTERN PetscErrorCode PetscFEGetTileSizes(PetscFE, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
116 PETSC_EXTERN PetscErrorCode PetscFESetTileSizes(PetscFE, PetscInt, PetscInt, PetscInt, PetscInt);
117 PETSC_EXTERN PetscErrorCode PetscFESetBasisSpace(PetscFE, PetscSpace);
118 PETSC_EXTERN PetscErrorCode PetscFEGetBasisSpace(PetscFE, PetscSpace *);
119 PETSC_EXTERN PetscErrorCode PetscFESetDualSpace(PetscFE, PetscDualSpace);
120 PETSC_EXTERN PetscErrorCode PetscFEGetDualSpace(PetscFE, PetscDualSpace *);
121 PETSC_EXTERN PetscErrorCode PetscFESetQuadrature(PetscFE, PetscQuadrature);
122 PETSC_EXTERN PetscErrorCode PetscFEGetQuadrature(PetscFE, PetscQuadrature *);
123 PETSC_EXTERN PetscErrorCode PetscFEGetNumDof(PetscFE, const PetscInt **);
124 PETSC_EXTERN PetscErrorCode PetscFEGetDefaultTabulation(PetscFE, PetscReal **, PetscReal **, PetscReal **);
125 PETSC_EXTERN PetscErrorCode PetscFEGetTabulation(PetscFE, PetscInt, const PetscReal[], PetscReal **, PetscReal **, PetscReal **);
126 PETSC_EXTERN PetscErrorCode PetscFERestoreTabulation(PetscFE, PetscInt, const PetscReal[], PetscReal **, PetscReal **, PetscReal **);
127 
128 PETSC_EXTERN PetscErrorCode PetscFEIntegrateResidual(PetscFE, PetscInt, PetscInt, PetscFE[], PetscInt, PetscCellGeometry, const PetscScalar[],
129                                                      PetscInt, PetscFE[], const PetscScalar[],
130                                                      void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
131                                                      void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
132                                                      PetscScalar[]);
133 PETSC_EXTERN PetscErrorCode PetscFEIntegrateBdResidual(PetscFE, PetscInt, PetscInt, PetscFE[], PetscInt, PetscCellGeometry, const PetscScalar[],
134                                                        PetscInt, PetscFE[], const PetscScalar[],
135                                                        void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]),
136                                                        void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]),
137                                                        PetscScalar[]);
138 PETSC_EXTERN PetscErrorCode PetscFEIntegrateJacobianAction(PetscFE, PetscInt, PetscInt, PetscFE[], PetscInt, PetscCellGeometry, const PetscScalar[], const PetscScalar[],
139                                                            void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
140                                                            void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
141                                                            void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
142                                                            void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
143                                                            PetscScalar[]);
144 PETSC_EXTERN PetscErrorCode PetscFEIntegrateJacobian(PetscFE, PetscInt, PetscInt, PetscFE[], PetscInt, PetscInt, PetscCellGeometry, const PetscScalar[],
145                                                      PetscInt, PetscFE[], const PetscScalar[],
146                                                      void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
147                                                      void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
148                                                      void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
149                                                      void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
150                                                      PetscScalar[]);
151 
152 PETSC_EXTERN PetscErrorCode PetscFEOpenCLSetRealType(PetscFE, PetscDataType);
153 PETSC_EXTERN PetscErrorCode PetscFEOpenCLGetRealType(PetscFE, PetscDataType *);
154 
155 /* TODO: Should be moved inside DM */
156 typedef struct {
157   PetscFE         *fe;
158   PetscFE         *feAux;
159   void (**f0Funcs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The f_0 functions for each field */
160   void (**f1Funcs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The f_1 functions for each field */
161   void (**g0Funcs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The g_0 functions for each field pair */
162   void (**g1Funcs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The g_1 functions for each field pair */
163   void (**g2Funcs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The g_2 functions for each field pair */
164   void (**g3Funcs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The g_3 functions for each field pair */
165   void (**bcFuncs)(const PetscReal[], PetscScalar *); /* The boundary condition function for each field component */
166   PetscFE         *feBd;
167   void (**f0BdFuncs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]); /* The f_0 functions for each field */
168   void (**f1BdFuncs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]); /* The f_1 functions for each field */
169 } PetscFEM;
170 
171 #endif
172