xref: /petsc/src/dm/dt/dualspace/interface/dualspace.c (revision e5a36eccef3d6b83a2c625c30d0dfd5adb4001f2)
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 /*@C
10   PetscDualSpaceRegister - Adds a new PetscDualSpace implementation
11 
12   Not Collective
13 
14   Input Parameters:
15 + name        - The name of a new user-defined creation routine
16 - create_func - The creation routine itself
17 
18   Notes:
19   PetscDualSpaceRegister() may be called multiple times to add several user-defined PetscDualSpaces
20 
21   Sample usage:
22 .vb
23     PetscDualSpaceRegister("my_space", MyPetscDualSpaceCreate);
24 .ve
25 
26   Then, your PetscDualSpace type can be chosen with the procedural interface via
27 .vb
28     PetscDualSpaceCreate(MPI_Comm, PetscDualSpace *);
29     PetscDualSpaceSetType(PetscDualSpace, "my_dual_space");
30 .ve
31    or at runtime via the option
32 .vb
33     -petscdualspace_type my_dual_space
34 .ve
35 
36   Level: advanced
37 
38 .keywords: PetscDualSpace, register
39 .seealso: PetscDualSpaceRegisterAll(), PetscDualSpaceRegisterDestroy()
40 
41 @*/
42 PetscErrorCode PetscDualSpaceRegister(const char sname[], PetscErrorCode (*function)(PetscDualSpace))
43 {
44   PetscErrorCode ierr;
45 
46   PetscFunctionBegin;
47   ierr = PetscFunctionListAdd(&PetscDualSpaceList, sname, function);CHKERRQ(ierr);
48   PetscFunctionReturn(0);
49 }
50 
51 /*@C
52   PetscDualSpaceSetType - Builds a particular PetscDualSpace
53 
54   Collective on PetscDualSpace
55 
56   Input Parameters:
57 + sp   - The PetscDualSpace object
58 - name - The kind of space
59 
60   Options Database Key:
61 . -petscdualspace_type <type> - Sets the PetscDualSpace type; use -help for a list of available types
62 
63   Level: intermediate
64 
65 .keywords: PetscDualSpace, set, type
66 .seealso: PetscDualSpaceGetType(), PetscDualSpaceCreate()
67 @*/
68 PetscErrorCode PetscDualSpaceSetType(PetscDualSpace sp, PetscDualSpaceType name)
69 {
70   PetscErrorCode (*r)(PetscDualSpace);
71   PetscBool      match;
72   PetscErrorCode ierr;
73 
74   PetscFunctionBegin;
75   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
76   ierr = PetscObjectTypeCompare((PetscObject) sp, name, &match);CHKERRQ(ierr);
77   if (match) PetscFunctionReturn(0);
78 
79   if (!PetscDualSpaceRegisterAllCalled) {ierr = PetscDualSpaceRegisterAll();CHKERRQ(ierr);}
80   ierr = PetscFunctionListFind(PetscDualSpaceList, name, &r);CHKERRQ(ierr);
81   if (!r) SETERRQ1(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscDualSpace type: %s", name);
82 
83   if (sp->ops->destroy) {
84     ierr             = (*sp->ops->destroy)(sp);CHKERRQ(ierr);
85     sp->ops->destroy = NULL;
86   }
87   ierr = (*r)(sp);CHKERRQ(ierr);
88   ierr = PetscObjectChangeTypeName((PetscObject) sp, name);CHKERRQ(ierr);
89   PetscFunctionReturn(0);
90 }
91 
92 /*@C
93   PetscDualSpaceGetType - Gets the PetscDualSpace type name (as a string) from the object.
94 
95   Not Collective
96 
97   Input Parameter:
98 . sp  - The PetscDualSpace
99 
100   Output Parameter:
101 . name - The PetscDualSpace type name
102 
103   Level: intermediate
104 
105 .keywords: PetscDualSpace, get, type, name
106 .seealso: PetscDualSpaceSetType(), PetscDualSpaceCreate()
107 @*/
108 PetscErrorCode PetscDualSpaceGetType(PetscDualSpace sp, PetscDualSpaceType *name)
109 {
110   PetscErrorCode ierr;
111 
112   PetscFunctionBegin;
113   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
114   PetscValidPointer(name, 2);
115   if (!PetscDualSpaceRegisterAllCalled) {
116     ierr = PetscDualSpaceRegisterAll();CHKERRQ(ierr);
117   }
118   *name = ((PetscObject) sp)->type_name;
119   PetscFunctionReturn(0);
120 }
121 
122 /*@
123   PetscDualSpaceView - Views a PetscDualSpace
124 
125   Collective on PetscDualSpace
126 
127   Input Parameter:
128 + sp - the PetscDualSpace object to view
129 - v  - the viewer
130 
131   Level: developer
132 
133 .seealso PetscDualSpaceDestroy()
134 @*/
135 PetscErrorCode PetscDualSpaceView(PetscDualSpace sp, PetscViewer v)
136 {
137   PetscBool      iascii;
138   PetscErrorCode ierr;
139 
140   PetscFunctionBegin;
141   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
142   if (v) PetscValidHeaderSpecific(v, PETSC_VIEWER_CLASSID, 2);
143   if (!v) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) sp), &v);CHKERRQ(ierr);}
144   ierr = PetscObjectPrintClassNamePrefixType((PetscObject)sp, v);CHKERRQ(ierr);
145   ierr = PetscObjectTypeCompare((PetscObject) v, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
146   ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr);
147   if (iascii) {ierr = PetscViewerASCIIPrintf(v, "Dual space of order %D with %D components\n", sp->order, sp->Nc);CHKERRQ(ierr);}
148   if (sp->ops->view) {ierr = (*sp->ops->view)(sp, v);CHKERRQ(ierr);}
149   ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr);
150   PetscFunctionReturn(0);
151 }
152 
153 /*@
154   PetscDualSpaceSetFromOptions - sets parameters in a PetscDualSpace from the options database
155 
156   Collective on PetscDualSpace
157 
158   Input Parameter:
159 . sp - the PetscDualSpace object to set options for
160 
161   Options Database:
162 . -petscspace_degree the approximation order of the space
163 
164   Level: developer
165 
166 .seealso PetscDualSpaceView()
167 @*/
168 PetscErrorCode PetscDualSpaceSetFromOptions(PetscDualSpace sp)
169 {
170   const char    *defaultType;
171   char           name[256];
172   PetscBool      flg;
173   PetscErrorCode ierr;
174 
175   PetscFunctionBegin;
176   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
177   if (!((PetscObject) sp)->type_name) {
178     defaultType = PETSCDUALSPACELAGRANGE;
179   } else {
180     defaultType = ((PetscObject) sp)->type_name;
181   }
182   if (!PetscSpaceRegisterAllCalled) {ierr = PetscSpaceRegisterAll();CHKERRQ(ierr);}
183 
184   ierr = PetscObjectOptionsBegin((PetscObject) sp);CHKERRQ(ierr);
185   ierr = PetscOptionsFList("-petscdualspace_type", "Dual space", "PetscDualSpaceSetType", PetscDualSpaceList, defaultType, name, 256, &flg);CHKERRQ(ierr);
186   if (flg) {
187     ierr = PetscDualSpaceSetType(sp, name);CHKERRQ(ierr);
188   } else if (!((PetscObject) sp)->type_name) {
189     ierr = PetscDualSpaceSetType(sp, defaultType);CHKERRQ(ierr);
190   }
191   ierr = PetscOptionsInt("-petscdualspace_degree", "The approximation order", "PetscDualSpaceSetOrder", sp->order, &sp->order, NULL);CHKERRQ(ierr);
192   ierr = PetscOptionsInt("-petscdualspace_components", "The number of components", "PetscDualSpaceSetNumComponents", sp->Nc, &sp->Nc, NULL);CHKERRQ(ierr);
193   if (sp->ops->setfromoptions) {
194     ierr = (*sp->ops->setfromoptions)(PetscOptionsObject,sp);CHKERRQ(ierr);
195   }
196   /* process any options handlers added with PetscObjectAddOptionsHandler() */
197   ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) sp);CHKERRQ(ierr);
198   ierr = PetscOptionsEnd();CHKERRQ(ierr);
199   ierr = PetscDualSpaceViewFromOptions(sp, NULL, "-petscdualspace_view");CHKERRQ(ierr);
200   PetscFunctionReturn(0);
201 }
202 
203 /*@
204   PetscDualSpaceSetUp - Construct a basis for the PetscDualSpace
205 
206   Collective on PetscDualSpace
207 
208   Input Parameter:
209 . sp - the PetscDualSpace object to setup
210 
211   Level: developer
212 
213 .seealso PetscDualSpaceView(), PetscDualSpaceDestroy()
214 @*/
215 PetscErrorCode PetscDualSpaceSetUp(PetscDualSpace sp)
216 {
217   PetscErrorCode ierr;
218 
219   PetscFunctionBegin;
220   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
221   if (sp->setupcalled) PetscFunctionReturn(0);
222   sp->setupcalled = PETSC_TRUE;
223   if (sp->ops->setup) {ierr = (*sp->ops->setup)(sp);CHKERRQ(ierr);}
224   PetscFunctionReturn(0);
225 }
226 
227 /*@
228   PetscDualSpaceDestroy - Destroys a PetscDualSpace object
229 
230   Collective on PetscDualSpace
231 
232   Input Parameter:
233 . sp - the PetscDualSpace object to destroy
234 
235   Level: developer
236 
237 .seealso PetscDualSpaceView()
238 @*/
239 PetscErrorCode PetscDualSpaceDestroy(PetscDualSpace *sp)
240 {
241   PetscInt       dim, f;
242   PetscErrorCode ierr;
243 
244   PetscFunctionBegin;
245   if (!*sp) PetscFunctionReturn(0);
246   PetscValidHeaderSpecific((*sp), PETSCDUALSPACE_CLASSID, 1);
247 
248   if (--((PetscObject)(*sp))->refct > 0) {*sp = 0; PetscFunctionReturn(0);}
249   ((PetscObject) (*sp))->refct = 0;
250 
251   ierr = PetscDualSpaceGetDimension(*sp, &dim);CHKERRQ(ierr);
252   for (f = 0; f < dim; ++f) {
253     ierr = PetscQuadratureDestroy(&(*sp)->functional[f]);CHKERRQ(ierr);
254   }
255   ierr = PetscFree((*sp)->functional);CHKERRQ(ierr);
256   ierr = PetscQuadratureDestroy(&(*sp)->allPoints);CHKERRQ(ierr);
257   ierr = DMDestroy(&(*sp)->dm);CHKERRQ(ierr);
258 
259   if ((*sp)->ops->destroy) {ierr = (*(*sp)->ops->destroy)(*sp);CHKERRQ(ierr);}
260   ierr = PetscHeaderDestroy(sp);CHKERRQ(ierr);
261   PetscFunctionReturn(0);
262 }
263 
264 /*@
265   PetscDualSpaceCreate - Creates an empty PetscDualSpace object. The type can then be set with PetscDualSpaceSetType().
266 
267   Collective on MPI_Comm
268 
269   Input Parameter:
270 . comm - The communicator for the PetscDualSpace object
271 
272   Output Parameter:
273 . sp - The PetscDualSpace object
274 
275   Level: beginner
276 
277 .seealso: PetscDualSpaceSetType(), PETSCDUALSPACELAGRANGE
278 @*/
279 PetscErrorCode PetscDualSpaceCreate(MPI_Comm comm, PetscDualSpace *sp)
280 {
281   PetscDualSpace s;
282   PetscErrorCode ierr;
283 
284   PetscFunctionBegin;
285   PetscValidPointer(sp, 2);
286   ierr = PetscCitationsRegister(FECitation,&FEcite);CHKERRQ(ierr);
287   *sp  = NULL;
288   ierr = PetscFEInitializePackage();CHKERRQ(ierr);
289 
290   ierr = PetscHeaderCreate(s, PETSCDUALSPACE_CLASSID, "PetscDualSpace", "Dual Space", "PetscDualSpace", comm, PetscDualSpaceDestroy, PetscDualSpaceView);CHKERRQ(ierr);
291 
292   s->order = 0;
293   s->Nc    = 1;
294   s->setupcalled = PETSC_FALSE;
295 
296   *sp = s;
297   PetscFunctionReturn(0);
298 }
299 
300 /*@
301   PetscDualSpaceDuplicate - Creates a duplicate PetscDualSpace object, however it is not setup.
302 
303   Collective on PetscDualSpace
304 
305   Input Parameter:
306 . sp - The original PetscDualSpace
307 
308   Output Parameter:
309 . spNew - The duplicate PetscDualSpace
310 
311   Level: beginner
312 
313 .seealso: PetscDualSpaceCreate(), PetscDualSpaceSetType()
314 @*/
315 PetscErrorCode PetscDualSpaceDuplicate(PetscDualSpace sp, PetscDualSpace *spNew)
316 {
317   PetscErrorCode ierr;
318 
319   PetscFunctionBegin;
320   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
321   PetscValidPointer(spNew, 2);
322   ierr = (*sp->ops->duplicate)(sp, spNew);CHKERRQ(ierr);
323   PetscFunctionReturn(0);
324 }
325 
326 /*@
327   PetscDualSpaceGetDM - Get the DM representing the reference cell
328 
329   Not collective
330 
331   Input Parameter:
332 . sp - The PetscDualSpace
333 
334   Output Parameter:
335 . dm - The reference cell
336 
337   Level: intermediate
338 
339 .seealso: PetscDualSpaceSetDM(), PetscDualSpaceCreate()
340 @*/
341 PetscErrorCode PetscDualSpaceGetDM(PetscDualSpace sp, DM *dm)
342 {
343   PetscFunctionBegin;
344   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
345   PetscValidPointer(dm, 2);
346   *dm = sp->dm;
347   PetscFunctionReturn(0);
348 }
349 
350 /*@
351   PetscDualSpaceSetDM - Get the DM representing the reference cell
352 
353   Not collective
354 
355   Input Parameters:
356 + sp - The PetscDualSpace
357 - dm - The reference cell
358 
359   Level: intermediate
360 
361 .seealso: PetscDualSpaceGetDM(), PetscDualSpaceCreate()
362 @*/
363 PetscErrorCode PetscDualSpaceSetDM(PetscDualSpace sp, DM dm)
364 {
365   PetscErrorCode ierr;
366 
367   PetscFunctionBegin;
368   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
369   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
370   ierr = DMDestroy(&sp->dm);CHKERRQ(ierr);
371   ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);
372   sp->dm = dm;
373   PetscFunctionReturn(0);
374 }
375 
376 /*@
377   PetscDualSpaceGetOrder - Get the order of the dual space
378 
379   Not collective
380 
381   Input Parameter:
382 . sp - The PetscDualSpace
383 
384   Output Parameter:
385 . order - The order
386 
387   Level: intermediate
388 
389 .seealso: PetscDualSpaceSetOrder(), PetscDualSpaceCreate()
390 @*/
391 PetscErrorCode PetscDualSpaceGetOrder(PetscDualSpace sp, PetscInt *order)
392 {
393   PetscFunctionBegin;
394   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
395   PetscValidPointer(order, 2);
396   *order = sp->order;
397   PetscFunctionReturn(0);
398 }
399 
400 /*@
401   PetscDualSpaceSetOrder - Set the order of the dual space
402 
403   Not collective
404 
405   Input Parameters:
406 + sp - The PetscDualSpace
407 - order - The order
408 
409   Level: intermediate
410 
411 .seealso: PetscDualSpaceGetOrder(), PetscDualSpaceCreate()
412 @*/
413 PetscErrorCode PetscDualSpaceSetOrder(PetscDualSpace sp, PetscInt order)
414 {
415   PetscFunctionBegin;
416   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
417   sp->order = order;
418   PetscFunctionReturn(0);
419 }
420 
421 /*@
422   PetscDualSpaceGetNumComponents - Return the number of components for this space
423 
424   Input Parameter:
425 . sp - The PetscDualSpace
426 
427   Output Parameter:
428 . Nc - The number of components
429 
430   Note: A vector space, for example, will have d components, where d is the spatial dimension
431 
432   Level: intermediate
433 
434 .seealso: PetscDualSpaceSetNumComponents(), PetscDualSpaceGetDimension(), PetscDualSpaceCreate(), PetscDualSpace
435 @*/
436 PetscErrorCode PetscDualSpaceGetNumComponents(PetscDualSpace sp, PetscInt *Nc)
437 {
438   PetscFunctionBegin;
439   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
440   PetscValidPointer(Nc, 2);
441   *Nc = sp->Nc;
442   PetscFunctionReturn(0);
443 }
444 
445 /*@
446   PetscDualSpaceSetNumComponents - Set the number of components for this space
447 
448   Input Parameters:
449 + sp - The PetscDualSpace
450 - order - The number of components
451 
452   Level: intermediate
453 
454 .seealso: PetscDualSpaceGetNumComponents(), PetscDualSpaceCreate(), PetscDualSpace
455 @*/
456 PetscErrorCode PetscDualSpaceSetNumComponents(PetscDualSpace sp, PetscInt Nc)
457 {
458   PetscFunctionBegin;
459   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
460   sp->Nc = Nc;
461   PetscFunctionReturn(0);
462 }
463 
464 /*@
465   PetscDualSpaceLagrangeGetTensor - Get the tensor nature of the dual space
466 
467   Not collective
468 
469   Input Parameter:
470 . sp - The PetscDualSpace
471 
472   Output Parameter:
473 . tensor - Whether the dual space has tensor layout (vs. simplicial)
474 
475   Level: intermediate
476 
477 .seealso: PetscDualSpaceLagrangeSetTensor(), PetscDualSpaceCreate()
478 @*/
479 PetscErrorCode PetscDualSpaceLagrangeGetTensor(PetscDualSpace sp, PetscBool *tensor)
480 {
481   PetscErrorCode ierr;
482 
483   PetscFunctionBegin;
484   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
485   PetscValidPointer(tensor, 2);
486   ierr = PetscTryMethod(sp,"PetscDualSpaceLagrangeGetTensor_C",(PetscDualSpace,PetscBool *),(sp,tensor));CHKERRQ(ierr);
487   PetscFunctionReturn(0);
488 }
489 
490 /*@
491   PetscDualSpaceLagrangeSetTensor - Set the tensor nature of the dual space
492 
493   Not collective
494 
495   Input Parameters:
496 + sp - The PetscDualSpace
497 - tensor - Whether the dual space has tensor layout (vs. simplicial)
498 
499   Level: intermediate
500 
501 .seealso: PetscDualSpaceLagrangeGetTensor(), PetscDualSpaceCreate()
502 @*/
503 PetscErrorCode PetscDualSpaceLagrangeSetTensor(PetscDualSpace sp, PetscBool tensor)
504 {
505   PetscErrorCode ierr;
506 
507   PetscFunctionBegin;
508   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
509   ierr = PetscTryMethod(sp,"PetscDualSpaceLagrangeSetTensor_C",(PetscDualSpace,PetscBool),(sp,tensor));CHKERRQ(ierr);
510   PetscFunctionReturn(0);
511 }
512 
513 /*@
514   PetscDualSpaceGetFunctional - Get the i-th basis functional in the dual space
515 
516   Not collective
517 
518   Input Parameters:
519 + sp - The PetscDualSpace
520 - i  - The basis number
521 
522   Output Parameter:
523 . functional - The basis functional
524 
525   Level: intermediate
526 
527 .seealso: PetscDualSpaceGetDimension(), PetscDualSpaceCreate()
528 @*/
529 PetscErrorCode PetscDualSpaceGetFunctional(PetscDualSpace sp, PetscInt i, PetscQuadrature *functional)
530 {
531   PetscInt       dim;
532   PetscErrorCode ierr;
533 
534   PetscFunctionBegin;
535   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
536   PetscValidPointer(functional, 3);
537   ierr = PetscDualSpaceGetDimension(sp, &dim);CHKERRQ(ierr);
538   if ((i < 0) || (i >= dim)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Functional index %d must be in [0, %d)", i, dim);
539   *functional = sp->functional[i];
540   PetscFunctionReturn(0);
541 }
542 
543 /*@
544   PetscDualSpaceGetDimension - Get the dimension of the dual space, i.e. the number of basis functionals
545 
546   Not collective
547 
548   Input Parameter:
549 . sp - The PetscDualSpace
550 
551   Output Parameter:
552 . dim - The dimension
553 
554   Level: intermediate
555 
556 .seealso: PetscDualSpaceGetFunctional(), PetscDualSpaceCreate()
557 @*/
558 PetscErrorCode PetscDualSpaceGetDimension(PetscDualSpace sp, PetscInt *dim)
559 {
560   PetscErrorCode ierr;
561 
562   PetscFunctionBegin;
563   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
564   PetscValidPointer(dim, 2);
565   *dim = 0;
566   if (sp->ops->getdimension) {ierr = (*sp->ops->getdimension)(sp, dim);CHKERRQ(ierr);}
567   PetscFunctionReturn(0);
568 }
569 
570 /*@C
571   PetscDualSpaceGetNumDof - Get the number of degrees of freedom for each spatial (topological) dimension
572 
573   Not collective
574 
575   Input Parameter:
576 . sp - The PetscDualSpace
577 
578   Output Parameter:
579 . numDof - An array of length dim+1 which holds the number of dofs for each dimension
580 
581   Level: intermediate
582 
583 .seealso: PetscDualSpaceGetFunctional(), PetscDualSpaceCreate()
584 @*/
585 PetscErrorCode PetscDualSpaceGetNumDof(PetscDualSpace sp, const PetscInt **numDof)
586 {
587   PetscErrorCode ierr;
588 
589   PetscFunctionBegin;
590   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
591   PetscValidPointer(numDof, 2);
592   ierr = (*sp->ops->getnumdof)(sp, numDof);CHKERRQ(ierr);
593   if (!*numDof) SETERRQ(PetscObjectComm((PetscObject) sp), PETSC_ERR_LIB, "Empty numDof[] returned from dual space implementation");
594   PetscFunctionReturn(0);
595 }
596 
597 PetscErrorCode PetscDualSpaceCreateSection(PetscDualSpace sp, PetscSection *section)
598 {
599   DM             dm;
600   PetscInt       pStart, pEnd, depth, h, offset;
601   const PetscInt *numDof;
602   PetscErrorCode ierr;
603 
604   PetscFunctionBegin;
605   ierr = PetscDualSpaceGetDM(sp,&dm);CHKERRQ(ierr);
606   ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr);
607   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)sp),section);CHKERRQ(ierr);
608   ierr = PetscSectionSetChart(*section,pStart,pEnd);CHKERRQ(ierr);
609   ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr);
610   ierr = PetscDualSpaceGetNumDof(sp,&numDof);CHKERRQ(ierr);
611   for (h = 0; h <= depth; h++) {
612     PetscInt hStart, hEnd, p, dof;
613 
614     ierr = DMPlexGetHeightStratum(dm,h,&hStart,&hEnd);CHKERRQ(ierr);
615     dof = numDof[depth - h];
616     for (p = hStart; p < hEnd; p++) {
617       ierr = PetscSectionSetDof(*section,p,dof);CHKERRQ(ierr);
618     }
619   }
620   ierr = PetscSectionSetUp(*section);CHKERRQ(ierr);
621   for (h = 0, offset = 0; h <= depth; h++) {
622     PetscInt hStart, hEnd, p, dof;
623 
624     ierr = DMPlexGetHeightStratum(dm,h,&hStart,&hEnd);CHKERRQ(ierr);
625     dof = numDof[depth - h];
626     for (p = hStart; p < hEnd; p++) {
627       ierr = PetscSectionGetDof(*section,p,&dof);CHKERRQ(ierr);
628       ierr = PetscSectionSetOffset(*section,p,offset);CHKERRQ(ierr);
629       offset += dof;
630     }
631   }
632   PetscFunctionReturn(0);
633 }
634 
635 /*@
636   PetscDualSpaceCreateReferenceCell - Create a DMPLEX with the appropriate FEM reference cell
637 
638   Collective on PetscDualSpace
639 
640   Input Parameters:
641 + sp      - The PetscDualSpace
642 . dim     - The spatial dimension
643 - simplex - Flag for simplex, otherwise use a tensor-product cell
644 
645   Output Parameter:
646 . refdm - The reference cell
647 
648   Level: advanced
649 
650 .keywords: PetscDualSpace, reference cell
651 .seealso: PetscDualSpaceCreate(), DMPLEX
652 @*/
653 PetscErrorCode PetscDualSpaceCreateReferenceCell(PetscDualSpace sp, PetscInt dim, PetscBool simplex, DM *refdm)
654 {
655   PetscErrorCode ierr;
656 
657   PetscFunctionBeginUser;
658   ierr = DMPlexCreateReferenceCell(PetscObjectComm((PetscObject) sp), dim, simplex, refdm);CHKERRQ(ierr);
659   PetscFunctionReturn(0);
660 }
661 
662 /*@C
663   PetscDualSpaceApply - Apply a functional from the dual space basis to an input function
664 
665   Input Parameters:
666 + sp      - The PetscDualSpace object
667 . f       - The basis functional index
668 . time    - The time
669 . 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)
670 . numComp - The number of components for the function
671 . func    - The input function
672 - ctx     - A context for the function
673 
674   Output Parameter:
675 . value   - numComp output values
676 
677   Note: The calling sequence for the callback func is given by:
678 
679 $ func(PetscInt dim, PetscReal time, const PetscReal x[],
680 $      PetscInt numComponents, PetscScalar values[], void *ctx)
681 
682   Level: developer
683 
684 .seealso: PetscDualSpaceCreate()
685 @*/
686 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)
687 {
688   PetscErrorCode ierr;
689 
690   PetscFunctionBegin;
691   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
692   PetscValidPointer(cgeom, 4);
693   PetscValidPointer(value, 8);
694   ierr = (*sp->ops->apply)(sp, f, time, cgeom, numComp, func, ctx, value);CHKERRQ(ierr);
695   PetscFunctionReturn(0);
696 }
697 
698 /*@C
699   PetscDualSpaceApplyAll - Apply all functionals from the dual space basis to the result of an evaluation at the points returned by PetscDualSpaceGetAllPoints()
700 
701   Input Parameters:
702 + sp        - The PetscDualSpace object
703 - pointEval - Evaluation at the points returned by PetscDualSpaceGetAllPoints()
704 
705   Output Parameter:
706 . spValue   - The values of all dual space functionals
707 
708   Level: developer
709 
710 .seealso: PetscDualSpaceCreate()
711 @*/
712 PetscErrorCode PetscDualSpaceApplyAll(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue)
713 {
714   PetscErrorCode ierr;
715 
716   PetscFunctionBegin;
717   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
718   ierr = (*sp->ops->applyall)(sp, pointEval, spValue);CHKERRQ(ierr);
719   PetscFunctionReturn(0);
720 }
721 
722 /*@C
723   PetscDualSpaceApplyDefault - Apply a functional from the dual space basis to an input function by assuming a point evaluation functional.
724 
725   Input Parameters:
726 + sp    - The PetscDualSpace object
727 . f     - The basis functional index
728 . time  - The time
729 . cgeom - A context with geometric information for this cell, we use v0 (the initial vertex) and J (the Jacobian)
730 . Nc    - The number of components for the function
731 . func  - The input function
732 - ctx   - A context for the function
733 
734   Output Parameter:
735 . value   - The output value
736 
737   Note: The calling sequence for the callback func is given by:
738 
739 $ func(PetscInt dim, PetscReal time, const PetscReal x[],
740 $      PetscInt numComponents, PetscScalar values[], void *ctx)
741 
742 and the idea is to evaluate the functional as an integral
743 
744 $ n(f) = int dx n(x) . f(x)
745 
746 where both n and f have Nc components.
747 
748   Level: developer
749 
750 .seealso: PetscDualSpaceCreate()
751 @*/
752 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)
753 {
754   DM               dm;
755   PetscQuadrature  n;
756   const PetscReal *points, *weights;
757   PetscReal        x[3];
758   PetscScalar     *val;
759   PetscInt         dim, dE, qNc, c, Nq, q;
760   PetscBool        isAffine;
761   PetscErrorCode   ierr;
762 
763   PetscFunctionBegin;
764   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
765   PetscValidPointer(value, 5);
766   ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr);
767   ierr = PetscDualSpaceGetFunctional(sp, f, &n);CHKERRQ(ierr);
768   ierr = PetscQuadratureGetData(n, &dim, &qNc, &Nq, &points, &weights);CHKERRQ(ierr);
769   if (dim != cgeom->dim) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature spatial dimension %D != cell geometry dimension %D", dim, cgeom->dim);
770   if (qNc != Nc) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature components %D != function components %D", qNc, Nc);
771   ierr = DMGetWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr);
772   *value = 0.0;
773   isAffine = cgeom->isAffine;
774   dE = cgeom->dimEmbed;
775   for (q = 0; q < Nq; ++q) {
776     if (isAffine) {
777       CoordinatesRefToReal(dE, cgeom->dim, cgeom->xi, cgeom->v, cgeom->J, &points[q*dim], x);
778       ierr = (*func)(dE, time, x, Nc, val, ctx);CHKERRQ(ierr);
779     } else {
780       ierr = (*func)(dE, time, &cgeom->v[dE*q], Nc, val, ctx);CHKERRQ(ierr);
781     }
782     for (c = 0; c < Nc; ++c) {
783       *value += val[c]*weights[q*Nc+c];
784     }
785   }
786   ierr = DMRestoreWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr);
787   PetscFunctionReturn(0);
788 }
789 
790 /*@C
791   PetscDualSpaceApplyAllDefault - Apply all functionals from the dual space basis to the result of an evaluation at the points returned by PetscDualSpaceGetAllPoints()
792 
793   Input Parameters:
794 + sp        - The PetscDualSpace object
795 - pointEval - Evaluation at the points returned by PetscDualSpaceGetAllPoints()
796 
797   Output Parameter:
798 . spValue   - The values of all dual space functionals
799 
800   Level: developer
801 
802 .seealso: PetscDualSpaceCreate()
803 @*/
804 PetscErrorCode PetscDualSpaceApplyAllDefault(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue)
805 {
806   PetscQuadrature  n;
807   const PetscReal *points, *weights;
808   PetscInt         qNc, c, Nq, q, f, spdim, Nc;
809   PetscInt         offset;
810   PetscErrorCode   ierr;
811 
812   PetscFunctionBegin;
813   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
814   PetscValidScalarPointer(pointEval, 2);
815   PetscValidScalarPointer(spValue, 5);
816   ierr = PetscDualSpaceGetDimension(sp, &spdim);CHKERRQ(ierr);
817   ierr = PetscDualSpaceGetNumComponents(sp, &Nc);CHKERRQ(ierr);
818   for (f = 0, offset = 0; f < spdim; f++) {
819     ierr = PetscDualSpaceGetFunctional(sp, f, &n);CHKERRQ(ierr);
820     ierr = PetscQuadratureGetData(n, NULL, &qNc, &Nq, &points, &weights);CHKERRQ(ierr);
821     if (qNc != Nc) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature components %D != function components %D", qNc, Nc);
822     spValue[f] = 0.0;
823     for (q = 0; q < Nq; ++q) {
824       for (c = 0; c < Nc; ++c) {
825         spValue[f] += pointEval[offset++]*weights[q*Nc+c];
826       }
827     }
828   }
829   PetscFunctionReturn(0);
830 }
831 
832 PetscErrorCode PetscDualSpaceGetAllPoints(PetscDualSpace sp, PetscQuadrature *allPoints)
833 {
834   PetscErrorCode ierr;
835 
836   PetscFunctionBegin;
837   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
838   PetscValidPointer(allPoints,2);
839   if (!sp->allPoints && sp->ops->createallpoints) {
840     ierr = (*sp->ops->createallpoints)(sp,&sp->allPoints);CHKERRQ(ierr);
841   }
842   *allPoints = sp->allPoints;
843   PetscFunctionReturn(0);
844 }
845 
846 PetscErrorCode PetscDualSpaceCreateAllPointsDefault(PetscDualSpace sp, PetscQuadrature *allPoints)
847 {
848   PetscInt        spdim;
849   PetscInt        numPoints, offset;
850   PetscReal       *points;
851   PetscInt        f, dim;
852   PetscQuadrature q;
853   PetscErrorCode  ierr;
854 
855   PetscFunctionBegin;
856   ierr = PetscDualSpaceGetDimension(sp,&spdim);CHKERRQ(ierr);
857   if (!spdim) {
858     ierr = PetscQuadratureCreate(PETSC_COMM_SELF,allPoints);CHKERRQ(ierr);
859     ierr = PetscQuadratureSetData(*allPoints,0,0,0,NULL,NULL);CHKERRQ(ierr);
860   }
861   ierr = PetscDualSpaceGetFunctional(sp,0,&q);CHKERRQ(ierr);
862   ierr = PetscQuadratureGetData(q,&dim,NULL,&numPoints,NULL,NULL);CHKERRQ(ierr);
863   for (f = 1; f < spdim; f++) {
864     PetscInt Np;
865 
866     ierr = PetscDualSpaceGetFunctional(sp,f,&q);CHKERRQ(ierr);
867     ierr = PetscQuadratureGetData(q,NULL,NULL,&Np,NULL,NULL);CHKERRQ(ierr);
868     numPoints += Np;
869   }
870   ierr = PetscMalloc1(dim*numPoints,&points);CHKERRQ(ierr);
871   for (f = 0, offset = 0; f < spdim; f++) {
872     const PetscReal *p;
873     PetscInt        Np, i;
874 
875     ierr = PetscDualSpaceGetFunctional(sp,f,&q);CHKERRQ(ierr);
876     ierr = PetscQuadratureGetData(q,NULL,NULL,&Np,&p,NULL);CHKERRQ(ierr);
877     for (i = 0; i < Np * dim; i++) {
878       points[offset + i] = p[i];
879     }
880     offset += Np * dim;
881   }
882   ierr = PetscQuadratureCreate(PETSC_COMM_SELF,allPoints);CHKERRQ(ierr);
883   ierr = PetscQuadratureSetData(*allPoints,dim,0,numPoints,points,NULL);CHKERRQ(ierr);
884   PetscFunctionReturn(0);
885 }
886 
887 /*@C
888   PetscDualSpaceApplyFVM - Apply a functional from the dual space basis to an input function by assuming a point evaluation functional at the cell centroid.
889 
890   Input Parameters:
891 + sp    - The PetscDualSpace object
892 . f     - The basis functional index
893 . time  - The time
894 . cgeom - A context with geometric information for this cell, we currently just use the centroid
895 . Nc    - The number of components for the function
896 . func  - The input function
897 - ctx   - A context for the function
898 
899   Output Parameter:
900 . value - The output value (scalar)
901 
902   Note: The calling sequence for the callback func is given by:
903 
904 $ func(PetscInt dim, PetscReal time, const PetscReal x[],
905 $      PetscInt numComponents, PetscScalar values[], void *ctx)
906 
907 and the idea is to evaluate the functional as an integral
908 
909 $ n(f) = int dx n(x) . f(x)
910 
911 where both n and f have Nc components.
912 
913   Level: developer
914 
915 .seealso: PetscDualSpaceCreate()
916 @*/
917 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)
918 {
919   DM               dm;
920   PetscQuadrature  n;
921   const PetscReal *points, *weights;
922   PetscScalar     *val;
923   PetscInt         dimEmbed, qNc, c, Nq, q;
924   PetscErrorCode   ierr;
925 
926   PetscFunctionBegin;
927   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
928   PetscValidPointer(value, 5);
929   ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr);
930   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
931   ierr = PetscDualSpaceGetFunctional(sp, f, &n);CHKERRQ(ierr);
932   ierr = PetscQuadratureGetData(n, NULL, &qNc, &Nq, &points, &weights);CHKERRQ(ierr);
933   if (qNc != Nc) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature components %D != function components %D", qNc, Nc);
934   ierr = DMGetWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr);
935   *value = 0.;
936   for (q = 0; q < Nq; ++q) {
937     ierr = (*func)(dimEmbed, time, cgeom->centroid, Nc, val, ctx);CHKERRQ(ierr);
938     for (c = 0; c < Nc; ++c) {
939       *value += val[c]*weights[q*Nc+c];
940     }
941   }
942   ierr = DMRestoreWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr);
943   PetscFunctionReturn(0);
944 }
945 
946 /*@
947   PetscDualSpaceGetHeightSubspace - Get the subset of the dual space basis that is supported on a mesh point of a
948   given height.  This assumes that the reference cell is symmetric over points of this height.
949 
950   If the dual space is not defined on mesh points of the given height (e.g. if the space is discontinuous and
951   pointwise values are not defined on the element boundaries), or if the implementation of PetscDualSpace does not
952   support extracting subspaces, then NULL is returned.
953 
954   This does not increment the reference count on the returned dual space, and the user should not destroy it.
955 
956   Not collective
957 
958   Input Parameters:
959 + sp - the PetscDualSpace object
960 - height - the height of the mesh point for which the subspace is desired
961 
962   Output Parameter:
963 . subsp - the subspace.  Note that the functionals in the subspace are with respect to the intrinsic geometry of the
964   point, which will be of lesser dimension if height > 0.
965 
966   Level: advanced
967 
968 .seealso: PetscSpaceGetHeightSubspace(), PetscDualSpace
969 @*/
970 PetscErrorCode PetscDualSpaceGetHeightSubspace(PetscDualSpace sp, PetscInt height, PetscDualSpace *subsp)
971 {
972   PetscErrorCode ierr;
973 
974   PetscFunctionBegin;
975   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
976   PetscValidPointer(subsp, 3);
977   *subsp = NULL;
978   if (sp->ops->getheightsubspace) {
979     ierr = (*sp->ops->getheightsubspace)(sp, height, subsp);CHKERRQ(ierr);
980   }
981   PetscFunctionReturn(0);
982 }
983 
984 /*@
985   PetscDualSpaceGetPointSubspace - Get the subset of the dual space basis that is supported on a particular mesh point.
986 
987   If the dual space is not defined on the mesh point (e.g. if the space is discontinuous and pointwise values are not
988   defined on the element boundaries), or if the implementation of PetscDualSpace does not support extracting
989   subspaces, then NULL is returned.
990 
991   This does not increment the reference count on the returned dual space, and the user should not destroy it.
992 
993   Not collective
994 
995   Input Parameters:
996 + sp - the PetscDualSpace object
997 - point - the point (in the dual space's DM) for which the subspace is desired
998 
999   Output Parameters:
1000   bdsp - the subspace.  Note that the functionals in the subspace are with respect to the intrinsic geometry of the
1001   point, which will be of lesser dimension if height > 0.
1002 
1003   Level: advanced
1004 
1005 .seealso: PetscDualSpace
1006 @*/
1007 PetscErrorCode PetscDualSpaceGetPointSubspace(PetscDualSpace sp, PetscInt point, PetscDualSpace *bdsp)
1008 {
1009   PetscErrorCode ierr;
1010 
1011   PetscFunctionBegin;
1012   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1013   PetscValidPointer(bdsp,2);
1014   *bdsp = NULL;
1015   if (sp->ops->getpointsubspace) {
1016     ierr = (*sp->ops->getpointsubspace)(sp,point,bdsp);CHKERRQ(ierr);
1017   } else if (sp->ops->getheightsubspace) {
1018     DM       dm;
1019     DMLabel  label;
1020     PetscInt dim, depth, height;
1021 
1022     ierr = PetscDualSpaceGetDM(sp,&dm);CHKERRQ(ierr);
1023     ierr = DMPlexGetDepth(dm,&dim);CHKERRQ(ierr);
1024     ierr = DMPlexGetDepthLabel(dm,&label);CHKERRQ(ierr);
1025     ierr = DMLabelGetValue(label,point,&depth);CHKERRQ(ierr);
1026     height = dim - depth;
1027     ierr = (*sp->ops->getheightsubspace)(sp,height,bdsp);CHKERRQ(ierr);
1028   }
1029   PetscFunctionReturn(0);
1030 }
1031 
1032