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