xref: /petsc/include/petscds.h (revision 4d0b9603208f03c394dfdee779423082bcb528af)
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 PetscDSGetFieldIndex(PetscDS, PetscObject, PetscInt *);
62 PETSC_EXTERN PetscErrorCode PetscDSGetFieldSize(PetscDS, PetscInt, PetscInt *);
63 PETSC_EXTERN PetscErrorCode PetscDSGetFieldOffset(PetscDS, PetscInt, PetscInt *);
64 PETSC_EXTERN PetscErrorCode PetscDSGetBdFieldOffset(PetscDS, PetscInt, PetscInt *);
65 PETSC_EXTERN PetscErrorCode PetscDSGetDimensions(PetscDS, PetscInt *[]);
66 PETSC_EXTERN PetscErrorCode PetscDSGetComponents(PetscDS, PetscInt *[]);
67 PETSC_EXTERN PetscErrorCode PetscDSGetComponentOffset(PetscDS, PetscInt, PetscInt *);
68 PETSC_EXTERN PetscErrorCode PetscDSGetComponentOffsets(PetscDS, PetscInt *[]);
69 PETSC_EXTERN PetscErrorCode PetscDSGetComponentBdOffsets(PetscDS, PetscInt *[]);
70 PETSC_EXTERN PetscErrorCode PetscDSGetComponentDerivativeOffsets(PetscDS, PetscInt *[]);
71 PETSC_EXTERN PetscErrorCode PetscDSGetComponentBdDerivativeOffsets(PetscDS, PetscInt *[]);
72 
73 PETSC_EXTERN PetscErrorCode PetscDSGetDiscretization(PetscDS, PetscInt, PetscObject *);
74 PETSC_EXTERN PetscErrorCode PetscDSSetDiscretization(PetscDS, PetscInt, PetscObject);
75 PETSC_EXTERN PetscErrorCode PetscDSAddDiscretization(PetscDS, PetscObject);
76 PETSC_EXTERN PetscErrorCode PetscDSGetBdDiscretization(PetscDS, PetscInt, PetscObject *);
77 PETSC_EXTERN PetscErrorCode PetscDSSetBdDiscretization(PetscDS, PetscInt, PetscObject);
78 PETSC_EXTERN PetscErrorCode PetscDSAddBdDiscretization(PetscDS, PetscObject);
79 PETSC_EXTERN PetscErrorCode PetscDSGetImplicit(PetscDS, PetscInt, PetscBool*);
80 PETSC_EXTERN PetscErrorCode PetscDSSetImplicit(PetscDS, PetscInt, PetscBool);
81 PETSC_EXTERN PetscErrorCode PetscDSGetAdjacency(PetscDS, PetscInt, PetscBool*, PetscBool*);
82 PETSC_EXTERN PetscErrorCode PetscDSSetAdjacency(PetscDS, PetscInt, PetscBool,  PetscBool);
83 PETSC_EXTERN PetscErrorCode PetscDSGetObjective(PetscDS, PetscInt,
84                                                 void (**)(PetscInt, PetscInt, PetscInt,
85                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
86                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
87                                                           PetscReal, const PetscReal[], PetscScalar[]));
88 PETSC_EXTERN PetscErrorCode PetscDSSetObjective(PetscDS, PetscInt,
89                                                 void (*)(PetscInt, PetscInt, PetscInt,
90                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
91                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
92                                                          PetscReal, const PetscReal[], PetscScalar[]));
93 PETSC_EXTERN PetscErrorCode PetscDSGetResidual(PetscDS, PetscInt,
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                                                void (**)(PetscInt, PetscInt, PetscInt,
99                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
100                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
101                                                          PetscReal, const PetscReal[], PetscScalar[]));
102 PETSC_EXTERN PetscErrorCode PetscDSSetResidual(PetscDS, PetscInt,
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                                                void (*)(PetscInt, PetscInt, PetscInt,
108                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
109                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
110                                                         PetscReal, const PetscReal[], PetscScalar[]));
111 PETSC_EXTERN PetscErrorCode PetscDSHasJacobian(PetscDS, PetscBool *);
112 PETSC_EXTERN PetscErrorCode PetscDSGetJacobian(PetscDS, PetscInt, PetscInt,
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                                                void (**)(PetscInt, PetscInt, PetscInt,
126                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
127                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
128                                                          PetscReal, PetscReal, const PetscReal[], PetscScalar[]));
129 PETSC_EXTERN PetscErrorCode PetscDSSetJacobian(PetscDS, PetscInt, PetscInt,
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                                                void (*)(PetscInt, PetscInt, PetscInt,
143                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
144                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
145                                                         PetscReal, PetscReal, const PetscReal[], PetscScalar[]));
146 PETSC_EXTERN PetscErrorCode PetscDSHasJacobianPreconditioner(PetscDS, PetscBool *);
147 PETSC_EXTERN PetscErrorCode PetscDSGetJacobianPreconditioner(PetscDS, PetscInt, PetscInt,
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                                                void (**)(PetscInt, PetscInt, PetscInt,
161                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
162                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
163                                                          PetscReal, PetscReal, const PetscReal[], PetscScalar[]));
164 PETSC_EXTERN PetscErrorCode PetscDSSetJacobianPreconditioner(PetscDS, PetscInt, PetscInt,
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                                                void (*)(PetscInt, PetscInt, PetscInt,
178                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
179                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
180                                                         PetscReal, PetscReal, const PetscReal[], PetscScalar[]));
181 PETSC_EXTERN PetscErrorCode PetscDSHasDynamicJacobian(PetscDS, PetscBool *);
182 PETSC_EXTERN PetscErrorCode PetscDSGetDynamicJacobian(PetscDS, PetscInt, PetscInt,
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                                                       void (**)(PetscInt, PetscInt, PetscInt,
196                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
197                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
198                                                                 PetscReal, PetscReal, const PetscReal[], PetscScalar[]));
199 PETSC_EXTERN PetscErrorCode PetscDSSetDynamicJacobian(PetscDS, PetscInt, PetscInt,
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                                                       void (*)(PetscInt, PetscInt, PetscInt,
213                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
214                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
215                                                                PetscReal, PetscReal, const PetscReal[], PetscScalar[]));
216 PETSC_EXTERN PetscErrorCode PetscDSGetRiemannSolver(PetscDS, PetscInt,
217                                                     void (**)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscScalar[], void *));
218 PETSC_EXTERN PetscErrorCode PetscDSSetRiemannSolver(PetscDS, PetscInt,
219                                                     void (*)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscScalar[], void *));
220 PETSC_EXTERN PetscErrorCode PetscDSGetContext(PetscDS, PetscInt, void **);
221 PETSC_EXTERN PetscErrorCode PetscDSSetContext(PetscDS, PetscInt, void *);
222 PETSC_EXTERN PetscErrorCode PetscDSGetBdResidual(PetscDS, PetscInt,
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                                                  void (**)(PetscInt, PetscInt, PetscInt,
228                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
229                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
230                                                            PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]));
231 PETSC_EXTERN PetscErrorCode PetscDSSetBdResidual(PetscDS, PetscInt,
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                                                  void (*)(PetscInt, PetscInt, PetscInt,
237                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
238                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
239                                                           PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]));
240 PETSC_EXTERN PetscErrorCode PetscDSGetBdJacobian(PetscDS, PetscInt, PetscInt,
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                                                  void (**)(PetscInt, PetscInt, PetscInt,
254                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
255                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
256                                                            PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]));
257 PETSC_EXTERN PetscErrorCode PetscDSSetBdJacobian(PetscDS, PetscInt, PetscInt,
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                                                  void (*)(PetscInt, PetscInt, PetscInt,
271                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
272                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
273                                                           PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]));
274 PETSC_EXTERN PetscErrorCode PetscDSGetTabulation(PetscDS, PetscReal ***, PetscReal ***);
275 PETSC_EXTERN PetscErrorCode PetscDSGetFaceTabulation(PetscDS, PetscReal ***, PetscReal ***);
276 PETSC_EXTERN PetscErrorCode PetscDSGetBdTabulation(PetscDS, PetscReal ***, PetscReal ***);
277 PETSC_EXTERN PetscErrorCode PetscDSGetEvaluationArrays(PetscDS, PetscScalar **, PetscScalar **, PetscScalar **);
278 PETSC_EXTERN PetscErrorCode PetscDSGetWeakFormArrays(PetscDS, PetscScalar **, PetscScalar **, PetscScalar **, PetscScalar **, PetscScalar **, PetscScalar **);
279 PETSC_EXTERN PetscErrorCode PetscDSGetRefCoordArrays(PetscDS, PetscReal **, PetscScalar **);
280 PETSC_EXTERN PetscErrorCode PetscDSCopyEquations(PetscDS, PetscDS);
281 PETSC_EXTERN PetscErrorCode PetscDSAddBoundary(PetscDS, PetscBool, const char[], const char[], PetscInt, PetscInt, const PetscInt *, void (*)(), PetscInt, const PetscInt *, void *);
282 PETSC_EXTERN PetscErrorCode PetscDSGetNumBoundary(PetscDS, PetscInt *);
283 PETSC_EXTERN PetscErrorCode PetscDSGetBoundary(PetscDS, PetscInt, PetscBool *, const char **, const char **, PetscInt *, PetscInt *, const PetscInt **, void (**)(), PetscInt *, const PetscInt **, void **);
284 PETSC_EXTERN PetscErrorCode PetscDSCopyBoundary(PetscDS, PetscDS);
285 
286 #endif
287