xref: /petsc/src/dm/dt/dualspace/interface/dualspace.c (revision 4bee2e389ac4efdf19d1420f70098a911b40ccd1)
1 #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
2 #include <petscdmplex.h>
3 
4 PetscClassId PETSCDUALSPACE_CLASSID = 0;
5 
6 PetscFunctionList PetscDualSpaceList              = NULL;
7 PetscBool         PetscDualSpaceRegisterAllCalled = PETSC_FALSE;
8 
9 const char *const PetscDualSpaceReferenceCells[] = {"SIMPLEX", "TENSOR", "PetscDualSpaceReferenceCell", "PETSCDUALSPACE_REFCELL_",0};
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}, {2,0}.
15 
16   Input Parameters:
17 + len - The length of the tuple
18 . max - The maximum sum
19 - tup - A tuple of length len+1: tup[len] > 0 indicates a stopping condition
20 
21   Output Parameter:
22 . tup - A tuple of len integers whos sum is at most 'max'
23 
24   Level: developer
25 
26 .seealso: PetscDualSpaceTensorPointLexicographic_Internal()
27 */
28 PetscErrorCode PetscDualSpaceLatticePointLexicographic_Internal(PetscInt len, PetscInt max, PetscInt tup[])
29 {
30   PetscFunctionBegin;
31   while (len--) {
32     max -= tup[len];
33     if (!max) {
34       tup[len] = 0;
35       break;
36     }
37   }
38   tup[++len]++;
39   PetscFunctionReturn(0);
40 }
41 
42 /*
43   PetscDualSpaceTensorPointLexicographic_Internal - Returns all tuples of size 'len' with nonnegative integers that are all less than or equal to 'max'.
44                                                     Ordering is lexicographic with lowest index as least significant in ordering.
45                                                     e.g. for len == 2 and max == 2, this will return, in order, {0,0}, {1,0}, {2,0}, {0,1}, {1,1}, {2,1}, {0,2}, {1,2}, {2,2}.
46 
47   Input Parameters:
48 + len - The length of the tuple
49 . max - The maximum value
50 - tup - A tuple of length len+1: tup[len] > 0 indicates a stopping condition
51 
52   Output Parameter:
53 . tup - A tuple of len integers whos sum is at most 'max'
54 
55   Level: developer
56 
57 .seealso: PetscDualSpaceLatticePointLexicographic_Internal()
58 */
59 PetscErrorCode PetscDualSpaceTensorPointLexicographic_Internal(PetscInt len, PetscInt max, PetscInt tup[])
60 {
61   PetscInt       i;
62 
63   PetscFunctionBegin;
64   for (i = 0; i < len; i++) {
65     if (tup[i] < max) {
66       break;
67     } else {
68       tup[i] = 0;
69     }
70   }
71   tup[i]++;
72   PetscFunctionReturn(0);
73 }
74 
75 /*@C
76   PetscDualSpaceRegister - Adds a new PetscDualSpace implementation
77 
78   Not Collective
79 
80   Input Parameters:
81 + name        - The name of a new user-defined creation routine
82 - create_func - The creation routine itself
83 
84   Notes:
85   PetscDualSpaceRegister() may be called multiple times to add several user-defined PetscDualSpaces
86 
87   Sample usage:
88 .vb
89     PetscDualSpaceRegister("my_space", MyPetscDualSpaceCreate);
90 .ve
91 
92   Then, your PetscDualSpace type can be chosen with the procedural interface via
93 .vb
94     PetscDualSpaceCreate(MPI_Comm, PetscDualSpace *);
95     PetscDualSpaceSetType(PetscDualSpace, "my_dual_space");
96 .ve
97    or at runtime via the option
98 .vb
99     -petscdualspace_type my_dual_space
100 .ve
101 
102   Level: advanced
103 
104 .keywords: PetscDualSpace, register
105 .seealso: PetscDualSpaceRegisterAll(), PetscDualSpaceRegisterDestroy()
106 
107 @*/
108 PetscErrorCode PetscDualSpaceRegister(const char sname[], PetscErrorCode (*function)(PetscDualSpace))
109 {
110   PetscErrorCode ierr;
111 
112   PetscFunctionBegin;
113   ierr = PetscFunctionListAdd(&PetscDualSpaceList, sname, function);CHKERRQ(ierr);
114   PetscFunctionReturn(0);
115 }
116 
117 /*@C
118   PetscDualSpaceSetType - Builds a particular PetscDualSpace
119 
120   Collective on PetscDualSpace
121 
122   Input Parameters:
123 + sp   - The PetscDualSpace object
124 - name - The kind of space
125 
126   Options Database Key:
127 . -petscdualspace_type <type> - Sets the PetscDualSpace type; use -help for a list of available types
128 
129   Level: intermediate
130 
131 .keywords: PetscDualSpace, set, type
132 .seealso: PetscDualSpaceGetType(), PetscDualSpaceCreate()
133 @*/
134 PetscErrorCode PetscDualSpaceSetType(PetscDualSpace sp, PetscDualSpaceType name)
135 {
136   PetscErrorCode (*r)(PetscDualSpace);
137   PetscBool      match;
138   PetscErrorCode ierr;
139 
140   PetscFunctionBegin;
141   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
142   ierr = PetscObjectTypeCompare((PetscObject) sp, name, &match);CHKERRQ(ierr);
143   if (match) PetscFunctionReturn(0);
144 
145   if (!PetscDualSpaceRegisterAllCalled) {ierr = PetscDualSpaceRegisterAll();CHKERRQ(ierr);}
146   ierr = PetscFunctionListFind(PetscDualSpaceList, name, &r);CHKERRQ(ierr);
147   if (!r) SETERRQ1(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscDualSpace type: %s", name);
148 
149   if (sp->ops->destroy) {
150     ierr             = (*sp->ops->destroy)(sp);CHKERRQ(ierr);
151     sp->ops->destroy = NULL;
152   }
153   ierr = (*r)(sp);CHKERRQ(ierr);
154   ierr = PetscObjectChangeTypeName((PetscObject) sp, name);CHKERRQ(ierr);
155   PetscFunctionReturn(0);
156 }
157 
158 /*@C
159   PetscDualSpaceGetType - Gets the PetscDualSpace type name (as a string) from the object.
160 
161   Not Collective
162 
163   Input Parameter:
164 . sp  - The PetscDualSpace
165 
166   Output Parameter:
167 . name - The PetscDualSpace type name
168 
169   Level: intermediate
170 
171 .keywords: PetscDualSpace, get, type, name
172 .seealso: PetscDualSpaceSetType(), PetscDualSpaceCreate()
173 @*/
174 PetscErrorCode PetscDualSpaceGetType(PetscDualSpace sp, PetscDualSpaceType *name)
175 {
176   PetscErrorCode ierr;
177 
178   PetscFunctionBegin;
179   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
180   PetscValidPointer(name, 2);
181   if (!PetscDualSpaceRegisterAllCalled) {
182     ierr = PetscDualSpaceRegisterAll();CHKERRQ(ierr);
183   }
184   *name = ((PetscObject) sp)->type_name;
185   PetscFunctionReturn(0);
186 }
187 
188 static PetscErrorCode PetscDualSpaceView_ASCII(PetscDualSpace sp, PetscViewer v)
189 {
190   PetscViewerFormat format;
191   PetscInt          pdim, f;
192   PetscErrorCode    ierr;
193 
194   PetscFunctionBegin;
195   ierr = PetscDualSpaceGetDimension(sp, &pdim);CHKERRQ(ierr);
196   ierr = PetscObjectPrintClassNamePrefixType((PetscObject) sp, v);CHKERRQ(ierr);
197   ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr);
198   ierr = PetscViewerASCIIPrintf(v, "Dual space with %D components, size %D\n", sp->Nc, pdim);CHKERRQ(ierr);
199   if (sp->ops->view) {ierr = (*sp->ops->view)(sp, v);CHKERRQ(ierr);}
200   ierr = PetscViewerGetFormat(v, &format);CHKERRQ(ierr);
201   if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
202     ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr);
203     for (f = 0; f < pdim; ++f) {
204       ierr = PetscViewerASCIIPrintf(v, "Dual basis vector %D\n", f);CHKERRQ(ierr);
205       ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr);
206       ierr = PetscQuadratureView(sp->functional[f], v);CHKERRQ(ierr);
207       ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr);
208     }
209     ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr);
210   }
211   ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr);
212   PetscFunctionReturn(0);
213 }
214 
215 /*@
216   PetscDualSpaceView - Views a PetscDualSpace
217 
218   Collective on PetscDualSpace
219 
220   Input Parameter:
221 + sp - the PetscDualSpace object to view
222 - v  - the viewer
223 
224   Level: developer
225 
226 .seealso PetscDualSpaceDestroy()
227 @*/
228 PetscErrorCode PetscDualSpaceView(PetscDualSpace sp, PetscViewer v)
229 {
230   PetscBool      iascii;
231   PetscErrorCode ierr;
232 
233   PetscFunctionBegin;
234   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
235   if (v) PetscValidHeaderSpecific(v, PETSC_VIEWER_CLASSID, 2);
236   if (!v) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) sp), &v);CHKERRQ(ierr);}
237   ierr = PetscObjectTypeCompare((PetscObject) v, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
238   if (iascii) {ierr = PetscDualSpaceView_ASCII(sp, v);CHKERRQ(ierr);}
239   PetscFunctionReturn(0);
240 }
241 
242 /*@
243   PetscDualSpaceSetFromOptions - sets parameters in a PetscDualSpace from the options database
244 
245   Collective on PetscDualSpace
246 
247   Input Parameter:
248 . sp - the PetscDualSpace object to set options for
249 
250   Options Database:
251 . -petscspace_degree the approximation order of the space
252 
253   Level: developer
254 
255 .seealso PetscDualSpaceView()
256 @*/
257 PetscErrorCode PetscDualSpaceSetFromOptions(PetscDualSpace sp)
258 {
259   PetscDualSpaceReferenceCell refCell = PETSCDUALSPACE_REFCELL_SIMPLEX;
260   PetscInt                    refDim  = 0;
261   PetscBool                   flg;
262   const char                 *defaultType;
263   char                        name[256];
264   PetscErrorCode              ierr;
265 
266   PetscFunctionBegin;
267   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
268   if (!((PetscObject) sp)->type_name) {
269     defaultType = PETSCDUALSPACELAGRANGE;
270   } else {
271     defaultType = ((PetscObject) sp)->type_name;
272   }
273   if (!PetscSpaceRegisterAllCalled) {ierr = PetscSpaceRegisterAll();CHKERRQ(ierr);}
274 
275   ierr = PetscObjectOptionsBegin((PetscObject) sp);CHKERRQ(ierr);
276   ierr = PetscOptionsFList("-petscdualspace_type", "Dual space", "PetscDualSpaceSetType", PetscDualSpaceList, defaultType, name, 256, &flg);CHKERRQ(ierr);
277   if (flg) {
278     ierr = PetscDualSpaceSetType(sp, name);CHKERRQ(ierr);
279   } else if (!((PetscObject) sp)->type_name) {
280     ierr = PetscDualSpaceSetType(sp, defaultType);CHKERRQ(ierr);
281   }
282   ierr = PetscOptionsInt("-petscdualspace_degree", "The approximation order", "PetscDualSpaceSetOrder", sp->order, &sp->order, NULL);CHKERRQ(ierr);
283   ierr = PetscOptionsInt("-petscdualspace_components", "The number of components", "PetscDualSpaceSetNumComponents", sp->Nc, &sp->Nc, NULL);CHKERRQ(ierr);
284   if (sp->ops->setfromoptions) {
285     ierr = (*sp->ops->setfromoptions)(PetscOptionsObject,sp);CHKERRQ(ierr);
286   }
287   ierr = PetscOptionsInt("-petscdualspace_refdim", "The spatial dimension of the reference cell", "PetscDualSpaceSetReferenceCell", refDim, &refDim, NULL);CHKERRQ(ierr);
288   ierr = PetscOptionsEnum("-petscdualspace_refcell", "Reference cell", "PetscDualSpaceSetReferenceCell", PetscDualSpaceReferenceCells, (PetscEnum) refCell, (PetscEnum *) &refCell, &flg);CHKERRQ(ierr);
289   if (flg) {
290     DM K;
291 
292     if (!refDim) SETERRQ(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_INCOMP, "Reference cell specified without a dimension. Use -petscdualspace_refdim.");
293     ierr = PetscDualSpaceCreateReferenceCell(sp, refDim, refCell == PETSCDUALSPACE_REFCELL_SIMPLEX ? PETSC_TRUE : PETSC_FALSE, &K);CHKERRQ(ierr);
294     ierr = PetscDualSpaceSetDM(sp, K);CHKERRQ(ierr);
295     ierr = DMDestroy(&K);CHKERRQ(ierr);
296   }
297 
298   /* process any options handlers added with PetscObjectAddOptionsHandler() */
299   ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) sp);CHKERRQ(ierr);
300   ierr = PetscOptionsEnd();CHKERRQ(ierr);
301   sp->setfromoptionscalled = PETSC_TRUE;
302   PetscFunctionReturn(0);
303 }
304 
305 /*@
306   PetscDualSpaceSetUp - Construct a basis for the PetscDualSpace
307 
308   Collective on PetscDualSpace
309 
310   Input Parameter:
311 . sp - the PetscDualSpace object to setup
312 
313   Level: developer
314 
315 .seealso PetscDualSpaceView(), PetscDualSpaceDestroy()
316 @*/
317 PetscErrorCode PetscDualSpaceSetUp(PetscDualSpace sp)
318 {
319   PetscErrorCode ierr;
320 
321   PetscFunctionBegin;
322   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
323   if (sp->setupcalled) PetscFunctionReturn(0);
324   sp->setupcalled = PETSC_TRUE;
325   if (sp->ops->setup) {ierr = (*sp->ops->setup)(sp);CHKERRQ(ierr);}
326   if (sp->setfromoptionscalled) {ierr = PetscDualSpaceViewFromOptions(sp, NULL, "-petscdualspace_view");CHKERRQ(ierr);}
327   PetscFunctionReturn(0);
328 }
329 
330 /*@
331   PetscDualSpaceDestroy - Destroys a PetscDualSpace object
332 
333   Collective on PetscDualSpace
334 
335   Input Parameter:
336 . sp - the PetscDualSpace object to destroy
337 
338   Level: developer
339 
340 .seealso PetscDualSpaceView()
341 @*/
342 PetscErrorCode PetscDualSpaceDestroy(PetscDualSpace *sp)
343 {
344   PetscInt       dim, f;
345   PetscErrorCode ierr;
346 
347   PetscFunctionBegin;
348   if (!*sp) PetscFunctionReturn(0);
349   PetscValidHeaderSpecific((*sp), PETSCDUALSPACE_CLASSID, 1);
350 
351   if (--((PetscObject)(*sp))->refct > 0) {*sp = 0; PetscFunctionReturn(0);}
352   ((PetscObject) (*sp))->refct = 0;
353 
354   ierr = PetscDualSpaceGetDimension(*sp, &dim);CHKERRQ(ierr);
355   for (f = 0; f < dim; ++f) {
356     ierr = PetscQuadratureDestroy(&(*sp)->functional[f]);CHKERRQ(ierr);
357   }
358   ierr = PetscFree((*sp)->functional);CHKERRQ(ierr);
359   ierr = PetscQuadratureDestroy(&(*sp)->allPoints);CHKERRQ(ierr);
360   ierr = DMDestroy(&(*sp)->dm);CHKERRQ(ierr);
361 
362   if ((*sp)->ops->destroy) {ierr = (*(*sp)->ops->destroy)(*sp);CHKERRQ(ierr);}
363   ierr = PetscHeaderDestroy(sp);CHKERRQ(ierr);
364   PetscFunctionReturn(0);
365 }
366 
367 /*@
368   PetscDualSpaceCreate - Creates an empty PetscDualSpace object. The type can then be set with PetscDualSpaceSetType().
369 
370   Collective on MPI_Comm
371 
372   Input Parameter:
373 . comm - The communicator for the PetscDualSpace object
374 
375   Output Parameter:
376 . sp - The PetscDualSpace object
377 
378   Level: beginner
379 
380 .seealso: PetscDualSpaceSetType(), PETSCDUALSPACELAGRANGE
381 @*/
382 PetscErrorCode PetscDualSpaceCreate(MPI_Comm comm, PetscDualSpace *sp)
383 {
384   PetscDualSpace s;
385   PetscErrorCode ierr;
386 
387   PetscFunctionBegin;
388   PetscValidPointer(sp, 2);
389   ierr = PetscCitationsRegister(FECitation,&FEcite);CHKERRQ(ierr);
390   *sp  = NULL;
391   ierr = PetscFEInitializePackage();CHKERRQ(ierr);
392 
393   ierr = PetscHeaderCreate(s, PETSCDUALSPACE_CLASSID, "PetscDualSpace", "Dual Space", "PetscDualSpace", comm, PetscDualSpaceDestroy, PetscDualSpaceView);CHKERRQ(ierr);
394 
395   s->order = 0;
396   s->Nc    = 1;
397   s->k     = 0;
398   s->setupcalled = PETSC_FALSE;
399 
400   *sp = s;
401   PetscFunctionReturn(0);
402 }
403 
404 /*@
405   PetscDualSpaceDuplicate - Creates a duplicate PetscDualSpace object, however it is not setup.
406 
407   Collective on PetscDualSpace
408 
409   Input Parameter:
410 . sp - The original PetscDualSpace
411 
412   Output Parameter:
413 . spNew - The duplicate PetscDualSpace
414 
415   Level: beginner
416 
417 .seealso: PetscDualSpaceCreate(), PetscDualSpaceSetType()
418 @*/
419 PetscErrorCode PetscDualSpaceDuplicate(PetscDualSpace sp, PetscDualSpace *spNew)
420 {
421   PetscErrorCode ierr;
422 
423   PetscFunctionBegin;
424   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
425   PetscValidPointer(spNew, 2);
426   ierr = (*sp->ops->duplicate)(sp, spNew);CHKERRQ(ierr);
427   PetscFunctionReturn(0);
428 }
429 
430 /*@
431   PetscDualSpaceGetDM - Get the DM representing the reference cell
432 
433   Not collective
434 
435   Input Parameter:
436 . sp - The PetscDualSpace
437 
438   Output Parameter:
439 . dm - The reference cell
440 
441   Level: intermediate
442 
443 .seealso: PetscDualSpaceSetDM(), PetscDualSpaceCreate()
444 @*/
445 PetscErrorCode PetscDualSpaceGetDM(PetscDualSpace sp, DM *dm)
446 {
447   PetscFunctionBegin;
448   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
449   PetscValidPointer(dm, 2);
450   *dm = sp->dm;
451   PetscFunctionReturn(0);
452 }
453 
454 /*@
455   PetscDualSpaceSetDM - Get the DM representing the reference cell
456 
457   Not collective
458 
459   Input Parameters:
460 + sp - The PetscDualSpace
461 - dm - The reference cell
462 
463   Level: intermediate
464 
465 .seealso: PetscDualSpaceGetDM(), PetscDualSpaceCreate()
466 @*/
467 PetscErrorCode PetscDualSpaceSetDM(PetscDualSpace sp, DM dm)
468 {
469   PetscErrorCode ierr;
470 
471   PetscFunctionBegin;
472   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
473   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
474   ierr = DMDestroy(&sp->dm);CHKERRQ(ierr);
475   ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);
476   sp->dm = dm;
477   PetscFunctionReturn(0);
478 }
479 
480 /*@
481   PetscDualSpaceGetOrder - Get the order of the dual space
482 
483   Not collective
484 
485   Input Parameter:
486 . sp - The PetscDualSpace
487 
488   Output Parameter:
489 . order - The order
490 
491   Level: intermediate
492 
493 .seealso: PetscDualSpaceSetOrder(), PetscDualSpaceCreate()
494 @*/
495 PetscErrorCode PetscDualSpaceGetOrder(PetscDualSpace sp, PetscInt *order)
496 {
497   PetscFunctionBegin;
498   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
499   PetscValidPointer(order, 2);
500   *order = sp->order;
501   PetscFunctionReturn(0);
502 }
503 
504 /*@
505   PetscDualSpaceSetOrder - Set the order of the dual space
506 
507   Not collective
508 
509   Input Parameters:
510 + sp - The PetscDualSpace
511 - order - The order
512 
513   Level: intermediate
514 
515 .seealso: PetscDualSpaceGetOrder(), PetscDualSpaceCreate()
516 @*/
517 PetscErrorCode PetscDualSpaceSetOrder(PetscDualSpace sp, PetscInt order)
518 {
519   PetscFunctionBegin;
520   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
521   sp->order = order;
522   PetscFunctionReturn(0);
523 }
524 
525 /*@
526   PetscDualSpaceGetNumComponents - Return the number of components for this space
527 
528   Input Parameter:
529 . sp - The PetscDualSpace
530 
531   Output Parameter:
532 . Nc - The number of components
533 
534   Note: A vector space, for example, will have d components, where d is the spatial dimension
535 
536   Level: intermediate
537 
538 .seealso: PetscDualSpaceSetNumComponents(), PetscDualSpaceGetDimension(), PetscDualSpaceCreate(), PetscDualSpace
539 @*/
540 PetscErrorCode PetscDualSpaceGetNumComponents(PetscDualSpace sp, PetscInt *Nc)
541 {
542   PetscFunctionBegin;
543   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
544   PetscValidPointer(Nc, 2);
545   *Nc = sp->Nc;
546   PetscFunctionReturn(0);
547 }
548 
549 /*@
550   PetscDualSpaceSetNumComponents - Set the number of components for this space
551 
552   Input Parameters:
553 + sp - The PetscDualSpace
554 - order - The number of components
555 
556   Level: intermediate
557 
558 .seealso: PetscDualSpaceGetNumComponents(), PetscDualSpaceCreate(), PetscDualSpace
559 @*/
560 PetscErrorCode PetscDualSpaceSetNumComponents(PetscDualSpace sp, PetscInt Nc)
561 {
562   PetscFunctionBegin;
563   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
564   sp->Nc = Nc;
565   PetscFunctionReturn(0);
566 }
567 
568 /*@
569   PetscDualSpaceGetFunctional - Get the i-th basis functional in the dual space
570 
571   Not collective
572 
573   Input Parameters:
574 + sp - The PetscDualSpace
575 - i  - The basis number
576 
577   Output Parameter:
578 . functional - The basis functional
579 
580   Level: intermediate
581 
582 .seealso: PetscDualSpaceGetDimension(), PetscDualSpaceCreate()
583 @*/
584 PetscErrorCode PetscDualSpaceGetFunctional(PetscDualSpace sp, PetscInt i, PetscQuadrature *functional)
585 {
586   PetscInt       dim;
587   PetscErrorCode ierr;
588 
589   PetscFunctionBegin;
590   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
591   PetscValidPointer(functional, 3);
592   ierr = PetscDualSpaceGetDimension(sp, &dim);CHKERRQ(ierr);
593   if ((i < 0) || (i >= dim)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Functional index %d must be in [0, %d)", i, dim);
594   *functional = sp->functional[i];
595   PetscFunctionReturn(0);
596 }
597 
598 /*@
599   PetscDualSpaceGetDimension - Get the dimension of the dual space, i.e. the number of basis functionals
600 
601   Not collective
602 
603   Input Parameter:
604 . sp - The PetscDualSpace
605 
606   Output Parameter:
607 . dim - The dimension
608 
609   Level: intermediate
610 
611 .seealso: PetscDualSpaceGetFunctional(), PetscDualSpaceCreate()
612 @*/
613 PetscErrorCode PetscDualSpaceGetDimension(PetscDualSpace sp, PetscInt *dim)
614 {
615   PetscErrorCode ierr;
616 
617   PetscFunctionBegin;
618   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
619   PetscValidPointer(dim, 2);
620   *dim = 0;
621   if (sp->ops->getdimension) {ierr = (*sp->ops->getdimension)(sp, dim);CHKERRQ(ierr);}
622   PetscFunctionReturn(0);
623 }
624 
625 /*@C
626   PetscDualSpaceGetNumDof - Get the number of degrees of freedom for each spatial (topological) dimension
627 
628   Not collective
629 
630   Input Parameter:
631 . sp - The PetscDualSpace
632 
633   Output Parameter:
634 . numDof - An array of length dim+1 which holds the number of dofs for each dimension
635 
636   Level: intermediate
637 
638 .seealso: PetscDualSpaceGetFunctional(), PetscDualSpaceCreate()
639 @*/
640 PetscErrorCode PetscDualSpaceGetNumDof(PetscDualSpace sp, const PetscInt **numDof)
641 {
642   PetscErrorCode ierr;
643 
644   PetscFunctionBegin;
645   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
646   PetscValidPointer(numDof, 2);
647   ierr = (*sp->ops->getnumdof)(sp, numDof);CHKERRQ(ierr);
648   if (!*numDof) SETERRQ(PetscObjectComm((PetscObject) sp), PETSC_ERR_LIB, "Empty numDof[] returned from dual space implementation");
649   PetscFunctionReturn(0);
650 }
651 
652 PetscErrorCode PetscDualSpaceCreateSection(PetscDualSpace sp, PetscSection *section)
653 {
654   DM             dm;
655   PetscInt       pStart, pEnd, depth, h, offset;
656   const PetscInt *numDof;
657   PetscErrorCode ierr;
658 
659   PetscFunctionBegin;
660   ierr = PetscDualSpaceGetDM(sp,&dm);CHKERRQ(ierr);
661   ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr);
662   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)sp),section);CHKERRQ(ierr);
663   ierr = PetscSectionSetChart(*section,pStart,pEnd);CHKERRQ(ierr);
664   ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr);
665   ierr = PetscDualSpaceGetNumDof(sp,&numDof);CHKERRQ(ierr);
666   for (h = 0; h <= depth; h++) {
667     PetscInt hStart, hEnd, p, dof;
668 
669     ierr = DMPlexGetHeightStratum(dm,h,&hStart,&hEnd);CHKERRQ(ierr);
670     dof = numDof[depth - h];
671     for (p = hStart; p < hEnd; p++) {
672       ierr = PetscSectionSetDof(*section,p,dof);CHKERRQ(ierr);
673     }
674   }
675   ierr = PetscSectionSetUp(*section);CHKERRQ(ierr);
676   for (h = 0, offset = 0; h <= depth; h++) {
677     PetscInt hStart, hEnd, p, dof;
678 
679     ierr = DMPlexGetHeightStratum(dm,h,&hStart,&hEnd);CHKERRQ(ierr);
680     dof = numDof[depth - h];
681     for (p = hStart; p < hEnd; p++) {
682       ierr = PetscSectionGetDof(*section,p,&dof);CHKERRQ(ierr);
683       ierr = PetscSectionSetOffset(*section,p,offset);CHKERRQ(ierr);
684       offset += dof;
685     }
686   }
687   PetscFunctionReturn(0);
688 }
689 
690 /*@
691   PetscDualSpaceCreateReferenceCell - Create a DMPLEX with the appropriate FEM reference cell
692 
693   Collective on PetscDualSpace
694 
695   Input Parameters:
696 + sp      - The PetscDualSpace
697 . dim     - The spatial dimension
698 - simplex - Flag for simplex, otherwise use a tensor-product cell
699 
700   Output Parameter:
701 . refdm - The reference cell
702 
703   Level: advanced
704 
705 .keywords: PetscDualSpace, reference cell
706 .seealso: PetscDualSpaceCreate(), DMPLEX
707 @*/
708 PetscErrorCode PetscDualSpaceCreateReferenceCell(PetscDualSpace sp, PetscInt dim, PetscBool simplex, DM *refdm)
709 {
710   PetscErrorCode ierr;
711 
712   PetscFunctionBeginUser;
713   ierr = DMPlexCreateReferenceCell(PetscObjectComm((PetscObject) sp), dim, simplex, refdm);CHKERRQ(ierr);
714   PetscFunctionReturn(0);
715 }
716 
717 /*@C
718   PetscDualSpaceApply - Apply a functional from the dual space basis to an input function
719 
720   Input Parameters:
721 + sp      - The PetscDualSpace object
722 . f       - The basis functional index
723 . time    - The time
724 . 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)
725 . numComp - The number of components for the function
726 . func    - The input function
727 - ctx     - A context for the function
728 
729   Output Parameter:
730 . value   - numComp output values
731 
732   Note: The calling sequence for the callback func is given by:
733 
734 $ func(PetscInt dim, PetscReal time, const PetscReal x[],
735 $      PetscInt numComponents, PetscScalar values[], void *ctx)
736 
737   Level: developer
738 
739 .seealso: PetscDualSpaceCreate()
740 @*/
741 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)
742 {
743   PetscErrorCode ierr;
744 
745   PetscFunctionBegin;
746   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
747   PetscValidPointer(cgeom, 4);
748   PetscValidPointer(value, 8);
749   ierr = (*sp->ops->apply)(sp, f, time, cgeom, numComp, func, ctx, value);CHKERRQ(ierr);
750   PetscFunctionReturn(0);
751 }
752 
753 /*@C
754   PetscDualSpaceApplyAll - Apply all functionals from the dual space basis to the result of an evaluation at the points returned by PetscDualSpaceGetAllPoints()
755 
756   Input Parameters:
757 + sp        - The PetscDualSpace object
758 - pointEval - Evaluation at the points returned by PetscDualSpaceGetAllPoints()
759 
760   Output Parameter:
761 . spValue   - The values of all dual space functionals
762 
763   Level: developer
764 
765 .seealso: PetscDualSpaceCreate()
766 @*/
767 PetscErrorCode PetscDualSpaceApplyAll(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue)
768 {
769   PetscErrorCode ierr;
770 
771   PetscFunctionBegin;
772   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
773   ierr = (*sp->ops->applyall)(sp, pointEval, spValue);CHKERRQ(ierr);
774   PetscFunctionReturn(0);
775 }
776 
777 /*@C
778   PetscDualSpaceApplyDefault - Apply a functional from the dual space basis to an input function by assuming a point evaluation functional.
779 
780   Input Parameters:
781 + sp    - The PetscDualSpace object
782 . f     - The basis functional index
783 . time  - The time
784 . cgeom - A context with geometric information for this cell, we use v0 (the initial vertex) and J (the Jacobian)
785 . Nc    - The number of components for the function
786 . func  - The input function
787 - ctx   - A context for the function
788 
789   Output Parameter:
790 . value   - The output value
791 
792   Note: The calling sequence for the callback func is given by:
793 
794 $ func(PetscInt dim, PetscReal time, const PetscReal x[],
795 $      PetscInt numComponents, PetscScalar values[], void *ctx)
796 
797 and the idea is to evaluate the functional as an integral
798 
799 $ n(f) = int dx n(x) . f(x)
800 
801 where both n and f have Nc components.
802 
803   Level: developer
804 
805 .seealso: PetscDualSpaceCreate()
806 @*/
807 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)
808 {
809   DM               dm;
810   PetscQuadrature  n;
811   const PetscReal *points, *weights;
812   PetscReal        x[3];
813   PetscScalar     *val;
814   PetscInt         dim, dE, qNc, c, Nq, q;
815   PetscBool        isAffine;
816   PetscErrorCode   ierr;
817 
818   PetscFunctionBegin;
819   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
820   PetscValidPointer(value, 5);
821   ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr);
822   ierr = PetscDualSpaceGetFunctional(sp, f, &n);CHKERRQ(ierr);
823   ierr = PetscQuadratureGetData(n, &dim, &qNc, &Nq, &points, &weights);CHKERRQ(ierr);
824   if (dim != cgeom->dim) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature spatial dimension %D != cell geometry dimension %D", dim, cgeom->dim);
825   if (qNc != Nc) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature components %D != function components %D", qNc, Nc);
826   ierr = DMGetWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr);
827   *value = 0.0;
828   isAffine = cgeom->isAffine;
829   dE = cgeom->dimEmbed;
830   for (q = 0; q < Nq; ++q) {
831     if (isAffine) {
832       CoordinatesRefToReal(dE, cgeom->dim, cgeom->xi, cgeom->v, cgeom->J, &points[q*dim], x);
833       ierr = (*func)(dE, time, x, Nc, val, ctx);CHKERRQ(ierr);
834     } else {
835       ierr = (*func)(dE, time, &cgeom->v[dE*q], Nc, val, ctx);CHKERRQ(ierr);
836     }
837     for (c = 0; c < Nc; ++c) {
838       *value += val[c]*weights[q*Nc+c];
839     }
840   }
841   ierr = DMRestoreWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr);
842   PetscFunctionReturn(0);
843 }
844 
845 /*@C
846   PetscDualSpaceApplyAllDefault - Apply all functionals from the dual space basis to the result of an evaluation at the points returned by PetscDualSpaceGetAllPoints()
847 
848   Input Parameters:
849 + sp        - The PetscDualSpace object
850 - pointEval - Evaluation at the points returned by PetscDualSpaceGetAllPoints()
851 
852   Output Parameter:
853 . spValue   - The values of all dual space functionals
854 
855   Level: developer
856 
857 .seealso: PetscDualSpaceCreate()
858 @*/
859 PetscErrorCode PetscDualSpaceApplyAllDefault(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue)
860 {
861   PetscQuadrature  n;
862   const PetscReal *points, *weights;
863   PetscInt         qNc, c, Nq, q, f, spdim, Nc;
864   PetscInt         offset;
865   PetscErrorCode   ierr;
866 
867   PetscFunctionBegin;
868   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
869   PetscValidScalarPointer(pointEval, 2);
870   PetscValidScalarPointer(spValue, 5);
871   ierr = PetscDualSpaceGetDimension(sp, &spdim);CHKERRQ(ierr);
872   ierr = PetscDualSpaceGetNumComponents(sp, &Nc);CHKERRQ(ierr);
873   for (f = 0, offset = 0; f < spdim; f++) {
874     ierr = PetscDualSpaceGetFunctional(sp, f, &n);CHKERRQ(ierr);
875     ierr = PetscQuadratureGetData(n, NULL, &qNc, &Nq, &points, &weights);CHKERRQ(ierr);
876     if (qNc != Nc) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature components %D != function components %D", qNc, Nc);
877     spValue[f] = 0.0;
878     for (q = 0; q < Nq; ++q) {
879       for (c = 0; c < Nc; ++c) {
880         spValue[f] += pointEval[offset++]*weights[q*Nc+c];
881       }
882     }
883   }
884   PetscFunctionReturn(0);
885 }
886 
887 PetscErrorCode PetscDualSpaceGetAllPoints(PetscDualSpace sp, PetscQuadrature *allPoints)
888 {
889   PetscErrorCode ierr;
890 
891   PetscFunctionBegin;
892   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
893   PetscValidPointer(allPoints,2);
894   if (!sp->allPoints && sp->ops->createallpoints) {
895     ierr = (*sp->ops->createallpoints)(sp,&sp->allPoints);CHKERRQ(ierr);
896   }
897   *allPoints = sp->allPoints;
898   PetscFunctionReturn(0);
899 }
900 
901 PetscErrorCode PetscDualSpaceCreateAllPointsDefault(PetscDualSpace sp, PetscQuadrature *allPoints)
902 {
903   PetscInt        spdim;
904   PetscInt        numPoints, offset;
905   PetscReal       *points;
906   PetscInt        f, dim;
907   PetscQuadrature q;
908   PetscErrorCode  ierr;
909 
910   PetscFunctionBegin;
911   ierr = PetscDualSpaceGetDimension(sp,&spdim);CHKERRQ(ierr);
912   if (!spdim) {
913     ierr = PetscQuadratureCreate(PETSC_COMM_SELF,allPoints);CHKERRQ(ierr);
914     ierr = PetscQuadratureSetData(*allPoints,0,0,0,NULL,NULL);CHKERRQ(ierr);
915   }
916   ierr = PetscDualSpaceGetFunctional(sp,0,&q);CHKERRQ(ierr);
917   ierr = PetscQuadratureGetData(q,&dim,NULL,&numPoints,NULL,NULL);CHKERRQ(ierr);
918   for (f = 1; f < spdim; f++) {
919     PetscInt Np;
920 
921     ierr = PetscDualSpaceGetFunctional(sp,f,&q);CHKERRQ(ierr);
922     ierr = PetscQuadratureGetData(q,NULL,NULL,&Np,NULL,NULL);CHKERRQ(ierr);
923     numPoints += Np;
924   }
925   ierr = PetscMalloc1(dim*numPoints,&points);CHKERRQ(ierr);
926   for (f = 0, offset = 0; f < spdim; f++) {
927     const PetscReal *p;
928     PetscInt        Np, i;
929 
930     ierr = PetscDualSpaceGetFunctional(sp,f,&q);CHKERRQ(ierr);
931     ierr = PetscQuadratureGetData(q,NULL,NULL,&Np,&p,NULL);CHKERRQ(ierr);
932     for (i = 0; i < Np * dim; i++) {
933       points[offset + i] = p[i];
934     }
935     offset += Np * dim;
936   }
937   ierr = PetscQuadratureCreate(PETSC_COMM_SELF,allPoints);CHKERRQ(ierr);
938   ierr = PetscQuadratureSetData(*allPoints,dim,0,numPoints,points,NULL);CHKERRQ(ierr);
939   PetscFunctionReturn(0);
940 }
941 
942 /*@C
943   PetscDualSpaceApplyFVM - Apply a functional from the dual space basis to an input function by assuming a point evaluation functional at the cell centroid.
944 
945   Input Parameters:
946 + sp    - The PetscDualSpace object
947 . f     - The basis functional index
948 . time  - The time
949 . cgeom - A context with geometric information for this cell, we currently just use the centroid
950 . Nc    - The number of components for the function
951 . func  - The input function
952 - ctx   - A context for the function
953 
954   Output Parameter:
955 . value - The output value (scalar)
956 
957   Note: The calling sequence for the callback func is given by:
958 
959 $ func(PetscInt dim, PetscReal time, const PetscReal x[],
960 $      PetscInt numComponents, PetscScalar values[], void *ctx)
961 
962 and the idea is to evaluate the functional as an integral
963 
964 $ n(f) = int dx n(x) . f(x)
965 
966 where both n and f have Nc components.
967 
968   Level: developer
969 
970 .seealso: PetscDualSpaceCreate()
971 @*/
972 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)
973 {
974   DM               dm;
975   PetscQuadrature  n;
976   const PetscReal *points, *weights;
977   PetscScalar     *val;
978   PetscInt         dimEmbed, qNc, c, Nq, q;
979   PetscErrorCode   ierr;
980 
981   PetscFunctionBegin;
982   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
983   PetscValidPointer(value, 5);
984   ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr);
985   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
986   ierr = PetscDualSpaceGetFunctional(sp, f, &n);CHKERRQ(ierr);
987   ierr = PetscQuadratureGetData(n, NULL, &qNc, &Nq, &points, &weights);CHKERRQ(ierr);
988   if (qNc != Nc) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature components %D != function components %D", qNc, Nc);
989   ierr = DMGetWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr);
990   *value = 0.;
991   for (q = 0; q < Nq; ++q) {
992     ierr = (*func)(dimEmbed, time, cgeom->centroid, Nc, val, ctx);CHKERRQ(ierr);
993     for (c = 0; c < Nc; ++c) {
994       *value += val[c]*weights[q*Nc+c];
995     }
996   }
997   ierr = DMRestoreWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr);
998   PetscFunctionReturn(0);
999 }
1000 
1001 /*@
1002   PetscDualSpaceGetHeightSubspace - Get the subset of the dual space basis that is supported on a mesh point of a
1003   given height.  This assumes that the reference cell is symmetric over points of this height.
1004 
1005   If the dual space is not defined on mesh points of the given height (e.g. if the space is discontinuous and
1006   pointwise values are not defined on the element boundaries), or if the implementation of PetscDualSpace does not
1007   support extracting subspaces, then NULL is returned.
1008 
1009   This does not increment the reference count on the returned dual space, and the user should not destroy it.
1010 
1011   Not collective
1012 
1013   Input Parameters:
1014 + sp - the PetscDualSpace object
1015 - height - the height of the mesh point for which the subspace is desired
1016 
1017   Output Parameter:
1018 . subsp - the subspace.  Note that the functionals in the subspace are with respect to the intrinsic geometry of the
1019   point, which will be of lesser dimension if height > 0.
1020 
1021   Level: advanced
1022 
1023 .seealso: PetscSpaceGetHeightSubspace(), PetscDualSpace
1024 @*/
1025 PetscErrorCode PetscDualSpaceGetHeightSubspace(PetscDualSpace sp, PetscInt height, PetscDualSpace *subsp)
1026 {
1027   PetscErrorCode ierr;
1028 
1029   PetscFunctionBegin;
1030   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1031   PetscValidPointer(subsp, 3);
1032   *subsp = NULL;
1033   if (sp->ops->getheightsubspace) {ierr = (*sp->ops->getheightsubspace)(sp, height, subsp);CHKERRQ(ierr);}
1034   PetscFunctionReturn(0);
1035 }
1036 
1037 /*@
1038   PetscDualSpaceGetPointSubspace - Get the subset of the dual space basis that is supported on a particular mesh point.
1039 
1040   If the dual space is not defined on the mesh point (e.g. if the space is discontinuous and pointwise values are not
1041   defined on the element boundaries), or if the implementation of PetscDualSpace does not support extracting
1042   subspaces, then NULL is returned.
1043 
1044   This does not increment the reference count on the returned dual space, and the user should not destroy it.
1045 
1046   Not collective
1047 
1048   Input Parameters:
1049 + sp - the PetscDualSpace object
1050 - point - the point (in the dual space's DM) for which the subspace is desired
1051 
1052   Output Parameters:
1053   bdsp - the subspace.  Note that the functionals in the subspace are with respect to the intrinsic geometry of the
1054   point, which will be of lesser dimension if height > 0.
1055 
1056   Level: advanced
1057 
1058 .seealso: PetscDualSpace
1059 @*/
1060 PetscErrorCode PetscDualSpaceGetPointSubspace(PetscDualSpace sp, PetscInt point, PetscDualSpace *bdsp)
1061 {
1062   PetscErrorCode ierr;
1063 
1064   PetscFunctionBegin;
1065   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1066   PetscValidPointer(bdsp,2);
1067   *bdsp = NULL;
1068   if (sp->ops->getpointsubspace) {
1069     ierr = (*sp->ops->getpointsubspace)(sp,point,bdsp);CHKERRQ(ierr);
1070   } else if (sp->ops->getheightsubspace) {
1071     DM       dm;
1072     DMLabel  label;
1073     PetscInt dim, depth, height;
1074 
1075     ierr = PetscDualSpaceGetDM(sp,&dm);CHKERRQ(ierr);
1076     ierr = DMPlexGetDepth(dm,&dim);CHKERRQ(ierr);
1077     ierr = DMPlexGetDepthLabel(dm,&label);CHKERRQ(ierr);
1078     ierr = DMLabelGetValue(label,point,&depth);CHKERRQ(ierr);
1079     height = dim - depth;
1080     ierr = (*sp->ops->getheightsubspace)(sp,height,bdsp);CHKERRQ(ierr);
1081   }
1082   PetscFunctionReturn(0);
1083 }
1084 
1085 /*@C
1086   PetscDualSpaceGetSymmetries - Returns a description of the symmetries of this basis
1087 
1088   Not collective
1089 
1090   Input Parameter:
1091 . sp - the PetscDualSpace object
1092 
1093   Output Parameters:
1094 + perms - Permutations of the local degrees of freedom, parameterized by the point orientation
1095 - flips - Sign reversal of the local degrees of freedom, parameterized by the point orientation
1096 
1097   Note: The permutation and flip arrays are organized in the following way
1098 $ perms[p][ornt][dof # on point] = new local dof #
1099 $ flips[p][ornt][dof # on point] = reversal or not
1100 
1101   Level: developer
1102 
1103 .seealso: PetscDualSpaceSetSymmetries()
1104 @*/
1105 PetscErrorCode PetscDualSpaceGetSymmetries(PetscDualSpace sp, const PetscInt ****perms, const PetscScalar ****flips)
1106 {
1107   PetscErrorCode ierr;
1108 
1109   PetscFunctionBegin;
1110   PetscValidHeaderSpecific(sp,PETSCDUALSPACE_CLASSID,1);
1111   if (perms) {PetscValidPointer(perms,2); *perms = NULL;}
1112   if (flips) {PetscValidPointer(flips,3); *flips = NULL;}
1113   if (sp->ops->getsymmetries) {ierr = (sp->ops->getsymmetries)(sp,perms,flips);CHKERRQ(ierr);}
1114   PetscFunctionReturn(0);
1115 }
1116 
1117 /*@
1118   PetscDualSpaceGetDeRahm - Get the k-simplex associated with the functionals in this dual space
1119 
1120   Input Parameter:
1121 . dsp - The PetscDualSpace
1122 
1123   Output Parameter:
1124 . k   - The simplex dimension
1125 
1126   Level: advanced
1127 
1128   Note: Currently supported values are
1129 $ 0: These are H_1 methods that only transform coordinates
1130 $ 1: These are Hcurl methods that transform functions using the covariant Piola transform (COVARIANT_PIOLA_TRANSFORM)
1131 $ 2: These are the same as 1
1132 $ 3: These are Hdiv methods that transform functions using the contravariant Piola transform (CONTRAVARIANT_PIOLA_TRANSFORM)
1133 
1134 .seealso: PetscDualSpacePullback(), PetscDualSpacePushforward(), PetscDualSpaceTransform(), PetscDualSpaceTransformType
1135 @*/
1136 PetscErrorCode PetscDualSpaceGetDeRahm(PetscDualSpace dsp, PetscInt *k)
1137 {
1138   PetscFunctionBeginHot;
1139   PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1);
1140   PetscValidPointer(k, 2);
1141   *k = dsp->k;
1142   PetscFunctionReturn(0);
1143 }
1144 
1145 /*@C
1146   PetscDualSpaceTransform - Transform the function values
1147 
1148   Input Parameters:
1149 + dsp       - The PetscDualSpace
1150 . trans     - The type of transform
1151 . isInverse - Flag to invert the transform
1152 . fegeom    - The cell geometry
1153 . Nv        - The number of function samples
1154 . Nc        - The number of function components
1155 - vals      - The function values
1156 
1157   Output Parameter:
1158 . vals      - The transformed function values
1159 
1160   Level: developer
1161 
1162 .seealso: PetscDualSpaceTransformGradient(), PetscDualSpacePullback(), PetscDualSpacePushforward(), PetscDualSpaceTransform
1163 @*/
1164 PetscErrorCode PetscDualSpaceTransform(PetscDualSpace dsp, PetscDualSpaceTransformType trans, PetscBool isInverse, PetscFEGeom *fegeom, PetscInt Nv, PetscInt Nc, PetscScalar vals[])
1165 {
1166   DM             dm;
1167   PetscInt       dim, v, c;
1168   PetscErrorCode ierr;
1169 
1170   PetscFunctionBeginHot;
1171   PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1);
1172   PetscValidPointer(fegeom, 4);
1173   PetscValidPointer(vals, 7);
1174   ierr = PetscDualSpaceGetDM(dsp, &dm);CHKERRQ(ierr);
1175   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1176   /* Assume its a vector, otherwise assume its a bunch of scalars */
1177   if (Nc == 1 || Nc != dim) PetscFunctionReturn(0);
1178   switch (trans) {
1179     case IDENTITY_TRANSFORM: break;
1180     case COVARIANT_PIOLA_TRANSFORM: /* Covariant Piola mapping $\sigma^*(F) = J^{-T} F \circ \phi^{-1)$ */
1181     if (isInverse) {
1182       for (v = 0; v < Nv; ++v) {
1183         switch (dim)
1184         {
1185           case 2: DMPlex_MultTranspose2D_Internal(fegeom->J, 1, &vals[v*Nc], &vals[v*Nc]);break;
1186           case 3: DMPlex_MultTranspose3D_Internal(fegeom->J, 1, &vals[v*Nc], &vals[v*Nc]);break;
1187           default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim);
1188         }
1189       }
1190     } else {
1191       for (v = 0; v < Nv; ++v) {
1192         switch (dim)
1193         {
1194           case 2: DMPlex_MultTranspose2D_Internal(fegeom->invJ, 1, &vals[v*Nc], &vals[v*Nc]);break;
1195           case 3: DMPlex_MultTranspose3D_Internal(fegeom->invJ, 1, &vals[v*Nc], &vals[v*Nc]);break;
1196           default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim);
1197         }
1198       }
1199     }
1200     break;
1201     case CONTRAVARIANT_PIOLA_TRANSFORM: /* Contravariant Piola mapping $\sigma^*(F) = \frac{1}{|\det J|} J F \circ \phi^{-1}$ */
1202     if (isInverse) {
1203       for (v = 0; v < Nv; ++v) {
1204         switch (dim)
1205         {
1206           case 2: DMPlex_Mult2D_Internal(fegeom->invJ, 1, &vals[v*Nc], &vals[v*Nc]);break;
1207           case 3: DMPlex_Mult3D_Internal(fegeom->invJ, 1, &vals[v*Nc], &vals[v*Nc]);break;
1208           default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim);
1209         }
1210         for (c = 0; c < Nc; ++c) vals[v*Nc+c] *= fegeom->detJ[0];
1211       }
1212     } else {
1213       for (v = 0; v < Nv; ++v) {
1214         switch (dim)
1215         {
1216           case 2: DMPlex_Mult2D_Internal(fegeom->J, 1, &vals[v*Nc], &vals[v*Nc]);break;
1217           case 3: DMPlex_Mult3D_Internal(fegeom->J, 1, &vals[v*Nc], &vals[v*Nc]);break;
1218           default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim);
1219         }
1220         for (c = 0; c < Nc; ++c) vals[v*Nc+c] /= fegeom->detJ[0];
1221       }
1222     }
1223     break;
1224   }
1225   PetscFunctionReturn(0);
1226 }
1227 /*@C
1228   PetscDualSpaceTransformGradient - Transform the function gradient values
1229 
1230   Input Parameters:
1231 + dsp       - The PetscDualSpace
1232 . trans     - The type of transform
1233 . isInverse - Flag to invert the transform
1234 . fegeom    - The cell geometry
1235 . Nv        - The number of function gradient samples
1236 . Nc        - The number of function components
1237 - vals      - The function gradient values
1238 
1239   Output Parameter:
1240 . vals      - The transformed function values
1241 
1242   Level: developer
1243 
1244 .seealso: PetscDualSpaceTransform(), PetscDualSpacePullback(), PetscDualSpacePushforward(), PetscDualSpaceTransform
1245 @*/
1246 PetscErrorCode PetscDualSpaceTransformGradient(PetscDualSpace dsp, PetscDualSpaceTransformType trans, PetscBool isInverse, PetscFEGeom *fegeom, PetscInt Nv, PetscInt Nc, PetscScalar vals[])
1247 {
1248   DM             dm;
1249   PetscInt       dim, v, c, d;
1250   PetscErrorCode ierr;
1251 
1252   PetscFunctionBeginHot;
1253   PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1);
1254   PetscValidPointer(fegeom, 4);
1255   PetscValidPointer(vals, 7);
1256   ierr = PetscDualSpaceGetDM(dsp, &dm);CHKERRQ(ierr);
1257   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1258   /* Transform gradient */
1259   for (v = 0; v < Nv; ++v) {
1260     for (c = 0; c < Nc; ++c) {
1261       switch (dim)
1262       {
1263         case 2: DMPlex_MultTranspose2D_Internal(fegeom->invJ, 1, &vals[(v*Nc+c)*dim], &vals[(v*Nc+c)*dim]);break;
1264         case 3: DMPlex_MultTranspose3D_Internal(fegeom->invJ, 1, &vals[(v*Nc+c)*dim], &vals[(v*Nc+c)*dim]);break;
1265         default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim);
1266       }
1267     }
1268   }
1269   /* Assume its a vector, otherwise assume its a bunch of scalars */
1270   if (Nc == 1 || Nc != dim) PetscFunctionReturn(0);
1271   switch (trans) {
1272     case IDENTITY_TRANSFORM: break;
1273     case COVARIANT_PIOLA_TRANSFORM: /* Covariant Piola mapping $\sigma^*(F) = J^{-T} F \circ \phi^{-1)$ */
1274     if (isInverse) {
1275       for (v = 0; v < Nv; ++v) {
1276         for (d = 0; d < dim; ++d) {
1277           switch (dim)
1278           {
1279             case 2: DMPlex_MultTranspose2D_Internal(fegeom->J, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
1280             case 3: DMPlex_MultTranspose3D_Internal(fegeom->J, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
1281             default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim);
1282           }
1283         }
1284       }
1285     } else {
1286       for (v = 0; v < Nv; ++v) {
1287         for (d = 0; d < dim; ++d) {
1288           switch (dim)
1289           {
1290             case 2: DMPlex_MultTranspose2D_Internal(fegeom->invJ, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
1291             case 3: DMPlex_MultTranspose3D_Internal(fegeom->invJ, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
1292             default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim);
1293           }
1294         }
1295       }
1296     }
1297     break;
1298     case CONTRAVARIANT_PIOLA_TRANSFORM: /* Contravariant Piola mapping $\sigma^*(F) = \frac{1}{|\det J|} J F \circ \phi^{-1}$ */
1299     if (isInverse) {
1300       for (v = 0; v < Nv; ++v) {
1301         for (d = 0; d < dim; ++d) {
1302           switch (dim)
1303           {
1304             case 2: DMPlex_Mult2D_Internal(fegeom->invJ, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
1305             case 3: DMPlex_Mult3D_Internal(fegeom->invJ, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
1306             default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim);
1307           }
1308           for (c = 0; c < Nc; ++c) vals[(v*Nc+c)*dim+d] *= fegeom->detJ[0];
1309         }
1310       }
1311     } else {
1312       for (v = 0; v < Nv; ++v) {
1313         for (d = 0; d < dim; ++d) {
1314           switch (dim)
1315           {
1316             case 2: DMPlex_Mult2D_Internal(fegeom->J, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
1317             case 3: DMPlex_Mult3D_Internal(fegeom->J, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
1318             default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim);
1319           }
1320           for (c = 0; c < Nc; ++c) vals[(v*Nc+c)*dim+d] /= fegeom->detJ[0];
1321         }
1322       }
1323     }
1324     break;
1325   }
1326   PetscFunctionReturn(0);
1327 }
1328 
1329 /*@C
1330   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.
1331 
1332   Input Parameters:
1333 + dsp        - The PetscDualSpace
1334 . fegeom     - The geometry for this cell
1335 . Nq         - The number of function samples
1336 . Nc         - The number of function components
1337 - pointEval  - The function values
1338 
1339   Output Parameter:
1340 . pointEval  - The transformed function values
1341 
1342   Level: advanced
1343 
1344   Note: Functions transform in a complementary way (pushforward) to functionals, so that the scalar product is invariant. The type of transform is dependent on the associated k-simplex from the DeRahm complex.
1345 
1346 .seealso: PetscDualSpacePushforward(), PetscDualSpaceTransform(), PetscDualSpaceGetDeRahm()
1347 @*/
1348 PetscErrorCode PetscDualSpacePullback(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, PetscScalar pointEval[])
1349 {
1350   PetscDualSpaceTransformType trans;
1351   PetscInt                    k;
1352   PetscErrorCode              ierr;
1353 
1354   PetscFunctionBeginHot;
1355   PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1);
1356   PetscValidPointer(fegeom, 2);
1357   PetscValidPointer(pointEval, 5);
1358   /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k.
1359      This determines their transformation properties. */
1360   ierr = PetscDualSpaceGetDeRahm(dsp, &k);CHKERRQ(ierr);
1361   switch (k)
1362   {
1363     case 0: /* H^1 point evaluations */
1364     trans = IDENTITY_TRANSFORM;break;
1365     case 1: /* Hcurl preserves tangential edge traces  */
1366     case 2:
1367     trans = COVARIANT_PIOLA_TRANSFORM;break;
1368     case 3: /* Hdiv preserve normal traces */
1369     trans = CONTRAVARIANT_PIOLA_TRANSFORM;break;
1370     default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %D for transformation", k);
1371   }
1372   ierr = PetscDualSpaceTransform(dsp, trans, PETSC_TRUE, fegeom, Nq, Nc, pointEval);CHKERRQ(ierr);
1373   PetscFunctionReturn(0);
1374 }
1375 
1376 /*@C
1377   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.
1378 
1379   Input Parameters:
1380 + dsp        - The PetscDualSpace
1381 . fegeom     - The geometry for this cell
1382 . Nq         - The number of function samples
1383 . Nc         - The number of function components
1384 - pointEval  - The function values
1385 
1386   Output Parameter:
1387 . pointEval  - The transformed function values
1388 
1389   Level: advanced
1390 
1391   Note: Functionals transform in a complementary way (pullback) to functions, so that the scalar product is invariant. The type of transform is dependent on the associated k-simplex from the DeRahm complex.
1392 
1393 .seealso: PetscDualSpacePullback(), PetscDualSpaceTransform(), PetscDualSpaceGetDeRahm()
1394 @*/
1395 PetscErrorCode PetscDualSpacePushforward(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, PetscScalar pointEval[])
1396 {
1397   PetscDualSpaceTransformType trans;
1398   PetscInt                    k;
1399   PetscErrorCode              ierr;
1400 
1401   PetscFunctionBeginHot;
1402   PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1);
1403   PetscValidPointer(fegeom, 2);
1404   PetscValidPointer(pointEval, 5);
1405   /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k.
1406      This determines their transformation properties. */
1407   ierr = PetscDualSpaceGetDeRahm(dsp, &k);CHKERRQ(ierr);
1408   switch (k)
1409   {
1410     case 0: /* H^1 point evaluations */
1411     trans = IDENTITY_TRANSFORM;break;
1412     case 1: /* Hcurl preserves tangential edge traces  */
1413     case 2:
1414     trans = COVARIANT_PIOLA_TRANSFORM;break;
1415     case 3: /* Hdiv preserve normal traces */
1416     trans = CONTRAVARIANT_PIOLA_TRANSFORM;break;
1417     default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %D for transformation", k);
1418   }
1419   ierr = PetscDualSpaceTransform(dsp, trans, PETSC_FALSE, fegeom, Nq, Nc, pointEval);CHKERRQ(ierr);
1420   PetscFunctionReturn(0);
1421 }
1422 
1423 /*@C
1424   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.
1425 
1426   Input Parameters:
1427 + dsp        - The PetscDualSpace
1428 . fegeom     - The geometry for this cell
1429 . Nq         - The number of function gradient samples
1430 . Nc         - The number of function components
1431 - pointEval  - The function gradient values
1432 
1433   Output Parameter:
1434 . pointEval  - The transformed function gradient values
1435 
1436   Level: advanced
1437 
1438   Note: Functionals transform in a complementary way (pullback) to functions, so that the scalar product is invariant. The type of transform is dependent on the associated k-simplex from the DeRahm complex.
1439 
1440 .seealso: PetscDualSpacePushforward(), PPetscDualSpacePullback(), PetscDualSpaceTransform(), PetscDualSpaceGetDeRahm()
1441 @*/PetscErrorCode PetscDualSpacePushforwardGradient(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, PetscScalar pointEval[])
1442 {
1443   PetscDualSpaceTransformType trans;
1444   PetscInt                    k;
1445   PetscErrorCode              ierr;
1446 
1447   PetscFunctionBeginHot;
1448   PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1);
1449   PetscValidPointer(fegeom, 2);
1450   PetscValidPointer(pointEval, 5);
1451   /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k.
1452      This determines their transformation properties. */
1453   ierr = PetscDualSpaceGetDeRahm(dsp, &k);CHKERRQ(ierr);
1454   switch (k)
1455   {
1456     case 0: /* H^1 point evaluations */
1457     trans = IDENTITY_TRANSFORM;break;
1458     case 1: /* Hcurl preserves tangential edge traces  */
1459     case 2:
1460     trans = COVARIANT_PIOLA_TRANSFORM;break;
1461     case 3: /* Hdiv preserve normal traces */
1462     trans = CONTRAVARIANT_PIOLA_TRANSFORM;break;
1463     default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %D for transformation", k);
1464   }
1465   ierr = PetscDualSpaceTransformGradient(dsp, trans, PETSC_FALSE, fegeom, Nq, Nc, pointEval);CHKERRQ(ierr);
1466   PetscFunctionReturn(0);
1467 }
1468