xref: /petsc/src/dm/dt/dualspace/interface/dualspace.c (revision 6f905325d484e237113fc7e4a9e8f803271bb6e3)
120cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
220cf1dd8SToby Isaac #include <petscdmplex.h>
320cf1dd8SToby Isaac 
420cf1dd8SToby Isaac PetscClassId PETSCDUALSPACE_CLASSID = 0;
520cf1dd8SToby Isaac 
620cf1dd8SToby Isaac PetscFunctionList PetscDualSpaceList              = NULL;
720cf1dd8SToby Isaac PetscBool         PetscDualSpaceRegisterAllCalled = PETSC_FALSE;
820cf1dd8SToby Isaac 
955cc6565SMatthew G. Knepley const char *const PetscDualSpaceReferenceCells[] = {"SIMPLEX", "TENSOR", "PetscDualSpaceReferenceCell", "PETSCDUALSPACE_REFCELL_",0};
1055cc6565SMatthew G. Knepley 
11*6f905325SMatthew G. Knepley /*
12*6f905325SMatthew G. Knepley   PetscDualSpaceLatticePointLexicographic_Internal - Returns all tuples of size 'len' with nonnegative integers that sum up to at most 'max'.
13*6f905325SMatthew G. Knepley                                                      Ordering is lexicographic with lowest index as least significant in ordering.
14*6f905325SMatthew G. Knepley                                                      e.g. for len == 2 and max == 2, this will return, in order, {0,0}, {1,0}, {2,0}, {0,1}, {1,1}, {2,0}.
15*6f905325SMatthew G. Knepley 
16*6f905325SMatthew G. Knepley   Input Parameters:
17*6f905325SMatthew G. Knepley + len - The length of the tuple
18*6f905325SMatthew G. Knepley . max - The maximum sum
19*6f905325SMatthew G. Knepley - tup - A tuple of length len+1: tup[len] > 0 indicates a stopping condition
20*6f905325SMatthew G. Knepley 
21*6f905325SMatthew G. Knepley   Output Parameter:
22*6f905325SMatthew G. Knepley . tup - A tuple of len integers whos sum is at most 'max'
23*6f905325SMatthew G. Knepley 
24*6f905325SMatthew G. Knepley   Level: developer
25*6f905325SMatthew G. Knepley 
26*6f905325SMatthew G. Knepley .seealso: PetscDualSpaceTensorPointLexicographic_Internal()
27*6f905325SMatthew G. Knepley */
28*6f905325SMatthew G. Knepley PetscErrorCode PetscDualSpaceLatticePointLexicographic_Internal(PetscInt len, PetscInt max, PetscInt tup[])
29*6f905325SMatthew G. Knepley {
30*6f905325SMatthew G. Knepley   PetscFunctionBegin;
31*6f905325SMatthew G. Knepley   while (len--) {
32*6f905325SMatthew G. Knepley     max -= tup[len];
33*6f905325SMatthew G. Knepley     if (!max) {
34*6f905325SMatthew G. Knepley       tup[len] = 0;
35*6f905325SMatthew G. Knepley       break;
36*6f905325SMatthew G. Knepley     }
37*6f905325SMatthew G. Knepley   }
38*6f905325SMatthew G. Knepley   tup[++len]++;
39*6f905325SMatthew G. Knepley   PetscFunctionReturn(0);
40*6f905325SMatthew G. Knepley }
41*6f905325SMatthew G. Knepley 
42*6f905325SMatthew G. Knepley /*
43*6f905325SMatthew G. Knepley   PetscDualSpaceTensorPointLexicographic_Internal - Returns all tuples of size 'len' with nonnegative integers that are all less than or equal to 'max'.
44*6f905325SMatthew G. Knepley                                                     Ordering is lexicographic with lowest index as least significant in ordering.
45*6f905325SMatthew G. Knepley                                                     e.g. for len == 2 and max == 2, this will return, in order, {0,0}, {1,0}, {2,0}, {0,1}, {1,1}, {2,1}, {0,2}, {1,2}, {2,2}.
46*6f905325SMatthew G. Knepley 
47*6f905325SMatthew G. Knepley   Input Parameters:
48*6f905325SMatthew G. Knepley + len - The length of the tuple
49*6f905325SMatthew G. Knepley . max - The maximum value
50*6f905325SMatthew G. Knepley - tup - A tuple of length len+1: tup[len] > 0 indicates a stopping condition
51*6f905325SMatthew G. Knepley 
52*6f905325SMatthew G. Knepley   Output Parameter:
53*6f905325SMatthew G. Knepley . tup - A tuple of len integers whos sum is at most 'max'
54*6f905325SMatthew G. Knepley 
55*6f905325SMatthew G. Knepley   Level: developer
56*6f905325SMatthew G. Knepley 
57*6f905325SMatthew G. Knepley .seealso: PetscDualSpaceLatticePointLexicographic_Internal()
58*6f905325SMatthew G. Knepley */
59*6f905325SMatthew G. Knepley PetscErrorCode PetscDualSpaceTensorPointLexicographic_Internal(PetscInt len, PetscInt max, PetscInt tup[])
60*6f905325SMatthew G. Knepley {
61*6f905325SMatthew G. Knepley   PetscInt       i;
62*6f905325SMatthew G. Knepley 
63*6f905325SMatthew G. Knepley   PetscFunctionBegin;
64*6f905325SMatthew G. Knepley   for (i = 0; i < len; i++) {
65*6f905325SMatthew G. Knepley     if (tup[i] < max) {
66*6f905325SMatthew G. Knepley       break;
67*6f905325SMatthew G. Knepley     } else {
68*6f905325SMatthew G. Knepley       tup[i] = 0;
69*6f905325SMatthew G. Knepley     }
70*6f905325SMatthew G. Knepley   }
71*6f905325SMatthew G. Knepley   tup[i]++;
72*6f905325SMatthew G. Knepley   PetscFunctionReturn(0);
73*6f905325SMatthew G. Knepley }
74*6f905325SMatthew G. Knepley 
7520cf1dd8SToby Isaac /*@C
7620cf1dd8SToby Isaac   PetscDualSpaceRegister - Adds a new PetscDualSpace implementation
7720cf1dd8SToby Isaac 
7820cf1dd8SToby Isaac   Not Collective
7920cf1dd8SToby Isaac 
8020cf1dd8SToby Isaac   Input Parameters:
8120cf1dd8SToby Isaac + name        - The name of a new user-defined creation routine
8220cf1dd8SToby Isaac - create_func - The creation routine itself
8320cf1dd8SToby Isaac 
8420cf1dd8SToby Isaac   Notes:
8520cf1dd8SToby Isaac   PetscDualSpaceRegister() may be called multiple times to add several user-defined PetscDualSpaces
8620cf1dd8SToby Isaac 
8720cf1dd8SToby Isaac   Sample usage:
8820cf1dd8SToby Isaac .vb
8920cf1dd8SToby Isaac     PetscDualSpaceRegister("my_space", MyPetscDualSpaceCreate);
9020cf1dd8SToby Isaac .ve
9120cf1dd8SToby Isaac 
9220cf1dd8SToby Isaac   Then, your PetscDualSpace type can be chosen with the procedural interface via
9320cf1dd8SToby Isaac .vb
9420cf1dd8SToby Isaac     PetscDualSpaceCreate(MPI_Comm, PetscDualSpace *);
9520cf1dd8SToby Isaac     PetscDualSpaceSetType(PetscDualSpace, "my_dual_space");
9620cf1dd8SToby Isaac .ve
9720cf1dd8SToby Isaac    or at runtime via the option
9820cf1dd8SToby Isaac .vb
9920cf1dd8SToby Isaac     -petscdualspace_type my_dual_space
10020cf1dd8SToby Isaac .ve
10120cf1dd8SToby Isaac 
10220cf1dd8SToby Isaac   Level: advanced
10320cf1dd8SToby Isaac 
10420cf1dd8SToby Isaac .keywords: PetscDualSpace, register
10520cf1dd8SToby Isaac .seealso: PetscDualSpaceRegisterAll(), PetscDualSpaceRegisterDestroy()
10620cf1dd8SToby Isaac 
10720cf1dd8SToby Isaac @*/
10820cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceRegister(const char sname[], PetscErrorCode (*function)(PetscDualSpace))
10920cf1dd8SToby Isaac {
11020cf1dd8SToby Isaac   PetscErrorCode ierr;
11120cf1dd8SToby Isaac 
11220cf1dd8SToby Isaac   PetscFunctionBegin;
11320cf1dd8SToby Isaac   ierr = PetscFunctionListAdd(&PetscDualSpaceList, sname, function);CHKERRQ(ierr);
11420cf1dd8SToby Isaac   PetscFunctionReturn(0);
11520cf1dd8SToby Isaac }
11620cf1dd8SToby Isaac 
11720cf1dd8SToby Isaac /*@C
11820cf1dd8SToby Isaac   PetscDualSpaceSetType - Builds a particular PetscDualSpace
11920cf1dd8SToby Isaac 
12020cf1dd8SToby Isaac   Collective on PetscDualSpace
12120cf1dd8SToby Isaac 
12220cf1dd8SToby Isaac   Input Parameters:
12320cf1dd8SToby Isaac + sp   - The PetscDualSpace object
12420cf1dd8SToby Isaac - name - The kind of space
12520cf1dd8SToby Isaac 
12620cf1dd8SToby Isaac   Options Database Key:
12720cf1dd8SToby Isaac . -petscdualspace_type <type> - Sets the PetscDualSpace type; use -help for a list of available types
12820cf1dd8SToby Isaac 
12920cf1dd8SToby Isaac   Level: intermediate
13020cf1dd8SToby Isaac 
13120cf1dd8SToby Isaac .keywords: PetscDualSpace, set, type
13220cf1dd8SToby Isaac .seealso: PetscDualSpaceGetType(), PetscDualSpaceCreate()
13320cf1dd8SToby Isaac @*/
13420cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetType(PetscDualSpace sp, PetscDualSpaceType name)
13520cf1dd8SToby Isaac {
13620cf1dd8SToby Isaac   PetscErrorCode (*r)(PetscDualSpace);
13720cf1dd8SToby Isaac   PetscBool      match;
13820cf1dd8SToby Isaac   PetscErrorCode ierr;
13920cf1dd8SToby Isaac 
14020cf1dd8SToby Isaac   PetscFunctionBegin;
14120cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
14220cf1dd8SToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) sp, name, &match);CHKERRQ(ierr);
14320cf1dd8SToby Isaac   if (match) PetscFunctionReturn(0);
14420cf1dd8SToby Isaac 
14520cf1dd8SToby Isaac   if (!PetscDualSpaceRegisterAllCalled) {ierr = PetscDualSpaceRegisterAll();CHKERRQ(ierr);}
14620cf1dd8SToby Isaac   ierr = PetscFunctionListFind(PetscDualSpaceList, name, &r);CHKERRQ(ierr);
14720cf1dd8SToby Isaac   if (!r) SETERRQ1(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscDualSpace type: %s", name);
14820cf1dd8SToby Isaac 
14920cf1dd8SToby Isaac   if (sp->ops->destroy) {
15020cf1dd8SToby Isaac     ierr             = (*sp->ops->destroy)(sp);CHKERRQ(ierr);
15120cf1dd8SToby Isaac     sp->ops->destroy = NULL;
15220cf1dd8SToby Isaac   }
15320cf1dd8SToby Isaac   ierr = (*r)(sp);CHKERRQ(ierr);
15420cf1dd8SToby Isaac   ierr = PetscObjectChangeTypeName((PetscObject) sp, name);CHKERRQ(ierr);
15520cf1dd8SToby Isaac   PetscFunctionReturn(0);
15620cf1dd8SToby Isaac }
15720cf1dd8SToby Isaac 
15820cf1dd8SToby Isaac /*@C
15920cf1dd8SToby Isaac   PetscDualSpaceGetType - Gets the PetscDualSpace type name (as a string) from the object.
16020cf1dd8SToby Isaac 
16120cf1dd8SToby Isaac   Not Collective
16220cf1dd8SToby Isaac 
16320cf1dd8SToby Isaac   Input Parameter:
16420cf1dd8SToby Isaac . sp  - The PetscDualSpace
16520cf1dd8SToby Isaac 
16620cf1dd8SToby Isaac   Output Parameter:
16720cf1dd8SToby Isaac . name - The PetscDualSpace type name
16820cf1dd8SToby Isaac 
16920cf1dd8SToby Isaac   Level: intermediate
17020cf1dd8SToby Isaac 
17120cf1dd8SToby Isaac .keywords: PetscDualSpace, get, type, name
17220cf1dd8SToby Isaac .seealso: PetscDualSpaceSetType(), PetscDualSpaceCreate()
17320cf1dd8SToby Isaac @*/
17420cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetType(PetscDualSpace sp, PetscDualSpaceType *name)
17520cf1dd8SToby Isaac {
17620cf1dd8SToby Isaac   PetscErrorCode ierr;
17720cf1dd8SToby Isaac 
17820cf1dd8SToby Isaac   PetscFunctionBegin;
17920cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
18020cf1dd8SToby Isaac   PetscValidPointer(name, 2);
18120cf1dd8SToby Isaac   if (!PetscDualSpaceRegisterAllCalled) {
18220cf1dd8SToby Isaac     ierr = PetscDualSpaceRegisterAll();CHKERRQ(ierr);
18320cf1dd8SToby Isaac   }
18420cf1dd8SToby Isaac   *name = ((PetscObject) sp)->type_name;
18520cf1dd8SToby Isaac   PetscFunctionReturn(0);
18620cf1dd8SToby Isaac }
18720cf1dd8SToby Isaac 
188221d6281SMatthew G. Knepley static PetscErrorCode PetscDualSpaceView_ASCII(PetscDualSpace sp, PetscViewer v)
189221d6281SMatthew G. Knepley {
190221d6281SMatthew G. Knepley   PetscViewerFormat format;
191221d6281SMatthew G. Knepley   PetscInt          pdim, f;
192221d6281SMatthew G. Knepley   PetscErrorCode    ierr;
193221d6281SMatthew G. Knepley 
194221d6281SMatthew G. Knepley   PetscFunctionBegin;
195221d6281SMatthew G. Knepley   ierr = PetscDualSpaceGetDimension(sp, &pdim);CHKERRQ(ierr);
196221d6281SMatthew G. Knepley   ierr = PetscObjectPrintClassNamePrefixType((PetscObject) sp, v);CHKERRQ(ierr);
197221d6281SMatthew G. Knepley   ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr);
198221d6281SMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(v, "Dual space with %D components, size %D\n", sp->Nc, pdim);CHKERRQ(ierr);
199221d6281SMatthew G. Knepley   if (sp->ops->view) {ierr = (*sp->ops->view)(sp, v);CHKERRQ(ierr);}
200221d6281SMatthew G. Knepley   ierr = PetscViewerGetFormat(v, &format);CHKERRQ(ierr);
201221d6281SMatthew G. Knepley   if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
202221d6281SMatthew G. Knepley     ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr);
203221d6281SMatthew G. Knepley     for (f = 0; f < pdim; ++f) {
204221d6281SMatthew G. Knepley       ierr = PetscViewerASCIIPrintf(v, "Dual basis vector %D\n", f);CHKERRQ(ierr);
205221d6281SMatthew G. Knepley       ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr);
206221d6281SMatthew G. Knepley       ierr = PetscQuadratureView(sp->functional[f], v);CHKERRQ(ierr);
207221d6281SMatthew G. Knepley       ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr);
208221d6281SMatthew G. Knepley     }
209221d6281SMatthew G. Knepley     ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr);
210221d6281SMatthew G. Knepley   }
211221d6281SMatthew G. Knepley   ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr);
212221d6281SMatthew G. Knepley   PetscFunctionReturn(0);
213221d6281SMatthew G. Knepley }
214221d6281SMatthew G. Knepley 
21520cf1dd8SToby Isaac /*@
21620cf1dd8SToby Isaac   PetscDualSpaceView - Views a PetscDualSpace
21720cf1dd8SToby Isaac 
21820cf1dd8SToby Isaac   Collective on PetscDualSpace
21920cf1dd8SToby Isaac 
22020cf1dd8SToby Isaac   Input Parameter:
22120cf1dd8SToby Isaac + sp - the PetscDualSpace object to view
22220cf1dd8SToby Isaac - v  - the viewer
22320cf1dd8SToby Isaac 
22420cf1dd8SToby Isaac   Level: developer
22520cf1dd8SToby Isaac 
22620cf1dd8SToby Isaac .seealso PetscDualSpaceDestroy()
22720cf1dd8SToby Isaac @*/
22820cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceView(PetscDualSpace sp, PetscViewer v)
22920cf1dd8SToby Isaac {
230d9bac1caSLisandro Dalcin   PetscBool      iascii;
23120cf1dd8SToby Isaac   PetscErrorCode ierr;
23220cf1dd8SToby Isaac 
23320cf1dd8SToby Isaac   PetscFunctionBegin;
23420cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
235d9bac1caSLisandro Dalcin   if (v) PetscValidHeaderSpecific(v, PETSC_VIEWER_CLASSID, 2);
23620cf1dd8SToby Isaac   if (!v) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) sp), &v);CHKERRQ(ierr);}
237d9bac1caSLisandro Dalcin   ierr = PetscObjectTypeCompare((PetscObject) v, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
238221d6281SMatthew G. Knepley   if (iascii) {ierr = PetscDualSpaceView_ASCII(sp, v);CHKERRQ(ierr);}
23920cf1dd8SToby Isaac   PetscFunctionReturn(0);
24020cf1dd8SToby Isaac }
24120cf1dd8SToby Isaac 
24220cf1dd8SToby Isaac /*@
24320cf1dd8SToby Isaac   PetscDualSpaceSetFromOptions - sets parameters in a PetscDualSpace from the options database
24420cf1dd8SToby Isaac 
24520cf1dd8SToby Isaac   Collective on PetscDualSpace
24620cf1dd8SToby Isaac 
24720cf1dd8SToby Isaac   Input Parameter:
24820cf1dd8SToby Isaac . sp - the PetscDualSpace object to set options for
24920cf1dd8SToby Isaac 
25020cf1dd8SToby Isaac   Options Database:
2517be5e748SToby Isaac . -petscspace_degree the approximation order of the space
25220cf1dd8SToby Isaac 
25320cf1dd8SToby Isaac   Level: developer
25420cf1dd8SToby Isaac 
25520cf1dd8SToby Isaac .seealso PetscDualSpaceView()
25620cf1dd8SToby Isaac @*/
25720cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetFromOptions(PetscDualSpace sp)
25820cf1dd8SToby Isaac {
259063ee4adSMatthew G. Knepley   PetscDualSpaceReferenceCell refCell = PETSCDUALSPACE_REFCELL_SIMPLEX;
260063ee4adSMatthew G. Knepley   PetscInt                    refDim  = 0;
261063ee4adSMatthew G. Knepley   PetscBool                   flg;
26220cf1dd8SToby Isaac   const char                 *defaultType;
26320cf1dd8SToby Isaac   char                        name[256];
26420cf1dd8SToby Isaac   PetscErrorCode              ierr;
26520cf1dd8SToby Isaac 
26620cf1dd8SToby Isaac   PetscFunctionBegin;
26720cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
26820cf1dd8SToby Isaac   if (!((PetscObject) sp)->type_name) {
26920cf1dd8SToby Isaac     defaultType = PETSCDUALSPACELAGRANGE;
27020cf1dd8SToby Isaac   } else {
27120cf1dd8SToby Isaac     defaultType = ((PetscObject) sp)->type_name;
27220cf1dd8SToby Isaac   }
27320cf1dd8SToby Isaac   if (!PetscSpaceRegisterAllCalled) {ierr = PetscSpaceRegisterAll();CHKERRQ(ierr);}
27420cf1dd8SToby Isaac 
27520cf1dd8SToby Isaac   ierr = PetscObjectOptionsBegin((PetscObject) sp);CHKERRQ(ierr);
27620cf1dd8SToby Isaac   ierr = PetscOptionsFList("-petscdualspace_type", "Dual space", "PetscDualSpaceSetType", PetscDualSpaceList, defaultType, name, 256, &flg);CHKERRQ(ierr);
27720cf1dd8SToby Isaac   if (flg) {
27820cf1dd8SToby Isaac     ierr = PetscDualSpaceSetType(sp, name);CHKERRQ(ierr);
27920cf1dd8SToby Isaac   } else if (!((PetscObject) sp)->type_name) {
28020cf1dd8SToby Isaac     ierr = PetscDualSpaceSetType(sp, defaultType);CHKERRQ(ierr);
28120cf1dd8SToby Isaac   }
2827be5e748SToby Isaac   ierr = PetscOptionsInt("-petscdualspace_degree", "The approximation order", "PetscDualSpaceSetOrder", sp->order, &sp->order, NULL);CHKERRQ(ierr);
28320cf1dd8SToby Isaac   ierr = PetscOptionsInt("-petscdualspace_components", "The number of components", "PetscDualSpaceSetNumComponents", sp->Nc, &sp->Nc, NULL);CHKERRQ(ierr);
28420cf1dd8SToby Isaac   if (sp->ops->setfromoptions) {
28520cf1dd8SToby Isaac     ierr = (*sp->ops->setfromoptions)(PetscOptionsObject,sp);CHKERRQ(ierr);
28620cf1dd8SToby Isaac   }
287063ee4adSMatthew G. Knepley   ierr = PetscOptionsInt("-petscdualspace_refdim", "The spatial dimension of the reference cell", "PetscDualSpaceSetReferenceCell", refDim, &refDim, NULL);CHKERRQ(ierr);
288063ee4adSMatthew G. Knepley   ierr = PetscOptionsEnum("-petscdualspace_refcell", "Reference cell", "PetscDualSpaceSetReferenceCell", PetscDualSpaceReferenceCells, (PetscEnum) refCell, (PetscEnum *) &refCell, &flg);CHKERRQ(ierr);
289063ee4adSMatthew G. Knepley   if (flg) {
290063ee4adSMatthew G. Knepley     DM K;
291063ee4adSMatthew G. Knepley 
292063ee4adSMatthew G. Knepley     if (!refDim) SETERRQ(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_INCOMP, "Reference cell specified without a dimension. Use -petscdualspace_refdim.");
293063ee4adSMatthew G. Knepley     ierr = PetscDualSpaceCreateReferenceCell(sp, refDim, refCell == PETSCDUALSPACE_REFCELL_SIMPLEX ? PETSC_TRUE : PETSC_FALSE, &K);CHKERRQ(ierr);
294063ee4adSMatthew G. Knepley     ierr = PetscDualSpaceSetDM(sp, K);CHKERRQ(ierr);
295063ee4adSMatthew G. Knepley     ierr = DMDestroy(&K);CHKERRQ(ierr);
296063ee4adSMatthew G. Knepley   }
297063ee4adSMatthew G. Knepley 
29820cf1dd8SToby Isaac   /* process any options handlers added with PetscObjectAddOptionsHandler() */
29920cf1dd8SToby Isaac   ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) sp);CHKERRQ(ierr);
30020cf1dd8SToby Isaac   ierr = PetscOptionsEnd();CHKERRQ(ierr);
301063ee4adSMatthew G. Knepley   sp->setfromoptionscalled = PETSC_TRUE;
30220cf1dd8SToby Isaac   PetscFunctionReturn(0);
30320cf1dd8SToby Isaac }
30420cf1dd8SToby Isaac 
30520cf1dd8SToby Isaac /*@
30620cf1dd8SToby Isaac   PetscDualSpaceSetUp - Construct a basis for the PetscDualSpace
30720cf1dd8SToby Isaac 
30820cf1dd8SToby Isaac   Collective on PetscDualSpace
30920cf1dd8SToby Isaac 
31020cf1dd8SToby Isaac   Input Parameter:
31120cf1dd8SToby Isaac . sp - the PetscDualSpace object to setup
31220cf1dd8SToby Isaac 
31320cf1dd8SToby Isaac   Level: developer
31420cf1dd8SToby Isaac 
31520cf1dd8SToby Isaac .seealso PetscDualSpaceView(), PetscDualSpaceDestroy()
31620cf1dd8SToby Isaac @*/
31720cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetUp(PetscDualSpace sp)
31820cf1dd8SToby Isaac {
31920cf1dd8SToby Isaac   PetscErrorCode ierr;
32020cf1dd8SToby Isaac 
32120cf1dd8SToby Isaac   PetscFunctionBegin;
32220cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
32320cf1dd8SToby Isaac   if (sp->setupcalled) PetscFunctionReturn(0);
32420cf1dd8SToby Isaac   sp->setupcalled = PETSC_TRUE;
32520cf1dd8SToby Isaac   if (sp->ops->setup) {ierr = (*sp->ops->setup)(sp);CHKERRQ(ierr);}
326063ee4adSMatthew G. Knepley   if (sp->setfromoptionscalled) {ierr = PetscDualSpaceViewFromOptions(sp, NULL, "-petscdualspace_view");CHKERRQ(ierr);}
32720cf1dd8SToby Isaac   PetscFunctionReturn(0);
32820cf1dd8SToby Isaac }
32920cf1dd8SToby Isaac 
33020cf1dd8SToby Isaac /*@
33120cf1dd8SToby Isaac   PetscDualSpaceDestroy - Destroys a PetscDualSpace object
33220cf1dd8SToby Isaac 
33320cf1dd8SToby Isaac   Collective on PetscDualSpace
33420cf1dd8SToby Isaac 
33520cf1dd8SToby Isaac   Input Parameter:
33620cf1dd8SToby Isaac . sp - the PetscDualSpace object to destroy
33720cf1dd8SToby Isaac 
33820cf1dd8SToby Isaac   Level: developer
33920cf1dd8SToby Isaac 
34020cf1dd8SToby Isaac .seealso PetscDualSpaceView()
34120cf1dd8SToby Isaac @*/
34220cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceDestroy(PetscDualSpace *sp)
34320cf1dd8SToby Isaac {
34420cf1dd8SToby Isaac   PetscInt       dim, f;
34520cf1dd8SToby Isaac   PetscErrorCode ierr;
34620cf1dd8SToby Isaac 
34720cf1dd8SToby Isaac   PetscFunctionBegin;
34820cf1dd8SToby Isaac   if (!*sp) PetscFunctionReturn(0);
34920cf1dd8SToby Isaac   PetscValidHeaderSpecific((*sp), PETSCDUALSPACE_CLASSID, 1);
35020cf1dd8SToby Isaac 
35120cf1dd8SToby Isaac   if (--((PetscObject)(*sp))->refct > 0) {*sp = 0; PetscFunctionReturn(0);}
35220cf1dd8SToby Isaac   ((PetscObject) (*sp))->refct = 0;
35320cf1dd8SToby Isaac 
35420cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDimension(*sp, &dim);CHKERRQ(ierr);
35520cf1dd8SToby Isaac   for (f = 0; f < dim; ++f) {
35620cf1dd8SToby Isaac     ierr = PetscQuadratureDestroy(&(*sp)->functional[f]);CHKERRQ(ierr);
35720cf1dd8SToby Isaac   }
35820cf1dd8SToby Isaac   ierr = PetscFree((*sp)->functional);CHKERRQ(ierr);
35920cf1dd8SToby Isaac   ierr = PetscQuadratureDestroy(&(*sp)->allPoints);CHKERRQ(ierr);
36020cf1dd8SToby Isaac   ierr = DMDestroy(&(*sp)->dm);CHKERRQ(ierr);
36120cf1dd8SToby Isaac 
36220cf1dd8SToby Isaac   if ((*sp)->ops->destroy) {ierr = (*(*sp)->ops->destroy)(*sp);CHKERRQ(ierr);}
36320cf1dd8SToby Isaac   ierr = PetscHeaderDestroy(sp);CHKERRQ(ierr);
36420cf1dd8SToby Isaac   PetscFunctionReturn(0);
36520cf1dd8SToby Isaac }
36620cf1dd8SToby Isaac 
36720cf1dd8SToby Isaac /*@
36820cf1dd8SToby Isaac   PetscDualSpaceCreate - Creates an empty PetscDualSpace object. The type can then be set with PetscDualSpaceSetType().
36920cf1dd8SToby Isaac 
37020cf1dd8SToby Isaac   Collective on MPI_Comm
37120cf1dd8SToby Isaac 
37220cf1dd8SToby Isaac   Input Parameter:
37320cf1dd8SToby Isaac . comm - The communicator for the PetscDualSpace object
37420cf1dd8SToby Isaac 
37520cf1dd8SToby Isaac   Output Parameter:
37620cf1dd8SToby Isaac . sp - The PetscDualSpace object
37720cf1dd8SToby Isaac 
37820cf1dd8SToby Isaac   Level: beginner
37920cf1dd8SToby Isaac 
38020cf1dd8SToby Isaac .seealso: PetscDualSpaceSetType(), PETSCDUALSPACELAGRANGE
38120cf1dd8SToby Isaac @*/
38220cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceCreate(MPI_Comm comm, PetscDualSpace *sp)
38320cf1dd8SToby Isaac {
38420cf1dd8SToby Isaac   PetscDualSpace s;
38520cf1dd8SToby Isaac   PetscErrorCode ierr;
38620cf1dd8SToby Isaac 
38720cf1dd8SToby Isaac   PetscFunctionBegin;
38820cf1dd8SToby Isaac   PetscValidPointer(sp, 2);
38920cf1dd8SToby Isaac   ierr = PetscCitationsRegister(FECitation,&FEcite);CHKERRQ(ierr);
39020cf1dd8SToby Isaac   *sp  = NULL;
39120cf1dd8SToby Isaac   ierr = PetscFEInitializePackage();CHKERRQ(ierr);
39220cf1dd8SToby Isaac 
39320cf1dd8SToby Isaac   ierr = PetscHeaderCreate(s, PETSCDUALSPACE_CLASSID, "PetscDualSpace", "Dual Space", "PetscDualSpace", comm, PetscDualSpaceDestroy, PetscDualSpaceView);CHKERRQ(ierr);
39420cf1dd8SToby Isaac 
39520cf1dd8SToby Isaac   s->order = 0;
39620cf1dd8SToby Isaac   s->Nc    = 1;
39720cf1dd8SToby Isaac   s->setupcalled = PETSC_FALSE;
39820cf1dd8SToby Isaac 
39920cf1dd8SToby Isaac   *sp = s;
40020cf1dd8SToby Isaac   PetscFunctionReturn(0);
40120cf1dd8SToby Isaac }
40220cf1dd8SToby Isaac 
40320cf1dd8SToby Isaac /*@
40420cf1dd8SToby Isaac   PetscDualSpaceDuplicate - Creates a duplicate PetscDualSpace object, however it is not setup.
40520cf1dd8SToby Isaac 
40620cf1dd8SToby Isaac   Collective on PetscDualSpace
40720cf1dd8SToby Isaac 
40820cf1dd8SToby Isaac   Input Parameter:
40920cf1dd8SToby Isaac . sp - The original PetscDualSpace
41020cf1dd8SToby Isaac 
41120cf1dd8SToby Isaac   Output Parameter:
41220cf1dd8SToby Isaac . spNew - The duplicate PetscDualSpace
41320cf1dd8SToby Isaac 
41420cf1dd8SToby Isaac   Level: beginner
41520cf1dd8SToby Isaac 
41620cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate(), PetscDualSpaceSetType()
41720cf1dd8SToby Isaac @*/
41820cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceDuplicate(PetscDualSpace sp, PetscDualSpace *spNew)
41920cf1dd8SToby Isaac {
42020cf1dd8SToby Isaac   PetscErrorCode ierr;
42120cf1dd8SToby Isaac 
42220cf1dd8SToby Isaac   PetscFunctionBegin;
42320cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
42420cf1dd8SToby Isaac   PetscValidPointer(spNew, 2);
42520cf1dd8SToby Isaac   ierr = (*sp->ops->duplicate)(sp, spNew);CHKERRQ(ierr);
42620cf1dd8SToby Isaac   PetscFunctionReturn(0);
42720cf1dd8SToby Isaac }
42820cf1dd8SToby Isaac 
42920cf1dd8SToby Isaac /*@
43020cf1dd8SToby Isaac   PetscDualSpaceGetDM - Get the DM representing the reference cell
43120cf1dd8SToby Isaac 
43220cf1dd8SToby Isaac   Not collective
43320cf1dd8SToby Isaac 
43420cf1dd8SToby Isaac   Input Parameter:
43520cf1dd8SToby Isaac . sp - The PetscDualSpace
43620cf1dd8SToby Isaac 
43720cf1dd8SToby Isaac   Output Parameter:
43820cf1dd8SToby Isaac . dm - The reference cell
43920cf1dd8SToby Isaac 
44020cf1dd8SToby Isaac   Level: intermediate
44120cf1dd8SToby Isaac 
44220cf1dd8SToby Isaac .seealso: PetscDualSpaceSetDM(), PetscDualSpaceCreate()
44320cf1dd8SToby Isaac @*/
44420cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetDM(PetscDualSpace sp, DM *dm)
44520cf1dd8SToby Isaac {
44620cf1dd8SToby Isaac   PetscFunctionBegin;
44720cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
44820cf1dd8SToby Isaac   PetscValidPointer(dm, 2);
44920cf1dd8SToby Isaac   *dm = sp->dm;
45020cf1dd8SToby Isaac   PetscFunctionReturn(0);
45120cf1dd8SToby Isaac }
45220cf1dd8SToby Isaac 
45320cf1dd8SToby Isaac /*@
45420cf1dd8SToby Isaac   PetscDualSpaceSetDM - Get the DM representing the reference cell
45520cf1dd8SToby Isaac 
45620cf1dd8SToby Isaac   Not collective
45720cf1dd8SToby Isaac 
45820cf1dd8SToby Isaac   Input Parameters:
45920cf1dd8SToby Isaac + sp - The PetscDualSpace
46020cf1dd8SToby Isaac - dm - The reference cell
46120cf1dd8SToby Isaac 
46220cf1dd8SToby Isaac   Level: intermediate
46320cf1dd8SToby Isaac 
46420cf1dd8SToby Isaac .seealso: PetscDualSpaceGetDM(), PetscDualSpaceCreate()
46520cf1dd8SToby Isaac @*/
46620cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetDM(PetscDualSpace sp, DM dm)
46720cf1dd8SToby Isaac {
46820cf1dd8SToby Isaac   PetscErrorCode ierr;
46920cf1dd8SToby Isaac 
47020cf1dd8SToby Isaac   PetscFunctionBegin;
47120cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
47220cf1dd8SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
47320cf1dd8SToby Isaac   ierr = DMDestroy(&sp->dm);CHKERRQ(ierr);
47420cf1dd8SToby Isaac   ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);
47520cf1dd8SToby Isaac   sp->dm = dm;
47620cf1dd8SToby Isaac   PetscFunctionReturn(0);
47720cf1dd8SToby Isaac }
47820cf1dd8SToby Isaac 
47920cf1dd8SToby Isaac /*@
48020cf1dd8SToby Isaac   PetscDualSpaceGetOrder - Get the order of the dual space
48120cf1dd8SToby Isaac 
48220cf1dd8SToby Isaac   Not collective
48320cf1dd8SToby Isaac 
48420cf1dd8SToby Isaac   Input Parameter:
48520cf1dd8SToby Isaac . sp - The PetscDualSpace
48620cf1dd8SToby Isaac 
48720cf1dd8SToby Isaac   Output Parameter:
48820cf1dd8SToby Isaac . order - The order
48920cf1dd8SToby Isaac 
49020cf1dd8SToby Isaac   Level: intermediate
49120cf1dd8SToby Isaac 
49220cf1dd8SToby Isaac .seealso: PetscDualSpaceSetOrder(), PetscDualSpaceCreate()
49320cf1dd8SToby Isaac @*/
49420cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetOrder(PetscDualSpace sp, PetscInt *order)
49520cf1dd8SToby Isaac {
49620cf1dd8SToby Isaac   PetscFunctionBegin;
49720cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
49820cf1dd8SToby Isaac   PetscValidPointer(order, 2);
49920cf1dd8SToby Isaac   *order = sp->order;
50020cf1dd8SToby Isaac   PetscFunctionReturn(0);
50120cf1dd8SToby Isaac }
50220cf1dd8SToby Isaac 
50320cf1dd8SToby Isaac /*@
50420cf1dd8SToby Isaac   PetscDualSpaceSetOrder - Set the order of the dual space
50520cf1dd8SToby Isaac 
50620cf1dd8SToby Isaac   Not collective
50720cf1dd8SToby Isaac 
50820cf1dd8SToby Isaac   Input Parameters:
50920cf1dd8SToby Isaac + sp - The PetscDualSpace
51020cf1dd8SToby Isaac - order - The order
51120cf1dd8SToby Isaac 
51220cf1dd8SToby Isaac   Level: intermediate
51320cf1dd8SToby Isaac 
51420cf1dd8SToby Isaac .seealso: PetscDualSpaceGetOrder(), PetscDualSpaceCreate()
51520cf1dd8SToby Isaac @*/
51620cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetOrder(PetscDualSpace sp, PetscInt order)
51720cf1dd8SToby Isaac {
51820cf1dd8SToby Isaac   PetscFunctionBegin;
51920cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
52020cf1dd8SToby Isaac   sp->order = order;
52120cf1dd8SToby Isaac   PetscFunctionReturn(0);
52220cf1dd8SToby Isaac }
52320cf1dd8SToby Isaac 
52420cf1dd8SToby Isaac /*@
52520cf1dd8SToby Isaac   PetscDualSpaceGetNumComponents - Return the number of components for this space
52620cf1dd8SToby Isaac 
52720cf1dd8SToby Isaac   Input Parameter:
52820cf1dd8SToby Isaac . sp - The PetscDualSpace
52920cf1dd8SToby Isaac 
53020cf1dd8SToby Isaac   Output Parameter:
53120cf1dd8SToby Isaac . Nc - The number of components
53220cf1dd8SToby Isaac 
53320cf1dd8SToby Isaac   Note: A vector space, for example, will have d components, where d is the spatial dimension
53420cf1dd8SToby Isaac 
53520cf1dd8SToby Isaac   Level: intermediate
53620cf1dd8SToby Isaac 
53720cf1dd8SToby Isaac .seealso: PetscDualSpaceSetNumComponents(), PetscDualSpaceGetDimension(), PetscDualSpaceCreate(), PetscDualSpace
53820cf1dd8SToby Isaac @*/
53920cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetNumComponents(PetscDualSpace sp, PetscInt *Nc)
54020cf1dd8SToby Isaac {
54120cf1dd8SToby Isaac   PetscFunctionBegin;
54220cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
54320cf1dd8SToby Isaac   PetscValidPointer(Nc, 2);
54420cf1dd8SToby Isaac   *Nc = sp->Nc;
54520cf1dd8SToby Isaac   PetscFunctionReturn(0);
54620cf1dd8SToby Isaac }
54720cf1dd8SToby Isaac 
54820cf1dd8SToby Isaac /*@
54920cf1dd8SToby Isaac   PetscDualSpaceSetNumComponents - Set the number of components for this space
55020cf1dd8SToby Isaac 
55120cf1dd8SToby Isaac   Input Parameters:
55220cf1dd8SToby Isaac + sp - The PetscDualSpace
55320cf1dd8SToby Isaac - order - The number of components
55420cf1dd8SToby Isaac 
55520cf1dd8SToby Isaac   Level: intermediate
55620cf1dd8SToby Isaac 
55720cf1dd8SToby Isaac .seealso: PetscDualSpaceGetNumComponents(), PetscDualSpaceCreate(), PetscDualSpace
55820cf1dd8SToby Isaac @*/
55920cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetNumComponents(PetscDualSpace sp, PetscInt Nc)
56020cf1dd8SToby Isaac {
56120cf1dd8SToby Isaac   PetscFunctionBegin;
56220cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
56320cf1dd8SToby Isaac   sp->Nc = Nc;
56420cf1dd8SToby Isaac   PetscFunctionReturn(0);
56520cf1dd8SToby Isaac }
56620cf1dd8SToby Isaac 
56720cf1dd8SToby Isaac /*@
56820cf1dd8SToby Isaac   PetscDualSpaceGetFunctional - Get the i-th basis functional in the dual space
56920cf1dd8SToby Isaac 
57020cf1dd8SToby Isaac   Not collective
57120cf1dd8SToby Isaac 
57220cf1dd8SToby Isaac   Input Parameters:
57320cf1dd8SToby Isaac + sp - The PetscDualSpace
57420cf1dd8SToby Isaac - i  - The basis number
57520cf1dd8SToby Isaac 
57620cf1dd8SToby Isaac   Output Parameter:
57720cf1dd8SToby Isaac . functional - The basis functional
57820cf1dd8SToby Isaac 
57920cf1dd8SToby Isaac   Level: intermediate
58020cf1dd8SToby Isaac 
58120cf1dd8SToby Isaac .seealso: PetscDualSpaceGetDimension(), PetscDualSpaceCreate()
58220cf1dd8SToby Isaac @*/
58320cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetFunctional(PetscDualSpace sp, PetscInt i, PetscQuadrature *functional)
58420cf1dd8SToby Isaac {
58520cf1dd8SToby Isaac   PetscInt       dim;
58620cf1dd8SToby Isaac   PetscErrorCode ierr;
58720cf1dd8SToby Isaac 
58820cf1dd8SToby Isaac   PetscFunctionBegin;
58920cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
59020cf1dd8SToby Isaac   PetscValidPointer(functional, 3);
59120cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDimension(sp, &dim);CHKERRQ(ierr);
59220cf1dd8SToby Isaac   if ((i < 0) || (i >= dim)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Functional index %d must be in [0, %d)", i, dim);
59320cf1dd8SToby Isaac   *functional = sp->functional[i];
59420cf1dd8SToby Isaac   PetscFunctionReturn(0);
59520cf1dd8SToby Isaac }
59620cf1dd8SToby Isaac 
59720cf1dd8SToby Isaac /*@
59820cf1dd8SToby Isaac   PetscDualSpaceGetDimension - Get the dimension of the dual space, i.e. the number of basis functionals
59920cf1dd8SToby Isaac 
60020cf1dd8SToby Isaac   Not collective
60120cf1dd8SToby Isaac 
60220cf1dd8SToby Isaac   Input Parameter:
60320cf1dd8SToby Isaac . sp - The PetscDualSpace
60420cf1dd8SToby Isaac 
60520cf1dd8SToby Isaac   Output Parameter:
60620cf1dd8SToby Isaac . dim - The dimension
60720cf1dd8SToby Isaac 
60820cf1dd8SToby Isaac   Level: intermediate
60920cf1dd8SToby Isaac 
61020cf1dd8SToby Isaac .seealso: PetscDualSpaceGetFunctional(), PetscDualSpaceCreate()
61120cf1dd8SToby Isaac @*/
61220cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetDimension(PetscDualSpace sp, PetscInt *dim)
61320cf1dd8SToby Isaac {
61420cf1dd8SToby Isaac   PetscErrorCode ierr;
61520cf1dd8SToby Isaac 
61620cf1dd8SToby Isaac   PetscFunctionBegin;
61720cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
61820cf1dd8SToby Isaac   PetscValidPointer(dim, 2);
61920cf1dd8SToby Isaac   *dim = 0;
62020cf1dd8SToby Isaac   if (sp->ops->getdimension) {ierr = (*sp->ops->getdimension)(sp, dim);CHKERRQ(ierr);}
62120cf1dd8SToby Isaac   PetscFunctionReturn(0);
62220cf1dd8SToby Isaac }
62320cf1dd8SToby Isaac 
62420cf1dd8SToby Isaac /*@C
62520cf1dd8SToby Isaac   PetscDualSpaceGetNumDof - Get the number of degrees of freedom for each spatial (topological) dimension
62620cf1dd8SToby Isaac 
62720cf1dd8SToby Isaac   Not collective
62820cf1dd8SToby Isaac 
62920cf1dd8SToby Isaac   Input Parameter:
63020cf1dd8SToby Isaac . sp - The PetscDualSpace
63120cf1dd8SToby Isaac 
63220cf1dd8SToby Isaac   Output Parameter:
63320cf1dd8SToby Isaac . numDof - An array of length dim+1 which holds the number of dofs for each dimension
63420cf1dd8SToby Isaac 
63520cf1dd8SToby Isaac   Level: intermediate
63620cf1dd8SToby Isaac 
63720cf1dd8SToby Isaac .seealso: PetscDualSpaceGetFunctional(), PetscDualSpaceCreate()
63820cf1dd8SToby Isaac @*/
63920cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetNumDof(PetscDualSpace sp, const PetscInt **numDof)
64020cf1dd8SToby Isaac {
64120cf1dd8SToby Isaac   PetscErrorCode ierr;
64220cf1dd8SToby Isaac 
64320cf1dd8SToby Isaac   PetscFunctionBegin;
64420cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
64520cf1dd8SToby Isaac   PetscValidPointer(numDof, 2);
64620cf1dd8SToby Isaac   ierr = (*sp->ops->getnumdof)(sp, numDof);CHKERRQ(ierr);
64720cf1dd8SToby Isaac   if (!*numDof) SETERRQ(PetscObjectComm((PetscObject) sp), PETSC_ERR_LIB, "Empty numDof[] returned from dual space implementation");
64820cf1dd8SToby Isaac   PetscFunctionReturn(0);
64920cf1dd8SToby Isaac }
65020cf1dd8SToby Isaac 
65120cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceCreateSection(PetscDualSpace sp, PetscSection *section)
65220cf1dd8SToby Isaac {
65320cf1dd8SToby Isaac   DM             dm;
65420cf1dd8SToby Isaac   PetscInt       pStart, pEnd, depth, h, offset;
65520cf1dd8SToby Isaac   const PetscInt *numDof;
65620cf1dd8SToby Isaac   PetscErrorCode ierr;
65720cf1dd8SToby Isaac 
65820cf1dd8SToby Isaac   PetscFunctionBegin;
65920cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDM(sp,&dm);CHKERRQ(ierr);
66020cf1dd8SToby Isaac   ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr);
66120cf1dd8SToby Isaac   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)sp),section);CHKERRQ(ierr);
66220cf1dd8SToby Isaac   ierr = PetscSectionSetChart(*section,pStart,pEnd);CHKERRQ(ierr);
66320cf1dd8SToby Isaac   ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr);
66420cf1dd8SToby Isaac   ierr = PetscDualSpaceGetNumDof(sp,&numDof);CHKERRQ(ierr);
66520cf1dd8SToby Isaac   for (h = 0; h <= depth; h++) {
66620cf1dd8SToby Isaac     PetscInt hStart, hEnd, p, dof;
66720cf1dd8SToby Isaac 
66820cf1dd8SToby Isaac     ierr = DMPlexGetHeightStratum(dm,h,&hStart,&hEnd);CHKERRQ(ierr);
66920cf1dd8SToby Isaac     dof = numDof[depth - h];
67020cf1dd8SToby Isaac     for (p = hStart; p < hEnd; p++) {
67120cf1dd8SToby Isaac       ierr = PetscSectionSetDof(*section,p,dof);CHKERRQ(ierr);
67220cf1dd8SToby Isaac     }
67320cf1dd8SToby Isaac   }
67420cf1dd8SToby Isaac   ierr = PetscSectionSetUp(*section);CHKERRQ(ierr);
67520cf1dd8SToby Isaac   for (h = 0, offset = 0; h <= depth; h++) {
67620cf1dd8SToby Isaac     PetscInt hStart, hEnd, p, dof;
67720cf1dd8SToby Isaac 
67820cf1dd8SToby Isaac     ierr = DMPlexGetHeightStratum(dm,h,&hStart,&hEnd);CHKERRQ(ierr);
67920cf1dd8SToby Isaac     dof = numDof[depth - h];
68020cf1dd8SToby Isaac     for (p = hStart; p < hEnd; p++) {
68120cf1dd8SToby Isaac       ierr = PetscSectionGetDof(*section,p,&dof);CHKERRQ(ierr);
68220cf1dd8SToby Isaac       ierr = PetscSectionSetOffset(*section,p,offset);CHKERRQ(ierr);
68320cf1dd8SToby Isaac       offset += dof;
68420cf1dd8SToby Isaac     }
68520cf1dd8SToby Isaac   }
68620cf1dd8SToby Isaac   PetscFunctionReturn(0);
68720cf1dd8SToby Isaac }
68820cf1dd8SToby Isaac 
68920cf1dd8SToby Isaac /*@
69020cf1dd8SToby Isaac   PetscDualSpaceCreateReferenceCell - Create a DMPLEX with the appropriate FEM reference cell
69120cf1dd8SToby Isaac 
69220cf1dd8SToby Isaac   Collective on PetscDualSpace
69320cf1dd8SToby Isaac 
69420cf1dd8SToby Isaac   Input Parameters:
69520cf1dd8SToby Isaac + sp      - The PetscDualSpace
69620cf1dd8SToby Isaac . dim     - The spatial dimension
69720cf1dd8SToby Isaac - simplex - Flag for simplex, otherwise use a tensor-product cell
69820cf1dd8SToby Isaac 
69920cf1dd8SToby Isaac   Output Parameter:
70020cf1dd8SToby Isaac . refdm - The reference cell
70120cf1dd8SToby Isaac 
70220cf1dd8SToby Isaac   Level: advanced
70320cf1dd8SToby Isaac 
70420cf1dd8SToby Isaac .keywords: PetscDualSpace, reference cell
70520cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate(), DMPLEX
70620cf1dd8SToby Isaac @*/
70720cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceCreateReferenceCell(PetscDualSpace sp, PetscInt dim, PetscBool simplex, DM *refdm)
70820cf1dd8SToby Isaac {
70920cf1dd8SToby Isaac   PetscErrorCode ierr;
71020cf1dd8SToby Isaac 
71120cf1dd8SToby Isaac   PetscFunctionBeginUser;
71220cf1dd8SToby Isaac   ierr = DMPlexCreateReferenceCell(PetscObjectComm((PetscObject) sp), dim, simplex, refdm);CHKERRQ(ierr);
71320cf1dd8SToby Isaac   PetscFunctionReturn(0);
71420cf1dd8SToby Isaac }
71520cf1dd8SToby Isaac 
71620cf1dd8SToby Isaac /*@C
71720cf1dd8SToby Isaac   PetscDualSpaceApply - Apply a functional from the dual space basis to an input function
71820cf1dd8SToby Isaac 
71920cf1dd8SToby Isaac   Input Parameters:
72020cf1dd8SToby Isaac + sp      - The PetscDualSpace object
72120cf1dd8SToby Isaac . f       - The basis functional index
72220cf1dd8SToby Isaac . time    - The time
72320cf1dd8SToby Isaac . 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)
72420cf1dd8SToby Isaac . numComp - The number of components for the function
72520cf1dd8SToby Isaac . func    - The input function
72620cf1dd8SToby Isaac - ctx     - A context for the function
72720cf1dd8SToby Isaac 
72820cf1dd8SToby Isaac   Output Parameter:
72920cf1dd8SToby Isaac . value   - numComp output values
73020cf1dd8SToby Isaac 
73120cf1dd8SToby Isaac   Note: The calling sequence for the callback func is given by:
73220cf1dd8SToby Isaac 
73320cf1dd8SToby Isaac $ func(PetscInt dim, PetscReal time, const PetscReal x[],
73420cf1dd8SToby Isaac $      PetscInt numComponents, PetscScalar values[], void *ctx)
73520cf1dd8SToby Isaac 
73620cf1dd8SToby Isaac   Level: developer
73720cf1dd8SToby Isaac 
73820cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate()
73920cf1dd8SToby Isaac @*/
74020cf1dd8SToby Isaac 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)
74120cf1dd8SToby Isaac {
74220cf1dd8SToby Isaac   PetscErrorCode ierr;
74320cf1dd8SToby Isaac 
74420cf1dd8SToby Isaac   PetscFunctionBegin;
74520cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
74620cf1dd8SToby Isaac   PetscValidPointer(cgeom, 4);
74720cf1dd8SToby Isaac   PetscValidPointer(value, 8);
74820cf1dd8SToby Isaac   ierr = (*sp->ops->apply)(sp, f, time, cgeom, numComp, func, ctx, value);CHKERRQ(ierr);
74920cf1dd8SToby Isaac   PetscFunctionReturn(0);
75020cf1dd8SToby Isaac }
75120cf1dd8SToby Isaac 
75220cf1dd8SToby Isaac /*@C
75320cf1dd8SToby Isaac   PetscDualSpaceApplyAll - Apply all functionals from the dual space basis to the result of an evaluation at the points returned by PetscDualSpaceGetAllPoints()
75420cf1dd8SToby Isaac 
75520cf1dd8SToby Isaac   Input Parameters:
75620cf1dd8SToby Isaac + sp        - The PetscDualSpace object
75720cf1dd8SToby Isaac - pointEval - Evaluation at the points returned by PetscDualSpaceGetAllPoints()
75820cf1dd8SToby Isaac 
75920cf1dd8SToby Isaac   Output Parameter:
76020cf1dd8SToby Isaac . spValue   - The values of all dual space functionals
76120cf1dd8SToby Isaac 
76220cf1dd8SToby Isaac   Level: developer
76320cf1dd8SToby Isaac 
76420cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate()
76520cf1dd8SToby Isaac @*/
76620cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceApplyAll(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue)
76720cf1dd8SToby Isaac {
76820cf1dd8SToby Isaac   PetscErrorCode ierr;
76920cf1dd8SToby Isaac 
77020cf1dd8SToby Isaac   PetscFunctionBegin;
77120cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
77220cf1dd8SToby Isaac   ierr = (*sp->ops->applyall)(sp, pointEval, spValue);CHKERRQ(ierr);
77320cf1dd8SToby Isaac   PetscFunctionReturn(0);
77420cf1dd8SToby Isaac }
77520cf1dd8SToby Isaac 
77620cf1dd8SToby Isaac /*@C
77720cf1dd8SToby Isaac   PetscDualSpaceApplyDefault - Apply a functional from the dual space basis to an input function by assuming a point evaluation functional.
77820cf1dd8SToby Isaac 
77920cf1dd8SToby Isaac   Input Parameters:
78020cf1dd8SToby Isaac + sp    - The PetscDualSpace object
78120cf1dd8SToby Isaac . f     - The basis functional index
78220cf1dd8SToby Isaac . time  - The time
78320cf1dd8SToby Isaac . cgeom - A context with geometric information for this cell, we use v0 (the initial vertex) and J (the Jacobian)
78420cf1dd8SToby Isaac . Nc    - The number of components for the function
78520cf1dd8SToby Isaac . func  - The input function
78620cf1dd8SToby Isaac - ctx   - A context for the function
78720cf1dd8SToby Isaac 
78820cf1dd8SToby Isaac   Output Parameter:
78920cf1dd8SToby Isaac . value   - The output value
79020cf1dd8SToby Isaac 
79120cf1dd8SToby Isaac   Note: The calling sequence for the callback func is given by:
79220cf1dd8SToby Isaac 
79320cf1dd8SToby Isaac $ func(PetscInt dim, PetscReal time, const PetscReal x[],
79420cf1dd8SToby Isaac $      PetscInt numComponents, PetscScalar values[], void *ctx)
79520cf1dd8SToby Isaac 
79620cf1dd8SToby Isaac and the idea is to evaluate the functional as an integral
79720cf1dd8SToby Isaac 
79820cf1dd8SToby Isaac $ n(f) = int dx n(x) . f(x)
79920cf1dd8SToby Isaac 
80020cf1dd8SToby Isaac where both n and f have Nc components.
80120cf1dd8SToby Isaac 
80220cf1dd8SToby Isaac   Level: developer
80320cf1dd8SToby Isaac 
80420cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate()
80520cf1dd8SToby Isaac @*/
80620cf1dd8SToby Isaac 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)
80720cf1dd8SToby Isaac {
80820cf1dd8SToby Isaac   DM               dm;
80920cf1dd8SToby Isaac   PetscQuadrature  n;
81020cf1dd8SToby Isaac   const PetscReal *points, *weights;
81120cf1dd8SToby Isaac   PetscReal        x[3];
81220cf1dd8SToby Isaac   PetscScalar     *val;
81320cf1dd8SToby Isaac   PetscInt         dim, dE, qNc, c, Nq, q;
81420cf1dd8SToby Isaac   PetscBool        isAffine;
81520cf1dd8SToby Isaac   PetscErrorCode   ierr;
81620cf1dd8SToby Isaac 
81720cf1dd8SToby Isaac   PetscFunctionBegin;
81820cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
81920cf1dd8SToby Isaac   PetscValidPointer(value, 5);
82020cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr);
82120cf1dd8SToby Isaac   ierr = PetscDualSpaceGetFunctional(sp, f, &n);CHKERRQ(ierr);
82220cf1dd8SToby Isaac   ierr = PetscQuadratureGetData(n, &dim, &qNc, &Nq, &points, &weights);CHKERRQ(ierr);
82320cf1dd8SToby Isaac   if (dim != cgeom->dim) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature spatial dimension %D != cell geometry dimension %D", dim, cgeom->dim);
82420cf1dd8SToby Isaac   if (qNc != Nc) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature components %D != function components %D", qNc, Nc);
82520cf1dd8SToby Isaac   ierr = DMGetWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr);
82620cf1dd8SToby Isaac   *value = 0.0;
82720cf1dd8SToby Isaac   isAffine = cgeom->isAffine;
82820cf1dd8SToby Isaac   dE = cgeom->dimEmbed;
82920cf1dd8SToby Isaac   for (q = 0; q < Nq; ++q) {
83020cf1dd8SToby Isaac     if (isAffine) {
83120cf1dd8SToby Isaac       CoordinatesRefToReal(dE, cgeom->dim, cgeom->xi, cgeom->v, cgeom->J, &points[q*dim], x);
83220cf1dd8SToby Isaac       ierr = (*func)(dE, time, x, Nc, val, ctx);CHKERRQ(ierr);
83320cf1dd8SToby Isaac     } else {
83420cf1dd8SToby Isaac       ierr = (*func)(dE, time, &cgeom->v[dE*q], Nc, val, ctx);CHKERRQ(ierr);
83520cf1dd8SToby Isaac     }
83620cf1dd8SToby Isaac     for (c = 0; c < Nc; ++c) {
83720cf1dd8SToby Isaac       *value += val[c]*weights[q*Nc+c];
83820cf1dd8SToby Isaac     }
83920cf1dd8SToby Isaac   }
84020cf1dd8SToby Isaac   ierr = DMRestoreWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr);
84120cf1dd8SToby Isaac   PetscFunctionReturn(0);
84220cf1dd8SToby Isaac }
84320cf1dd8SToby Isaac 
84420cf1dd8SToby Isaac /*@C
84520cf1dd8SToby Isaac   PetscDualSpaceApplyAllDefault - Apply all functionals from the dual space basis to the result of an evaluation at the points returned by PetscDualSpaceGetAllPoints()
84620cf1dd8SToby Isaac 
84720cf1dd8SToby Isaac   Input Parameters:
84820cf1dd8SToby Isaac + sp        - The PetscDualSpace object
84920cf1dd8SToby Isaac - pointEval - Evaluation at the points returned by PetscDualSpaceGetAllPoints()
85020cf1dd8SToby Isaac 
85120cf1dd8SToby Isaac   Output Parameter:
85220cf1dd8SToby Isaac . spValue   - The values of all dual space functionals
85320cf1dd8SToby Isaac 
85420cf1dd8SToby Isaac   Level: developer
85520cf1dd8SToby Isaac 
85620cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate()
85720cf1dd8SToby Isaac @*/
85820cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceApplyAllDefault(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue)
85920cf1dd8SToby Isaac {
86020cf1dd8SToby Isaac   PetscQuadrature  n;
86120cf1dd8SToby Isaac   const PetscReal *points, *weights;
86220cf1dd8SToby Isaac   PetscInt         qNc, c, Nq, q, f, spdim, Nc;
86320cf1dd8SToby Isaac   PetscInt         offset;
86420cf1dd8SToby Isaac   PetscErrorCode   ierr;
86520cf1dd8SToby Isaac 
86620cf1dd8SToby Isaac   PetscFunctionBegin;
86720cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
86820cf1dd8SToby Isaac   PetscValidScalarPointer(pointEval, 2);
86920cf1dd8SToby Isaac   PetscValidScalarPointer(spValue, 5);
87020cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDimension(sp, &spdim);CHKERRQ(ierr);
87120cf1dd8SToby Isaac   ierr = PetscDualSpaceGetNumComponents(sp, &Nc);CHKERRQ(ierr);
87220cf1dd8SToby Isaac   for (f = 0, offset = 0; f < spdim; f++) {
87320cf1dd8SToby Isaac     ierr = PetscDualSpaceGetFunctional(sp, f, &n);CHKERRQ(ierr);
87420cf1dd8SToby Isaac     ierr = PetscQuadratureGetData(n, NULL, &qNc, &Nq, &points, &weights);CHKERRQ(ierr);
87520cf1dd8SToby Isaac     if (qNc != Nc) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature components %D != function components %D", qNc, Nc);
87620cf1dd8SToby Isaac     spValue[f] = 0.0;
87720cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
87820cf1dd8SToby Isaac       for (c = 0; c < Nc; ++c) {
87920cf1dd8SToby Isaac         spValue[f] += pointEval[offset++]*weights[q*Nc+c];
88020cf1dd8SToby Isaac       }
88120cf1dd8SToby Isaac     }
88220cf1dd8SToby Isaac   }
88320cf1dd8SToby Isaac   PetscFunctionReturn(0);
88420cf1dd8SToby Isaac }
88520cf1dd8SToby Isaac 
88620cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetAllPoints(PetscDualSpace sp, PetscQuadrature *allPoints)
88720cf1dd8SToby Isaac {
88820cf1dd8SToby Isaac   PetscErrorCode ierr;
88920cf1dd8SToby Isaac 
89020cf1dd8SToby Isaac   PetscFunctionBegin;
89120cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
89220cf1dd8SToby Isaac   PetscValidPointer(allPoints,2);
89320cf1dd8SToby Isaac   if (!sp->allPoints && sp->ops->createallpoints) {
89420cf1dd8SToby Isaac     ierr = (*sp->ops->createallpoints)(sp,&sp->allPoints);CHKERRQ(ierr);
89520cf1dd8SToby Isaac   }
89620cf1dd8SToby Isaac   *allPoints = sp->allPoints;
89720cf1dd8SToby Isaac   PetscFunctionReturn(0);
89820cf1dd8SToby Isaac }
89920cf1dd8SToby Isaac 
90020cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceCreateAllPointsDefault(PetscDualSpace sp, PetscQuadrature *allPoints)
90120cf1dd8SToby Isaac {
90220cf1dd8SToby Isaac   PetscInt        spdim;
90320cf1dd8SToby Isaac   PetscInt        numPoints, offset;
90420cf1dd8SToby Isaac   PetscReal       *points;
90520cf1dd8SToby Isaac   PetscInt        f, dim;
90620cf1dd8SToby Isaac   PetscQuadrature q;
90720cf1dd8SToby Isaac   PetscErrorCode  ierr;
90820cf1dd8SToby Isaac 
90920cf1dd8SToby Isaac   PetscFunctionBegin;
91020cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDimension(sp,&spdim);CHKERRQ(ierr);
91120cf1dd8SToby Isaac   if (!spdim) {
91220cf1dd8SToby Isaac     ierr = PetscQuadratureCreate(PETSC_COMM_SELF,allPoints);CHKERRQ(ierr);
91320cf1dd8SToby Isaac     ierr = PetscQuadratureSetData(*allPoints,0,0,0,NULL,NULL);CHKERRQ(ierr);
91420cf1dd8SToby Isaac   }
91520cf1dd8SToby Isaac   ierr = PetscDualSpaceGetFunctional(sp,0,&q);CHKERRQ(ierr);
91620cf1dd8SToby Isaac   ierr = PetscQuadratureGetData(q,&dim,NULL,&numPoints,NULL,NULL);CHKERRQ(ierr);
91720cf1dd8SToby Isaac   for (f = 1; f < spdim; f++) {
91820cf1dd8SToby Isaac     PetscInt Np;
91920cf1dd8SToby Isaac 
92020cf1dd8SToby Isaac     ierr = PetscDualSpaceGetFunctional(sp,f,&q);CHKERRQ(ierr);
92120cf1dd8SToby Isaac     ierr = PetscQuadratureGetData(q,NULL,NULL,&Np,NULL,NULL);CHKERRQ(ierr);
92220cf1dd8SToby Isaac     numPoints += Np;
92320cf1dd8SToby Isaac   }
92420cf1dd8SToby Isaac   ierr = PetscMalloc1(dim*numPoints,&points);CHKERRQ(ierr);
92520cf1dd8SToby Isaac   for (f = 0, offset = 0; f < spdim; f++) {
92620cf1dd8SToby Isaac     const PetscReal *p;
92720cf1dd8SToby Isaac     PetscInt        Np, i;
92820cf1dd8SToby Isaac 
92920cf1dd8SToby Isaac     ierr = PetscDualSpaceGetFunctional(sp,f,&q);CHKERRQ(ierr);
93020cf1dd8SToby Isaac     ierr = PetscQuadratureGetData(q,NULL,NULL,&Np,&p,NULL);CHKERRQ(ierr);
93120cf1dd8SToby Isaac     for (i = 0; i < Np * dim; i++) {
93220cf1dd8SToby Isaac       points[offset + i] = p[i];
93320cf1dd8SToby Isaac     }
93420cf1dd8SToby Isaac     offset += Np * dim;
93520cf1dd8SToby Isaac   }
93620cf1dd8SToby Isaac   ierr = PetscQuadratureCreate(PETSC_COMM_SELF,allPoints);CHKERRQ(ierr);
93720cf1dd8SToby Isaac   ierr = PetscQuadratureSetData(*allPoints,dim,0,numPoints,points,NULL);CHKERRQ(ierr);
93820cf1dd8SToby Isaac   PetscFunctionReturn(0);
93920cf1dd8SToby Isaac }
94020cf1dd8SToby Isaac 
94120cf1dd8SToby Isaac /*@C
94220cf1dd8SToby Isaac   PetscDualSpaceApplyFVM - Apply a functional from the dual space basis to an input function by assuming a point evaluation functional at the cell centroid.
94320cf1dd8SToby Isaac 
94420cf1dd8SToby Isaac   Input Parameters:
94520cf1dd8SToby Isaac + sp    - The PetscDualSpace object
94620cf1dd8SToby Isaac . f     - The basis functional index
94720cf1dd8SToby Isaac . time  - The time
94820cf1dd8SToby Isaac . cgeom - A context with geometric information for this cell, we currently just use the centroid
94920cf1dd8SToby Isaac . Nc    - The number of components for the function
95020cf1dd8SToby Isaac . func  - The input function
95120cf1dd8SToby Isaac - ctx   - A context for the function
95220cf1dd8SToby Isaac 
95320cf1dd8SToby Isaac   Output Parameter:
95420cf1dd8SToby Isaac . value - The output value (scalar)
95520cf1dd8SToby Isaac 
95620cf1dd8SToby Isaac   Note: The calling sequence for the callback func is given by:
95720cf1dd8SToby Isaac 
95820cf1dd8SToby Isaac $ func(PetscInt dim, PetscReal time, const PetscReal x[],
95920cf1dd8SToby Isaac $      PetscInt numComponents, PetscScalar values[], void *ctx)
96020cf1dd8SToby Isaac 
96120cf1dd8SToby Isaac and the idea is to evaluate the functional as an integral
96220cf1dd8SToby Isaac 
96320cf1dd8SToby Isaac $ n(f) = int dx n(x) . f(x)
96420cf1dd8SToby Isaac 
96520cf1dd8SToby Isaac where both n and f have Nc components.
96620cf1dd8SToby Isaac 
96720cf1dd8SToby Isaac   Level: developer
96820cf1dd8SToby Isaac 
96920cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate()
97020cf1dd8SToby Isaac @*/
97120cf1dd8SToby Isaac 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)
97220cf1dd8SToby Isaac {
97320cf1dd8SToby Isaac   DM               dm;
97420cf1dd8SToby Isaac   PetscQuadrature  n;
97520cf1dd8SToby Isaac   const PetscReal *points, *weights;
97620cf1dd8SToby Isaac   PetscScalar     *val;
97720cf1dd8SToby Isaac   PetscInt         dimEmbed, qNc, c, Nq, q;
97820cf1dd8SToby Isaac   PetscErrorCode   ierr;
97920cf1dd8SToby Isaac 
98020cf1dd8SToby Isaac   PetscFunctionBegin;
98120cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
98220cf1dd8SToby Isaac   PetscValidPointer(value, 5);
98320cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr);
98420cf1dd8SToby Isaac   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
98520cf1dd8SToby Isaac   ierr = PetscDualSpaceGetFunctional(sp, f, &n);CHKERRQ(ierr);
98620cf1dd8SToby Isaac   ierr = PetscQuadratureGetData(n, NULL, &qNc, &Nq, &points, &weights);CHKERRQ(ierr);
98720cf1dd8SToby Isaac   if (qNc != Nc) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature components %D != function components %D", qNc, Nc);
98820cf1dd8SToby Isaac   ierr = DMGetWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr);
98920cf1dd8SToby Isaac   *value = 0.;
99020cf1dd8SToby Isaac   for (q = 0; q < Nq; ++q) {
99120cf1dd8SToby Isaac     ierr = (*func)(dimEmbed, time, cgeom->centroid, Nc, val, ctx);CHKERRQ(ierr);
99220cf1dd8SToby Isaac     for (c = 0; c < Nc; ++c) {
99320cf1dd8SToby Isaac       *value += val[c]*weights[q*Nc+c];
99420cf1dd8SToby Isaac     }
99520cf1dd8SToby Isaac   }
99620cf1dd8SToby Isaac   ierr = DMRestoreWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr);
99720cf1dd8SToby Isaac   PetscFunctionReturn(0);
99820cf1dd8SToby Isaac }
99920cf1dd8SToby Isaac 
100020cf1dd8SToby Isaac /*@
100120cf1dd8SToby Isaac   PetscDualSpaceGetHeightSubspace - Get the subset of the dual space basis that is supported on a mesh point of a
100220cf1dd8SToby Isaac   given height.  This assumes that the reference cell is symmetric over points of this height.
100320cf1dd8SToby Isaac 
100420cf1dd8SToby Isaac   If the dual space is not defined on mesh points of the given height (e.g. if the space is discontinuous and
100520cf1dd8SToby Isaac   pointwise values are not defined on the element boundaries), or if the implementation of PetscDualSpace does not
100620cf1dd8SToby Isaac   support extracting subspaces, then NULL is returned.
100720cf1dd8SToby Isaac 
100820cf1dd8SToby Isaac   This does not increment the reference count on the returned dual space, and the user should not destroy it.
100920cf1dd8SToby Isaac 
101020cf1dd8SToby Isaac   Not collective
101120cf1dd8SToby Isaac 
101220cf1dd8SToby Isaac   Input Parameters:
101320cf1dd8SToby Isaac + sp - the PetscDualSpace object
101420cf1dd8SToby Isaac - height - the height of the mesh point for which the subspace is desired
101520cf1dd8SToby Isaac 
101620cf1dd8SToby Isaac   Output Parameter:
101720cf1dd8SToby Isaac . subsp - the subspace.  Note that the functionals in the subspace are with respect to the intrinsic geometry of the
101820cf1dd8SToby Isaac   point, which will be of lesser dimension if height > 0.
101920cf1dd8SToby Isaac 
102020cf1dd8SToby Isaac   Level: advanced
102120cf1dd8SToby Isaac 
102220cf1dd8SToby Isaac .seealso: PetscSpaceGetHeightSubspace(), PetscDualSpace
102320cf1dd8SToby Isaac @*/
102420cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetHeightSubspace(PetscDualSpace sp, PetscInt height, PetscDualSpace *subsp)
102520cf1dd8SToby Isaac {
102620cf1dd8SToby Isaac   PetscErrorCode ierr;
102720cf1dd8SToby Isaac 
102820cf1dd8SToby Isaac   PetscFunctionBegin;
102920cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
103020cf1dd8SToby Isaac   PetscValidPointer(subsp, 3);
103120cf1dd8SToby Isaac   *subsp = NULL;
1032aa545514SMatthew G. Knepley   if (sp->ops->getheightsubspace) {ierr = (*sp->ops->getheightsubspace)(sp, height, subsp);CHKERRQ(ierr);}
103320cf1dd8SToby Isaac   PetscFunctionReturn(0);
103420cf1dd8SToby Isaac }
103520cf1dd8SToby Isaac 
103620cf1dd8SToby Isaac /*@
103720cf1dd8SToby Isaac   PetscDualSpaceGetPointSubspace - Get the subset of the dual space basis that is supported on a particular mesh point.
103820cf1dd8SToby Isaac 
103920cf1dd8SToby Isaac   If the dual space is not defined on the mesh point (e.g. if the space is discontinuous and pointwise values are not
104020cf1dd8SToby Isaac   defined on the element boundaries), or if the implementation of PetscDualSpace does not support extracting
104120cf1dd8SToby Isaac   subspaces, then NULL is returned.
104220cf1dd8SToby Isaac 
104320cf1dd8SToby Isaac   This does not increment the reference count on the returned dual space, and the user should not destroy it.
104420cf1dd8SToby Isaac 
104520cf1dd8SToby Isaac   Not collective
104620cf1dd8SToby Isaac 
104720cf1dd8SToby Isaac   Input Parameters:
104820cf1dd8SToby Isaac + sp - the PetscDualSpace object
104920cf1dd8SToby Isaac - point - the point (in the dual space's DM) for which the subspace is desired
105020cf1dd8SToby Isaac 
105120cf1dd8SToby Isaac   Output Parameters:
105220cf1dd8SToby Isaac   bdsp - the subspace.  Note that the functionals in the subspace are with respect to the intrinsic geometry of the
105320cf1dd8SToby Isaac   point, which will be of lesser dimension if height > 0.
105420cf1dd8SToby Isaac 
105520cf1dd8SToby Isaac   Level: advanced
105620cf1dd8SToby Isaac 
105720cf1dd8SToby Isaac .seealso: PetscDualSpace
105820cf1dd8SToby Isaac @*/
105920cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetPointSubspace(PetscDualSpace sp, PetscInt point, PetscDualSpace *bdsp)
106020cf1dd8SToby Isaac {
106120cf1dd8SToby Isaac   PetscErrorCode ierr;
106220cf1dd8SToby Isaac 
106320cf1dd8SToby Isaac   PetscFunctionBegin;
106420cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
106520cf1dd8SToby Isaac   PetscValidPointer(bdsp,2);
106620cf1dd8SToby Isaac   *bdsp = NULL;
106720cf1dd8SToby Isaac   if (sp->ops->getpointsubspace) {
106820cf1dd8SToby Isaac     ierr = (*sp->ops->getpointsubspace)(sp,point,bdsp);CHKERRQ(ierr);
106920cf1dd8SToby Isaac   } else if (sp->ops->getheightsubspace) {
107020cf1dd8SToby Isaac     DM       dm;
107120cf1dd8SToby Isaac     DMLabel  label;
107220cf1dd8SToby Isaac     PetscInt dim, depth, height;
107320cf1dd8SToby Isaac 
107420cf1dd8SToby Isaac     ierr = PetscDualSpaceGetDM(sp,&dm);CHKERRQ(ierr);
107520cf1dd8SToby Isaac     ierr = DMPlexGetDepth(dm,&dim);CHKERRQ(ierr);
107620cf1dd8SToby Isaac     ierr = DMPlexGetDepthLabel(dm,&label);CHKERRQ(ierr);
107720cf1dd8SToby Isaac     ierr = DMLabelGetValue(label,point,&depth);CHKERRQ(ierr);
107820cf1dd8SToby Isaac     height = dim - depth;
107920cf1dd8SToby Isaac     ierr = (*sp->ops->getheightsubspace)(sp,height,bdsp);CHKERRQ(ierr);
108020cf1dd8SToby Isaac   }
108120cf1dd8SToby Isaac   PetscFunctionReturn(0);
108220cf1dd8SToby Isaac }
108320cf1dd8SToby Isaac 
1084*6f905325SMatthew G. Knepley /*@C
1085*6f905325SMatthew G. Knepley   PetscDualSpaceGetSymmetries - Returns a description of the symmetries of this basis
1086*6f905325SMatthew G. Knepley 
1087*6f905325SMatthew G. Knepley   Not collective
1088*6f905325SMatthew G. Knepley 
1089*6f905325SMatthew G. Knepley   Input Parameter:
1090*6f905325SMatthew G. Knepley . sp - the PetscDualSpace object
1091*6f905325SMatthew G. Knepley 
1092*6f905325SMatthew G. Knepley   Output Parameters:
1093*6f905325SMatthew G. Knepley + perms - Permutations of the local degrees of freedom, parameterized by the point orientation
1094*6f905325SMatthew G. Knepley - flips - Sign reversal of the local degrees of freedom, parameterized by the point orientation
1095*6f905325SMatthew G. Knepley 
1096*6f905325SMatthew G. Knepley   Note: The permutation and flip arrays are organized in the following way
1097*6f905325SMatthew G. Knepley $ perms[p][ornt][dof # on point] = new local dof #
1098*6f905325SMatthew G. Knepley $ flips[p][ornt][dof # on point] = reversal or not
1099*6f905325SMatthew G. Knepley 
1100*6f905325SMatthew G. Knepley   Level: developer
1101*6f905325SMatthew G. Knepley 
1102*6f905325SMatthew G. Knepley .seealso: PetscDualSpaceSetSymmetries()
1103*6f905325SMatthew G. Knepley @*/
1104*6f905325SMatthew G. Knepley PetscErrorCode PetscDualSpaceGetSymmetries(PetscDualSpace sp, const PetscInt ****perms, const PetscScalar ****flips)
1105*6f905325SMatthew G. Knepley {
1106*6f905325SMatthew G. Knepley   PetscErrorCode ierr;
1107*6f905325SMatthew G. Knepley 
1108*6f905325SMatthew G. Knepley   PetscFunctionBegin;
1109*6f905325SMatthew G. Knepley   PetscValidHeaderSpecific(sp,PETSCDUALSPACE_CLASSID,1);
1110*6f905325SMatthew G. Knepley   if (perms) {PetscValidPointer(perms,2); *perms = NULL;}
1111*6f905325SMatthew G. Knepley   if (flips) {PetscValidPointer(flips,3); *flips = NULL;}
1112*6f905325SMatthew G. Knepley   if (sp->ops->getsymmetries) {ierr = (sp->ops->getsymmetries)(sp,perms,flips);CHKERRQ(ierr);}
1113*6f905325SMatthew G. Knepley   PetscFunctionReturn(0);
1114*6f905325SMatthew G. Knepley }
1115