xref: /petsc/src/dm/dt/interface/dtweakform.c (revision 2fb1d9daa8d289b4e375fe89ef72d82069342f85)
1 #include <petsc/private/petscdsimpl.h> /*I "petscds.h" I*/
2 
3 PetscClassId PETSCWEAKFORM_CLASSID = 0;
4 
5 static PetscErrorCode PetscChunkBufferCreate(size_t unitbytes, size_t expected, PetscChunkBuffer **buffer)
6 {
7   PetscErrorCode ierr;
8 
9   PetscFunctionBegin;
10   ierr = PetscNew(buffer);CHKERRQ(ierr);
11   ierr = PetscCalloc1(expected*unitbytes, &(*buffer)->array);CHKERRQ(ierr);
12   (*buffer)->size      = expected;
13   (*buffer)->unitbytes = unitbytes;
14   (*buffer)->alloc     = expected*unitbytes;
15   PetscFunctionReturn(0);
16 }
17 
18 static PetscErrorCode PetscChunkBufferDestroy(PetscChunkBuffer **buffer)
19 {
20   PetscErrorCode ierr;
21 
22   PetscFunctionBegin;
23   ierr = PetscFree((*buffer)->array);CHKERRQ(ierr);
24   ierr = PetscFree(*buffer);CHKERRQ(ierr);
25   PetscFunctionReturn(0);
26 }
27 
28 static PetscErrorCode PetscChunkBufferCreateChunk(PetscChunkBuffer *buffer, PetscInt size, PetscChunk *chunk)
29 {
30   PetscErrorCode ierr;
31 
32   PetscFunctionBegin;
33   if ((buffer->size + size)*buffer->unitbytes > buffer->alloc) {
34     char *tmp;
35 
36     if (!buffer->alloc) buffer->alloc = (buffer->size + size)*buffer->unitbytes;
37     while ((buffer->size + size)*buffer->unitbytes > buffer->alloc) buffer->alloc *= 2;
38     ierr = PetscMalloc(buffer->alloc, &tmp);CHKERRQ(ierr);
39     ierr = PetscMemcpy(tmp, buffer->array, buffer->size*buffer->unitbytes);CHKERRQ(ierr);
40     ierr = PetscFree(buffer->array);CHKERRQ(ierr);
41     buffer->array = tmp;
42   }
43   chunk->start    = buffer->size*buffer->unitbytes;
44   chunk->size     = size;
45   chunk->reserved = size;
46   buffer->size   += size;
47   PetscFunctionReturn(0);
48 }
49 
50 static PetscErrorCode PetscChunkBufferEnlargeChunk(PetscChunkBuffer *buffer, PetscInt size, PetscChunk *chunk)
51 {
52   size_t         siz = size;
53   PetscErrorCode ierr;
54 
55   PetscFunctionBegin;
56   if (chunk->size + size > chunk->reserved) {
57     PetscChunk newchunk;
58     PetscInt   reserved = chunk->size;
59 
60     /* TODO Here if we had a chunk list, we could update them all to reclaim unused space */
61     while (reserved < chunk->size+size) reserved *= 2;
62     ierr = PetscChunkBufferCreateChunk(buffer, (size_t) reserved, &newchunk);CHKERRQ(ierr);
63     newchunk.size = chunk->size+size;
64     ierr = PetscMemcpy(&buffer->array[newchunk.start], &buffer->array[chunk->start], chunk->size * buffer->unitbytes);CHKERRQ(ierr);
65     *chunk = newchunk;
66   } else {
67     chunk->size += siz;
68   }
69   PetscFunctionReturn(0);
70 }
71 
72 /*@C
73   PetscHashFormKeySort - Sorts an array of PetscHashFormKey in place in increasing order.
74 
75   Not Collective
76 
77   Input Parameters:
78 + n - number of values
79 - X - array of PetscHashFormKey
80 
81   Level: intermediate
82 
83 .seealso: PetscIntSortSemiOrdered(), PetscSortInt()
84 @*/
85 PetscErrorCode PetscHashFormKeySort(PetscInt n, PetscHashFormKey arr[])
86 {
87   PetscErrorCode ierr;
88 
89   PetscFunctionBegin;
90   if (n <= 1) PetscFunctionReturn(0);
91   PetscValidPointer(arr, 2);
92   ierr = PetscTimSort(n, arr, sizeof(PetscHashFormKey), Compare_PetscHashFormKey_Private, NULL);CHKERRQ(ierr);
93   PetscFunctionReturn(0);
94 }
95 
96 PetscErrorCode PetscWeakFormGetFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt *n, void (***func)())
97 {
98   PetscHashFormKey key;
99   PetscChunk       chunk;
100   PetscErrorCode   ierr;
101 
102   PetscFunctionBegin;
103   key.label = label; key.value = value; key.field = f;
104   ierr = PetscHMapFormGet(ht, key, &chunk);CHKERRQ(ierr);
105   if (chunk.size < 0) {*n = 0;          *func = NULL;}
106   else                {*n = chunk.size; *func = (void (**)()) &wf->funcs->array[chunk.start];}
107   PetscFunctionReturn(0);
108 }
109 
110 /* A NULL argument for func causes this to clear the key */
111 PetscErrorCode PetscWeakFormSetFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt n, void (**func)())
112 {
113   PetscHashFormKey key;
114   PetscChunk       chunk;
115   PetscInt         i;
116   PetscErrorCode   ierr;
117 
118   PetscFunctionBegin;
119   key.label = label; key.value = value; key.field = f;
120   if (!func) {
121     ierr = PetscHMapFormDel(ht, key);CHKERRQ(ierr);
122     PetscFunctionReturn(0);
123   } else {
124     ierr = PetscHMapFormGet(ht, key, &chunk);CHKERRQ(ierr);
125   }
126   if (chunk.size < 0) {
127     ierr = PetscChunkBufferCreateChunk(wf->funcs, n, &chunk);CHKERRQ(ierr);
128     ierr = PetscHMapFormSet(ht, key, chunk);CHKERRQ(ierr);
129   } else if (chunk.size <= n) {
130     ierr = PetscChunkBufferEnlargeChunk(wf->funcs, n - chunk.size, &chunk);CHKERRQ(ierr);
131     ierr = PetscHMapFormSet(ht, key, chunk);CHKERRQ(ierr);
132   }
133   for (i = 0; i < n; ++i) ((void (**)()) &wf->funcs->array[chunk.start])[i] = func[i];
134   PetscFunctionReturn(0);
135 }
136 
137 PetscErrorCode PetscWeakFormAddFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, void (*func)())
138 {
139   PetscHashFormKey key;
140   PetscChunk       chunk;
141   PetscErrorCode   ierr;
142 
143   PetscFunctionBegin;
144   if (!func) PetscFunctionReturn(0);
145   key.label = label; key.value = value; key.field = f;
146   ierr = PetscHMapFormGet(ht, key, &chunk);CHKERRQ(ierr);
147   if (chunk.size < 0) {
148     ierr = PetscChunkBufferCreateChunk(wf->funcs, 1, &chunk);CHKERRQ(ierr);
149     ierr = PetscHMapFormSet(ht, key, chunk);CHKERRQ(ierr);
150     ((void (**)()) &wf->funcs->array[chunk.start])[0] = func;
151   } else {
152     ierr = PetscChunkBufferEnlargeChunk(wf->funcs, 1, &chunk);CHKERRQ(ierr);
153     ierr = PetscHMapFormSet(ht, key, chunk);CHKERRQ(ierr);
154     ((void (**)()) &wf->funcs->array[chunk.start])[chunk.size-1] = func;
155   }
156   PetscFunctionReturn(0);
157 }
158 
159 PetscErrorCode PetscWeakFormGetIndexFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt ind, void (**func)())
160 {
161   PetscHashFormKey key;
162   PetscChunk       chunk;
163   PetscErrorCode   ierr;
164 
165   PetscFunctionBegin;
166   key.label = label; key.value = value; key.field = f;
167   ierr = PetscHMapFormGet(ht, key, &chunk);CHKERRQ(ierr);
168   if (chunk.size < 0) {*func = NULL;}
169   else {
170     if (ind >= chunk.size) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Index %D not in [0, %D)", ind, chunk.size);
171     *func = ((void (**)()) &wf->funcs->array[chunk.start])[ind];
172   }
173   PetscFunctionReturn(0);
174 }
175 
176 /* A NULL argument for func causes this to clear the slot, and if there is nothing else, clear the key */
177 PetscErrorCode PetscWeakFormSetIndexFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt ind, void (*func)())
178 {
179   PetscHashFormKey key;
180   PetscChunk       chunk;
181   PetscErrorCode   ierr;
182 
183   PetscFunctionBegin;
184   key.label = label; key.value = value; key.field = f;
185   ierr = PetscHMapFormGet(ht, key, &chunk);CHKERRQ(ierr);
186   if (chunk.size < 0) {
187     if (!func) PetscFunctionReturn(0);
188     ierr = PetscChunkBufferCreateChunk(wf->funcs, ind+1, &chunk);CHKERRQ(ierr);
189     ierr = PetscHMapFormSet(ht, key, chunk);CHKERRQ(ierr);
190   } else if (!func && !ind && chunk.size == 1) {
191     ierr = PetscHMapFormDel(ht, key);CHKERRQ(ierr);
192     PetscFunctionReturn(0);
193   } else if (chunk.size <= ind) {
194     ierr = PetscChunkBufferEnlargeChunk(wf->funcs, ind - chunk.size + 1, &chunk);CHKERRQ(ierr);
195     ierr = PetscHMapFormSet(ht, key, chunk);CHKERRQ(ierr);
196   }
197   ((void (**)()) &wf->funcs->array[chunk.start])[ind] = func;
198   PetscFunctionReturn(0);
199 }
200 
201 PetscErrorCode PetscWeakFormGetObjective(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt *n,
202                                          void (***obj)(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, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
206 {
207   PetscErrorCode ierr;
208 
209   PetscFunctionBegin;
210   ierr = PetscWeakFormGetFunction_Private(wf, wf->obj, label, val, f, n, (void (***)(void)) obj);CHKERRQ(ierr);
211   PetscFunctionReturn(0);
212 }
213 
214 PetscErrorCode PetscWeakFormSetObjective(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt n,
215                                          void (**obj)(PetscInt, PetscInt, PetscInt,
216                                                       const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const  PetscScalar[],
217                                                       const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
218                                                       PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
219 {
220   PetscErrorCode ierr;
221 
222   PetscFunctionBegin;
223   ierr = PetscWeakFormSetFunction_Private(wf, wf->obj, label, val, f, n, (void (**)(void)) obj);CHKERRQ(ierr);
224   PetscFunctionReturn(0);
225 }
226 
227 PetscErrorCode PetscWeakFormAddObjective(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f,
228                                          void (*obj)(PetscInt, PetscInt, PetscInt,
229                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
230                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
231                                                      PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
232 {
233   PetscErrorCode ierr;
234 
235   PetscFunctionBegin;
236   ierr = PetscWeakFormAddFunction_Private(wf, wf->obj, label, val, f, (void (*)(void)) obj);CHKERRQ(ierr);
237   PetscFunctionReturn(0);
238 }
239 
240 PetscErrorCode PetscWeakFormGetIndexObjective(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt ind,
241                                               void (**obj)(PetscInt, PetscInt, PetscInt,
242                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
243                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
244                                                            PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
245 {
246   PetscErrorCode ierr;
247 
248   PetscFunctionBegin;
249   ierr = PetscWeakFormGetIndexFunction_Private(wf, wf->obj, label, val, f, ind, (void (**)(void)) obj);CHKERRQ(ierr);
250   PetscFunctionReturn(0);
251 }
252 
253 PetscErrorCode PetscWeakFormSetIndexObjective(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt ind,
254                                               void (*obj)(PetscInt, PetscInt, PetscInt,
255                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
256                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
257                                                           PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
258 {
259   PetscErrorCode ierr;
260 
261   PetscFunctionBegin;
262   ierr = PetscWeakFormSetIndexFunction_Private(wf, wf->obj, label, val, f, ind, (void (*)(void)) obj);CHKERRQ(ierr);
263   PetscFunctionReturn(0);
264 }
265 
266 PetscErrorCode PetscWeakFormGetResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f,
267                                         PetscInt *n0,
268                                         void (***f0)(PetscInt, PetscInt, PetscInt,
269                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
270                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
271                                                      PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
272                                         PetscInt *n1,
273                                         void (***f1)(PetscInt, PetscInt, PetscInt,
274                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
275                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
276                                                      PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
277 {
278   PetscErrorCode ierr;
279 
280   PetscFunctionBegin;
281   ierr = PetscWeakFormGetFunction_Private(wf, wf->f0, label, val, f, n0, (void (***)(void)) f0);CHKERRQ(ierr);
282   ierr = PetscWeakFormGetFunction_Private(wf, wf->f1, label, val, f, n1, (void (***)(void)) f1);CHKERRQ(ierr);
283   PetscFunctionReturn(0);
284 }
285 
286 PetscErrorCode PetscWeakFormAddResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f,
287                                         void (*f0)(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, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
291                                         void (*f1)(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, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
295 {
296   PetscErrorCode ierr;
297 
298   PetscFunctionBegin;
299   ierr = PetscWeakFormAddFunction_Private(wf, wf->f0, label, val, f, (void (*)(void)) f0);CHKERRQ(ierr);
300   ierr = PetscWeakFormAddFunction_Private(wf, wf->f1, label, val, f, (void (*)(void)) f1);CHKERRQ(ierr);
301   PetscFunctionReturn(0);
302 }
303 
304 PetscErrorCode PetscWeakFormSetResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f,
305                                         PetscInt n0,
306                                         void (**f0)(PetscInt, PetscInt, PetscInt,
307                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
308                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
309                                                    PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
310                                         PetscInt n1,
311                                         void (**f1)(PetscInt, PetscInt, PetscInt,
312                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
313                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
314                                                    PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
315 {
316   PetscErrorCode ierr;
317 
318   PetscFunctionBegin;
319   ierr = PetscWeakFormSetFunction_Private(wf, wf->f0, label, val, f, n0, (void (**)(void)) f0);CHKERRQ(ierr);
320   ierr = PetscWeakFormSetFunction_Private(wf, wf->f1, label, val, f, n1, (void (**)(void)) f1);CHKERRQ(ierr);
321   PetscFunctionReturn(0);
322 }
323 
324 PetscErrorCode PetscWeakFormSetIndexResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f,
325                                         PetscInt i0,
326                                         void (*f0)(PetscInt, PetscInt, PetscInt,
327                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
328                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
329                                                    PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
330                                         PetscInt i1,
331                                         void (*f1)(PetscInt, PetscInt, PetscInt,
332                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
333                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
334                                                    PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
335 {
336   PetscErrorCode ierr;
337 
338   PetscFunctionBegin;
339   ierr = PetscWeakFormSetIndexFunction_Private(wf, wf->f0, label, val, f, i0, (void (*)(void)) f0);CHKERRQ(ierr);
340   ierr = PetscWeakFormSetIndexFunction_Private(wf, wf->f1, label, val, f, i1, (void (*)(void)) f1);CHKERRQ(ierr);
341   PetscFunctionReturn(0);
342 }
343 
344 PetscErrorCode PetscWeakFormGetBdResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f,
345                                           PetscInt *n0,
346                                         void (***f0)(PetscInt, PetscInt, PetscInt,
347                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
348                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
349                                                      PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
350                                         PetscInt *n1,
351                                         void (***f1)(PetscInt, PetscInt, PetscInt,
352                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
353                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
354                                                      PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
355 {
356   PetscErrorCode ierr;
357 
358   PetscFunctionBegin;
359   ierr = PetscWeakFormGetFunction_Private(wf, wf->bdf0, label, val, f, n0, (void (***)(void)) f0);CHKERRQ(ierr);
360   ierr = PetscWeakFormGetFunction_Private(wf, wf->bdf1, label, val, f, n1, (void (***)(void)) f1);CHKERRQ(ierr);
361   PetscFunctionReturn(0);
362 }
363 
364 PetscErrorCode PetscWeakFormAddBdResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f,
365                                           void (*f0)(PetscInt, PetscInt, PetscInt,
366                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
367                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
368                                                      PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
369                                           void (*f1)(PetscInt, PetscInt, PetscInt,
370                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
371                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
372                                                      PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
373 {
374   PetscErrorCode ierr;
375 
376   PetscFunctionBegin;
377   ierr = PetscWeakFormAddFunction_Private(wf, wf->bdf0, label, val, f, (void (*)(void)) f0);CHKERRQ(ierr);
378   ierr = PetscWeakFormAddFunction_Private(wf, wf->bdf1, label, val, f, (void (*)(void)) f1);CHKERRQ(ierr);
379   PetscFunctionReturn(0);
380 }
381 
382 PetscErrorCode PetscWeakFormSetBdResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f,
383                                           PetscInt n0,
384                                           void (**f0)(PetscInt, PetscInt, PetscInt,
385                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
386                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
387                                                      PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
388                                           PetscInt n1,
389                                           void (**f1)(PetscInt, PetscInt, PetscInt,
390                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
391                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
392                                                      PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
393 {
394   PetscErrorCode ierr;
395 
396   PetscFunctionBegin;
397   ierr = PetscWeakFormSetFunction_Private(wf, wf->bdf0, label, val, f, n0, (void (**)(void)) f0);CHKERRQ(ierr);
398   ierr = PetscWeakFormSetFunction_Private(wf, wf->bdf1, label, val, f, n1, (void (**)(void)) f1);CHKERRQ(ierr);
399   PetscFunctionReturn(0);
400 }
401 
402 PetscErrorCode PetscWeakFormSetIndexBdResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f,
403                                           PetscInt i0,
404                                           void (*f0)(PetscInt, PetscInt, PetscInt,
405                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
406                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
407                                                      PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
408                                           PetscInt i1,
409                                           void (*f1)(PetscInt, PetscInt, PetscInt,
410                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
411                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
412                                                      PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
413 {
414   PetscErrorCode ierr;
415 
416   PetscFunctionBegin;
417   ierr = PetscWeakFormSetIndexFunction_Private(wf, wf->bdf0, label, val, f, i0, (void (*)(void)) f0);CHKERRQ(ierr);
418   ierr = PetscWeakFormSetIndexFunction_Private(wf, wf->bdf1, label, val, f, i1, (void (*)(void)) f1);CHKERRQ(ierr);
419   PetscFunctionReturn(0);
420 }
421 
422 PetscErrorCode PetscWeakFormHasJacobian(PetscWeakForm wf, PetscBool *hasJac)
423 {
424   PetscInt       n0, n1, n2, n3;
425   PetscErrorCode ierr;
426 
427   PetscFunctionBegin;
428   PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1);
429   PetscValidBoolPointer(hasJac, 2);
430   ierr = PetscHMapFormGetSize(wf->g0, &n0);CHKERRQ(ierr);
431   ierr = PetscHMapFormGetSize(wf->g1, &n1);CHKERRQ(ierr);
432   ierr = PetscHMapFormGetSize(wf->g2, &n2);CHKERRQ(ierr);
433   ierr = PetscHMapFormGetSize(wf->g3, &n3);CHKERRQ(ierr);
434   *hasJac = n0+n1+n2+n3 ? PETSC_TRUE : PETSC_FALSE;
435   PetscFunctionReturn(0);
436 }
437 
438 PetscErrorCode PetscWeakFormGetJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g,
439                                         PetscInt *n0,
440                                         void (***g0)(PetscInt, PetscInt, PetscInt,
441                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
442                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
443                                                      PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
444                                         PetscInt *n1,
445                                         void (***g1)(PetscInt, PetscInt, PetscInt,
446                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
447                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
448                                                      PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
449                                         PetscInt *n2,
450                                         void (***g2)(PetscInt, PetscInt, PetscInt,
451                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
452                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
453                                                      PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
454                                         PetscInt *n3,
455                                         void (***g3)(PetscInt, PetscInt, PetscInt,
456                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
457                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
458                                                      PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
459 {
460   PetscInt       find = f*wf->Nf + g;
461   PetscErrorCode ierr;
462 
463   PetscFunctionBegin;
464   ierr = PetscWeakFormGetFunction_Private(wf, wf->g0, label, val, find, n0, (void (***)(void)) g0);CHKERRQ(ierr);
465   ierr = PetscWeakFormGetFunction_Private(wf, wf->g1, label, val, find, n1, (void (***)(void)) g1);CHKERRQ(ierr);
466   ierr = PetscWeakFormGetFunction_Private(wf, wf->g2, label, val, find, n2, (void (***)(void)) g2);CHKERRQ(ierr);
467   ierr = PetscWeakFormGetFunction_Private(wf, wf->g3, label, val, find, n3, (void (***)(void)) g3);CHKERRQ(ierr);
468   PetscFunctionReturn(0);
469 }
470 
471 PetscErrorCode PetscWeakFormAddJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g,
472                                         void (*g0)(PetscInt, PetscInt, PetscInt,
473                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
474                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
475                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
476                                         void (*g1)(PetscInt, PetscInt, PetscInt,
477                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
478                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
479                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
480                                         void (*g2)(PetscInt, PetscInt, PetscInt,
481                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
482                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
483                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
484                                         void (*g3)(PetscInt, PetscInt, PetscInt,
485                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
486                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
487                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
488 {
489   PetscInt       find = f*wf->Nf + g;
490   PetscErrorCode ierr;
491 
492   PetscFunctionBegin;
493   ierr = PetscWeakFormAddFunction_Private(wf, wf->g0, label, val, find, (void (*)(void)) g0);CHKERRQ(ierr);
494   ierr = PetscWeakFormAddFunction_Private(wf, wf->g1, label, val, find, (void (*)(void)) g1);CHKERRQ(ierr);
495   ierr = PetscWeakFormAddFunction_Private(wf, wf->g2, label, val, find, (void (*)(void)) g2);CHKERRQ(ierr);
496   ierr = PetscWeakFormAddFunction_Private(wf, wf->g3, label, val, find, (void (*)(void)) g3);CHKERRQ(ierr);
497   PetscFunctionReturn(0);
498 }
499 
500 PetscErrorCode PetscWeakFormSetJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g,
501                                         PetscInt n0,
502                                         void (**g0)(PetscInt, PetscInt, PetscInt,
503                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
504                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
505                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
506                                         PetscInt n1,
507                                         void (**g1)(PetscInt, PetscInt, PetscInt,
508                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
509                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
510                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
511                                         PetscInt n2,
512                                         void (**g2)(PetscInt, PetscInt, PetscInt,
513                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
514                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
515                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
516                                         PetscInt n3,
517                                         void (**g3)(PetscInt, PetscInt, PetscInt,
518                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
519                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
520                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
521 {
522   PetscInt       find = f*wf->Nf + g;
523   PetscErrorCode ierr;
524 
525   PetscFunctionBegin;
526   ierr = PetscWeakFormSetFunction_Private(wf, wf->g0, label, val, find, n0, (void (**)(void)) g0);CHKERRQ(ierr);
527   ierr = PetscWeakFormSetFunction_Private(wf, wf->g1, label, val, find, n1, (void (**)(void)) g1);CHKERRQ(ierr);
528   ierr = PetscWeakFormSetFunction_Private(wf, wf->g2, label, val, find, n2, (void (**)(void)) g2);CHKERRQ(ierr);
529   ierr = PetscWeakFormSetFunction_Private(wf, wf->g3, label, val, find, n3, (void (**)(void)) g3);CHKERRQ(ierr);
530   PetscFunctionReturn(0);
531 }
532 
533 PetscErrorCode PetscWeakFormSetIndexJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g,
534                                         PetscInt i0,
535                                         void (*g0)(PetscInt, PetscInt, PetscInt,
536                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
537                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
538                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
539                                         PetscInt i1,
540                                         void (*g1)(PetscInt, PetscInt, PetscInt,
541                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
542                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
543                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
544                                         PetscInt i2,
545                                         void (*g2)(PetscInt, PetscInt, PetscInt,
546                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
547                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
548                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
549                                         PetscInt i3,
550                                         void (*g3)(PetscInt, PetscInt, PetscInt,
551                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
552                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
553                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
554 {
555   PetscInt       find = f*wf->Nf + g;
556   PetscErrorCode ierr;
557 
558   PetscFunctionBegin;
559   ierr = PetscWeakFormSetIndexFunction_Private(wf, wf->g0, label, val, find, i0, (void (*)(void)) g0);CHKERRQ(ierr);
560   ierr = PetscWeakFormSetIndexFunction_Private(wf, wf->g1, label, val, find, i1, (void (*)(void)) g1);CHKERRQ(ierr);
561   ierr = PetscWeakFormSetIndexFunction_Private(wf, wf->g2, label, val, find, i2, (void (*)(void)) g2);CHKERRQ(ierr);
562   ierr = PetscWeakFormSetIndexFunction_Private(wf, wf->g3, label, val, find, i3, (void (*)(void)) g3);CHKERRQ(ierr);
563   PetscFunctionReturn(0);
564 }
565 
566 PetscErrorCode PetscWeakFormHasJacobianPreconditioner(PetscWeakForm wf, PetscBool *hasJacPre)
567 {
568   PetscInt       n0, n1, n2, n3;
569   PetscErrorCode ierr;
570 
571   PetscFunctionBegin;
572   PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1);
573   PetscValidBoolPointer(hasJacPre, 2);
574   ierr = PetscHMapFormGetSize(wf->gp0, &n0);CHKERRQ(ierr);
575   ierr = PetscHMapFormGetSize(wf->gp1, &n1);CHKERRQ(ierr);
576   ierr = PetscHMapFormGetSize(wf->gp2, &n2);CHKERRQ(ierr);
577   ierr = PetscHMapFormGetSize(wf->gp3, &n3);CHKERRQ(ierr);
578   *hasJacPre = n0+n1+n2+n3 ? PETSC_TRUE : PETSC_FALSE;
579   PetscFunctionReturn(0);
580 }
581 
582 PetscErrorCode PetscWeakFormGetJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g,
583                                                       PetscInt *n0,
584                                                       void (***g0)(PetscInt, PetscInt, PetscInt,
585                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
586                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
587                                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
588                                                       PetscInt *n1,
589                                                       void (***g1)(PetscInt, PetscInt, PetscInt,
590                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
591                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
592                                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
593                                                       PetscInt *n2,
594                                                       void (***g2)(PetscInt, PetscInt, PetscInt,
595                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
596                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
597                                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
598                                                       PetscInt *n3,
599                                                       void (***g3)(PetscInt, PetscInt, PetscInt,
600                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
601                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
602                                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
603 {
604   PetscInt       find = f*wf->Nf + g;
605   PetscErrorCode ierr;
606 
607   PetscFunctionBegin;
608   ierr = PetscWeakFormGetFunction_Private(wf, wf->gp0, label, val, find, n0, (void (***)(void)) g0);CHKERRQ(ierr);
609   ierr = PetscWeakFormGetFunction_Private(wf, wf->gp1, label, val, find, n1, (void (***)(void)) g1);CHKERRQ(ierr);
610   ierr = PetscWeakFormGetFunction_Private(wf, wf->gp2, label, val, find, n2, (void (***)(void)) g2);CHKERRQ(ierr);
611   ierr = PetscWeakFormGetFunction_Private(wf, wf->gp3, label, val, find, n3, (void (***)(void)) g3);CHKERRQ(ierr);
612   PetscFunctionReturn(0);
613 }
614 
615 PetscErrorCode PetscWeakFormAddJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g,
616                                         void (*g0)(PetscInt, PetscInt, PetscInt,
617                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
618                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
619                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
620                                         void (*g1)(PetscInt, PetscInt, PetscInt,
621                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
622                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
623                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
624                                         void (*g2)(PetscInt, PetscInt, PetscInt,
625                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
626                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
627                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
628                                         void (*g3)(PetscInt, PetscInt, PetscInt,
629                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
630                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
631                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
632 {
633   PetscInt       find = f*wf->Nf + g;
634   PetscErrorCode ierr;
635 
636   PetscFunctionBegin;
637   ierr = PetscWeakFormAddFunction_Private(wf, wf->gp0, label, val, find, (void (*)(void)) g0);CHKERRQ(ierr);
638   ierr = PetscWeakFormAddFunction_Private(wf, wf->gp1, label, val, find, (void (*)(void)) g1);CHKERRQ(ierr);
639   ierr = PetscWeakFormAddFunction_Private(wf, wf->gp2, label, val, find, (void (*)(void)) g2);CHKERRQ(ierr);
640   ierr = PetscWeakFormAddFunction_Private(wf, wf->gp3, label, val, find, (void (*)(void)) g3);CHKERRQ(ierr);
641   PetscFunctionReturn(0);
642 }
643 
644 PetscErrorCode PetscWeakFormSetJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g,
645                                                       PetscInt n0,
646                                                       void (**g0)(PetscInt, PetscInt, PetscInt,
647                                                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
648                                                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
649                                                                  PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
650                                                       PetscInt n1,
651                                                       void (**g1)(PetscInt, PetscInt, PetscInt,
652                                                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
653                                                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
654                                                                  PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
655                                                       PetscInt n2,
656                                                       void (**g2)(PetscInt, PetscInt, PetscInt,
657                                                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
658                                                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
659                                                                  PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
660                                                       PetscInt n3,
661                                                       void (**g3)(PetscInt, PetscInt, PetscInt,
662                                                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
663                                                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
664                                                                  PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
665 {
666   PetscInt       find = f*wf->Nf + g;
667   PetscErrorCode ierr;
668 
669   PetscFunctionBegin;
670   ierr = PetscWeakFormSetFunction_Private(wf, wf->gp0, label, val, find, n0, (void (**)(void)) g0);CHKERRQ(ierr);
671   ierr = PetscWeakFormSetFunction_Private(wf, wf->gp1, label, val, find, n1, (void (**)(void)) g1);CHKERRQ(ierr);
672   ierr = PetscWeakFormSetFunction_Private(wf, wf->gp2, label, val, find, n2, (void (**)(void)) g2);CHKERRQ(ierr);
673   ierr = PetscWeakFormSetFunction_Private(wf, wf->gp3, label, val, find, n3, (void (**)(void)) g3);CHKERRQ(ierr);
674   PetscFunctionReturn(0);
675 }
676 
677 PetscErrorCode PetscWeakFormSetIndexJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g,
678                                                       PetscInt i0,
679                                                       void (*g0)(PetscInt, PetscInt, PetscInt,
680                                                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
681                                                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
682                                                                  PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
683                                                       PetscInt i1,
684                                                       void (*g1)(PetscInt, PetscInt, PetscInt,
685                                                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
686                                                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
687                                                                  PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
688                                                       PetscInt i2,
689                                                       void (*g2)(PetscInt, PetscInt, PetscInt,
690                                                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
691                                                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
692                                                                  PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
693                                                       PetscInt i3,
694                                                       void (*g3)(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 {
699   PetscInt       find = f*wf->Nf + g;
700   PetscErrorCode ierr;
701 
702   PetscFunctionBegin;
703   ierr = PetscWeakFormSetIndexFunction_Private(wf, wf->gp0, label, val, find, i0, (void (*)(void)) g0);CHKERRQ(ierr);
704   ierr = PetscWeakFormSetIndexFunction_Private(wf, wf->gp1, label, val, find, i1, (void (*)(void)) g1);CHKERRQ(ierr);
705   ierr = PetscWeakFormSetIndexFunction_Private(wf, wf->gp2, label, val, find, i2, (void (*)(void)) g2);CHKERRQ(ierr);
706   ierr = PetscWeakFormSetIndexFunction_Private(wf, wf->gp3, label, val, find, i3, (void (*)(void)) g3);CHKERRQ(ierr);
707   PetscFunctionReturn(0);
708 }
709 
710 PetscErrorCode PetscWeakFormHasBdJacobian(PetscWeakForm wf, PetscBool *hasJac)
711 {
712   PetscInt       n0, n1, n2, n3;
713   PetscErrorCode ierr;
714 
715   PetscFunctionBegin;
716   PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1);
717   PetscValidBoolPointer(hasJac, 2);
718   ierr = PetscHMapFormGetSize(wf->bdg0, &n0);CHKERRQ(ierr);
719   ierr = PetscHMapFormGetSize(wf->bdg1, &n1);CHKERRQ(ierr);
720   ierr = PetscHMapFormGetSize(wf->bdg2, &n2);CHKERRQ(ierr);
721   ierr = PetscHMapFormGetSize(wf->bdg3, &n3);CHKERRQ(ierr);
722   *hasJac = n0+n1+n2+n3 ? PETSC_TRUE : PETSC_FALSE;
723   PetscFunctionReturn(0);
724 }
725 
726 PetscErrorCode PetscWeakFormGetBdJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g,
727                                           PetscInt *n0,
728                                           void (***g0)(PetscInt, PetscInt, PetscInt,
729                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
730                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
731                                                        PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
732                                           PetscInt *n1,
733                                           void (***g1)(PetscInt, PetscInt, PetscInt,
734                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
735                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
736                                                        PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
737                                           PetscInt *n2,
738                                           void (***g2)(PetscInt, PetscInt, PetscInt,
739                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
740                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
741                                                        PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
742                                           PetscInt *n3,
743                                           void (***g3)(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[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
747 {
748   PetscInt       find = f*wf->Nf + g;
749   PetscErrorCode ierr;
750 
751   PetscFunctionBegin;
752   ierr = PetscWeakFormGetFunction_Private(wf, wf->bdg0, label, val, find, n0, (void (***)(void)) g0);CHKERRQ(ierr);
753   ierr = PetscWeakFormGetFunction_Private(wf, wf->bdg1, label, val, find, n1, (void (***)(void)) g1);CHKERRQ(ierr);
754   ierr = PetscWeakFormGetFunction_Private(wf, wf->bdg2, label, val, find, n2, (void (***)(void)) g2);CHKERRQ(ierr);
755   ierr = PetscWeakFormGetFunction_Private(wf, wf->bdg3, label, val, find, n3, (void (***)(void)) g3);CHKERRQ(ierr);
756   PetscFunctionReturn(0);
757 }
758 
759 PetscErrorCode PetscWeakFormAddBdJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g,
760                                           void (*g0)(PetscInt, PetscInt, PetscInt,
761                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
762                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
763                                                      PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
764                                           void (*g1)(PetscInt, PetscInt, PetscInt,
765                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
766                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
767                                                      PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
768                                           void (*g2)(PetscInt, PetscInt, PetscInt,
769                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
770                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
771                                                      PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
772                                           void (*g3)(PetscInt, PetscInt, PetscInt,
773                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
774                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
775                                                      PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
776 {
777   PetscInt       find = f*wf->Nf + g;
778   PetscErrorCode ierr;
779 
780   PetscFunctionBegin;
781   ierr = PetscWeakFormAddFunction_Private(wf, wf->bdg0, label, val, find, (void (*)(void)) g0);CHKERRQ(ierr);
782   ierr = PetscWeakFormAddFunction_Private(wf, wf->bdg1, label, val, find, (void (*)(void)) g1);CHKERRQ(ierr);
783   ierr = PetscWeakFormAddFunction_Private(wf, wf->bdg2, label, val, find, (void (*)(void)) g2);CHKERRQ(ierr);
784   ierr = PetscWeakFormAddFunction_Private(wf, wf->bdg3, label, val, find, (void (*)(void)) g3);CHKERRQ(ierr);
785   PetscFunctionReturn(0);
786 }
787 
788 PetscErrorCode PetscWeakFormSetBdJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g,
789                                           PetscInt n0,
790                                           void (**g0)(PetscInt, PetscInt, PetscInt,
791                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
792                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
793                                                      PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
794                                           PetscInt n1,
795                                           void (**g1)(PetscInt, PetscInt, PetscInt,
796                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
797                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
798                                                      PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
799                                           PetscInt n2,
800                                           void (**g2)(PetscInt, PetscInt, PetscInt,
801                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
802                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
803                                                      PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
804                                           PetscInt n3,
805                                           void (**g3)(PetscInt, PetscInt, PetscInt,
806                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
807                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
808                                                      PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
809 {
810   PetscInt       find = f*wf->Nf + g;
811   PetscErrorCode ierr;
812 
813   PetscFunctionBegin;
814   ierr = PetscWeakFormSetFunction_Private(wf, wf->bdg0, label, val, find, n0, (void (**)(void)) g0);CHKERRQ(ierr);
815   ierr = PetscWeakFormSetFunction_Private(wf, wf->bdg1, label, val, find, n1, (void (**)(void)) g1);CHKERRQ(ierr);
816   ierr = PetscWeakFormSetFunction_Private(wf, wf->bdg2, label, val, find, n2, (void (**)(void)) g2);CHKERRQ(ierr);
817   ierr = PetscWeakFormSetFunction_Private(wf, wf->bdg3, label, val, find, n3, (void (**)(void)) g3);CHKERRQ(ierr);
818   PetscFunctionReturn(0);
819 }
820 
821 PetscErrorCode PetscWeakFormSetIndexBdJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g,
822                                           PetscInt i0,
823                                           void (*g0)(PetscInt, PetscInt, PetscInt,
824                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
825                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
826                                                      PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
827                                           PetscInt i1,
828                                           void (*g1)(PetscInt, PetscInt, PetscInt,
829                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
830                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
831                                                      PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
832                                           PetscInt i2,
833                                           void (*g2)(PetscInt, PetscInt, PetscInt,
834                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
835                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
836                                                      PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
837                                           PetscInt i3,
838                                           void (*g3)(PetscInt, PetscInt, PetscInt,
839                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
840                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
841                                                      PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
842 {
843   PetscInt       find = f*wf->Nf + g;
844   PetscErrorCode ierr;
845 
846   PetscFunctionBegin;
847   ierr = PetscWeakFormSetIndexFunction_Private(wf, wf->bdg0, label, val, find, i0, (void (*)(void)) g0);CHKERRQ(ierr);
848   ierr = PetscWeakFormSetIndexFunction_Private(wf, wf->bdg1, label, val, find, i1, (void (*)(void)) g1);CHKERRQ(ierr);
849   ierr = PetscWeakFormSetIndexFunction_Private(wf, wf->bdg2, label, val, find, i2, (void (*)(void)) g2);CHKERRQ(ierr);
850   ierr = PetscWeakFormSetIndexFunction_Private(wf, wf->bdg3, label, val, find, i3, (void (*)(void)) g3);CHKERRQ(ierr);
851   PetscFunctionReturn(0);
852 }
853 
854 PetscErrorCode PetscWeakFormHasBdJacobianPreconditioner(PetscWeakForm wf, PetscBool *hasJacPre)
855 {
856   PetscInt       n0, n1, n2, n3;
857   PetscErrorCode ierr;
858 
859   PetscFunctionBegin;
860   PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1);
861   PetscValidBoolPointer(hasJacPre, 2);
862   ierr = PetscHMapFormGetSize(wf->bdgp0, &n0);CHKERRQ(ierr);
863   ierr = PetscHMapFormGetSize(wf->bdgp1, &n1);CHKERRQ(ierr);
864   ierr = PetscHMapFormGetSize(wf->bdgp2, &n2);CHKERRQ(ierr);
865   ierr = PetscHMapFormGetSize(wf->bdgp3, &n3);CHKERRQ(ierr);
866   *hasJacPre = n0+n1+n2+n3 ? PETSC_TRUE : PETSC_FALSE;
867   PetscFunctionReturn(0);
868 }
869 
870 PetscErrorCode PetscWeakFormGetBdJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g,
871                                                         PetscInt *n0,
872                                                         void (***g0)(PetscInt, PetscInt, PetscInt,
873                                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
874                                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
875                                                                      PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
876                                                         PetscInt *n1,
877                                                         void (***g1)(PetscInt, PetscInt, PetscInt,
878                                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
879                                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
880                                                                      PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
881                                                         PetscInt *n2,
882                                                         void (***g2)(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                                                         PetscInt *n3,
887                                                         void (***g3)(PetscInt, PetscInt, PetscInt,
888                                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
889                                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
890                                                                      PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
891 {
892   PetscInt       find = f*wf->Nf + g;
893   PetscErrorCode ierr;
894 
895   PetscFunctionBegin;
896   ierr = PetscWeakFormGetFunction_Private(wf, wf->bdgp0, label, val, find, n0, (void (***)(void)) g0);CHKERRQ(ierr);
897   ierr = PetscWeakFormGetFunction_Private(wf, wf->bdgp1, label, val, find, n1, (void (***)(void)) g1);CHKERRQ(ierr);
898   ierr = PetscWeakFormGetFunction_Private(wf, wf->bdgp2, label, val, find, n2, (void (***)(void)) g2);CHKERRQ(ierr);
899   ierr = PetscWeakFormGetFunction_Private(wf, wf->bdgp3, label, val, find, n3, (void (***)(void)) g3);CHKERRQ(ierr);
900   PetscFunctionReturn(0);
901 }
902 
903 PetscErrorCode PetscWeakFormAddBdJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g,
904                                                         void (*g0)(PetscInt, PetscInt, PetscInt,
905                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
906                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
907                                                                    PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
908                                                         void (*g1)(PetscInt, PetscInt, PetscInt,
909                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
910                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
911                                                                    PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
912                                                         void (*g2)(PetscInt, PetscInt, PetscInt,
913                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
914                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
915                                                                    PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
916                                                         void (*g3)(PetscInt, PetscInt, PetscInt,
917                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
918                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
919                                                                    PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
920 {
921   PetscInt       find = f*wf->Nf + g;
922   PetscErrorCode ierr;
923 
924   PetscFunctionBegin;
925   ierr = PetscWeakFormAddFunction_Private(wf, wf->bdgp0, label, val, find, (void (*)(void)) g0);CHKERRQ(ierr);
926   ierr = PetscWeakFormAddFunction_Private(wf, wf->bdgp1, label, val, find, (void (*)(void)) g1);CHKERRQ(ierr);
927   ierr = PetscWeakFormAddFunction_Private(wf, wf->bdgp2, label, val, find, (void (*)(void)) g2);CHKERRQ(ierr);
928   ierr = PetscWeakFormAddFunction_Private(wf, wf->bdgp3, label, val, find, (void (*)(void)) g3);CHKERRQ(ierr);
929   PetscFunctionReturn(0);
930 }
931 
932 PetscErrorCode PetscWeakFormSetBdJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g,
933                                                         PetscInt n0,
934                                                         void (**g0)(PetscInt, PetscInt, PetscInt,
935                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
936                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
937                                                                    PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
938                                                         PetscInt n1,
939                                                         void (**g1)(PetscInt, PetscInt, PetscInt,
940                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
941                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
942                                                                    PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
943                                                         PetscInt n2,
944                                                         void (**g2)(PetscInt, PetscInt, PetscInt,
945                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
946                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
947                                                                    PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
948                                                         PetscInt n3,
949                                                         void (**g3)(PetscInt, PetscInt, PetscInt,
950                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
951                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
952                                                                    PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
953 {
954   PetscInt       find = f*wf->Nf + g;
955   PetscErrorCode ierr;
956 
957   PetscFunctionBegin;
958   ierr = PetscWeakFormSetFunction_Private(wf, wf->bdgp0, label, val, find, n0, (void (**)(void)) g0);CHKERRQ(ierr);
959   ierr = PetscWeakFormSetFunction_Private(wf, wf->bdgp1, label, val, find, n1, (void (**)(void)) g1);CHKERRQ(ierr);
960   ierr = PetscWeakFormSetFunction_Private(wf, wf->bdgp2, label, val, find, n2, (void (**)(void)) g2);CHKERRQ(ierr);
961   ierr = PetscWeakFormSetFunction_Private(wf, wf->bdgp3, label, val, find, n3, (void (**)(void)) g3);CHKERRQ(ierr);
962   PetscFunctionReturn(0);
963 }
964 
965 PetscErrorCode PetscWeakFormSetIndexBdJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g,
966                                                         PetscInt i0,
967                                                         void (*g0)(PetscInt, PetscInt, PetscInt,
968                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
969                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
970                                                                    PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
971                                                         PetscInt i1,
972                                                         void (*g1)(PetscInt, PetscInt, PetscInt,
973                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
974                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
975                                                                    PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
976                                                         PetscInt i2,
977                                                         void (*g2)(PetscInt, PetscInt, PetscInt,
978                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
979                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
980                                                                    PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
981                                                         PetscInt i3,
982                                                         void (*g3)(PetscInt, PetscInt, PetscInt,
983                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
984                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
985                                                                    PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
986 {
987   PetscInt       find = f*wf->Nf + g;
988   PetscErrorCode ierr;
989 
990   PetscFunctionBegin;
991   ierr = PetscWeakFormSetIndexFunction_Private(wf, wf->bdgp0, label, val, find, i0, (void (*)(void)) g0);CHKERRQ(ierr);
992   ierr = PetscWeakFormSetIndexFunction_Private(wf, wf->bdgp1, label, val, find, i1, (void (*)(void)) g1);CHKERRQ(ierr);
993   ierr = PetscWeakFormSetIndexFunction_Private(wf, wf->bdgp2, label, val, find, i2, (void (*)(void)) g2);CHKERRQ(ierr);
994   ierr = PetscWeakFormSetIndexFunction_Private(wf, wf->bdgp3, label, val, find, i3, (void (*)(void)) g3);CHKERRQ(ierr);
995   PetscFunctionReturn(0);
996 }
997 
998 PetscErrorCode PetscWeakFormHasDynamicJacobian(PetscWeakForm wf, PetscBool *hasDynJac)
999 {
1000   PetscInt       n0, n1, n2, n3;
1001   PetscErrorCode ierr;
1002 
1003   PetscFunctionBegin;
1004   PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1);
1005   PetscValidBoolPointer(hasDynJac, 2);
1006   ierr = PetscHMapFormGetSize(wf->gt0, &n0);CHKERRQ(ierr);
1007   ierr = PetscHMapFormGetSize(wf->gt1, &n1);CHKERRQ(ierr);
1008   ierr = PetscHMapFormGetSize(wf->gt2, &n2);CHKERRQ(ierr);
1009   ierr = PetscHMapFormGetSize(wf->gt3, &n3);CHKERRQ(ierr);
1010   *hasDynJac = n0+n1+n2+n3 ? PETSC_TRUE : PETSC_FALSE;
1011   PetscFunctionReturn(0);
1012 }
1013 
1014 PetscErrorCode PetscWeakFormGetDynamicJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g,
1015                                         PetscInt *n0,
1016                                         void (***g0)(PetscInt, PetscInt, PetscInt,
1017                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1018                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1019                                                      PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1020                                         PetscInt *n1,
1021                                         void (***g1)(PetscInt, PetscInt, PetscInt,
1022                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1023                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1024                                                      PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1025                                         PetscInt *n2,
1026                                         void (***g2)(PetscInt, PetscInt, PetscInt,
1027                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1028                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1029                                                      PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1030                                         PetscInt *n3,
1031                                         void (***g3)(PetscInt, PetscInt, PetscInt,
1032                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1033                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1034                                                      PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
1035 {
1036   PetscInt       find = f*wf->Nf + g;
1037   PetscErrorCode ierr;
1038 
1039   PetscFunctionBegin;
1040   ierr = PetscWeakFormGetFunction_Private(wf, wf->gt0, label, val, find, n0, (void (***)(void)) g0);CHKERRQ(ierr);
1041   ierr = PetscWeakFormGetFunction_Private(wf, wf->gt1, label, val, find, n1, (void (***)(void)) g1);CHKERRQ(ierr);
1042   ierr = PetscWeakFormGetFunction_Private(wf, wf->gt2, label, val, find, n2, (void (***)(void)) g2);CHKERRQ(ierr);
1043   ierr = PetscWeakFormGetFunction_Private(wf, wf->gt3, label, val, find, n3, (void (***)(void)) g3);CHKERRQ(ierr);
1044   PetscFunctionReturn(0);
1045 }
1046 
1047 PetscErrorCode PetscWeakFormAddDynamicJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g,
1048                                         void (*g0)(PetscInt, PetscInt, PetscInt,
1049                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1050                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1051                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1052                                         void (*g1)(PetscInt, PetscInt, PetscInt,
1053                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1054                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1055                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1056                                         void (*g2)(PetscInt, PetscInt, PetscInt,
1057                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1058                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1059                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1060                                         void (*g3)(PetscInt, PetscInt, PetscInt,
1061                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1062                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1063                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
1064 {
1065   PetscInt       find = f*wf->Nf + g;
1066   PetscErrorCode ierr;
1067 
1068   PetscFunctionBegin;
1069   ierr = PetscWeakFormAddFunction_Private(wf, wf->gt0, label, val, find, (void (*)(void)) g0);CHKERRQ(ierr);
1070   ierr = PetscWeakFormAddFunction_Private(wf, wf->gt1, label, val, find, (void (*)(void)) g1);CHKERRQ(ierr);
1071   ierr = PetscWeakFormAddFunction_Private(wf, wf->gt2, label, val, find, (void (*)(void)) g2);CHKERRQ(ierr);
1072   ierr = PetscWeakFormAddFunction_Private(wf, wf->gt3, label, val, find, (void (*)(void)) g3);CHKERRQ(ierr);
1073   PetscFunctionReturn(0);
1074 }
1075 
1076 PetscErrorCode PetscWeakFormSetDynamicJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g,
1077                                                PetscInt n0,
1078                                                void (**g0)(PetscInt, PetscInt, PetscInt,
1079                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1080                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1081                                                           PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1082                                                PetscInt n1,
1083                                                void (**g1)(PetscInt, PetscInt, PetscInt,
1084                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1085                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1086                                                           PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1087                                                PetscInt n2,
1088                                                void (**g2)(PetscInt, PetscInt, PetscInt,
1089                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1090                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1091                                                           PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1092                                                PetscInt n3,
1093                                                void (**g3)(PetscInt, PetscInt, PetscInt,
1094                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1095                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1096                                                           PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
1097 {
1098   PetscInt       find = f*wf->Nf + g;
1099   PetscErrorCode ierr;
1100 
1101   PetscFunctionBegin;
1102   ierr = PetscWeakFormSetFunction_Private(wf, wf->gt0, label, val, find, n0, (void (**)(void)) g0);CHKERRQ(ierr);
1103   ierr = PetscWeakFormSetFunction_Private(wf, wf->gt1, label, val, find, n1, (void (**)(void)) g1);CHKERRQ(ierr);
1104   ierr = PetscWeakFormSetFunction_Private(wf, wf->gt2, label, val, find, n2, (void (**)(void)) g2);CHKERRQ(ierr);
1105   ierr = PetscWeakFormSetFunction_Private(wf, wf->gt3, label, val, find, n3, (void (**)(void)) g3);CHKERRQ(ierr);
1106   PetscFunctionReturn(0);
1107 }
1108 
1109 PetscErrorCode PetscWeakFormSetIndexDynamicJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g,
1110                                                PetscInt i0,
1111                                                void (*g0)(PetscInt, PetscInt, PetscInt,
1112                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1113                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1114                                                           PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1115                                                PetscInt i1,
1116                                                void (*g1)(PetscInt, PetscInt, PetscInt,
1117                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1118                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1119                                                           PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1120                                                PetscInt i2,
1121                                                void (*g2)(PetscInt, PetscInt, PetscInt,
1122                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1123                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1124                                                           PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1125                                                PetscInt i3,
1126                                                void (*g3)(PetscInt, PetscInt, PetscInt,
1127                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1128                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1129                                                           PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
1130 {
1131   PetscInt       find = f*wf->Nf + g;
1132   PetscErrorCode ierr;
1133 
1134   PetscFunctionBegin;
1135   ierr = PetscWeakFormSetIndexFunction_Private(wf, wf->gt0, label, val, find, i0, (void (*)(void)) g0);CHKERRQ(ierr);
1136   ierr = PetscWeakFormSetIndexFunction_Private(wf, wf->gt1, label, val, find, i1, (void (*)(void)) g1);CHKERRQ(ierr);
1137   ierr = PetscWeakFormSetIndexFunction_Private(wf, wf->gt2, label, val, find, i2, (void (*)(void)) g2);CHKERRQ(ierr);
1138   ierr = PetscWeakFormSetIndexFunction_Private(wf, wf->gt3, label, val, find, i3, (void (*)(void)) g3);CHKERRQ(ierr);
1139   PetscFunctionReturn(0);
1140 }
1141 
1142 PetscErrorCode PetscWeakFormGetRiemannSolver(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt *n,
1143                                              void (***r)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *))
1144 {
1145   PetscErrorCode ierr;
1146 
1147   PetscFunctionBegin;
1148   ierr = PetscWeakFormGetFunction_Private(wf, wf->r, label, val, f, n, (void (***)(void)) r);CHKERRQ(ierr);
1149   PetscFunctionReturn(0);
1150 }
1151 
1152 PetscErrorCode PetscWeakFormSetRiemannSolver(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f,
1153                                              PetscInt n,
1154                                              void (**r)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *))
1155 {
1156   PetscErrorCode ierr;
1157 
1158   PetscFunctionBegin;
1159   ierr = PetscWeakFormSetFunction_Private(wf, wf->r, label, val, f, n, (void (**)(void)) r);CHKERRQ(ierr);
1160   PetscFunctionReturn(0);
1161 }
1162 
1163 PetscErrorCode PetscWeakFormSetIndexRiemannSolver(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f,
1164                                                   PetscInt i,
1165                                                   void (*r)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *))
1166 {
1167   PetscErrorCode ierr;
1168 
1169   PetscFunctionBegin;
1170   ierr = PetscWeakFormSetIndexFunction_Private(wf, wf->r, label, val, f, i, (void (*)(void)) r);CHKERRQ(ierr);
1171   PetscFunctionReturn(0);
1172 }
1173 
1174 /*@
1175   PetscWeakFormGetNumFields - Returns the number of fields
1176 
1177   Not collective
1178 
1179   Input Parameter:
1180 . wf - The PetscWeakForm object
1181 
1182   Output Parameter:
1183 . Nf - The nubmer of fields
1184 
1185   Level: beginner
1186 
1187 .seealso: PetscWeakFormSetNumFields(), PetscWeakFormCreate()
1188 @*/
1189 PetscErrorCode PetscWeakFormGetNumFields(PetscWeakForm wf, PetscInt *Nf)
1190 {
1191   PetscFunctionBegin;
1192   PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1);
1193   PetscValidPointer(Nf, 2);
1194   *Nf = wf->Nf;
1195   PetscFunctionReturn(0);
1196 }
1197 
1198 /*@
1199   PetscWeakFormSetNumFields - Sets the number of fields
1200 
1201   Not collective
1202 
1203   Input Parameters:
1204 + wf - The PetscWeakForm object
1205 - Nf - The number of fields
1206 
1207   Level: beginner
1208 
1209 .seealso: PetscWeakFormGetNumFields(), PetscWeakFormCreate()
1210 @*/
1211 PetscErrorCode PetscWeakFormSetNumFields(PetscWeakForm wf, PetscInt Nf)
1212 {
1213   PetscFunctionBegin;
1214   PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1);
1215   wf->Nf = Nf;
1216   PetscFunctionReturn(0);
1217 }
1218 
1219 /*@
1220   PetscWeakFormDestroy - Destroys a PetscWeakForm object
1221 
1222   Collective on wf
1223 
1224   Input Parameter:
1225 . wf - the PetscWeakForm object to destroy
1226 
1227   Level: developer
1228 
1229 .seealso PetscWeakFormCreate(), PetscWeakFormView()
1230 @*/
1231 PetscErrorCode PetscWeakFormDestroy(PetscWeakForm *wf)
1232 {
1233   PetscErrorCode ierr;
1234 
1235   PetscFunctionBegin;
1236   if (!*wf) PetscFunctionReturn(0);
1237   PetscValidHeaderSpecific((*wf), PETSCWEAKFORM_CLASSID, 1);
1238 
1239   if (--((PetscObject)(*wf))->refct > 0) {*wf = NULL; PetscFunctionReturn(0);}
1240   ((PetscObject) (*wf))->refct = 0;
1241   ierr = PetscChunkBufferDestroy(&(*wf)->funcs);CHKERRQ(ierr);
1242   ierr = PetscHMapFormDestroy(&(*wf)->obj);CHKERRQ(ierr);
1243   ierr = PetscHMapFormDestroy(&(*wf)->f0);CHKERRQ(ierr);
1244   ierr = PetscHMapFormDestroy(&(*wf)->f1);CHKERRQ(ierr);
1245   ierr = PetscHMapFormDestroy(&(*wf)->g0);CHKERRQ(ierr);
1246   ierr = PetscHMapFormDestroy(&(*wf)->g1);CHKERRQ(ierr);
1247   ierr = PetscHMapFormDestroy(&(*wf)->g2);CHKERRQ(ierr);
1248   ierr = PetscHMapFormDestroy(&(*wf)->g3);CHKERRQ(ierr);
1249   ierr = PetscHMapFormDestroy(&(*wf)->gp0);CHKERRQ(ierr);
1250   ierr = PetscHMapFormDestroy(&(*wf)->gp1);CHKERRQ(ierr);
1251   ierr = PetscHMapFormDestroy(&(*wf)->gp2);CHKERRQ(ierr);
1252   ierr = PetscHMapFormDestroy(&(*wf)->gp3);CHKERRQ(ierr);
1253   ierr = PetscHMapFormDestroy(&(*wf)->gt0);CHKERRQ(ierr);
1254   ierr = PetscHMapFormDestroy(&(*wf)->gt1);CHKERRQ(ierr);
1255   ierr = PetscHMapFormDestroy(&(*wf)->gt2);CHKERRQ(ierr);
1256   ierr = PetscHMapFormDestroy(&(*wf)->gt3);CHKERRQ(ierr);
1257   ierr = PetscHMapFormDestroy(&(*wf)->bdf0);CHKERRQ(ierr);
1258   ierr = PetscHMapFormDestroy(&(*wf)->bdf1);CHKERRQ(ierr);
1259   ierr = PetscHMapFormDestroy(&(*wf)->bdg0);CHKERRQ(ierr);
1260   ierr = PetscHMapFormDestroy(&(*wf)->bdg1);CHKERRQ(ierr);
1261   ierr = PetscHMapFormDestroy(&(*wf)->bdg2);CHKERRQ(ierr);
1262   ierr = PetscHMapFormDestroy(&(*wf)->bdg3);CHKERRQ(ierr);
1263   ierr = PetscHMapFormDestroy(&(*wf)->bdgp0);CHKERRQ(ierr);
1264   ierr = PetscHMapFormDestroy(&(*wf)->bdgp1);CHKERRQ(ierr);
1265   ierr = PetscHMapFormDestroy(&(*wf)->bdgp2);CHKERRQ(ierr);
1266   ierr = PetscHMapFormDestroy(&(*wf)->bdgp3);CHKERRQ(ierr);
1267   ierr = PetscHMapFormDestroy(&(*wf)->r);CHKERRQ(ierr);
1268   ierr = PetscHeaderDestroy(wf);CHKERRQ(ierr);
1269   PetscFunctionReturn(0);
1270 }
1271 
1272 static PetscErrorCode PetscWeakFormViewTable_Ascii(PetscWeakForm wf, PetscViewer viewer, const char tableName[], PetscHMapForm map)
1273 {
1274   PetscInt       Nk, k;
1275   PetscErrorCode ierr;
1276 
1277   PetscFunctionBegin;
1278   ierr = PetscHMapFormGetSize(map, &Nk);CHKERRQ(ierr);
1279   if (Nk) {
1280     PetscHashFormKey *keys;
1281     void           (**funcs)(void);
1282     const char       *name;
1283     PetscInt          off = 0, n, i;
1284 
1285     ierr = PetscMalloc1(Nk, &keys);CHKERRQ(ierr);
1286     ierr = PetscHMapFormGetKeys(map, &off, keys);CHKERRQ(ierr);
1287     ierr = PetscViewerASCIIPrintf(viewer, "%s\n", tableName);CHKERRQ(ierr);
1288     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1289     for (k = 0; k < Nk; ++k) {
1290       if (keys[k].label) {ierr = PetscObjectGetName((PetscObject) keys[k].label, &name);CHKERRQ(ierr);}
1291       ierr = PetscViewerASCIIPrintf(viewer, "Key (%s, %D, %D) ", keys[k].label ? name : "None", keys[k].value, keys[k].field);CHKERRQ(ierr);
1292       ierr = PetscViewerASCIIUseTabs(viewer, PETSC_FALSE);CHKERRQ(ierr);
1293       ierr = PetscWeakFormGetFunction_Private(wf, map, keys[k].label, keys[k].value, keys[k].field, &n, &funcs);CHKERRQ(ierr);
1294       for (i = 0; i < n; ++i) {
1295         if (i > 0) {ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr);}
1296         ierr = PetscDLAddr(funcs[i], &name);CHKERRQ(ierr);
1297         if (name) {ierr = PetscViewerASCIIPrintf(viewer, "%s", name);CHKERRQ(ierr);}
1298         else      {ierr = PetscViewerASCIIPrintf(viewer, "%p", funcs[i]);CHKERRQ(ierr);}
1299       }
1300       ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr);
1301       ierr = PetscViewerASCIIUseTabs(viewer, PETSC_TRUE);CHKERRQ(ierr);
1302     }
1303     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1304     ierr = PetscFree(keys);CHKERRQ(ierr);
1305   }
1306   PetscFunctionReturn(0);
1307 }
1308 
1309 static PetscErrorCode PetscWeakFormView_Ascii(PetscWeakForm wf, PetscViewer viewer)
1310 {
1311   PetscViewerFormat format;
1312   PetscErrorCode    ierr;
1313 
1314   PetscFunctionBegin;
1315   ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
1316   ierr = PetscViewerASCIIPrintf(viewer, "Weak Form System with %d fields\n", wf->Nf);CHKERRQ(ierr);
1317   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1318   ierr = PetscWeakFormViewTable_Ascii(wf, viewer, "Objective", wf->obj);CHKERRQ(ierr);
1319   ierr = PetscWeakFormViewTable_Ascii(wf, viewer, "Residual f0", wf->f0);CHKERRQ(ierr);
1320   ierr = PetscWeakFormViewTable_Ascii(wf, viewer, "Residual f1", wf->f1);CHKERRQ(ierr);
1321   ierr = PetscWeakFormViewTable_Ascii(wf, viewer, "Boudnary Residual f0", wf->bdf0);CHKERRQ(ierr);
1322   ierr = PetscWeakFormViewTable_Ascii(wf, viewer, "Boundary Residual f1", wf->bdf1);CHKERRQ(ierr);
1323   ierr = PetscWeakFormViewTable_Ascii(wf, viewer, "Jacobian g0", wf->g0);CHKERRQ(ierr);
1324   ierr = PetscWeakFormViewTable_Ascii(wf, viewer, "Jacobian g1", wf->g1);CHKERRQ(ierr);
1325   ierr = PetscWeakFormViewTable_Ascii(wf, viewer, "Jacobian g2", wf->g2);CHKERRQ(ierr);
1326   ierr = PetscWeakFormViewTable_Ascii(wf, viewer, "Jacobian g3", wf->g3);CHKERRQ(ierr);
1327   ierr = PetscWeakFormViewTable_Ascii(wf, viewer, "Jacobian Preconditioner g0", wf->gp0);CHKERRQ(ierr);
1328   ierr = PetscWeakFormViewTable_Ascii(wf, viewer, "Jacobian Preconditioner g1", wf->gp1);CHKERRQ(ierr);
1329   ierr = PetscWeakFormViewTable_Ascii(wf, viewer, "Jacobian Preconditioner g2", wf->gp2);CHKERRQ(ierr);
1330   ierr = PetscWeakFormViewTable_Ascii(wf, viewer, "Jacobian Preconditioner g3", wf->gp3);CHKERRQ(ierr);
1331   ierr = PetscWeakFormViewTable_Ascii(wf, viewer, "Dynamic Jacobian g0", wf->gt0);CHKERRQ(ierr);
1332   ierr = PetscWeakFormViewTable_Ascii(wf, viewer, "Dynamic Jacobian g1", wf->gt1);CHKERRQ(ierr);
1333   ierr = PetscWeakFormViewTable_Ascii(wf, viewer, "Dynamic Jacobian g2", wf->gt2);CHKERRQ(ierr);
1334   ierr = PetscWeakFormViewTable_Ascii(wf, viewer, "Dynamic Jacobian g3", wf->gt3);CHKERRQ(ierr);
1335   ierr = PetscWeakFormViewTable_Ascii(wf, viewer, "Boundary Jacobian g0", wf->bdg0);CHKERRQ(ierr);
1336   ierr = PetscWeakFormViewTable_Ascii(wf, viewer, "Boundary Jacobian g1", wf->bdg1);CHKERRQ(ierr);
1337   ierr = PetscWeakFormViewTable_Ascii(wf, viewer, "Boundary Jacobian g2", wf->bdg2);CHKERRQ(ierr);
1338   ierr = PetscWeakFormViewTable_Ascii(wf, viewer, "Boundary Jacobian g3", wf->bdg3);CHKERRQ(ierr);
1339   ierr = PetscWeakFormViewTable_Ascii(wf, viewer, "Boundary Jacobian Preconditioner g0", wf->bdgp0);CHKERRQ(ierr);
1340   ierr = PetscWeakFormViewTable_Ascii(wf, viewer, "Boundary Jacobian Preconditioner g1", wf->bdgp1);CHKERRQ(ierr);
1341   ierr = PetscWeakFormViewTable_Ascii(wf, viewer, "Boundary Jacobian Preconditioner g2", wf->bdgp2);CHKERRQ(ierr);
1342   ierr = PetscWeakFormViewTable_Ascii(wf, viewer, "Boundary Jacobian Preconditioner g3", wf->bdgp3);CHKERRQ(ierr);
1343   ierr = PetscWeakFormViewTable_Ascii(wf, viewer, "Riemann Solver", wf->r);CHKERRQ(ierr);
1344   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1345   PetscFunctionReturn(0);
1346 }
1347 
1348 /*@C
1349   PetscWeakFormView - Views a PetscWeakForm
1350 
1351   Collective on wf
1352 
1353   Input Parameter:
1354 + wf - the PetscWeakForm object to view
1355 - v  - the viewer
1356 
1357   Level: developer
1358 
1359 .seealso PetscWeakFormDestroy(), PetscWeakFormCreate()
1360 @*/
1361 PetscErrorCode PetscWeakFormView(PetscWeakForm wf, PetscViewer v)
1362 {
1363   PetscBool      iascii;
1364   PetscErrorCode ierr;
1365 
1366   PetscFunctionBegin;
1367   PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1);
1368   if (!v) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) wf), &v);CHKERRQ(ierr);}
1369   else    {PetscValidHeaderSpecific(v, PETSC_VIEWER_CLASSID, 2);}
1370   ierr = PetscObjectTypeCompare((PetscObject) v, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
1371   if (iascii) {ierr = PetscWeakFormView_Ascii(wf, v);CHKERRQ(ierr);}
1372   if (wf->ops->view) {ierr = (*wf->ops->view)(wf, v);CHKERRQ(ierr);}
1373   PetscFunctionReturn(0);
1374 }
1375 
1376 /*@
1377   PetscWeakFormCreate - Creates an empty PetscWeakForm object.
1378 
1379   Collective
1380 
1381   Input Parameter:
1382 . comm - The communicator for the PetscWeakForm object
1383 
1384   Output Parameter:
1385 . wf - The PetscWeakForm object
1386 
1387   Level: beginner
1388 
1389 .seealso: PetscDS, PetscWeakFormDestroy()
1390 @*/
1391 PetscErrorCode PetscWeakFormCreate(MPI_Comm comm, PetscWeakForm *wf)
1392 {
1393   PetscWeakForm  p;
1394   PetscErrorCode ierr;
1395 
1396   PetscFunctionBegin;
1397   PetscValidPointer(wf, 2);
1398   *wf  = NULL;
1399   ierr = PetscDSInitializePackage();CHKERRQ(ierr);
1400 
1401   ierr = PetscHeaderCreate(p, PETSCWEAKFORM_CLASSID, "PetscWeakForm", "Weak Form System", "PetscWeakForm", comm, PetscWeakFormDestroy, PetscWeakFormView);CHKERRQ(ierr);
1402 
1403   p->Nf = 0;
1404   ierr = PetscChunkBufferCreate(sizeof(&PetscWeakFormCreate), 2, &p->funcs);CHKERRQ(ierr);
1405   ierr = PetscHMapFormCreate(&p->obj);CHKERRQ(ierr);
1406   ierr = PetscHMapFormCreate(&p->f0);CHKERRQ(ierr);
1407   ierr = PetscHMapFormCreate(&p->f1);CHKERRQ(ierr);
1408   ierr = PetscHMapFormCreate(&p->g0);CHKERRQ(ierr);
1409   ierr = PetscHMapFormCreate(&p->g1);CHKERRQ(ierr);
1410   ierr = PetscHMapFormCreate(&p->g2);CHKERRQ(ierr);
1411   ierr = PetscHMapFormCreate(&p->g3);CHKERRQ(ierr);
1412   ierr = PetscHMapFormCreate(&p->gp0);CHKERRQ(ierr);
1413   ierr = PetscHMapFormCreate(&p->gp1);CHKERRQ(ierr);
1414   ierr = PetscHMapFormCreate(&p->gp2);CHKERRQ(ierr);
1415   ierr = PetscHMapFormCreate(&p->gp3);CHKERRQ(ierr);
1416   ierr = PetscHMapFormCreate(&p->gt0);CHKERRQ(ierr);
1417   ierr = PetscHMapFormCreate(&p->gt1);CHKERRQ(ierr);
1418   ierr = PetscHMapFormCreate(&p->gt2);CHKERRQ(ierr);
1419   ierr = PetscHMapFormCreate(&p->gt3);CHKERRQ(ierr);
1420   ierr = PetscHMapFormCreate(&p->bdf0);CHKERRQ(ierr);
1421   ierr = PetscHMapFormCreate(&p->bdf1);CHKERRQ(ierr);
1422   ierr = PetscHMapFormCreate(&p->bdg0);CHKERRQ(ierr);
1423   ierr = PetscHMapFormCreate(&p->bdg1);CHKERRQ(ierr);
1424   ierr = PetscHMapFormCreate(&p->bdg2);CHKERRQ(ierr);
1425   ierr = PetscHMapFormCreate(&p->bdg3);CHKERRQ(ierr);
1426   ierr = PetscHMapFormCreate(&p->bdgp0);CHKERRQ(ierr);
1427   ierr = PetscHMapFormCreate(&p->bdgp1);CHKERRQ(ierr);
1428   ierr = PetscHMapFormCreate(&p->bdgp2);CHKERRQ(ierr);
1429   ierr = PetscHMapFormCreate(&p->bdgp3);CHKERRQ(ierr);
1430   ierr = PetscHMapFormCreate(&p->r);CHKERRQ(ierr);
1431   *wf = p;
1432   PetscFunctionReturn(0);
1433 }
1434