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