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