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