xref: /petsc/src/dm/dt/dualspace/interface/dualspace.c (revision 5bd1e5768a3b141382a3ebb2e780cd2d3b3bfbf0)
1 #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
2 #include <petscdmplex.h>
3 
4 PetscClassId PETSCDUALSPACE_CLASSID = 0;
5 
6 PetscFunctionList PetscDualSpaceList              = NULL;
7 PetscBool         PetscDualSpaceRegisterAllCalled = PETSC_FALSE;
8 
9 const char *const PetscDualSpaceReferenceCells[] = {"SIMPLEX", "TENSOR", "PetscDualSpaceReferenceCell", "PETSCDUALSPACE_REFCELL_",0};
10 
11 /*
12   PetscDualSpaceLatticePointLexicographic_Internal - Returns all tuples of size 'len' with nonnegative integers that sum up to at most 'max'.
13                                                      Ordering is lexicographic with lowest index as least significant in ordering.
14                                                      e.g. for len == 2 and max == 2, this will return, in order, {0,0}, {1,0}, {2,0}, {0,1}, {1,1}, {2,0}.
15 
16   Input Parameters:
17 + len - The length of the tuple
18 . max - The maximum sum
19 - tup - A tuple of length len+1: tup[len] > 0 indicates a stopping condition
20 
21   Output Parameter:
22 . tup - A tuple of len integers whos sum is at most 'max'
23 
24   Level: developer
25 
26 .seealso: PetscDualSpaceTensorPointLexicographic_Internal()
27 */
28 PetscErrorCode PetscDualSpaceLatticePointLexicographic_Internal(PetscInt len, PetscInt max, PetscInt tup[])
29 {
30   PetscFunctionBegin;
31   while (len--) {
32     max -= tup[len];
33     if (!max) {
34       tup[len] = 0;
35       break;
36     }
37   }
38   tup[++len]++;
39   PetscFunctionReturn(0);
40 }
41 
42 /*
43   PetscDualSpaceTensorPointLexicographic_Internal - Returns all tuples of size 'len' with nonnegative integers that are all less than or equal to 'max'.
44                                                     Ordering is lexicographic with lowest index as least significant in ordering.
45                                                     e.g. for len == 2 and max == 2, this will return, in order, {0,0}, {1,0}, {2,0}, {0,1}, {1,1}, {2,1}, {0,2}, {1,2}, {2,2}.
46 
47   Input Parameters:
48 + len - The length of the tuple
49 . max - The maximum value
50 - tup - A tuple of length len+1: tup[len] > 0 indicates a stopping condition
51 
52   Output Parameter:
53 . tup - A tuple of len integers whos sum is at most 'max'
54 
55   Level: developer
56 
57 .seealso: PetscDualSpaceLatticePointLexicographic_Internal()
58 */
59 PetscErrorCode PetscDualSpaceTensorPointLexicographic_Internal(PetscInt len, PetscInt max, PetscInt tup[])
60 {
61   PetscInt       i;
62 
63   PetscFunctionBegin;
64   for (i = 0; i < len; i++) {
65     if (tup[i] < max) {
66       break;
67     } else {
68       tup[i] = 0;
69     }
70   }
71   tup[i]++;
72   PetscFunctionReturn(0);
73 }
74 
75 /*@C
76   PetscDualSpaceRegister - Adds a new PetscDualSpace implementation
77 
78   Not Collective
79 
80   Input Parameters:
81 + name        - The name of a new user-defined creation routine
82 - create_func - The creation routine itself
83 
84   Notes:
85   PetscDualSpaceRegister() may be called multiple times to add several user-defined PetscDualSpaces
86 
87   Sample usage:
88 .vb
89     PetscDualSpaceRegister("my_space", MyPetscDualSpaceCreate);
90 .ve
91 
92   Then, your PetscDualSpace type can be chosen with the procedural interface via
93 .vb
94     PetscDualSpaceCreate(MPI_Comm, PetscDualSpace *);
95     PetscDualSpaceSetType(PetscDualSpace, "my_dual_space");
96 .ve
97    or at runtime via the option
98 .vb
99     -petscdualspace_type my_dual_space
100 .ve
101 
102   Level: advanced
103 
104 .seealso: PetscDualSpaceRegisterAll(), PetscDualSpaceRegisterDestroy()
105 
106 @*/
107 PetscErrorCode PetscDualSpaceRegister(const char sname[], PetscErrorCode (*function)(PetscDualSpace))
108 {
109   PetscErrorCode ierr;
110 
111   PetscFunctionBegin;
112   ierr = PetscFunctionListAdd(&PetscDualSpaceList, sname, function);CHKERRQ(ierr);
113   PetscFunctionReturn(0);
114 }
115 
116 /*@C
117   PetscDualSpaceSetType - Builds a particular PetscDualSpace
118 
119   Collective on sp
120 
121   Input Parameters:
122 + sp   - The PetscDualSpace object
123 - name - The kind of space
124 
125   Options Database Key:
126 . -petscdualspace_type <type> - Sets the PetscDualSpace type; use -help for a list of available types
127 
128   Level: intermediate
129 
130 .seealso: PetscDualSpaceGetType(), PetscDualSpaceCreate()
131 @*/
132 PetscErrorCode PetscDualSpaceSetType(PetscDualSpace sp, PetscDualSpaceType name)
133 {
134   PetscErrorCode (*r)(PetscDualSpace);
135   PetscBool      match;
136   PetscErrorCode ierr;
137 
138   PetscFunctionBegin;
139   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
140   ierr = PetscObjectTypeCompare((PetscObject) sp, name, &match);CHKERRQ(ierr);
141   if (match) PetscFunctionReturn(0);
142 
143   if (!PetscDualSpaceRegisterAllCalled) {ierr = PetscDualSpaceRegisterAll();CHKERRQ(ierr);}
144   ierr = PetscFunctionListFind(PetscDualSpaceList, name, &r);CHKERRQ(ierr);
145   if (!r) SETERRQ1(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscDualSpace type: %s", name);
146 
147   if (sp->ops->destroy) {
148     ierr             = (*sp->ops->destroy)(sp);CHKERRQ(ierr);
149     sp->ops->destroy = NULL;
150   }
151   ierr = (*r)(sp);CHKERRQ(ierr);
152   ierr = PetscObjectChangeTypeName((PetscObject) sp, name);CHKERRQ(ierr);
153   PetscFunctionReturn(0);
154 }
155 
156 /*@C
157   PetscDualSpaceGetType - Gets the PetscDualSpace type name (as a string) from the object.
158 
159   Not Collective
160 
161   Input Parameter:
162 . sp  - The PetscDualSpace
163 
164   Output Parameter:
165 . name - The PetscDualSpace type name
166 
167   Level: intermediate
168 
169 .seealso: PetscDualSpaceSetType(), PetscDualSpaceCreate()
170 @*/
171 PetscErrorCode PetscDualSpaceGetType(PetscDualSpace sp, PetscDualSpaceType *name)
172 {
173   PetscErrorCode ierr;
174 
175   PetscFunctionBegin;
176   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
177   PetscValidPointer(name, 2);
178   if (!PetscDualSpaceRegisterAllCalled) {
179     ierr = PetscDualSpaceRegisterAll();CHKERRQ(ierr);
180   }
181   *name = ((PetscObject) sp)->type_name;
182   PetscFunctionReturn(0);
183 }
184 
185 static PetscErrorCode PetscDualSpaceView_ASCII(PetscDualSpace sp, PetscViewer v)
186 {
187   PetscViewerFormat format;
188   PetscInt          pdim, f;
189   PetscErrorCode    ierr;
190 
191   PetscFunctionBegin;
192   ierr = PetscDualSpaceGetDimension(sp, &pdim);CHKERRQ(ierr);
193   ierr = PetscObjectPrintClassNamePrefixType((PetscObject) sp, v);CHKERRQ(ierr);
194   ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr);
195   ierr = PetscViewerASCIIPrintf(v, "Dual space with %D components, size %D\n", sp->Nc, pdim);CHKERRQ(ierr);
196   if (sp->ops->view) {ierr = (*sp->ops->view)(sp, v);CHKERRQ(ierr);}
197   ierr = PetscViewerGetFormat(v, &format);CHKERRQ(ierr);
198   if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
199     ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr);
200     for (f = 0; f < pdim; ++f) {
201       ierr = PetscViewerASCIIPrintf(v, "Dual basis vector %D\n", f);CHKERRQ(ierr);
202       ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr);
203       ierr = PetscQuadratureView(sp->functional[f], v);CHKERRQ(ierr);
204       ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr);
205     }
206     ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr);
207   }
208   ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr);
209   PetscFunctionReturn(0);
210 }
211 
212 /*@
213   PetscDualSpaceView - Views a PetscDualSpace
214 
215   Collective on sp
216 
217   Input Parameter:
218 + sp - the PetscDualSpace object to view
219 - v  - the viewer
220 
221   Level: beginner
222 
223 .seealso PetscDualSpaceDestroy()
224 @*/
225 PetscErrorCode PetscDualSpaceView(PetscDualSpace sp, PetscViewer v)
226 {
227   PetscBool      iascii;
228   PetscErrorCode ierr;
229 
230   PetscFunctionBegin;
231   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
232   if (v) PetscValidHeaderSpecific(v, PETSC_VIEWER_CLASSID, 2);
233   if (!v) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) sp), &v);CHKERRQ(ierr);}
234   ierr = PetscObjectTypeCompare((PetscObject) v, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
235   if (iascii) {ierr = PetscDualSpaceView_ASCII(sp, v);CHKERRQ(ierr);}
236   PetscFunctionReturn(0);
237 }
238 
239 /*@
240   PetscDualSpaceSetFromOptions - sets parameters in a PetscDualSpace from the options database
241 
242   Collective on sp
243 
244   Input Parameter:
245 . sp - the PetscDualSpace object to set options for
246 
247   Options Database:
248 . -petscspace_degree the approximation order of the space
249 
250   Level: intermediate
251 
252 .seealso PetscDualSpaceView()
253 @*/
254 PetscErrorCode PetscDualSpaceSetFromOptions(PetscDualSpace sp)
255 {
256   PetscDualSpaceReferenceCell refCell = PETSCDUALSPACE_REFCELL_SIMPLEX;
257   PetscInt                    refDim  = 0;
258   PetscBool                   flg;
259   const char                 *defaultType;
260   char                        name[256];
261   PetscErrorCode              ierr;
262 
263   PetscFunctionBegin;
264   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
265   if (!((PetscObject) sp)->type_name) {
266     defaultType = PETSCDUALSPACELAGRANGE;
267   } else {
268     defaultType = ((PetscObject) sp)->type_name;
269   }
270   if (!PetscSpaceRegisterAllCalled) {ierr = PetscSpaceRegisterAll();CHKERRQ(ierr);}
271 
272   ierr = PetscObjectOptionsBegin((PetscObject) sp);CHKERRQ(ierr);
273   ierr = PetscOptionsFList("-petscdualspace_type", "Dual space", "PetscDualSpaceSetType", PetscDualSpaceList, defaultType, name, 256, &flg);CHKERRQ(ierr);
274   if (flg) {
275     ierr = PetscDualSpaceSetType(sp, name);CHKERRQ(ierr);
276   } else if (!((PetscObject) sp)->type_name) {
277     ierr = PetscDualSpaceSetType(sp, defaultType);CHKERRQ(ierr);
278   }
279   ierr = PetscOptionsBoundedInt("-petscdualspace_degree", "The approximation order", "PetscDualSpaceSetOrder", sp->order, &sp->order, NULL,0);CHKERRQ(ierr);
280   ierr = PetscOptionsBoundedInt("-petscdualspace_components", "The number of components", "PetscDualSpaceSetNumComponents", sp->Nc, &sp->Nc, NULL,1);CHKERRQ(ierr);
281   if (sp->ops->setfromoptions) {
282     ierr = (*sp->ops->setfromoptions)(PetscOptionsObject,sp);CHKERRQ(ierr);
283   }
284   ierr = PetscOptionsBoundedInt("-petscdualspace_refdim", "The spatial dimension of the reference cell", "PetscDualSpaceSetReferenceCell", refDim, &refDim, NULL,0);CHKERRQ(ierr);
285   ierr = PetscOptionsEnum("-petscdualspace_refcell", "Reference cell", "PetscDualSpaceSetReferenceCell", PetscDualSpaceReferenceCells, (PetscEnum) refCell, (PetscEnum *) &refCell, &flg);CHKERRQ(ierr);
286   if (flg) {
287     DM K;
288 
289     if (!refDim) SETERRQ(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_INCOMP, "Reference cell specified without a dimension. Use -petscdualspace_refdim.");
290     ierr = PetscDualSpaceCreateReferenceCell(sp, refDim, refCell == PETSCDUALSPACE_REFCELL_SIMPLEX ? PETSC_TRUE : PETSC_FALSE, &K);CHKERRQ(ierr);
291     ierr = PetscDualSpaceSetDM(sp, K);CHKERRQ(ierr);
292     ierr = DMDestroy(&K);CHKERRQ(ierr);
293   }
294 
295   /* process any options handlers added with PetscObjectAddOptionsHandler() */
296   ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) sp);CHKERRQ(ierr);
297   ierr = PetscOptionsEnd();CHKERRQ(ierr);
298   sp->setfromoptionscalled = PETSC_TRUE;
299   PetscFunctionReturn(0);
300 }
301 
302 /*@
303   PetscDualSpaceSetUp - Construct a basis for the PetscDualSpace
304 
305   Collective on sp
306 
307   Input Parameter:
308 . sp - the PetscDualSpace object to setup
309 
310   Level: intermediate
311 
312 .seealso PetscDualSpaceView(), PetscDualSpaceDestroy()
313 @*/
314 PetscErrorCode PetscDualSpaceSetUp(PetscDualSpace sp)
315 {
316   PetscErrorCode ierr;
317 
318   PetscFunctionBegin;
319   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
320   if (sp->setupcalled) PetscFunctionReturn(0);
321   sp->setupcalled = PETSC_TRUE;
322   if (sp->ops->setup) {ierr = (*sp->ops->setup)(sp);CHKERRQ(ierr);}
323   if (sp->setfromoptionscalled) {ierr = PetscDualSpaceViewFromOptions(sp, NULL, "-petscdualspace_view");CHKERRQ(ierr);}
324   PetscFunctionReturn(0);
325 }
326 
327 /*@
328   PetscDualSpaceDestroy - Destroys a PetscDualSpace object
329 
330   Collective on sp
331 
332   Input Parameter:
333 . sp - the PetscDualSpace object to destroy
334 
335   Level: beginner
336 
337 .seealso PetscDualSpaceView()
338 @*/
339 PetscErrorCode PetscDualSpaceDestroy(PetscDualSpace *sp)
340 {
341   PetscInt       dim, f;
342   PetscErrorCode ierr;
343 
344   PetscFunctionBegin;
345   if (!*sp) PetscFunctionReturn(0);
346   PetscValidHeaderSpecific((*sp), PETSCDUALSPACE_CLASSID, 1);
347 
348   if (--((PetscObject)(*sp))->refct > 0) {*sp = 0; PetscFunctionReturn(0);}
349   ((PetscObject) (*sp))->refct = 0;
350 
351   ierr = PetscDualSpaceGetDimension(*sp, &dim);CHKERRQ(ierr);
352   for (f = 0; f < dim; ++f) {
353     ierr = PetscQuadratureDestroy(&(*sp)->functional[f]);CHKERRQ(ierr);
354   }
355   ierr = PetscFree((*sp)->functional);CHKERRQ(ierr);
356   ierr = PetscQuadratureDestroy(&(*sp)->allPoints);CHKERRQ(ierr);
357   ierr = DMDestroy(&(*sp)->dm);CHKERRQ(ierr);
358 
359   if ((*sp)->ops->destroy) {ierr = (*(*sp)->ops->destroy)(*sp);CHKERRQ(ierr);}
360   ierr = PetscHeaderDestroy(sp);CHKERRQ(ierr);
361   PetscFunctionReturn(0);
362 }
363 
364 /*@
365   PetscDualSpaceCreate - Creates an empty PetscDualSpace object. The type can then be set with PetscDualSpaceSetType().
366 
367   Collective
368 
369   Input Parameter:
370 . comm - The communicator for the PetscDualSpace object
371 
372   Output Parameter:
373 . sp - The PetscDualSpace object
374 
375   Level: beginner
376 
377 .seealso: PetscDualSpaceSetType(), PETSCDUALSPACELAGRANGE
378 @*/
379 PetscErrorCode PetscDualSpaceCreate(MPI_Comm comm, PetscDualSpace *sp)
380 {
381   PetscDualSpace s;
382   PetscErrorCode ierr;
383 
384   PetscFunctionBegin;
385   PetscValidPointer(sp, 2);
386   ierr = PetscCitationsRegister(FECitation,&FEcite);CHKERRQ(ierr);
387   *sp  = NULL;
388   ierr = PetscFEInitializePackage();CHKERRQ(ierr);
389 
390   ierr = PetscHeaderCreate(s, PETSCDUALSPACE_CLASSID, "PetscDualSpace", "Dual Space", "PetscDualSpace", comm, PetscDualSpaceDestroy, PetscDualSpaceView);CHKERRQ(ierr);
391 
392   s->order = 0;
393   s->Nc    = 1;
394   s->k     = 0;
395   s->setupcalled = PETSC_FALSE;
396 
397   *sp = s;
398   PetscFunctionReturn(0);
399 }
400 
401 /*@
402   PetscDualSpaceDuplicate - Creates a duplicate PetscDualSpace object, however it is not setup.
403 
404   Collective on sp
405 
406   Input Parameter:
407 . sp - The original PetscDualSpace
408 
409   Output Parameter:
410 . spNew - The duplicate PetscDualSpace
411 
412   Level: beginner
413 
414 .seealso: PetscDualSpaceCreate(), PetscDualSpaceSetType()
415 @*/
416 PetscErrorCode PetscDualSpaceDuplicate(PetscDualSpace sp, PetscDualSpace *spNew)
417 {
418   PetscErrorCode ierr;
419 
420   PetscFunctionBegin;
421   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
422   PetscValidPointer(spNew, 2);
423   ierr = (*sp->ops->duplicate)(sp, spNew);CHKERRQ(ierr);
424   PetscFunctionReturn(0);
425 }
426 
427 /*@
428   PetscDualSpaceGetDM - Get the DM representing the reference cell
429 
430   Not collective
431 
432   Input Parameter:
433 . sp - The PetscDualSpace
434 
435   Output Parameter:
436 . dm - The reference cell
437 
438   Level: intermediate
439 
440 .seealso: PetscDualSpaceSetDM(), PetscDualSpaceCreate()
441 @*/
442 PetscErrorCode PetscDualSpaceGetDM(PetscDualSpace sp, DM *dm)
443 {
444   PetscFunctionBegin;
445   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
446   PetscValidPointer(dm, 2);
447   *dm = sp->dm;
448   PetscFunctionReturn(0);
449 }
450 
451 /*@
452   PetscDualSpaceSetDM - Get the DM representing the reference cell
453 
454   Not collective
455 
456   Input Parameters:
457 + sp - The PetscDualSpace
458 - dm - The reference cell
459 
460   Level: intermediate
461 
462 .seealso: PetscDualSpaceGetDM(), PetscDualSpaceCreate()
463 @*/
464 PetscErrorCode PetscDualSpaceSetDM(PetscDualSpace sp, DM dm)
465 {
466   PetscErrorCode ierr;
467 
468   PetscFunctionBegin;
469   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
470   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
471   ierr = DMDestroy(&sp->dm);CHKERRQ(ierr);
472   ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);
473   sp->dm = dm;
474   PetscFunctionReturn(0);
475 }
476 
477 /*@
478   PetscDualSpaceGetOrder - Get the order of the dual space
479 
480   Not collective
481 
482   Input Parameter:
483 . sp - The PetscDualSpace
484 
485   Output Parameter:
486 . order - The order
487 
488   Level: intermediate
489 
490 .seealso: PetscDualSpaceSetOrder(), PetscDualSpaceCreate()
491 @*/
492 PetscErrorCode PetscDualSpaceGetOrder(PetscDualSpace sp, PetscInt *order)
493 {
494   PetscFunctionBegin;
495   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
496   PetscValidPointer(order, 2);
497   *order = sp->order;
498   PetscFunctionReturn(0);
499 }
500 
501 /*@
502   PetscDualSpaceSetOrder - Set the order of the dual space
503 
504   Not collective
505 
506   Input Parameters:
507 + sp - The PetscDualSpace
508 - order - The order
509 
510   Level: intermediate
511 
512 .seealso: PetscDualSpaceGetOrder(), PetscDualSpaceCreate()
513 @*/
514 PetscErrorCode PetscDualSpaceSetOrder(PetscDualSpace sp, PetscInt order)
515 {
516   PetscFunctionBegin;
517   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
518   sp->order = order;
519   PetscFunctionReturn(0);
520 }
521 
522 /*@
523   PetscDualSpaceGetNumComponents - Return the number of components for this space
524 
525   Input Parameter:
526 . sp - The PetscDualSpace
527 
528   Output Parameter:
529 . Nc - The number of components
530 
531   Note: A vector space, for example, will have d components, where d is the spatial dimension
532 
533   Level: intermediate
534 
535 .seealso: PetscDualSpaceSetNumComponents(), PetscDualSpaceGetDimension(), PetscDualSpaceCreate(), PetscDualSpace
536 @*/
537 PetscErrorCode PetscDualSpaceGetNumComponents(PetscDualSpace sp, PetscInt *Nc)
538 {
539   PetscFunctionBegin;
540   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
541   PetscValidPointer(Nc, 2);
542   *Nc = sp->Nc;
543   PetscFunctionReturn(0);
544 }
545 
546 /*@
547   PetscDualSpaceSetNumComponents - Set the number of components for this space
548 
549   Input Parameters:
550 + sp - The PetscDualSpace
551 - order - The number of components
552 
553   Level: intermediate
554 
555 .seealso: PetscDualSpaceGetNumComponents(), PetscDualSpaceCreate(), PetscDualSpace
556 @*/
557 PetscErrorCode PetscDualSpaceSetNumComponents(PetscDualSpace sp, PetscInt Nc)
558 {
559   PetscFunctionBegin;
560   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
561   sp->Nc = Nc;
562   PetscFunctionReturn(0);
563 }
564 
565 /*@
566   PetscDualSpaceGetFunctional - Get the i-th basis functional in the dual space
567 
568   Not collective
569 
570   Input Parameters:
571 + sp - The PetscDualSpace
572 - i  - The basis number
573 
574   Output Parameter:
575 . functional - The basis functional
576 
577   Level: intermediate
578 
579 .seealso: PetscDualSpaceGetDimension(), PetscDualSpaceCreate()
580 @*/
581 PetscErrorCode PetscDualSpaceGetFunctional(PetscDualSpace sp, PetscInt i, PetscQuadrature *functional)
582 {
583   PetscInt       dim;
584   PetscErrorCode ierr;
585 
586   PetscFunctionBegin;
587   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
588   PetscValidPointer(functional, 3);
589   ierr = PetscDualSpaceGetDimension(sp, &dim);CHKERRQ(ierr);
590   if ((i < 0) || (i >= dim)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Functional index %d must be in [0, %d)", i, dim);
591   *functional = sp->functional[i];
592   PetscFunctionReturn(0);
593 }
594 
595 /*@
596   PetscDualSpaceGetDimension - Get the dimension of the dual space, i.e. the number of basis functionals
597 
598   Not collective
599 
600   Input Parameter:
601 . sp - The PetscDualSpace
602 
603   Output Parameter:
604 . dim - The dimension
605 
606   Level: intermediate
607 
608 .seealso: PetscDualSpaceGetFunctional(), PetscDualSpaceCreate()
609 @*/
610 PetscErrorCode PetscDualSpaceGetDimension(PetscDualSpace sp, PetscInt *dim)
611 {
612   PetscErrorCode ierr;
613 
614   PetscFunctionBegin;
615   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
616   PetscValidPointer(dim, 2);
617   *dim = 0;
618   if (sp->ops->getdimension) {ierr = (*sp->ops->getdimension)(sp, dim);CHKERRQ(ierr);}
619   PetscFunctionReturn(0);
620 }
621 
622 /*@C
623   PetscDualSpaceGetNumDof - Get the number of degrees of freedom for each spatial (topological) dimension
624 
625   Not collective
626 
627   Input Parameter:
628 . sp - The PetscDualSpace
629 
630   Output Parameter:
631 . numDof - An array of length dim+1 which holds the number of dofs for each dimension
632 
633   Level: intermediate
634 
635 .seealso: PetscDualSpaceGetFunctional(), PetscDualSpaceCreate()
636 @*/
637 PetscErrorCode PetscDualSpaceGetNumDof(PetscDualSpace sp, const PetscInt **numDof)
638 {
639   PetscErrorCode ierr;
640 
641   PetscFunctionBegin;
642   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
643   PetscValidPointer(numDof, 2);
644   ierr = (*sp->ops->getnumdof)(sp, numDof);CHKERRQ(ierr);
645   if (!*numDof) SETERRQ(PetscObjectComm((PetscObject) sp), PETSC_ERR_LIB, "Empty numDof[] returned from dual space implementation");
646   PetscFunctionReturn(0);
647 }
648 
649 /*@
650   PetscDualSpaceCreateSection - Create a PetscSection over the reference cell with the layout from this space
651 
652   Collective on sp
653 
654   Input Parameters:
655 + sp      - The PetscDualSpace
656 
657   Output Parameter:
658 . section - The section
659 
660   Level: advanced
661 
662 .seealso: PetscDualSpaceCreate(), DMPLEX
663 @*/
664 PetscErrorCode PetscDualSpaceCreateSection(PetscDualSpace sp, PetscSection *section)
665 {
666   DM             dm;
667   PetscInt       pStart, pEnd, depth, h, offset;
668   const PetscInt *numDof;
669   PetscErrorCode ierr;
670 
671   PetscFunctionBegin;
672   ierr = PetscDualSpaceGetDM(sp,&dm);CHKERRQ(ierr);
673   ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr);
674   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)sp),section);CHKERRQ(ierr);
675   ierr = PetscSectionSetChart(*section,pStart,pEnd);CHKERRQ(ierr);
676   ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr);
677   ierr = PetscDualSpaceGetNumDof(sp,&numDof);CHKERRQ(ierr);
678   for (h = 0; h <= depth; h++) {
679     PetscInt hStart, hEnd, p, dof;
680 
681     ierr = DMPlexGetHeightStratum(dm,h,&hStart,&hEnd);CHKERRQ(ierr);
682     dof = numDof[depth - h];
683     for (p = hStart; p < hEnd; p++) {
684       ierr = PetscSectionSetDof(*section,p,dof);CHKERRQ(ierr);
685     }
686   }
687   ierr = PetscSectionSetUp(*section);CHKERRQ(ierr);
688   for (h = 0, offset = 0; h <= depth; h++) {
689     PetscInt hStart, hEnd, p, dof;
690 
691     ierr = DMPlexGetHeightStratum(dm,h,&hStart,&hEnd);CHKERRQ(ierr);
692     dof = numDof[depth - h];
693     for (p = hStart; p < hEnd; p++) {
694       ierr = PetscSectionGetDof(*section,p,&dof);CHKERRQ(ierr);
695       ierr = PetscSectionSetOffset(*section,p,offset);CHKERRQ(ierr);
696       offset += dof;
697     }
698   }
699   PetscFunctionReturn(0);
700 }
701 
702 /*@
703   PetscDualSpaceCreateReferenceCell - Create a DMPLEX with the appropriate FEM reference cell
704 
705   Collective on sp
706 
707   Input Parameters:
708 + sp      - The PetscDualSpace
709 . dim     - The spatial dimension
710 - simplex - Flag for simplex, otherwise use a tensor-product cell
711 
712   Output Parameter:
713 . refdm - The reference cell
714 
715   Level: intermediate
716 
717 .seealso: PetscDualSpaceCreate(), DMPLEX
718 @*/
719 PetscErrorCode PetscDualSpaceCreateReferenceCell(PetscDualSpace sp, PetscInt dim, PetscBool simplex, DM *refdm)
720 {
721   PetscErrorCode ierr;
722 
723   PetscFunctionBeginUser;
724   ierr = DMPlexCreateReferenceCell(PetscObjectComm((PetscObject) sp), dim, simplex, refdm);CHKERRQ(ierr);
725   PetscFunctionReturn(0);
726 }
727 
728 /*@C
729   PetscDualSpaceApply - Apply a functional from the dual space basis to an input function
730 
731   Input Parameters:
732 + sp      - The PetscDualSpace object
733 . f       - The basis functional index
734 . time    - The time
735 . 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)
736 . numComp - The number of components for the function
737 . func    - The input function
738 - ctx     - A context for the function
739 
740   Output Parameter:
741 . value   - numComp output values
742 
743   Note: The calling sequence for the callback func is given by:
744 
745 $ func(PetscInt dim, PetscReal time, const PetscReal x[],
746 $      PetscInt numComponents, PetscScalar values[], void *ctx)
747 
748   Level: beginner
749 
750 .seealso: PetscDualSpaceCreate()
751 @*/
752 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)
753 {
754   PetscErrorCode ierr;
755 
756   PetscFunctionBegin;
757   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
758   PetscValidPointer(cgeom, 4);
759   PetscValidPointer(value, 8);
760   ierr = (*sp->ops->apply)(sp, f, time, cgeom, numComp, func, ctx, value);CHKERRQ(ierr);
761   PetscFunctionReturn(0);
762 }
763 
764 /*@C
765   PetscDualSpaceApplyAll - Apply all functionals from the dual space basis to the result of an evaluation at the points returned by PetscDualSpaceGetAllPoints()
766 
767   Input Parameters:
768 + sp        - The PetscDualSpace object
769 - pointEval - Evaluation at the points returned by PetscDualSpaceGetAllPoints()
770 
771   Output Parameter:
772 . spValue   - The values of all dual space functionals
773 
774   Level: beginner
775 
776 .seealso: PetscDualSpaceCreate()
777 @*/
778 PetscErrorCode PetscDualSpaceApplyAll(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue)
779 {
780   PetscErrorCode ierr;
781 
782   PetscFunctionBegin;
783   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
784   ierr = (*sp->ops->applyall)(sp, pointEval, spValue);CHKERRQ(ierr);
785   PetscFunctionReturn(0);
786 }
787 
788 /*@C
789   PetscDualSpaceApplyDefault - Apply a functional from the dual space basis to an input function by assuming a point evaluation functional.
790 
791   Input Parameters:
792 + sp    - The PetscDualSpace object
793 . f     - The basis functional index
794 . time  - The time
795 . cgeom - A context with geometric information for this cell, we use v0 (the initial vertex) and J (the Jacobian)
796 . Nc    - The number of components for the function
797 . func  - The input function
798 - ctx   - A context for the function
799 
800   Output Parameter:
801 . value   - The output value
802 
803   Note: The calling sequence for the callback func is given by:
804 
805 $ func(PetscInt dim, PetscReal time, const PetscReal x[],
806 $      PetscInt numComponents, PetscScalar values[], void *ctx)
807 
808 and the idea is to evaluate the functional as an integral
809 
810 $ n(f) = int dx n(x) . f(x)
811 
812 where both n and f have Nc components.
813 
814   Level: beginner
815 
816 .seealso: PetscDualSpaceCreate()
817 @*/
818 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)
819 {
820   DM               dm;
821   PetscQuadrature  n;
822   const PetscReal *points, *weights;
823   PetscReal        x[3];
824   PetscScalar     *val;
825   PetscInt         dim, dE, qNc, c, Nq, q;
826   PetscBool        isAffine;
827   PetscErrorCode   ierr;
828 
829   PetscFunctionBegin;
830   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
831   PetscValidPointer(value, 5);
832   ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr);
833   ierr = PetscDualSpaceGetFunctional(sp, f, &n);CHKERRQ(ierr);
834   ierr = PetscQuadratureGetData(n, &dim, &qNc, &Nq, &points, &weights);CHKERRQ(ierr);
835   if (dim != cgeom->dim) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature spatial dimension %D != cell geometry dimension %D", dim, cgeom->dim);
836   if (qNc != Nc) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature components %D != function components %D", qNc, Nc);
837   ierr = DMGetWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr);
838   *value = 0.0;
839   isAffine = cgeom->isAffine;
840   dE = cgeom->dimEmbed;
841   for (q = 0; q < Nq; ++q) {
842     if (isAffine) {
843       CoordinatesRefToReal(dE, cgeom->dim, cgeom->xi, cgeom->v, cgeom->J, &points[q*dim], x);
844       ierr = (*func)(dE, time, x, Nc, val, ctx);CHKERRQ(ierr);
845     } else {
846       ierr = (*func)(dE, time, &cgeom->v[dE*q], Nc, val, ctx);CHKERRQ(ierr);
847     }
848     for (c = 0; c < Nc; ++c) {
849       *value += val[c]*weights[q*Nc+c];
850     }
851   }
852   ierr = DMRestoreWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr);
853   PetscFunctionReturn(0);
854 }
855 
856 /*@C
857   PetscDualSpaceApplyAllDefault - Apply all functionals from the dual space basis to the result of an evaluation at the points returned by PetscDualSpaceGetAllPoints()
858 
859   Input Parameters:
860 + sp        - The PetscDualSpace object
861 - pointEval - Evaluation at the points returned by PetscDualSpaceGetAllPoints()
862 
863   Output Parameter:
864 . spValue   - The values of all dual space functionals
865 
866   Level: beginner
867 
868 .seealso: PetscDualSpaceCreate()
869 @*/
870 PetscErrorCode PetscDualSpaceApplyAllDefault(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue)
871 {
872   PetscQuadrature  n;
873   const PetscReal *points, *weights;
874   PetscInt         qNc, c, Nq, q, f, spdim, Nc;
875   PetscInt         offset;
876   PetscErrorCode   ierr;
877 
878   PetscFunctionBegin;
879   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
880   PetscValidScalarPointer(pointEval, 2);
881   PetscValidScalarPointer(spValue, 5);
882   ierr = PetscDualSpaceGetDimension(sp, &spdim);CHKERRQ(ierr);
883   ierr = PetscDualSpaceGetNumComponents(sp, &Nc);CHKERRQ(ierr);
884   for (f = 0, offset = 0; f < spdim; f++) {
885     ierr = PetscDualSpaceGetFunctional(sp, f, &n);CHKERRQ(ierr);
886     ierr = PetscQuadratureGetData(n, NULL, &qNc, &Nq, &points, &weights);CHKERRQ(ierr);
887     if (qNc != Nc) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature components %D != function components %D", qNc, Nc);
888     spValue[f] = 0.0;
889     for (q = 0; q < Nq; ++q) {
890       for (c = 0; c < Nc; ++c) {
891         spValue[f] += pointEval[offset++]*weights[q*Nc+c];
892       }
893     }
894   }
895   PetscFunctionReturn(0);
896 }
897 
898 /*@
899   PetscDualSpaceGetAllPoints - Get all quadrature points from this space
900 
901   Input Parameter:
902 . sp - The dualspace
903 
904   Output Parameter:
905 . allPoints - A PetscQuadrature object containing all evaluation points
906 
907   Level: advanced
908 
909 .seealso: PetscDualSpaceCreate()
910 @*/
911 PetscErrorCode PetscDualSpaceGetAllPoints(PetscDualSpace sp, PetscQuadrature *allPoints)
912 {
913   PetscErrorCode ierr;
914 
915   PetscFunctionBegin;
916   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
917   PetscValidPointer(allPoints,2);
918   if (!sp->allPoints && sp->ops->createallpoints) {
919     ierr = (*sp->ops->createallpoints)(sp,&sp->allPoints);CHKERRQ(ierr);
920   }
921   *allPoints = sp->allPoints;
922   PetscFunctionReturn(0);
923 }
924 
925 /*@
926   PetscDualSpaceCreateAllPointsDefault - Create all evaluation points by examining functionals
927 
928   Input Parameter:
929 . sp - The dualspace
930 
931   Output Parameter:
932 . allPoints - A PetscQuadrature object containing all evaluation points
933 
934   Level: advanced
935 
936 .seealso: PetscDualSpaceCreate()
937 @*/
938 PetscErrorCode PetscDualSpaceCreateAllPointsDefault(PetscDualSpace sp, PetscQuadrature *allPoints)
939 {
940   PetscInt        spdim;
941   PetscInt        numPoints, offset;
942   PetscReal       *points;
943   PetscInt        f, dim;
944   PetscQuadrature q;
945   PetscErrorCode  ierr;
946 
947   PetscFunctionBegin;
948   ierr = PetscDualSpaceGetDimension(sp,&spdim);CHKERRQ(ierr);
949   if (!spdim) {
950     ierr = PetscQuadratureCreate(PETSC_COMM_SELF,allPoints);CHKERRQ(ierr);
951     ierr = PetscQuadratureSetData(*allPoints,0,0,0,NULL,NULL);CHKERRQ(ierr);
952   }
953   ierr = PetscDualSpaceGetFunctional(sp,0,&q);CHKERRQ(ierr);
954   ierr = PetscQuadratureGetData(q,&dim,NULL,&numPoints,NULL,NULL);CHKERRQ(ierr);
955   for (f = 1; f < spdim; f++) {
956     PetscInt Np;
957 
958     ierr = PetscDualSpaceGetFunctional(sp,f,&q);CHKERRQ(ierr);
959     ierr = PetscQuadratureGetData(q,NULL,NULL,&Np,NULL,NULL);CHKERRQ(ierr);
960     numPoints += Np;
961   }
962   ierr = PetscMalloc1(dim*numPoints,&points);CHKERRQ(ierr);
963   for (f = 0, offset = 0; f < spdim; f++) {
964     const PetscReal *p;
965     PetscInt        Np, i;
966 
967     ierr = PetscDualSpaceGetFunctional(sp,f,&q);CHKERRQ(ierr);
968     ierr = PetscQuadratureGetData(q,NULL,NULL,&Np,&p,NULL);CHKERRQ(ierr);
969     for (i = 0; i < Np * dim; i++) {
970       points[offset + i] = p[i];
971     }
972     offset += Np * dim;
973   }
974   ierr = PetscQuadratureCreate(PETSC_COMM_SELF,allPoints);CHKERRQ(ierr);
975   ierr = PetscQuadratureSetData(*allPoints,dim,0,numPoints,points,NULL);CHKERRQ(ierr);
976   PetscFunctionReturn(0);
977 }
978 
979 /*@C
980   PetscDualSpaceApplyFVM - Apply a functional from the dual space basis to an input function by assuming a point evaluation functional at the cell centroid.
981 
982   Input Parameters:
983 + sp    - The PetscDualSpace object
984 . f     - The basis functional index
985 . time  - The time
986 . cgeom - A context with geometric information for this cell, we currently just use the centroid
987 . Nc    - The number of components for the function
988 . func  - The input function
989 - ctx   - A context for the function
990 
991   Output Parameter:
992 . value - The output value (scalar)
993 
994   Note: The calling sequence for the callback func is given by:
995 
996 $ func(PetscInt dim, PetscReal time, const PetscReal x[],
997 $      PetscInt numComponents, PetscScalar values[], void *ctx)
998 
999 and the idea is to evaluate the functional as an integral
1000 
1001 $ n(f) = int dx n(x) . f(x)
1002 
1003 where both n and f have Nc components.
1004 
1005   Level: beginner
1006 
1007 .seealso: PetscDualSpaceCreate()
1008 @*/
1009 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)
1010 {
1011   DM               dm;
1012   PetscQuadrature  n;
1013   const PetscReal *points, *weights;
1014   PetscScalar     *val;
1015   PetscInt         dimEmbed, qNc, c, Nq, q;
1016   PetscErrorCode   ierr;
1017 
1018   PetscFunctionBegin;
1019   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1020   PetscValidPointer(value, 5);
1021   ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr);
1022   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
1023   ierr = PetscDualSpaceGetFunctional(sp, f, &n);CHKERRQ(ierr);
1024   ierr = PetscQuadratureGetData(n, NULL, &qNc, &Nq, &points, &weights);CHKERRQ(ierr);
1025   if (qNc != Nc) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature components %D != function components %D", qNc, Nc);
1026   ierr = DMGetWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr);
1027   *value = 0.;
1028   for (q = 0; q < Nq; ++q) {
1029     ierr = (*func)(dimEmbed, time, cgeom->centroid, Nc, val, ctx);CHKERRQ(ierr);
1030     for (c = 0; c < Nc; ++c) {
1031       *value += val[c]*weights[q*Nc+c];
1032     }
1033   }
1034   ierr = DMRestoreWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr);
1035   PetscFunctionReturn(0);
1036 }
1037 
1038 /*@
1039   PetscDualSpaceGetHeightSubspace - Get the subset of the dual space basis that is supported on a mesh point of a
1040   given height.  This assumes that the reference cell is symmetric over points of this height.
1041 
1042   If the dual space is not defined on mesh points of the given height (e.g. if the space is discontinuous and
1043   pointwise values are not defined on the element boundaries), or if the implementation of PetscDualSpace does not
1044   support extracting subspaces, then NULL is returned.
1045 
1046   This does not increment the reference count on the returned dual space, and the user should not destroy it.
1047 
1048   Not collective
1049 
1050   Input Parameters:
1051 + sp - the PetscDualSpace object
1052 - height - the height of the mesh point for which the subspace is desired
1053 
1054   Output Parameter:
1055 . subsp - the subspace.  Note that the functionals in the subspace are with respect to the intrinsic geometry of the
1056   point, which will be of lesser dimension if height > 0.
1057 
1058   Level: advanced
1059 
1060 .seealso: PetscSpaceGetHeightSubspace(), PetscDualSpace
1061 @*/
1062 PetscErrorCode PetscDualSpaceGetHeightSubspace(PetscDualSpace sp, PetscInt height, PetscDualSpace *subsp)
1063 {
1064   PetscErrorCode ierr;
1065 
1066   PetscFunctionBegin;
1067   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1068   PetscValidPointer(subsp, 3);
1069   *subsp = NULL;
1070   if (sp->ops->getheightsubspace) {ierr = (*sp->ops->getheightsubspace)(sp, height, subsp);CHKERRQ(ierr);}
1071   PetscFunctionReturn(0);
1072 }
1073 
1074 /*@
1075   PetscDualSpaceGetPointSubspace - Get the subset of the dual space basis that is supported on a particular mesh point.
1076 
1077   If the dual space is not defined on the mesh point (e.g. if the space is discontinuous and pointwise values are not
1078   defined on the element boundaries), or if the implementation of PetscDualSpace does not support extracting
1079   subspaces, then NULL is returned.
1080 
1081   This does not increment the reference count on the returned dual space, and the user should not destroy it.
1082 
1083   Not collective
1084 
1085   Input Parameters:
1086 + sp - the PetscDualSpace object
1087 - point - the point (in the dual space's DM) for which the subspace is desired
1088 
1089   Output Parameters:
1090   bdsp - the subspace.  Note that the functionals in the subspace are with respect to the intrinsic geometry of the
1091   point, which will be of lesser dimension if height > 0.
1092 
1093   Level: advanced
1094 
1095 .seealso: PetscDualSpace
1096 @*/
1097 PetscErrorCode PetscDualSpaceGetPointSubspace(PetscDualSpace sp, PetscInt point, PetscDualSpace *bdsp)
1098 {
1099   PetscErrorCode ierr;
1100 
1101   PetscFunctionBegin;
1102   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1103   PetscValidPointer(bdsp,2);
1104   *bdsp = NULL;
1105   if (sp->ops->getpointsubspace) {
1106     ierr = (*sp->ops->getpointsubspace)(sp,point,bdsp);CHKERRQ(ierr);
1107   } else if (sp->ops->getheightsubspace) {
1108     DM       dm;
1109     DMLabel  label;
1110     PetscInt dim, depth, height;
1111 
1112     ierr = PetscDualSpaceGetDM(sp,&dm);CHKERRQ(ierr);
1113     ierr = DMPlexGetDepth(dm,&dim);CHKERRQ(ierr);
1114     ierr = DMPlexGetDepthLabel(dm,&label);CHKERRQ(ierr);
1115     ierr = DMLabelGetValue(label,point,&depth);CHKERRQ(ierr);
1116     height = dim - depth;
1117     ierr = (*sp->ops->getheightsubspace)(sp,height,bdsp);CHKERRQ(ierr);
1118   }
1119   PetscFunctionReturn(0);
1120 }
1121 
1122 /*@C
1123   PetscDualSpaceGetSymmetries - Returns a description of the symmetries of this basis
1124 
1125   Not collective
1126 
1127   Input Parameter:
1128 . sp - the PetscDualSpace object
1129 
1130   Output Parameters:
1131 + perms - Permutations of the local degrees of freedom, parameterized by the point orientation
1132 - flips - Sign reversal of the local degrees of freedom, parameterized by the point orientation
1133 
1134   Note: The permutation and flip arrays are organized in the following way
1135 $ perms[p][ornt][dof # on point] = new local dof #
1136 $ flips[p][ornt][dof # on point] = reversal or not
1137 
1138   Level: developer
1139 
1140 .seealso: PetscDualSpaceSetSymmetries()
1141 @*/
1142 PetscErrorCode PetscDualSpaceGetSymmetries(PetscDualSpace sp, const PetscInt ****perms, const PetscScalar ****flips)
1143 {
1144   PetscErrorCode ierr;
1145 
1146   PetscFunctionBegin;
1147   PetscValidHeaderSpecific(sp,PETSCDUALSPACE_CLASSID,1);
1148   if (perms) {PetscValidPointer(perms,2); *perms = NULL;}
1149   if (flips) {PetscValidPointer(flips,3); *flips = NULL;}
1150   if (sp->ops->getsymmetries) {ierr = (sp->ops->getsymmetries)(sp,perms,flips);CHKERRQ(ierr);}
1151   PetscFunctionReturn(0);
1152 }
1153 
1154 /*@
1155   PetscDualSpaceGetDeRahm - Get the k-simplex associated with the functionals in this dual space
1156 
1157   Input Parameter:
1158 . dsp - The PetscDualSpace
1159 
1160   Output Parameter:
1161 . k   - The simplex dimension
1162 
1163   Level: developer
1164 
1165   Note: Currently supported values are
1166 $ 0: These are H_1 methods that only transform coordinates
1167 $ 1: These are Hcurl methods that transform functions using the covariant Piola transform (COVARIANT_PIOLA_TRANSFORM)
1168 $ 2: These are the same as 1
1169 $ 3: These are Hdiv methods that transform functions using the contravariant Piola transform (CONTRAVARIANT_PIOLA_TRANSFORM)
1170 
1171 .seealso: PetscDualSpacePullback(), PetscDualSpacePushforward(), PetscDualSpaceTransform(), PetscDualSpaceTransformType
1172 @*/
1173 PetscErrorCode PetscDualSpaceGetDeRahm(PetscDualSpace dsp, PetscInt *k)
1174 {
1175   PetscFunctionBeginHot;
1176   PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1);
1177   PetscValidPointer(k, 2);
1178   *k = dsp->k;
1179   PetscFunctionReturn(0);
1180 }
1181 
1182 /*@C
1183   PetscDualSpaceTransform - Transform the function values
1184 
1185   Input Parameters:
1186 + dsp       - The PetscDualSpace
1187 . trans     - The type of transform
1188 . isInverse - Flag to invert the transform
1189 . fegeom    - The cell geometry
1190 . Nv        - The number of function samples
1191 . Nc        - The number of function components
1192 - vals      - The function values
1193 
1194   Output Parameter:
1195 . vals      - The transformed function values
1196 
1197   Level: intermediate
1198 
1199 .seealso: PetscDualSpaceTransformGradient(), PetscDualSpacePullback(), PetscDualSpacePushforward(), PetscDualSpaceTransformType
1200 @*/
1201 PetscErrorCode PetscDualSpaceTransform(PetscDualSpace dsp, PetscDualSpaceTransformType trans, PetscBool isInverse, PetscFEGeom *fegeom, PetscInt Nv, PetscInt Nc, PetscScalar vals[])
1202 {
1203   PetscInt dim, v, c;
1204 
1205   PetscFunctionBeginHot;
1206   PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1);
1207   PetscValidPointer(fegeom, 4);
1208   PetscValidPointer(vals, 7);
1209   dim = dsp->dm->dim;
1210   /* Assume its a vector, otherwise assume its a bunch of scalars */
1211   if (Nc == 1 || Nc != dim) PetscFunctionReturn(0);
1212   switch (trans) {
1213     case IDENTITY_TRANSFORM: break;
1214     case COVARIANT_PIOLA_TRANSFORM: /* Covariant Piola mapping $\sigma^*(F) = J^{-T} F \circ \phi^{-1)$ */
1215     if (isInverse) {
1216       for (v = 0; v < Nv; ++v) {
1217         switch (dim)
1218         {
1219           case 2: DMPlex_MultTranspose2DReal_Internal(fegeom->J, 1, &vals[v*Nc], &vals[v*Nc]);break;
1220           case 3: DMPlex_MultTranspose3DReal_Internal(fegeom->J, 1, &vals[v*Nc], &vals[v*Nc]);break;
1221           default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim);
1222         }
1223       }
1224     } else {
1225       for (v = 0; v < Nv; ++v) {
1226         switch (dim)
1227         {
1228           case 2: DMPlex_MultTranspose2DReal_Internal(fegeom->invJ, 1, &vals[v*Nc], &vals[v*Nc]);break;
1229           case 3: DMPlex_MultTranspose3DReal_Internal(fegeom->invJ, 1, &vals[v*Nc], &vals[v*Nc]);break;
1230           default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim);
1231         }
1232       }
1233     }
1234     break;
1235     case CONTRAVARIANT_PIOLA_TRANSFORM: /* Contravariant Piola mapping $\sigma^*(F) = \frac{1}{|\det J|} J F \circ \phi^{-1}$ */
1236     if (isInverse) {
1237       for (v = 0; v < Nv; ++v) {
1238         switch (dim)
1239         {
1240           case 2: DMPlex_Mult2DReal_Internal(fegeom->invJ, 1, &vals[v*Nc], &vals[v*Nc]);break;
1241           case 3: DMPlex_Mult3DReal_Internal(fegeom->invJ, 1, &vals[v*Nc], &vals[v*Nc]);break;
1242           default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim);
1243         }
1244         for (c = 0; c < Nc; ++c) vals[v*Nc+c] *= fegeom->detJ[0];
1245       }
1246     } else {
1247       for (v = 0; v < Nv; ++v) {
1248         switch (dim)
1249         {
1250           case 2: DMPlex_Mult2DReal_Internal(fegeom->J, 1, &vals[v*Nc], &vals[v*Nc]);break;
1251           case 3: DMPlex_Mult3DReal_Internal(fegeom->J, 1, &vals[v*Nc], &vals[v*Nc]);break;
1252           default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim);
1253         }
1254         for (c = 0; c < Nc; ++c) vals[v*Nc+c] /= fegeom->detJ[0];
1255       }
1256     }
1257     break;
1258   }
1259   PetscFunctionReturn(0);
1260 }
1261 /*@C
1262   PetscDualSpaceTransformGradient - Transform the function gradient values
1263 
1264   Input Parameters:
1265 + dsp       - The PetscDualSpace
1266 . trans     - The type of transform
1267 . isInverse - Flag to invert the transform
1268 . fegeom    - The cell geometry
1269 . Nv        - The number of function gradient samples
1270 . Nc        - The number of function components
1271 - vals      - The function gradient values
1272 
1273   Output Parameter:
1274 . vals      - The transformed function values
1275 
1276   Level: intermediate
1277 
1278 .seealso: PetscDualSpaceTransform(), PetscDualSpacePullback(), PetscDualSpacePushforward(), PetscDualSpaceTransformType
1279 @*/
1280 PetscErrorCode PetscDualSpaceTransformGradient(PetscDualSpace dsp, PetscDualSpaceTransformType trans, PetscBool isInverse, PetscFEGeom *fegeom, PetscInt Nv, PetscInt Nc, PetscScalar vals[])
1281 {
1282   PetscInt dim, v, c, d;
1283 
1284   PetscFunctionBeginHot;
1285   PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1);
1286   PetscValidPointer(fegeom, 4);
1287   PetscValidPointer(vals, 7);
1288   dim = dsp->dm->dim;
1289   /* Transform gradient */
1290   for (v = 0; v < Nv; ++v) {
1291     for (c = 0; c < Nc; ++c) {
1292       switch (dim)
1293       {
1294         case 1: vals[(v*Nc+c)*dim] *= fegeom->invJ[0];
1295         case 2: DMPlex_MultTranspose2DReal_Internal(fegeom->invJ, 1, &vals[(v*Nc+c)*dim], &vals[(v*Nc+c)*dim]);break;
1296         case 3: DMPlex_MultTranspose3DReal_Internal(fegeom->invJ, 1, &vals[(v*Nc+c)*dim], &vals[(v*Nc+c)*dim]);break;
1297         default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim);
1298       }
1299     }
1300   }
1301   /* Assume its a vector, otherwise assume its a bunch of scalars */
1302   if (Nc == 1 || Nc != dim) PetscFunctionReturn(0);
1303   switch (trans) {
1304     case IDENTITY_TRANSFORM: break;
1305     case COVARIANT_PIOLA_TRANSFORM: /* Covariant Piola mapping $\sigma^*(F) = J^{-T} F \circ \phi^{-1)$ */
1306     if (isInverse) {
1307       for (v = 0; v < Nv; ++v) {
1308         for (d = 0; d < dim; ++d) {
1309           switch (dim)
1310           {
1311             case 2: DMPlex_MultTranspose2DReal_Internal(fegeom->J, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
1312             case 3: DMPlex_MultTranspose3DReal_Internal(fegeom->J, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
1313             default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim);
1314           }
1315         }
1316       }
1317     } else {
1318       for (v = 0; v < Nv; ++v) {
1319         for (d = 0; d < dim; ++d) {
1320           switch (dim)
1321           {
1322             case 2: DMPlex_MultTranspose2DReal_Internal(fegeom->invJ, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
1323             case 3: DMPlex_MultTranspose3DReal_Internal(fegeom->invJ, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
1324             default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim);
1325           }
1326         }
1327       }
1328     }
1329     break;
1330     case CONTRAVARIANT_PIOLA_TRANSFORM: /* Contravariant Piola mapping $\sigma^*(F) = \frac{1}{|\det J|} J F \circ \phi^{-1}$ */
1331     if (isInverse) {
1332       for (v = 0; v < Nv; ++v) {
1333         for (d = 0; d < dim; ++d) {
1334           switch (dim)
1335           {
1336             case 2: DMPlex_Mult2DReal_Internal(fegeom->invJ, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
1337             case 3: DMPlex_Mult3DReal_Internal(fegeom->invJ, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
1338             default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim);
1339           }
1340           for (c = 0; c < Nc; ++c) vals[(v*Nc+c)*dim+d] *= fegeom->detJ[0];
1341         }
1342       }
1343     } else {
1344       for (v = 0; v < Nv; ++v) {
1345         for (d = 0; d < dim; ++d) {
1346           switch (dim)
1347           {
1348             case 2: DMPlex_Mult2DReal_Internal(fegeom->J, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
1349             case 3: DMPlex_Mult3DReal_Internal(fegeom->J, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
1350             default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim);
1351           }
1352           for (c = 0; c < Nc; ++c) vals[(v*Nc+c)*dim+d] /= fegeom->detJ[0];
1353         }
1354       }
1355     }
1356     break;
1357   }
1358   PetscFunctionReturn(0);
1359 }
1360 
1361 /*@C
1362   PetscDualSpacePullback - Transform the given functional so that it operates on real space, rather than the reference element. Operationally, this means that we map the function evaluations depending on continuity requirements of our finite element method.
1363 
1364   Input Parameters:
1365 + dsp        - The PetscDualSpace
1366 . fegeom     - The geometry for this cell
1367 . Nq         - The number of function samples
1368 . Nc         - The number of function components
1369 - pointEval  - The function values
1370 
1371   Output Parameter:
1372 . pointEval  - The transformed function values
1373 
1374   Level: advanced
1375 
1376   Note: Functions transform in a complementary way (pushforward) to functionals, so that the scalar product is invariant. The type of transform is dependent on the associated k-simplex from the DeRahm complex.
1377 
1378 .seealso: PetscDualSpacePushforward(), PetscDualSpaceTransform(), PetscDualSpaceGetDeRahm()
1379 @*/
1380 PetscErrorCode PetscDualSpacePullback(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, PetscScalar pointEval[])
1381 {
1382   PetscDualSpaceTransformType trans;
1383   PetscErrorCode              ierr;
1384 
1385   PetscFunctionBeginHot;
1386   PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1);
1387   PetscValidPointer(fegeom, 2);
1388   PetscValidPointer(pointEval, 5);
1389   /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k.
1390      This determines their transformation properties. */
1391   switch (dsp->k)
1392   {
1393     case 0: /* H^1 point evaluations */
1394     trans = IDENTITY_TRANSFORM;break;
1395     case 1: /* Hcurl preserves tangential edge traces  */
1396     case 2:
1397     trans = COVARIANT_PIOLA_TRANSFORM;break;
1398     case 3: /* Hdiv preserve normal traces */
1399     trans = CONTRAVARIANT_PIOLA_TRANSFORM;break;
1400     default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %D for transformation", dsp->k);
1401   }
1402   ierr = PetscDualSpaceTransform(dsp, trans, PETSC_TRUE, fegeom, Nq, Nc, pointEval);CHKERRQ(ierr);
1403   PetscFunctionReturn(0);
1404 }
1405 
1406 /*@C
1407   PetscDualSpacePushforward - Transform the given function so that it operates on real space, rather than the reference element. Operationally, this means that we map the function evaluations depending on continuity requirements of our finite element method.
1408 
1409   Input Parameters:
1410 + dsp        - The PetscDualSpace
1411 . fegeom     - The geometry for this cell
1412 . Nq         - The number of function samples
1413 . Nc         - The number of function components
1414 - pointEval  - The function values
1415 
1416   Output Parameter:
1417 . pointEval  - The transformed function values
1418 
1419   Level: advanced
1420 
1421   Note: Functionals transform in a complementary way (pullback) to functions, so that the scalar product is invariant. The type of transform is dependent on the associated k-simplex from the DeRahm complex.
1422 
1423 .seealso: PetscDualSpacePullback(), PetscDualSpaceTransform(), PetscDualSpaceGetDeRahm()
1424 @*/
1425 PetscErrorCode PetscDualSpacePushforward(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, PetscScalar pointEval[])
1426 {
1427   PetscDualSpaceTransformType trans;
1428   PetscErrorCode              ierr;
1429 
1430   PetscFunctionBeginHot;
1431   PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1);
1432   PetscValidPointer(fegeom, 2);
1433   PetscValidPointer(pointEval, 5);
1434   /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k.
1435      This determines their transformation properties. */
1436   switch (dsp->k)
1437   {
1438     case 0: /* H^1 point evaluations */
1439     trans = IDENTITY_TRANSFORM;break;
1440     case 1: /* Hcurl preserves tangential edge traces  */
1441     case 2:
1442     trans = COVARIANT_PIOLA_TRANSFORM;break;
1443     case 3: /* Hdiv preserve normal traces */
1444     trans = CONTRAVARIANT_PIOLA_TRANSFORM;break;
1445     default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %D for transformation", dsp->k);
1446   }
1447   ierr = PetscDualSpaceTransform(dsp, trans, PETSC_FALSE, fegeom, Nq, Nc, pointEval);CHKERRQ(ierr);
1448   PetscFunctionReturn(0);
1449 }
1450 
1451 /*@C
1452   PetscDualSpacePushforwardGradient - Transform the given function gradient so that it operates on real space, rather than the reference element. Operationally, this means that we map the function evaluations depending on continuity requirements of our finite element method.
1453 
1454   Input Parameters:
1455 + dsp        - The PetscDualSpace
1456 . fegeom     - The geometry for this cell
1457 . Nq         - The number of function gradient samples
1458 . Nc         - The number of function components
1459 - pointEval  - The function gradient values
1460 
1461   Output Parameter:
1462 . pointEval  - The transformed function gradient values
1463 
1464   Level: advanced
1465 
1466   Note: Functionals transform in a complementary way (pullback) to functions, so that the scalar product is invariant. The type of transform is dependent on the associated k-simplex from the DeRahm complex.
1467 
1468 .seealso: PetscDualSpacePushforward(), PPetscDualSpacePullback(), PetscDualSpaceTransform(), PetscDualSpaceGetDeRahm()
1469 @*/
1470 PetscErrorCode PetscDualSpacePushforwardGradient(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, PetscScalar pointEval[])
1471 {
1472   PetscDualSpaceTransformType trans;
1473   PetscErrorCode              ierr;
1474 
1475   PetscFunctionBeginHot;
1476   PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1);
1477   PetscValidPointer(fegeom, 2);
1478   PetscValidPointer(pointEval, 5);
1479   /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k.
1480      This determines their transformation properties. */
1481   switch (dsp->k)
1482   {
1483     case 0: /* H^1 point evaluations */
1484     trans = IDENTITY_TRANSFORM;break;
1485     case 1: /* Hcurl preserves tangential edge traces  */
1486     case 2:
1487     trans = COVARIANT_PIOLA_TRANSFORM;break;
1488     case 3: /* Hdiv preserve normal traces */
1489     trans = CONTRAVARIANT_PIOLA_TRANSFORM;break;
1490     default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %D for transformation", dsp->k);
1491   }
1492   ierr = PetscDualSpaceTransformGradient(dsp, trans, PETSC_FALSE, fegeom, Nq, Nc, pointEval);CHKERRQ(ierr);
1493   PetscFunctionReturn(0);
1494 }
1495