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