Actual source code: dtweakform.c

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

  3: PetscClassId PETSCWEAKFORM_CLASSID = 0;

  5: const char *const PetscWeakFormKinds[] = {"objective", "residual_f0", "residual_f1", "jacobian_g0", "jacobian_g1", "jacobian_g2", "jacobian_g3", "jacobian_preconditioner_g0", "jacobian_preconditioner_g1", "jacobian_preconditioner_g2", "jacobian_preconditioner_g3", "dynamic_jacobian_g0", "dynamic_jacobian_g1", "dynamic_jacobian_g2", "dynamic_jacobian_g3", "boundary_residual_f0", "boundary_residual_f1", "boundary_jacobian_g0", "boundary_jacobian_g1", "boundary_jacobian_g2", "boundary_jacobian_g3", "boundary_jacobian_preconditioner_g0", "boundary_jacobian_preconditioner_g1", "boundary_jacobian_preconditioner_g2", "boundary_jacobian_preconditioner_g3", "riemann_solver", "PetscWeakFormKind", "PETSC_WF_", NULL};

  7: static PetscErrorCode PetscChunkBufferCreate(size_t unitbytes, size_t expected, PetscChunkBuffer **buffer)
  8: {

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

 20: static PetscErrorCode PetscChunkBufferDuplicate(PetscChunkBuffer *buffer, PetscChunkBuffer **bufferNew)
 21: {

 25:   PetscNew(bufferNew);
 26:   PetscCalloc1(buffer->size*buffer->unitbytes, &(*bufferNew)->array);
 27:   PetscMemcpy((*bufferNew)->array, buffer->array, buffer->size*buffer->unitbytes);
 28:   (*bufferNew)->size      = buffer->size;
 29:   (*bufferNew)->unitbytes = buffer->unitbytes;
 30:   (*bufferNew)->alloc     = buffer->size*buffer->unitbytes;
 31:   return(0);
 32: }

 34: static PetscErrorCode PetscChunkBufferDestroy(PetscChunkBuffer **buffer)
 35: {

 39:   PetscFree((*buffer)->array);
 40:   PetscFree(*buffer);
 41:   return(0);
 42: }

 44: static PetscErrorCode PetscChunkBufferCreateChunk(PetscChunkBuffer *buffer, PetscInt size, PetscChunk *chunk)
 45: {

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

 52:     if (!buffer->alloc) buffer->alloc = (buffer->size + size)*buffer->unitbytes;
 53:     while ((buffer->size + size)*buffer->unitbytes > buffer->alloc) buffer->alloc *= 2;
 54:     PetscMalloc(buffer->alloc, &tmp);
 55:     PetscMemcpy(tmp, buffer->array, buffer->size*buffer->unitbytes);
 56:     PetscFree(buffer->array);
 57:     buffer->array = tmp;
 58:   }
 59:   chunk->start    = buffer->size*buffer->unitbytes;
 60:   chunk->size     = size;
 61:   chunk->reserved = size;
 62:   buffer->size   += size;
 63:   return(0);
 64: }

 66: static PetscErrorCode PetscChunkBufferEnlargeChunk(PetscChunkBuffer *buffer, PetscInt size, PetscChunk *chunk)
 67: {
 68:   size_t         siz = size;

 72:   if (chunk->size + size > chunk->reserved) {
 73:     PetscChunk newchunk;
 74:     PetscInt   reserved = chunk->size;

 76:     /* TODO Here if we had a chunk list, we could update them all to reclaim unused space */
 77:     while (reserved < chunk->size+size) reserved *= 2;
 78:     PetscChunkBufferCreateChunk(buffer, (size_t) reserved, &newchunk);
 79:     newchunk.size = chunk->size+size;
 80:     PetscMemcpy(&buffer->array[newchunk.start], &buffer->array[chunk->start], chunk->size * buffer->unitbytes);
 81:     *chunk = newchunk;
 82:   } else {
 83:     chunk->size += siz;
 84:   }
 85:   return(0);
 86: }

 88: /*@C
 89:   PetscFormKeySort - Sorts an array of PetscFormKey in place in increasing order.

 91:   Not Collective

 93:   Input Parameters:
 94: + n - number of values
 95: - X - array of PetscFormKey

 97:   Level: intermediate

 99: .seealso: PetscIntSortSemiOrdered(), PetscSortInt()
100: @*/
101: PetscErrorCode PetscFormKeySort(PetscInt n, PetscFormKey arr[])
102: {

106:   if (n <= 1) return(0);
108:   PetscTimSort(n, arr, sizeof(PetscFormKey), Compare_PetscFormKey_Private, NULL);
109:   return(0);
110: }

112: PetscErrorCode PetscWeakFormGetFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt part, PetscInt *n, void (***func)())
113: {
114:   PetscFormKey   key;
115:   PetscChunk     chunk;

119:   key.label = label; key.value = value; key.field = f; key.part = part;
120:   PetscHMapFormGet(ht, key, &chunk);
121:   if (chunk.size < 0) {*n = 0;          *func = NULL;}
122:   else                {*n = chunk.size; *func = (void (**)()) &wf->funcs->array[chunk.start];}
123:   return(0);
124: }

126: /* A NULL argument for func causes this to clear the key */
127: PetscErrorCode PetscWeakFormSetFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt part, PetscInt n, void (**func)())
128: {
129:   PetscFormKey   key;
130:   PetscChunk     chunk;
131:   PetscInt       i;

135:   key.label = label; key.value = value; key.field = f; key.part = part;
136:   if (!func) {
137:     PetscHMapFormDel(ht, key);
138:     return(0);
139:   } else {
140:     PetscHMapFormGet(ht, key, &chunk);
141:   }
142:   if (chunk.size < 0) {
143:     PetscChunkBufferCreateChunk(wf->funcs, n, &chunk);
144:     PetscHMapFormSet(ht, key, chunk);
145:   } else if (chunk.size <= n) {
146:     PetscChunkBufferEnlargeChunk(wf->funcs, n - chunk.size, &chunk);
147:     PetscHMapFormSet(ht, key, chunk);
148:   }
149:   for (i = 0; i < n; ++i) ((void (**)()) &wf->funcs->array[chunk.start])[i] = func[i];
150:   return(0);
151: }

153: PetscErrorCode PetscWeakFormAddFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt part, void (*func)())
154: {
155:   PetscFormKey   key;
156:   PetscChunk     chunk;

160:   if (!func) return(0);
161:   key.label = label; key.value = value; key.field = f; key.part = part;
162:   PetscHMapFormGet(ht, key, &chunk);
163:   if (chunk.size < 0) {
164:     PetscChunkBufferCreateChunk(wf->funcs, 1, &chunk);
165:     PetscHMapFormSet(ht, key, chunk);
166:     ((void (**)()) &wf->funcs->array[chunk.start])[0] = func;
167:   } else {
168:     PetscChunkBufferEnlargeChunk(wf->funcs, 1, &chunk);
169:     PetscHMapFormSet(ht, key, chunk);
170:     ((void (**)()) &wf->funcs->array[chunk.start])[chunk.size-1] = func;
171:   }
172:   return(0);
173: }

