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