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