175: PetscErrorCode PetscWeakFormGetIndexFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt part, PetscInt ind, void (**func)())
176: {
177:   PetscFormKey key;
178:   PetscChunk       chunk;
179:   PetscErrorCode   ierr;

182:   key.label = label; key.value = value; key.field = f; key.part = part;
183:   PetscHMapFormGet(ht, key, &chunk);
184:   if (chunk.size < 0) {*func = NULL;}
185:   else {
186:     if (ind >= chunk.size) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Index %D not in [0, %D)", ind, chunk.size);
187:     *func = ((void (**)()) &wf->funcs->array[chunk.start])[ind];
188:   }
189:   return(0);
190: }

192: /* Ignore a NULL func */
193: PetscErrorCode PetscWeakFormSetIndexFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt part, PetscInt ind, void (*func)())
194: {
195:   PetscFormKey key;
196:   PetscChunk       chunk;
197:   PetscErrorCode   ierr;

200:   if (!func) return(0);
201:   key.label = label; key.value = value; key.field = f; key.part = part;
202:   PetscHMapFormGet(ht, key, &chunk);
203:   if (chunk.size < 0) {
204:     PetscChunkBufferCreateChunk(wf->funcs, ind+1, &chunk);
205:     PetscHMapFormSet(ht, key, chunk);
206:   } else if (chunk.size <= ind) {
207:     PetscChunkBufferEnlargeChunk(wf->funcs, ind - chunk.size + 1, &chunk);
208:     PetscHMapFormSet(ht, key, chunk);
209:   }
210:   ((void (**)()) &wf->funcs->array[chunk.start])[ind] = func;
211:   return(0);
212: }

214: PetscErrorCode PetscWeakFormClearIndexFunction_Private(PetscWeakForm wf, PetscHMapForm ht, DMLabel label, PetscInt value, PetscInt f, PetscInt part, PetscInt ind)
215: {
216:   PetscFormKey key;
217:   PetscChunk       chunk;
218:   PetscErrorCode   ierr;

221:   key.label = label; key.value = value; key.field = f; key.part = part;
222:   PetscHMapFormGet(ht, key, &chunk);
223:   if (chunk.size < 0) {
224:     return(0);
225:   } else if (!ind && chunk.size == 1) {
226:     PetscHMapFormDel(ht, key);
227:     return(0);
228:   } else if (chunk.size <= ind) {
229:     return(0);
230:   }
231:   ((void (**)()) &wf->funcs->array[chunk.start])[ind] = NULL;
232:   return(0);
233: }

235: /*@
236:   PetscWeakFormCopy - Copy the pointwise functions to another PetscWeakForm

238:   Not Collective

240:   Input Parameter:
241: . wf - The original PetscWeakForm

243:   Output Parameter:
244: . wfNew - The copy PetscWeakForm

246:   Level: intermediate

248: .seealso: PetscWeakFormCreate(), PetscWeakFormDestroy()
249: @*/
250: PetscErrorCode PetscWeakFormCopy(PetscWeakForm wf, PetscWeakForm wfNew)
251: {
252:   PetscInt       f;

256:   wfNew->Nf = wf->Nf;
257:   PetscChunkBufferDestroy(&wfNew->funcs);
258:   PetscChunkBufferDuplicate(wf->funcs, &wfNew->funcs);
259:   for (f = 0; f < PETSC_NUM_WF; ++f) {
260:     PetscHMapFormDestroy(&wfNew->form[f]);
261:     PetscHMapFormDuplicate(wf->form[f], &wfNew->form[f]);
262:   }
263:   return(0);
264: }

266: /*@
267:   PetscWeakFormClear - Clear all functions from the PetscWeakForm

269:   Not Collective

271:   Input Parameter:
272: . wf - The original PetscWeakForm

274:   Level: intermediate

276: .seealso: PetscWeakFormCopy(), PetscWeakFormCreate(), PetscWeakFormDestroy()
277: @*/
278: PetscErrorCode PetscWeakFormClear(PetscWeakForm wf)
279: {
280:   PetscInt       f;

284:   for (f = 0; f < PETSC_NUM_WF; ++f) {PetscHMapFormClear(wf->form[f]);}
285:   return(0);
286: }

288: static PetscErrorCode PetscWeakFormRewriteKeys_Internal(PetscWeakForm wf, PetscHMapForm hmap, DMLabel label, PetscInt Nv, const PetscInt values[])
289: {
290:   PetscFormKey  *keys;
291:   PetscInt       n, i, v, off = 0;

295:   PetscHMapFormGetSize(hmap, &n);
296:   PetscMalloc1(n, &keys);
297:   PetscHMapFormGetKeys(hmap, &off, keys);
298:   for (i = 0; i < n; ++i) {
299:     if (keys[i].label == label) {
300:       PetscBool clear = PETSC_TRUE;
301:       void   (**funcs)();
302:       PetscInt  Nf;

304:       PetscWeakFormGetFunction_Private(wf, hmap, keys[i].label, keys[i].value, keys[i].field, keys[i].part, &Nf, &funcs);
305:       for (v = 0; v < Nv; ++v) {
306:         PetscWeakFormSetFunction_Private(wf, hmap, keys[i].label, values[v], keys[i].field, keys[i].part, Nf, funcs);
307:         if (values[v] == keys[i].value) clear = PETSC_FALSE;
308:       }
309:       if (clear) {PetscWeakFormSetFunction_Private(wf, hmap, keys[i].label, keys[i].value, keys[i].field, keys[i].part, 0, NULL);}
310:     }
311:   }
312:   PetscFree(keys);
313:   return(0);
314: }

316: /*@C
317:   PetscWeakFormRewriteKeys - Change any key on the given label to use the new set of label values

319:   Not Collective

321:   Input Parameters:
322: + wf     - The original PetscWeakForm
323: . label  - The label to change keys for
324: . Nv     - The number of new label values
325: - values - The set of new values to relabel keys with

327:   Note: This is used internally when boundary label values are specified from the command line.

329:   Level: intermediate

331: .seealso: PetscWeakFormReplaceLabel(), PetscWeakFormCreate(), PetscWeakFormDestroy()
332: @*/
333: PetscErrorCode PetscWeakFormRewriteKeys(PetscWeakForm wf, DMLabel label, PetscInt Nv, const PetscInt values[])
334: {
335:   PetscInt       f;

339:   for (f = 0; f < PETSC_NUM_WF; ++f) {PetscWeakFormRewriteKeys_Internal(wf, wf->form[f], label, Nv, values);}
340:   return(0);
341: }

