xref: /petsc/src/dm/dt/interface/dtds.c (revision bbd3f3368baff017a6ea7962313571cf23b2b4ec)
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   PetscValidCharPointer(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         fe   = (PetscFE) prob->disc[f];
263     PetscFE         feBd = (PetscFE) prob->discBd[f];
264     PetscQuadrature q;
265     PetscInt        Nq, Nb, Nc;
266 
267     /* TODO Dispatch on discretization type*/
268     ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
269     ierr = PetscQuadratureGetData(q, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);
270     ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
271     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
272     ierr = PetscFEGetDefaultTabulation(fe, &prob->basis[f], &prob->basisDer[f], NULL);CHKERRQ(ierr);
273     NqMax          = PetscMax(NqMax, Nq);
274     NcMax          = PetscMax(NcMax, Nc);
275     prob->totDim  += Nb*Nc;
276     prob->totComp += Nc;
277     if (feBd) {
278       ierr = PetscFEGetDimension(feBd, &Nb);CHKERRQ(ierr);
279       ierr = PetscFEGetNumComponents(feBd, &Nc);CHKERRQ(ierr);
280       ierr = PetscFEGetDefaultTabulation(feBd, &prob->basisBd[f], &prob->basisDerBd[f], NULL);CHKERRQ(ierr);
281       prob->totDimBd += Nb*Nc;
282     }
283   }
284   work = PetscMax(prob->totComp*dim, PetscSqr(NcMax*dim));
285   /* Allocate works space */
286   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);
287   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);
288   if (prob->ops->setup) {ierr = (*prob->ops->setup)(prob);CHKERRQ(ierr);}
289   prob->setup = PETSC_TRUE;
290   PetscFunctionReturn(0);
291 }
292 
293 #undef __FUNCT__
294 #define __FUNCT__ "PetscDSDestroyStructs_Static"
295 static PetscErrorCode PetscDSDestroyStructs_Static(PetscDS prob)
296 {
297   PetscErrorCode ierr;
298 
299   PetscFunctionBegin;
300   ierr = PetscFree4(prob->basis,prob->basisDer,prob->basisBd,prob->basisDerBd);CHKERRQ(ierr);
301   ierr = PetscFree5(prob->u,prob->u_t,prob->u_x,prob->x,prob->refSpaceDer);CHKERRQ(ierr);
302   ierr = PetscFree6(prob->f0,prob->f1,prob->g0,prob->g1,prob->g2,prob->g3);CHKERRQ(ierr);
303   PetscFunctionReturn(0);
304 }
305 
306 #undef __FUNCT__
307 #define __FUNCT__ "PetscDSEnlarge_Static"
308 static PetscErrorCode PetscDSEnlarge_Static(PetscDS prob, PetscInt NfNew)
309 {
310   PetscObject   *tmpd, *tmpdbd;
311   PointFunc     *tmpobj, *tmpf, *tmpg;
312   BdPointFunc   *tmpfbd, *tmpgbd;
313   PetscInt       Nf = prob->Nf, f;
314   PetscErrorCode ierr;
315 
316   PetscFunctionBegin;
317   if (Nf >= NfNew) PetscFunctionReturn(0);
318   prob->setup = PETSC_FALSE;
319   ierr = PetscDSDestroyStructs_Static(prob);CHKERRQ(ierr);
320   ierr = PetscMalloc1(NfNew, &tmpd);CHKERRQ(ierr);
321   for (f = 0; f < Nf; ++f) tmpd[f] = prob->disc[f];
322   for (f = Nf; f < NfNew; ++f) {ierr = PetscContainerCreate(PetscObjectComm((PetscObject) prob), (PetscContainer *) &tmpd[f]);CHKERRQ(ierr);}
323   ierr = PetscFree(prob->disc);CHKERRQ(ierr);
324   prob->Nf   = NfNew;
325   prob->disc = tmpd;
326   ierr = PetscCalloc3(NfNew, &tmpobj, NfNew*2, &tmpf, NfNew*NfNew*4, &tmpg);CHKERRQ(ierr);
327   for (f = 0; f < Nf; ++f) tmpobj[f] = prob->obj[f];
328   for (f = 0; f < Nf*2; ++f) tmpf[f] = prob->f[f];
329   for (f = 0; f < Nf*Nf*4; ++f) tmpg[f] = prob->g[f];
330   for (f = Nf; f < NfNew; ++f) tmpobj[f] = NULL;
331   for (f = Nf*2; f < NfNew*2; ++f) tmpf[f] = NULL;
332   for (f = Nf*Nf*4; f < NfNew*NfNew*4; ++f) tmpg[f] = NULL;
333   ierr = PetscFree3(prob->obj, prob->f, prob->g);CHKERRQ(ierr);
334   prob->obj = tmpobj;
335   prob->f   = tmpf;
336   prob->g   = tmpg;
337   ierr = PetscMalloc1(NfNew, &tmpdbd);CHKERRQ(ierr);
338   for (f = 0; f < Nf; ++f) tmpdbd[f] = prob->discBd[f];
339   for (f = Nf; f < NfNew; ++f) tmpdbd[f] = NULL;
340   ierr = PetscFree(prob->discBd);CHKERRQ(ierr);
341   prob->discBd = tmpdbd;
342   ierr = PetscCalloc2(NfNew*2, &tmpfbd, NfNew*NfNew*4, &tmpgbd);CHKERRQ(ierr);
343   for (f = 0; f < Nf*2; ++f) tmpfbd[f] = prob->fBd[f];
344   for (f = 0; f < Nf*Nf*4; ++f) tmpgbd[f] = prob->gBd[f];
345   for (f = Nf*2; f < NfNew*2; ++f) tmpfbd[f] = NULL;
346   for (f = Nf*Nf*4; f < NfNew*NfNew*4; ++f) tmpgbd[f] = NULL;
347   ierr = PetscFree2(prob->fBd, prob->gBd);CHKERRQ(ierr);
348   prob->fBd = tmpfbd;
349   prob->gBd = tmpgbd;
350   PetscFunctionReturn(0);
351 }
352 
353 #undef __FUNCT__
354 #define __FUNCT__ "PetscDSDestroy"
355 /*@
356   PetscDSDestroy - Destroys a PetscDS object
357 
358   Collective on PetscDS
359 
360   Input Parameter:
361 . prob - the PetscDS object to destroy
362 
363   Level: developer
364 
365 .seealso PetscDSView()
366 @*/
367 PetscErrorCode PetscDSDestroy(PetscDS *prob)
368 {
369   PetscInt       f;
370   PetscErrorCode ierr;
371 
372   PetscFunctionBegin;
373   if (!*prob) PetscFunctionReturn(0);
374   PetscValidHeaderSpecific((*prob), PETSCDS_CLASSID, 1);
375 
376   if (--((PetscObject)(*prob))->refct > 0) {*prob = 0; PetscFunctionReturn(0);}
377   ((PetscObject) (*prob))->refct = 0;
378   ierr = PetscDSDestroyStructs_Static(*prob);CHKERRQ(ierr);
379   for (f = 0; f < (*prob)->Nf; ++f) {
380     ierr = PetscObjectDereference((*prob)->disc[f]);CHKERRQ(ierr);
381     ierr = PetscObjectDereference((*prob)->discBd[f]);CHKERRQ(ierr);
382   }
383   ierr = PetscFree((*prob)->disc);CHKERRQ(ierr);
384   ierr = PetscFree((*prob)->discBd);CHKERRQ(ierr);
385   ierr = PetscFree3((*prob)->obj,(*prob)->f,(*prob)->g);CHKERRQ(ierr);
386   ierr = PetscFree2((*prob)->fBd,(*prob)->gBd);CHKERRQ(ierr);
387   if ((*prob)->ops->destroy) {ierr = (*(*prob)->ops->destroy)(*prob);CHKERRQ(ierr);}
388   ierr = PetscHeaderDestroy(prob);CHKERRQ(ierr);
389   PetscFunctionReturn(0);
390 }
391 
392 #undef __FUNCT__
393 #define __FUNCT__ "PetscDSCreate"
394 /*@
395   PetscDSCreate - Creates an empty PetscDS object. The type can then be set with PetscDSSetType().
396 
397   Collective on MPI_Comm
398 
399   Input Parameter:
400 . comm - The communicator for the PetscDS object
401 
402   Output Parameter:
403 . prob - The PetscDS object
404 
405   Level: beginner
406 
407 .seealso: PetscDSSetType(), PETSCDSBASIC
408 @*/
409 PetscErrorCode PetscDSCreate(MPI_Comm comm, PetscDS *prob)
410 {
411   PetscDS   p;
412   PetscErrorCode ierr;
413 
414   PetscFunctionBegin;
415   PetscValidPointer(prob, 2);
416   *prob  = NULL;
417   ierr = PetscDSInitializePackage();CHKERRQ(ierr);
418 
419   ierr = PetscHeaderCreate(p, _p_PetscDS, struct _PetscDSOps, PETSCDS_CLASSID, "PetscDS", "Discrete System", "PetscDS", comm, PetscDSDestroy, PetscDSView);CHKERRQ(ierr);
420   ierr = PetscMemzero(p->ops, sizeof(struct _PetscDSOps));CHKERRQ(ierr);
421 
422   p->Nf    = 0;
423   p->setup = PETSC_FALSE;
424 
425   *prob = p;
426   PetscFunctionReturn(0);
427 }
428 
429 #undef __FUNCT__
430 #define __FUNCT__ "PetscDSGetNumFields"
431 PetscErrorCode PetscDSGetNumFields(PetscDS prob, PetscInt *Nf)
432 {
433   PetscFunctionBegin;
434   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
435   PetscValidPointer(Nf, 2);
436   *Nf = prob->Nf;
437   PetscFunctionReturn(0);
438 }
439 
440 #undef __FUNCT__
441 #define __FUNCT__ "PetscDSGetSpatialDimension"
442 PetscErrorCode PetscDSGetSpatialDimension(PetscDS prob, PetscInt *dim)
443 {
444   PetscErrorCode ierr;
445 
446   PetscFunctionBegin;
447   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
448   PetscValidPointer(dim, 2);
449   *dim = 0;
450   if (prob->Nf) {ierr = PetscFEGetSpatialDimension((PetscFE) prob->disc[0], dim);CHKERRQ(ierr);}
451   PetscFunctionReturn(0);
452 }
453 
454 #undef __FUNCT__
455 #define __FUNCT__ "PetscDSGetTotalDimension"
456 PetscErrorCode PetscDSGetTotalDimension(PetscDS prob, PetscInt *dim)
457 {
458   PetscErrorCode ierr;
459 
460   PetscFunctionBegin;
461   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
462   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
463   PetscValidPointer(dim, 2);
464   *dim = prob->totDim;
465   PetscFunctionReturn(0);
466 }
467 
468 #undef __FUNCT__
469 #define __FUNCT__ "PetscDSGetTotalBdDimension"
470 PetscErrorCode PetscDSGetTotalBdDimension(PetscDS prob, PetscInt *dim)
471 {
472   PetscErrorCode ierr;
473 
474   PetscFunctionBegin;
475   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
476   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
477   PetscValidPointer(dim, 2);
478   *dim = prob->totDimBd;
479   PetscFunctionReturn(0);
480 }
481 
482 #undef __FUNCT__
483 #define __FUNCT__ "PetscDSGetTotalComponents"
484 PetscErrorCode PetscDSGetTotalComponents(PetscDS prob, PetscInt *Nc)
485 {
486   PetscErrorCode ierr;
487 
488   PetscFunctionBegin;
489   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
490   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
491   PetscValidPointer(Nc, 2);
492   *Nc = prob->totComp;
493   PetscFunctionReturn(0);
494 }
495 
496 #undef __FUNCT__
497 #define __FUNCT__ "PetscDSGetDiscretization"
498 PetscErrorCode PetscDSGetDiscretization(PetscDS prob, PetscInt f, PetscObject *disc)
499 {
500   PetscFunctionBegin;
501   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
502   PetscValidPointer(disc, 3);
503   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);
504   *disc = prob->disc[f];
505   PetscFunctionReturn(0);
506 }
507 
508 #undef __FUNCT__
509 #define __FUNCT__ "PetscDSGetBdDiscretization"
510 PetscErrorCode PetscDSGetBdDiscretization(PetscDS prob, PetscInt f, PetscObject *disc)
511 {
512   PetscFunctionBegin;
513   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
514   PetscValidPointer(disc, 3);
515   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);
516   *disc = prob->discBd[f];
517   PetscFunctionReturn(0);
518 }
519 
520 #undef __FUNCT__
521 #define __FUNCT__ "PetscDSSetDiscretization"
522 PetscErrorCode PetscDSSetDiscretization(PetscDS prob, PetscInt f, PetscObject disc)
523 {
524   PetscErrorCode ierr;
525 
526   PetscFunctionBegin;
527   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
528   PetscValidPointer(disc, 3);
529   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
530   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
531   if (prob->disc[f]) {ierr = PetscObjectDereference(prob->disc[f]);CHKERRQ(ierr);}
532   prob->disc[f] = disc;
533   ierr = PetscObjectReference(disc);CHKERRQ(ierr);
534   PetscFunctionReturn(0);
535 }
536 
537 #undef __FUNCT__
538 #define __FUNCT__ "PetscDSSetBdDiscretization"
539 PetscErrorCode PetscDSSetBdDiscretization(PetscDS prob, PetscInt f, PetscObject disc)
540 {
541   PetscErrorCode ierr;
542 
543   PetscFunctionBegin;
544   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
545   if (disc) PetscValidPointer(disc, 3);
546   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
547   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
548   if (prob->discBd[f]) {ierr = PetscObjectDereference(prob->discBd[f]);CHKERRQ(ierr);}
549   prob->discBd[f] = disc;
550   ierr = PetscObjectReference(disc);CHKERRQ(ierr);
551   PetscFunctionReturn(0);
552 }
553 
554 #undef __FUNCT__
555 #define __FUNCT__ "PetscDSAddDiscretization"
556 PetscErrorCode PetscDSAddDiscretization(PetscDS prob, PetscObject disc)
557 {
558   PetscErrorCode ierr;
559 
560   PetscFunctionBegin;
561   ierr = PetscDSSetDiscretization(prob, prob->Nf, disc);CHKERRQ(ierr);
562   PetscFunctionReturn(0);
563 }
564 
565 #undef __FUNCT__
566 #define __FUNCT__ "PetscDSAddBdDiscretization"
567 PetscErrorCode PetscDSAddBdDiscretization(PetscDS prob, PetscObject disc)
568 {
569   PetscErrorCode ierr;
570 
571   PetscFunctionBegin;
572   ierr = PetscDSSetBdDiscretization(prob, prob->Nf, disc);CHKERRQ(ierr);
573   PetscFunctionReturn(0);
574 }
575 
576 #undef __FUNCT__
577 #define __FUNCT__ "PetscDSGetObjective"
578 PetscErrorCode PetscDSGetObjective(PetscDS prob, PetscInt f,
579                                         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[]))
580 {
581   PetscFunctionBegin;
582   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
583   PetscValidPointer(obj, 2);
584   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);
585   *obj = prob->obj[f];
586   PetscFunctionReturn(0);
587 }
588 
589 #undef __FUNCT__
590 #define __FUNCT__ "PetscDSSetObjective"
591 PetscErrorCode PetscDSSetObjective(PetscDS prob, PetscInt f,
592                                         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[]))
593 {
594   PetscErrorCode ierr;
595 
596   PetscFunctionBegin;
597   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
598   PetscValidFunction(obj, 2);
599   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
600   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
601   prob->obj[f] = obj;
602   PetscFunctionReturn(0);
603 }
604 
605 #undef __FUNCT__
606 #define __FUNCT__ "PetscDSGetResidual"
607 PetscErrorCode PetscDSGetResidual(PetscDS prob, PetscInt f,
608                                        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[]),
609                                        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[]))
610 {
611   PetscFunctionBegin;
612   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
613   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);
614   if (f0) {PetscValidPointer(f0, 3); *f0 = prob->f[f*2+0];}
615   if (f1) {PetscValidPointer(f1, 4); *f1 = prob->f[f*2+1];}
616   PetscFunctionReturn(0);
617 }
618 
619 #undef __FUNCT__
620 #define __FUNCT__ "PetscDSSetResidual"
621 PetscErrorCode PetscDSSetResidual(PetscDS prob, PetscInt f,
622                                        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[]),
623                                        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[]))
624 {
625   PetscErrorCode ierr;
626 
627   PetscFunctionBegin;
628   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
629   PetscValidFunction(f0, 3);
630   PetscValidFunction(f1, 4);
631   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
632   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
633   prob->f[f*2+0] = f0;
634   prob->f[f*2+1] = f1;
635   PetscFunctionReturn(0);
636 }
637 
638 #undef __FUNCT__
639 #define __FUNCT__ "PetscDSGetJacobian"
640 PetscErrorCode PetscDSGetJacobian(PetscDS prob, PetscInt f, PetscInt g,
641                                        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[]),
642                                        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[]),
643                                        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[]),
644                                        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[]))
645 {
646   PetscFunctionBegin;
647   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
648   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);
649   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);
650   if (g0) {PetscValidPointer(g0, 4); *g0 = prob->g[(f*prob->Nf + g)*4+0];}
651   if (g1) {PetscValidPointer(g1, 5); *g1 = prob->g[(f*prob->Nf + g)*4+1];}
652   if (g2) {PetscValidPointer(g2, 6); *g2 = prob->g[(f*prob->Nf + g)*4+2];}
653   if (g3) {PetscValidPointer(g3, 7); *g3 = prob->g[(f*prob->Nf + g)*4+3];}
654   PetscFunctionReturn(0);
655 }
656 
657 #undef __FUNCT__
658 #define __FUNCT__ "PetscDSSetJacobian"
659 PetscErrorCode PetscDSSetJacobian(PetscDS prob, PetscInt f, PetscInt g,
660                                        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[]),
661                                        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[]),
662                                        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[]),
663                                        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[]))
664 {
665   PetscErrorCode ierr;
666 
667   PetscFunctionBegin;
668   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
669   if (g0) PetscValidFunction(g0, 4);
670   if (g1) PetscValidFunction(g1, 5);
671   if (g2) PetscValidFunction(g2, 6);
672   if (g3) PetscValidFunction(g3, 7);
673   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
674   if (g < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", g);
675   ierr = PetscDSEnlarge_Static(prob, PetscMax(f, g)+1);CHKERRQ(ierr);
676   prob->g[(f*prob->Nf + g)*4+0] = g0;
677   prob->g[(f*prob->Nf + g)*4+1] = g1;
678   prob->g[(f*prob->Nf + g)*4+2] = g2;
679   prob->g[(f*prob->Nf + g)*4+3] = g3;
680   PetscFunctionReturn(0);
681 }
682 
683 #undef __FUNCT__
684 #define __FUNCT__ "PetscDSGetBdResidual"
685 PetscErrorCode PetscDSGetBdResidual(PetscDS prob, PetscInt f,
686                                          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[]),
687                                          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[]))
688 {
689   PetscFunctionBegin;
690   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
691   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);
692   if (f0) {PetscValidPointer(f0, 3); *f0 = prob->fBd[f*2+0];}
693   if (f1) {PetscValidPointer(f1, 4); *f1 = prob->fBd[f*2+1];}
694   PetscFunctionReturn(0);
695 }
696 
697 #undef __FUNCT__
698 #define __FUNCT__ "PetscDSSetBdResidual"
699 PetscErrorCode PetscDSSetBdResidual(PetscDS prob, PetscInt f,
700                                          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[]),
701                                          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[]))
702 {
703   PetscErrorCode ierr;
704 
705   PetscFunctionBegin;
706   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
707   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
708   ierr = PetscDSEnlarge_Static(prob, f+1);CHKERRQ(ierr);
709   if (f0) {PetscValidFunction(f0, 3); prob->fBd[f*2+0] = f0;}
710   if (f1) {PetscValidFunction(f1, 4); prob->fBd[f*2+1] = f1;}
711   PetscFunctionReturn(0);
712 }
713 
714 #undef __FUNCT__
715 #define __FUNCT__ "PetscDSGetBdJacobian"
716 PetscErrorCode PetscDSGetBdJacobian(PetscDS prob, PetscInt f, PetscInt g,
717                                          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[]),
718                                          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[]),
719                                          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[]),
720                                          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[]))
721 {
722   PetscFunctionBegin;
723   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
724   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);
725   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);
726   if (g0) {PetscValidPointer(g0, 4); *g0 = prob->gBd[(f*prob->Nf + g)*4+0];}
727   if (g1) {PetscValidPointer(g1, 5); *g1 = prob->gBd[(f*prob->Nf + g)*4+1];}
728   if (g2) {PetscValidPointer(g2, 6); *g2 = prob->gBd[(f*prob->Nf + g)*4+2];}
729   if (g3) {PetscValidPointer(g3, 7); *g3 = prob->gBd[(f*prob->Nf + g)*4+3];}
730   PetscFunctionReturn(0);
731 }
732 
733 #undef __FUNCT__
734 #define __FUNCT__ "PetscDSSetBdJacobian"
735 PetscErrorCode PetscDSSetBdJacobian(PetscDS prob, PetscInt f, PetscInt g,
736                                          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[]),
737                                          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[]),
738                                          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[]),
739                                          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[]))
740 {
741   PetscErrorCode ierr;
742 
743   PetscFunctionBegin;
744   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
745   if (g0) PetscValidFunction(g0, 4);
746   if (g1) PetscValidFunction(g1, 5);
747   if (g2) PetscValidFunction(g2, 6);
748   if (g3) PetscValidFunction(g3, 7);
749   if (f < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", f);
750   if (g < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Field number %d must be non-negative", g);
751   ierr = PetscDSEnlarge_Static(prob, PetscMax(f, g)+1);CHKERRQ(ierr);
752   prob->gBd[(f*prob->Nf + g)*4+0] = g0;
753   prob->gBd[(f*prob->Nf + g)*4+1] = g1;
754   prob->gBd[(f*prob->Nf + g)*4+2] = g2;
755   prob->gBd[(f*prob->Nf + g)*4+3] = g3;
756   PetscFunctionReturn(0);
757 }
758 
759 #undef __FUNCT__
760 #define __FUNCT__ "PetscDSGetFieldOffset"
761 PetscErrorCode PetscDSGetFieldOffset(PetscDS prob, PetscInt f, PetscInt *off)
762 {
763   PetscInt       g;
764   PetscErrorCode ierr;
765 
766   PetscFunctionBegin;
767   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
768   PetscValidPointer(off, 3);
769   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);
770   *off = 0;
771   for (g = 0; g < f; ++g) {
772     PetscFE  fe = (PetscFE) prob->disc[g];
773     PetscInt Nb, Nc;
774 
775     ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
776     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
777     *off += Nb*Nc;
778   }
779   PetscFunctionReturn(0);
780 }
781 
782 #undef __FUNCT__
783 #define __FUNCT__ "PetscDSGetBdFieldOffset"
784 PetscErrorCode PetscDSGetBdFieldOffset(PetscDS prob, PetscInt f, PetscInt *off)
785 {
786   PetscInt       g;
787   PetscErrorCode ierr;
788 
789   PetscFunctionBegin;
790   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
791   PetscValidPointer(off, 3);
792   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);
793   *off = 0;
794   for (g = 0; g < f; ++g) {
795     PetscFE  fe = (PetscFE) prob->discBd[g];
796     PetscInt Nb, Nc;
797 
798     ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
799     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
800     *off += Nb*Nc;
801   }
802   PetscFunctionReturn(0);
803 }
804 
805 #undef __FUNCT__
806 #define __FUNCT__ "PetscDSGetTabulation"
807 PetscErrorCode PetscDSGetTabulation(PetscDS prob, PetscReal ***basis, PetscReal ***basisDer)
808 {
809   PetscErrorCode ierr;
810 
811   PetscFunctionBegin;
812   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
813   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
814   if (basis)    {PetscValidPointer(basis, 2);    *basis    = prob->basis;}
815   if (basisDer) {PetscValidPointer(basisDer, 3); *basisDer = prob->basisDer;}
816   PetscFunctionReturn(0);
817 }
818 
819 #undef __FUNCT__
820 #define __FUNCT__ "PetscDSGetBdTabulation"
821 PetscErrorCode PetscDSGetBdTabulation(PetscDS prob, PetscReal ***basis, PetscReal ***basisDer)
822 {
823   PetscErrorCode ierr;
824 
825   PetscFunctionBegin;
826   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
827   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
828   if (basis)    {PetscValidPointer(basis, 2);    *basis    = prob->basisBd;}
829   if (basisDer) {PetscValidPointer(basisDer, 3); *basisDer = prob->basisDerBd;}
830   PetscFunctionReturn(0);
831 }
832 
833 #undef __FUNCT__
834 #define __FUNCT__ "PetscDSGetEvaluationArrays"
835 PetscErrorCode PetscDSGetEvaluationArrays(PetscDS prob, PetscScalar **u, PetscScalar **u_t, PetscScalar **u_x)
836 {
837   PetscErrorCode ierr;
838 
839   PetscFunctionBegin;
840   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
841   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
842   if (u)   {PetscValidPointer(u, 2);   *u   = prob->u;}
843   if (u_t) {PetscValidPointer(u_t, 3); *u_t = prob->u_t;}
844   if (u_x) {PetscValidPointer(u_x, 4); *u_x = prob->u_x;}
845   PetscFunctionReturn(0);
846 }
847 
848 #undef __FUNCT__
849 #define __FUNCT__ "PetscDSGetWeakFormArrays"
850 PetscErrorCode PetscDSGetWeakFormArrays(PetscDS prob, PetscScalar **f0, PetscScalar **f1, PetscScalar **g0, PetscScalar **g1, PetscScalar **g2, PetscScalar **g3)
851 {
852   PetscErrorCode ierr;
853 
854   PetscFunctionBegin;
855   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
856   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
857   if (f0) {PetscValidPointer(f0, 2); *f0 = prob->f0;}
858   if (f1) {PetscValidPointer(f1, 3); *f1 = prob->f1;}
859   if (g0) {PetscValidPointer(g0, 4); *g0 = prob->g0;}
860   if (g1) {PetscValidPointer(g1, 5); *g1 = prob->g1;}
861   if (g2) {PetscValidPointer(g2, 6); *g2 = prob->g2;}
862   if (g3) {PetscValidPointer(g3, 7); *g3 = prob->g3;}
863   PetscFunctionReturn(0);
864 }
865 
866 #undef __FUNCT__
867 #define __FUNCT__ "PetscDSGetRefCoordArrays"
868 PetscErrorCode PetscDSGetRefCoordArrays(PetscDS prob, PetscReal **x, PetscScalar **refSpaceDer)
869 {
870   PetscErrorCode ierr;
871 
872   PetscFunctionBegin;
873   PetscValidHeaderSpecific(prob, PETSCDS_CLASSID, 1);
874   ierr = PetscDSSetUp(prob);CHKERRQ(ierr);
875   if (x)           {PetscValidPointer(x, 2);           *x           = prob->x;}
876   if (refSpaceDer) {PetscValidPointer(refSpaceDer, 3); *refSpaceDer = prob->refSpaceDer;}
877   PetscFunctionReturn(0);
878 }
879 
880 #undef __FUNCT__
881 #define __FUNCT__ "PetscDSDestroy_Basic"
882 PetscErrorCode PetscDSDestroy_Basic(PetscDS prob)
883 {
884   PetscFunctionBegin;
885   PetscFunctionReturn(0);
886 }
887 
888 #undef __FUNCT__
889 #define __FUNCT__ "PetscDSInitialize_Basic"
890 PetscErrorCode PetscDSInitialize_Basic(PetscDS prob)
891 {
892   PetscFunctionBegin;
893   prob->ops->setfromoptions = NULL;
894   prob->ops->setup          = NULL;
895   prob->ops->view           = NULL;
896   prob->ops->destroy        = PetscDSDestroy_Basic;
897   PetscFunctionReturn(0);
898 }
899 
900 /*MC
901   PETSCDSBASIC = "basic" - A discrete system with pointwise residual and boundary residual functions
902 
903   Level: intermediate
904 
905 .seealso: PetscDSType, PetscDSCreate(), PetscDSSetType()
906 M*/
907 
908 #undef __FUNCT__
909 #define __FUNCT__ "PetscDSCreate_Basic"
910 PETSC_EXTERN PetscErrorCode PetscDSCreate_Basic(PetscDS prob)
911 {
912   PetscDS_Basic *b;
913   PetscErrorCode      ierr;
914 
915   PetscFunctionBegin;
916   PetscValidHeaderSpecific(prob, PETSCSPACE_CLASSID, 1);
917   ierr       = PetscNewLog(prob, &b);CHKERRQ(ierr);
918   prob->data = b;
919 
920   ierr = PetscDSInitialize_Basic(prob);CHKERRQ(ierr);
921   PetscFunctionReturn(0);
922 }
923