Actual source code: dtweakform.c

  1: #include <petsc/private/petscdsimpl.h>

  3: PetscClassId PETSCWEAKFORM_CLASSID = 0;

  5: static PetscErrorCode PetscChunkBufferCreate(size_t unitbytes, size_t expected, PetscChunkBuffer **buffer)
  6: {

 10:   PetscNew(buffer);
 11:   PetscCalloc1(expected*unitbytes, &(*buffer)->array);
 12:   (*buffer)->size      = expected;
 13:   (*buffer)->unitbytes = unitbytes;
 14:   (*buffer)->alloc     = expected*unitbytes;
 15:   return(0);
 16: }

 18: static PetscErrorCode PetscChunkBufferDestroy(PetscChunkBuffer **buffer)
 19: {

 23:   PetscFree((*buffer)->array);
 24:   PetscFree(*buffer);
 25:   return(0);
 26: }

 28: static PetscErrorCode PetscChunkBufferCreateChunk(PetscChunkBuffer *buffer, PetscInt size, PetscChunk *chunk)
 29: {

 33:   if ((buffer->size + size)*buffer->unitbytes > buffer->alloc) {
 34:     char *tmp;

 36:     if (!buffer->alloc) buffer->alloc = (buffer->size + size)*buffer->unitbytes;
 37:     while ((buffer->size + size)*buffer->unitbytes > buffer->alloc) buffer->alloc *= 2;
 38:     PetscMalloc(buffer->alloc, &tmp);
 39:     PetscMemcpy(tmp, buffer->array, buffer->size*buffer->unitbytes);
 40:     PetscFree(buffer->array);
 41:     buffer->array = tmp;
 42:   }
 43:   chunk->start    = buffer->size*buffer->unitbytes;
 44:   chunk->size     = size;
 45:   chunk->reserved = size;
 46:   buffer->size   += size;
 47:   return(0);
 48: }

 50: static PetscErrorCode PetscChunkBufferEnlargeChunk(PetscChunkBuffer *buffer, PetscInt size, PetscChunk *chunk)
 51: {
 52:   size_t         siz = size;

 56:   if (chunk->size + size > chunk->reserved) {
 57:     PetscChunk newchunk;
 58:     PetscInt   reserved = chunk->size;

 60:     /* TODO Here if we had a chunk list, we could update them all to reclaim unused space */
 61:     while (reserved < chunk->size+size) reserved *= 2;
 62:     PetscChunkBufferCreateChunk(buffer, (size_t) reserved, &newchunk);
 63:     newchunk.size = chunk->size+size;
 64:     PetscMemcpy(&buffer->array[newchunk.start], &buffer->array[chunk->start], chunk->size * buffer->unitbytes);
 65:     *chunk = newchunk;
 66:   } else {
 67:     chunk->size += siz;
 68:   }
 69:   return(0);
 70: }

 72: /*@C
 73:   PetscHashFormKeySort - Sorts an array of PetscHashFormKey in place in increasing order.

 75:   Not Collective

 77:   Input Parameters:
 78: + n - number of values
 79: - X - array of PetscHashFormKey

 81:   Level: intermediate

 83: .seealso: PetscIntSortSemiOrdered(), PetscSortInt()
 84: @*/
 85: PetscErrorCode PetscHashFormKeySort(PetscInt n, PetscHashFormKey arr[])
 86: {

 90:   if (n <= 1) return(0);
 92:   PetscTimSort(n, arr, sizeof(PetscHashFormKey), Compare_PetscHashFormKey_Private, NULL);
 93:   return(0);
 94: }

 96: PetscErrorCode PetscWeakFormGetFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt *n, void (***func)())
 97: {
 98:   PetscHashFormKey key;
 99:   PetscChunk       chunk;
100:   PetscErrorCode   ierr;

103:   key.label = label; key.value = value; key.field = f;
104:   PetscHMapFormGet(ht, key, &chunk);
105:   if (chunk.size < 0) {*n = 0;          *func = NULL;}
106:   else                {*n = chunk.size; *func = (void (**)()) &wf->funcs->array[chunk.start];}
107:   return(0);
108: }

110: /* A NULL argument for func causes this to clear the key */
111: PetscErrorCode PetscWeakFormSetFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt n, void (**func)())
112: {
113:   PetscHashFormKey key;
114:   PetscChunk       chunk;
115:   PetscInt         i;
116:   PetscErrorCode   ierr;

119:   key.label = label; key.value = value; key.field = f;
120:   if (!func) {
121:     PetscHMapFormDel(ht, key);
122:     return(0);
123:   } else {
124:     PetscHMapFormGet(ht, key, &chunk);
125:   }
126:   if (chunk.size < 0) {
127:     PetscChunkBufferCreateChunk(wf->funcs, n, &chunk);
128:     PetscHMapFormSet(ht, key, chunk);
129:   } else if (chunk.size <= n) {
130:     PetscChunkBufferEnlargeChunk(wf->funcs, n - chunk.size, &chunk);
131:     PetscHMapFormSet(ht, key, chunk);
132:   }
133:   for (i = 0; i < n; ++i) ((void (**)()) &wf->funcs->array[chunk.start])[i] = func[i];
134:   return(0);
135: }

137: PetscErrorCode PetscWeakFormAddFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, void (*func)())
138: {
139:   PetscHashFormKey key;
140:   PetscChunk       chunk;
141:   PetscErrorCode   ierr;

144:   if (!func) return(0);
145:   key.label = label; key.value = value; key.field = f;
146:   PetscHMapFormGet(ht, key, &chunk);
147:   if (chunk.size < 0) {
148:     PetscChunkBufferCreateChunk(wf->funcs, 1, &chunk);
149:     PetscHMapFormSet(ht, key, chunk);
150:     ((void (**)()) &wf->funcs->array[chunk.start])[0] = func;
151:   } else {
152:     PetscChunkBufferEnlargeChunk(wf->funcs, 1, &chunk);
153:     PetscHMapFormSet(ht, key, chunk);
154:     ((void (**)()) &wf->funcs->array[chunk.start])[chunk.size-1] = func;
155:   }
156:   return(0);
157: }

159: PetscErrorCode PetscWeakFormGetIndexFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt ind, void (**func)())
160: {
161:   PetscHashFormKey key;
162:   PetscChunk       chunk;
163:   PetscErrorCode   ierr;

166:   key.label = label; key.value = value; key.field = f;
167:   PetscHMapFormGet(ht, key, &chunk);
168:   if (chunk.size < 0) {*func = NULL;}
169:   else {
170:     if (ind >= chunk.size) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Index %D not in [0, %D)", ind, chunk.size);
171:     *func = ((void (**)()) &wf->funcs->array[chunk.start])[ind];
172:   }
173:   return(0);
174: }