343: static PetscErrorCode PetscWeakFormReplaceLabel_Internal(PetscWeakForm wf, PetscHMapForm hmap, DMLabel label)
344: {
345:   PetscFormKey  *keys;
346:   PetscInt       n, i, off = 0, maxFuncs = 0;
347:   void       (**tmpf)();
348:   const char    *name = NULL;

352:   if (label) {PetscObjectGetName((PetscObject) label, &name);}
353:   PetscHMapFormGetSize(hmap, &n);
354:   PetscMalloc1(n, &keys);
355:   PetscHMapFormGetKeys(hmap, &off, keys);
356:   for (i = 0; i < n; ++i) {
357:     PetscBool   match = PETSC_FALSE;
358:     const char *lname = NULL;

360:     if (label == keys[i].label) continue;
361:     if (keys[i].label) {PetscObjectGetName((PetscObject) keys[i].label, &lname);}
362:     PetscStrcmp(name, lname, &match);
363:     if ((!name && !lname) || match) {
364:       void  (**funcs)();
365:       PetscInt Nf;

367:       PetscWeakFormGetFunction_Private(wf, hmap, keys[i].label, keys[i].value, keys[i].field, keys[i].part, &Nf, &funcs);
368:       maxFuncs = PetscMax(maxFuncs, Nf);
369:     }
370:   }
371:   /* Need temp space because chunk buffer can be reallocated in SetFunction() call */
372:   PetscMalloc1(maxFuncs, &tmpf);
373:   for (i = 0; i < n; ++i) {
374:     PetscBool   match = PETSC_FALSE;
375:     const char *lname = NULL;

377:     if (label == keys[i].label) continue;
378:     if (keys[i].label) {PetscObjectGetName((PetscObject) keys[i].label, &lname);}
379:     PetscStrcmp(name, lname, &match);
380:     if ((!name && !lname) || match) {
381:       void  (**funcs)();
382:       PetscInt Nf, j;

384:       PetscWeakFormGetFunction_Private(wf, hmap, keys[i].label, keys[i].value, keys[i].field, keys[i].part, &Nf, &funcs);
385:       for (j = 0; j < Nf; ++j) tmpf[j] = funcs[j];
386:       PetscWeakFormSetFunction_Private(wf, hmap, label,         keys[i].value, keys[i].field, keys[i].part,  Nf,  tmpf);
387:       PetscWeakFormSetFunction_Private(wf, hmap, keys[i].label, keys[i].value, keys[i].field, keys[i].part,  0,   NULL);
388:     }
389:   }
390:   PetscFree(tmpf);
391:   PetscFree(keys);
392:   return(0);
393: }

395: /*@C
396:   PetscWeakFormReplaceLabel - Change any key on a label of the same name to use the new label

398:   Not Collective

400:   Input Parameters:
401: + wf    - The original PetscWeakForm
402: - label - The label to change keys for

404:   Note: This is used internally when meshes are modified

406:   Level: intermediate

408: .seealso: PetscWeakFormRewriteKeys(), PetscWeakFormCreate(), PetscWeakFormDestroy()
409: @*/
410: PetscErrorCode PetscWeakFormReplaceLabel(PetscWeakForm wf, DMLabel label)
411: {
412:   PetscInt       f;

416:   for (f = 0; f < PETSC_NUM_WF; ++f) {PetscWeakFormReplaceLabel_Internal(wf, wf->form[f], label);}
417:   return(0);
418: }

420: PetscErrorCode PetscWeakFormClearIndex(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscWeakFormKind kind, PetscInt ind)
421: {

425:   PetscWeakFormClearIndexFunction_Private(wf, wf->form[kind], label, val, f, part, ind);
426:   return(0);
427: }

429: PetscErrorCode PetscWeakFormGetObjective(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt *n,
430:                                          void (***obj)(PetscInt, PetscInt, PetscInt,
431:                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
432:                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
433:                                                        PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
434: {

438:   PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_OBJECTIVE], label, val, f, part, n, (void (***)(void)) obj);
439:   return(0);
440: }

442: PetscErrorCode PetscWeakFormSetObjective(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt n,
443:                                          void (**obj)(PetscInt, PetscInt, PetscInt,
444:                                                       const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const  PetscScalar[],
445:                                                       const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
446:                                                       PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
447: {

451:   PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_OBJECTIVE], label, val, f, part, n, (void (**)(void)) obj);
452:   return(0);
453: }

455: PetscErrorCode PetscWeakFormAddObjective(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part,
456:                                          void (*obj)(PetscInt, PetscInt, PetscInt,
457:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
458:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
459:                                                      PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
460: {

464:   PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_OBJECTIVE], label, val, f, part, (void (*)(void)) obj);
465:   return(0);
466: }

468: PetscErrorCode PetscWeakFormGetIndexObjective(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt ind,
469:                                               void (**obj)(PetscInt, PetscInt, PetscInt,
470:                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
471:                                                            const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
472:                                                            PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
473: {

477:   PetscWeakFormGetIndexFunction_Private(wf, wf->form[PETSC_WF_OBJECTIVE], label, val, f, part, ind, (void (**)(void)) obj);
478:   return(0);
479: }

481: PetscErrorCode PetscWeakFormSetIndexObjective(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt ind,
482:                                               void (*obj)(PetscInt, PetscInt, PetscInt,
483:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
484:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
485:                                                           PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
486: {

490:   PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_OBJECTIVE], label, val, f, part, ind, (void (*)(void)) obj);
491:   return(0);
492: }

494: PetscErrorCode PetscWeakFormGetResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part,
495:                                         PetscInt *n0,
496:                                         void (***f0)(PetscInt, PetscInt, PetscInt,
497:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
498:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
499:                                                      PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
500:                                         PetscInt *n1,
501:                                         void (***f1)(PetscInt, PetscInt, PetscInt,
502:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
503:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
504:                                                      PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
505: {

509:   PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_F0], label, val, f, part, n0, (void (***)(void)) f0);
510:   PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_F1], label, val, f, part, n1, (void (***)(void)) f1);
511:   return(0);
512: }

514: PetscErrorCode PetscWeakFormAddResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part,
515:                                         void (*f0)(PetscInt, PetscInt, PetscInt,
516:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
517:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
518:                                                    PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
519:                                         void (*f1)(PetscInt, PetscInt, PetscInt,
520:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
521:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
522:                                                    PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
523: {

527:   PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_F0], label, val, f, part, (void (*)(void)) f0);
528:   PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_F1], label, val, f, part, (void (*)(void)) f1);
529:   return(0);
530: }

532: PetscErrorCode PetscWeakFormSetResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part,
533:                                         PetscInt n0,
534:                                         void (**f0)(PetscInt, PetscInt, PetscInt,
535:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
536:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
537:                                                    PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
538:                                         PetscInt n1,
539:                                         void (**f1)(PetscInt, PetscInt, PetscInt,
540:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
541:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
542:                                                    PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
543: {

547:   PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_F0], label, val, f, part, n0, (void (**)(void)) f0);
548:   PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_F1], label, val, f, part, n1, (void (**)(void)) f1);
549:   return(0);
550: }

