xref: /petsc/include/petscfe.h (revision 772ec989c7718bb40cfbe6762601ac78951cc001)
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 PetscSpaceViewFromOptions(PetscSpace,const char[],const char[]);
34 PETSC_EXTERN PetscErrorCode PetscSpaceView(PetscSpace,PetscViewer);
35 PETSC_EXTERN PetscErrorCode PetscSpaceRegister(const char [], PetscErrorCode (*)(PetscSpace));
36 PETSC_EXTERN PetscErrorCode PetscSpaceRegisterAll(void);
37 PETSC_EXTERN PetscErrorCode PetscSpaceRegisterDestroy(void);
38 
39 PETSC_EXTERN PetscErrorCode PetscSpaceGetDimension(PetscSpace, PetscInt *);
40 PETSC_EXTERN PetscErrorCode PetscSpaceSetOrder(PetscSpace, PetscInt);
41 PETSC_EXTERN PetscErrorCode PetscSpaceGetOrder(PetscSpace, PetscInt *);
42 PETSC_EXTERN PetscErrorCode PetscSpaceEvaluate(PetscSpace, PetscInt, const PetscReal[], PetscReal[], PetscReal[], PetscReal[]);
43 
44 PETSC_EXTERN PetscErrorCode PetscSpacePolynomialSetNumVariables(PetscSpace, PetscInt);
45 PETSC_EXTERN PetscErrorCode PetscSpacePolynomialGetNumVariables(PetscSpace, PetscInt *);
46 PETSC_EXTERN PetscErrorCode PetscSpacePolynomialSetSymmetric(PetscSpace, PetscBool);
47 PETSC_EXTERN PetscErrorCode PetscSpacePolynomialGetSymmetric(PetscSpace, PetscBool *);
48 
49 PETSC_EXTERN PetscErrorCode PetscSpaceDGSetQuadrature(PetscSpace, PetscQuadrature);
50 PETSC_EXTERN PetscErrorCode PetscSpaceDGGetQuadrature(PetscSpace, PetscQuadrature *);
51 
52 PETSC_EXTERN PetscClassId PETSCDUALSPACE_CLASSID;
53 
54 /*J
55   PetscDualSpaceType - String with the name of a PETSc dual space
56 
57   Level: beginner
58 
59 .seealso: PetscDualSpaceSetType(), PetscDualSpace
60 J*/
61 typedef const char *PetscDualSpaceType;
62 #define PETSCDUALSPACELAGRANGE "lagrange"
63 
64 PETSC_EXTERN PetscFunctionList PetscDualSpaceList;
65 PETSC_EXTERN PetscBool         PetscDualSpaceRegisterAllCalled;
66 PETSC_EXTERN PetscErrorCode PetscDualSpaceCreate(MPI_Comm, PetscDualSpace *);
67 PETSC_EXTERN PetscErrorCode PetscDualSpaceDestroy(PetscDualSpace *);
68 PETSC_EXTERN PetscErrorCode PetscDualSpaceSetType(PetscDualSpace, PetscDualSpaceType);
69 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetType(PetscDualSpace, PetscDualSpaceType *);
70 PETSC_EXTERN PetscErrorCode PetscDualSpaceSetUp(PetscDualSpace);
71 PETSC_EXTERN PetscErrorCode PetscDualSpaceSetFromOptions(PetscDualSpace);
72 PETSC_EXTERN PetscErrorCode PetscDualSpaceViewFromOptions(PetscDualSpace,const char[],const char[]);
73 PETSC_EXTERN PetscErrorCode PetscDualSpaceView(PetscDualSpace,PetscViewer);
74 PETSC_EXTERN PetscErrorCode PetscDualSpaceRegister(const char [], PetscErrorCode (*)(PetscDualSpace));
75 PETSC_EXTERN PetscErrorCode PetscDualSpaceRegisterAll(void);
76 PETSC_EXTERN PetscErrorCode PetscDualSpaceRegisterDestroy(void);
77 
78 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetDimension(PetscDualSpace, PetscInt *);
79 PETSC_EXTERN PetscErrorCode PetscDualSpaceSetOrder(PetscDualSpace, PetscInt);
80 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetOrder(PetscDualSpace, PetscInt *);
81 PETSC_EXTERN PetscErrorCode PetscDualSpaceSetDM(PetscDualSpace, DM);
82 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetDM(PetscDualSpace, DM *);
83 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetFunctional(PetscDualSpace, PetscInt, PetscQuadrature *);
84 PETSC_EXTERN PetscErrorCode PetscDualSpaceCreateReferenceCell(PetscDualSpace, PetscInt, PetscBool, DM *);
85 
86 PETSC_EXTERN PetscErrorCode PetscDualSpaceApply(PetscDualSpace, PetscInt, PetscCellGeometry, PetscInt, void (*)(const PetscReal [], PetscScalar *, void *), void *, PetscScalar *);
87 
88 PETSC_EXTERN PetscClassId PETSCFE_CLASSID;
89 
90 /*J
91   PetscFEType - String with the name of a PETSc finite element space
92 
93   Level: beginner
94 
95   Note: Currently, the classes are concerned with the implementation of element integration
96 
97 .seealso: PetscFESetType(), PetscFE
98 J*/
99 typedef const char *PetscFEType;
100 #define PETSCFEBASIC     "basic"
101 #define PETSCFENONAFFINE "nonaffine"
102 #define PETSCFEOPENCL    "opencl"
103 #define PETSCFECOMPOSITE "composite"
104 
105 PETSC_EXTERN PetscFunctionList PetscFEList;
106 PETSC_EXTERN PetscBool         PetscFERegisterAllCalled;
107 PETSC_EXTERN PetscErrorCode PetscFECreate(MPI_Comm, PetscFE *);
108 PETSC_EXTERN PetscErrorCode PetscFEDestroy(PetscFE *);
109 PETSC_EXTERN PetscErrorCode PetscFESetType(PetscFE, PetscFEType);
110 PETSC_EXTERN PetscErrorCode PetscFEGetType(PetscFE, PetscFEType *);
111 PETSC_EXTERN PetscErrorCode PetscFESetUp(PetscFE);
112 PETSC_EXTERN PetscErrorCode PetscFESetFromOptions(PetscFE);
113 PETSC_EXTERN PetscErrorCode PetscFEViewFromOptions(PetscFE,const char[],const char[]);
114 PETSC_EXTERN PetscErrorCode PetscFEView(PetscFE,PetscViewer);
115 PETSC_EXTERN PetscErrorCode PetscFERegister(const char [], PetscErrorCode (*)(PetscFE));
116 PETSC_EXTERN PetscErrorCode PetscFERegisterAll(void);
117 PETSC_EXTERN PetscErrorCode PetscFERegisterDestroy(void);
118 PETSC_EXTERN PetscErrorCode PetscFECreateDefault(DM, PetscInt, PetscInt, PetscBool, const char [], PetscInt, PetscFE *);
119 
120 PETSC_EXTERN PetscErrorCode PetscFEGetDimension(PetscFE, PetscInt *);
121 PETSC_EXTERN PetscErrorCode PetscFEGetSpatialDimension(PetscFE, PetscInt *);
122 PETSC_EXTERN PetscErrorCode PetscFESetNumComponents(PetscFE, PetscInt);
123 PETSC_EXTERN PetscErrorCode PetscFEGetNumComponents(PetscFE, PetscInt *);
124 PETSC_EXTERN PetscErrorCode PetscFEGetTileSizes(PetscFE, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
125 PETSC_EXTERN PetscErrorCode PetscFESetTileSizes(PetscFE, PetscInt, PetscInt, PetscInt, PetscInt);
126 PETSC_EXTERN PetscErrorCode PetscFESetBasisSpace(PetscFE, PetscSpace);
127 PETSC_EXTERN PetscErrorCode PetscFEGetBasisSpace(PetscFE, PetscSpace *);
128 PETSC_EXTERN PetscErrorCode PetscFESetDualSpace(PetscFE, PetscDualSpace);
129 PETSC_EXTERN PetscErrorCode PetscFEGetDualSpace(PetscFE, PetscDualSpace *);
130 PETSC_EXTERN PetscErrorCode PetscFESetQuadrature(PetscFE, PetscQuadrature);
131 PETSC_EXTERN PetscErrorCode PetscFEGetQuadrature(PetscFE, PetscQuadrature *);
132 PETSC_EXTERN PetscErrorCode PetscFEGetNumDof(PetscFE, const PetscInt **);
133 PETSC_EXTERN PetscErrorCode PetscFEGetDefaultTabulation(PetscFE, PetscReal **, PetscReal **, PetscReal **);
134 PETSC_EXTERN PetscErrorCode PetscFEGetTabulation(PetscFE, PetscInt, const PetscReal[], PetscReal **, PetscReal **, PetscReal **);
135 PETSC_EXTERN PetscErrorCode PetscFERestoreTabulation(PetscFE, PetscInt, const PetscReal[], PetscReal **, PetscReal **, PetscReal **);
136 PETSC_EXTERN PetscErrorCode PetscFERefine(PetscFE, PetscFE *);
137 
138 PETSC_EXTERN PetscErrorCode PetscFEIntegrateResidual(PetscFE, PetscInt, PetscInt, PetscFE[], PetscInt, PetscCellGeometry, const PetscScalar[],
139                                                      PetscInt, PetscFE[], const 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                                                      PetscScalar[]);
143 PETSC_EXTERN PetscErrorCode PetscFEIntegrateBdResidual(PetscFE, PetscInt, PetscInt, PetscFE[], PetscInt, PetscCellGeometry, const PetscScalar[],
144                                                        PetscInt, PetscFE[], const PetscScalar[],
145                                                        void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]),
146                                                        void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]),
147                                                        PetscScalar[]);
148 PETSC_EXTERN PetscErrorCode PetscFEIntegrateJacobianAction(PetscFE, PetscInt, PetscInt, PetscFE[], PetscInt, PetscCellGeometry, const PetscScalar[], const PetscScalar[],
149                                                            void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
150                                                            void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
151                                                            void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
152                                                            void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
153                                                            PetscScalar[]);
154 PETSC_EXTERN PetscErrorCode PetscFEIntegrateJacobian(PetscFE, PetscInt, PetscInt, PetscFE[], PetscInt, PetscInt, PetscCellGeometry, const PetscScalar[],
155                                                      PetscInt, PetscFE[], const PetscScalar[],
156                                                      void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
157                                                      void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
158                                                      void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
159                                                      void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
160                                                      PetscScalar[]);
161 PETSC_EXTERN PetscErrorCode PetscFEIntegrateBdJacobian(PetscFE, PetscInt, PetscInt, PetscFE[], PetscInt, PetscInt, PetscCellGeometry, const PetscScalar[],
162                                                        PetscInt, PetscFE[], const PetscScalar[],
163                                                        void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]),
164                                                        void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]),
165                                                        void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]),
166                                                        void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]),
167                                                        PetscScalar[]);
168 PETSC_EXTERN PetscErrorCode PetscFEIntegrateIFunction(PetscFE, PetscInt, PetscInt, PetscFE[], PetscInt, PetscCellGeometry, const PetscScalar[], const PetscScalar[],
169                                                       PetscInt, PetscFE[], const PetscScalar[],
170                                                       void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
171                                                       void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]),
172                                                       PetscScalar[]);
173 PETSC_EXTERN PetscErrorCode PetscFEIntegrateBdIFunction(PetscFE, PetscInt, PetscInt, PetscFE[], PetscInt, PetscCellGeometry, const PetscScalar[], const PetscScalar[],
174                                                         PetscInt, PetscFE[], const PetscScalar[],
175                                                         void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]),
176                                                         void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]),
177                                                         PetscScalar[]);
178 
179 PETSC_EXTERN PetscErrorCode PetscFEOpenCLSetRealType(PetscFE, PetscDataType);
180 PETSC_EXTERN PetscErrorCode PetscFEOpenCLGetRealType(PetscFE, PetscDataType *);
181 
182 /* TODO: Should be moved inside DM */
183 typedef struct {
184   PetscFE         *fe;
185   PetscFE         *feAux;
186   void (**f0Funcs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The f_0 functions for each field */
187   void (**f1Funcs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The f_1 functions for each field */
188   void (**f0IFuncs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The f_0 functions for each field */
189   void (**f1IFuncs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The f_1 functions for each field */
190   void (**g0Funcs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The g_0 functions for each field pair */
191   void (**g1Funcs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The g_1 functions for each field pair */
192   void (**g2Funcs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The g_2 functions for each field pair */
193   void (**g3Funcs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The g_3 functions for each field pair */
194   void (**bcFuncs)(const PetscReal[], PetscScalar *, void *); /* The boundary condition function for each field component */
195   void **bcCtxs; /* Contexts for each boundary condition function, or null if all contexts are null */
196   PetscFE         *feBd;
197   void (**f0BdFuncs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]); /* The f_0 functions for each field */
198   void (**f1BdFuncs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]); /* The f_1 functions for each field */
199   void (**f0BdIFuncs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]); /* The f_0 functions for each field */
200   void (**f1BdIFuncs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]); /* The f_1 functions for each field */
201   void (**g0BdFuncs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]); /* The g_0 functions for each field pair */
202   void (**g1BdFuncs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]); /* The g_1 functions for each field pair */
203   void (**g2BdFuncs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]); /* The g_2 functions for each field pair */
204   void (**g3BdFuncs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]); /* The g_3 functions for each field pair */
205 } PetscFEM;
206 
207 #endif
208