xref: /petsc/include/petscds.h (revision 6e2a8d7422d05df32a0fe3bd7da52877537e6496)
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 PetscDSGetTotalBdDimension(PetscDS, PetscInt *);
60 PETSC_EXTERN PetscErrorCode PetscDSGetTotalComponents(PetscDS, PetscInt *);
61 PETSC_EXTERN PetscErrorCode PetscDSGetFieldOffset(PetscDS, PetscInt, PetscInt *);
62 PETSC_EXTERN PetscErrorCode PetscDSGetBdFieldOffset(PetscDS, PetscInt, PetscInt *);
63 PETSC_EXTERN PetscErrorCode PetscDSGetComponentOffset(PetscDS, PetscInt, PetscInt *);
64 PETSC_EXTERN PetscErrorCode PetscDSGetComponentOffsets(PetscDS, PetscInt *[]);
65 PETSC_EXTERN PetscErrorCode PetscDSGetComponentBdOffsets(PetscDS, PetscInt *[]);
66 PETSC_EXTERN PetscErrorCode PetscDSGetComponentDerivativeOffsets(PetscDS, PetscInt *[]);
67 PETSC_EXTERN PetscErrorCode PetscDSGetComponentBdDerivativeOffsets(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 PetscDSGetBdDiscretization(PetscDS, PetscInt, PetscObject *);
73 PETSC_EXTERN PetscErrorCode PetscDSSetBdDiscretization(PetscDS, PetscInt, PetscObject);
74 PETSC_EXTERN PetscErrorCode PetscDSAddBdDiscretization(PetscDS, PetscObject);
75 PETSC_EXTERN PetscErrorCode PetscDSGetImplicit(PetscDS, PetscInt, PetscBool*);
76 PETSC_EXTERN PetscErrorCode PetscDSSetImplicit(PetscDS, PetscInt, PetscBool);
77 PETSC_EXTERN PetscErrorCode PetscDSGetAdjacency(PetscDS, PetscInt, PetscBool*, PetscBool*);
78 PETSC_EXTERN PetscErrorCode PetscDSSetAdjacency(PetscDS, PetscInt, PetscBool,  PetscBool);
79 PETSC_EXTERN PetscErrorCode PetscDSGetObjective(PetscDS, PetscInt,
80                                                 void (**)(PetscInt, PetscInt, PetscInt,
81                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
82                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
83                                                           PetscReal, const PetscReal[], PetscScalar[]));
84 PETSC_EXTERN PetscErrorCode PetscDSSetObjective(PetscDS, PetscInt,
85                                                 void (*)(PetscInt, PetscInt, PetscInt,
86                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
87                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
88                                                          PetscReal, const PetscReal[], PetscScalar[]));
89 PETSC_EXTERN PetscErrorCode PetscDSGetResidual(PetscDS, PetscInt,
90                                                void (**)(PetscInt, PetscInt, PetscInt,
91                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
92                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
93                                                          PetscReal, const PetscReal[], PetscScalar[]),
94                                                void (**)(PetscInt, PetscInt, PetscInt,
95                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
96                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
97                                                          PetscReal, const PetscReal[], PetscScalar[]));
98 PETSC_EXTERN PetscErrorCode PetscDSSetResidual(PetscDS, PetscInt,
99                                                void (*)(PetscInt, PetscInt, PetscInt,
100                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
101                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
102                                                         PetscReal, const PetscReal[], PetscScalar[]),
103                                                void (*)(PetscInt, PetscInt, PetscInt,
104                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
105                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
106                                                         PetscReal, const PetscReal[], PetscScalar[]));
107 PETSC_EXTERN PetscErrorCode PetscDSHasJacobian(PetscDS, PetscBool *);
108 PETSC_EXTERN PetscErrorCode PetscDSGetJacobian(PetscDS, PetscInt, PetscInt,
109                                                void (**)(PetscInt, 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                                                          PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
113                                                void (**)(PetscInt, PetscInt, PetscInt,
114                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
115                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
116                                                          PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
117                                                void (**)(PetscInt, PetscInt, PetscInt,
118                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
119                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
120                                                          PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
121                                                void (**)(PetscInt, PetscInt, PetscInt,
122                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
123                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
124                                                          PetscReal, PetscReal, const PetscReal[], PetscScalar[]));
125 PETSC_EXTERN PetscErrorCode PetscDSSetJacobian(PetscDS, PetscInt, PetscInt,
126                                                void (*)(PetscInt, 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                                                         PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
130                                                void (*)(PetscInt, PetscInt, PetscInt,
131                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
132                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
133                                                         PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
134                                                void (*)(PetscInt, PetscInt, PetscInt,
135                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
136                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
137                                                         PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
138                                                void (*)(PetscInt, PetscInt, PetscInt,
139                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
140                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
141                                                         PetscReal, PetscReal, const PetscReal[], PetscScalar[]));
142 PETSC_EXTERN PetscErrorCode PetscDSHasJacobianPreconditioner(PetscDS, PetscBool *);
143 PETSC_EXTERN PetscErrorCode PetscDSGetJacobianPreconditioner(PetscDS, PetscInt, PetscInt,
144                                                void (**)(PetscInt, PetscInt, PetscInt,
145                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
146                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
147                                                          PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
148                                                void (**)(PetscInt, PetscInt, PetscInt,
149                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
150                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
151                                                          PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
152                                                void (**)(PetscInt, PetscInt, PetscInt,
153                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
154                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
155                                                          PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
156                                                void (**)(PetscInt, PetscInt, PetscInt,
157                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
158                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
159                                                          PetscReal, PetscReal, const PetscReal[], PetscScalar[]));
160 PETSC_EXTERN PetscErrorCode PetscDSSetJacobianPreconditioner(PetscDS, PetscInt, PetscInt,
161                                                void (*)(PetscInt, PetscInt, PetscInt,
162                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
163                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
164                                                         PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
165                                                void (*)(PetscInt, PetscInt, PetscInt,
166                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
167                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
168                                                         PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
169                                                void (*)(PetscInt, PetscInt, PetscInt,
170                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
171                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
172                                                         PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
173                                                void (*)(PetscInt, PetscInt, PetscInt,
174                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
175                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
176                                                         PetscReal, PetscReal, const PetscReal[], PetscScalar[]));
177 PETSC_EXTERN PetscErrorCode PetscDSHasDynamicJacobian(PetscDS, PetscBool *);
178 PETSC_EXTERN PetscErrorCode PetscDSGetDynamicJacobian(PetscDS, PetscInt, PetscInt,
179                                                       void (**)(PetscInt, PetscInt, PetscInt,
180                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
181                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
182                                                                 PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
183                                                       void (**)(PetscInt, PetscInt, PetscInt,
184                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
185                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
186                                                                 PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
187                                                       void (**)(PetscInt, PetscInt, PetscInt,
188                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
189                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
190                                                                 PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
191                                                       void (**)(PetscInt, PetscInt, PetscInt,
192                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
193                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
194                                                                 PetscReal, PetscReal, const PetscReal[], PetscScalar[]));
195 PETSC_EXTERN PetscErrorCode PetscDSSetDynamicJacobian(PetscDS, PetscInt, PetscInt,
196                                                       void (*)(PetscInt, PetscInt, PetscInt,
197                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
198                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
199                                                                PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
200                                                       void (*)(PetscInt, PetscInt, PetscInt,
201                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
202                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
203                                                                PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
204                                                       void (*)(PetscInt, PetscInt, PetscInt,
205                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
206                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
207                                                                PetscReal, PetscReal, const PetscReal[], PetscScalar[]),
208                                                       void (*)(PetscInt, PetscInt, PetscInt,
209                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
210                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
211                                                                PetscReal, PetscReal, const PetscReal[], PetscScalar[]));
212 PETSC_EXTERN PetscErrorCode PetscDSGetRiemannSolver(PetscDS, PetscInt,
213                                                     void (**)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscScalar[], void *));
214 PETSC_EXTERN PetscErrorCode PetscDSSetRiemannSolver(PetscDS, PetscInt,
215                                                     void (*)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscScalar[], void *));
216 PETSC_EXTERN PetscErrorCode PetscDSGetContext(PetscDS, PetscInt, void **);
217 PETSC_EXTERN PetscErrorCode PetscDSSetContext(PetscDS, PetscInt, void *);
218 PETSC_EXTERN PetscErrorCode PetscDSGetBdResidual(PetscDS, PetscInt,
219                                                  void (**)(PetscInt, PetscInt, PetscInt,
220                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
221                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
222                                                            PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]),
223                                                  void (**)(PetscInt, PetscInt, PetscInt,
224                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
225                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
226                                                            PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]));
227 PETSC_EXTERN PetscErrorCode PetscDSSetBdResidual(PetscDS, PetscInt,
228                                                  void (*)(PetscInt, PetscInt, PetscInt,
229                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
230                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
231                                                           PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]),
232                                                  void (*)(PetscInt, PetscInt, PetscInt,
233                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
234                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
235                                                           PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]));
236 PETSC_EXTERN PetscErrorCode PetscDSGetBdJacobian(PetscDS, PetscInt, PetscInt,
237                                                  void (**)(PetscInt, PetscInt, PetscInt,
238                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
239                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
240                                                            PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]),
241                                                  void (**)(PetscInt, PetscInt, PetscInt,
242                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
243                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
244                                                            PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]),
245                                                  void (**)(PetscInt, PetscInt, PetscInt,
246                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
247                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
248                                                            PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]),
249                                                  void (**)(PetscInt, PetscInt, PetscInt,
250                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
251                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
252                                                            PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]));
253 PETSC_EXTERN PetscErrorCode PetscDSSetBdJacobian(PetscDS, PetscInt, PetscInt,
254                                                  void (*)(PetscInt, PetscInt, PetscInt,
255                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
256                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
257                                                           PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]),
258                                                  void (*)(PetscInt, PetscInt, PetscInt,
259                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
260                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
261                                                           PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]),
262                                                  void (*)(PetscInt, PetscInt, PetscInt,
263                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
264                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
265                                                           PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]),
266                                                  void (*)(PetscInt, PetscInt, PetscInt,
267                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
268                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
269                                                           PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]));
270 PETSC_EXTERN PetscErrorCode PetscDSGetTabulation(PetscDS, PetscReal ***, PetscReal ***);
271 PETSC_EXTERN PetscErrorCode PetscDSGetBdTabulation(PetscDS, PetscReal ***, PetscReal ***);
272 PETSC_EXTERN PetscErrorCode PetscDSGetEvaluationArrays(PetscDS, PetscScalar **, PetscScalar **, PetscScalar **);
273 PETSC_EXTERN PetscErrorCode PetscDSGetWeakFormArrays(PetscDS, PetscScalar **, PetscScalar **, PetscScalar **, PetscScalar **, PetscScalar **, PetscScalar **);
274 PETSC_EXTERN PetscErrorCode PetscDSGetRefCoordArrays(PetscDS, PetscReal **, PetscScalar **);
275 PETSC_EXTERN PetscErrorCode PetscDSCopyEquations(PetscDS, PetscDS);
276 
277 #endif
278