176: /* A NULL argument for func causes this to clear the slot, and if there is nothing else, clear the key */
177: PetscErrorCode PetscWeakFormSetIndexFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt ind, void (*func)())
178: {
179:   PetscHashFormKey key;
180:   PetscChunk       chunk;
181:   PetscErrorCode   ierr;

184:   key.label = label; key.value = value; key.field = f;
185:   PetscHMapFormGet(ht, key, &chunk);
186:   if (chunk.size < 0) {
187:     if (!func) return(0);
188:     PetscChunkBufferCreateChunk(wf->funcs, ind+1, &chunk);
189:     PetscHMapFormSet(ht, key, chunk);
190:   } else if (!func && !ind && chunk.size == 1) {
191:     PetscHMapFormDel(ht, key);
192:     return(0);
193:   } else if (chunk.size <= ind) {
194:     PetscChunkBufferEnlargeChunk(wf->funcs, ind - chunk.size + 1, &chunk);
195:     PetscHMapFormSet(ht, key, chunk);
196:   }
197:   ((void (**)()) &wf->funcs->array[chunk.start])[ind] = func;
198:   return(0);
199: }

201: PetscErrorCode PetscWeakFormGetObjective(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt *n,
202:                                          void (***obj)(PetscInt, PetscInt, PetscInt,
203:                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
204:                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
205:                                                        PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
206: {

210:   PetscWeakFormGetFunction_Private(wf, wf->obj, label, val, f, n, (void (***)(void)) obj);
211:   return(0);
212: }

214: PetscErrorCode PetscWeakFormSetObjective(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt n,
215:                                          void (**obj)(PetscInt, PetscInt, PetscInt,
216:                                                       const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const  PetscScalar[],
217:                                                       const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
218:                                                       PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
219: {

223:   PetscWeakFormSetFunction_Private(wf, wf->obj, label, val, f, n, (void (**)(void)) obj);
224:   return(0);
225: }

227: PetscErrorCode PetscWeakFormAddObjective(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f,
228:                                          void (*obj)(PetscInt, PetscInt, PetscInt,
229:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
230:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
231:                                                      PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
232: {

236:   PetscWeakFormAddFunction_Private(wf, wf->obj, label, val, f, (void (*)(void)) obj);
237:   return(0);
238: }

240: PetscErrorCode PetscWeakFormGetIndexObjective(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt ind,
241:                                               void (**obj)(PetscInt, PetscInt, PetscInt,
242:                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
243:                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
244:                                                            PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
245: {

249:   PetscWeakFormGetIndexFunction_Private(wf, wf->obj, label, val, f, ind, (void (**)(void)) obj);
250:   return(0);
251: }

253: PetscErrorCode PetscWeakFormSetIndexObjective(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt ind,
254:                                               void (*obj)(PetscInt, PetscInt, PetscInt,
255:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
256:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
257:                                                           PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
258: {

262:   PetscWeakFormSetIndexFunction_Private(wf, wf->obj, label, val, f, ind, (void (*)(void)) obj);
263:   return(0);
264: }

266: PetscErrorCode PetscWeakFormGetResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f,
267:                                         PetscInt *n0,
268:                                         void (***f0)(PetscInt, PetscInt, PetscInt,
269:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
270:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
271:                                                      PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
272:                                         PetscInt *n1,
273:                                         void (***f1)(PetscInt, PetscInt, PetscInt,
274:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
275:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
276:                                                      PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
277: {

281:   PetscWeakFormGetFunction_Private(wf, wf->f0, label, val, f, n0, (void (***)(void)) f0);
282:   PetscWeakFormGetFunction_Private(wf, wf->f1, label, val, f, n1, (void (***)(void)) f1);
283:   return(0);
284: }

286: PetscErrorCode PetscWeakFormAddResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f,
287:                                         void (*f0)(PetscInt, PetscInt, PetscInt,
288:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
289:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
290:                                                    PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
291:                                         void (*f1)(PetscInt, PetscInt, PetscInt,
292:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
293:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
294:                                                    PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
295: {

299:   PetscWeakFormAddFunction_Private(wf, wf->f0, label, val, f, (void (*)(void)) f0);
300:   PetscWeakFormAddFunction_Private(wf, wf->f1, label, val, f, (void (*)(void)) f1);
301:   return(0);
302: }

304: PetscErrorCode PetscWeakFormSetResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f,
305:                                         PetscInt n0,
306:                                         void (**f0)(PetscInt, PetscInt, PetscInt,
307:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
308:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
309:                                                    PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
310:                                         PetscInt n1,
311:                                         void (**f1)(PetscInt, PetscInt, PetscInt,
312:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
313:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
314:                                                    PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
315: {

319:   PetscWeakFormSetFunction_Private(wf, wf->f0, label, val, f, n0, (void (**)(void)) f0);
320:   PetscWeakFormSetFunction_Private(wf, wf->f1, label, val, f, n1, (void (**)(void)) f1);
321:   return(0);
322: }

324: PetscErrorCode PetscWeakFormSetIndexResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f,
325:                                         PetscInt i0,
326:                                         void (*f0)(PetscInt, PetscInt, PetscInt,
327:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
328:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
329:                                                    PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
330:                                         PetscInt i1,
331:                                         void (*f1)(PetscInt, PetscInt, PetscInt,
332:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
333:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
334:                                                    PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
335: {

339:   PetscWeakFormSetIndexFunction_Private(wf, wf->f0, label, val, f, i0, (void (*)(void)) f0);
340:   PetscWeakFormSetIndexFunction_Private(wf, wf->f1, label, val, f, i1, (void (*)(void)) f1);
341:   return(0);
342: }

344: PetscErrorCode PetscWeakFormGetBdResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f,
345:                                           PetscInt *n0,
346:                                         void (***f0)(PetscInt, PetscInt, PetscInt,
347:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
348:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
349:                                                      PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
350:                                         PetscInt *n1,
351:                                         void (***f1)(PetscInt, PetscInt, PetscInt,
352:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
353:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
354:                                                      PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
355: {

359:   PetscWeakFormGetFunction_Private(wf, wf->bdf0, label, val, f, n0, (void (***)(void)) f0);
360:   PetscWeakFormGetFunction_Private(wf, wf->bdf1, label, val, f, n1, (void (***)(void)) f1);
361:   return(0);
362: }

364: PetscErrorCode PetscWeakFormAddBdResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f,
365:                                           void (*f0)(PetscInt, PetscInt, PetscInt,
366:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
367:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
368:                                                      PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
369:                                           void (*f1)(PetscInt, PetscInt, PetscInt,
370:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
371:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
372:                                                      PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
373: {

377:   PetscWeakFormAddFunction_Private(wf, wf->bdf0, label, val, f, (void (*)(void)) f0);
378:   PetscWeakFormAddFunction_Private(wf, wf->bdf1, label, val, f, (void (*)(void)) f1);
379:   return(0);
380: }

382: PetscErrorCode PetscWeakFormSetBdResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f,
383:                                           PetscInt n0,
384:                                           void (**f0)(PetscInt, PetscInt, PetscInt,
385:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
386:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
387:                                                      PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
388:                                           PetscInt n1,
389:                                           void (**f1)(PetscInt, PetscInt, PetscInt,
390:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
391:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
392:                                                      PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
393: {

397:   PetscWeakFormSetFunction_Private(wf, wf->bdf0, label, val, f, n0, (void (**)(void)) f0);
398:   PetscWeakFormSetFunction_Private(wf, wf->bdf1, label, val, f, n1, (void (**)(void)) f1);
399:   return(0);
400: }

402: PetscErrorCode PetscWeakFormSetIndexBdResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f,
403:                                           PetscInt i0,
404:                                           void (*f0)(PetscInt, PetscInt, PetscInt,
405:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
406:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
407:                                                      PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
408:                                           PetscInt i1,
409:                                           void (*f1)(PetscInt, PetscInt, PetscInt,
410:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
411:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
412:                                                      PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
413: {

417:   PetscWeakFormSetIndexFunction_Private(wf, wf->bdf0, label, val, f, i0, (void (*)(void)) f0);
418:   PetscWeakFormSetIndexFunction_Private(wf, wf->bdf1, label, val, f, i1, (void (*)(void)) f1);
419:   return(0);
420: }

422: PetscErrorCode PetscWeakFormHasJacobian(PetscWeakForm wf, PetscBool *hasJac)
423: {
424:   PetscInt       n0, n1, n2, n3;

430:   PetscHMapFormGetSize(wf->g0, &n0);
431:   PetscHMapFormGetSize(wf->g1, &n1);
432:   PetscHMapFormGetSize(wf->g2, &n2);
433:   PetscHMapFormGetSize(wf->g3, &n3);
434:   *hasJac = n0+n1+n2+n3 ? PETSC_TRUE : PETSC_FALSE;
435:   return(0);
436: }

438: PetscErrorCode PetscWeakFormGetJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g,
439:                                         PetscInt *n0,
440:                                         void (***g0)(PetscInt, PetscInt, PetscInt,
441:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
442:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
443:                                                      PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
444:                                         PetscInt *n1,
445:                                         void (***g1)(PetscInt, PetscInt, PetscInt,
446:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
447:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
448:                                                      PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
449:                                         PetscInt *n2,
450:                                         void (***g2)(PetscInt, PetscInt, PetscInt,
451:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
452:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
453:                                                      PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
454:                                         PetscInt *n3,
455:                                         void (***g3)(PetscInt, PetscInt, PetscInt,
456:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
457:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
458:                                                      PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
459: {
460:   PetscInt       find = f*wf->Nf + g;

464:   PetscWeakFormGetFunction_Private(wf, wf->g0, label, val, find, n0, (void (***)(void)) g0);
465:   PetscWeakFormGetFunction_Private(wf, wf->g1, label, val, find, n1, (void (***)(void)) g1);
466:   PetscWeakFormGetFunction_Private(wf, wf->g2, label, val, find, n2, (void (***)(void)) g2);
467:   PetscWeakFormGetFunction_Private(wf, wf->g3, label, val, find, n3, (void (***)(void)) g3);
468:   return(0);
469: }

471: PetscErrorCode PetscWeakFormAddJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g,
472:                                         void (*g0)(PetscInt, PetscInt, PetscInt,
473:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
474:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
475:                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
476:                                         void (*g1)(PetscInt, PetscInt, PetscInt,
477:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
478:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
479:                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
480:                                         void (*g2)(PetscInt, PetscInt, PetscInt,
481:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
482:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
483:                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
484:                                         void (*g3)(PetscInt, PetscInt, PetscInt,
485:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
486:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
487:                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
488: {
489:   PetscInt       find = f*wf->Nf + g;

493:   PetscWeakFormAddFunction_Private(wf, wf->g0, label, val, find, (void (*)(void)) g0);
494:   PetscWeakFormAddFunction_Private(wf, wf->g1, label, val, find, (void (*)(void)) g1);
495:   PetscWeakFormAddFunction_Private(wf, wf->g2, label, val, find, (void (*)(void)) g2);
496:   PetscWeakFormAddFunction_Private(wf, wf->g3, label, val, find, (void (*)(void)) g3);
497:   return(0);
498: }

500: PetscErrorCode PetscWeakFormSetJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g,
501:                                         PetscInt n0,
502:                                         void (**g0)(PetscInt, PetscInt, PetscInt,
503:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
504:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
505:                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
506:                                         PetscInt n1,
507:                                         void (**g1)(PetscInt, PetscInt, PetscInt,
508:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
509:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
510:                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
511:                                         PetscInt n2,
512:                                         void (**g2)(PetscInt, PetscInt, PetscInt,
513:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
514:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
515:                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
516:                                         PetscInt n3,
517:                                         void (**g3)(PetscInt, PetscInt, PetscInt,
518:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
519:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
520:                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
521: {
522:   PetscInt       find = f*wf->Nf + g;

526:   PetscWeakFormSetFunction_Private(wf, wf->g0, label, val, find, n0, (void (**)(void)) g0);
527:   PetscWeakFormSetFunction_Private(wf, wf->g1, label, val, find, n1, (void (**)(void)) g1);
528:   PetscWeakFormSetFunction_Private(wf, wf->g2, label, val, find, n2, (void (**)(void)) g2);
529:   PetscWeakFormSetFunction_Private(wf, wf->g3, label, val, find, n3, (void (**)(void)) g3);
530:   return(0);
531: }

533: PetscErrorCode PetscWeakFormSetIndexJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g,
534:                                         PetscInt i0,
535:                                         void (*g0)(PetscInt, PetscInt, PetscInt,
536:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
537:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
538:                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
539:                                         PetscInt i1,
540:                                         void (*g1)(PetscInt, PetscInt, PetscInt,
541:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
542:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
543:                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
544:                                         PetscInt i2,
545:                                         void (*g2)(PetscInt, PetscInt, PetscInt,
546:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
547:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
548:                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
549:                                         PetscInt i3,
550:                                         void (*g3)(PetscInt, PetscInt, PetscInt,
551:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
552:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
553:                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
554: {
555:   PetscInt       find = f*wf->Nf + g;

559:   PetscWeakFormSetIndexFunction_Private(wf, wf->g0, label, val, find, i0, (void (*)(void)) g0);
560:   PetscWeakFormSetIndexFunction_Private(wf, wf->g1, label, val, find, i1, (void (*)(void)) g1);
561:   PetscWeakFormSetIndexFunction_Private(wf, wf->g2, label, val, find, i2, (void (*)(void)) g2);
562:   PetscWeakFormSetIndexFunction_Private(wf, wf->g3, label, val, find, i3, (void (*)(void)) g3);
563:   return(0);
564: }

566: PetscErrorCode PetscWeakFormHasJacobianPreconditioner(PetscWeakForm wf, PetscBool *hasJacPre)
567: {
568:   PetscInt       n0, n1, n2, n3;

574:   PetscHMapFormGetSize(wf->gp0, &n0);
575:   PetscHMapFormGetSize(wf->gp1, &n1);
576:   PetscHMapFormGetSize(wf->gp2, &n2);
577:   PetscHMapFormGetSize(wf->gp3, &n3);
578:   *hasJacPre = n0+n1+n2+n3 ? PETSC_TRUE : PETSC_FALSE;
579:   return(0);
580: }

582: PetscErrorCode PetscWeakFormGetJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g,
583:                                                       PetscInt *n0,
584:                                                       void (***g0)(PetscInt, PetscInt, PetscInt,
585:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
586:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
587:                                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
588:                                                       PetscInt *n1,
589:                                                       void (***g1)(PetscInt, PetscInt, PetscInt,
590:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
591:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
592:                                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
593:                                                       PetscInt *n2,
594:                                                       void (***g2)(PetscInt, PetscInt, PetscInt,
595:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
596:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
597:                                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
598:                                                       PetscInt *n3,
599:                                                       void (***g3)(PetscInt, PetscInt, PetscInt,
600:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
601:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
602:                                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
603: {
604:   PetscInt       find = f*wf->Nf + g;

608:   PetscWeakFormGetFunction_Private(wf, wf->gp0, label, val, find, n0, (void (***)(void)) g0);
609:   PetscWeakFormGetFunction_Private(wf, wf->gp1, label, val, find, n1, (void (***)(void)) g1);
610:   PetscWeakFormGetFunction_Private(wf, wf->gp2, label, val, find, n2, (void (***)(void)) g2);
611:   PetscWeakFormGetFunction_Private(wf, wf->gp3, label, val, find, n3, (void (***)(void)) g3);
612:   return(0);
613: }

615: PetscErrorCode PetscWeakFormAddJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g,
616:                                         void (*g0)(PetscInt, PetscInt, PetscInt,
617:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
618:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
619:                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
620:                                         void (*g1)(PetscInt, PetscInt, PetscInt,
621:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
622:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
623:                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
624:                                         void (*g2)(PetscInt, PetscInt, PetscInt,
625:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
626:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
627:                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
628:                                         void (*g3)(PetscInt, PetscInt, PetscInt,
629:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
630:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
631:                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
632: {
633:   PetscInt       find = f*wf->Nf + g;

637:   PetscWeakFormAddFunction_Private(wf, wf->gp0, label, val, find, (void (*)(void)) g0);
638:   PetscWeakFormAddFunction_Private(wf, wf->gp1, label, val, find, (void (*)(void)) g1);
639:   PetscWeakFormAddFunction_Private(wf, wf->gp2, label, val, find, (void (*)(void)) g2);
640:   PetscWeakFormAddFunction_Private(wf, wf->gp3, label, val, find, (void (*)(void)) g3);
641:   return(0);
642: }

644: PetscErrorCode PetscWeakFormSetJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g,
645:                                                       PetscInt n0,
646:                                                       void (**g0)(PetscInt, PetscInt, PetscInt,
647:                                                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
648:                                                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
649:                                                                  PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
650:                                                       PetscInt n1,
651:                                                       void (**g1)(PetscInt, PetscInt, PetscInt,
652:                                                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
653:                                                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
654:                                                                  PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
655:                                                       PetscInt n2,
656:                                                       void (**g2)(PetscInt, PetscInt, PetscInt,
657:                                                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
658:                                                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
659:                                                                  PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
660:                                                       PetscInt n3,
661:                                                       void (**g3)(PetscInt, PetscInt, PetscInt,
662:                                                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
663:                                                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
664:                                                                  PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
665: {
666:   PetscInt       find = f*wf->Nf + g;

670:   PetscWeakFormSetFunction_Private(wf, wf->gp0, label, val, find, n0, (void (**)(void)) g0);
671:   PetscWeakFormSetFunction_Private(wf, wf->gp1, label, val, find, n1, (void (**)(void)) g1);
672:   PetscWeakFormSetFunction_Private(wf, wf->gp2, label, val, find, n2, (void (**)(void)) g2);
673:   PetscWeakFormSetFunction_Private(wf, wf->gp3, label, val, find, n3, (void (**)(void)) g3);
674:   return(0);
675: }

677: PetscErrorCode PetscWeakFormSetIndexJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g,
678:                                                       PetscInt i0,
679:                                                       void (*g0)(PetscInt, PetscInt, PetscInt,
680:                                                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
681:                                                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
682:                                                                  PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
683:                                                       PetscInt i1,
684:                                                       void (*g1)(PetscInt, PetscInt, PetscInt,
685:                                                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
686:                                                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
687:                                                                  PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
688:                                                       PetscInt i2,
689:                                                       void (*g2)(PetscInt, PetscInt, PetscInt,
690:                                                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
691:                                                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
692:                                                                  PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
693:                                                       PetscInt i3,
694:                                                       void (*g3)(PetscInt, PetscInt, PetscInt,
695:                                                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
696:                                                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
697:                                                                  PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
698: {
699:   PetscInt       find = f*wf->Nf + g;

703:   PetscWeakFormSetIndexFunction_Private(wf, wf->gp0, label, val, find, i0, (void (*)(void)) g0);
704:   PetscWeakFormSetIndexFunction_Private(wf, wf->gp1, label, val, find, i1, (void (*)(void)) g1);
705:   PetscWeakFormSetIndexFunction_Private(wf, wf->gp2, label, val, find, i2, (void (*)(void)) g2);
706:   PetscWeakFormSetIndexFunction_Private(wf, wf->gp3, label, val, find, i3, (void (*)(void)) g3);
707:   return(0);
708: }

710: PetscErrorCode PetscWeakFormHasBdJacobian(PetscWeakForm wf, PetscBool *hasJac)
711: {
712:   PetscInt       n0, n1, n2, n3;

718:   PetscHMapFormGetSize(wf->bdg0, &n0);
719:   PetscHMapFormGetSize(wf->bdg1, &n1);
720:   PetscHMapFormGetSize(wf->bdg2, &n2);
721:   PetscHMapFormGetSize(wf->bdg3, &n3);
722:   *hasJac = n0+n1+n2+n3 ? PETSC_TRUE : PETSC_FALSE;
723:   return(0);
724: }

726: PetscErrorCode PetscWeakFormGetBdJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g,
727:                                           PetscInt *n0,
728:                                           void (***g0)(PetscInt, PetscInt, PetscInt,
729:                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
730:                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
731:                                                        PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
732:                                           PetscInt *n1,
733:                                           void (***g1)(PetscInt, PetscInt, PetscInt,
734:                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
735:                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
736:                                                        PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
737:                                           PetscInt *n2,
738:                                           void (***g2)(PetscInt, PetscInt, PetscInt,
739:                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
740:                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
741:                                                        PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
742:                                           PetscInt *n3,
743:                                           void (***g3)(PetscInt, PetscInt, PetscInt,
744:                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
745:                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
746:                                                        PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
747: {
748:   PetscInt       find = f*wf->Nf + g;

752:   PetscWeakFormGetFunction_Private(wf, wf->bdg0, label, val, find, n0, (void (***)(void)) g0);
753:   PetscWeakFormGetFunction_Private(wf, wf->bdg1, label, val, find, n1, (void (***)(void)) g1);
754:   PetscWeakFormGetFunction_Private(wf, wf->bdg2, label, val, find, n2, (void (***)(void)) g2);
755:   PetscWeakFormGetFunction_Private(wf, wf->bdg3, label, val, find, n3, (void (***)(void)) g3);
756:   return(0);
757: }

759: PetscErrorCode PetscWeakFormAddBdJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g,
760:                                           void (*g0)(PetscInt, PetscInt, PetscInt,
761:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
762:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
763:                                                      PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
764:                                           void (*g1)(PetscInt, PetscInt, PetscInt,
765:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
766:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
767:                                                      PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
768:                                           void (*g2)(PetscInt, PetscInt, PetscInt,
769:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
770:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
771:                                                      PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
772:                                           void (*g3)(PetscInt, PetscInt, PetscInt,
773:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
774:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
775:                                                      PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
776: {
777:   PetscInt       find = f*wf->Nf + g;

781:   PetscWeakFormAddFunction_Private(wf, wf->bdg0, label, val, find, (void (*)(void)) g0);
782:   PetscWeakFormAddFunction_Private(wf, wf->bdg1, label, val, find, (void (*)(void)) g1);
783:   PetscWeakFormAddFunction_Private(wf, wf->bdg2, label, val, find, (void (*)(void)) g2);
784:   PetscWeakFormAddFunction_Private(wf, wf->bdg3, label, val, find, (void (*)(void)) g3);
785:   return(0);
786: }

788: PetscErrorCode PetscWeakFormSetBdJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g,
789:                                           PetscInt n0,
790:                                           void (**g0)(PetscInt, PetscInt, PetscInt,
791:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
792:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
793:                                                      PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
794:                                           PetscInt n1,
795:                                           void (**g1)(PetscInt, PetscInt, PetscInt,
796:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
797:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
798:                                                      PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
799:                                           PetscInt n2,
800:                                           void (**g2)(PetscInt, PetscInt, PetscInt,
801:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
802:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
803:                                                      PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
804:                                           PetscInt n3,
805:                                           void (**g3)(PetscInt, PetscInt, PetscInt,
806:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
807:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
808:                                                      PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
809: {
810:   PetscInt       find = f*wf->Nf + g;

814:   PetscWeakFormSetFunction_Private(wf, wf->bdg0, label, val, find, n0, (void (**)(void)) g0);
815:   PetscWeakFormSetFunction_Private(wf, wf->bdg1, label, val, find, n1, (void (**)(void)) g1);
816:   PetscWeakFormSetFunction_Private(wf, wf->bdg2, label, val, find, n2, (void (**)(void)) g2);
817:   PetscWeakFormSetFunction_Private(wf, wf->bdg3, label, val, find, n3, (void (**)(void)) g3);
818:   return(0);
819: }

821: PetscErrorCode PetscWeakFormSetIndexBdJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g,
822:                                           PetscInt i0,
823:                                           void (*g0)(PetscInt, PetscInt, PetscInt,
824:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
825:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
826:                                                      PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
827:                                           PetscInt i1,
828:                                           void (*g1)(PetscInt, PetscInt, PetscInt,
829:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
830:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
831:                                                      PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
832:                                           PetscInt i2,
833:                                           void (*g2)(PetscInt, PetscInt, PetscInt,
834:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
835:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
836:                                                      PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
837:                                           PetscInt i3,
838:                                           void (*g3)(PetscInt, PetscInt, PetscInt,
839:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
840:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
841:                                                      PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
842: {
843:   PetscInt       find = f*wf->Nf + g;

847:   PetscWeakFormSetIndexFunction_Private(wf, wf->bdg0, label, val, find, i0, (void (*)(void)) g0);
848:   PetscWeakFormSetIndexFunction_Private(wf, wf->bdg1, label, val, find, i1, (void (*)(void)) g1);
849:   PetscWeakFormSetIndexFunction_Private(wf, wf->bdg2, label, val, find, i2, (void (*)(void)) g2);
850:   PetscWeakFormSetIndexFunction_Private(wf, wf->bdg3, label, val, find, i3, (void (*)(void)) g3);
851:   return(0);
852: }

854: PetscErrorCode PetscWeakFormHasBdJacobianPreconditioner(PetscWeakForm wf, PetscBool *hasJacPre)
855: {
856:   PetscInt       n0, n1, n2, n3;

862:   PetscHMapFormGetSize(wf->bdgp0, &n0);
863:   PetscHMapFormGetSize(wf->bdgp1, &n1);
864:   PetscHMapFormGetSize(wf->bdgp2, &n2);
865:   PetscHMapFormGetSize(wf->bdgp3, &n3);
866:   *hasJacPre = n0+n1+n2+n3 ? PETSC_TRUE : PETSC_FALSE;
867:   return(0);
868: }

870: PetscErrorCode PetscWeakFormGetBdJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g,
871:                                                         PetscInt *n0,
872:                                                         void (***g0)(PetscInt, PetscInt, PetscInt,
873:                                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
874:                                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
875:                                                                      PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
876:                                                         PetscInt *n1,
877:                                                         void (***g1)(PetscInt, PetscInt, PetscInt,
878:                                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
879:                                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
880:                                                                      PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
881:                                                         PetscInt *n2,
882:                                                         void (***g2)(PetscInt, PetscInt, PetscInt,
883:                                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
884:                                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
885:                                                                      PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
886:                                                         PetscInt *n3,
887:                                                         void (***g3)(PetscInt, PetscInt, PetscInt,
888:                                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
889:                                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
890:                                                                      PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
891: {
892:   PetscInt       find = f*wf->Nf + g;

896:   PetscWeakFormGetFunction_Private(wf, wf->bdgp0, label, val, find, n0, (void (***)(void)) g0);
897:   PetscWeakFormGetFunction_Private(wf, wf->bdgp1, label, val, find, n1, (void (***)(void)) g1);
898:   PetscWeakFormGetFunction_Private(wf, wf->bdgp2, label, val, find, n2, (void (***)(void)) g2);
899:   PetscWeakFormGetFunction_Private(wf, wf->bdgp3, label, val, find, n3, (void (***)(void)) g3);
900:   return(0);
901: }

903: PetscErrorCode PetscWeakFormAddBdJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g,
904:                                                         void (*g0)(PetscInt, PetscInt, PetscInt,
905:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
906:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
907:                                                                    PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
908:                                                         void (*g1)(PetscInt, PetscInt, PetscInt,
909:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
910:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
911:                                                                    PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
912:                                                         void (*g2)(PetscInt, PetscInt, PetscInt,
913:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
914:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
915:                                                                    PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
916:                                                         void (*g3)(PetscInt, PetscInt, PetscInt,
917:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
918:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
919:                                                                    PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
920: {
921:   PetscInt       find = f*wf->Nf + g;

925:   PetscWeakFormAddFunction_Private(wf, wf->bdgp0, label, val, find, (void (*)(void)) g0);
926:   PetscWeakFormAddFunction_Private(wf, wf->bdgp1, label, val, find, (void (*)(void)) g1);
927:   PetscWeakFormAddFunction_Private(wf, wf->bdgp2, label, val, find, (void (*)(void)) g2);
928:   PetscWeakFormAddFunction_Private(wf, wf->bdgp3, label, val, find, (void (*)(void)) g3);
929:   return(0);
930: }

932: PetscErrorCode PetscWeakFormSetBdJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g,
933:                                                         PetscInt n0,
934:                                                         void (**g0)(PetscInt, PetscInt, PetscInt,
935:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
936:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
937:                                                                    PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
938:                                                         PetscInt n1,
939:                                                         void (**g1)(PetscInt, PetscInt, PetscInt,
940:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
941:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
942:                                                                    PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
943:                                                         PetscInt n2,
944:                                                         void (**g2)(PetscInt, PetscInt, PetscInt,
945:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
946:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
947:                                                                    PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
948:                                                         PetscInt n3,
949:                                                         void (**g3)(PetscInt, PetscInt, PetscInt,
950:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
951:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
952:                                                                    PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
953: {
954:   PetscInt       find = f*wf->Nf + g;

958:   PetscWeakFormSetFunction_Private(wf, wf->bdgp0, label, val, find, n0, (void (**)(void)) g0);
959:   PetscWeakFormSetFunction_Private(wf, wf->bdgp1, label, val, find, n1, (void (**)(void)) g1);
960:   PetscWeakFormSetFunction_Private(wf, wf->bdgp2, label, val, find, n2, (void (**)(void)) g2);
961:   PetscWeakFormSetFunction_Private(wf, wf->bdgp3, label, val, find, n3, (void (**)(void)) g3);
962:   return(0);
963: }

965: PetscErrorCode PetscWeakFormSetIndexBdJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g,
966:                                                         PetscInt i0,
967:                                                         void (*g0)(PetscInt, PetscInt, PetscInt,
968:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
969:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
970:                                                                    PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
971:                                                         PetscInt i1,
972:                                                         void (*g1)(PetscInt, PetscInt, PetscInt,
973:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
974:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
975:                                                                    PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
976:                                                         PetscInt i2,
977:                                                         void (*g2)(PetscInt, PetscInt, PetscInt,
978:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
979:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
980:                                                                    PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
981:                                                         PetscInt i3,
982:                                                         void (*g3)(PetscInt, PetscInt, PetscInt,
983:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
984:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
985:                                                                    PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
986: {
987:   PetscInt       find = f*wf->Nf + g;

991:   PetscWeakFormSetIndexFunction_Private(wf, wf->bdgp0, label, val, find, i0, (void (*)(void)) g0);
992:   PetscWeakFormSetIndexFunction_Private(wf, wf->bdgp1, label, val, find, i1, (void (*)(void)) g1);
993:   PetscWeakFormSetIndexFunction_Private(wf, wf->bdgp2, label, val, find, i2, (void (*)(void)) g2);
994:   PetscWeakFormSetIndexFunction_Private(wf, wf->bdgp3, label, val, find, i3, (void (*)(void)) g3);
995:   return(0);
996: }

998: PetscErrorCode PetscWeakFormHasDynamicJacobian(PetscWeakForm wf, PetscBool *hasDynJac)
999: {
1000:   PetscInt       n0, n1, n2, n3;

1006:   PetscHMapFormGetSize(wf->gt0, &n0);
1007:   PetscHMapFormGetSize(wf->gt1, &n1);
1008:   PetscHMapFormGetSize(wf->gt2, &n2);
1009:   PetscHMapFormGetSize(wf->gt3, &n3);
1010:   *hasDynJac = n0+n1+n2+n3 ? PETSC_TRUE : PETSC_FALSE;
1011:   return(0);
1012: }

1014: PetscErrorCode PetscWeakFormGetDynamicJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g,
1015:                                         PetscInt *n0,
1016:                                         void (***g0)(PetscInt, PetscInt, PetscInt,
1017:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1018:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1019:                                                      PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1020:                                         PetscInt *n1,
1021:                                         void (***g1)(PetscInt, PetscInt, PetscInt,
1022:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1023:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1024:                                                      PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1025:                                         PetscInt *n2,
1026:                                         void (***g2)(PetscInt, PetscInt, PetscInt,
1027:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1028:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1029:                                                      PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1030:                                         PetscInt *n3,
1031:                                         void (***g3)(PetscInt, PetscInt, PetscInt,
1032:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1033:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1034:                                                      PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
1035: {
1036:   PetscInt       find = f*wf->Nf + g;

1040:   PetscWeakFormGetFunction_Private(wf, wf->gt0, label, val, find, n0, (void (***)(void)) g0);
1041:   PetscWeakFormGetFunction_Private(wf, wf->gt1, label, val, find, n1, (void (***)(void)) g1);
1042:   PetscWeakFormGetFunction_Private(wf, wf->gt2, label, val, find, n2, (void (***)(void)) g2);
1043:   PetscWeakFormGetFunction_Private(wf, wf->gt3, label, val, find, n3, (void (***)(void)) g3);
1044:   return(0);
1045: }

1047: PetscErrorCode PetscWeakFormAddDynamicJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g,
1048:                                         void (*g0)(PetscInt, PetscInt, PetscInt,
1049:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1050:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1051:                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1052:                                         void (*g1)(PetscInt, PetscInt, PetscInt,
1053:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1054:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1055:                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1056:                                         void (*g2)(PetscInt, PetscInt, PetscInt,
1057:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1058:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1059:                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1060:                                         void (*g3)(PetscInt, PetscInt, PetscInt,
1061:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1062:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1063:                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
1064: {
1065:   PetscInt       find = f*wf->Nf + g;

1069:   PetscWeakFormAddFunction_Private(wf, wf->gt0, label, val, find, (void (*)(void)) g0);
1070:   PetscWeakFormAddFunction_Private(wf, wf->gt1, label, val, find, (void (*)(void)) g1);
1071:   PetscWeakFormAddFunction_Private(wf, wf->gt2, label, val, find, (void (*)(void)) g2);
1072:   PetscWeakFormAddFunction_Private(wf, wf->gt3, label, val, find, (void (*)(void)) g3);
1073:   return(0);
1074: }

1076: PetscErrorCode PetscWeakFormSetDynamicJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g,
1077:                                                PetscInt n0,
1078:                                                void (**g0)(PetscInt, PetscInt, PetscInt,
1079:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1080:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1081:                                                           PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1082:                                                PetscInt n1,
1083:                                                void (**g1)(PetscInt, PetscInt, PetscInt,
1084:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1085:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1086:                                                           PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1087:                                                PetscInt n2,
1088:                                                void (**g2)(PetscInt, PetscInt, PetscInt,
1089:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1090:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1091:                                                           PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1092:                                                PetscInt n3,
1093:                                                void (**g3)(PetscInt, PetscInt, PetscInt,
1094:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1095:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1096:                                                           PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
1097: {
1098:   PetscInt       find = f*wf->Nf + g;

1102:   PetscWeakFormSetFunction_Private(wf, wf->gt0, label, val, find, n0, (void (**)(void)) g0);
1103:   PetscWeakFormSetFunction_Private(wf, wf->gt1, label, val, find, n1, (void (**)(void)) g1);
1104:   PetscWeakFormSetFunction_Private(wf, wf->gt2, label, val, find, n2, (void (**)(void)) g2);
1105:   PetscWeakFormSetFunction_Private(wf, wf->gt3, label, val, find, n3, (void (**)(void)) g3);
1106:   return(0);
1107: }

1109: PetscErrorCode PetscWeakFormSetIndexDynamicJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g,
1110:                                                PetscInt i0,
1111:                                                void (*g0)(PetscInt, PetscInt, PetscInt,
1112:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1113:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1114:                                                           PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1115:                                                PetscInt i1,
1116:                                                void (*g1)(PetscInt, PetscInt, PetscInt,
1117:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1118:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1119:                                                           PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1120:                                                PetscInt i2,
1121:                                                void (*g2)(PetscInt, PetscInt, PetscInt,
1122:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1123:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1124:                                                           PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1125:                                                PetscInt i3,
1126:                                                void (*g3)(PetscInt, PetscInt, PetscInt,
1127:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1128:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1129:                                                           PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
1130: {
1131:   PetscInt       find = f*wf->Nf + g;

1135:   PetscWeakFormSetIndexFunction_Private(wf, wf->gt0, label, val, find, i0, (void (*)(void)) g0);
1136:   PetscWeakFormSetIndexFunction_Private(wf, wf->gt1, label, val, find, i1, (void (*)(void)) g1);
1137:   PetscWeakFormSetIndexFunction_Private(wf, wf->gt2, label, val, find, i2, (void (*)(void)) g2);
1138:   PetscWeakFormSetIndexFunction_Private(wf, wf->gt3, label, val, find, i3, (void (*)(void)) g3);
1139:   return(0);
1140: }

1142: PetscErrorCode PetscWeakFormGetRiemannSolver(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt *n,
1143:                                              void (***r)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *))
1144: {

1148:   PetscWeakFormGetFunction_Private(wf, wf->r, label, val, f, n, (void (***)(void)) r);
1149:   return(0);
1150: }

1152: PetscErrorCode PetscWeakFormSetRiemannSolver(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f,
1153:                                              PetscInt n,
1154:                                              void (**r)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *))
1155: {

1159:   PetscWeakFormSetFunction_Private(wf, wf->r, label, val, f, n, (void (**)(void)) r);
1160:   return(0);
1161: }

1163: PetscErrorCode PetscWeakFormSetIndexRiemannSolver(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f,
1164:                                                   PetscInt i,
1165:                                                   void (*r)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *))
1166: {

1170:   PetscWeakFormSetIndexFunction_Private(wf, wf->r, label, val, f, i, (void (*)(void)) r);
1171:   return(0);
1172: }

1174: /*@
1175:   PetscWeakFormGetNumFields - Returns the number of fields

1177:   Not collective

1179:   Input Parameter:
1180: . wf - The PetscWeakForm object

1182:   Output Parameter:
1183: . Nf - The nubmer of fields

1185:   Level: beginner

1187: .seealso: PetscWeakFormSetNumFields(), PetscWeakFormCreate()
1188: @*/
1189: PetscErrorCode PetscWeakFormGetNumFields(PetscWeakForm wf, PetscInt *Nf)
1190: {
1194:   *Nf = wf->Nf;
1195:   return(0);
1196: }

1198: /*@
1199:   PetscWeakFormSetNumFields - Sets the number of fields

1201:   Not collective

1203:   Input Parameters:
1204: + wf - The PetscWeakForm object
1205: - Nf - The number of fields

1207:   Level: beginner

1209: .seealso: PetscWeakFormGetNumFields(), PetscWeakFormCreate()
1210: @*/
1211: PetscErrorCode PetscWeakFormSetNumFields(PetscWeakForm wf, PetscInt Nf)
1212: {
1215:   wf->Nf = Nf;
1216:   return(0);
1217: }

1219: /*@
1220:   PetscWeakFormDestroy - Destroys a PetscWeakForm object

1222:   Collective on wf

1224:   Input Parameter:
1225: . wf - the PetscWeakForm object to destroy

1227:   Level: developer

1229: .seealso PetscWeakFormCreate(), PetscWeakFormView()
1230: @*/
1231: PetscErrorCode PetscWeakFormDestroy(PetscWeakForm *wf)
1232: {

1236:   if (!*wf) return(0);

1239:   if (--((PetscObject)(*wf))->refct > 0) {*wf = NULL; return(0);}
1240:   ((PetscObject) (*wf))->refct = 0;
1241:   PetscChunkBufferDestroy(&(*wf)->funcs);
1242:   PetscHMapFormDestroy(&(*wf)->obj);
1243:   PetscHMapFormDestroy(&(*wf)->f0);
1244:   PetscHMapFormDestroy(&(*wf)->f1);
1245:   PetscHMapFormDestroy(&(*wf)->g0);
1246:   PetscHMapFormDestroy(&(*wf)->g1);
1247:   PetscHMapFormDestroy(&(*wf)->g2);
1248:   PetscHMapFormDestroy(&(*wf)->g3);
1249:   PetscHMapFormDestroy(&(*wf)->gp0);
1250:   PetscHMapFormDestroy(&(*wf)->gp1);
1251:   PetscHMapFormDestroy(&(*wf)->gp2);
1252:   PetscHMapFormDestroy(&(*wf)->gp3);
1253:   PetscHMapFormDestroy(&(*wf)->gt0);
1254:   PetscHMapFormDestroy(&(*wf)->gt1);
1255:   PetscHMapFormDestroy(&(*wf)->gt2);
1256:   PetscHMapFormDestroy(&(*wf)->gt3);
1257:   PetscHMapFormDestroy(&(*wf)->bdf0);
1258:   PetscHMapFormDestroy(&(*wf)->bdf1);
1259:   PetscHMapFormDestroy(&(*wf)->bdg0);
1260:   PetscHMapFormDestroy(&(*wf)->bdg1);
1261:   PetscHMapFormDestroy(&(*wf)->bdg2);
1262:   PetscHMapFormDestroy(&(*wf)->bdg3);
1263:   PetscHMapFormDestroy(&(*wf)->bdgp0);
1264:   PetscHMapFormDestroy(&(*wf)->bdgp1);
1265:   PetscHMapFormDestroy(&(*wf)->bdgp2);
1266:   PetscHMapFormDestroy(&(*wf)->bdgp3);
1267:   PetscHMapFormDestroy(&(*wf)->r);
1268:   PetscHeaderDestroy(wf);
1269:   return(0);
1270: }

1272: static PetscErrorCode PetscWeakFormViewTable_Ascii(PetscWeakForm wf, PetscViewer viewer, const char tableName[], PetscHMapForm map)
1273: {
1274:   PetscInt       Nk, k;

1278:   PetscHMapFormGetSize(map, &Nk);
1279:   if (Nk) {
1280:     PetscHashFormKey *keys;
1281:     void           (**funcs)(void);
1282:     const char       *name;
1283:     PetscInt          off = 0, n, i;

1285:     PetscMalloc1(Nk, &keys);
1286:     PetscHMapFormGetKeys(map, &off, keys);
1287:     PetscViewerASCIIPrintf(viewer, "%s\n", tableName);
1288:     PetscViewerASCIIPushTab(viewer);
1289:     for (k = 0; k < Nk; ++k) {
1290:       if (keys[k].label) {PetscObjectGetName((PetscObject) keys[k].label, &name);}
1291:       PetscViewerASCIIPrintf(viewer, "Key (%s, %D, %D) ", keys[k].label ? name : "None", keys[k].value, keys[k].field);
1292:       PetscViewerASCIIUseTabs(viewer, PETSC_FALSE);
1293:       PetscWeakFormGetFunction_Private(wf, map, keys[k].label, keys[k].value, keys[k].field, &n, &funcs);
1294:       for (i = 0; i < n; ++i) {
1295:         if (i > 0) {PetscViewerASCIIPrintf(viewer, ", ");}
1296:         PetscDLAddr(funcs[i], &name);
1297:         if (name) {PetscViewerASCIIPrintf(viewer, "%s", name);}
1298:         else      {PetscViewerASCIIPrintf(viewer, "%p", funcs[i]);}
1299:       }
1300:       PetscViewerASCIIPrintf(viewer, "\n");
1301:       PetscViewerASCIIUseTabs(viewer, PETSC_TRUE);
1302:     }
1303:     PetscViewerASCIIPopTab(viewer);
1304:     PetscFree(keys);
1305:   }
1306:   return(0);
1307: }

1309: static PetscErrorCode PetscWeakFormView_Ascii(PetscWeakForm wf, PetscViewer viewer)
1310: {
1311:   PetscViewerFormat format;
1312:   PetscErrorCode    ierr;

1315:   PetscViewerGetFormat(viewer, &format);
1316:   PetscViewerASCIIPrintf(viewer, "Weak Form System with %d fields\n", wf->Nf);
1317:   PetscViewerASCIIPushTab(viewer);
1318:   PetscWeakFormViewTable_Ascii(wf, viewer, "Objective", wf->obj);
1319:   PetscWeakFormViewTable_Ascii(wf, viewer, "Residual f0", wf->f0);
1320:   PetscWeakFormViewTable_Ascii(wf, viewer, "Residual f1", wf->f1);
1321:   PetscWeakFormViewTable_Ascii(wf, viewer, "Boudnary Residual f0", wf->bdf0);
1322:   PetscWeakFormViewTable_Ascii(wf, viewer, "Boundary Residual f1", wf->bdf1);
1323:   PetscWeakFormViewTable_Ascii(wf, viewer, "Jacobian g0", wf->g0);
1324:   PetscWeakFormViewTable_Ascii(wf, viewer, "Jacobian g1", wf->g1);
1325:   PetscWeakFormViewTable_Ascii(wf, viewer, "Jacobian g2", wf->g2);
1326:   PetscWeakFormViewTable_Ascii(wf, viewer, "Jacobian g3", wf->g3);
1327:   PetscWeakFormViewTable_Ascii(wf, viewer, "Jacobian Preconditioner g0", wf->gp0);
1328:   PetscWeakFormViewTable_Ascii(wf, viewer, "Jacobian Preconditioner g1", wf->gp1);
1329:   PetscWeakFormViewTable_Ascii(wf, viewer, "Jacobian Preconditioner g2", wf->gp2);
1330:   PetscWeakFormViewTable_Ascii(wf, viewer, "Jacobian Preconditioner g3", wf->gp3);
1331:   PetscWeakFormViewTable_Ascii(wf, viewer, "Dynamic Jacobian g0", wf->gt0);
1332:   PetscWeakFormViewTable_Ascii(wf, viewer, "Dynamic Jacobian g1", wf->gt1);
1333:   PetscWeakFormViewTable_Ascii(wf, viewer, "Dynamic Jacobian g2", wf->gt2);
1334:   PetscWeakFormViewTable_Ascii(wf, viewer, "Dynamic Jacobian g3", wf->gt3);
1335:   PetscWeakFormViewTable_Ascii(wf, viewer, "Boundary Jacobian g0", wf->bdg0);
1336:   PetscWeakFormViewTable_Ascii(wf, viewer, "Boundary Jacobian g1", wf->bdg1);
1337:   PetscWeakFormViewTable_Ascii(wf, viewer, "Boundary Jacobian g2", wf->bdg2);
1338:   PetscWeakFormViewTable_Ascii(wf, viewer, "Boundary Jacobian g3", wf->bdg3);
1339:   PetscWeakFormViewTable_Ascii(wf, viewer, "Boundary Jacobian Preconditioner g0", wf->bdgp0);
1340:   PetscWeakFormViewTable_Ascii(wf, viewer, "Boundary Jacobian Preconditioner g1", wf->bdgp1);
1341:   PetscWeakFormViewTable_Ascii(wf, viewer, "Boundary Jacobian Preconditioner g2", wf->bdgp2);
1342:   PetscWeakFormViewTable_Ascii(wf, viewer, "Boundary Jacobian Preconditioner g3", wf->bdgp3);
1343:   PetscWeakFormViewTable_Ascii(wf, viewer, "Riemann Solver", wf->r);
1344:   PetscViewerASCIIPopTab(viewer);
1345:   return(0);
1346: }

1348: /*@C
1349:   PetscWeakFormView - Views a PetscWeakForm

1351:   Collective on wf

1353:   Input Parameter:
1354: + wf - the PetscWeakForm object to view
1355: - v  - the viewer

1357:   Level: developer

1359: .seealso PetscWeakFormDestroy(), PetscWeakFormCreate()
1360: @*/
1361: PetscErrorCode PetscWeakFormView(PetscWeakForm wf, PetscViewer v)
1362: {
1363:   PetscBool      iascii;

1368:   if (!v) {PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) wf), &v);}
1370:   PetscObjectTypeCompare((PetscObject) v, PETSCVIEWERASCII, &iascii);
1371:   if (iascii) {PetscWeakFormView_Ascii(wf, v);}
1372:   if (wf->ops->view) {(*wf->ops->view)(wf, v);}
1373:   return(0);
1374: }

1376: /*@
1377:   PetscWeakFormCreate - Creates an empty PetscWeakForm object.

1379:   Collective

1381:   Input Parameter:
1382: . comm - The communicator for the PetscWeakForm object

1384:   Output Parameter:
1385: . wf - The PetscWeakForm object

1387:   Level: beginner

1389: .seealso: PetscDS, PetscWeakFormDestroy()
1390: @*/
1391: PetscErrorCode PetscWeakFormCreate(MPI_Comm comm, PetscWeakForm *wf)
1392: {
1393:   PetscWeakForm  p;

1398:   *wf  = NULL;
1399:   PetscDSInitializePackage();

1401:   PetscHeaderCreate(p, PETSCWEAKFORM_CLASSID, "PetscWeakForm", "Weak Form System", "PetscWeakForm", comm, PetscWeakFormDestroy, PetscWeakFormView);

1403:   p->Nf = 0;
1404:   PetscChunkBufferCreate(sizeof(&PetscWeakFormCreate), 2, &p->funcs);
1405:   PetscHMapFormCreate(&p->obj);
1406:   PetscHMapFormCreate(&p->f0);
1407:   PetscHMapFormCreate(&p->f1);
1408:   PetscHMapFormCreate(&p->g0);
1409:   PetscHMapFormCreate(&p->g1);
1410:   PetscHMapFormCreate(&p->g2);
1411:   PetscHMapFormCreate(&p->g3);
1412:   PetscHMapFormCreate(&p->gp0);
1413:   PetscHMapFormCreate(&p->gp1);
1414:   PetscHMapFormCreate(&p->gp2);
1415:   PetscHMapFormCreate(&p->gp3);
1416:   PetscHMapFormCreate(&p->gt0);
1417:   PetscHMapFormCreate(&p->gt1);
1418:   PetscHMapFormCreate(&p->gt2);
1419:   PetscHMapFormCreate(&p->gt3);
1420:   PetscHMapFormCreate(&p->bdf0);
1421:   PetscHMapFormCreate(&p->bdf1);
1422:   PetscHMapFormCreate(&p->bdg0);
1423:   PetscHMapFormCreate(&p->bdg1);
1424:   PetscHMapFormCreate(&p->bdg2);
1425:   PetscHMapFormCreate(&p->bdg3);
1426:   PetscHMapFormCreate(&p->bdgp0);
1427:   PetscHMapFormCreate(&p->bdgp1);
1428:   PetscHMapFormCreate(&p->bdgp2);
1429:   PetscHMapFormCreate(&p->bdgp3);
1430:   PetscHMapFormCreate(&p->r);
1431:   *wf = p;
1432:   return(0);
1433: }