552: PetscErrorCode PetscWeakFormSetIndexResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part,
553:                                         PetscInt i0,
554:                                         void (*f0)(PetscInt, PetscInt, PetscInt,
555:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
556:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
557:                                                    PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
558:                                         PetscInt i1,
559:                                         void (*f1)(PetscInt, PetscInt, PetscInt,
560:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
561:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
562:                                                    PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
563: {

567:   PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_F0], label, val, f, part, i0, (void (*)(void)) f0);
568:   PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_F1], label, val, f, part, i1, (void (*)(void)) f1);
569:   return(0);
570: }

572: PetscErrorCode PetscWeakFormGetBdResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part,
573:                                           PetscInt *n0,
574:                                         void (***f0)(PetscInt, PetscInt, PetscInt,
575:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
576:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
577:                                                      PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
578:                                         PetscInt *n1,
579:                                         void (***f1)(PetscInt, PetscInt, PetscInt,
580:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
581:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
582:                                                      PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
583: {

587:   PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDF0], label, val, f, part, n0, (void (***)(void)) f0);
588:   PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDF1], label, val, f, part, n1, (void (***)(void)) f1);
589:   return(0);
590: }

592: PetscErrorCode PetscWeakFormAddBdResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part,
593:                                           void (*f0)(PetscInt, PetscInt, PetscInt,
594:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
595:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
596:                                                      PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
597:                                           void (*f1)(PetscInt, PetscInt, PetscInt,
598:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
599:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
600:                                                      PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
601: {

605:   PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDF0], label, val, f, part, (void (*)(void)) f0);
606:   PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDF1], label, val, f, part, (void (*)(void)) f1);
607:   return(0);
608: }

610: PetscErrorCode PetscWeakFormSetBdResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part,
611:                                           PetscInt n0,
612:                                           void (**f0)(PetscInt, PetscInt, PetscInt,
613:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
614:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
615:                                                      PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
616:                                           PetscInt n1,
617:                                           void (**f1)(PetscInt, PetscInt, PetscInt,
618:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
619:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
620:                                                      PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
621: {

625:   PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDF0], label, val, f, part, n0, (void (**)(void)) f0);
626:   PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDF1], label, val, f, part, n1, (void (**)(void)) f1);
627:   return(0);
628: }

630: PetscErrorCode PetscWeakFormSetIndexBdResidual(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part,
631:                                           PetscInt i0,
632:                                           void (*f0)(PetscInt, PetscInt, PetscInt,
633:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
634:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
635:                                                      PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
636:                                           PetscInt i1,
637:                                           void (*f1)(PetscInt, PetscInt, PetscInt,
638:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
639:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
640:                                                      PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
641: {

645:   PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDF0], label, val, f, part, i0, (void (*)(void)) f0);
646:   PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDF1], label, val, f, part, i1, (void (*)(void)) f1);
647:   return(0);
648: }

650: PetscErrorCode PetscWeakFormHasJacobian(PetscWeakForm wf, PetscBool *hasJac)
651: {
652:   PetscInt       n0, n1, n2, n3;

658:   PetscHMapFormGetSize(wf->form[PETSC_WF_G0], &n0);
659:   PetscHMapFormGetSize(wf->form[PETSC_WF_G1], &n1);
660:   PetscHMapFormGetSize(wf->form[PETSC_WF_G2], &n2);
661:   PetscHMapFormGetSize(wf->form[PETSC_WF_G3], &n3);
662:   *hasJac = n0+n1+n2+n3 ? PETSC_TRUE : PETSC_FALSE;
663:   return(0);
664: }

666: PetscErrorCode PetscWeakFormGetJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part,
667:                                         PetscInt *n0,
668:                                         void (***g0)(PetscInt, PetscInt, PetscInt,
669:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
670:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
671:                                                      PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
672:                                         PetscInt *n1,
673:                                         void (***g1)(PetscInt, PetscInt, PetscInt,
674:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
675:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
676:                                                      PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
677:                                         PetscInt *n2,
678:                                         void (***g2)(PetscInt, PetscInt, PetscInt,
679:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
680:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
681:                                                      PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
682:                                         PetscInt *n3,
683:                                         void (***g3)(PetscInt, PetscInt, PetscInt,
684:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
685:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
686:                                                      PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
687: {
688:   PetscInt       find = f*wf->Nf + g;

692:   PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_G0], label, val, find, part, n0, (void (***)(void)) g0);
693:   PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_G1], label, val, find, part, n1, (void (***)(void)) g1);
694:   PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_G2], label, val, find, part, n2, (void (***)(void)) g2);
695:   PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_G3], label, val, find, part, n3, (void (***)(void)) g3);
696:   return(0);
697: }

699: PetscErrorCode PetscWeakFormAddJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part,
700:                                         void (*g0)(PetscInt, PetscInt, PetscInt,
701:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
702:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
703:                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
704:                                         void (*g1)(PetscInt, PetscInt, PetscInt,
705:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
706:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
707:                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
708:                                         void (*g2)(PetscInt, PetscInt, PetscInt,
709:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
710:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
711:                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
712:                                         void (*g3)(PetscInt, PetscInt, PetscInt,
713:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
714:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
715:                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
716: {
717:   PetscInt       find = f*wf->Nf + g;

721:   PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_G0], label, val, find, part, (void (*)(void)) g0);
722:   PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_G1], label, val, find, part, (void (*)(void)) g1);
723:   PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_G2], label, val, find, part, (void (*)(void)) g2);
724:   PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_G3], label, val, find, part, (void (*)(void)) g3);
725:   return(0);
726: }

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

754:   PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_G0], label, val, find, part, n0, (void (**)(void)) g0);
755:   PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_G1], label, val, find, part, n1, (void (**)(void)) g1);
756:   PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_G2], label, val, find, part, n2, (void (**)(void)) g2);
757:   PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_G3], label, val, find, part, n3, (void (**)(void)) g3);
758:   return(0);
759: }

761: PetscErrorCode PetscWeakFormSetIndexJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part,
762:                                         PetscInt i0,
763:                                         void (*g0)(PetscInt, PetscInt, PetscInt,
764:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
765:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
766:                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
767:                                         PetscInt i1,
768:                                         void (*g1)(PetscInt, PetscInt, PetscInt,
769:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
770:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
771:                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
772:                                         PetscInt i2,
773:                                         void (*g2)(PetscInt, PetscInt, PetscInt,
774:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
775:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
776:                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
777:                                         PetscInt i3,
778:                                         void (*g3)(PetscInt, PetscInt, PetscInt,
779:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
780:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
781:                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
782: {
783:   PetscInt       find = f*wf->Nf + g;

787:   PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_G0], label, val, find, part, i0, (void (*)(void)) g0);
788:   PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_G1], label, val, find, part, i1, (void (*)(void)) g1);
789:   PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_G2], label, val, find, part, i2, (void (*)(void)) g2);
790:   PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_G3], label, val, find, part, i3, (void (*)(void)) g3);
791:   return(0);
792: }

