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