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