xref: /petsc/src/dm/dt/interface/dtds.c (revision bebe2cf65d55febe21a5af8db2bd2e168caaa2e7)
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;
387   PetscPointJac    *tmpg;
388   PetscBdPointFunc *tmpfbd;
389   PetscBdPointJac  *tmpgbd;
390   PetscRiemannFunc *tmpr;
391   void            **tmpctx;
392   PetscInt          Nf = prob->Nf, f, i;
393   PetscErrorCode    ierr;
394 
395   PetscFunctionBegin;
396   if (Nf >= NfNew) PetscFunctionReturn(0);
397   prob->setup = PETSC_FALSE;
398   ierr = PetscDSDestroyStructs_Static(prob);CHKERRQ(ierr);
399   ierr = PetscMalloc4(NfNew, &tmpd, NfNew, &tmpdbd, NfNew, &tmpi, NfNew*2, &tmpa);CHKERRQ(ierr);
400   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];}
401   for (f = Nf; f < NfNew; ++f) {ierr = PetscContainerCreate(PetscObjectComm((PetscObject) prob), (PetscContainer *) &tmpd[f]);CHKERRQ(ierr);
402     tmpdbd[f] = NULL; tmpi[f] = PETSC_TRUE; tmpa[f*2+0] = PETSC_FALSE; tmpa[f*2+1] = PETSC_TRUE;}
403   ierr = PetscFree4(prob->disc, prob->discBd, prob->implicit, prob->adjacency);CHKERRQ(ierr);
404   prob->Nf        = NfNew;
405   prob->disc      = tmpd;
406   prob->discBd    = tmpdbd;
407   prob->implicit  = tmpi;
408   prob->adjacency = tmpa;
409   ierr = PetscCalloc5(NfNew, &tmpobj, NfNew*2, &tmpf, NfNew*NfNew*4, &tmpg, NfNew, &tmpr, NfNew, &tmpctx);CHKERRQ(ierr);
410   for (f = 0; f < Nf; ++f) tmpobj[f] = prob->obj[f];
411   for (f = 0; f < Nf*2; ++f) tmpf[f] = prob->f[f];
412   for (f = 0; f < Nf*Nf*4; ++f) tmpg[f] = prob->g[f];
413   for (f = 0; f < Nf; ++f) tmpr[f] = prob->r[f];
414   for (f = 0; f < Nf; ++f) tmpctx[f] = prob->ctx[f];
415   for (f = Nf; f < NfNew; ++f) tmpobj[f] = NULL;
416   for (f = Nf*2; f < NfNew*2; ++f) tmpf[f] = NULL;
417   for (f = Nf*Nf*4; f < NfNew*NfNew*4; ++f) tmpg[f] = NULL;
418   for (f = Nf; f < NfNew; ++f) tmpr[f] = NULL;
419   for (f = Nf; f < NfNew; ++f) tmpctx[f] = NULL;
420   ierr = PetscFree5(prob->obj, prob->f, prob->g, prob->r, prob->ctx);CHKERRQ(ierr);
421   prob->obj = tmpobj;
422   prob->f   = tmpf;
423   prob->g   = tmpg;
424   prob->r   = tmpr;
425   prob->ctx = tmpctx;
426   ierr = PetscCalloc2(NfNew*2, &tmpfbd, NfNew*NfNew*4, &tmpgbd);CHKERRQ(ierr);
427   for (f = 0; f < Nf*2; ++f) tmpfbd[f] = prob->fBd[f];
428   for (f = 0; f < Nf*Nf*4; ++f) tmpgbd[f] = prob->gBd[f];
429   for (f = Nf*2; f < NfNew*2; ++f) tmpfbd[f] = NULL;
430   for (f = Nf*Nf*4; f < NfNew*NfNew*4; ++f) tmpgbd[f] = NULL;
431   ierr = PetscFree2(prob->fBd, prob->gBd);CHKERRQ(ierr);
432   prob->fBd = tmpfbd;
433   prob->gBd = tmpgbd;
434   PetscFunctionReturn(0);
435 }
436 
437 #undef __FUNCT__
438 #define __FUNCT__ "PetscDSDestroy"
439 /*@
440   PetscDSDestroy - Destroys a PetscDS object
441 
442   Collective on PetscDS
443 
444   Input Parameter:
445 . prob - the PetscDS object to destroy
446 
447   Level: developer
448 
449 .seealso PetscDSView()
450 @*/
451 PetscErrorCode PetscDSDestroy(PetscDS *prob)
452 {
453   PetscInt       f;
454   PetscErrorCode ierr;
455 
456   PetscFunctionBegin;
457   if (!*prob) PetscFunctionReturn(0);
458   PetscValidHeaderSpecific((*prob), PETSCDS_CLASSID, 1);
459 
460   if (--((PetscObject)(*prob))->refct > 0) {*prob = 0; PetscFunctionReturn(0);}
461   ((PetscObject) (*prob))->refct = 0;
462   ierr = PetscDSDestroyStructs_Static(*prob);CHKERRQ(ierr);
463   for (f = 0; f < (*prob)->Nf; ++f) {
464     ierr = PetscObjectDereference((*prob)->disc[f]);CHKERRQ(ierr);
465     ierr = PetscObjectDereference((*prob)->discBd[f]);CHKERRQ(ierr);
466   }
467   ierr = PetscFree4((*prob)->disc, (*prob)->discBd, (*prob)->implicit, (*prob)->adjacency);CHKERRQ(ierr);
468   ierr = PetscFree5((*prob)->obj,(*prob)->f,(*prob)->g,(*prob)->r,(*prob)->ctx);CHKERRQ(ierr);
469   ierr = PetscFree2((*prob)->fBd,(*prob)->gBd);CHKERRQ(ierr);
470   if ((*prob)->ops->destroy) {ierr = (*(*prob)->ops->destroy)(*prob);CHKERRQ(ierr);}
471   ierr = PetscHeaderDestroy(prob);CHKERRQ(ierr);
472   PetscFunctionReturn(0);
473 }
474 
475 #undef __FUNCT__
476 #define __FUNCT__ "PetscDSCreate"
477 /*@
478   PetscDSCreate - Creates an empty PetscDS object. The type can then be set with PetscDSSetType().
479 
480   Collective on MPI_Comm
481 
482   Input Parameter:
483 . comm - The communicator for the PetscDS object
484 
485   Output Parameter:
486 . prob - The PetscDS object
487 
488   Level: beginner
489 
490 .seealso: PetscDSSetType(), PETSCDSBASIC
491 @*/
492 PetscErrorCode PetscDSCreate(MPI_Comm comm, PetscDS *prob)
493 {
494   PetscDS   p;
495   PetscErrorCode ierr;
496 
497   PetscFunctionBegin;
498   PetscValidPointer(prob, 2);
499   *prob  = NULL;
500   ierr = PetscDSInitializePackage();CHKERRQ(ierr);
501 
502   ierr = PetscHeaderCreate(p, PETSCDS_CLASSID, "PetscDS", "Discrete System", "PetscDS", comm, PetscDSDestroy, PetscDSView);CHKERRQ(ierr);
503 
504   p->Nf    = 0;
505   p->setup = PETSC_FALSE;
506 
507   *prob = p;
508   PetscFunctionReturn(0);
509 }
510 
511 #undef __FUNCT__
512 #define __FUNCT__ "PetscDSGetNumFields"
513 /*@
514   PetscDSGetNumFields - Returns the number of fields in the DS
515 
516   Not collective
517 
518   Input Parameter:
519 . prob - The PetscDS object
520 
521   Output Parameter:
522 . Nf - The number of fields
523 
524   Level: beginner
525 
526 .seealso: PetscDSGetSpatialDimension(), PetscDSCreate()
527 @*/
528 PetscErrorCode PetscDSGetNumFields(PetscDS prob, PetscInt *Nf)
529 {
530   PetscFunctionBegin;
531   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
532   PetscValidPointer(Nf, 2);
533   *Nf = prob->Nf;
534   PetscFunctionReturn(0);
535 }
536 
537 #undef __FUNCT__
538 #define __FUNCT__ "PetscDSGetSpatialDimension"
539 /*@
540   PetscDSGetSpatialDimension - Returns the spatial dimension of the DS
541 
542   Not collective
543 
544   Input Parameter:
545 . prob - The PetscDS object
546 
547   Output Parameter:
548 . dim - The spatial dimension
549 
550   Level: beginner
551 
552 .seealso: PetscDSGetNumFields(), PetscDSCreate()
553 @*/
554 PetscErrorCode PetscDSGetSpatialDimension(PetscDS prob, PetscInt *dim)
555 {
556   PetscErrorCode ierr;
557 
558   PetscFunctionBegin;
559   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
560   PetscValidPointer(dim, 2);
561   *dim = 0;
562   if (prob->Nf) {
563     PetscObject  obj;
564     PetscClassId id;
565 
566     ierr = PetscDSGetDiscretization(prob, 0, &obj);CHKERRQ(ierr);
567     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
568     if (id == PETSCFE_CLASSID)      {ierr = PetscFEGetSpatialDimension((PetscFE) obj, dim);CHKERRQ(ierr);}
569     else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetSpatialDimension((PetscFV) obj, dim);CHKERRQ(ierr);}
570     else SETERRQ1(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", 0);
571   }
572   PetscFunctionReturn(0);
573 }
574 
575 #undef __FUNCT__
576 #define __FUNCT__ "PetscDSGetTotalDimension"
577 /*@
578   PetscDSGetTotalDimension - Returns the total size of the approximation space for this system
579 
580   Not collective
581 
582   Input Parameter:
583 . prob - The PetscDS object
584 
585   Output Parameter:
586 . dim - The total problem dimension
587 
588   Level: beginner
589 
590 .seealso: PetscDSGetNumFields(), PetscDSCreate()
591 @*/
592 PetscErrorCode PetscDSGetTotalDimension(PetscDS prob, PetscInt *dim)
593 {
594   PetscErrorCode ierr;
595 
596   PetscFunctionBegin;
597   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
598   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
599   PetscValidPointer(dim, 2);
600   *dim = prob->totDim;
601   PetscFunctionReturn(0);
602 }
603 
604 #undef __FUNCT__
605 #define __FUNCT__ "PetscDSGetTotalBdDimension"
606 /*@
607   PetscDSGetTotalBdDimension - Returns the total size of the boundary approximation space for this system
608 
609   Not collective
610 
611   Input Parameter:
612 . prob - The PetscDS object
613 
614   Output Parameter:
615 . dim - The total boundary problem dimension
616 
617   Level: beginner
618 
619 .seealso: PetscDSGetNumFields(), PetscDSCreate()
620 @*/
621 PetscErrorCode PetscDSGetTotalBdDimension(PetscDS prob, PetscInt *dim)
622 {
623   PetscErrorCode ierr;
624 
625   PetscFunctionBegin;
626   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
627   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
628   PetscValidPointer(dim, 2);
629   *dim = prob->totDimBd;
630   PetscFunctionReturn(0);
631 }
632 
633 #undef __FUNCT__
634 #define __FUNCT__ "PetscDSGetTotalComponents"
635 /*@
636   PetscDSGetTotalComponents - Returns the total number of components in this system
637 
638   Not collective
639 
640   Input Parameter:
641 . prob - The PetscDS object
642 
643   Output Parameter:
644 . dim - The total number of components
645 
646   Level: beginner
647 
648 .seealso: PetscDSGetNumFields(), PetscDSCreate()
649 @*/
650 PetscErrorCode PetscDSGetTotalComponents(PetscDS prob, PetscInt *Nc)
651 {
652   PetscErrorCode ierr;
653 
654   PetscFunctionBegin;
655   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
656   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
657   PetscValidPointer(Nc, 2);
658   *Nc = prob->totComp;
659   PetscFunctionReturn(0);
660 }
661 
662 #undef __FUNCT__
663 #define __FUNCT__ "PetscDSGetDiscretization"
664 /*@
665   PetscDSGetDiscretization - Returns the discretization object for the given field
666 
667   Not collective
668 
669   Input Parameters:
670 + prob - The PetscDS object
671 - f - The field number
672 
673   Output Parameter:
674 . disc - The discretization object
675 
676   Level: beginner
677 
678 .seealso: PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
679 @*/
680 PetscErrorCode PetscDSGetDiscretization(PetscDS prob, PetscInt f, PetscObject *disc)
681 {
682   PetscFunctionBegin;
683   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
684   PetscValidPointer(disc, 3);
685   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);
686   *disc = prob->disc[f];
687   PetscFunctionReturn(0);
688 }
689 
690 #undef __FUNCT__
691 #define __FUNCT__ "PetscDSGetBdDiscretization"
692 /*@
693   PetscDSGetBdDiscretization - Returns the boundary discretization object for the given field
694 
695   Not collective
696 
697   Input Parameters:
698 + prob - The PetscDS object
699 - f - The field number
700 
701   Output Parameter:
702 . disc - The boundary discretization object
703 
704   Level: beginner
705 
706 .seealso: PetscDSSetBdDiscretization(), PetscDSAddBdDiscretization(), PetscDSGetDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
707 @*/
708 PetscErrorCode PetscDSGetBdDiscretization(PetscDS prob, PetscInt f, PetscObject *disc)
709 {
710   PetscFunctionBegin;
711   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
712   PetscValidPointer(disc, 3);
713   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);
714   *disc = prob->discBd[f];
715   PetscFunctionReturn(0);
716 }
717 
718 #undef __FUNCT__
719 #define __FUNCT__ "PetscDSSetDiscretization"
720 /*@
721   PetscDSSetDiscretization - Sets the discretization object for the given field
722 
723   Not collective
724 
725   Input Parameters:
726 + prob - The PetscDS object
727 . f - The field number
728 - disc - The discretization object
729 
730   Level: beginner
731 
732 .seealso: PetscDSGetDiscretization(), PetscDSAddDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
733 @*/
734 PetscErrorCode PetscDSSetDiscretization(PetscDS prob, PetscInt f, PetscObject disc)
735 {
736   PetscErrorCode ierr;
737 
738   PetscFunctionBegin;
739   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
740   PetscValidPointer(disc, 3);
741   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
742   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
743   if (prob->disc[f]) {ierr = PetscObjectDereference(prob->disc[f]);CHKERRQ(ierr);}
744   prob->disc[f] = disc;
745   ierr = PetscObjectReference(disc);CHKERRQ(ierr);
746   {
747     PetscClassId id;
748 
749     ierr = PetscObjectGetClassId(disc, &id);CHKERRQ(ierr);
750     if (id == PETSCFV_CLASSID) {
751       prob->implicit[f]      = PETSC_FALSE;
752       prob->adjacency[f*2+0] = PETSC_TRUE;
753       prob->adjacency[f*2+1] = PETSC_FALSE;
754     }
755   }
756   PetscFunctionReturn(0);
757 }
758 
759 #undef __FUNCT__
760 #define __FUNCT__ "PetscDSSetBdDiscretization"
761 /*@
762   PetscDSSetBdDiscretization - Sets the boundary discretization object for the given field
763 
764   Not collective
765 
766   Input Parameters:
767 + prob - The PetscDS object
768 . f - The field number
769 - disc - The boundary discretization object
770 
771   Level: beginner
772 
773 .seealso: PetscDSGetBdDiscretization(), PetscDSAddBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
774 @*/
775 PetscErrorCode PetscDSSetBdDiscretization(PetscDS prob, PetscInt f, PetscObject disc)
776 {
777   PetscErrorCode ierr;
778 
779   PetscFunctionBegin;
780   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
781   if (disc) PetscValidPointer(disc, 3);
782   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
783   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
784   if (prob->discBd[f]) {ierr = PetscObjectDereference(prob->discBd[f]);CHKERRQ(ierr);}
785   prob->discBd[f] = disc;
786   ierr = PetscObjectReference(disc);CHKERRQ(ierr);
787   PetscFunctionReturn(0);
788 }
789 
790 #undef __FUNCT__
791 #define __FUNCT__ "PetscDSAddDiscretization"
792 /*@
793   PetscDSAddDiscretization - Adds a discretization object
794 
795   Not collective
796 
797   Input Parameters:
798 + prob - The PetscDS object
799 - disc - The boundary discretization object
800 
801   Level: beginner
802 
803 .seealso: PetscDSGetDiscretization(), PetscDSSetDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
804 @*/
805 PetscErrorCode PetscDSAddDiscretization(PetscDS prob, PetscObject disc)
806 {
807   PetscErrorCode ierr;
808 
809   PetscFunctionBegin;
810   ierr = PetscDSSetDiscretization(prob, prob->Nf, disc);CHKERRQ(ierr);
811   PetscFunctionReturn(0);
812 }
813 
814 #undef __FUNCT__
815 #define __FUNCT__ "PetscDSAddBdDiscretization"
816 /*@
817   PetscDSAddBdDiscretization - Adds a boundary discretization object
818 
819   Not collective
820 
821   Input Parameters:
822 + prob - The PetscDS object
823 - disc - The boundary discretization object
824 
825   Level: beginner
826 
827 .seealso: PetscDSGetBdDiscretization(), PetscDSSetBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
828 @*/
829 PetscErrorCode PetscDSAddBdDiscretization(PetscDS prob, PetscObject disc)
830 {
831   PetscErrorCode ierr;
832 
833   PetscFunctionBegin;
834   ierr = PetscDSSetBdDiscretization(prob, prob->Nf, disc);CHKERRQ(ierr);
835   PetscFunctionReturn(0);
836 }
837 
838 #undef __FUNCT__
839 #define __FUNCT__ "PetscDSGetImplicit"
840 /*@
841   PetscDSGetImplicit - Returns the flag for implicit solve for this field. This is just a guide for IMEX
842 
843   Not collective
844 
845   Input Parameters:
846 + prob - The PetscDS object
847 - f - The field number
848 
849   Output Parameter:
850 . implicit - The flag indicating what kind of solve to use for this field
851 
852   Level: developer
853 
854 .seealso: PetscDSSetImplicit(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
855 @*/
856 PetscErrorCode PetscDSGetImplicit(PetscDS prob, PetscInt f, PetscBool *implicit)
857 {
858   PetscFunctionBegin;
859   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
860   PetscValidPointer(implicit, 3);
861   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);
862   *implicit = prob->implicit[f];
863   PetscFunctionReturn(0);
864 }
865 
866 #undef __FUNCT__
867 #define __FUNCT__ "PetscDSSetImplicit"
868 /*@
869   PetscDSSetImplicit - Set the flag for implicit solve for this field. This is just a guide for IMEX
870 
871   Not collective
872 
873   Input Parameters:
874 + prob - The PetscDS object
875 . f - The field number
876 - implicit - The flag indicating what kind of solve to use for this field
877 
878   Level: developer
879 
880 .seealso: PetscDSGetImplicit(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
881 @*/
882 PetscErrorCode PetscDSSetImplicit(PetscDS prob, PetscInt f, PetscBool implicit)
883 {
884   PetscFunctionBegin;
885   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
886   if ((f < 0) || (f >= prob->Nf)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be in [0, %d)", f, prob->Nf);
887   prob->implicit[f] = implicit;
888   PetscFunctionReturn(0);
889 }
890 
891 #undef __FUNCT__
892 #define __FUNCT__ "PetscDSGetAdjacency"
893 /*@
894   PetscDSGetAdjacency - Returns the flags for determining variable influence
895 
896   Not collective
897 
898   Input Parameters:
899 + prob - The PetscDS object
900 - f - The field number
901 
902   Output Parameter:
903 + useCone    - Flag for variable influence starting with the cone operation
904 - useClosure - Flag for variable influence using transitive closure
905 
906   Note: See the discussion in DMPlexGetAdjacencyUseCone() and DMPlexGetAdjacencyUseClosure()
907 
908   Level: developer
909 
910 .seealso: PetscDSSetAdjacency(), DMPlexGetAdjacencyUseCone(), DMPlexGetAdjacencyUseClosure(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
911 @*/
912 PetscErrorCode PetscDSGetAdjacency(PetscDS prob, PetscInt f, PetscBool *useCone, PetscBool *useClosure)
913 {
914   PetscFunctionBegin;
915   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
916   PetscValidPointer(useCone, 3);
917   PetscValidPointer(useClosure, 4);
918   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);
919   *useCone    = prob->adjacency[f*2+0];
920   *useClosure = prob->adjacency[f*2+1];
921   PetscFunctionReturn(0);
922 }
923 
924 #undef __FUNCT__
925 #define __FUNCT__ "PetscDSSetAdjacency"
926 /*@
927   PetscDSSetAdjacency - Set the flags for determining variable influence
928 
929   Not collective
930 
931   Input Parameters:
932 + prob - The PetscDS object
933 . f - The field number
934 . useCone    - Flag for variable influence starting with the cone operation
935 - useClosure - Flag for variable influence using transitive closure
936 
937   Note: See the discussion in DMPlexGetAdjacencyUseCone() and DMPlexGetAdjacencyUseClosure()
938 
939   Level: developer
940 
941 .seealso: PetscDSGetAdjacency(), DMPlexGetAdjacencyUseCone(), DMPlexGetAdjacencyUseClosure(), PetscDSSetDiscretization(), PetscDSAddDiscretization(), PetscDSGetBdDiscretization(), PetscDSGetNumFields(), PetscDSCreate()
942 @*/
943 PetscErrorCode PetscDSSetAdjacency(PetscDS prob, PetscInt f, PetscBool useCone, PetscBool useClosure)
944 {
945   PetscFunctionBegin;
946   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
947   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);
948   prob->adjacency[f*2+0] = useCone;
949   prob->adjacency[f*2+1] = useClosure;
950   PetscFunctionReturn(0);
951 }
952 
953 #undef __FUNCT__
954 #define __FUNCT__ "PetscDSGetObjective"
955 PetscErrorCode PetscDSGetObjective(PetscDS prob, PetscInt f,
956                                    void (**obj)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
957                                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
958                                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
959                                                 PetscReal t, const PetscReal x[], PetscScalar obj[]))
960 {
961   PetscFunctionBegin;
962   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
963   PetscValidPointer(obj, 2);
964   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);
965   *obj = prob->obj[f];
966   PetscFunctionReturn(0);
967 }
968 
969 #undef __FUNCT__
970 #define __FUNCT__ "PetscDSSetObjective"
971 PetscErrorCode PetscDSSetObjective(PetscDS prob, PetscInt f,
972                                    void (*obj)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
973                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
974                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
975                                                PetscReal t, const PetscReal x[], PetscScalar obj[]))
976 {
977   PetscErrorCode ierr;
978 
979   PetscFunctionBegin;
980   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
981   PetscValidFunction(obj, 2);
982   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
983   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
984   prob->obj[f] = obj;
985   PetscFunctionReturn(0);
986 }
987 
988 #undef __FUNCT__
989 #define __FUNCT__ "PetscDSGetResidual"
990 /*@C
991   PetscDSGetResidual - Get the pointwise residual function for a given test field
992 
993   Not collective
994 
995   Input Parameters:
996 + prob - The PetscDS
997 - f    - The test field number
998 
999   Output Parameters:
1000 + f0 - integrand for the test function term
1001 - f1 - integrand for the test function gradient term
1002 
1003   Note: We are using a first order FEM model for the weak form:
1004 
1005   \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)
1006 
1007 The calling sequence for the callbacks f0 and f1 is given by:
1008 
1009 $ f0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1010 $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1011 $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1012 $    PetscReal t, const PetscReal x[], PetscScalar f0[])
1013 
1014 + dim - the spatial dimension
1015 . Nf - the number of fields
1016 . uOff - the offset into u[] and u_t[] for each field
1017 . uOff_x - the offset into u_x[] for each field
1018 . u - each field evaluated at the current point
1019 . u_t - the time derivative of each field evaluated at the current point
1020 . u_x - the gradient of each field evaluated at the current point
1021 . aOff - the offset into a[] and a_t[] for each auxiliary field
1022 . aOff_x - the offset into a_x[] for each auxiliary field
1023 . a - each auxiliary field evaluated at the current point
1024 . a_t - the time derivative of each auxiliary field evaluated at the current point
1025 . a_x - the gradient of auxiliary each field evaluated at the current point
1026 . t - current time
1027 . x - coordinates of the current point
1028 - f0 - output values at the current point
1029 
1030   Level: intermediate
1031 
1032 .seealso: PetscDSSetResidual()
1033 @*/
1034 PetscErrorCode PetscDSGetResidual(PetscDS prob, PetscInt f,
1035                                   void (**f0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1036                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1037                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1038                                               PetscReal t, const PetscReal x[], PetscScalar f0[]),
1039                                   void (**f1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1040                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1041                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1042                                               PetscReal t, const PetscReal x[], PetscScalar f1[]))
1043 {
1044   PetscFunctionBegin;
1045   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1046   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);
1047   if (f0) {PetscValidPointer(f0, 3); *f0 = prob->f[f*2+0];}
1048   if (f1) {PetscValidPointer(f1, 4); *f1 = prob->f[f*2+1];}
1049   PetscFunctionReturn(0);
1050 }
1051 
1052 #undef __FUNCT__
1053 #define __FUNCT__ "PetscDSSetResidual"
1054 /*@C
1055   PetscDSSetResidual - Set the pointwise residual function for a given test field
1056 
1057   Not collective
1058 
1059   Input Parameters:
1060 + prob - The PetscDS
1061 . f    - The test field number
1062 . f0 - integrand for the test function term
1063 - f1 - integrand for the test function gradient term
1064 
1065   Note: We are using a first order FEM model for the weak form:
1066 
1067   \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)
1068 
1069 The calling sequence for the callbacks f0 and f1 is given by:
1070 
1071 $ f0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1072 $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1073 $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1074 $    PetscReal t, const PetscReal x[], PetscScalar f0[])
1075 
1076 + dim - the spatial dimension
1077 . Nf - the number of fields
1078 . uOff - the offset into u[] and u_t[] for each field
1079 . uOff_x - the offset into u_x[] for each field
1080 . u - each field evaluated at the current point
1081 . u_t - the time derivative of each field evaluated at the current point
1082 . u_x - the gradient of each field evaluated at the current point
1083 . aOff - the offset into a[] and a_t[] for each auxiliary field
1084 . aOff_x - the offset into a_x[] for each auxiliary field
1085 . a - each auxiliary field evaluated at the current point
1086 . a_t - the time derivative of each auxiliary field evaluated at the current point
1087 . a_x - the gradient of auxiliary each field evaluated at the current point
1088 . t - current time
1089 . x - coordinates of the current point
1090 - f0 - output values at the current point
1091 
1092   Level: intermediate
1093 
1094 .seealso: PetscDSGetResidual()
1095 @*/
1096 PetscErrorCode PetscDSSetResidual(PetscDS prob, PetscInt f,
1097                                   void (*f0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1098                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1099                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1100                                              PetscReal t, const PetscReal x[], PetscScalar f0[]),
1101                                   void (*f1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1102                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1103                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1104                                              PetscReal t, const PetscReal x[], PetscScalar f1[]))
1105 {
1106   PetscErrorCode ierr;
1107 
1108   PetscFunctionBegin;
1109   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1110   if (f0) PetscValidFunction(f0, 3);
1111   if (f1) PetscValidFunction(f1, 4);
1112   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
1113   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
1114   prob->f[f*2+0] = f0;
1115   prob->f[f*2+1] = f1;
1116   PetscFunctionReturn(0);
1117 }
1118 
1119 #undef __FUNCT__
1120 #define __FUNCT__ "PetscDSGetJacobian"
1121 /*@C
1122   PetscDSGetJacobian - Get the pointwise Jacobian function for given test and basis field
1123 
1124   Not collective
1125 
1126   Input Parameters:
1127 + prob - The PetscDS
1128 . f    - The test field number
1129 - g    - The field number
1130 
1131   Output Parameters:
1132 + g0 - integrand for the test and basis function term
1133 . g1 - integrand for the test function and basis function gradient term
1134 . g2 - integrand for the test function gradient and basis function term
1135 - g3 - integrand for the test function gradient and basis function gradient term
1136 
1137   Note: We are using a first order FEM model for the weak form:
1138 
1139   \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
1140 
1141 The calling sequence for the callbacks g0, g1, g2 and g3 is given by:
1142 
1143 $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1144 $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1145 $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1146 $    PetscReal t, const PetscReal u_tShift, const PetscReal x[], PetscScalar g0[])
1147 
1148 + dim - the spatial dimension
1149 . Nf - the number of fields
1150 . uOff - the offset into u[] and u_t[] for each field
1151 . uOff_x - the offset into u_x[] for each field
1152 . u - each field evaluated at the current point
1153 . u_t - the time derivative of each field evaluated at the current point
1154 . u_x - the gradient of each field evaluated at the current point
1155 . aOff - the offset into a[] and a_t[] for each auxiliary field
1156 . aOff_x - the offset into a_x[] for each auxiliary field
1157 . a - each auxiliary field evaluated at the current point
1158 . a_t - the time derivative of each auxiliary field evaluated at the current point
1159 . a_x - the gradient of auxiliary each field evaluated at the current point
1160 . t - current time
1161 . u_tShift - the multiplier a for dF/dU_t
1162 . x - coordinates of the current point
1163 - g0 - output values at the current point
1164 
1165   Level: intermediate
1166 
1167 .seealso: PetscDSSetJacobian()
1168 @*/
1169 PetscErrorCode PetscDSGetJacobian(PetscDS prob, PetscInt f, PetscInt g,
1170                                   void (**g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1171                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1172                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1173                                               PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g0[]),
1174                                   void (**g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1175                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1176                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1177                                               PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g1[]),
1178                                   void (**g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1179                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1180                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1181                                               PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g2[]),
1182                                   void (**g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1183                                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1184                                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1185                                               PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g3[]))
1186 {
1187   PetscFunctionBegin;
1188   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1189   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);
1190   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);
1191   if (g0) {PetscValidPointer(g0, 4); *g0 = prob->g[(f*prob->Nf + g)*4+0];}
1192   if (g1) {PetscValidPointer(g1, 5); *g1 = prob->g[(f*prob->Nf + g)*4+1];}
1193   if (g2) {PetscValidPointer(g2, 6); *g2 = prob->g[(f*prob->Nf + g)*4+2];}
1194   if (g3) {PetscValidPointer(g3, 7); *g3 = prob->g[(f*prob->Nf + g)*4+3];}
1195   PetscFunctionReturn(0);
1196 }
1197 
1198 #undef __FUNCT__
1199 #define __FUNCT__ "PetscDSSetJacobian"
1200 /*@C
1201   PetscDSSetJacobian - Set the pointwise Jacobian function for given test and basis fields
1202 
1203   Not collective
1204 
1205   Input Parameters:
1206 + prob - The PetscDS
1207 . f    - The test field number
1208 . g    - The field number
1209 . g0 - integrand for the test and basis function term
1210 . g1 - integrand for the test function and basis function gradient term
1211 . g2 - integrand for the test function gradient and basis function term
1212 - g3 - integrand for the test function gradient and basis function gradient term
1213 
1214   Note: We are using a first order FEM model for the weak form:
1215 
1216   \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
1217 
1218 The calling sequence for the callbacks g0, g1, g2 and g3 is given by:
1219 
1220 $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1221 $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1222 $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1223 $    PetscReal t, const PetscReal x[], PetscScalar g0[])
1224 
1225 + dim - the spatial dimension
1226 . Nf - the number of fields
1227 . uOff - the offset into u[] and u_t[] for each field
1228 . uOff_x - the offset into u_x[] for each field
1229 . u - each field evaluated at the current point
1230 . u_t - the time derivative of each field evaluated at the current point
1231 . u_x - the gradient of each field evaluated at the current point
1232 . aOff - the offset into a[] and a_t[] for each auxiliary field
1233 . aOff_x - the offset into a_x[] for each auxiliary field
1234 . a - each auxiliary field evaluated at the current point
1235 . a_t - the time derivative of each auxiliary field evaluated at the current point
1236 . a_x - the gradient of auxiliary each field evaluated at the current point
1237 . t - current time
1238 . u_tShift - the multiplier a for dF/dU_t
1239 . x - coordinates of the current point
1240 - g0 - output values at the current point
1241 
1242   Level: intermediate
1243 
1244 .seealso: PetscDSGetJacobian()
1245 @*/
1246 PetscErrorCode PetscDSSetJacobian(PetscDS prob, PetscInt f, PetscInt g,
1247                                   void (*g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1248                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1249                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1250                                              PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g0[]),
1251                                   void (*g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1252                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1253                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1254                                              PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g1[]),
1255                                   void (*g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1256                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1257                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1258                                              PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g2[]),
1259                                   void (*g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1260                                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1261                                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1262                                              PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscScalar g3[]))
1263 {
1264   PetscErrorCode ierr;
1265 
1266   PetscFunctionBegin;
1267   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1268   if (g0) PetscValidFunction(g0, 4);
1269   if (g1) PetscValidFunction(g1, 5);
1270   if (g2) PetscValidFunction(g2, 6);
1271   if (g3) PetscValidFunction(g3, 7);
1272   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
1273   if (g < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", g);
1274   ierr = PetscDSEnlarge_Static(prob, PetscMax(f, g)+1);CHKERRQ(ierr);
1275   prob->g[(f*prob->Nf + g)*4+0] = g0;
1276   prob->g[(f*prob->Nf + g)*4+1] = g1;
1277   prob->g[(f*prob->Nf + g)*4+2] = g2;
1278   prob->g[(f*prob->Nf + g)*4+3] = g3;
1279   PetscFunctionReturn(0);
1280 }
1281 
1282 #undef __FUNCT__
1283 #define __FUNCT__ "PetscDSGetRiemannSolver"
1284 /*@C
1285   PetscDSGetRiemannSolver - Returns the Riemann solver for the given field
1286 
1287   Not collective
1288 
1289   Input Arguments:
1290 + prob - The PetscDS object
1291 - f    - The field number
1292 
1293   Output Argument:
1294 . r    - Riemann solver
1295 
1296   Calling sequence for r:
1297 
1298 $ r(PetscInt dim, PetscInt Nf, const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx)
1299 
1300 + dim  - The spatial dimension
1301 . Nf   - The number of fields
1302 . x    - The coordinates at a point on the interface
1303 . n    - The normal vector to the interface
1304 . uL   - The state vector to the left of the interface
1305 . uR   - The state vector to the right of the interface
1306 . flux - output array of flux through the interface
1307 - ctx  - optional user context
1308 
1309   Level: intermediate
1310 
1311 .seealso: PetscDSSetRiemannSolver()
1312 @*/
1313 PetscErrorCode PetscDSGetRiemannSolver(PetscDS prob, PetscInt f,
1314                                        void (**r)(PetscInt dim, PetscInt Nf, const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx))
1315 {
1316   PetscFunctionBegin;
1317   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1318   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);
1319   PetscValidPointer(r, 3);
1320   *r = prob->r[f];
1321   PetscFunctionReturn(0);
1322 }
1323 
1324 #undef __FUNCT__
1325 #define __FUNCT__ "PetscDSSetRiemannSolver"
1326 /*@C
1327   PetscDSSetRiemannSolver - Sets the Riemann solver for the given field
1328 
1329   Not collective
1330 
1331   Input Arguments:
1332 + prob - The PetscDS object
1333 . f    - The field number
1334 - r    - Riemann solver
1335 
1336   Calling sequence for r:
1337 
1338 $ r(PetscInt dim, PetscInt Nf, const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx)
1339 
1340 + dim  - The spatial dimension
1341 . Nf   - The number of fields
1342 . x    - The coordinates at a point on the interface
1343 . n    - The normal vector to the interface
1344 . uL   - The state vector to the left of the interface
1345 . uR   - The state vector to the right of the interface
1346 . flux - output array of flux through the interface
1347 - ctx  - optional user context
1348 
1349   Level: intermediate
1350 
1351 .seealso: PetscDSGetRiemannSolver()
1352 @*/
1353 PetscErrorCode PetscDSSetRiemannSolver(PetscDS prob, PetscInt f,
1354                                        void (*r)(PetscInt dim, PetscInt Nf, const PetscReal x[], const PetscReal n[], const PetscScalar uL[], const PetscScalar uR[], PetscScalar flux[], void *ctx))
1355 {
1356   PetscErrorCode ierr;
1357 
1358   PetscFunctionBegin;
1359   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1360   PetscValidFunction(r, 3);
1361   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
1362   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
1363   prob->r[f] = r;
1364   PetscFunctionReturn(0);
1365 }
1366 
1367 #undef __FUNCT__
1368 #define __FUNCT__ "PetscDSGetContext"
1369 PetscErrorCode PetscDSGetContext(PetscDS prob, PetscInt f, void **ctx)
1370 {
1371   PetscFunctionBegin;
1372   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1373   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);
1374   PetscValidPointer(ctx, 3);
1375   *ctx = prob->ctx[f];
1376   PetscFunctionReturn(0);
1377 }
1378 
1379 #undef __FUNCT__
1380 #define __FUNCT__ "PetscDSSetContext"
1381 PetscErrorCode PetscDSSetContext(PetscDS prob, PetscInt f, void *ctx)
1382 {
1383   PetscErrorCode ierr;
1384 
1385   PetscFunctionBegin;
1386   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1387   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
1388   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
1389   prob->ctx[f] = ctx;
1390   PetscFunctionReturn(0);
1391 }
1392 
1393 #undef __FUNCT__
1394 #define __FUNCT__ "PetscDSGetBdResidual"
1395 /*@C
1396   PetscDSGetBdResidual - Get the pointwise boundary residual function for a given test field
1397 
1398   Not collective
1399 
1400   Input Parameters:
1401 + prob - The PetscDS
1402 - f    - The test field number
1403 
1404   Output Parameters:
1405 + f0 - boundary integrand for the test function term
1406 - f1 - boundary integrand for the test function gradient term
1407 
1408   Note: We are using a first order FEM model for the weak form:
1409 
1410   \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
1411 
1412 The calling sequence for the callbacks f0 and f1 is given by:
1413 
1414 $ f0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1415 $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1416 $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1417 $    PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f0[])
1418 
1419 + dim - the spatial dimension
1420 . Nf - the number of fields
1421 . uOff - the offset into u[] and u_t[] for each field
1422 . uOff_x - the offset into u_x[] for each field
1423 . u - each field evaluated at the current point
1424 . u_t - the time derivative of each field evaluated at the current point
1425 . u_x - the gradient of each field evaluated at the current point
1426 . aOff - the offset into a[] and a_t[] for each auxiliary field
1427 . aOff_x - the offset into a_x[] for each auxiliary field
1428 . a - each auxiliary field evaluated at the current point
1429 . a_t - the time derivative of each auxiliary field evaluated at the current point
1430 . a_x - the gradient of auxiliary each field evaluated at the current point
1431 . t - current time
1432 . x - coordinates of the current point
1433 . n - unit normal at the current point
1434 - f0 - output values at the current point
1435 
1436   Level: intermediate
1437 
1438 .seealso: PetscDSSetBdResidual()
1439 @*/
1440 PetscErrorCode PetscDSGetBdResidual(PetscDS prob, PetscInt f,
1441                                     void (**f0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1442                                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1443                                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1444                                                 PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f0[]),
1445                                     void (**f1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1446                                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1447                                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1448                                                 PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f1[]))
1449 {
1450   PetscFunctionBegin;
1451   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1452   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);
1453   if (f0) {PetscValidPointer(f0, 3); *f0 = prob->fBd[f*2+0];}
1454   if (f1) {PetscValidPointer(f1, 4); *f1 = prob->fBd[f*2+1];}
1455   PetscFunctionReturn(0);
1456 }
1457 
1458 #undef __FUNCT__
1459 #define __FUNCT__ "PetscDSSetBdResidual"
1460 /*@C
1461   PetscDSSetBdResidual - Get the pointwise boundary residual function for a given test field
1462 
1463   Not collective
1464 
1465   Input Parameters:
1466 + prob - The PetscDS
1467 . f    - The test field number
1468 . f0 - boundary integrand for the test function term
1469 - f1 - boundary integrand for the test function gradient term
1470 
1471   Note: We are using a first order FEM model for the weak form:
1472 
1473   \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
1474 
1475 The calling sequence for the callbacks f0 and f1 is given by:
1476 
1477 $ f0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1478 $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1479 $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1480 $    PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f0[])
1481 
1482 + dim - the spatial dimension
1483 . Nf - the number of fields
1484 . uOff - the offset into u[] and u_t[] for each field
1485 . uOff_x - the offset into u_x[] for each field
1486 . u - each field evaluated at the current point
1487 . u_t - the time derivative of each field evaluated at the current point
1488 . u_x - the gradient of each field evaluated at the current point
1489 . aOff - the offset into a[] and a_t[] for each auxiliary field
1490 . aOff_x - the offset into a_x[] for each auxiliary field
1491 . a - each auxiliary field evaluated at the current point
1492 . a_t - the time derivative of each auxiliary field evaluated at the current point
1493 . a_x - the gradient of auxiliary each field evaluated at the current point
1494 . t - current time
1495 . x - coordinates of the current point
1496 . n - unit normal at the current point
1497 - f0 - output values at the current point
1498 
1499   Level: intermediate
1500 
1501 .seealso: PetscDSGetBdResidual()
1502 @*/
1503 PetscErrorCode PetscDSSetBdResidual(PetscDS prob, PetscInt f,
1504                                     void (*f0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1505                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1506                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1507                                                PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f0[]),
1508                                     void (*f1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1509                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1510                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1511                                                PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar f1[]))
1512 {
1513   PetscErrorCode ierr;
1514 
1515   PetscFunctionBegin;
1516   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1517   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
1518   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
1519   if (f0) {PetscValidFunction(f0, 3); prob->fBd[f*2+0] = f0;}
1520   if (f1) {PetscValidFunction(f1, 4); prob->fBd[f*2+1] = f1;}
1521   PetscFunctionReturn(0);
1522 }
1523 
1524 #undef __FUNCT__
1525 #define __FUNCT__ "PetscDSGetBdJacobian"
1526 /*@C
1527   PetscDSGetBdJacobian - Get the pointwise boundary Jacobian function for given test and basis field
1528 
1529   Not collective
1530 
1531   Input Parameters:
1532 + prob - The PetscDS
1533 . f    - The test field number
1534 - g    - The field number
1535 
1536   Output Parameters:
1537 + g0 - integrand for the test and basis function term
1538 . g1 - integrand for the test function and basis function gradient term
1539 . g2 - integrand for the test function gradient and basis function term
1540 - g3 - integrand for the test function gradient and basis function gradient term
1541 
1542   Note: We are using a first order FEM model for the weak form:
1543 
1544   \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
1545 
1546 The calling sequence for the callbacks g0, g1, g2 and g3 is given by:
1547 
1548 $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1549 $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1550 $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1551 $    PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar g0[])
1552 
1553 + dim - the spatial dimension
1554 . Nf - the number of fields
1555 . uOff - the offset into u[] and u_t[] for each field
1556 . uOff_x - the offset into u_x[] for each field
1557 . u - each field evaluated at the current point
1558 . u_t - the time derivative of each field evaluated at the current point
1559 . u_x - the gradient of each field evaluated at the current point
1560 . aOff - the offset into a[] and a_t[] for each auxiliary field
1561 . aOff_x - the offset into a_x[] for each auxiliary field
1562 . a - each auxiliary field evaluated at the current point
1563 . a_t - the time derivative of each auxiliary field evaluated at the current point
1564 . a_x - the gradient of auxiliary each field evaluated at the current point
1565 . t - current time
1566 . u_tShift - the multiplier a for dF/dU_t
1567 . x - coordinates of the current point
1568 . n - normal at the current point
1569 - g0 - output values at the current point
1570 
1571   Level: intermediate
1572 
1573 .seealso: PetscDSSetBdJacobian()
1574 @*/
1575 PetscErrorCode PetscDSGetBdJacobian(PetscDS prob, PetscInt f, PetscInt g,
1576                                     void (**g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
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, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscScalar g0[]),
1580                                     void (**g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
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, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscScalar g1[]),
1584                                     void (**g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1585                                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1586                                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1587                                                 PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscScalar g2[]),
1588                                     void (**g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1589                                                 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1590                                                 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1591                                                 PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscScalar g3[]))
1592 {
1593   PetscFunctionBegin;
1594   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1595   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);
1596   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);
1597   if (g0) {PetscValidPointer(g0, 4); *g0 = prob->gBd[(f*prob->Nf + g)*4+0];}
1598   if (g1) {PetscValidPointer(g1, 5); *g1 = prob->gBd[(f*prob->Nf + g)*4+1];}
1599   if (g2) {PetscValidPointer(g2, 6); *g2 = prob->gBd[(f*prob->Nf + g)*4+2];}
1600   if (g3) {PetscValidPointer(g3, 7); *g3 = prob->gBd[(f*prob->Nf + g)*4+3];}
1601   PetscFunctionReturn(0);
1602 }
1603 
1604 #undef __FUNCT__
1605 #define __FUNCT__ "PetscDSSetBdJacobian"
1606 /*@C
1607   PetscDSSetBdJacobian - Set the pointwise boundary Jacobian function for given test and basis field
1608 
1609   Not collective
1610 
1611   Input Parameters:
1612 + prob - The PetscDS
1613 . f    - The test field number
1614 . g    - The field number
1615 . g0 - integrand for the test and basis function term
1616 . g1 - integrand for the test function and basis function gradient term
1617 . g2 - integrand for the test function gradient and basis function term
1618 - g3 - integrand for the test function gradient and basis function gradient term
1619 
1620   Note: We are using a first order FEM model for the weak form:
1621 
1622   \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
1623 
1624 The calling sequence for the callbacks g0, g1, g2 and g3 is given by:
1625 
1626 $ g0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1627 $    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1628 $    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1629 $    PetscReal t, const PetscReal x[], const PetscReal n[], PetscScalar g0[])
1630 
1631 + dim - the spatial dimension
1632 . Nf - the number of fields
1633 . uOff - the offset into u[] and u_t[] for each field
1634 . uOff_x - the offset into u_x[] for each field
1635 . u - each field evaluated at the current point
1636 . u_t - the time derivative of each field evaluated at the current point
1637 . u_x - the gradient of each field evaluated at the current point
1638 . aOff - the offset into a[] and a_t[] for each auxiliary field
1639 . aOff_x - the offset into a_x[] for each auxiliary field
1640 . a - each auxiliary field evaluated at the current point
1641 . a_t - the time derivative of each auxiliary field evaluated at the current point
1642 . a_x - the gradient of auxiliary each field evaluated at the current point
1643 . t - current time
1644 . u_tShift - the multiplier a for dF/dU_t
1645 . x - coordinates of the current point
1646 . n - normal at the current point
1647 - g0 - output values at the current point
1648 
1649   Level: intermediate
1650 
1651 .seealso: PetscDSGetBdJacobian()
1652 @*/
1653 PetscErrorCode PetscDSSetBdJacobian(PetscDS prob, PetscInt f, PetscInt g,
1654                                     void (*g0)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1655                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1656                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1657                                                PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscScalar g0[]),
1658                                     void (*g1)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1659                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1660                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1661                                                PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscScalar g1[]),
1662                                     void (*g2)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1663                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1664                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1665                                                PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscScalar g2[]),
1666                                     void (*g3)(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1667                                                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1668                                                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1669                                                PetscReal t, PetscReal u_tShift, const PetscReal x[], const PetscReal n[], PetscScalar g3[]))
1670 {
1671   PetscErrorCode ierr;
1672 
1673   PetscFunctionBegin;
1674   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1675   if (g0) PetscValidFunction(g0, 4);
1676   if (g1) PetscValidFunction(g1, 5);
1677   if (g2) PetscValidFunction(g2, 6);
1678   if (g3) PetscValidFunction(g3, 7);
1679   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
1680   if (g < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", g);
1681   ierr = PetscDSEnlarge_Static(prob, PetscMax(f, g)+1);CHKERRQ(ierr);
1682   prob->gBd[(f*prob->Nf + g)*4+0] = g0;
1683   prob->gBd[(f*prob->Nf + g)*4+1] = g1;
1684   prob->gBd[(f*prob->Nf + g)*4+2] = g2;
1685   prob->gBd[(f*prob->Nf + g)*4+3] = g3;
1686   PetscFunctionReturn(0);
1687 }
1688 
1689 #undef __FUNCT__
1690 #define __FUNCT__ "PetscDSGetFieldOffset"
1691 /*@
1692   PetscDSGetFieldOffset - Returns the offset of the given field in the full space basis
1693 
1694   Not collective
1695 
1696   Input Parameters:
1697 + prob - The PetscDS object
1698 - f - The field number
1699 
1700   Output Parameter:
1701 . off - The offset
1702 
1703   Level: beginner
1704 
1705 .seealso: PetscDSGetBdFieldOffset(), PetscDSGetNumFields(), PetscDSCreate()
1706 @*/
1707 PetscErrorCode PetscDSGetFieldOffset(PetscDS prob, PetscInt f, PetscInt *off)
1708 {
1709   PetscInt       g;
1710   PetscErrorCode ierr;
1711 
1712   PetscFunctionBegin;
1713   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1714   PetscValidPointer(off, 3);
1715   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);
1716   *off = 0;
1717   for (g = 0; g < f; ++g) {
1718     PetscFE  fe = (PetscFE) prob->disc[g];
1719     PetscInt Nb, Nc;
1720 
1721     ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
1722     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
1723     *off += Nb*Nc;
1724   }
1725   PetscFunctionReturn(0);
1726 }
1727 
1728 #undef __FUNCT__
1729 #define __FUNCT__ "PetscDSGetBdFieldOffset"
1730 /*@
1731   PetscDSGetBdFieldOffset - Returns the offset of the given field in the full space boundary basis
1732 
1733   Not collective
1734 
1735   Input Parameters:
1736 + prob - The PetscDS object
1737 - f - The field number
1738 
1739   Output Parameter:
1740 . off - The boundary offset
1741 
1742   Level: beginner
1743 
1744 .seealso: PetscDSGetFieldOffset(), PetscDSGetNumFields(), PetscDSCreate()
1745 @*/
1746 PetscErrorCode PetscDSGetBdFieldOffset(PetscDS prob, PetscInt f, PetscInt *off)
1747 {
1748   PetscInt       g;
1749   PetscErrorCode ierr;
1750 
1751   PetscFunctionBegin;
1752   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1753   PetscValidPointer(off, 3);
1754   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);
1755   *off = 0;
1756   for (g = 0; g < f; ++g) {
1757     PetscFE  fe = (PetscFE) prob->discBd[g];
1758     PetscInt Nb, Nc;
1759 
1760     ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
1761     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
1762     *off += Nb*Nc;
1763   }
1764   PetscFunctionReturn(0);
1765 }
1766 
1767 #undef __FUNCT__
1768 #define __FUNCT__ "PetscDSGetComponentOffset"
1769 /*@
1770   PetscDSGetComponentOffset - Returns the offset of the given field on an evaluation point
1771 
1772   Not collective
1773 
1774   Input Parameters:
1775 + prob - The PetscDS object
1776 - f - The field number
1777 
1778   Output Parameter:
1779 . off - The offset
1780 
1781   Level: beginner
1782 
1783 .seealso: PetscDSGetBdFieldOffset(), PetscDSGetNumFields(), PetscDSCreate()
1784 @*/
1785 PetscErrorCode PetscDSGetComponentOffset(PetscDS prob, PetscInt f, PetscInt *off)
1786 {
1787   PetscInt       g;
1788   PetscErrorCode ierr;
1789 
1790   PetscFunctionBegin;
1791   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1792   PetscValidPointer(off, 3);
1793   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);
1794   *off = 0;
1795   for (g = 0; g < f; ++g) {
1796     PetscFE  fe = (PetscFE) prob->disc[g];
1797     PetscInt Nc;
1798 
1799     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
1800     *off += Nc;
1801   }
1802   PetscFunctionReturn(0);
1803 }
1804 
1805 #undef __FUNCT__
1806 #define __FUNCT__ "PetscDSGetComponentOffsets"
1807 /*@
1808   PetscDSGetComponentOffsets - Returns the offset of each field on an evaluation point
1809 
1810   Not collective
1811 
1812   Input Parameter:
1813 . prob - The PetscDS object
1814 
1815   Output Parameter:
1816 . offsets - The offsets
1817 
1818   Level: beginner
1819 
1820 .seealso: PetscDSGetBdFieldOffset(), PetscDSGetNumFields(), PetscDSCreate()
1821 @*/
1822 PetscErrorCode PetscDSGetComponentOffsets(PetscDS prob, PetscInt *offsets[])
1823 {
1824   PetscFunctionBegin;
1825   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1826   PetscValidPointer(offsets, 2);
1827   *offsets = prob->off;
1828   PetscFunctionReturn(0);
1829 }
1830 
1831 #undef __FUNCT__
1832 #define __FUNCT__ "PetscDSGetComponentDerivativeOffsets"
1833 /*@
1834   PetscDSGetComponentDerivativeOffsets - Returns the offset of each field derivative on an evaluation point
1835 
1836   Not collective
1837 
1838   Input Parameter:
1839 . prob - The PetscDS object
1840 
1841   Output Parameter:
1842 . offsets - The offsets
1843 
1844   Level: beginner
1845 
1846 .seealso: PetscDSGetBdFieldOffset(), PetscDSGetNumFields(), PetscDSCreate()
1847 @*/
1848 PetscErrorCode PetscDSGetComponentDerivativeOffsets(PetscDS prob, PetscInt *offsets[])
1849 {
1850   PetscFunctionBegin;
1851   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1852   PetscValidPointer(offsets, 2);
1853   *offsets = prob->offDer;
1854   PetscFunctionReturn(0);
1855 }
1856 
1857 #undef __FUNCT__
1858 #define __FUNCT__ "PetscDSGetComponentBdOffsets"
1859 /*@
1860   PetscDSGetComponentBdOffsets - Returns the offset of each field on a boundary evaluation point
1861 
1862   Not collective
1863 
1864   Input Parameter:
1865 . prob - The PetscDS object
1866 
1867   Output Parameter:
1868 . offsets - The offsets
1869 
1870   Level: beginner
1871 
1872 .seealso: PetscDSGetBdFieldOffset(), PetscDSGetNumFields(), PetscDSCreate()
1873 @*/
1874 PetscErrorCode PetscDSGetComponentBdOffsets(PetscDS prob, PetscInt *offsets[])
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   PetscFunctionBegin;
1903   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1904   PetscValidPointer(offsets, 2);
1905   *offsets = prob->offDerBd;
1906   PetscFunctionReturn(0);
1907 }
1908 
1909 #undef __FUNCT__
1910 #define __FUNCT__ "PetscDSGetTabulation"
1911 /*@C
1912   PetscDSGetTabulation - Return the basis tabulation at quadrature points for the volume discretization
1913 
1914   Not collective
1915 
1916   Input Parameter:
1917 . prob - The PetscDS object
1918 
1919   Output Parameters:
1920 + basis - The basis function tabulation at quadrature points
1921 - basisDer - The basis function derivative tabulation at quadrature points
1922 
1923   Level: intermediate
1924 
1925 .seealso: PetscDSGetBdTabulation(), PetscDSCreate()
1926 @*/
1927 PetscErrorCode PetscDSGetTabulation(PetscDS prob, PetscReal ***basis, PetscReal ***basisDer)
1928 {
1929   PetscErrorCode ierr;
1930 
1931   PetscFunctionBegin;
1932   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1933   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
1934   if (basis)    {PetscValidPointer(basis, 2);    *basis    = prob->basis;}
1935   if (basisDer) {PetscValidPointer(basisDer, 3); *basisDer = prob->basisDer;}
1936   PetscFunctionReturn(0);
1937 }
1938 
1939 #undef __FUNCT__
1940 #define __FUNCT__ "PetscDSGetBdTabulation"
1941 /*@C
1942   PetscDSGetBdTabulation - Return the basis tabulation at quadrature points for the boundary discretization
1943 
1944   Not collective
1945 
1946   Input Parameter:
1947 . prob - The PetscDS object
1948 
1949   Output Parameters:
1950 + basis - The basis function tabulation at quadrature points
1951 - basisDer - The basis function derivative tabulation at quadrature points
1952 
1953   Level: intermediate
1954 
1955 .seealso: PetscDSGetTabulation(), PetscDSCreate()
1956 @*/
1957 PetscErrorCode PetscDSGetBdTabulation(PetscDS prob, PetscReal ***basis, PetscReal ***basisDer)
1958 {
1959   PetscErrorCode ierr;
1960 
1961   PetscFunctionBegin;
1962   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1963   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
1964   if (basis)    {PetscValidPointer(basis, 2);    *basis    = prob->basisBd;}
1965   if (basisDer) {PetscValidPointer(basisDer, 3); *basisDer = prob->basisDerBd;}
1966   PetscFunctionReturn(0);
1967 }
1968 
1969 #undef __FUNCT__
1970 #define __FUNCT__ "PetscDSGetEvaluationArrays"
1971 PetscErrorCode PetscDSGetEvaluationArrays(PetscDS prob, PetscScalar **u, PetscScalar **u_t, PetscScalar **u_x)
1972 {
1973   PetscErrorCode ierr;
1974 
1975   PetscFunctionBegin;
1976   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1977   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
1978   if (u)   {PetscValidPointer(u, 2);   *u   = prob->u;}
1979   if (u_t) {PetscValidPointer(u_t, 3); *u_t = prob->u_t;}
1980   if (u_x) {PetscValidPointer(u_x, 4); *u_x = prob->u_x;}
1981   PetscFunctionReturn(0);
1982 }
1983 
1984 #undef __FUNCT__
1985 #define __FUNCT__ "PetscDSGetWeakFormArrays"
1986 PetscErrorCode PetscDSGetWeakFormArrays(PetscDS prob, PetscScalar **f0, PetscScalar **f1, PetscScalar **g0, PetscScalar **g1, PetscScalar **g2, PetscScalar **g3)
1987 {
1988   PetscErrorCode ierr;
1989 
1990   PetscFunctionBegin;
1991   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
1992   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
1993   if (f0) {PetscValidPointer(f0, 2); *f0 = prob->f0;}
1994   if (f1) {PetscValidPointer(f1, 3); *f1 = prob->f1;}
1995   if (g0) {PetscValidPointer(g0, 4); *g0 = prob->g0;}
1996   if (g1) {PetscValidPointer(g1, 5); *g1 = prob->g1;}
1997   if (g2) {PetscValidPointer(g2, 6); *g2 = prob->g2;}
1998   if (g3) {PetscValidPointer(g3, 7); *g3 = prob->g3;}
1999   PetscFunctionReturn(0);
2000 }
2001 
2002 #undef __FUNCT__
2003 #define __FUNCT__ "PetscDSGetRefCoordArrays"
2004 PetscErrorCode PetscDSGetRefCoordArrays(PetscDS prob, PetscReal **x, PetscScalar **refSpaceDer)
2005 {
2006   PetscErrorCode ierr;
2007 
2008   PetscFunctionBegin;
2009   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
2010   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
2011   if (x)           {PetscValidPointer(x, 2);           *x           = prob->x;}
2012   if (refSpaceDer) {PetscValidPointer(refSpaceDer, 3); *refSpaceDer = prob->refSpaceDer;}
2013   PetscFunctionReturn(0);
2014 }
2015 
2016 #undef __FUNCT__
2017 #define __FUNCT__ "PetscDSDestroy_Basic"
2018 static PetscErrorCode PetscDSDestroy_Basic(PetscDS prob)
2019 {
2020   PetscFunctionBegin;
2021   PetscFunctionReturn(0);
2022 }
2023 
2024 #undef __FUNCT__
2025 #define __FUNCT__ "PetscDSInitialize_Basic"
2026 static PetscErrorCode PetscDSInitialize_Basic(PetscDS prob)
2027 {
2028   PetscFunctionBegin;
2029   prob->ops->setfromoptions = NULL;
2030   prob->ops->setup          = NULL;
2031   prob->ops->view           = NULL;
2032   prob->ops->destroy        = PetscDSDestroy_Basic;
2033   PetscFunctionReturn(0);
2034 }
2035 
2036 /*MC
2037   PETSCDSBASIC = "basic" - A discrete system with pointwise residual and boundary residual functions
2038 
2039   Level: intermediate
2040 
2041 .seealso: PetscDSType, PetscDSCreate(), PetscDSSetType()
2042 M*/
2043 
2044 #undef __FUNCT__
2045 #define __FUNCT__ "PetscDSCreate_Basic"
2046 PETSC_EXTERN PetscErrorCode PetscDSCreate_Basic(PetscDS prob)
2047 {
2048   PetscDS_Basic *b;
2049   PetscErrorCode      ierr;
2050 
2051   PetscFunctionBegin;
2052   PetscValidHeaderSpecific(prob, PETSCSPACE_CLASSID, 1);
2053   ierr       = PetscNewLog(prob, &b);CHKERRQ(ierr);
2054   prob->data = b;
2055 
2056   ierr = PetscDSInitialize_Basic(prob);CHKERRQ(ierr);
2057   PetscFunctionReturn(0);
2058 }
2059