794: PetscErrorCode PetscWeakFormHasJacobianPreconditioner(PetscWeakForm wf, PetscBool *hasJacPre)
795: {
796:   PetscInt       n0, n1, n2, n3;

802:   PetscHMapFormGetSize(wf->form[PETSC_WF_GP0], &n0);
803:   PetscHMapFormGetSize(wf->form[PETSC_WF_GP1], &n1);
804:   PetscHMapFormGetSize(wf->form[PETSC_WF_GP2], &n2);
805:   PetscHMapFormGetSize(wf->form[PETSC_WF_GP3], &n3);
806:   *hasJacPre = n0+n1+n2+n3 ? PETSC_TRUE : PETSC_FALSE;
807:   return(0);
808: }

810: PetscErrorCode PetscWeakFormGetJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part,
811:                                                       PetscInt *n0,
812:                                                       void (***g0)(PetscInt, PetscInt, PetscInt,
813:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
814:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
815:                                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
816:                                                       PetscInt *n1,
817:                                                       void (***g1)(PetscInt, PetscInt, PetscInt,
818:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
819:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
820:                                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
821:                                                       PetscInt *n2,
822:                                                       void (***g2)(PetscInt, PetscInt, PetscInt,
823:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
824:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
825:                                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
826:                                                       PetscInt *n3,
827:                                                       void (***g3)(PetscInt, PetscInt, PetscInt,
828:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
829:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
830:                                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
831: {
832:   PetscInt       find = f*wf->Nf + g;

836:   PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GP0], label, val, find, part, n0, (void (***)(void)) g0);
837:   PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GP1], label, val, find, part, n1, (void (***)(void)) g1);
838:   PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GP2], label, val, find, part, n2, (void (***)(void)) g2);
839:   PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GP3], label, val, find, part, n3, (void (***)(void)) g3);
840:   return(0);
841: }

843: PetscErrorCode PetscWeakFormAddJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part,
844:                                         void (*g0)(PetscInt, PetscInt, PetscInt,
845:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
846:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
847:                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
848:                                         void (*g1)(PetscInt, PetscInt, PetscInt,
849:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
850:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
851:                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
852:                                         void (*g2)(PetscInt, PetscInt, PetscInt,
853:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
854:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
855:                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
856:                                         void (*g3)(PetscInt, PetscInt, PetscInt,
857:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
858:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
859:                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
860: {
861:   PetscInt       find = f*wf->Nf + g;

865:   PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GP0], label, val, find, part, (void (*)(void)) g0);
866:   PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GP1], label, val, find, part, (void (*)(void)) g1);
867:   PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GP2], label, val, find, part, (void (*)(void)) g2);
868:   PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GP3], label, val, find, part, (void (*)(void)) g3);
869:   return(0);
870: }

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

898:   PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GP0], label, val, find, part, n0, (void (**)(void)) g0);
899:   PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GP1], label, val, find, part, n1, (void (**)(void)) g1);
900:   PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GP2], label, val, find, part, n2, (void (**)(void)) g2);
901:   PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GP3], label, val, find, part, n3, (void (**)(void)) g3);
902:   return(0);
903: }

905: PetscErrorCode PetscWeakFormSetIndexJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part,
906:                                                       PetscInt i0,
907:                                                       void (*g0)(PetscInt, PetscInt, PetscInt,
908:                                                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
909:                                                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
910:                                                                  PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
911:                                                       PetscInt i1,
912:                                                       void (*g1)(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[], PetscInt, const PetscScalar[], PetscScalar[]),
916:                                                       PetscInt i2,
917:                                                       void (*g2)(PetscInt, PetscInt, PetscInt,
918:                                                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
919:                                                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
920:                                                                  PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
921:                                                       PetscInt i3,
922:                                                       void (*g3)(PetscInt, PetscInt, PetscInt,
923:                                                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
924:                                                                  const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
925:                                                                  PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
926: {
927:   PetscInt       find = f*wf->Nf + g;

931:   PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GP0], label, val, find, part, i0, (void (*)(void)) g0);
932:   PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GP1], label, val, find, part, i1, (void (*)(void)) g1);
933:   PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GP2], label, val, find, part, i2, (void (*)(void)) g2);
934:   PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GP3], label, val, find, part, i3, (void (*)(void)) g3);
935:   return(0);
936: }

938: PetscErrorCode PetscWeakFormHasBdJacobian(PetscWeakForm wf, PetscBool *hasJac)
939: {
940:   PetscInt       n0, n1, n2, n3;

946:   PetscHMapFormGetSize(wf->form[PETSC_WF_BDG0], &n0);
947:   PetscHMapFormGetSize(wf->form[PETSC_WF_BDG1], &n1);
948:   PetscHMapFormGetSize(wf->form[PETSC_WF_BDG2], &n2);
949:   PetscHMapFormGetSize(wf->form[PETSC_WF_BDG3], &n3);
950:   *hasJac = n0+n1+n2+n3 ? PETSC_TRUE : PETSC_FALSE;
951:   return(0);
952: }

954: PetscErrorCode PetscWeakFormGetBdJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part,
955:                                           PetscInt *n0,
956:                                           void (***g0)(PetscInt, PetscInt, PetscInt,
957:                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
958:                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
959:                                                        PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
960:                                           PetscInt *n1,
961:                                           void (***g1)(PetscInt, PetscInt, PetscInt,
962:                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
963:                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
964:                                                        PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
965:                                           PetscInt *n2,
966:                                           void (***g2)(PetscInt, PetscInt, PetscInt,
967:                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
968:                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
969:                                                        PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
970:                                           PetscInt *n3,
971:                                           void (***g3)(PetscInt, PetscInt, PetscInt,
972:                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
973:                                                        const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
974:                                                        PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
975: {
976:   PetscInt       find = f*wf->Nf + g;

980:   PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDG0], label, val, find, part, n0, (void (***)(void)) g0);
981:   PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDG1], label, val, find, part, n1, (void (***)(void)) g1);
982:   PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDG2], label, val, find, part, n2, (void (***)(void)) g2);
983:   PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDG3], label, val, find, part, n3, (void (***)(void)) g3);
984:   return(0);
985: }

987: PetscErrorCode PetscWeakFormAddBdJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part,
988:                                           void (*g0)(PetscInt, PetscInt, PetscInt,
989:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
990:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
991:                                                      PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
992:                                           void (*g1)(PetscInt, PetscInt, PetscInt,
993:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
994:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
995:                                                      PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
996:                                           void (*g2)(PetscInt, PetscInt, PetscInt,
997:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
998:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
999:                                                      PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1000:                                           void (*g3)(PetscInt, PetscInt, PetscInt,
1001:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1002:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1003:                                                      PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
1004: {
1005:   PetscInt       find = f*wf->Nf + g;

1009:   PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDG0], label, val, find, part, (void (*)(void)) g0);
1010:   PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDG1], label, val, find, part, (void (*)(void)) g1);
1011:   PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDG2], label, val, find, part, (void (*)(void)) g2);
1012:   PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDG3], label, val, find, part, (void (*)(void)) g3);
1013:   return(0);
1014: }

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

