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