xref: /petsc/src/dm/dt/interface/dtds.c (revision 7d8a60ead99f2456a3d64dee7a317e4840463142)
1 #include <petsc-private/petscdsimpl.h> /*I "petscds.h" I*/
2 
3 PetscClassId PETSCDS_CLASSID = 0;
4 
5 PetscFunctionList PetscDSList              = NULL;
6 PetscBool         PetscDSRegisterAllCalled = PETSC_FALSE;
7 
8 #undef __FUNCT__
9 #define __FUNCT__ "PetscDSRegister"
10 /*@C
11   PetscDSRegister - Adds a new PetscDS implementation
12 
13   Not Collective
14 
15   Input Parameters:
16 + name        - The name of a new user-defined creation routine
17 - create_func - The creation routine itself
18 
19   Notes:
20   PetscDSRegister() may be called multiple times to add several user-defined PetscDSs
21 
22   Sample usage:
23 .vb
24     PetscDSRegister("my_ds", MyPetscDSCreate);
25 .ve
26 
27   Then, your PetscDS type can be chosen with the procedural interface via
28 .vb
29     PetscDSCreate(MPI_Comm, PetscDS *);
30     PetscDSSetType(PetscDS, "my_ds");
31 .ve
32    or at runtime via the option
33 .vb
34     -petscds_type my_ds
35 .ve
36 
37   Level: advanced
38 
39 .keywords: PetscDS, register
40 .seealso: PetscDSRegisterAll(), PetscDSRegisterDestroy()
41 
42 @*/
43 PetscErrorCode PetscDSRegister(const char sname[], PetscErrorCode (*function)(PetscDS))
44 {
45   PetscErrorCode ierr;
46 
47   PetscFunctionBegin;
48   ierr = PetscFunctionListAdd(&PetscDSList, sname, function);CHKERRQ(ierr);
49   PetscFunctionReturn(0);
50 }
51 
52 #undef __FUNCT__
53 #define __FUNCT__ "PetscDSSetType"
54 /*@C
55   PetscDSSetType - Builds a particular PetscDS
56 
57   Collective on PetscDS
58 
59   Input Parameters:
60 + prob - The PetscDS object
61 - name - The kind of system
62 
63   Options Database Key:
64 . -petscds_type <type> - Sets the PetscDS type; use -help for a list of available types
65 
66   Level: intermediate
67 
68 .keywords: PetscDS, set, type
69 .seealso: PetscDSGetType(), PetscDSCreate()
70 @*/
71 PetscErrorCode PetscDSSetType(PetscDS prob, PetscDSType name)
72 {
73   PetscErrorCode (*r)(PetscDS);
74   PetscBool      match;
75   PetscErrorCode ierr;
76 
77   PetscFunctionBegin;
78   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
79   ierr = PetscObjectTypeCompare((PetscObject) prob, name, &match);CHKERRQ(ierr);
80   if (match) PetscFunctionReturn(0);
81 
82   if (!PetscDSRegisterAllCalled) {ierr = PetscDSRegisterAll();CHKERRQ(ierr);}
83   ierr = PetscFunctionListFind(PetscDSList, name, &r);CHKERRQ(ierr);
84   if (!r) SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscDS type: %s", name);
85 
86   if (prob->ops->destroy) {
87     ierr             = (*prob->ops->destroy)(prob);CHKERRQ(ierr);
88     prob->ops->destroy = NULL;
89   }
90   ierr = (*r)(prob);CHKERRQ(ierr);
91   ierr = PetscObjectChangeTypeName((PetscObject) prob, name);CHKERRQ(ierr);
92   PetscFunctionReturn(0);
93 }
94 
95 #undef __FUNCT__
96 #define __FUNCT__ "PetscDSGetType"
97 /*@C
98   PetscDSGetType - Gets the PetscDS type name (as a string) from the object.
99 
100   Not Collective
101 
102   Input Parameter:
103 . prob  - The PetscDS
104 
105   Output Parameter:
106 . name - The PetscDS type name
107 
108   Level: intermediate
109 
110 .keywords: PetscDS, get, type, name
111 .seealso: PetscDSSetType(), PetscDSCreate()
112 @*/
113 PetscErrorCode PetscDSGetType(PetscDS prob, PetscDSType *name)
114 {
115   PetscErrorCode ierr;
116 
117   PetscFunctionBegin;
118   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
119   PetscValidPointer(name, 2);
120   if (!PetscDSRegisterAllCalled) {ierr = PetscDSRegisterAll();CHKERRQ(ierr);}
121   *name = ((PetscObject) prob)->type_name;
122   PetscFunctionReturn(0);
123 }
124 
125 #undef __FUNCT__
126 #define __FUNCT__ "PetscDSView_Ascii"
127 static PetscErrorCode PetscDSView_Ascii(PetscDS prob, PetscViewer viewer)
128 {
129   PetscViewerFormat format;
130   PetscInt          f;
131   PetscErrorCode    ierr;
132 
133   PetscFunctionBegin;
134   ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
135   ierr = PetscViewerASCIIPrintf(viewer, "Discrete System with %d fields\n", prob->Nf);CHKERRQ(ierr);
136   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
137   for (f = 0; f < prob->Nf; ++f) {
138     PetscObject  obj;
139     PetscClassId id;
140     const char  *name;
141     PetscInt     Nc;
142 
143     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
144     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
145     ierr = PetscObjectGetName(obj, &name);CHKERRQ(ierr);
146     ierr = PetscViewerASCIIPrintf(viewer, "Field %s", name ? name : "<unknown>");CHKERRQ(ierr);
147     if (id == PETSCFE_CLASSID)      {
148       ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);
149       ierr = PetscViewerASCIIPrintf(viewer, " FEM");CHKERRQ(ierr);
150     } else if (id == PETSCFV_CLASSID) {
151       ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);
152       ierr = PetscViewerASCIIPrintf(viewer, " FVM");CHKERRQ(ierr);
153     }
154     else SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
155     if (Nc > 1) {ierr = PetscViewerASCIIPrintf(viewer, "%d components", Nc);CHKERRQ(ierr);}
156     else        {ierr = PetscViewerASCIIPrintf(viewer, "%d component", Nc);CHKERRQ(ierr);}
157     ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr);
158     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
159       if (id == PETSCFE_CLASSID)      {ierr = PetscFEView((PetscFE) obj, viewer);CHKERRQ(ierr);}
160       else if (id == PETSCFV_CLASSID) {ierr = PetscFVView((PetscFV) obj, viewer);CHKERRQ(ierr);}
161     }
162   }
163   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
164   PetscFunctionReturn(0);
165 }
166 
167 #undef __FUNCT__
168 #define __FUNCT__ "PetscDSView"
169 /*@C
170   PetscDSView - Views a PetscDS
171 
172   Collective on PetscDS
173 
174   Input Parameter:
175 + prob - the PetscDS object to view
176 - v  - the viewer
177 
178   Level: developer
179 
180 .seealso PetscDSDestroy()
181 @*/
182 PetscErrorCode PetscDSView(PetscDS prob, PetscViewer v)
183 {
184   PetscBool      iascii;
185   PetscErrorCode ierr;
186 
187   PetscFunctionBegin;
188   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
189   if (!v) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) prob), &v);CHKERRQ(ierr);}
190   else    {PetscValidHeaderSpecific(v, PETSC_VIEWER_CLASSID, 2);}
191   ierr = PetscObjectTypeCompare((PetscObject) v, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
192   if (iascii) {ierr = PetscDSView_Ascii(prob, v);CHKERRQ(ierr);}
193   if (prob->ops->view) {ierr = (*prob->ops->view)(prob, v);CHKERRQ(ierr);}
194   PetscFunctionReturn(0);
195 }
196 
197 #undef __FUNCT__
198 #define __FUNCT__ "PetscDSViewFromOptions"
199 /*
200   PetscDSViewFromOptions - Processes command line options to determine if/how a PetscDS is to be viewed.
201 
202   Collective on PetscDS
203 
204   Input Parameters:
205 + prob   - the PetscDS
206 . prefix - prefix to use for viewing, or NULL to use prefix of 'rnd'
207 - optionname - option to activate viewing
208 
209   Level: intermediate
210 
211 .keywords: PetscDS, view, options, database
212 .seealso: VecViewFromOptions(), MatViewFromOptions()
213 */
214 PetscErrorCode PetscDSViewFromOptions(PetscDS prob, const char prefix[], const char optionname[])
215 {
216   PetscViewer       viewer;
217   PetscViewerFormat format;
218   PetscBool         flg;
219   PetscErrorCode    ierr;
220 
221   PetscFunctionBegin;
222   if (prefix) {ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject) prob), prefix, optionname, &viewer, &format, &flg);CHKERRQ(ierr);}
223   else        {ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject) prob), ((PetscObject) prob)->prefix, optionname, &viewer, &format, &flg);CHKERRQ(ierr);}
224   if (flg) {
225     ierr = PetscViewerPushFormat(viewer, format);CHKERRQ(ierr);
226     ierr = PetscDSView(prob, viewer);CHKERRQ(ierr);
227     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
228     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
229   }
230   PetscFunctionReturn(0);
231 }
232 
233 #undef __FUNCT__
234 #define __FUNCT__ "PetscDSSetFromOptions"
235 /*@
236   PetscDSSetFromOptions - sets parameters in a PetscDS from the options database
237 
238   Collective on PetscDS
239 
240   Input Parameter:
241 . prob - the PetscDS object to set options for
242 
243   Options Database:
244 
245   Level: developer
246 
247 .seealso PetscDSView()
248 @*/
249 PetscErrorCode PetscDSSetFromOptions(PetscDS prob)
250 {
251   const char    *defaultType;
252   char           name[256];
253   PetscBool      flg;
254   PetscErrorCode ierr;
255 
256   PetscFunctionBegin;
257   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
258   if (!((PetscObject) prob)->type_name) {
259     defaultType = PETSCDSBASIC;
260   } else {
261     defaultType = ((PetscObject) prob)->type_name;
262   }
263   if (!PetscDSRegisterAllCalled) {ierr = PetscDSRegisterAll();CHKERRQ(ierr);}
264 
265   ierr = PetscObjectOptionsBegin((PetscObject) prob);CHKERRQ(ierr);
266   ierr = PetscOptionsFList("-petscds_type", "Discrete System", "PetscDSSetType", PetscDSList, defaultType, name, 256, &flg);CHKERRQ(ierr);
267   if (flg) {
268     ierr = PetscDSSetType(prob, name);CHKERRQ(ierr);
269   } else if (!((PetscObject) prob)->type_name) {
270     ierr = PetscDSSetType(prob, defaultType);CHKERRQ(ierr);
271   }
272   if (prob->ops->setfromoptions) {ierr = (*prob->ops->setfromoptions)(prob);CHKERRQ(ierr);}
273   /* process any options handlers added with PetscObjectAddOptionsHandler() */
274   ierr = PetscObjectProcessOptionsHandlers((PetscObject) prob);CHKERRQ(ierr);
275   ierr = PetscOptionsEnd();CHKERRQ(ierr);
276   ierr = PetscDSViewFromOptions(prob, NULL, "-petscds_view");CHKERRQ(ierr);
277   PetscFunctionReturn(0);
278 }
279 
280 #undef __FUNCT__
281 #define __FUNCT__ "PetscDSSetUp"
282 /*@C
283   PetscDSSetUp - Construct data structures for the PetscDS
284 
285   Collective on PetscDS
286 
287   Input Parameter:
288 . prob - the PetscDS object to setup
289 
290   Level: developer
291 
292 .seealso PetscDSView(), PetscDSDestroy()
293 @*/
294 PetscErrorCode PetscDSSetUp(PetscDS prob)
295 {
296   const PetscInt Nf = prob->Nf;
297   PetscInt       dim, work, NcMax = 0, NqMax = 0, f;
298   PetscErrorCode ierr;
299 
300   PetscFunctionBegin;
301   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
302   if (prob->setup) PetscFunctionReturn(0);
303   /* Calculate sizes */
304   ierr = PetscDSGetSpatialDimension(prob, &dim);CHKERRQ(ierr);
305   prob->totDim = prob->totDimBd = prob->totComp = 0;
306   ierr = PetscMalloc4(Nf,&prob->basis,Nf,&prob->basisDer,Nf,&prob->basisBd,Nf,&prob->basisDerBd);CHKERRQ(ierr);
307   for (f = 0; f < Nf; ++f) {
308     PetscFE         feBd = (PetscFE) prob->discBd[f];
309     PetscObject     obj;
310     PetscClassId    id;
311     PetscQuadrature q;
312     PetscInt        Nq = 0, Nb, Nc;
313 
314     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
315     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
316     if (id == PETSCFE_CLASSID)      {
317       PetscFE fe = (PetscFE) obj;
318 
319       ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
320       ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
321       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
322       ierr = PetscFEGetDefaultTabulation(fe, &prob->basis[f], &prob->basisDer[f], NULL);CHKERRQ(ierr);
323     } else if (id == PETSCFV_CLASSID) {
324       PetscFV fv = (PetscFV) obj;
325 
326       ierr = PetscFVGetQuadrature(fv, &q);CHKERRQ(ierr);
327       Nb   = 1;
328       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
329     } else SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", 0);
330     if (q) {ierr = PetscQuadratureGetData(q, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);}
331     NqMax          = PetscMax(NqMax, Nq);
332     NcMax          = PetscMax(NcMax, Nc);
333     prob->totDim  += Nb*Nc;
334     prob->totComp += Nc;
335     if (feBd) {
336       ierr = PetscFEGetDimension(feBd, &Nb);CHKERRQ(ierr);
337       ierr = PetscFEGetNumComponents(feBd, &Nc);CHKERRQ(ierr);
338       ierr = PetscFEGetDefaultTabulation(feBd, &prob->basisBd[f], &prob->basisDerBd[f], NULL);CHKERRQ(ierr);
339       prob->totDimBd += Nb*Nc;
340     }
341   }
342   work = PetscMax(prob->totComp*dim, PetscSqr(NcMax*dim));
343   /* Allocate works space */
344   ierr = PetscMalloc5(prob->totComp,&prob->u,prob->totComp,&prob->u_t,prob->totComp*dim,&prob->u_x,dim,&prob->x,work,&prob->refSpaceDer);CHKERRQ(ierr);
345   ierr = PetscMalloc6(NqMax*NcMax,&prob->f0,NqMax*NcMax*dim,&prob->f1,NqMax*NcMax*NcMax,&prob->g0,NqMax*NcMax*NcMax*dim,&prob->g1,NqMax*NcMax*NcMax*dim,&prob->g2,NqMax*NcMax*NcMax*dim*dim,&prob->g3);CHKERRQ(ierr);
346   if (prob->ops->setup) {ierr = (*prob->ops->setup)(prob);CHKERRQ(ierr);}
347   prob->setup = PETSC_TRUE;
348   PetscFunctionReturn(0);
349 }
350 
351 #undef __FUNCT__
352 #define __FUNCT__ "PetscDSDestroyStructs_Static"
353 static PetscErrorCode PetscDSDestroyStructs_Static(PetscDS prob)
354 {
355   PetscErrorCode ierr;
356 
357   PetscFunctionBegin;
358   ierr = PetscFree4(prob->basis,prob->basisDer,prob->basisBd,prob->basisDerBd);CHKERRQ(ierr);
359   ierr = PetscFree5(prob->u,prob->u_t,prob->u_x,prob->x,prob->refSpaceDer);CHKERRQ(ierr);
360   ierr = PetscFree6(prob->f0,prob->f1,prob->g0,prob->g1,prob->g2,prob->g3);CHKERRQ(ierr);
361   PetscFunctionReturn(0);
362 }
363 
364 #undef __FUNCT__
365 #define __FUNCT__ "PetscDSEnlarge_Static"
366 static PetscErrorCode PetscDSEnlarge_Static(PetscDS prob, PetscInt NfNew)
367 {
368   PetscObject   *tmpd, *tmpdbd;
369   PointFunc     *tmpobj, *tmpf, *tmpg;
370   BdPointFunc   *tmpfbd, *tmpgbd;
371   RiemannFunc   *tmpr;
372   void         **tmpctx;
373   PetscInt       Nf = prob->Nf, f;
374   PetscErrorCode ierr;
375 
376   PetscFunctionBegin;
377   if (Nf >= NfNew) PetscFunctionReturn(0);
378   prob->setup = PETSC_FALSE;
379   ierr = PetscDSDestroyStructs_Static(prob);CHKERRQ(ierr);
380   ierr = PetscMalloc1(NfNew, &tmpd);CHKERRQ(ierr);
381   for (f = 0; f < Nf; ++f) tmpd[f] = prob->disc[f];
382   for (f = Nf; f < NfNew; ++f) {ierr = PetscContainerCreate(PetscObjectComm((PetscObject) prob), (PetscContainer *) &tmpd[f]);CHKERRQ(ierr);}
383   ierr = PetscFree(prob->disc);CHKERRQ(ierr);
384   prob->Nf   = NfNew;
385   prob->disc = tmpd;
386   ierr = PetscCalloc5(NfNew, &tmpobj, NfNew*2, &tmpf, NfNew*NfNew*4, &tmpg, NfNew, &tmpr, NfNew, &tmpctx);CHKERRQ(ierr);
387   for (f = 0; f < Nf; ++f) tmpobj[f] = prob->obj[f];
388   for (f = 0; f < Nf*2; ++f) tmpf[f] = prob->f[f];
389   for (f = 0; f < Nf*Nf*4; ++f) tmpg[f] = prob->g[f];
390   for (f = 0; f < Nf; ++f) tmpr[f] = prob->r[f];
391   for (f = 0; f < Nf; ++f) tmpctx[f] = prob->ctx[f];
392   for (f = Nf; f < NfNew; ++f) tmpobj[f] = NULL;
393   for (f = Nf*2; f < NfNew*2; ++f) tmpf[f] = NULL;
394   for (f = Nf*Nf*4; f < NfNew*NfNew*4; ++f) tmpg[f] = NULL;
395   for (f = Nf; f < NfNew; ++f) tmpr[f] = NULL;
396   for (f = Nf; f < NfNew; ++f) tmpctx[f] = NULL;
397   ierr = PetscFree5(prob->obj, prob->f, prob->g, prob->r, prob->ctx);CHKERRQ(ierr);
398   prob->obj = tmpobj;
399   prob->f   = tmpf;
400   prob->g   = tmpg;
401   prob->r   = tmpr;
402   prob->ctx = tmpctx;
403   ierr = PetscMalloc1(NfNew, &tmpdbd);CHKERRQ(ierr);
404   for (f = 0; f < Nf; ++f) tmpdbd[f] = prob->discBd[f];
405   for (f = Nf; f < NfNew; ++f) tmpdbd[f] = NULL;
406   ierr = PetscFree(prob->discBd);CHKERRQ(ierr);
407   prob->discBd = tmpdbd;
408   ierr = PetscCalloc2(NfNew*2, &tmpfbd, NfNew*NfNew*4, &tmpgbd);CHKERRQ(ierr);
409   for (f = 0; f < Nf*2; ++f) tmpfbd[f] = prob->fBd[f];
410   for (f = 0; f < Nf*Nf*4; ++f) tmpgbd[f] = prob->gBd[f];
411   for (f = Nf*2; f < NfNew*2; ++f) tmpfbd[f] = NULL;
412   for (f = Nf*Nf*4; f < NfNew*NfNew*4; ++f) tmpgbd[f] = NULL;
413   ierr = PetscFree2(prob->fBd, prob->gBd);CHKERRQ(ierr);
414   prob->fBd = tmpfbd;
415   prob->gBd = tmpgbd;
416   PetscFunctionReturn(0);
417 }
418 
419 #undef __FUNCT__
420 #define __FUNCT__ "PetscDSDestroy"
421 /*@
422   PetscDSDestroy - Destroys a PetscDS object
423 
424   Collective on PetscDS
425 
426   Input Parameter:
427 . prob - the PetscDS object to destroy
428 
429   Level: developer
430 
431 .seealso PetscDSView()
432 @*/
433 PetscErrorCode PetscDSDestroy(PetscDS *prob)
434 {
435   PetscInt       f;
436   PetscErrorCode ierr;
437 
438   PetscFunctionBegin;
439   if (!*prob) PetscFunctionReturn(0);
440   PetscValidHeaderSpecific((*prob), PETSCDS_CLASSID, 1);
441 
442   if (--((PetscObject)(*prob))->refct > 0) {*prob = 0; PetscFunctionReturn(0);}
443   ((PetscObject) (*prob))->refct = 0;
444   ierr = PetscDSDestroyStructs_Static(*prob);CHKERRQ(ierr);
445   for (f = 0; f < (*prob)->Nf; ++f) {
446     ierr = PetscObjectDereference((*prob)->disc[f]);CHKERRQ(ierr);
447     ierr = PetscObjectDereference((*prob)->discBd[f]);CHKERRQ(ierr);
448   }
449   ierr = PetscFree((*prob)->disc);CHKERRQ(ierr);
450   ierr = PetscFree((*prob)->discBd);CHKERRQ(ierr);
451   ierr = PetscFree5((*prob)->obj,(*prob)->f,(*prob)->g,(*prob)->r,(*prob)->ctx);CHKERRQ(ierr);
452   ierr = PetscFree2((*prob)->fBd,(*prob)->gBd);CHKERRQ(ierr);
453   if ((*prob)->ops->destroy) {ierr = (*(*prob)->ops->destroy)(*prob);CHKERRQ(ierr);}
454   ierr = PetscHeaderDestroy(prob);CHKERRQ(ierr);
455   PetscFunctionReturn(0);
456 }
457 
458 #undef __FUNCT__
459 #define __FUNCT__ "PetscDSCreate"
460 /*@
461   PetscDSCreate - Creates an empty PetscDS object. The type can then be set with PetscDSSetType().
462 
463   Collective on MPI_Comm
464 
465   Input Parameter:
466 . comm - The communicator for the PetscDS object
467 
468   Output Parameter:
469 . prob - The PetscDS object
470 
471   Level: beginner
472 
473 .seealso: PetscDSSetType(), PETSCDSBASIC
474 @*/
475 PetscErrorCode PetscDSCreate(MPI_Comm comm, PetscDS *prob)
476 {
477   PetscDS   p;
478   PetscErrorCode ierr;
479 
480   PetscFunctionBegin;
481   PetscValidPointer(prob, 2);
482   *prob  = NULL;
483   ierr = PetscDSInitializePackage();CHKERRQ(ierr);
484 
485   ierr = PetscHeaderCreate(p, _p_PetscDS, struct _PetscDSOps, PETSCDS_CLASSID, "PetscDS", "Discrete System", "PetscDS", comm, PetscDSDestroy, PetscDSView);CHKERRQ(ierr);
486   ierr = PetscMemzero(p->ops, sizeof(struct _PetscDSOps));CHKERRQ(ierr);
487 
488   p->Nf    = 0;
489   p->setup = PETSC_FALSE;
490 
491   *prob = p;
492   PetscFunctionReturn(0);
493 }
494 
495 #undef __FUNCT__
496 #define __FUNCT__ "PetscDSGetNumFields"
497 /*@
498   PetscDSGetNumFields - Returns the number of fields in the DS
499 
500   Not collective
501 
502   Input Parameter:
503 . prob - The PetscDS object
504 
505   Output Parameter:
506 . Nf - The number of fields
507 
508   Level: beginner
509 
510 .seealso: PetscDSGetSpatialDimension(), PetscDSCreate()
511 @*/
512 PetscErrorCode PetscDSGetNumFields(PetscDS prob, PetscInt *Nf)
513 {
514   PetscFunctionBegin;
515   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
516   PetscValidPointer(Nf, 2);
517   *Nf = prob->Nf;
518   PetscFunctionReturn(0);
519 }
520 
521 #undef __FUNCT__
522 #define __FUNCT__ "PetscDSGetSpatialDimension"
523 /*@
524   PetscDSGetSpatialDimension - Returns the spatial dimension of the DS
525 
526   Not collective
527 
528   Input Parameter:
529 . prob - The PetscDS object
530 
531   Output Parameter:
532 . dim - The spatial dimension
533 
534   Level: beginner
535 
536 .seealso: PetscDSGetNumFields(), PetscDSCreate()
537 @*/
538 PetscErrorCode PetscDSGetSpatialDimension(PetscDS prob, PetscInt *dim)
539 {
540   PetscErrorCode ierr;
541 
542   PetscFunctionBegin;
543   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
544   PetscValidPointer(dim, 2);
545   *dim = 0;
546   if (prob->Nf) {
547     PetscObject  obj;
548     PetscClassId id;
549 
550     ierr = PetscDSGetDiscretization(prob, 0, &obj);CHKERRQ(ierr);
551     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
552     if (id == PETSCFE_CLASSID)      {ierr = PetscFEGetSpatialDimension((PetscFE) obj, dim);CHKERRQ(ierr);}
553     else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetSpatialDimension((PetscFV) obj, dim);CHKERRQ(ierr);}
554     else SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", 0);
555   }
556   PetscFunctionReturn(0);
557 }
558 
559 #undef __FUNCT__
560 #define __FUNCT__ "PetscDSGetTotalDimension"
561 /*@
562   PetscDSGetTotalDimension - Returns the total size of the approximation space for this system
563 
564   Not collective
565 
566   Input Parameter:
567 . prob - The PetscDS object
568 
569   Output Parameter:
570 . dim - The total problem dimension
571 
572   Level: beginner
573 
574 .seealso: PetscDSGetNumFields(), PetscDSCreate()
575 @*/
576 PetscErrorCode PetscDSGetTotalDimension(PetscDS prob, PetscInt *dim)
577 {
578   PetscErrorCode ierr;
579 
580   PetscFunctionBegin;
581   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
582   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
583   PetscValidPointer(dim, 2);
584   *dim = prob->totDim;
585   PetscFunctionReturn(0);
586 }
587 
588 #undef __FUNCT__
589 #define __FUNCT__ "PetscDSGetTotalBdDimension"
590 /*@
591   PetscDSGetTotalBdDimension - Returns the total size of the boundary approximation space for this system
592 
593   Not collective
594 
595   Input Parameter:
596 . prob - The PetscDS object
597 
598   Output Parameter:
599 . dim - The total boundary problem dimension
600 
601   Level: beginner
602 
603 .seealso: PetscDSGetNumFields(), PetscDSCreate()
604 @*/
605 PetscErrorCode PetscDSGetTotalBdDimension(PetscDS prob, PetscInt *dim)
606 {
607   PetscErrorCode ierr;
608 
609   PetscFunctionBegin;
610   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
611   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
612   PetscValidPointer(dim, 2);
613   *dim = prob->totDimBd;
614   PetscFunctionReturn(0);
615 }
616 
617 #undef __FUNCT__
618 #define __FUNCT__ "PetscDSGetTotalComponents"
619 /*@
620   PetscDSGetTotalComponents - Returns the total number of components in this system
621 
622   Not collective
623 
624   Input Parameter:
625 . prob - The PetscDS object
626 
627   Output Parameter:
628 . dim - The total number of components
629 
630   Level: beginner
631 
632 .seealso: PetscDSGetNumFields(), PetscDSCreate()
633 @*/
634 PetscErrorCode PetscDSGetTotalComponents(PetscDS prob, PetscInt *Nc)
635 {
636   PetscErrorCode ierr;
637 
638   PetscFunctionBegin;
639   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
640   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
641   PetscValidPointer(Nc, 2);
642   *Nc = prob->totComp;
643   PetscFunctionReturn(0);
644 }
645 
646 #undef __FUNCT__
647 #define __FUNCT__ "PetscDSGetDiscretization"
648 /*@
649   PetscDSGetDiscretization - Returns the discretization object for the given field
650 
651   Not collective
652 
653   Input Parameters:
654 + prob - The PetscDS object
655 - f - The field number
656 
657   Output Parameter:
658 . disc - The discretization object
659 
660   Level: beginner
661 
662 .seealso: PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
663 @*/
664 PetscErrorCode PetscDSGetDiscretization(PetscDS prob, PetscInt f, PetscObject *disc)
665 {
666   PetscFunctionBegin;
667   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
668   PetscValidPointer(disc, 3);
669   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);
670   *disc = prob->disc[f];
671   PetscFunctionReturn(0);
672 }
673 
674 #undef __FUNCT__
675 #define __FUNCT__ "PetscDSGetBdDiscretization"
676 /*@
677   PetscDSGetBdDiscretization - Returns the boundary discretization object for the given field
678 
679   Not collective
680 
681   Input Parameters:
682 + prob - The PetscDS object
683 - f - The field number
684 
685   Output Parameter:
686 . disc - The boundary discretization object
687 
688   Level: beginner
689 
690 .seealso: PetscDSSetBdDiscretization(), PetscDSAddBdDiscretization(), PetscDSGetDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
691 @*/
692 PetscErrorCode PetscDSGetBdDiscretization(PetscDS prob, PetscInt f, PetscObject *disc)
693 {
694   PetscFunctionBegin;
695   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
696   PetscValidPointer(disc, 3);
697   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);
698   *disc = prob->discBd[f];
699   PetscFunctionReturn(0);
700 }
701 
702 #undef __FUNCT__
703 #define __FUNCT__ "PetscDSSetDiscretization"
704 /*@
705   PetscDSSetDiscretization - Sets the discretization object for the given field
706 
707   Not collective
708 
709   Input Parameters:
710 + prob - The PetscDS object
711 . f - The field number
712 - disc - The discretization object
713 
714   Level: beginner
715 
716 .seealso: PetscDSGetDiscretization(), PetscDSAddDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
717 @*/
718 PetscErrorCode PetscDSSetDiscretization(PetscDS prob, PetscInt f, PetscObject disc)
719 {
720   PetscErrorCode ierr;
721 
722   PetscFunctionBegin;
723   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
724   PetscValidPointer(disc, 3);
725   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
726   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
727   if (prob->disc[f]) {ierr = PetscObjectDereference(prob->disc[f]);CHKERRQ(ierr);}
728   prob->disc[f] = disc;
729   ierr = PetscObjectReference(disc);CHKERRQ(ierr);
730   PetscFunctionReturn(0);
731 }
732 
733 #undef __FUNCT__
734 #define __FUNCT__ "PetscDSSetBdDiscretization"
735 /*@
736   PetscDSSetBdDiscretization - Sets the boundary discretization object for the given field
737 
738   Not collective
739 
740   Input Parameters:
741 + prob - The PetscDS object
742 . f - The field number
743 - disc - The boundary discretization object
744 
745   Level: beginner
746 
747 .seealso: PetscDSGetBdDiscretization(), PetscDSAddBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
748 @*/
749 PetscErrorCode PetscDSSetBdDiscretization(PetscDS prob, PetscInt f, PetscObject disc)
750 {
751   PetscErrorCode ierr;
752 
753   PetscFunctionBegin;
754   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
755   if (disc) PetscValidPointer(disc, 3);
756   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
757   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
758   if (prob->discBd[f]) {ierr = PetscObjectDereference(prob->discBd[f]);CHKERRQ(ierr);}
759   prob->discBd[f] = disc;
760   ierr = PetscObjectReference(disc);CHKERRQ(ierr);
761   PetscFunctionReturn(0);
762 }
763 
764 #undef __FUNCT__
765 #define __FUNCT__ "PetscDSAddDiscretization"
766 /*@
767   PetscDSAddDiscretization - Adds a discretization object
768 
769   Not collective
770 
771   Input Parameters:
772 + prob - The PetscDS object
773 - disc - The boundary discretization object
774 
775   Level: beginner
776 
777 .seealso: PetscDSGetDiscretization(), PetscDSSetDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
778 @*/
779 PetscErrorCode PetscDSAddDiscretization(PetscDS prob, PetscObject disc)
780 {
781   PetscErrorCode ierr;
782 
783   PetscFunctionBegin;
784   ierr = PetscDSSetDiscretization(prob, prob->Nf, disc);CHKERRQ(ierr);
785   PetscFunctionReturn(0);
786 }
787 
788 #undef __FUNCT__
789 #define __FUNCT__ "PetscDSAddBdDiscretization"
790 /*@
791   PetscDSAddBdDiscretization - Adds a boundary discretization object
792 
793   Not collective
794 
795   Input Parameters:
796 + prob - The PetscDS object
797 - disc - The boundary discretization object
798 
799   Level: beginner
800 
801 .seealso: PetscDSGetBdDiscretization(), PetscDSSetBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
802 @*/
803 PetscErrorCode PetscDSAddBdDiscretization(PetscDS prob, PetscObject disc)
804 {
805   PetscErrorCode ierr;
806 
807   PetscFunctionBegin;
808   ierr = PetscDSSetBdDiscretization(prob, prob->Nf, disc);CHKERRQ(ierr);
809   PetscFunctionReturn(0);
810 }
811 
812 #undef __FUNCT__
813 #define __FUNCT__ "PetscDSGetObjective"
814 PetscErrorCode PetscDSGetObjective(PetscDS prob, PetscInt f,
815                                         void (**obj)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar obj[]))
816 {
817   PetscFunctionBegin;
818   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
819   PetscValidPointer(obj, 2);
820   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);
821   *obj = prob->obj[f];
822   PetscFunctionReturn(0);
823 }
824 
825 #undef __FUNCT__
826 #define __FUNCT__ "PetscDSSetObjective"
827 PetscErrorCode PetscDSSetObjective(PetscDS prob, PetscInt f,
828                                         void (*obj)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar obj[]))
829 {
830   PetscErrorCode ierr;
831 
832   PetscFunctionBegin;
833   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
834   PetscValidFunction(obj, 2);
835   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
836   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
837   prob->obj[f] = obj;
838   PetscFunctionReturn(0);
839 }
840 
841 #undef __FUNCT__
842 #define __FUNCT__ "PetscDSGetResidual"
843 PetscErrorCode PetscDSGetResidual(PetscDS prob, PetscInt f,
844                                        void (**f0)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar f0[]),
845                                        void (**f1)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar f1[]))
846 {
847   PetscFunctionBegin;
848   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
849   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);
850   if (f0) {PetscValidPointer(f0, 3); *f0 = prob->f[f*2+0];}
851   if (f1) {PetscValidPointer(f1, 4); *f1 = prob->f[f*2+1];}
852   PetscFunctionReturn(0);
853 }
854 
855 #undef __FUNCT__
856 #define __FUNCT__ "PetscDSSetResidual"
857 PetscErrorCode PetscDSSetResidual(PetscDS prob, PetscInt f,
858                                        void (*f0)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar f0[]),
859                                        void (*f1)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar f1[]))
860 {
861   PetscErrorCode ierr;
862 
863   PetscFunctionBegin;
864   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
865   PetscValidFunction(f0, 3);
866   PetscValidFunction(f1, 4);
867   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
868   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
869   prob->f[f*2+0] = f0;
870   prob->f[f*2+1] = f1;
871   PetscFunctionReturn(0);
872 }
873 
874 #undef __FUNCT__
875 #define __FUNCT__ "PetscDSGetJacobian"
876 PetscErrorCode PetscDSGetJacobian(PetscDS prob, PetscInt f, PetscInt g,
877                                        void (**g0)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar g0[]),
878                                        void (**g1)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar g1[]),
879                                        void (**g2)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar g2[]),
880                                        void (**g3)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar g3[]))
881 {
882   PetscFunctionBegin;
883   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
884   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);
885   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);
886   if (g0) {PetscValidPointer(g0, 4); *g0 = prob->g[(f*prob->Nf + g)*4+0];}
887   if (g1) {PetscValidPointer(g1, 5); *g1 = prob->g[(f*prob->Nf + g)*4+1];}
888   if (g2) {PetscValidPointer(g2, 6); *g2 = prob->g[(f*prob->Nf + g)*4+2];}
889   if (g3) {PetscValidPointer(g3, 7); *g3 = prob->g[(f*prob->Nf + g)*4+3];}
890   PetscFunctionReturn(0);
891 }
892 
893 #undef __FUNCT__
894 #define __FUNCT__ "PetscDSSetJacobian"
895 PetscErrorCode PetscDSSetJacobian(PetscDS prob, PetscInt f, PetscInt g,
896                                        void (*g0)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar g0[]),
897                                        void (*g1)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar g1[]),
898                                        void (*g2)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar g2[]),
899                                        void (*g3)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], PetscScalar g3[]))
900 {
901   PetscErrorCode ierr;
902 
903   PetscFunctionBegin;
904   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
905   if (g0) PetscValidFunction(g0, 4);
906   if (g1) PetscValidFunction(g1, 5);
907   if (g2) PetscValidFunction(g2, 6);
908   if (g3) PetscValidFunction(g3, 7);
909   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
910   if (g < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", g);
911   ierr = PetscDSEnlarge_Static(prob, PetscMax(f, g)+1);CHKERRQ(ierr);
912   prob->g[(f*prob->Nf + g)*4+0] = g0;
913   prob->g[(f*prob->Nf + g)*4+1] = g1;
914   prob->g[(f*prob->Nf + g)*4+2] = g2;
915   prob->g[(f*prob->Nf + g)*4+3] = g3;
916   PetscFunctionReturn(0);
917 }
918 
919 #undef __FUNCT__
920 #define __FUNCT__ "PetscDSGetRiemannSolver"
921 /*@C
922   PetscDSGetRiemannSolver - Returns the Riemann solver for the given field
923 
924   Not collective
925 
926   Input Arguments:
927 + prob - The PetscDS object
928 - f    - The field number
929 
930   Output Argument:
931 . r    - Riemann solver
932 
933   Calling sequence for r:
934 
935 $ r(const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx)
936 
937 + x    - The coordinates at a point on the interface
938 . n    - The normal vector to the interface
939 . uL   - The state vector to the left of the interface
940 . uR   - The state vector to the right of the interface
941 . flux - output array of flux through the interface
942 - ctx  - optional user context
943 
944   Level: intermediate
945 
946 .seealso: PetscDSSetRiemannSolver()
947 @*/
948 PetscErrorCode PetscDSGetRiemannSolver(PetscDS prob, PetscInt f,
949                                        void (**r)(const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx))
950 {
951   PetscFunctionBegin;
952   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
953   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);
954   PetscValidPointer(r, 3);
955   *r = prob->r[f];
956   PetscFunctionReturn(0);
957 }
958 
959 #undef __FUNCT__
960 #define __FUNCT__ "PetscDSSetRiemannSolver"
961 /*@C
962   PetscDSSetRiemannSolver - Sets the Riemann solver for the given field
963 
964   Not collective
965 
966   Input Arguments:
967 + prob - The PetscDS object
968 . f    - The field number
969 - r    - Riemann solver
970 
971   Calling sequence for r:
972 
973 $ r(const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx)
974 
975 + x    - The coordinates at a point on the interface
976 . n    - The normal vector to the interface
977 . uL   - The state vector to the left of the interface
978 . uR   - The state vector to the right of the interface
979 . flux - output array of flux through the interface
980 - ctx  - optional user context
981 
982   Level: intermediate
983 
984 .seealso: PetscDSGetRiemannSolver()
985 @*/
986 PetscErrorCode PetscDSSetRiemannSolver(PetscDS prob, PetscInt f,
987                                        void (*r)(const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx))
988 {
989   PetscErrorCode ierr;
990 
991   PetscFunctionBegin;
992   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
993   PetscValidFunction(r, 3);
994   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
995   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
996   prob->r[f] = r;
997   PetscFunctionReturn(0);
998 }
999 
1000 #undef __FUNCT__
1001 #define __FUNCT__ "PetscDSGetContext"
1002 PetscErrorCode PetscDSGetContext(PetscDS prob, PetscInt f, void **ctx)
1003 {
1004   PetscFunctionBegin;
1005   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1006   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);
1007   PetscValidPointer(ctx, 3);
1008   *ctx = prob->ctx[f];
1009   PetscFunctionReturn(0);
1010 }
1011 
1012 #undef __FUNCT__
1013 #define __FUNCT__ "PetscDSSetContext"
1014 PetscErrorCode PetscDSSetContext(PetscDS prob, PetscInt f, void *ctx)
1015 {
1016   PetscErrorCode ierr;
1017 
1018   PetscFunctionBegin;
1019   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1020   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
1021   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
1022   prob->ctx[f] = ctx;
1023   PetscFunctionReturn(0);
1024 }
1025 
1026 #undef __FUNCT__
1027 #define __FUNCT__ "PetscDSGetBdResidual"
1028 PetscErrorCode PetscDSGetBdResidual(PetscDS prob, PetscInt f,
1029                                          void (**f0)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar f0[]),
1030                                          void (**f1)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar f1[]))
1031 {
1032   PetscFunctionBegin;
1033   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1034   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);
1035   if (f0) {PetscValidPointer(f0, 3); *f0 = prob->fBd[f*2+0];}
1036   if (f1) {PetscValidPointer(f1, 4); *f1 = prob->fBd[f*2+1];}
1037   PetscFunctionReturn(0);
1038 }
1039 
1040 #undef __FUNCT__
1041 #define __FUNCT__ "PetscDSSetBdResidual"
1042 PetscErrorCode PetscDSSetBdResidual(PetscDS prob, PetscInt f,
1043                                          void (*f0)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar f0[]),
1044                                          void (*f1)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar f1[]))
1045 {
1046   PetscErrorCode ierr;
1047 
1048   PetscFunctionBegin;
1049   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1050   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
1051   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
1052   if (f0) {PetscValidFunction(f0, 3); prob->fBd[f*2+0] = f0;}
1053   if (f1) {PetscValidFunction(f1, 4); prob->fBd[f*2+1] = f1;}
1054   PetscFunctionReturn(0);
1055 }
1056 
1057 #undef __FUNCT__
1058 #define __FUNCT__ "PetscDSGetBdJacobian"
1059 PetscErrorCode PetscDSGetBdJacobian(PetscDS prob, PetscInt f, PetscInt g,
1060                                          void (**g0)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar g0[]),
1061                                          void (**g1)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar g1[]),
1062                                          void (**g2)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar g2[]),
1063                                          void (**g3)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar g3[]))
1064 {
1065   PetscFunctionBegin;
1066   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1067   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);
1068   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);
1069   if (g0) {PetscValidPointer(g0, 4); *g0 = prob->gBd[(f*prob->Nf + g)*4+0];}
1070   if (g1) {PetscValidPointer(g1, 5); *g1 = prob->gBd[(f*prob->Nf + g)*4+1];}
1071   if (g2) {PetscValidPointer(g2, 6); *g2 = prob->gBd[(f*prob->Nf + g)*4+2];}
1072   if (g3) {PetscValidPointer(g3, 7); *g3 = prob->gBd[(f*prob->Nf + g)*4+3];}
1073   PetscFunctionReturn(0);
1074 }
1075 
1076 #undef __FUNCT__
1077 #define __FUNCT__ "PetscDSSetBdJacobian"
1078 PetscErrorCode PetscDSSetBdJacobian(PetscDS prob, PetscInt f, PetscInt g,
1079                                          void (*g0)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar g0[]),
1080                                          void (*g1)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar g1[]),
1081                                          void (*g2)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar g2[]),
1082                                          void (*g3)(const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], const PetscReal x[], const PetscReal n[], PetscScalar g3[]))
1083 {
1084   PetscErrorCode ierr;
1085 
1086   PetscFunctionBegin;
1087   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1088   if (g0) PetscValidFunction(g0, 4);
1089   if (g1) PetscValidFunction(g1, 5);
1090   if (g2) PetscValidFunction(g2, 6);
1091   if (g3) PetscValidFunction(g3, 7);
1092   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
1093   if (g < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", g);
1094   ierr = PetscDSEnlarge_Static(prob, PetscMax(f, g)+1);CHKERRQ(ierr);
1095   prob->gBd[(f*prob->Nf + g)*4+0] = g0;
1096   prob->gBd[(f*prob->Nf + g)*4+1] = g1;
1097   prob->gBd[(f*prob->Nf + g)*4+2] = g2;
1098   prob->gBd[(f*prob->Nf + g)*4+3] = g3;
1099   PetscFunctionReturn(0);
1100 }
1101 
1102 #undef __FUNCT__
1103 #define __FUNCT__ "PetscDSGetFieldOffset"
1104 /*@
1105   PetscDSGetFieldOffset - Returns the offset of the given field in the full space basis
1106 
1107   Not collective
1108 
1109   Input Parameters:
1110 + prob - The PetscDS object
1111 - f - The field number
1112 
1113   Output Parameter:
1114 . off - The offset
1115 
1116   Level: beginner
1117 
1118 .seealso: PetscDSGetBdFieldOffset(), PetscDSGetNumFields(), PetscDSCreate()
1119 @*/
1120 PetscErrorCode PetscDSGetFieldOffset(PetscDS prob, PetscInt f, PetscInt *off)
1121 {
1122   PetscInt       g;
1123   PetscErrorCode ierr;
1124 
1125   PetscFunctionBegin;
1126   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1127   PetscValidPointer(off, 3);
1128   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);
1129   *off = 0;
1130   for (g = 0; g < f; ++g) {
1131     PetscFE  fe = (PetscFE) prob->disc[g];
1132     PetscInt Nb, Nc;
1133 
1134     ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
1135     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
1136     *off += Nb*Nc;
1137   }
1138   PetscFunctionReturn(0);
1139 }
1140 
1141 #undef __FUNCT__
1142 #define __FUNCT__ "PetscDSGetBdFieldOffset"
1143 /*@
1144   PetscDSGetBdFieldOffset - Returns the offset of the given field in the full space boundary basis
1145 
1146   Not collective
1147 
1148   Input Parameters:
1149 + prob - The PetscDS object
1150 - f - The field number
1151 
1152   Output Parameter:
1153 . off - The boundary offset
1154 
1155   Level: beginner
1156 
1157 .seealso: PetscDSGetFieldOffset(), PetscDSGetNumFields(), PetscDSCreate()
1158 @*/
1159 PetscErrorCode PetscDSGetBdFieldOffset(PetscDS prob, PetscInt f, PetscInt *off)
1160 {
1161   PetscInt       g;
1162   PetscErrorCode ierr;
1163 
1164   PetscFunctionBegin;
1165   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1166   PetscValidPointer(off, 3);
1167   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);
1168   *off = 0;
1169   for (g = 0; g < f; ++g) {
1170     PetscFE  fe = (PetscFE) prob->discBd[g];
1171     PetscInt Nb, Nc;
1172 
1173     ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
1174     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
1175     *off += Nb*Nc;
1176   }
1177   PetscFunctionReturn(0);
1178 }
1179 
1180 #undef __FUNCT__
1181 #define __FUNCT__ "PetscDSGetTabulation"
1182 /*@C
1183   PetscDSGetTabulation - Return the basis tabulation at quadrature points for the volume discretization
1184 
1185   Not collective
1186 
1187   Input Parameter:
1188 . prob - The PetscDS object
1189 
1190   Output Parameters:
1191 + basis - The basis function tabulation at quadrature points
1192 - basisDer - The basis function derivative tabulation at quadrature points
1193 
1194   Level: intermediate
1195 
1196 .seealso: PetscDSGetBdTabulation(), PetscDSCreate()
1197 @*/
1198 PetscErrorCode PetscDSGetTabulation(PetscDS prob, PetscReal ***basis, PetscReal ***basisDer)
1199 {
1200   PetscErrorCode ierr;
1201 
1202   PetscFunctionBegin;
1203   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1204   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
1205   if (basis)    {PetscValidPointer(basis, 2);    *basis    = prob->basis;}
1206   if (basisDer) {PetscValidPointer(basisDer, 3); *basisDer = prob->basisDer;}
1207   PetscFunctionReturn(0);
1208 }
1209 
1210 #undef __FUNCT__
1211 #define __FUNCT__ "PetscDSGetBdTabulation"
1212 /*@C
1213   PetscDSGetBdTabulation - Return the basis tabulation at quadrature points for the boundary discretization
1214 
1215   Not collective
1216 
1217   Input Parameter:
1218 . prob - The PetscDS object
1219 
1220   Output Parameters:
1221 + basis - The basis function tabulation at quadrature points
1222 - basisDer - The basis function derivative tabulation at quadrature points
1223 
1224   Level: intermediate
1225 
1226 .seealso: PetscDSGetTabulation(), PetscDSCreate()
1227 @*/
1228 PetscErrorCode PetscDSGetBdTabulation(PetscDS prob, PetscReal ***basis, PetscReal ***basisDer)
1229 {
1230   PetscErrorCode ierr;
1231 
1232   PetscFunctionBegin;
1233   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1234   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
1235   if (basis)    {PetscValidPointer(basis, 2);    *basis    = prob->basisBd;}
1236   if (basisDer) {PetscValidPointer(basisDer, 3); *basisDer = prob->basisDerBd;}
1237   PetscFunctionReturn(0);
1238 }
1239 
1240 #undef __FUNCT__
1241 #define __FUNCT__ "PetscDSGetEvaluationArrays"
1242 PetscErrorCode PetscDSGetEvaluationArrays(PetscDS prob, PetscScalar **u, PetscScalar **u_t, PetscScalar **u_x)
1243 {
1244   PetscErrorCode ierr;
1245 
1246   PetscFunctionBegin;
1247   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1248   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
1249   if (u)   {PetscValidPointer(u, 2);   *u   = prob->u;}
1250   if (u_t) {PetscValidPointer(u_t, 3); *u_t = prob->u_t;}
1251   if (u_x) {PetscValidPointer(u_x, 4); *u_x = prob->u_x;}
1252   PetscFunctionReturn(0);
1253 }
1254 
1255 #undef __FUNCT__
1256 #define __FUNCT__ "PetscDSGetWeakFormArrays"
1257 PetscErrorCode PetscDSGetWeakFormArrays(PetscDS prob, PetscScalar **f0, PetscScalar **f1, PetscScalar **g0, PetscScalar **g1, PetscScalar **g2, PetscScalar **g3)
1258 {
1259   PetscErrorCode ierr;
1260 
1261   PetscFunctionBegin;
1262   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1263   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
1264   if (f0) {PetscValidPointer(f0, 2); *f0 = prob->f0;}
1265   if (f1) {PetscValidPointer(f1, 3); *f1 = prob->f1;}
1266   if (g0) {PetscValidPointer(g0, 4); *g0 = prob->g0;}
1267   if (g1) {PetscValidPointer(g1, 5); *g1 = prob->g1;}
1268   if (g2) {PetscValidPointer(g2, 6); *g2 = prob->g2;}
1269   if (g3) {PetscValidPointer(g3, 7); *g3 = prob->g3;}
1270   PetscFunctionReturn(0);
1271 }
1272 
1273 #undef __FUNCT__
1274 #define __FUNCT__ "PetscDSGetRefCoordArrays"
1275 PetscErrorCode PetscDSGetRefCoordArrays(PetscDS prob, PetscReal **x, PetscScalar **refSpaceDer)
1276 {
1277   PetscErrorCode ierr;
1278 
1279   PetscFunctionBegin;
1280   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1281   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
1282   if (x)           {PetscValidPointer(x, 2);           *x           = prob->x;}
1283   if (refSpaceDer) {PetscValidPointer(refSpaceDer, 3); *refSpaceDer = prob->refSpaceDer;}
1284   PetscFunctionReturn(0);
1285 }
1286 
1287 #undef __FUNCT__
1288 #define __FUNCT__ "PetscDSDestroy_Basic"
1289 static PetscErrorCode PetscDSDestroy_Basic(PetscDS prob)
1290 {
1291   PetscFunctionBegin;
1292   PetscFunctionReturn(0);
1293 }
1294 
1295 #undef __FUNCT__
1296 #define __FUNCT__ "PetscDSInitialize_Basic"
1297 static PetscErrorCode PetscDSInitialize_Basic(PetscDS prob)
1298 {
1299   PetscFunctionBegin;
1300   prob->ops->setfromoptions = NULL;
1301   prob->ops->setup          = NULL;
1302   prob->ops->view           = NULL;
1303   prob->ops->destroy        = PetscDSDestroy_Basic;
1304   PetscFunctionReturn(0);
1305 }
1306 
1307 /*MC
1308   PETSCDSBASIC = "basic" - A discrete system with pointwise residual and boundary residual functions
1309 
1310   Level: intermediate
1311 
1312 .seealso: PetscDSType, PetscDSCreate(), PetscDSSetType()
1313 M*/
1314 
1315 #undef __FUNCT__
1316 #define __FUNCT__ "PetscDSCreate_Basic"
1317 PETSC_EXTERN PetscErrorCode PetscDSCreate_Basic(PetscDS prob)
1318 {
1319   PetscDS_Basic *b;
1320   PetscErrorCode      ierr;
1321 
1322   PetscFunctionBegin;
1323   PetscValidHeaderSpecific(prob, PETSCSPACE_CLASSID, 1);
1324   ierr       = PetscNewLog(prob, &b);CHKERRQ(ierr);
1325   prob->data = b;
1326 
1327   ierr = PetscDSInitialize_Basic(prob);CHKERRQ(ierr);
1328   PetscFunctionReturn(0);
1329 }
1330