1042:   PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDG0], label, val, find, part, n0, (void (**)(void)) g0);
1043:   PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDG1], label, val, find, part, n1, (void (**)(void)) g1);
1044:   PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDG2], label, val, find, part, n2, (void (**)(void)) g2);
1045:   PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDG3], label, val, find, part, n3, (void (**)(void)) g3);
1046:   return(0);
1047: }

1049: PetscErrorCode PetscWeakFormSetIndexBdJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part,
1050:                                           PetscInt i0,
1051:                                           void (*g0)(PetscInt, PetscInt, PetscInt,
1052:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1053:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1054:                                                      PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1055:                                           PetscInt i1,
1056:                                           void (*g1)(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[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1060:                                           PetscInt i2,
1061:                                           void (*g2)(PetscInt, PetscInt, PetscInt,
1062:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1063:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1064:                                                      PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1065:                                           PetscInt i3,
1066:                                           void (*g3)(PetscInt, PetscInt, PetscInt,
1067:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1068:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1069:                                                      PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
1070: {
1071:   PetscInt       find = f*wf->Nf + g;

1075:   PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDG0], label, val, find, part, i0, (void (*)(void)) g0);
1076:   PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDG1], label, val, find, part, i1, (void (*)(void)) g1);
1077:   PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDG2], label, val, find, part, i2, (void (*)(void)) g2);
1078:   PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDG3], label, val, find, part, i3, (void (*)(void)) g3);
1079:   return(0);
1080: }

1082: PetscErrorCode PetscWeakFormHasBdJacobianPreconditioner(PetscWeakForm wf, PetscBool *hasJacPre)
1083: {
1084:   PetscInt       n0, n1, n2, n3;

1090:   PetscHMapFormGetSize(wf->form[PETSC_WF_BDGP0], &n0);
1091:   PetscHMapFormGetSize(wf->form[PETSC_WF_BDGP1], &n1);
1092:   PetscHMapFormGetSize(wf->form[PETSC_WF_BDGP2], &n2);
1093:   PetscHMapFormGetSize(wf->form[PETSC_WF_BDGP3], &n3);
1094:   *hasJacPre = n0+n1+n2+n3 ? PETSC_TRUE : PETSC_FALSE;
1095:   return(0);
1096: }

1098: PetscErrorCode PetscWeakFormGetBdJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part,
1099:                                                         PetscInt *n0,
1100:                                                         void (***g0)(PetscInt, PetscInt, PetscInt,
1101:                                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1102:                                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1103:                                                                      PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1104:                                                         PetscInt *n1,
1105:                                                         void (***g1)(PetscInt, PetscInt, PetscInt,
1106:                                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1107:                                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1108:                                                                      PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1109:                                                         PetscInt *n2,
1110:                                                         void (***g2)(PetscInt, PetscInt, PetscInt,
1111:                                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1112:                                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1113:                                                                      PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1114:                                                         PetscInt *n3,
1115:                                                         void (***g3)(PetscInt, PetscInt, PetscInt,
1116:                                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1117:                                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1118:                                                                      PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
1119: {
1120:   PetscInt       find = f*wf->Nf + g;

1124:   PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDGP0], label, val, find, part, n0, (void (***)(void)) g0);
1125:   PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDGP1], label, val, find, part, n1, (void (***)(void)) g1);
1126:   PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDGP2], label, val, find, part, n2, (void (***)(void)) g2);
1127:   PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_BDGP3], label, val, find, part, n3, (void (***)(void)) g3);
1128:   return(0);
1129: }

1131: PetscErrorCode PetscWeakFormAddBdJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part,
1132:                                                         void (*g0)(PetscInt, PetscInt, PetscInt,
1133:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1134:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1135:                                                                    PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1136:                                                         void (*g1)(PetscInt, PetscInt, PetscInt,
1137:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1138:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1139:                                                                    PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1140:                                                         void (*g2)(PetscInt, PetscInt, PetscInt,
1141:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1142:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1143:                                                                    PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1144:                                                         void (*g3)(PetscInt, PetscInt, PetscInt,
1145:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1146:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1147:                                                                    PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
1148: {
1149:   PetscInt       find = f*wf->Nf + g;

1153:   PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDGP0], label, val, find, part, (void (*)(void)) g0);
1154:   PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDGP1], label, val, find, part, (void (*)(void)) g1);
1155:   PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDGP2], label, val, find, part, (void (*)(void)) g2);
1156:   PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_BDGP3], label, val, find, part, (void (*)(void)) g3);
1157:   return(0);
1158: }

1160: PetscErrorCode PetscWeakFormSetBdJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part,
1161:                                                         PetscInt n0,
1162:                                                         void (**g0)(PetscInt, PetscInt, PetscInt,
1163:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1164:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1165:                                                                    PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1166:                                                         PetscInt n1,
1167:                                                         void (**g1)(PetscInt, PetscInt, PetscInt,
1168:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1169:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1170:                                                                    PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1171:                                                         PetscInt n2,
1172:                                                         void (**g2)(PetscInt, PetscInt, PetscInt,
1173:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1174:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1175:                                                                    PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1176:                                                         PetscInt n3,
1177:                                                         void (**g3)(PetscInt, PetscInt, PetscInt,
1178:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1179:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1180:                                                                    PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
1181: {
1182:   PetscInt       find = f*wf->Nf + g;

1186:   PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDGP0], label, val, find, part, n0, (void (**)(void)) g0);
1187:   PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDGP1], label, val, find, part, n1, (void (**)(void)) g1);
1188:   PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDGP2], label, val, find, part, n2, (void (**)(void)) g2);
1189:   PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_BDGP3], label, val, find, part, n3, (void (**)(void)) g3);
1190:   return(0);
1191: }

1193: PetscErrorCode PetscWeakFormSetIndexBdJacobianPreconditioner(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part,
1194:                                                         PetscInt i0,
1195:                                                         void (*g0)(PetscInt, PetscInt, PetscInt,
1196:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1197:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1198:                                                                    PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1199:                                                         PetscInt i1,
1200:                                                         void (*g1)(PetscInt, PetscInt, PetscInt,
1201:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1202:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1203:                                                                    PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1204:                                                         PetscInt i2,
1205:                                                         void (*g2)(PetscInt, PetscInt, PetscInt,
1206:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1207:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1208:                                                                    PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1209:                                                         PetscInt i3,
1210:                                                         void (*g3)(PetscInt, PetscInt, PetscInt,
1211:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1212:                                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1213:                                                                    PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
1214: {
1215:   PetscInt       find = f*wf->Nf + g;

1219:   PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDGP0], label, val, find, part, i0, (void (*)(void)) g0);
1220:   PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDGP1], label, val, find, part, i1, (void (*)(void)) g1);
1221:   PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDGP2], label, val, find, part, i2, (void (*)(void)) g2);
1222:   PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_BDGP3], label, val, find, part, i3, (void (*)(void)) g3);
1223:   return(0);
1224: }

