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