xref: /petsc/include/petscds.h (revision ac09b9214d23ea9ad238aa607de9fa447fd4e91b)
1 /*
2       Objects which encapsulate discretizations+continuum residuals
3 */
4 #if !defined(PETSCDS_H)
5 #define PETSCDS_H
6 
7 #include <petscfe.h>
8 #include <petscfv.h>
9 #include <petscdstypes.h>
10 
11 /* SUBMANSEC = DT */
12 
13 PETSC_EXTERN PetscClassId PETSCWEAKFORM_CLASSID;
14 
15 PETSC_EXTERN PetscErrorCode PetscWeakFormCreate(MPI_Comm, PetscWeakForm *);
16 PETSC_EXTERN PetscErrorCode PetscWeakFormDestroy(PetscWeakForm *);
17 PETSC_EXTERN PetscErrorCode PetscWeakFormView(PetscWeakForm, PetscViewer);
18 PETSC_EXTERN PetscErrorCode PetscWeakFormCopy(PetscWeakForm, PetscWeakForm);
19 PETSC_EXTERN PetscErrorCode PetscWeakFormClear(PetscWeakForm);
20 PETSC_EXTERN PetscErrorCode PetscWeakFormGetNumFields(PetscWeakForm, PetscInt *);
21 PETSC_EXTERN PetscErrorCode PetscWeakFormSetNumFields(PetscWeakForm, PetscInt);
22 PETSC_EXTERN PetscErrorCode PetscFormKeySort(PetscInt, PetscFormKey[]);
23 PETSC_EXTERN PetscErrorCode PetscWeakFormRewriteKeys(PetscWeakForm, DMLabel, PetscInt, const PetscInt[]);
24 PETSC_EXTERN PetscErrorCode PetscWeakFormReplaceLabel(PetscWeakForm, DMLabel);
25 PETSC_EXTERN PetscErrorCode PetscWeakFormClearIndex(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscWeakFormKind, PetscInt);
26 
27 PETSC_EXTERN PetscErrorCode PetscWeakFormGetObjective(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt *,
28                                                       void (***)(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, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
32 PETSC_EXTERN PetscErrorCode PetscWeakFormAddObjective(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt,
33                                                       void (*)(PetscInt, PetscInt, PetscInt,
34                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
35                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
36                                                                PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
37 PETSC_EXTERN PetscErrorCode PetscWeakFormSetObjective(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
38                                                       void (**)(PetscInt, PetscInt, PetscInt,
39                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
40                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
41                                                                PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
42 PETSC_EXTERN PetscErrorCode PetscWeakFormGetIndexObjective(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
43                                                       void (**)(PetscInt, PetscInt, PetscInt,
44                                                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
45                                                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
46                                                                  PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
47 PETSC_EXTERN PetscErrorCode PetscWeakFormSetIndexObjective(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
48                                                       void (*)(PetscInt, PetscInt, PetscInt,
49                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
50                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
51                                                                PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
52 PETSC_EXTERN PetscErrorCode PetscWeakFormGetResidual(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt,
53                                                      PetscInt *,
54                                                      void (***)(PetscInt, PetscInt, PetscInt,
55                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
56                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
57                                                                 PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
58                                                      PetscInt *,
59                                                      void (***)(PetscInt, PetscInt, PetscInt,
60                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
61                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
62                                                                 PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
63 PETSC_EXTERN PetscErrorCode PetscWeakFormAddResidual(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt,
64                                                      void (*)(PetscInt, PetscInt, PetscInt,
65                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
66                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
67                                                               PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
68                                                       void (*)(PetscInt, PetscInt, PetscInt,
69                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
70                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
71                                                                PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
72 PETSC_EXTERN PetscErrorCode PetscWeakFormSetResidual(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt,
73                                                      PetscInt,
74                                                      void (**)(PetscInt, PetscInt, PetscInt,
75                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
76                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
77                                                               PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
78                                                      PetscInt,
79                                                      void (**)(PetscInt, PetscInt, PetscInt,
80                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
81                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
82                                                               PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
83 PETSC_EXTERN PetscErrorCode PetscWeakFormSetIndexResidual(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt,
84                                                      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[], PetscInt, const PetscScalar[], PetscScalar[]),
89                                                      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[], PetscInt, const PetscScalar[], PetscScalar[]));
94 PETSC_EXTERN PetscErrorCode PetscWeakFormHasJacobian(PetscWeakForm, PetscBool *);
95 PETSC_EXTERN PetscErrorCode PetscWeakFormGetJacobian(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
96                                                      PetscInt *,
97                                                      void (***)(PetscInt, PetscInt, PetscInt,
98                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
99                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
100                                                                 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
101                                                      PetscInt *,
102                                                      void (***)(PetscInt, PetscInt, PetscInt,
103                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
104                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
105                                                                 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
106                                                      PetscInt *,
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, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
111                                                      PetscInt *,
112                                                      void (***)(PetscInt, PetscInt, PetscInt,
113                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
114                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
115                                                                 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
116 PETSC_EXTERN PetscErrorCode PetscWeakFormAddJacobian(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
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[], PetscInt, const PetscScalar[], 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[], PetscInt, const PetscScalar[], 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[], PetscInt, const PetscScalar[], PetscScalar[]),
129                                                       void (*)(PetscInt, PetscInt, PetscInt,
130                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
131                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
132                                                                PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
133 PETSC_EXTERN PetscErrorCode PetscWeakFormSetJacobian(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
134                                                      PetscInt,
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[], PetscInt, const PetscScalar[], PetscScalar[]),
139                                                      PetscInt,
140                                                      void (**)(PetscInt, PetscInt, PetscInt,
141                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
142                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
143                                                               PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
144                                                      PetscInt,
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[], PetscInt, const PetscScalar[], PetscScalar[]),
149                                                      PetscInt,
150                                                      void (**)(PetscInt, PetscInt, PetscInt,
151                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
152                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
153                                                               PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
154 PETSC_EXTERN PetscErrorCode PetscWeakFormSetIndexJacobian(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
155                                                      PetscInt,
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[], PetscInt, const PetscScalar[], PetscScalar[]),
160                                                      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[], PetscInt, const PetscScalar[], PetscScalar[]),
165                                                      PetscInt,
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[], PetscInt, const PetscScalar[], PetscScalar[]),
170                                                      PetscInt,
171                                                      void (*)(PetscInt, PetscInt, PetscInt,
172                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
173                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
174                                                               PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
175 PETSC_EXTERN PetscErrorCode PetscWeakFormHasJacobianPreconditioner(PetscWeakForm, PetscBool *);
176 PETSC_EXTERN PetscErrorCode PetscWeakFormGetJacobianPreconditioner(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
177                                                      PetscInt *,
178                                                      void (***)(PetscInt, PetscInt, PetscInt,
179                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
180                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
181                                                                 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
182                                                      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[], PetscInt, const PetscScalar[], PetscScalar[]),
187                                                      PetscInt *,
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[], PetscInt, const PetscScalar[], PetscScalar[]),
192                                                      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[], PetscInt, const PetscScalar[], PetscScalar[]));
197 PETSC_EXTERN PetscErrorCode PetscWeakFormAddJacobianPreconditioner(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
198                                                      void (*)(PetscInt, PetscInt, PetscInt,
199                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
200                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
201                                                               PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
202                                                       void (*)(PetscInt, PetscInt, PetscInt,
203                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
204                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
205                                                                PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
206                                                       void (*)(PetscInt, PetscInt, PetscInt,
207                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
208                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
209                                                                PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
210                                                       void (*)(PetscInt, PetscInt, PetscInt,
211                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
212                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
213                                                                PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
214 PETSC_EXTERN PetscErrorCode PetscWeakFormSetJacobianPreconditioner(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
215                                                      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, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
220                                                      PetscInt,
221                                                      void (**)(PetscInt, PetscInt, PetscInt,
222                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
223                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
224                                                               PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
225                                                      PetscInt,
226                                                      void (**)(PetscInt, PetscInt, PetscInt,
227                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
228                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
229                                                               PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
230                                                      PetscInt,
231                                                      void (**)(PetscInt, PetscInt, PetscInt,
232                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
233                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
234                                                               PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
235 PETSC_EXTERN PetscErrorCode PetscWeakFormSetIndexJacobianPreconditioner(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
236                                                      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[], PetscInt, const PetscScalar[], PetscScalar[]),
241                                                      PetscInt,
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[], PetscInt, const PetscScalar[], PetscScalar[]),
246                                                      PetscInt,
247                                                      void (*)(PetscInt, PetscInt, PetscInt,
248                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
249                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
250                                                               PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
251                                                      PetscInt,
252                                                      void (*)(PetscInt, PetscInt, PetscInt,
253                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
254                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
255                                                               PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
256 PETSC_EXTERN PetscErrorCode PetscWeakFormHasDynamicJacobian(PetscWeakForm, PetscBool *);
257 PETSC_EXTERN PetscErrorCode PetscWeakFormGetDynamicJacobian(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
258                                                      PetscInt *,
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[], PetscInt, const PetscScalar[], PetscScalar[]),
263                                                      PetscInt *,
264                                                      void (***)(PetscInt, PetscInt, PetscInt,
265                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
266                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
267                                                                 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
268                                                      PetscInt *,
269                                                      void (***)(PetscInt, PetscInt, PetscInt,
270                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
271                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
272                                                                 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
273                                                      PetscInt *,
274                                                      void (***)(PetscInt, PetscInt, PetscInt,
275                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
276                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
277                                                                 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
278 PETSC_EXTERN PetscErrorCode PetscWeakFormAddDynamicJacobian(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
279                                                      void (*)(PetscInt, PetscInt, PetscInt,
280                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
281                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
282                                                               PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
283                                                       void (*)(PetscInt, PetscInt, PetscInt,
284                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
285                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
286                                                                PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
287                                                       void (*)(PetscInt, PetscInt, PetscInt,
288                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
289                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
290                                                                PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
291                                                       void (*)(PetscInt, PetscInt, PetscInt,
292                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
293                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
294                                                                PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
295 PETSC_EXTERN PetscErrorCode PetscWeakFormSetDynamicJacobian(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
296                                                      PetscInt,
297                                                      void (**)(PetscInt, PetscInt, PetscInt,
298                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
299                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
300                                                               PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
301                                                      PetscInt,
302                                                      void (**)(PetscInt, PetscInt, PetscInt,
303                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
304                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
305                                                               PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
306                                                      PetscInt,
307                                                      void (**)(PetscInt, PetscInt, PetscInt,
308                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
309                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
310                                                               PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
311                                                      PetscInt,
312                                                      void (**)(PetscInt, PetscInt, PetscInt,
313                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
314                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
315                                                               PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
316 PETSC_EXTERN PetscErrorCode PetscWeakFormSetIndexDynamicJacobian(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
317                                                      PetscInt,
318                                                      void (*)(PetscInt, PetscInt, PetscInt,
319                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
320                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
321                                                               PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
322                                                      PetscInt,
323                                                      void (*)(PetscInt, PetscInt, PetscInt,
324                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
325                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
326                                                               PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
327                                                      PetscInt,
328                                                      void (*)(PetscInt, PetscInt, PetscInt,
329                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
330                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
331                                                               PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
332                                                      PetscInt,
333                                                      void (*)(PetscInt, PetscInt, PetscInt,
334                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
335                                                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
336                                                               PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
337 PETSC_EXTERN PetscErrorCode PetscWeakFormGetBdResidual(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt,
338                                                        PetscInt *,
339                                                        void (***)(PetscInt, PetscInt, PetscInt,
340                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
341                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
342                                                                   PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
343                                                        PetscInt *,
344                                                        void (***)(PetscInt, PetscInt, PetscInt,
345                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
346                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
347                                                                   PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
348 PETSC_EXTERN PetscErrorCode PetscWeakFormAddBdResidual(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt,
349                                                        void (*)(PetscInt, PetscInt, PetscInt,
350                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
351                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
352                                                                 PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
353                                                        void (*)(PetscInt, PetscInt, PetscInt,
354                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
355                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
356                                                                 PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
357 PETSC_EXTERN PetscErrorCode PetscWeakFormSetBdResidual(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt,
358                                                        PetscInt,
359                                                        void (**)(PetscInt, PetscInt, PetscInt,
360                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
361                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
362                                                                 PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
363                                                        PetscInt,
364                                                        void (**)(PetscInt, PetscInt, PetscInt,
365                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
366                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
367                                                                 PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
368 PETSC_EXTERN PetscErrorCode PetscWeakFormSetIndexBdResidual(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt,
369                                                        PetscInt,
370                                                        void (*)(PetscInt, PetscInt, PetscInt,
371                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
372                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
373                                                                 PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
374                                                        PetscInt,
375                                                        void (*)(PetscInt, PetscInt, PetscInt,
376                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
377                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
378                                                                 PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
379 PETSC_EXTERN PetscErrorCode PetscWeakFormHasBdJacobian(PetscWeakForm, PetscBool *);
380 PETSC_EXTERN PetscErrorCode PetscWeakFormGetBdJacobian(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
381                                                        PetscInt *,
382                                                        void (***)(PetscInt, PetscInt, PetscInt,
383                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
384                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
385                                                                   PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
386                                                        PetscInt *,
387                                                        void (***)(PetscInt, PetscInt, PetscInt,
388                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
389                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
390                                                                   PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
391                                                        PetscInt *,
392                                                        void (***)(PetscInt, PetscInt, PetscInt,
393                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
394                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
395                                                                   PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
396                                                        PetscInt *,
397                                                        void (***)(PetscInt, PetscInt, PetscInt,
398                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
399                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
400                                                                   PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
401 PETSC_EXTERN PetscErrorCode PetscWeakFormAddBdJacobian(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
402                                                        void (*)(PetscInt, PetscInt, PetscInt,
403                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
404                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
405                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
406                                                        void (*)(PetscInt, PetscInt, PetscInt,
407                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
408                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
409                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
410                                                        void (*)(PetscInt, PetscInt, PetscInt,
411                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
412                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
413                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
414                                                        void (*)(PetscInt, PetscInt, PetscInt,
415                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
416                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
417                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
418 PETSC_EXTERN PetscErrorCode PetscWeakFormSetBdJacobian(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
419                                                        PetscInt,
420                                                        void (**)(PetscInt, PetscInt, PetscInt,
421                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
422                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
423                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
424                                                        PetscInt,
425                                                        void (**)(PetscInt, PetscInt, PetscInt,
426                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
427                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
428                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
429                                                        PetscInt,
430                                                        void (**)(PetscInt, PetscInt, PetscInt,
431                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
432                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
433                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
434                                                        PetscInt,
435                                                        void (**)(PetscInt, PetscInt, PetscInt,
436                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
437                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
438                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
439 PETSC_EXTERN PetscErrorCode PetscWeakFormSetIndexBdJacobian(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
440                                                        PetscInt,
441                                                        void (*)(PetscInt, PetscInt, PetscInt,
442                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
443                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
444                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
445                                                        PetscInt,
446                                                        void (*)(PetscInt, PetscInt, PetscInt,
447                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
448                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
449                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
450                                                        PetscInt,
451                                                        void (*)(PetscInt, PetscInt, PetscInt,
452                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
453                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
454                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
455                                                        PetscInt,
456                                                        void (*)(PetscInt, PetscInt, PetscInt,
457                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
458                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
459                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
460 PETSC_EXTERN PetscErrorCode PetscWeakFormHasBdJacobianPreconditioner(PetscWeakForm, PetscBool *);
461 PETSC_EXTERN PetscErrorCode PetscWeakFormGetBdJacobianPreconditioner(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
462                                                        PetscInt *,
463                                                        void (***)(PetscInt, PetscInt, PetscInt,
464                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
465                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
466                                                                   PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
467                                                        PetscInt *,
468                                                        void (***)(PetscInt, PetscInt, PetscInt,
469                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
470                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
471                                                                   PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
472                                                        PetscInt *,
473                                                        void (***)(PetscInt, PetscInt, PetscInt,
474                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
475                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
476                                                                   PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
477                                                        PetscInt *,
478                                                        void (***)(PetscInt, PetscInt, PetscInt,
479                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
480                                                                   const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
481                                                                   PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
482 PETSC_EXTERN PetscErrorCode PetscWeakFormAddBdJacobianPreconditioner(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
483                                                        void (*)(PetscInt, PetscInt, PetscInt,
484                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
485                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
486                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
487                                                        void (*)(PetscInt, PetscInt, PetscInt,
488                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
489                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
490                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
491                                                        void (*)(PetscInt, PetscInt, PetscInt,
492                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
493                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
494                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
495                                                        void (*)(PetscInt, PetscInt, PetscInt,
496                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
497                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
498                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
499 PETSC_EXTERN PetscErrorCode PetscWeakFormSetBdJacobianPreconditioner(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
500                                                        PetscInt,
501                                                        void (**)(PetscInt, PetscInt, PetscInt,
502                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
503                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
504                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
505                                                        PetscInt,
506                                                        void (**)(PetscInt, PetscInt, PetscInt,
507                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
508                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
509                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
510                                                        PetscInt,
511                                                        void (**)(PetscInt, PetscInt, PetscInt,
512                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
513                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
514                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
515                                                        PetscInt,
516                                                        void (**)(PetscInt, PetscInt, PetscInt,
517                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
518                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
519                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
520 PETSC_EXTERN PetscErrorCode PetscWeakFormSetIndexBdJacobianPreconditioner(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt, PetscInt,
521                                                        PetscInt,
522                                                        void (*)(PetscInt, PetscInt, PetscInt,
523                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
524                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
525                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
526                                                        PetscInt,
527                                                        void (*)(PetscInt, PetscInt, PetscInt,
528                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
529                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
530                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
531                                                        PetscInt,
532                                                        void (*)(PetscInt, PetscInt, PetscInt,
533                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
534                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
535                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
536                                                        PetscInt,
537                                                        void (*)(PetscInt, PetscInt, PetscInt,
538                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
539                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
540                                                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
541 PETSC_EXTERN PetscErrorCode PetscWeakFormGetRiemannSolver(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt,
542                                                           PetscInt *,
543                                                           void (***)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *));
544 PETSC_EXTERN PetscErrorCode PetscWeakFormSetRiemannSolver(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt,
545                                                           PetscInt,
546                                                           void (**)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *));
547 PETSC_EXTERN PetscErrorCode PetscWeakFormSetIndexRiemannSolver(PetscWeakForm, DMLabel, PetscInt, PetscInt, PetscInt,
548                                                           PetscInt,
549                                                           void (*)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *));
550 
551 PETSC_EXTERN PetscErrorCode PetscDSInitializePackage(void);
552 
553 PETSC_EXTERN PetscClassId PETSCDS_CLASSID;
554 
555 /*J
556   PetscDSType - String with the name of a PETSc discrete system
557 
558   Level: beginner
559 
560 .seealso: `PetscDSSetType()`, `PetscDS`
561 J*/
562 typedef const char *PetscDSType;
563 #define PETSCDSBASIC "basic"
564 
565 typedef enum {PETSC_DISC_NONE, PETSC_DISC_FE, PETSC_DISC_FV} PetscDiscType;
566 
567 typedef void (*PetscPointFunc)(PetscInt, PetscInt, PetscInt,
568                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
569                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
570                                PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]);
571 typedef void (*PetscPointJac)(PetscInt, PetscInt, PetscInt,
572                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
573                               const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
574                               PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]);
575 typedef void (*PetscBdPointFunc)(PetscInt, PetscInt, PetscInt,
576                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
577                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
578                                  PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]);
579 typedef void (*PetscBdPointJac)(PetscInt, PetscInt, PetscInt,
580                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
581                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
582                                 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]);
583 typedef void (*PetscRiemannFunc)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *);
584 typedef PetscErrorCode (*PetscSimplePointFunc)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *);
585 
586 PETSC_EXTERN PetscFunctionList PetscDSList;
587 PETSC_EXTERN PetscErrorCode PetscDSCreate(MPI_Comm, PetscDS *);
588 PETSC_EXTERN PetscErrorCode PetscDSDestroy(PetscDS *);
589 PETSC_EXTERN PetscErrorCode PetscDSSetType(PetscDS, PetscDSType);
590 PETSC_EXTERN PetscErrorCode PetscDSGetType(PetscDS, PetscDSType *);
591 PETSC_EXTERN PetscErrorCode PetscDSSetUp(PetscDS);
592 PETSC_EXTERN PetscErrorCode PetscDSSetFromOptions(PetscDS);
593 PETSC_EXTERN PetscErrorCode PetscDSViewFromOptions(PetscDS,PetscObject,const char[]);
594 
595 PETSC_EXTERN PetscErrorCode PetscDSView(PetscDS,PetscViewer);
596 PETSC_EXTERN PetscErrorCode PetscDSRegister(const char [], PetscErrorCode (*)(PetscDS));
597 PETSC_EXTERN PetscErrorCode PetscDSRegisterDestroy(void);
598 
599 PETSC_EXTERN PetscErrorCode PetscDSGetHeightSubspace(PetscDS, PetscInt, PetscDS *);
600 PETSC_EXTERN PetscErrorCode PetscDSGetSpatialDimension(PetscDS, PetscInt *);
601 PETSC_EXTERN PetscErrorCode PetscDSGetCoordinateDimension(PetscDS, PetscInt *);
602 PETSC_EXTERN PetscErrorCode PetscDSSetCoordinateDimension(PetscDS, PetscInt);
603 PETSC_EXTERN PetscErrorCode PetscDSIsCohesive(PetscDS, PetscBool *);
604 PETSC_EXTERN PetscErrorCode PetscDSGetNumCohesive(PetscDS, PetscInt *);
605 PETSC_EXTERN PetscErrorCode PetscDSGetCohesive(PetscDS, PetscInt, PetscBool *);
606 PETSC_EXTERN PetscErrorCode PetscDSSetCohesive(PetscDS, PetscInt, PetscBool);
607 PETSC_EXTERN PetscErrorCode PetscDSGetNumFields(PetscDS, PetscInt *);
608 PETSC_EXTERN PetscErrorCode PetscDSGetTotalDimension(PetscDS, PetscInt *);
609 PETSC_EXTERN PetscErrorCode PetscDSGetTotalComponents(PetscDS, PetscInt *);
610 PETSC_EXTERN PetscErrorCode PetscDSGetFieldIndex(PetscDS, PetscObject, PetscInt *);
611 PETSC_EXTERN PetscErrorCode PetscDSGetFieldSize(PetscDS, PetscInt, PetscInt *);
612 PETSC_EXTERN PetscErrorCode PetscDSGetFieldOffset(PetscDS, PetscInt, PetscInt *);
613 PETSC_EXTERN PetscErrorCode PetscDSGetFieldOffsetCohesive(PetscDS, PetscInt, PetscInt *);
614 PETSC_EXTERN PetscErrorCode PetscDSGetDimensions(PetscDS, PetscInt *[]);
615 PETSC_EXTERN PetscErrorCode PetscDSGetComponents(PetscDS, PetscInt *[]);
616 PETSC_EXTERN PetscErrorCode PetscDSGetComponentOffset(PetscDS, PetscInt, PetscInt *);
617 PETSC_EXTERN PetscErrorCode PetscDSGetComponentOffsets(PetscDS, PetscInt *[]);
618 PETSC_EXTERN PetscErrorCode PetscDSGetComponentDerivativeOffsets(PetscDS, PetscInt *[]);
619 PETSC_EXTERN PetscErrorCode PetscDSGetComponentOffsetsCohesive(PetscDS, PetscInt, PetscInt *[]);
620 PETSC_EXTERN PetscErrorCode PetscDSGetComponentDerivativeOffsetsCohesive(PetscDS, PetscInt, PetscInt *[]);
621 
622 PETSC_EXTERN PetscErrorCode PetscDSGetWeakForm(PetscDS, PetscWeakForm *);
623 PETSC_EXTERN PetscErrorCode PetscDSSetWeakForm(PetscDS, PetscWeakForm);
624 PETSC_EXTERN PetscErrorCode PetscDSGetDiscretization(PetscDS, PetscInt, PetscObject *);
625 PETSC_EXTERN PetscErrorCode PetscDSSetDiscretization(PetscDS, PetscInt, PetscObject);
626 PETSC_EXTERN PetscErrorCode PetscDSAddDiscretization(PetscDS, PetscObject);
627 PETSC_EXTERN PetscErrorCode PetscDSGetQuadrature(PetscDS, PetscQuadrature*);
628 PETSC_EXTERN PetscErrorCode PetscDSGetImplicit(PetscDS, PetscInt, PetscBool*);
629 PETSC_EXTERN PetscErrorCode PetscDSSetImplicit(PetscDS, PetscInt, PetscBool);
630 PETSC_EXTERN PetscErrorCode PetscDSGetJetDegree(PetscDS, PetscInt, PetscInt*);
631 PETSC_EXTERN PetscErrorCode PetscDSSetJetDegree(PetscDS, PetscInt, PetscInt);
632 PETSC_EXTERN PetscErrorCode PetscDSGetConstants(PetscDS, PetscInt *, const PetscScalar *[]);
633 PETSC_EXTERN PetscErrorCode PetscDSSetConstants(PetscDS, PetscInt, PetscScalar[]);
634 PETSC_EXTERN PetscErrorCode PetscDSGetObjective(PetscDS, PetscInt,
635                                                 void (**)(PetscInt, PetscInt, PetscInt,
636                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
637                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
638                                                           PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
639 PETSC_EXTERN PetscErrorCode PetscDSSetObjective(PetscDS, PetscInt,
640                                                 void (*)(PetscInt, PetscInt, PetscInt,
641                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
642                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
643                                                          PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
644 PETSC_EXTERN PetscErrorCode PetscDSGetResidual(PetscDS, PetscInt,
645                                                void (**)(PetscInt, PetscInt, PetscInt,
646                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
647                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
648                                                          PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
649                                                void (**)(PetscInt, PetscInt, PetscInt,
650                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
651                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
652                                                          PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
653 PETSC_EXTERN PetscErrorCode PetscDSSetResidual(PetscDS, PetscInt,
654                                                void (*)(PetscInt, PetscInt, PetscInt,
655                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
656                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
657                                                         PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
658                                                void (*)(PetscInt, PetscInt, PetscInt,
659                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
660                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
661                                                         PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
662 PETSC_EXTERN PetscErrorCode PetscDSGetRHSResidual(PetscDS, PetscInt,
663                                                   void (**)(PetscInt, PetscInt, PetscInt,
664                                                             const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
665                                                             const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
666                                                             PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
667                                                   void (**)(PetscInt, PetscInt, PetscInt,
668                                                             const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
669                                                             const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
670                                                             PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
671 PETSC_EXTERN PetscErrorCode PetscDSSetRHSResidual(PetscDS, PetscInt,
672                                                   void (*)(PetscInt, PetscInt, PetscInt,
673                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
674                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
675                                                            PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
676                                                   void (*)(PetscInt, PetscInt, PetscInt,
677                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
678                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
679                                                            PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
680 PETSC_EXTERN PetscErrorCode PetscDSHasJacobian(PetscDS, PetscBool *);
681 PETSC_EXTERN PetscErrorCode PetscDSGetJacobian(PetscDS, PetscInt, PetscInt,
682                                                void (**)(PetscInt, PetscInt, PetscInt,
683                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
684                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
685                                                          PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
686                                                void (**)(PetscInt, PetscInt, PetscInt,
687                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
688                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
689                                                          PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
690                                                void (**)(PetscInt, PetscInt, PetscInt,
691                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
692                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
693                                                          PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
694                                                void (**)(PetscInt, PetscInt, PetscInt,
695                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
696                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
697                                                          PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
698 PETSC_EXTERN PetscErrorCode PetscDSSetJacobian(PetscDS, PetscInt, PetscInt,
699                                                void (*)(PetscInt, PetscInt, PetscInt,
700                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
701                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
702                                                         PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
703                                                void (*)(PetscInt, PetscInt, PetscInt,
704                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
705                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
706                                                         PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
707                                                void (*)(PetscInt, PetscInt, PetscInt,
708                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
709                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
710                                                         PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
711                                                void (*)(PetscInt, PetscInt, PetscInt,
712                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
713                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
714                                                         PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
715 PETSC_EXTERN PetscErrorCode PetscDSUseJacobianPreconditioner(PetscDS, PetscBool);
716 PETSC_EXTERN PetscErrorCode PetscDSHasJacobianPreconditioner(PetscDS, PetscBool *);
717 PETSC_EXTERN PetscErrorCode PetscDSGetJacobianPreconditioner(PetscDS, PetscInt, PetscInt,
718                                                void (**)(PetscInt, PetscInt, PetscInt,
719                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
720                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
721                                                          PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
722                                                void (**)(PetscInt, PetscInt, PetscInt,
723                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
724                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
725                                                          PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
726                                                void (**)(PetscInt, PetscInt, PetscInt,
727                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
728                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
729                                                          PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
730                                                void (**)(PetscInt, PetscInt, PetscInt,
731                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
732                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
733                                                          PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
734 PETSC_EXTERN PetscErrorCode PetscDSSetJacobianPreconditioner(PetscDS, PetscInt, PetscInt,
735                                                void (*)(PetscInt, PetscInt, PetscInt,
736                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
737                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
738                                                         PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
739                                                void (*)(PetscInt, PetscInt, PetscInt,
740                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
741                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
742                                                         PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
743                                                void (*)(PetscInt, PetscInt, PetscInt,
744                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
745                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
746                                                         PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
747                                                void (*)(PetscInt, PetscInt, PetscInt,
748                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
749                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
750                                                         PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
751 PETSC_EXTERN PetscErrorCode PetscDSHasDynamicJacobian(PetscDS, PetscBool *);
752 PETSC_EXTERN PetscErrorCode PetscDSGetDynamicJacobian(PetscDS, PetscInt, PetscInt,
753                                                       void (**)(PetscInt, PetscInt, PetscInt,
754                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
755                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
756                                                                 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
757                                                       void (**)(PetscInt, PetscInt, PetscInt,
758                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
759                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
760                                                                 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
761                                                       void (**)(PetscInt, PetscInt, PetscInt,
762                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
763                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
764                                                                 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
765                                                       void (**)(PetscInt, PetscInt, PetscInt,
766                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
767                                                                 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
768                                                                 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
769 PETSC_EXTERN PetscErrorCode PetscDSSetDynamicJacobian(PetscDS, PetscInt, PetscInt,
770                                                       void (*)(PetscInt, PetscInt, PetscInt,
771                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
772                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
773                                                                PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
774                                                       void (*)(PetscInt, PetscInt, PetscInt,
775                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
776                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
777                                                                PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
778                                                       void (*)(PetscInt, PetscInt, PetscInt,
779                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
780                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
781                                                                PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
782                                                       void (*)(PetscInt, PetscInt, PetscInt,
783                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
784                                                                const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
785                                                                PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
786 PETSC_EXTERN PetscErrorCode PetscDSGetRiemannSolver(PetscDS, PetscInt,
787                                                     void (**)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *));
788 PETSC_EXTERN PetscErrorCode PetscDSSetRiemannSolver(PetscDS, PetscInt,
789                                                     void (*)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *));
790 PETSC_EXTERN PetscErrorCode PetscDSGetUpdate(PetscDS, PetscInt,
791                                              void (**)(PetscInt, PetscInt, PetscInt,
792                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
793                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
794                                                        PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
795 PETSC_EXTERN PetscErrorCode PetscDSSetUpdate(PetscDS, PetscInt,
796                                              void (*)(PetscInt, PetscInt, PetscInt,
797                                                       const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
798                                                       const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
799                                                       PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
800 PETSC_EXTERN PetscErrorCode PetscDSGetContext(PetscDS, PetscInt, void *);
801 PETSC_EXTERN PetscErrorCode PetscDSSetContext(PetscDS, PetscInt, void *);
802 PETSC_EXTERN PetscErrorCode PetscDSGetBdResidual(PetscDS, PetscInt,
803                                                  void (**)(PetscInt, PetscInt, PetscInt,
804                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
805                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
806                                                            PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
807                                                  void (**)(PetscInt, PetscInt, PetscInt,
808                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
809                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
810                                                            PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
811 PETSC_EXTERN PetscErrorCode PetscDSSetBdResidual(PetscDS, PetscInt,
812                                                  void (*)(PetscInt, PetscInt, PetscInt,
813                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
814                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
815                                                           PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
816                                                  void (*)(PetscInt, PetscInt, PetscInt,
817                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
818                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
819                                                           PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
820 PETSC_EXTERN PetscErrorCode PetscDSHasBdJacobian(PetscDS, PetscBool *);
821 PETSC_EXTERN PetscErrorCode PetscDSGetBdJacobian(PetscDS, PetscInt, PetscInt,
822                                                  void (**)(PetscInt, PetscInt, PetscInt,
823                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
824                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
825                                                            PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
826                                                  void (**)(PetscInt, PetscInt, PetscInt,
827                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
828                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
829                                                            PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
830                                                  void (**)(PetscInt, PetscInt, PetscInt,
831                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
832                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
833                                                            PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
834                                                  void (**)(PetscInt, PetscInt, PetscInt,
835                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
836                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
837                                                            PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
838 PETSC_EXTERN PetscErrorCode PetscDSSetBdJacobian(PetscDS, PetscInt, PetscInt,
839                                                  void (*)(PetscInt, PetscInt, PetscInt,
840                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
841                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
842                                                           PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
843                                                  void (*)(PetscInt, PetscInt, PetscInt,
844                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
845                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
846                                                           PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
847                                                  void (*)(PetscInt, PetscInt, PetscInt,
848                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
849                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
850                                                           PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
851                                                  void (*)(PetscInt, PetscInt, PetscInt,
852                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
853                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
854                                                           PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
855 PETSC_EXTERN PetscErrorCode PetscDSHasBdJacobianPreconditioner(PetscDS, PetscBool *);
856 PETSC_EXTERN PetscErrorCode PetscDSGetBdJacobianPreconditioner(PetscDS, PetscInt, PetscInt,
857                                                                void (**)(PetscInt, PetscInt, PetscInt,
858                                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
859                                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
860                                                                          PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
861                                                                void (**)(PetscInt, PetscInt, PetscInt,
862                                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
863                                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
864                                                                          PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
865                                                                void (**)(PetscInt, PetscInt, PetscInt,
866                                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
867                                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
868                                                                          PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
869                                                                void (**)(PetscInt, PetscInt, PetscInt,
870                                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
871                                                                          const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
872                                                                          PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
873 PETSC_EXTERN PetscErrorCode PetscDSSetBdJacobianPreconditioner(PetscDS, PetscInt, PetscInt,
874                                                                void (*)(PetscInt, PetscInt, PetscInt,
875                                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
876                                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
877                                                                         PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
878                                                                void (*)(PetscInt, PetscInt, PetscInt,
879                                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
880                                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
881                                                                         PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
882                                                                void (*)(PetscInt, PetscInt, PetscInt,
883                                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
884                                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
885                                                                         PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
886                                                                void (*)(PetscInt, PetscInt, PetscInt,
887                                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
888                                                                         const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
889                                                                         PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]));
890 PETSC_EXTERN PetscErrorCode PetscDSGetExactSolution(PetscDS, PetscInt, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *), void **);
891 PETSC_EXTERN PetscErrorCode PetscDSSetExactSolution(PetscDS, PetscInt, PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *), void *);
892 PETSC_EXTERN PetscErrorCode PetscDSGetExactSolutionTimeDerivative(PetscDS, PetscInt, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *), void **);
893 PETSC_EXTERN PetscErrorCode PetscDSSetExactSolutionTimeDerivative(PetscDS, PetscInt, PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *), void *);
894 PETSC_EXTERN PetscErrorCode PetscDSGetTabulation(PetscDS, PetscTabulation *[]);
895 PETSC_EXTERN PetscErrorCode PetscDSGetFaceTabulation(PetscDS, PetscTabulation *[]);
896 PETSC_EXTERN PetscErrorCode PetscDSGetEvaluationArrays(PetscDS, PetscScalar **, PetscScalar **, PetscScalar **);
897 PETSC_EXTERN PetscErrorCode PetscDSGetWeakFormArrays(PetscDS, PetscScalar **, PetscScalar **, PetscScalar **, PetscScalar **, PetscScalar **, PetscScalar **);
898 PETSC_EXTERN PetscErrorCode PetscDSGetWorkspace(PetscDS, PetscReal **, PetscScalar **, PetscScalar **, PetscScalar **, PetscScalar **);
899 PETSC_EXTERN PetscErrorCode PetscDSCopyConstants(PetscDS, PetscDS);
900 PETSC_EXTERN PetscErrorCode PetscDSCopyExactSolutions(PetscDS, PetscDS);
901 PETSC_EXTERN PetscErrorCode PetscDSCopyEquations(PetscDS, PetscDS);
902 PETSC_EXTERN PetscErrorCode PetscDSSelectDiscretizations(PetscDS, PetscInt, const PetscInt[], PetscDS);
903 PETSC_EXTERN PetscErrorCode PetscDSSelectEquations(PetscDS, PetscInt, const PetscInt[], PetscDS);
904 PETSC_EXTERN PetscErrorCode PetscDSAddBoundary(PetscDS, DMBoundaryConditionType, const char[], DMLabel, PetscInt, const PetscInt[], PetscInt, PetscInt, const PetscInt[], void (*)(void), void (*)(void), void *, PetscInt *);
905 PETSC_EXTERN PetscErrorCode PetscDSAddBoundaryByName(PetscDS, DMBoundaryConditionType, const char[], const char[], PetscInt, const PetscInt[], PetscInt, PetscInt, const PetscInt[], void (*)(void), void (*)(void), void *, PetscInt *);
906 PETSC_EXTERN PetscErrorCode PetscDSUpdateBoundary(PetscDS, PetscInt, DMBoundaryConditionType, const char[], DMLabel, PetscInt, const PetscInt[], PetscInt, PetscInt, const PetscInt[], void (*)(void), void (*)(void), void *);
907 PETSC_EXTERN PetscErrorCode PetscDSGetNumBoundary(PetscDS, PetscInt *);
908 PETSC_EXTERN PetscErrorCode PetscDSGetBoundary(PetscDS, PetscInt, PetscWeakForm *, DMBoundaryConditionType *, const char *[], DMLabel *, PetscInt *, const PetscInt *[], PetscInt *, PetscInt *, const PetscInt *[], void (**)(void), void (**)(void), void **);
909 PETSC_EXTERN PetscErrorCode PetscDSCopyBoundary(PetscDS, PetscInt, const PetscInt[], PetscDS);
910 PETSC_EXTERN PetscErrorCode PetscDSDestroyBoundary(PetscDS);
911 
912 #endif
913