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