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