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