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