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