xref: /petsc/include/petscds.h (revision 420e96eda417c4a9936f208d8e146b6072205ea6)
1 /*
2       Objects which encapsulate discretizations+continuum residuals
3 */
4 #if !defined(__PETSCDS_H)
5 #define __PETSCDS_H
6 #include <petscfe.h>
7 #include <petscfv.h>
8 #include <petscdstypes.h>
9 
10 PETSC_EXTERN PetscErrorCode PetscDSInitializePackage(void);
11 
12 PETSC_EXTERN PetscClassId PETSCDS_CLASSID;
13 
14 /*J
15   PetscDSType - String with the name of a PETSc discrete system
16 
17   Level: beginner
18 
19 .seealso: PetscDSSetType(), PetscDS
20 J*/
21 typedef const char *PetscDSType;
22 #define PETSCDSBASIC "basic"
23 
24 typedef void (*PetscPointFunc)(PetscInt, PetscInt,
25                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
26                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
27                                PetscReal, const PetscReal[], PetscScalar[]);
28 typedef void (*PetscBdPointFunc)(PetscInt, PetscInt,
29                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
30                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
31                                  PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]);
32 typedef void (*PetscRiemannFunc)(const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscScalar[], void *);
33 
34 
35 PETSC_EXTERN PetscFunctionList PetscDSList;
36 PETSC_EXTERN PetscErrorCode PetscDSCreate(MPI_Comm, PetscDS *);
37 PETSC_EXTERN PetscErrorCode PetscDSDestroy(PetscDS *);
38 PETSC_EXTERN PetscErrorCode PetscDSSetType(PetscDS, PetscDSType);
39 PETSC_EXTERN PetscErrorCode PetscDSGetType(PetscDS, PetscDSType *);
40 PETSC_EXTERN PetscErrorCode PetscDSSetUp(PetscDS);
41 PETSC_EXTERN PetscErrorCode PetscDSSetFromOptions(PetscDS);
42 PETSC_EXTERN PetscErrorCode PetscDSViewFromOptions(PetscDS,const char[],const char[]);
43 PETSC_EXTERN PetscErrorCode PetscDSView(PetscDS,PetscViewer);
44 PETSC_EXTERN PetscErrorCode PetscDSRegister(const char [], PetscErrorCode (*)(PetscDS));
45 PETSC_EXTERN PetscErrorCode PetscDSRegisterDestroy(void);
46 
47 PETSC_EXTERN PetscErrorCode PetscDSGetSpatialDimension(PetscDS, PetscInt *);
48 PETSC_EXTERN PetscErrorCode PetscDSGetNumFields(PetscDS, PetscInt *);
49 PETSC_EXTERN PetscErrorCode PetscDSGetTotalDimension(PetscDS, PetscInt *);
50 PETSC_EXTERN PetscErrorCode PetscDSGetTotalBdDimension(PetscDS, PetscInt *);
51 PETSC_EXTERN PetscErrorCode PetscDSGetTotalComponents(PetscDS, PetscInt *);
52 PETSC_EXTERN PetscErrorCode PetscDSGetFieldOffset(PetscDS, PetscInt, PetscInt *);
53 PETSC_EXTERN PetscErrorCode PetscDSGetBdFieldOffset(PetscDS, PetscInt, PetscInt *);
54 PETSC_EXTERN PetscErrorCode PetscDSGetComponentOffset(PetscDS, PetscInt, PetscInt *);
55 PETSC_EXTERN PetscErrorCode PetscDSGetComponentOffsets(PetscDS, PetscInt *[]);
56 PETSC_EXTERN PetscErrorCode PetscDSGetComponentDerivativeOffsets(PetscDS, PetscInt *[]);
57 
58 PETSC_EXTERN PetscErrorCode PetscDSGetDiscretization(PetscDS, PetscInt, PetscObject *);
59 PETSC_EXTERN PetscErrorCode PetscDSSetDiscretization(PetscDS, PetscInt, PetscObject);
60 PETSC_EXTERN PetscErrorCode PetscDSAddDiscretization(PetscDS, PetscObject);
61 PETSC_EXTERN PetscErrorCode PetscDSGetBdDiscretization(PetscDS, PetscInt, PetscObject *);
62 PETSC_EXTERN PetscErrorCode PetscDSSetBdDiscretization(PetscDS, PetscInt, PetscObject);
63 PETSC_EXTERN PetscErrorCode PetscDSAddBdDiscretization(PetscDS, PetscObject);
64 PETSC_EXTERN PetscErrorCode PetscDSGetImplicit(PetscDS, PetscInt, PetscBool*);
65 PETSC_EXTERN PetscErrorCode PetscDSSetImplicit(PetscDS, PetscInt, PetscBool);
66 PETSC_EXTERN PetscErrorCode PetscDSGetAdjacency(PetscDS, PetscInt, PetscBool*, PetscBool*);
67 PETSC_EXTERN PetscErrorCode PetscDSSetAdjacency(PetscDS, PetscInt, PetscBool,  PetscBool);
68 PETSC_EXTERN PetscErrorCode PetscDSGetObjective(PetscDS, PetscInt,
69                                                 void (**)(PetscInt, PetscInt,
70                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
71                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
72                                                           const PetscReal, const PetscReal[], PetscScalar[]));
73 PETSC_EXTERN PetscErrorCode PetscDSSetObjective(PetscDS, PetscInt,
74                                                 void (*)(PetscInt, PetscInt,
75                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
76                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
77                                                          const PetscReal, const PetscReal[], PetscScalar[]));
78 PETSC_EXTERN PetscErrorCode PetscDSGetResidual(PetscDS, PetscInt,
79                                                void (**)(PetscInt, PetscInt,
80                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
81                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
82                                                          const PetscReal, const PetscReal[], PetscScalar[]),
83                                                void (**)(PetscInt, PetscInt,
84                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
85                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
86                                                          const PetscReal, const PetscReal[], PetscScalar[]));
87 PETSC_EXTERN PetscErrorCode PetscDSSetResidual(PetscDS, PetscInt,
88                                                void (*)(PetscInt, PetscInt,
89                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
90                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
91                                                         const PetscReal, const PetscReal[], PetscScalar[]),
92                                                void (*)(PetscInt, PetscInt,
93                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
94                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
95                                                         const PetscReal, const PetscReal[], PetscScalar[]));
96 PETSC_EXTERN PetscErrorCode PetscDSGetJacobian(PetscDS, PetscInt, PetscInt,
97                                                void (**)(PetscInt, PetscInt,
98                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
99                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
100                                                          const PetscReal, const PetscReal[], PetscScalar[]),
101                                                void (**)(PetscInt, PetscInt,
102                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
103                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
104                                                          const PetscReal, const PetscReal[], PetscScalar[]),
105                                                void (**)(PetscInt, PetscInt,
106                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
107                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
108                                                          const PetscReal, const PetscReal[], PetscScalar[]),
109                                                void (**)(PetscInt, PetscInt,
110                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
111                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
112                                                          const PetscReal, const PetscReal[], PetscScalar[]));
113 PETSC_EXTERN PetscErrorCode PetscDSSetJacobian(PetscDS, PetscInt, PetscInt,
114                                                void (*)(PetscInt, PetscInt,
115                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
116                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
117                                                         const PetscReal, const PetscReal[], PetscScalar[]),
118                                                void (*)(PetscInt, PetscInt,
119                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
120                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
121                                                         const PetscReal, const PetscReal[], PetscScalar[]),
122                                                void (*)(PetscInt, PetscInt,
123                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
124                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
125                                                         const PetscReal, const PetscReal[], PetscScalar[]),
126                                                void (*)(PetscInt, PetscInt,
127                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
128                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
129                                                         const PetscReal, const PetscReal[], PetscScalar[]));
130 PETSC_EXTERN PetscErrorCode PetscDSGetRiemannSolver(PetscDS, PetscInt,
131                                                     void (**)(const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscScalar[], void *));
132 PETSC_EXTERN PetscErrorCode PetscDSSetRiemannSolver(PetscDS, PetscInt,
133                                                     void (*)(const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscScalar[], void *));
134 PETSC_EXTERN PetscErrorCode PetscDSGetContext(PetscDS, PetscInt, void **);
135 PETSC_EXTERN PetscErrorCode PetscDSSetContext(PetscDS, PetscInt, void *);
136 PETSC_EXTERN PetscErrorCode PetscDSGetBdResidual(PetscDS, PetscInt,
137                                                  void (**)(PetscInt, PetscInt,
138                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
139                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
140                                                            const PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]),
141                                                  void (**)(PetscInt, PetscInt,
142                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
143                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
144                                                            const PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]));
145 PETSC_EXTERN PetscErrorCode PetscDSSetBdResidual(PetscDS, PetscInt,
146                                                  void (*)(PetscInt, PetscInt,
147                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
148                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
149                                                           const PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]),
150                                                  void (*)(PetscInt, PetscInt,
151                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
152                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
153                                                           const PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]));
154 PETSC_EXTERN PetscErrorCode PetscDSGetBdJacobian(PetscDS, PetscInt, PetscInt,
155                                                  void (**)(PetscInt, PetscInt,
156                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
157                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
158                                                            const PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]),
159                                                  void (**)(PetscInt, PetscInt,
160                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
161                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
162                                                            const PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]),
163                                                  void (**)(PetscInt, PetscInt,
164                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
165                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
166                                                            const PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]),
167                                                  void (**)(PetscInt, PetscInt,
168                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
169                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
170                                                            const PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]));
171 PETSC_EXTERN PetscErrorCode PetscDSSetBdJacobian(PetscDS, PetscInt, PetscInt,
172                                                  void (*)(PetscInt, PetscInt,
173                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
174                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
175                                                           const PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]),
176                                                  void (*)(PetscInt, PetscInt,
177                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
178                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
179                                                           const PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]),
180                                                  void (*)(PetscInt, PetscInt,
181                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
182                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
183                                                           const PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]),
184                                                  void (*)(PetscInt, PetscInt,
185                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
186                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
187                                                           const PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]));
188 PETSC_EXTERN PetscErrorCode PetscDSGetTabulation(PetscDS, PetscReal ***, PetscReal ***);
189 PETSC_EXTERN PetscErrorCode PetscDSGetBdTabulation(PetscDS, PetscReal ***, PetscReal ***);
190 PETSC_EXTERN PetscErrorCode PetscDSGetEvaluationArrays(PetscDS, PetscScalar **, PetscScalar **, PetscScalar **);
191 PETSC_EXTERN PetscErrorCode PetscDSGetWeakFormArrays(PetscDS, PetscScalar **, PetscScalar **, PetscScalar **, PetscScalar **, PetscScalar **, PetscScalar **);
192 PETSC_EXTERN PetscErrorCode PetscDSGetRefCoordArrays(PetscDS, PetscReal **, PetscScalar **);
193 
194 #endif
195