xref: /petsc/src/dm/dt/dualspace/interface/dualspace.c (revision 20cf1dd8cd656e27510885a71ea4be028bf48813)
1*20cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
2*20cf1dd8SToby Isaac #include <petscdmplex.h>
3*20cf1dd8SToby Isaac 
4*20cf1dd8SToby Isaac PetscClassId PETSCDUALSPACE_CLASSID = 0;
5*20cf1dd8SToby Isaac 
6*20cf1dd8SToby Isaac PetscFunctionList PetscDualSpaceList              = NULL;
7*20cf1dd8SToby Isaac PetscBool         PetscDualSpaceRegisterAllCalled = PETSC_FALSE;
8*20cf1dd8SToby Isaac 
9*20cf1dd8SToby Isaac /*@C
10*20cf1dd8SToby Isaac   PetscDualSpaceRegister - Adds a new PetscDualSpace implementation
11*20cf1dd8SToby Isaac 
12*20cf1dd8SToby Isaac   Not Collective
13*20cf1dd8SToby Isaac 
14*20cf1dd8SToby Isaac   Input Parameters:
15*20cf1dd8SToby Isaac + name        - The name of a new user-defined creation routine
16*20cf1dd8SToby Isaac - create_func - The creation routine itself
17*20cf1dd8SToby Isaac 
18*20cf1dd8SToby Isaac   Notes:
19*20cf1dd8SToby Isaac   PetscDualSpaceRegister() may be called multiple times to add several user-defined PetscDualSpaces
20*20cf1dd8SToby Isaac 
21*20cf1dd8SToby Isaac   Sample usage:
22*20cf1dd8SToby Isaac .vb
23*20cf1dd8SToby Isaac     PetscDualSpaceRegister("my_space", MyPetscDualSpaceCreate);
24*20cf1dd8SToby Isaac .ve
25*20cf1dd8SToby Isaac 
26*20cf1dd8SToby Isaac   Then, your PetscDualSpace type can be chosen with the procedural interface via
27*20cf1dd8SToby Isaac .vb
28*20cf1dd8SToby Isaac     PetscDualSpaceCreate(MPI_Comm, PetscDualSpace *);
29*20cf1dd8SToby Isaac     PetscDualSpaceSetType(PetscDualSpace, "my_dual_space");
30*20cf1dd8SToby Isaac .ve
31*20cf1dd8SToby Isaac    or at runtime via the option
32*20cf1dd8SToby Isaac .vb
33*20cf1dd8SToby Isaac     -petscdualspace_type my_dual_space
34*20cf1dd8SToby Isaac .ve
35*20cf1dd8SToby Isaac 
36*20cf1dd8SToby Isaac   Level: advanced
37*20cf1dd8SToby Isaac 
38*20cf1dd8SToby Isaac .keywords: PetscDualSpace, register
39*20cf1dd8SToby Isaac .seealso: PetscDualSpaceRegisterAll(), PetscDualSpaceRegisterDestroy()
40*20cf1dd8SToby Isaac 
41*20cf1dd8SToby Isaac @*/
42*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceRegister(const char sname[], PetscErrorCode (*function)(PetscDualSpace))
43*20cf1dd8SToby Isaac {
44*20cf1dd8SToby Isaac   PetscErrorCode ierr;
45*20cf1dd8SToby Isaac 
46*20cf1dd8SToby Isaac   PetscFunctionBegin;
47*20cf1dd8SToby Isaac   ierr = PetscFunctionListAdd(&PetscDualSpaceList, sname, function);CHKERRQ(ierr);
48*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
49*20cf1dd8SToby Isaac }
50*20cf1dd8SToby Isaac 
51*20cf1dd8SToby Isaac /*@C
52*20cf1dd8SToby Isaac   PetscDualSpaceSetType - Builds a particular PetscDualSpace
53*20cf1dd8SToby Isaac 
54*20cf1dd8SToby Isaac   Collective on PetscDualSpace
55*20cf1dd8SToby Isaac 
56*20cf1dd8SToby Isaac   Input Parameters:
57*20cf1dd8SToby Isaac + sp   - The PetscDualSpace object
58*20cf1dd8SToby Isaac - name - The kind of space
59*20cf1dd8SToby Isaac 
60*20cf1dd8SToby Isaac   Options Database Key:
61*20cf1dd8SToby Isaac . -petscdualspace_type <type> - Sets the PetscDualSpace type; use -help for a list of available types
62*20cf1dd8SToby Isaac 
63*20cf1dd8SToby Isaac   Level: intermediate
64*20cf1dd8SToby Isaac 
65*20cf1dd8SToby Isaac .keywords: PetscDualSpace, set, type
66*20cf1dd8SToby Isaac .seealso: PetscDualSpaceGetType(), PetscDualSpaceCreate()
67*20cf1dd8SToby Isaac @*/
68*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetType(PetscDualSpace sp, PetscDualSpaceType name)
69*20cf1dd8SToby Isaac {
70*20cf1dd8SToby Isaac   PetscErrorCode (*r)(PetscDualSpace);
71*20cf1dd8SToby Isaac   PetscBool      match;
72*20cf1dd8SToby Isaac   PetscErrorCode ierr;
73*20cf1dd8SToby Isaac 
74*20cf1dd8SToby Isaac   PetscFunctionBegin;
75*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
76*20cf1dd8SToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) sp, name, &match);CHKERRQ(ierr);
77*20cf1dd8SToby Isaac   if (match) PetscFunctionReturn(0);
78*20cf1dd8SToby Isaac 
79*20cf1dd8SToby Isaac   if (!PetscDualSpaceRegisterAllCalled) {ierr = PetscDualSpaceRegisterAll();CHKERRQ(ierr);}
80*20cf1dd8SToby Isaac   ierr = PetscFunctionListFind(PetscDualSpaceList, name, &r);CHKERRQ(ierr);
81*20cf1dd8SToby Isaac   if (!r) SETERRQ1(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscDualSpace type: %s", name);
82*20cf1dd8SToby Isaac 
83*20cf1dd8SToby Isaac   if (sp->ops->destroy) {
84*20cf1dd8SToby Isaac     ierr             = (*sp->ops->destroy)(sp);CHKERRQ(ierr);
85*20cf1dd8SToby Isaac     sp->ops->destroy = NULL;
86*20cf1dd8SToby Isaac   }
87*20cf1dd8SToby Isaac   ierr = (*r)(sp);CHKERRQ(ierr);
88*20cf1dd8SToby Isaac   ierr = PetscObjectChangeTypeName((PetscObject) sp, name);CHKERRQ(ierr);
89*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
90*20cf1dd8SToby Isaac }
91*20cf1dd8SToby Isaac 
92*20cf1dd8SToby Isaac /*@C
93*20cf1dd8SToby Isaac   PetscDualSpaceGetType - Gets the PetscDualSpace type name (as a string) from the object.
94*20cf1dd8SToby Isaac 
95*20cf1dd8SToby Isaac   Not Collective
96*20cf1dd8SToby Isaac 
97*20cf1dd8SToby Isaac   Input Parameter:
98*20cf1dd8SToby Isaac . sp  - The PetscDualSpace
99*20cf1dd8SToby Isaac 
100*20cf1dd8SToby Isaac   Output Parameter:
101*20cf1dd8SToby Isaac . name - The PetscDualSpace type name
102*20cf1dd8SToby Isaac 
103*20cf1dd8SToby Isaac   Level: intermediate
104*20cf1dd8SToby Isaac 
105*20cf1dd8SToby Isaac .keywords: PetscDualSpace, get, type, name
106*20cf1dd8SToby Isaac .seealso: PetscDualSpaceSetType(), PetscDualSpaceCreate()
107*20cf1dd8SToby Isaac @*/
108*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetType(PetscDualSpace sp, PetscDualSpaceType *name)
109*20cf1dd8SToby Isaac {
110*20cf1dd8SToby Isaac   PetscErrorCode ierr;
111*20cf1dd8SToby Isaac 
112*20cf1dd8SToby Isaac   PetscFunctionBegin;
113*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
114*20cf1dd8SToby Isaac   PetscValidPointer(name, 2);
115*20cf1dd8SToby Isaac   if (!PetscDualSpaceRegisterAllCalled) {
116*20cf1dd8SToby Isaac     ierr = PetscDualSpaceRegisterAll();CHKERRQ(ierr);
117*20cf1dd8SToby Isaac   }
118*20cf1dd8SToby Isaac   *name = ((PetscObject) sp)->type_name;
119*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
120*20cf1dd8SToby Isaac }
121*20cf1dd8SToby Isaac 
122*20cf1dd8SToby Isaac /*@
123*20cf1dd8SToby Isaac   PetscDualSpaceView - Views a PetscDualSpace
124*20cf1dd8SToby Isaac 
125*20cf1dd8SToby Isaac   Collective on PetscDualSpace
126*20cf1dd8SToby Isaac 
127*20cf1dd8SToby Isaac   Input Parameter:
128*20cf1dd8SToby Isaac + sp - the PetscDualSpace object to view
129*20cf1dd8SToby Isaac - v  - the viewer
130*20cf1dd8SToby Isaac 
131*20cf1dd8SToby Isaac   Level: developer
132*20cf1dd8SToby Isaac 
133*20cf1dd8SToby Isaac .seealso PetscDualSpaceDestroy()
134*20cf1dd8SToby Isaac @*/
135*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceView(PetscDualSpace sp, PetscViewer v)
136*20cf1dd8SToby Isaac {
137*20cf1dd8SToby Isaac   PetscErrorCode ierr;
138*20cf1dd8SToby Isaac 
139*20cf1dd8SToby Isaac   PetscFunctionBegin;
140*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
141*20cf1dd8SToby Isaac   if (!v) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) sp), &v);CHKERRQ(ierr);}
142*20cf1dd8SToby Isaac   if (sp->ops->view) {ierr = (*sp->ops->view)(sp, v);CHKERRQ(ierr);}
143*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
144*20cf1dd8SToby Isaac }
145*20cf1dd8SToby Isaac 
146*20cf1dd8SToby Isaac /*@
147*20cf1dd8SToby Isaac   PetscDualSpaceSetFromOptions - sets parameters in a PetscDualSpace from the options database
148*20cf1dd8SToby Isaac 
149*20cf1dd8SToby Isaac   Collective on PetscDualSpace
150*20cf1dd8SToby Isaac 
151*20cf1dd8SToby Isaac   Input Parameter:
152*20cf1dd8SToby Isaac . sp - the PetscDualSpace object to set options for
153*20cf1dd8SToby Isaac 
154*20cf1dd8SToby Isaac   Options Database:
155*20cf1dd8SToby Isaac . -petscspace_order the approximation order of the space
156*20cf1dd8SToby Isaac 
157*20cf1dd8SToby Isaac   Level: developer
158*20cf1dd8SToby Isaac 
159*20cf1dd8SToby Isaac .seealso PetscDualSpaceView()
160*20cf1dd8SToby Isaac @*/
161*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetFromOptions(PetscDualSpace sp)
162*20cf1dd8SToby Isaac {
163*20cf1dd8SToby Isaac   const char    *defaultType;
164*20cf1dd8SToby Isaac   char           name[256];
165*20cf1dd8SToby Isaac   PetscBool      flg;
166*20cf1dd8SToby Isaac   PetscErrorCode ierr;
167*20cf1dd8SToby Isaac 
168*20cf1dd8SToby Isaac   PetscFunctionBegin;
169*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
170*20cf1dd8SToby Isaac   if (!((PetscObject) sp)->type_name) {
171*20cf1dd8SToby Isaac     defaultType = PETSCDUALSPACELAGRANGE;
172*20cf1dd8SToby Isaac   } else {
173*20cf1dd8SToby Isaac     defaultType = ((PetscObject) sp)->type_name;
174*20cf1dd8SToby Isaac   }
175*20cf1dd8SToby Isaac   if (!PetscSpaceRegisterAllCalled) {ierr = PetscSpaceRegisterAll();CHKERRQ(ierr);}
176*20cf1dd8SToby Isaac 
177*20cf1dd8SToby Isaac   ierr = PetscObjectOptionsBegin((PetscObject) sp);CHKERRQ(ierr);
178*20cf1dd8SToby Isaac   ierr = PetscOptionsFList("-petscdualspace_type", "Dual space", "PetscDualSpaceSetType", PetscDualSpaceList, defaultType, name, 256, &flg);CHKERRQ(ierr);
179*20cf1dd8SToby Isaac   if (flg) {
180*20cf1dd8SToby Isaac     ierr = PetscDualSpaceSetType(sp, name);CHKERRQ(ierr);
181*20cf1dd8SToby Isaac   } else if (!((PetscObject) sp)->type_name) {
182*20cf1dd8SToby Isaac     ierr = PetscDualSpaceSetType(sp, defaultType);CHKERRQ(ierr);
183*20cf1dd8SToby Isaac   }
184*20cf1dd8SToby Isaac   ierr = PetscOptionsInt("-petscdualspace_order", "The approximation order", "PetscDualSpaceSetOrder", sp->order, &sp->order, NULL);CHKERRQ(ierr);
185*20cf1dd8SToby Isaac   ierr = PetscOptionsInt("-petscdualspace_components", "The number of components", "PetscDualSpaceSetNumComponents", sp->Nc, &sp->Nc, NULL);CHKERRQ(ierr);
186*20cf1dd8SToby Isaac   if (sp->ops->setfromoptions) {
187*20cf1dd8SToby Isaac     ierr = (*sp->ops->setfromoptions)(PetscOptionsObject,sp);CHKERRQ(ierr);
188*20cf1dd8SToby Isaac   }
189*20cf1dd8SToby Isaac   /* process any options handlers added with PetscObjectAddOptionsHandler() */
190*20cf1dd8SToby Isaac   ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) sp);CHKERRQ(ierr);
191*20cf1dd8SToby Isaac   ierr = PetscOptionsEnd();CHKERRQ(ierr);
192*20cf1dd8SToby Isaac   ierr = PetscDualSpaceViewFromOptions(sp, NULL, "-petscdualspace_view");CHKERRQ(ierr);
193*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
194*20cf1dd8SToby Isaac }
195*20cf1dd8SToby Isaac 
196*20cf1dd8SToby Isaac /*@
197*20cf1dd8SToby Isaac   PetscDualSpaceSetUp - Construct a basis for the PetscDualSpace
198*20cf1dd8SToby Isaac 
199*20cf1dd8SToby Isaac   Collective on PetscDualSpace
200*20cf1dd8SToby Isaac 
201*20cf1dd8SToby Isaac   Input Parameter:
202*20cf1dd8SToby Isaac . sp - the PetscDualSpace object to setup
203*20cf1dd8SToby Isaac 
204*20cf1dd8SToby Isaac   Level: developer
205*20cf1dd8SToby Isaac 
206*20cf1dd8SToby Isaac .seealso PetscDualSpaceView(), PetscDualSpaceDestroy()
207*20cf1dd8SToby Isaac @*/
208*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetUp(PetscDualSpace sp)
209*20cf1dd8SToby Isaac {
210*20cf1dd8SToby Isaac   PetscErrorCode ierr;
211*20cf1dd8SToby Isaac 
212*20cf1dd8SToby Isaac   PetscFunctionBegin;
213*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
214*20cf1dd8SToby Isaac   if (sp->setupcalled) PetscFunctionReturn(0);
215*20cf1dd8SToby Isaac   sp->setupcalled = PETSC_TRUE;
216*20cf1dd8SToby Isaac   if (sp->ops->setup) {ierr = (*sp->ops->setup)(sp);CHKERRQ(ierr);}
217*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
218*20cf1dd8SToby Isaac }
219*20cf1dd8SToby Isaac 
220*20cf1dd8SToby Isaac /*@
221*20cf1dd8SToby Isaac   PetscDualSpaceDestroy - Destroys a PetscDualSpace object
222*20cf1dd8SToby Isaac 
223*20cf1dd8SToby Isaac   Collective on PetscDualSpace
224*20cf1dd8SToby Isaac 
225*20cf1dd8SToby Isaac   Input Parameter:
226*20cf1dd8SToby Isaac . sp - the PetscDualSpace object to destroy
227*20cf1dd8SToby Isaac 
228*20cf1dd8SToby Isaac   Level: developer
229*20cf1dd8SToby Isaac 
230*20cf1dd8SToby Isaac .seealso PetscDualSpaceView()
231*20cf1dd8SToby Isaac @*/
232*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceDestroy(PetscDualSpace *sp)
233*20cf1dd8SToby Isaac {
234*20cf1dd8SToby Isaac   PetscInt       dim, f;
235*20cf1dd8SToby Isaac   PetscErrorCode ierr;
236*20cf1dd8SToby Isaac 
237*20cf1dd8SToby Isaac   PetscFunctionBegin;
238*20cf1dd8SToby Isaac   if (!*sp) PetscFunctionReturn(0);
239*20cf1dd8SToby Isaac   PetscValidHeaderSpecific((*sp), PETSCDUALSPACE_CLASSID, 1);
240*20cf1dd8SToby Isaac 
241*20cf1dd8SToby Isaac   if (--((PetscObject)(*sp))->refct > 0) {*sp = 0; PetscFunctionReturn(0);}
242*20cf1dd8SToby Isaac   ((PetscObject) (*sp))->refct = 0;
243*20cf1dd8SToby Isaac 
244*20cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDimension(*sp, &dim);CHKERRQ(ierr);
245*20cf1dd8SToby Isaac   for (f = 0; f < dim; ++f) {
246*20cf1dd8SToby Isaac     ierr = PetscQuadratureDestroy(&(*sp)->functional[f]);CHKERRQ(ierr);
247*20cf1dd8SToby Isaac   }
248*20cf1dd8SToby Isaac   ierr = PetscFree((*sp)->functional);CHKERRQ(ierr);
249*20cf1dd8SToby Isaac   ierr = PetscQuadratureDestroy(&(*sp)->allPoints);CHKERRQ(ierr);
250*20cf1dd8SToby Isaac   ierr = DMDestroy(&(*sp)->dm);CHKERRQ(ierr);
251*20cf1dd8SToby Isaac 
252*20cf1dd8SToby Isaac   if ((*sp)->ops->destroy) {ierr = (*(*sp)->ops->destroy)(*sp);CHKERRQ(ierr);}
253*20cf1dd8SToby Isaac   ierr = PetscHeaderDestroy(sp);CHKERRQ(ierr);
254*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
255*20cf1dd8SToby Isaac }
256*20cf1dd8SToby Isaac 
257*20cf1dd8SToby Isaac /*@
258*20cf1dd8SToby Isaac   PetscDualSpaceCreate - Creates an empty PetscDualSpace object. The type can then be set with PetscDualSpaceSetType().
259*20cf1dd8SToby Isaac 
260*20cf1dd8SToby Isaac   Collective on MPI_Comm
261*20cf1dd8SToby Isaac 
262*20cf1dd8SToby Isaac   Input Parameter:
263*20cf1dd8SToby Isaac . comm - The communicator for the PetscDualSpace object
264*20cf1dd8SToby Isaac 
265*20cf1dd8SToby Isaac   Output Parameter:
266*20cf1dd8SToby Isaac . sp - The PetscDualSpace object
267*20cf1dd8SToby Isaac 
268*20cf1dd8SToby Isaac   Level: beginner
269*20cf1dd8SToby Isaac 
270*20cf1dd8SToby Isaac .seealso: PetscDualSpaceSetType(), PETSCDUALSPACELAGRANGE
271*20cf1dd8SToby Isaac @*/
272*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceCreate(MPI_Comm comm, PetscDualSpace *sp)
273*20cf1dd8SToby Isaac {
274*20cf1dd8SToby Isaac   PetscDualSpace s;
275*20cf1dd8SToby Isaac   PetscErrorCode ierr;
276*20cf1dd8SToby Isaac 
277*20cf1dd8SToby Isaac   PetscFunctionBegin;
278*20cf1dd8SToby Isaac   PetscValidPointer(sp, 2);
279*20cf1dd8SToby Isaac   ierr = PetscCitationsRegister(FECitation,&FEcite);CHKERRQ(ierr);
280*20cf1dd8SToby Isaac   *sp  = NULL;
281*20cf1dd8SToby Isaac   ierr = PetscFEInitializePackage();CHKERRQ(ierr);
282*20cf1dd8SToby Isaac 
283*20cf1dd8SToby Isaac   ierr = PetscHeaderCreate(s, PETSCDUALSPACE_CLASSID, "PetscDualSpace", "Dual Space", "PetscDualSpace", comm, PetscDualSpaceDestroy, PetscDualSpaceView);CHKERRQ(ierr);
284*20cf1dd8SToby Isaac 
285*20cf1dd8SToby Isaac   s->order = 0;
286*20cf1dd8SToby Isaac   s->Nc    = 1;
287*20cf1dd8SToby Isaac   s->setupcalled = PETSC_FALSE;
288*20cf1dd8SToby Isaac 
289*20cf1dd8SToby Isaac   *sp = s;
290*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
291*20cf1dd8SToby Isaac }
292*20cf1dd8SToby Isaac 
293*20cf1dd8SToby Isaac /*@
294*20cf1dd8SToby Isaac   PetscDualSpaceDuplicate - Creates a duplicate PetscDualSpace object, however it is not setup.
295*20cf1dd8SToby Isaac 
296*20cf1dd8SToby Isaac   Collective on PetscDualSpace
297*20cf1dd8SToby Isaac 
298*20cf1dd8SToby Isaac   Input Parameter:
299*20cf1dd8SToby Isaac . sp - The original PetscDualSpace
300*20cf1dd8SToby Isaac 
301*20cf1dd8SToby Isaac   Output Parameter:
302*20cf1dd8SToby Isaac . spNew - The duplicate PetscDualSpace
303*20cf1dd8SToby Isaac 
304*20cf1dd8SToby Isaac   Level: beginner
305*20cf1dd8SToby Isaac 
306*20cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate(), PetscDualSpaceSetType()
307*20cf1dd8SToby Isaac @*/
308*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceDuplicate(PetscDualSpace sp, PetscDualSpace *spNew)
309*20cf1dd8SToby Isaac {
310*20cf1dd8SToby Isaac   PetscErrorCode ierr;
311*20cf1dd8SToby Isaac 
312*20cf1dd8SToby Isaac   PetscFunctionBegin;
313*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
314*20cf1dd8SToby Isaac   PetscValidPointer(spNew, 2);
315*20cf1dd8SToby Isaac   ierr = (*sp->ops->duplicate)(sp, spNew);CHKERRQ(ierr);
316*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
317*20cf1dd8SToby Isaac }
318*20cf1dd8SToby Isaac 
319*20cf1dd8SToby Isaac /*@
320*20cf1dd8SToby Isaac   PetscDualSpaceGetDM - Get the DM representing the reference cell
321*20cf1dd8SToby Isaac 
322*20cf1dd8SToby Isaac   Not collective
323*20cf1dd8SToby Isaac 
324*20cf1dd8SToby Isaac   Input Parameter:
325*20cf1dd8SToby Isaac . sp - The PetscDualSpace
326*20cf1dd8SToby Isaac 
327*20cf1dd8SToby Isaac   Output Parameter:
328*20cf1dd8SToby Isaac . dm - The reference cell
329*20cf1dd8SToby Isaac 
330*20cf1dd8SToby Isaac   Level: intermediate
331*20cf1dd8SToby Isaac 
332*20cf1dd8SToby Isaac .seealso: PetscDualSpaceSetDM(), PetscDualSpaceCreate()
333*20cf1dd8SToby Isaac @*/
334*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetDM(PetscDualSpace sp, DM *dm)
335*20cf1dd8SToby Isaac {
336*20cf1dd8SToby Isaac   PetscFunctionBegin;
337*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
338*20cf1dd8SToby Isaac   PetscValidPointer(dm, 2);
339*20cf1dd8SToby Isaac   *dm = sp->dm;
340*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
341*20cf1dd8SToby Isaac }
342*20cf1dd8SToby Isaac 
343*20cf1dd8SToby Isaac /*@
344*20cf1dd8SToby Isaac   PetscDualSpaceSetDM - Get the DM representing the reference cell
345*20cf1dd8SToby Isaac 
346*20cf1dd8SToby Isaac   Not collective
347*20cf1dd8SToby Isaac 
348*20cf1dd8SToby Isaac   Input Parameters:
349*20cf1dd8SToby Isaac + sp - The PetscDualSpace
350*20cf1dd8SToby Isaac - dm - The reference cell
351*20cf1dd8SToby Isaac 
352*20cf1dd8SToby Isaac   Level: intermediate
353*20cf1dd8SToby Isaac 
354*20cf1dd8SToby Isaac .seealso: PetscDualSpaceGetDM(), PetscDualSpaceCreate()
355*20cf1dd8SToby Isaac @*/
356*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetDM(PetscDualSpace sp, DM dm)
357*20cf1dd8SToby Isaac {
358*20cf1dd8SToby Isaac   PetscErrorCode ierr;
359*20cf1dd8SToby Isaac 
360*20cf1dd8SToby Isaac   PetscFunctionBegin;
361*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
362*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
363*20cf1dd8SToby Isaac   ierr = DMDestroy(&sp->dm);CHKERRQ(ierr);
364*20cf1dd8SToby Isaac   ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);
365*20cf1dd8SToby Isaac   sp->dm = dm;
366*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
367*20cf1dd8SToby Isaac }
368*20cf1dd8SToby Isaac 
369*20cf1dd8SToby Isaac /*@
370*20cf1dd8SToby Isaac   PetscDualSpaceGetOrder - Get the order of the dual space
371*20cf1dd8SToby Isaac 
372*20cf1dd8SToby Isaac   Not collective
373*20cf1dd8SToby Isaac 
374*20cf1dd8SToby Isaac   Input Parameter:
375*20cf1dd8SToby Isaac . sp - The PetscDualSpace
376*20cf1dd8SToby Isaac 
377*20cf1dd8SToby Isaac   Output Parameter:
378*20cf1dd8SToby Isaac . order - The order
379*20cf1dd8SToby Isaac 
380*20cf1dd8SToby Isaac   Level: intermediate
381*20cf1dd8SToby Isaac 
382*20cf1dd8SToby Isaac .seealso: PetscDualSpaceSetOrder(), PetscDualSpaceCreate()
383*20cf1dd8SToby Isaac @*/
384*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetOrder(PetscDualSpace sp, PetscInt *order)
385*20cf1dd8SToby Isaac {
386*20cf1dd8SToby Isaac   PetscFunctionBegin;
387*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
388*20cf1dd8SToby Isaac   PetscValidPointer(order, 2);
389*20cf1dd8SToby Isaac   *order = sp->order;
390*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
391*20cf1dd8SToby Isaac }
392*20cf1dd8SToby Isaac 
393*20cf1dd8SToby Isaac /*@
394*20cf1dd8SToby Isaac   PetscDualSpaceSetOrder - Set the order of the dual space
395*20cf1dd8SToby Isaac 
396*20cf1dd8SToby Isaac   Not collective
397*20cf1dd8SToby Isaac 
398*20cf1dd8SToby Isaac   Input Parameters:
399*20cf1dd8SToby Isaac + sp - The PetscDualSpace
400*20cf1dd8SToby Isaac - order - The order
401*20cf1dd8SToby Isaac 
402*20cf1dd8SToby Isaac   Level: intermediate
403*20cf1dd8SToby Isaac 
404*20cf1dd8SToby Isaac .seealso: PetscDualSpaceGetOrder(), PetscDualSpaceCreate()
405*20cf1dd8SToby Isaac @*/
406*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetOrder(PetscDualSpace sp, PetscInt order)
407*20cf1dd8SToby Isaac {
408*20cf1dd8SToby Isaac   PetscFunctionBegin;
409*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
410*20cf1dd8SToby Isaac   sp->order = order;
411*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
412*20cf1dd8SToby Isaac }
413*20cf1dd8SToby Isaac 
414*20cf1dd8SToby Isaac /*@
415*20cf1dd8SToby Isaac   PetscDualSpaceGetNumComponents - Return the number of components for this space
416*20cf1dd8SToby Isaac 
417*20cf1dd8SToby Isaac   Input Parameter:
418*20cf1dd8SToby Isaac . sp - The PetscDualSpace
419*20cf1dd8SToby Isaac 
420*20cf1dd8SToby Isaac   Output Parameter:
421*20cf1dd8SToby Isaac . Nc - The number of components
422*20cf1dd8SToby Isaac 
423*20cf1dd8SToby Isaac   Note: A vector space, for example, will have d components, where d is the spatial dimension
424*20cf1dd8SToby Isaac 
425*20cf1dd8SToby Isaac   Level: intermediate
426*20cf1dd8SToby Isaac 
427*20cf1dd8SToby Isaac .seealso: PetscDualSpaceSetNumComponents(), PetscDualSpaceGetDimension(), PetscDualSpaceCreate(), PetscDualSpace
428*20cf1dd8SToby Isaac @*/
429*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetNumComponents(PetscDualSpace sp, PetscInt *Nc)
430*20cf1dd8SToby Isaac {
431*20cf1dd8SToby Isaac   PetscFunctionBegin;
432*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
433*20cf1dd8SToby Isaac   PetscValidPointer(Nc, 2);
434*20cf1dd8SToby Isaac   *Nc = sp->Nc;
435*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
436*20cf1dd8SToby Isaac }
437*20cf1dd8SToby Isaac 
438*20cf1dd8SToby Isaac /*@
439*20cf1dd8SToby Isaac   PetscDualSpaceSetNumComponents - Set the number of components for this space
440*20cf1dd8SToby Isaac 
441*20cf1dd8SToby Isaac   Input Parameters:
442*20cf1dd8SToby Isaac + sp - The PetscDualSpace
443*20cf1dd8SToby Isaac - order - The number of components
444*20cf1dd8SToby Isaac 
445*20cf1dd8SToby Isaac   Level: intermediate
446*20cf1dd8SToby Isaac 
447*20cf1dd8SToby Isaac .seealso: PetscDualSpaceGetNumComponents(), PetscDualSpaceCreate(), PetscDualSpace
448*20cf1dd8SToby Isaac @*/
449*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetNumComponents(PetscDualSpace sp, PetscInt Nc)
450*20cf1dd8SToby Isaac {
451*20cf1dd8SToby Isaac   PetscFunctionBegin;
452*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
453*20cf1dd8SToby Isaac   sp->Nc = Nc;
454*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
455*20cf1dd8SToby Isaac }
456*20cf1dd8SToby Isaac 
457*20cf1dd8SToby Isaac /*@
458*20cf1dd8SToby Isaac   PetscDualSpaceLagrangeGetTensor - Get the tensor nature of the dual space
459*20cf1dd8SToby Isaac 
460*20cf1dd8SToby Isaac   Not collective
461*20cf1dd8SToby Isaac 
462*20cf1dd8SToby Isaac   Input Parameter:
463*20cf1dd8SToby Isaac . sp - The PetscDualSpace
464*20cf1dd8SToby Isaac 
465*20cf1dd8SToby Isaac   Output Parameter:
466*20cf1dd8SToby Isaac . tensor - Whether the dual space has tensor layout (vs. simplicial)
467*20cf1dd8SToby Isaac 
468*20cf1dd8SToby Isaac   Level: intermediate
469*20cf1dd8SToby Isaac 
470*20cf1dd8SToby Isaac .seealso: PetscDualSpaceLagrangeSetTensor(), PetscDualSpaceCreate()
471*20cf1dd8SToby Isaac @*/
472*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceLagrangeGetTensor(PetscDualSpace sp, PetscBool *tensor)
473*20cf1dd8SToby Isaac {
474*20cf1dd8SToby Isaac   PetscErrorCode ierr;
475*20cf1dd8SToby Isaac 
476*20cf1dd8SToby Isaac   PetscFunctionBegin;
477*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
478*20cf1dd8SToby Isaac   PetscValidPointer(tensor, 2);
479*20cf1dd8SToby Isaac   ierr = PetscTryMethod(sp,"PetscDualSpaceLagrangeGetTensor_C",(PetscDualSpace,PetscBool *),(sp,tensor));CHKERRQ(ierr);
480*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
481*20cf1dd8SToby Isaac }
482*20cf1dd8SToby Isaac 
483*20cf1dd8SToby Isaac /*@
484*20cf1dd8SToby Isaac   PetscDualSpaceLagrangeSetTensor - Set the tensor nature of the dual space
485*20cf1dd8SToby Isaac 
486*20cf1dd8SToby Isaac   Not collective
487*20cf1dd8SToby Isaac 
488*20cf1dd8SToby Isaac   Input Parameters:
489*20cf1dd8SToby Isaac + sp - The PetscDualSpace
490*20cf1dd8SToby Isaac - tensor - Whether the dual space has tensor layout (vs. simplicial)
491*20cf1dd8SToby Isaac 
492*20cf1dd8SToby Isaac   Level: intermediate
493*20cf1dd8SToby Isaac 
494*20cf1dd8SToby Isaac .seealso: PetscDualSpaceLagrangeGetTensor(), PetscDualSpaceCreate()
495*20cf1dd8SToby Isaac @*/
496*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceLagrangeSetTensor(PetscDualSpace sp, PetscBool tensor)
497*20cf1dd8SToby Isaac {
498*20cf1dd8SToby Isaac   PetscErrorCode ierr;
499*20cf1dd8SToby Isaac 
500*20cf1dd8SToby Isaac   PetscFunctionBegin;
501*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
502*20cf1dd8SToby Isaac   ierr = PetscTryMethod(sp,"PetscDualSpaceLagrangeSetTensor_C",(PetscDualSpace,PetscBool),(sp,tensor));CHKERRQ(ierr);
503*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
504*20cf1dd8SToby Isaac }
505*20cf1dd8SToby Isaac 
506*20cf1dd8SToby Isaac /*@
507*20cf1dd8SToby Isaac   PetscDualSpaceGetFunctional - Get the i-th basis functional in the dual space
508*20cf1dd8SToby Isaac 
509*20cf1dd8SToby Isaac   Not collective
510*20cf1dd8SToby Isaac 
511*20cf1dd8SToby Isaac   Input Parameters:
512*20cf1dd8SToby Isaac + sp - The PetscDualSpace
513*20cf1dd8SToby Isaac - i  - The basis number
514*20cf1dd8SToby Isaac 
515*20cf1dd8SToby Isaac   Output Parameter:
516*20cf1dd8SToby Isaac . functional - The basis functional
517*20cf1dd8SToby Isaac 
518*20cf1dd8SToby Isaac   Level: intermediate
519*20cf1dd8SToby Isaac 
520*20cf1dd8SToby Isaac .seealso: PetscDualSpaceGetDimension(), PetscDualSpaceCreate()
521*20cf1dd8SToby Isaac @*/
522*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetFunctional(PetscDualSpace sp, PetscInt i, PetscQuadrature *functional)
523*20cf1dd8SToby Isaac {
524*20cf1dd8SToby Isaac   PetscInt       dim;
525*20cf1dd8SToby Isaac   PetscErrorCode ierr;
526*20cf1dd8SToby Isaac 
527*20cf1dd8SToby Isaac   PetscFunctionBegin;
528*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
529*20cf1dd8SToby Isaac   PetscValidPointer(functional, 3);
530*20cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDimension(sp, &dim);CHKERRQ(ierr);
531*20cf1dd8SToby Isaac   if ((i < 0) || (i >= dim)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Functional index %d must be in [0, %d)", i, dim);
532*20cf1dd8SToby Isaac   *functional = sp->functional[i];
533*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
534*20cf1dd8SToby Isaac }
535*20cf1dd8SToby Isaac 
536*20cf1dd8SToby Isaac /*@
537*20cf1dd8SToby Isaac   PetscDualSpaceGetDimension - Get the dimension of the dual space, i.e. the number of basis functionals
538*20cf1dd8SToby Isaac 
539*20cf1dd8SToby Isaac   Not collective
540*20cf1dd8SToby Isaac 
541*20cf1dd8SToby Isaac   Input Parameter:
542*20cf1dd8SToby Isaac . sp - The PetscDualSpace
543*20cf1dd8SToby Isaac 
544*20cf1dd8SToby Isaac   Output Parameter:
545*20cf1dd8SToby Isaac . dim - The dimension
546*20cf1dd8SToby Isaac 
547*20cf1dd8SToby Isaac   Level: intermediate
548*20cf1dd8SToby Isaac 
549*20cf1dd8SToby Isaac .seealso: PetscDualSpaceGetFunctional(), PetscDualSpaceCreate()
550*20cf1dd8SToby Isaac @*/
551*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetDimension(PetscDualSpace sp, PetscInt *dim)
552*20cf1dd8SToby Isaac {
553*20cf1dd8SToby Isaac   PetscErrorCode ierr;
554*20cf1dd8SToby Isaac 
555*20cf1dd8SToby Isaac   PetscFunctionBegin;
556*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
557*20cf1dd8SToby Isaac   PetscValidPointer(dim, 2);
558*20cf1dd8SToby Isaac   *dim = 0;
559*20cf1dd8SToby Isaac   if (sp->ops->getdimension) {ierr = (*sp->ops->getdimension)(sp, dim);CHKERRQ(ierr);}
560*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
561*20cf1dd8SToby Isaac }
562*20cf1dd8SToby Isaac 
563*20cf1dd8SToby Isaac /*@C
564*20cf1dd8SToby Isaac   PetscDualSpaceGetNumDof - Get the number of degrees of freedom for each spatial (topological) dimension
565*20cf1dd8SToby Isaac 
566*20cf1dd8SToby Isaac   Not collective
567*20cf1dd8SToby Isaac 
568*20cf1dd8SToby Isaac   Input Parameter:
569*20cf1dd8SToby Isaac . sp - The PetscDualSpace
570*20cf1dd8SToby Isaac 
571*20cf1dd8SToby Isaac   Output Parameter:
572*20cf1dd8SToby Isaac . numDof - An array of length dim+1 which holds the number of dofs for each dimension
573*20cf1dd8SToby Isaac 
574*20cf1dd8SToby Isaac   Level: intermediate
575*20cf1dd8SToby Isaac 
576*20cf1dd8SToby Isaac .seealso: PetscDualSpaceGetFunctional(), PetscDualSpaceCreate()
577*20cf1dd8SToby Isaac @*/
578*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetNumDof(PetscDualSpace sp, const PetscInt **numDof)
579*20cf1dd8SToby Isaac {
580*20cf1dd8SToby Isaac   PetscErrorCode ierr;
581*20cf1dd8SToby Isaac 
582*20cf1dd8SToby Isaac   PetscFunctionBegin;
583*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
584*20cf1dd8SToby Isaac   PetscValidPointer(numDof, 2);
585*20cf1dd8SToby Isaac   ierr = (*sp->ops->getnumdof)(sp, numDof);CHKERRQ(ierr);
586*20cf1dd8SToby Isaac   if (!*numDof) SETERRQ(PetscObjectComm((PetscObject) sp), PETSC_ERR_LIB, "Empty numDof[] returned from dual space implementation");
587*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
588*20cf1dd8SToby Isaac }
589*20cf1dd8SToby Isaac 
590*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceCreateSection(PetscDualSpace sp, PetscSection *section)
591*20cf1dd8SToby Isaac {
592*20cf1dd8SToby Isaac   DM             dm;
593*20cf1dd8SToby Isaac   PetscInt       pStart, pEnd, depth, h, offset;
594*20cf1dd8SToby Isaac   const PetscInt *numDof;
595*20cf1dd8SToby Isaac   PetscErrorCode ierr;
596*20cf1dd8SToby Isaac 
597*20cf1dd8SToby Isaac   PetscFunctionBegin;
598*20cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDM(sp,&dm);CHKERRQ(ierr);
599*20cf1dd8SToby Isaac   ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr);
600*20cf1dd8SToby Isaac   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)sp),section);CHKERRQ(ierr);
601*20cf1dd8SToby Isaac   ierr = PetscSectionSetChart(*section,pStart,pEnd);CHKERRQ(ierr);
602*20cf1dd8SToby Isaac   ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr);
603*20cf1dd8SToby Isaac   ierr = PetscDualSpaceGetNumDof(sp,&numDof);CHKERRQ(ierr);
604*20cf1dd8SToby Isaac   for (h = 0; h <= depth; h++) {
605*20cf1dd8SToby Isaac     PetscInt hStart, hEnd, p, dof;
606*20cf1dd8SToby Isaac 
607*20cf1dd8SToby Isaac     ierr = DMPlexGetHeightStratum(dm,h,&hStart,&hEnd);CHKERRQ(ierr);
608*20cf1dd8SToby Isaac     dof = numDof[depth - h];
609*20cf1dd8SToby Isaac     for (p = hStart; p < hEnd; p++) {
610*20cf1dd8SToby Isaac       ierr = PetscSectionSetDof(*section,p,dof);CHKERRQ(ierr);
611*20cf1dd8SToby Isaac     }
612*20cf1dd8SToby Isaac   }
613*20cf1dd8SToby Isaac   ierr = PetscSectionSetUp(*section);CHKERRQ(ierr);
614*20cf1dd8SToby Isaac   for (h = 0, offset = 0; h <= depth; h++) {
615*20cf1dd8SToby Isaac     PetscInt hStart, hEnd, p, dof;
616*20cf1dd8SToby Isaac 
617*20cf1dd8SToby Isaac     ierr = DMPlexGetHeightStratum(dm,h,&hStart,&hEnd);CHKERRQ(ierr);
618*20cf1dd8SToby Isaac     dof = numDof[depth - h];
619*20cf1dd8SToby Isaac     for (p = hStart; p < hEnd; p++) {
620*20cf1dd8SToby Isaac       ierr = PetscSectionGetDof(*section,p,&dof);CHKERRQ(ierr);
621*20cf1dd8SToby Isaac       ierr = PetscSectionSetOffset(*section,p,offset);CHKERRQ(ierr);
622*20cf1dd8SToby Isaac       offset += dof;
623*20cf1dd8SToby Isaac     }
624*20cf1dd8SToby Isaac   }
625*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
626*20cf1dd8SToby Isaac }
627*20cf1dd8SToby Isaac 
628*20cf1dd8SToby Isaac /*@
629*20cf1dd8SToby Isaac   PetscDualSpaceCreateReferenceCell - Create a DMPLEX with the appropriate FEM reference cell
630*20cf1dd8SToby Isaac 
631*20cf1dd8SToby Isaac   Collective on PetscDualSpace
632*20cf1dd8SToby Isaac 
633*20cf1dd8SToby Isaac   Input Parameters:
634*20cf1dd8SToby Isaac + sp      - The PetscDualSpace
635*20cf1dd8SToby Isaac . dim     - The spatial dimension
636*20cf1dd8SToby Isaac - simplex - Flag for simplex, otherwise use a tensor-product cell
637*20cf1dd8SToby Isaac 
638*20cf1dd8SToby Isaac   Output Parameter:
639*20cf1dd8SToby Isaac . refdm - The reference cell
640*20cf1dd8SToby Isaac 
641*20cf1dd8SToby Isaac   Level: advanced
642*20cf1dd8SToby Isaac 
643*20cf1dd8SToby Isaac .keywords: PetscDualSpace, reference cell
644*20cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate(), DMPLEX
645*20cf1dd8SToby Isaac @*/
646*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceCreateReferenceCell(PetscDualSpace sp, PetscInt dim, PetscBool simplex, DM *refdm)
647*20cf1dd8SToby Isaac {
648*20cf1dd8SToby Isaac   PetscErrorCode ierr;
649*20cf1dd8SToby Isaac 
650*20cf1dd8SToby Isaac   PetscFunctionBeginUser;
651*20cf1dd8SToby Isaac   ierr = DMPlexCreateReferenceCell(PetscObjectComm((PetscObject) sp), dim, simplex, refdm);CHKERRQ(ierr);
652*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
653*20cf1dd8SToby Isaac }
654*20cf1dd8SToby Isaac 
655*20cf1dd8SToby Isaac /*@C
656*20cf1dd8SToby Isaac   PetscDualSpaceApply - Apply a functional from the dual space basis to an input function
657*20cf1dd8SToby Isaac 
658*20cf1dd8SToby Isaac   Input Parameters:
659*20cf1dd8SToby Isaac + sp      - The PetscDualSpace object
660*20cf1dd8SToby Isaac . f       - The basis functional index
661*20cf1dd8SToby Isaac . time    - The time
662*20cf1dd8SToby Isaac . cgeom   - A context with geometric information for this cell, we use v0 (the initial vertex) and J (the Jacobian) (or evaluated at the coordinates of the functional)
663*20cf1dd8SToby Isaac . numComp - The number of components for the function
664*20cf1dd8SToby Isaac . func    - The input function
665*20cf1dd8SToby Isaac - ctx     - A context for the function
666*20cf1dd8SToby Isaac 
667*20cf1dd8SToby Isaac   Output Parameter:
668*20cf1dd8SToby Isaac . value   - numComp output values
669*20cf1dd8SToby Isaac 
670*20cf1dd8SToby Isaac   Note: The calling sequence for the callback func is given by:
671*20cf1dd8SToby Isaac 
672*20cf1dd8SToby Isaac $ func(PetscInt dim, PetscReal time, const PetscReal x[],
673*20cf1dd8SToby Isaac $      PetscInt numComponents, PetscScalar values[], void *ctx)
674*20cf1dd8SToby Isaac 
675*20cf1dd8SToby Isaac   Level: developer
676*20cf1dd8SToby Isaac 
677*20cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate()
678*20cf1dd8SToby Isaac @*/
679*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceApply(PetscDualSpace sp, PetscInt f, PetscReal time, PetscFEGeom *cgeom, PetscInt numComp, PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void *ctx, PetscScalar *value)
680*20cf1dd8SToby Isaac {
681*20cf1dd8SToby Isaac   PetscErrorCode ierr;
682*20cf1dd8SToby Isaac 
683*20cf1dd8SToby Isaac   PetscFunctionBegin;
684*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
685*20cf1dd8SToby Isaac   PetscValidPointer(cgeom, 4);
686*20cf1dd8SToby Isaac   PetscValidPointer(value, 8);
687*20cf1dd8SToby Isaac   ierr = (*sp->ops->apply)(sp, f, time, cgeom, numComp, func, ctx, value);CHKERRQ(ierr);
688*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
689*20cf1dd8SToby Isaac }
690*20cf1dd8SToby Isaac 
691*20cf1dd8SToby Isaac /*@C
692*20cf1dd8SToby Isaac   PetscDualSpaceApplyAll - Apply all functionals from the dual space basis to the result of an evaluation at the points returned by PetscDualSpaceGetAllPoints()
693*20cf1dd8SToby Isaac 
694*20cf1dd8SToby Isaac   Input Parameters:
695*20cf1dd8SToby Isaac + sp        - The PetscDualSpace object
696*20cf1dd8SToby Isaac - pointEval - Evaluation at the points returned by PetscDualSpaceGetAllPoints()
697*20cf1dd8SToby Isaac 
698*20cf1dd8SToby Isaac   Output Parameter:
699*20cf1dd8SToby Isaac . spValue   - The values of all dual space functionals
700*20cf1dd8SToby Isaac 
701*20cf1dd8SToby Isaac   Level: developer
702*20cf1dd8SToby Isaac 
703*20cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate()
704*20cf1dd8SToby Isaac @*/
705*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceApplyAll(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue)
706*20cf1dd8SToby Isaac {
707*20cf1dd8SToby Isaac   PetscErrorCode ierr;
708*20cf1dd8SToby Isaac 
709*20cf1dd8SToby Isaac   PetscFunctionBegin;
710*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
711*20cf1dd8SToby Isaac   ierr = (*sp->ops->applyall)(sp, pointEval, spValue);CHKERRQ(ierr);
712*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
713*20cf1dd8SToby Isaac }
714*20cf1dd8SToby Isaac 
715*20cf1dd8SToby Isaac /*@C
716*20cf1dd8SToby Isaac   PetscDualSpaceApplyDefault - Apply a functional from the dual space basis to an input function by assuming a point evaluation functional.
717*20cf1dd8SToby Isaac 
718*20cf1dd8SToby Isaac   Input Parameters:
719*20cf1dd8SToby Isaac + sp    - The PetscDualSpace object
720*20cf1dd8SToby Isaac . f     - The basis functional index
721*20cf1dd8SToby Isaac . time  - The time
722*20cf1dd8SToby Isaac . cgeom - A context with geometric information for this cell, we use v0 (the initial vertex) and J (the Jacobian)
723*20cf1dd8SToby Isaac . Nc    - The number of components for the function
724*20cf1dd8SToby Isaac . func  - The input function
725*20cf1dd8SToby Isaac - ctx   - A context for the function
726*20cf1dd8SToby Isaac 
727*20cf1dd8SToby Isaac   Output Parameter:
728*20cf1dd8SToby Isaac . value   - The output value
729*20cf1dd8SToby Isaac 
730*20cf1dd8SToby Isaac   Note: The calling sequence for the callback func is given by:
731*20cf1dd8SToby Isaac 
732*20cf1dd8SToby Isaac $ func(PetscInt dim, PetscReal time, const PetscReal x[],
733*20cf1dd8SToby Isaac $      PetscInt numComponents, PetscScalar values[], void *ctx)
734*20cf1dd8SToby Isaac 
735*20cf1dd8SToby Isaac and the idea is to evaluate the functional as an integral
736*20cf1dd8SToby Isaac 
737*20cf1dd8SToby Isaac $ n(f) = int dx n(x) . f(x)
738*20cf1dd8SToby Isaac 
739*20cf1dd8SToby Isaac where both n and f have Nc components.
740*20cf1dd8SToby Isaac 
741*20cf1dd8SToby Isaac   Level: developer
742*20cf1dd8SToby Isaac 
743*20cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate()
744*20cf1dd8SToby Isaac @*/
745*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceApplyDefault(PetscDualSpace sp, PetscInt f, PetscReal time, PetscFEGeom *cgeom, PetscInt Nc, PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void *ctx, PetscScalar *value)
746*20cf1dd8SToby Isaac {
747*20cf1dd8SToby Isaac   DM               dm;
748*20cf1dd8SToby Isaac   PetscQuadrature  n;
749*20cf1dd8SToby Isaac   const PetscReal *points, *weights;
750*20cf1dd8SToby Isaac   PetscReal        x[3];
751*20cf1dd8SToby Isaac   PetscScalar     *val;
752*20cf1dd8SToby Isaac   PetscInt         dim, dE, qNc, c, Nq, q;
753*20cf1dd8SToby Isaac   PetscBool        isAffine;
754*20cf1dd8SToby Isaac   PetscErrorCode   ierr;
755*20cf1dd8SToby Isaac 
756*20cf1dd8SToby Isaac   PetscFunctionBegin;
757*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
758*20cf1dd8SToby Isaac   PetscValidPointer(value, 5);
759*20cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr);
760*20cf1dd8SToby Isaac   ierr = PetscDualSpaceGetFunctional(sp, f, &n);CHKERRQ(ierr);
761*20cf1dd8SToby Isaac   ierr = PetscQuadratureGetData(n, &dim, &qNc, &Nq, &points, &weights);CHKERRQ(ierr);
762*20cf1dd8SToby Isaac   if (dim != cgeom->dim) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature spatial dimension %D != cell geometry dimension %D", dim, cgeom->dim);
763*20cf1dd8SToby Isaac   if (qNc != Nc) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature components %D != function components %D", qNc, Nc);
764*20cf1dd8SToby Isaac   ierr = DMGetWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr);
765*20cf1dd8SToby Isaac   *value = 0.0;
766*20cf1dd8SToby Isaac   isAffine = cgeom->isAffine;
767*20cf1dd8SToby Isaac   dE = cgeom->dimEmbed;
768*20cf1dd8SToby Isaac   for (q = 0; q < Nq; ++q) {
769*20cf1dd8SToby Isaac     if (isAffine) {
770*20cf1dd8SToby Isaac       CoordinatesRefToReal(dE, cgeom->dim, cgeom->xi, cgeom->v, cgeom->J, &points[q*dim], x);
771*20cf1dd8SToby Isaac       ierr = (*func)(dE, time, x, Nc, val, ctx);CHKERRQ(ierr);
772*20cf1dd8SToby Isaac     } else {
773*20cf1dd8SToby Isaac       ierr = (*func)(dE, time, &cgeom->v[dE*q], Nc, val, ctx);CHKERRQ(ierr);
774*20cf1dd8SToby Isaac     }
775*20cf1dd8SToby Isaac     for (c = 0; c < Nc; ++c) {
776*20cf1dd8SToby Isaac       *value += val[c]*weights[q*Nc+c];
777*20cf1dd8SToby Isaac     }
778*20cf1dd8SToby Isaac   }
779*20cf1dd8SToby Isaac   ierr = DMRestoreWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr);
780*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
781*20cf1dd8SToby Isaac }
782*20cf1dd8SToby Isaac 
783*20cf1dd8SToby Isaac /*@C
784*20cf1dd8SToby Isaac   PetscDualSpaceApplyAllDefault - Apply all functionals from the dual space basis to the result of an evaluation at the points returned by PetscDualSpaceGetAllPoints()
785*20cf1dd8SToby Isaac 
786*20cf1dd8SToby Isaac   Input Parameters:
787*20cf1dd8SToby Isaac + sp        - The PetscDualSpace object
788*20cf1dd8SToby Isaac - pointEval - Evaluation at the points returned by PetscDualSpaceGetAllPoints()
789*20cf1dd8SToby Isaac 
790*20cf1dd8SToby Isaac   Output Parameter:
791*20cf1dd8SToby Isaac . spValue   - The values of all dual space functionals
792*20cf1dd8SToby Isaac 
793*20cf1dd8SToby Isaac   Level: developer
794*20cf1dd8SToby Isaac 
795*20cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate()
796*20cf1dd8SToby Isaac @*/
797*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceApplyAllDefault(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue)
798*20cf1dd8SToby Isaac {
799*20cf1dd8SToby Isaac   PetscQuadrature  n;
800*20cf1dd8SToby Isaac   const PetscReal *points, *weights;
801*20cf1dd8SToby Isaac   PetscInt         qNc, c, Nq, q, f, spdim, Nc;
802*20cf1dd8SToby Isaac   PetscInt         offset;
803*20cf1dd8SToby Isaac   PetscErrorCode   ierr;
804*20cf1dd8SToby Isaac 
805*20cf1dd8SToby Isaac   PetscFunctionBegin;
806*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
807*20cf1dd8SToby Isaac   PetscValidScalarPointer(pointEval, 2);
808*20cf1dd8SToby Isaac   PetscValidScalarPointer(spValue, 5);
809*20cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDimension(sp, &spdim);CHKERRQ(ierr);
810*20cf1dd8SToby Isaac   ierr = PetscDualSpaceGetNumComponents(sp, &Nc);CHKERRQ(ierr);
811*20cf1dd8SToby Isaac   for (f = 0, offset = 0; f < spdim; f++) {
812*20cf1dd8SToby Isaac     ierr = PetscDualSpaceGetFunctional(sp, f, &n);CHKERRQ(ierr);
813*20cf1dd8SToby Isaac     ierr = PetscQuadratureGetData(n, NULL, &qNc, &Nq, &points, &weights);CHKERRQ(ierr);
814*20cf1dd8SToby Isaac     if (qNc != Nc) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature components %D != function components %D", qNc, Nc);
815*20cf1dd8SToby Isaac     spValue[f] = 0.0;
816*20cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
817*20cf1dd8SToby Isaac       for (c = 0; c < Nc; ++c) {
818*20cf1dd8SToby Isaac         spValue[f] += pointEval[offset++]*weights[q*Nc+c];
819*20cf1dd8SToby Isaac       }
820*20cf1dd8SToby Isaac     }
821*20cf1dd8SToby Isaac   }
822*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
823*20cf1dd8SToby Isaac }
824*20cf1dd8SToby Isaac 
825*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetAllPoints(PetscDualSpace sp, PetscQuadrature *allPoints)
826*20cf1dd8SToby Isaac {
827*20cf1dd8SToby Isaac   PetscErrorCode ierr;
828*20cf1dd8SToby Isaac 
829*20cf1dd8SToby Isaac   PetscFunctionBegin;
830*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
831*20cf1dd8SToby Isaac   PetscValidPointer(allPoints,2);
832*20cf1dd8SToby Isaac   if (!sp->allPoints && sp->ops->createallpoints) {
833*20cf1dd8SToby Isaac     ierr = (*sp->ops->createallpoints)(sp,&sp->allPoints);CHKERRQ(ierr);
834*20cf1dd8SToby Isaac   }
835*20cf1dd8SToby Isaac   *allPoints = sp->allPoints;
836*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
837*20cf1dd8SToby Isaac }
838*20cf1dd8SToby Isaac 
839*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceCreateAllPointsDefault(PetscDualSpace sp, PetscQuadrature *allPoints)
840*20cf1dd8SToby Isaac {
841*20cf1dd8SToby Isaac   PetscInt        spdim;
842*20cf1dd8SToby Isaac   PetscInt        numPoints, offset;
843*20cf1dd8SToby Isaac   PetscReal       *points;
844*20cf1dd8SToby Isaac   PetscInt        f, dim;
845*20cf1dd8SToby Isaac   PetscQuadrature q;
846*20cf1dd8SToby Isaac   PetscErrorCode  ierr;
847*20cf1dd8SToby Isaac 
848*20cf1dd8SToby Isaac   PetscFunctionBegin;
849*20cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDimension(sp,&spdim);CHKERRQ(ierr);
850*20cf1dd8SToby Isaac   if (!spdim) {
851*20cf1dd8SToby Isaac     ierr = PetscQuadratureCreate(PETSC_COMM_SELF,allPoints);CHKERRQ(ierr);
852*20cf1dd8SToby Isaac     ierr = PetscQuadratureSetData(*allPoints,0,0,0,NULL,NULL);CHKERRQ(ierr);
853*20cf1dd8SToby Isaac   }
854*20cf1dd8SToby Isaac   ierr = PetscDualSpaceGetFunctional(sp,0,&q);CHKERRQ(ierr);
855*20cf1dd8SToby Isaac   ierr = PetscQuadratureGetData(q,&dim,NULL,&numPoints,NULL,NULL);CHKERRQ(ierr);
856*20cf1dd8SToby Isaac   for (f = 1; f < spdim; f++) {
857*20cf1dd8SToby Isaac     PetscInt Np;
858*20cf1dd8SToby Isaac 
859*20cf1dd8SToby Isaac     ierr = PetscDualSpaceGetFunctional(sp,f,&q);CHKERRQ(ierr);
860*20cf1dd8SToby Isaac     ierr = PetscQuadratureGetData(q,NULL,NULL,&Np,NULL,NULL);CHKERRQ(ierr);
861*20cf1dd8SToby Isaac     numPoints += Np;
862*20cf1dd8SToby Isaac   }
863*20cf1dd8SToby Isaac   ierr = PetscMalloc1(dim*numPoints,&points);CHKERRQ(ierr);
864*20cf1dd8SToby Isaac   for (f = 0, offset = 0; f < spdim; f++) {
865*20cf1dd8SToby Isaac     const PetscReal *p;
866*20cf1dd8SToby Isaac     PetscInt        Np, i;
867*20cf1dd8SToby Isaac 
868*20cf1dd8SToby Isaac     ierr = PetscDualSpaceGetFunctional(sp,f,&q);CHKERRQ(ierr);
869*20cf1dd8SToby Isaac     ierr = PetscQuadratureGetData(q,NULL,NULL,&Np,&p,NULL);CHKERRQ(ierr);
870*20cf1dd8SToby Isaac     for (i = 0; i < Np * dim; i++) {
871*20cf1dd8SToby Isaac       points[offset + i] = p[i];
872*20cf1dd8SToby Isaac     }
873*20cf1dd8SToby Isaac     offset += Np * dim;
874*20cf1dd8SToby Isaac   }
875*20cf1dd8SToby Isaac   ierr = PetscQuadratureCreate(PETSC_COMM_SELF,allPoints);CHKERRQ(ierr);
876*20cf1dd8SToby Isaac   ierr = PetscQuadratureSetData(*allPoints,dim,0,numPoints,points,NULL);CHKERRQ(ierr);
877*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
878*20cf1dd8SToby Isaac }
879*20cf1dd8SToby Isaac 
880*20cf1dd8SToby Isaac /*@C
881*20cf1dd8SToby Isaac   PetscDualSpaceApplyFVM - Apply a functional from the dual space basis to an input function by assuming a point evaluation functional at the cell centroid.
882*20cf1dd8SToby Isaac 
883*20cf1dd8SToby Isaac   Input Parameters:
884*20cf1dd8SToby Isaac + sp    - The PetscDualSpace object
885*20cf1dd8SToby Isaac . f     - The basis functional index
886*20cf1dd8SToby Isaac . time  - The time
887*20cf1dd8SToby Isaac . cgeom - A context with geometric information for this cell, we currently just use the centroid
888*20cf1dd8SToby Isaac . Nc    - The number of components for the function
889*20cf1dd8SToby Isaac . func  - The input function
890*20cf1dd8SToby Isaac - ctx   - A context for the function
891*20cf1dd8SToby Isaac 
892*20cf1dd8SToby Isaac   Output Parameter:
893*20cf1dd8SToby Isaac . value - The output value (scalar)
894*20cf1dd8SToby Isaac 
895*20cf1dd8SToby Isaac   Note: The calling sequence for the callback func is given by:
896*20cf1dd8SToby Isaac 
897*20cf1dd8SToby Isaac $ func(PetscInt dim, PetscReal time, const PetscReal x[],
898*20cf1dd8SToby Isaac $      PetscInt numComponents, PetscScalar values[], void *ctx)
899*20cf1dd8SToby Isaac 
900*20cf1dd8SToby Isaac and the idea is to evaluate the functional as an integral
901*20cf1dd8SToby Isaac 
902*20cf1dd8SToby Isaac $ n(f) = int dx n(x) . f(x)
903*20cf1dd8SToby Isaac 
904*20cf1dd8SToby Isaac where both n and f have Nc components.
905*20cf1dd8SToby Isaac 
906*20cf1dd8SToby Isaac   Level: developer
907*20cf1dd8SToby Isaac 
908*20cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate()
909*20cf1dd8SToby Isaac @*/
910*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceApplyFVM(PetscDualSpace sp, PetscInt f, PetscReal time, PetscFVCellGeom *cgeom, PetscInt Nc, PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void *ctx, PetscScalar *value)
911*20cf1dd8SToby Isaac {
912*20cf1dd8SToby Isaac   DM               dm;
913*20cf1dd8SToby Isaac   PetscQuadrature  n;
914*20cf1dd8SToby Isaac   const PetscReal *points, *weights;
915*20cf1dd8SToby Isaac   PetscScalar     *val;
916*20cf1dd8SToby Isaac   PetscInt         dimEmbed, qNc, c, Nq, q;
917*20cf1dd8SToby Isaac   PetscErrorCode   ierr;
918*20cf1dd8SToby Isaac 
919*20cf1dd8SToby Isaac   PetscFunctionBegin;
920*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
921*20cf1dd8SToby Isaac   PetscValidPointer(value, 5);
922*20cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr);
923*20cf1dd8SToby Isaac   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
924*20cf1dd8SToby Isaac   ierr = PetscDualSpaceGetFunctional(sp, f, &n);CHKERRQ(ierr);
925*20cf1dd8SToby Isaac   ierr = PetscQuadratureGetData(n, NULL, &qNc, &Nq, &points, &weights);CHKERRQ(ierr);
926*20cf1dd8SToby Isaac   if (qNc != Nc) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature components %D != function components %D", qNc, Nc);
927*20cf1dd8SToby Isaac   ierr = DMGetWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr);
928*20cf1dd8SToby Isaac   *value = 0.;
929*20cf1dd8SToby Isaac   for (q = 0; q < Nq; ++q) {
930*20cf1dd8SToby Isaac     ierr = (*func)(dimEmbed, time, cgeom->centroid, Nc, val, ctx);CHKERRQ(ierr);
931*20cf1dd8SToby Isaac     for (c = 0; c < Nc; ++c) {
932*20cf1dd8SToby Isaac       *value += val[c]*weights[q*Nc+c];
933*20cf1dd8SToby Isaac     }
934*20cf1dd8SToby Isaac   }
935*20cf1dd8SToby Isaac   ierr = DMRestoreWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr);
936*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
937*20cf1dd8SToby Isaac }
938*20cf1dd8SToby Isaac 
939*20cf1dd8SToby Isaac /*@
940*20cf1dd8SToby Isaac   PetscDualSpaceGetHeightSubspace - Get the subset of the dual space basis that is supported on a mesh point of a
941*20cf1dd8SToby Isaac   given height.  This assumes that the reference cell is symmetric over points of this height.
942*20cf1dd8SToby Isaac 
943*20cf1dd8SToby Isaac   If the dual space is not defined on mesh points of the given height (e.g. if the space is discontinuous and
944*20cf1dd8SToby Isaac   pointwise values are not defined on the element boundaries), or if the implementation of PetscDualSpace does not
945*20cf1dd8SToby Isaac   support extracting subspaces, then NULL is returned.
946*20cf1dd8SToby Isaac 
947*20cf1dd8SToby Isaac   This does not increment the reference count on the returned dual space, and the user should not destroy it.
948*20cf1dd8SToby Isaac 
949*20cf1dd8SToby Isaac   Not collective
950*20cf1dd8SToby Isaac 
951*20cf1dd8SToby Isaac   Input Parameters:
952*20cf1dd8SToby Isaac + sp - the PetscDualSpace object
953*20cf1dd8SToby Isaac - height - the height of the mesh point for which the subspace is desired
954*20cf1dd8SToby Isaac 
955*20cf1dd8SToby Isaac   Output Parameter:
956*20cf1dd8SToby Isaac . subsp - the subspace.  Note that the functionals in the subspace are with respect to the intrinsic geometry of the
957*20cf1dd8SToby Isaac   point, which will be of lesser dimension if height > 0.
958*20cf1dd8SToby Isaac 
959*20cf1dd8SToby Isaac   Level: advanced
960*20cf1dd8SToby Isaac 
961*20cf1dd8SToby Isaac .seealso: PetscSpaceGetHeightSubspace(), PetscDualSpace
962*20cf1dd8SToby Isaac @*/
963*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetHeightSubspace(PetscDualSpace sp, PetscInt height, PetscDualSpace *subsp)
964*20cf1dd8SToby Isaac {
965*20cf1dd8SToby Isaac   PetscErrorCode ierr;
966*20cf1dd8SToby Isaac 
967*20cf1dd8SToby Isaac   PetscFunctionBegin;
968*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
969*20cf1dd8SToby Isaac   PetscValidPointer(subsp, 3);
970*20cf1dd8SToby Isaac   *subsp = NULL;
971*20cf1dd8SToby Isaac   if (sp->ops->getheightsubspace) {
972*20cf1dd8SToby Isaac     ierr = (*sp->ops->getheightsubspace)(sp, height, subsp);CHKERRQ(ierr);
973*20cf1dd8SToby Isaac   }
974*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
975*20cf1dd8SToby Isaac }
976*20cf1dd8SToby Isaac 
977*20cf1dd8SToby Isaac /*@
978*20cf1dd8SToby Isaac   PetscDualSpaceGetPointSubspace - Get the subset of the dual space basis that is supported on a particular mesh point.
979*20cf1dd8SToby Isaac 
980*20cf1dd8SToby Isaac   If the dual space is not defined on the mesh point (e.g. if the space is discontinuous and pointwise values are not
981*20cf1dd8SToby Isaac   defined on the element boundaries), or if the implementation of PetscDualSpace does not support extracting
982*20cf1dd8SToby Isaac   subspaces, then NULL is returned.
983*20cf1dd8SToby Isaac 
984*20cf1dd8SToby Isaac   This does not increment the reference count on the returned dual space, and the user should not destroy it.
985*20cf1dd8SToby Isaac 
986*20cf1dd8SToby Isaac   Not collective
987*20cf1dd8SToby Isaac 
988*20cf1dd8SToby Isaac   Input Parameters:
989*20cf1dd8SToby Isaac + sp - the PetscDualSpace object
990*20cf1dd8SToby Isaac - point - the point (in the dual space's DM) for which the subspace is desired
991*20cf1dd8SToby Isaac 
992*20cf1dd8SToby Isaac   Output Parameters:
993*20cf1dd8SToby Isaac   bdsp - the subspace.  Note that the functionals in the subspace are with respect to the intrinsic geometry of the
994*20cf1dd8SToby Isaac   point, which will be of lesser dimension if height > 0.
995*20cf1dd8SToby Isaac 
996*20cf1dd8SToby Isaac   Level: advanced
997*20cf1dd8SToby Isaac 
998*20cf1dd8SToby Isaac .seealso: PetscDualSpace
999*20cf1dd8SToby Isaac @*/
1000*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetPointSubspace(PetscDualSpace sp, PetscInt point, PetscDualSpace *bdsp)
1001*20cf1dd8SToby Isaac {
1002*20cf1dd8SToby Isaac   PetscErrorCode ierr;
1003*20cf1dd8SToby Isaac 
1004*20cf1dd8SToby Isaac   PetscFunctionBegin;
1005*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1006*20cf1dd8SToby Isaac   PetscValidPointer(bdsp,2);
1007*20cf1dd8SToby Isaac   *bdsp = NULL;
1008*20cf1dd8SToby Isaac   if (sp->ops->getpointsubspace) {
1009*20cf1dd8SToby Isaac     ierr = (*sp->ops->getpointsubspace)(sp,point,bdsp);CHKERRQ(ierr);
1010*20cf1dd8SToby Isaac   } else if (sp->ops->getheightsubspace) {
1011*20cf1dd8SToby Isaac     DM       dm;
1012*20cf1dd8SToby Isaac     DMLabel  label;
1013*20cf1dd8SToby Isaac     PetscInt dim, depth, height;
1014*20cf1dd8SToby Isaac 
1015*20cf1dd8SToby Isaac     ierr = PetscDualSpaceGetDM(sp,&dm);CHKERRQ(ierr);
1016*20cf1dd8SToby Isaac     ierr = DMPlexGetDepth(dm,&dim);CHKERRQ(ierr);
1017*20cf1dd8SToby Isaac     ierr = DMPlexGetDepthLabel(dm,&label);CHKERRQ(ierr);
1018*20cf1dd8SToby Isaac     ierr = DMLabelGetValue(label,point,&depth);CHKERRQ(ierr);
1019*20cf1dd8SToby Isaac     height = dim - depth;
1020*20cf1dd8SToby Isaac     ierr = (*sp->ops->getheightsubspace)(sp,height,bdsp);CHKERRQ(ierr);
1021*20cf1dd8SToby Isaac   }
1022*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
1023*20cf1dd8SToby Isaac }
1024*20cf1dd8SToby Isaac 
1025