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