xref: /petsc/src/dm/dt/dualspace/interface/dualspace.c (revision ebead697dbf761eb322f829370bbe90b3bd93fa3)
1 #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
2 #include <petscdmplex.h>
3 
4 PetscClassId PETSCDUALSPACE_CLASSID = 0;
5 
6 PetscLogEvent PETSCDUALSPACE_SetUp;
7 
8 PetscFunctionList PetscDualSpaceList              = NULL;
9 PetscBool         PetscDualSpaceRegisterAllCalled = PETSC_FALSE;
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}, {0,2}.
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   PetscFunctionBegin;
110   PetscCall(PetscFunctionListAdd(&PetscDualSpaceList, sname, function));
111   PetscFunctionReturn(0);
112 }
113 
114 /*@C
115   PetscDualSpaceSetType - Builds a particular PetscDualSpace
116 
117   Collective on sp
118 
119   Input Parameters:
120 + sp   - The PetscDualSpace object
121 - name - The kind of space
122 
123   Options Database Key:
124 . -petscdualspace_type <type> - Sets the PetscDualSpace type; use -help for a list of available types
125 
126   Level: intermediate
127 
128 .seealso: `PetscDualSpaceGetType()`, `PetscDualSpaceCreate()`
129 @*/
130 PetscErrorCode PetscDualSpaceSetType(PetscDualSpace sp, PetscDualSpaceType name)
131 {
132   PetscErrorCode (*r)(PetscDualSpace);
133   PetscBool      match;
134 
135   PetscFunctionBegin;
136   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
137   PetscCall(PetscObjectTypeCompare((PetscObject) sp, name, &match));
138   if (match) PetscFunctionReturn(0);
139 
140   if (!PetscDualSpaceRegisterAllCalled) PetscCall(PetscDualSpaceRegisterAll());
141   PetscCall(PetscFunctionListFind(PetscDualSpaceList, name, &r));
142   PetscCheck(r,PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscDualSpace type: %s", name);
143 
144   if (sp->ops->destroy) {
145     PetscCall((*sp->ops->destroy)(sp));
146     sp->ops->destroy = NULL;
147   }
148   PetscCall((*r)(sp));
149   PetscCall(PetscObjectChangeTypeName((PetscObject) sp, name));
150   PetscFunctionReturn(0);
151 }
152 
153 /*@C
154   PetscDualSpaceGetType - Gets the PetscDualSpace type name (as a string) from the object.
155 
156   Not Collective
157 
158   Input Parameter:
159 . sp  - The PetscDualSpace
160 
161   Output Parameter:
162 . name - The PetscDualSpace type name
163 
164   Level: intermediate
165 
166 .seealso: `PetscDualSpaceSetType()`, `PetscDualSpaceCreate()`
167 @*/
168 PetscErrorCode PetscDualSpaceGetType(PetscDualSpace sp, PetscDualSpaceType *name)
169 {
170   PetscFunctionBegin;
171   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
172   PetscValidPointer(name, 2);
173   if (!PetscDualSpaceRegisterAllCalled) {
174     PetscCall(PetscDualSpaceRegisterAll());
175   }
176   *name = ((PetscObject) sp)->type_name;
177   PetscFunctionReturn(0);
178 }
179 
180 static PetscErrorCode PetscDualSpaceView_ASCII(PetscDualSpace sp, PetscViewer v)
181 {
182   PetscViewerFormat format;
183   PetscInt          pdim, f;
184 
185   PetscFunctionBegin;
186   PetscCall(PetscDualSpaceGetDimension(sp, &pdim));
187   PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject) sp, v));
188   PetscCall(PetscViewerASCIIPushTab(v));
189   if (sp->k) {
190     PetscCall(PetscViewerASCIIPrintf(v, "Dual space for %" PetscInt_FMT "-forms %swith %" PetscInt_FMT " components, size %" PetscInt_FMT "\n", PetscAbsInt(sp->k), sp->k < 0 ? "(stored in dual form) ": "", sp->Nc, pdim));
191   } else {
192     PetscCall(PetscViewerASCIIPrintf(v, "Dual space with %" PetscInt_FMT " components, size %" PetscInt_FMT "\n", sp->Nc, pdim));
193   }
194   if (sp->ops->view) PetscCall((*sp->ops->view)(sp, v));
195   PetscCall(PetscViewerGetFormat(v, &format));
196   if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
197     PetscCall(PetscViewerASCIIPushTab(v));
198     for (f = 0; f < pdim; ++f) {
199       PetscCall(PetscViewerASCIIPrintf(v, "Dual basis vector %" PetscInt_FMT "\n", f));
200       PetscCall(PetscViewerASCIIPushTab(v));
201       PetscCall(PetscQuadratureView(sp->functional[f], v));
202       PetscCall(PetscViewerASCIIPopTab(v));
203     }
204     PetscCall(PetscViewerASCIIPopTab(v));
205   }
206   PetscCall(PetscViewerASCIIPopTab(v));
207   PetscFunctionReturn(0);
208 }
209 
210 /*@C
211    PetscDualSpaceViewFromOptions - View from Options
212 
213    Collective on PetscDualSpace
214 
215    Input Parameters:
216 +  A - the PetscDualSpace object
217 .  obj - Optional object, proivides prefix
218 -  name - command line option
219 
220    Level: intermediate
221 .seealso: `PetscDualSpace`, `PetscDualSpaceView()`, `PetscObjectViewFromOptions()`, `PetscDualSpaceCreate()`
222 @*/
223 PetscErrorCode  PetscDualSpaceViewFromOptions(PetscDualSpace A,PetscObject obj,const char name[])
224 {
225   PetscFunctionBegin;
226   PetscValidHeaderSpecific(A,PETSCDUALSPACE_CLASSID,1);
227   PetscCall(PetscObjectViewFromOptions((PetscObject)A,obj,name));
228   PetscFunctionReturn(0);
229 }
230 
231 /*@
232   PetscDualSpaceView - Views a PetscDualSpace
233 
234   Collective on sp
235 
236   Input Parameters:
237 + sp - the PetscDualSpace object to view
238 - v  - the viewer
239 
240   Level: beginner
241 
242 .seealso `PetscDualSpaceDestroy()`, `PetscDualSpace`
243 @*/
244 PetscErrorCode PetscDualSpaceView(PetscDualSpace sp, PetscViewer v)
245 {
246   PetscBool      iascii;
247 
248   PetscFunctionBegin;
249   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
250   if (v) PetscValidHeaderSpecific(v, PETSC_VIEWER_CLASSID, 2);
251   if (!v) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) sp), &v));
252   PetscCall(PetscObjectTypeCompare((PetscObject) v, PETSCVIEWERASCII, &iascii));
253   if (iascii) PetscCall(PetscDualSpaceView_ASCII(sp, v));
254   PetscFunctionReturn(0);
255 }
256 
257 /*@
258   PetscDualSpaceSetFromOptions - sets parameters in a PetscDualSpace from the options database
259 
260   Collective on sp
261 
262   Input Parameter:
263 . sp - the PetscDualSpace object to set options for
264 
265   Options Database:
266 + -petscdualspace_order <order>      - the approximation order of the space
267 . -petscdualspace_form_degree <deg>  - the form degree, say 0 for point evaluations, or 2 for area integrals
268 . -petscdualspace_components <c>     - the number of components, say d for a vector field
269 - -petscdualspace_refcell <celltype> - Reference cell type name
270 
271   Level: intermediate
272 
273 .seealso `PetscDualSpaceView()`, `PetscDualSpace`, `PetscObjectSetFromOptions()`
274 @*/
275 PetscErrorCode PetscDualSpaceSetFromOptions(PetscDualSpace sp)
276 {
277   DMPolytopeType refCell = DM_POLYTOPE_TRIANGLE;
278   const char    *defaultType;
279   char           name[256];
280   PetscBool      flg;
281 
282   PetscFunctionBegin;
283   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
284   if (!((PetscObject) sp)->type_name) {
285     defaultType = PETSCDUALSPACELAGRANGE;
286   } else {
287     defaultType = ((PetscObject) sp)->type_name;
288   }
289   if (!PetscSpaceRegisterAllCalled) PetscCall(PetscSpaceRegisterAll());
290 
291   PetscObjectOptionsBegin((PetscObject) sp);
292   PetscCall(PetscOptionsFList("-petscdualspace_type", "Dual space", "PetscDualSpaceSetType", PetscDualSpaceList, defaultType, name, 256, &flg));
293   if (flg) {
294     PetscCall(PetscDualSpaceSetType(sp, name));
295   } else if (!((PetscObject) sp)->type_name) {
296     PetscCall(PetscDualSpaceSetType(sp, defaultType));
297   }
298   PetscCall(PetscOptionsBoundedInt("-petscdualspace_order", "The approximation order", "PetscDualSpaceSetOrder", sp->order, &sp->order, NULL,0));
299   PetscCall(PetscOptionsInt("-petscdualspace_form_degree", "The form degree of the dofs", "PetscDualSpaceSetFormDegree", sp->k, &sp->k, NULL));
300   PetscCall(PetscOptionsBoundedInt("-petscdualspace_components", "The number of components", "PetscDualSpaceSetNumComponents", sp->Nc, &sp->Nc, NULL,1));
301   if (sp->ops->setfromoptions) PetscCall((*sp->ops->setfromoptions)(PetscOptionsObject,sp));
302   PetscCall(PetscOptionsEnum("-petscdualspace_refcell", "Reference cell shape", "PetscDualSpaceSetReferenceCell", DMPolytopeTypes, (PetscEnum) refCell, (PetscEnum *) &refCell, &flg));
303   if (flg) {
304     DM K;
305 
306     PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, refCell, &K));
307     PetscCall(PetscDualSpaceSetDM(sp, K));
308     PetscCall(DMDestroy(&K));
309   }
310 
311   /* process any options handlers added with PetscObjectAddOptionsHandler() */
312   PetscCall(PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) sp));
313   PetscOptionsEnd();
314   sp->setfromoptionscalled = PETSC_TRUE;
315   PetscFunctionReturn(0);
316 }
317 
318 /*@
319   PetscDualSpaceSetUp - Construct a basis for the PetscDualSpace
320 
321   Collective on sp
322 
323   Input Parameter:
324 . sp - the PetscDualSpace object to setup
325 
326   Level: intermediate
327 
328 .seealso `PetscDualSpaceView()`, `PetscDualSpaceDestroy()`, `PetscDualSpace`
329 @*/
330 PetscErrorCode PetscDualSpaceSetUp(PetscDualSpace sp)
331 {
332   PetscFunctionBegin;
333   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
334   if (sp->setupcalled) PetscFunctionReturn(0);
335   PetscCall(PetscLogEventBegin(PETSCDUALSPACE_SetUp, sp, 0, 0, 0));
336   sp->setupcalled = PETSC_TRUE;
337   if (sp->ops->setup) PetscCall((*sp->ops->setup)(sp));
338   PetscCall(PetscLogEventEnd(PETSCDUALSPACE_SetUp, sp, 0, 0, 0));
339   if (sp->setfromoptionscalled) PetscCall(PetscDualSpaceViewFromOptions(sp, NULL, "-petscdualspace_view"));
340   PetscFunctionReturn(0);
341 }
342 
343 static PetscErrorCode PetscDualSpaceClearDMData_Internal(PetscDualSpace sp, DM dm)
344 {
345   PetscInt       pStart = -1, pEnd = -1, depth = -1;
346 
347   PetscFunctionBegin;
348   if (!dm) PetscFunctionReturn(0);
349   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
350   PetscCall(DMPlexGetDepth(dm, &depth));
351 
352   if (sp->pointSpaces) {
353     PetscInt i;
354 
355     for (i = 0; i < pEnd - pStart; i++) {
356       PetscCall(PetscDualSpaceDestroy(&(sp->pointSpaces[i])));
357     }
358   }
359   PetscCall(PetscFree(sp->pointSpaces));
360 
361   if (sp->heightSpaces) {
362     PetscInt i;
363 
364     for (i = 0; i <= depth; i++) {
365       PetscCall(PetscDualSpaceDestroy(&(sp->heightSpaces[i])));
366     }
367   }
368   PetscCall(PetscFree(sp->heightSpaces));
369 
370   PetscCall(PetscSectionDestroy(&(sp->pointSection)));
371   PetscCall(PetscQuadratureDestroy(&(sp->intNodes)));
372   PetscCall(VecDestroy(&(sp->intDofValues)));
373   PetscCall(VecDestroy(&(sp->intNodeValues)));
374   PetscCall(MatDestroy(&(sp->intMat)));
375   PetscCall(PetscQuadratureDestroy(&(sp->allNodes)));
376   PetscCall(VecDestroy(&(sp->allDofValues)));
377   PetscCall(VecDestroy(&(sp->allNodeValues)));
378   PetscCall(MatDestroy(&(sp->allMat)));
379   PetscCall(PetscFree(sp->numDof));
380   PetscFunctionReturn(0);
381 }
382 
383 /*@
384   PetscDualSpaceDestroy - Destroys a PetscDualSpace object
385 
386   Collective on sp
387 
388   Input Parameter:
389 . sp - the PetscDualSpace object to destroy
390 
391   Level: beginner
392 
393 .seealso `PetscDualSpaceView()`, `PetscDualSpace()`, `PetscDualSpaceCreate()`
394 @*/
395 PetscErrorCode PetscDualSpaceDestroy(PetscDualSpace *sp)
396 {
397   PetscInt       dim, f;
398   DM             dm;
399 
400   PetscFunctionBegin;
401   if (!*sp) PetscFunctionReturn(0);
402   PetscValidHeaderSpecific((*sp), PETSCDUALSPACE_CLASSID, 1);
403 
404   if (--((PetscObject)(*sp))->refct > 0) {*sp = NULL; PetscFunctionReturn(0);}
405   ((PetscObject) (*sp))->refct = 0;
406 
407   PetscCall(PetscDualSpaceGetDimension(*sp, &dim));
408   dm = (*sp)->dm;
409 
410   if ((*sp)->ops->destroy) PetscCall((*(*sp)->ops->destroy)(*sp));
411   PetscCall(PetscDualSpaceClearDMData_Internal(*sp, dm));
412 
413   for (f = 0; f < dim; ++f) {
414     PetscCall(PetscQuadratureDestroy(&(*sp)->functional[f]));
415   }
416   PetscCall(PetscFree((*sp)->functional));
417   PetscCall(DMDestroy(&(*sp)->dm));
418   PetscCall(PetscHeaderDestroy(sp));
419   PetscFunctionReturn(0);
420 }
421 
422 /*@
423   PetscDualSpaceCreate - Creates an empty PetscDualSpace object. The type can then be set with PetscDualSpaceSetType().
424 
425   Collective
426 
427   Input Parameter:
428 . comm - The communicator for the PetscDualSpace object
429 
430   Output Parameter:
431 . sp - The PetscDualSpace object
432 
433   Level: beginner
434 
435 .seealso: `PetscDualSpaceSetType()`, `PETSCDUALSPACELAGRANGE`
436 @*/
437 PetscErrorCode PetscDualSpaceCreate(MPI_Comm comm, PetscDualSpace *sp)
438 {
439   PetscDualSpace s;
440 
441   PetscFunctionBegin;
442   PetscValidPointer(sp, 2);
443   PetscCall(PetscCitationsRegister(FECitation,&FEcite));
444   *sp  = NULL;
445   PetscCall(PetscFEInitializePackage());
446 
447   PetscCall(PetscHeaderCreate(s, PETSCDUALSPACE_CLASSID, "PetscDualSpace", "Dual Space", "PetscDualSpace", comm, PetscDualSpaceDestroy, PetscDualSpaceView));
448 
449   s->order       = 0;
450   s->Nc          = 1;
451   s->k           = 0;
452   s->spdim       = -1;
453   s->spintdim    = -1;
454   s->uniform     = PETSC_TRUE;
455   s->setupcalled = PETSC_FALSE;
456 
457   *sp = s;
458   PetscFunctionReturn(0);
459 }
460 
461 /*@
462   PetscDualSpaceDuplicate - Creates a duplicate PetscDualSpace object, however it is not setup.
463 
464   Collective on sp
465 
466   Input Parameter:
467 . sp - The original PetscDualSpace
468 
469   Output Parameter:
470 . spNew - The duplicate PetscDualSpace
471 
472   Level: beginner
473 
474 .seealso: `PetscDualSpaceCreate()`, `PetscDualSpaceSetType()`
475 @*/
476 PetscErrorCode PetscDualSpaceDuplicate(PetscDualSpace sp, PetscDualSpace *spNew)
477 {
478   DM             dm;
479   PetscDualSpaceType type;
480   const char     *name;
481 
482   PetscFunctionBegin;
483   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
484   PetscValidPointer(spNew, 2);
485   PetscCall(PetscDualSpaceCreate(PetscObjectComm((PetscObject)sp), spNew));
486   PetscCall(PetscObjectGetName((PetscObject) sp,     &name));
487   PetscCall(PetscObjectSetName((PetscObject) *spNew,  name));
488   PetscCall(PetscDualSpaceGetType(sp, &type));
489   PetscCall(PetscDualSpaceSetType(*spNew, type));
490   PetscCall(PetscDualSpaceGetDM(sp, &dm));
491   PetscCall(PetscDualSpaceSetDM(*spNew, dm));
492 
493   (*spNew)->order   = sp->order;
494   (*spNew)->k       = sp->k;
495   (*spNew)->Nc      = sp->Nc;
496   (*spNew)->uniform = sp->uniform;
497   if (sp->ops->duplicate) PetscCall((*sp->ops->duplicate)(sp, *spNew));
498   PetscFunctionReturn(0);
499 }
500 
501 /*@
502   PetscDualSpaceGetDM - Get the DM representing the reference cell
503 
504   Not collective
505 
506   Input Parameter:
507 . sp - The PetscDualSpace
508 
509   Output Parameter:
510 . dm - The reference cell
511 
512   Level: intermediate
513 
514 .seealso: `PetscDualSpaceSetDM()`, `PetscDualSpaceCreate()`
515 @*/
516 PetscErrorCode PetscDualSpaceGetDM(PetscDualSpace sp, DM *dm)
517 {
518   PetscFunctionBegin;
519   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
520   PetscValidPointer(dm, 2);
521   *dm = sp->dm;
522   PetscFunctionReturn(0);
523 }
524 
525 /*@
526   PetscDualSpaceSetDM - Get the DM representing the reference cell
527 
528   Not collective
529 
530   Input Parameters:
531 + sp - The PetscDualSpace
532 - dm - The reference cell
533 
534   Level: intermediate
535 
536 .seealso: `PetscDualSpaceGetDM()`, `PetscDualSpaceCreate()`
537 @*/
538 PetscErrorCode PetscDualSpaceSetDM(PetscDualSpace sp, DM dm)
539 {
540   PetscFunctionBegin;
541   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
542   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
543   PetscCheck(!sp->setupcalled,PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change DM after dualspace is set up");
544   PetscCall(PetscObjectReference((PetscObject) dm));
545   if (sp->dm && sp->dm != dm) {
546     PetscCall(PetscDualSpaceClearDMData_Internal(sp, sp->dm));
547   }
548   PetscCall(DMDestroy(&sp->dm));
549   sp->dm = dm;
550   PetscFunctionReturn(0);
551 }
552 
553 /*@
554   PetscDualSpaceGetOrder - Get the order of the dual space
555 
556   Not collective
557 
558   Input Parameter:
559 . sp - The PetscDualSpace
560 
561   Output Parameter:
562 . order - The order
563 
564   Level: intermediate
565 
566 .seealso: `PetscDualSpaceSetOrder()`, `PetscDualSpaceCreate()`
567 @*/
568 PetscErrorCode PetscDualSpaceGetOrder(PetscDualSpace sp, PetscInt *order)
569 {
570   PetscFunctionBegin;
571   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
572   PetscValidIntPointer(order, 2);
573   *order = sp->order;
574   PetscFunctionReturn(0);
575 }
576 
577 /*@
578   PetscDualSpaceSetOrder - Set the order of the dual space
579 
580   Not collective
581 
582   Input Parameters:
583 + sp - The PetscDualSpace
584 - order - The order
585 
586   Level: intermediate
587 
588 .seealso: `PetscDualSpaceGetOrder()`, `PetscDualSpaceCreate()`
589 @*/
590 PetscErrorCode PetscDualSpaceSetOrder(PetscDualSpace sp, PetscInt order)
591 {
592   PetscFunctionBegin;
593   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
594   PetscCheck(!sp->setupcalled,PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change order after dualspace is set up");
595   sp->order = order;
596   PetscFunctionReturn(0);
597 }
598 
599 /*@
600   PetscDualSpaceGetNumComponents - Return the number of components for this space
601 
602   Input Parameter:
603 . sp - The PetscDualSpace
604 
605   Output Parameter:
606 . Nc - The number of components
607 
608   Note: A vector space, for example, will have d components, where d is the spatial dimension
609 
610   Level: intermediate
611 
612 .seealso: `PetscDualSpaceSetNumComponents()`, `PetscDualSpaceGetDimension()`, `PetscDualSpaceCreate()`, `PetscDualSpace`
613 @*/
614 PetscErrorCode PetscDualSpaceGetNumComponents(PetscDualSpace sp, PetscInt *Nc)
615 {
616   PetscFunctionBegin;
617   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
618   PetscValidIntPointer(Nc, 2);
619   *Nc = sp->Nc;
620   PetscFunctionReturn(0);
621 }
622 
623 /*@
624   PetscDualSpaceSetNumComponents - Set the number of components for this space
625 
626   Input Parameters:
627 + sp - The PetscDualSpace
628 - order - The number of components
629 
630   Level: intermediate
631 
632 .seealso: `PetscDualSpaceGetNumComponents()`, `PetscDualSpaceCreate()`, `PetscDualSpace`
633 @*/
634 PetscErrorCode PetscDualSpaceSetNumComponents(PetscDualSpace sp, PetscInt Nc)
635 {
636   PetscFunctionBegin;
637   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
638   PetscCheck(!sp->setupcalled,PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change number of components after dualspace is set up");
639   sp->Nc = Nc;
640   PetscFunctionReturn(0);
641 }
642 
643 /*@
644   PetscDualSpaceGetFunctional - Get the i-th basis functional in the dual space
645 
646   Not collective
647 
648   Input Parameters:
649 + sp - The PetscDualSpace
650 - i  - The basis number
651 
652   Output Parameter:
653 . functional - The basis functional
654 
655   Level: intermediate
656 
657 .seealso: `PetscDualSpaceGetDimension()`, `PetscDualSpaceCreate()`
658 @*/
659 PetscErrorCode PetscDualSpaceGetFunctional(PetscDualSpace sp, PetscInt i, PetscQuadrature *functional)
660 {
661   PetscInt       dim;
662 
663   PetscFunctionBegin;
664   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
665   PetscValidPointer(functional, 3);
666   PetscCall(PetscDualSpaceGetDimension(sp, &dim));
667   PetscCheck(!(i < 0) && !(i >= dim),PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Functional index %" PetscInt_FMT " must be in [0, %" PetscInt_FMT ")", i, dim);
668   *functional = sp->functional[i];
669   PetscFunctionReturn(0);
670 }
671 
672 /*@
673   PetscDualSpaceGetDimension - Get the dimension of the dual space, i.e. the number of basis functionals
674 
675   Not collective
676 
677   Input Parameter:
678 . sp - The PetscDualSpace
679 
680   Output Parameter:
681 . dim - The dimension
682 
683   Level: intermediate
684 
685 .seealso: `PetscDualSpaceGetFunctional()`, `PetscDualSpaceCreate()`
686 @*/
687 PetscErrorCode PetscDualSpaceGetDimension(PetscDualSpace sp, PetscInt *dim)
688 {
689   PetscFunctionBegin;
690   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
691   PetscValidIntPointer(dim, 2);
692   if (sp->spdim < 0) {
693     PetscSection section;
694 
695     PetscCall(PetscDualSpaceGetSection(sp, &section));
696     if (section) {
697       PetscCall(PetscSectionGetStorageSize(section, &(sp->spdim)));
698     } else sp->spdim = 0;
699   }
700   *dim = sp->spdim;
701   PetscFunctionReturn(0);
702 }
703 
704 /*@
705   PetscDualSpaceGetInteriorDimension - Get the interior dimension of the dual space, i.e. the number of basis functionals assigned to the interior of the reference domain
706 
707   Not collective
708 
709   Input Parameter:
710 . sp - The PetscDualSpace
711 
712   Output Parameter:
713 . dim - The dimension
714 
715   Level: intermediate
716 
717 .seealso: `PetscDualSpaceGetFunctional()`, `PetscDualSpaceCreate()`
718 @*/
719 PetscErrorCode PetscDualSpaceGetInteriorDimension(PetscDualSpace sp, PetscInt *intdim)
720 {
721   PetscFunctionBegin;
722   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
723   PetscValidIntPointer(intdim, 2);
724   if (sp->spintdim < 0) {
725     PetscSection section;
726 
727     PetscCall(PetscDualSpaceGetSection(sp, &section));
728     if (section) {
729       PetscCall(PetscSectionGetConstrainedStorageSize(section, &(sp->spintdim)));
730     } else sp->spintdim = 0;
731   }
732   *intdim = sp->spintdim;
733   PetscFunctionReturn(0);
734 }
735 
736 /*@
737    PetscDualSpaceGetUniform - Whether this dual space is uniform
738 
739    Not collective
740 
741    Input Parameters:
742 .  sp - A dual space
743 
744    Output Parameters:
745 .  uniform - PETSC_TRUE if (a) the dual space is the same for each point in a stratum of the reference DMPlex, and
746              (b) every symmetry of each point in the reference DMPlex is also a symmetry of the point's dual space.
747 
748    Level: advanced
749 
750    Note: all of the usual spaces on simplex or tensor-product elements will be uniform, only reference cells
751    with non-uniform strata (like trianguar-prisms) or anisotropic hp dual spaces will not be uniform.
752 
753 .seealso: `PetscDualSpaceGetPointSubspace()`, `PetscDualSpaceGetSymmetries()`
754 @*/
755 PetscErrorCode PetscDualSpaceGetUniform(PetscDualSpace sp, PetscBool *uniform)
756 {
757   PetscFunctionBegin;
758   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
759   PetscValidBoolPointer(uniform, 2);
760   *uniform = sp->uniform;
761   PetscFunctionReturn(0);
762 }
763 
764 /*@C
765   PetscDualSpaceGetNumDof - Get the number of degrees of freedom for each spatial (topological) dimension
766 
767   Not collective
768 
769   Input Parameter:
770 . sp - The PetscDualSpace
771 
772   Output Parameter:
773 . numDof - An array of length dim+1 which holds the number of dofs for each dimension
774 
775   Level: intermediate
776 
777 .seealso: `PetscDualSpaceGetFunctional()`, `PetscDualSpaceCreate()`
778 @*/
779 PetscErrorCode PetscDualSpaceGetNumDof(PetscDualSpace sp, const PetscInt **numDof)
780 {
781   PetscFunctionBegin;
782   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
783   PetscValidPointer(numDof, 2);
784   PetscCheck(sp->uniform,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "A non-uniform space does not have a fixed number of dofs for each height");
785   if (!sp->numDof) {
786     DM       dm;
787     PetscInt depth, d;
788     PetscSection section;
789 
790     PetscCall(PetscDualSpaceGetDM(sp, &dm));
791     PetscCall(DMPlexGetDepth(dm, &depth));
792     PetscCall(PetscCalloc1(depth+1,&(sp->numDof)));
793     PetscCall(PetscDualSpaceGetSection(sp, &section));
794     for (d = 0; d <= depth; d++) {
795       PetscInt dStart, dEnd;
796 
797       PetscCall(DMPlexGetDepthStratum(dm, d, &dStart, &dEnd));
798       if (dEnd <= dStart) continue;
799       PetscCall(PetscSectionGetDof(section, dStart, &(sp->numDof[d])));
800 
801     }
802   }
803   *numDof = sp->numDof;
804   PetscCheck(*numDof,PetscObjectComm((PetscObject) sp), PETSC_ERR_LIB, "Empty numDof[] returned from dual space implementation");
805   PetscFunctionReturn(0);
806 }
807 
808 /* create the section of the right size and set a permutation for topological ordering */
809 PetscErrorCode PetscDualSpaceSectionCreate_Internal(PetscDualSpace sp, PetscSection *topSection)
810 {
811   DM             dm;
812   PetscInt       pStart, pEnd, cStart, cEnd, c, depth, count, i;
813   PetscInt       *seen, *perm;
814   PetscSection   section;
815 
816   PetscFunctionBegin;
817   dm = sp->dm;
818   PetscCall(PetscSectionCreate(PETSC_COMM_SELF, &section));
819   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
820   PetscCall(PetscSectionSetChart(section, pStart, pEnd));
821   PetscCall(PetscCalloc1(pEnd - pStart, &seen));
822   PetscCall(PetscMalloc1(pEnd - pStart, &perm));
823   PetscCall(DMPlexGetDepth(dm, &depth));
824   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
825   for (c = cStart, count = 0; c < cEnd; c++) {
826     PetscInt closureSize = -1, e;
827     PetscInt *closure = NULL;
828 
829     perm[count++] = c;
830     seen[c-pStart] = 1;
831     PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
832     for (e = 0; e < closureSize; e++) {
833       PetscInt point = closure[2*e];
834 
835       if (seen[point-pStart]) continue;
836       perm[count++] = point;
837       seen[point-pStart] = 1;
838     }
839     PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure));
840   }
841   PetscCheck(count == pEnd - pStart,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Bad topological ordering");
842   for (i = 0; i < pEnd - pStart; i++) if (perm[i] != i) break;
843   if (i < pEnd - pStart) {
844     IS permIS;
845 
846     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, pEnd - pStart, perm, PETSC_OWN_POINTER, &permIS));
847     PetscCall(ISSetPermutation(permIS));
848     PetscCall(PetscSectionSetPermutation(section, permIS));
849     PetscCall(ISDestroy(&permIS));
850   } else {
851     PetscCall(PetscFree(perm));
852   }
853   PetscCall(PetscFree(seen));
854   *topSection = section;
855   PetscFunctionReturn(0);
856 }
857 
858 /* mark boundary points and set up */
859 PetscErrorCode PetscDualSpaceSectionSetUp_Internal(PetscDualSpace sp, PetscSection section)
860 {
861   DM             dm;
862   DMLabel        boundary;
863   PetscInt       pStart, pEnd, p;
864 
865   PetscFunctionBegin;
866   dm = sp->dm;
867   PetscCall(DMLabelCreate(PETSC_COMM_SELF,"boundary",&boundary));
868   PetscCall(PetscDualSpaceGetDM(sp,&dm));
869   PetscCall(DMPlexMarkBoundaryFaces(dm,1,boundary));
870   PetscCall(DMPlexLabelComplete(dm,boundary));
871   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
872   for (p = pStart; p < pEnd; p++) {
873     PetscInt bval;
874 
875     PetscCall(DMLabelGetValue(boundary, p, &bval));
876     if (bval == 1) {
877       PetscInt dof;
878 
879       PetscCall(PetscSectionGetDof(section, p, &dof));
880       PetscCall(PetscSectionSetConstraintDof(section, p, dof));
881     }
882   }
883   PetscCall(DMLabelDestroy(&boundary));
884   PetscCall(PetscSectionSetUp(section));
885   PetscFunctionReturn(0);
886 }
887 
888 /*@
889   PetscDualSpaceGetSection - Create a PetscSection over the reference cell with the layout from this space
890 
891   Collective on sp
892 
893   Input Parameters:
894 . sp      - The PetscDualSpace
895 
896   Output Parameter:
897 . section - The section
898 
899   Level: advanced
900 
901 .seealso: `PetscDualSpaceCreate()`, `DMPLEX`
902 @*/
903 PetscErrorCode PetscDualSpaceGetSection(PetscDualSpace sp, PetscSection *section)
904 {
905   PetscInt       pStart, pEnd, p;
906 
907   PetscFunctionBegin;
908   if (!sp->dm) {
909     *section = NULL;
910     PetscFunctionReturn(0);
911   }
912   if (!sp->pointSection) {
913     /* mark the boundary */
914     PetscCall(PetscDualSpaceSectionCreate_Internal(sp, &(sp->pointSection)));
915     PetscCall(DMPlexGetChart(sp->dm,&pStart,&pEnd));
916     for (p = pStart; p < pEnd; p++) {
917       PetscDualSpace psp;
918 
919       PetscCall(PetscDualSpaceGetPointSubspace(sp, p, &psp));
920       if (psp) {
921         PetscInt dof;
922 
923         PetscCall(PetscDualSpaceGetInteriorDimension(psp, &dof));
924         PetscCall(PetscSectionSetDof(sp->pointSection,p,dof));
925       }
926     }
927     PetscCall(PetscDualSpaceSectionSetUp_Internal(sp,sp->pointSection));
928   }
929   *section = sp->pointSection;
930   PetscFunctionReturn(0);
931 }
932 
933 /* this assumes that all of the point dual spaces store their interior dofs first, which is true when the point DMs
934  * have one cell */
935 PetscErrorCode PetscDualSpacePushForwardSubspaces_Internal(PetscDualSpace sp, PetscInt sStart, PetscInt sEnd)
936 {
937   PetscReal *sv0, *v0, *J;
938   PetscSection section;
939   PetscInt     dim, s, k;
940   DM             dm;
941 
942   PetscFunctionBegin;
943   PetscCall(PetscDualSpaceGetDM(sp, &dm));
944   PetscCall(DMGetDimension(dm, &dim));
945   PetscCall(PetscDualSpaceGetSection(sp, &section));
946   PetscCall(PetscMalloc3(dim, &v0, dim, &sv0, dim*dim, &J));
947   PetscCall(PetscDualSpaceGetFormDegree(sp, &k));
948   for (s = sStart; s < sEnd; s++) {
949     PetscReal detJ, hdetJ;
950     PetscDualSpace ssp;
951     PetscInt dof, off, f, sdim;
952     PetscInt i, j;
953     DM sdm;
954 
955     PetscCall(PetscDualSpaceGetPointSubspace(sp, s, &ssp));
956     if (!ssp) continue;
957     PetscCall(PetscSectionGetDof(section, s, &dof));
958     PetscCall(PetscSectionGetOffset(section, s, &off));
959     /* get the first vertex of the reference cell */
960     PetscCall(PetscDualSpaceGetDM(ssp, &sdm));
961     PetscCall(DMGetDimension(sdm, &sdim));
962     PetscCall(DMPlexComputeCellGeometryAffineFEM(sdm, 0, sv0, NULL, NULL, &hdetJ));
963     PetscCall(DMPlexComputeCellGeometryAffineFEM(dm, s, v0, J, NULL, &detJ));
964     /* compactify Jacobian */
965     for (i = 0; i < dim; i++) for (j = 0; j < sdim; j++) J[i* sdim + j] = J[i * dim + j];
966     for (f = 0; f < dof; f++) {
967       PetscQuadrature fn;
968 
969       PetscCall(PetscDualSpaceGetFunctional(ssp, f, &fn));
970       PetscCall(PetscQuadraturePushForward(fn, dim, sv0, v0, J, k, &(sp->functional[off+f])));
971     }
972   }
973   PetscCall(PetscFree3(v0, sv0, J));
974   PetscFunctionReturn(0);
975 }
976 
977 /*@C
978   PetscDualSpaceApply - Apply a functional from the dual space basis to an input function
979 
980   Input Parameters:
981 + sp      - The PetscDualSpace object
982 . f       - The basis functional index
983 . time    - The time
984 . 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)
985 . numComp - The number of components for the function
986 . func    - The input function
987 - ctx     - A context for the function
988 
989   Output Parameter:
990 . value   - numComp output values
991 
992   Note: The calling sequence for the callback func is given by:
993 
994 $ func(PetscInt dim, PetscReal time, const PetscReal x[],
995 $      PetscInt numComponents, PetscScalar values[], void *ctx)
996 
997   Level: beginner
998 
999 .seealso: `PetscDualSpaceCreate()`
1000 @*/
1001 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)
1002 {
1003   PetscFunctionBegin;
1004   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1005   PetscValidPointer(cgeom, 4);
1006   PetscValidScalarPointer(value, 8);
1007   PetscCall((*sp->ops->apply)(sp, f, time, cgeom, numComp, func, ctx, value));
1008   PetscFunctionReturn(0);
1009 }
1010 
1011 /*@C
1012   PetscDualSpaceApplyAll - Apply all functionals from the dual space basis to the result of an evaluation at the points returned by PetscDualSpaceGetAllData()
1013 
1014   Input Parameters:
1015 + sp        - The PetscDualSpace object
1016 - pointEval - Evaluation at the points returned by PetscDualSpaceGetAllData()
1017 
1018   Output Parameter:
1019 . spValue   - The values of all dual space functionals
1020 
1021   Level: beginner
1022 
1023 .seealso: `PetscDualSpaceCreate()`
1024 @*/
1025 PetscErrorCode PetscDualSpaceApplyAll(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue)
1026 {
1027   PetscFunctionBegin;
1028   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1029   PetscCall((*sp->ops->applyall)(sp, pointEval, spValue));
1030   PetscFunctionReturn(0);
1031 }
1032 
1033 /*@C
1034   PetscDualSpaceApplyInterior - Apply interior functionals from the dual space basis to the result of an evaluation at the points returned by PetscDualSpaceGetInteriorData()
1035 
1036   Input Parameters:
1037 + sp        - The PetscDualSpace object
1038 - pointEval - Evaluation at the points returned by PetscDualSpaceGetInteriorData()
1039 
1040   Output Parameter:
1041 . spValue   - The values of interior dual space functionals
1042 
1043   Level: beginner
1044 
1045 .seealso: `PetscDualSpaceCreate()`
1046 @*/
1047 PetscErrorCode PetscDualSpaceApplyInterior(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue)
1048 {
1049   PetscFunctionBegin;
1050   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1051   PetscCall((*sp->ops->applyint)(sp, pointEval, spValue));
1052   PetscFunctionReturn(0);
1053 }
1054 
1055 /*@C
1056   PetscDualSpaceApplyDefault - Apply a functional from the dual space basis to an input function by assuming a point evaluation functional.
1057 
1058   Input Parameters:
1059 + sp    - The PetscDualSpace object
1060 . f     - The basis functional index
1061 . time  - The time
1062 . cgeom - A context with geometric information for this cell, we use v0 (the initial vertex) and J (the Jacobian)
1063 . Nc    - The number of components for the function
1064 . func  - The input function
1065 - ctx   - A context for the function
1066 
1067   Output Parameter:
1068 . value   - The output value
1069 
1070   Note: The calling sequence for the callback func is given by:
1071 
1072 $ func(PetscInt dim, PetscReal time, const PetscReal x[],
1073 $      PetscInt numComponents, PetscScalar values[], void *ctx)
1074 
1075 and the idea is to evaluate the functional as an integral
1076 
1077 $ n(f) = int dx n(x) . f(x)
1078 
1079 where both n and f have Nc components.
1080 
1081   Level: beginner
1082 
1083 .seealso: `PetscDualSpaceCreate()`
1084 @*/
1085 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)
1086 {
1087   DM               dm;
1088   PetscQuadrature  n;
1089   const PetscReal *points, *weights;
1090   PetscReal        x[3];
1091   PetscScalar     *val;
1092   PetscInt         dim, dE, qNc, c, Nq, q;
1093   PetscBool        isAffine;
1094 
1095   PetscFunctionBegin;
1096   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1097   PetscValidScalarPointer(value, 8);
1098   PetscCall(PetscDualSpaceGetDM(sp, &dm));
1099   PetscCall(PetscDualSpaceGetFunctional(sp, f, &n));
1100   PetscCall(PetscQuadratureGetData(n, &dim, &qNc, &Nq, &points, &weights));
1101   PetscCheck(dim == cgeom->dim,PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature spatial dimension %" PetscInt_FMT " != cell geometry dimension %" PetscInt_FMT, dim, cgeom->dim);
1102   PetscCheck(qNc == Nc,PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature components %" PetscInt_FMT " != function components %" PetscInt_FMT, qNc, Nc);
1103   PetscCall(DMGetWorkArray(dm, Nc, MPIU_SCALAR, &val));
1104   *value = 0.0;
1105   isAffine = cgeom->isAffine;
1106   dE = cgeom->dimEmbed;
1107   for (q = 0; q < Nq; ++q) {
1108     if (isAffine) {
1109       CoordinatesRefToReal(dE, cgeom->dim, cgeom->xi, cgeom->v, cgeom->J, &points[q*dim], x);
1110       PetscCall((*func)(dE, time, x, Nc, val, ctx));
1111     } else {
1112       PetscCall((*func)(dE, time, &cgeom->v[dE*q], Nc, val, ctx));
1113     }
1114     for (c = 0; c < Nc; ++c) {
1115       *value += val[c]*weights[q*Nc+c];
1116     }
1117   }
1118   PetscCall(DMRestoreWorkArray(dm, Nc, MPIU_SCALAR, &val));
1119   PetscFunctionReturn(0);
1120 }
1121 
1122 /*@C
1123   PetscDualSpaceApplyAllDefault - Apply all functionals from the dual space basis to the result of an evaluation at the points returned by PetscDualSpaceGetAllData()
1124 
1125   Input Parameters:
1126 + sp        - The PetscDualSpace object
1127 - pointEval - Evaluation at the points returned by PetscDualSpaceGetAllData()
1128 
1129   Output Parameter:
1130 . spValue   - The values of all dual space functionals
1131 
1132   Level: beginner
1133 
1134 .seealso: `PetscDualSpaceCreate()`
1135 @*/
1136 PetscErrorCode PetscDualSpaceApplyAllDefault(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue)
1137 {
1138   Vec              pointValues, dofValues;
1139   Mat              allMat;
1140 
1141   PetscFunctionBegin;
1142   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1143   PetscValidScalarPointer(pointEval, 2);
1144   PetscValidScalarPointer(spValue, 3);
1145   PetscCall(PetscDualSpaceGetAllData(sp, NULL, &allMat));
1146   if (!(sp->allNodeValues)) {
1147     PetscCall(MatCreateVecs(allMat, &(sp->allNodeValues), NULL));
1148   }
1149   pointValues = sp->allNodeValues;
1150   if (!(sp->allDofValues)) {
1151     PetscCall(MatCreateVecs(allMat, NULL, &(sp->allDofValues)));
1152   }
1153   dofValues = sp->allDofValues;
1154   PetscCall(VecPlaceArray(pointValues, pointEval));
1155   PetscCall(VecPlaceArray(dofValues, spValue));
1156   PetscCall(MatMult(allMat, pointValues, dofValues));
1157   PetscCall(VecResetArray(dofValues));
1158   PetscCall(VecResetArray(pointValues));
1159   PetscFunctionReturn(0);
1160 }
1161 
1162 /*@C
1163   PetscDualSpaceApplyInteriorDefault - Apply interior functionals from the dual space basis to the result of an evaluation at the points returned by PetscDualSpaceGetInteriorData()
1164 
1165   Input Parameters:
1166 + sp        - The PetscDualSpace object
1167 - pointEval - Evaluation at the points returned by PetscDualSpaceGetInteriorData()
1168 
1169   Output Parameter:
1170 . spValue   - The values of interior dual space functionals
1171 
1172   Level: beginner
1173 
1174 .seealso: `PetscDualSpaceCreate()`
1175 @*/
1176 PetscErrorCode PetscDualSpaceApplyInteriorDefault(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue)
1177 {
1178   Vec              pointValues, dofValues;
1179   Mat              intMat;
1180 
1181   PetscFunctionBegin;
1182   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1183   PetscValidScalarPointer(pointEval, 2);
1184   PetscValidScalarPointer(spValue, 3);
1185   PetscCall(PetscDualSpaceGetInteriorData(sp, NULL, &intMat));
1186   if (!(sp->intNodeValues)) {
1187     PetscCall(MatCreateVecs(intMat, &(sp->intNodeValues), NULL));
1188   }
1189   pointValues = sp->intNodeValues;
1190   if (!(sp->intDofValues)) {
1191     PetscCall(MatCreateVecs(intMat, NULL, &(sp->intDofValues)));
1192   }
1193   dofValues = sp->intDofValues;
1194   PetscCall(VecPlaceArray(pointValues, pointEval));
1195   PetscCall(VecPlaceArray(dofValues, spValue));
1196   PetscCall(MatMult(intMat, pointValues, dofValues));
1197   PetscCall(VecResetArray(dofValues));
1198   PetscCall(VecResetArray(pointValues));
1199   PetscFunctionReturn(0);
1200 }
1201 
1202 /*@
1203   PetscDualSpaceGetAllData - Get all quadrature nodes from this space, and the matrix that sends quadrature node values to degree-of-freedom values
1204 
1205   Input Parameter:
1206 . sp - The dualspace
1207 
1208   Output Parameters:
1209 + allNodes - A PetscQuadrature object containing all evaluation nodes
1210 - allMat - A Mat for the node-to-dof transformation
1211 
1212   Level: advanced
1213 
1214 .seealso: `PetscDualSpaceCreate()`
1215 @*/
1216 PetscErrorCode PetscDualSpaceGetAllData(PetscDualSpace sp, PetscQuadrature *allNodes, Mat *allMat)
1217 {
1218   PetscFunctionBegin;
1219   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1220   if (allNodes) PetscValidPointer(allNodes,2);
1221   if (allMat) PetscValidPointer(allMat,3);
1222   if ((!sp->allNodes || !sp->allMat) && sp->ops->createalldata) {
1223     PetscQuadrature qpoints;
1224     Mat amat;
1225 
1226     PetscCall((*sp->ops->createalldata)(sp,&qpoints,&amat));
1227     PetscCall(PetscQuadratureDestroy(&(sp->allNodes)));
1228     PetscCall(MatDestroy(&(sp->allMat)));
1229     sp->allNodes = qpoints;
1230     sp->allMat = amat;
1231   }
1232   if (allNodes) *allNodes = sp->allNodes;
1233   if (allMat) *allMat = sp->allMat;
1234   PetscFunctionReturn(0);
1235 }
1236 
1237 /*@
1238   PetscDualSpaceCreateAllDataDefault - Create all evaluation nodes and the node-to-dof matrix by examining functionals
1239 
1240   Input Parameter:
1241 . sp - The dualspace
1242 
1243   Output Parameters:
1244 + allNodes - A PetscQuadrature object containing all evaluation nodes
1245 - allMat - A Mat for the node-to-dof transformation
1246 
1247   Level: advanced
1248 
1249 .seealso: `PetscDualSpaceCreate()`
1250 @*/
1251 PetscErrorCode PetscDualSpaceCreateAllDataDefault(PetscDualSpace sp, PetscQuadrature *allNodes, Mat *allMat)
1252 {
1253   PetscInt        spdim;
1254   PetscInt        numPoints, offset;
1255   PetscReal       *points;
1256   PetscInt        f, dim;
1257   PetscInt        Nc, nrows, ncols;
1258   PetscInt        maxNumPoints;
1259   PetscQuadrature q;
1260   Mat             A;
1261 
1262   PetscFunctionBegin;
1263   PetscCall(PetscDualSpaceGetNumComponents(sp, &Nc));
1264   PetscCall(PetscDualSpaceGetDimension(sp,&spdim));
1265   if (!spdim) {
1266     PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF,allNodes));
1267     PetscCall(PetscQuadratureSetData(*allNodes,0,0,0,NULL,NULL));
1268   }
1269   nrows = spdim;
1270   PetscCall(PetscDualSpaceGetFunctional(sp,0,&q));
1271   PetscCall(PetscQuadratureGetData(q,&dim,NULL,&numPoints,NULL,NULL));
1272   maxNumPoints = numPoints;
1273   for (f = 1; f < spdim; f++) {
1274     PetscInt Np;
1275 
1276     PetscCall(PetscDualSpaceGetFunctional(sp,f,&q));
1277     PetscCall(PetscQuadratureGetData(q,NULL,NULL,&Np,NULL,NULL));
1278     numPoints += Np;
1279     maxNumPoints = PetscMax(maxNumPoints,Np);
1280   }
1281   ncols = numPoints * Nc;
1282   PetscCall(PetscMalloc1(dim*numPoints,&points));
1283   PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, nrows, ncols, maxNumPoints * Nc, NULL, &A));
1284   for (f = 0, offset = 0; f < spdim; f++) {
1285     const PetscReal *p, *w;
1286     PetscInt        Np, i;
1287     PetscInt        fnc;
1288 
1289     PetscCall(PetscDualSpaceGetFunctional(sp,f,&q));
1290     PetscCall(PetscQuadratureGetData(q,NULL,&fnc,&Np,&p,&w));
1291     PetscCheck(fnc == Nc,PETSC_COMM_SELF, PETSC_ERR_PLIB, "functional component mismatch");
1292     for (i = 0; i < Np * dim; i++) {
1293       points[offset* dim + i] = p[i];
1294     }
1295     for (i = 0; i < Np * Nc; i++) {
1296       PetscCall(MatSetValue(A, f, offset * Nc, w[i], INSERT_VALUES));
1297     }
1298     offset += Np;
1299   }
1300   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1301   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
1302   PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF,allNodes));
1303   PetscCall(PetscQuadratureSetData(*allNodes,dim,0,numPoints,points,NULL));
1304   *allMat = A;
1305   PetscFunctionReturn(0);
1306 }
1307 
1308 /*@
1309   PetscDualSpaceGetInteriorData - Get all quadrature points necessary to compute the interior degrees of freedom from
1310   this space, as well as the matrix that computes the degrees of freedom from the quadrature values.  Degrees of
1311   freedom are interior degrees of freedom if they belong (by PetscDualSpaceGetSection()) to interior points in the
1312   reference DMPlex: complementary boundary degrees of freedom are marked as constrained in the section returned by
1313   PetscDualSpaceGetSection()).
1314 
1315   Input Parameter:
1316 . sp - The dualspace
1317 
1318   Output Parameters:
1319 + intNodes - A PetscQuadrature object containing all evaluation points needed to evaluate interior degrees of freedom
1320 - intMat   - A matrix that computes dual space values from point values: size [spdim0 x (npoints * nc)], where spdim0 is
1321              the size of the constrained layout (PetscSectionGetConstrainStorageSize()) of the dual space section,
1322              npoints is the number of points in intNodes and nc is PetscDualSpaceGetNumComponents().
1323 
1324   Level: advanced
1325 
1326 .seealso: `PetscDualSpaceCreate()`, `PetscDualSpaceGetDimension()`, `PetscDualSpaceGetNumComponents()`, `PetscQuadratureGetData()`
1327 @*/
1328 PetscErrorCode PetscDualSpaceGetInteriorData(PetscDualSpace sp, PetscQuadrature *intNodes, Mat *intMat)
1329 {
1330   PetscFunctionBegin;
1331   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1332   if (intNodes) PetscValidPointer(intNodes,2);
1333   if (intMat) PetscValidPointer(intMat,3);
1334   if ((!sp->intNodes || !sp->intMat) && sp->ops->createintdata) {
1335     PetscQuadrature qpoints;
1336     Mat imat;
1337 
1338     PetscCall((*sp->ops->createintdata)(sp,&qpoints,&imat));
1339     PetscCall(PetscQuadratureDestroy(&(sp->intNodes)));
1340     PetscCall(MatDestroy(&(sp->intMat)));
1341     sp->intNodes = qpoints;
1342     sp->intMat = imat;
1343   }
1344   if (intNodes) *intNodes = sp->intNodes;
1345   if (intMat) *intMat = sp->intMat;
1346   PetscFunctionReturn(0);
1347 }
1348 
1349 /*@
1350   PetscDualSpaceCreateInteriorDataDefault - Create quadrature points by examining interior functionals and create the matrix mapping quadrature point values to interior dual space values
1351 
1352   Input Parameter:
1353 . sp - The dualspace
1354 
1355   Output Parameters:
1356 + intNodes - A PetscQuadrature object containing all evaluation points needed to evaluate interior degrees of freedom
1357 - intMat    - A matrix that computes dual space values from point values: size [spdim0 x (npoints * nc)], where spdim0 is
1358               the size of the constrained layout (PetscSectionGetConstrainStorageSize()) of the dual space section,
1359               npoints is the number of points in allNodes and nc is PetscDualSpaceGetNumComponents().
1360 
1361   Level: advanced
1362 
1363 .seealso: `PetscDualSpaceCreate()`, `PetscDualSpaceGetInteriorData()`
1364 @*/
1365 PetscErrorCode PetscDualSpaceCreateInteriorDataDefault(PetscDualSpace sp, PetscQuadrature *intNodes, Mat *intMat)
1366 {
1367   DM              dm;
1368   PetscInt        spdim0;
1369   PetscInt        Nc;
1370   PetscInt        pStart, pEnd, p, f;
1371   PetscSection    section;
1372   PetscInt        numPoints, offset, matoffset;
1373   PetscReal       *points;
1374   PetscInt        dim;
1375   PetscInt        *nnz;
1376   PetscQuadrature q;
1377   Mat             imat;
1378 
1379   PetscFunctionBegin;
1380   PetscValidHeaderSpecific(sp,PETSCDUALSPACE_CLASSID,1);
1381   PetscCall(PetscDualSpaceGetSection(sp, &section));
1382   PetscCall(PetscSectionGetConstrainedStorageSize(section, &spdim0));
1383   if (!spdim0) {
1384     *intNodes = NULL;
1385     *intMat = NULL;
1386     PetscFunctionReturn(0);
1387   }
1388   PetscCall(PetscDualSpaceGetNumComponents(sp, &Nc));
1389   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
1390   PetscCall(PetscDualSpaceGetDM(sp, &dm));
1391   PetscCall(DMGetDimension(dm, &dim));
1392   PetscCall(PetscMalloc1(spdim0, &nnz));
1393   for (p = pStart, f = 0, numPoints = 0; p < pEnd; p++) {
1394     PetscInt dof, cdof, off, d;
1395 
1396     PetscCall(PetscSectionGetDof(section, p, &dof));
1397     PetscCall(PetscSectionGetConstraintDof(section, p, &cdof));
1398     if (!(dof - cdof)) continue;
1399     PetscCall(PetscSectionGetOffset(section, p, &off));
1400     for (d = 0; d < dof; d++, off++, f++) {
1401       PetscInt Np;
1402 
1403       PetscCall(PetscDualSpaceGetFunctional(sp,off,&q));
1404       PetscCall(PetscQuadratureGetData(q,NULL,NULL,&Np,NULL,NULL));
1405       nnz[f] = Np * Nc;
1406       numPoints += Np;
1407     }
1408   }
1409   PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, spdim0, numPoints * Nc, 0, nnz, &imat));
1410   PetscCall(PetscFree(nnz));
1411   PetscCall(PetscMalloc1(dim*numPoints,&points));
1412   for (p = pStart, f = 0, offset = 0, matoffset = 0; p < pEnd; p++) {
1413     PetscInt dof, cdof, off, d;
1414 
1415     PetscCall(PetscSectionGetDof(section, p, &dof));
1416     PetscCall(PetscSectionGetConstraintDof(section, p, &cdof));
1417     if (!(dof - cdof)) continue;
1418     PetscCall(PetscSectionGetOffset(section, p, &off));
1419     for (d = 0; d < dof; d++, off++, f++) {
1420       const PetscReal *p;
1421       const PetscReal *w;
1422       PetscInt        Np, i;
1423 
1424       PetscCall(PetscDualSpaceGetFunctional(sp,off,&q));
1425       PetscCall(PetscQuadratureGetData(q,NULL,NULL,&Np,&p,&w));
1426       for (i = 0; i < Np * dim; i++) {
1427         points[offset + i] = p[i];
1428       }
1429       for (i = 0; i < Np * Nc; i++) {
1430         PetscCall(MatSetValue(imat, f, matoffset + i, w[i],INSERT_VALUES));
1431       }
1432       offset += Np * dim;
1433       matoffset += Np * Nc;
1434     }
1435   }
1436   PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF,intNodes));
1437   PetscCall(PetscQuadratureSetData(*intNodes,dim,0,numPoints,points,NULL));
1438   PetscCall(MatAssemblyBegin(imat, MAT_FINAL_ASSEMBLY));
1439   PetscCall(MatAssemblyEnd(imat, MAT_FINAL_ASSEMBLY));
1440   *intMat = imat;
1441   PetscFunctionReturn(0);
1442 }
1443 
1444 /*@
1445   PetscDualSpaceEqual - Determine if a dual space is equivalent
1446 
1447   Input Parameters:
1448 + A    - A PetscDualSpace object
1449 - B    - Another PetscDualSpace object
1450 
1451   Output Parameter:
1452 . equal - PETSC_TRUE if the dual spaces are equivalent
1453 
1454   Level: advanced
1455 
1456 .seealso: `PetscDualSpaceCreate()`
1457 @*/
1458 PetscErrorCode PetscDualSpaceEqual(PetscDualSpace A, PetscDualSpace B, PetscBool *equal)
1459 {
1460   PetscInt sizeA, sizeB, dimA, dimB;
1461   const PetscInt *dofA, *dofB;
1462   PetscQuadrature quadA, quadB;
1463   Mat matA, matB;
1464 
1465   PetscFunctionBegin;
1466   PetscValidHeaderSpecific(A,PETSCDUALSPACE_CLASSID,1);
1467   PetscValidHeaderSpecific(B,PETSCDUALSPACE_CLASSID,2);
1468   PetscValidBoolPointer(equal,3);
1469   *equal = PETSC_FALSE;
1470   PetscCall(PetscDualSpaceGetDimension(A, &sizeA));
1471   PetscCall(PetscDualSpaceGetDimension(B, &sizeB));
1472   if (sizeB != sizeA) {
1473     PetscFunctionReturn(0);
1474   }
1475   PetscCall(DMGetDimension(A->dm, &dimA));
1476   PetscCall(DMGetDimension(B->dm, &dimB));
1477   if (dimA != dimB) {
1478     PetscFunctionReturn(0);
1479   }
1480 
1481   PetscCall(PetscDualSpaceGetNumDof(A, &dofA));
1482   PetscCall(PetscDualSpaceGetNumDof(B, &dofB));
1483   for (PetscInt d=0; d<dimA; d++) {
1484     if (dofA[d] != dofB[d]) {
1485       PetscFunctionReturn(0);
1486     }
1487   }
1488 
1489   PetscCall(PetscDualSpaceGetInteriorData(A, &quadA, &matA));
1490   PetscCall(PetscDualSpaceGetInteriorData(B, &quadB, &matB));
1491   if (!quadA && !quadB) {
1492     *equal = PETSC_TRUE;
1493   } else if (quadA && quadB) {
1494     PetscCall(PetscQuadratureEqual(quadA, quadB, equal));
1495     if (*equal == PETSC_FALSE) PetscFunctionReturn(0);
1496     if (!matA && !matB) PetscFunctionReturn(0);
1497     if (matA && matB) PetscCall(MatEqual(matA, matB, equal));
1498     else *equal = PETSC_FALSE;
1499   }
1500   PetscFunctionReturn(0);
1501 }
1502 
1503 /*@C
1504   PetscDualSpaceApplyFVM - Apply a functional from the dual space basis to an input function by assuming a point evaluation functional at the cell centroid.
1505 
1506   Input Parameters:
1507 + sp    - The PetscDualSpace object
1508 . f     - The basis functional index
1509 . time  - The time
1510 . cgeom - A context with geometric information for this cell, we currently just use the centroid
1511 . Nc    - The number of components for the function
1512 . func  - The input function
1513 - ctx   - A context for the function
1514 
1515   Output Parameter:
1516 . value - The output value (scalar)
1517 
1518   Note: The calling sequence for the callback func is given by:
1519 
1520 $ func(PetscInt dim, PetscReal time, const PetscReal x[],
1521 $      PetscInt numComponents, PetscScalar values[], void *ctx)
1522 
1523 and the idea is to evaluate the functional as an integral
1524 
1525 $ n(f) = int dx n(x) . f(x)
1526 
1527 where both n and f have Nc components.
1528 
1529   Level: beginner
1530 
1531 .seealso: `PetscDualSpaceCreate()`
1532 @*/
1533 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)
1534 {
1535   DM               dm;
1536   PetscQuadrature  n;
1537   const PetscReal *points, *weights;
1538   PetscScalar     *val;
1539   PetscInt         dimEmbed, qNc, c, Nq, q;
1540 
1541   PetscFunctionBegin;
1542   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1543   PetscValidScalarPointer(value, 8);
1544   PetscCall(PetscDualSpaceGetDM(sp, &dm));
1545   PetscCall(DMGetCoordinateDim(dm, &dimEmbed));
1546   PetscCall(PetscDualSpaceGetFunctional(sp, f, &n));
1547   PetscCall(PetscQuadratureGetData(n, NULL, &qNc, &Nq, &points, &weights));
1548   PetscCheck(qNc == Nc,PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature components %" PetscInt_FMT " != function components %" PetscInt_FMT, qNc, Nc);
1549   PetscCall(DMGetWorkArray(dm, Nc, MPIU_SCALAR, &val));
1550   *value = 0.;
1551   for (q = 0; q < Nq; ++q) {
1552     PetscCall((*func)(dimEmbed, time, cgeom->centroid, Nc, val, ctx));
1553     for (c = 0; c < Nc; ++c) {
1554       *value += val[c]*weights[q*Nc+c];
1555     }
1556   }
1557   PetscCall(DMRestoreWorkArray(dm, Nc, MPIU_SCALAR, &val));
1558   PetscFunctionReturn(0);
1559 }
1560 
1561 /*@
1562   PetscDualSpaceGetHeightSubspace - Get the subset of the dual space basis that is supported on a mesh point of a
1563   given height.  This assumes that the reference cell is symmetric over points of this height.
1564 
1565   If the dual space is not defined on mesh points of the given height (e.g. if the space is discontinuous and
1566   pointwise values are not defined on the element boundaries), or if the implementation of PetscDualSpace does not
1567   support extracting subspaces, then NULL is returned.
1568 
1569   This does not increment the reference count on the returned dual space, and the user should not destroy it.
1570 
1571   Not collective
1572 
1573   Input Parameters:
1574 + sp - the PetscDualSpace object
1575 - height - the height of the mesh point for which the subspace is desired
1576 
1577   Output Parameter:
1578 . subsp - the subspace.  Note that the functionals in the subspace are with respect to the intrinsic geometry of the
1579   point, which will be of lesser dimension if height > 0.
1580 
1581   Level: advanced
1582 
1583 .seealso: `PetscSpaceGetHeightSubspace()`, `PetscDualSpace`
1584 @*/
1585 PetscErrorCode PetscDualSpaceGetHeightSubspace(PetscDualSpace sp, PetscInt height, PetscDualSpace *subsp)
1586 {
1587   PetscInt       depth = -1, cStart, cEnd;
1588   DM             dm;
1589 
1590   PetscFunctionBegin;
1591   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1592   PetscValidPointer(subsp,3);
1593   PetscCheck((sp->uniform),PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "A non-uniform dual space does not have a single dual space at each height");
1594   *subsp = NULL;
1595   dm = sp->dm;
1596   PetscCall(DMPlexGetDepth(dm, &depth));
1597   PetscCheck(height >= 0 && height <= depth,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid height");
1598   PetscCall(DMPlexGetHeightStratum(dm,0,&cStart,&cEnd));
1599   if (height == 0 && cEnd == cStart + 1) {
1600     *subsp = sp;
1601     PetscFunctionReturn(0);
1602   }
1603   if (!sp->heightSpaces) {
1604     PetscInt h;
1605     PetscCall(PetscCalloc1(depth+1, &(sp->heightSpaces)));
1606 
1607     for (h = 0; h <= depth; h++) {
1608       if (h == 0 && cEnd == cStart + 1) continue;
1609       if (sp->ops->createheightsubspace) PetscCall((*sp->ops->createheightsubspace)(sp,height,&(sp->heightSpaces[h])));
1610       else if (sp->pointSpaces) {
1611         PetscInt hStart, hEnd;
1612 
1613         PetscCall(DMPlexGetHeightStratum(dm,h,&hStart,&hEnd));
1614         if (hEnd > hStart) {
1615           const char *name;
1616 
1617           PetscCall(PetscObjectReference((PetscObject)(sp->pointSpaces[hStart])));
1618           if (sp->pointSpaces[hStart]) {
1619             PetscCall(PetscObjectGetName((PetscObject) sp,                     &name));
1620             PetscCall(PetscObjectSetName((PetscObject) sp->pointSpaces[hStart], name));
1621           }
1622           sp->heightSpaces[h] = sp->pointSpaces[hStart];
1623         }
1624       }
1625     }
1626   }
1627   *subsp = sp->heightSpaces[height];
1628   PetscFunctionReturn(0);
1629 }
1630 
1631 /*@
1632   PetscDualSpaceGetPointSubspace - Get the subset of the dual space basis that is supported on a particular mesh point.
1633 
1634   If the dual space is not defined on the mesh point (e.g. if the space is discontinuous and pointwise values are not
1635   defined on the element boundaries), or if the implementation of PetscDualSpace does not support extracting
1636   subspaces, then NULL is returned.
1637 
1638   This does not increment the reference count on the returned dual space, and the user should not destroy it.
1639 
1640   Not collective
1641 
1642   Input Parameters:
1643 + sp - the PetscDualSpace object
1644 - point - the point (in the dual space's DM) for which the subspace is desired
1645 
1646   Output Parameters:
1647   bdsp - the subspace.  Note that the functionals in the subspace are with respect to the intrinsic geometry of the
1648   point, which will be of lesser dimension if height > 0.
1649 
1650   Level: advanced
1651 
1652 .seealso: `PetscDualSpace`
1653 @*/
1654 PetscErrorCode PetscDualSpaceGetPointSubspace(PetscDualSpace sp, PetscInt point, PetscDualSpace *bdsp)
1655 {
1656   PetscInt       pStart = 0, pEnd = 0, cStart, cEnd;
1657   DM             dm;
1658 
1659   PetscFunctionBegin;
1660   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1661   PetscValidPointer(bdsp,3);
1662   *bdsp = NULL;
1663   dm = sp->dm;
1664   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
1665   PetscCheck(point >= pStart && point <= pEnd,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid point");
1666   PetscCall(DMPlexGetHeightStratum(dm,0,&cStart,&cEnd));
1667   if (point == cStart && cEnd == cStart + 1) { /* the dual space is only equivalent to the dual space on a cell if the reference mesh has just one cell */
1668     *bdsp = sp;
1669     PetscFunctionReturn(0);
1670   }
1671   if (!sp->pointSpaces) {
1672     PetscInt p;
1673     PetscCall(PetscCalloc1(pEnd - pStart, &(sp->pointSpaces)));
1674 
1675     for (p = 0; p < pEnd - pStart; p++) {
1676       if (p + pStart == cStart && cEnd == cStart + 1) continue;
1677       if (sp->ops->createpointsubspace) PetscCall((*sp->ops->createpointsubspace)(sp,p+pStart,&(sp->pointSpaces[p])));
1678       else if (sp->heightSpaces || sp->ops->createheightsubspace) {
1679         PetscInt dim, depth, height;
1680         DMLabel  label;
1681 
1682         PetscCall(DMPlexGetDepth(dm,&dim));
1683         PetscCall(DMPlexGetDepthLabel(dm,&label));
1684         PetscCall(DMLabelGetValue(label,p+pStart,&depth));
1685         height = dim - depth;
1686         PetscCall(PetscDualSpaceGetHeightSubspace(sp, height, &(sp->pointSpaces[p])));
1687         PetscCall(PetscObjectReference((PetscObject)sp->pointSpaces[p]));
1688       }
1689     }
1690   }
1691   *bdsp = sp->pointSpaces[point - pStart];
1692   PetscFunctionReturn(0);
1693 }
1694 
1695 /*@C
1696   PetscDualSpaceGetSymmetries - Returns a description of the symmetries of this basis
1697 
1698   Not collective
1699 
1700   Input Parameter:
1701 . sp - the PetscDualSpace object
1702 
1703   Output Parameters:
1704 + perms - Permutations of the interior degrees of freedom, parameterized by the point orientation
1705 - flips - Sign reversal of the interior degrees of freedom, parameterized by the point orientation
1706 
1707   Note: The permutation and flip arrays are organized in the following way
1708 $ perms[p][ornt][dof # on point] = new local dof #
1709 $ flips[p][ornt][dof # on point] = reversal or not
1710 
1711   Level: developer
1712 
1713 @*/
1714 PetscErrorCode PetscDualSpaceGetSymmetries(PetscDualSpace sp, const PetscInt ****perms, const PetscScalar ****flips)
1715 {
1716   PetscFunctionBegin;
1717   PetscValidHeaderSpecific(sp,PETSCDUALSPACE_CLASSID,1);
1718   if (perms) {PetscValidPointer(perms,2); *perms = NULL;}
1719   if (flips) {PetscValidPointer(flips,3); *flips = NULL;}
1720   if (sp->ops->getsymmetries) PetscCall((sp->ops->getsymmetries)(sp,perms,flips));
1721   PetscFunctionReturn(0);
1722 }
1723 
1724 /*@
1725   PetscDualSpaceGetFormDegree - Get the form degree k for the k-form the describes the pushforwards/pullbacks of this
1726   dual space's functionals.
1727 
1728   Input Parameter:
1729 . dsp - The PetscDualSpace
1730 
1731   Output Parameter:
1732 . k   - The *signed* degree k of the k.  If k >= 0, this means that the degrees of freedom are k-forms, and are stored
1733         in lexicographic order according to the basis of k-forms constructed from the wedge product of 1-forms.  So for example,
1734         the 1-form basis in 3-D is (dx, dy, dz), and the 2-form basis in 3-D is (dx wedge dy, dx wedge dz, dy wedge dz).
1735         If k < 0, this means that the degrees transform as k-forms, but are stored as (N-k) forms according to the
1736         Hodge star map.  So for example if k = -2 and N = 3, this means that the degrees of freedom transform as 2-forms
1737         but are stored as 1-forms.
1738 
1739   Level: developer
1740 
1741 .seealso: `PetscDTAltV`, `PetscDualSpacePullback()`, `PetscDualSpacePushforward()`, `PetscDualSpaceTransform()`, `PetscDualSpaceTransformType`
1742 @*/
1743 PetscErrorCode PetscDualSpaceGetFormDegree(PetscDualSpace dsp, PetscInt *k)
1744 {
1745   PetscFunctionBeginHot;
1746   PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1);
1747   PetscValidIntPointer(k, 2);
1748   *k = dsp->k;
1749   PetscFunctionReturn(0);
1750 }
1751 
1752 /*@
1753   PetscDualSpaceSetFormDegree - Set the form degree k for the k-form the describes the pushforwards/pullbacks of this
1754   dual space's functionals.
1755 
1756   Input Parameters:
1757 + dsp - The PetscDualSpace
1758 - k   - The *signed* degree k of the k.  If k >= 0, this means that the degrees of freedom are k-forms, and are stored
1759         in lexicographic order according to the basis of k-forms constructed from the wedge product of 1-forms.  So for example,
1760         the 1-form basis in 3-D is (dx, dy, dz), and the 2-form basis in 3-D is (dx wedge dy, dx wedge dz, dy wedge dz).
1761         If k < 0, this means that the degrees transform as k-forms, but are stored as (N-k) forms according to the
1762         Hodge star map.  So for example if k = -2 and N = 3, this means that the degrees of freedom transform as 2-forms
1763         but are stored as 1-forms.
1764 
1765   Level: developer
1766 
1767 .seealso: `PetscDTAltV`, `PetscDualSpacePullback()`, `PetscDualSpacePushforward()`, `PetscDualSpaceTransform()`, `PetscDualSpaceTransformType`
1768 @*/
1769 PetscErrorCode PetscDualSpaceSetFormDegree(PetscDualSpace dsp, PetscInt k)
1770 {
1771   PetscInt dim;
1772 
1773   PetscFunctionBeginHot;
1774   PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1);
1775   PetscCheck(!dsp->setupcalled,PetscObjectComm((PetscObject)dsp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change number of components after dualspace is set up");
1776   dim = dsp->dm->dim;
1777   PetscCheck(k >= -dim && k <= dim,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unsupported %" PetscInt_FMT "-form on %" PetscInt_FMT "-dimensional reference cell", PetscAbsInt(k), dim);
1778   dsp->k = k;
1779   PetscFunctionReturn(0);
1780 }
1781 
1782 /*@
1783   PetscDualSpaceGetDeRahm - Get the k-simplex associated with the functionals in this dual space
1784 
1785   Input Parameter:
1786 . dsp - The PetscDualSpace
1787 
1788   Output Parameter:
1789 . k   - The simplex dimension
1790 
1791   Level: developer
1792 
1793   Note: Currently supported values are
1794 $ 0: These are H_1 methods that only transform coordinates
1795 $ 1: These are Hcurl methods that transform functions using the covariant Piola transform (COVARIANT_PIOLA_TRANSFORM)
1796 $ 2: These are the same as 1
1797 $ 3: These are Hdiv methods that transform functions using the contravariant Piola transform (CONTRAVARIANT_PIOLA_TRANSFORM)
1798 
1799 .seealso: `PetscDualSpacePullback()`, `PetscDualSpacePushforward()`, `PetscDualSpaceTransform()`, `PetscDualSpaceTransformType`
1800 @*/
1801 PetscErrorCode PetscDualSpaceGetDeRahm(PetscDualSpace dsp, PetscInt *k)
1802 {
1803   PetscInt dim;
1804 
1805   PetscFunctionBeginHot;
1806   PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1);
1807   PetscValidIntPointer(k, 2);
1808   dim = dsp->dm->dim;
1809   if (!dsp->k) *k = IDENTITY_TRANSFORM;
1810   else if (dsp->k == 1) *k = COVARIANT_PIOLA_TRANSFORM;
1811   else if (dsp->k == -(dim - 1)) *k = CONTRAVARIANT_PIOLA_TRANSFORM;
1812   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unsupported transformation");
1813   PetscFunctionReturn(0);
1814 }
1815 
1816 /*@C
1817   PetscDualSpaceTransform - Transform the function values
1818 
1819   Input Parameters:
1820 + dsp       - The PetscDualSpace
1821 . trans     - The type of transform
1822 . isInverse - Flag to invert the transform
1823 . fegeom    - The cell geometry
1824 . Nv        - The number of function samples
1825 . Nc        - The number of function components
1826 - vals      - The function values
1827 
1828   Output Parameter:
1829 . vals      - The transformed function values
1830 
1831   Level: intermediate
1832 
1833   Note: This only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension.
1834 
1835 .seealso: `PetscDualSpaceTransformGradient()`, `PetscDualSpaceTransformHessian()`, `PetscDualSpacePullback()`, `PetscDualSpacePushforward()`, `PetscDualSpaceTransformType`
1836 @*/
1837 PetscErrorCode PetscDualSpaceTransform(PetscDualSpace dsp, PetscDualSpaceTransformType trans, PetscBool isInverse, PetscFEGeom *fegeom, PetscInt Nv, PetscInt Nc, PetscScalar vals[])
1838 {
1839   PetscReal Jstar[9] = {0};
1840   PetscInt dim, v, c, Nk;
1841 
1842   PetscFunctionBeginHot;
1843   PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1);
1844   PetscValidPointer(fegeom, 4);
1845   PetscValidScalarPointer(vals, 7);
1846   /* TODO: not handling dimEmbed != dim right now */
1847   dim = dsp->dm->dim;
1848   /* No change needed for 0-forms */
1849   if (!dsp->k) PetscFunctionReturn(0);
1850   PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(dsp->k), &Nk));
1851   /* TODO: use fegeom->isAffine */
1852   PetscCall(PetscDTAltVPullbackMatrix(dim, dim, isInverse ? fegeom->J : fegeom->invJ, dsp->k, Jstar));
1853   for (v = 0; v < Nv; ++v) {
1854     switch (Nk) {
1855     case 1:
1856       for (c = 0; c < Nc; c++) vals[v*Nc + c] *= Jstar[0];
1857       break;
1858     case 2:
1859       for (c = 0; c < Nc; c += 2) DMPlex_Mult2DReal_Internal(Jstar, 1, &vals[v*Nc + c], &vals[v*Nc + c]);
1860       break;
1861     case 3:
1862       for (c = 0; c < Nc; c += 3) DMPlex_Mult3DReal_Internal(Jstar, 1, &vals[v*Nc + c], &vals[v*Nc + c]);
1863       break;
1864     default: SETERRQ(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported form size %" PetscInt_FMT " for transformation", Nk);
1865     }
1866   }
1867   PetscFunctionReturn(0);
1868 }
1869 
1870 /*@C
1871   PetscDualSpaceTransformGradient - Transform the function gradient values
1872 
1873   Input Parameters:
1874 + dsp       - The PetscDualSpace
1875 . trans     - The type of transform
1876 . isInverse - Flag to invert the transform
1877 . fegeom    - The cell geometry
1878 . Nv        - The number of function gradient samples
1879 . Nc        - The number of function components
1880 - vals      - The function gradient values
1881 
1882   Output Parameter:
1883 . vals      - The transformed function gradient values
1884 
1885   Level: intermediate
1886 
1887   Note: This only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension.
1888 
1889 .seealso: `PetscDualSpaceTransform()`, `PetscDualSpacePullback()`, `PetscDualSpacePushforward()`, `PetscDualSpaceTransformType`
1890 @*/
1891 PetscErrorCode PetscDualSpaceTransformGradient(PetscDualSpace dsp, PetscDualSpaceTransformType trans, PetscBool isInverse, PetscFEGeom *fegeom, PetscInt Nv, PetscInt Nc, PetscScalar vals[])
1892 {
1893   const PetscInt dim = dsp->dm->dim, dE = fegeom->dimEmbed;
1894   PetscInt       v, c, d;
1895 
1896   PetscFunctionBeginHot;
1897   PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1);
1898   PetscValidPointer(fegeom, 4);
1899   PetscValidScalarPointer(vals, 7);
1900 #ifdef PETSC_USE_DEBUG
1901   PetscCheck(dE > 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid embedding dimension %" PetscInt_FMT, dE);
1902 #endif
1903   /* Transform gradient */
1904   if (dim == dE) {
1905     for (v = 0; v < Nv; ++v) {
1906       for (c = 0; c < Nc; ++c) {
1907         switch (dim)
1908         {
1909           case 1: vals[(v*Nc+c)*dim] *= fegeom->invJ[0];break;
1910           case 2: DMPlex_MultTranspose2DReal_Internal(fegeom->invJ, 1, &vals[(v*Nc+c)*dim], &vals[(v*Nc+c)*dim]);break;
1911           case 3: DMPlex_MultTranspose3DReal_Internal(fegeom->invJ, 1, &vals[(v*Nc+c)*dim], &vals[(v*Nc+c)*dim]);break;
1912           default: SETERRQ(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %" PetscInt_FMT " for transformation", dim);
1913         }
1914       }
1915     }
1916   } else {
1917     for (v = 0; v < Nv; ++v) {
1918       for (c = 0; c < Nc; ++c) {
1919         DMPlex_MultTransposeReal_Internal(fegeom->invJ, dim, dE, 1, &vals[(v*Nc+c)*dE], &vals[(v*Nc+c)*dE]);
1920       }
1921     }
1922   }
1923   /* Assume its a vector, otherwise assume its a bunch of scalars */
1924   if (Nc == 1 || Nc != dim) PetscFunctionReturn(0);
1925   switch (trans) {
1926     case IDENTITY_TRANSFORM: break;
1927     case COVARIANT_PIOLA_TRANSFORM: /* Covariant Piola mapping $\sigma^*(F) = J^{-T} F \circ \phi^{-1)$ */
1928     if (isInverse) {
1929       for (v = 0; v < Nv; ++v) {
1930         for (d = 0; d < dim; ++d) {
1931           switch (dim)
1932           {
1933             case 2: DMPlex_MultTranspose2DReal_Internal(fegeom->J, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
1934             case 3: DMPlex_MultTranspose3DReal_Internal(fegeom->J, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
1935             default: SETERRQ(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %" PetscInt_FMT " for transformation", dim);
1936           }
1937         }
1938       }
1939     } else {
1940       for (v = 0; v < Nv; ++v) {
1941         for (d = 0; d < dim; ++d) {
1942           switch (dim)
1943           {
1944             case 2: DMPlex_MultTranspose2DReal_Internal(fegeom->invJ, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
1945             case 3: DMPlex_MultTranspose3DReal_Internal(fegeom->invJ, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
1946             default: SETERRQ(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %" PetscInt_FMT " for transformation", dim);
1947           }
1948         }
1949       }
1950     }
1951     break;
1952     case CONTRAVARIANT_PIOLA_TRANSFORM: /* Contravariant Piola mapping $\sigma^*(F) = \frac{1}{|\det J|} J F \circ \phi^{-1}$ */
1953     if (isInverse) {
1954       for (v = 0; v < Nv; ++v) {
1955         for (d = 0; d < dim; ++d) {
1956           switch (dim)
1957           {
1958             case 2: DMPlex_Mult2DReal_Internal(fegeom->invJ, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
1959             case 3: DMPlex_Mult3DReal_Internal(fegeom->invJ, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
1960             default: SETERRQ(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %" PetscInt_FMT " for transformation", dim);
1961           }
1962           for (c = 0; c < Nc; ++c) vals[(v*Nc+c)*dim+d] *= fegeom->detJ[0];
1963         }
1964       }
1965     } else {
1966       for (v = 0; v < Nv; ++v) {
1967         for (d = 0; d < dim; ++d) {
1968           switch (dim)
1969           {
1970             case 2: DMPlex_Mult2DReal_Internal(fegeom->J, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
1971             case 3: DMPlex_Mult3DReal_Internal(fegeom->J, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
1972             default: SETERRQ(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %" PetscInt_FMT " for transformation", dim);
1973           }
1974           for (c = 0; c < Nc; ++c) vals[(v*Nc+c)*dim+d] /= fegeom->detJ[0];
1975         }
1976       }
1977     }
1978     break;
1979   }
1980   PetscFunctionReturn(0);
1981 }
1982 
1983 /*@C
1984   PetscDualSpaceTransformHessian - Transform the function Hessian values
1985 
1986   Input Parameters:
1987 + dsp       - The PetscDualSpace
1988 . trans     - The type of transform
1989 . isInverse - Flag to invert the transform
1990 . fegeom    - The cell geometry
1991 . Nv        - The number of function Hessian samples
1992 . Nc        - The number of function components
1993 - vals      - The function gradient values
1994 
1995   Output Parameter:
1996 . vals      - The transformed function Hessian values
1997 
1998   Level: intermediate
1999 
2000   Note: This only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension.
2001 
2002 .seealso: `PetscDualSpaceTransform()`, `PetscDualSpacePullback()`, `PetscDualSpacePushforward()`, `PetscDualSpaceTransformType`
2003 @*/
2004 PetscErrorCode PetscDualSpaceTransformHessian(PetscDualSpace dsp, PetscDualSpaceTransformType trans, PetscBool isInverse, PetscFEGeom *fegeom, PetscInt Nv, PetscInt Nc, PetscScalar vals[])
2005 {
2006   const PetscInt dim = dsp->dm->dim, dE = fegeom->dimEmbed;
2007   PetscInt       v, c;
2008 
2009   PetscFunctionBeginHot;
2010   PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1);
2011   PetscValidPointer(fegeom, 4);
2012   PetscValidScalarPointer(vals, 7);
2013 #ifdef PETSC_USE_DEBUG
2014   PetscCheck(dE > 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid embedding dimension %" PetscInt_FMT, dE);
2015 #endif
2016   /* Transform Hessian: J^{-T}_{ik} J^{-T}_{jl} H(f)_{kl} = J^{-T}_{ik} H(f)_{kl} J^{-1}_{lj} */
2017   if (dim == dE) {
2018     for (v = 0; v < Nv; ++v) {
2019       for (c = 0; c < Nc; ++c) {
2020         switch (dim)
2021         {
2022           case 1: vals[(v*Nc+c)*dim*dim] *= PetscSqr(fegeom->invJ[0]);break;
2023           case 2: DMPlex_PTAP2DReal_Internal(fegeom->invJ, &vals[(v*Nc+c)*dim*dim], &vals[(v*Nc+c)*dim*dim]);break;
2024           case 3: DMPlex_PTAP3DReal_Internal(fegeom->invJ, &vals[(v*Nc+c)*dim*dim], &vals[(v*Nc+c)*dim*dim]);break;
2025           default: SETERRQ(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %" PetscInt_FMT " for transformation", dim);
2026         }
2027       }
2028     }
2029   } else {
2030     for (v = 0; v < Nv; ++v) {
2031       for (c = 0; c < Nc; ++c) {
2032         DMPlex_PTAPReal_Internal(fegeom->invJ, dim, dE, &vals[(v*Nc+c)*dE*dE], &vals[(v*Nc+c)*dE*dE]);
2033       }
2034     }
2035   }
2036   /* Assume its a vector, otherwise assume its a bunch of scalars */
2037   if (Nc == 1 || Nc != dim) PetscFunctionReturn(0);
2038   switch (trans) {
2039     case IDENTITY_TRANSFORM: break;
2040     case COVARIANT_PIOLA_TRANSFORM: /* Covariant Piola mapping $\sigma^*(F) = J^{-T} F \circ \phi^{-1)$ */
2041     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Piola mapping for Hessians not yet supported");
2042     case CONTRAVARIANT_PIOLA_TRANSFORM: /* Contravariant Piola mapping $\sigma^*(F) = \frac{1}{|\det J|} J F \circ \phi^{-1}$ */
2043     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Piola mapping for Hessians not yet supported");
2044   }
2045   PetscFunctionReturn(0);
2046 }
2047 
2048 /*@C
2049   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.
2050 
2051   Input Parameters:
2052 + dsp        - The PetscDualSpace
2053 . fegeom     - The geometry for this cell
2054 . Nq         - The number of function samples
2055 . Nc         - The number of function components
2056 - pointEval  - The function values
2057 
2058   Output Parameter:
2059 . pointEval  - The transformed function values
2060 
2061   Level: advanced
2062 
2063   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.
2064 
2065   Note: This only handles tranformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension.
2066 
2067 .seealso: `PetscDualSpacePushforward()`, `PetscDualSpaceTransform()`, `PetscDualSpaceGetDeRahm()`
2068 @*/
2069 PetscErrorCode PetscDualSpacePullback(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, PetscScalar pointEval[])
2070 {
2071   PetscDualSpaceTransformType trans;
2072   PetscInt                    k;
2073 
2074   PetscFunctionBeginHot;
2075   PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1);
2076   PetscValidPointer(fegeom, 2);
2077   PetscValidScalarPointer(pointEval, 5);
2078   /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k.
2079      This determines their transformation properties. */
2080   PetscCall(PetscDualSpaceGetDeRahm(dsp, &k));
2081   switch (k)
2082   {
2083     case 0: /* H^1 point evaluations */
2084     trans = IDENTITY_TRANSFORM;break;
2085     case 1: /* Hcurl preserves tangential edge traces  */
2086     trans = COVARIANT_PIOLA_TRANSFORM;break;
2087     case 2:
2088     case 3: /* Hdiv preserve normal traces */
2089     trans = CONTRAVARIANT_PIOLA_TRANSFORM;break;
2090     default: SETERRQ(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %" PetscInt_FMT " for transformation", k);
2091   }
2092   PetscCall(PetscDualSpaceTransform(dsp, trans, PETSC_TRUE, fegeom, Nq, Nc, pointEval));
2093   PetscFunctionReturn(0);
2094 }
2095 
2096 /*@C
2097   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.
2098 
2099   Input Parameters:
2100 + dsp        - The PetscDualSpace
2101 . fegeom     - The geometry for this cell
2102 . Nq         - The number of function samples
2103 . Nc         - The number of function components
2104 - pointEval  - The function values
2105 
2106   Output Parameter:
2107 . pointEval  - The transformed function values
2108 
2109   Level: advanced
2110 
2111   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.
2112 
2113   Note: This only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension.
2114 
2115 .seealso: `PetscDualSpacePullback()`, `PetscDualSpaceTransform()`, `PetscDualSpaceGetDeRahm()`
2116 @*/
2117 PetscErrorCode PetscDualSpacePushforward(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, PetscScalar pointEval[])
2118 {
2119   PetscDualSpaceTransformType trans;
2120   PetscInt                    k;
2121 
2122   PetscFunctionBeginHot;
2123   PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1);
2124   PetscValidPointer(fegeom, 2);
2125   PetscValidScalarPointer(pointEval, 5);
2126   /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k.
2127      This determines their transformation properties. */
2128   PetscCall(PetscDualSpaceGetDeRahm(dsp, &k));
2129   switch (k)
2130   {
2131     case 0: /* H^1 point evaluations */
2132     trans = IDENTITY_TRANSFORM;break;
2133     case 1: /* Hcurl preserves tangential edge traces  */
2134     trans = COVARIANT_PIOLA_TRANSFORM;break;
2135     case 2:
2136     case 3: /* Hdiv preserve normal traces */
2137     trans = CONTRAVARIANT_PIOLA_TRANSFORM;break;
2138     default: SETERRQ(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %" PetscInt_FMT " for transformation", k);
2139   }
2140   PetscCall(PetscDualSpaceTransform(dsp, trans, PETSC_FALSE, fegeom, Nq, Nc, pointEval));
2141   PetscFunctionReturn(0);
2142 }
2143 
2144 /*@C
2145   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.
2146 
2147   Input Parameters:
2148 + dsp        - The PetscDualSpace
2149 . fegeom     - The geometry for this cell
2150 . Nq         - The number of function gradient samples
2151 . Nc         - The number of function components
2152 - pointEval  - The function gradient values
2153 
2154   Output Parameter:
2155 . pointEval  - The transformed function gradient values
2156 
2157   Level: advanced
2158 
2159   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.
2160 
2161   Note: This only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension.
2162 
2163 .seealso: `PetscDualSpacePushforward()`, `PPetscDualSpacePullback()`, `PetscDualSpaceTransform()`, `PetscDualSpaceGetDeRahm()`
2164 @*/
2165 PetscErrorCode PetscDualSpacePushforwardGradient(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, PetscScalar pointEval[])
2166 {
2167   PetscDualSpaceTransformType trans;
2168   PetscInt                    k;
2169 
2170   PetscFunctionBeginHot;
2171   PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1);
2172   PetscValidPointer(fegeom, 2);
2173   PetscValidScalarPointer(pointEval, 5);
2174   /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k.
2175      This determines their transformation properties. */
2176   PetscCall(PetscDualSpaceGetDeRahm(dsp, &k));
2177   switch (k)
2178   {
2179     case 0: /* H^1 point evaluations */
2180     trans = IDENTITY_TRANSFORM;break;
2181     case 1: /* Hcurl preserves tangential edge traces  */
2182     trans = COVARIANT_PIOLA_TRANSFORM;break;
2183     case 2:
2184     case 3: /* Hdiv preserve normal traces */
2185     trans = CONTRAVARIANT_PIOLA_TRANSFORM;break;
2186     default: SETERRQ(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %" PetscInt_FMT " for transformation", k);
2187   }
2188   PetscCall(PetscDualSpaceTransformGradient(dsp, trans, PETSC_FALSE, fegeom, Nq, Nc, pointEval));
2189   PetscFunctionReturn(0);
2190 }
2191 
2192 /*@C
2193   PetscDualSpacePushforwardHessian - Transform the given function Hessian 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.
2194 
2195   Input Parameters:
2196 + dsp        - The PetscDualSpace
2197 . fegeom     - The geometry for this cell
2198 . Nq         - The number of function Hessian samples
2199 . Nc         - The number of function components
2200 - pointEval  - The function gradient values
2201 
2202   Output Parameter:
2203 . pointEval  - The transformed function Hessian values
2204 
2205   Level: advanced
2206 
2207   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.
2208 
2209   Note: This only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension.
2210 
2211 .seealso: `PetscDualSpacePushforward()`, `PPetscDualSpacePullback()`, `PetscDualSpaceTransform()`, `PetscDualSpaceGetDeRahm()`
2212 @*/
2213 PetscErrorCode PetscDualSpacePushforwardHessian(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, PetscScalar pointEval[])
2214 {
2215   PetscDualSpaceTransformType trans;
2216   PetscInt                    k;
2217 
2218   PetscFunctionBeginHot;
2219   PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1);
2220   PetscValidPointer(fegeom, 2);
2221   PetscValidScalarPointer(pointEval, 5);
2222   /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k.
2223      This determines their transformation properties. */
2224   PetscCall(PetscDualSpaceGetDeRahm(dsp, &k));
2225   switch (k)
2226   {
2227     case 0: /* H^1 point evaluations */
2228     trans = IDENTITY_TRANSFORM;break;
2229     case 1: /* Hcurl preserves tangential edge traces  */
2230     trans = COVARIANT_PIOLA_TRANSFORM;break;
2231     case 2:
2232     case 3: /* Hdiv preserve normal traces */
2233     trans = CONTRAVARIANT_PIOLA_TRANSFORM;break;
2234     default: SETERRQ(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %" PetscInt_FMT " for transformation", k);
2235   }
2236   PetscCall(PetscDualSpaceTransformHessian(dsp, trans, PETSC_FALSE, fegeom, Nq, Nc, pointEval));
2237   PetscFunctionReturn(0);
2238 }
2239