xref: /petsc/include/petscds.h (revision 609bdbee21ea3be08735c64dbe00a9ab27759925)
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, 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 (*PetscPointJac)(PetscInt, 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, PetscReal, const PetscReal[], PetscScalar[]);
32 typedef void (*PetscBdPointFunc)(PetscInt, PetscInt, PetscInt,
33                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
34                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
35                                  PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]);
36 typedef void (*PetscBdPointJac)(PetscInt, PetscInt, PetscInt,
37                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
38                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
39                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]);
40 typedef void (*PetscRiemannFunc)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscScalar[], void *);
41 
42 
43 PETSC_EXTERN PetscFunctionList PetscDSList;
44 PETSC_EXTERN PetscErrorCode PetscDSCreate(MPI_Comm, PetscDS *);
45 PETSC_EXTERN PetscErrorCode PetscDSDestroy(PetscDS *);
46 PETSC_EXTERN PetscErrorCode PetscDSSetType(PetscDS, PetscDSType);
47 PETSC_EXTERN PetscErrorCode PetscDSGetType(PetscDS, PetscDSType *);
48 PETSC_EXTERN PetscErrorCode PetscDSSetUp(PetscDS);
49 PETSC_EXTERN PetscErrorCode PetscDSSetFromOptions(PetscDS);
50 PETSC_STATIC_INLINE PetscErrorCode PetscDSViewFromOptions(PetscDS A,PetscObject B,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,B,name);}
51 
52 PETSC_EXTERN PetscErrorCode PetscDSView(PetscDS,PetscViewer);
53 PETSC_EXTERN PetscErrorCode PetscDSRegister(const char [], PetscErrorCode (*)(PetscDS));
54 PETSC_EXTERN PetscErrorCode PetscDSRegisterDestroy(void);
55 
56 PETSC_EXTERN PetscErrorCode PetscDSGetSpatialDimension(PetscDS, PetscInt *);
57 PETSC_EXTERN PetscErrorCode PetscDSGetNumFields(PetscDS, PetscInt *);
58 PETSC_EXTERN PetscErrorCode PetscDSGetTotalDimension(PetscDS, PetscInt *);
59 PETSC_EXTERN PetscErrorCode PetscDSGetTotalComponents(PetscDS, PetscInt *);
60 PETSC_EXTERN PetscErrorCode PetscDSGetFieldIndex(PetscDS, PetscObject, PetscInt *);
61 PETSC_EXTERN PetscErrorCode PetscDSGetFieldSize(PetscDS, PetscInt, PetscInt *);
62 PETSC_EXTERN PetscErrorCode PetscDSGetFieldOffset(PetscDS, PetscInt, PetscInt *);
63 PETSC_EXTERN PetscErrorCode PetscDSGetDimensions(PetscDS, PetscInt *[]);
64 PETSC_EXTERN PetscErrorCode PetscDSGetComponents(PetscDS, PetscInt *[]);
65 PETSC_EXTERN PetscErrorCode PetscDSGetComponentOffset(PetscDS, PetscInt, PetscInt *);
66 PETSC_EXTERN PetscErrorCode PetscDSGetComponentOffsets(PetscDS, PetscInt *[]);
67 PETSC_EXTERN PetscErrorCode PetscDSGetComponentDerivativeOffsets(PetscDS, PetscInt *[]);
68 
69 PETSC_EXTERN PetscErrorCode PetscDSGetDiscretization(PetscDS, PetscInt, PetscObject *);
70 PETSC_EXTERN PetscErrorCode PetscDSSetDiscretization(PetscDS, PetscInt, PetscObject);
71 PETSC_EXTERN PetscErrorCode PetscDSAddDiscretization(PetscDS, PetscObject);
72 PETSC_EXTERN PetscErrorCode PetscDSGetImplicit(PetscDS, PetscInt, PetscBool*);
73 PETSC_EXTERN PetscErrorCode PetscDSSetImplicit(PetscDS, PetscInt, PetscBool);
74 PETSC_EXTERN PetscErrorCode PetscDSGetAdjacency(PetscDS, PetscInt, PetscBool*, PetscBool*);
75 PETSC_EXTERN PetscErrorCode PetscDSSetAdjacency(PetscDS, PetscInt, PetscBool,  PetscBool);
76 PETSC_EXTERN PetscErrorCode PetscDSGetObjective(PetscDS, PetscInt,
77                                                 void (**)(PetscInt, PetscInt, PetscInt,
78                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
79                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
80                                                           PetscReal, const PetscReal[], PetscScalar[]));
81 PETSC_EXTERN PetscErrorCode PetscDSSetObjective(PetscDS, PetscInt,
82                                                 void (*)(PetscInt, PetscInt, PetscInt,
83                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
84                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
85                                                          PetscReal, const PetscReal[], PetscScalar[]));
86 PETSC_EXTERN PetscErrorCode PetscDSGetResidual(PetscDS, PetscInt,
87                                                void (**)(PetscInt, PetscInt, PetscInt,
88                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
89                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
90                                                          PetscReal, const PetscReal[], PetscScalar[]),
91                                                void (**)(PetscInt, PetscInt, PetscInt,
92                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
93                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
94                                                          PetscReal, const PetscReal[], PetscScalar[]));
95 PETSC_EXTERN PetscErrorCode PetscDSSetResidual(PetscDS, PetscInt,
96                                                void (*)(PetscInt, PetscInt, PetscInt,
97                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
98                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
99                                                         PetscReal, const PetscReal[], PetscScalar[]),
100                                                void (*)(PetscInt, PetscInt, PetscInt,
101                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
102                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
103                                                         PetscReal, const PetscReal[], PetscScalar[]));
104 PETSC_EXTERN PetscErrorCode PetscDSHasJacobian(PetscDS, PetscBool *);
105 PETSC_EXTERN PetscErrorCode PetscDSGetJacobian(PetscDS, PetscInt, PetscInt,
106                                                void (**)(PetscInt, PetscInt, PetscInt,
107                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
108                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
109                                                          PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
110                                                void (**)(PetscInt, PetscInt, PetscInt,
111                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
112                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
113                                                          PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
114                                                void (**)(PetscInt, 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                                                          PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
118                                                void (**)(PetscInt, 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                                                          PetscReal, PetscReal, const PetscReal[], PetscScalar[]));
122 PETSC_EXTERN PetscErrorCode PetscDSSetJacobian(PetscDS, PetscInt, PetscInt,
123                                                void (*)(PetscInt, PetscInt, PetscInt,
124                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
125                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
126                                                         PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
127                                                void (*)(PetscInt, PetscInt, PetscInt,
128                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
129                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
130                                                         PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
131                                                void (*)(PetscInt, PetscInt, PetscInt,
132                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
133                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
134                                                         PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
135                                                void (*)(PetscInt, PetscInt, PetscInt,
136                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
137                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
138                                                         PetscReal, PetscReal, const PetscReal[], PetscScalar[]));
139 PETSC_EXTERN PetscErrorCode PetscDSHasJacobianPreconditioner(PetscDS, PetscBool *);
140 PETSC_EXTERN PetscErrorCode PetscDSGetJacobianPreconditioner(PetscDS, PetscInt, PetscInt,
141                                                void (**)(PetscInt, 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                                                          PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
145                                                void (**)(PetscInt, PetscInt, PetscInt,
146                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
147                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
148                                                          PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
149                                                void (**)(PetscInt, PetscInt, PetscInt,
150                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
151                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
152                                                          PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
153                                                void (**)(PetscInt, PetscInt, PetscInt,
154                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
155                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
156                                                          PetscReal, PetscReal, const PetscReal[], PetscScalar[]));
157 PETSC_EXTERN PetscErrorCode PetscDSSetJacobianPreconditioner(PetscDS, PetscInt, PetscInt,
158                                                void (*)(PetscInt, PetscInt, PetscInt,
159                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
160                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
161                                                         PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
162                                                void (*)(PetscInt, PetscInt, PetscInt,
163                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
164                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
165                                                         PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
166                                                void (*)(PetscInt, PetscInt, PetscInt,
167                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
168                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
169                                                         PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
170                                                void (*)(PetscInt, PetscInt, PetscInt,
171                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
172                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
173                                                         PetscReal, PetscReal, const PetscReal[], PetscScalar[]));
174 PETSC_EXTERN PetscErrorCode PetscDSHasDynamicJacobian(PetscDS, PetscBool *);
175 PETSC_EXTERN PetscErrorCode PetscDSGetDynamicJacobian(PetscDS, PetscInt, PetscInt,
176                                                       void (**)(PetscInt, 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                                                                 PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
180                                                       void (**)(PetscInt, 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                                                                 PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
184                                                       void (**)(PetscInt, 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                                                                 PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
188                                                       void (**)(PetscInt, PetscInt, PetscInt,
189                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
190                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
191                                                                 PetscReal, PetscReal, const PetscReal[], PetscScalar[]));
192 PETSC_EXTERN PetscErrorCode PetscDSSetDynamicJacobian(PetscDS, PetscInt, PetscInt,
193                                                       void (*)(PetscInt, PetscInt, PetscInt,
194                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
195                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
196                                                                PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
197                                                       void (*)(PetscInt, PetscInt, PetscInt,
198                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
199                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
200                                                                PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
201                                                       void (*)(PetscInt, PetscInt, PetscInt,
202                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
203                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
204                                                                PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
205                                                       void (*)(PetscInt, PetscInt, PetscInt,
206                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
207                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
208                                                                PetscReal, PetscReal, const PetscReal[], PetscScalar[]));
209 PETSC_EXTERN PetscErrorCode PetscDSGetRiemannSolver(PetscDS, PetscInt,
210                                                     void (**)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscScalar[], void *));
211 PETSC_EXTERN PetscErrorCode PetscDSSetRiemannSolver(PetscDS, PetscInt,
212                                                     void (*)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscScalar[], void *));
213 PETSC_EXTERN PetscErrorCode PetscDSGetContext(PetscDS, PetscInt, void **);
214 PETSC_EXTERN PetscErrorCode PetscDSSetContext(PetscDS, PetscInt, void *);
215 PETSC_EXTERN PetscErrorCode PetscDSGetBdResidual(PetscDS, PetscInt,
216                                                  void (**)(PetscInt, PetscInt, PetscInt,
217                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
218                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
219                                                            PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]),
220                                                  void (**)(PetscInt, PetscInt, PetscInt,
221                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
222                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
223                                                            PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]));
224 PETSC_EXTERN PetscErrorCode PetscDSSetBdResidual(PetscDS, PetscInt,
225                                                  void (*)(PetscInt, PetscInt, PetscInt,
226                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
227                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
228                                                           PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]),
229                                                  void (*)(PetscInt, PetscInt, PetscInt,
230                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
231                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
232                                                           PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]));
233 PETSC_EXTERN PetscErrorCode PetscDSGetBdJacobian(PetscDS, PetscInt, PetscInt,
234                                                  void (**)(PetscInt, PetscInt, PetscInt,
235                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
236                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
237                                                            PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]),
238                                                  void (**)(PetscInt, PetscInt, PetscInt,
239                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
240                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
241                                                            PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]),
242                                                  void (**)(PetscInt, PetscInt, PetscInt,
243                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
244                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
245                                                            PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]),
246                                                  void (**)(PetscInt, PetscInt, PetscInt,
247                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
248                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
249                                                            PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]));
250 PETSC_EXTERN PetscErrorCode PetscDSSetBdJacobian(PetscDS, PetscInt, PetscInt,
251                                                  void (*)(PetscInt, PetscInt, PetscInt,
252                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
253                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
254                                                           PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]),
255                                                  void (*)(PetscInt, PetscInt, PetscInt,
256                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
257                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
258                                                           PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]),
259                                                  void (*)(PetscInt, PetscInt, PetscInt,
260                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
261                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
262                                                           PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]),
263                                                  void (*)(PetscInt, PetscInt, PetscInt,
264                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
265                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
266                                                           PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]));
267 PETSC_EXTERN PetscErrorCode PetscDSGetTabulation(PetscDS, PetscReal ***, PetscReal ***);
268 PETSC_EXTERN PetscErrorCode PetscDSGetFaceTabulation(PetscDS, PetscReal ***, PetscReal ***);
269 PETSC_EXTERN PetscErrorCode PetscDSGetEvaluationArrays(PetscDS, PetscScalar **, PetscScalar **, PetscScalar **);
270 PETSC_EXTERN PetscErrorCode PetscDSGetWeakFormArrays(PetscDS, PetscScalar **, PetscScalar **, PetscScalar **, PetscScalar **, PetscScalar **, PetscScalar **);
271 PETSC_EXTERN PetscErrorCode PetscDSGetRefCoordArrays(PetscDS, PetscReal **, PetscScalar **);
272 PETSC_EXTERN PetscErrorCode PetscDSCopyEquations(PetscDS, PetscDS);
273 PETSC_EXTERN PetscErrorCode PetscDSAddBoundary(PetscDS, DMBoundaryConditionType, const char[], const char[], PetscInt, PetscInt, const PetscInt *, void (*)(), PetscInt, const PetscInt *, void *);
274 PETSC_EXTERN PetscErrorCode PetscDSGetNumBoundary(PetscDS, PetscInt *);
275 PETSC_EXTERN PetscErrorCode PetscDSGetBoundary(PetscDS, PetscInt, DMBoundaryConditionType *, const char **, const char **, PetscInt *, PetscInt *, const PetscInt **, void (**)(), PetscInt *, const PetscInt **, void **);
276 PETSC_EXTERN PetscErrorCode PetscDSCopyBoundary(PetscDS, PetscDS);
277 
278 #endif
279