1226: PetscErrorCode PetscWeakFormHasDynamicJacobian(PetscWeakForm wf, PetscBool *hasDynJac)
1227: {
1228:   PetscInt       n0, n1, n2, n3;

1234:   PetscHMapFormGetSize(wf->form[PETSC_WF_GT0], &n0);
1235:   PetscHMapFormGetSize(wf->form[PETSC_WF_GT1], &n1);
1236:   PetscHMapFormGetSize(wf->form[PETSC_WF_GT2], &n2);
1237:   PetscHMapFormGetSize(wf->form[PETSC_WF_GT3], &n3);
1238:   *hasDynJac = n0+n1+n2+n3 ? PETSC_TRUE : PETSC_FALSE;
1239:   return(0);
1240: }

1242: PetscErrorCode PetscWeakFormGetDynamicJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part,
1243:                                         PetscInt *n0,
1244:                                         void (***g0)(PetscInt, PetscInt, PetscInt,
1245:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1246:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1247:                                                      PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1248:                                         PetscInt *n1,
1249:                                         void (***g1)(PetscInt, PetscInt, PetscInt,
1250:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1251:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1252:                                                      PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1253:                                         PetscInt *n2,
1254:                                         void (***g2)(PetscInt, PetscInt, PetscInt,
1255:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1256:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1257:                                                      PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1258:                                         PetscInt *n3,
1259:                                         void (***g3)(PetscInt, PetscInt, PetscInt,
1260:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1261:                                                      const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1262:                                                      PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
1263: {
1264:   PetscInt       find = f*wf->Nf + g;

1268:   PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GT0], label, val, find, part, n0, (void (***)(void)) g0);
1269:   PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GT1], label, val, find, part, n1, (void (***)(void)) g1);
1270:   PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GT2], label, val, find, part, n2, (void (***)(void)) g2);
1271:   PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_GT3], label, val, find, part, n3, (void (***)(void)) g3);
1272:   return(0);
1273: }

1275: PetscErrorCode PetscWeakFormAddDynamicJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part,
1276:                                         void (*g0)(PetscInt, PetscInt, PetscInt,
1277:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1278:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1279:                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1280:                                         void (*g1)(PetscInt, PetscInt, PetscInt,
1281:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1282:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1283:                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1284:                                         void (*g2)(PetscInt, PetscInt, PetscInt,
1285:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1286:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1287:                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1288:                                         void (*g3)(PetscInt, PetscInt, PetscInt,
1289:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1290:                                                    const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1291:                                                    PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
1292: {
1293:   PetscInt       find = f*wf->Nf + g;

1297:   PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GT0], label, val, find, part, (void (*)(void)) g0);
1298:   PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GT1], label, val, find, part, (void (*)(void)) g1);
1299:   PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GT2], label, val, find, part, (void (*)(void)) g2);
1300:   PetscWeakFormAddFunction_Private(wf, wf->form[PETSC_WF_GT3], label, val, find, part, (void (*)(void)) g3);
1301:   return(0);
1302: }

1304: PetscErrorCode PetscWeakFormSetDynamicJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part,
1305:                                                PetscInt n0,
1306:                                                void (**g0)(PetscInt, PetscInt, PetscInt,
1307:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1308:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1309:                                                           PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1310:                                                PetscInt n1,
1311:                                                void (**g1)(PetscInt, PetscInt, PetscInt,
1312:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1313:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1314:                                                           PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1315:                                                PetscInt n2,
1316:                                                void (**g2)(PetscInt, PetscInt, PetscInt,
1317:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1318:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1319:                                                           PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1320:                                                PetscInt n3,
1321:                                                void (**g3)(PetscInt, PetscInt, PetscInt,
1322:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1323:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1324:                                                           PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
1325: {
1326:   PetscInt       find = f*wf->Nf + g;

1330:   PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GT0], label, val, find, part, n0, (void (**)(void)) g0);
1331:   PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GT1], label, val, find, part, n1, (void (**)(void)) g1);
1332:   PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GT2], label, val, find, part, n2, (void (**)(void)) g2);
1333:   PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_GT3], label, val, find, part, n3, (void (**)(void)) g3);
1334:   return(0);
1335: }

1337: PetscErrorCode PetscWeakFormSetIndexDynamicJacobian(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt g, PetscInt part,
1338:                                                PetscInt i0,
1339:                                                void (*g0)(PetscInt, PetscInt, PetscInt,
1340:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1341:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1342:                                                           PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1343:                                                PetscInt i1,
1344:                                                void (*g1)(PetscInt, PetscInt, PetscInt,
1345:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1346:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1347:                                                           PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1348:                                                PetscInt i2,
1349:                                                void (*g2)(PetscInt, PetscInt, PetscInt,
1350:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1351:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1352:                                                           PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]),
1353:                                                PetscInt i3,
1354:                                                void (*g3)(PetscInt, PetscInt, PetscInt,
1355:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1356:                                                           const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[],
1357:                                                           PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]))
1358: {
1359:   PetscInt       find = f*wf->Nf + g;

1363:   PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GT0], label, val, find, part, i0, (void (*)(void)) g0);
1364:   PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GT1], label, val, find, part, i1, (void (*)(void)) g1);
1365:   PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GT2], label, val, find, part, i2, (void (*)(void)) g2);
1366:   PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_GT3], label, val, find, part, i3, (void (*)(void)) g3);
1367:   return(0);
1368: }

1370: PetscErrorCode PetscWeakFormGetRiemannSolver(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part, PetscInt *n,
1371:                                              void (***r)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *))
1372: {

1376:   PetscWeakFormGetFunction_Private(wf, wf->form[PETSC_WF_R], label, val, f, part, n, (void (***)(void)) r);
1377:   return(0);
1378: }

1380: PetscErrorCode PetscWeakFormSetRiemannSolver(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part,
1381:                                              PetscInt n,
1382:                                              void (**r)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *))
1383: {

1387:   PetscWeakFormSetFunction_Private(wf, wf->form[PETSC_WF_R], label, val, f, part, n, (void (**)(void)) r);
1388:   return(0);
1389: }

1391: PetscErrorCode PetscWeakFormSetIndexRiemannSolver(PetscWeakForm wf, DMLabel label, PetscInt val, PetscInt f, PetscInt part,
1392:                                                   PetscInt i,
1393:                                                   void (*r)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *))
1394: {

1398:   PetscWeakFormSetIndexFunction_Private(wf, wf->form[PETSC_WF_R], label, val, f, part, i, (void (*)(void)) r);
1399:   return(0);
1400: }

1402: /*@
1403:   PetscWeakFormGetNumFields - Returns the number of fields

1405:   Not collective

1407:   Input Parameter:
1408: . wf - The PetscWeakForm object

1410:   Output Parameter:
1411: . Nf - The number of fields

1413:   Level: beginner

1415: .seealso: PetscWeakFormSetNumFields(), PetscWeakFormCreate()
1416: @*/
1417: PetscErrorCode PetscWeakFormGetNumFields(PetscWeakForm wf, PetscInt *Nf)
1418: {
1422:   *Nf = wf->Nf;
1423:   return(0);
1424: }

