1 #include <petsc/private/petscdsimpl.h> /*I "petscds.h" I*/
2
3 PetscClassId PETSCWEAKFORM_CLASSID = 0;
4
5 const char *const PetscWeakFormKinds[] = {"objective", "residual_f0", "residual_f1", "jacobian_g0", "jacobian_g1", "jacobian_g2", "jacobian_g3", "jacobian_preconditioner_g0", "jacobian_preconditioner_g1", "jacobian_preconditioner_g2", "jacobian_preconditioner_g3", "dynamic_jacobian_g0", "dynamic_jacobian_g1", "dynamic_jacobian_g2", "dynamic_jacobian_g3", "boundary_residual_f0", "boundary_residual_f1", "boundary_jacobian_g0", "boundary_jacobian_g1", "boundary_jacobian_g2", "boundary_jacobian_g3", "boundary_jacobian_preconditioner_g0", "boundary_jacobian_preconditioner_g1", "boundary_jacobian_preconditioner_g2", "boundary_jacobian_preconditioner_g3", "riemann_solver", "PetscWeakFormKind", "PETSC_WF_", NULL};
6
PetscChunkBufferCreate(size_t unitbytes,PetscCount expected,PetscChunkBuffer * buffer[])7 static PetscErrorCode PetscChunkBufferCreate(size_t unitbytes, PetscCount expected, PetscChunkBuffer *buffer[])
8 {
9 PetscFunctionBegin;
10 PetscCall(PetscNew(buffer));
11 PetscCall(PetscCalloc1(expected * unitbytes, &(*buffer)->array));
12 (*buffer)->size = expected;
13 (*buffer)->unitbytes = unitbytes;
14 (*buffer)->alloc = expected * unitbytes;
15 PetscFunctionReturn(PETSC_SUCCESS);
16 }
17
PetscChunkBufferDuplicate(PetscChunkBuffer * buffer,PetscChunkBuffer * bufferNew[])18 static PetscErrorCode PetscChunkBufferDuplicate(PetscChunkBuffer *buffer, PetscChunkBuffer *bufferNew[])
19 {
20 PetscFunctionBegin;
21 PetscCall(PetscNew(bufferNew));
22 PetscCall(PetscCalloc1(buffer->size * buffer->unitbytes, &(*bufferNew)->array));
23 PetscCall(PetscMemcpy((*bufferNew)->array, buffer->array, buffer->size * buffer->unitbytes));
24 (*bufferNew)->size = buffer->size;
25 (*bufferNew)->unitbytes = buffer->unitbytes;
26 (*bufferNew)->alloc = buffer->size * buffer->unitbytes;
27 PetscFunctionReturn(PETSC_SUCCESS);
28 }
29
PetscChunkBufferDestroy(PetscChunkBuffer ** buffer)30 static PetscErrorCode PetscChunkBufferDestroy(PetscChunkBuffer **buffer)
31 {
32 PetscFunctionBegin;
33 PetscCall(PetscFree((*buffer)->array));
34 PetscCall(PetscFree(*buffer));
35 PetscFunctionReturn(PETSC_SUCCESS);
36 }
37
PetscChunkBufferCreateChunk(PetscChunkBuffer * buffer,PetscCount size,PetscChunk * chunk)38 static PetscErrorCode PetscChunkBufferCreateChunk(PetscChunkBuffer *buffer, PetscCount size, PetscChunk *chunk)
39 {
40 PetscFunctionBegin;
41 if ((buffer->size + size) * buffer->unitbytes > buffer->alloc) {
42 char *tmp;
43
44 if (!buffer->alloc) buffer->alloc = (buffer->size + size) * buffer->unitbytes;
45 while ((buffer->size + size) * buffer->unitbytes > buffer->alloc) buffer->alloc *= 2;
46 PetscCall(PetscMalloc(buffer->alloc, &tmp));
47 PetscCall(PetscMemcpy(tmp, buffer->array, buffer->size * buffer->unitbytes));
48 PetscCall(PetscFree(buffer->array));
49 buffer->array = tmp;
50 }
51 chunk->start = buffer->size * buffer->unitbytes;
52 chunk->size = size;
53 chunk->reserved = size;
54 buffer->size += size;
55 PetscFunctionReturn(PETSC_SUCCESS);
56 }
57
PetscChunkBufferEnlargeChunk(PetscChunkBuffer * buffer,PetscCount size,PetscChunk * chunk)58 static PetscErrorCode PetscChunkBufferEnlargeChunk(PetscChunkBuffer *buffer, PetscCount size, PetscChunk *chunk)
59 {
60 size_t siz = size;
61
62 PetscFunctionBegin;
63 if (chunk->size + size > chunk->reserved) {
64 PetscChunk newchunk;
65 PetscCount reserved = chunk->size;
66
67 /* TODO Here if we had a chunk list, we could update them all to reclaim unused space */
68 while (reserved < chunk->size + size) reserved *= 2;
69 PetscCall(PetscChunkBufferCreateChunk(buffer, (size_t)reserved, &newchunk));
70 newchunk.size = chunk->size + size;
71 PetscCall(PetscMemcpy(&buffer->array[newchunk.start], &buffer->array[chunk->start], chunk->size * buffer->unitbytes));
72 *chunk = newchunk;
73 } else {
74 chunk->size += siz;
75 }
76 PetscFunctionReturn(PETSC_SUCCESS);
77 }
78
79 /*@C
80 PetscFormKeySort - Sorts an array of `PetscFormKey` in place in increasing order.
81
82 Not Collective
83
84 Input Parameters:
85 + n - number of values
86 - arr - array of `PetscFormKey`
87
88 Level: intermediate
89
90 .seealso: `PetscFormKey`, `PetscIntSortSemiOrdered()`, `PetscSortInt()`
91 @*/
PetscFormKeySort(PetscInt n,PetscFormKey arr[])92 PetscErrorCode PetscFormKeySort(PetscInt n, PetscFormKey arr[])
93 {
94 PetscFunctionBegin;
95 if (n <= 1) PetscFunctionReturn(PETSC_SUCCESS);
96 PetscAssertPointer(arr, 2);
97 PetscCall(PetscTimSort(n, arr, sizeof(PetscFormKey), Compare_PetscFormKey_Private, NULL));
98 PetscFunctionReturn(PETSC_SUCCESS);
99 }
100
PetscWeakFormGetFunction_Private(PetscWeakForm wf,PetscHMapForm ht,DMLabel label,PetscInt value,PetscInt f,PetscInt part,PetscInt * n,void (*** func)(void))101 static PetscErrorCode PetscWeakFormGetFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt part, PetscInt *n, void (***func)(void))
102 {
103 PetscFormKey key;
104 PetscChunk chunk;
105
106 PetscFunctionBegin;
107 key.label = label;
108 key.value = value;
109 key.field = f;
110 key.part = part;
111 PetscCall(PetscHMapFormGet(ht, key, &chunk));
112 if (chunk.size < 0) {
113 *n = 0;
114 *func = NULL;
115 } else {
116 PetscCall(PetscIntCast(chunk.size, n));
117 *func = (PetscVoidFn **)&wf->funcs->array[chunk.start];
118 }
119 PetscFunctionReturn(PETSC_SUCCESS);
120 }
121
122 /* A NULL argument for func causes this to clear the key */
PetscWeakFormSetFunction_Private(PetscWeakForm wf,PetscHMapForm ht,DMLabel label,PetscInt value,PetscInt f,PetscInt part,PetscInt n,PetscVoidFn ** func)123 static PetscErrorCode PetscWeakFormSetFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt part, PetscInt n, PetscVoidFn **func)
124 {
125 PetscFormKey key;
126 PetscChunk chunk;
127 PetscInt i;
128
129 PetscFunctionBegin;
130 key.label = label;
131 key.value = value;
132 key.field = f;
133 key.part = part;
134 if (!func) {
135 PetscCall(PetscHMapFormDel(ht, key));
136 PetscFunctionReturn(PETSC_SUCCESS);
137 } else PetscCall(PetscHMapFormGet(ht, key, &chunk));
138 if (chunk.size < 0) {
139 PetscCall(PetscChunkBufferCreateChunk(wf->funcs, n, &chunk));
140 PetscCall(PetscHMapFormSet(ht, key, chunk));
141 } else if (chunk.size <= n) {
142 PetscCall(PetscChunkBufferEnlargeChunk(wf->funcs, n - chunk.size, &chunk));
143 PetscCall(PetscHMapFormSet(ht, key, chunk));
144 }
145 for (i = 0; i < n; ++i) ((PetscVoidFn **)&wf->funcs->array[chunk.start])[i] = func[i];
146 PetscFunctionReturn(PETSC_SUCCESS);
147 }
148
PetscWeakFormAddFunction_Private(PetscWeakForm wf,PetscHMapForm ht,DMLabel label,PetscInt value,PetscInt f,PetscInt part,PetscVoidFn * func)149 static PetscErrorCode PetscWeakFormAddFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt part, PetscVoidFn *func)
150 {
151 PetscFormKey key;
152 PetscChunk chunk;
153
154 PetscFunctionBegin;
155 if (!func) PetscFunctionReturn(PETSC_SUCCESS);
156 key.label = label;
157 key.value = value;
158 key.field = f;
159 key.part = part;
160 PetscCall(PetscHMapFormGet(ht, key, &chunk));
161 if (chunk.size < 0) {
162 PetscCall(PetscChunkBufferCreateChunk(wf->funcs, 1, &chunk));
163 PetscCall(PetscHMapFormSet(ht, key, chunk));
164 ((PetscVoidFn **)&wf->funcs->array[chunk.start])[0] = func;
165 } else {
166 PetscCall(PetscChunkBufferEnlargeChunk(wf->funcs, 1, &chunk));
167 PetscCall(PetscHMapFormSet(ht, key, chunk));
168 ((PetscVoidFn **)&wf->funcs->array[chunk.start])[chunk.size - 1] = func;
169 }
170 PetscFunctionReturn(PETSC_SUCCESS);
171 }
172
PetscWeakFormGetIndexFunction_Private(PetscWeakForm wf,PetscHMapForm ht,DMLabel label,PetscInt value,PetscInt f,PetscInt part,PetscInt ind,PetscVoidFn ** func)173 static PetscErrorCode PetscWeakFormGetIndexFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt part, PetscInt ind, PetscVoidFn **func)
174 {
175 PetscFormKey key;
176 PetscChunk chunk;
177
178 PetscFunctionBegin;
179 key.label = label;
180 key.value = value;
181 key.field = f;
182 key.part = part;
183 PetscCall(PetscHMapFormGet(ht, key, &chunk));
184 if (chunk.size < 0) {
185 *func = NULL;
186 } else {
187 PetscCheck(ind < chunk.size, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Index %" PetscInt_FMT " not in [0, %" PetscCount_FMT ")", ind, chunk.size);
188 *func = ((PetscVoidFn **)&wf->funcs->array[chunk.start])[ind];
189 }
190 PetscFunctionReturn(PETSC_SUCCESS);
191 }
192
193 /* Ignore a NULL func */
PetscWeakFormSetIndexFunction_Private(PetscWeakForm wf,PetscHMapForm ht,DMLabel label,PetscInt value,PetscInt f,PetscInt part,PetscInt ind,PetscVoidFn * func)194 static PetscErrorCode PetscWeakFormSetIndexFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt part, PetscInt ind, PetscVoidFn *func)
195 {
196 PetscFormKey key;
197 PetscChunk chunk;
198
199 PetscFunctionBegin;
200 if (!func) PetscFunctionReturn(PETSC_SUCCESS);
201 key.label = label;
202 key.value = value;
203 key.field = f;
204 key.part = part;
205 PetscCall(PetscHMapFormGet(ht, key, &chunk));
206 if (chunk.size < 0) {
207 PetscCall(PetscChunkBufferCreateChunk(wf->funcs, ind + 1, &chunk));
208 PetscCall(PetscHMapFormSet(ht, key, chunk));
209 } else if (chunk.size <= ind) {
210 PetscCall(PetscChunkBufferEnlargeChunk(wf->funcs, ind - chunk.size + 1, &chunk));
211 PetscCall(PetscHMapFormSet(ht, key, chunk));
212 }
213 ((PetscVoidFn **)&wf->funcs->array[chunk.start])[ind] = func;
214 PetscFunctionReturn(PETSC_SUCCESS);
215 }
216
PetscWeakFormClearIndexFunction_Private(PetscWeakForm wf,PetscHMapForm ht,DMLabel label,PetscInt value,PetscInt f,PetscInt part,PetscInt ind)217 static PetscErrorCode PetscWeakFormClearIndexFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt part, PetscInt ind)
218 {
219 PetscFormKey key;
220 PetscChunk chunk;
221
222 PetscFunctionBegin;
223 key.label = label;
224 key.value = value;
225 key.field = f;
226 key.part = part;
227 PetscCall(PetscHMapFormGet(ht, key, &chunk));
228 if (chunk.size < 0) {
229 PetscFunctionReturn(PETSC_SUCCESS);
230 } else if (!ind && chunk.size == 1) {
231 PetscCall(PetscHMapFormDel(ht, key));
232 PetscFunctionReturn(PETSC_SUCCESS);
233 } else if (chunk.size <= ind) {
234 PetscFunctionReturn(PETSC_SUCCESS);
235 }
236 ((PetscVoidFn **)&wf->funcs->array[chunk.start])[ind] = NULL;
237 PetscFunctionReturn(PETSC_SUCCESS);
238 }
239
240 /*@
241 PetscWeakFormCopy - Copy the pointwise functions to another `PetscWeakForm`
242
243 Not Collective
244
245 Input Parameter:
246 . wf - The original `PetscWeakForm`
247
248 Output Parameter:
249 . wfNew - The copy of the `PetscWeakForm`
250
251 Level: intermediate
252
253 .seealso: `PetscWeakForm`, `PetscWeakFormCreate()`, `PetscWeakFormDestroy()`
254 @*/
PetscWeakFormCopy(PetscWeakForm wf,PetscWeakForm wfNew)255 PetscErrorCode PetscWeakFormCopy(PetscWeakForm wf, PetscWeakForm wfNew)
256 {
257 PetscInt f;
258
259 PetscFunctionBegin;
260 wfNew->Nf = wf->Nf;
261 PetscCall(PetscChunkBufferDestroy(&wfNew->funcs));
262 PetscCall(PetscChunkBufferDuplicate(wf->funcs, &wfNew->funcs));
263 for (f = 0; f < PETSC_NUM_WF; ++f) {
264 PetscCall(PetscHMapFormDestroy(&wfNew->form[f]));
265 PetscCall(PetscHMapFormDuplicate(wf->form[f], &wfNew->form[f]));
266 }
267 PetscFunctionReturn(PETSC_SUCCESS);
268 }
269
270 /*@
271 PetscWeakFormClear - Clear all functions from the `PetscWeakForm`
272
273 Not Collective
274
275 Input Parameter:
276 . wf - The original `PetscWeakForm`
277
278 Level: intermediate
279
280 .seealso: `PetscWeakForm`, `PetscWeakFormCopy()`, `PetscWeakFormCreate()`, `PetscWeakFormDestroy()`
281 @*/
PetscWeakFormClear(PetscWeakForm wf)282 PetscErrorCode PetscWeakFormClear(PetscWeakForm wf)
283 {
284 PetscInt f;
285
286 PetscFunctionBegin;
287 for (f = 0; f < PETSC_NUM_WF; ++f) PetscCall(PetscHMapFormClear(wf->form[f]));
288 PetscFunctionReturn(PETSC_SUCCESS);
289 }
290
PetscWeakFormRewriteKeys_Internal(PetscWeakForm wf,PetscHMapForm hmap,DMLabel label,PetscInt Nv,const PetscInt values[])291 static PetscErrorCode PetscWeakFormRewriteKeys_Internal(PetscWeakForm wf, PetscHMapForm hmap, DMLabel label, PetscInt Nv, const PetscInt values[])
292 {
293 PetscFormKey *keys;
294 PetscVoidFn **tmpfuncs;
295 PetscInt n, off = 0, maxNf = 0;
296
297 PetscFunctionBegin;
298 PetscCall(PetscHMapFormGetSize(hmap, &n));
299 PetscCall(PetscMalloc1(n, &keys));
300 PetscCall(PetscHMapFormGetKeys(hmap, &off, keys));
301 // Need to make a copy since SetFunction() can invalidate the storage
302 for (PetscInt i = 0; i < n; ++i) {
303 if (keys[i].label == label) {
304 PetscVoidFn **funcs;
305 PetscInt Nf;
306
307 PetscCall(PetscWeakFormGetFunction_Private(wf, hmap, keys[i].label, keys[i].value, keys[i].field, keys[i].part, &Nf, &funcs));
308 maxNf = PetscMax(maxNf, Nf);
309 }
310 }
311 PetscCall(PetscMalloc1(maxNf, &tmpfuncs));
312 for (PetscInt i = 0; i < n; ++i) {
313 if (keys[i].label == label) {
314 PetscBool clear = PETSC_TRUE;
315 PetscVoidFn **funcs;
316 PetscInt Nf;
317
318 PetscCall(PetscWeakFormGetFunction_Private(wf, hmap, keys[i].label, keys[i].value, keys[i].field, keys[i].part, &Nf, &funcs));
319 for (PetscInt f = 0; f < Nf; ++f) tmpfuncs[f] = funcs[f];
320 for (PetscInt v = 0; v < Nv; ++v) {
321 PetscCall(PetscWeakFormSetFunction_Private(wf, hmap, keys[i].label, values[v], keys[i].field, keys[i].part, Nf, tmpfuncs));
322 if (values[v] == keys[i].value) clear = PETSC_FALSE;
323 }
324 if (clear) PetscCall(PetscWeakFormSetFunction_Private(wf, hmap, keys[i].label, keys[i].value, keys[i].field, keys[i].part, 0, NULL));
325 }
326 }
327 PetscCall(PetscFree(tmpfuncs));
328 PetscCall(PetscFree(keys));
329 PetscFunctionReturn(PETSC_SUCCESS);
330 }
331
332 /*@
333 PetscWeakFormRewriteKeys - Change any key on the given label to use the new set of label values
334
335 Not Collective
336
337 Input Parameters:
338 + wf - The original `PetscWeakForm`
339 . label - The label to change keys for
340 . Nv - The number of new label values
341 - values - The set of new values to relabel keys with
342
343 Level: intermediate
344
345 Note:
346 This is used internally when boundary label values are specified from the command line.
347
348 .seealso: `PetscWeakForm`, `DMLabel`, `PetscWeakFormReplaceLabel()`, `PetscWeakFormCreate()`, `PetscWeakFormDestroy()`
349 @*/
PetscWeakFormRewriteKeys(PetscWeakForm wf,DMLabel label,PetscInt Nv,const PetscInt values[])350 PetscErrorCode PetscWeakFormRewriteKeys(PetscWeakForm wf, DMLabel label, PetscInt Nv, const PetscInt values[])
351 {
352 PetscInt f;
353
354 PetscFunctionBegin;
355 for (f = 0; f < PETSC_NUM_WF; ++f) PetscCall(PetscWeakFormRewriteKeys_Internal(wf, wf->form[f], label, Nv, values));
356 PetscFunctionReturn(PETSC_SUCCESS);
357 }
358
PetscWeakFormReplaceLabel_Internal(PetscWeakForm wf,PetscHMapForm hmap,DMLabel label)359 static PetscErrorCode PetscWeakFormReplaceLabel_Internal(PetscWeakForm wf, PetscHMapForm hmap, DMLabel label)
360 {
361 PetscFormKey *keys;
362 PetscInt n, i, off = 0, maxFuncs = 0;
363 PetscVoidFn **tmpf;
364 const char *name = NULL;
365
366 PetscFunctionBegin;
367 if (label) PetscCall(PetscObjectGetName((PetscObject)label, &name));
368 PetscCall(PetscHMapFormGetSize(hmap, &n));
369 PetscCall(PetscMalloc1(n, &keys));
370 PetscCall(PetscHMapFormGetKeys(hmap, &off, keys));
371 for (i = 0; i < n; ++i) {
372 PetscBool match = PETSC_FALSE;
373 const char *lname = NULL;
374
375 if (label == keys[i].label) continue;
376 if (keys[i].label) PetscCall(PetscObjectGetName((PetscObject)keys[i].label, &lname));
377 PetscCall(PetscStrcmp(name, lname, &match));
378 if ((!name && !lname) || match) {
379 PetscVoidFn **funcs;
380 PetscInt Nf;
381
382 PetscCall(PetscWeakFormGetFunction_Private(wf, hmap, keys[i].label, keys[i].value, keys[i].field, keys[i].part, &Nf, &funcs));
383 maxFuncs = PetscMax(maxFuncs, Nf);
384 }
385 }
386 /* Need temp space because chunk buffer can be reallocated in SetFunction() call */
387 PetscCall(PetscMalloc1(maxFuncs, &tmpf));
388 for (i = 0; i < n; ++i) {
389 PetscBool match = PETSC_FALSE;
390 const char *lname = NULL;
391
392 if (label == keys[i].label) continue;
393 if (keys[i].label) PetscCall(PetscObjectGetName((PetscObject)keys[i].label, &lname));
394 PetscCall(PetscStrcmp(name, lname, &match));
395 if ((!name && !lname) || match) {
396 PetscVoidFn **funcs;
397 PetscInt Nf, j;
398
399 PetscCall(PetscWeakFormGetFunction_Private(wf, hmap, keys[i].label, keys[i].value, keys[i].field, keys[i].part, &Nf, &funcs));
400 for (j = 0; j < Nf; ++j) tmpf[j] = funcs[j];
401 PetscCall(PetscWeakFormSetFunction_Private(wf, hmap, label, keys[i].value, keys[i].field, keys[i].part, Nf, tmpf));
402 PetscCall(PetscWeakFormSetFunction_Private(wf, hmap, keys[i].label, keys[i].value, keys[i].field, keys[i].part, 0, NULL));
403 }
404 }
405 PetscCall(PetscFree(tmpf));
406 PetscCall(PetscFree(keys));
407 PetscFunctionReturn(PETSC_SUCCESS);
408 }
409
410 /*@
411 PetscWeakFormReplaceLabel - Change any key on a label of the same name to use the new label
412
413 Not Collective
414
415 Input Parameters:
416 + wf - The original `PetscWeakForm`
417 - label - The label to change keys for
418
419 Level: intermediate
420
421 Note:
422 This is used internally when meshes are modified
423
424 .seealso: `PetscWeakForm`, `DMLabel`, `PetscWeakFormRewriteKeys()`, `PetscWeakFormCreate()`, `PetscWeakFormDestroy()`
425 @*/
PetscWeakFormReplaceLabel(PetscWeakForm wf,DMLabel label)426 PetscErrorCode PetscWeakFormReplaceLabel(PetscWeakForm wf, DMLabel label)
427 {
428 PetscInt f;
429
430 PetscFunctionBegin;
431 for (f = 0; f < PETSC_NUM_WF; ++f) PetscCall(PetscWeakFormReplaceLabel_Internal(wf, wf->form[f], label));
432 PetscFunctionReturn(PETSC_SUCCESS);
433 }
434
PetscWeakFormClearIndex(PetscWeakForm wf,DMLabel label,PetscInt val,PetscInt f,PetscInt part,PetscWeakFormKind kind,PetscInt ind)435 PetscErrorCode PetscWeakFormClearIndex(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscWeakFormKind kind, PetscInt ind)
436 {
437 PetscFunctionBegin;
438 PetscCall(PetscWeakFormClearIndexFunction_Private(wf, wf->form[kind], label, val, f, part, ind));
439 PetscFunctionReturn(PETSC_SUCCESS);
440 }
441
PetscWeakFormGetObjective(PetscWeakForm wf,DMLabel label,PetscInt val,PetscInt f,PetscInt part,PetscInt * n,void (*** obj)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]))442 PetscErrorCode PetscWeakFormGetObjective(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt *n, void (***obj)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
443 {
444 PetscFunctionBegin;
445 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_OBJECTIVE], label, val, f, part, n, (void (***)(void))obj));
446 PetscFunctionReturn(PETSC_SUCCESS);
447 }
448
PetscWeakFormSetObjective(PetscWeakForm wf,DMLabel label,PetscInt val,PetscInt f,PetscInt part,PetscInt n,void (** obj)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]))449 PetscErrorCode PetscWeakFormSetObjective(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt n, void (**obj)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
450 {
451 PetscFunctionBegin;
452 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_OBJECTIVE], label, val, f, part, n, (PetscVoidFn **)obj));
453 PetscFunctionReturn(PETSC_SUCCESS);
454 }
455
PetscWeakFormAddObjective(PetscWeakForm wf,DMLabel label,PetscInt val,PetscInt f,PetscInt part,void (* obj)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]))456 PetscErrorCode PetscWeakFormAddObjective(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, void (*obj)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
457 {
458 PetscFunctionBegin;
459 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_OBJECTIVE], label, val, f, part, (PetscVoidFn *)obj));
460 PetscFunctionReturn(PETSC_SUCCESS);
461 }
462
PetscWeakFormGetIndexObjective(PetscWeakForm wf,DMLabel label,PetscInt val,PetscInt f,PetscInt part,PetscInt ind,void (** obj)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]))463 PetscErrorCode PetscWeakFormGetIndexObjective(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt ind, void (**obj)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
464 {
465 PetscFunctionBegin;
466 PetscCall(PetscWeakFormGetIndexFunction_Private(wf, wf->form[PETSC_WF_OBJECTIVE], label, val, f, part, ind, (PetscVoidFn **)obj));
467 PetscFunctionReturn(PETSC_SUCCESS);
468 }
469
PetscWeakFormSetIndexObjective(PetscWeakForm wf,DMLabel label,PetscInt val,PetscInt f,PetscInt part,PetscInt ind,void (* obj)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]))470 PetscErrorCode PetscWeakFormSetIndexObjective(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt ind, void (*obj)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
471 {
472 PetscFunctionBegin;
473 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_OBJECTIVE], label, val, f, part, ind, (PetscVoidFn *)obj));
474 PetscFunctionReturn(PETSC_SUCCESS);
475 }
476
PetscWeakFormGetResidual(PetscWeakForm wf,DMLabel label,PetscInt val,PetscInt f,PetscInt part,PetscInt * n0,void (*** f0)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt * n1,void (*** f1)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]))477 PetscErrorCode PetscWeakFormGetResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt *n0, void (***f0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n1, void (***f1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
478 {
479 PetscFunctionBegin;
480 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_F0], label, val, f, part, n0, (void (***)(void))f0));
481 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_F1], label, val, f, part, n1, (void (***)(void))f1));
482 PetscFunctionReturn(PETSC_SUCCESS);
483 }
484
PetscWeakFormAddResidual(PetscWeakForm wf,DMLabel label,PetscInt val,PetscInt f,PetscInt part,void (* f0)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),void (* f1)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]))485 PetscErrorCode PetscWeakFormAddResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, void (*f0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*f1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
486 {
487 PetscFunctionBegin;
488 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_F0], label, val, f, part, (PetscVoidFn *)f0));
489 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_F1], label, val, f, part, (PetscVoidFn *)f1));
490 PetscFunctionReturn(PETSC_SUCCESS);
491 }
492
PetscWeakFormSetResidual(PetscWeakForm wf,DMLabel label,PetscInt val,PetscInt f,PetscInt part,PetscInt n0,void (** f0)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt n1,void (** f1)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]))493 PetscErrorCode PetscWeakFormSetResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt n0, void (**f0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n1, void (**f1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
494 {
495 PetscFunctionBegin;
496 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_F0], label, val, f, part, n0, (PetscVoidFn **)f0));
497 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_F1], label, val, f, part, n1, (PetscVoidFn **)f1));
498 PetscFunctionReturn(PETSC_SUCCESS);
499 }
500
PetscWeakFormSetIndexResidual(PetscWeakForm wf,DMLabel label,PetscInt val,PetscInt f,PetscInt part,PetscInt i0,void (* f0)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt i1,void (* f1)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]))501 PetscErrorCode PetscWeakFormSetIndexResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt i0, void (*f0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i1, void (*f1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
502 {
503 PetscFunctionBegin;
504 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_F0], label, val, f, part, i0, (PetscVoidFn *)f0));
505 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_F1], label, val, f, part, i1, (PetscVoidFn *)f1));
506 PetscFunctionReturn(PETSC_SUCCESS);
507 }
508
PetscWeakFormGetBdResidual(PetscWeakForm wf,DMLabel label,PetscInt val,PetscInt f,PetscInt part,PetscInt * n0,void (*** f0)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,const PetscReal[],const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt * n1,void (*** f1)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,const PetscReal[],const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]))509 PetscErrorCode PetscWeakFormGetBdResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt *n0, void (***f0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n1, void (***f1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
510 {
511 PetscFunctionBegin;
512 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDF0], label, val, f, part, n0, (void (***)(void))f0));
513 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDF1], label, val, f, part, n1, (void (***)(void))f1));
514 PetscFunctionReturn(PETSC_SUCCESS);
515 }
516
PetscWeakFormAddBdResidual(PetscWeakForm wf,DMLabel label,PetscInt val,PetscInt f,PetscInt part,void (* f0)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,const PetscReal[],const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),void (* f1)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,const PetscReal[],const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]))517 PetscErrorCode PetscWeakFormAddBdResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, void (*f0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*f1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
518 {
519 PetscFunctionBegin;
520 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDF0], label, val, f, part, (PetscVoidFn *)f0));
521 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDF1], label, val, f, part, (PetscVoidFn *)f1));
522 PetscFunctionReturn(PETSC_SUCCESS);
523 }
524
PetscWeakFormSetBdResidual(PetscWeakForm wf,DMLabel label,PetscInt val,PetscInt f,PetscInt part,PetscInt n0,void (** f0)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,const PetscReal[],const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt n1,void (** f1)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,const PetscReal[],const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]))525 PetscErrorCode PetscWeakFormSetBdResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt n0, void (**f0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n1, void (**f1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
526 {
527 PetscFunctionBegin;
528 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDF0], label, val, f, part, n0, (PetscVoidFn **)f0));
529 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDF1], label, val, f, part, n1, (PetscVoidFn **)f1));
530 PetscFunctionReturn(PETSC_SUCCESS);
531 }
532
PetscWeakFormSetIndexBdResidual(PetscWeakForm wf,DMLabel label,PetscInt val,PetscInt f,PetscInt part,PetscInt i0,void (* f0)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,const PetscReal[],const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt i1,void (* f1)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,const PetscReal[],const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]))533 PetscErrorCode PetscWeakFormSetIndexBdResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt i0, void (*f0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i1, void (*f1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
534 {
535 PetscFunctionBegin;
536 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDF0], label, val, f, part, i0, (PetscVoidFn *)f0));
537 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDF1], label, val, f, part, i1, (PetscVoidFn *)f1));
538 PetscFunctionReturn(PETSC_SUCCESS);
539 }
540
PetscWeakFormHasJacobian(PetscWeakForm wf,PetscBool * hasJac)541 PetscErrorCode PetscWeakFormHasJacobian(PetscWeakForm wf, PetscBool *hasJac)
542 {
543 PetscInt n0, n1, n2, n3;
544
545 PetscFunctionBegin;
546 PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1);
547 PetscAssertPointer(hasJac, 2);
548 PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_G0], &n0));
549 PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_G1], &n1));
550 PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_G2], &n2));
551 PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_G3], &n3));
552 *hasJac = n0 + n1 + n2 + n3 ? PETSC_TRUE : PETSC_FALSE;
553 PetscFunctionReturn(PETSC_SUCCESS);
554 }
555
PetscWeakFormGetJacobian(PetscWeakForm wf,DMLabel label,PetscInt val,PetscInt f,PetscInt g,PetscInt part,PetscInt * n0,void (*** g0)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt * n1,void (*** g1)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt * n2,void (*** g2)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt * n3,void (*** g3)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]))556 PetscErrorCode PetscWeakFormGetJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt *n0, void (***g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n1, void (***g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n2, void (***g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n3, void (***g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
557 {
558 PetscInt find = f * wf->Nf + g;
559
560 PetscFunctionBegin;
561 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_G0], label, val, find, part, n0, (void (***)(void))g0));
562 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_G1], label, val, find, part, n1, (void (***)(void))g1));
563 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_G2], label, val, find, part, n2, (void (***)(void))g2));
564 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_G3], label, val, find, part, n3, (void (***)(void))g3));
565 PetscFunctionReturn(PETSC_SUCCESS);
566 }
567
PetscWeakFormAddJacobian(PetscWeakForm wf,DMLabel label,PetscInt val,PetscInt f,PetscInt g,PetscInt part,void (* g0)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),void (* g1)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),void (* g2)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),void (* g3)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]))568 PetscErrorCode PetscWeakFormAddJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, void (*g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
569 {
570 PetscInt find = f * wf->Nf + g;
571
572 PetscFunctionBegin;
573 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_G0], label, val, find, part, (PetscVoidFn *)g0));
574 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_G1], label, val, find, part, (PetscVoidFn *)g1));
575 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_G2], label, val, find, part, (PetscVoidFn *)g2));
576 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_G3], label, val, find, part, (PetscVoidFn *)g3));
577 PetscFunctionReturn(PETSC_SUCCESS);
578 }
579
PetscWeakFormSetJacobian(PetscWeakForm wf,DMLabel label,PetscInt val,PetscInt f,PetscInt g,PetscInt part,PetscInt n0,void (** g0)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt n1,void (** g1)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt n2,void (** g2)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt n3,void (** g3)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]))580 PetscErrorCode PetscWeakFormSetJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt n0, void (**g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n1, void (**g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n2, void (**g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n3, void (**g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
581 {
582 PetscInt find = f * wf->Nf + g;
583
584 PetscFunctionBegin;
585 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_G0], label, val, find, part, n0, (PetscVoidFn **)g0));
586 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_G1], label, val, find, part, n1, (PetscVoidFn **)g1));
587 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_G2], label, val, find, part, n2, (PetscVoidFn **)g2));
588 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_G3], label, val, find, part, n3, (PetscVoidFn **)g3));
589 PetscFunctionReturn(PETSC_SUCCESS);
590 }
591
PetscWeakFormSetIndexJacobian(PetscWeakForm wf,DMLabel label,PetscInt val,PetscInt f,PetscInt g,PetscInt part,PetscInt i0,void (* g0)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt i1,void (* g1)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt i2,void (* g2)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt i3,void (* g3)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]))592 PetscErrorCode PetscWeakFormSetIndexJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt i0, void (*g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i1, void (*g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i2, void (*g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i3, void (*g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
593 {
594 PetscInt find = f * wf->Nf + g;
595
596 PetscFunctionBegin;
597 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_G0], label, val, find, part, i0, (PetscVoidFn *)g0));
598 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_G1], label, val, find, part, i1, (PetscVoidFn *)g1));
599 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_G2], label, val, find, part, i2, (PetscVoidFn *)g2));
600 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_G3], label, val, find, part, i3, (PetscVoidFn *)g3));
601 PetscFunctionReturn(PETSC_SUCCESS);
602 }
603
PetscWeakFormHasJacobianPreconditioner(PetscWeakForm wf,PetscBool * hasJacPre)604 PetscErrorCode PetscWeakFormHasJacobianPreconditioner(PetscWeakForm wf, PetscBool *hasJacPre)
605 {
606 PetscInt n0, n1, n2, n3;
607
608 PetscFunctionBegin;
609 PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1);
610 PetscAssertPointer(hasJacPre, 2);
611 PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_GP0], &n0));
612 PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_GP1], &n1));
613 PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_GP2], &n2));
614 PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_GP3], &n3));
615 *hasJacPre = n0 + n1 + n2 + n3 ? PETSC_TRUE : PETSC_FALSE;
616 PetscFunctionReturn(PETSC_SUCCESS);
617 }
618
PetscWeakFormGetJacobianPreconditioner(PetscWeakForm wf,DMLabel label,PetscInt val,PetscInt f,PetscInt g,PetscInt part,PetscInt * n0,void (*** g0)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt * n1,void (*** g1)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt * n2,void (*** g2)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt * n3,void (*** g3)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]))619 PetscErrorCode PetscWeakFormGetJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt *n0, void (***g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n1, void (***g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n2, void (***g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n3, void (***g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
620 {
621 PetscInt find = f * wf->Nf + g;
622
623 PetscFunctionBegin;
624 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GP0], label, val, find, part, n0, (void (***)(void))g0));
625 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GP1], label, val, find, part, n1, (void (***)(void))g1));
626 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GP2], label, val, find, part, n2, (void (***)(void))g2));
627 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GP3], label, val, find, part, n3, (void (***)(void))g3));
628 PetscFunctionReturn(PETSC_SUCCESS);
629 }
630
PetscWeakFormAddJacobianPreconditioner(PetscWeakForm wf,DMLabel label,PetscInt val,PetscInt f,PetscInt g,PetscInt part,void (* g0)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),void (* g1)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),void (* g2)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),void (* g3)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]))631 PetscErrorCode PetscWeakFormAddJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, void (*g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
632 {
633 PetscInt find = f * wf->Nf + g;
634
635 PetscFunctionBegin;
636 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GP0], label, val, find, part, (PetscVoidFn *)g0));
637 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GP1], label, val, find, part, (PetscVoidFn *)g1));
638 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GP2], label, val, find, part, (PetscVoidFn *)g2));
639 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GP3], label, val, find, part, (PetscVoidFn *)g3));
640 PetscFunctionReturn(PETSC_SUCCESS);
641 }
642
PetscWeakFormSetJacobianPreconditioner(PetscWeakForm wf,DMLabel label,PetscInt val,PetscInt f,PetscInt g,PetscInt part,PetscInt n0,void (** g0)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt n1,void (** g1)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt n2,void (** g2)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt n3,void (** g3)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]))643 PetscErrorCode PetscWeakFormSetJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt n0, void (**g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n1, void (**g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n2, void (**g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n3, void (**g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
644 {
645 PetscInt find = f * wf->Nf + g;
646
647 PetscFunctionBegin;
648 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GP0], label, val, find, part, n0, (PetscVoidFn **)g0));
649 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GP1], label, val, find, part, n1, (PetscVoidFn **)g1));
650 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GP2], label, val, find, part, n2, (PetscVoidFn **)g2));
651 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GP3], label, val, find, part, n3, (PetscVoidFn **)g3));
652 PetscFunctionReturn(PETSC_SUCCESS);
653 }
654
PetscWeakFormSetIndexJacobianPreconditioner(PetscWeakForm wf,DMLabel label,PetscInt val,PetscInt f,PetscInt g,PetscInt part,PetscInt i0,void (* g0)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt i1,void (* g1)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt i2,void (* g2)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt i3,void (* g3)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]))655 PetscErrorCode PetscWeakFormSetIndexJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt i0, void (*g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i1, void (*g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i2, void (*g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i3, void (*g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
656 {
657 PetscInt find = f * wf->Nf + g;
658
659 PetscFunctionBegin;
660 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GP0], label, val, find, part, i0, (PetscVoidFn *)g0));
661 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GP1], label, val, find, part, i1, (PetscVoidFn *)g1));
662 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GP2], label, val, find, part, i2, (PetscVoidFn *)g2));
663 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GP3], label, val, find, part, i3, (PetscVoidFn *)g3));
664 PetscFunctionReturn(PETSC_SUCCESS);
665 }
666
PetscWeakFormHasBdJacobian(PetscWeakForm wf,PetscBool * hasJac)667 PetscErrorCode PetscWeakFormHasBdJacobian(PetscWeakForm wf, PetscBool *hasJac)
668 {
669 PetscInt n0, n1, n2, n3;
670
671 PetscFunctionBegin;
672 PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1);
673 PetscAssertPointer(hasJac, 2);
674 PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_BDG0], &n0));
675 PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_BDG1], &n1));
676 PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_BDG2], &n2));
677 PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_BDG3], &n3));
678 *hasJac = n0 + n1 + n2 + n3 ? PETSC_TRUE : PETSC_FALSE;
679 PetscFunctionReturn(PETSC_SUCCESS);
680 }
681
PetscWeakFormGetBdJacobian(PetscWeakForm wf,DMLabel label,PetscInt val,PetscInt f,PetscInt g,PetscInt part,PetscInt * n0,void (*** g0)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt * n1,void (*** g1)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt * n2,void (*** g2)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt * n3,void (*** g3)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]))682 PetscErrorCode PetscWeakFormGetBdJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt *n0, void (***g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n1, void (***g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n2, void (***g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n3, void (***g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
683 {
684 PetscInt find = f * wf->Nf + g;
685
686 PetscFunctionBegin;
687 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDG0], label, val, find, part, n0, (void (***)(void))g0));
688 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDG1], label, val, find, part, n1, (void (***)(void))g1));
689 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDG2], label, val, find, part, n2, (void (***)(void))g2));
690 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDG3], label, val, find, part, n3, (void (***)(void))g3));
691 PetscFunctionReturn(PETSC_SUCCESS);
692 }
693
PetscWeakFormAddBdJacobian(PetscWeakForm wf,DMLabel label,PetscInt val,PetscInt f,PetscInt g,PetscInt part,void (* g0)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),void (* g1)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),void (* g2)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),void (* g3)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]))694 PetscErrorCode PetscWeakFormAddBdJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, void (*g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
695 {
696 PetscInt find = f * wf->Nf + g;
697
698 PetscFunctionBegin;
699 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDG0], label, val, find, part, (PetscVoidFn *)g0));
700 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDG1], label, val, find, part, (PetscVoidFn *)g1));
701 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDG2], label, val, find, part, (PetscVoidFn *)g2));
702 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDG3], label, val, find, part, (PetscVoidFn *)g3));
703 PetscFunctionReturn(PETSC_SUCCESS);
704 }
705
PetscWeakFormSetBdJacobian(PetscWeakForm wf,DMLabel label,PetscInt val,PetscInt f,PetscInt g,PetscInt part,PetscInt n0,void (** g0)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt n1,void (** g1)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt n2,void (** g2)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt n3,void (** g3)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]))706 PetscErrorCode PetscWeakFormSetBdJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt n0, void (**g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n1, void (**g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n2, void (**g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n3, void (**g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
707 {
708 PetscInt find = f * wf->Nf + g;
709
710 PetscFunctionBegin;
711 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDG0], label, val, find, part, n0, (PetscVoidFn **)g0));
712 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDG1], label, val, find, part, n1, (PetscVoidFn **)g1));
713 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDG2], label, val, find, part, n2, (PetscVoidFn **)g2));
714 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDG3], label, val, find, part, n3, (PetscVoidFn **)g3));
715 PetscFunctionReturn(PETSC_SUCCESS);
716 }
717
PetscWeakFormSetIndexBdJacobian(PetscWeakForm wf,DMLabel label,PetscInt val,PetscInt f,PetscInt g,PetscInt part,PetscInt i0,void (* g0)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt i1,void (* g1)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt i2,void (* g2)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt i3,void (* g3)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]))718 PetscErrorCode PetscWeakFormSetIndexBdJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt i0, void (*g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i1, void (*g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i2, void (*g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i3, void (*g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
719 {
720 PetscInt find = f * wf->Nf + g;
721
722 PetscFunctionBegin;
723 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDG0], label, val, find, part, i0, (PetscVoidFn *)g0));
724 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDG1], label, val, find, part, i1, (PetscVoidFn *)g1));
725 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDG2], label, val, find, part, i2, (PetscVoidFn *)g2));
726 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDG3], label, val, find, part, i3, (PetscVoidFn *)g3));
727 PetscFunctionReturn(PETSC_SUCCESS);
728 }
729
PetscWeakFormHasBdJacobianPreconditioner(PetscWeakForm wf,PetscBool * hasJacPre)730 PetscErrorCode PetscWeakFormHasBdJacobianPreconditioner(PetscWeakForm wf, PetscBool *hasJacPre)
731 {
732 PetscInt n0, n1, n2, n3;
733
734 PetscFunctionBegin;
735 PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1);
736 PetscAssertPointer(hasJacPre, 2);
737 PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_BDGP0], &n0));
738 PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_BDGP1], &n1));
739 PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_BDGP2], &n2));
740 PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_BDGP3], &n3));
741 *hasJacPre = n0 + n1 + n2 + n3 ? PETSC_TRUE : PETSC_FALSE;
742 PetscFunctionReturn(PETSC_SUCCESS);
743 }
744
PetscWeakFormGetBdJacobianPreconditioner(PetscWeakForm wf,DMLabel label,PetscInt val,PetscInt f,PetscInt g,PetscInt part,PetscInt * n0,void (*** g0)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt * n1,void (*** g1)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt * n2,void (*** g2)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt * n3,void (*** g3)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]))745 PetscErrorCode PetscWeakFormGetBdJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt *n0, void (***g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n1, void (***g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n2, void (***g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n3, void (***g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
746 {
747 PetscInt find = f * wf->Nf + g;
748
749 PetscFunctionBegin;
750 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDGP0], label, val, find, part, n0, (void (***)(void))g0));
751 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDGP1], label, val, find, part, n1, (void (***)(void))g1));
752 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDGP2], label, val, find, part, n2, (void (***)(void))g2));
753 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDGP3], label, val, find, part, n3, (void (***)(void))g3));
754 PetscFunctionReturn(PETSC_SUCCESS);
755 }
756
PetscWeakFormAddBdJacobianPreconditioner(PetscWeakForm wf,DMLabel label,PetscInt val,PetscInt f,PetscInt g,PetscInt part,void (* g0)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),void (* g1)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),void (* g2)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),void (* g3)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]))757 PetscErrorCode PetscWeakFormAddBdJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, void (*g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
758 {
759 PetscInt find = f * wf->Nf + g;
760
761 PetscFunctionBegin;
762 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDGP0], label, val, find, part, (PetscVoidFn *)g0));
763 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDGP1], label, val, find, part, (PetscVoidFn *)g1));
764 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDGP2], label, val, find, part, (PetscVoidFn *)g2));
765 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDGP3], label, val, find, part, (PetscVoidFn *)g3));
766 PetscFunctionReturn(PETSC_SUCCESS);
767 }
768
PetscWeakFormSetBdJacobianPreconditioner(PetscWeakForm wf,DMLabel label,PetscInt val,PetscInt f,PetscInt g,PetscInt part,PetscInt n0,void (** g0)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt n1,void (** g1)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt n2,void (** g2)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt n3,void (** g3)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]))769 PetscErrorCode PetscWeakFormSetBdJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt n0, void (**g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n1, void (**g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n2, void (**g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n3, void (**g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
770 {
771 PetscInt find = f * wf->Nf + g;
772
773 PetscFunctionBegin;
774 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDGP0], label, val, find, part, n0, (PetscVoidFn **)g0));
775 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDGP1], label, val, find, part, n1, (PetscVoidFn **)g1));
776 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDGP2], label, val, find, part, n2, (PetscVoidFn **)g2));
777 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDGP3], label, val, find, part, n3, (PetscVoidFn **)g3));
778 PetscFunctionReturn(PETSC_SUCCESS);
779 }
780
PetscWeakFormSetIndexBdJacobianPreconditioner(PetscWeakForm wf,DMLabel label,PetscInt val,PetscInt f,PetscInt g,PetscInt part,PetscInt i0,void (* g0)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt i1,void (* g1)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt i2,void (* g2)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt i3,void (* g3)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]))781 PetscErrorCode PetscWeakFormSetIndexBdJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt i0, void (*g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i1, void (*g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i2, void (*g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i3, void (*g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
782 {
783 PetscInt find = f * wf->Nf + g;
784
785 PetscFunctionBegin;
786 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDGP0], label, val, find, part, i0, (PetscVoidFn *)g0));
787 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDGP1], label, val, find, part, i1, (PetscVoidFn *)g1));
788 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDGP2], label, val, find, part, i2, (PetscVoidFn *)g2));
789 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDGP3], label, val, find, part, i3, (PetscVoidFn *)g3));
790 PetscFunctionReturn(PETSC_SUCCESS);
791 }
792
PetscWeakFormHasDynamicJacobian(PetscWeakForm wf,PetscBool * hasDynJac)793 PetscErrorCode PetscWeakFormHasDynamicJacobian(PetscWeakForm wf, PetscBool *hasDynJac)
794 {
795 PetscInt n0, n1, n2, n3;
796
797 PetscFunctionBegin;
798 PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1);
799 PetscAssertPointer(hasDynJac, 2);
800 PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_GT0], &n0));
801 PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_GT1], &n1));
802 PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_GT2], &n2));
803 PetscCall(PetscHMapFormGetSize(wf->form[PETSC_WF_GT3], &n3));
804 *hasDynJac = n0 + n1 + n2 + n3 ? PETSC_TRUE : PETSC_FALSE;
805 PetscFunctionReturn(PETSC_SUCCESS);
806 }
807
PetscWeakFormGetDynamicJacobian(PetscWeakForm wf,DMLabel label,PetscInt val,PetscInt f,PetscInt g,PetscInt part,PetscInt * n0,void (*** g0)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt * n1,void (*** g1)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt * n2,void (*** g2)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt * n3,void (*** g3)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]))808 PetscErrorCode PetscWeakFormGetDynamicJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt *n0, void (***g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n1, void (***g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n2, void (***g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt *n3, void (***g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
809 {
810 PetscInt find = f * wf->Nf + g;
811
812 PetscFunctionBegin;
813 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GT0], label, val, find, part, n0, (void (***)(void))g0));
814 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GT1], label, val, find, part, n1, (void (***)(void))g1));
815 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GT2], label, val, find, part, n2, (void (***)(void))g2));
816 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GT3], label, val, find, part, n3, (void (***)(void))g3));
817 PetscFunctionReturn(PETSC_SUCCESS);
818 }
819
PetscWeakFormAddDynamicJacobian(PetscWeakForm wf,DMLabel label,PetscInt val,PetscInt f,PetscInt g,PetscInt part,void (* g0)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),void (* g1)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),void (* g2)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),void (* g3)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]))820 PetscErrorCode PetscWeakFormAddDynamicJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, void (*g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), void (*g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
821 {
822 PetscInt find = f * wf->Nf + g;
823
824 PetscFunctionBegin;
825 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GT0], label, val, find, part, (PetscVoidFn *)g0));
826 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GT1], label, val, find, part, (PetscVoidFn *)g1));
827 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GT2], label, val, find, part, (PetscVoidFn *)g2));
828 PetscCall(PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GT3], label, val, find, part, (PetscVoidFn *)g3));
829 PetscFunctionReturn(PETSC_SUCCESS);
830 }
831
PetscWeakFormSetDynamicJacobian(PetscWeakForm wf,DMLabel label,PetscInt val,PetscInt f,PetscInt g,PetscInt part,PetscInt n0,void (** g0)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt n1,void (** g1)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt n2,void (** g2)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt n3,void (** g3)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]))832 PetscErrorCode PetscWeakFormSetDynamicJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt n0, void (**g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n1, void (**g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n2, void (**g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt n3, void (**g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
833 {
834 PetscInt find = f * wf->Nf + g;
835
836 PetscFunctionBegin;
837 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GT0], label, val, find, part, n0, (PetscVoidFn **)g0));
838 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GT1], label, val, find, part, n1, (PetscVoidFn **)g1));
839 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GT2], label, val, find, part, n2, (PetscVoidFn **)g2));
840 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GT3], label, val, find, part, n3, (PetscVoidFn **)g3));
841 PetscFunctionReturn(PETSC_SUCCESS);
842 }
843
PetscWeakFormSetIndexDynamicJacobian(PetscWeakForm wf,DMLabel label,PetscInt val,PetscInt f,PetscInt g,PetscInt part,PetscInt i0,void (* g0)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt i1,void (* g1)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt i2,void (* g2)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]),PetscInt i3,void (* g3)(PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],const PetscInt[],const PetscInt[],const PetscScalar[],const PetscScalar[],const PetscScalar[],PetscReal,PetscReal,const PetscReal[],PetscInt,const PetscScalar[],PetscScalar[]))844 PetscErrorCode PetscWeakFormSetIndexDynamicJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part, PetscInt i0, void (*g0)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i1, void (*g1)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i2, void (*g2)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt i3, void (*g3)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
845 {
846 PetscInt find = f * wf->Nf + g;
847
848 PetscFunctionBegin;
849 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GT0], label, val, find, part, i0, (PetscVoidFn *)g0));
850 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GT1], label, val, find, part, i1, (PetscVoidFn *)g1));
851 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GT2], label, val, find, part, i2, (PetscVoidFn *)g2));
852 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GT3], label, val, find, part, i3, (PetscVoidFn *)g3));
853 PetscFunctionReturn(PETSC_SUCCESS);
854 }
855
PetscWeakFormGetRiemannSolver(PetscWeakForm wf,DMLabel label,PetscInt val,PetscInt f,PetscInt part,PetscInt * n,void (*** r)(PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscScalar[],const PetscScalar[],PetscInt,const PetscScalar[],PetscScalar[],void *))856 PetscErrorCode PetscWeakFormGetRiemannSolver(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt *n, void (***r)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *))
857 {
858 PetscFunctionBegin;
859 PetscCall(PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_R], label, val, f, part, n, (void (***)(void))r));
860 PetscFunctionReturn(PETSC_SUCCESS);
861 }
862
PetscWeakFormSetRiemannSolver(PetscWeakForm wf,DMLabel label,PetscInt val,PetscInt f,PetscInt part,PetscInt n,void (** r)(PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscScalar[],const PetscScalar[],PetscInt,const PetscScalar[],PetscScalar[],void *))863 PetscErrorCode PetscWeakFormSetRiemannSolver(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt n, void (**r)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *))
864 {
865 PetscFunctionBegin;
866 PetscCall(PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_R], label, val, f, part, n, (PetscVoidFn **)r));
867 PetscFunctionReturn(PETSC_SUCCESS);
868 }
869
PetscWeakFormSetIndexRiemannSolver(PetscWeakForm wf,DMLabel label,PetscInt val,PetscInt f,PetscInt part,PetscInt i,void (* r)(PetscInt,PetscInt,const PetscReal[],const PetscReal[],const PetscScalar[],const PetscScalar[],PetscInt,const PetscScalar[],PetscScalar[],void *))870 PetscErrorCode PetscWeakFormSetIndexRiemannSolver(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt i, void (*r)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *))
871 {
872 PetscFunctionBegin;
873 PetscCall(PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_R], label, val, f, part, i, (PetscVoidFn *)r));
874 PetscFunctionReturn(PETSC_SUCCESS);
875 }
876
877 /*@
878 PetscWeakFormGetNumFields - Returns the number of fields in a `PetscWeakForm`
879
880 Not Collective
881
882 Input Parameter:
883 . wf - The `PetscWeakForm` object
884
885 Output Parameter:
886 . Nf - The number of fields
887
888 Level: beginner
889
890 .seealso: `PetscWeakForm`, `PetscWeakFormSetNumFields()`, `PetscWeakFormCreate()`
891 @*/
PetscWeakFormGetNumFields(PetscWeakForm wf,PetscInt * Nf)892 PetscErrorCode PetscWeakFormGetNumFields(PetscWeakForm wf, PetscInt *Nf)
893 {
894 PetscFunctionBegin;
895 PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1);
896 PetscAssertPointer(Nf, 2);
897 *Nf = wf->Nf;
898 PetscFunctionReturn(PETSC_SUCCESS);
899 }
900
901 /*@
902 PetscWeakFormSetNumFields - Sets the number of fields
903
904 Not Collective
905
906 Input Parameters:
907 + wf - The `PetscWeakForm` object
908 - Nf - The number of fields
909
910 Level: beginner
911
912 .seealso: `PetscWeakForm`, `PetscWeakFormGetNumFields()`, `PetscWeakFormCreate()`
913 @*/
PetscWeakFormSetNumFields(PetscWeakForm wf,PetscInt Nf)914 PetscErrorCode PetscWeakFormSetNumFields(PetscWeakForm wf, PetscInt Nf)
915 {
916 PetscFunctionBegin;
917 PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1);
918 wf->Nf = Nf;
919 PetscFunctionReturn(PETSC_SUCCESS);
920 }
921
922 /*@
923 PetscWeakFormDestroy - Destroys a `PetscWeakForm` object
924
925 Collective
926
927 Input Parameter:
928 . wf - the `PetscWeakForm` object to destroy
929
930 Level: developer
931
932 .seealso: `PetscWeakForm`, `PetscWeakFormCreate()`, `PetscWeakFormView()`
933 @*/
PetscWeakFormDestroy(PetscWeakForm * wf)934 PetscErrorCode PetscWeakFormDestroy(PetscWeakForm *wf)
935 {
936 PetscInt f;
937
938 PetscFunctionBegin;
939 if (!*wf) PetscFunctionReturn(PETSC_SUCCESS);
940 PetscValidHeaderSpecific(*wf, PETSCWEAKFORM_CLASSID, 1);
941
942 if (--((PetscObject)*wf)->refct > 0) {
943 *wf = NULL;
944 PetscFunctionReturn(PETSC_SUCCESS);
945 }
946 ((PetscObject)*wf)->refct = 0;
947 PetscCall(PetscChunkBufferDestroy(&(*wf)->funcs));
948 for (f = 0; f < PETSC_NUM_WF; ++f) PetscCall(PetscHMapFormDestroy(&(*wf)->form[f]));
949 PetscCall(PetscFree((*wf)->form));
950 PetscCall(PetscHeaderDestroy(wf));
951 PetscFunctionReturn(PETSC_SUCCESS);
952 }
953
PetscWeakFormViewTable_Ascii(PetscWeakForm wf,PetscViewer viewer,PetscBool splitField,const char tableName[],PetscHMapForm map)954 static PetscErrorCode PetscWeakFormViewTable_Ascii(PetscWeakForm wf, PetscViewer viewer, PetscBool splitField, const char tableName[], PetscHMapForm map)
955 {
956 PetscInt Nf = wf->Nf, Nk, k;
957
958 PetscFunctionBegin;
959 PetscCall(PetscHMapFormGetSize(map, &Nk));
960 if (Nk) {
961 PetscFormKey *keys;
962 PetscVoidFn **funcs = NULL;
963 const char **names;
964 PetscInt *values, *idx1, *idx2, *idx;
965 PetscBool showPart = PETSC_FALSE, showPointer = PETSC_FALSE;
966 PetscInt off = 0;
967
968 PetscCall(PetscMalloc6(Nk, &keys, Nk, &names, Nk, &values, Nk, &idx1, Nk, &idx2, Nk, &idx));
969 PetscCall(PetscHMapFormGetKeys(map, &off, keys));
970 /* Sort keys by label name and value */
971 {
972 /* First sort values */
973 for (k = 0; k < Nk; ++k) {
974 values[k] = keys[k].value;
975 idx1[k] = k;
976 }
977 PetscCall(PetscSortIntWithPermutation(Nk, values, idx1));
978 /* If the string sort is stable, it will be sorted correctly overall */
979 for (k = 0; k < Nk; ++k) {
980 if (keys[idx1[k]].label) PetscCall(PetscObjectGetName((PetscObject)keys[idx1[k]].label, &names[k]));
981 else names[k] = "";
982 idx2[k] = k;
983 }
984 PetscCall(PetscSortStrWithPermutation(Nk, names, idx2));
985 for (k = 0; k < Nk; ++k) {
986 if (keys[k].label) PetscCall(PetscObjectGetName((PetscObject)keys[k].label, &names[k]));
987 else names[k] = "";
988 idx[k] = idx1[idx2[k]];
989 }
990 }
991 PetscCall(PetscViewerASCIIPrintf(viewer, "%s\n", tableName));
992 PetscCall(PetscViewerASCIIPushTab(viewer));
993 for (k = 0; k < Nk; ++k) {
994 if (keys[k].part != 0) showPart = PETSC_TRUE;
995 }
996 for (k = 0; k < Nk; ++k) {
997 const PetscInt i = idx[k];
998 PetscInt n, f;
999
1000 if (keys[i].label) {
1001 if (showPointer) PetscCall(PetscViewerASCIIPrintf(viewer, "(%s:%p, %" PetscInt_FMT ") ", names[i], (void *)keys[i].label, keys[i].value));
1002 else PetscCall(PetscViewerASCIIPrintf(viewer, "(%s, %" PetscInt_FMT ") ", names[i], keys[i].value));
1003 }
1004 PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
1005 if (splitField) PetscCall(PetscViewerASCIIPrintf(viewer, "(%" PetscInt_FMT ", %" PetscInt_FMT ") ", keys[i].field / Nf, keys[i].field % Nf));
1006 else PetscCall(PetscViewerASCIIPrintf(viewer, "(%" PetscInt_FMT ") ", keys[i].field));
1007 if (showPart) PetscCall(PetscViewerASCIIPrintf(viewer, "(%" PetscInt_FMT ") ", keys[i].part));
1008 PetscCall(PetscWeakFormGetFunction_Private(wf, map, keys[i].label, keys[i].value, keys[i].field, keys[i].part, &n, &funcs));
1009 for (f = 0; f < n; ++f) {
1010 char *fname;
1011 size_t len, l;
1012
1013 if (f > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ", "));
1014 PetscCall(PetscDLAddr(funcs[f], &fname));
1015 if (fname) {
1016 /* Eliminate argument types */
1017 PetscCall(PetscStrlen(fname, &len));
1018 for (l = 0; l < len; ++l)
1019 if (fname[l] == '(') {
1020 fname[l] = '\0';
1021 break;
1022 }
1023 PetscCall(PetscViewerASCIIPrintf(viewer, "%s", fname));
1024 } else if (showPointer) {
1025 #if defined(__clang__)
1026 PETSC_PRAGMA_DIAGNOSTIC_IGNORED_BEGIN("-Wformat-pedantic")
1027 #elif defined(__GNUC__) || defined(__GNUG__)
1028 PETSC_PRAGMA_DIAGNOSTIC_IGNORED_BEGIN("-Wformat")
1029 #endif
1030 PetscCall(PetscViewerASCIIPrintf(viewer, "%p", funcs[f]));
1031 PETSC_PRAGMA_DIAGNOSTIC_IGNORED_END()
1032 }
1033 PetscCall(PetscFree(fname));
1034 }
1035 PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
1036 PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
1037 }
1038 PetscCall(PetscViewerASCIIPopTab(viewer));
1039 PetscCall(PetscFree6(keys, names, values, idx1, idx2, idx));
1040 }
1041 PetscFunctionReturn(PETSC_SUCCESS);
1042 }
1043
PetscWeakFormView_Ascii(PetscWeakForm wf,PetscViewer viewer)1044 static PetscErrorCode PetscWeakFormView_Ascii(PetscWeakForm wf, PetscViewer viewer)
1045 {
1046 PetscViewerFormat format;
1047 PetscInt f;
1048
1049 PetscFunctionBegin;
1050 PetscCall(PetscViewerGetFormat(viewer, &format));
1051 PetscCall(PetscViewerASCIIPrintf(viewer, "Weak Form System with %" PetscInt_FMT " fields\n", wf->Nf));
1052 PetscCall(PetscViewerASCIIPushTab(viewer));
1053 for (f = 0; f < PETSC_NUM_WF; ++f) PetscCall(PetscWeakFormViewTable_Ascii(wf, viewer, PETSC_TRUE, PetscWeakFormKinds[f], wf->form[f]));
1054 PetscCall(PetscViewerASCIIPopTab(viewer));
1055 PetscFunctionReturn(PETSC_SUCCESS);
1056 }
1057
1058 /*@
1059 PetscWeakFormView - Views a `PetscWeakForm`
1060
1061 Collective
1062
1063 Input Parameters:
1064 + wf - the `PetscWeakForm` object to view
1065 - v - the viewer
1066
1067 Level: developer
1068
1069 .seealso: `PetscViewer`, `PetscWeakForm`, `PetscWeakFormDestroy()`, `PetscWeakFormCreate()`
1070 @*/
PetscWeakFormView(PetscWeakForm wf,PetscViewer v)1071 PetscErrorCode PetscWeakFormView(PetscWeakForm wf, PetscViewer v)
1072 {
1073 PetscBool isascii;
1074
1075 PetscFunctionBegin;
1076 PetscValidHeaderSpecific(wf, PETSCWEAKFORM_CLASSID, 1);
1077 if (!v) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)wf), &v));
1078 else PetscValidHeaderSpecific(v, PETSC_VIEWER_CLASSID, 2);
1079 PetscCall(PetscObjectTypeCompare((PetscObject)v, PETSCVIEWERASCII, &isascii));
1080 if (isascii) PetscCall(PetscWeakFormView_Ascii(wf, v));
1081 PetscTryTypeMethod(wf, view, v);
1082 PetscFunctionReturn(PETSC_SUCCESS);
1083 }
1084
1085 /*@
1086 PetscWeakFormCreate - Creates an empty `PetscWeakForm` object.
1087
1088 Collective
1089
1090 Input Parameter:
1091 . comm - The communicator for the `PetscWeakForm` object
1092
1093 Output Parameter:
1094 . wf - The `PetscWeakForm` object
1095
1096 Level: beginner
1097
1098 .seealso: `PetscWeakForm`, `PetscDS`, `PetscWeakFormDestroy()`
1099 @*/
PetscWeakFormCreate(MPI_Comm comm,PetscWeakForm * wf)1100 PetscErrorCode PetscWeakFormCreate(MPI_Comm comm, PetscWeakForm *wf)
1101 {
1102 PetscWeakForm p;
1103 PetscInt f;
1104
1105 PetscFunctionBegin;
1106 PetscAssertPointer(wf, 2);
1107 PetscCall(PetscDSInitializePackage());
1108
1109 PetscCall(PetscHeaderCreate(p, PETSCWEAKFORM_CLASSID, "PetscWeakForm", "Weak Form System", "PetscWeakForm", comm, PetscWeakFormDestroy, PetscWeakFormView));
1110 p->Nf = 0;
1111 PetscCall(PetscChunkBufferCreate(sizeof(&PetscWeakFormCreate), 2, &p->funcs));
1112 PetscCall(PetscMalloc1(PETSC_NUM_WF, &p->form));
1113 for (f = 0; f < PETSC_NUM_WF; ++f) PetscCall(PetscHMapFormCreate(&p->form[f]));
1114 *wf = p;
1115 PetscFunctionReturn(PETSC_SUCCESS);
1116 }
1117