xref: /petsc/src/dm/dt/dualspace/interface/dualspace.c (revision b44575274871f84dc4018b7bb4f8b6493c17eb57)
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 
116f905325SMatthew G. Knepley /*
126f905325SMatthew G. Knepley   PetscDualSpaceLatticePointLexicographic_Internal - Returns all tuples of size 'len' with nonnegative integers that sum up to at most 'max'.
136f905325SMatthew G. Knepley                                                      Ordering is lexicographic with lowest index as least significant in ordering.
14*b4457527SToby Isaac                                                      e.g. for len == 2 and max == 2, this will return, in order, {0,0}, {1,0}, {2,0}, {0,1}, {1,1}, {0,2}.
156f905325SMatthew G. Knepley 
166f905325SMatthew G. Knepley   Input Parameters:
176f905325SMatthew G. Knepley + len - The length of the tuple
186f905325SMatthew G. Knepley . max - The maximum sum
196f905325SMatthew G. Knepley - tup - A tuple of length len+1: tup[len] > 0 indicates a stopping condition
206f905325SMatthew G. Knepley 
216f905325SMatthew G. Knepley   Output Parameter:
226f905325SMatthew G. Knepley . tup - A tuple of len integers whos sum is at most 'max'
236f905325SMatthew G. Knepley 
246f905325SMatthew G. Knepley   Level: developer
256f905325SMatthew G. Knepley 
266f905325SMatthew G. Knepley .seealso: PetscDualSpaceTensorPointLexicographic_Internal()
276f905325SMatthew G. Knepley */
286f905325SMatthew G. Knepley PetscErrorCode PetscDualSpaceLatticePointLexicographic_Internal(PetscInt len, PetscInt max, PetscInt tup[])
296f905325SMatthew G. Knepley {
306f905325SMatthew G. Knepley   PetscFunctionBegin;
316f905325SMatthew G. Knepley   while (len--) {
326f905325SMatthew G. Knepley     max -= tup[len];
336f905325SMatthew G. Knepley     if (!max) {
346f905325SMatthew G. Knepley       tup[len] = 0;
356f905325SMatthew G. Knepley       break;
366f905325SMatthew G. Knepley     }
376f905325SMatthew G. Knepley   }
386f905325SMatthew G. Knepley   tup[++len]++;
396f905325SMatthew G. Knepley   PetscFunctionReturn(0);
406f905325SMatthew G. Knepley }
416f905325SMatthew G. Knepley 
426f905325SMatthew G. Knepley /*
436f905325SMatthew G. Knepley   PetscDualSpaceTensorPointLexicographic_Internal - Returns all tuples of size 'len' with nonnegative integers that are all less than or equal to 'max'.
446f905325SMatthew G. Knepley                                                     Ordering is lexicographic with lowest index as least significant in ordering.
456f905325SMatthew 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}.
466f905325SMatthew G. Knepley 
476f905325SMatthew G. Knepley   Input Parameters:
486f905325SMatthew G. Knepley + len - The length of the tuple
496f905325SMatthew G. Knepley . max - The maximum value
506f905325SMatthew G. Knepley - tup - A tuple of length len+1: tup[len] > 0 indicates a stopping condition
516f905325SMatthew G. Knepley 
526f905325SMatthew G. Knepley   Output Parameter:
536f905325SMatthew G. Knepley . tup - A tuple of len integers whos sum is at most 'max'
546f905325SMatthew G. Knepley 
556f905325SMatthew G. Knepley   Level: developer
566f905325SMatthew G. Knepley 
576f905325SMatthew G. Knepley .seealso: PetscDualSpaceLatticePointLexicographic_Internal()
586f905325SMatthew G. Knepley */
596f905325SMatthew G. Knepley PetscErrorCode PetscDualSpaceTensorPointLexicographic_Internal(PetscInt len, PetscInt max, PetscInt tup[])
606f905325SMatthew G. Knepley {
616f905325SMatthew G. Knepley   PetscInt       i;
626f905325SMatthew G. Knepley 
636f905325SMatthew G. Knepley   PetscFunctionBegin;
646f905325SMatthew G. Knepley   for (i = 0; i < len; i++) {
656f905325SMatthew G. Knepley     if (tup[i] < max) {
666f905325SMatthew G. Knepley       break;
676f905325SMatthew G. Knepley     } else {
686f905325SMatthew G. Knepley       tup[i] = 0;
696f905325SMatthew G. Knepley     }
706f905325SMatthew G. Knepley   }
716f905325SMatthew G. Knepley   tup[i]++;
726f905325SMatthew G. Knepley   PetscFunctionReturn(0);
736f905325SMatthew G. Knepley }
746f905325SMatthew 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 .seealso: PetscDualSpaceRegisterAll(), PetscDualSpaceRegisterDestroy()
10520cf1dd8SToby Isaac 
10620cf1dd8SToby Isaac @*/
10720cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceRegister(const char sname[], PetscErrorCode (*function)(PetscDualSpace))
10820cf1dd8SToby Isaac {
10920cf1dd8SToby Isaac   PetscErrorCode ierr;
11020cf1dd8SToby Isaac 
11120cf1dd8SToby Isaac   PetscFunctionBegin;
11220cf1dd8SToby Isaac   ierr = PetscFunctionListAdd(&PetscDualSpaceList, sname, function);CHKERRQ(ierr);
11320cf1dd8SToby Isaac   PetscFunctionReturn(0);
11420cf1dd8SToby Isaac }
11520cf1dd8SToby Isaac 
11620cf1dd8SToby Isaac /*@C
11720cf1dd8SToby Isaac   PetscDualSpaceSetType - Builds a particular PetscDualSpace
11820cf1dd8SToby Isaac 
119d083f849SBarry Smith   Collective on sp
12020cf1dd8SToby Isaac 
12120cf1dd8SToby Isaac   Input Parameters:
12220cf1dd8SToby Isaac + sp   - The PetscDualSpace object
12320cf1dd8SToby Isaac - name - The kind of space
12420cf1dd8SToby Isaac 
12520cf1dd8SToby Isaac   Options Database Key:
12620cf1dd8SToby Isaac . -petscdualspace_type <type> - Sets the PetscDualSpace type; use -help for a list of available types
12720cf1dd8SToby Isaac 
12820cf1dd8SToby Isaac   Level: intermediate
12920cf1dd8SToby Isaac 
13020cf1dd8SToby Isaac .seealso: PetscDualSpaceGetType(), PetscDualSpaceCreate()
13120cf1dd8SToby Isaac @*/
13220cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetType(PetscDualSpace sp, PetscDualSpaceType name)
13320cf1dd8SToby Isaac {
13420cf1dd8SToby Isaac   PetscErrorCode (*r)(PetscDualSpace);
13520cf1dd8SToby Isaac   PetscBool      match;
13620cf1dd8SToby Isaac   PetscErrorCode ierr;
13720cf1dd8SToby Isaac 
13820cf1dd8SToby Isaac   PetscFunctionBegin;
13920cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
14020cf1dd8SToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) sp, name, &match);CHKERRQ(ierr);
14120cf1dd8SToby Isaac   if (match) PetscFunctionReturn(0);
14220cf1dd8SToby Isaac 
14320cf1dd8SToby Isaac   if (!PetscDualSpaceRegisterAllCalled) {ierr = PetscDualSpaceRegisterAll();CHKERRQ(ierr);}
14420cf1dd8SToby Isaac   ierr = PetscFunctionListFind(PetscDualSpaceList, name, &r);CHKERRQ(ierr);
14520cf1dd8SToby Isaac   if (!r) SETERRQ1(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscDualSpace type: %s", name);
14620cf1dd8SToby Isaac 
14720cf1dd8SToby Isaac   if (sp->ops->destroy) {
14820cf1dd8SToby Isaac     ierr             = (*sp->ops->destroy)(sp);CHKERRQ(ierr);
14920cf1dd8SToby Isaac     sp->ops->destroy = NULL;
15020cf1dd8SToby Isaac   }
15120cf1dd8SToby Isaac   ierr = (*r)(sp);CHKERRQ(ierr);
15220cf1dd8SToby Isaac   ierr = PetscObjectChangeTypeName((PetscObject) sp, name);CHKERRQ(ierr);
15320cf1dd8SToby Isaac   PetscFunctionReturn(0);
15420cf1dd8SToby Isaac }
15520cf1dd8SToby Isaac 
15620cf1dd8SToby Isaac /*@C
15720cf1dd8SToby Isaac   PetscDualSpaceGetType - Gets the PetscDualSpace type name (as a string) from the object.
15820cf1dd8SToby Isaac 
15920cf1dd8SToby Isaac   Not Collective
16020cf1dd8SToby Isaac 
16120cf1dd8SToby Isaac   Input Parameter:
16220cf1dd8SToby Isaac . sp  - The PetscDualSpace
16320cf1dd8SToby Isaac 
16420cf1dd8SToby Isaac   Output Parameter:
16520cf1dd8SToby Isaac . name - The PetscDualSpace type name
16620cf1dd8SToby Isaac 
16720cf1dd8SToby Isaac   Level: intermediate
16820cf1dd8SToby Isaac 
16920cf1dd8SToby Isaac .seealso: PetscDualSpaceSetType(), PetscDualSpaceCreate()
17020cf1dd8SToby Isaac @*/
17120cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetType(PetscDualSpace sp, PetscDualSpaceType *name)
17220cf1dd8SToby Isaac {
17320cf1dd8SToby Isaac   PetscErrorCode ierr;
17420cf1dd8SToby Isaac 
17520cf1dd8SToby Isaac   PetscFunctionBegin;
17620cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
17720cf1dd8SToby Isaac   PetscValidPointer(name, 2);
17820cf1dd8SToby Isaac   if (!PetscDualSpaceRegisterAllCalled) {
17920cf1dd8SToby Isaac     ierr = PetscDualSpaceRegisterAll();CHKERRQ(ierr);
18020cf1dd8SToby Isaac   }
18120cf1dd8SToby Isaac   *name = ((PetscObject) sp)->type_name;
18220cf1dd8SToby Isaac   PetscFunctionReturn(0);
18320cf1dd8SToby Isaac }
18420cf1dd8SToby Isaac 
185221d6281SMatthew G. Knepley static PetscErrorCode PetscDualSpaceView_ASCII(PetscDualSpace sp, PetscViewer v)
186221d6281SMatthew G. Knepley {
187221d6281SMatthew G. Knepley   PetscViewerFormat format;
188221d6281SMatthew G. Knepley   PetscInt          pdim, f;
189221d6281SMatthew G. Knepley   PetscErrorCode    ierr;
190221d6281SMatthew G. Knepley 
191221d6281SMatthew G. Knepley   PetscFunctionBegin;
192221d6281SMatthew G. Knepley   ierr = PetscDualSpaceGetDimension(sp, &pdim);CHKERRQ(ierr);
193221d6281SMatthew G. Knepley   ierr = PetscObjectPrintClassNamePrefixType((PetscObject) sp, v);CHKERRQ(ierr);
194221d6281SMatthew G. Knepley   ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr);
195*b4457527SToby Isaac   if (sp->k) {
196*b4457527SToby Isaac     ierr = PetscViewerASCIIPrintf(v, "Dual space for %D-forms %swith %D components, size %D\n", PetscAbsInt(sp->k), sp->k < 0 ? "(stored in dual form) ": "", sp->Nc, pdim);CHKERRQ(ierr);
197*b4457527SToby Isaac   } else {
198221d6281SMatthew G. Knepley     ierr = PetscViewerASCIIPrintf(v, "Dual space with %D components, size %D\n", sp->Nc, pdim);CHKERRQ(ierr);
199*b4457527SToby Isaac   }
200221d6281SMatthew G. Knepley   if (sp->ops->view) {ierr = (*sp->ops->view)(sp, v);CHKERRQ(ierr);}
201221d6281SMatthew G. Knepley   ierr = PetscViewerGetFormat(v, &format);CHKERRQ(ierr);
202221d6281SMatthew G. Knepley   if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
203221d6281SMatthew G. Knepley     ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr);
204221d6281SMatthew G. Knepley     for (f = 0; f < pdim; ++f) {
205221d6281SMatthew G. Knepley       ierr = PetscViewerASCIIPrintf(v, "Dual basis vector %D\n", f);CHKERRQ(ierr);
206221d6281SMatthew G. Knepley       ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr);
207221d6281SMatthew G. Knepley       ierr = PetscQuadratureView(sp->functional[f], v);CHKERRQ(ierr);
208221d6281SMatthew G. Knepley       ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr);
209221d6281SMatthew G. Knepley     }
210221d6281SMatthew G. Knepley     ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr);
211221d6281SMatthew G. Knepley   }
212221d6281SMatthew G. Knepley   ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr);
213221d6281SMatthew G. Knepley   PetscFunctionReturn(0);
214221d6281SMatthew G. Knepley }
215221d6281SMatthew G. Knepley 
216fe2efc57SMark /*@C
217fe2efc57SMark    PetscDualSpaceViewFromOptions - View from Options
218fe2efc57SMark 
219fe2efc57SMark    Collective on PetscDualSpace
220fe2efc57SMark 
221fe2efc57SMark    Input Parameters:
222fe2efc57SMark +  A - the PetscDualSpace object
223736c3998SJose E. Roman .  obj - Optional object, proivides prefix
224736c3998SJose E. Roman -  name - command line option
225fe2efc57SMark 
226fe2efc57SMark    Level: intermediate
227fe2efc57SMark .seealso:  PetscDualSpace, PetscDualSpaceView(), PetscObjectViewFromOptions(), PetscDualSpaceCreate()
228fe2efc57SMark @*/
229fe2efc57SMark PetscErrorCode  PetscDualSpaceViewFromOptions(PetscDualSpace A,PetscObject obj,const char name[])
230fe2efc57SMark {
231fe2efc57SMark   PetscErrorCode ierr;
232fe2efc57SMark 
233fe2efc57SMark   PetscFunctionBegin;
234fe2efc57SMark   PetscValidHeaderSpecific(A,PETSCDUALSPACE_CLASSID,1);
235fe2efc57SMark   ierr = PetscObjectViewFromOptions((PetscObject)A,obj,name);CHKERRQ(ierr);
236fe2efc57SMark   PetscFunctionReturn(0);
237fe2efc57SMark }
238fe2efc57SMark 
23920cf1dd8SToby Isaac /*@
24020cf1dd8SToby Isaac   PetscDualSpaceView - Views a PetscDualSpace
24120cf1dd8SToby Isaac 
242d083f849SBarry Smith   Collective on sp
24320cf1dd8SToby Isaac 
24420cf1dd8SToby Isaac   Input Parameter:
24520cf1dd8SToby Isaac + sp - the PetscDualSpace object to view
24620cf1dd8SToby Isaac - v  - the viewer
24720cf1dd8SToby Isaac 
248a4ce7ad1SMatthew G. Knepley   Level: beginner
24920cf1dd8SToby Isaac 
250fe2efc57SMark .seealso PetscDualSpaceDestroy(), PetscDualSpace
25120cf1dd8SToby Isaac @*/
25220cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceView(PetscDualSpace sp, PetscViewer v)
25320cf1dd8SToby Isaac {
254d9bac1caSLisandro Dalcin   PetscBool      iascii;
25520cf1dd8SToby Isaac   PetscErrorCode ierr;
25620cf1dd8SToby Isaac 
25720cf1dd8SToby Isaac   PetscFunctionBegin;
25820cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
259d9bac1caSLisandro Dalcin   if (v) PetscValidHeaderSpecific(v, PETSC_VIEWER_CLASSID, 2);
26020cf1dd8SToby Isaac   if (!v) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) sp), &v);CHKERRQ(ierr);}
261d9bac1caSLisandro Dalcin   ierr = PetscObjectTypeCompare((PetscObject) v, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
262221d6281SMatthew G. Knepley   if (iascii) {ierr = PetscDualSpaceView_ASCII(sp, v);CHKERRQ(ierr);}
26320cf1dd8SToby Isaac   PetscFunctionReturn(0);
26420cf1dd8SToby Isaac }
26520cf1dd8SToby Isaac 
26620cf1dd8SToby Isaac /*@
26720cf1dd8SToby Isaac   PetscDualSpaceSetFromOptions - sets parameters in a PetscDualSpace from the options database
26820cf1dd8SToby Isaac 
269d083f849SBarry Smith   Collective on sp
27020cf1dd8SToby Isaac 
27120cf1dd8SToby Isaac   Input Parameter:
27220cf1dd8SToby Isaac . sp - the PetscDualSpace object to set options for
27320cf1dd8SToby Isaac 
27420cf1dd8SToby Isaac   Options Database:
2757be5e748SToby Isaac . -petscspace_degree the approximation order of the space
27620cf1dd8SToby Isaac 
277a4ce7ad1SMatthew G. Knepley   Level: intermediate
27820cf1dd8SToby Isaac 
279fe2efc57SMark .seealso PetscDualSpaceView(), PetscDualSpace, PetscObjectSetFromOptions()
28020cf1dd8SToby Isaac @*/
28120cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetFromOptions(PetscDualSpace sp)
28220cf1dd8SToby Isaac {
283063ee4adSMatthew G. Knepley   PetscDualSpaceReferenceCell refCell = PETSCDUALSPACE_REFCELL_SIMPLEX;
284063ee4adSMatthew G. Knepley   PetscInt                    refDim  = 0;
285063ee4adSMatthew G. Knepley   PetscBool                   flg;
28620cf1dd8SToby Isaac   const char                 *defaultType;
28720cf1dd8SToby Isaac   char                        name[256];
28820cf1dd8SToby Isaac   PetscErrorCode              ierr;
28920cf1dd8SToby Isaac 
29020cf1dd8SToby Isaac   PetscFunctionBegin;
29120cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
29220cf1dd8SToby Isaac   if (!((PetscObject) sp)->type_name) {
29320cf1dd8SToby Isaac     defaultType = PETSCDUALSPACELAGRANGE;
29420cf1dd8SToby Isaac   } else {
29520cf1dd8SToby Isaac     defaultType = ((PetscObject) sp)->type_name;
29620cf1dd8SToby Isaac   }
29720cf1dd8SToby Isaac   if (!PetscSpaceRegisterAllCalled) {ierr = PetscSpaceRegisterAll();CHKERRQ(ierr);}
29820cf1dd8SToby Isaac 
29920cf1dd8SToby Isaac   ierr = PetscObjectOptionsBegin((PetscObject) sp);CHKERRQ(ierr);
30020cf1dd8SToby Isaac   ierr = PetscOptionsFList("-petscdualspace_type", "Dual space", "PetscDualSpaceSetType", PetscDualSpaceList, defaultType, name, 256, &flg);CHKERRQ(ierr);
30120cf1dd8SToby Isaac   if (flg) {
30220cf1dd8SToby Isaac     ierr = PetscDualSpaceSetType(sp, name);CHKERRQ(ierr);
30320cf1dd8SToby Isaac   } else if (!((PetscObject) sp)->type_name) {
30420cf1dd8SToby Isaac     ierr = PetscDualSpaceSetType(sp, defaultType);CHKERRQ(ierr);
30520cf1dd8SToby Isaac   }
306*b4457527SToby Isaac   ierr = PetscOptionsBoundedInt("-petscdualspace_order", "The approximation order", "PetscDualSpaceSetOrder", sp->order, &sp->order, NULL,0);CHKERRQ(ierr);
307*b4457527SToby Isaac   ierr = PetscOptionsInt("-petscdualspace_form_degree", "The form degree of the dofs", "PetscDualSpaceSetFormDegree", sp->k, &sp->k, NULL);CHKERRQ(ierr);
3085a856986SBarry Smith   ierr = PetscOptionsBoundedInt("-petscdualspace_components", "The number of components", "PetscDualSpaceSetNumComponents", sp->Nc, &sp->Nc, NULL,1);CHKERRQ(ierr);
30920cf1dd8SToby Isaac   if (sp->ops->setfromoptions) {
31020cf1dd8SToby Isaac     ierr = (*sp->ops->setfromoptions)(PetscOptionsObject,sp);CHKERRQ(ierr);
31120cf1dd8SToby Isaac   }
3125a856986SBarry Smith   ierr = PetscOptionsBoundedInt("-petscdualspace_refdim", "The spatial dimension of the reference cell", "PetscDualSpaceSetReferenceCell", refDim, &refDim, NULL,0);CHKERRQ(ierr);
313063ee4adSMatthew G. Knepley   ierr = PetscOptionsEnum("-petscdualspace_refcell", "Reference cell", "PetscDualSpaceSetReferenceCell", PetscDualSpaceReferenceCells, (PetscEnum) refCell, (PetscEnum *) &refCell, &flg);CHKERRQ(ierr);
314063ee4adSMatthew G. Knepley   if (flg) {
315063ee4adSMatthew G. Knepley     DM K;
316063ee4adSMatthew G. Knepley 
317063ee4adSMatthew G. Knepley     if (!refDim) SETERRQ(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_INCOMP, "Reference cell specified without a dimension. Use -petscdualspace_refdim.");
318063ee4adSMatthew G. Knepley     ierr = PetscDualSpaceCreateReferenceCell(sp, refDim, refCell == PETSCDUALSPACE_REFCELL_SIMPLEX ? PETSC_TRUE : PETSC_FALSE, &K);CHKERRQ(ierr);
319063ee4adSMatthew G. Knepley     ierr = PetscDualSpaceSetDM(sp, K);CHKERRQ(ierr);
320063ee4adSMatthew G. Knepley     ierr = DMDestroy(&K);CHKERRQ(ierr);
321063ee4adSMatthew G. Knepley   }
322063ee4adSMatthew G. Knepley 
32320cf1dd8SToby Isaac   /* process any options handlers added with PetscObjectAddOptionsHandler() */
32420cf1dd8SToby Isaac   ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) sp);CHKERRQ(ierr);
32520cf1dd8SToby Isaac   ierr = PetscOptionsEnd();CHKERRQ(ierr);
326063ee4adSMatthew G. Knepley   sp->setfromoptionscalled = PETSC_TRUE;
32720cf1dd8SToby Isaac   PetscFunctionReturn(0);
32820cf1dd8SToby Isaac }
32920cf1dd8SToby Isaac 
33020cf1dd8SToby Isaac /*@
33120cf1dd8SToby Isaac   PetscDualSpaceSetUp - Construct a basis for the PetscDualSpace
33220cf1dd8SToby Isaac 
333d083f849SBarry Smith   Collective on sp
33420cf1dd8SToby Isaac 
33520cf1dd8SToby Isaac   Input Parameter:
33620cf1dd8SToby Isaac . sp - the PetscDualSpace object to setup
33720cf1dd8SToby Isaac 
338a4ce7ad1SMatthew G. Knepley   Level: intermediate
33920cf1dd8SToby Isaac 
340fe2efc57SMark .seealso PetscDualSpaceView(), PetscDualSpaceDestroy(), PetscDualSpace
34120cf1dd8SToby Isaac @*/
34220cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetUp(PetscDualSpace sp)
34320cf1dd8SToby Isaac {
34420cf1dd8SToby Isaac   PetscErrorCode ierr;
34520cf1dd8SToby Isaac 
34620cf1dd8SToby Isaac   PetscFunctionBegin;
34720cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
34820cf1dd8SToby Isaac   if (sp->setupcalled) PetscFunctionReturn(0);
34920cf1dd8SToby Isaac   sp->setupcalled = PETSC_TRUE;
35020cf1dd8SToby Isaac   if (sp->ops->setup) {ierr = (*sp->ops->setup)(sp);CHKERRQ(ierr);}
351063ee4adSMatthew G. Knepley   if (sp->setfromoptionscalled) {ierr = PetscDualSpaceViewFromOptions(sp, NULL, "-petscdualspace_view");CHKERRQ(ierr);}
35220cf1dd8SToby Isaac   PetscFunctionReturn(0);
35320cf1dd8SToby Isaac }
35420cf1dd8SToby Isaac 
355*b4457527SToby Isaac static PetscErrorCode PetscDualSpaceClearDMData_Internal(PetscDualSpace sp, DM dm)
356*b4457527SToby Isaac {
357*b4457527SToby Isaac   PetscInt       pStart = -1, pEnd = -1, depth = -1;
358*b4457527SToby Isaac   PetscErrorCode ierr;
359*b4457527SToby Isaac 
360*b4457527SToby Isaac   PetscFunctionBegin;
361*b4457527SToby Isaac   if (!dm) PetscFunctionReturn(0);
362*b4457527SToby Isaac   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
363*b4457527SToby Isaac   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
364*b4457527SToby Isaac 
365*b4457527SToby Isaac   if (sp->pointSpaces) {
366*b4457527SToby Isaac     PetscInt i;
367*b4457527SToby Isaac 
368*b4457527SToby Isaac     for (i = 0; i < pEnd - pStart; i++) {
369*b4457527SToby Isaac       ierr = PetscDualSpaceDestroy(&(sp->pointSpaces[i]));CHKERRQ(ierr);
370*b4457527SToby Isaac     }
371*b4457527SToby Isaac   }
372*b4457527SToby Isaac   ierr = PetscFree(sp->pointSpaces);CHKERRQ(ierr);
373*b4457527SToby Isaac 
374*b4457527SToby Isaac   if (sp->heightSpaces) {
375*b4457527SToby Isaac     PetscInt i;
376*b4457527SToby Isaac 
377*b4457527SToby Isaac     for (i = 0; i <= depth; i++) {
378*b4457527SToby Isaac       ierr = PetscDualSpaceDestroy(&(sp->heightSpaces[i]));CHKERRQ(ierr);
379*b4457527SToby Isaac     }
380*b4457527SToby Isaac   }
381*b4457527SToby Isaac   ierr = PetscFree(sp->heightSpaces);CHKERRQ(ierr);
382*b4457527SToby Isaac 
383*b4457527SToby Isaac   ierr = PetscSectionDestroy(&(sp->pointSection));CHKERRQ(ierr);
384*b4457527SToby Isaac   ierr = PetscQuadratureDestroy(&(sp->intNodes));CHKERRQ(ierr);
385*b4457527SToby Isaac   ierr = VecDestroy(&(sp->intDofValues));CHKERRQ(ierr);
386*b4457527SToby Isaac   ierr = VecDestroy(&(sp->intNodeValues));CHKERRQ(ierr);
387*b4457527SToby Isaac   ierr = MatDestroy(&(sp->intMat));CHKERRQ(ierr);
388*b4457527SToby Isaac   ierr = PetscQuadratureDestroy(&(sp->allNodes));CHKERRQ(ierr);
389*b4457527SToby Isaac   ierr = VecDestroy(&(sp->allDofValues));CHKERRQ(ierr);
390*b4457527SToby Isaac   ierr = VecDestroy(&(sp->allNodeValues));CHKERRQ(ierr);
391*b4457527SToby Isaac   ierr = MatDestroy(&(sp->allMat));CHKERRQ(ierr);
392*b4457527SToby Isaac   ierr = PetscFree(sp->numDof);CHKERRQ(ierr);
393*b4457527SToby Isaac   PetscFunctionReturn(0);
394*b4457527SToby Isaac }
395*b4457527SToby Isaac 
396*b4457527SToby Isaac 
39720cf1dd8SToby Isaac /*@
39820cf1dd8SToby Isaac   PetscDualSpaceDestroy - Destroys a PetscDualSpace object
39920cf1dd8SToby Isaac 
400d083f849SBarry Smith   Collective on sp
40120cf1dd8SToby Isaac 
40220cf1dd8SToby Isaac   Input Parameter:
40320cf1dd8SToby Isaac . sp - the PetscDualSpace object to destroy
40420cf1dd8SToby Isaac 
405a4ce7ad1SMatthew G. Knepley   Level: beginner
40620cf1dd8SToby Isaac 
407fe2efc57SMark .seealso PetscDualSpaceView(), PetscDualSpace(), PetscDualSpaceCreate()
40820cf1dd8SToby Isaac @*/
40920cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceDestroy(PetscDualSpace *sp)
41020cf1dd8SToby Isaac {
41120cf1dd8SToby Isaac   PetscInt       dim, f;
412*b4457527SToby Isaac   DM             dm;
41320cf1dd8SToby Isaac   PetscErrorCode ierr;
41420cf1dd8SToby Isaac 
41520cf1dd8SToby Isaac   PetscFunctionBegin;
41620cf1dd8SToby Isaac   if (!*sp) PetscFunctionReturn(0);
41720cf1dd8SToby Isaac   PetscValidHeaderSpecific((*sp), PETSCDUALSPACE_CLASSID, 1);
41820cf1dd8SToby Isaac 
41920cf1dd8SToby Isaac   if (--((PetscObject)(*sp))->refct > 0) {*sp = 0; PetscFunctionReturn(0);}
42020cf1dd8SToby Isaac   ((PetscObject) (*sp))->refct = 0;
42120cf1dd8SToby Isaac 
42220cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDimension(*sp, &dim);CHKERRQ(ierr);
423*b4457527SToby Isaac   dm = (*sp)->dm;
424*b4457527SToby Isaac 
425*b4457527SToby Isaac   if ((*sp)->ops->destroy) {ierr = (*(*sp)->ops->destroy)(*sp);CHKERRQ(ierr);}
426*b4457527SToby Isaac   ierr = PetscDualSpaceClearDMData_Internal(*sp, dm);CHKERRQ(ierr);
427*b4457527SToby Isaac 
42820cf1dd8SToby Isaac   for (f = 0; f < dim; ++f) {
42920cf1dd8SToby Isaac     ierr = PetscQuadratureDestroy(&(*sp)->functional[f]);CHKERRQ(ierr);
43020cf1dd8SToby Isaac   }
43120cf1dd8SToby Isaac   ierr = PetscFree((*sp)->functional);CHKERRQ(ierr);
43220cf1dd8SToby Isaac   ierr = DMDestroy(&(*sp)->dm);CHKERRQ(ierr);
43320cf1dd8SToby Isaac   ierr = PetscHeaderDestroy(sp);CHKERRQ(ierr);
43420cf1dd8SToby Isaac   PetscFunctionReturn(0);
43520cf1dd8SToby Isaac }
43620cf1dd8SToby Isaac 
43720cf1dd8SToby Isaac /*@
43820cf1dd8SToby Isaac   PetscDualSpaceCreate - Creates an empty PetscDualSpace object. The type can then be set with PetscDualSpaceSetType().
43920cf1dd8SToby Isaac 
440d083f849SBarry Smith   Collective
44120cf1dd8SToby Isaac 
44220cf1dd8SToby Isaac   Input Parameter:
44320cf1dd8SToby Isaac . comm - The communicator for the PetscDualSpace object
44420cf1dd8SToby Isaac 
44520cf1dd8SToby Isaac   Output Parameter:
44620cf1dd8SToby Isaac . sp - The PetscDualSpace object
44720cf1dd8SToby Isaac 
44820cf1dd8SToby Isaac   Level: beginner
44920cf1dd8SToby Isaac 
45020cf1dd8SToby Isaac .seealso: PetscDualSpaceSetType(), PETSCDUALSPACELAGRANGE
45120cf1dd8SToby Isaac @*/
45220cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceCreate(MPI_Comm comm, PetscDualSpace *sp)
45320cf1dd8SToby Isaac {
45420cf1dd8SToby Isaac   PetscDualSpace s;
45520cf1dd8SToby Isaac   PetscErrorCode ierr;
45620cf1dd8SToby Isaac 
45720cf1dd8SToby Isaac   PetscFunctionBegin;
45820cf1dd8SToby Isaac   PetscValidPointer(sp, 2);
45920cf1dd8SToby Isaac   ierr = PetscCitationsRegister(FECitation,&FEcite);CHKERRQ(ierr);
46020cf1dd8SToby Isaac   *sp  = NULL;
46120cf1dd8SToby Isaac   ierr = PetscFEInitializePackage();CHKERRQ(ierr);
46220cf1dd8SToby Isaac 
46320cf1dd8SToby Isaac   ierr = PetscHeaderCreate(s, PETSCDUALSPACE_CLASSID, "PetscDualSpace", "Dual Space", "PetscDualSpace", comm, PetscDualSpaceDestroy, PetscDualSpaceView);CHKERRQ(ierr);
46420cf1dd8SToby Isaac 
46520cf1dd8SToby Isaac   s->order       = 0;
46620cf1dd8SToby Isaac   s->Nc          = 1;
4674bee2e38SMatthew G. Knepley   s->k           = 0;
468*b4457527SToby Isaac   s->spdim       = -1;
469*b4457527SToby Isaac   s->spintdim    = -1;
470*b4457527SToby Isaac   s->uniform     = PETSC_TRUE;
47120cf1dd8SToby Isaac   s->setupcalled = PETSC_FALSE;
47220cf1dd8SToby Isaac 
47320cf1dd8SToby Isaac   *sp = s;
47420cf1dd8SToby Isaac   PetscFunctionReturn(0);
47520cf1dd8SToby Isaac }
47620cf1dd8SToby Isaac 
47720cf1dd8SToby Isaac /*@
47820cf1dd8SToby Isaac   PetscDualSpaceDuplicate - Creates a duplicate PetscDualSpace object, however it is not setup.
47920cf1dd8SToby Isaac 
480d083f849SBarry Smith   Collective on sp
48120cf1dd8SToby Isaac 
48220cf1dd8SToby Isaac   Input Parameter:
48320cf1dd8SToby Isaac . sp - The original PetscDualSpace
48420cf1dd8SToby Isaac 
48520cf1dd8SToby Isaac   Output Parameter:
48620cf1dd8SToby Isaac . spNew - The duplicate PetscDualSpace
48720cf1dd8SToby Isaac 
48820cf1dd8SToby Isaac   Level: beginner
48920cf1dd8SToby Isaac 
49020cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate(), PetscDualSpaceSetType()
49120cf1dd8SToby Isaac @*/
49220cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceDuplicate(PetscDualSpace sp, PetscDualSpace *spNew)
49320cf1dd8SToby Isaac {
494*b4457527SToby Isaac   DM             dm;
495*b4457527SToby Isaac   PetscDualSpaceType type;
496*b4457527SToby Isaac   const char     *name;
49720cf1dd8SToby Isaac   PetscErrorCode ierr;
49820cf1dd8SToby Isaac 
49920cf1dd8SToby Isaac   PetscFunctionBegin;
50020cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
50120cf1dd8SToby Isaac   PetscValidPointer(spNew, 2);
502*b4457527SToby Isaac   ierr = PetscDualSpaceCreate(PetscObjectComm((PetscObject)sp), spNew);CHKERRQ(ierr);
503*b4457527SToby Isaac   ierr = PetscObjectGetName((PetscObject) sp,     &name);CHKERRQ(ierr);
504*b4457527SToby Isaac   ierr = PetscObjectSetName((PetscObject) *spNew,  name);CHKERRQ(ierr);
505*b4457527SToby Isaac   ierr = PetscDualSpaceGetType(sp, &type);CHKERRQ(ierr);
506*b4457527SToby Isaac   ierr = PetscDualSpaceSetType(*spNew, type);CHKERRQ(ierr);
507*b4457527SToby Isaac   ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr);
508*b4457527SToby Isaac   ierr = PetscDualSpaceSetDM(*spNew, dm);CHKERRQ(ierr);
509*b4457527SToby Isaac 
510*b4457527SToby Isaac   (*spNew)->order   = sp->order;
511*b4457527SToby Isaac   (*spNew)->k       = sp->k;
512*b4457527SToby Isaac   (*spNew)->Nc      = sp->Nc;
513*b4457527SToby Isaac   (*spNew)->uniform = sp->uniform;
514*b4457527SToby Isaac   if (sp->ops->duplicate) {ierr = (*sp->ops->duplicate)(sp, *spNew);CHKERRQ(ierr);}
51520cf1dd8SToby Isaac   PetscFunctionReturn(0);
51620cf1dd8SToby Isaac }
51720cf1dd8SToby Isaac 
51820cf1dd8SToby Isaac /*@
51920cf1dd8SToby Isaac   PetscDualSpaceGetDM - Get the DM representing the reference cell
52020cf1dd8SToby Isaac 
52120cf1dd8SToby Isaac   Not collective
52220cf1dd8SToby Isaac 
52320cf1dd8SToby Isaac   Input Parameter:
52420cf1dd8SToby Isaac . sp - The PetscDualSpace
52520cf1dd8SToby Isaac 
52620cf1dd8SToby Isaac   Output Parameter:
52720cf1dd8SToby Isaac . dm - The reference cell
52820cf1dd8SToby Isaac 
52920cf1dd8SToby Isaac   Level: intermediate
53020cf1dd8SToby Isaac 
53120cf1dd8SToby Isaac .seealso: PetscDualSpaceSetDM(), PetscDualSpaceCreate()
53220cf1dd8SToby Isaac @*/
53320cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetDM(PetscDualSpace sp, DM *dm)
53420cf1dd8SToby Isaac {
53520cf1dd8SToby Isaac   PetscFunctionBegin;
53620cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
53720cf1dd8SToby Isaac   PetscValidPointer(dm, 2);
53820cf1dd8SToby Isaac   *dm = sp->dm;
53920cf1dd8SToby Isaac   PetscFunctionReturn(0);
54020cf1dd8SToby Isaac }
54120cf1dd8SToby Isaac 
54220cf1dd8SToby Isaac /*@
54320cf1dd8SToby Isaac   PetscDualSpaceSetDM - Get the DM representing the reference cell
54420cf1dd8SToby Isaac 
54520cf1dd8SToby Isaac   Not collective
54620cf1dd8SToby Isaac 
54720cf1dd8SToby Isaac   Input Parameters:
54820cf1dd8SToby Isaac + sp - The PetscDualSpace
54920cf1dd8SToby Isaac - dm - The reference cell
55020cf1dd8SToby Isaac 
55120cf1dd8SToby Isaac   Level: intermediate
55220cf1dd8SToby Isaac 
55320cf1dd8SToby Isaac .seealso: PetscDualSpaceGetDM(), PetscDualSpaceCreate()
55420cf1dd8SToby Isaac @*/
55520cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetDM(PetscDualSpace sp, DM dm)
55620cf1dd8SToby Isaac {
55720cf1dd8SToby Isaac   PetscErrorCode ierr;
55820cf1dd8SToby Isaac 
55920cf1dd8SToby Isaac   PetscFunctionBegin;
56020cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
56120cf1dd8SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
562*b4457527SToby Isaac   if (sp->setupcalled) SETERRQ(PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change DM after dualspace is set up");
56320cf1dd8SToby Isaac   ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);
564*b4457527SToby Isaac   if (sp->dm && sp->dm != dm) {
565*b4457527SToby Isaac     ierr = PetscDualSpaceClearDMData_Internal(sp, sp->dm);CHKERRQ(ierr);
566*b4457527SToby Isaac   }
567*b4457527SToby Isaac   ierr = DMDestroy(&sp->dm);CHKERRQ(ierr);
56820cf1dd8SToby Isaac   sp->dm = dm;
56920cf1dd8SToby Isaac   PetscFunctionReturn(0);
57020cf1dd8SToby Isaac }
57120cf1dd8SToby Isaac 
57220cf1dd8SToby Isaac /*@
57320cf1dd8SToby Isaac   PetscDualSpaceGetOrder - Get the order of the dual space
57420cf1dd8SToby Isaac 
57520cf1dd8SToby Isaac   Not collective
57620cf1dd8SToby Isaac 
57720cf1dd8SToby Isaac   Input Parameter:
57820cf1dd8SToby Isaac . sp - The PetscDualSpace
57920cf1dd8SToby Isaac 
58020cf1dd8SToby Isaac   Output Parameter:
58120cf1dd8SToby Isaac . order - The order
58220cf1dd8SToby Isaac 
58320cf1dd8SToby Isaac   Level: intermediate
58420cf1dd8SToby Isaac 
58520cf1dd8SToby Isaac .seealso: PetscDualSpaceSetOrder(), PetscDualSpaceCreate()
58620cf1dd8SToby Isaac @*/
58720cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetOrder(PetscDualSpace sp, PetscInt *order)
58820cf1dd8SToby Isaac {
58920cf1dd8SToby Isaac   PetscFunctionBegin;
59020cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
59120cf1dd8SToby Isaac   PetscValidPointer(order, 2);
59220cf1dd8SToby Isaac   *order = sp->order;
59320cf1dd8SToby Isaac   PetscFunctionReturn(0);
59420cf1dd8SToby Isaac }
59520cf1dd8SToby Isaac 
59620cf1dd8SToby Isaac /*@
59720cf1dd8SToby Isaac   PetscDualSpaceSetOrder - Set the order of the dual space
59820cf1dd8SToby Isaac 
59920cf1dd8SToby Isaac   Not collective
60020cf1dd8SToby Isaac 
60120cf1dd8SToby Isaac   Input Parameters:
60220cf1dd8SToby Isaac + sp - The PetscDualSpace
60320cf1dd8SToby Isaac - order - The order
60420cf1dd8SToby Isaac 
60520cf1dd8SToby Isaac   Level: intermediate
60620cf1dd8SToby Isaac 
60720cf1dd8SToby Isaac .seealso: PetscDualSpaceGetOrder(), PetscDualSpaceCreate()
60820cf1dd8SToby Isaac @*/
60920cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetOrder(PetscDualSpace sp, PetscInt order)
61020cf1dd8SToby Isaac {
61120cf1dd8SToby Isaac   PetscFunctionBegin;
61220cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
613*b4457527SToby Isaac   if (sp->setupcalled) SETERRQ(PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change order after dualspace is set up");
61420cf1dd8SToby Isaac   sp->order = order;
61520cf1dd8SToby Isaac   PetscFunctionReturn(0);
61620cf1dd8SToby Isaac }
61720cf1dd8SToby Isaac 
61820cf1dd8SToby Isaac /*@
61920cf1dd8SToby Isaac   PetscDualSpaceGetNumComponents - Return the number of components for this space
62020cf1dd8SToby Isaac 
62120cf1dd8SToby Isaac   Input Parameter:
62220cf1dd8SToby Isaac . sp - The PetscDualSpace
62320cf1dd8SToby Isaac 
62420cf1dd8SToby Isaac   Output Parameter:
62520cf1dd8SToby Isaac . Nc - The number of components
62620cf1dd8SToby Isaac 
62720cf1dd8SToby Isaac   Note: A vector space, for example, will have d components, where d is the spatial dimension
62820cf1dd8SToby Isaac 
62920cf1dd8SToby Isaac   Level: intermediate
63020cf1dd8SToby Isaac 
63120cf1dd8SToby Isaac .seealso: PetscDualSpaceSetNumComponents(), PetscDualSpaceGetDimension(), PetscDualSpaceCreate(), PetscDualSpace
63220cf1dd8SToby Isaac @*/
63320cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetNumComponents(PetscDualSpace sp, PetscInt *Nc)
63420cf1dd8SToby Isaac {
63520cf1dd8SToby Isaac   PetscFunctionBegin;
63620cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
63720cf1dd8SToby Isaac   PetscValidPointer(Nc, 2);
63820cf1dd8SToby Isaac   *Nc = sp->Nc;
63920cf1dd8SToby Isaac   PetscFunctionReturn(0);
64020cf1dd8SToby Isaac }
64120cf1dd8SToby Isaac 
64220cf1dd8SToby Isaac /*@
64320cf1dd8SToby Isaac   PetscDualSpaceSetNumComponents - Set the number of components for this space
64420cf1dd8SToby Isaac 
64520cf1dd8SToby Isaac   Input Parameters:
64620cf1dd8SToby Isaac + sp - The PetscDualSpace
64720cf1dd8SToby Isaac - order - The number of components
64820cf1dd8SToby Isaac 
64920cf1dd8SToby Isaac   Level: intermediate
65020cf1dd8SToby Isaac 
65120cf1dd8SToby Isaac .seealso: PetscDualSpaceGetNumComponents(), PetscDualSpaceCreate(), PetscDualSpace
65220cf1dd8SToby Isaac @*/
65320cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetNumComponents(PetscDualSpace sp, PetscInt Nc)
65420cf1dd8SToby Isaac {
65520cf1dd8SToby Isaac   PetscFunctionBegin;
65620cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
657*b4457527SToby Isaac   if (sp->setupcalled) SETERRQ(PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change number of components after dualspace is set up");
65820cf1dd8SToby Isaac   sp->Nc = Nc;
65920cf1dd8SToby Isaac   PetscFunctionReturn(0);
66020cf1dd8SToby Isaac }
66120cf1dd8SToby Isaac 
66220cf1dd8SToby Isaac /*@
66320cf1dd8SToby Isaac   PetscDualSpaceGetFunctional - Get the i-th basis functional in the dual space
66420cf1dd8SToby Isaac 
66520cf1dd8SToby Isaac   Not collective
66620cf1dd8SToby Isaac 
66720cf1dd8SToby Isaac   Input Parameters:
66820cf1dd8SToby Isaac + sp - The PetscDualSpace
66920cf1dd8SToby Isaac - i  - The basis number
67020cf1dd8SToby Isaac 
67120cf1dd8SToby Isaac   Output Parameter:
67220cf1dd8SToby Isaac . functional - The basis functional
67320cf1dd8SToby Isaac 
67420cf1dd8SToby Isaac   Level: intermediate
67520cf1dd8SToby Isaac 
67620cf1dd8SToby Isaac .seealso: PetscDualSpaceGetDimension(), PetscDualSpaceCreate()
67720cf1dd8SToby Isaac @*/
67820cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetFunctional(PetscDualSpace sp, PetscInt i, PetscQuadrature *functional)
67920cf1dd8SToby Isaac {
68020cf1dd8SToby Isaac   PetscInt       dim;
68120cf1dd8SToby Isaac   PetscErrorCode ierr;
68220cf1dd8SToby Isaac 
68320cf1dd8SToby Isaac   PetscFunctionBegin;
68420cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
68520cf1dd8SToby Isaac   PetscValidPointer(functional, 3);
68620cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDimension(sp, &dim);CHKERRQ(ierr);
68720cf1dd8SToby Isaac   if ((i < 0) || (i >= dim)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Functional index %d must be in [0, %d)", i, dim);
68820cf1dd8SToby Isaac   *functional = sp->functional[i];
68920cf1dd8SToby Isaac   PetscFunctionReturn(0);
69020cf1dd8SToby Isaac }
69120cf1dd8SToby Isaac 
69220cf1dd8SToby Isaac /*@
69320cf1dd8SToby Isaac   PetscDualSpaceGetDimension - Get the dimension of the dual space, i.e. the number of basis functionals
69420cf1dd8SToby Isaac 
69520cf1dd8SToby Isaac   Not collective
69620cf1dd8SToby Isaac 
69720cf1dd8SToby Isaac   Input Parameter:
69820cf1dd8SToby Isaac . sp - The PetscDualSpace
69920cf1dd8SToby Isaac 
70020cf1dd8SToby Isaac   Output Parameter:
70120cf1dd8SToby Isaac . dim - The dimension
70220cf1dd8SToby Isaac 
70320cf1dd8SToby Isaac   Level: intermediate
70420cf1dd8SToby Isaac 
70520cf1dd8SToby Isaac .seealso: PetscDualSpaceGetFunctional(), PetscDualSpaceCreate()
70620cf1dd8SToby Isaac @*/
70720cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetDimension(PetscDualSpace sp, PetscInt *dim)
70820cf1dd8SToby Isaac {
70920cf1dd8SToby Isaac   PetscErrorCode ierr;
71020cf1dd8SToby Isaac 
71120cf1dd8SToby Isaac   PetscFunctionBegin;
71220cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
71320cf1dd8SToby Isaac   PetscValidPointer(dim, 2);
714*b4457527SToby Isaac   if (sp->spdim < 0) {
715*b4457527SToby Isaac     PetscSection section;
716*b4457527SToby Isaac 
717*b4457527SToby Isaac     ierr = PetscDualSpaceGetSection(sp, &section);CHKERRQ(ierr);
718*b4457527SToby Isaac     if (section) {
719*b4457527SToby Isaac       ierr = PetscSectionGetStorageSize(section, &(sp->spdim));CHKERRQ(ierr);
720*b4457527SToby Isaac     } else sp->spdim = 0;
721*b4457527SToby Isaac   }
722*b4457527SToby Isaac   *dim = sp->spdim;
72320cf1dd8SToby Isaac   PetscFunctionReturn(0);
72420cf1dd8SToby Isaac }
72520cf1dd8SToby Isaac 
726*b4457527SToby Isaac /*@
727*b4457527SToby Isaac   PetscDualSpaceGetInteriorDimension - Get the interior dimension of the dual space, i.e. the number of basis functionals assigned to the interior of the reference domain
728*b4457527SToby Isaac 
729*b4457527SToby Isaac   Not collective
730*b4457527SToby Isaac 
731*b4457527SToby Isaac   Input Parameter:
732*b4457527SToby Isaac . sp - The PetscDualSpace
733*b4457527SToby Isaac 
734*b4457527SToby Isaac   Output Parameter:
735*b4457527SToby Isaac . dim - The dimension
736*b4457527SToby Isaac 
737*b4457527SToby Isaac   Level: intermediate
738*b4457527SToby Isaac 
739*b4457527SToby Isaac .seealso: PetscDualSpaceGetFunctional(), PetscDualSpaceCreate()
740*b4457527SToby Isaac @*/
741*b4457527SToby Isaac PetscErrorCode PetscDualSpaceGetInteriorDimension(PetscDualSpace sp, PetscInt *intdim)
742*b4457527SToby Isaac {
743*b4457527SToby Isaac   PetscErrorCode ierr;
744*b4457527SToby Isaac 
745*b4457527SToby Isaac   PetscFunctionBegin;
746*b4457527SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
747*b4457527SToby Isaac   PetscValidPointer(intdim, 2);
748*b4457527SToby Isaac   if (sp->spintdim < 0) {
749*b4457527SToby Isaac     PetscSection section;
750*b4457527SToby Isaac 
751*b4457527SToby Isaac     ierr = PetscDualSpaceGetSection(sp, &section);CHKERRQ(ierr);
752*b4457527SToby Isaac     if (section) {
753*b4457527SToby Isaac       ierr = PetscSectionGetConstrainedStorageSize(section, &(sp->spintdim));CHKERRQ(ierr);
754*b4457527SToby Isaac     } else sp->spintdim = 0;
755*b4457527SToby Isaac   }
756*b4457527SToby Isaac   *intdim = sp->spintdim;
757*b4457527SToby Isaac   PetscFunctionReturn(0);
758*b4457527SToby Isaac }
759*b4457527SToby Isaac 
760*b4457527SToby Isaac /*@
761*b4457527SToby Isaac    PetscDualSpaceGetUniform - Whether this dual space is uniform
762*b4457527SToby Isaac 
763*b4457527SToby Isaac    Not collective
764*b4457527SToby Isaac 
765*b4457527SToby Isaac    Input Parameters:
766*b4457527SToby Isaac .  sp - A dual space
767*b4457527SToby Isaac 
768*b4457527SToby Isaac    Output Parameters:
769*b4457527SToby Isaac .  uniform - PETSC_TRUE if (a) the dual space is the same for each point in a stratum of the reference DMPlex, and
770*b4457527SToby Isaac              (b) every symmetry of each point in the reference DMPlex is also a symmetry of the point's dual space.
771*b4457527SToby Isaac 
772*b4457527SToby Isaac 
773*b4457527SToby Isaac    Level: advanced
774*b4457527SToby Isaac 
775*b4457527SToby Isaac    Note: all of the usual spaces on simplex or tensor-product elements will be uniform, only reference cells
776*b4457527SToby Isaac    with non-uniform strata (like trianguar-prisms) or anisotropic hp dual spaces will not be uniform.
777*b4457527SToby Isaac 
778*b4457527SToby Isaac .seealso: PetscDualSpaceGetPointSubspace(), PetscDualSpaceGetSymmetries()
779*b4457527SToby Isaac @*/
780*b4457527SToby Isaac PetscErrorCode PetscDualSpaceGetUniform(PetscDualSpace sp, PetscBool *uniform)
781*b4457527SToby Isaac {
782*b4457527SToby Isaac   PetscFunctionBegin;
783*b4457527SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
784*b4457527SToby Isaac   PetscValidPointer(uniform, 2);
785*b4457527SToby Isaac   *uniform = sp->uniform;
786*b4457527SToby Isaac   PetscFunctionReturn(0);
787*b4457527SToby Isaac }
788*b4457527SToby Isaac 
789*b4457527SToby Isaac 
79020cf1dd8SToby Isaac /*@C
79120cf1dd8SToby Isaac   PetscDualSpaceGetNumDof - Get the number of degrees of freedom for each spatial (topological) dimension
79220cf1dd8SToby Isaac 
79320cf1dd8SToby Isaac   Not collective
79420cf1dd8SToby Isaac 
79520cf1dd8SToby Isaac   Input Parameter:
79620cf1dd8SToby Isaac . sp - The PetscDualSpace
79720cf1dd8SToby Isaac 
79820cf1dd8SToby Isaac   Output Parameter:
79920cf1dd8SToby Isaac . numDof - An array of length dim+1 which holds the number of dofs for each dimension
80020cf1dd8SToby Isaac 
80120cf1dd8SToby Isaac   Level: intermediate
80220cf1dd8SToby Isaac 
80320cf1dd8SToby Isaac .seealso: PetscDualSpaceGetFunctional(), PetscDualSpaceCreate()
80420cf1dd8SToby Isaac @*/
80520cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetNumDof(PetscDualSpace sp, const PetscInt **numDof)
80620cf1dd8SToby Isaac {
80720cf1dd8SToby Isaac   PetscErrorCode ierr;
80820cf1dd8SToby Isaac 
80920cf1dd8SToby Isaac   PetscFunctionBegin;
81020cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
81120cf1dd8SToby Isaac   PetscValidPointer(numDof, 2);
812*b4457527SToby Isaac   if (!sp->uniform) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "A non-uniform space does not have a fixed number of dofs for each height");
813*b4457527SToby Isaac   if (!sp->numDof) {
814*b4457527SToby Isaac     DM       dm;
815*b4457527SToby Isaac     PetscInt depth, d;
816*b4457527SToby Isaac     PetscSection section;
817*b4457527SToby Isaac 
818*b4457527SToby Isaac     ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr);
819*b4457527SToby Isaac     ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
820*b4457527SToby Isaac     ierr = PetscCalloc1(depth+1,&(sp->numDof));CHKERRQ(ierr);
821*b4457527SToby Isaac     ierr = PetscDualSpaceGetSection(sp, &section);CHKERRQ(ierr);
822*b4457527SToby Isaac     for (d = 0; d <= depth; d++) {
823*b4457527SToby Isaac       PetscInt dStart, dEnd;
824*b4457527SToby Isaac 
825*b4457527SToby Isaac       ierr = DMPlexGetDepthStratum(dm, d, &dStart, &dEnd);CHKERRQ(ierr);
826*b4457527SToby Isaac       if (dEnd <= dStart) continue;
827*b4457527SToby Isaac       ierr = PetscSectionGetDof(section, dStart, &(sp->numDof[d]));CHKERRQ(ierr);
828*b4457527SToby Isaac 
829*b4457527SToby Isaac     }
830*b4457527SToby Isaac   }
831*b4457527SToby Isaac   *numDof = sp->numDof;
83220cf1dd8SToby Isaac   if (!*numDof) SETERRQ(PetscObjectComm((PetscObject) sp), PETSC_ERR_LIB, "Empty numDof[] returned from dual space implementation");
83320cf1dd8SToby Isaac   PetscFunctionReturn(0);
83420cf1dd8SToby Isaac }
83520cf1dd8SToby Isaac 
836*b4457527SToby Isaac /* create the section of the right size and set a permutation for topological ordering */
837*b4457527SToby Isaac PetscErrorCode PetscDualSpaceSectionCreate_Internal(PetscDualSpace sp, PetscSection *topSection)
838*b4457527SToby Isaac {
839*b4457527SToby Isaac   DM             dm;
840*b4457527SToby Isaac   PetscInt       pStart, pEnd, cStart, cEnd, c, depth, count, i;
841*b4457527SToby Isaac   PetscInt       *seen, *perm;
842*b4457527SToby Isaac   PetscSection   section;
843*b4457527SToby Isaac   PetscErrorCode ierr;
844*b4457527SToby Isaac 
845*b4457527SToby Isaac   PetscFunctionBegin;
846*b4457527SToby Isaac   dm = sp->dm;
847*b4457527SToby Isaac   ierr = PetscSectionCreate(PETSC_COMM_SELF, &section);CHKERRQ(ierr);
848*b4457527SToby Isaac   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
849*b4457527SToby Isaac   ierr = PetscSectionSetChart(section, pStart, pEnd);CHKERRQ(ierr);
850*b4457527SToby Isaac   ierr = PetscCalloc1(pEnd - pStart, &seen);CHKERRQ(ierr);
851*b4457527SToby Isaac   ierr = PetscMalloc1(pEnd - pStart, &perm);CHKERRQ(ierr);
852*b4457527SToby Isaac   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
853*b4457527SToby Isaac   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
854*b4457527SToby Isaac   for (c = cStart, count = 0; c < cEnd; c++) {
855*b4457527SToby Isaac     PetscInt closureSize = -1, e;
856*b4457527SToby Isaac     PetscInt *closure = NULL;
857*b4457527SToby Isaac 
858*b4457527SToby Isaac     perm[count++] = c;
859*b4457527SToby Isaac     seen[c-pStart] = 1;
860*b4457527SToby Isaac     ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
861*b4457527SToby Isaac     for (e = 0; e < closureSize; e++) {
862*b4457527SToby Isaac       PetscInt point = closure[2*e];
863*b4457527SToby Isaac 
864*b4457527SToby Isaac       if (seen[point-pStart]) continue;
865*b4457527SToby Isaac       perm[count++] = point;
866*b4457527SToby Isaac       seen[point-pStart] = 1;
867*b4457527SToby Isaac     }
868*b4457527SToby Isaac     ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
869*b4457527SToby Isaac   }
870*b4457527SToby Isaac   if (count != pEnd - pStart) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Bad topological ordering");
871*b4457527SToby Isaac   for (i = 0; i < pEnd - pStart; i++) if (perm[i] != i) break;
872*b4457527SToby Isaac   if (i < pEnd - pStart) {
873*b4457527SToby Isaac     IS permIS;
874*b4457527SToby Isaac 
875*b4457527SToby Isaac     ierr = ISCreateGeneral(PETSC_COMM_SELF, pEnd - pStart, perm, PETSC_OWN_POINTER, &permIS);CHKERRQ(ierr);
876*b4457527SToby Isaac     ierr = ISSetPermutation(permIS);CHKERRQ(ierr);
877*b4457527SToby Isaac     ierr = PetscSectionSetPermutation(section, permIS);CHKERRQ(ierr);
878*b4457527SToby Isaac     ierr = ISDestroy(&permIS);CHKERRQ(ierr);
879*b4457527SToby Isaac   } else {
880*b4457527SToby Isaac     ierr = PetscFree(perm);CHKERRQ(ierr);
881*b4457527SToby Isaac   }
882*b4457527SToby Isaac   ierr = PetscFree(seen);CHKERRQ(ierr);
883*b4457527SToby Isaac   *topSection = section;
884*b4457527SToby Isaac   PetscFunctionReturn(0);
885*b4457527SToby Isaac }
886*b4457527SToby Isaac 
887*b4457527SToby Isaac /* mark boundary points and set up */
888*b4457527SToby Isaac PetscErrorCode PetscDualSpaceSectionSetUp_Internal(PetscDualSpace sp, PetscSection section)
889*b4457527SToby Isaac {
890*b4457527SToby Isaac   DM             dm;
891*b4457527SToby Isaac   DMLabel        boundary;
892*b4457527SToby Isaac   PetscInt       pStart, pEnd, p;
893*b4457527SToby Isaac   PetscErrorCode ierr;
894*b4457527SToby Isaac 
895*b4457527SToby Isaac   PetscFunctionBegin;
896*b4457527SToby Isaac   dm = sp->dm;
897*b4457527SToby Isaac   ierr = DMLabelCreate(PETSC_COMM_SELF,"boundary",&boundary);CHKERRQ(ierr);
898*b4457527SToby Isaac   ierr = PetscDualSpaceGetDM(sp,&dm);CHKERRQ(ierr);
899*b4457527SToby Isaac   ierr = DMPlexMarkBoundaryFaces(dm,1,boundary);CHKERRQ(ierr);
900*b4457527SToby Isaac   ierr = DMPlexLabelComplete(dm,boundary);CHKERRQ(ierr);
901*b4457527SToby Isaac   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
902*b4457527SToby Isaac   for (p = pStart; p < pEnd; p++) {
903*b4457527SToby Isaac     PetscInt bval;
904*b4457527SToby Isaac 
905*b4457527SToby Isaac     ierr = DMLabelGetValue(boundary, p, &bval);CHKERRQ(ierr);
906*b4457527SToby Isaac     if (bval == 1) {
907*b4457527SToby Isaac       PetscInt dof;
908*b4457527SToby Isaac 
909*b4457527SToby Isaac       ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
910*b4457527SToby Isaac       ierr = PetscSectionSetConstraintDof(section, p, dof);CHKERRQ(ierr);
911*b4457527SToby Isaac     }
912*b4457527SToby Isaac   }
913*b4457527SToby Isaac   ierr = DMLabelDestroy(&boundary);CHKERRQ(ierr);
914*b4457527SToby Isaac   ierr = PetscSectionSetUp(section);
915*b4457527SToby Isaac   PetscFunctionReturn(0);
916*b4457527SToby Isaac }
917*b4457527SToby Isaac 
918a4ce7ad1SMatthew G. Knepley /*@
919*b4457527SToby Isaac   PetscDualSpaceGetSection - Create a PetscSection over the reference cell with the layout from this space
920a4ce7ad1SMatthew G. Knepley 
921a4ce7ad1SMatthew G. Knepley   Collective on sp
922a4ce7ad1SMatthew G. Knepley 
923a4ce7ad1SMatthew G. Knepley   Input Parameters:
924f0fc11ceSJed Brown . sp      - The PetscDualSpace
925a4ce7ad1SMatthew G. Knepley 
926a4ce7ad1SMatthew G. Knepley   Output Parameter:
927a4ce7ad1SMatthew G. Knepley . section - The section
928a4ce7ad1SMatthew G. Knepley 
929a4ce7ad1SMatthew G. Knepley   Level: advanced
930a4ce7ad1SMatthew G. Knepley 
931a4ce7ad1SMatthew G. Knepley .seealso: PetscDualSpaceCreate(), DMPLEX
932a4ce7ad1SMatthew G. Knepley @*/
933*b4457527SToby Isaac PetscErrorCode PetscDualSpaceGetSection(PetscDualSpace sp, PetscSection *section)
93420cf1dd8SToby Isaac {
935*b4457527SToby Isaac   PetscInt       pStart, pEnd, p;
936*b4457527SToby Isaac   PetscErrorCode ierr;
937*b4457527SToby Isaac 
938*b4457527SToby Isaac   PetscFunctionBegin;
939*b4457527SToby Isaac   if (!sp->pointSection) {
940*b4457527SToby Isaac     /* mark the boundary */
941*b4457527SToby Isaac     ierr = PetscDualSpaceSectionCreate_Internal(sp, &(sp->pointSection));CHKERRQ(ierr);
942*b4457527SToby Isaac     ierr = DMPlexGetChart(sp->dm,&pStart,&pEnd);CHKERRQ(ierr);
943*b4457527SToby Isaac     for (p = pStart; p < pEnd; p++) {
944*b4457527SToby Isaac       PetscDualSpace psp;
945*b4457527SToby Isaac 
946*b4457527SToby Isaac       ierr = PetscDualSpaceGetPointSubspace(sp, p, &psp);CHKERRQ(ierr);
947*b4457527SToby Isaac       if (psp) {
948*b4457527SToby Isaac         PetscInt dof;
949*b4457527SToby Isaac 
950*b4457527SToby Isaac         ierr = PetscDualSpaceGetInteriorDimension(psp, &dof);CHKERRQ(ierr);
951*b4457527SToby Isaac         ierr = PetscSectionSetDof(sp->pointSection,p,dof);CHKERRQ(ierr);
952*b4457527SToby Isaac       }
953*b4457527SToby Isaac     }
954*b4457527SToby Isaac     ierr = PetscDualSpaceSectionSetUp_Internal(sp,sp->pointSection);CHKERRQ(ierr);
955*b4457527SToby Isaac   }
956*b4457527SToby Isaac   *section = sp->pointSection;
957*b4457527SToby Isaac   PetscFunctionReturn(0);
958*b4457527SToby Isaac }
959*b4457527SToby Isaac 
960*b4457527SToby Isaac /* this assumes that all of the point dual spaces store their interior dofs first, which is true when the point DMs
961*b4457527SToby Isaac  * have one cell */
962*b4457527SToby Isaac PetscErrorCode PetscDualSpacePushForwardSubspaces_Internal(PetscDualSpace sp, PetscInt sStart, PetscInt sEnd)
963*b4457527SToby Isaac {
964*b4457527SToby Isaac   PetscReal *sv0, *v0, *J;
965*b4457527SToby Isaac   PetscSection section;
966*b4457527SToby Isaac   PetscInt     dim, s, k;
96720cf1dd8SToby Isaac   DM             dm;
96820cf1dd8SToby Isaac   PetscErrorCode ierr;
96920cf1dd8SToby Isaac 
97020cf1dd8SToby Isaac   PetscFunctionBegin;
97120cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr);
972*b4457527SToby Isaac   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
973*b4457527SToby Isaac   ierr = PetscDualSpaceGetSection(sp, &section);CHKERRQ(ierr);
974*b4457527SToby Isaac   ierr = PetscMalloc3(dim, &v0, dim, &sv0, dim*dim, &J);CHKERRQ(ierr);
975*b4457527SToby Isaac   ierr = PetscDualSpaceGetFormDegree(sp, &k);CHKERRQ(ierr);
976*b4457527SToby Isaac   for (s = sStart; s < sEnd; s++) {
977*b4457527SToby Isaac     PetscReal detJ, hdetJ;
978*b4457527SToby Isaac     PetscDualSpace ssp;
979*b4457527SToby Isaac     PetscInt dof, off, f, sdim;
980*b4457527SToby Isaac     PetscInt i, j;
981*b4457527SToby Isaac     DM sdm;
98220cf1dd8SToby Isaac 
983*b4457527SToby Isaac     ierr = PetscDualSpaceGetPointSubspace(sp, s, &ssp);CHKERRQ(ierr);
984*b4457527SToby Isaac     if (!ssp) continue;
985*b4457527SToby Isaac     ierr = PetscSectionGetDof(section, s, &dof);CHKERRQ(ierr);
986*b4457527SToby Isaac     ierr = PetscSectionGetOffset(section, s, &off);CHKERRQ(ierr);
987*b4457527SToby Isaac     /* get the first vertex of the reference cell */
988*b4457527SToby Isaac     ierr = PetscDualSpaceGetDM(ssp, &sdm);CHKERRQ(ierr);
989*b4457527SToby Isaac     ierr = DMGetDimension(sdm, &sdim);CHKERRQ(ierr);
990*b4457527SToby Isaac     ierr = DMPlexComputeCellGeometryAffineFEM(sdm, 0, sv0, NULL, NULL, &hdetJ);CHKERRQ(ierr);
991*b4457527SToby Isaac     ierr = DMPlexComputeCellGeometryAffineFEM(dm, s, v0, J, NULL, &detJ);CHKERRQ(ierr);
992*b4457527SToby Isaac     /* compactify Jacobian */
993*b4457527SToby Isaac     for (i = 0; i < dim; i++) for (j = 0; j < sdim; j++) J[i* sdim + j] = J[i * dim + j];
994*b4457527SToby Isaac     for (f = 0; f < dof; f++) {
995*b4457527SToby Isaac       PetscQuadrature fn;
99620cf1dd8SToby Isaac 
997*b4457527SToby Isaac       ierr = PetscDualSpaceGetFunctional(ssp, f, &fn);CHKERRQ(ierr);
998*b4457527SToby Isaac       ierr = PetscQuadraturePushForward(fn, dim, sv0, v0, J, k, &(sp->functional[off+f]));CHKERRQ(ierr);
99920cf1dd8SToby Isaac     }
100020cf1dd8SToby Isaac   }
1001*b4457527SToby Isaac   ierr = PetscFree3(v0, sv0, J);CHKERRQ(ierr);
100220cf1dd8SToby Isaac   PetscFunctionReturn(0);
100320cf1dd8SToby Isaac }
100420cf1dd8SToby Isaac 
100520cf1dd8SToby Isaac /*@
100620cf1dd8SToby Isaac   PetscDualSpaceCreateReferenceCell - Create a DMPLEX with the appropriate FEM reference cell
100720cf1dd8SToby Isaac 
1008d083f849SBarry Smith   Collective on sp
100920cf1dd8SToby Isaac 
101020cf1dd8SToby Isaac   Input Parameters:
101120cf1dd8SToby Isaac + sp      - The PetscDualSpace
101220cf1dd8SToby Isaac . dim     - The spatial dimension
101320cf1dd8SToby Isaac - simplex - Flag for simplex, otherwise use a tensor-product cell
101420cf1dd8SToby Isaac 
101520cf1dd8SToby Isaac   Output Parameter:
101620cf1dd8SToby Isaac . refdm - The reference cell
101720cf1dd8SToby Isaac 
1018a4ce7ad1SMatthew G. Knepley   Level: intermediate
101920cf1dd8SToby Isaac 
102020cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate(), DMPLEX
102120cf1dd8SToby Isaac @*/
102220cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceCreateReferenceCell(PetscDualSpace sp, PetscInt dim, PetscBool simplex, DM *refdm)
102320cf1dd8SToby Isaac {
102420cf1dd8SToby Isaac   PetscErrorCode ierr;
102520cf1dd8SToby Isaac 
102620cf1dd8SToby Isaac   PetscFunctionBeginUser;
102720cf1dd8SToby Isaac   ierr = DMPlexCreateReferenceCell(PetscObjectComm((PetscObject) sp), dim, simplex, refdm);CHKERRQ(ierr);
102820cf1dd8SToby Isaac   PetscFunctionReturn(0);
102920cf1dd8SToby Isaac }
103020cf1dd8SToby Isaac 
103120cf1dd8SToby Isaac /*@C
103220cf1dd8SToby Isaac   PetscDualSpaceApply - Apply a functional from the dual space basis to an input function
103320cf1dd8SToby Isaac 
103420cf1dd8SToby Isaac   Input Parameters:
103520cf1dd8SToby Isaac + sp      - The PetscDualSpace object
103620cf1dd8SToby Isaac . f       - The basis functional index
103720cf1dd8SToby Isaac . time    - The time
103820cf1dd8SToby 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)
103920cf1dd8SToby Isaac . numComp - The number of components for the function
104020cf1dd8SToby Isaac . func    - The input function
104120cf1dd8SToby Isaac - ctx     - A context for the function
104220cf1dd8SToby Isaac 
104320cf1dd8SToby Isaac   Output Parameter:
104420cf1dd8SToby Isaac . value   - numComp output values
104520cf1dd8SToby Isaac 
104620cf1dd8SToby Isaac   Note: The calling sequence for the callback func is given by:
104720cf1dd8SToby Isaac 
104820cf1dd8SToby Isaac $ func(PetscInt dim, PetscReal time, const PetscReal x[],
104920cf1dd8SToby Isaac $      PetscInt numComponents, PetscScalar values[], void *ctx)
105020cf1dd8SToby Isaac 
1051a4ce7ad1SMatthew G. Knepley   Level: beginner
105220cf1dd8SToby Isaac 
105320cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate()
105420cf1dd8SToby Isaac @*/
105520cf1dd8SToby 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)
105620cf1dd8SToby Isaac {
105720cf1dd8SToby Isaac   PetscErrorCode ierr;
105820cf1dd8SToby Isaac 
105920cf1dd8SToby Isaac   PetscFunctionBegin;
106020cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
106120cf1dd8SToby Isaac   PetscValidPointer(cgeom, 4);
106220cf1dd8SToby Isaac   PetscValidPointer(value, 8);
106320cf1dd8SToby Isaac   ierr = (*sp->ops->apply)(sp, f, time, cgeom, numComp, func, ctx, value);CHKERRQ(ierr);
106420cf1dd8SToby Isaac   PetscFunctionReturn(0);
106520cf1dd8SToby Isaac }
106620cf1dd8SToby Isaac 
106720cf1dd8SToby Isaac /*@C
1068*b4457527SToby Isaac   PetscDualSpaceApplyAll - Apply all functionals from the dual space basis to the result of an evaluation at the points returned by PetscDualSpaceGetAllData()
106920cf1dd8SToby Isaac 
107020cf1dd8SToby Isaac   Input Parameters:
107120cf1dd8SToby Isaac + sp        - The PetscDualSpace object
1072*b4457527SToby Isaac - pointEval - Evaluation at the points returned by PetscDualSpaceGetAllData()
107320cf1dd8SToby Isaac 
107420cf1dd8SToby Isaac   Output Parameter:
107520cf1dd8SToby Isaac . spValue   - The values of all dual space functionals
107620cf1dd8SToby Isaac 
1077a4ce7ad1SMatthew G. Knepley   Level: beginner
107820cf1dd8SToby Isaac 
107920cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate()
108020cf1dd8SToby Isaac @*/
108120cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceApplyAll(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue)
108220cf1dd8SToby Isaac {
108320cf1dd8SToby Isaac   PetscErrorCode ierr;
108420cf1dd8SToby Isaac 
108520cf1dd8SToby Isaac   PetscFunctionBegin;
108620cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
108720cf1dd8SToby Isaac   ierr = (*sp->ops->applyall)(sp, pointEval, spValue);CHKERRQ(ierr);
108820cf1dd8SToby Isaac   PetscFunctionReturn(0);
108920cf1dd8SToby Isaac }
109020cf1dd8SToby Isaac 
109120cf1dd8SToby Isaac /*@C
1092*b4457527SToby Isaac   PetscDualSpaceApplyInterior - Apply interior functionals from the dual space basis to the result of an evaluation at the points returned by PetscDualSpaceGetInteriorData()
1093*b4457527SToby Isaac 
1094*b4457527SToby Isaac   Input Parameters:
1095*b4457527SToby Isaac + sp        - The PetscDualSpace object
1096*b4457527SToby Isaac - pointEval - Evaluation at the points returned by PetscDualSpaceGetInteriorData()
1097*b4457527SToby Isaac 
1098*b4457527SToby Isaac   Output Parameter:
1099*b4457527SToby Isaac . spValue   - The values of interior dual space functionals
1100*b4457527SToby Isaac 
1101*b4457527SToby Isaac   Level: beginner
1102*b4457527SToby Isaac 
1103*b4457527SToby Isaac .seealso: PetscDualSpaceCreate()
1104*b4457527SToby Isaac @*/
1105*b4457527SToby Isaac PetscErrorCode PetscDualSpaceApplyInterior(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue)
1106*b4457527SToby Isaac {
1107*b4457527SToby Isaac   PetscErrorCode ierr;
1108*b4457527SToby Isaac 
1109*b4457527SToby Isaac   PetscFunctionBegin;
1110*b4457527SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1111*b4457527SToby Isaac   ierr = (*sp->ops->applyint)(sp, pointEval, spValue);CHKERRQ(ierr);
1112*b4457527SToby Isaac   PetscFunctionReturn(0);
1113*b4457527SToby Isaac }
1114*b4457527SToby Isaac 
1115*b4457527SToby Isaac /*@C
111620cf1dd8SToby Isaac   PetscDualSpaceApplyDefault - Apply a functional from the dual space basis to an input function by assuming a point evaluation functional.
111720cf1dd8SToby Isaac 
111820cf1dd8SToby Isaac   Input Parameters:
111920cf1dd8SToby Isaac + sp    - The PetscDualSpace object
112020cf1dd8SToby Isaac . f     - The basis functional index
112120cf1dd8SToby Isaac . time  - The time
112220cf1dd8SToby Isaac . cgeom - A context with geometric information for this cell, we use v0 (the initial vertex) and J (the Jacobian)
112320cf1dd8SToby Isaac . Nc    - The number of components for the function
112420cf1dd8SToby Isaac . func  - The input function
112520cf1dd8SToby Isaac - ctx   - A context for the function
112620cf1dd8SToby Isaac 
112720cf1dd8SToby Isaac   Output Parameter:
112820cf1dd8SToby Isaac . value   - The output value
112920cf1dd8SToby Isaac 
113020cf1dd8SToby Isaac   Note: The calling sequence for the callback func is given by:
113120cf1dd8SToby Isaac 
113220cf1dd8SToby Isaac $ func(PetscInt dim, PetscReal time, const PetscReal x[],
113320cf1dd8SToby Isaac $      PetscInt numComponents, PetscScalar values[], void *ctx)
113420cf1dd8SToby Isaac 
113520cf1dd8SToby Isaac and the idea is to evaluate the functional as an integral
113620cf1dd8SToby Isaac 
113720cf1dd8SToby Isaac $ n(f) = int dx n(x) . f(x)
113820cf1dd8SToby Isaac 
113920cf1dd8SToby Isaac where both n and f have Nc components.
114020cf1dd8SToby Isaac 
1141a4ce7ad1SMatthew G. Knepley   Level: beginner
114220cf1dd8SToby Isaac 
114320cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate()
114420cf1dd8SToby Isaac @*/
114520cf1dd8SToby 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)
114620cf1dd8SToby Isaac {
114720cf1dd8SToby Isaac   DM               dm;
114820cf1dd8SToby Isaac   PetscQuadrature  n;
114920cf1dd8SToby Isaac   const PetscReal *points, *weights;
115020cf1dd8SToby Isaac   PetscReal        x[3];
115120cf1dd8SToby Isaac   PetscScalar     *val;
115220cf1dd8SToby Isaac   PetscInt         dim, dE, qNc, c, Nq, q;
115320cf1dd8SToby Isaac   PetscBool        isAffine;
115420cf1dd8SToby Isaac   PetscErrorCode   ierr;
115520cf1dd8SToby Isaac 
115620cf1dd8SToby Isaac   PetscFunctionBegin;
115720cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
115820cf1dd8SToby Isaac   PetscValidPointer(value, 5);
115920cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr);
116020cf1dd8SToby Isaac   ierr = PetscDualSpaceGetFunctional(sp, f, &n);CHKERRQ(ierr);
116120cf1dd8SToby Isaac   ierr = PetscQuadratureGetData(n, &dim, &qNc, &Nq, &points, &weights);CHKERRQ(ierr);
116220cf1dd8SToby 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);
116320cf1dd8SToby Isaac   if (qNc != Nc) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature components %D != function components %D", qNc, Nc);
116420cf1dd8SToby Isaac   ierr = DMGetWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr);
116520cf1dd8SToby Isaac   *value = 0.0;
116620cf1dd8SToby Isaac   isAffine = cgeom->isAffine;
116720cf1dd8SToby Isaac   dE = cgeom->dimEmbed;
116820cf1dd8SToby Isaac   for (q = 0; q < Nq; ++q) {
116920cf1dd8SToby Isaac     if (isAffine) {
117020cf1dd8SToby Isaac       CoordinatesRefToReal(dE, cgeom->dim, cgeom->xi, cgeom->v, cgeom->J, &points[q*dim], x);
117120cf1dd8SToby Isaac       ierr = (*func)(dE, time, x, Nc, val, ctx);CHKERRQ(ierr);
117220cf1dd8SToby Isaac     } else {
117320cf1dd8SToby Isaac       ierr = (*func)(dE, time, &cgeom->v[dE*q], Nc, val, ctx);CHKERRQ(ierr);
117420cf1dd8SToby Isaac     }
117520cf1dd8SToby Isaac     for (c = 0; c < Nc; ++c) {
117620cf1dd8SToby Isaac       *value += val[c]*weights[q*Nc+c];
117720cf1dd8SToby Isaac     }
117820cf1dd8SToby Isaac   }
117920cf1dd8SToby Isaac   ierr = DMRestoreWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr);
118020cf1dd8SToby Isaac   PetscFunctionReturn(0);
118120cf1dd8SToby Isaac }
118220cf1dd8SToby Isaac 
118320cf1dd8SToby Isaac /*@C
1184*b4457527SToby Isaac   PetscDualSpaceApplyAllDefault - Apply all functionals from the dual space basis to the result of an evaluation at the points returned by PetscDualSpaceGetAllData()
118520cf1dd8SToby Isaac 
118620cf1dd8SToby Isaac   Input Parameters:
118720cf1dd8SToby Isaac + sp        - The PetscDualSpace object
1188*b4457527SToby Isaac - pointEval - Evaluation at the points returned by PetscDualSpaceGetAllData()
118920cf1dd8SToby Isaac 
119020cf1dd8SToby Isaac   Output Parameter:
119120cf1dd8SToby Isaac . spValue   - The values of all dual space functionals
119220cf1dd8SToby Isaac 
1193a4ce7ad1SMatthew G. Knepley   Level: beginner
119420cf1dd8SToby Isaac 
119520cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate()
119620cf1dd8SToby Isaac @*/
119720cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceApplyAllDefault(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue)
119820cf1dd8SToby Isaac {
1199*b4457527SToby Isaac   Vec              pointValues, dofValues;
1200*b4457527SToby Isaac   Mat              allMat;
120120cf1dd8SToby Isaac   PetscErrorCode   ierr;
120220cf1dd8SToby Isaac 
120320cf1dd8SToby Isaac   PetscFunctionBegin;
120420cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
120520cf1dd8SToby Isaac   PetscValidScalarPointer(pointEval, 2);
120620cf1dd8SToby Isaac   PetscValidScalarPointer(spValue, 5);
1207*b4457527SToby Isaac   ierr = PetscDualSpaceGetAllData(sp, NULL, &allMat);CHKERRQ(ierr);
1208*b4457527SToby Isaac   if (!(sp->allNodeValues)) {
1209*b4457527SToby Isaac     ierr = MatCreateVecs(allMat, &(sp->allNodeValues), NULL);CHKERRQ(ierr);
121020cf1dd8SToby Isaac   }
1211*b4457527SToby Isaac   pointValues = sp->allNodeValues;
1212*b4457527SToby Isaac   if (!(sp->allDofValues)) {
1213*b4457527SToby Isaac     ierr = MatCreateVecs(allMat, NULL, &(sp->allDofValues));CHKERRQ(ierr);
121420cf1dd8SToby Isaac   }
1215*b4457527SToby Isaac   dofValues = sp->allDofValues;
1216*b4457527SToby Isaac   ierr = VecPlaceArray(pointValues, pointEval);CHKERRQ(ierr);
1217*b4457527SToby Isaac   ierr = VecPlaceArray(dofValues, spValue);CHKERRQ(ierr);
1218*b4457527SToby Isaac   ierr = MatMult(allMat, pointValues, dofValues);CHKERRQ(ierr);
1219*b4457527SToby Isaac   ierr = VecResetArray(dofValues);CHKERRQ(ierr);
1220*b4457527SToby Isaac   ierr = VecResetArray(pointValues);CHKERRQ(ierr);
1221*b4457527SToby Isaac   PetscFunctionReturn(0);
122220cf1dd8SToby Isaac }
1223*b4457527SToby Isaac 
1224*b4457527SToby Isaac /*@C
1225*b4457527SToby Isaac   PetscDualSpaceApplyInteriorDefault - Apply interior functionals from the dual space basis to the result of an evaluation at the points returned by PetscDualSpaceGetInteriorData()
1226*b4457527SToby Isaac 
1227*b4457527SToby Isaac   Input Parameters:
1228*b4457527SToby Isaac + sp        - The PetscDualSpace object
1229*b4457527SToby Isaac - pointEval - Evaluation at the points returned by PetscDualSpaceGetInteriorData()
1230*b4457527SToby Isaac 
1231*b4457527SToby Isaac   Output Parameter:
1232*b4457527SToby Isaac . spValue   - The values of interior dual space functionals
1233*b4457527SToby Isaac 
1234*b4457527SToby Isaac   Level: beginner
1235*b4457527SToby Isaac 
1236*b4457527SToby Isaac .seealso: PetscDualSpaceCreate()
1237*b4457527SToby Isaac @*/
1238*b4457527SToby Isaac PetscErrorCode PetscDualSpaceApplyInteriorDefault(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue)
1239*b4457527SToby Isaac {
1240*b4457527SToby Isaac   Vec              pointValues, dofValues;
1241*b4457527SToby Isaac   Mat              intMat;
1242*b4457527SToby Isaac   PetscErrorCode   ierr;
1243*b4457527SToby Isaac 
1244*b4457527SToby Isaac   PetscFunctionBegin;
1245*b4457527SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1246*b4457527SToby Isaac   PetscValidScalarPointer(pointEval, 2);
1247*b4457527SToby Isaac   PetscValidScalarPointer(spValue, 5);
1248*b4457527SToby Isaac   ierr = PetscDualSpaceGetInteriorData(sp, NULL, &intMat);CHKERRQ(ierr);
1249*b4457527SToby Isaac   if (!(sp->intNodeValues)) {
1250*b4457527SToby Isaac     ierr = MatCreateVecs(intMat, &(sp->intNodeValues), NULL);CHKERRQ(ierr);
1251*b4457527SToby Isaac   }
1252*b4457527SToby Isaac   pointValues = sp->intNodeValues;
1253*b4457527SToby Isaac   if (!(sp->intDofValues)) {
1254*b4457527SToby Isaac     ierr = MatCreateVecs(intMat, NULL, &(sp->intDofValues));CHKERRQ(ierr);
1255*b4457527SToby Isaac   }
1256*b4457527SToby Isaac   dofValues = sp->intDofValues;
1257*b4457527SToby Isaac   ierr = VecPlaceArray(pointValues, pointEval);CHKERRQ(ierr);
1258*b4457527SToby Isaac   ierr = VecPlaceArray(dofValues, spValue);CHKERRQ(ierr);
1259*b4457527SToby Isaac   ierr = MatMult(intMat, pointValues, dofValues);CHKERRQ(ierr);
1260*b4457527SToby Isaac   ierr = VecResetArray(dofValues);CHKERRQ(ierr);
1261*b4457527SToby Isaac   ierr = VecResetArray(pointValues);CHKERRQ(ierr);
126220cf1dd8SToby Isaac   PetscFunctionReturn(0);
126320cf1dd8SToby Isaac }
126420cf1dd8SToby Isaac 
1265a4ce7ad1SMatthew G. Knepley /*@
1266*b4457527SToby Isaac   PetscDualSpaceGetAllData - Get all quadrature nodes from this space, and the matrix that sends quadrature node values to degree-of-freedom values
1267a4ce7ad1SMatthew G. Knepley 
1268a4ce7ad1SMatthew G. Knepley   Input Parameter:
1269a4ce7ad1SMatthew G. Knepley . sp - The dualspace
1270a4ce7ad1SMatthew G. Knepley 
1271a4ce7ad1SMatthew G. Knepley   Output Parameter:
1272*b4457527SToby Isaac + allNodes - A PetscQuadrature object containing all evaluation nodes
1273*b4457527SToby Isaac - allMat - A Mat for the node-to-dof transformation
1274a4ce7ad1SMatthew G. Knepley 
1275a4ce7ad1SMatthew G. Knepley   Level: advanced
1276a4ce7ad1SMatthew G. Knepley 
1277a4ce7ad1SMatthew G. Knepley .seealso: PetscDualSpaceCreate()
1278a4ce7ad1SMatthew G. Knepley @*/
1279*b4457527SToby Isaac PetscErrorCode PetscDualSpaceGetAllData(PetscDualSpace sp, PetscQuadrature *allNodes, Mat *allMat)
128020cf1dd8SToby Isaac {
128120cf1dd8SToby Isaac   PetscErrorCode ierr;
128220cf1dd8SToby Isaac 
128320cf1dd8SToby Isaac   PetscFunctionBegin;
128420cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1285*b4457527SToby Isaac   if (allNodes) PetscValidPointer(allNodes,2);
1286*b4457527SToby Isaac   if (allMat) PetscValidPointer(allMat,3);
1287*b4457527SToby Isaac   if ((!sp->allNodes || !sp->allMat) && sp->ops->createalldata) {
1288*b4457527SToby Isaac     PetscQuadrature qpoints;
1289*b4457527SToby Isaac     Mat amat;
1290*b4457527SToby Isaac 
1291*b4457527SToby Isaac     ierr = (*sp->ops->createalldata)(sp,&qpoints,&amat);CHKERRQ(ierr);
1292*b4457527SToby Isaac     ierr = PetscQuadratureDestroy(&(sp->allNodes));CHKERRQ(ierr);
1293*b4457527SToby Isaac     ierr = MatDestroy(&(sp->allMat));CHKERRQ(ierr);
1294*b4457527SToby Isaac     sp->allNodes = qpoints;
1295*b4457527SToby Isaac     sp->allMat = amat;
129620cf1dd8SToby Isaac   }
1297*b4457527SToby Isaac   if (allNodes) *allNodes = sp->allNodes;
1298*b4457527SToby Isaac   if (allMat) *allMat = sp->allMat;
129920cf1dd8SToby Isaac   PetscFunctionReturn(0);
130020cf1dd8SToby Isaac }
130120cf1dd8SToby Isaac 
1302a4ce7ad1SMatthew G. Knepley /*@
1303*b4457527SToby Isaac   PetscDualSpaceCreateAllDataDefault - Create all evaluation nodes and the node-to-dof matrix by examining functionals
1304a4ce7ad1SMatthew G. Knepley 
1305a4ce7ad1SMatthew G. Knepley   Input Parameter:
1306a4ce7ad1SMatthew G. Knepley . sp - The dualspace
1307a4ce7ad1SMatthew G. Knepley 
1308a4ce7ad1SMatthew G. Knepley   Output Parameter:
1309*b4457527SToby Isaac + allNodes - A PetscQuadrature object containing all evaluation nodes
1310*b4457527SToby Isaac - allMat - A Mat for the node-to-dof transformation
1311a4ce7ad1SMatthew G. Knepley 
1312a4ce7ad1SMatthew G. Knepley   Level: advanced
1313a4ce7ad1SMatthew G. Knepley 
1314a4ce7ad1SMatthew G. Knepley .seealso: PetscDualSpaceCreate()
1315a4ce7ad1SMatthew G. Knepley @*/
1316*b4457527SToby Isaac PetscErrorCode PetscDualSpaceCreateAllDataDefault(PetscDualSpace sp, PetscQuadrature *allNodes, Mat *allMat)
131720cf1dd8SToby Isaac {
131820cf1dd8SToby Isaac   PetscInt        spdim;
131920cf1dd8SToby Isaac   PetscInt        numPoints, offset;
132020cf1dd8SToby Isaac   PetscReal       *points;
132120cf1dd8SToby Isaac   PetscInt        f, dim;
1322*b4457527SToby Isaac   PetscInt        Nc, nrows, ncols;
1323*b4457527SToby Isaac   PetscInt        maxNumPoints;
132420cf1dd8SToby Isaac   PetscQuadrature q;
1325*b4457527SToby Isaac   Mat             A;
132620cf1dd8SToby Isaac   PetscErrorCode  ierr;
132720cf1dd8SToby Isaac 
132820cf1dd8SToby Isaac   PetscFunctionBegin;
1329*b4457527SToby Isaac   ierr = PetscDualSpaceGetNumComponents(sp, &Nc);CHKERRQ(ierr);
133020cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDimension(sp,&spdim);CHKERRQ(ierr);
133120cf1dd8SToby Isaac   if (!spdim) {
1332*b4457527SToby Isaac     ierr = PetscQuadratureCreate(PETSC_COMM_SELF,allNodes);CHKERRQ(ierr);
1333*b4457527SToby Isaac     ierr = PetscQuadratureSetData(*allNodes,0,0,0,NULL,NULL);CHKERRQ(ierr);
133420cf1dd8SToby Isaac   }
1335*b4457527SToby Isaac   nrows = spdim;
133620cf1dd8SToby Isaac   ierr = PetscDualSpaceGetFunctional(sp,0,&q);CHKERRQ(ierr);
133720cf1dd8SToby Isaac   ierr = PetscQuadratureGetData(q,&dim,NULL,&numPoints,NULL,NULL);CHKERRQ(ierr);
1338*b4457527SToby Isaac   maxNumPoints = numPoints;
133920cf1dd8SToby Isaac   for (f = 1; f < spdim; f++) {
134020cf1dd8SToby Isaac     PetscInt Np;
134120cf1dd8SToby Isaac 
134220cf1dd8SToby Isaac     ierr = PetscDualSpaceGetFunctional(sp,f,&q);CHKERRQ(ierr);
134320cf1dd8SToby Isaac     ierr = PetscQuadratureGetData(q,NULL,NULL,&Np,NULL,NULL);CHKERRQ(ierr);
134420cf1dd8SToby Isaac     numPoints += Np;
1345*b4457527SToby Isaac     maxNumPoints = PetscMax(maxNumPoints,Np);
134620cf1dd8SToby Isaac   }
1347*b4457527SToby Isaac   ncols = numPoints * Nc;
134820cf1dd8SToby Isaac   ierr = PetscMalloc1(dim*numPoints,&points);CHKERRQ(ierr);
1349*b4457527SToby Isaac   ierr = MatCreateSeqAIJ(PETSC_COMM_SELF, nrows, ncols, maxNumPoints * Nc, NULL, &A);CHKERRQ(ierr);
135020cf1dd8SToby Isaac   for (f = 0, offset = 0; f < spdim; f++) {
1351*b4457527SToby Isaac     const PetscReal *p, *w;
135220cf1dd8SToby Isaac     PetscInt        Np, i;
1353*b4457527SToby Isaac     PetscInt        fnc;
135420cf1dd8SToby Isaac 
135520cf1dd8SToby Isaac     ierr = PetscDualSpaceGetFunctional(sp,f,&q);CHKERRQ(ierr);
1356*b4457527SToby Isaac     ierr = PetscQuadratureGetData(q,NULL,&fnc,&Np,&p,&w);CHKERRQ(ierr);
1357*b4457527SToby Isaac     if (fnc != Nc) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "functional component mismatch");
1358*b4457527SToby Isaac     for (i = 0; i < Np * dim; i++) {
1359*b4457527SToby Isaac       points[offset* dim + i] = p[i];
1360*b4457527SToby Isaac     }
1361*b4457527SToby Isaac     for (i = 0; i < Np * Nc; i++) {
1362*b4457527SToby Isaac       ierr = MatSetValue(A, f, offset * Nc, w[i], INSERT_VALUES);CHKERRQ(ierr);
1363*b4457527SToby Isaac     }
1364*b4457527SToby Isaac     offset += Np;
1365*b4457527SToby Isaac   }
1366*b4457527SToby Isaac   ierr = MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1367*b4457527SToby Isaac   ierr = MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1368*b4457527SToby Isaac   ierr = PetscQuadratureCreate(PETSC_COMM_SELF,allNodes);CHKERRQ(ierr);
1369*b4457527SToby Isaac   ierr = PetscQuadratureSetData(*allNodes,dim,0,numPoints,points,NULL);CHKERRQ(ierr);
1370*b4457527SToby Isaac   *allMat = A;
1371*b4457527SToby Isaac   PetscFunctionReturn(0);
1372*b4457527SToby Isaac }
1373*b4457527SToby Isaac 
1374*b4457527SToby Isaac /*@
1375*b4457527SToby Isaac   PetscDualSpaceGetInteriorData - Get all quadrature points necessary to compute the interior degrees of freedom from
1376*b4457527SToby Isaac   this space, as well as the matrix that computes the degrees of freedom from the quadrature values.  Degrees of
1377*b4457527SToby Isaac   freedom are interior degrees of freedom if they belong (by PetscDualSpaceGetSection()) to interior points in the
1378*b4457527SToby Isaac   reference DMPlex: complementary boundary degrees of freedom are marked as constrained in the section returned by
1379*b4457527SToby Isaac   PetscDualSpaceGetSection()).
1380*b4457527SToby Isaac 
1381*b4457527SToby Isaac   Input Parameter:
1382*b4457527SToby Isaac . sp - The dualspace
1383*b4457527SToby Isaac 
1384*b4457527SToby Isaac   Output Parameter:
1385*b4457527SToby Isaac + intNodes - A PetscQuadrature object containing all evaluation points needed to evaluate interior degrees of freedom
1386*b4457527SToby Isaac - intMat   - A matrix that computes dual space values from point values: size [spdim0 x (npoints * nc)], where spdim0 is
1387*b4457527SToby Isaac              the size of the constrained layout (PetscSectionGetConstrainStorageSize()) of the dual space section,
1388*b4457527SToby Isaac              npoints is the number of points in intNodes and nc is PetscDualSpaceGetNumComponents().
1389*b4457527SToby Isaac 
1390*b4457527SToby Isaac   Level: advanced
1391*b4457527SToby Isaac 
1392*b4457527SToby Isaac .seealso: PetscDualSpaceCreate(), PetscDualSpaceGetDimension(), PetscDualSpaceGetNumComponents(), PetscQuadratureGetData()
1393*b4457527SToby Isaac @*/
1394*b4457527SToby Isaac PetscErrorCode PetscDualSpaceGetInteriorData(PetscDualSpace sp, PetscQuadrature *intNodes, Mat *intMat)
1395*b4457527SToby Isaac {
1396*b4457527SToby Isaac   PetscErrorCode ierr;
1397*b4457527SToby Isaac 
1398*b4457527SToby Isaac   PetscFunctionBegin;
1399*b4457527SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1400*b4457527SToby Isaac   if (intNodes) PetscValidPointer(intNodes,2);
1401*b4457527SToby Isaac   if (intMat) PetscValidPointer(intMat,3);
1402*b4457527SToby Isaac   if ((!sp->intNodes || !sp->intMat) && sp->ops->createintdata) {
1403*b4457527SToby Isaac     PetscQuadrature qpoints;
1404*b4457527SToby Isaac     Mat imat;
1405*b4457527SToby Isaac 
1406*b4457527SToby Isaac     ierr = (*sp->ops->createintdata)(sp,&qpoints,&imat);CHKERRQ(ierr);
1407*b4457527SToby Isaac     ierr = PetscQuadratureDestroy(&(sp->intNodes));CHKERRQ(ierr);
1408*b4457527SToby Isaac     ierr = MatDestroy(&(sp->intMat));CHKERRQ(ierr);
1409*b4457527SToby Isaac     sp->intNodes = qpoints;
1410*b4457527SToby Isaac     sp->intMat = imat;
1411*b4457527SToby Isaac   }
1412*b4457527SToby Isaac   if (intNodes) *intNodes = sp->intNodes;
1413*b4457527SToby Isaac   if (intMat) *intMat = sp->intMat;
1414*b4457527SToby Isaac   PetscFunctionReturn(0);
1415*b4457527SToby Isaac }
1416*b4457527SToby Isaac 
1417*b4457527SToby Isaac /*@
1418*b4457527SToby Isaac   PetscDualSpaceCreateInteriorDataDefault - Create quadrature points by examining interior functionals and create the matrix mapping quadrature point values to interior dual space values
1419*b4457527SToby Isaac 
1420*b4457527SToby Isaac   Input Parameter:
1421*b4457527SToby Isaac . sp - The dualspace
1422*b4457527SToby Isaac 
1423*b4457527SToby Isaac   Output Parameter:
1424*b4457527SToby Isaac + intNodes - A PetscQuadrature object containing all evaluation points needed to evaluate interior degrees of freedom
1425*b4457527SToby Isaac - intMat    - A matrix that computes dual space values from point values: size [spdim0 x (npoints * nc)], where spdim0 is
1426*b4457527SToby Isaac               the size of the constrained layout (PetscSectionGetConstrainStorageSize()) of the dual space section,
1427*b4457527SToby Isaac               npoints is the number of points in allNodes and nc is PetscDualSpaceGetNumComponents().
1428*b4457527SToby Isaac 
1429*b4457527SToby Isaac   Level: advanced
1430*b4457527SToby Isaac 
1431*b4457527SToby Isaac .seealso: PetscDualSpaceCreate(), PetscDualSpaceGetInteriorData()
1432*b4457527SToby Isaac @*/
1433*b4457527SToby Isaac PetscErrorCode PetscDualSpaceCreateInteriorDataDefault(PetscDualSpace sp, PetscQuadrature *intNodes, Mat *intMat)
1434*b4457527SToby Isaac {
1435*b4457527SToby Isaac   DM              dm;
1436*b4457527SToby Isaac   PetscInt        spdim0;
1437*b4457527SToby Isaac   PetscInt        Nc;
1438*b4457527SToby Isaac   PetscInt        pStart, pEnd, p, f;
1439*b4457527SToby Isaac   PetscSection    section;
1440*b4457527SToby Isaac   PetscInt        numPoints, offset, matoffset;
1441*b4457527SToby Isaac   PetscReal       *points;
1442*b4457527SToby Isaac   PetscInt        dim;
1443*b4457527SToby Isaac   PetscInt        *nnz;
1444*b4457527SToby Isaac   PetscQuadrature q;
1445*b4457527SToby Isaac   Mat             imat;
1446*b4457527SToby Isaac   PetscErrorCode  ierr;
1447*b4457527SToby Isaac 
1448*b4457527SToby Isaac   PetscFunctionBegin;
1449*b4457527SToby Isaac   PetscValidHeaderSpecific(sp,PETSCDUALSPACE_CLASSID,1);
1450*b4457527SToby Isaac   ierr = PetscDualSpaceGetSection(sp, &section);CHKERRQ(ierr);
1451*b4457527SToby Isaac   ierr = PetscSectionGetConstrainedStorageSize(section, &spdim0);CHKERRQ(ierr);
1452*b4457527SToby Isaac   if (!spdim0) {
1453*b4457527SToby Isaac     *intNodes = NULL;
1454*b4457527SToby Isaac     *intMat = NULL;
1455*b4457527SToby Isaac     PetscFunctionReturn(0);
1456*b4457527SToby Isaac   }
1457*b4457527SToby Isaac   ierr = PetscDualSpaceGetNumComponents(sp, &Nc);CHKERRQ(ierr);
1458*b4457527SToby Isaac   ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
1459*b4457527SToby Isaac   ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr);
1460*b4457527SToby Isaac   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1461*b4457527SToby Isaac   ierr = PetscMalloc1(spdim0, &nnz);CHKERRQ(ierr);
1462*b4457527SToby Isaac   for (p = pStart, f = 0, numPoints = 0; p < pEnd; p++) {
1463*b4457527SToby Isaac     PetscInt dof, cdof, off, d;
1464*b4457527SToby Isaac 
1465*b4457527SToby Isaac     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
1466*b4457527SToby Isaac     ierr = PetscSectionGetConstraintDof(section, p, &cdof);CHKERRQ(ierr);
1467*b4457527SToby Isaac     if (!(dof - cdof)) continue;
1468*b4457527SToby Isaac     ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
1469*b4457527SToby Isaac     for (d = 0; d < dof; d++, off++, f++) {
1470*b4457527SToby Isaac       PetscInt Np;
1471*b4457527SToby Isaac 
1472*b4457527SToby Isaac       ierr = PetscDualSpaceGetFunctional(sp,off,&q);CHKERRQ(ierr);
1473*b4457527SToby Isaac       ierr = PetscQuadratureGetData(q,NULL,NULL,&Np,NULL,NULL);CHKERRQ(ierr);
1474*b4457527SToby Isaac       nnz[f] = Np * Nc;
1475*b4457527SToby Isaac       numPoints += Np;
1476*b4457527SToby Isaac     }
1477*b4457527SToby Isaac   }
1478*b4457527SToby Isaac   ierr = MatCreateSeqAIJ(PETSC_COMM_SELF, spdim0, numPoints * Nc, 0, nnz, &imat);CHKERRQ(ierr);
1479*b4457527SToby Isaac   ierr = PetscFree(nnz);CHKERRQ(ierr);
1480*b4457527SToby Isaac   ierr = PetscMalloc1(dim*numPoints,&points);CHKERRQ(ierr);
1481*b4457527SToby Isaac   for (p = pStart, f = 0, offset = 0, matoffset = 0; p < pEnd; p++) {
1482*b4457527SToby Isaac     PetscInt dof, cdof, off, d;
1483*b4457527SToby Isaac 
1484*b4457527SToby Isaac     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
1485*b4457527SToby Isaac     ierr = PetscSectionGetConstraintDof(section, p, &cdof);CHKERRQ(ierr);
1486*b4457527SToby Isaac     if (!(dof - cdof)) continue;
1487*b4457527SToby Isaac     ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
1488*b4457527SToby Isaac     for (d = 0; d < dof; d++, off++, f++) {
1489*b4457527SToby Isaac       const PetscReal *p;
1490*b4457527SToby Isaac       const PetscReal *w;
1491*b4457527SToby Isaac       PetscInt        Np, i;
1492*b4457527SToby Isaac 
1493*b4457527SToby Isaac       ierr = PetscDualSpaceGetFunctional(sp,off,&q);CHKERRQ(ierr);
1494*b4457527SToby Isaac       ierr = PetscQuadratureGetData(q,NULL,NULL,&Np,&p,&w);CHKERRQ(ierr);
149520cf1dd8SToby Isaac       for (i = 0; i < Np * dim; i++) {
149620cf1dd8SToby Isaac         points[offset + i] = p[i];
149720cf1dd8SToby Isaac       }
1498*b4457527SToby Isaac       for (i = 0; i < Np * Nc; i++) {
1499*b4457527SToby Isaac         ierr = MatSetValue(imat, f, matoffset + i, w[i],INSERT_VALUES);CHKERRQ(ierr);
150020cf1dd8SToby Isaac       }
1501*b4457527SToby Isaac       offset += Np * dim;
1502*b4457527SToby Isaac       matoffset += Np * Nc;
1503*b4457527SToby Isaac     }
1504*b4457527SToby Isaac   }
1505*b4457527SToby Isaac   ierr = PetscQuadratureCreate(PETSC_COMM_SELF,intNodes);CHKERRQ(ierr);
1506*b4457527SToby Isaac   ierr = PetscQuadratureSetData(*intNodes,dim,0,numPoints,points,NULL);CHKERRQ(ierr);
1507*b4457527SToby Isaac   ierr = MatAssemblyBegin(imat, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1508*b4457527SToby Isaac   ierr = MatAssemblyEnd(imat, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1509*b4457527SToby Isaac   *intMat = imat;
151020cf1dd8SToby Isaac   PetscFunctionReturn(0);
151120cf1dd8SToby Isaac }
151220cf1dd8SToby Isaac 
151320cf1dd8SToby Isaac /*@C
151420cf1dd8SToby Isaac   PetscDualSpaceApplyFVM - Apply a functional from the dual space basis to an input function by assuming a point evaluation functional at the cell centroid.
151520cf1dd8SToby Isaac 
151620cf1dd8SToby Isaac   Input Parameters:
151720cf1dd8SToby Isaac + sp    - The PetscDualSpace object
151820cf1dd8SToby Isaac . f     - The basis functional index
151920cf1dd8SToby Isaac . time  - The time
152020cf1dd8SToby Isaac . cgeom - A context with geometric information for this cell, we currently just use the centroid
152120cf1dd8SToby Isaac . Nc    - The number of components for the function
152220cf1dd8SToby Isaac . func  - The input function
152320cf1dd8SToby Isaac - ctx   - A context for the function
152420cf1dd8SToby Isaac 
152520cf1dd8SToby Isaac   Output Parameter:
152620cf1dd8SToby Isaac . value - The output value (scalar)
152720cf1dd8SToby Isaac 
152820cf1dd8SToby Isaac   Note: The calling sequence for the callback func is given by:
152920cf1dd8SToby Isaac 
153020cf1dd8SToby Isaac $ func(PetscInt dim, PetscReal time, const PetscReal x[],
153120cf1dd8SToby Isaac $      PetscInt numComponents, PetscScalar values[], void *ctx)
153220cf1dd8SToby Isaac 
153320cf1dd8SToby Isaac and the idea is to evaluate the functional as an integral
153420cf1dd8SToby Isaac 
153520cf1dd8SToby Isaac $ n(f) = int dx n(x) . f(x)
153620cf1dd8SToby Isaac 
153720cf1dd8SToby Isaac where both n and f have Nc components.
153820cf1dd8SToby Isaac 
1539a4ce7ad1SMatthew G. Knepley   Level: beginner
154020cf1dd8SToby Isaac 
154120cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate()
154220cf1dd8SToby Isaac @*/
154320cf1dd8SToby 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)
154420cf1dd8SToby Isaac {
154520cf1dd8SToby Isaac   DM               dm;
154620cf1dd8SToby Isaac   PetscQuadrature  n;
154720cf1dd8SToby Isaac   const PetscReal *points, *weights;
154820cf1dd8SToby Isaac   PetscScalar     *val;
154920cf1dd8SToby Isaac   PetscInt         dimEmbed, qNc, c, Nq, q;
155020cf1dd8SToby Isaac   PetscErrorCode   ierr;
155120cf1dd8SToby Isaac 
155220cf1dd8SToby Isaac   PetscFunctionBegin;
155320cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
155420cf1dd8SToby Isaac   PetscValidPointer(value, 5);
155520cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr);
155620cf1dd8SToby Isaac   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
155720cf1dd8SToby Isaac   ierr = PetscDualSpaceGetFunctional(sp, f, &n);CHKERRQ(ierr);
155820cf1dd8SToby Isaac   ierr = PetscQuadratureGetData(n, NULL, &qNc, &Nq, &points, &weights);CHKERRQ(ierr);
155920cf1dd8SToby Isaac   if (qNc != Nc) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature components %D != function components %D", qNc, Nc);
156020cf1dd8SToby Isaac   ierr = DMGetWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr);
156120cf1dd8SToby Isaac   *value = 0.;
156220cf1dd8SToby Isaac   for (q = 0; q < Nq; ++q) {
156320cf1dd8SToby Isaac     ierr = (*func)(dimEmbed, time, cgeom->centroid, Nc, val, ctx);CHKERRQ(ierr);
156420cf1dd8SToby Isaac     for (c = 0; c < Nc; ++c) {
156520cf1dd8SToby Isaac       *value += val[c]*weights[q*Nc+c];
156620cf1dd8SToby Isaac     }
156720cf1dd8SToby Isaac   }
156820cf1dd8SToby Isaac   ierr = DMRestoreWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr);
156920cf1dd8SToby Isaac   PetscFunctionReturn(0);
157020cf1dd8SToby Isaac }
157120cf1dd8SToby Isaac 
157220cf1dd8SToby Isaac /*@
157320cf1dd8SToby Isaac   PetscDualSpaceGetHeightSubspace - Get the subset of the dual space basis that is supported on a mesh point of a
157420cf1dd8SToby Isaac   given height.  This assumes that the reference cell is symmetric over points of this height.
157520cf1dd8SToby Isaac 
157620cf1dd8SToby Isaac   If the dual space is not defined on mesh points of the given height (e.g. if the space is discontinuous and
157720cf1dd8SToby Isaac   pointwise values are not defined on the element boundaries), or if the implementation of PetscDualSpace does not
157820cf1dd8SToby Isaac   support extracting subspaces, then NULL is returned.
157920cf1dd8SToby Isaac 
158020cf1dd8SToby Isaac   This does not increment the reference count on the returned dual space, and the user should not destroy it.
158120cf1dd8SToby Isaac 
158220cf1dd8SToby Isaac   Not collective
158320cf1dd8SToby Isaac 
158420cf1dd8SToby Isaac   Input Parameters:
158520cf1dd8SToby Isaac + sp - the PetscDualSpace object
158620cf1dd8SToby Isaac - height - the height of the mesh point for which the subspace is desired
158720cf1dd8SToby Isaac 
158820cf1dd8SToby Isaac   Output Parameter:
158920cf1dd8SToby Isaac . subsp - the subspace.  Note that the functionals in the subspace are with respect to the intrinsic geometry of the
159020cf1dd8SToby Isaac   point, which will be of lesser dimension if height > 0.
159120cf1dd8SToby Isaac 
159220cf1dd8SToby Isaac   Level: advanced
159320cf1dd8SToby Isaac 
159420cf1dd8SToby Isaac .seealso: PetscSpaceGetHeightSubspace(), PetscDualSpace
159520cf1dd8SToby Isaac @*/
159620cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetHeightSubspace(PetscDualSpace sp, PetscInt height, PetscDualSpace *subsp)
159720cf1dd8SToby Isaac {
1598*b4457527SToby Isaac   PetscInt       depth = -1, cStart, cEnd;
1599*b4457527SToby Isaac   DM             dm;
160020cf1dd8SToby Isaac   PetscErrorCode ierr;
160120cf1dd8SToby Isaac 
160220cf1dd8SToby Isaac   PetscFunctionBegin;
160320cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1604*b4457527SToby Isaac   PetscValidPointer(subsp,2);
1605*b4457527SToby Isaac   if (!(sp->uniform)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "A non-uniform dual space does not have a single dual space at each height");
160620cf1dd8SToby Isaac   *subsp = NULL;
1607*b4457527SToby Isaac   dm = sp->dm;
1608*b4457527SToby Isaac   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
1609*b4457527SToby Isaac   if (height < 0 || height > depth) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid height");
1610*b4457527SToby Isaac   ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr);
1611*b4457527SToby Isaac   if (height == 0 && cEnd == cStart + 1) {
1612*b4457527SToby Isaac     *subsp = sp;
1613*b4457527SToby Isaac     PetscFunctionReturn(0);
1614*b4457527SToby Isaac   }
1615*b4457527SToby Isaac   if (!sp->heightSpaces) {
1616*b4457527SToby Isaac     PetscInt h;
1617*b4457527SToby Isaac     ierr = PetscCalloc1(depth+1, &(sp->heightSpaces));CHKERRQ(ierr);
1618*b4457527SToby Isaac 
1619*b4457527SToby Isaac     for (h = 0; h <= depth; h++) {
1620*b4457527SToby Isaac       if (h == 0 && cEnd == cStart + 1) continue;
1621*b4457527SToby Isaac       if (sp->ops->createheightsubspace) {ierr = (*sp->ops->createheightsubspace)(sp,height,&(sp->heightSpaces[h]));CHKERRQ(ierr);}
1622*b4457527SToby Isaac       else if (sp->pointSpaces) {
1623*b4457527SToby Isaac         PetscInt hStart, hEnd;
1624*b4457527SToby Isaac 
1625*b4457527SToby Isaac         ierr = DMPlexGetHeightStratum(dm,h,&hStart,&hEnd);CHKERRQ(ierr);
1626*b4457527SToby Isaac         if (hEnd > hStart) {
1627*b4457527SToby Isaac           ierr = PetscObjectReference((PetscObject)(sp->pointSpaces[hStart]));CHKERRQ(ierr);
1628*b4457527SToby Isaac           sp->heightSpaces[h] = sp->pointSpaces[hStart];
1629*b4457527SToby Isaac         }
1630*b4457527SToby Isaac       }
1631*b4457527SToby Isaac     }
1632*b4457527SToby Isaac   }
1633*b4457527SToby Isaac   *subsp = sp->heightSpaces[height];
163420cf1dd8SToby Isaac   PetscFunctionReturn(0);
163520cf1dd8SToby Isaac }
163620cf1dd8SToby Isaac 
163720cf1dd8SToby Isaac /*@
163820cf1dd8SToby Isaac   PetscDualSpaceGetPointSubspace - Get the subset of the dual space basis that is supported on a particular mesh point.
163920cf1dd8SToby Isaac 
164020cf1dd8SToby Isaac   If the dual space is not defined on the mesh point (e.g. if the space is discontinuous and pointwise values are not
164120cf1dd8SToby Isaac   defined on the element boundaries), or if the implementation of PetscDualSpace does not support extracting
164220cf1dd8SToby Isaac   subspaces, then NULL is returned.
164320cf1dd8SToby Isaac 
164420cf1dd8SToby Isaac   This does not increment the reference count on the returned dual space, and the user should not destroy it.
164520cf1dd8SToby Isaac 
164620cf1dd8SToby Isaac   Not collective
164720cf1dd8SToby Isaac 
164820cf1dd8SToby Isaac   Input Parameters:
164920cf1dd8SToby Isaac + sp - the PetscDualSpace object
165020cf1dd8SToby Isaac - point - the point (in the dual space's DM) for which the subspace is desired
165120cf1dd8SToby Isaac 
165220cf1dd8SToby Isaac   Output Parameters:
165320cf1dd8SToby Isaac   bdsp - the subspace.  Note that the functionals in the subspace are with respect to the intrinsic geometry of the
165420cf1dd8SToby Isaac   point, which will be of lesser dimension if height > 0.
165520cf1dd8SToby Isaac 
165620cf1dd8SToby Isaac   Level: advanced
165720cf1dd8SToby Isaac 
165820cf1dd8SToby Isaac .seealso: PetscDualSpace
165920cf1dd8SToby Isaac @*/
166020cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetPointSubspace(PetscDualSpace sp, PetscInt point, PetscDualSpace *bdsp)
166120cf1dd8SToby Isaac {
1662*b4457527SToby Isaac   PetscInt       pStart = 0, pEnd = 0, cStart, cEnd;
1663*b4457527SToby Isaac   DM             dm;
166420cf1dd8SToby Isaac   PetscErrorCode ierr;
166520cf1dd8SToby Isaac 
166620cf1dd8SToby Isaac   PetscFunctionBegin;
166720cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
166820cf1dd8SToby Isaac   PetscValidPointer(bdsp,2);
166920cf1dd8SToby Isaac   *bdsp = NULL;
1670*b4457527SToby Isaac   dm = sp->dm;
1671*b4457527SToby Isaac   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
1672*b4457527SToby Isaac   if (point < pStart || point > pEnd) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid point");
1673*b4457527SToby Isaac   ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr);
1674*b4457527SToby Isaac   if (point == cStart && cEnd == cStart + 1) { /* the dual space is only equivalent to the dual space on a cell if the reference mesh has just one cell */
1675*b4457527SToby Isaac     *bdsp = sp;
1676*b4457527SToby Isaac     PetscFunctionReturn(0);
1677*b4457527SToby Isaac   }
1678*b4457527SToby Isaac   if (!sp->pointSpaces) {
1679*b4457527SToby Isaac     PetscInt p;
1680*b4457527SToby Isaac     ierr = PetscCalloc1(pEnd - pStart, &(sp->pointSpaces));CHKERRQ(ierr);
168120cf1dd8SToby Isaac 
1682*b4457527SToby Isaac     for (p = 0; p < pEnd - pStart; p++) {
1683*b4457527SToby Isaac       if (p + pStart == cStart && cEnd == cStart + 1) continue;
1684*b4457527SToby Isaac       if (sp->ops->createpointsubspace) {ierr = (*sp->ops->createpointsubspace)(sp,p+pStart,&(sp->pointSpaces[p]));CHKERRQ(ierr);}
1685*b4457527SToby Isaac       else if (sp->heightSpaces || sp->ops->createheightsubspace) {
1686*b4457527SToby Isaac         PetscInt dim, depth, height;
1687*b4457527SToby Isaac         DMLabel  label;
1688*b4457527SToby Isaac 
168920cf1dd8SToby Isaac         ierr = DMPlexGetDepth(dm,&dim);CHKERRQ(ierr);
169020cf1dd8SToby Isaac         ierr = DMPlexGetDepthLabel(dm,&label);CHKERRQ(ierr);
1691*b4457527SToby Isaac         ierr = DMLabelGetValue(label,p+pStart,&depth);CHKERRQ(ierr);
169220cf1dd8SToby Isaac         height = dim - depth;
1693*b4457527SToby Isaac         ierr = PetscDualSpaceGetHeightSubspace(sp, height, &(sp->pointSpaces[p]));CHKERRQ(ierr);
1694*b4457527SToby Isaac         ierr = PetscObjectReference((PetscObject)sp->pointSpaces[p]);CHKERRQ(ierr);
169520cf1dd8SToby Isaac       }
1696*b4457527SToby Isaac     }
1697*b4457527SToby Isaac   }
1698*b4457527SToby Isaac   *bdsp = sp->pointSpaces[point - pStart];
169920cf1dd8SToby Isaac   PetscFunctionReturn(0);
170020cf1dd8SToby Isaac }
170120cf1dd8SToby Isaac 
17026f905325SMatthew G. Knepley /*@C
17036f905325SMatthew G. Knepley   PetscDualSpaceGetSymmetries - Returns a description of the symmetries of this basis
17046f905325SMatthew G. Knepley 
17056f905325SMatthew G. Knepley   Not collective
17066f905325SMatthew G. Knepley 
17076f905325SMatthew G. Knepley   Input Parameter:
17086f905325SMatthew G. Knepley . sp - the PetscDualSpace object
17096f905325SMatthew G. Knepley 
17106f905325SMatthew G. Knepley   Output Parameters:
1711*b4457527SToby Isaac + perms - Permutations of the interior degrees of freedom, parameterized by the point orientation
1712*b4457527SToby Isaac - flips - Sign reversal of the interior degrees of freedom, parameterized by the point orientation
17136f905325SMatthew G. Knepley 
17146f905325SMatthew G. Knepley   Note: The permutation and flip arrays are organized in the following way
17156f905325SMatthew G. Knepley $ perms[p][ornt][dof # on point] = new local dof #
17166f905325SMatthew G. Knepley $ flips[p][ornt][dof # on point] = reversal or not
17176f905325SMatthew G. Knepley 
17186f905325SMatthew G. Knepley   Level: developer
17196f905325SMatthew G. Knepley 
17206f905325SMatthew G. Knepley @*/
17216f905325SMatthew G. Knepley PetscErrorCode PetscDualSpaceGetSymmetries(PetscDualSpace sp, const PetscInt ****perms, const PetscScalar ****flips)
17226f905325SMatthew G. Knepley {
17236f905325SMatthew G. Knepley   PetscErrorCode ierr;
17246f905325SMatthew G. Knepley 
17256f905325SMatthew G. Knepley   PetscFunctionBegin;
17266f905325SMatthew G. Knepley   PetscValidHeaderSpecific(sp,PETSCDUALSPACE_CLASSID,1);
17276f905325SMatthew G. Knepley   if (perms) {PetscValidPointer(perms,2); *perms = NULL;}
17286f905325SMatthew G. Knepley   if (flips) {PetscValidPointer(flips,3); *flips = NULL;}
17296f905325SMatthew G. Knepley   if (sp->ops->getsymmetries) {ierr = (sp->ops->getsymmetries)(sp,perms,flips);CHKERRQ(ierr);}
17306f905325SMatthew G. Knepley   PetscFunctionReturn(0);
17316f905325SMatthew G. Knepley }
17324bee2e38SMatthew G. Knepley 
17334bee2e38SMatthew G. Knepley /*@
1734*b4457527SToby Isaac   PetscDualSpaceGetFormDegree - Get the form degree k for the k-form the describes the pushforwards/pullbacks of this
1735*b4457527SToby Isaac   dual space's functionals.
1736*b4457527SToby Isaac 
1737*b4457527SToby Isaac   Input Parameter:
1738*b4457527SToby Isaac . dsp - The PetscDualSpace
1739*b4457527SToby Isaac 
1740*b4457527SToby Isaac   Output Parameter:
1741*b4457527SToby Isaac . k   - The *signed* degree k of the k.  If k >= 0, this means that the degrees of freedom are k-forms, and are stored
1742*b4457527SToby Isaac         in lexicographic order according to the basis of k-forms constructed from the wedge product of 1-forms.  So for example,
1743*b4457527SToby Isaac         the 1-form basis in 3-D is (dx, dy, dz), and the 2-form basis in 3-D is (dx wedge dy, dx wedge dz, dy wedge dz).
1744*b4457527SToby Isaac         If k < 0, this means that the degrees transform as k-forms, but are stored as (N-k) forms according to the
1745*b4457527SToby Isaac         Hodge star map.  So for example if k = -2 and N = 3, this means that the degrees of freedom transform as 2-forms
1746*b4457527SToby Isaac         but are stored as 1-forms.
1747*b4457527SToby Isaac 
1748*b4457527SToby Isaac   Level: developer
1749*b4457527SToby Isaac 
1750*b4457527SToby Isaac .seealso: PetscDTAltV, PetscDualSpacePullback(), PetscDualSpacePushforward(), PetscDualSpaceTransform(), PetscDualSpaceTransformType
1751*b4457527SToby Isaac @*/
1752*b4457527SToby Isaac PetscErrorCode PetscDualSpaceGetFormDegree(PetscDualSpace dsp, PetscInt *k)
1753*b4457527SToby Isaac {
1754*b4457527SToby Isaac   PetscFunctionBeginHot;
1755*b4457527SToby Isaac   PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1);
1756*b4457527SToby Isaac   PetscValidPointer(k, 2);
1757*b4457527SToby Isaac   *k = dsp->k;
1758*b4457527SToby Isaac   PetscFunctionReturn(0);
1759*b4457527SToby Isaac }
1760*b4457527SToby Isaac 
1761*b4457527SToby Isaac /*@
1762*b4457527SToby Isaac   PetscDualSpaceSetFormDegree - Set the form degree k for the k-form the describes the pushforwards/pullbacks of this
1763*b4457527SToby Isaac   dual space's functionals.
1764*b4457527SToby Isaac 
1765*b4457527SToby Isaac   Input Parameter:
1766*b4457527SToby Isaac + dsp - The PetscDualSpace
1767*b4457527SToby Isaac - k   - The *signed* degree k of the k.  If k >= 0, this means that the degrees of freedom are k-forms, and are stored
1768*b4457527SToby Isaac         in lexicographic order according to the basis of k-forms constructed from the wedge product of 1-forms.  So for example,
1769*b4457527SToby Isaac         the 1-form basis in 3-D is (dx, dy, dz), and the 2-form basis in 3-D is (dx wedge dy, dx wedge dz, dy wedge dz).
1770*b4457527SToby Isaac         If k < 0, this means that the degrees transform as k-forms, but are stored as (N-k) forms according to the
1771*b4457527SToby Isaac         Hodge star map.  So for example if k = -2 and N = 3, this means that the degrees of freedom transform as 2-forms
1772*b4457527SToby Isaac         but are stored as 1-forms.
1773*b4457527SToby Isaac 
1774*b4457527SToby Isaac   Level: developer
1775*b4457527SToby Isaac 
1776*b4457527SToby Isaac .seealso: PetscDTAltV, PetscDualSpacePullback(), PetscDualSpacePushforward(), PetscDualSpaceTransform(), PetscDualSpaceTransformType
1777*b4457527SToby Isaac @*/
1778*b4457527SToby Isaac PetscErrorCode PetscDualSpaceSetFormDegree(PetscDualSpace dsp, PetscInt k)
1779*b4457527SToby Isaac {
1780*b4457527SToby Isaac   PetscInt dim;
1781*b4457527SToby Isaac 
1782*b4457527SToby Isaac   PetscFunctionBeginHot;
1783*b4457527SToby Isaac   PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1);
1784*b4457527SToby Isaac   if (dsp->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dsp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change number of components after dualspace is set up");
1785*b4457527SToby Isaac   dim = dsp->dm->dim;
1786*b4457527SToby Isaac   if (k < -dim || k > dim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unsupported %D-form on %D-dimensional reference cell", PetscAbsInt(k), dim);
1787*b4457527SToby Isaac   dsp->k = k;
1788*b4457527SToby Isaac   PetscFunctionReturn(0);
1789*b4457527SToby Isaac }
1790*b4457527SToby Isaac 
1791*b4457527SToby Isaac /*@
17924bee2e38SMatthew G. Knepley   PetscDualSpaceGetDeRahm - Get the k-simplex associated with the functionals in this dual space
17934bee2e38SMatthew G. Knepley 
17944bee2e38SMatthew G. Knepley   Input Parameter:
17954bee2e38SMatthew G. Knepley . dsp - The PetscDualSpace
17964bee2e38SMatthew G. Knepley 
17974bee2e38SMatthew G. Knepley   Output Parameter:
17984bee2e38SMatthew G. Knepley . k   - The simplex dimension
17994bee2e38SMatthew G. Knepley 
1800a4ce7ad1SMatthew G. Knepley   Level: developer
18014bee2e38SMatthew G. Knepley 
18024bee2e38SMatthew G. Knepley   Note: Currently supported values are
18034bee2e38SMatthew G. Knepley $ 0: These are H_1 methods that only transform coordinates
18044bee2e38SMatthew G. Knepley $ 1: These are Hcurl methods that transform functions using the covariant Piola transform (COVARIANT_PIOLA_TRANSFORM)
18054bee2e38SMatthew G. Knepley $ 2: These are the same as 1
18064bee2e38SMatthew G. Knepley $ 3: These are Hdiv methods that transform functions using the contravariant Piola transform (CONTRAVARIANT_PIOLA_TRANSFORM)
18074bee2e38SMatthew G. Knepley 
18084bee2e38SMatthew G. Knepley .seealso: PetscDualSpacePullback(), PetscDualSpacePushforward(), PetscDualSpaceTransform(), PetscDualSpaceTransformType
18094bee2e38SMatthew G. Knepley @*/
18104bee2e38SMatthew G. Knepley PetscErrorCode PetscDualSpaceGetDeRahm(PetscDualSpace dsp, PetscInt *k)
18114bee2e38SMatthew G. Knepley {
1812*b4457527SToby Isaac   PetscInt dim;
1813*b4457527SToby Isaac 
18144bee2e38SMatthew G. Knepley   PetscFunctionBeginHot;
18154bee2e38SMatthew G. Knepley   PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1);
18164bee2e38SMatthew G. Knepley   PetscValidPointer(k, 2);
1817*b4457527SToby Isaac   dim = dsp->dm->dim;
1818*b4457527SToby Isaac   if (!dsp->k) *k = IDENTITY_TRANSFORM;
1819*b4457527SToby Isaac   else if (dsp->k == 1) *k = COVARIANT_PIOLA_TRANSFORM;
1820*b4457527SToby Isaac   else if (dsp->k == -(dim - 1)) *k = CONTRAVARIANT_PIOLA_TRANSFORM;
1821*b4457527SToby Isaac   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unsupported transformation");
18224bee2e38SMatthew G. Knepley   PetscFunctionReturn(0);
18234bee2e38SMatthew G. Knepley }
18244bee2e38SMatthew G. Knepley 
18254bee2e38SMatthew G. Knepley /*@C
18264bee2e38SMatthew G. Knepley   PetscDualSpaceTransform - Transform the function values
18274bee2e38SMatthew G. Knepley 
18284bee2e38SMatthew G. Knepley   Input Parameters:
18294bee2e38SMatthew G. Knepley + dsp       - The PetscDualSpace
18304bee2e38SMatthew G. Knepley . trans     - The type of transform
18314bee2e38SMatthew G. Knepley . isInverse - Flag to invert the transform
18324bee2e38SMatthew G. Knepley . fegeom    - The cell geometry
18334bee2e38SMatthew G. Knepley . Nv        - The number of function samples
18344bee2e38SMatthew G. Knepley . Nc        - The number of function components
18354bee2e38SMatthew G. Knepley - vals      - The function values
18364bee2e38SMatthew G. Knepley 
18374bee2e38SMatthew G. Knepley   Output Parameter:
18384bee2e38SMatthew G. Knepley . vals      - The transformed function values
18394bee2e38SMatthew G. Knepley 
1840a4ce7ad1SMatthew G. Knepley   Level: intermediate
18414bee2e38SMatthew G. Knepley 
1842625e0966SMatthew G. Knepley .seealso: PetscDualSpaceTransformGradient(), PetscDualSpacePullback(), PetscDualSpacePushforward(), PetscDualSpaceTransformType
18434bee2e38SMatthew G. Knepley @*/
18444bee2e38SMatthew G. Knepley PetscErrorCode PetscDualSpaceTransform(PetscDualSpace dsp, PetscDualSpaceTransformType trans, PetscBool isInverse, PetscFEGeom *fegeom, PetscInt Nv, PetscInt Nc, PetscScalar vals[])
18454bee2e38SMatthew G. Knepley {
1846*b4457527SToby Isaac   PetscReal Jstar[9] = {0};
1847*b4457527SToby Isaac   PetscInt dim, v, c, Nk;
1848*b4457527SToby Isaac   PetscErrorCode ierr;
18494bee2e38SMatthew G. Knepley 
18504bee2e38SMatthew G. Knepley   PetscFunctionBeginHot;
18514bee2e38SMatthew G. Knepley   PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1);
18524bee2e38SMatthew G. Knepley   PetscValidPointer(fegeom, 4);
18534bee2e38SMatthew G. Knepley   PetscValidPointer(vals, 7);
1854*b4457527SToby Isaac   /* TODO: not handling dimEmbed != dim right now */
18552ae266adSMatthew G. Knepley   dim = dsp->dm->dim;
1856*b4457527SToby Isaac   /* No change needed for 0-forms */
1857*b4457527SToby Isaac   if (!dsp->k) PetscFunctionReturn(0);
1858*b4457527SToby Isaac   ierr = PetscDTBinomialInt(dim, PetscAbsInt(dsp->k), &Nk);CHKERRQ(ierr);
1859*b4457527SToby Isaac   /* TODO: use fegeom->isAffine */
1860*b4457527SToby Isaac   ierr = PetscDTAltVPullbackMatrix(dim, dim, isInverse ? fegeom->J : fegeom->invJ, dsp->k, Jstar);CHKERRQ(ierr);
18614bee2e38SMatthew G. Knepley   for (v = 0; v < Nv; ++v) {
1862*b4457527SToby Isaac     switch (Nk) {
1863*b4457527SToby Isaac     case 1:
1864*b4457527SToby Isaac       for (c = 0; c < Nc; c++) vals[v*Nc + c] *= Jstar[0];
18654bee2e38SMatthew G. Knepley       break;
1866*b4457527SToby Isaac     case 2:
1867*b4457527SToby Isaac       for (c = 0; c < Nc; c += 2) DMPlex_Mult2DReal_Internal(Jstar, 1, &vals[v*Nc + c], &vals[v*Nc + c]);
18684bee2e38SMatthew G. Knepley       break;
1869*b4457527SToby Isaac     case 3:
1870*b4457527SToby Isaac       for (c = 0; c < Nc; c += 3) DMPlex_Mult3DReal_Internal(Jstar, 1, &vals[v*Nc + c], &vals[v*Nc + c]);
1871*b4457527SToby Isaac       break;
1872*b4457527SToby Isaac     default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported form size %D for transformation", Nk);
1873*b4457527SToby Isaac     }
18744bee2e38SMatthew G. Knepley   }
18754bee2e38SMatthew G. Knepley   PetscFunctionReturn(0);
18764bee2e38SMatthew G. Knepley }
1877*b4457527SToby Isaac 
18784bee2e38SMatthew G. Knepley /*@C
18794bee2e38SMatthew G. Knepley   PetscDualSpaceTransformGradient - Transform the function gradient values
18804bee2e38SMatthew G. Knepley 
18814bee2e38SMatthew G. Knepley   Input Parameters:
18824bee2e38SMatthew G. Knepley + dsp       - The PetscDualSpace
18834bee2e38SMatthew G. Knepley . trans     - The type of transform
18844bee2e38SMatthew G. Knepley . isInverse - Flag to invert the transform
18854bee2e38SMatthew G. Knepley . fegeom    - The cell geometry
18864bee2e38SMatthew G. Knepley . Nv        - The number of function gradient samples
18874bee2e38SMatthew G. Knepley . Nc        - The number of function components
18884bee2e38SMatthew G. Knepley - vals      - The function gradient values
18894bee2e38SMatthew G. Knepley 
18904bee2e38SMatthew G. Knepley   Output Parameter:
18914bee2e38SMatthew G. Knepley . vals      - The transformed function values
18924bee2e38SMatthew G. Knepley 
1893a4ce7ad1SMatthew G. Knepley   Level: intermediate
18944bee2e38SMatthew G. Knepley 
1895625e0966SMatthew G. Knepley .seealso: PetscDualSpaceTransform(), PetscDualSpacePullback(), PetscDualSpacePushforward(), PetscDualSpaceTransformType
18964bee2e38SMatthew G. Knepley @*/
18974bee2e38SMatthew G. Knepley PetscErrorCode PetscDualSpaceTransformGradient(PetscDualSpace dsp, PetscDualSpaceTransformType trans, PetscBool isInverse, PetscFEGeom *fegeom, PetscInt Nv, PetscInt Nc, PetscScalar vals[])
18984bee2e38SMatthew G. Knepley {
18994bee2e38SMatthew G. Knepley   PetscInt dim, v, c, d;
19004bee2e38SMatthew G. Knepley 
19014bee2e38SMatthew G. Knepley   PetscFunctionBeginHot;
19024bee2e38SMatthew G. Knepley   PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1);
19034bee2e38SMatthew G. Knepley   PetscValidPointer(fegeom, 4);
19044bee2e38SMatthew G. Knepley   PetscValidPointer(vals, 7);
19052ae266adSMatthew G. Knepley   dim = dsp->dm->dim;
19064bee2e38SMatthew G. Knepley   /* Transform gradient */
19074bee2e38SMatthew G. Knepley   for (v = 0; v < Nv; ++v) {
19084bee2e38SMatthew G. Knepley     for (c = 0; c < Nc; ++c) {
19094bee2e38SMatthew G. Knepley       switch (dim)
19104bee2e38SMatthew G. Knepley       {
1911100a78e1SStefano Zampini         case 1: vals[(v*Nc+c)*dim] *= fegeom->invJ[0];break;
19126142fa51SMatthew G. Knepley         case 2: DMPlex_MultTranspose2DReal_Internal(fegeom->invJ, 1, &vals[(v*Nc+c)*dim], &vals[(v*Nc+c)*dim]);break;
19136142fa51SMatthew G. Knepley         case 3: DMPlex_MultTranspose3DReal_Internal(fegeom->invJ, 1, &vals[(v*Nc+c)*dim], &vals[(v*Nc+c)*dim]);break;
19144bee2e38SMatthew G. Knepley         default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim);
19154bee2e38SMatthew G. Knepley       }
19164bee2e38SMatthew G. Knepley     }
19174bee2e38SMatthew G. Knepley   }
19184bee2e38SMatthew G. Knepley   /* Assume its a vector, otherwise assume its a bunch of scalars */
19194bee2e38SMatthew G. Knepley   if (Nc == 1 || Nc != dim) PetscFunctionReturn(0);
19204bee2e38SMatthew G. Knepley   switch (trans) {
19214bee2e38SMatthew G. Knepley     case IDENTITY_TRANSFORM: break;
19224bee2e38SMatthew G. Knepley     case COVARIANT_PIOLA_TRANSFORM: /* Covariant Piola mapping $\sigma^*(F) = J^{-T} F \circ \phi^{-1)$ */
19234bee2e38SMatthew G. Knepley     if (isInverse) {
19244bee2e38SMatthew G. Knepley       for (v = 0; v < Nv; ++v) {
19254bee2e38SMatthew G. Knepley         for (d = 0; d < dim; ++d) {
19264bee2e38SMatthew G. Knepley           switch (dim)
19274bee2e38SMatthew G. Knepley           {
19286142fa51SMatthew G. Knepley             case 2: DMPlex_MultTranspose2DReal_Internal(fegeom->J, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
19296142fa51SMatthew G. Knepley             case 3: DMPlex_MultTranspose3DReal_Internal(fegeom->J, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
19304bee2e38SMatthew G. Knepley             default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim);
19314bee2e38SMatthew G. Knepley           }
19324bee2e38SMatthew G. Knepley         }
19334bee2e38SMatthew G. Knepley       }
19344bee2e38SMatthew G. Knepley     } else {
19354bee2e38SMatthew G. Knepley       for (v = 0; v < Nv; ++v) {
19364bee2e38SMatthew G. Knepley         for (d = 0; d < dim; ++d) {
19374bee2e38SMatthew G. Knepley           switch (dim)
19384bee2e38SMatthew G. Knepley           {
19396142fa51SMatthew G. Knepley             case 2: DMPlex_MultTranspose2DReal_Internal(fegeom->invJ, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
19406142fa51SMatthew G. Knepley             case 3: DMPlex_MultTranspose3DReal_Internal(fegeom->invJ, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
19414bee2e38SMatthew G. Knepley             default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim);
19424bee2e38SMatthew G. Knepley           }
19434bee2e38SMatthew G. Knepley         }
19444bee2e38SMatthew G. Knepley       }
19454bee2e38SMatthew G. Knepley     }
19464bee2e38SMatthew G. Knepley     break;
19474bee2e38SMatthew G. Knepley     case CONTRAVARIANT_PIOLA_TRANSFORM: /* Contravariant Piola mapping $\sigma^*(F) = \frac{1}{|\det J|} J F \circ \phi^{-1}$ */
19484bee2e38SMatthew G. Knepley     if (isInverse) {
19494bee2e38SMatthew G. Knepley       for (v = 0; v < Nv; ++v) {
19504bee2e38SMatthew G. Knepley         for (d = 0; d < dim; ++d) {
19514bee2e38SMatthew G. Knepley           switch (dim)
19524bee2e38SMatthew G. Knepley           {
19536142fa51SMatthew G. Knepley             case 2: DMPlex_Mult2DReal_Internal(fegeom->invJ, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
19546142fa51SMatthew G. Knepley             case 3: DMPlex_Mult3DReal_Internal(fegeom->invJ, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
19554bee2e38SMatthew G. Knepley             default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim);
19564bee2e38SMatthew G. Knepley           }
19574bee2e38SMatthew G. Knepley           for (c = 0; c < Nc; ++c) vals[(v*Nc+c)*dim+d] *= fegeom->detJ[0];
19584bee2e38SMatthew G. Knepley         }
19594bee2e38SMatthew G. Knepley       }
19604bee2e38SMatthew G. Knepley     } else {
19614bee2e38SMatthew G. Knepley       for (v = 0; v < Nv; ++v) {
19624bee2e38SMatthew G. Knepley         for (d = 0; d < dim; ++d) {
19634bee2e38SMatthew G. Knepley           switch (dim)
19644bee2e38SMatthew G. Knepley           {
19656142fa51SMatthew G. Knepley             case 2: DMPlex_Mult2DReal_Internal(fegeom->J, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
19666142fa51SMatthew G. Knepley             case 3: DMPlex_Mult3DReal_Internal(fegeom->J, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
19674bee2e38SMatthew G. Knepley             default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim);
19684bee2e38SMatthew G. Knepley           }
19694bee2e38SMatthew G. Knepley           for (c = 0; c < Nc; ++c) vals[(v*Nc+c)*dim+d] /= fegeom->detJ[0];
19704bee2e38SMatthew G. Knepley         }
19714bee2e38SMatthew G. Knepley       }
19724bee2e38SMatthew G. Knepley     }
19734bee2e38SMatthew G. Knepley     break;
19744bee2e38SMatthew G. Knepley   }
19754bee2e38SMatthew G. Knepley   PetscFunctionReturn(0);
19764bee2e38SMatthew G. Knepley }
19774bee2e38SMatthew G. Knepley 
19784bee2e38SMatthew G. Knepley /*@C
19794bee2e38SMatthew G. Knepley   PetscDualSpacePullback - Transform the given functional so that it operates on real space, rather than the reference element. Operationally, this means that we map the function evaluations depending on continuity requirements of our finite element method.
19804bee2e38SMatthew G. Knepley 
19814bee2e38SMatthew G. Knepley   Input Parameters:
19824bee2e38SMatthew G. Knepley + dsp        - The PetscDualSpace
19834bee2e38SMatthew G. Knepley . fegeom     - The geometry for this cell
19844bee2e38SMatthew G. Knepley . Nq         - The number of function samples
19854bee2e38SMatthew G. Knepley . Nc         - The number of function components
19864bee2e38SMatthew G. Knepley - pointEval  - The function values
19874bee2e38SMatthew G. Knepley 
19884bee2e38SMatthew G. Knepley   Output Parameter:
19894bee2e38SMatthew G. Knepley . pointEval  - The transformed function values
19904bee2e38SMatthew G. Knepley 
19914bee2e38SMatthew G. Knepley   Level: advanced
19924bee2e38SMatthew G. Knepley 
19934bee2e38SMatthew G. Knepley   Note: Functions transform in a complementary way (pushforward) to functionals, so that the scalar product is invariant. The type of transform is dependent on the associated k-simplex from the DeRahm complex.
19944bee2e38SMatthew G. Knepley 
19954bee2e38SMatthew G. Knepley .seealso: PetscDualSpacePushforward(), PetscDualSpaceTransform(), PetscDualSpaceGetDeRahm()
19964bee2e38SMatthew G. Knepley @*/
1997*b4457527SToby Isaac PetscErrorCode PetscDualSpacePullback(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, const PetscScalar realEval[], PetscScalar refEval[])
19984bee2e38SMatthew G. Knepley {
19994bee2e38SMatthew G. Knepley   PetscDualSpaceTransformType trans;
2000*b4457527SToby Isaac   PetscInt                    k;
20014bee2e38SMatthew G. Knepley   PetscErrorCode              ierr;
20024bee2e38SMatthew G. Knepley 
20034bee2e38SMatthew G. Knepley   PetscFunctionBeginHot;
20044bee2e38SMatthew G. Knepley   PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1);
20054bee2e38SMatthew G. Knepley   PetscValidPointer(fegeom, 2);
2006*b4457527SToby Isaac   PetscValidPointer(refEval, 6);
20074bee2e38SMatthew G. Knepley   /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k.
20084bee2e38SMatthew G. Knepley      This determines their transformation properties. */
2009*b4457527SToby Isaac   ierr = PetscDualSpaceGetDeRahm(dsp, &k);CHKERRQ(ierr);
2010*b4457527SToby Isaac   switch (k)
20114bee2e38SMatthew G. Knepley   {
20124bee2e38SMatthew G. Knepley     case 0: /* H^1 point evaluations */
20134bee2e38SMatthew G. Knepley     trans = IDENTITY_TRANSFORM;break;
20144bee2e38SMatthew G. Knepley     case 1: /* Hcurl preserves tangential edge traces  */
20154bee2e38SMatthew G. Knepley     trans = COVARIANT_PIOLA_TRANSFORM;break;
2016*b4457527SToby Isaac     case 2:
20174bee2e38SMatthew G. Knepley     case 3: /* Hdiv preserve normal traces */
20184bee2e38SMatthew G. Knepley     trans = CONTRAVARIANT_PIOLA_TRANSFORM;break;
2019*b4457527SToby Isaac     default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %D for transformation", k);
20204bee2e38SMatthew G. Knepley   }
2021*b4457527SToby Isaac   ierr = PetscDualSpaceTransform(dsp, trans, PETSC_TRUE, fegeom, Nq, Nc, refEval);CHKERRQ(ierr);
20224bee2e38SMatthew G. Knepley   PetscFunctionReturn(0);
20234bee2e38SMatthew G. Knepley }
20244bee2e38SMatthew G. Knepley 
20254bee2e38SMatthew G. Knepley /*@C
20264bee2e38SMatthew G. Knepley   PetscDualSpacePushforward - Transform the given function so that it operates on real space, rather than the reference element. Operationally, this means that we map the function evaluations depending on continuity requirements of our finite element method.
20274bee2e38SMatthew G. Knepley 
20284bee2e38SMatthew G. Knepley   Input Parameters:
20294bee2e38SMatthew G. Knepley + dsp        - The PetscDualSpace
20304bee2e38SMatthew G. Knepley . fegeom     - The geometry for this cell
20314bee2e38SMatthew G. Knepley . Nq         - The number of function samples
20324bee2e38SMatthew G. Knepley . Nc         - The number of function components
20334bee2e38SMatthew G. Knepley - pointEval  - The function values
20344bee2e38SMatthew G. Knepley 
20354bee2e38SMatthew G. Knepley   Output Parameter:
20364bee2e38SMatthew G. Knepley . pointEval  - The transformed function values
20374bee2e38SMatthew G. Knepley 
20384bee2e38SMatthew G. Knepley   Level: advanced
20394bee2e38SMatthew G. Knepley 
20404bee2e38SMatthew G. Knepley   Note: Functionals transform in a complementary way (pullback) to functions, so that the scalar product is invariant. The type of transform is dependent on the associated k-simplex from the DeRahm complex.
20414bee2e38SMatthew G. Knepley 
20424bee2e38SMatthew G. Knepley .seealso: PetscDualSpacePullback(), PetscDualSpaceTransform(), PetscDualSpaceGetDeRahm()
20434bee2e38SMatthew G. Knepley @*/
2044*b4457527SToby Isaac PetscErrorCode PetscDualSpacePushforward(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, const PetscScalar refEval[], PetscScalar realEval[])
20454bee2e38SMatthew G. Knepley {
20464bee2e38SMatthew G. Knepley   PetscDualSpaceTransformType trans;
2047*b4457527SToby Isaac   PetscInt                    k;
20484bee2e38SMatthew G. Knepley   PetscErrorCode              ierr;
20494bee2e38SMatthew G. Knepley 
20504bee2e38SMatthew G. Knepley   PetscFunctionBeginHot;
20514bee2e38SMatthew G. Knepley   PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1);
20524bee2e38SMatthew G. Knepley   PetscValidPointer(fegeom, 2);
2053*b4457527SToby Isaac   PetscValidPointer(realEval, 6);
20544bee2e38SMatthew G. Knepley   /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k.
20554bee2e38SMatthew G. Knepley      This determines their transformation properties. */
2056*b4457527SToby Isaac   ierr = PetscDualSpaceGetDeRahm(dsp, &k);CHKERRQ(ierr);
2057*b4457527SToby Isaac   switch (k)
20584bee2e38SMatthew G. Knepley   {
20594bee2e38SMatthew G. Knepley     case 0: /* H^1 point evaluations */
20604bee2e38SMatthew G. Knepley     trans = IDENTITY_TRANSFORM;break;
20614bee2e38SMatthew G. Knepley     case 1: /* Hcurl preserves tangential edge traces  */
20624bee2e38SMatthew G. Knepley     trans = COVARIANT_PIOLA_TRANSFORM;break;
2063*b4457527SToby Isaac     case 2:
20644bee2e38SMatthew G. Knepley     case 3: /* Hdiv preserve normal traces */
20654bee2e38SMatthew G. Knepley     trans = CONTRAVARIANT_PIOLA_TRANSFORM;break;
2066*b4457527SToby Isaac     default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %D for transformation", k);
20674bee2e38SMatthew G. Knepley   }
2068*b4457527SToby Isaac   ierr = PetscDualSpaceTransform(dsp, trans, PETSC_FALSE, fegeom, Nq, Nc, realEval);CHKERRQ(ierr);
20694bee2e38SMatthew G. Knepley   PetscFunctionReturn(0);
20704bee2e38SMatthew G. Knepley }
20714bee2e38SMatthew G. Knepley 
20724bee2e38SMatthew G. Knepley /*@C
20734bee2e38SMatthew G. Knepley   PetscDualSpacePushforwardGradient - Transform the given function gradient so that it operates on real space, rather than the reference element. Operationally, this means that we map the function evaluations depending on continuity requirements of our finite element method.
20744bee2e38SMatthew G. Knepley 
20754bee2e38SMatthew G. Knepley   Input Parameters:
20764bee2e38SMatthew G. Knepley + dsp        - The PetscDualSpace
20774bee2e38SMatthew G. Knepley . fegeom     - The geometry for this cell
20784bee2e38SMatthew G. Knepley . Nq         - The number of function gradient samples
20794bee2e38SMatthew G. Knepley . Nc         - The number of function components
20804bee2e38SMatthew G. Knepley - pointEval  - The function gradient values
20814bee2e38SMatthew G. Knepley 
20824bee2e38SMatthew G. Knepley   Output Parameter:
20834bee2e38SMatthew G. Knepley . pointEval  - The transformed function gradient values
20844bee2e38SMatthew G. Knepley 
20854bee2e38SMatthew G. Knepley   Level: advanced
20864bee2e38SMatthew G. Knepley 
20874bee2e38SMatthew G. Knepley   Note: Functionals transform in a complementary way (pullback) to functions, so that the scalar product is invariant. The type of transform is dependent on the associated k-simplex from the DeRahm complex.
20884bee2e38SMatthew G. Knepley 
20894bee2e38SMatthew G. Knepley .seealso: PetscDualSpacePushforward(), PPetscDualSpacePullback(), PetscDualSpaceTransform(), PetscDualSpaceGetDeRahm()
2090dc0529c6SBarry Smith @*/
2091*b4457527SToby Isaac PetscErrorCode PetscDualSpacePushforwardGradient(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, const PetscScalar refEval[], PetscScalar realEval[])
20924bee2e38SMatthew G. Knepley {
20934bee2e38SMatthew G. Knepley   PetscDualSpaceTransformType trans;
2094*b4457527SToby Isaac   PetscInt                    k;
20954bee2e38SMatthew G. Knepley   PetscErrorCode              ierr;
20964bee2e38SMatthew G. Knepley 
20974bee2e38SMatthew G. Knepley   PetscFunctionBeginHot;
20984bee2e38SMatthew G. Knepley   PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1);
20994bee2e38SMatthew G. Knepley   PetscValidPointer(fegeom, 2);
2100*b4457527SToby Isaac   PetscValidPointer(realEval, 6);
21014bee2e38SMatthew G. Knepley   /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k.
21024bee2e38SMatthew G. Knepley      This determines their transformation properties. */
2103*b4457527SToby Isaac   ierr = PetscDualSpaceGetDeRahm(dsp, &k);CHKERRQ(ierr);
2104*b4457527SToby Isaac   switch (k)
21054bee2e38SMatthew G. Knepley   {
21064bee2e38SMatthew G. Knepley     case 0: /* H^1 point evaluations */
21074bee2e38SMatthew G. Knepley     trans = IDENTITY_TRANSFORM;break;
21084bee2e38SMatthew G. Knepley     case 1: /* Hcurl preserves tangential edge traces  */
21094bee2e38SMatthew G. Knepley     trans = COVARIANT_PIOLA_TRANSFORM;break;
2110*b4457527SToby Isaac     case 2:
21114bee2e38SMatthew G. Knepley     case 3: /* Hdiv preserve normal traces */
21124bee2e38SMatthew G. Knepley     trans = CONTRAVARIANT_PIOLA_TRANSFORM;break;
2113*b4457527SToby Isaac     default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %D for transformation", k);
21144bee2e38SMatthew G. Knepley   }
2115*b4457527SToby Isaac   ierr = PetscDualSpaceTransformGradient(dsp, trans, PETSC_FALSE, fegeom, Nq, Nc, realEval);CHKERRQ(ierr);
21164bee2e38SMatthew G. Knepley   PetscFunctionReturn(0);
21174bee2e38SMatthew G. Knepley }
2118