xref: /petsc/src/dm/dt/interface/dtds.c (revision f744cafa0d2dce869fb9012d006ed1e50ca28605)
1 #include <petsc/private/petscdsimpl.h> /*I "petscds.h" I*/
2 
3 PetscClassId PETSCDS_CLASSID = 0;
4 
5 PetscFunctionList PetscDSList              = NULL;
6 PetscBool         PetscDSRegisterAllCalled = PETSC_FALSE;
7 
8 #undef __FUNCT__
9 #define __FUNCT__ "PetscDSRegister"
10 /*@C
11   PetscDSRegister - Adds a new PetscDS implementation
12 
13   Not Collective
14 
15   Input Parameters:
16 + name        - The name of a new user-defined creation routine
17 - create_func - The creation routine itself
18 
19   Notes:
20   PetscDSRegister() may be called multiple times to add several user-defined PetscDSs
21 
22   Sample usage:
23 .vb
24     PetscDSRegister("my_ds", MyPetscDSCreate);
25 .ve
26 
27   Then, your PetscDS type can be chosen with the procedural interface via
28 .vb
29     PetscDSCreate(MPI_Comm, PetscDS *);
30     PetscDSSetType(PetscDS, "my_ds");
31 .ve
32    or at runtime via the option
33 .vb
34     -petscds_type my_ds
35 .ve
36 
37   Level: advanced
38 
39 .keywords: PetscDS, register
40 .seealso: PetscDSRegisterAll(), PetscDSRegisterDestroy()
41 
42 @*/
43 PetscErrorCode PetscDSRegister(const char sname[], PetscErrorCode (*function)(PetscDS))
44 {
45   PetscErrorCode ierr;
46 
47   PetscFunctionBegin;
48   ierr = PetscFunctionListAdd(&PetscDSList, sname, function);CHKERRQ(ierr);
49   PetscFunctionReturn(0);
50 }
51 
52 #undef __FUNCT__
53 #define __FUNCT__ "PetscDSSetType"
54 /*@C
55   PetscDSSetType - Builds a particular PetscDS
56 
57   Collective on PetscDS
58 
59   Input Parameters:
60 + prob - The PetscDS object
61 - name - The kind of system
62 
63   Options Database Key:
64 . -petscds_type <type> - Sets the PetscDS type; use -help for a list of available types
65 
66   Level: intermediate
67 
68 .keywords: PetscDS, set, type
69 .seealso: PetscDSGetType(), PetscDSCreate()
70 @*/
71 PetscErrorCode PetscDSSetType(PetscDS prob, PetscDSType name)
72 {
73   PetscErrorCode (*r)(PetscDS);
74   PetscBool      match;
75   PetscErrorCode ierr;
76 
77   PetscFunctionBegin;
78   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
79   ierr = PetscObjectTypeCompare((PetscObject) prob, name, &match);CHKERRQ(ierr);
80   if (match) PetscFunctionReturn(0);
81 
82   ierr = PetscDSRegisterAll();CHKERRQ(ierr);
83   ierr = PetscFunctionListFind(PetscDSList, name, &r);CHKERRQ(ierr);
84   if (!r) SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscDS type: %s", name);
85 
86   if (prob->ops->destroy) {
87     ierr             = (*prob->ops->destroy)(prob);CHKERRQ(ierr);
88     prob->ops->destroy = NULL;
89   }
90   ierr = (*r)(prob);CHKERRQ(ierr);
91   ierr = PetscObjectChangeTypeName((PetscObject) prob, name);CHKERRQ(ierr);
92   PetscFunctionReturn(0);
93 }
94 
95 #undef __FUNCT__
96 #define __FUNCT__ "PetscDSGetType"
97 /*@C
98   PetscDSGetType - Gets the PetscDS type name (as a string) from the object.
99 
100   Not Collective
101 
102   Input Parameter:
103 . prob  - The PetscDS
104 
105   Output Parameter:
106 . name - The PetscDS type name
107 
108   Level: intermediate
109 
110 .keywords: PetscDS, get, type, name
111 .seealso: PetscDSSetType(), PetscDSCreate()
112 @*/
113 PetscErrorCode PetscDSGetType(PetscDS prob, PetscDSType *name)
114 {
115   PetscErrorCode ierr;
116 
117   PetscFunctionBegin;
118   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
119   PetscValidPointer(name, 2);
120   ierr = PetscDSRegisterAll();CHKERRQ(ierr);
121   *name = ((PetscObject) prob)->type_name;
122   PetscFunctionReturn(0);
123 }
124 
125 #undef __FUNCT__
126 #define __FUNCT__ "PetscDSView_Ascii"
127 static PetscErrorCode PetscDSView_Ascii(PetscDS prob, PetscViewer viewer)
128 {
129   PetscViewerFormat format;
130   PetscInt          f;
131   PetscErrorCode    ierr;
132 
133   PetscFunctionBegin;
134   ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
135   ierr = PetscViewerASCIIPrintf(viewer, "Discrete System with %d fields\n", prob->Nf);CHKERRQ(ierr);
136   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
137   for (f = 0; f < prob->Nf; ++f) {
138     PetscObject  obj;
139     PetscClassId id;
140     const char  *name;
141     PetscInt     Nc;
142 
143     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
144     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
145     ierr = PetscObjectGetName(obj, &name);CHKERRQ(ierr);
146     ierr = PetscViewerASCIIPrintf(viewer, "Field %s", name ? name : "<unknown>");CHKERRQ(ierr);
147     if (id == PETSCFE_CLASSID)      {
148       ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);
149       ierr = PetscViewerASCIIPrintf(viewer, " FEM");CHKERRQ(ierr);
150     } else if (id == PETSCFV_CLASSID) {
151       ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);
152       ierr = PetscViewerASCIIPrintf(viewer, " FVM");CHKERRQ(ierr);
153     }
154     else SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
155     if (Nc > 1) {ierr = PetscViewerASCIIPrintf(viewer, "%d components", Nc);CHKERRQ(ierr);}
156     else        {ierr = PetscViewerASCIIPrintf(viewer, "%d component ", Nc);CHKERRQ(ierr);}
157     if (prob->implicit[f]) {ierr = PetscViewerASCIIPrintf(viewer, " (implicit)");CHKERRQ(ierr);}
158     else                   {ierr = PetscViewerASCIIPrintf(viewer, " (explicit)");CHKERRQ(ierr);}
159     if (prob->adjacency[f*2+0]) {
160       if (prob->adjacency[f*2+1]) {ierr = PetscViewerASCIIPrintf(viewer, " (adj FVM++)");CHKERRQ(ierr);}
161       else                        {ierr = PetscViewerASCIIPrintf(viewer, " (adj FVM)");CHKERRQ(ierr);}
162     } else {
163       if (prob->adjacency[f*2+1]) {ierr = PetscViewerASCIIPrintf(viewer, " (adj FEM)");CHKERRQ(ierr);}
164       else                        {ierr = PetscViewerASCIIPrintf(viewer, " (adj FUNKY)");CHKERRQ(ierr);}
165     }
166     ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr);
167     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
168       if (id == PETSCFE_CLASSID)      {ierr = PetscFEView((PetscFE) obj, viewer);CHKERRQ(ierr);}
169       else if (id == PETSCFV_CLASSID) {ierr = PetscFVView((PetscFV) obj, viewer);CHKERRQ(ierr);}
170     }
171   }
172   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
173   PetscFunctionReturn(0);
174 }
175 
176 #undef __FUNCT__
177 #define __FUNCT__ "PetscDSView"
178 /*@C
179   PetscDSView - Views a PetscDS
180 
181   Collective on PetscDS
182 
183   Input Parameter:
184 + prob - the PetscDS object to view
185 - v  - the viewer
186 
187   Level: developer
188 
189 .seealso PetscDSDestroy()
190 @*/
191 PetscErrorCode PetscDSView(PetscDS prob, PetscViewer v)
192 {
193   PetscBool      iascii;
194   PetscErrorCode ierr;
195 
196   PetscFunctionBegin;
197   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
198   if (!v) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) prob), &v);CHKERRQ(ierr);}
199   else    {PetscValidHeaderSpecific(v, PETSC_VIEWER_CLASSID, 2);}
200   ierr = PetscObjectTypeCompare((PetscObject) v, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
201   if (iascii) {ierr = PetscDSView_Ascii(prob, v);CHKERRQ(ierr);}
202   if (prob->ops->view) {ierr = (*prob->ops->view)(prob, v);CHKERRQ(ierr);}
203   PetscFunctionReturn(0);
204 }
205 
206 #undef __FUNCT__
207 #define __FUNCT__ "PetscDSSetFromOptions"
208 /*@
209   PetscDSSetFromOptions - sets parameters in a PetscDS from the options database
210 
211   Collective on PetscDS
212 
213   Input Parameter:
214 . prob - the PetscDS object to set options for
215 
216   Options Database:
217 
218   Level: developer
219 
220 .seealso PetscDSView()
221 @*/
222 PetscErrorCode PetscDSSetFromOptions(PetscDS prob)
223 {
224   DSBoundary     b;
225   const char    *defaultType;
226   char           name[256];
227   PetscBool      flg;
228   PetscErrorCode ierr;
229 
230   PetscFunctionBegin;
231   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
232   if (!((PetscObject) prob)->type_name) {
233     defaultType = PETSCDSBASIC;
234   } else {
235     defaultType = ((PetscObject) prob)->type_name;
236   }
237   ierr = PetscDSRegisterAll();CHKERRQ(ierr);
238 
239   ierr = PetscObjectOptionsBegin((PetscObject) prob);CHKERRQ(ierr);
240   for (b = prob->boundary; b; b = b->next) {
241     char       optname[1024];
242     PetscInt   ids[1024], len = 1024;
243     PetscBool  flg;
244 
245     ierr = PetscSNPrintf(optname, sizeof(optname), "-bc_%s", b->name);CHKERRQ(ierr);
246     ierr = PetscMemzero(ids, sizeof(ids));CHKERRQ(ierr);
247     ierr = PetscOptionsIntArray(optname, "List of boundary IDs", "", ids, &len, &flg);CHKERRQ(ierr);
248     if (flg) {
249       b->numids = len;
250       ierr = PetscFree(b->ids);CHKERRQ(ierr);
251       ierr = PetscMalloc1(len, &b->ids);CHKERRQ(ierr);
252       ierr = PetscMemcpy(b->ids, ids, len*sizeof(PetscInt));CHKERRQ(ierr);
253     }
254     ierr = PetscSNPrintf(optname, sizeof(optname), "-bc_%s_comp", b->name);CHKERRQ(ierr);
255     ierr = PetscMemzero(ids, sizeof(ids));CHKERRQ(ierr);
256     ierr = PetscOptionsIntArray(optname, "List of boundary field components", "", ids, &len, &flg);CHKERRQ(ierr);
257     if (flg) {
258       b->numcomps = len;
259       ierr = PetscFree(b->comps);CHKERRQ(ierr);
260       ierr = PetscMalloc1(len, &b->comps);CHKERRQ(ierr);
261       ierr = PetscMemcpy(b->comps, ids, len*sizeof(PetscInt));CHKERRQ(ierr);
262     }
263   }
264   ierr = PetscOptionsFList("-petscds_type", "Discrete System", "PetscDSSetType", PetscDSList, defaultType, name, 256, &flg);CHKERRQ(ierr);
265   if (flg) {
266     ierr = PetscDSSetType(prob, name);CHKERRQ(ierr);
267   } else if (!((PetscObject) prob)->type_name) {
268     ierr = PetscDSSetType(prob, defaultType);CHKERRQ(ierr);
269   }
270   if (prob->ops->setfromoptions) {ierr = (*prob->ops->setfromoptions)(prob);CHKERRQ(ierr);}
271   /* process any options handlers added with PetscObjectAddOptionsHandler() */
272   ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) prob);CHKERRQ(ierr);
273   ierr = PetscOptionsEnd();CHKERRQ(ierr);
274   ierr = PetscDSViewFromOptions(prob, NULL, "-petscds_view");CHKERRQ(ierr);
275   PetscFunctionReturn(0);
276 }
277 
278 #undef __FUNCT__
279 #define __FUNCT__ "PetscDSSetUp"
280 /*@C
281   PetscDSSetUp - Construct data structures for the PetscDS
282 
283   Collective on PetscDS
284 
285   Input Parameter:
286 . prob - the PetscDS object to setup
287 
288   Level: developer
289 
290 .seealso PetscDSView(), PetscDSDestroy()
291 @*/
292 PetscErrorCode PetscDSSetUp(PetscDS prob)
293 {
294   const PetscInt Nf = prob->Nf;
295   PetscInt       dim, work, NcMax = 0, NqMax = 0, f;
296   PetscErrorCode ierr;
297 
298   PetscFunctionBegin;
299   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
300   if (prob->setup) PetscFunctionReturn(0);
301   /* Calculate sizes */
302   ierr = PetscDSGetSpatialDimension(prob, &dim);CHKERRQ(ierr);
303   prob->totDim = prob->totComp = 0;
304   ierr = PetscMalloc2(Nf,&prob->Nc,Nf,&prob->Nb);CHKERRQ(ierr);
305   ierr = PetscCalloc2(Nf+1,&prob->off,Nf+1,&prob->offDer);CHKERRQ(ierr);
306   ierr = PetscMalloc4(Nf,&prob->basis,Nf,&prob->basisDer,Nf,&prob->basisFace,Nf,&prob->basisDerFace);CHKERRQ(ierr);
307   for (f = 0; f < Nf; ++f) {
308     PetscObject     obj;
309     PetscClassId    id;
310     PetscQuadrature q;
311     PetscInt        Nq = 0, Nb, Nc;
312 
313     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
314     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
315     if (id == PETSCFE_CLASSID)      {
316       PetscFE fe = (PetscFE) obj;
317 
318       ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
319       ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
320       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
321       ierr = PetscFEGetDefaultTabulation(fe, &prob->basis[f], &prob->basisDer[f], NULL);CHKERRQ(ierr);
322       ierr = PetscFEGetFaceTabulation(fe, &prob->basisFace[f], &prob->basisDerFace[f], NULL);CHKERRQ(ierr);
323     } else if (id == PETSCFV_CLASSID) {
324       PetscFV fv = (PetscFV) obj;
325 
326       ierr = PetscFVGetQuadrature(fv, &q);CHKERRQ(ierr);
327       Nb   = 1;
328       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
329       ierr = PetscFVGetDefaultTabulation(fv, &prob->basis[f], &prob->basisDer[f], NULL);CHKERRQ(ierr);
330       /* TODO: should PetscFV also have face tabulation? Otherwise there will be a null pointer in prob->basisFace */
331     } else SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
332     prob->Nc[f]       = Nc;
333     prob->Nb[f]       = Nb;
334     prob->off[f+1]    = Nc     + prob->off[f];
335     prob->offDer[f+1] = Nc*dim + prob->offDer[f];
336     if (q) {ierr = PetscQuadratureGetData(q, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);}
337     NqMax          = PetscMax(NqMax, Nq);
338     NcMax          = PetscMax(NcMax, Nc);
339     prob->totDim  += Nb*Nc;
340     prob->totComp += Nc;
341   }
342   work = PetscMax(prob->totComp*dim, PetscSqr(NcMax*dim));
343   /* Allocate works space */
344   ierr = PetscMalloc5(prob->totComp,&prob->u,prob->totComp,&prob->u_t,prob->totComp*dim,&prob->u_x,dim,&prob->x,work,&prob->refSpaceDer);CHKERRQ(ierr);
345   ierr = PetscMalloc6(NqMax*NcMax,&prob->f0,NqMax*NcMax*dim,&prob->f1,NqMax*NcMax*NcMax,&prob->g0,NqMax*NcMax*NcMax*dim,&prob->g1,NqMax*NcMax*NcMax*dim,&prob->g2,NqMax*NcMax*NcMax*dim*dim,&prob->g3);CHKERRQ(ierr);
346   if (prob->ops->setup) {ierr = (*prob->ops->setup)(prob);CHKERRQ(ierr);}
347   prob->setup = PETSC_TRUE;
348   PetscFunctionReturn(0);
349 }
350 
351 #undef __FUNCT__
352 #define __FUNCT__ "PetscDSDestroyStructs_Static"
353 static PetscErrorCode PetscDSDestroyStructs_Static(PetscDS prob)
354 {
355   PetscErrorCode ierr;
356 
357   PetscFunctionBegin;
358   ierr = PetscFree2(prob->Nc,prob->Nb);CHKERRQ(ierr);
359   ierr = PetscFree2(prob->off,prob->offDer);CHKERRQ(ierr);
360   ierr = PetscFree4(prob->basis,prob->basisDer,prob->basisFace,prob->basisDerFace);CHKERRQ(ierr);
361   ierr = PetscFree5(prob->u,prob->u_t,prob->u_x,prob->x,prob->refSpaceDer);CHKERRQ(ierr);
362   ierr = PetscFree6(prob->f0,prob->f1,prob->g0,prob->g1,prob->g2,prob->g3);CHKERRQ(ierr);
363   PetscFunctionReturn(0);
364 }
365 
366 #undef __FUNCT__
367 #define __FUNCT__ "PetscDSEnlarge_Static"
368 static PetscErrorCode PetscDSEnlarge_Static(PetscDS prob, PetscInt NfNew)
369 {
370   PetscObject      *tmpd;
371   PetscBool        *tmpi, *tmpa;
372   PetscPointFunc   *tmpobj, *tmpf;
373   PetscPointJac    *tmpg, *tmpgp, *tmpgt;
374   PetscBdPointFunc *tmpfbd;
375   PetscBdPointJac  *tmpgbd;
376   PetscRiemannFunc *tmpr;
377   void            **tmpctx;
378   PetscInt          Nf = prob->Nf, f, i;
379   PetscErrorCode    ierr;
380 
381   PetscFunctionBegin;
382   if (Nf >= NfNew) PetscFunctionReturn(0);
383   prob->setup = PETSC_FALSE;
384   ierr = PetscDSDestroyStructs_Static(prob);CHKERRQ(ierr);
385   ierr = PetscMalloc3(NfNew, &tmpd, NfNew, &tmpi, NfNew*2, &tmpa);CHKERRQ(ierr);
386   for (f = 0; f < Nf; ++f) {tmpd[f] = prob->disc[f]; tmpi[f] = prob->implicit[f]; for (i = 0; i < 2; ++i) tmpa[f*2+i] = prob->adjacency[f*2+i];}
387   for (f = Nf; f < NfNew; ++f) {tmpd[f] = NULL; tmpi[f] = PETSC_TRUE; tmpa[f*2+0] = PETSC_FALSE; tmpa[f*2+1] = PETSC_TRUE;}
388   ierr = PetscFree3(prob->disc, prob->implicit, prob->adjacency);CHKERRQ(ierr);
389   prob->Nf        = NfNew;
390   prob->disc      = tmpd;
391   prob->implicit  = tmpi;
392   prob->adjacency = tmpa;
393   ierr = PetscCalloc7(NfNew, &tmpobj, NfNew*2, &tmpf, NfNew*NfNew*4, &tmpg, NfNew*NfNew*4, &tmpgp, NfNew*NfNew*4, &tmpgt, NfNew, &tmpr, NfNew, &tmpctx);CHKERRQ(ierr);
394   for (f = 0; f < Nf; ++f) tmpobj[f] = prob->obj[f];
395   for (f = 0; f < Nf*2; ++f) tmpf[f] = prob->f[f];
396   for (f = 0; f < Nf*Nf*4; ++f) tmpg[f] = prob->g[f];
397   for (f = 0; f < Nf*Nf*4; ++f) tmpgp[f] = prob->gp[f];
398   for (f = 0; f < Nf; ++f) tmpr[f] = prob->r[f];
399   for (f = 0; f < Nf; ++f) tmpctx[f] = prob->ctx[f];
400   for (f = Nf; f < NfNew; ++f) tmpobj[f] = NULL;
401   for (f = Nf*2; f < NfNew*2; ++f) tmpf[f] = NULL;
402   for (f = Nf*Nf*4; f < NfNew*NfNew*4; ++f) tmpg[f] = NULL;
403   for (f = Nf*Nf*4; f < NfNew*NfNew*4; ++f) tmpgp[f] = NULL;
404   for (f = Nf*Nf*4; f < NfNew*NfNew*4; ++f) tmpgt[f] = NULL;
405   for (f = Nf; f < NfNew; ++f) tmpr[f] = NULL;
406   for (f = Nf; f < NfNew; ++f) tmpctx[f] = NULL;
407   ierr = PetscFree7(prob->obj, prob->f, prob->g, prob->gp, prob->gt, prob->r, prob->ctx);CHKERRQ(ierr);
408   prob->obj = tmpobj;
409   prob->f   = tmpf;
410   prob->g   = tmpg;
411   prob->gp  = tmpgp;
412   prob->gt  = tmpgt;
413   prob->r   = tmpr;
414   prob->ctx = tmpctx;
415   ierr = PetscCalloc2(NfNew*2, &tmpfbd, NfNew*NfNew*4, &tmpgbd);CHKERRQ(ierr);
416   for (f = 0; f < Nf*2; ++f) tmpfbd[f] = prob->fBd[f];
417   for (f = 0; f < Nf*Nf*4; ++f) tmpgbd[f] = prob->gBd[f];
418   for (f = Nf*2; f < NfNew*2; ++f) tmpfbd[f] = NULL;
419   for (f = Nf*Nf*4; f < NfNew*NfNew*4; ++f) tmpgbd[f] = NULL;
420   ierr = PetscFree2(prob->fBd, prob->gBd);CHKERRQ(ierr);
421   prob->fBd = tmpfbd;
422   prob->gBd = tmpgbd;
423   PetscFunctionReturn(0);
424 }
425 
426 #undef __FUNCT__
427 #define __FUNCT__ "PetscDSDestroy"
428 /*@
429   PetscDSDestroy - Destroys a PetscDS object
430 
431   Collective on PetscDS
432 
433   Input Parameter:
434 . prob - the PetscDS object to destroy
435 
436   Level: developer
437 
438 .seealso PetscDSView()
439 @*/
440 PetscErrorCode PetscDSDestroy(PetscDS *prob)
441 {
442   PetscInt       f;
443   DSBoundary     next;
444   PetscErrorCode ierr;
445 
446   PetscFunctionBegin;
447   if (!*prob) PetscFunctionReturn(0);
448   PetscValidHeaderSpecific((*prob), PETSCDS_CLASSID, 1);
449 
450   if (--((PetscObject)(*prob))->refct > 0) {*prob = 0; PetscFunctionReturn(0);}
451   ((PetscObject) (*prob))->refct = 0;
452   ierr = PetscDSDestroyStructs_Static(*prob);CHKERRQ(ierr);
453   for (f = 0; f < (*prob)->Nf; ++f) {
454     ierr = PetscObjectDereference((*prob)->disc[f]);CHKERRQ(ierr);
455   }
456   ierr = PetscFree3((*prob)->disc, (*prob)->implicit, (*prob)->adjacency);CHKERRQ(ierr);
457   ierr = PetscFree7((*prob)->obj,(*prob)->f,(*prob)->g,(*prob)->gp,(*prob)->gt,(*prob)->r,(*prob)->ctx);CHKERRQ(ierr);
458   ierr = PetscFree2((*prob)->fBd,(*prob)->gBd);CHKERRQ(ierr);
459   if ((*prob)->ops->destroy) {ierr = (*(*prob)->ops->destroy)(*prob);CHKERRQ(ierr);}
460   next = (*prob)->boundary;
461   while (next) {
462     DSBoundary b = next;
463 
464     next = b->next;
465     ierr = PetscFree(b->comps);CHKERRQ(ierr);
466     ierr = PetscFree(b->ids);CHKERRQ(ierr);
467     ierr = PetscFree(b->name);CHKERRQ(ierr);
468     ierr = PetscFree(b->labelname);CHKERRQ(ierr);
469     ierr = PetscFree(b);CHKERRQ(ierr);
470   }
471   ierr = PetscHeaderDestroy(prob);CHKERRQ(ierr);
472   PetscFunctionReturn(0);
473 }
474 
475 #undef __FUNCT__
476 #define __FUNCT__ "PetscDSCreate"
477 /*@
478   PetscDSCreate - Creates an empty PetscDS object. The type can then be set with PetscDSSetType().
479 
480   Collective on MPI_Comm
481 
482   Input Parameter:
483 . comm - The communicator for the PetscDS object
484 
485   Output Parameter:
486 . prob - The PetscDS object
487 
488   Level: beginner
489 
490 .seealso: PetscDSSetType(), PETSCDSBASIC
491 @*/
492 PetscErrorCode PetscDSCreate(MPI_Comm comm, PetscDS *prob)
493 {
494   PetscDS   p;
495   PetscErrorCode ierr;
496 
497   PetscFunctionBegin;
498   PetscValidPointer(prob, 2);
499   *prob  = NULL;
500   ierr = PetscDSInitializePackage();CHKERRQ(ierr);
501 
502   ierr = PetscHeaderCreate(p, PETSCDS_CLASSID, "PetscDS", "Discrete System", "PetscDS", comm, PetscDSDestroy, PetscDSView);CHKERRQ(ierr);
503 
504   p->Nf    = 0;
505   p->setup = PETSC_FALSE;
506 
507   *prob = p;
508   PetscFunctionReturn(0);
509 }
510 
511 #undef __FUNCT__
512 #define __FUNCT__ "PetscDSGetNumFields"
513 /*@
514   PetscDSGetNumFields - Returns the number of fields in the DS
515 
516   Not collective
517 
518   Input Parameter:
519 . prob - The PetscDS object
520 
521   Output Parameter:
522 . Nf - The number of fields
523 
524   Level: beginner
525 
526 .seealso: PetscDSGetSpatialDimension(), PetscDSCreate()
527 @*/
528 PetscErrorCode PetscDSGetNumFields(PetscDS prob, PetscInt *Nf)
529 {
530   PetscFunctionBegin;
531   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
532   PetscValidPointer(Nf, 2);
533   *Nf = prob->Nf;
534   PetscFunctionReturn(0);
535 }
536 
537 #undef __FUNCT__
538 #define __FUNCT__ "PetscDSGetSpatialDimension"
539 /*@
540   PetscDSGetSpatialDimension - Returns the spatial dimension of the DS
541 
542   Not collective
543 
544   Input Parameter:
545 . prob - The PetscDS object
546 
547   Output Parameter:
548 . dim - The spatial dimension
549 
550   Level: beginner
551 
552 .seealso: PetscDSGetNumFields(), PetscDSCreate()
553 @*/
554 PetscErrorCode PetscDSGetSpatialDimension(PetscDS prob, PetscInt *dim)
555 {
556   PetscErrorCode ierr;
557 
558   PetscFunctionBegin;
559   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
560   PetscValidPointer(dim, 2);
561   *dim = 0;
562   if (prob->Nf) {
563     PetscObject  obj;
564     PetscClassId id;
565 
566     ierr = PetscDSGetDiscretization(prob, 0, &obj);CHKERRQ(ierr);
567     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
568     if (id == PETSCFE_CLASSID)      {ierr = PetscFEGetSpatialDimension((PetscFE) obj, dim);CHKERRQ(ierr);}
569     else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetSpatialDimension((PetscFV) obj, dim);CHKERRQ(ierr);}
570     else SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", 0);
571   }
572   PetscFunctionReturn(0);
573 }
574 
575 #undef __FUNCT__
576 #define __FUNCT__ "PetscDSGetTotalDimension"
577 /*@
578   PetscDSGetTotalDimension - Returns the total size of the approximation space for this system
579 
580   Not collective
581 
582   Input Parameter:
583 . prob - The PetscDS object
584 
585   Output Parameter:
586 . dim - The total problem dimension
587 
588   Level: beginner
589 
590 .seealso: PetscDSGetNumFields(), PetscDSCreate()
591 @*/
592 PetscErrorCode PetscDSGetTotalDimension(PetscDS prob, PetscInt *dim)
593 {
594   PetscErrorCode ierr;
595 
596   PetscFunctionBegin;
597   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
598   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
599   PetscValidPointer(dim, 2);
600   *dim = prob->totDim;
601   PetscFunctionReturn(0);
602 }
603 
604 #undef __FUNCT__
605 #define __FUNCT__ "PetscDSGetTotalComponents"
606 /*@
607   PetscDSGetTotalComponents - Returns the total number of components in this system
608 
609   Not collective
610 
611   Input Parameter:
612 . prob - The PetscDS object
613 
614   Output Parameter:
615 . dim - The total number of components
616 
617   Level: beginner
618 
619 .seealso: PetscDSGetNumFields(), PetscDSCreate()
620 @*/
621 PetscErrorCode PetscDSGetTotalComponents(PetscDS prob, PetscInt *Nc)
622 {
623   PetscErrorCode ierr;
624 
625   PetscFunctionBegin;
626   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
627   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
628   PetscValidPointer(Nc, 2);
629   *Nc = prob->totComp;
630   PetscFunctionReturn(0);
631 }
632 
633 #undef __FUNCT__
634 #define __FUNCT__ "PetscDSGetDiscretization"
635 /*@
636   PetscDSGetDiscretization - Returns the discretization object for the given field
637 
638   Not collective
639 
640   Input Parameters:
641 + prob - The PetscDS object
642 - f - The field number
643 
644   Output Parameter:
645 . disc - The discretization object
646 
647   Level: beginner
648 
649 .seealso: PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
650 @*/
651 PetscErrorCode PetscDSGetDiscretization(PetscDS prob, PetscInt f, PetscObject *disc)
652 {
653   PetscFunctionBegin;
654   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
655   PetscValidPointer(disc, 3);
656   if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf);
657   *disc = prob->disc[f];
658   PetscFunctionReturn(0);
659 }
660 
661 #undef __FUNCT__
662 #define __FUNCT__ "PetscDSSetDiscretization"
663 /*@
664   PetscDSSetDiscretization - Sets the discretization object for the given field
665 
666   Not collective
667 
668   Input Parameters:
669 + prob - The PetscDS object
670 . f - The field number
671 - disc - The discretization object
672 
673   Level: beginner
674 
675 .seealso: PetscDSGetDiscretization(), PetscDSAddDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
676 @*/
677 PetscErrorCode PetscDSSetDiscretization(PetscDS prob, PetscInt f, PetscObject disc)
678 {
679   PetscErrorCode ierr;
680 
681   PetscFunctionBegin;
682   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
683   PetscValidPointer(disc, 3);
684   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
685   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
686   if (prob->disc[f]) {ierr = PetscObjectDereference(prob->disc[f]);CHKERRQ(ierr);}
687   prob->disc[f] = disc;
688   ierr = PetscObjectReference(disc);CHKERRQ(ierr);
689   {
690     PetscClassId id;
691 
692     ierr = PetscObjectGetClassId(disc, &id);CHKERRQ(ierr);
693     if (id == PETSCFV_CLASSID) {
694       prob->implicit[f]      = PETSC_FALSE;
695       prob->adjacency[f*2+0] = PETSC_TRUE;
696       prob->adjacency[f*2+1] = PETSC_FALSE;
697     }
698   }
699   PetscFunctionReturn(0);
700 }
701 
702 #undef __FUNCT__
703 #define __FUNCT__ "PetscDSAddDiscretization"
704 /*@
705   PetscDSAddDiscretization - Adds a discretization object
706 
707   Not collective
708 
709   Input Parameters:
710 + prob - The PetscDS object
711 - disc - The boundary discretization object
712 
713   Level: beginner
714 
715 .seealso: PetscDSGetDiscretization(), PetscDSSetDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
716 @*/
717 PetscErrorCode PetscDSAddDiscretization(PetscDS prob, PetscObject disc)
718 {
719   PetscErrorCode ierr;
720 
721   PetscFunctionBegin;
722   ierr = PetscDSSetDiscretization(prob, prob->Nf, disc);CHKERRQ(ierr);
723   PetscFunctionReturn(0);
724 }
725 
726 #undef __FUNCT__
727 #define __FUNCT__ "PetscDSGetImplicit"
728 /*@
729   PetscDSGetImplicit - Returns the flag for implicit solve for this field. This is just a guide for IMEX
730 
731   Not collective
732 
733   Input Parameters:
734 + prob - The PetscDS object
735 - f - The field number
736 
737   Output Parameter:
738 . implicit - The flag indicating what kind of solve to use for this field
739 
740   Level: developer
741 
742 .seealso: PetscDSSetImplicit(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
743 @*/
744 PetscErrorCode PetscDSGetImplicit(PetscDS prob, PetscInt f, PetscBool *implicit)
745 {
746   PetscFunctionBegin;
747   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
748   PetscValidPointer(implicit, 3);
749   if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf);
750   *implicit = prob->implicit[f];
751   PetscFunctionReturn(0);
752 }
753 
754 #undef __FUNCT__
755 #define __FUNCT__ "PetscDSSetImplicit"
756 /*@
757   PetscDSSetImplicit - Set the flag for implicit solve for this field. This is just a guide for IMEX
758 
759   Not collective
760 
761   Input Parameters:
762 + prob - The PetscDS object
763 . f - The field number
764 - implicit - The flag indicating what kind of solve to use for this field
765 
766   Level: developer
767 
768 .seealso: PetscDSGetImplicit(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
769 @*/
770 PetscErrorCode PetscDSSetImplicit(PetscDS prob, PetscInt f, PetscBool implicit)
771 {
772   PetscFunctionBegin;
773   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
774   if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf);
775   prob->implicit[f] = implicit;
776   PetscFunctionReturn(0);
777 }
778 
779 #undef __FUNCT__
780 #define __FUNCT__ "PetscDSGetAdjacency"
781 /*@
782   PetscDSGetAdjacency - Returns the flags for determining variable influence
783 
784   Not collective
785 
786   Input Parameters:
787 + prob - The PetscDS object
788 - f - The field number
789 
790   Output Parameter:
791 + useCone    - Flag for variable influence starting with the cone operation
792 - useClosure - Flag for variable influence using transitive closure
793 
794   Note: See the discussion in DMPlexGetAdjacencyUseCone() and DMPlexGetAdjacencyUseClosure()
795 
796   Level: developer
797 
798 .seealso: PetscDSSetAdjacency(), DMPlexGetAdjacencyUseCone(), DMPlexGetAdjacencyUseClosure(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
799 @*/
800 PetscErrorCode PetscDSGetAdjacency(PetscDS prob, PetscInt f, PetscBool *useCone, PetscBool *useClosure)
801 {
802   PetscFunctionBegin;
803   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
804   PetscValidPointer(useCone, 3);
805   PetscValidPointer(useClosure, 4);
806   if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf);
807   *useCone    = prob->adjacency[f*2+0];
808   *useClosure = prob->adjacency[f*2+1];
809   PetscFunctionReturn(0);
810 }
811 
812 #undef __FUNCT__
813 #define __FUNCT__ "PetscDSSetAdjacency"
814 /*@
815   PetscDSSetAdjacency - Set the flags for determining variable influence
816 
817   Not collective
818 
819   Input Parameters:
820 + prob - The PetscDS object
821 . f - The field number
822 . useCone    - Flag for variable influence starting with the cone operation
823 - useClosure - Flag for variable influence using transitive closure
824 
825   Note: See the discussion in DMPlexGetAdjacencyUseCone() and DMPlexGetAdjacencyUseClosure()
826 
827   Level: developer
828 
829 .seealso: PetscDSGetAdjacency(), DMPlexGetAdjacencyUseCone(), DMPlexGetAdjacencyUseClosure(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
830 @*/
831 PetscErrorCode PetscDSSetAdjacency(PetscDS prob, PetscInt f, PetscBool useCone, PetscBool useClosure)
832 {
833   PetscFunctionBegin;
834   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
835   if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf);
836   prob->adjacency[f*2+0] = useCone;
837   prob->adjacency[f*2+1] = useClosure;
838   PetscFunctionReturn(0);
839 }
840 
841 #undef __FUNCT__
842 #define __FUNCT__ "PetscDSGetObjective"
843 PetscErrorCode PetscDSGetObjective(PetscDS prob, PetscInt f,
844                                    void (**obj)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
845                                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
846                                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
847                                                 PetscReal t, const PetscReal x[], PetscScalar obj[]))
848 {
849   PetscFunctionBegin;
850   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
851   PetscValidPointer(obj, 2);
852   if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf);
853   *obj = prob->obj[f];
854   PetscFunctionReturn(0);
855 }
856 
857 #undef __FUNCT__
858 #define __FUNCT__ "PetscDSSetObjective"
859 PetscErrorCode PetscDSSetObjective(PetscDS prob, PetscInt f,
860                                    void (*obj)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
861                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
862                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
863                                                PetscReal t, const PetscReal x[], PetscScalar obj[]))
864 {
865   PetscErrorCode ierr;
866 
867   PetscFunctionBegin;
868   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
869   if (obj) PetscValidFunction(obj, 2);
870   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
871   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
872   prob->obj[f] = obj;
873   PetscFunctionReturn(0);
874 }
875 
876 #undef __FUNCT__
877 #define __FUNCT__ "PetscDSGetResidual"
878 /*@C
879   PetscDSGetResidual - Get the pointwise residual function for a given test field
880 
881   Not collective
882 
883   Input Parameters:
884 + prob - The PetscDS
885 - f    - The test field number
886 
887   Output Parameters:
888 + f0 - integrand for the test function term
889 - f1 - integrand for the test function gradient term
890 
891   Note: We are using a first order FEM model for the weak form:
892 
893   \int_\Omega \phi f_0(u, u_t, \nabla u, x, t) + \nabla\phi \cdot {\vec f}_1(u, u_t, \nabla u, x, t)
894 
895 The calling sequence for the callbacks f0 and f1 is given by:
896 
897 $ f0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
898 $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
899 $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
900 $    PetscReal t, const PetscReal x[], PetscScalar f0[])
901 
902 + dim - the spatial dimension
903 . Nf - the number of fields
904 . uOff - the offset into u[] and u_t[] for each field
905 . uOff_x - the offset into u_x[] for each field
906 . u - each field evaluated at the current point
907 . u_t - the time derivative of each field evaluated at the current point
908 . u_x - the gradient of each field evaluated at the current point
909 . aOff - the offset into a[] and a_t[] for each auxiliary field
910 . aOff_x - the offset into a_x[] for each auxiliary field
911 . a - each auxiliary field evaluated at the current point
912 . a_t - the time derivative of each auxiliary field evaluated at the current point
913 . a_x - the gradient of auxiliary each field evaluated at the current point
914 . t - current time
915 . x - coordinates of the current point
916 - f0 - output values at the current point
917 
918   Level: intermediate
919 
920 .seealso: PetscDSSetResidual()
921 @*/
922 PetscErrorCode PetscDSGetResidual(PetscDS prob, PetscInt f,
923                                   void (**f0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
924                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
925                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
926                                               PetscReal t, const PetscReal x[], PetscScalar f0[]),
927                                   void (**f1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
928                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
929                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
930                                               PetscReal t, const PetscReal x[], PetscScalar f1[]))
931 {
932   PetscFunctionBegin;
933   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
934   if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf);
935   if (f0) {PetscValidPointer(f0, 3); *f0 = prob->f[f*2+0];}
936   if (f1) {PetscValidPointer(f1, 4); *f1 = prob->f[f*2+1];}
937   PetscFunctionReturn(0);
938 }
939 
940 #undef __FUNCT__
941 #define __FUNCT__ "PetscDSSetResidual"
942 /*@C
943   PetscDSSetResidual - Set the pointwise residual function for a given test field
944 
945   Not collective
946 
947   Input Parameters:
948 + prob - The PetscDS
949 . f    - The test field number
950 . f0 - integrand for the test function term
951 - f1 - integrand for the test function gradient term
952 
953   Note: We are using a first order FEM model for the weak form:
954 
955   \int_\Omega \phi f_0(u, u_t, \nabla u, x, t) + \nabla\phi \cdot {\vec f}_1(u, u_t, \nabla u, x, t)
956 
957 The calling sequence for the callbacks f0 and f1 is given by:
958 
959 $ f0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
960 $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
961 $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
962 $    PetscReal t, const PetscReal x[], PetscScalar f0[])
963 
964 + dim - the spatial dimension
965 . Nf - the number of fields
966 . uOff - the offset into u[] and u_t[] for each field
967 . uOff_x - the offset into u_x[] for each field
968 . u - each field evaluated at the current point
969 . u_t - the time derivative of each field evaluated at the current point
970 . u_x - the gradient of each field evaluated at the current point
971 . aOff - the offset into a[] and a_t[] for each auxiliary field
972 . aOff_x - the offset into a_x[] for each auxiliary field
973 . a - each auxiliary field evaluated at the current point
974 . a_t - the time derivative of each auxiliary field evaluated at the current point
975 . a_x - the gradient of auxiliary each field evaluated at the current point
976 . t - current time
977 . x - coordinates of the current point
978 - f0 - output values at the current point
979 
980   Level: intermediate
981 
982 .seealso: PetscDSGetResidual()
983 @*/
984 PetscErrorCode PetscDSSetResidual(PetscDS prob, PetscInt f,
985                                   void (*f0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
986                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
987                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
988                                              PetscReal t, const PetscReal x[], PetscScalar f0[]),
989                                   void (*f1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
990                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
991                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
992                                              PetscReal t, const PetscReal x[], PetscScalar f1[]))
993 {
994   PetscErrorCode ierr;
995 
996   PetscFunctionBegin;
997   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
998   if (f0) PetscValidFunction(f0, 3);
999   if (f1) PetscValidFunction(f1, 4);
1000   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
1001   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
1002   prob->f[f*2+0] = f0;
1003   prob->f[f*2+1] = f1;
1004   PetscFunctionReturn(0);
1005 }
1006 
1007 #undef __FUNCT__
1008 #define __FUNCT__ "PetscDSHasJacobian"
1009 /*@C
1010   PetscDSHasJacobian - Signals that Jacobian functions have been set
1011 
1012   Not collective
1013 
1014   Input Parameter:
1015 . prob - The PetscDS
1016 
1017   Output Parameter:
1018 . hasJac - flag that pointwise function for the Jacobian has been set
1019 
1020   Level: intermediate
1021 
1022 .seealso: PetscDSGetJacobianPreconditioner(), PetscDSSetJacobianPreconditioner(), PetscDSGetJacobian()
1023 @*/
1024 PetscErrorCode PetscDSHasJacobian(PetscDS prob, PetscBool *hasJac)
1025 {
1026   PetscInt f, g, h;
1027 
1028   PetscFunctionBegin;
1029   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1030   *hasJac = PETSC_FALSE;
1031   for (f = 0; f < prob->Nf; ++f) {
1032     for (g = 0; g < prob->Nf; ++g) {
1033       for (h = 0; h < 4; ++h) {
1034         if (prob->g[(f*prob->Nf + g)*4+h]) *hasJac = PETSC_TRUE;
1035       }
1036     }
1037   }
1038   PetscFunctionReturn(0);
1039 }
1040 
1041 #undef __FUNCT__
1042 #define __FUNCT__ "PetscDSGetJacobian"
1043 /*@C
1044   PetscDSGetJacobian - Get the pointwise Jacobian function for given test and basis field
1045 
1046   Not collective
1047 
1048   Input Parameters:
1049 + prob - The PetscDS
1050 . f    - The test field number
1051 - g    - The field number
1052 
1053   Output Parameters:
1054 + g0 - integrand for the test and basis function term
1055 . g1 - integrand for the test function and basis function gradient term
1056 . g2 - integrand for the test function gradient and basis function term
1057 - g3 - integrand for the test function gradient and basis function gradient term
1058 
1059   Note: We are using a first order FEM model for the weak form:
1060 
1061   \int_\Omega \phi g_0(u, u_t, \nabla u, x, t) \psi + \phi {\vec g}_1(u, u_t, \nabla u, x, t) \nabla \psi + \nabla\phi \cdot {\vec g}_2(u, u_t, \nabla u, x, t) \psi + \nabla\phi \cdot {\overleftrightarrow g}_3(u, u_t, \nabla u, x, t) \cdot \nabla \psi
1062 
1063 The calling sequence for the callbacks g0, g1, g2 and g3 is given by:
1064 
1065 $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1066 $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1067 $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1068 $    PetscReal t, const PetscReal u_tShift, const PetscReal x[], PetscScalar g0[])
1069 
1070 + dim - the spatial dimension
1071 . Nf - the number of fields
1072 . uOff - the offset into u[] and u_t[] for each field
1073 . uOff_x - the offset into u_x[] for each field
1074 . u - each field evaluated at the current point
1075 . u_t - the time derivative of each field evaluated at the current point
1076 . u_x - the gradient of each field evaluated at the current point
1077 . aOff - the offset into a[] and a_t[] for each auxiliary field
1078 . aOff_x - the offset into a_x[] for each auxiliary field
1079 . a - each auxiliary field evaluated at the current point
1080 . a_t - the time derivative of each auxiliary field evaluated at the current point
1081 . a_x - the gradient of auxiliary each field evaluated at the current point
1082 . t - current time
1083 . u_tShift - the multiplier a for dF/dU_t
1084 . x - coordinates of the current point
1085 - g0 - output values at the current point
1086 
1087   Level: intermediate
1088 
1089 .seealso: PetscDSSetJacobian()
1090 @*/
1091 PetscErrorCode PetscDSGetJacobian(PetscDS prob, PetscInt f, PetscInt g,
1092                                   void (**g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1093                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1094                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1095                                               PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g0[]),
1096                                   void (**g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1097                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1098                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1099                                               PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g1[]),
1100                                   void (**g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1101                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1102                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1103                                               PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g2[]),
1104                                   void (**g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1105                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1106                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1107                                               PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g3[]))
1108 {
1109   PetscFunctionBegin;
1110   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1111   if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf);
1112   if ((g < 0) || (g >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", g, prob->Nf);
1113   if (g0) {PetscValidPointer(g0, 4); *g0 = prob->g[(f*prob->Nf + g)*4+0];}
1114   if (g1) {PetscValidPointer(g1, 5); *g1 = prob->g[(f*prob->Nf + g)*4+1];}
1115   if (g2) {PetscValidPointer(g2, 6); *g2 = prob->g[(f*prob->Nf + g)*4+2];}
1116   if (g3) {PetscValidPointer(g3, 7); *g3 = prob->g[(f*prob->Nf + g)*4+3];}
1117   PetscFunctionReturn(0);
1118 }
1119 
1120 #undef __FUNCT__
1121 #define __FUNCT__ "PetscDSSetJacobian"
1122 /*@C
1123   PetscDSSetJacobian - Set the pointwise Jacobian function for given test and basis fields
1124 
1125   Not collective
1126 
1127   Input Parameters:
1128 + prob - The PetscDS
1129 . f    - The test field number
1130 . g    - The field number
1131 . g0 - integrand for the test and basis function term
1132 . g1 - integrand for the test function and basis function gradient term
1133 . g2 - integrand for the test function gradient and basis function term
1134 - g3 - integrand for the test function gradient and basis function gradient term
1135 
1136   Note: We are using a first order FEM model for the weak form:
1137 
1138   \int_\Omega \phi g_0(u, u_t, \nabla u, x, t) \psi + \phi {\vec g}_1(u, u_t, \nabla u, x, t) \nabla \psi + \nabla\phi \cdot {\vec g}_2(u, u_t, \nabla u, x, t) \psi + \nabla\phi \cdot {\overleftrightarrow g}_3(u, u_t, \nabla u, x, t) \cdot \nabla \psi
1139 
1140 The calling sequence for the callbacks g0, g1, g2 and g3 is given by:
1141 
1142 $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1143 $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1144 $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1145 $    PetscReal t, const PetscReal x[], PetscScalar g0[])
1146 
1147 + dim - the spatial dimension
1148 . Nf - the number of fields
1149 . uOff - the offset into u[] and u_t[] for each field
1150 . uOff_x - the offset into u_x[] for each field
1151 . u - each field evaluated at the current point
1152 . u_t - the time derivative of each field evaluated at the current point
1153 . u_x - the gradient of each field evaluated at the current point
1154 . aOff - the offset into a[] and a_t[] for each auxiliary field
1155 . aOff_x - the offset into a_x[] for each auxiliary field
1156 . a - each auxiliary field evaluated at the current point
1157 . a_t - the time derivative of each auxiliary field evaluated at the current point
1158 . a_x - the gradient of auxiliary each field evaluated at the current point
1159 . t - current time
1160 . u_tShift - the multiplier a for dF/dU_t
1161 . x - coordinates of the current point
1162 - g0 - output values at the current point
1163 
1164   Level: intermediate
1165 
1166 .seealso: PetscDSGetJacobian()
1167 @*/
1168 PetscErrorCode PetscDSSetJacobian(PetscDS prob, PetscInt f, PetscInt g,
1169                                   void (*g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1170                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1171                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1172                                              PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g0[]),
1173                                   void (*g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1174                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1175                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1176                                              PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g1[]),
1177                                   void (*g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1178                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1179                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1180                                              PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g2[]),
1181                                   void (*g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1182                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1183                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1184                                              PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g3[]))
1185 {
1186   PetscErrorCode ierr;
1187 
1188   PetscFunctionBegin;
1189   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1190   if (g0) PetscValidFunction(g0, 4);
1191   if (g1) PetscValidFunction(g1, 5);
1192   if (g2) PetscValidFunction(g2, 6);
1193   if (g3) PetscValidFunction(g3, 7);
1194   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
1195   if (g < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", g);
1196   ierr = PetscDSEnlarge_Static(prob, PetscMax(f, g)+1);CHKERRQ(ierr);
1197   prob->g[(f*prob->Nf + g)*4+0] = g0;
1198   prob->g[(f*prob->Nf + g)*4+1] = g1;
1199   prob->g[(f*prob->Nf + g)*4+2] = g2;
1200   prob->g[(f*prob->Nf + g)*4+3] = g3;
1201   PetscFunctionReturn(0);
1202 }
1203 
1204 #undef __FUNCT__
1205 #define __FUNCT__ "PetscDSHasJacobianPreconditioner"
1206 /*@C
1207   PetscDSHasJacobianPreconditioner - Signals that a Jacobian preconditioner matrix has been set
1208 
1209   Not collective
1210 
1211   Input Parameter:
1212 . prob - The PetscDS
1213 
1214   Output Parameter:
1215 . hasJacPre - flag that pointwise function for Jacobian preconditioner matrix has been set
1216 
1217   Level: intermediate
1218 
1219 .seealso: PetscDSGetJacobianPreconditioner(), PetscDSSetJacobianPreconditioner(), PetscDSGetJacobian()
1220 @*/
1221 PetscErrorCode PetscDSHasJacobianPreconditioner(PetscDS prob, PetscBool *hasJacPre)
1222 {
1223   PetscInt f, g, h;
1224 
1225   PetscFunctionBegin;
1226   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1227   *hasJacPre = PETSC_FALSE;
1228   for (f = 0; f < prob->Nf; ++f) {
1229     for (g = 0; g < prob->Nf; ++g) {
1230       for (h = 0; h < 4; ++h) {
1231         if (prob->gp[(f*prob->Nf + g)*4+h]) *hasJacPre = PETSC_TRUE;
1232       }
1233     }
1234   }
1235   PetscFunctionReturn(0);
1236 }
1237 
1238 #undef __FUNCT__
1239 #define __FUNCT__ "PetscDSGetJacobianPreconditioner"
1240 /*@C
1241   PetscDSGetJacobianPreconditioner - Get the pointwise Jacobian preconditioner function for given test and basis field. If this is missing, the system matrix is used to build the preconditioner.
1242 
1243   Not collective
1244 
1245   Input Parameters:
1246 + prob - The PetscDS
1247 . f    - The test field number
1248 - g    - The field number
1249 
1250   Output Parameters:
1251 + g0 - integrand for the test and basis function term
1252 . g1 - integrand for the test function and basis function gradient term
1253 . g2 - integrand for the test function gradient and basis function term
1254 - g3 - integrand for the test function gradient and basis function gradient term
1255 
1256   Note: We are using a first order FEM model for the weak form:
1257 
1258   \int_\Omega \phi g_0(u, u_t, \nabla u, x, t) \psi + \phi {\vec g}_1(u, u_t, \nabla u, x, t) \nabla \psi + \nabla\phi \cdot {\vec g}_2(u, u_t, \nabla u, x, t) \psi + \nabla\phi \cdot {\overleftrightarrow g}_3(u, u_t, \nabla u, x, t) \cdot \nabla \psi
1259 
1260 The calling sequence for the callbacks g0, g1, g2 and g3 is given by:
1261 
1262 $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1263 $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1264 $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1265 $    PetscReal t, const PetscReal u_tShift, const PetscReal x[], PetscScalar g0[])
1266 
1267 + dim - the spatial dimension
1268 . Nf - the number of fields
1269 . uOff - the offset into u[] and u_t[] for each field
1270 . uOff_x - the offset into u_x[] for each field
1271 . u - each field evaluated at the current point
1272 . u_t - the time derivative of each field evaluated at the current point
1273 . u_x - the gradient of each field evaluated at the current point
1274 . aOff - the offset into a[] and a_t[] for each auxiliary field
1275 . aOff_x - the offset into a_x[] for each auxiliary field
1276 . a - each auxiliary field evaluated at the current point
1277 . a_t - the time derivative of each auxiliary field evaluated at the current point
1278 . a_x - the gradient of auxiliary each field evaluated at the current point
1279 . t - current time
1280 . u_tShift - the multiplier a for dF/dU_t
1281 . x - coordinates of the current point
1282 - g0 - output values at the current point
1283 
1284   Level: intermediate
1285 
1286 .seealso: PetscDSSetJacobianPreconditioner(), PetscDSGetJacobian()
1287 @*/
1288 PetscErrorCode PetscDSGetJacobianPreconditioner(PetscDS prob, PetscInt f, PetscInt g,
1289                                   void (**g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1290                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1291                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1292                                               PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g0[]),
1293                                   void (**g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1294                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1295                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1296                                               PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g1[]),
1297                                   void (**g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1298                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1299                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1300                                               PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g2[]),
1301                                   void (**g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1302                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1303                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1304                                               PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g3[]))
1305 {
1306   PetscFunctionBegin;
1307   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1308   if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf);
1309   if ((g < 0) || (g >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", g, prob->Nf);
1310   if (g0) {PetscValidPointer(g0, 4); *g0 = prob->gp[(f*prob->Nf + g)*4+0];}
1311   if (g1) {PetscValidPointer(g1, 5); *g1 = prob->gp[(f*prob->Nf + g)*4+1];}
1312   if (g2) {PetscValidPointer(g2, 6); *g2 = prob->gp[(f*prob->Nf + g)*4+2];}
1313   if (g3) {PetscValidPointer(g3, 7); *g3 = prob->gp[(f*prob->Nf + g)*4+3];}
1314   PetscFunctionReturn(0);
1315 }
1316 
1317 #undef __FUNCT__
1318 #define __FUNCT__ "PetscDSSetJacobianPreconditioner"
1319 /*@C
1320   PetscDSSetJacobianPreconditioner - Set the pointwise Jacobian preconditioner function for given test and basis fields. If this is missing, the system matrix is used to build the preconditioner.
1321 
1322   Not collective
1323 
1324   Input Parameters:
1325 + prob - The PetscDS
1326 . f    - The test field number
1327 . g    - The field number
1328 . g0 - integrand for the test and basis function term
1329 . g1 - integrand for the test function and basis function gradient term
1330 . g2 - integrand for the test function gradient and basis function term
1331 - g3 - integrand for the test function gradient and basis function gradient term
1332 
1333   Note: We are using a first order FEM model for the weak form:
1334 
1335   \int_\Omega \phi g_0(u, u_t, \nabla u, x, t) \psi + \phi {\vec g}_1(u, u_t, \nabla u, x, t) \nabla \psi + \nabla\phi \cdot {\vec g}_2(u, u_t, \nabla u, x, t) \psi + \nabla\phi \cdot {\overleftrightarrow g}_3(u, u_t, \nabla u, x, t) \cdot \nabla \psi
1336 
1337 The calling sequence for the callbacks g0, g1, g2 and g3 is given by:
1338 
1339 $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1340 $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1341 $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1342 $    PetscReal t, const PetscReal x[], PetscScalar g0[])
1343 
1344 + dim - the spatial dimension
1345 . Nf - the number of fields
1346 . uOff - the offset into u[] and u_t[] for each field
1347 . uOff_x - the offset into u_x[] for each field
1348 . u - each field evaluated at the current point
1349 . u_t - the time derivative of each field evaluated at the current point
1350 . u_x - the gradient of each field evaluated at the current point
1351 . aOff - the offset into a[] and a_t[] for each auxiliary field
1352 . aOff_x - the offset into a_x[] for each auxiliary field
1353 . a - each auxiliary field evaluated at the current point
1354 . a_t - the time derivative of each auxiliary field evaluated at the current point
1355 . a_x - the gradient of auxiliary each field evaluated at the current point
1356 . t - current time
1357 . u_tShift - the multiplier a for dF/dU_t
1358 . x - coordinates of the current point
1359 - g0 - output values at the current point
1360 
1361   Level: intermediate
1362 
1363 .seealso: PetscDSGetJacobianPreconditioner(), PetscDSSetJacobian()
1364 @*/
1365 PetscErrorCode PetscDSSetJacobianPreconditioner(PetscDS prob, PetscInt f, PetscInt g,
1366                                   void (*g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1367                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1368                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1369                                              PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g0[]),
1370                                   void (*g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1371                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1372                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1373                                              PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g1[]),
1374                                   void (*g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1375                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1376                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1377                                              PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g2[]),
1378                                   void (*g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1379                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1380                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1381                                              PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g3[]))
1382 {
1383   PetscErrorCode ierr;
1384 
1385   PetscFunctionBegin;
1386   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1387   if (g0) PetscValidFunction(g0, 4);
1388   if (g1) PetscValidFunction(g1, 5);
1389   if (g2) PetscValidFunction(g2, 6);
1390   if (g3) PetscValidFunction(g3, 7);
1391   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
1392   if (g < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", g);
1393   ierr = PetscDSEnlarge_Static(prob, PetscMax(f, g)+1);CHKERRQ(ierr);
1394   prob->gp[(f*prob->Nf + g)*4+0] = g0;
1395   prob->gp[(f*prob->Nf + g)*4+1] = g1;
1396   prob->gp[(f*prob->Nf + g)*4+2] = g2;
1397   prob->gp[(f*prob->Nf + g)*4+3] = g3;
1398   PetscFunctionReturn(0);
1399 }
1400 
1401 #undef __FUNCT__
1402 #define __FUNCT__ "PetscDSHasDynamicJacobian"
1403 /*@C
1404   PetscDSHasDynamicJacobian - Signals that a dynamic Jacobian, dF/du_t, has been set
1405 
1406   Not collective
1407 
1408   Input Parameter:
1409 . prob - The PetscDS
1410 
1411   Output Parameter:
1412 . hasDynJac - flag that pointwise function for dynamic Jacobian has been set
1413 
1414   Level: intermediate
1415 
1416 .seealso: PetscDSGetDynamicJacobian(), PetscDSSetDynamicJacobian(), PetscDSGetJacobian()
1417 @*/
1418 PetscErrorCode PetscDSHasDynamicJacobian(PetscDS prob, PetscBool *hasDynJac)
1419 {
1420   PetscInt f, g, h;
1421 
1422   PetscFunctionBegin;
1423   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1424   *hasDynJac = PETSC_FALSE;
1425   for (f = 0; f < prob->Nf; ++f) {
1426     for (g = 0; g < prob->Nf; ++g) {
1427       for (h = 0; h < 4; ++h) {
1428         if (prob->gt[(f*prob->Nf + g)*4+h]) *hasDynJac = PETSC_TRUE;
1429       }
1430     }
1431   }
1432   PetscFunctionReturn(0);
1433 }
1434 
1435 #undef __FUNCT__
1436 #define __FUNCT__ "PetscDSGetDynamicJacobian"
1437 /*@C
1438   PetscDSGetDynamicJacobian - Get the pointwise dynamic Jacobian, dF/du_t, function for given test and basis field
1439 
1440   Not collective
1441 
1442   Input Parameters:
1443 + prob - The PetscDS
1444 . f    - The test field number
1445 - g    - The field number
1446 
1447   Output Parameters:
1448 + g0 - integrand for the test and basis function term
1449 . g1 - integrand for the test function and basis function gradient term
1450 . g2 - integrand for the test function gradient and basis function term
1451 - g3 - integrand for the test function gradient and basis function gradient term
1452 
1453   Note: We are using a first order FEM model for the weak form:
1454 
1455   \int_\Omega \phi g_0(u, u_t, \nabla u, x, t) \psi + \phi {\vec g}_1(u, u_t, \nabla u, x, t) \nabla \psi + \nabla\phi \cdot {\vec g}_2(u, u_t, \nabla u, x, t) \psi + \nabla\phi \cdot {\overleftrightarrow g}_3(u, u_t, \nabla u, x, t) \cdot \nabla \psi
1456 
1457 The calling sequence for the callbacks g0, g1, g2 and g3 is given by:
1458 
1459 $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1460 $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1461 $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1462 $    PetscReal t, const PetscReal u_tShift, const PetscReal x[], PetscScalar g0[])
1463 
1464 + dim - the spatial dimension
1465 . Nf - the number of fields
1466 . uOff - the offset into u[] and u_t[] for each field
1467 . uOff_x - the offset into u_x[] for each field
1468 . u - each field evaluated at the current point
1469 . u_t - the time derivative of each field evaluated at the current point
1470 . u_x - the gradient of each field evaluated at the current point
1471 . aOff - the offset into a[] and a_t[] for each auxiliary field
1472 . aOff_x - the offset into a_x[] for each auxiliary field
1473 . a - each auxiliary field evaluated at the current point
1474 . a_t - the time derivative of each auxiliary field evaluated at the current point
1475 . a_x - the gradient of auxiliary each field evaluated at the current point
1476 . t - current time
1477 . u_tShift - the multiplier a for dF/dU_t
1478 . x - coordinates of the current point
1479 - g0 - output values at the current point
1480 
1481   Level: intermediate
1482 
1483 .seealso: PetscDSSetJacobian()
1484 @*/
1485 PetscErrorCode PetscDSGetDynamicJacobian(PetscDS prob, PetscInt f, PetscInt g,
1486                                          void (**g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1487                                                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1488                                                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1489                                                      PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g0[]),
1490                                          void (**g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1491                                                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1492                                                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1493                                                      PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g1[]),
1494                                          void (**g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1495                                                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1496                                                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1497                                               PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g2[]),
1498                                          void (**g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1499                                                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1500                                                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1501                                                      PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g3[]))
1502 {
1503   PetscFunctionBegin;
1504   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1505   if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf);
1506   if ((g < 0) || (g >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", g, prob->Nf);
1507   if (g0) {PetscValidPointer(g0, 4); *g0 = prob->gt[(f*prob->Nf + g)*4+0];}
1508   if (g1) {PetscValidPointer(g1, 5); *g1 = prob->gt[(f*prob->Nf + g)*4+1];}
1509   if (g2) {PetscValidPointer(g2, 6); *g2 = prob->gt[(f*prob->Nf + g)*4+2];}
1510   if (g3) {PetscValidPointer(g3, 7); *g3 = prob->gt[(f*prob->Nf + g)*4+3];}
1511   PetscFunctionReturn(0);
1512 }
1513 
1514 #undef __FUNCT__
1515 #define __FUNCT__ "PetscDSSetDynamicJacobian"
1516 /*@C
1517   PetscDSSetDynamicJacobian - Set the pointwise dynamic Jacobian, dF/du_t, function for given test and basis fields
1518 
1519   Not collective
1520 
1521   Input Parameters:
1522 + prob - The PetscDS
1523 . f    - The test field number
1524 . g    - The field number
1525 . g0 - integrand for the test and basis function term
1526 . g1 - integrand for the test function and basis function gradient term
1527 . g2 - integrand for the test function gradient and basis function term
1528 - g3 - integrand for the test function gradient and basis function gradient term
1529 
1530   Note: We are using a first order FEM model for the weak form:
1531 
1532   \int_\Omega \phi g_0(u, u_t, \nabla u, x, t) \psi + \phi {\vec g}_1(u, u_t, \nabla u, x, t) \nabla \psi + \nabla\phi \cdot {\vec g}_2(u, u_t, \nabla u, x, t) \psi + \nabla\phi \cdot {\overleftrightarrow g}_3(u, u_t, \nabla u, x, t) \cdot \nabla \psi
1533 
1534 The calling sequence for the callbacks g0, g1, g2 and g3 is given by:
1535 
1536 $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1537 $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1538 $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1539 $    PetscReal t, const PetscReal x[], PetscScalar g0[])
1540 
1541 + dim - the spatial dimension
1542 . Nf - the number of fields
1543 . uOff - the offset into u[] and u_t[] for each field
1544 . uOff_x - the offset into u_x[] for each field
1545 . u - each field evaluated at the current point
1546 . u_t - the time derivative of each field evaluated at the current point
1547 . u_x - the gradient of each field evaluated at the current point
1548 . aOff - the offset into a[] and a_t[] for each auxiliary field
1549 . aOff_x - the offset into a_x[] for each auxiliary field
1550 . a - each auxiliary field evaluated at the current point
1551 . a_t - the time derivative of each auxiliary field evaluated at the current point
1552 . a_x - the gradient of auxiliary each field evaluated at the current point
1553 . t - current time
1554 . u_tShift - the multiplier a for dF/dU_t
1555 . x - coordinates of the current point
1556 - g0 - output values at the current point
1557 
1558   Level: intermediate
1559 
1560 .seealso: PetscDSGetJacobian()
1561 @*/
1562 PetscErrorCode PetscDSSetDynamicJacobian(PetscDS prob, PetscInt f, PetscInt g,
1563                                          void (*g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1564                                                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1565                                                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1566                                                     PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g0[]),
1567                                          void (*g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1568                                                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1569                                                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1570                                                     PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g1[]),
1571                                          void (*g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1572                                                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1573                                                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1574                                                     PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g2[]),
1575                                          void (*g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1576                                                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1577                                                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1578                                                     PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g3[]))
1579 {
1580   PetscErrorCode ierr;
1581 
1582   PetscFunctionBegin;
1583   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1584   if (g0) PetscValidFunction(g0, 4);
1585   if (g1) PetscValidFunction(g1, 5);
1586   if (g2) PetscValidFunction(g2, 6);
1587   if (g3) PetscValidFunction(g3, 7);
1588   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
1589   if (g < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", g);
1590   ierr = PetscDSEnlarge_Static(prob, PetscMax(f, g)+1);CHKERRQ(ierr);
1591   prob->gt[(f*prob->Nf + g)*4+0] = g0;
1592   prob->gt[(f*prob->Nf + g)*4+1] = g1;
1593   prob->gt[(f*prob->Nf + g)*4+2] = g2;
1594   prob->gt[(f*prob->Nf + g)*4+3] = g3;
1595   PetscFunctionReturn(0);
1596 }
1597 
1598 #undef __FUNCT__
1599 #define __FUNCT__ "PetscDSGetRiemannSolver"
1600 /*@C
1601   PetscDSGetRiemannSolver - Returns the Riemann solver for the given field
1602 
1603   Not collective
1604 
1605   Input Arguments:
1606 + prob - The PetscDS object
1607 - f    - The field number
1608 
1609   Output Argument:
1610 . r    - Riemann solver
1611 
1612   Calling sequence for r:
1613 
1614 $ r(PetscInt dim, PetscInt Nf, const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx)
1615 
1616 + dim  - The spatial dimension
1617 . Nf   - The number of fields
1618 . x    - The coordinates at a point on the interface
1619 . n    - The normal vector to the interface
1620 . uL   - The state vector to the left of the interface
1621 . uR   - The state vector to the right of the interface
1622 . flux - output array of flux through the interface
1623 - ctx  - optional user context
1624 
1625   Level: intermediate
1626 
1627 .seealso: PetscDSSetRiemannSolver()
1628 @*/
1629 PetscErrorCode PetscDSGetRiemannSolver(PetscDS prob, PetscInt f,
1630                                        void (**r)(PetscInt dim, PetscInt Nf, const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx))
1631 {
1632   PetscFunctionBegin;
1633   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1634   if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf);
1635   PetscValidPointer(r, 3);
1636   *r = prob->r[f];
1637   PetscFunctionReturn(0);
1638 }
1639 
1640 #undef __FUNCT__
1641 #define __FUNCT__ "PetscDSSetRiemannSolver"
1642 /*@C
1643   PetscDSSetRiemannSolver - Sets the Riemann solver for the given field
1644 
1645   Not collective
1646 
1647   Input Arguments:
1648 + prob - The PetscDS object
1649 . f    - The field number
1650 - r    - Riemann solver
1651 
1652   Calling sequence for r:
1653 
1654 $ r(PetscInt dim, PetscInt Nf, const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx)
1655 
1656 + dim  - The spatial dimension
1657 . Nf   - The number of fields
1658 . x    - The coordinates at a point on the interface
1659 . n    - The normal vector to the interface
1660 . uL   - The state vector to the left of the interface
1661 . uR   - The state vector to the right of the interface
1662 . flux - output array of flux through the interface
1663 - ctx  - optional user context
1664 
1665   Level: intermediate
1666 
1667 .seealso: PetscDSGetRiemannSolver()
1668 @*/
1669 PetscErrorCode PetscDSSetRiemannSolver(PetscDS prob, PetscInt f,
1670                                        void (*r)(PetscInt dim, PetscInt Nf, const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx))
1671 {
1672   PetscErrorCode ierr;
1673 
1674   PetscFunctionBegin;
1675   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1676   if (r) PetscValidFunction(r, 3);
1677   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
1678   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
1679   prob->r[f] = r;
1680   PetscFunctionReturn(0);
1681 }
1682 
1683 #undef __FUNCT__
1684 #define __FUNCT__ "PetscDSGetContext"
1685 PetscErrorCode PetscDSGetContext(PetscDS prob, PetscInt f, void **ctx)
1686 {
1687   PetscFunctionBegin;
1688   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1689   if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf);
1690   PetscValidPointer(ctx, 3);
1691   *ctx = prob->ctx[f];
1692   PetscFunctionReturn(0);
1693 }
1694 
1695 #undef __FUNCT__
1696 #define __FUNCT__ "PetscDSSetContext"
1697 PetscErrorCode PetscDSSetContext(PetscDS prob, PetscInt f, void *ctx)
1698 {
1699   PetscErrorCode ierr;
1700 
1701   PetscFunctionBegin;
1702   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1703   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
1704   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
1705   prob->ctx[f] = ctx;
1706   PetscFunctionReturn(0);
1707 }
1708 
1709 #undef __FUNCT__
1710 #define __FUNCT__ "PetscDSGetBdResidual"
1711 /*@C
1712   PetscDSGetBdResidual - Get the pointwise boundary residual function for a given test field
1713 
1714   Not collective
1715 
1716   Input Parameters:
1717 + prob - The PetscDS
1718 - f    - The test field number
1719 
1720   Output Parameters:
1721 + f0 - boundary integrand for the test function term
1722 - f1 - boundary integrand for the test function gradient term
1723 
1724   Note: We are using a first order FEM model for the weak form:
1725 
1726   \int_\Gamma \phi {\vec f}_0(u, u_t, \nabla u, x, t) \cdot \hat n + \nabla\phi \cdot {\overleftrightarrow f}_1(u, u_t, \nabla u, x, t) \cdot \hat n
1727 
1728 The calling sequence for the callbacks f0 and f1 is given by:
1729 
1730 $ f0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1731 $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1732 $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1733 $    PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f0[])
1734 
1735 + dim - the spatial dimension
1736 . Nf - the number of fields
1737 . uOff - the offset into u[] and u_t[] for each field
1738 . uOff_x - the offset into u_x[] for each field
1739 . u - each field evaluated at the current point
1740 . u_t - the time derivative of each field evaluated at the current point
1741 . u_x - the gradient of each field evaluated at the current point
1742 . aOff - the offset into a[] and a_t[] for each auxiliary field
1743 . aOff_x - the offset into a_x[] for each auxiliary field
1744 . a - each auxiliary field evaluated at the current point
1745 . a_t - the time derivative of each auxiliary field evaluated at the current point
1746 . a_x - the gradient of auxiliary each field evaluated at the current point
1747 . t - current time
1748 . x - coordinates of the current point
1749 . n - unit normal at the current point
1750 - f0 - output values at the current point
1751 
1752   Level: intermediate
1753 
1754 .seealso: PetscDSSetBdResidual()
1755 @*/
1756 PetscErrorCode PetscDSGetBdResidual(PetscDS prob, PetscInt f,
1757                                     void (**f0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1758                                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1759                                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1760                                                 PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f0[]),
1761                                     void (**f1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1762                                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1763                                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1764                                                 PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f1[]))
1765 {
1766   PetscFunctionBegin;
1767   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1768   if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf);
1769   if (f0) {PetscValidPointer(f0, 3); *f0 = prob->fBd[f*2+0];}
1770   if (f1) {PetscValidPointer(f1, 4); *f1 = prob->fBd[f*2+1];}
1771   PetscFunctionReturn(0);
1772 }
1773 
1774 #undef __FUNCT__
1775 #define __FUNCT__ "PetscDSSetBdResidual"
1776 /*@C
1777   PetscDSSetBdResidual - Get the pointwise boundary residual function for a given test field
1778 
1779   Not collective
1780 
1781   Input Parameters:
1782 + prob - The PetscDS
1783 . f    - The test field number
1784 . f0 - boundary integrand for the test function term
1785 - f1 - boundary integrand for the test function gradient term
1786 
1787   Note: We are using a first order FEM model for the weak form:
1788 
1789   \int_\Gamma \phi {\vec f}_0(u, u_t, \nabla u, x, t) \cdot \hat n + \nabla\phi \cdot {\overleftrightarrow f}_1(u, u_t, \nabla u, x, t) \cdot \hat n
1790 
1791 The calling sequence for the callbacks f0 and f1 is given by:
1792 
1793 $ f0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1794 $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1795 $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1796 $    PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f0[])
1797 
1798 + dim - the spatial dimension
1799 . Nf - the number of fields
1800 . uOff - the offset into u[] and u_t[] for each field
1801 . uOff_x - the offset into u_x[] for each field
1802 . u - each field evaluated at the current point
1803 . u_t - the time derivative of each field evaluated at the current point
1804 . u_x - the gradient of each field evaluated at the current point
1805 . aOff - the offset into a[] and a_t[] for each auxiliary field
1806 . aOff_x - the offset into a_x[] for each auxiliary field
1807 . a - each auxiliary field evaluated at the current point
1808 . a_t - the time derivative of each auxiliary field evaluated at the current point
1809 . a_x - the gradient of auxiliary each field evaluated at the current point
1810 . t - current time
1811 . x - coordinates of the current point
1812 . n - unit normal at the current point
1813 - f0 - output values at the current point
1814 
1815   Level: intermediate
1816 
1817 .seealso: PetscDSGetBdResidual()
1818 @*/
1819 PetscErrorCode PetscDSSetBdResidual(PetscDS prob, PetscInt f,
1820                                     void (*f0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1821                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1822                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1823                                                PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f0[]),
1824                                     void (*f1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1825                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1826                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1827                                                PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f1[]))
1828 {
1829   PetscErrorCode ierr;
1830 
1831   PetscFunctionBegin;
1832   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1833   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
1834   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
1835   if (f0) {PetscValidFunction(f0, 3); prob->fBd[f*2+0] = f0;}
1836   if (f1) {PetscValidFunction(f1, 4); prob->fBd[f*2+1] = f1;}
1837   PetscFunctionReturn(0);
1838 }
1839 
1840 #undef __FUNCT__
1841 #define __FUNCT__ "PetscDSGetBdJacobian"
1842 /*@C
1843   PetscDSGetBdJacobian - Get the pointwise boundary Jacobian function for given test and basis field
1844 
1845   Not collective
1846 
1847   Input Parameters:
1848 + prob - The PetscDS
1849 . f    - The test field number
1850 - g    - The field number
1851 
1852   Output Parameters:
1853 + g0 - integrand for the test and basis function term
1854 . g1 - integrand for the test function and basis function gradient term
1855 . g2 - integrand for the test function gradient and basis function term
1856 - g3 - integrand for the test function gradient and basis function gradient term
1857 
1858   Note: We are using a first order FEM model for the weak form:
1859 
1860   \int_\Gamma \phi {\vec g}_0(u, u_t, \nabla u, x, t) \cdot \hat n \psi + \phi {\vec g}_1(u, u_t, \nabla u, x, t) \cdot \hat n \nabla \psi + \nabla\phi \cdot {\vec g}_2(u, u_t, \nabla u, x, t) \cdot \hat n \psi + \nabla\phi \cdot {\overleftrightarrow g}_3(u, u_t, \nabla u, x, t) \cdot \hat n \cdot \nabla \psi
1861 
1862 The calling sequence for the callbacks g0, g1, g2 and g3 is given by:
1863 
1864 $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1865 $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1866 $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1867 $    PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar g0[])
1868 
1869 + dim - the spatial dimension
1870 . Nf - the number of fields
1871 . uOff - the offset into u[] and u_t[] for each field
1872 . uOff_x - the offset into u_x[] for each field
1873 . u - each field evaluated at the current point
1874 . u_t - the time derivative of each field evaluated at the current point
1875 . u_x - the gradient of each field evaluated at the current point
1876 . aOff - the offset into a[] and a_t[] for each auxiliary field
1877 . aOff_x - the offset into a_x[] for each auxiliary field
1878 . a - each auxiliary field evaluated at the current point
1879 . a_t - the time derivative of each auxiliary field evaluated at the current point
1880 . a_x - the gradient of auxiliary each field evaluated at the current point
1881 . t - current time
1882 . u_tShift - the multiplier a for dF/dU_t
1883 . x - coordinates of the current point
1884 . n - normal at the current point
1885 - g0 - output values at the current point
1886 
1887   Level: intermediate
1888 
1889 .seealso: PetscDSSetBdJacobian()
1890 @*/
1891 PetscErrorCode PetscDSGetBdJacobian(PetscDS prob, PetscInt f, PetscInt g,
1892                                     void (**g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1893                                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1894                                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1895                                                 PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscScalar g0[]),
1896                                     void (**g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1897                                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1898                                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1899                                                 PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscScalar g1[]),
1900                                     void (**g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1901                                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1902                                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1903                                                 PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscScalar g2[]),
1904                                     void (**g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1905                                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1906                                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1907                                                 PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscScalar g3[]))
1908 {
1909   PetscFunctionBegin;
1910   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1911   if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf);
1912   if ((g < 0) || (g >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", g, prob->Nf);
1913   if (g0) {PetscValidPointer(g0, 4); *g0 = prob->gBd[(f*prob->Nf + g)*4+0];}
1914   if (g1) {PetscValidPointer(g1, 5); *g1 = prob->gBd[(f*prob->Nf + g)*4+1];}
1915   if (g2) {PetscValidPointer(g2, 6); *g2 = prob->gBd[(f*prob->Nf + g)*4+2];}
1916   if (g3) {PetscValidPointer(g3, 7); *g3 = prob->gBd[(f*prob->Nf + g)*4+3];}
1917   PetscFunctionReturn(0);
1918 }
1919 
1920 #undef __FUNCT__
1921 #define __FUNCT__ "PetscDSSetBdJacobian"
1922 /*@C
1923   PetscDSSetBdJacobian - Set the pointwise boundary Jacobian function for given test and basis field
1924 
1925   Not collective
1926 
1927   Input Parameters:
1928 + prob - The PetscDS
1929 . f    - The test field number
1930 . g    - The field number
1931 . g0 - integrand for the test and basis function term
1932 . g1 - integrand for the test function and basis function gradient term
1933 . g2 - integrand for the test function gradient and basis function term
1934 - g3 - integrand for the test function gradient and basis function gradient term
1935 
1936   Note: We are using a first order FEM model for the weak form:
1937 
1938   \int_\Gamma \phi {\vec g}_0(u, u_t, \nabla u, x, t) \cdot \hat n \psi + \phi {\vec g}_1(u, u_t, \nabla u, x, t) \cdot \hat n \nabla \psi + \nabla\phi \cdot {\vec g}_2(u, u_t, \nabla u, x, t) \cdot \hat n \psi + \nabla\phi \cdot {\overleftrightarrow g}_3(u, u_t, \nabla u, x, t) \cdot \hat n \cdot \nabla \psi
1939 
1940 The calling sequence for the callbacks g0, g1, g2 and g3 is given by:
1941 
1942 $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1943 $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1944 $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1945 $    PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar g0[])
1946 
1947 + dim - the spatial dimension
1948 . Nf - the number of fields
1949 . uOff - the offset into u[] and u_t[] for each field
1950 . uOff_x - the offset into u_x[] for each field
1951 . u - each field evaluated at the current point
1952 . u_t - the time derivative of each field evaluated at the current point
1953 . u_x - the gradient of each field evaluated at the current point
1954 . aOff - the offset into a[] and a_t[] for each auxiliary field
1955 . aOff_x - the offset into a_x[] for each auxiliary field
1956 . a - each auxiliary field evaluated at the current point
1957 . a_t - the time derivative of each auxiliary field evaluated at the current point
1958 . a_x - the gradient of auxiliary each field evaluated at the current point
1959 . t - current time
1960 . u_tShift - the multiplier a for dF/dU_t
1961 . x - coordinates of the current point
1962 . n - normal at the current point
1963 - g0 - output values at the current point
1964 
1965   Level: intermediate
1966 
1967 .seealso: PetscDSGetBdJacobian()
1968 @*/
1969 PetscErrorCode PetscDSSetBdJacobian(PetscDS prob, PetscInt f, PetscInt g,
1970                                     void (*g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1971                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1972                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1973                                                PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscScalar g0[]),
1974                                     void (*g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1975                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1976                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1977                                                PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscScalar g1[]),
1978                                     void (*g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1979                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1980                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1981                                                PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscScalar g2[]),
1982                                     void (*g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1983                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1984                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1985                                                PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscScalar g3[]))
1986 {
1987   PetscErrorCode ierr;
1988 
1989   PetscFunctionBegin;
1990   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1991   if (g0) PetscValidFunction(g0, 4);
1992   if (g1) PetscValidFunction(g1, 5);
1993   if (g2) PetscValidFunction(g2, 6);
1994   if (g3) PetscValidFunction(g3, 7);
1995   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
1996   if (g < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", g);
1997   ierr = PetscDSEnlarge_Static(prob, PetscMax(f, g)+1);CHKERRQ(ierr);
1998   prob->gBd[(f*prob->Nf + g)*4+0] = g0;
1999   prob->gBd[(f*prob->Nf + g)*4+1] = g1;
2000   prob->gBd[(f*prob->Nf + g)*4+2] = g2;
2001   prob->gBd[(f*prob->Nf + g)*4+3] = g3;
2002   PetscFunctionReturn(0);
2003 }
2004 
2005 #undef __FUNCT__
2006 #define __FUNCT__ "PetscDSGetFieldIndex"
2007 /*@
2008   PetscDSGetFieldIndex - Returns the index of the given field
2009 
2010   Not collective
2011 
2012   Input Parameters:
2013 + prob - The PetscDS object
2014 - disc - The discretization object
2015 
2016   Output Parameter:
2017 . f - The field number
2018 
2019   Level: beginner
2020 
2021 .seealso: PetscGetDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
2022 @*/
2023 PetscErrorCode PetscDSGetFieldIndex(PetscDS prob, PetscObject disc, PetscInt *f)
2024 {
2025   PetscInt g;
2026 
2027   PetscFunctionBegin;
2028   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
2029   PetscValidPointer(f, 3);
2030   *f = -1;
2031   for (g = 0; g < prob->Nf; ++g) {if (disc == prob->disc[g]) break;}
2032   if (g == prob->Nf) SETERRQ(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Field not found in PetscDS.");
2033   *f = g;
2034   PetscFunctionReturn(0);
2035 }
2036 
2037 #undef __FUNCT__
2038 #define __FUNCT__ "PetscDSGetFieldSize"
2039 /*@
2040   PetscDSGetFieldSize - Returns the size of the given field in the full space basis
2041 
2042   Not collective
2043 
2044   Input Parameters:
2045 + prob - The PetscDS object
2046 - f - The field number
2047 
2048   Output Parameter:
2049 . size - The size
2050 
2051   Level: beginner
2052 
2053 .seealso: PetscDSGetFieldOffset(), PetscDSGetNumFields(), PetscDSCreate()
2054 @*/
2055 PetscErrorCode PetscDSGetFieldSize(PetscDS prob, PetscInt f, PetscInt *size)
2056 {
2057   PetscInt       Nb, Nc;
2058   PetscErrorCode ierr;
2059 
2060   PetscFunctionBegin;
2061   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
2062   PetscValidPointer(size, 3);
2063   if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf);
2064   *size = prob->Nc[f] * prob->Nb[f];
2065   PetscFunctionReturn(0);
2066 }
2067 
2068 #undef __FUNCT__
2069 #define __FUNCT__ "PetscDSGetFieldOffset"
2070 /*@
2071   PetscDSGetFieldOffset - Returns the offset of the given field in the full space basis
2072 
2073   Not collective
2074 
2075   Input Parameters:
2076 + prob - The PetscDS object
2077 - f - The field number
2078 
2079   Output Parameter:
2080 . off - The offset
2081 
2082   Level: beginner
2083 
2084 .seealso: PetscDSGetFieldSize(), PetscDSGetNumFields(), PetscDSCreate()
2085 @*/
2086 PetscErrorCode PetscDSGetFieldOffset(PetscDS prob, PetscInt f, PetscInt *off)
2087 {
2088   PetscInt       size, g;
2089   PetscErrorCode ierr;
2090 
2091   PetscFunctionBegin;
2092   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
2093   PetscValidPointer(off, 3);
2094   if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf);
2095   *off = 0;
2096   for (g = 0; g < f; ++g) {
2097     ierr = PetscDSGetFieldSize(prob, g, &size);CHKERRQ(ierr);
2098     *off += size;
2099   }
2100   PetscFunctionReturn(0);
2101 }
2102 
2103 #undef __FUNCT__
2104 #define __FUNCT__ "PetscDSGetDimensions"
2105 /*@
2106   PetscDSGetDimensions - Returns the size of the approximation space for each field on an evaluation point
2107 
2108   Not collective
2109 
2110   Input Parameter:
2111 . prob - The PetscDS object
2112 
2113   Output Parameter:
2114 . dimensions - The number of dimensions
2115 
2116   Level: beginner
2117 
2118 .seealso: PetscDSGetComponentOffsets(), PetscDSGetNumFields(), PetscDSCreate()
2119 @*/
2120 PetscErrorCode PetscDSGetDimensions(PetscDS prob, PetscInt *dimensions[])
2121 {
2122   PetscErrorCode ierr;
2123 
2124   PetscFunctionBegin;
2125   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
2126   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
2127   PetscValidPointer(dimensions, 2);
2128   *dimensions = prob->Nb;
2129   PetscFunctionReturn(0);
2130 }
2131 
2132 #undef __FUNCT__
2133 #define __FUNCT__ "PetscDSGetComponents"
2134 /*@
2135   PetscDSGetComponents - Returns the number of components for each field on an evaluation point
2136 
2137   Not collective
2138 
2139   Input Parameter:
2140 . prob - The PetscDS object
2141 
2142   Output Parameter:
2143 . components - The number of components
2144 
2145   Level: beginner
2146 
2147 .seealso: PetscDSGetComponentOffsets(), PetscDSGetNumFields(), PetscDSCreate()
2148 @*/
2149 PetscErrorCode PetscDSGetComponents(PetscDS prob, PetscInt *components[])
2150 {
2151   PetscErrorCode ierr;
2152 
2153   PetscFunctionBegin;
2154   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
2155   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
2156   PetscValidPointer(components, 2);
2157   *components = prob->Nc;
2158   PetscFunctionReturn(0);
2159 }
2160 
2161 #undef __FUNCT__
2162 #define __FUNCT__ "PetscDSGetComponentOffset"
2163 /*@
2164   PetscDSGetComponentOffset - Returns the offset of the given field on an evaluation point
2165 
2166   Not collective
2167 
2168   Input Parameters:
2169 + prob - The PetscDS object
2170 - f - The field number
2171 
2172   Output Parameter:
2173 . off - The offset
2174 
2175   Level: beginner
2176 
2177 .seealso: PetscDSGetNumFields(), PetscDSCreate()
2178 @*/
2179 PetscErrorCode PetscDSGetComponentOffset(PetscDS prob, PetscInt f, PetscInt *off)
2180 {
2181   PetscInt       g;
2182   PetscErrorCode ierr;
2183 
2184   PetscFunctionBegin;
2185   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
2186   PetscValidPointer(off, 3);
2187   if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf);
2188   *off = prob->off[f];
2189   PetscFunctionReturn(0);
2190 }
2191 
2192 #undef __FUNCT__
2193 #define __FUNCT__ "PetscDSGetComponentOffsets"
2194 /*@
2195   PetscDSGetComponentOffsets - Returns the offset of each field on an evaluation point
2196 
2197   Not collective
2198 
2199   Input Parameter:
2200 . prob - The PetscDS object
2201 
2202   Output Parameter:
2203 . offsets - The offsets
2204 
2205   Level: beginner
2206 
2207 .seealso: PetscDSGetNumFields(), PetscDSCreate()
2208 @*/
2209 PetscErrorCode PetscDSGetComponentOffsets(PetscDS prob, PetscInt *offsets[])
2210 {
2211   PetscFunctionBegin;
2212   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
2213   PetscValidPointer(offsets, 2);
2214   *offsets = prob->off;
2215   PetscFunctionReturn(0);
2216 }
2217 
2218 #undef __FUNCT__
2219 #define __FUNCT__ "PetscDSGetComponentDerivativeOffsets"
2220 /*@
2221   PetscDSGetComponentDerivativeOffsets - Returns the offset of each field derivative on an evaluation point
2222 
2223   Not collective
2224 
2225   Input Parameter:
2226 . prob - The PetscDS object
2227 
2228   Output Parameter:
2229 . offsets - The offsets
2230 
2231   Level: beginner
2232 
2233 .seealso: PetscDSGetNumFields(), PetscDSCreate()
2234 @*/
2235 PetscErrorCode PetscDSGetComponentDerivativeOffsets(PetscDS prob, PetscInt *offsets[])
2236 {
2237   PetscFunctionBegin;
2238   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
2239   PetscValidPointer(offsets, 2);
2240   *offsets = prob->offDer;
2241   PetscFunctionReturn(0);
2242 }
2243 
2244 #undef __FUNCT__
2245 #define __FUNCT__ "PetscDSGetTabulation"
2246 /*@C
2247   PetscDSGetTabulation - Return the basis tabulation at quadrature points for the volume discretization
2248 
2249   Not collective
2250 
2251   Input Parameter:
2252 . prob - The PetscDS object
2253 
2254   Output Parameters:
2255 + basis - The basis function tabulation at quadrature points
2256 - basisDer - The basis function derivative tabulation at quadrature points
2257 
2258   Level: intermediate
2259 
2260 .seealso: PetscDSCreate()
2261 @*/
2262 PetscErrorCode PetscDSGetTabulation(PetscDS prob, PetscReal ***basis, PetscReal ***basisDer)
2263 {
2264   PetscErrorCode ierr;
2265 
2266   PetscFunctionBegin;
2267   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
2268   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
2269   if (basis)    {PetscValidPointer(basis, 2);    *basis    = prob->basis;}
2270   if (basisDer) {PetscValidPointer(basisDer, 3); *basisDer = prob->basisDer;}
2271   PetscFunctionReturn(0);
2272 }
2273 
2274 #undef __FUNCT__
2275 #define __FUNCT__ "PetscDSGetFaceTabulation"
2276 /*@C
2277   PetscDSGetFaceTabulation - Return the basis tabulation at quadrature points on the faces
2278 
2279   Not collective
2280 
2281   Input Parameter:
2282 . prob - The PetscDS object
2283 
2284   Output Parameters:
2285 + basisFace - The basis function tabulation at quadrature points
2286 - basisDerFace - The basis function derivative tabulation at quadrature points
2287 
2288   Level: intermediate
2289 
2290 .seealso: PetscDSGetTabulation(), PetscDSCreate()
2291 @*/
2292 PetscErrorCode PetscDSGetFaceTabulation(PetscDS prob, PetscReal ***basis, PetscReal ***basisDer)
2293 {
2294   PetscErrorCode ierr;
2295 
2296   PetscFunctionBegin;
2297   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
2298   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
2299   if (basis)    {PetscValidPointer(basis, 2);    *basis    = prob->basisFace;}
2300   if (basisDer) {PetscValidPointer(basisDer, 3); *basisDer = prob->basisDerFace;}
2301   PetscFunctionReturn(0);
2302 }
2303 
2304 #undef __FUNCT__
2305 #define __FUNCT__ "PetscDSGetEvaluationArrays"
2306 PetscErrorCode PetscDSGetEvaluationArrays(PetscDS prob, PetscScalar **u, PetscScalar **u_t, PetscScalar **u_x)
2307 {
2308   PetscErrorCode ierr;
2309 
2310   PetscFunctionBegin;
2311   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
2312   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
2313   if (u)   {PetscValidPointer(u, 2);   *u   = prob->u;}
2314   if (u_t) {PetscValidPointer(u_t, 3); *u_t = prob->u_t;}
2315   if (u_x) {PetscValidPointer(u_x, 4); *u_x = prob->u_x;}
2316   PetscFunctionReturn(0);
2317 }
2318 
2319 #undef __FUNCT__
2320 #define __FUNCT__ "PetscDSGetWeakFormArrays"
2321 PetscErrorCode PetscDSGetWeakFormArrays(PetscDS prob, PetscScalar **f0, PetscScalar **f1, PetscScalar **g0, PetscScalar **g1, PetscScalar **g2, PetscScalar **g3)
2322 {
2323   PetscErrorCode ierr;
2324 
2325   PetscFunctionBegin;
2326   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
2327   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
2328   if (f0) {PetscValidPointer(f0, 2); *f0 = prob->f0;}
2329   if (f1) {PetscValidPointer(f1, 3); *f1 = prob->f1;}
2330   if (g0) {PetscValidPointer(g0, 4); *g0 = prob->g0;}
2331   if (g1) {PetscValidPointer(g1, 5); *g1 = prob->g1;}
2332   if (g2) {PetscValidPointer(g2, 6); *g2 = prob->g2;}
2333   if (g3) {PetscValidPointer(g3, 7); *g3 = prob->g3;}
2334   PetscFunctionReturn(0);
2335 }
2336 
2337 #undef __FUNCT__
2338 #define __FUNCT__ "PetscDSGetRefCoordArrays"
2339 PetscErrorCode PetscDSGetRefCoordArrays(PetscDS prob, PetscReal **x, PetscScalar **refSpaceDer)
2340 {
2341   PetscErrorCode ierr;
2342 
2343   PetscFunctionBegin;
2344   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
2345   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
2346   if (x)           {PetscValidPointer(x, 2);           *x           = prob->x;}
2347   if (refSpaceDer) {PetscValidPointer(refSpaceDer, 3); *refSpaceDer = prob->refSpaceDer;}
2348   PetscFunctionReturn(0);
2349 }
2350 
2351 #undef __FUNCT__
2352 #define __FUNCT__ "PetscDSAddBoundary"
2353 /*@C
2354   PetscDSAddBoundary - Add a boundary condition to the model
2355 
2356   Input Parameters:
2357 + ds          - The PetscDS object
2358 . isEssential - Flag for an essential (Dirichlet) condition, as opposed to a natural (Neumann) condition
2359 . name        - The BC name
2360 . labelname   - The label defining constrained points
2361 . field       - The field to constrain
2362 . numcomps    - The number of constrained field components
2363 . comps       - An array of constrained component numbers
2364 . bcFunc      - A pointwise function giving boundary values
2365 . numids      - The number of DMLabel ids for constrained points
2366 . ids         - An array of ids for constrained points
2367 - ctx         - An optional user context for bcFunc
2368 
2369   Options Database Keys:
2370 + -bc_<boundary name> <num> - Overrides the boundary ids
2371 - -bc_<boundary name>_comp <num> - Overrides the boundary components
2372 
2373   Level: developer
2374 
2375 .seealso: PetscDSGetBoundary()
2376 @*/
2377 PetscErrorCode PetscDSAddBoundary(PetscDS ds, PetscBool isEssential, const char name[], const char labelname[], PetscInt field, PetscInt numcomps, const PetscInt *comps, void (*bcFunc)(), PetscInt numids, const PetscInt *ids, void *ctx)
2378 {
2379   DSBoundary     b;
2380   PetscErrorCode ierr;
2381 
2382   PetscFunctionBegin;
2383   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
2384   ierr = PetscNew(&b);CHKERRQ(ierr);
2385   ierr = PetscStrallocpy(name, (char **) &b->name);CHKERRQ(ierr);
2386   ierr = PetscStrallocpy(labelname, (char **) &b->labelname);CHKERRQ(ierr);
2387   ierr = PetscMalloc1(numcomps, &b->comps);CHKERRQ(ierr);
2388   if (numcomps) {ierr = PetscMemcpy(b->comps, comps, numcomps*sizeof(PetscInt));CHKERRQ(ierr);}
2389   ierr = PetscMalloc1(numids, &b->ids);CHKERRQ(ierr);
2390   if (numids) {ierr = PetscMemcpy(b->ids, ids, numids*sizeof(PetscInt));CHKERRQ(ierr);}
2391   b->essential       = isEssential;
2392   b->field           = field;
2393   b->numcomps        = numcomps;
2394   b->func            = bcFunc;
2395   b->numids          = numids;
2396   b->ctx             = ctx;
2397   b->next            = ds->boundary;
2398   ds->boundary       = b;
2399   PetscFunctionReturn(0);
2400 }
2401 
2402 #undef __FUNCT__
2403 #define __FUNCT__ "PetscDSGetNumBoundary"
2404 /*@
2405   PetscDSGetNumBoundary - Get the number of registered BC
2406 
2407   Input Parameters:
2408 . ds - The PetscDS object
2409 
2410   Output Parameters:
2411 . numBd - The number of BC
2412 
2413   Level: intermediate
2414 
2415 .seealso: PetscDSAddBoundary(), PetscDSGetBoundary()
2416 @*/
2417 PetscErrorCode PetscDSGetNumBoundary(PetscDS ds, PetscInt *numBd)
2418 {
2419   DSBoundary b = ds->boundary;
2420 
2421   PetscFunctionBegin;
2422   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
2423   PetscValidPointer(numBd, 2);
2424   *numBd = 0;
2425   while (b) {++(*numBd); b = b->next;}
2426   PetscFunctionReturn(0);
2427 }
2428 
2429 #undef __FUNCT__
2430 #define __FUNCT__ "PetscDSGetBoundary"
2431 /*@C
2432   PetscDSGetBoundary - Add a boundary condition to the model
2433 
2434   Input Parameters:
2435 + ds          - The PetscDS object
2436 - bd          - The BC number
2437 
2438   Output Parameters:
2439 + isEssential - Flag for an essential (Dirichlet) condition, as opposed to a natural (Neumann) condition
2440 . name        - The BC name
2441 . labelname   - The label defining constrained points
2442 . field       - The field to constrain
2443 . numcomps    - The number of constrained field components
2444 . comps       - An array of constrained component numbers
2445 . bcFunc      - A pointwise function giving boundary values
2446 . numids      - The number of DMLabel ids for constrained points
2447 . ids         - An array of ids for constrained points
2448 - ctx         - An optional user context for bcFunc
2449 
2450   Options Database Keys:
2451 + -bc_<boundary name> <num> - Overrides the boundary ids
2452 - -bc_<boundary name>_comp <num> - Overrides the boundary components
2453 
2454   Level: developer
2455 
2456 .seealso: PetscDSAddBoundary()
2457 @*/
2458 PetscErrorCode PetscDSGetBoundary(PetscDS ds, PetscInt bd, PetscBool *isEssential, const char **name, const char **labelname, PetscInt *field, PetscInt *numcomps, const PetscInt **comps, void (**func)(), PetscInt *numids, const PetscInt **ids, void **ctx)
2459 {
2460   DSBoundary b    = ds->boundary;
2461   PetscInt   n    = 0;
2462 
2463   PetscFunctionBegin;
2464   PetscValidHeaderSpecific(ds, PETSCDS_CLASSID, 1);
2465   while (b) {
2466     if (n == bd) break;
2467     b = b->next;
2468     ++n;
2469   }
2470   if (!b) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Boundary %d is not in [0, %d)", bd, n);
2471   if (isEssential) {
2472     PetscValidPointer(isEssential, 3);
2473     *isEssential = b->essential;
2474   }
2475   if (name) {
2476     PetscValidPointer(name, 4);
2477     *name = b->name;
2478   }
2479   if (labelname) {
2480     PetscValidPointer(labelname, 5);
2481     *labelname = b->labelname;
2482   }
2483   if (field) {
2484     PetscValidPointer(field, 6);
2485     *field = b->field;
2486   }
2487   if (numcomps) {
2488     PetscValidPointer(numcomps, 7);
2489     *numcomps = b->numcomps;
2490   }
2491   if (comps) {
2492     PetscValidPointer(comps, 8);
2493     *comps = b->comps;
2494   }
2495   if (func) {
2496     PetscValidPointer(func, 9);
2497     *func = b->func;
2498   }
2499   if (numids) {
2500     PetscValidPointer(numids, 10);
2501     *numids = b->numids;
2502   }
2503   if (ids) {
2504     PetscValidPointer(ids, 11);
2505     *ids = b->ids;
2506   }
2507   if (ctx) {
2508     PetscValidPointer(ctx, 12);
2509     *ctx = b->ctx;
2510   }
2511   PetscFunctionReturn(0);
2512 }
2513 
2514 #undef __FUNCT__
2515 #define __FUNCT__ "PetscDSCopyBoundary"
2516 PetscErrorCode PetscDSCopyBoundary(PetscDS probA, PetscDS probB)
2517 {
2518   DSBoundary     b, next, *lastnext;
2519   PetscErrorCode ierr;
2520 
2521   PetscFunctionBegin;
2522   PetscValidHeaderSpecific(probA, PETSCDS_CLASSID, 1);
2523   PetscValidHeaderSpecific(probB, PETSCDS_CLASSID, 2);
2524   if (probA == probB) PetscFunctionReturn(0);
2525   next = probB->boundary;
2526   while (next) {
2527     DSBoundary b = next;
2528 
2529     next = b->next;
2530     ierr = PetscFree(b->comps);CHKERRQ(ierr);
2531     ierr = PetscFree(b->ids);CHKERRQ(ierr);
2532     ierr = PetscFree(b->name);CHKERRQ(ierr);
2533     ierr = PetscFree(b->labelname);CHKERRQ(ierr);
2534     ierr = PetscFree(b);CHKERRQ(ierr);
2535   }
2536   lastnext = &(probB->boundary);
2537   for (b = probA->boundary; b; b = b->next) {
2538     DSBoundary bNew;
2539 
2540     ierr = PetscNew(&bNew);
2541     bNew->numcomps = b->numcomps;
2542     ierr = PetscMalloc1(bNew->numcomps, &bNew->comps);CHKERRQ(ierr);
2543     ierr = PetscMemcpy(bNew->comps, b->comps, bNew->numcomps*sizeof(PetscInt));CHKERRQ(ierr);
2544     bNew->numids = b->numids;
2545     ierr = PetscMalloc1(bNew->numids, &bNew->ids);CHKERRQ(ierr);
2546     ierr = PetscMemcpy(bNew->ids, b->ids, bNew->numids*sizeof(PetscInt));CHKERRQ(ierr);
2547     ierr = PetscStrallocpy(b->labelname,(char **) &(bNew->labelname));CHKERRQ(ierr);
2548     ierr = PetscStrallocpy(b->name,(char **) &(bNew->name));CHKERRQ(ierr);
2549     bNew->ctx       = b->ctx;
2550     bNew->essential = b->essential;
2551     bNew->field     = b->field;
2552     bNew->func      = b->func;
2553 
2554     *lastnext = bNew;
2555     lastnext = &(bNew->next);
2556   }
2557   PetscFunctionReturn(0);
2558 }
2559 
2560 #undef __FUNCT__
2561 #define __FUNCT__ "PetscDSCopyEquations"
2562 /*@
2563   PetscDSCopyEquations - Copy all pointwise function pointers to the new problem
2564 
2565   Not collective
2566 
2567   Input Parameter:
2568 . prob - The PetscDS object
2569 
2570   Output Parameter:
2571 . newprob - The PetscDS copy
2572 
2573   Level: intermediate
2574 
2575 .seealso: PetscDSSetResidual(), PetscDSSetJacobian(), PetscDSSetRiemannSolver(), PetscDSSetBdResidual(), PetscDSSetBdJacobian(), PetscDSCreate()
2576 @*/
2577 PetscErrorCode PetscDSCopyEquations(PetscDS prob, PetscDS newprob)
2578 {
2579   PetscInt       Nf, Ng, f, g;
2580   PetscErrorCode ierr;
2581 
2582   PetscFunctionBegin;
2583   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
2584   PetscValidHeaderSpecific(newprob, PETSCDS_CLASSID, 2);
2585   ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr);
2586   ierr = PetscDSGetNumFields(newprob, &Ng);CHKERRQ(ierr);
2587   if (Nf != Ng) SETERRQ2(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_SIZ, "Number of fields must match %D != %D", Nf, Ng);CHKERRQ(ierr);
2588   for (f = 0; f < Nf; ++f) {
2589     PetscPointFunc   obj;
2590     PetscPointFunc   f0, f1;
2591     PetscPointJac    g0, g1, g2, g3;
2592     PetscBdPointFunc f0Bd, f1Bd;
2593     PetscBdPointJac  g0Bd, g1Bd, g2Bd, g3Bd;
2594     PetscRiemannFunc r;
2595 
2596     ierr = PetscDSGetObjective(prob, f, &obj);CHKERRQ(ierr);
2597     ierr = PetscDSGetResidual(prob, f, &f0, &f1);CHKERRQ(ierr);
2598     ierr = PetscDSGetBdResidual(prob, f, &f0Bd, &f1Bd);CHKERRQ(ierr);
2599     ierr = PetscDSGetRiemannSolver(prob, f, &r);CHKERRQ(ierr);
2600     ierr = PetscDSSetObjective(newprob, f, obj);CHKERRQ(ierr);
2601     ierr = PetscDSSetResidual(newprob, f, f0, f1);CHKERRQ(ierr);
2602     ierr = PetscDSSetBdResidual(newprob, f, f0Bd, f1Bd);CHKERRQ(ierr);
2603     ierr = PetscDSSetRiemannSolver(newprob, f, r);CHKERRQ(ierr);
2604     for (g = 0; g < Nf; ++g) {
2605       ierr = PetscDSGetJacobian(prob, f, g, &g0, &g1, &g2, &g3);CHKERRQ(ierr);
2606       ierr = PetscDSGetBdJacobian(prob, f, g, &g0Bd, &g1Bd, &g2Bd, &g3Bd);CHKERRQ(ierr);
2607       ierr = PetscDSSetJacobian(newprob, f, g, g0, g1, g2, g3);CHKERRQ(ierr);
2608       ierr = PetscDSSetBdJacobian(newprob, f, g, g0Bd, g1Bd, g2Bd, g3Bd);CHKERRQ(ierr);
2609     }
2610   }
2611   PetscFunctionReturn(0);
2612 }
2613 
2614 #undef __FUNCT__
2615 #define __FUNCT__ "PetscDSDestroy_Basic"
2616 static PetscErrorCode PetscDSDestroy_Basic(PetscDS prob)
2617 {
2618   PetscErrorCode      ierr;
2619 
2620   PetscFunctionBegin;
2621   ierr = PetscFree(prob->data);CHKERRQ(ierr);
2622   PetscFunctionReturn(0);
2623 }
2624 
2625 #undef __FUNCT__
2626 #define __FUNCT__ "PetscDSInitialize_Basic"
2627 static PetscErrorCode PetscDSInitialize_Basic(PetscDS prob)
2628 {
2629   PetscFunctionBegin;
2630   prob->ops->setfromoptions = NULL;
2631   prob->ops->setup          = NULL;
2632   prob->ops->view           = NULL;
2633   prob->ops->destroy        = PetscDSDestroy_Basic;
2634   PetscFunctionReturn(0);
2635 }
2636 
2637 /*MC
2638   PETSCDSBASIC = "basic" - A discrete system with pointwise residual and boundary residual functions
2639 
2640   Level: intermediate
2641 
2642 .seealso: PetscDSType, PetscDSCreate(), PetscDSSetType()
2643 M*/
2644 
2645 #undef __FUNCT__
2646 #define __FUNCT__ "PetscDSCreate_Basic"
2647 PETSC_EXTERN PetscErrorCode PetscDSCreate_Basic(PetscDS prob)
2648 {
2649   PetscDS_Basic *b;
2650   PetscErrorCode      ierr;
2651 
2652   PetscFunctionBegin;
2653   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
2654   ierr       = PetscNewLog(prob, &b);CHKERRQ(ierr);
2655   prob->data = b;
2656 
2657   ierr = PetscDSInitialize_Basic(prob);CHKERRQ(ierr);
2658   PetscFunctionReturn(0);
2659 }
2660