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