1426: /*@
1427:   PetscWeakFormSetNumFields - Sets the number of fields

1429:   Not collective

1431:   Input Parameters:
1432: + wf - The PetscWeakForm object
1433: - Nf - The number of fields

1435:   Level: beginner

1437: .seealso: PetscWeakFormGetNumFields(), PetscWeakFormCreate()
1438: @*/
1439: PetscErrorCode PetscWeakFormSetNumFields(PetscWeakForm wf, PetscInt Nf)
1440: {
1443:   wf->Nf = Nf;
1444:   return(0);
1445: }

1447: /*@
1448:   PetscWeakFormDestroy - Destroys a PetscWeakForm object

1450:   Collective on wf

1452:   Input Parameter:
1453: . wf - the PetscWeakForm object to destroy

1455:   Level: developer

1457: .seealso PetscWeakFormCreate(), PetscWeakFormView()
1458: @*/
1459: PetscErrorCode PetscWeakFormDestroy(PetscWeakForm *wf)
1460: {
1461:   PetscInt       f;

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

1468:   if (--((PetscObject)(*wf))->refct > 0) {*wf = NULL; return(0);}
1469:   ((PetscObject) (*wf))->refct = 0;
1470:   PetscChunkBufferDestroy(&(*wf)->funcs);
1471:   for (f = 0; f < PETSC_NUM_WF; ++f) {PetscHMapFormDestroy(&(*wf)->form[f]);}
1472:   PetscFree((*wf)->form);
1473:   PetscHeaderDestroy(wf);
1474:   return(0);
1475: }

1477: static PetscErrorCode PetscWeakFormViewTable_Ascii(PetscWeakForm wf, PetscViewer viewer, PetscBool splitField, const char tableName[], PetscHMapForm map)
1478: {
1479:   PetscInt       Nf = wf->Nf, Nk, k;

1483:   PetscHMapFormGetSize(map, &Nk);
1484:   if (Nk) {
1485:     PetscFormKey *keys;
1486:     void           (**funcs)(void);
1487:     const char       *name;
1488:     PetscBool         showPart = PETSC_FALSE;
1489:     PetscInt          off = 0, n, i;

1491:     PetscMalloc1(Nk, &keys);
1492:     PetscHMapFormGetKeys(map, &off, keys);
1493:     PetscViewerASCIIPrintf(viewer, "%s\n", tableName);
1494:     PetscViewerASCIIPushTab(viewer);
1495:     for (k = 0; k < Nk; ++k) {
1496:       if (keys[k].part != 0) showPart = PETSC_TRUE;
1497:     }
1498:     for (k = 0; k < Nk; ++k) {
1499:       if (keys[k].label) {
1500:         PetscObjectGetName((PetscObject) keys[k].label, &name);
1501:         PetscViewerASCIIPrintf(viewer, "(%s:%p, %D) ", name, keys[k].label, keys[k].value);
1502:       } else {PetscViewerASCIIPrintf(viewer, "");}
1503:       PetscViewerASCIIUseTabs(viewer, PETSC_FALSE);
1504:       if (splitField) {PetscViewerASCIIPrintf(viewer, "(%D, %D) ", keys[k].field/Nf, keys[k].field%Nf);}
1505:       else            {PetscViewerASCIIPrintf(viewer, "(%D) ", keys[k].field);}
1506:       if (showPart)   {PetscViewerASCIIPrintf(viewer, "(%D) ", keys[k].part);}
1507:       PetscWeakFormGetFunction_Private(wf, map, keys[k].label, keys[k].value, keys[k].field, keys[k].part, &n, &funcs);
1508:       for (i = 0; i < n; ++i) {
1509:         char *fname;

1511:         if (i > 0) {PetscViewerASCIIPrintf(viewer, ", ");}
1512:         PetscDLAddr(funcs[i], &fname);
1513:         if (fname) {PetscViewerASCIIPrintf(viewer, "%s", fname);}
1514:         else       {PetscViewerASCIIPrintf(viewer, "%p", funcs[i]);}
1515:         PetscFree(fname);
1516:       }
1517:       PetscViewerASCIIPrintf(viewer, "\n");
1518:       PetscViewerASCIIUseTabs(viewer, PETSC_TRUE);
1519:     }
1520:     PetscViewerASCIIPopTab(viewer);
1521:     PetscFree(keys);
1522:   }
1523:   return(0);
1524: }

1526: static PetscErrorCode PetscWeakFormView_Ascii(PetscWeakForm wf, PetscViewer viewer)
1527: {
1528:   PetscViewerFormat format;
1529:   PetscInt          f;
1530:   PetscErrorCode    ierr;

1533:   PetscViewerGetFormat(viewer, &format);
1534:   PetscViewerASCIIPrintf(viewer, "Weak Form System with %d fields\n", wf->Nf);
1535:   PetscViewerASCIIPushTab(viewer);
1536:   for (f = 0; f < PETSC_NUM_WF; ++f) {
1537:     PetscWeakFormViewTable_Ascii(wf, viewer, PETSC_TRUE, PetscWeakFormKinds[f], wf->form[f]);
1538:   }
1539:   PetscViewerASCIIPopTab(viewer);
1540:   return(0);
1541: }

1543: /*@C
1544:   PetscWeakFormView - Views a PetscWeakForm

1546:   Collective on wf

1548:   Input Parameters:
1549: + wf - the PetscWeakForm object to view
1550: - v  - the viewer

1552:   Level: developer

1554: .seealso PetscWeakFormDestroy(), PetscWeakFormCreate()
1555: @*/
1556: PetscErrorCode PetscWeakFormView(PetscWeakForm wf, PetscViewer v)
1557: {
1558:   PetscBool      iascii;

1563:   if (!v) {PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) wf), &v);}
1565:   PetscObjectTypeCompare((PetscObject) v, PETSCVIEWERASCII, &iascii);
1566:   if (iascii) {PetscWeakFormView_Ascii(wf, v);}
1567:   if (wf->ops->view) {(*wf->ops->view)(wf, v);}
1568:   return(0);
1569: }

1571: /*@
1572:   PetscWeakFormCreate - Creates an empty PetscWeakForm object.

1574:   Collective

1576:   Input Parameter:
1577: . comm - The communicator for the PetscWeakForm object

1579:   Output Parameter:
1580: . wf - The PetscWeakForm object

1582:   Level: beginner

1584: .seealso: PetscDS, PetscWeakFormDestroy()
1585: @*/
1586: PetscErrorCode PetscWeakFormCreate(MPI_Comm comm, PetscWeakForm *wf)
1587: {
1588:   PetscWeakForm  p;
1589:   PetscInt       f;

1594:   *wf  = NULL;
1595:   PetscDSInitializePackage();

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

1599:   p->Nf = 0;
1600:   PetscChunkBufferCreate(sizeof(&PetscWeakFormCreate), 2, &p->funcs);
1601:   PetscMalloc1(PETSC_NUM_WF, &p->form);
1602:   for (f = 0; f < PETSC_NUM_WF; ++f) {PetscHMapFormCreate(&p->form[f]);}
1603:   *wf = p;
1604:   return(0);
1605: }