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