xref: /petsc/src/dm/dt/dualspace/interface/dualspace.c (revision 063ee4adbc11f743c8903782eb67073d0c6dd346)
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 
1120cf1dd8SToby Isaac /*@C
1220cf1dd8SToby Isaac   PetscDualSpaceRegister - Adds a new PetscDualSpace implementation
1320cf1dd8SToby Isaac 
1420cf1dd8SToby Isaac   Not Collective
1520cf1dd8SToby Isaac 
1620cf1dd8SToby Isaac   Input Parameters:
1720cf1dd8SToby Isaac + name        - The name of a new user-defined creation routine
1820cf1dd8SToby Isaac - create_func - The creation routine itself
1920cf1dd8SToby Isaac 
2020cf1dd8SToby Isaac   Notes:
2120cf1dd8SToby Isaac   PetscDualSpaceRegister() may be called multiple times to add several user-defined PetscDualSpaces
2220cf1dd8SToby Isaac 
2320cf1dd8SToby Isaac   Sample usage:
2420cf1dd8SToby Isaac .vb
2520cf1dd8SToby Isaac     PetscDualSpaceRegister("my_space", MyPetscDualSpaceCreate);
2620cf1dd8SToby Isaac .ve
2720cf1dd8SToby Isaac 
2820cf1dd8SToby Isaac   Then, your PetscDualSpace type can be chosen with the procedural interface via
2920cf1dd8SToby Isaac .vb
3020cf1dd8SToby Isaac     PetscDualSpaceCreate(MPI_Comm, PetscDualSpace *);
3120cf1dd8SToby Isaac     PetscDualSpaceSetType(PetscDualSpace, "my_dual_space");
3220cf1dd8SToby Isaac .ve
3320cf1dd8SToby Isaac    or at runtime via the option
3420cf1dd8SToby Isaac .vb
3520cf1dd8SToby Isaac     -petscdualspace_type my_dual_space
3620cf1dd8SToby Isaac .ve
3720cf1dd8SToby Isaac 
3820cf1dd8SToby Isaac   Level: advanced
3920cf1dd8SToby Isaac 
4020cf1dd8SToby Isaac .keywords: PetscDualSpace, register
4120cf1dd8SToby Isaac .seealso: PetscDualSpaceRegisterAll(), PetscDualSpaceRegisterDestroy()
4220cf1dd8SToby Isaac 
4320cf1dd8SToby Isaac @*/
4420cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceRegister(const char sname[], PetscErrorCode (*function)(PetscDualSpace))
4520cf1dd8SToby Isaac {
4620cf1dd8SToby Isaac   PetscErrorCode ierr;
4720cf1dd8SToby Isaac 
4820cf1dd8SToby Isaac   PetscFunctionBegin;
4920cf1dd8SToby Isaac   ierr = PetscFunctionListAdd(&PetscDualSpaceList, sname, function);CHKERRQ(ierr);
5020cf1dd8SToby Isaac   PetscFunctionReturn(0);
5120cf1dd8SToby Isaac }
5220cf1dd8SToby Isaac 
5320cf1dd8SToby Isaac /*@C
5420cf1dd8SToby Isaac   PetscDualSpaceSetType - Builds a particular PetscDualSpace
5520cf1dd8SToby Isaac 
5620cf1dd8SToby Isaac   Collective on PetscDualSpace
5720cf1dd8SToby Isaac 
5820cf1dd8SToby Isaac   Input Parameters:
5920cf1dd8SToby Isaac + sp   - The PetscDualSpace object
6020cf1dd8SToby Isaac - name - The kind of space
6120cf1dd8SToby Isaac 
6220cf1dd8SToby Isaac   Options Database Key:
6320cf1dd8SToby Isaac . -petscdualspace_type <type> - Sets the PetscDualSpace type; use -help for a list of available types
6420cf1dd8SToby Isaac 
6520cf1dd8SToby Isaac   Level: intermediate
6620cf1dd8SToby Isaac 
6720cf1dd8SToby Isaac .keywords: PetscDualSpace, set, type
6820cf1dd8SToby Isaac .seealso: PetscDualSpaceGetType(), PetscDualSpaceCreate()
6920cf1dd8SToby Isaac @*/
7020cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetType(PetscDualSpace sp, PetscDualSpaceType name)
7120cf1dd8SToby Isaac {
7220cf1dd8SToby Isaac   PetscErrorCode (*r)(PetscDualSpace);
7320cf1dd8SToby Isaac   PetscBool      match;
7420cf1dd8SToby Isaac   PetscErrorCode ierr;
7520cf1dd8SToby Isaac 
7620cf1dd8SToby Isaac   PetscFunctionBegin;
7720cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
7820cf1dd8SToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) sp, name, &match);CHKERRQ(ierr);
7920cf1dd8SToby Isaac   if (match) PetscFunctionReturn(0);
8020cf1dd8SToby Isaac 
8120cf1dd8SToby Isaac   if (!PetscDualSpaceRegisterAllCalled) {ierr = PetscDualSpaceRegisterAll();CHKERRQ(ierr);}
8220cf1dd8SToby Isaac   ierr = PetscFunctionListFind(PetscDualSpaceList, name, &r);CHKERRQ(ierr);
8320cf1dd8SToby Isaac   if (!r) SETERRQ1(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscDualSpace type: %s", name);
8420cf1dd8SToby Isaac 
8520cf1dd8SToby Isaac   if (sp->ops->destroy) {
8620cf1dd8SToby Isaac     ierr             = (*sp->ops->destroy)(sp);CHKERRQ(ierr);
8720cf1dd8SToby Isaac     sp->ops->destroy = NULL;
8820cf1dd8SToby Isaac   }
8920cf1dd8SToby Isaac   ierr = (*r)(sp);CHKERRQ(ierr);
9020cf1dd8SToby Isaac   ierr = PetscObjectChangeTypeName((PetscObject) sp, name);CHKERRQ(ierr);
9120cf1dd8SToby Isaac   PetscFunctionReturn(0);
9220cf1dd8SToby Isaac }
9320cf1dd8SToby Isaac 
9420cf1dd8SToby Isaac /*@C
9520cf1dd8SToby Isaac   PetscDualSpaceGetType - Gets the PetscDualSpace type name (as a string) from the object.
9620cf1dd8SToby Isaac 
9720cf1dd8SToby Isaac   Not Collective
9820cf1dd8SToby Isaac 
9920cf1dd8SToby Isaac   Input Parameter:
10020cf1dd8SToby Isaac . sp  - The PetscDualSpace
10120cf1dd8SToby Isaac 
10220cf1dd8SToby Isaac   Output Parameter:
10320cf1dd8SToby Isaac . name - The PetscDualSpace type name
10420cf1dd8SToby Isaac 
10520cf1dd8SToby Isaac   Level: intermediate
10620cf1dd8SToby Isaac 
10720cf1dd8SToby Isaac .keywords: PetscDualSpace, get, type, name
10820cf1dd8SToby Isaac .seealso: PetscDualSpaceSetType(), PetscDualSpaceCreate()
10920cf1dd8SToby Isaac @*/
11020cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetType(PetscDualSpace sp, PetscDualSpaceType *name)
11120cf1dd8SToby Isaac {
11220cf1dd8SToby Isaac   PetscErrorCode ierr;
11320cf1dd8SToby Isaac 
11420cf1dd8SToby Isaac   PetscFunctionBegin;
11520cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
11620cf1dd8SToby Isaac   PetscValidPointer(name, 2);
11720cf1dd8SToby Isaac   if (!PetscDualSpaceRegisterAllCalled) {
11820cf1dd8SToby Isaac     ierr = PetscDualSpaceRegisterAll();CHKERRQ(ierr);
11920cf1dd8SToby Isaac   }
12020cf1dd8SToby Isaac   *name = ((PetscObject) sp)->type_name;
12120cf1dd8SToby Isaac   PetscFunctionReturn(0);
12220cf1dd8SToby Isaac }
12320cf1dd8SToby Isaac 
12420cf1dd8SToby Isaac /*@
12520cf1dd8SToby Isaac   PetscDualSpaceView - Views a PetscDualSpace
12620cf1dd8SToby Isaac 
12720cf1dd8SToby Isaac   Collective on PetscDualSpace
12820cf1dd8SToby Isaac 
12920cf1dd8SToby Isaac   Input Parameter:
13020cf1dd8SToby Isaac + sp - the PetscDualSpace object to view
13120cf1dd8SToby Isaac - v  - the viewer
13220cf1dd8SToby Isaac 
13320cf1dd8SToby Isaac   Level: developer
13420cf1dd8SToby Isaac 
13520cf1dd8SToby Isaac .seealso PetscDualSpaceDestroy()
13620cf1dd8SToby Isaac @*/
13720cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceView(PetscDualSpace sp, PetscViewer v)
13820cf1dd8SToby Isaac {
139d9bac1caSLisandro Dalcin   PetscBool      iascii;
14020cf1dd8SToby Isaac   PetscErrorCode ierr;
14120cf1dd8SToby Isaac 
14220cf1dd8SToby Isaac   PetscFunctionBegin;
14320cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
144d9bac1caSLisandro Dalcin   if (v) PetscValidHeaderSpecific(v, PETSC_VIEWER_CLASSID, 2);
14520cf1dd8SToby Isaac   if (!v) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) sp), &v);CHKERRQ(ierr);}
146d9bac1caSLisandro Dalcin   ierr = PetscObjectPrintClassNamePrefixType((PetscObject)sp, v);CHKERRQ(ierr);
147d9bac1caSLisandro Dalcin   ierr = PetscObjectTypeCompare((PetscObject) v, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
148d9bac1caSLisandro Dalcin   ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr);
149d9bac1caSLisandro Dalcin   if (iascii) {ierr = PetscViewerASCIIPrintf(v, "Dual space of order %D with %D components\n", sp->order, sp->Nc);CHKERRQ(ierr);}
15020cf1dd8SToby Isaac   if (sp->ops->view) {ierr = (*sp->ops->view)(sp, v);CHKERRQ(ierr);}
151d9bac1caSLisandro Dalcin   ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr);
15220cf1dd8SToby Isaac   PetscFunctionReturn(0);
15320cf1dd8SToby Isaac }
15420cf1dd8SToby Isaac 
15520cf1dd8SToby Isaac /*@
15620cf1dd8SToby Isaac   PetscDualSpaceSetFromOptions - sets parameters in a PetscDualSpace from the options database
15720cf1dd8SToby Isaac 
15820cf1dd8SToby Isaac   Collective on PetscDualSpace
15920cf1dd8SToby Isaac 
16020cf1dd8SToby Isaac   Input Parameter:
16120cf1dd8SToby Isaac . sp - the PetscDualSpace object to set options for
16220cf1dd8SToby Isaac 
16320cf1dd8SToby Isaac   Options Database:
1647be5e748SToby Isaac . -petscspace_degree the approximation order of the space
16520cf1dd8SToby Isaac 
16620cf1dd8SToby Isaac   Level: developer
16720cf1dd8SToby Isaac 
16820cf1dd8SToby Isaac .seealso PetscDualSpaceView()
16920cf1dd8SToby Isaac @*/
17020cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetFromOptions(PetscDualSpace sp)
17120cf1dd8SToby Isaac {
172*063ee4adSMatthew G. Knepley   PetscDualSpaceReferenceCell refCell = PETSCDUALSPACE_REFCELL_SIMPLEX;
173*063ee4adSMatthew G. Knepley   PetscInt                    refDim  = 0;
174*063ee4adSMatthew G. Knepley   PetscBool                   flg;
17520cf1dd8SToby Isaac   const char                 *defaultType;
17620cf1dd8SToby Isaac   char                        name[256];
17720cf1dd8SToby Isaac   PetscErrorCode              ierr;
17820cf1dd8SToby Isaac 
17920cf1dd8SToby Isaac   PetscFunctionBegin;
18020cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
18120cf1dd8SToby Isaac   if (!((PetscObject) sp)->type_name) {
18220cf1dd8SToby Isaac     defaultType = PETSCDUALSPACELAGRANGE;
18320cf1dd8SToby Isaac   } else {
18420cf1dd8SToby Isaac     defaultType = ((PetscObject) sp)->type_name;
18520cf1dd8SToby Isaac   }
18620cf1dd8SToby Isaac   if (!PetscSpaceRegisterAllCalled) {ierr = PetscSpaceRegisterAll();CHKERRQ(ierr);}
18720cf1dd8SToby Isaac 
18820cf1dd8SToby Isaac   ierr = PetscObjectOptionsBegin((PetscObject) sp);CHKERRQ(ierr);
18920cf1dd8SToby Isaac   ierr = PetscOptionsFList("-petscdualspace_type", "Dual space", "PetscDualSpaceSetType", PetscDualSpaceList, defaultType, name, 256, &flg);CHKERRQ(ierr);
19020cf1dd8SToby Isaac   if (flg) {
19120cf1dd8SToby Isaac     ierr = PetscDualSpaceSetType(sp, name);CHKERRQ(ierr);
19220cf1dd8SToby Isaac   } else if (!((PetscObject) sp)->type_name) {
19320cf1dd8SToby Isaac     ierr = PetscDualSpaceSetType(sp, defaultType);CHKERRQ(ierr);
19420cf1dd8SToby Isaac   }
1957be5e748SToby Isaac   ierr = PetscOptionsInt("-petscdualspace_degree", "The approximation order", "PetscDualSpaceSetOrder", sp->order, &sp->order, NULL);CHKERRQ(ierr);
19620cf1dd8SToby Isaac   ierr = PetscOptionsInt("-petscdualspace_components", "The number of components", "PetscDualSpaceSetNumComponents", sp->Nc, &sp->Nc, NULL);CHKERRQ(ierr);
19720cf1dd8SToby Isaac   if (sp->ops->setfromoptions) {
19820cf1dd8SToby Isaac     ierr = (*sp->ops->setfromoptions)(PetscOptionsObject,sp);CHKERRQ(ierr);
19920cf1dd8SToby Isaac   }
200*063ee4adSMatthew G. Knepley   ierr = PetscOptionsInt("-petscdualspace_refdim", "The spatial dimension of the reference cell", "PetscDualSpaceSetReferenceCell", refDim, &refDim, NULL);CHKERRQ(ierr);
201*063ee4adSMatthew G. Knepley   ierr = PetscOptionsEnum("-petscdualspace_refcell", "Reference cell", "PetscDualSpaceSetReferenceCell", PetscDualSpaceReferenceCells, (PetscEnum) refCell, (PetscEnum *) &refCell, &flg);CHKERRQ(ierr);
202*063ee4adSMatthew G. Knepley   if (flg) {
203*063ee4adSMatthew G. Knepley     DM K;
204*063ee4adSMatthew G. Knepley 
205*063ee4adSMatthew G. Knepley     if (!refDim) SETERRQ(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_INCOMP, "Reference cell specified without a dimension. Use -petscdualspace_refdim.");
206*063ee4adSMatthew G. Knepley     ierr = PetscDualSpaceCreateReferenceCell(sp, refDim, refCell == PETSCDUALSPACE_REFCELL_SIMPLEX ? PETSC_TRUE : PETSC_FALSE, &K);CHKERRQ(ierr);
207*063ee4adSMatthew G. Knepley     ierr = PetscDualSpaceSetDM(sp, K);CHKERRQ(ierr);
208*063ee4adSMatthew G. Knepley     ierr = DMDestroy(&K);CHKERRQ(ierr);
209*063ee4adSMatthew G. Knepley   }
210*063ee4adSMatthew G. Knepley 
21120cf1dd8SToby Isaac   /* process any options handlers added with PetscObjectAddOptionsHandler() */
21220cf1dd8SToby Isaac   ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) sp);CHKERRQ(ierr);
21320cf1dd8SToby Isaac   ierr = PetscOptionsEnd();CHKERRQ(ierr);
214*063ee4adSMatthew G. Knepley   sp->setfromoptionscalled = PETSC_TRUE;
21520cf1dd8SToby Isaac   PetscFunctionReturn(0);
21620cf1dd8SToby Isaac }
21720cf1dd8SToby Isaac 
21820cf1dd8SToby Isaac /*@
21920cf1dd8SToby Isaac   PetscDualSpaceSetUp - Construct a basis for the PetscDualSpace
22020cf1dd8SToby Isaac 
22120cf1dd8SToby Isaac   Collective on PetscDualSpace
22220cf1dd8SToby Isaac 
22320cf1dd8SToby Isaac   Input Parameter:
22420cf1dd8SToby Isaac . sp - the PetscDualSpace object to setup
22520cf1dd8SToby Isaac 
22620cf1dd8SToby Isaac   Level: developer
22720cf1dd8SToby Isaac 
22820cf1dd8SToby Isaac .seealso PetscDualSpaceView(), PetscDualSpaceDestroy()
22920cf1dd8SToby Isaac @*/
23020cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetUp(PetscDualSpace sp)
23120cf1dd8SToby Isaac {
23220cf1dd8SToby Isaac   PetscErrorCode ierr;
23320cf1dd8SToby Isaac 
23420cf1dd8SToby Isaac   PetscFunctionBegin;
23520cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
23620cf1dd8SToby Isaac   if (sp->setupcalled) PetscFunctionReturn(0);
23720cf1dd8SToby Isaac   sp->setupcalled = PETSC_TRUE;
23820cf1dd8SToby Isaac   if (sp->ops->setup) {ierr = (*sp->ops->setup)(sp);CHKERRQ(ierr);}
239*063ee4adSMatthew G. Knepley   if (sp->setfromoptionscalled) {ierr = PetscDualSpaceViewFromOptions(sp, NULL, "-petscdualspace_view");CHKERRQ(ierr);}
24020cf1dd8SToby Isaac   PetscFunctionReturn(0);
24120cf1dd8SToby Isaac }
24220cf1dd8SToby Isaac 
24320cf1dd8SToby Isaac /*@
24420cf1dd8SToby Isaac   PetscDualSpaceDestroy - Destroys a PetscDualSpace object
24520cf1dd8SToby Isaac 
24620cf1dd8SToby Isaac   Collective on PetscDualSpace
24720cf1dd8SToby Isaac 
24820cf1dd8SToby Isaac   Input Parameter:
24920cf1dd8SToby Isaac . sp - the PetscDualSpace object to destroy
25020cf1dd8SToby Isaac 
25120cf1dd8SToby Isaac   Level: developer
25220cf1dd8SToby Isaac 
25320cf1dd8SToby Isaac .seealso PetscDualSpaceView()
25420cf1dd8SToby Isaac @*/
25520cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceDestroy(PetscDualSpace *sp)
25620cf1dd8SToby Isaac {
25720cf1dd8SToby Isaac   PetscInt       dim, f;
25820cf1dd8SToby Isaac   PetscErrorCode ierr;
25920cf1dd8SToby Isaac 
26020cf1dd8SToby Isaac   PetscFunctionBegin;
26120cf1dd8SToby Isaac   if (!*sp) PetscFunctionReturn(0);
26220cf1dd8SToby Isaac   PetscValidHeaderSpecific((*sp), PETSCDUALSPACE_CLASSID, 1);
26320cf1dd8SToby Isaac 
26420cf1dd8SToby Isaac   if (--((PetscObject)(*sp))->refct > 0) {*sp = 0; PetscFunctionReturn(0);}
26520cf1dd8SToby Isaac   ((PetscObject) (*sp))->refct = 0;
26620cf1dd8SToby Isaac 
26720cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDimension(*sp, &dim);CHKERRQ(ierr);
26820cf1dd8SToby Isaac   for (f = 0; f < dim; ++f) {
26920cf1dd8SToby Isaac     ierr = PetscQuadratureDestroy(&(*sp)->functional[f]);CHKERRQ(ierr);
27020cf1dd8SToby Isaac   }
27120cf1dd8SToby Isaac   ierr = PetscFree((*sp)->functional);CHKERRQ(ierr);
27220cf1dd8SToby Isaac   ierr = PetscQuadratureDestroy(&(*sp)->allPoints);CHKERRQ(ierr);
27320cf1dd8SToby Isaac   ierr = DMDestroy(&(*sp)->dm);CHKERRQ(ierr);
27420cf1dd8SToby Isaac 
27520cf1dd8SToby Isaac   if ((*sp)->ops->destroy) {ierr = (*(*sp)->ops->destroy)(*sp);CHKERRQ(ierr);}
27620cf1dd8SToby Isaac   ierr = PetscHeaderDestroy(sp);CHKERRQ(ierr);
27720cf1dd8SToby Isaac   PetscFunctionReturn(0);
27820cf1dd8SToby Isaac }
27920cf1dd8SToby Isaac 
28020cf1dd8SToby Isaac /*@
28120cf1dd8SToby Isaac   PetscDualSpaceCreate - Creates an empty PetscDualSpace object. The type can then be set with PetscDualSpaceSetType().
28220cf1dd8SToby Isaac 
28320cf1dd8SToby Isaac   Collective on MPI_Comm
28420cf1dd8SToby Isaac 
28520cf1dd8SToby Isaac   Input Parameter:
28620cf1dd8SToby Isaac . comm - The communicator for the PetscDualSpace object
28720cf1dd8SToby Isaac 
28820cf1dd8SToby Isaac   Output Parameter:
28920cf1dd8SToby Isaac . sp - The PetscDualSpace object
29020cf1dd8SToby Isaac 
29120cf1dd8SToby Isaac   Level: beginner
29220cf1dd8SToby Isaac 
29320cf1dd8SToby Isaac .seealso: PetscDualSpaceSetType(), PETSCDUALSPACELAGRANGE
29420cf1dd8SToby Isaac @*/
29520cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceCreate(MPI_Comm comm, PetscDualSpace *sp)
29620cf1dd8SToby Isaac {
29720cf1dd8SToby Isaac   PetscDualSpace s;
29820cf1dd8SToby Isaac   PetscErrorCode ierr;
29920cf1dd8SToby Isaac 
30020cf1dd8SToby Isaac   PetscFunctionBegin;
30120cf1dd8SToby Isaac   PetscValidPointer(sp, 2);
30220cf1dd8SToby Isaac   ierr = PetscCitationsRegister(FECitation,&FEcite);CHKERRQ(ierr);
30320cf1dd8SToby Isaac   *sp  = NULL;
30420cf1dd8SToby Isaac   ierr = PetscFEInitializePackage();CHKERRQ(ierr);
30520cf1dd8SToby Isaac 
30620cf1dd8SToby Isaac   ierr = PetscHeaderCreate(s, PETSCDUALSPACE_CLASSID, "PetscDualSpace", "Dual Space", "PetscDualSpace", comm, PetscDualSpaceDestroy, PetscDualSpaceView);CHKERRQ(ierr);
30720cf1dd8SToby Isaac 
30820cf1dd8SToby Isaac   s->order = 0;
30920cf1dd8SToby Isaac   s->Nc    = 1;
31020cf1dd8SToby Isaac   s->setupcalled = PETSC_FALSE;
31120cf1dd8SToby Isaac 
31220cf1dd8SToby Isaac   *sp = s;
31320cf1dd8SToby Isaac   PetscFunctionReturn(0);
31420cf1dd8SToby Isaac }
31520cf1dd8SToby Isaac 
31620cf1dd8SToby Isaac /*@
31720cf1dd8SToby Isaac   PetscDualSpaceDuplicate - Creates a duplicate PetscDualSpace object, however it is not setup.
31820cf1dd8SToby Isaac 
31920cf1dd8SToby Isaac   Collective on PetscDualSpace
32020cf1dd8SToby Isaac 
32120cf1dd8SToby Isaac   Input Parameter:
32220cf1dd8SToby Isaac . sp - The original PetscDualSpace
32320cf1dd8SToby Isaac 
32420cf1dd8SToby Isaac   Output Parameter:
32520cf1dd8SToby Isaac . spNew - The duplicate PetscDualSpace
32620cf1dd8SToby Isaac 
32720cf1dd8SToby Isaac   Level: beginner
32820cf1dd8SToby Isaac 
32920cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate(), PetscDualSpaceSetType()
33020cf1dd8SToby Isaac @*/
33120cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceDuplicate(PetscDualSpace sp, PetscDualSpace *spNew)
33220cf1dd8SToby Isaac {
33320cf1dd8SToby Isaac   PetscErrorCode ierr;
33420cf1dd8SToby Isaac 
33520cf1dd8SToby Isaac   PetscFunctionBegin;
33620cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
33720cf1dd8SToby Isaac   PetscValidPointer(spNew, 2);
33820cf1dd8SToby Isaac   ierr = (*sp->ops->duplicate)(sp, spNew);CHKERRQ(ierr);
33920cf1dd8SToby Isaac   PetscFunctionReturn(0);
34020cf1dd8SToby Isaac }
34120cf1dd8SToby Isaac 
34220cf1dd8SToby Isaac /*@
34320cf1dd8SToby Isaac   PetscDualSpaceGetDM - Get the DM representing the reference cell
34420cf1dd8SToby Isaac 
34520cf1dd8SToby Isaac   Not collective
34620cf1dd8SToby Isaac 
34720cf1dd8SToby Isaac   Input Parameter:
34820cf1dd8SToby Isaac . sp - The PetscDualSpace
34920cf1dd8SToby Isaac 
35020cf1dd8SToby Isaac   Output Parameter:
35120cf1dd8SToby Isaac . dm - The reference cell
35220cf1dd8SToby Isaac 
35320cf1dd8SToby Isaac   Level: intermediate
35420cf1dd8SToby Isaac 
35520cf1dd8SToby Isaac .seealso: PetscDualSpaceSetDM(), PetscDualSpaceCreate()
35620cf1dd8SToby Isaac @*/
35720cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetDM(PetscDualSpace sp, DM *dm)
35820cf1dd8SToby Isaac {
35920cf1dd8SToby Isaac   PetscFunctionBegin;
36020cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
36120cf1dd8SToby Isaac   PetscValidPointer(dm, 2);
36220cf1dd8SToby Isaac   *dm = sp->dm;
36320cf1dd8SToby Isaac   PetscFunctionReturn(0);
36420cf1dd8SToby Isaac }
36520cf1dd8SToby Isaac 
36620cf1dd8SToby Isaac /*@
36720cf1dd8SToby Isaac   PetscDualSpaceSetDM - Get the DM representing the reference cell
36820cf1dd8SToby Isaac 
36920cf1dd8SToby Isaac   Not collective
37020cf1dd8SToby Isaac 
37120cf1dd8SToby Isaac   Input Parameters:
37220cf1dd8SToby Isaac + sp - The PetscDualSpace
37320cf1dd8SToby Isaac - dm - The reference cell
37420cf1dd8SToby Isaac 
37520cf1dd8SToby Isaac   Level: intermediate
37620cf1dd8SToby Isaac 
37720cf1dd8SToby Isaac .seealso: PetscDualSpaceGetDM(), PetscDualSpaceCreate()
37820cf1dd8SToby Isaac @*/
37920cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetDM(PetscDualSpace sp, DM dm)
38020cf1dd8SToby Isaac {
38120cf1dd8SToby Isaac   PetscErrorCode ierr;
38220cf1dd8SToby Isaac 
38320cf1dd8SToby Isaac   PetscFunctionBegin;
38420cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
38520cf1dd8SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
38620cf1dd8SToby Isaac   ierr = DMDestroy(&sp->dm);CHKERRQ(ierr);
38720cf1dd8SToby Isaac   ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);
38820cf1dd8SToby Isaac   sp->dm = dm;
38920cf1dd8SToby Isaac   PetscFunctionReturn(0);
39020cf1dd8SToby Isaac }
39120cf1dd8SToby Isaac 
39220cf1dd8SToby Isaac /*@
39320cf1dd8SToby Isaac   PetscDualSpaceGetOrder - Get the order of the dual space
39420cf1dd8SToby Isaac 
39520cf1dd8SToby Isaac   Not collective
39620cf1dd8SToby Isaac 
39720cf1dd8SToby Isaac   Input Parameter:
39820cf1dd8SToby Isaac . sp - The PetscDualSpace
39920cf1dd8SToby Isaac 
40020cf1dd8SToby Isaac   Output Parameter:
40120cf1dd8SToby Isaac . order - The order
40220cf1dd8SToby Isaac 
40320cf1dd8SToby Isaac   Level: intermediate
40420cf1dd8SToby Isaac 
40520cf1dd8SToby Isaac .seealso: PetscDualSpaceSetOrder(), PetscDualSpaceCreate()
40620cf1dd8SToby Isaac @*/
40720cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetOrder(PetscDualSpace sp, PetscInt *order)
40820cf1dd8SToby Isaac {
40920cf1dd8SToby Isaac   PetscFunctionBegin;
41020cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
41120cf1dd8SToby Isaac   PetscValidPointer(order, 2);
41220cf1dd8SToby Isaac   *order = sp->order;
41320cf1dd8SToby Isaac   PetscFunctionReturn(0);
41420cf1dd8SToby Isaac }
41520cf1dd8SToby Isaac 
41620cf1dd8SToby Isaac /*@
41720cf1dd8SToby Isaac   PetscDualSpaceSetOrder - Set the order of the dual space
41820cf1dd8SToby Isaac 
41920cf1dd8SToby Isaac   Not collective
42020cf1dd8SToby Isaac 
42120cf1dd8SToby Isaac   Input Parameters:
42220cf1dd8SToby Isaac + sp - The PetscDualSpace
42320cf1dd8SToby Isaac - order - The order
42420cf1dd8SToby Isaac 
42520cf1dd8SToby Isaac   Level: intermediate
42620cf1dd8SToby Isaac 
42720cf1dd8SToby Isaac .seealso: PetscDualSpaceGetOrder(), PetscDualSpaceCreate()
42820cf1dd8SToby Isaac @*/
42920cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetOrder(PetscDualSpace sp, PetscInt order)
43020cf1dd8SToby Isaac {
43120cf1dd8SToby Isaac   PetscFunctionBegin;
43220cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
43320cf1dd8SToby Isaac   sp->order = order;
43420cf1dd8SToby Isaac   PetscFunctionReturn(0);
43520cf1dd8SToby Isaac }
43620cf1dd8SToby Isaac 
43720cf1dd8SToby Isaac /*@
43820cf1dd8SToby Isaac   PetscDualSpaceGetNumComponents - Return the number of components for this space
43920cf1dd8SToby Isaac 
44020cf1dd8SToby Isaac   Input Parameter:
44120cf1dd8SToby Isaac . sp - The PetscDualSpace
44220cf1dd8SToby Isaac 
44320cf1dd8SToby Isaac   Output Parameter:
44420cf1dd8SToby Isaac . Nc - The number of components
44520cf1dd8SToby Isaac 
44620cf1dd8SToby Isaac   Note: A vector space, for example, will have d components, where d is the spatial dimension
44720cf1dd8SToby Isaac 
44820cf1dd8SToby Isaac   Level: intermediate
44920cf1dd8SToby Isaac 
45020cf1dd8SToby Isaac .seealso: PetscDualSpaceSetNumComponents(), PetscDualSpaceGetDimension(), PetscDualSpaceCreate(), PetscDualSpace
45120cf1dd8SToby Isaac @*/
45220cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetNumComponents(PetscDualSpace sp, PetscInt *Nc)
45320cf1dd8SToby Isaac {
45420cf1dd8SToby Isaac   PetscFunctionBegin;
45520cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
45620cf1dd8SToby Isaac   PetscValidPointer(Nc, 2);
45720cf1dd8SToby Isaac   *Nc = sp->Nc;
45820cf1dd8SToby Isaac   PetscFunctionReturn(0);
45920cf1dd8SToby Isaac }
46020cf1dd8SToby Isaac 
46120cf1dd8SToby Isaac /*@
46220cf1dd8SToby Isaac   PetscDualSpaceSetNumComponents - Set the number of components for this space
46320cf1dd8SToby Isaac 
46420cf1dd8SToby Isaac   Input Parameters:
46520cf1dd8SToby Isaac + sp - The PetscDualSpace
46620cf1dd8SToby Isaac - order - The number of components
46720cf1dd8SToby Isaac 
46820cf1dd8SToby Isaac   Level: intermediate
46920cf1dd8SToby Isaac 
47020cf1dd8SToby Isaac .seealso: PetscDualSpaceGetNumComponents(), PetscDualSpaceCreate(), PetscDualSpace
47120cf1dd8SToby Isaac @*/
47220cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetNumComponents(PetscDualSpace sp, PetscInt Nc)
47320cf1dd8SToby Isaac {
47420cf1dd8SToby Isaac   PetscFunctionBegin;
47520cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
47620cf1dd8SToby Isaac   sp->Nc = Nc;
47720cf1dd8SToby Isaac   PetscFunctionReturn(0);
47820cf1dd8SToby Isaac }
47920cf1dd8SToby Isaac 
48020cf1dd8SToby Isaac /*@
48120cf1dd8SToby Isaac   PetscDualSpaceLagrangeGetTensor - Get the tensor nature of the dual space
48220cf1dd8SToby Isaac 
48320cf1dd8SToby Isaac   Not collective
48420cf1dd8SToby Isaac 
48520cf1dd8SToby Isaac   Input Parameter:
48620cf1dd8SToby Isaac . sp - The PetscDualSpace
48720cf1dd8SToby Isaac 
48820cf1dd8SToby Isaac   Output Parameter:
48920cf1dd8SToby Isaac . tensor - Whether the dual space has tensor layout (vs. simplicial)
49020cf1dd8SToby Isaac 
49120cf1dd8SToby Isaac   Level: intermediate
49220cf1dd8SToby Isaac 
49320cf1dd8SToby Isaac .seealso: PetscDualSpaceLagrangeSetTensor(), PetscDualSpaceCreate()
49420cf1dd8SToby Isaac @*/
49520cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceLagrangeGetTensor(PetscDualSpace sp, PetscBool *tensor)
49620cf1dd8SToby Isaac {
49720cf1dd8SToby Isaac   PetscErrorCode ierr;
49820cf1dd8SToby Isaac 
49920cf1dd8SToby Isaac   PetscFunctionBegin;
50020cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
50120cf1dd8SToby Isaac   PetscValidPointer(tensor, 2);
50220cf1dd8SToby Isaac   ierr = PetscTryMethod(sp,"PetscDualSpaceLagrangeGetTensor_C",(PetscDualSpace,PetscBool *),(sp,tensor));CHKERRQ(ierr);
50320cf1dd8SToby Isaac   PetscFunctionReturn(0);
50420cf1dd8SToby Isaac }
50520cf1dd8SToby Isaac 
50620cf1dd8SToby Isaac /*@
50720cf1dd8SToby Isaac   PetscDualSpaceLagrangeSetTensor - Set the tensor nature of the dual space
50820cf1dd8SToby Isaac 
50920cf1dd8SToby Isaac   Not collective
51020cf1dd8SToby Isaac 
51120cf1dd8SToby Isaac   Input Parameters:
51220cf1dd8SToby Isaac + sp - The PetscDualSpace
51320cf1dd8SToby Isaac - tensor - Whether the dual space has tensor layout (vs. simplicial)
51420cf1dd8SToby Isaac 
51520cf1dd8SToby Isaac   Level: intermediate
51620cf1dd8SToby Isaac 
51720cf1dd8SToby Isaac .seealso: PetscDualSpaceLagrangeGetTensor(), PetscDualSpaceCreate()
51820cf1dd8SToby Isaac @*/
51920cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceLagrangeSetTensor(PetscDualSpace sp, PetscBool tensor)
52020cf1dd8SToby Isaac {
52120cf1dd8SToby Isaac   PetscErrorCode ierr;
52220cf1dd8SToby Isaac 
52320cf1dd8SToby Isaac   PetscFunctionBegin;
52420cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
52520cf1dd8SToby Isaac   ierr = PetscTryMethod(sp,"PetscDualSpaceLagrangeSetTensor_C",(PetscDualSpace,PetscBool),(sp,tensor));CHKERRQ(ierr);
52620cf1dd8SToby Isaac   PetscFunctionReturn(0);
52720cf1dd8SToby Isaac }
52820cf1dd8SToby Isaac 
52920cf1dd8SToby Isaac /*@
53020cf1dd8SToby Isaac   PetscDualSpaceGetFunctional - Get the i-th basis functional in the dual space
53120cf1dd8SToby Isaac 
53220cf1dd8SToby Isaac   Not collective
53320cf1dd8SToby Isaac 
53420cf1dd8SToby Isaac   Input Parameters:
53520cf1dd8SToby Isaac + sp - The PetscDualSpace
53620cf1dd8SToby Isaac - i  - The basis number
53720cf1dd8SToby Isaac 
53820cf1dd8SToby Isaac   Output Parameter:
53920cf1dd8SToby Isaac . functional - The basis functional
54020cf1dd8SToby Isaac 
54120cf1dd8SToby Isaac   Level: intermediate
54220cf1dd8SToby Isaac 
54320cf1dd8SToby Isaac .seealso: PetscDualSpaceGetDimension(), PetscDualSpaceCreate()
54420cf1dd8SToby Isaac @*/
54520cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetFunctional(PetscDualSpace sp, PetscInt i, PetscQuadrature *functional)
54620cf1dd8SToby Isaac {
54720cf1dd8SToby Isaac   PetscInt       dim;
54820cf1dd8SToby Isaac   PetscErrorCode ierr;
54920cf1dd8SToby Isaac 
55020cf1dd8SToby Isaac   PetscFunctionBegin;
55120cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
55220cf1dd8SToby Isaac   PetscValidPointer(functional, 3);
55320cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDimension(sp, &dim);CHKERRQ(ierr);
55420cf1dd8SToby Isaac   if ((i < 0) || (i >= dim)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Functional index %d must be in [0, %d)", i, dim);
55520cf1dd8SToby Isaac   *functional = sp->functional[i];
55620cf1dd8SToby Isaac   PetscFunctionReturn(0);
55720cf1dd8SToby Isaac }
55820cf1dd8SToby Isaac 
55920cf1dd8SToby Isaac /*@
56020cf1dd8SToby Isaac   PetscDualSpaceGetDimension - Get the dimension of the dual space, i.e. the number of basis functionals
56120cf1dd8SToby Isaac 
56220cf1dd8SToby Isaac   Not collective
56320cf1dd8SToby Isaac 
56420cf1dd8SToby Isaac   Input Parameter:
56520cf1dd8SToby Isaac . sp - The PetscDualSpace
56620cf1dd8SToby Isaac 
56720cf1dd8SToby Isaac   Output Parameter:
56820cf1dd8SToby Isaac . dim - The dimension
56920cf1dd8SToby Isaac 
57020cf1dd8SToby Isaac   Level: intermediate
57120cf1dd8SToby Isaac 
57220cf1dd8SToby Isaac .seealso: PetscDualSpaceGetFunctional(), PetscDualSpaceCreate()
57320cf1dd8SToby Isaac @*/
57420cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetDimension(PetscDualSpace sp, PetscInt *dim)
57520cf1dd8SToby Isaac {
57620cf1dd8SToby Isaac   PetscErrorCode ierr;
57720cf1dd8SToby Isaac 
57820cf1dd8SToby Isaac   PetscFunctionBegin;
57920cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
58020cf1dd8SToby Isaac   PetscValidPointer(dim, 2);
58120cf1dd8SToby Isaac   *dim = 0;
58220cf1dd8SToby Isaac   if (sp->ops->getdimension) {ierr = (*sp->ops->getdimension)(sp, dim);CHKERRQ(ierr);}
58320cf1dd8SToby Isaac   PetscFunctionReturn(0);
58420cf1dd8SToby Isaac }
58520cf1dd8SToby Isaac 
58620cf1dd8SToby Isaac /*@C
58720cf1dd8SToby Isaac   PetscDualSpaceGetNumDof - Get the number of degrees of freedom for each spatial (topological) dimension
58820cf1dd8SToby Isaac 
58920cf1dd8SToby Isaac   Not collective
59020cf1dd8SToby Isaac 
59120cf1dd8SToby Isaac   Input Parameter:
59220cf1dd8SToby Isaac . sp - The PetscDualSpace
59320cf1dd8SToby Isaac 
59420cf1dd8SToby Isaac   Output Parameter:
59520cf1dd8SToby Isaac . numDof - An array of length dim+1 which holds the number of dofs for each dimension
59620cf1dd8SToby Isaac 
59720cf1dd8SToby Isaac   Level: intermediate
59820cf1dd8SToby Isaac 
59920cf1dd8SToby Isaac .seealso: PetscDualSpaceGetFunctional(), PetscDualSpaceCreate()
60020cf1dd8SToby Isaac @*/
60120cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetNumDof(PetscDualSpace sp, const PetscInt **numDof)
60220cf1dd8SToby Isaac {
60320cf1dd8SToby Isaac   PetscErrorCode ierr;
60420cf1dd8SToby Isaac 
60520cf1dd8SToby Isaac   PetscFunctionBegin;
60620cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
60720cf1dd8SToby Isaac   PetscValidPointer(numDof, 2);
60820cf1dd8SToby Isaac   ierr = (*sp->ops->getnumdof)(sp, numDof);CHKERRQ(ierr);
60920cf1dd8SToby Isaac   if (!*numDof) SETERRQ(PetscObjectComm((PetscObject) sp), PETSC_ERR_LIB, "Empty numDof[] returned from dual space implementation");
61020cf1dd8SToby Isaac   PetscFunctionReturn(0);
61120cf1dd8SToby Isaac }
61220cf1dd8SToby Isaac 
61320cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceCreateSection(PetscDualSpace sp, PetscSection *section)
61420cf1dd8SToby Isaac {
61520cf1dd8SToby Isaac   DM             dm;
61620cf1dd8SToby Isaac   PetscInt       pStart, pEnd, depth, h, offset;
61720cf1dd8SToby Isaac   const PetscInt *numDof;
61820cf1dd8SToby Isaac   PetscErrorCode ierr;
61920cf1dd8SToby Isaac 
62020cf1dd8SToby Isaac   PetscFunctionBegin;
62120cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDM(sp,&dm);CHKERRQ(ierr);
62220cf1dd8SToby Isaac   ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr);
62320cf1dd8SToby Isaac   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)sp),section);CHKERRQ(ierr);
62420cf1dd8SToby Isaac   ierr = PetscSectionSetChart(*section,pStart,pEnd);CHKERRQ(ierr);
62520cf1dd8SToby Isaac   ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr);
62620cf1dd8SToby Isaac   ierr = PetscDualSpaceGetNumDof(sp,&numDof);CHKERRQ(ierr);
62720cf1dd8SToby Isaac   for (h = 0; h <= depth; h++) {
62820cf1dd8SToby Isaac     PetscInt hStart, hEnd, p, dof;
62920cf1dd8SToby Isaac 
63020cf1dd8SToby Isaac     ierr = DMPlexGetHeightStratum(dm,h,&hStart,&hEnd);CHKERRQ(ierr);
63120cf1dd8SToby Isaac     dof = numDof[depth - h];
63220cf1dd8SToby Isaac     for (p = hStart; p < hEnd; p++) {
63320cf1dd8SToby Isaac       ierr = PetscSectionSetDof(*section,p,dof);CHKERRQ(ierr);
63420cf1dd8SToby Isaac     }
63520cf1dd8SToby Isaac   }
63620cf1dd8SToby Isaac   ierr = PetscSectionSetUp(*section);CHKERRQ(ierr);
63720cf1dd8SToby Isaac   for (h = 0, offset = 0; h <= depth; h++) {
63820cf1dd8SToby Isaac     PetscInt hStart, hEnd, p, dof;
63920cf1dd8SToby Isaac 
64020cf1dd8SToby Isaac     ierr = DMPlexGetHeightStratum(dm,h,&hStart,&hEnd);CHKERRQ(ierr);
64120cf1dd8SToby Isaac     dof = numDof[depth - h];
64220cf1dd8SToby Isaac     for (p = hStart; p < hEnd; p++) {
64320cf1dd8SToby Isaac       ierr = PetscSectionGetDof(*section,p,&dof);CHKERRQ(ierr);
64420cf1dd8SToby Isaac       ierr = PetscSectionSetOffset(*section,p,offset);CHKERRQ(ierr);
64520cf1dd8SToby Isaac       offset += dof;
64620cf1dd8SToby Isaac     }
64720cf1dd8SToby Isaac   }
64820cf1dd8SToby Isaac   PetscFunctionReturn(0);
64920cf1dd8SToby Isaac }
65020cf1dd8SToby Isaac 
65120cf1dd8SToby Isaac /*@
65220cf1dd8SToby Isaac   PetscDualSpaceCreateReferenceCell - Create a DMPLEX with the appropriate FEM reference cell
65320cf1dd8SToby Isaac 
65420cf1dd8SToby Isaac   Collective on PetscDualSpace
65520cf1dd8SToby Isaac 
65620cf1dd8SToby Isaac   Input Parameters:
65720cf1dd8SToby Isaac + sp      - The PetscDualSpace
65820cf1dd8SToby Isaac . dim     - The spatial dimension
65920cf1dd8SToby Isaac - simplex - Flag for simplex, otherwise use a tensor-product cell
66020cf1dd8SToby Isaac 
66120cf1dd8SToby Isaac   Output Parameter:
66220cf1dd8SToby Isaac . refdm - The reference cell
66320cf1dd8SToby Isaac 
66420cf1dd8SToby Isaac   Level: advanced
66520cf1dd8SToby Isaac 
66620cf1dd8SToby Isaac .keywords: PetscDualSpace, reference cell
66720cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate(), DMPLEX
66820cf1dd8SToby Isaac @*/
66920cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceCreateReferenceCell(PetscDualSpace sp, PetscInt dim, PetscBool simplex, DM *refdm)
67020cf1dd8SToby Isaac {
67120cf1dd8SToby Isaac   PetscErrorCode ierr;
67220cf1dd8SToby Isaac 
67320cf1dd8SToby Isaac   PetscFunctionBeginUser;
67420cf1dd8SToby Isaac   ierr = DMPlexCreateReferenceCell(PetscObjectComm((PetscObject) sp), dim, simplex, refdm);CHKERRQ(ierr);
67520cf1dd8SToby Isaac   PetscFunctionReturn(0);
67620cf1dd8SToby Isaac }
67720cf1dd8SToby Isaac 
67820cf1dd8SToby Isaac /*@C
67920cf1dd8SToby Isaac   PetscDualSpaceApply - Apply a functional from the dual space basis to an input function
68020cf1dd8SToby Isaac 
68120cf1dd8SToby Isaac   Input Parameters:
68220cf1dd8SToby Isaac + sp      - The PetscDualSpace object
68320cf1dd8SToby Isaac . f       - The basis functional index
68420cf1dd8SToby Isaac . time    - The time
68520cf1dd8SToby 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)
68620cf1dd8SToby Isaac . numComp - The number of components for the function
68720cf1dd8SToby Isaac . func    - The input function
68820cf1dd8SToby Isaac - ctx     - A context for the function
68920cf1dd8SToby Isaac 
69020cf1dd8SToby Isaac   Output Parameter:
69120cf1dd8SToby Isaac . value   - numComp output values
69220cf1dd8SToby Isaac 
69320cf1dd8SToby Isaac   Note: The calling sequence for the callback func is given by:
69420cf1dd8SToby Isaac 
69520cf1dd8SToby Isaac $ func(PetscInt dim, PetscReal time, const PetscReal x[],
69620cf1dd8SToby Isaac $      PetscInt numComponents, PetscScalar values[], void *ctx)
69720cf1dd8SToby Isaac 
69820cf1dd8SToby Isaac   Level: developer
69920cf1dd8SToby Isaac 
70020cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate()
70120cf1dd8SToby Isaac @*/
70220cf1dd8SToby 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)
70320cf1dd8SToby Isaac {
70420cf1dd8SToby Isaac   PetscErrorCode ierr;
70520cf1dd8SToby Isaac 
70620cf1dd8SToby Isaac   PetscFunctionBegin;
70720cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
70820cf1dd8SToby Isaac   PetscValidPointer(cgeom, 4);
70920cf1dd8SToby Isaac   PetscValidPointer(value, 8);
71020cf1dd8SToby Isaac   ierr = (*sp->ops->apply)(sp, f, time, cgeom, numComp, func, ctx, value);CHKERRQ(ierr);
71120cf1dd8SToby Isaac   PetscFunctionReturn(0);
71220cf1dd8SToby Isaac }
71320cf1dd8SToby Isaac 
71420cf1dd8SToby Isaac /*@C
71520cf1dd8SToby Isaac   PetscDualSpaceApplyAll - Apply all functionals from the dual space basis to the result of an evaluation at the points returned by PetscDualSpaceGetAllPoints()
71620cf1dd8SToby Isaac 
71720cf1dd8SToby Isaac   Input Parameters:
71820cf1dd8SToby Isaac + sp        - The PetscDualSpace object
71920cf1dd8SToby Isaac - pointEval - Evaluation at the points returned by PetscDualSpaceGetAllPoints()
72020cf1dd8SToby Isaac 
72120cf1dd8SToby Isaac   Output Parameter:
72220cf1dd8SToby Isaac . spValue   - The values of all dual space functionals
72320cf1dd8SToby Isaac 
72420cf1dd8SToby Isaac   Level: developer
72520cf1dd8SToby Isaac 
72620cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate()
72720cf1dd8SToby Isaac @*/
72820cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceApplyAll(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue)
72920cf1dd8SToby Isaac {
73020cf1dd8SToby Isaac   PetscErrorCode ierr;
73120cf1dd8SToby Isaac 
73220cf1dd8SToby Isaac   PetscFunctionBegin;
73320cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
73420cf1dd8SToby Isaac   ierr = (*sp->ops->applyall)(sp, pointEval, spValue);CHKERRQ(ierr);
73520cf1dd8SToby Isaac   PetscFunctionReturn(0);
73620cf1dd8SToby Isaac }
73720cf1dd8SToby Isaac 
73820cf1dd8SToby Isaac /*@C
73920cf1dd8SToby Isaac   PetscDualSpaceApplyDefault - Apply a functional from the dual space basis to an input function by assuming a point evaluation functional.
74020cf1dd8SToby Isaac 
74120cf1dd8SToby Isaac   Input Parameters:
74220cf1dd8SToby Isaac + sp    - The PetscDualSpace object
74320cf1dd8SToby Isaac . f     - The basis functional index
74420cf1dd8SToby Isaac . time  - The time
74520cf1dd8SToby Isaac . cgeom - A context with geometric information for this cell, we use v0 (the initial vertex) and J (the Jacobian)
74620cf1dd8SToby Isaac . Nc    - The number of components for the function
74720cf1dd8SToby Isaac . func  - The input function
74820cf1dd8SToby Isaac - ctx   - A context for the function
74920cf1dd8SToby Isaac 
75020cf1dd8SToby Isaac   Output Parameter:
75120cf1dd8SToby Isaac . value   - The output value
75220cf1dd8SToby Isaac 
75320cf1dd8SToby Isaac   Note: The calling sequence for the callback func is given by:
75420cf1dd8SToby Isaac 
75520cf1dd8SToby Isaac $ func(PetscInt dim, PetscReal time, const PetscReal x[],
75620cf1dd8SToby Isaac $      PetscInt numComponents, PetscScalar values[], void *ctx)
75720cf1dd8SToby Isaac 
75820cf1dd8SToby Isaac and the idea is to evaluate the functional as an integral
75920cf1dd8SToby Isaac 
76020cf1dd8SToby Isaac $ n(f) = int dx n(x) . f(x)
76120cf1dd8SToby Isaac 
76220cf1dd8SToby Isaac where both n and f have Nc components.
76320cf1dd8SToby Isaac 
76420cf1dd8SToby Isaac   Level: developer
76520cf1dd8SToby Isaac 
76620cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate()
76720cf1dd8SToby Isaac @*/
76820cf1dd8SToby 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)
76920cf1dd8SToby Isaac {
77020cf1dd8SToby Isaac   DM               dm;
77120cf1dd8SToby Isaac   PetscQuadrature  n;
77220cf1dd8SToby Isaac   const PetscReal *points, *weights;
77320cf1dd8SToby Isaac   PetscReal        x[3];
77420cf1dd8SToby Isaac   PetscScalar     *val;
77520cf1dd8SToby Isaac   PetscInt         dim, dE, qNc, c, Nq, q;
77620cf1dd8SToby Isaac   PetscBool        isAffine;
77720cf1dd8SToby Isaac   PetscErrorCode   ierr;
77820cf1dd8SToby Isaac 
77920cf1dd8SToby Isaac   PetscFunctionBegin;
78020cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
78120cf1dd8SToby Isaac   PetscValidPointer(value, 5);
78220cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr);
78320cf1dd8SToby Isaac   ierr = PetscDualSpaceGetFunctional(sp, f, &n);CHKERRQ(ierr);
78420cf1dd8SToby Isaac   ierr = PetscQuadratureGetData(n, &dim, &qNc, &Nq, &points, &weights);CHKERRQ(ierr);
78520cf1dd8SToby 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);
78620cf1dd8SToby Isaac   if (qNc != Nc) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature components %D != function components %D", qNc, Nc);
78720cf1dd8SToby Isaac   ierr = DMGetWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr);
78820cf1dd8SToby Isaac   *value = 0.0;
78920cf1dd8SToby Isaac   isAffine = cgeom->isAffine;
79020cf1dd8SToby Isaac   dE = cgeom->dimEmbed;
79120cf1dd8SToby Isaac   for (q = 0; q < Nq; ++q) {
79220cf1dd8SToby Isaac     if (isAffine) {
79320cf1dd8SToby Isaac       CoordinatesRefToReal(dE, cgeom->dim, cgeom->xi, cgeom->v, cgeom->J, &points[q*dim], x);
79420cf1dd8SToby Isaac       ierr = (*func)(dE, time, x, Nc, val, ctx);CHKERRQ(ierr);
79520cf1dd8SToby Isaac     } else {
79620cf1dd8SToby Isaac       ierr = (*func)(dE, time, &cgeom->v[dE*q], Nc, val, ctx);CHKERRQ(ierr);
79720cf1dd8SToby Isaac     }
79820cf1dd8SToby Isaac     for (c = 0; c < Nc; ++c) {
79920cf1dd8SToby Isaac       *value += val[c]*weights[q*Nc+c];
80020cf1dd8SToby Isaac     }
80120cf1dd8SToby Isaac   }
80220cf1dd8SToby Isaac   ierr = DMRestoreWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr);
80320cf1dd8SToby Isaac   PetscFunctionReturn(0);
80420cf1dd8SToby Isaac }
80520cf1dd8SToby Isaac 
80620cf1dd8SToby Isaac /*@C
80720cf1dd8SToby Isaac   PetscDualSpaceApplyAllDefault - Apply all functionals from the dual space basis to the result of an evaluation at the points returned by PetscDualSpaceGetAllPoints()
80820cf1dd8SToby Isaac 
80920cf1dd8SToby Isaac   Input Parameters:
81020cf1dd8SToby Isaac + sp        - The PetscDualSpace object
81120cf1dd8SToby Isaac - pointEval - Evaluation at the points returned by PetscDualSpaceGetAllPoints()
81220cf1dd8SToby Isaac 
81320cf1dd8SToby Isaac   Output Parameter:
81420cf1dd8SToby Isaac . spValue   - The values of all dual space functionals
81520cf1dd8SToby Isaac 
81620cf1dd8SToby Isaac   Level: developer
81720cf1dd8SToby Isaac 
81820cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate()
81920cf1dd8SToby Isaac @*/
82020cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceApplyAllDefault(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue)
82120cf1dd8SToby Isaac {
82220cf1dd8SToby Isaac   PetscQuadrature  n;
82320cf1dd8SToby Isaac   const PetscReal *points, *weights;
82420cf1dd8SToby Isaac   PetscInt         qNc, c, Nq, q, f, spdim, Nc;
82520cf1dd8SToby Isaac   PetscInt         offset;
82620cf1dd8SToby Isaac   PetscErrorCode   ierr;
82720cf1dd8SToby Isaac 
82820cf1dd8SToby Isaac   PetscFunctionBegin;
82920cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
83020cf1dd8SToby Isaac   PetscValidScalarPointer(pointEval, 2);
83120cf1dd8SToby Isaac   PetscValidScalarPointer(spValue, 5);
83220cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDimension(sp, &spdim);CHKERRQ(ierr);
83320cf1dd8SToby Isaac   ierr = PetscDualSpaceGetNumComponents(sp, &Nc);CHKERRQ(ierr);
83420cf1dd8SToby Isaac   for (f = 0, offset = 0; f < spdim; f++) {
83520cf1dd8SToby Isaac     ierr = PetscDualSpaceGetFunctional(sp, f, &n);CHKERRQ(ierr);
83620cf1dd8SToby Isaac     ierr = PetscQuadratureGetData(n, NULL, &qNc, &Nq, &points, &weights);CHKERRQ(ierr);
83720cf1dd8SToby Isaac     if (qNc != Nc) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature components %D != function components %D", qNc, Nc);
83820cf1dd8SToby Isaac     spValue[f] = 0.0;
83920cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
84020cf1dd8SToby Isaac       for (c = 0; c < Nc; ++c) {
84120cf1dd8SToby Isaac         spValue[f] += pointEval[offset++]*weights[q*Nc+c];
84220cf1dd8SToby Isaac       }
84320cf1dd8SToby Isaac     }
84420cf1dd8SToby Isaac   }
84520cf1dd8SToby Isaac   PetscFunctionReturn(0);
84620cf1dd8SToby Isaac }
84720cf1dd8SToby Isaac 
84820cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetAllPoints(PetscDualSpace sp, PetscQuadrature *allPoints)
84920cf1dd8SToby Isaac {
85020cf1dd8SToby Isaac   PetscErrorCode ierr;
85120cf1dd8SToby Isaac 
85220cf1dd8SToby Isaac   PetscFunctionBegin;
85320cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
85420cf1dd8SToby Isaac   PetscValidPointer(allPoints,2);
85520cf1dd8SToby Isaac   if (!sp->allPoints && sp->ops->createallpoints) {
85620cf1dd8SToby Isaac     ierr = (*sp->ops->createallpoints)(sp,&sp->allPoints);CHKERRQ(ierr);
85720cf1dd8SToby Isaac   }
85820cf1dd8SToby Isaac   *allPoints = sp->allPoints;
85920cf1dd8SToby Isaac   PetscFunctionReturn(0);
86020cf1dd8SToby Isaac }
86120cf1dd8SToby Isaac 
86220cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceCreateAllPointsDefault(PetscDualSpace sp, PetscQuadrature *allPoints)
86320cf1dd8SToby Isaac {
86420cf1dd8SToby Isaac   PetscInt        spdim;
86520cf1dd8SToby Isaac   PetscInt        numPoints, offset;
86620cf1dd8SToby Isaac   PetscReal       *points;
86720cf1dd8SToby Isaac   PetscInt        f, dim;
86820cf1dd8SToby Isaac   PetscQuadrature q;
86920cf1dd8SToby Isaac   PetscErrorCode  ierr;
87020cf1dd8SToby Isaac 
87120cf1dd8SToby Isaac   PetscFunctionBegin;
87220cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDimension(sp,&spdim);CHKERRQ(ierr);
87320cf1dd8SToby Isaac   if (!spdim) {
87420cf1dd8SToby Isaac     ierr = PetscQuadratureCreate(PETSC_COMM_SELF,allPoints);CHKERRQ(ierr);
87520cf1dd8SToby Isaac     ierr = PetscQuadratureSetData(*allPoints,0,0,0,NULL,NULL);CHKERRQ(ierr);
87620cf1dd8SToby Isaac   }
87720cf1dd8SToby Isaac   ierr = PetscDualSpaceGetFunctional(sp,0,&q);CHKERRQ(ierr);
87820cf1dd8SToby Isaac   ierr = PetscQuadratureGetData(q,&dim,NULL,&numPoints,NULL,NULL);CHKERRQ(ierr);
87920cf1dd8SToby Isaac   for (f = 1; f < spdim; f++) {
88020cf1dd8SToby Isaac     PetscInt Np;
88120cf1dd8SToby Isaac 
88220cf1dd8SToby Isaac     ierr = PetscDualSpaceGetFunctional(sp,f,&q);CHKERRQ(ierr);
88320cf1dd8SToby Isaac     ierr = PetscQuadratureGetData(q,NULL,NULL,&Np,NULL,NULL);CHKERRQ(ierr);
88420cf1dd8SToby Isaac     numPoints += Np;
88520cf1dd8SToby Isaac   }
88620cf1dd8SToby Isaac   ierr = PetscMalloc1(dim*numPoints,&points);CHKERRQ(ierr);
88720cf1dd8SToby Isaac   for (f = 0, offset = 0; f < spdim; f++) {
88820cf1dd8SToby Isaac     const PetscReal *p;
88920cf1dd8SToby Isaac     PetscInt        Np, i;
89020cf1dd8SToby Isaac 
89120cf1dd8SToby Isaac     ierr = PetscDualSpaceGetFunctional(sp,f,&q);CHKERRQ(ierr);
89220cf1dd8SToby Isaac     ierr = PetscQuadratureGetData(q,NULL,NULL,&Np,&p,NULL);CHKERRQ(ierr);
89320cf1dd8SToby Isaac     for (i = 0; i < Np * dim; i++) {
89420cf1dd8SToby Isaac       points[offset + i] = p[i];
89520cf1dd8SToby Isaac     }
89620cf1dd8SToby Isaac     offset += Np * dim;
89720cf1dd8SToby Isaac   }
89820cf1dd8SToby Isaac   ierr = PetscQuadratureCreate(PETSC_COMM_SELF,allPoints);CHKERRQ(ierr);
89920cf1dd8SToby Isaac   ierr = PetscQuadratureSetData(*allPoints,dim,0,numPoints,points,NULL);CHKERRQ(ierr);
90020cf1dd8SToby Isaac   PetscFunctionReturn(0);
90120cf1dd8SToby Isaac }
90220cf1dd8SToby Isaac 
90320cf1dd8SToby Isaac /*@C
90420cf1dd8SToby Isaac   PetscDualSpaceApplyFVM - Apply a functional from the dual space basis to an input function by assuming a point evaluation functional at the cell centroid.
90520cf1dd8SToby Isaac 
90620cf1dd8SToby Isaac   Input Parameters:
90720cf1dd8SToby Isaac + sp    - The PetscDualSpace object
90820cf1dd8SToby Isaac . f     - The basis functional index
90920cf1dd8SToby Isaac . time  - The time
91020cf1dd8SToby Isaac . cgeom - A context with geometric information for this cell, we currently just use the centroid
91120cf1dd8SToby Isaac . Nc    - The number of components for the function
91220cf1dd8SToby Isaac . func  - The input function
91320cf1dd8SToby Isaac - ctx   - A context for the function
91420cf1dd8SToby Isaac 
91520cf1dd8SToby Isaac   Output Parameter:
91620cf1dd8SToby Isaac . value - The output value (scalar)
91720cf1dd8SToby Isaac 
91820cf1dd8SToby Isaac   Note: The calling sequence for the callback func is given by:
91920cf1dd8SToby Isaac 
92020cf1dd8SToby Isaac $ func(PetscInt dim, PetscReal time, const PetscReal x[],
92120cf1dd8SToby Isaac $      PetscInt numComponents, PetscScalar values[], void *ctx)
92220cf1dd8SToby Isaac 
92320cf1dd8SToby Isaac and the idea is to evaluate the functional as an integral
92420cf1dd8SToby Isaac 
92520cf1dd8SToby Isaac $ n(f) = int dx n(x) . f(x)
92620cf1dd8SToby Isaac 
92720cf1dd8SToby Isaac where both n and f have Nc components.
92820cf1dd8SToby Isaac 
92920cf1dd8SToby Isaac   Level: developer
93020cf1dd8SToby Isaac 
93120cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate()
93220cf1dd8SToby Isaac @*/
93320cf1dd8SToby 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)
93420cf1dd8SToby Isaac {
93520cf1dd8SToby Isaac   DM               dm;
93620cf1dd8SToby Isaac   PetscQuadrature  n;
93720cf1dd8SToby Isaac   const PetscReal *points, *weights;
93820cf1dd8SToby Isaac   PetscScalar     *val;
93920cf1dd8SToby Isaac   PetscInt         dimEmbed, qNc, c, Nq, q;
94020cf1dd8SToby Isaac   PetscErrorCode   ierr;
94120cf1dd8SToby Isaac 
94220cf1dd8SToby Isaac   PetscFunctionBegin;
94320cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
94420cf1dd8SToby Isaac   PetscValidPointer(value, 5);
94520cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr);
94620cf1dd8SToby Isaac   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
94720cf1dd8SToby Isaac   ierr = PetscDualSpaceGetFunctional(sp, f, &n);CHKERRQ(ierr);
94820cf1dd8SToby Isaac   ierr = PetscQuadratureGetData(n, NULL, &qNc, &Nq, &points, &weights);CHKERRQ(ierr);
94920cf1dd8SToby Isaac   if (qNc != Nc) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature components %D != function components %D", qNc, Nc);
95020cf1dd8SToby Isaac   ierr = DMGetWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr);
95120cf1dd8SToby Isaac   *value = 0.;
95220cf1dd8SToby Isaac   for (q = 0; q < Nq; ++q) {
95320cf1dd8SToby Isaac     ierr = (*func)(dimEmbed, time, cgeom->centroid, Nc, val, ctx);CHKERRQ(ierr);
95420cf1dd8SToby Isaac     for (c = 0; c < Nc; ++c) {
95520cf1dd8SToby Isaac       *value += val[c]*weights[q*Nc+c];
95620cf1dd8SToby Isaac     }
95720cf1dd8SToby Isaac   }
95820cf1dd8SToby Isaac   ierr = DMRestoreWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr);
95920cf1dd8SToby Isaac   PetscFunctionReturn(0);
96020cf1dd8SToby Isaac }
96120cf1dd8SToby Isaac 
96220cf1dd8SToby Isaac /*@
96320cf1dd8SToby Isaac   PetscDualSpaceGetHeightSubspace - Get the subset of the dual space basis that is supported on a mesh point of a
96420cf1dd8SToby Isaac   given height.  This assumes that the reference cell is symmetric over points of this height.
96520cf1dd8SToby Isaac 
96620cf1dd8SToby Isaac   If the dual space is not defined on mesh points of the given height (e.g. if the space is discontinuous and
96720cf1dd8SToby Isaac   pointwise values are not defined on the element boundaries), or if the implementation of PetscDualSpace does not
96820cf1dd8SToby Isaac   support extracting subspaces, then NULL is returned.
96920cf1dd8SToby Isaac 
97020cf1dd8SToby Isaac   This does not increment the reference count on the returned dual space, and the user should not destroy it.
97120cf1dd8SToby Isaac 
97220cf1dd8SToby Isaac   Not collective
97320cf1dd8SToby Isaac 
97420cf1dd8SToby Isaac   Input Parameters:
97520cf1dd8SToby Isaac + sp - the PetscDualSpace object
97620cf1dd8SToby Isaac - height - the height of the mesh point for which the subspace is desired
97720cf1dd8SToby Isaac 
97820cf1dd8SToby Isaac   Output Parameter:
97920cf1dd8SToby Isaac . subsp - the subspace.  Note that the functionals in the subspace are with respect to the intrinsic geometry of the
98020cf1dd8SToby Isaac   point, which will be of lesser dimension if height > 0.
98120cf1dd8SToby Isaac 
98220cf1dd8SToby Isaac   Level: advanced
98320cf1dd8SToby Isaac 
98420cf1dd8SToby Isaac .seealso: PetscSpaceGetHeightSubspace(), PetscDualSpace
98520cf1dd8SToby Isaac @*/
98620cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetHeightSubspace(PetscDualSpace sp, PetscInt height, PetscDualSpace *subsp)
98720cf1dd8SToby Isaac {
98820cf1dd8SToby Isaac   PetscErrorCode ierr;
98920cf1dd8SToby Isaac 
99020cf1dd8SToby Isaac   PetscFunctionBegin;
99120cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
99220cf1dd8SToby Isaac   PetscValidPointer(subsp, 3);
99320cf1dd8SToby Isaac   *subsp = NULL;
99420cf1dd8SToby Isaac   if (sp->ops->getheightsubspace) {
99520cf1dd8SToby Isaac     ierr = (*sp->ops->getheightsubspace)(sp, height, subsp);CHKERRQ(ierr);
99620cf1dd8SToby Isaac   }
99720cf1dd8SToby Isaac   PetscFunctionReturn(0);
99820cf1dd8SToby Isaac }
99920cf1dd8SToby Isaac 
100020cf1dd8SToby Isaac /*@
100120cf1dd8SToby Isaac   PetscDualSpaceGetPointSubspace - Get the subset of the dual space basis that is supported on a particular mesh point.
100220cf1dd8SToby Isaac 
100320cf1dd8SToby Isaac   If the dual space is not defined on the mesh point (e.g. if the space is discontinuous and pointwise values are not
100420cf1dd8SToby Isaac   defined on the element boundaries), or if the implementation of PetscDualSpace does not support extracting
100520cf1dd8SToby Isaac   subspaces, then NULL is returned.
100620cf1dd8SToby Isaac 
100720cf1dd8SToby Isaac   This does not increment the reference count on the returned dual space, and the user should not destroy it.
100820cf1dd8SToby Isaac 
100920cf1dd8SToby Isaac   Not collective
101020cf1dd8SToby Isaac 
101120cf1dd8SToby Isaac   Input Parameters:
101220cf1dd8SToby Isaac + sp - the PetscDualSpace object
101320cf1dd8SToby Isaac - point - the point (in the dual space's DM) for which the subspace is desired
101420cf1dd8SToby Isaac 
101520cf1dd8SToby Isaac   Output Parameters:
101620cf1dd8SToby Isaac   bdsp - the subspace.  Note that the functionals in the subspace are with respect to the intrinsic geometry of the
101720cf1dd8SToby Isaac   point, which will be of lesser dimension if height > 0.
101820cf1dd8SToby Isaac 
101920cf1dd8SToby Isaac   Level: advanced
102020cf1dd8SToby Isaac 
102120cf1dd8SToby Isaac .seealso: PetscDualSpace
102220cf1dd8SToby Isaac @*/
102320cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetPointSubspace(PetscDualSpace sp, PetscInt point, PetscDualSpace *bdsp)
102420cf1dd8SToby Isaac {
102520cf1dd8SToby Isaac   PetscErrorCode ierr;
102620cf1dd8SToby Isaac 
102720cf1dd8SToby Isaac   PetscFunctionBegin;
102820cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
102920cf1dd8SToby Isaac   PetscValidPointer(bdsp,2);
103020cf1dd8SToby Isaac   *bdsp = NULL;
103120cf1dd8SToby Isaac   if (sp->ops->getpointsubspace) {
103220cf1dd8SToby Isaac     ierr = (*sp->ops->getpointsubspace)(sp,point,bdsp);CHKERRQ(ierr);
103320cf1dd8SToby Isaac   } else if (sp->ops->getheightsubspace) {
103420cf1dd8SToby Isaac     DM       dm;
103520cf1dd8SToby Isaac     DMLabel  label;
103620cf1dd8SToby Isaac     PetscInt dim, depth, height;
103720cf1dd8SToby Isaac 
103820cf1dd8SToby Isaac     ierr = PetscDualSpaceGetDM(sp,&dm);CHKERRQ(ierr);
103920cf1dd8SToby Isaac     ierr = DMPlexGetDepth(dm,&dim);CHKERRQ(ierr);
104020cf1dd8SToby Isaac     ierr = DMPlexGetDepthLabel(dm,&label);CHKERRQ(ierr);
104120cf1dd8SToby Isaac     ierr = DMLabelGetValue(label,point,&depth);CHKERRQ(ierr);
104220cf1dd8SToby Isaac     height = dim - depth;
104320cf1dd8SToby Isaac     ierr = (*sp->ops->getheightsubspace)(sp,height,bdsp);CHKERRQ(ierr);
104420cf1dd8SToby Isaac   }
104520cf1dd8SToby Isaac   PetscFunctionReturn(0);
104620cf1dd8SToby Isaac }
104720cf1dd8SToby Isaac 
1048