xref: /petsc/src/dm/dt/interface/dtds.c (revision 420e96eda417c4a9936f208d8e146b6072205ea6)
1 #include <petsc-private/petscdsimpl.h> /*I "petscds.h" I*/
2 
3 PetscClassId PETSCDS_CLASSID = 0;
4 
5 PetscFunctionList PetscDSList              = NULL;
6 PetscBool         PetscDSRegisterAllCalled = PETSC_FALSE;
7 
8 #undef __FUNCT__
9 #define __FUNCT__ "PetscDSRegister"
10 /*@C
11   PetscDSRegister - Adds a new PetscDS implementation
12 
13   Not Collective
14 
15   Input Parameters:
16 + name        - The name of a new user-defined creation routine
17 - create_func - The creation routine itself
18 
19   Notes:
20   PetscDSRegister() may be called multiple times to add several user-defined PetscDSs
21 
22   Sample usage:
23 .vb
24     PetscDSRegister("my_ds", MyPetscDSCreate);
25 .ve
26 
27   Then, your PetscDS type can be chosen with the procedural interface via
28 .vb
29     PetscDSCreate(MPI_Comm, PetscDS *);
30     PetscDSSetType(PetscDS, "my_ds");
31 .ve
32    or at runtime via the option
33 .vb
34     -petscds_type my_ds
35 .ve
36 
37   Level: advanced
38 
39 .keywords: PetscDS, register
40 .seealso: PetscDSRegisterAll(), PetscDSRegisterDestroy()
41 
42 @*/
43 PetscErrorCode PetscDSRegister(const char sname[], PetscErrorCode (*function)(PetscDS))
44 {
45   PetscErrorCode ierr;
46 
47   PetscFunctionBegin;
48   ierr = PetscFunctionListAdd(&PetscDSList, sname, function);CHKERRQ(ierr);
49   PetscFunctionReturn(0);
50 }
51 
52 #undef __FUNCT__
53 #define __FUNCT__ "PetscDSSetType"
54 /*@C
55   PetscDSSetType - Builds a particular PetscDS
56 
57   Collective on PetscDS
58 
59   Input Parameters:
60 + prob - The PetscDS object
61 - name - The kind of system
62 
63   Options Database Key:
64 . -petscds_type <type> - Sets the PetscDS type; use -help for a list of available types
65 
66   Level: intermediate
67 
68 .keywords: PetscDS, set, type
69 .seealso: PetscDSGetType(), PetscDSCreate()
70 @*/
71 PetscErrorCode PetscDSSetType(PetscDS prob, PetscDSType name)
72 {
73   PetscErrorCode (*r)(PetscDS);
74   PetscBool      match;
75   PetscErrorCode ierr;
76 
77   PetscFunctionBegin;
78   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
79   ierr = PetscObjectTypeCompare((PetscObject) prob, name, &match);CHKERRQ(ierr);
80   if (match) PetscFunctionReturn(0);
81 
82   ierr = PetscDSRegisterAll();CHKERRQ(ierr);
83   ierr = PetscFunctionListFind(PetscDSList, name, &r);CHKERRQ(ierr);
84   if (!r) SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscDS type: %s", name);
85 
86   if (prob->ops->destroy) {
87     ierr             = (*prob->ops->destroy)(prob);CHKERRQ(ierr);
88     prob->ops->destroy = NULL;
89   }
90   ierr = (*r)(prob);CHKERRQ(ierr);
91   ierr = PetscObjectChangeTypeName((PetscObject) prob, name);CHKERRQ(ierr);
92   PetscFunctionReturn(0);
93 }
94 
95 #undef __FUNCT__
96 #define __FUNCT__ "PetscDSGetType"
97 /*@C
98   PetscDSGetType - Gets the PetscDS type name (as a string) from the object.
99 
100   Not Collective
101 
102   Input Parameter:
103 . prob  - The PetscDS
104 
105   Output Parameter:
106 . name - The PetscDS type name
107 
108   Level: intermediate
109 
110 .keywords: PetscDS, get, type, name
111 .seealso: PetscDSSetType(), PetscDSCreate()
112 @*/
113 PetscErrorCode PetscDSGetType(PetscDS prob, PetscDSType *name)
114 {
115   PetscErrorCode ierr;
116 
117   PetscFunctionBegin;
118   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
119   PetscValidPointer(name, 2);
120   ierr = PetscDSRegisterAll();CHKERRQ(ierr);
121   *name = ((PetscObject) prob)->type_name;
122   PetscFunctionReturn(0);
123 }
124 
125 #undef __FUNCT__
126 #define __FUNCT__ "PetscDSView_Ascii"
127 static PetscErrorCode PetscDSView_Ascii(PetscDS prob, PetscViewer viewer)
128 {
129   PetscViewerFormat format;
130   PetscInt          f;
131   PetscErrorCode    ierr;
132 
133   PetscFunctionBegin;
134   ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
135   ierr = PetscViewerASCIIPrintf(viewer, "Discrete System with %d fields\n", prob->Nf);CHKERRQ(ierr);
136   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
137   for (f = 0; f < prob->Nf; ++f) {
138     PetscObject  obj;
139     PetscClassId id;
140     const char  *name;
141     PetscInt     Nc;
142 
143     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
144     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
145     ierr = PetscObjectGetName(obj, &name);CHKERRQ(ierr);
146     ierr = PetscViewerASCIIPrintf(viewer, "Field %s", name ? name : "<unknown>");CHKERRQ(ierr);
147     if (id == PETSCFE_CLASSID)      {
148       ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);
149       ierr = PetscViewerASCIIPrintf(viewer, " FEM");CHKERRQ(ierr);
150     } else if (id == PETSCFV_CLASSID) {
151       ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);
152       ierr = PetscViewerASCIIPrintf(viewer, " FVM");CHKERRQ(ierr);
153     }
154     else SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
155     if (Nc > 1) {ierr = PetscViewerASCIIPrintf(viewer, "%d components", Nc);CHKERRQ(ierr);}
156     else        {ierr = PetscViewerASCIIPrintf(viewer, "%d component ", Nc);CHKERRQ(ierr);}
157     if (prob->implicit[f]) {ierr = PetscViewerASCIIPrintf(viewer, " (implicit)");CHKERRQ(ierr);}
158     else                   {ierr = PetscViewerASCIIPrintf(viewer, " (explicit)");CHKERRQ(ierr);}
159     if (prob->adjacency[f*2+0]) {
160       if (prob->adjacency[f*2+1]) {ierr = PetscViewerASCIIPrintf(viewer, " (adj FVM++)");CHKERRQ(ierr);}
161       else                        {ierr = PetscViewerASCIIPrintf(viewer, " (adj FVM)");CHKERRQ(ierr);}
162     } else {
163       if (prob->adjacency[f*2+1]) {ierr = PetscViewerASCIIPrintf(viewer, " (adj FEM)");CHKERRQ(ierr);}
164       else                        {ierr = PetscViewerASCIIPrintf(viewer, " (adj FUNKY)");CHKERRQ(ierr);}
165     }
166     ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr);
167     if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
168       if (id == PETSCFE_CLASSID)      {ierr = PetscFEView((PetscFE) obj, viewer);CHKERRQ(ierr);}
169       else if (id == PETSCFV_CLASSID) {ierr = PetscFVView((PetscFV) obj, viewer);CHKERRQ(ierr);}
170     }
171   }
172   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
173   PetscFunctionReturn(0);
174 }
175 
176 #undef __FUNCT__
177 #define __FUNCT__ "PetscDSView"
178 /*@C
179   PetscDSView - Views a PetscDS
180 
181   Collective on PetscDS
182 
183   Input Parameter:
184 + prob - the PetscDS object to view
185 - v  - the viewer
186 
187   Level: developer
188 
189 .seealso PetscDSDestroy()
190 @*/
191 PetscErrorCode PetscDSView(PetscDS prob, PetscViewer v)
192 {
193   PetscBool      iascii;
194   PetscErrorCode ierr;
195 
196   PetscFunctionBegin;
197   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
198   if (!v) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) prob), &v);CHKERRQ(ierr);}
199   else    {PetscValidHeaderSpecific(v, PETSC_VIEWER_CLASSID, 2);}
200   ierr = PetscObjectTypeCompare((PetscObject) v, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
201   if (iascii) {ierr = PetscDSView_Ascii(prob, v);CHKERRQ(ierr);}
202   if (prob->ops->view) {ierr = (*prob->ops->view)(prob, v);CHKERRQ(ierr);}
203   PetscFunctionReturn(0);
204 }
205 
206 #undef __FUNCT__
207 #define __FUNCT__ "PetscDSViewFromOptions"
208 /*
209   PetscDSViewFromOptions - Processes command line options to determine if/how a PetscDS is to be viewed.
210 
211   Collective on PetscDS
212 
213   Input Parameters:
214 + prob   - the PetscDS
215 . prefix - prefix to use for viewing, or NULL to use prefix of 'rnd'
216 - optionname - option to activate viewing
217 
218   Level: intermediate
219 
220 .keywords: PetscDS, view, options, database
221 .seealso: VecViewFromOptions(), MatViewFromOptions()
222 */
223 PetscErrorCode PetscDSViewFromOptions(PetscDS prob, const char prefix[], const char optionname[])
224 {
225   PetscViewer       viewer;
226   PetscViewerFormat format;
227   PetscBool         flg;
228   PetscErrorCode    ierr;
229 
230   PetscFunctionBegin;
231   if (prefix) {ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject) prob), prefix, optionname, &viewer, &format, &flg);CHKERRQ(ierr);}
232   else        {ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject) prob), ((PetscObject) prob)->prefix, optionname, &viewer, &format, &flg);CHKERRQ(ierr);}
233   if (flg) {
234     ierr = PetscViewerPushFormat(viewer, format);CHKERRQ(ierr);
235     ierr = PetscDSView(prob, viewer);CHKERRQ(ierr);
236     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
237     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
238   }
239   PetscFunctionReturn(0);
240 }
241 
242 #undef __FUNCT__
243 #define __FUNCT__ "PetscDSSetFromOptions"
244 /*@
245   PetscDSSetFromOptions - sets parameters in a PetscDS from the options database
246 
247   Collective on PetscDS
248 
249   Input Parameter:
250 . prob - the PetscDS object to set options for
251 
252   Options Database:
253 
254   Level: developer
255 
256 .seealso PetscDSView()
257 @*/
258 PetscErrorCode PetscDSSetFromOptions(PetscDS prob)
259 {
260   const char    *defaultType;
261   char           name[256];
262   PetscBool      flg;
263   PetscErrorCode ierr;
264 
265   PetscFunctionBegin;
266   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
267   if (!((PetscObject) prob)->type_name) {
268     defaultType = PETSCDSBASIC;
269   } else {
270     defaultType = ((PetscObject) prob)->type_name;
271   }
272   ierr = PetscDSRegisterAll();CHKERRQ(ierr);
273 
274   ierr = PetscObjectOptionsBegin((PetscObject) prob);CHKERRQ(ierr);
275   ierr = PetscOptionsFList("-petscds_type", "Discrete System", "PetscDSSetType", PetscDSList, defaultType, name, 256, &flg);CHKERRQ(ierr);
276   if (flg) {
277     ierr = PetscDSSetType(prob, name);CHKERRQ(ierr);
278   } else if (!((PetscObject) prob)->type_name) {
279     ierr = PetscDSSetType(prob, defaultType);CHKERRQ(ierr);
280   }
281   if (prob->ops->setfromoptions) {ierr = (*prob->ops->setfromoptions)(prob);CHKERRQ(ierr);}
282   /* process any options handlers added with PetscObjectAddOptionsHandler() */
283   ierr = PetscObjectProcessOptionsHandlers((PetscObject) prob);CHKERRQ(ierr);
284   ierr = PetscOptionsEnd();CHKERRQ(ierr);
285   ierr = PetscDSViewFromOptions(prob, NULL, "-petscds_view");CHKERRQ(ierr);
286   PetscFunctionReturn(0);
287 }
288 
289 #undef __FUNCT__
290 #define __FUNCT__ "PetscDSSetUp"
291 /*@C
292   PetscDSSetUp - Construct data structures for the PetscDS
293 
294   Collective on PetscDS
295 
296   Input Parameter:
297 . prob - the PetscDS object to setup
298 
299   Level: developer
300 
301 .seealso PetscDSView(), PetscDSDestroy()
302 @*/
303 PetscErrorCode PetscDSSetUp(PetscDS prob)
304 {
305   const PetscInt Nf = prob->Nf;
306   PetscInt       dim, work, NcMax = 0, NqMax = 0, f;
307   PetscErrorCode ierr;
308 
309   PetscFunctionBegin;
310   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
311   if (prob->setup) PetscFunctionReturn(0);
312   /* Calculate sizes */
313   ierr = PetscDSGetSpatialDimension(prob, &dim);CHKERRQ(ierr);
314   prob->totDim = prob->totDimBd = prob->totComp = 0;
315   ierr = PetscCalloc4(Nf+1,&prob->off,Nf+1,&prob->offDer,Nf+1,&prob->offBd,Nf+1,&prob->offDerBd);CHKERRQ(ierr);
316   ierr = PetscMalloc4(Nf,&prob->basis,Nf,&prob->basisDer,Nf,&prob->basisBd,Nf,&prob->basisDerBd);CHKERRQ(ierr);
317   for (f = 0; f < Nf; ++f) {
318     PetscFE         feBd = (PetscFE) prob->discBd[f];
319     PetscObject     obj;
320     PetscClassId    id;
321     PetscQuadrature q;
322     PetscInt        Nq = 0, Nb, Nc;
323 
324     ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr);
325     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
326     if (id == PETSCFE_CLASSID)      {
327       PetscFE fe = (PetscFE) obj;
328 
329       ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
330       ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
331       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
332       ierr = PetscFEGetDefaultTabulation(fe, &prob->basis[f], &prob->basisDer[f], NULL);CHKERRQ(ierr);
333     } else if (id == PETSCFV_CLASSID) {
334       PetscFV fv = (PetscFV) obj;
335 
336       ierr = PetscFVGetQuadrature(fv, &q);CHKERRQ(ierr);
337       Nb   = 1;
338       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
339       ierr = PetscFVGetDefaultTabulation(fv, &prob->basis[f], &prob->basisDer[f], NULL);CHKERRQ(ierr);
340     } else SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
341     prob->off[f+1]    = Nc     + prob->off[f];
342     prob->offDer[f+1] = Nc*dim + prob->offDer[f];
343     if (q) {ierr = PetscQuadratureGetData(q, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);}
344     NqMax          = PetscMax(NqMax, Nq);
345     NcMax          = PetscMax(NcMax, Nc);
346     prob->totDim  += Nb*Nc;
347     prob->totComp += Nc;
348     if (feBd) {
349       ierr = PetscFEGetDimension(feBd, &Nb);CHKERRQ(ierr);
350       ierr = PetscFEGetNumComponents(feBd, &Nc);CHKERRQ(ierr);
351       ierr = PetscFEGetDefaultTabulation(feBd, &prob->basisBd[f], &prob->basisDerBd[f], NULL);CHKERRQ(ierr);
352       prob->totDimBd += Nb*Nc;
353       prob->offBd[f+1]    = Nc     + prob->offBd[f];
354       prob->offDerBd[f+1] = Nc*dim + prob->offDerBd[f];
355     }
356   }
357   work = PetscMax(prob->totComp*dim, PetscSqr(NcMax*dim));
358   /* Allocate works space */
359   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);
360   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);
361   if (prob->ops->setup) {ierr = (*prob->ops->setup)(prob);CHKERRQ(ierr);}
362   prob->setup = PETSC_TRUE;
363   PetscFunctionReturn(0);
364 }
365 
366 #undef __FUNCT__
367 #define __FUNCT__ "PetscDSDestroyStructs_Static"
368 static PetscErrorCode PetscDSDestroyStructs_Static(PetscDS prob)
369 {
370   PetscErrorCode ierr;
371 
372   PetscFunctionBegin;
373   ierr = PetscFree4(prob->off,prob->offDer,prob->offBd,prob->offDerBd);CHKERRQ(ierr);
374   ierr = PetscFree4(prob->basis,prob->basisDer,prob->basisBd,prob->basisDerBd);CHKERRQ(ierr);
375   ierr = PetscFree5(prob->u,prob->u_t,prob->u_x,prob->x,prob->refSpaceDer);CHKERRQ(ierr);
376   ierr = PetscFree6(prob->f0,prob->f1,prob->g0,prob->g1,prob->g2,prob->g3);CHKERRQ(ierr);
377   PetscFunctionReturn(0);
378 }
379 
380 #undef __FUNCT__
381 #define __FUNCT__ "PetscDSEnlarge_Static"
382 static PetscErrorCode PetscDSEnlarge_Static(PetscDS prob, PetscInt NfNew)
383 {
384   PetscObject      *tmpd, *tmpdbd;
385   PetscBool        *tmpi, *tmpa;
386   PetscPointFunc   *tmpobj, *tmpf, *tmpg;
387   PetscBdPointFunc *tmpfbd, *tmpgbd;
388   PetscRiemannFunc *tmpr;
389   void            **tmpctx;
390   PetscInt          Nf = prob->Nf, f, i;
391   PetscErrorCode    ierr;
392 
393   PetscFunctionBegin;
394   if (Nf >= NfNew) PetscFunctionReturn(0);
395   prob->setup = PETSC_FALSE;
396   ierr = PetscDSDestroyStructs_Static(prob);CHKERRQ(ierr);
397   ierr = PetscMalloc4(NfNew, &tmpd, NfNew, &tmpdbd, NfNew, &tmpi, NfNew*2, &tmpa);CHKERRQ(ierr);
398   for (f = 0; f < Nf; ++f) {tmpd[f] = prob->disc[f]; tmpdbd[f] = prob->discBd[f]; tmpi[f] = prob->implicit[f]; for (i = 0; i < 2; ++i) tmpa[f*2+i] = prob->adjacency[f*2+i];}
399   for (f = Nf; f < NfNew; ++f) {ierr = PetscContainerCreate(PetscObjectComm((PetscObject) prob), (PetscContainer *) &tmpd[f]);CHKERRQ(ierr);
400     tmpdbd[f] = NULL; tmpi[f] = PETSC_TRUE; tmpa[f*2+0] = PETSC_FALSE; tmpa[f*2+1] = PETSC_TRUE;}
401   ierr = PetscFree4(prob->disc, prob->discBd, prob->implicit, prob->adjacency);CHKERRQ(ierr);
402   prob->Nf        = NfNew;
403   prob->disc      = tmpd;
404   prob->discBd    = tmpdbd;
405   prob->implicit  = tmpi;
406   prob->adjacency = tmpa;
407   ierr = PetscCalloc5(NfNew, &tmpobj, NfNew*2, &tmpf, NfNew*NfNew*4, &tmpg, NfNew, &tmpr, NfNew, &tmpctx);CHKERRQ(ierr);
408   for (f = 0; f < Nf; ++f) tmpobj[f] = prob->obj[f];
409   for (f = 0; f < Nf*2; ++f) tmpf[f] = prob->f[f];
410   for (f = 0; f < Nf*Nf*4; ++f) tmpg[f] = prob->g[f];
411   for (f = 0; f < Nf; ++f) tmpr[f] = prob->r[f];
412   for (f = 0; f < Nf; ++f) tmpctx[f] = prob->ctx[f];
413   for (f = Nf; f < NfNew; ++f) tmpobj[f] = NULL;
414   for (f = Nf*2; f < NfNew*2; ++f) tmpf[f] = NULL;
415   for (f = Nf*Nf*4; f < NfNew*NfNew*4; ++f) tmpg[f] = NULL;
416   for (f = Nf; f < NfNew; ++f) tmpr[f] = NULL;
417   for (f = Nf; f < NfNew; ++f) tmpctx[f] = NULL;
418   ierr = PetscFree5(prob->obj, prob->f, prob->g, prob->r, prob->ctx);CHKERRQ(ierr);
419   prob->obj = tmpobj;
420   prob->f   = tmpf;
421   prob->g   = tmpg;
422   prob->r   = tmpr;
423   prob->ctx = tmpctx;
424   ierr = PetscCalloc2(NfNew*2, &tmpfbd, NfNew*NfNew*4, &tmpgbd);CHKERRQ(ierr);
425   for (f = 0; f < Nf*2; ++f) tmpfbd[f] = prob->fBd[f];
426   for (f = 0; f < Nf*Nf*4; ++f) tmpgbd[f] = prob->gBd[f];
427   for (f = Nf*2; f < NfNew*2; ++f) tmpfbd[f] = NULL;
428   for (f = Nf*Nf*4; f < NfNew*NfNew*4; ++f) tmpgbd[f] = NULL;
429   ierr = PetscFree2(prob->fBd, prob->gBd);CHKERRQ(ierr);
430   prob->fBd = tmpfbd;
431   prob->gBd = tmpgbd;
432   PetscFunctionReturn(0);
433 }
434 
435 #undef __FUNCT__
436 #define __FUNCT__ "PetscDSDestroy"
437 /*@
438   PetscDSDestroy - Destroys a PetscDS object
439 
440   Collective on PetscDS
441 
442   Input Parameter:
443 . prob - the PetscDS object to destroy
444 
445   Level: developer
446 
447 .seealso PetscDSView()
448 @*/
449 PetscErrorCode PetscDSDestroy(PetscDS *prob)
450 {
451   PetscInt       f;
452   PetscErrorCode ierr;
453 
454   PetscFunctionBegin;
455   if (!*prob) PetscFunctionReturn(0);
456   PetscValidHeaderSpecific((*prob), PETSCDS_CLASSID, 1);
457 
458   if (--((PetscObject)(*prob))->refct > 0) {*prob = 0; PetscFunctionReturn(0);}
459   ((PetscObject) (*prob))->refct = 0;
460   ierr = PetscDSDestroyStructs_Static(*prob);CHKERRQ(ierr);
461   for (f = 0; f < (*prob)->Nf; ++f) {
462     ierr = PetscObjectDereference((*prob)->disc[f]);CHKERRQ(ierr);
463     ierr = PetscObjectDereference((*prob)->discBd[f]);CHKERRQ(ierr);
464   }
465   ierr = PetscFree4((*prob)->disc, (*prob)->discBd, (*prob)->implicit, (*prob)->adjacency);CHKERRQ(ierr);
466   ierr = PetscFree5((*prob)->obj,(*prob)->f,(*prob)->g,(*prob)->r,(*prob)->ctx);CHKERRQ(ierr);
467   ierr = PetscFree2((*prob)->fBd,(*prob)->gBd);CHKERRQ(ierr);
468   if ((*prob)->ops->destroy) {ierr = (*(*prob)->ops->destroy)(*prob);CHKERRQ(ierr);}
469   ierr = PetscHeaderDestroy(prob);CHKERRQ(ierr);
470   PetscFunctionReturn(0);
471 }
472 
473 #undef __FUNCT__
474 #define __FUNCT__ "PetscDSCreate"
475 /*@
476   PetscDSCreate - Creates an empty PetscDS object. The type can then be set with PetscDSSetType().
477 
478   Collective on MPI_Comm
479 
480   Input Parameter:
481 . comm - The communicator for the PetscDS object
482 
483   Output Parameter:
484 . prob - The PetscDS object
485 
486   Level: beginner
487 
488 .seealso: PetscDSSetType(), PETSCDSBASIC
489 @*/
490 PetscErrorCode PetscDSCreate(MPI_Comm comm, PetscDS *prob)
491 {
492   PetscDS   p;
493   PetscErrorCode ierr;
494 
495   PetscFunctionBegin;
496   PetscValidPointer(prob, 2);
497   *prob  = NULL;
498   ierr = PetscDSInitializePackage();CHKERRQ(ierr);
499 
500   ierr = PetscHeaderCreate(p, _p_PetscDS, struct _PetscDSOps, PETSCDS_CLASSID, "PetscDS", "Discrete System", "PetscDS", comm, PetscDSDestroy, PetscDSView);CHKERRQ(ierr);
501   ierr = PetscMemzero(p->ops, sizeof(struct _PetscDSOps));CHKERRQ(ierr);
502 
503   p->Nf    = 0;
504   p->setup = PETSC_FALSE;
505 
506   *prob = p;
507   PetscFunctionReturn(0);
508 }
509 
510 #undef __FUNCT__
511 #define __FUNCT__ "PetscDSGetNumFields"
512 /*@
513   PetscDSGetNumFields - Returns the number of fields in the DS
514 
515   Not collective
516 
517   Input Parameter:
518 . prob - The PetscDS object
519 
520   Output Parameter:
521 . Nf - The number of fields
522 
523   Level: beginner
524 
525 .seealso: PetscDSGetSpatialDimension(), PetscDSCreate()
526 @*/
527 PetscErrorCode PetscDSGetNumFields(PetscDS prob, PetscInt *Nf)
528 {
529   PetscFunctionBegin;
530   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
531   PetscValidPointer(Nf, 2);
532   *Nf = prob->Nf;
533   PetscFunctionReturn(0);
534 }
535 
536 #undef __FUNCT__
537 #define __FUNCT__ "PetscDSGetSpatialDimension"
538 /*@
539   PetscDSGetSpatialDimension - Returns the spatial dimension of the DS
540 
541   Not collective
542 
543   Input Parameter:
544 . prob - The PetscDS object
545 
546   Output Parameter:
547 . dim - The spatial dimension
548 
549   Level: beginner
550 
551 .seealso: PetscDSGetNumFields(), PetscDSCreate()
552 @*/
553 PetscErrorCode PetscDSGetSpatialDimension(PetscDS prob, PetscInt *dim)
554 {
555   PetscErrorCode ierr;
556 
557   PetscFunctionBegin;
558   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
559   PetscValidPointer(dim, 2);
560   *dim = 0;
561   if (prob->Nf) {
562     PetscObject  obj;
563     PetscClassId id;
564 
565     ierr = PetscDSGetDiscretization(prob, 0, &obj);CHKERRQ(ierr);
566     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
567     if (id == PETSCFE_CLASSID)      {ierr = PetscFEGetSpatialDimension((PetscFE) obj, dim);CHKERRQ(ierr);}
568     else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetSpatialDimension((PetscFV) obj, dim);CHKERRQ(ierr);}
569     else SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", 0);
570   }
571   PetscFunctionReturn(0);
572 }
573 
574 #undef __FUNCT__
575 #define __FUNCT__ "PetscDSGetTotalDimension"
576 /*@
577   PetscDSGetTotalDimension - Returns the total size of the approximation space for this system
578 
579   Not collective
580 
581   Input Parameter:
582 . prob - The PetscDS object
583 
584   Output Parameter:
585 . dim - The total problem dimension
586 
587   Level: beginner
588 
589 .seealso: PetscDSGetNumFields(), PetscDSCreate()
590 @*/
591 PetscErrorCode PetscDSGetTotalDimension(PetscDS prob, PetscInt *dim)
592 {
593   PetscErrorCode ierr;
594 
595   PetscFunctionBegin;
596   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
597   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
598   PetscValidPointer(dim, 2);
599   *dim = prob->totDim;
600   PetscFunctionReturn(0);
601 }
602 
603 #undef __FUNCT__
604 #define __FUNCT__ "PetscDSGetTotalBdDimension"
605 /*@
606   PetscDSGetTotalBdDimension - Returns the total size of the boundary approximation space for this system
607 
608   Not collective
609 
610   Input Parameter:
611 . prob - The PetscDS object
612 
613   Output Parameter:
614 . dim - The total boundary problem dimension
615 
616   Level: beginner
617 
618 .seealso: PetscDSGetNumFields(), PetscDSCreate()
619 @*/
620 PetscErrorCode PetscDSGetTotalBdDimension(PetscDS prob, PetscInt *dim)
621 {
622   PetscErrorCode ierr;
623 
624   PetscFunctionBegin;
625   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
626   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
627   PetscValidPointer(dim, 2);
628   *dim = prob->totDimBd;
629   PetscFunctionReturn(0);
630 }
631 
632 #undef __FUNCT__
633 #define __FUNCT__ "PetscDSGetTotalComponents"
634 /*@
635   PetscDSGetTotalComponents - Returns the total number of components in this system
636 
637   Not collective
638 
639   Input Parameter:
640 . prob - The PetscDS object
641 
642   Output Parameter:
643 . dim - The total number of components
644 
645   Level: beginner
646 
647 .seealso: PetscDSGetNumFields(), PetscDSCreate()
648 @*/
649 PetscErrorCode PetscDSGetTotalComponents(PetscDS prob, PetscInt *Nc)
650 {
651   PetscErrorCode ierr;
652 
653   PetscFunctionBegin;
654   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
655   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
656   PetscValidPointer(Nc, 2);
657   *Nc = prob->totComp;
658   PetscFunctionReturn(0);
659 }
660 
661 #undef __FUNCT__
662 #define __FUNCT__ "PetscDSGetDiscretization"
663 /*@
664   PetscDSGetDiscretization - Returns the discretization object for the given field
665 
666   Not collective
667 
668   Input Parameters:
669 + prob - The PetscDS object
670 - f - The field number
671 
672   Output Parameter:
673 . disc - The discretization object
674 
675   Level: beginner
676 
677 .seealso: PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
678 @*/
679 PetscErrorCode PetscDSGetDiscretization(PetscDS prob, PetscInt f, PetscObject *disc)
680 {
681   PetscFunctionBegin;
682   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
683   PetscValidPointer(disc, 3);
684   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);
685   *disc = prob->disc[f];
686   PetscFunctionReturn(0);
687 }
688 
689 #undef __FUNCT__
690 #define __FUNCT__ "PetscDSGetBdDiscretization"
691 /*@
692   PetscDSGetBdDiscretization - Returns the boundary discretization object for the given field
693 
694   Not collective
695 
696   Input Parameters:
697 + prob - The PetscDS object
698 - f - The field number
699 
700   Output Parameter:
701 . disc - The boundary discretization object
702 
703   Level: beginner
704 
705 .seealso: PetscDSSetBdDiscretization(), PetscDSAddBdDiscretization(), PetscDSGetDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
706 @*/
707 PetscErrorCode PetscDSGetBdDiscretization(PetscDS prob, PetscInt f, PetscObject *disc)
708 {
709   PetscFunctionBegin;
710   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
711   PetscValidPointer(disc, 3);
712   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);
713   *disc = prob->discBd[f];
714   PetscFunctionReturn(0);
715 }
716 
717 #undef __FUNCT__
718 #define __FUNCT__ "PetscDSSetDiscretization"
719 /*@
720   PetscDSSetDiscretization - Sets the discretization object for the given field
721 
722   Not collective
723 
724   Input Parameters:
725 + prob - The PetscDS object
726 . f - The field number
727 - disc - The discretization object
728 
729   Level: beginner
730 
731 .seealso: PetscDSGetDiscretization(), PetscDSAddDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
732 @*/
733 PetscErrorCode PetscDSSetDiscretization(PetscDS prob, PetscInt f, PetscObject disc)
734 {
735   PetscErrorCode ierr;
736 
737   PetscFunctionBegin;
738   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
739   PetscValidPointer(disc, 3);
740   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
741   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
742   if (prob->disc[f]) {ierr = PetscObjectDereference(prob->disc[f]);CHKERRQ(ierr);}
743   prob->disc[f] = disc;
744   ierr = PetscObjectReference(disc);CHKERRQ(ierr);
745   {
746     PetscClassId id;
747 
748     ierr = PetscObjectGetClassId(disc, &id);CHKERRQ(ierr);
749     if (id == PETSCFV_CLASSID) {
750       prob->implicit[f]      = PETSC_FALSE;
751       prob->adjacency[f*2+0] = PETSC_TRUE;
752       prob->adjacency[f*2+1] = PETSC_FALSE;
753     }
754   }
755   PetscFunctionReturn(0);
756 }
757 
758 #undef __FUNCT__
759 #define __FUNCT__ "PetscDSSetBdDiscretization"
760 /*@
761   PetscDSSetBdDiscretization - Sets the boundary discretization object for the given field
762 
763   Not collective
764 
765   Input Parameters:
766 + prob - The PetscDS object
767 . f - The field number
768 - disc - The boundary discretization object
769 
770   Level: beginner
771 
772 .seealso: PetscDSGetBdDiscretization(), PetscDSAddBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
773 @*/
774 PetscErrorCode PetscDSSetBdDiscretization(PetscDS prob, PetscInt f, PetscObject disc)
775 {
776   PetscErrorCode ierr;
777 
778   PetscFunctionBegin;
779   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
780   if (disc) PetscValidPointer(disc, 3);
781   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
782   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
783   if (prob->discBd[f]) {ierr = PetscObjectDereference(prob->discBd[f]);CHKERRQ(ierr);}
784   prob->discBd[f] = disc;
785   ierr = PetscObjectReference(disc);CHKERRQ(ierr);
786   PetscFunctionReturn(0);
787 }
788 
789 #undef __FUNCT__
790 #define __FUNCT__ "PetscDSAddDiscretization"
791 /*@
792   PetscDSAddDiscretization - Adds a 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: PetscDSGetDiscretization(), PetscDSSetDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
803 @*/
804 PetscErrorCode PetscDSAddDiscretization(PetscDS prob, PetscObject disc)
805 {
806   PetscErrorCode ierr;
807 
808   PetscFunctionBegin;
809   ierr = PetscDSSetDiscretization(prob, prob->Nf, disc);CHKERRQ(ierr);
810   PetscFunctionReturn(0);
811 }
812 
813 #undef __FUNCT__
814 #define __FUNCT__ "PetscDSAddBdDiscretization"
815 /*@
816   PetscDSAddBdDiscretization - Adds a boundary discretization object
817 
818   Not collective
819 
820   Input Parameters:
821 + prob - The PetscDS object
822 - disc - The boundary discretization object
823 
824   Level: beginner
825 
826 .seealso: PetscDSGetBdDiscretization(), PetscDSSetBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
827 @*/
828 PetscErrorCode PetscDSAddBdDiscretization(PetscDS prob, PetscObject disc)
829 {
830   PetscErrorCode ierr;
831 
832   PetscFunctionBegin;
833   ierr = PetscDSSetBdDiscretization(prob, prob->Nf, disc);CHKERRQ(ierr);
834   PetscFunctionReturn(0);
835 }
836 
837 #undef __FUNCT__
838 #define __FUNCT__ "PetscDSGetImplicit"
839 /*@
840   PetscDSGetImplicit - Returns the flag for implicit solve for this field. This is just a guide for IMEX
841 
842   Not collective
843 
844   Input Parameters:
845 + prob - The PetscDS object
846 - f - The field number
847 
848   Output Parameter:
849 . implicit - The flag indicating what kind of solve to use for this field
850 
851   Level: developer
852 
853 .seealso: PetscDSSetImplicit(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
854 @*/
855 PetscErrorCode PetscDSGetImplicit(PetscDS prob, PetscInt f, PetscBool *implicit)
856 {
857   PetscFunctionBegin;
858   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
859   PetscValidPointer(implicit, 3);
860   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);
861   *implicit = prob->implicit[f];
862   PetscFunctionReturn(0);
863 }
864 
865 #undef __FUNCT__
866 #define __FUNCT__ "PetscDSSetImplicit"
867 /*@
868   PetscDSSetImplicit - Set the flag for implicit solve for this field. This is just a guide for IMEX
869 
870   Not collective
871 
872   Input Parameters:
873 + prob - The PetscDS object
874 . f - The field number
875 - implicit - The flag indicating what kind of solve to use for this field
876 
877   Level: developer
878 
879 .seealso: PetscDSGetImplicit(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
880 @*/
881 PetscErrorCode PetscDSSetImplicit(PetscDS prob, PetscInt f, PetscBool implicit)
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   prob->implicit[f] = implicit;
887   PetscFunctionReturn(0);
888 }
889 
890 #undef __FUNCT__
891 #define __FUNCT__ "PetscDSGetAdjacency"
892 /*@
893   PetscDSGetAdjacency - Returns the flags for determining variable influence
894 
895   Not collective
896 
897   Input Parameters:
898 + prob - The PetscDS object
899 - f - The field number
900 
901   Output Parameter:
902 + useCone    - Flag for variable influence starting with the cone operation
903 - useClosure - Flag for variable influence using transitive closure
904 
905   Note: See the discussion in DMPlexGetAdjacencyUseCone() and DMPlexGetAdjacencyUseClosure()
906 
907   Level: developer
908 
909 .seealso: PetscDSSetAdjacency(), DMPlexGetAdjacencyUseCone(), DMPlexGetAdjacencyUseClosure(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
910 @*/
911 PetscErrorCode PetscDSGetAdjacency(PetscDS prob, PetscInt f, PetscBool *useCone, PetscBool *useClosure)
912 {
913   PetscFunctionBegin;
914   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
915   PetscValidPointer(useCone, 3);
916   PetscValidPointer(useClosure, 4);
917   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);
918   *useCone    = prob->adjacency[f*2+0];
919   *useClosure = prob->adjacency[f*2+1];
920   PetscFunctionReturn(0);
921 }
922 
923 #undef __FUNCT__
924 #define __FUNCT__ "PetscDSSetAdjacency"
925 /*@
926   PetscDSSetAdjacency - Set the flags for determining variable influence
927 
928   Not collective
929 
930   Input Parameters:
931 + prob - The PetscDS object
932 . f - The field number
933 . useCone    - Flag for variable influence starting with the cone operation
934 - useClosure - Flag for variable influence using transitive closure
935 
936   Note: See the discussion in DMPlexGetAdjacencyUseCone() and DMPlexGetAdjacencyUseClosure()
937 
938   Level: developer
939 
940 .seealso: PetscDSGetAdjacency(), DMPlexGetAdjacencyUseCone(), DMPlexGetAdjacencyUseClosure(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
941 @*/
942 PetscErrorCode PetscDSSetAdjacency(PetscDS prob, PetscInt f, PetscBool useCone, PetscBool useClosure)
943 {
944   PetscFunctionBegin;
945   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
946   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);
947   prob->adjacency[f*2+0] = useCone;
948   prob->adjacency[f*2+1] = useClosure;
949   PetscFunctionReturn(0);
950 }
951 
952 #undef __FUNCT__
953 #define __FUNCT__ "PetscDSGetObjective"
954 PetscErrorCode PetscDSGetObjective(PetscDS prob, PetscInt f,
955                                    void (**obj)(PetscInt dim, PetscInt Nf,
956                                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
957                                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
958                                                 const PetscReal t, const PetscReal x[], PetscScalar obj[]))
959 {
960   PetscFunctionBegin;
961   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
962   PetscValidPointer(obj, 2);
963   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);
964   *obj = prob->obj[f];
965   PetscFunctionReturn(0);
966 }
967 
968 #undef __FUNCT__
969 #define __FUNCT__ "PetscDSSetObjective"
970 PetscErrorCode PetscDSSetObjective(PetscDS prob, PetscInt f,
971                                    void (*obj)(PetscInt dim, PetscInt Nf,
972                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
973                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
974                                                const PetscReal t, const PetscReal x[], PetscScalar obj[]))
975 {
976   PetscErrorCode ierr;
977 
978   PetscFunctionBegin;
979   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
980   PetscValidFunction(obj, 2);
981   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
982   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
983   prob->obj[f] = obj;
984   PetscFunctionReturn(0);
985 }
986 
987 #undef __FUNCT__
988 #define __FUNCT__ "PetscDSGetResidual"
989 /*@C
990   PetscDSGetResidual - Get the pointwise residual function for a given test field
991 
992   Not collective
993 
994   Input Parameters:
995 + prob - The PetscDS
996 - f    - The test field number
997 
998   Output Parameters:
999 + f0 - integrand for the test function term
1000 - f1 - integrand for the test function gradient term
1001 
1002   Note: We are using a first order FEM model for the weak form:
1003 
1004   \int_\Omega \phi f_0(u, u_t, \nabla u, x, t) + \nabla\phi \cdot {\vec f}_1(u, u_t, \nabla u, x, t)
1005 
1006 The calling sequence for the callbacks f0 and f1 is given by:
1007 
1008 $ f0(PetscInt dim, PetscInt Nf,
1009 $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1010 $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1011 $    const PetscReal t, const PetscReal x[], PetscScalar f0[])
1012 
1013 + dim - the spatial dimension
1014 . Nf - the number of fields
1015 . uOff - the offset into u[] and u_t[] for each field
1016 . uOff_x - the offset into u_x[] for each field
1017 . u - each field evaluated at the current point
1018 . u_t - the time derivative of each field evaluated at the current point
1019 . u_x - the gradient of each field evaluated at the current point
1020 . aOff - the offset into a[] and a_t[] for each auxiliary field
1021 . aOff_x - the offset into a_x[] for each auxiliary field
1022 . a - each auxiliary field evaluated at the current point
1023 . a_t - the time derivative of each auxiliary field evaluated at the current point
1024 . a_x - the gradient of auxiliary each field evaluated at the current point
1025 . t - current time
1026 . x - coordinates of the current point
1027 - f0 - output values at the current point
1028 
1029   Level: intermediate
1030 
1031 .seealso: PetscDSSetResidual()
1032 @*/
1033 PetscErrorCode PetscDSGetResidual(PetscDS prob, PetscInt f,
1034                                   void (**f0)(PetscInt dim, PetscInt Nf,
1035                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1036                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1037                                               const PetscReal t, const PetscReal x[], PetscScalar f0[]),
1038                                   void (**f1)(PetscInt dim, PetscInt Nf,
1039                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1040                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1041                                               const PetscReal t, const PetscReal x[], PetscScalar f1[]))
1042 {
1043   PetscFunctionBegin;
1044   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1045   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);
1046   if (f0) {PetscValidPointer(f0, 3); *f0 = prob->f[f*2+0];}
1047   if (f1) {PetscValidPointer(f1, 4); *f1 = prob->f[f*2+1];}
1048   PetscFunctionReturn(0);
1049 }
1050 
1051 #undef __FUNCT__
1052 #define __FUNCT__ "PetscDSSetResidual"
1053 /*@C
1054   PetscDSSetResidual - Set the pointwise residual function for a given test field
1055 
1056   Not collective
1057 
1058   Input Parameters:
1059 + prob - The PetscDS
1060 . f    - The test field number
1061 . f0 - integrand for the test function term
1062 - f1 - integrand for the test function gradient term
1063 
1064   Note: We are using a first order FEM model for the weak form:
1065 
1066   \int_\Omega \phi f_0(u, u_t, \nabla u, x, t) + \nabla\phi \cdot {\vec f}_1(u, u_t, \nabla u, x, t)
1067 
1068 The calling sequence for the callbacks f0 and f1 is given by:
1069 
1070 $ f0(PetscInt dim, PetscInt Nf,
1071 $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1072 $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1073 $    const PetscReal t, const PetscReal x[], PetscScalar f0[])
1074 
1075 + dim - the spatial dimension
1076 . Nf - the number of fields
1077 . uOff - the offset into u[] and u_t[] for each field
1078 . uOff_x - the offset into u_x[] for each field
1079 . u - each field evaluated at the current point
1080 . u_t - the time derivative of each field evaluated at the current point
1081 . u_x - the gradient of each field evaluated at the current point
1082 . aOff - the offset into a[] and a_t[] for each auxiliary field
1083 . aOff_x - the offset into a_x[] for each auxiliary field
1084 . a - each auxiliary field evaluated at the current point
1085 . a_t - the time derivative of each auxiliary field evaluated at the current point
1086 . a_x - the gradient of auxiliary each field evaluated at the current point
1087 . t - current time
1088 . x - coordinates of the current point
1089 - f0 - output values at the current point
1090 
1091   Level: intermediate
1092 
1093 .seealso: PetscDSGetResidual()
1094 @*/
1095 PetscErrorCode PetscDSSetResidual(PetscDS prob, PetscInt f,
1096                                   void (*f0)(PetscInt dim, PetscInt Nf,
1097                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1098                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1099                                              const PetscReal t, const PetscReal x[], PetscScalar f0[]),
1100                                   void (*f1)(PetscInt dim, PetscInt Nf,
1101                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1102                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1103                                              const PetscReal t, const PetscReal x[], PetscScalar f1[]))
1104 {
1105   PetscErrorCode ierr;
1106 
1107   PetscFunctionBegin;
1108   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1109   PetscValidFunction(f0, 3);
1110   PetscValidFunction(f1, 4);
1111   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
1112   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
1113   prob->f[f*2+0] = f0;
1114   prob->f[f*2+1] = f1;
1115   PetscFunctionReturn(0);
1116 }
1117 
1118 #undef __FUNCT__
1119 #define __FUNCT__ "PetscDSGetJacobian"
1120 /*@C
1121   PetscDSGetJacobian - Get the pointwise Jacobian function for given test and basis field
1122 
1123   Not collective
1124 
1125   Input Parameters:
1126 + prob - The PetscDS
1127 . f    - The test field number
1128 - g    - The field number
1129 
1130   Output Parameters:
1131 + g0 - integrand for the test and basis function term
1132 . g1 - integrand for the test function and basis function gradient term
1133 . g2 - integrand for the test function gradient and basis function term
1134 - g3 - integrand for the test function gradient and basis function gradient term
1135 
1136   Note: We are using a first order FEM model for the weak form:
1137 
1138   \int_\Omega \phi g_0(u, u_t, \nabla u, x, t) \psi + \phi {\vec g}_1(u, u_t, \nabla u, x, t) \nabla \psi + \nabla\phi \cdot {\vec g}_2(u, u_t, \nabla u, x, t) \psi + \nabla\phi \cdot {\overleftrightarrow g}_3(u, u_t, \nabla u, x, t) \cdot \nabla \psi
1139 
1140 The calling sequence for the callbacks g0, g1, g2 and g3 is given by:
1141 
1142 $ g0(PetscInt dim, PetscInt Nf,
1143 $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1144 $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1145 $    const PetscReal t, const PetscReal x[], PetscScalar g0[])
1146 
1147 + dim - the spatial dimension
1148 . Nf - the number of fields
1149 . uOff - the offset into u[] and u_t[] for each field
1150 . uOff_x - the offset into u_x[] for each field
1151 . u - each field evaluated at the current point
1152 . u_t - the time derivative of each field evaluated at the current point
1153 . u_x - the gradient of each field evaluated at the current point
1154 . aOff - the offset into a[] and a_t[] for each auxiliary field
1155 . aOff_x - the offset into a_x[] for each auxiliary field
1156 . a - each auxiliary field evaluated at the current point
1157 . a_t - the time derivative of each auxiliary field evaluated at the current point
1158 . a_x - the gradient of auxiliary each field evaluated at the current point
1159 . t - current time
1160 . x - coordinates of the current point
1161 - g0 - output values at the current point
1162 
1163   Level: intermediate
1164 
1165 .seealso: PetscDSSetJacobian()
1166 @*/
1167 PetscErrorCode PetscDSGetJacobian(PetscDS prob, PetscInt f, PetscInt g,
1168                                   void (**g0)(PetscInt dim, PetscInt Nf,
1169                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1170                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1171                                               PetscReal t, const PetscReal x[], PetscScalar g0[]),
1172                                   void (**g1)(PetscInt dim, PetscInt Nf,
1173                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1174                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1175                                               PetscReal t, const PetscReal x[], PetscScalar g1[]),
1176                                   void (**g2)(PetscInt dim, PetscInt Nf,
1177                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1178                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1179                                               PetscReal t, const PetscReal x[], PetscScalar g2[]),
1180                                   void (**g3)(PetscInt dim, PetscInt Nf,
1181                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1182                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1183                                               PetscReal t, const PetscReal x[], PetscScalar g3[]))
1184 {
1185   PetscFunctionBegin;
1186   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1187   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);
1188   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);
1189   if (g0) {PetscValidPointer(g0, 4); *g0 = prob->g[(f*prob->Nf + g)*4+0];}
1190   if (g1) {PetscValidPointer(g1, 5); *g1 = prob->g[(f*prob->Nf + g)*4+1];}
1191   if (g2) {PetscValidPointer(g2, 6); *g2 = prob->g[(f*prob->Nf + g)*4+2];}
1192   if (g3) {PetscValidPointer(g3, 7); *g3 = prob->g[(f*prob->Nf + g)*4+3];}
1193   PetscFunctionReturn(0);
1194 }
1195 
1196 #undef __FUNCT__
1197 #define __FUNCT__ "PetscDSSetJacobian"
1198 /*@C
1199   PetscDSSetJacobian - Set the pointwise Jacobian function for given test and basis fields
1200 
1201   Not collective
1202 
1203   Input Parameters:
1204 + prob - The PetscDS
1205 . f    - The test field number
1206 . g    - The field number
1207 . g0 - integrand for the test and basis function term
1208 . g1 - integrand for the test function and basis function gradient term
1209 . g2 - integrand for the test function gradient and basis function term
1210 - g3 - integrand for the test function gradient and basis function gradient term
1211 
1212   Note: We are using a first order FEM model for the weak form:
1213 
1214   \int_\Omega \phi g_0(u, u_t, \nabla u, x, t) \psi + \phi {\vec g}_1(u, u_t, \nabla u, x, t) \nabla \psi + \nabla\phi \cdot {\vec g}_2(u, u_t, \nabla u, x, t) \psi + \nabla\phi \cdot {\overleftrightarrow g}_3(u, u_t, \nabla u, x, t) \cdot \nabla \psi
1215 
1216 The calling sequence for the callbacks g0, g1, g2 and g3 is given by:
1217 
1218 $ g0(PetscInt dim, PetscInt Nf,
1219 $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1220 $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1221 $    const PetscReal t, const PetscReal x[], PetscScalar g0[])
1222 
1223 + dim - the spatial dimension
1224 . Nf - the number of fields
1225 . uOff - the offset into u[] and u_t[] for each field
1226 . uOff_x - the offset into u_x[] for each field
1227 . u - each field evaluated at the current point
1228 . u_t - the time derivative of each field evaluated at the current point
1229 . u_x - the gradient of each field evaluated at the current point
1230 . aOff - the offset into a[] and a_t[] for each auxiliary field
1231 . aOff_x - the offset into a_x[] for each auxiliary field
1232 . a - each auxiliary field evaluated at the current point
1233 . a_t - the time derivative of each auxiliary field evaluated at the current point
1234 . a_x - the gradient of auxiliary each field evaluated at the current point
1235 . t - current time
1236 . x - coordinates of the current point
1237 - g0 - output values at the current point
1238 
1239   Level: intermediate
1240 
1241 .seealso: PetscDSGetJacobian()
1242 @*/
1243 PetscErrorCode PetscDSSetJacobian(PetscDS prob, PetscInt f, PetscInt g,
1244                                   void (*g0)(PetscInt dim, PetscInt Nf,
1245                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1246                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1247                                              const PetscReal t, const PetscReal x[], PetscScalar g0[]),
1248                                   void (*g1)(PetscInt dim, PetscInt Nf,
1249                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1250                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1251                                              const PetscReal t, const PetscReal x[], PetscScalar g1[]),
1252                                   void (*g2)(PetscInt dim, PetscInt Nf,
1253                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1254                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1255                                              const PetscReal t, const PetscReal x[], PetscScalar g2[]),
1256                                   void (*g3)(PetscInt dim, PetscInt Nf,
1257                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1258                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1259                                              const PetscReal t, const PetscReal x[], PetscScalar g3[]))
1260 {
1261   PetscErrorCode ierr;
1262 
1263   PetscFunctionBegin;
1264   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1265   if (g0) PetscValidFunction(g0, 4);
1266   if (g1) PetscValidFunction(g1, 5);
1267   if (g2) PetscValidFunction(g2, 6);
1268   if (g3) PetscValidFunction(g3, 7);
1269   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
1270   if (g < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", g);
1271   ierr = PetscDSEnlarge_Static(prob, PetscMax(f, g)+1);CHKERRQ(ierr);
1272   prob->g[(f*prob->Nf + g)*4+0] = g0;
1273   prob->g[(f*prob->Nf + g)*4+1] = g1;
1274   prob->g[(f*prob->Nf + g)*4+2] = g2;
1275   prob->g[(f*prob->Nf + g)*4+3] = g3;
1276   PetscFunctionReturn(0);
1277 }
1278 
1279 #undef __FUNCT__
1280 #define __FUNCT__ "PetscDSGetRiemannSolver"
1281 /*@C
1282   PetscDSGetRiemannSolver - Returns the Riemann solver for the given field
1283 
1284   Not collective
1285 
1286   Input Arguments:
1287 + prob - The PetscDS object
1288 - f    - The field number
1289 
1290   Output Argument:
1291 . r    - Riemann solver
1292 
1293   Calling sequence for r:
1294 
1295 $ r(const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx)
1296 
1297 + x    - The coordinates at a point on the interface
1298 . n    - The normal vector to the interface
1299 . uL   - The state vector to the left of the interface
1300 . uR   - The state vector to the right of the interface
1301 . flux - output array of flux through the interface
1302 - ctx  - optional user context
1303 
1304   Level: intermediate
1305 
1306 .seealso: PetscDSSetRiemannSolver()
1307 @*/
1308 PetscErrorCode PetscDSGetRiemannSolver(PetscDS prob, PetscInt f,
1309                                        void (**r)(const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx))
1310 {
1311   PetscFunctionBegin;
1312   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1313   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);
1314   PetscValidPointer(r, 3);
1315   *r = prob->r[f];
1316   PetscFunctionReturn(0);
1317 }
1318 
1319 #undef __FUNCT__
1320 #define __FUNCT__ "PetscDSSetRiemannSolver"
1321 /*@C
1322   PetscDSSetRiemannSolver - Sets the Riemann solver for the given field
1323 
1324   Not collective
1325 
1326   Input Arguments:
1327 + prob - The PetscDS object
1328 . f    - The field number
1329 - r    - Riemann solver
1330 
1331   Calling sequence for r:
1332 
1333 $ r(const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx)
1334 
1335 + x    - The coordinates at a point on the interface
1336 . n    - The normal vector to the interface
1337 . uL   - The state vector to the left of the interface
1338 . uR   - The state vector to the right of the interface
1339 . flux - output array of flux through the interface
1340 - ctx  - optional user context
1341 
1342   Level: intermediate
1343 
1344 .seealso: PetscDSGetRiemannSolver()
1345 @*/
1346 PetscErrorCode PetscDSSetRiemannSolver(PetscDS prob, PetscInt f,
1347                                        void (*r)(const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx))
1348 {
1349   PetscErrorCode ierr;
1350 
1351   PetscFunctionBegin;
1352   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1353   PetscValidFunction(r, 3);
1354   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
1355   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
1356   prob->r[f] = r;
1357   PetscFunctionReturn(0);
1358 }
1359 
1360 #undef __FUNCT__
1361 #define __FUNCT__ "PetscDSGetContext"
1362 PetscErrorCode PetscDSGetContext(PetscDS prob, PetscInt f, void **ctx)
1363 {
1364   PetscFunctionBegin;
1365   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1366   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);
1367   PetscValidPointer(ctx, 3);
1368   *ctx = prob->ctx[f];
1369   PetscFunctionReturn(0);
1370 }
1371 
1372 #undef __FUNCT__
1373 #define __FUNCT__ "PetscDSSetContext"
1374 PetscErrorCode PetscDSSetContext(PetscDS prob, PetscInt f, void *ctx)
1375 {
1376   PetscErrorCode ierr;
1377 
1378   PetscFunctionBegin;
1379   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1380   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
1381   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
1382   prob->ctx[f] = ctx;
1383   PetscFunctionReturn(0);
1384 }
1385 
1386 #undef __FUNCT__
1387 #define __FUNCT__ "PetscDSGetBdResidual"
1388 /*@C
1389   PetscDSGetBdResidual - Get the pointwise boundary residual function for a given test field
1390 
1391   Not collective
1392 
1393   Input Parameters:
1394 + prob - The PetscDS
1395 - f    - The test field number
1396 
1397   Output Parameters:
1398 + f0 - boundary integrand for the test function term
1399 - f1 - boundary integrand for the test function gradient term
1400 
1401   Note: We are using a first order FEM model for the weak form:
1402 
1403   \int_\Gamma \phi {\vec f}_0(u, u_t, \nabla u, x, t) \cdot \hat n + \nabla\phi \cdot {\overleftrightarrow f}_1(u, u_t, \nabla u, x, t) \cdot \hat n
1404 
1405 The calling sequence for the callbacks f0 and f1 is given by:
1406 
1407 $ f0(PetscInt dim, PetscInt Nf,
1408 $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1409 $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1410 $    const PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f0[])
1411 
1412 + dim - the spatial dimension
1413 . Nf - the number of fields
1414 . uOff - the offset into u[] and u_t[] for each field
1415 . uOff_x - the offset into u_x[] for each field
1416 . u - each field evaluated at the current point
1417 . u_t - the time derivative of each field evaluated at the current point
1418 . u_x - the gradient of each field evaluated at the current point
1419 . aOff - the offset into a[] and a_t[] for each auxiliary field
1420 . aOff_x - the offset into a_x[] for each auxiliary field
1421 . a - each auxiliary field evaluated at the current point
1422 . a_t - the time derivative of each auxiliary field evaluated at the current point
1423 . a_x - the gradient of auxiliary each field evaluated at the current point
1424 . t - current time
1425 . x - coordinates of the current point
1426 . n - unit normal at the current point
1427 - f0 - output values at the current point
1428 
1429   Level: intermediate
1430 
1431 .seealso: PetscDSSetBdResidual()
1432 @*/
1433 PetscErrorCode PetscDSGetBdResidual(PetscDS prob, PetscInt f,
1434                                     void (**f0)(PetscInt dim, PetscInt Nf,
1435                                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1436                                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1437                                                 const PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f0[]),
1438                                     void (**f1)(PetscInt dim, PetscInt Nf,
1439                                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1440                                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1441                                                 const PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f1[]))
1442 {
1443   PetscFunctionBegin;
1444   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1445   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);
1446   if (f0) {PetscValidPointer(f0, 3); *f0 = prob->fBd[f*2+0];}
1447   if (f1) {PetscValidPointer(f1, 4); *f1 = prob->fBd[f*2+1];}
1448   PetscFunctionReturn(0);
1449 }
1450 
1451 #undef __FUNCT__
1452 #define __FUNCT__ "PetscDSSetBdResidual"
1453 /*@C
1454   PetscDSSetBdResidual - Get the pointwise boundary residual function for a given test field
1455 
1456   Not collective
1457 
1458   Input Parameters:
1459 + prob - The PetscDS
1460 . f    - The test field number
1461 . f0 - boundary integrand for the test function term
1462 - f1 - boundary integrand for the test function gradient term
1463 
1464   Note: We are using a first order FEM model for the weak form:
1465 
1466   \int_\Gamma \phi {\vec f}_0(u, u_t, \nabla u, x, t) \cdot \hat n + \nabla\phi \cdot {\overleftrightarrow f}_1(u, u_t, \nabla u, x, t) \cdot \hat n
1467 
1468 The calling sequence for the callbacks f0 and f1 is given by:
1469 
1470 $ f0(PetscInt dim, PetscInt Nf,
1471 $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1472 $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1473 $    const PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f0[])
1474 
1475 + dim - the spatial dimension
1476 . Nf - the number of fields
1477 . uOff - the offset into u[] and u_t[] for each field
1478 . uOff_x - the offset into u_x[] for each field
1479 . u - each field evaluated at the current point
1480 . u_t - the time derivative of each field evaluated at the current point
1481 . u_x - the gradient of each field evaluated at the current point
1482 . aOff - the offset into a[] and a_t[] for each auxiliary field
1483 . aOff_x - the offset into a_x[] for each auxiliary field
1484 . a - each auxiliary field evaluated at the current point
1485 . a_t - the time derivative of each auxiliary field evaluated at the current point
1486 . a_x - the gradient of auxiliary each field evaluated at the current point
1487 . t - current time
1488 . x - coordinates of the current point
1489 . n - unit normal at the current point
1490 - f0 - output values at the current point
1491 
1492   Level: intermediate
1493 
1494 .seealso: PetscDSGetBdResidual()
1495 @*/
1496 PetscErrorCode PetscDSSetBdResidual(PetscDS prob, PetscInt f,
1497                                     void (*f0)(PetscInt dim, PetscInt Nf,
1498                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1499                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1500                                                const PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f0[]),
1501                                     void (*f1)(PetscInt dim, PetscInt Nf,
1502                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1503                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1504                                                const PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f1[]))
1505 {
1506   PetscErrorCode ierr;
1507 
1508   PetscFunctionBegin;
1509   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1510   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
1511   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
1512   if (f0) {PetscValidFunction(f0, 3); prob->fBd[f*2+0] = f0;}
1513   if (f1) {PetscValidFunction(f1, 4); prob->fBd[f*2+1] = f1;}
1514   PetscFunctionReturn(0);
1515 }
1516 
1517 #undef __FUNCT__
1518 #define __FUNCT__ "PetscDSGetBdJacobian"
1519 /*@C
1520   PetscDSGetBdJacobian - Get the pointwise boundary Jacobian function for given test and basis field
1521 
1522   Not collective
1523 
1524   Input Parameters:
1525 + prob - The PetscDS
1526 . f    - The test field number
1527 - g    - The field number
1528 
1529   Output Parameters:
1530 + g0 - integrand for the test and basis function term
1531 . g1 - integrand for the test function and basis function gradient term
1532 . g2 - integrand for the test function gradient and basis function term
1533 - g3 - integrand for the test function gradient and basis function gradient term
1534 
1535   Note: We are using a first order FEM model for the weak form:
1536 
1537   \int_\Gamma \phi {\vec g}_0(u, u_t, \nabla u, x, t) \cdot \hat n \psi + \phi {\vec g}_1(u, u_t, \nabla u, x, t) \cdot \hat n \nabla \psi + \nabla\phi \cdot {\vec g}_2(u, u_t, \nabla u, x, t) \cdot \hat n \psi + \nabla\phi \cdot {\overleftrightarrow g}_3(u, u_t, \nabla u, x, t) \cdot \hat n \cdot \nabla \psi
1538 
1539 The calling sequence for the callbacks g0, g1, g2 and g3 is given by:
1540 
1541 $ g0(PetscInt dim, PetscInt Nf,
1542 $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1543 $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1544 $    const PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar g0[])
1545 
1546 + dim - the spatial dimension
1547 . Nf - the number of fields
1548 . uOff - the offset into u[] and u_t[] for each field
1549 . uOff_x - the offset into u_x[] for each field
1550 . u - each field evaluated at the current point
1551 . u_t - the time derivative of each field evaluated at the current point
1552 . u_x - the gradient of each field evaluated at the current point
1553 . aOff - the offset into a[] and a_t[] for each auxiliary field
1554 . aOff_x - the offset into a_x[] for each auxiliary field
1555 . a - each auxiliary field evaluated at the current point
1556 . a_t - the time derivative of each auxiliary field evaluated at the current point
1557 . a_x - the gradient of auxiliary each field evaluated at the current point
1558 . t - current time
1559 . x - coordinates of the current point
1560 . n - normal at the current point
1561 - g0 - output values at the current point
1562 
1563   Level: intermediate
1564 
1565 .seealso: PetscDSSetBdJacobian()
1566 @*/
1567 PetscErrorCode PetscDSGetBdJacobian(PetscDS prob, PetscInt f, PetscInt g,
1568                                     void (**g0)(PetscInt dim, PetscInt Nf,
1569                                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1570                                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1571                                                 PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar g0[]),
1572                                     void (**g1)(PetscInt dim, PetscInt Nf,
1573                                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1574                                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1575                                                 PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar g1[]),
1576                                     void (**g2)(PetscInt dim, PetscInt Nf,
1577                                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1578                                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1579                                                 PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar g2[]),
1580                                     void (**g3)(PetscInt dim, PetscInt Nf,
1581                                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1582                                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1583                                                 PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar g3[]))
1584 {
1585   PetscFunctionBegin;
1586   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1587   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);
1588   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);
1589   if (g0) {PetscValidPointer(g0, 4); *g0 = prob->gBd[(f*prob->Nf + g)*4+0];}
1590   if (g1) {PetscValidPointer(g1, 5); *g1 = prob->gBd[(f*prob->Nf + g)*4+1];}
1591   if (g2) {PetscValidPointer(g2, 6); *g2 = prob->gBd[(f*prob->Nf + g)*4+2];}
1592   if (g3) {PetscValidPointer(g3, 7); *g3 = prob->gBd[(f*prob->Nf + g)*4+3];}
1593   PetscFunctionReturn(0);
1594 }
1595 
1596 #undef __FUNCT__
1597 #define __FUNCT__ "PetscDSSetBdJacobian"
1598 /*@C
1599   PetscDSSetBdJacobian - Set the pointwise boundary Jacobian function for given test and basis field
1600 
1601   Not collective
1602 
1603   Input Parameters:
1604 + prob - The PetscDS
1605 . f    - The test field number
1606 . g    - The field number
1607 . g0 - integrand for the test and basis function term
1608 . g1 - integrand for the test function and basis function gradient term
1609 . g2 - integrand for the test function gradient and basis function term
1610 - g3 - integrand for the test function gradient and basis function gradient term
1611 
1612   Note: We are using a first order FEM model for the weak form:
1613 
1614   \int_\Gamma \phi {\vec g}_0(u, u_t, \nabla u, x, t) \cdot \hat n \psi + \phi {\vec g}_1(u, u_t, \nabla u, x, t) \cdot \hat n \nabla \psi + \nabla\phi \cdot {\vec g}_2(u, u_t, \nabla u, x, t) \cdot \hat n \psi + \nabla\phi \cdot {\overleftrightarrow g}_3(u, u_t, \nabla u, x, t) \cdot \hat n \cdot \nabla \psi
1615 
1616 The calling sequence for the callbacks g0, g1, g2 and g3 is given by:
1617 
1618 $ g0(PetscInt dim, PetscInt Nf,
1619 $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1620 $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1621 $    const PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar g0[])
1622 
1623 + dim - the spatial dimension
1624 . Nf - the number of fields
1625 . uOff - the offset into u[] and u_t[] for each field
1626 . uOff_x - the offset into u_x[] for each field
1627 . u - each field evaluated at the current point
1628 . u_t - the time derivative of each field evaluated at the current point
1629 . u_x - the gradient of each field evaluated at the current point
1630 . aOff - the offset into a[] and a_t[] for each auxiliary field
1631 . aOff_x - the offset into a_x[] for each auxiliary field
1632 . a - each auxiliary field evaluated at the current point
1633 . a_t - the time derivative of each auxiliary field evaluated at the current point
1634 . a_x - the gradient of auxiliary each field evaluated at the current point
1635 . t - current time
1636 . x - coordinates of the current point
1637 . n - normal at the current point
1638 - g0 - output values at the current point
1639 
1640   Level: intermediate
1641 
1642 .seealso: PetscDSGetBdJacobian()
1643 @*/
1644 PetscErrorCode PetscDSSetBdJacobian(PetscDS prob, PetscInt f, PetscInt g,
1645                                     void (*g0)(PetscInt dim, PetscInt Nf,
1646                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1647                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1648                                                PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar g0[]),
1649                                     void (*g1)(PetscInt dim, PetscInt Nf,
1650                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1651                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1652                                                PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar g1[]),
1653                                     void (*g2)(PetscInt dim, PetscInt Nf,
1654                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1655                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1656                                                PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar g2[]),
1657                                     void (*g3)(PetscInt dim, PetscInt Nf,
1658                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1659                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1660                                                PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar g3[]))
1661 {
1662   PetscErrorCode ierr;
1663 
1664   PetscFunctionBegin;
1665   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1666   if (g0) PetscValidFunction(g0, 4);
1667   if (g1) PetscValidFunction(g1, 5);
1668   if (g2) PetscValidFunction(g2, 6);
1669   if (g3) PetscValidFunction(g3, 7);
1670   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
1671   if (g < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", g);
1672   ierr = PetscDSEnlarge_Static(prob, PetscMax(f, g)+1);CHKERRQ(ierr);
1673   prob->gBd[(f*prob->Nf + g)*4+0] = g0;
1674   prob->gBd[(f*prob->Nf + g)*4+1] = g1;
1675   prob->gBd[(f*prob->Nf + g)*4+2] = g2;
1676   prob->gBd[(f*prob->Nf + g)*4+3] = g3;
1677   PetscFunctionReturn(0);
1678 }
1679 
1680 #undef __FUNCT__
1681 #define __FUNCT__ "PetscDSGetFieldOffset"
1682 /*@
1683   PetscDSGetFieldOffset - Returns the offset of the given field in the full space basis
1684 
1685   Not collective
1686 
1687   Input Parameters:
1688 + prob - The PetscDS object
1689 - f - The field number
1690 
1691   Output Parameter:
1692 . off - The offset
1693 
1694   Level: beginner
1695 
1696 .seealso: PetscDSGetBdFieldOffset(), PetscDSGetNumFields(), PetscDSCreate()
1697 @*/
1698 PetscErrorCode PetscDSGetFieldOffset(PetscDS prob, PetscInt f, PetscInt *off)
1699 {
1700   PetscInt       g;
1701   PetscErrorCode ierr;
1702 
1703   PetscFunctionBegin;
1704   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1705   PetscValidPointer(off, 3);
1706   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);
1707   *off = 0;
1708   for (g = 0; g < f; ++g) {
1709     PetscFE  fe = (PetscFE) prob->disc[g];
1710     PetscInt Nb, Nc;
1711 
1712     ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
1713     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
1714     *off += Nb*Nc;
1715   }
1716   PetscFunctionReturn(0);
1717 }
1718 
1719 #undef __FUNCT__
1720 #define __FUNCT__ "PetscDSGetBdFieldOffset"
1721 /*@
1722   PetscDSGetBdFieldOffset - Returns the offset of the given field in the full space boundary basis
1723 
1724   Not collective
1725 
1726   Input Parameters:
1727 + prob - The PetscDS object
1728 - f - The field number
1729 
1730   Output Parameter:
1731 . off - The boundary offset
1732 
1733   Level: beginner
1734 
1735 .seealso: PetscDSGetFieldOffset(), PetscDSGetNumFields(), PetscDSCreate()
1736 @*/
1737 PetscErrorCode PetscDSGetBdFieldOffset(PetscDS prob, PetscInt f, PetscInt *off)
1738 {
1739   PetscInt       g;
1740   PetscErrorCode ierr;
1741 
1742   PetscFunctionBegin;
1743   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1744   PetscValidPointer(off, 3);
1745   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);
1746   *off = 0;
1747   for (g = 0; g < f; ++g) {
1748     PetscFE  fe = (PetscFE) prob->discBd[g];
1749     PetscInt Nb, Nc;
1750 
1751     ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
1752     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
1753     *off += Nb*Nc;
1754   }
1755   PetscFunctionReturn(0);
1756 }
1757 
1758 #undef __FUNCT__
1759 #define __FUNCT__ "PetscDSGetComponentOffset"
1760 /*@
1761   PetscDSGetComponentOffset - Returns the offset of the given field on an evaluation point
1762 
1763   Not collective
1764 
1765   Input Parameters:
1766 + prob - The PetscDS object
1767 - f - The field number
1768 
1769   Output Parameter:
1770 . off - The offset
1771 
1772   Level: beginner
1773 
1774 .seealso: PetscDSGetBdFieldOffset(), PetscDSGetNumFields(), PetscDSCreate()
1775 @*/
1776 PetscErrorCode PetscDSGetComponentOffset(PetscDS prob, PetscInt f, PetscInt *off)
1777 {
1778   PetscInt       g;
1779   PetscErrorCode ierr;
1780 
1781   PetscFunctionBegin;
1782   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1783   PetscValidPointer(off, 3);
1784   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);
1785   *off = 0;
1786   for (g = 0; g < f; ++g) {
1787     PetscFE  fe = (PetscFE) prob->disc[g];
1788     PetscInt Nc;
1789 
1790     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
1791     *off += Nc;
1792   }
1793   PetscFunctionReturn(0);
1794 }
1795 
1796 #undef __FUNCT__
1797 #define __FUNCT__ "PetscDSGetComponentOffsets"
1798 /*@
1799   PetscDSGetComponentOffsets - Returns the offset of each field on an evaluation point
1800 
1801   Not collective
1802 
1803   Input Parameter:
1804 . prob - The PetscDS object
1805 
1806   Output Parameter:
1807 . offsets - The offsets
1808 
1809   Level: beginner
1810 
1811 .seealso: PetscDSGetBdFieldOffset(), PetscDSGetNumFields(), PetscDSCreate()
1812 @*/
1813 PetscErrorCode PetscDSGetComponentOffsets(PetscDS prob, PetscInt *offsets[])
1814 {
1815   PetscInt       g;
1816   PetscErrorCode ierr;
1817 
1818   PetscFunctionBegin;
1819   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1820   PetscValidPointer(offsets, 2);
1821   *offsets = prob->off;
1822   PetscFunctionReturn(0);
1823 }
1824 
1825 #undef __FUNCT__
1826 #define __FUNCT__ "PetscDSGetComponentDerivativeOffsets"
1827 /*@
1828   PetscDSGetComponentDerivativeOffsets - Returns the offset of each field derivative on an evaluation point
1829 
1830   Not collective
1831 
1832   Input Parameter:
1833 . prob - The PetscDS object
1834 
1835   Output Parameter:
1836 . offsets - The offsets
1837 
1838   Level: beginner
1839 
1840 .seealso: PetscDSGetBdFieldOffset(), PetscDSGetNumFields(), PetscDSCreate()
1841 @*/
1842 PetscErrorCode PetscDSGetComponentDerivativeOffsets(PetscDS prob, PetscInt *offsets[])
1843 {
1844   PetscInt       g;
1845   PetscErrorCode ierr;
1846 
1847   PetscFunctionBegin;
1848   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1849   PetscValidPointer(offsets, 2);
1850   *offsets = prob->offDer;
1851   PetscFunctionReturn(0);
1852 }
1853 
1854 #undef __FUNCT__
1855 #define __FUNCT__ "PetscDSGetComponentBdOffsets"
1856 /*@
1857   PetscDSGetComponentBdOffsets - Returns the offset of each field on a boundary evaluation point
1858 
1859   Not collective
1860 
1861   Input Parameter:
1862 . prob - The PetscDS object
1863 
1864   Output Parameter:
1865 . offsets - The offsets
1866 
1867   Level: beginner
1868 
1869 .seealso: PetscDSGetBdFieldOffset(), PetscDSGetNumFields(), PetscDSCreate()
1870 @*/
1871 PetscErrorCode PetscDSGetComponentBdOffsets(PetscDS prob, PetscInt *offsets[])
1872 {
1873   PetscInt       g;
1874   PetscErrorCode ierr;
1875 
1876   PetscFunctionBegin;
1877   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1878   PetscValidPointer(offsets, 2);
1879   *offsets = prob->offBd;
1880   PetscFunctionReturn(0);
1881 }
1882 
1883 #undef __FUNCT__
1884 #define __FUNCT__ "PetscDSGetComponentBdDerivativeOffsets"
1885 /*@
1886   PetscDSGetComponentBdDerivativeOffsets - Returns the offset of each field derivative on a boundary evaluation point
1887 
1888   Not collective
1889 
1890   Input Parameter:
1891 . prob - The PetscDS object
1892 
1893   Output Parameter:
1894 . offsets - The offsets
1895 
1896   Level: beginner
1897 
1898 .seealso: PetscDSGetBdFieldOffset(), PetscDSGetNumFields(), PetscDSCreate()
1899 @*/
1900 PetscErrorCode PetscDSGetComponentBdDerivativeOffsets(PetscDS prob, PetscInt *offsets[])
1901 {
1902   PetscInt       g;
1903   PetscErrorCode ierr;
1904 
1905   PetscFunctionBegin;
1906   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1907   PetscValidPointer(offsets, 2);
1908   *offsets = prob->offDerBd;
1909   PetscFunctionReturn(0);
1910 }
1911 
1912 #undef __FUNCT__
1913 #define __FUNCT__ "PetscDSGetTabulation"
1914 /*@C
1915   PetscDSGetTabulation - Return the basis tabulation at quadrature points for the volume discretization
1916 
1917   Not collective
1918 
1919   Input Parameter:
1920 . prob - The PetscDS object
1921 
1922   Output Parameters:
1923 + basis - The basis function tabulation at quadrature points
1924 - basisDer - The basis function derivative tabulation at quadrature points
1925 
1926   Level: intermediate
1927 
1928 .seealso: PetscDSGetBdTabulation(), PetscDSCreate()
1929 @*/
1930 PetscErrorCode PetscDSGetTabulation(PetscDS prob, PetscReal ***basis, PetscReal ***basisDer)
1931 {
1932   PetscErrorCode ierr;
1933 
1934   PetscFunctionBegin;
1935   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1936   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
1937   if (basis)    {PetscValidPointer(basis, 2);    *basis    = prob->basis;}
1938   if (basisDer) {PetscValidPointer(basisDer, 3); *basisDer = prob->basisDer;}
1939   PetscFunctionReturn(0);
1940 }
1941 
1942 #undef __FUNCT__
1943 #define __FUNCT__ "PetscDSGetBdTabulation"
1944 /*@C
1945   PetscDSGetBdTabulation - Return the basis tabulation at quadrature points for the boundary discretization
1946 
1947   Not collective
1948 
1949   Input Parameter:
1950 . prob - The PetscDS object
1951 
1952   Output Parameters:
1953 + basis - The basis function tabulation at quadrature points
1954 - basisDer - The basis function derivative tabulation at quadrature points
1955 
1956   Level: intermediate
1957 
1958 .seealso: PetscDSGetTabulation(), PetscDSCreate()
1959 @*/
1960 PetscErrorCode PetscDSGetBdTabulation(PetscDS prob, PetscReal ***basis, PetscReal ***basisDer)
1961 {
1962   PetscErrorCode ierr;
1963 
1964   PetscFunctionBegin;
1965   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1966   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
1967   if (basis)    {PetscValidPointer(basis, 2);    *basis    = prob->basisBd;}
1968   if (basisDer) {PetscValidPointer(basisDer, 3); *basisDer = prob->basisDerBd;}
1969   PetscFunctionReturn(0);
1970 }
1971 
1972 #undef __FUNCT__
1973 #define __FUNCT__ "PetscDSGetEvaluationArrays"
1974 PetscErrorCode PetscDSGetEvaluationArrays(PetscDS prob, PetscScalar **u, PetscScalar **u_t, PetscScalar **u_x)
1975 {
1976   PetscErrorCode ierr;
1977 
1978   PetscFunctionBegin;
1979   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1980   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
1981   if (u)   {PetscValidPointer(u, 2);   *u   = prob->u;}
1982   if (u_t) {PetscValidPointer(u_t, 3); *u_t = prob->u_t;}
1983   if (u_x) {PetscValidPointer(u_x, 4); *u_x = prob->u_x;}
1984   PetscFunctionReturn(0);
1985 }
1986 
1987 #undef __FUNCT__
1988 #define __FUNCT__ "PetscDSGetWeakFormArrays"
1989 PetscErrorCode PetscDSGetWeakFormArrays(PetscDS prob, PetscScalar **f0, PetscScalar **f1, PetscScalar **g0, PetscScalar **g1, PetscScalar **g2, PetscScalar **g3)
1990 {
1991   PetscErrorCode ierr;
1992 
1993   PetscFunctionBegin;
1994   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1995   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
1996   if (f0) {PetscValidPointer(f0, 2); *f0 = prob->f0;}
1997   if (f1) {PetscValidPointer(f1, 3); *f1 = prob->f1;}
1998   if (g0) {PetscValidPointer(g0, 4); *g0 = prob->g0;}
1999   if (g1) {PetscValidPointer(g1, 5); *g1 = prob->g1;}
2000   if (g2) {PetscValidPointer(g2, 6); *g2 = prob->g2;}
2001   if (g3) {PetscValidPointer(g3, 7); *g3 = prob->g3;}
2002   PetscFunctionReturn(0);
2003 }
2004 
2005 #undef __FUNCT__
2006 #define __FUNCT__ "PetscDSGetRefCoordArrays"
2007 PetscErrorCode PetscDSGetRefCoordArrays(PetscDS prob, PetscReal **x, PetscScalar **refSpaceDer)
2008 {
2009   PetscErrorCode ierr;
2010 
2011   PetscFunctionBegin;
2012   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
2013   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
2014   if (x)           {PetscValidPointer(x, 2);           *x           = prob->x;}
2015   if (refSpaceDer) {PetscValidPointer(refSpaceDer, 3); *refSpaceDer = prob->refSpaceDer;}
2016   PetscFunctionReturn(0);
2017 }
2018 
2019 #undef __FUNCT__
2020 #define __FUNCT__ "PetscDSDestroy_Basic"
2021 static PetscErrorCode PetscDSDestroy_Basic(PetscDS prob)
2022 {
2023   PetscFunctionBegin;
2024   PetscFunctionReturn(0);
2025 }
2026 
2027 #undef __FUNCT__
2028 #define __FUNCT__ "PetscDSInitialize_Basic"
2029 static PetscErrorCode PetscDSInitialize_Basic(PetscDS prob)
2030 {
2031   PetscFunctionBegin;
2032   prob->ops->setfromoptions = NULL;
2033   prob->ops->setup          = NULL;
2034   prob->ops->view           = NULL;
2035   prob->ops->destroy        = PetscDSDestroy_Basic;
2036   PetscFunctionReturn(0);
2037 }
2038 
2039 /*MC
2040   PETSCDSBASIC = "basic" - A discrete system with pointwise residual and boundary residual functions
2041 
2042   Level: intermediate
2043 
2044 .seealso: PetscDSType, PetscDSCreate(), PetscDSSetType()
2045 M*/
2046 
2047 #undef __FUNCT__
2048 #define __FUNCT__ "PetscDSCreate_Basic"
2049 PETSC_EXTERN PetscErrorCode PetscDSCreate_Basic(PetscDS prob)
2050 {
2051   PetscDS_Basic *b;
2052   PetscErrorCode      ierr;
2053 
2054   PetscFunctionBegin;
2055   PetscValidHeaderSpecific(prob, PETSCSPACE_CLASSID, 1);
2056   ierr       = PetscNewLog(prob, &b);CHKERRQ(ierr);
2057   prob->data = b;
2058 
2059   ierr = PetscDSInitialize_Basic(prob);CHKERRQ(ierr);
2060   PetscFunctionReturn(0);
2061 }
2062