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