xref: /petsc/src/dm/dt/interface/dtds.c (revision 8cc346eb00edc34721c4641da4e067bd7f773099)
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     ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr);
160     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
161       if (id == PETSCFE_CLASSID)      {ierr = PetscFEView((PetscFE) obj, viewer);CHKERRQ(ierr);}
162       else if (id == PETSCFV_CLASSID) {ierr = PetscFVView((PetscFV) obj, viewer);CHKERRQ(ierr);}
163     }
164   }
165   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
166   PetscFunctionReturn(0);
167 }
168 
169 #undef __FUNCT__
170 #define __FUNCT__ "PetscDSView"
171 /*@C
172   PetscDSView - Views a PetscDS
173 
174   Collective on PetscDS
175 
176   Input Parameter:
177 + prob - the PetscDS object to view
178 - v  - the viewer
179 
180   Level: developer
181 
182 .seealso PetscDSDestroy()
183 @*/
184 PetscErrorCode PetscDSView(PetscDS prob, PetscViewer v)
185 {
186   PetscBool      iascii;
187   PetscErrorCode ierr;
188 
189   PetscFunctionBegin;
190   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
191   if (!v) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) prob), &v);CHKERRQ(ierr);}
192   else    {PetscValidHeaderSpecific(v, PETSC_VIEWER_CLASSID, 2);}
193   ierr = PetscObjectTypeCompare((PetscObject) v, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
194   if (iascii) {ierr = PetscDSView_Ascii(prob, v);CHKERRQ(ierr);}
195   if (prob->ops->view) {ierr = (*prob->ops->view)(prob, v);CHKERRQ(ierr);}
196   PetscFunctionReturn(0);
197 }
198 
199 #undef __FUNCT__
200 #define __FUNCT__ "PetscDSViewFromOptions"
201 /*
202   PetscDSViewFromOptions - Processes command line options to determine if/how a PetscDS is to be viewed.
203 
204   Collective on PetscDS
205 
206   Input Parameters:
207 + prob   - the PetscDS
208 . prefix - prefix to use for viewing, or NULL to use prefix of 'rnd'
209 - optionname - option to activate viewing
210 
211   Level: intermediate
212 
213 .keywords: PetscDS, view, options, database
214 .seealso: VecViewFromOptions(), MatViewFromOptions()
215 */
216 PetscErrorCode PetscDSViewFromOptions(PetscDS prob, const char prefix[], const char optionname[])
217 {
218   PetscViewer       viewer;
219   PetscViewerFormat format;
220   PetscBool         flg;
221   PetscErrorCode    ierr;
222 
223   PetscFunctionBegin;
224   if (prefix) {ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject) prob), prefix, optionname, &viewer, &format, &flg);CHKERRQ(ierr);}
225   else        {ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject) prob), ((PetscObject) prob)->prefix, optionname, &viewer, &format, &flg);CHKERRQ(ierr);}
226   if (flg) {
227     ierr = PetscViewerPushFormat(viewer, format);CHKERRQ(ierr);
228     ierr = PetscDSView(prob, viewer);CHKERRQ(ierr);
229     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
230     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
231   }
232   PetscFunctionReturn(0);
233 }
234 
235 #undef __FUNCT__
236 #define __FUNCT__ "PetscDSSetFromOptions"
237 /*@
238   PetscDSSetFromOptions - sets parameters in a PetscDS from the options database
239 
240   Collective on PetscDS
241 
242   Input Parameter:
243 . prob - the PetscDS object to set options for
244 
245   Options Database:
246 
247   Level: developer
248 
249 .seealso PetscDSView()
250 @*/
251 PetscErrorCode PetscDSSetFromOptions(PetscDS prob)
252 {
253   const char    *defaultType;
254   char           name[256];
255   PetscBool      flg;
256   PetscErrorCode ierr;
257 
258   PetscFunctionBegin;
259   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
260   if (!((PetscObject) prob)->type_name) {
261     defaultType = PETSCDSBASIC;
262   } else {
263     defaultType = ((PetscObject) prob)->type_name;
264   }
265   ierr = PetscDSRegisterAll();CHKERRQ(ierr);
266 
267   ierr = PetscObjectOptionsBegin((PetscObject) prob);CHKERRQ(ierr);
268   ierr = PetscOptionsFList("-petscds_type", "Discrete System", "PetscDSSetType", PetscDSList, defaultType, name, 256, &flg);CHKERRQ(ierr);
269   if (flg) {
270     ierr = PetscDSSetType(prob, name);CHKERRQ(ierr);
271   } else if (!((PetscObject) prob)->type_name) {
272     ierr = PetscDSSetType(prob, defaultType);CHKERRQ(ierr);
273   }
274   if (prob->ops->setfromoptions) {ierr = (*prob->ops->setfromoptions)(prob);CHKERRQ(ierr);}
275   /* process any options handlers added with PetscObjectAddOptionsHandler() */
276   ierr = PetscObjectProcessOptionsHandlers((PetscObject) prob);CHKERRQ(ierr);
277   ierr = PetscOptionsEnd();CHKERRQ(ierr);
278   ierr = PetscDSViewFromOptions(prob, NULL, "-petscds_view");CHKERRQ(ierr);
279   PetscFunctionReturn(0);
280 }
281 
282 #undef __FUNCT__
283 #define __FUNCT__ "PetscDSSetUp"
284 /*@C
285   PetscDSSetUp - Construct data structures for the PetscDS
286 
287   Collective on PetscDS
288 
289   Input Parameter:
290 . prob - the PetscDS object to setup
291 
292   Level: developer
293 
294 .seealso PetscDSView(), PetscDSDestroy()
295 @*/
296 PetscErrorCode PetscDSSetUp(PetscDS prob)
297 {
298   const PetscInt Nf = prob->Nf;
299   PetscInt       dim, work, NcMax = 0, NqMax = 0, f;
300   PetscErrorCode ierr;
301 
302   PetscFunctionBegin;
303   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
304   if (prob->setup) PetscFunctionReturn(0);
305   /* Calculate sizes */
306   ierr = PetscDSGetSpatialDimension(prob, &dim);CHKERRQ(ierr);
307   prob->totDim = prob->totDimBd = prob->totComp = 0;
308   ierr = PetscMalloc4(Nf,&prob->basis,Nf,&prob->basisDer,Nf,&prob->basisBd,Nf,&prob->basisDerBd);CHKERRQ(ierr);
309   for (f = 0; f < Nf; ++f) {
310     PetscFE         feBd = (PetscFE) prob->discBd[f];
311     PetscObject     obj;
312     PetscClassId    id;
313     PetscQuadrature q;
314     PetscInt        Nq = 0, Nb, Nc;
315 
316     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
317     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
318     if (id == PETSCFE_CLASSID)      {
319       PetscFE fe = (PetscFE) obj;
320 
321       ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
322       ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
323       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
324       ierr = PetscFEGetDefaultTabulation(fe, &prob->basis[f], &prob->basisDer[f], NULL);CHKERRQ(ierr);
325     } else if (id == PETSCFV_CLASSID) {
326       PetscFV fv = (PetscFV) obj;
327 
328       ierr = PetscFVGetQuadrature(fv, &q);CHKERRQ(ierr);
329       Nb   = 1;
330       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
331       ierr = PetscFVGetDefaultTabulation(fv, &prob->basis[f], &prob->basisDer[f], NULL);CHKERRQ(ierr);
332     } else SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
333     if (q) {ierr = PetscQuadratureGetData(q, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);}
334     NqMax          = PetscMax(NqMax, Nq);
335     NcMax          = PetscMax(NcMax, Nc);
336     prob->totDim  += Nb*Nc;
337     prob->totComp += Nc;
338     if (feBd) {
339       ierr = PetscFEGetDimension(feBd, &Nb);CHKERRQ(ierr);
340       ierr = PetscFEGetNumComponents(feBd, &Nc);CHKERRQ(ierr);
341       ierr = PetscFEGetDefaultTabulation(feBd, &prob->basisBd[f], &prob->basisDerBd[f], NULL);CHKERRQ(ierr);
342       prob->totDimBd += Nb*Nc;
343     }
344   }
345   work = PetscMax(prob->totComp*dim, PetscSqr(NcMax*dim));
346   /* Allocate works space */
347   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);
348   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);
349   if (prob->ops->setup) {ierr = (*prob->ops->setup)(prob);CHKERRQ(ierr);}
350   prob->setup = PETSC_TRUE;
351   PetscFunctionReturn(0);
352 }
353 
354 #undef __FUNCT__
355 #define __FUNCT__ "PetscDSDestroyStructs_Static"
356 static PetscErrorCode PetscDSDestroyStructs_Static(PetscDS prob)
357 {
358   PetscErrorCode ierr;
359 
360   PetscFunctionBegin;
361   ierr = PetscFree4(prob->basis,prob->basisDer,prob->basisBd,prob->basisDerBd);CHKERRQ(ierr);
362   ierr = PetscFree5(prob->u,prob->u_t,prob->u_x,prob->x,prob->refSpaceDer);CHKERRQ(ierr);
363   ierr = PetscFree6(prob->f0,prob->f1,prob->g0,prob->g1,prob->g2,prob->g3);CHKERRQ(ierr);
364   PetscFunctionReturn(0);
365 }
366 
367 #undef __FUNCT__
368 #define __FUNCT__ "PetscDSEnlarge_Static"
369 static PetscErrorCode PetscDSEnlarge_Static(PetscDS prob, PetscInt NfNew)
370 {
371   PetscObject   *tmpd, *tmpdbd;
372   PetscBool     *tmpi;
373   PointFunc     *tmpobj, *tmpf, *tmpg;
374   BdPointFunc   *tmpfbd, *tmpgbd;
375   RiemannFunc   *tmpr;
376   void         **tmpctx;
377   PetscInt       Nf = prob->Nf, f;
378   PetscErrorCode ierr;
379 
380   PetscFunctionBegin;
381   if (Nf >= NfNew) PetscFunctionReturn(0);
382   prob->setup = PETSC_FALSE;
383   ierr = PetscDSDestroyStructs_Static(prob);CHKERRQ(ierr);
384   ierr = PetscMalloc1(NfNew, &tmpd);CHKERRQ(ierr);
385   ierr = PetscMalloc1(NfNew, &tmpi);CHKERRQ(ierr);
386   for (f = 0; f < Nf; ++f) {tmpd[f] = prob->disc[f]; tmpi[f] = prob->implicit[f];}
387   for (f = Nf; f < NfNew; ++f) {ierr = PetscContainerCreate(PetscObjectComm((PetscObject) prob), (PetscContainer *) &tmpd[f]);CHKERRQ(ierr); tmpi[f] = PETSC_TRUE;}
388   ierr = PetscFree(prob->disc);CHKERRQ(ierr);
389   ierr = PetscFree(prob->implicit);CHKERRQ(ierr);
390   prob->Nf       = NfNew;
391   prob->disc     = tmpd;
392   prob->implicit = tmpi;
393   ierr = PetscCalloc5(NfNew, &tmpobj, NfNew*2, &tmpf, NfNew*NfNew*4, &tmpg, NfNew, &tmpr, NfNew, &tmpctx);CHKERRQ(ierr);
394   for (f = 0; f < Nf; ++f) tmpobj[f] = prob->obj[f];
395   for (f = 0; f < Nf*2; ++f) tmpf[f] = prob->f[f];
396   for (f = 0; f < Nf*Nf*4; ++f) tmpg[f] = prob->g[f];
397   for (f = 0; f < Nf; ++f) tmpr[f] = prob->r[f];
398   for (f = 0; f < Nf; ++f) tmpctx[f] = prob->ctx[f];
399   for (f = Nf; f < NfNew; ++f) tmpobj[f] = NULL;
400   for (f = Nf*2; f < NfNew*2; ++f) tmpf[f] = NULL;
401   for (f = Nf*Nf*4; f < NfNew*NfNew*4; ++f) tmpg[f] = NULL;
402   for (f = Nf; f < NfNew; ++f) tmpr[f] = NULL;
403   for (f = Nf; f < NfNew; ++f) tmpctx[f] = NULL;
404   ierr = PetscFree5(prob->obj, prob->f, prob->g, prob->r, prob->ctx);CHKERRQ(ierr);
405   prob->obj = tmpobj;
406   prob->f   = tmpf;
407   prob->g   = tmpg;
408   prob->r   = tmpr;
409   prob->ctx = tmpctx;
410   ierr = PetscMalloc1(NfNew, &tmpdbd);CHKERRQ(ierr);
411   for (f = 0; f < Nf; ++f) tmpdbd[f] = prob->discBd[f];
412   for (f = Nf; f < NfNew; ++f) tmpdbd[f] = NULL;
413   ierr = PetscFree(prob->discBd);CHKERRQ(ierr);
414   prob->discBd = tmpdbd;
415   ierr = PetscCalloc2(NfNew*2, &tmpfbd, NfNew*NfNew*4, &tmpgbd);CHKERRQ(ierr);
416   for (f = 0; f < Nf*2; ++f) tmpfbd[f] = prob->fBd[f];
417   for (f = 0; f < Nf*Nf*4; ++f) tmpgbd[f] = prob->gBd[f];
418   for (f = Nf*2; f < NfNew*2; ++f) tmpfbd[f] = NULL;
419   for (f = Nf*Nf*4; f < NfNew*NfNew*4; ++f) tmpgbd[f] = NULL;
420   ierr = PetscFree2(prob->fBd, prob->gBd);CHKERRQ(ierr);
421   prob->fBd = tmpfbd;
422   prob->gBd = tmpgbd;
423   PetscFunctionReturn(0);
424 }
425 
426 #undef __FUNCT__
427 #define __FUNCT__ "PetscDSDestroy"
428 /*@
429   PetscDSDestroy - Destroys a PetscDS object
430 
431   Collective on PetscDS
432 
433   Input Parameter:
434 . prob - the PetscDS object to destroy
435 
436   Level: developer
437 
438 .seealso PetscDSView()
439 @*/
440 PetscErrorCode PetscDSDestroy(PetscDS *prob)
441 {
442   PetscInt       f;
443   PetscErrorCode ierr;
444 
445   PetscFunctionBegin;
446   if (!*prob) PetscFunctionReturn(0);
447   PetscValidHeaderSpecific((*prob), PETSCDS_CLASSID, 1);
448 
449   if (--((PetscObject)(*prob))->refct > 0) {*prob = 0; PetscFunctionReturn(0);}
450   ((PetscObject) (*prob))->refct = 0;
451   ierr = PetscDSDestroyStructs_Static(*prob);CHKERRQ(ierr);
452   for (f = 0; f < (*prob)->Nf; ++f) {
453     ierr = PetscObjectDereference((*prob)->disc[f]);CHKERRQ(ierr);
454     ierr = PetscObjectDereference((*prob)->discBd[f]);CHKERRQ(ierr);
455   }
456   ierr = PetscFree((*prob)->implicit);CHKERRQ(ierr);
457   ierr = PetscFree((*prob)->disc);CHKERRQ(ierr);
458   ierr = PetscFree((*prob)->discBd);CHKERRQ(ierr);
459   ierr = PetscFree5((*prob)->obj,(*prob)->f,(*prob)->g,(*prob)->r,(*prob)->ctx);CHKERRQ(ierr);
460   ierr = PetscFree2((*prob)->fBd,(*prob)->gBd);CHKERRQ(ierr);
461   if ((*prob)->ops->destroy) {ierr = (*(*prob)->ops->destroy)(*prob);CHKERRQ(ierr);}
462   ierr = PetscHeaderDestroy(prob);CHKERRQ(ierr);
463   PetscFunctionReturn(0);
464 }
465 
466 #undef __FUNCT__
467 #define __FUNCT__ "PetscDSCreate"
468 /*@
469   PetscDSCreate - Creates an empty PetscDS object. The type can then be set with PetscDSSetType().
470 
471   Collective on MPI_Comm
472 
473   Input Parameter:
474 . comm - The communicator for the PetscDS object
475 
476   Output Parameter:
477 . prob - The PetscDS object
478 
479   Level: beginner
480 
481 .seealso: PetscDSSetType(), PETSCDSBASIC
482 @*/
483 PetscErrorCode PetscDSCreate(MPI_Comm comm, PetscDS *prob)
484 {
485   PetscDS   p;
486   PetscErrorCode ierr;
487 
488   PetscFunctionBegin;
489   PetscValidPointer(prob, 2);
490   *prob  = NULL;
491   ierr = PetscDSInitializePackage();CHKERRQ(ierr);
492 
493   ierr = PetscHeaderCreate(p, _p_PetscDS, struct _PetscDSOps, PETSCDS_CLASSID, "PetscDS", "Discrete System", "PetscDS", comm, PetscDSDestroy, PetscDSView);CHKERRQ(ierr);
494   ierr = PetscMemzero(p->ops, sizeof(struct _PetscDSOps));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) prob->implicit[f] = PETSC_FALSE;
743   }
744   PetscFunctionReturn(0);
745 }
746 
747 #undef __FUNCT__
748 #define __FUNCT__ "PetscDSSetBdDiscretization"
749 /*@
750   PetscDSSetBdDiscretization - Sets the boundary discretization object for the given field
751 
752   Not collective
753 
754   Input Parameters:
755 + prob - The PetscDS object
756 . f - The field number
757 - disc - The boundary discretization object
758 
759   Level: beginner
760 
761 .seealso: PetscDSGetBdDiscretization(), PetscDSAddBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
762 @*/
763 PetscErrorCode PetscDSSetBdDiscretization(PetscDS prob, PetscInt f, PetscObject disc)
764 {
765   PetscErrorCode ierr;
766 
767   PetscFunctionBegin;
768   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
769   if (disc) PetscValidPointer(disc, 3);
770   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
771   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
772   if (prob->discBd[f]) {ierr = PetscObjectDereference(prob->discBd[f]);CHKERRQ(ierr);}
773   prob->discBd[f] = disc;
774   ierr = PetscObjectReference(disc);CHKERRQ(ierr);
775   PetscFunctionReturn(0);
776 }
777 
778 #undef __FUNCT__
779 #define __FUNCT__ "PetscDSAddDiscretization"
780 /*@
781   PetscDSAddDiscretization - Adds a discretization object
782 
783   Not collective
784 
785   Input Parameters:
786 + prob - The PetscDS object
787 - disc - The boundary discretization object
788 
789   Level: beginner
790 
791 .seealso: PetscDSGetDiscretization(), PetscDSSetDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
792 @*/
793 PetscErrorCode PetscDSAddDiscretization(PetscDS prob, PetscObject disc)
794 {
795   PetscErrorCode ierr;
796 
797   PetscFunctionBegin;
798   ierr = PetscDSSetDiscretization(prob, prob->Nf, disc);CHKERRQ(ierr);
799   PetscFunctionReturn(0);
800 }
801 
802 #undef __FUNCT__
803 #define __FUNCT__ "PetscDSAddBdDiscretization"
804 /*@
805   PetscDSAddBdDiscretization - Adds a boundary discretization object
806 
807   Not collective
808 
809   Input Parameters:
810 + prob - The PetscDS object
811 - disc - The boundary discretization object
812 
813   Level: beginner
814 
815 .seealso: PetscDSGetBdDiscretization(), PetscDSSetBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
816 @*/
817 PetscErrorCode PetscDSAddBdDiscretization(PetscDS prob, PetscObject disc)
818 {
819   PetscErrorCode ierr;
820 
821   PetscFunctionBegin;
822   ierr = PetscDSSetBdDiscretization(prob, prob->Nf, disc);CHKERRQ(ierr);
823   PetscFunctionReturn(0);
824 }
825 
826 #undef __FUNCT__
827 #define __FUNCT__ "PetscDSGetImplicit"
828 /*@
829   PetscDSGetImplicit - Returns the flag for implicit solve for this field. This is just a guide for IMEX
830 
831   Not collective
832 
833   Input Parameters:
834 + prob - The PetscDS object
835 - f - The field number
836 
837   Output Parameter:
838 . implicit - The flag indicating what kind of solve to use for this field
839 
840   Level: developer
841 
842 .seealso: PetscDSSetImplicit(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
843 @*/
844 PetscErrorCode PetscDSGetImplicit(PetscDS prob, PetscInt f, PetscBool *implicit)
845 {
846   PetscFunctionBegin;
847   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
848   PetscValidPointer(implicit, 3);
849   if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf);
850   *implicit = prob->implicit[f];
851   PetscFunctionReturn(0);
852 }
853 
854 #undef __FUNCT__
855 #define __FUNCT__ "PetscDSSetImplicit"
856 /*@
857   PetscDSSetImplicit - Set the flag for implicit solve for this field. This is just a guide for IMEX
858 
859   Not collective
860 
861   Input Parameters:
862 + prob - The PetscDS object
863 . f - The field number
864 - implicit - The flag indicating what kind of solve to use for this field
865 
866   Level: developer
867 
868 .seealso: PetscDSGetImplicit(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
869 @*/
870 PetscErrorCode PetscDSSetImplicit(PetscDS prob, PetscInt f, PetscBool implicit)
871 {
872   PetscFunctionBegin;
873   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
874   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);
875   prob->implicit[f] = implicit;
876   PetscFunctionReturn(0);
877 }
878 
879 #undef __FUNCT__
880 #define __FUNCT__ "PetscDSGetObjective"
881 PetscErrorCode PetscDSGetObjective(PetscDS prob, PetscInt f,
882                                         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[]))
883 {
884   PetscFunctionBegin;
885   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
886   PetscValidPointer(obj, 2);
887   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);
888   *obj = prob->obj[f];
889   PetscFunctionReturn(0);
890 }
891 
892 #undef __FUNCT__
893 #define __FUNCT__ "PetscDSSetObjective"
894 PetscErrorCode PetscDSSetObjective(PetscDS prob, PetscInt f,
895                                         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[]))
896 {
897   PetscErrorCode ierr;
898 
899   PetscFunctionBegin;
900   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
901   PetscValidFunction(obj, 2);
902   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
903   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
904   prob->obj[f] = obj;
905   PetscFunctionReturn(0);
906 }
907 
908 #undef __FUNCT__
909 #define __FUNCT__ "PetscDSGetResidual"
910 PetscErrorCode PetscDSGetResidual(PetscDS prob, PetscInt f,
911                                        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[]),
912                                        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[]))
913 {
914   PetscFunctionBegin;
915   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
916   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);
917   if (f0) {PetscValidPointer(f0, 3); *f0 = prob->f[f*2+0];}
918   if (f1) {PetscValidPointer(f1, 4); *f1 = prob->f[f*2+1];}
919   PetscFunctionReturn(0);
920 }
921 
922 #undef __FUNCT__
923 #define __FUNCT__ "PetscDSSetResidual"
924 PetscErrorCode PetscDSSetResidual(PetscDS prob, PetscInt f,
925                                   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[]),
926                                   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[]))
927 {
928   PetscErrorCode ierr;
929 
930   PetscFunctionBegin;
931   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
932   PetscValidFunction(f0, 3);
933   PetscValidFunction(f1, 4);
934   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
935   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
936   prob->f[f*2+0] = f0;
937   prob->f[f*2+1] = f1;
938   PetscFunctionReturn(0);
939 }
940 
941 #undef __FUNCT__
942 #define __FUNCT__ "PetscDSGetJacobian"
943 PetscErrorCode PetscDSGetJacobian(PetscDS prob, PetscInt f, PetscInt g,
944                                        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[]),
945                                        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[]),
946                                        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[]),
947                                        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[]))
948 {
949   PetscFunctionBegin;
950   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
951   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);
952   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);
953   if (g0) {PetscValidPointer(g0, 4); *g0 = prob->g[(f*prob->Nf + g)*4+0];}
954   if (g1) {PetscValidPointer(g1, 5); *g1 = prob->g[(f*prob->Nf + g)*4+1];}
955   if (g2) {PetscValidPointer(g2, 6); *g2 = prob->g[(f*prob->Nf + g)*4+2];}
956   if (g3) {PetscValidPointer(g3, 7); *g3 = prob->g[(f*prob->Nf + g)*4+3];}
957   PetscFunctionReturn(0);
958 }
959 
960 #undef __FUNCT__
961 #define __FUNCT__ "PetscDSSetJacobian"
962 PetscErrorCode PetscDSSetJacobian(PetscDS prob, PetscInt f, PetscInt g,
963                                        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[]),
964                                        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[]),
965                                        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[]),
966                                        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[]))
967 {
968   PetscErrorCode ierr;
969 
970   PetscFunctionBegin;
971   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
972   if (g0) PetscValidFunction(g0, 4);
973   if (g1) PetscValidFunction(g1, 5);
974   if (g2) PetscValidFunction(g2, 6);
975   if (g3) PetscValidFunction(g3, 7);
976   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
977   if (g < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", g);
978   ierr = PetscDSEnlarge_Static(prob, PetscMax(f, g)+1);CHKERRQ(ierr);
979   prob->g[(f*prob->Nf + g)*4+0] = g0;
980   prob->g[(f*prob->Nf + g)*4+1] = g1;
981   prob->g[(f*prob->Nf + g)*4+2] = g2;
982   prob->g[(f*prob->Nf + g)*4+3] = g3;
983   PetscFunctionReturn(0);
984 }
985 
986 #undef __FUNCT__
987 #define __FUNCT__ "PetscDSGetRiemannSolver"
988 /*@C
989   PetscDSGetRiemannSolver - Returns the Riemann solver for the given field
990 
991   Not collective
992 
993   Input Arguments:
994 + prob - The PetscDS object
995 - f    - The field number
996 
997   Output Argument:
998 . r    - Riemann solver
999 
1000   Calling sequence for r:
1001 
1002 $ r(const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx)
1003 
1004 + x    - The coordinates at a point on the interface
1005 . n    - The normal vector to the interface
1006 . uL   - The state vector to the left of the interface
1007 . uR   - The state vector to the right of the interface
1008 . flux - output array of flux through the interface
1009 - ctx  - optional user context
1010 
1011   Level: intermediate
1012 
1013 .seealso: PetscDSSetRiemannSolver()
1014 @*/
1015 PetscErrorCode PetscDSGetRiemannSolver(PetscDS prob, PetscInt f,
1016                                        void (**r)(const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx))
1017 {
1018   PetscFunctionBegin;
1019   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1020   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);
1021   PetscValidPointer(r, 3);
1022   *r = prob->r[f];
1023   PetscFunctionReturn(0);
1024 }
1025 
1026 #undef __FUNCT__
1027 #define __FUNCT__ "PetscDSSetRiemannSolver"
1028 /*@C
1029   PetscDSSetRiemannSolver - Sets the Riemann solver for the given field
1030 
1031   Not collective
1032 
1033   Input Arguments:
1034 + prob - The PetscDS object
1035 . f    - The field number
1036 - r    - Riemann solver
1037 
1038   Calling sequence for r:
1039 
1040 $ r(const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx)
1041 
1042 + x    - The coordinates at a point on the interface
1043 . n    - The normal vector to the interface
1044 . uL   - The state vector to the left of the interface
1045 . uR   - The state vector to the right of the interface
1046 . flux - output array of flux through the interface
1047 - ctx  - optional user context
1048 
1049   Level: intermediate
1050 
1051 .seealso: PetscDSGetRiemannSolver()
1052 @*/
1053 PetscErrorCode PetscDSSetRiemannSolver(PetscDS prob, PetscInt f,
1054                                        void (*r)(const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx))
1055 {
1056   PetscErrorCode ierr;
1057 
1058   PetscFunctionBegin;
1059   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1060   PetscValidFunction(r, 3);
1061   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
1062   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
1063   prob->r[f] = r;
1064   PetscFunctionReturn(0);
1065 }
1066 
1067 #undef __FUNCT__
1068 #define __FUNCT__ "PetscDSGetContext"
1069 PetscErrorCode PetscDSGetContext(PetscDS prob, PetscInt f, void **ctx)
1070 {
1071   PetscFunctionBegin;
1072   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1073   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);
1074   PetscValidPointer(ctx, 3);
1075   *ctx = prob->ctx[f];
1076   PetscFunctionReturn(0);
1077 }
1078 
1079 #undef __FUNCT__
1080 #define __FUNCT__ "PetscDSSetContext"
1081 PetscErrorCode PetscDSSetContext(PetscDS prob, PetscInt f, void *ctx)
1082 {
1083   PetscErrorCode ierr;
1084 
1085   PetscFunctionBegin;
1086   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1087   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
1088   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
1089   prob->ctx[f] = ctx;
1090   PetscFunctionReturn(0);
1091 }
1092 
1093 #undef __FUNCT__
1094 #define __FUNCT__ "PetscDSGetBdResidual"
1095 PetscErrorCode PetscDSGetBdResidual(PetscDS prob, PetscInt f,
1096                                          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[]),
1097                                          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[]))
1098 {
1099   PetscFunctionBegin;
1100   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1101   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);
1102   if (f0) {PetscValidPointer(f0, 3); *f0 = prob->fBd[f*2+0];}
1103   if (f1) {PetscValidPointer(f1, 4); *f1 = prob->fBd[f*2+1];}
1104   PetscFunctionReturn(0);
1105 }
1106 
1107 #undef __FUNCT__
1108 #define __FUNCT__ "PetscDSSetBdResidual"
1109 PetscErrorCode PetscDSSetBdResidual(PetscDS prob, PetscInt f,
1110                                          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[]),
1111                                          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[]))
1112 {
1113   PetscErrorCode ierr;
1114 
1115   PetscFunctionBegin;
1116   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1117   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
1118   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
1119   if (f0) {PetscValidFunction(f0, 3); prob->fBd[f*2+0] = f0;}
1120   if (f1) {PetscValidFunction(f1, 4); prob->fBd[f*2+1] = f1;}
1121   PetscFunctionReturn(0);
1122 }
1123 
1124 #undef __FUNCT__
1125 #define __FUNCT__ "PetscDSGetBdJacobian"
1126 PetscErrorCode PetscDSGetBdJacobian(PetscDS prob, PetscInt f, PetscInt g,
1127                                          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[]),
1128                                          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[]),
1129                                          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[]),
1130                                          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[]))
1131 {
1132   PetscFunctionBegin;
1133   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1134   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);
1135   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);
1136   if (g0) {PetscValidPointer(g0, 4); *g0 = prob->gBd[(f*prob->Nf + g)*4+0];}
1137   if (g1) {PetscValidPointer(g1, 5); *g1 = prob->gBd[(f*prob->Nf + g)*4+1];}
1138   if (g2) {PetscValidPointer(g2, 6); *g2 = prob->gBd[(f*prob->Nf + g)*4+2];}
1139   if (g3) {PetscValidPointer(g3, 7); *g3 = prob->gBd[(f*prob->Nf + g)*4+3];}
1140   PetscFunctionReturn(0);
1141 }
1142 
1143 #undef __FUNCT__
1144 #define __FUNCT__ "PetscDSSetBdJacobian"
1145 PetscErrorCode PetscDSSetBdJacobian(PetscDS prob, PetscInt f, PetscInt g,
1146                                          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[]),
1147                                          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[]),
1148                                          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[]),
1149                                          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[]))
1150 {
1151   PetscErrorCode ierr;
1152 
1153   PetscFunctionBegin;
1154   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1155   if (g0) PetscValidFunction(g0, 4);
1156   if (g1) PetscValidFunction(g1, 5);
1157   if (g2) PetscValidFunction(g2, 6);
1158   if (g3) PetscValidFunction(g3, 7);
1159   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
1160   if (g < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", g);
1161   ierr = PetscDSEnlarge_Static(prob, PetscMax(f, g)+1);CHKERRQ(ierr);
1162   prob->gBd[(f*prob->Nf + g)*4+0] = g0;
1163   prob->gBd[(f*prob->Nf + g)*4+1] = g1;
1164   prob->gBd[(f*prob->Nf + g)*4+2] = g2;
1165   prob->gBd[(f*prob->Nf + g)*4+3] = g3;
1166   PetscFunctionReturn(0);
1167 }
1168 
1169 #undef __FUNCT__
1170 #define __FUNCT__ "PetscDSGetFieldOffset"
1171 /*@
1172   PetscDSGetFieldOffset - Returns the offset of the given field in the full space basis
1173 
1174   Not collective
1175 
1176   Input Parameters:
1177 + prob - The PetscDS object
1178 - f - The field number
1179 
1180   Output Parameter:
1181 . off - The offset
1182 
1183   Level: beginner
1184 
1185 .seealso: PetscDSGetBdFieldOffset(), PetscDSGetNumFields(), PetscDSCreate()
1186 @*/
1187 PetscErrorCode PetscDSGetFieldOffset(PetscDS prob, PetscInt f, PetscInt *off)
1188 {
1189   PetscInt       g;
1190   PetscErrorCode ierr;
1191 
1192   PetscFunctionBegin;
1193   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1194   PetscValidPointer(off, 3);
1195   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);
1196   *off = 0;
1197   for (g = 0; g < f; ++g) {
1198     PetscFE  fe = (PetscFE) prob->disc[g];
1199     PetscInt Nb, Nc;
1200 
1201     ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
1202     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
1203     *off += Nb*Nc;
1204   }
1205   PetscFunctionReturn(0);
1206 }
1207 
1208 #undef __FUNCT__
1209 #define __FUNCT__ "PetscDSGetBdFieldOffset"
1210 /*@
1211   PetscDSGetBdFieldOffset - Returns the offset of the given field in the full space boundary basis
1212 
1213   Not collective
1214 
1215   Input Parameters:
1216 + prob - The PetscDS object
1217 - f - The field number
1218 
1219   Output Parameter:
1220 . off - The boundary offset
1221 
1222   Level: beginner
1223 
1224 .seealso: PetscDSGetFieldOffset(), PetscDSGetNumFields(), PetscDSCreate()
1225 @*/
1226 PetscErrorCode PetscDSGetBdFieldOffset(PetscDS prob, PetscInt f, PetscInt *off)
1227 {
1228   PetscInt       g;
1229   PetscErrorCode ierr;
1230 
1231   PetscFunctionBegin;
1232   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1233   PetscValidPointer(off, 3);
1234   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);
1235   *off = 0;
1236   for (g = 0; g < f; ++g) {
1237     PetscFE  fe = (PetscFE) prob->discBd[g];
1238     PetscInt Nb, Nc;
1239 
1240     ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
1241     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
1242     *off += Nb*Nc;
1243   }
1244   PetscFunctionReturn(0);
1245 }
1246 
1247 #undef __FUNCT__
1248 #define __FUNCT__ "PetscDSGetComponentOffset"
1249 /*@
1250   PetscDSGetComponentOffset - Returns the offset of the given field on an evaluation point
1251 
1252   Not collective
1253 
1254   Input Parameters:
1255 + prob - The PetscDS object
1256 - f - The field number
1257 
1258   Output Parameter:
1259 . off - The offset
1260 
1261   Level: beginner
1262 
1263 .seealso: PetscDSGetBdFieldOffset(), PetscDSGetNumFields(), PetscDSCreate()
1264 @*/
1265 PetscErrorCode PetscDSGetComponentOffset(PetscDS prob, PetscInt f, PetscInt *off)
1266 {
1267   PetscInt       g;
1268   PetscErrorCode ierr;
1269 
1270   PetscFunctionBegin;
1271   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1272   PetscValidPointer(off, 3);
1273   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);
1274   *off = 0;
1275   for (g = 0; g < f; ++g) {
1276     PetscFE  fe = (PetscFE) prob->disc[g];
1277     PetscInt Nc;
1278 
1279     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
1280     *off += Nc;
1281   }
1282   PetscFunctionReturn(0);
1283 }
1284 
1285 #undef __FUNCT__
1286 #define __FUNCT__ "PetscDSGetTabulation"
1287 /*@C
1288   PetscDSGetTabulation - Return the basis tabulation at quadrature points for the volume discretization
1289 
1290   Not collective
1291 
1292   Input Parameter:
1293 . prob - The PetscDS object
1294 
1295   Output Parameters:
1296 + basis - The basis function tabulation at quadrature points
1297 - basisDer - The basis function derivative tabulation at quadrature points
1298 
1299   Level: intermediate
1300 
1301 .seealso: PetscDSGetBdTabulation(), PetscDSCreate()
1302 @*/
1303 PetscErrorCode PetscDSGetTabulation(PetscDS prob, PetscReal ***basis, PetscReal ***basisDer)
1304 {
1305   PetscErrorCode ierr;
1306 
1307   PetscFunctionBegin;
1308   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1309   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
1310   if (basis)    {PetscValidPointer(basis, 2);    *basis    = prob->basis;}
1311   if (basisDer) {PetscValidPointer(basisDer, 3); *basisDer = prob->basisDer;}
1312   PetscFunctionReturn(0);
1313 }
1314 
1315 #undef __FUNCT__
1316 #define __FUNCT__ "PetscDSGetBdTabulation"
1317 /*@C
1318   PetscDSGetBdTabulation - Return the basis tabulation at quadrature points for the boundary discretization
1319 
1320   Not collective
1321 
1322   Input Parameter:
1323 . prob - The PetscDS object
1324 
1325   Output Parameters:
1326 + basis - The basis function tabulation at quadrature points
1327 - basisDer - The basis function derivative tabulation at quadrature points
1328 
1329   Level: intermediate
1330 
1331 .seealso: PetscDSGetTabulation(), PetscDSCreate()
1332 @*/
1333 PetscErrorCode PetscDSGetBdTabulation(PetscDS prob, PetscReal ***basis, PetscReal ***basisDer)
1334 {
1335   PetscErrorCode ierr;
1336 
1337   PetscFunctionBegin;
1338   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1339   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
1340   if (basis)    {PetscValidPointer(basis, 2);    *basis    = prob->basisBd;}
1341   if (basisDer) {PetscValidPointer(basisDer, 3); *basisDer = prob->basisDerBd;}
1342   PetscFunctionReturn(0);
1343 }
1344 
1345 #undef __FUNCT__
1346 #define __FUNCT__ "PetscDSGetEvaluationArrays"
1347 PetscErrorCode PetscDSGetEvaluationArrays(PetscDS prob, PetscScalar **u, PetscScalar **u_t, PetscScalar **u_x)
1348 {
1349   PetscErrorCode ierr;
1350 
1351   PetscFunctionBegin;
1352   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1353   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
1354   if (u)   {PetscValidPointer(u, 2);   *u   = prob->u;}
1355   if (u_t) {PetscValidPointer(u_t, 3); *u_t = prob->u_t;}
1356   if (u_x) {PetscValidPointer(u_x, 4); *u_x = prob->u_x;}
1357   PetscFunctionReturn(0);
1358 }
1359 
1360 #undef __FUNCT__
1361 #define __FUNCT__ "PetscDSGetWeakFormArrays"
1362 PetscErrorCode PetscDSGetWeakFormArrays(PetscDS prob, PetscScalar **f0, PetscScalar **f1, PetscScalar **g0, PetscScalar **g1, PetscScalar **g2, PetscScalar **g3)
1363 {
1364   PetscErrorCode ierr;
1365 
1366   PetscFunctionBegin;
1367   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1368   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
1369   if (f0) {PetscValidPointer(f0, 2); *f0 = prob->f0;}
1370   if (f1) {PetscValidPointer(f1, 3); *f1 = prob->f1;}
1371   if (g0) {PetscValidPointer(g0, 4); *g0 = prob->g0;}
1372   if (g1) {PetscValidPointer(g1, 5); *g1 = prob->g1;}
1373   if (g2) {PetscValidPointer(g2, 6); *g2 = prob->g2;}
1374   if (g3) {PetscValidPointer(g3, 7); *g3 = prob->g3;}
1375   PetscFunctionReturn(0);
1376 }
1377 
1378 #undef __FUNCT__
1379 #define __FUNCT__ "PetscDSGetRefCoordArrays"
1380 PetscErrorCode PetscDSGetRefCoordArrays(PetscDS prob, PetscReal **x, PetscScalar **refSpaceDer)
1381 {
1382   PetscErrorCode ierr;
1383 
1384   PetscFunctionBegin;
1385   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1386   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
1387   if (x)           {PetscValidPointer(x, 2);           *x           = prob->x;}
1388   if (refSpaceDer) {PetscValidPointer(refSpaceDer, 3); *refSpaceDer = prob->refSpaceDer;}
1389   PetscFunctionReturn(0);
1390 }
1391 
1392 #undef __FUNCT__
1393 #define __FUNCT__ "PetscDSDestroy_Basic"
1394 static PetscErrorCode PetscDSDestroy_Basic(PetscDS prob)
1395 {
1396   PetscFunctionBegin;
1397   PetscFunctionReturn(0);
1398 }
1399 
1400 #undef __FUNCT__
1401 #define __FUNCT__ "PetscDSInitialize_Basic"
1402 static PetscErrorCode PetscDSInitialize_Basic(PetscDS prob)
1403 {
1404   PetscFunctionBegin;
1405   prob->ops->setfromoptions = NULL;
1406   prob->ops->setup          = NULL;
1407   prob->ops->view           = NULL;
1408   prob->ops->destroy        = PetscDSDestroy_Basic;
1409   PetscFunctionReturn(0);
1410 }
1411 
1412 /*MC
1413   PETSCDSBASIC = "basic" - A discrete system with pointwise residual and boundary residual functions
1414 
1415   Level: intermediate
1416 
1417 .seealso: PetscDSType, PetscDSCreate(), PetscDSSetType()
1418 M*/
1419 
1420 #undef __FUNCT__
1421 #define __FUNCT__ "PetscDSCreate_Basic"
1422 PETSC_EXTERN PetscErrorCode PetscDSCreate_Basic(PetscDS prob)
1423 {
1424   PetscDS_Basic *b;
1425   PetscErrorCode      ierr;
1426 
1427   PetscFunctionBegin;
1428   PetscValidHeaderSpecific(prob, PETSCSPACE_CLASSID, 1);
1429   ierr       = PetscNewLog(prob, &b);CHKERRQ(ierr);
1430   prob->data = b;
1431 
1432   ierr = PetscDSInitialize_Basic(prob);CHKERRQ(ierr);
1433   PetscFunctionReturn(0);
1434 }
1435