xref: /petsc/src/dm/dt/dualspace/interface/dualspace.c (revision a4ce7ad1fcbd2ddf3c416448e98dcc9c4d96ab2c)
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.
146f905325SMatthew G. Knepley                                                      e.g. for len == 2 and max == 2, this will return, in order, {0,0}, {1,0}, {2,0}, {0,1}, {1,1}, {2,0}.
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);
195221d6281SMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(v, "Dual space with %D components, size %D\n", sp->Nc, pdim);CHKERRQ(ierr);
196221d6281SMatthew G. Knepley   if (sp->ops->view) {ierr = (*sp->ops->view)(sp, v);CHKERRQ(ierr);}
197221d6281SMatthew G. Knepley   ierr = PetscViewerGetFormat(v, &format);CHKERRQ(ierr);
198221d6281SMatthew G. Knepley   if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
199221d6281SMatthew G. Knepley     ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr);
200221d6281SMatthew G. Knepley     for (f = 0; f < pdim; ++f) {
201221d6281SMatthew G. Knepley       ierr = PetscViewerASCIIPrintf(v, "Dual basis vector %D\n", f);CHKERRQ(ierr);
202221d6281SMatthew G. Knepley       ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr);
203221d6281SMatthew G. Knepley       ierr = PetscQuadratureView(sp->functional[f], v);CHKERRQ(ierr);
204221d6281SMatthew G. Knepley       ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr);
205221d6281SMatthew G. Knepley     }
206221d6281SMatthew G. Knepley     ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr);
207221d6281SMatthew G. Knepley   }
208221d6281SMatthew G. Knepley   ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr);
209221d6281SMatthew G. Knepley   PetscFunctionReturn(0);
210221d6281SMatthew G. Knepley }
211221d6281SMatthew G. Knepley 
21220cf1dd8SToby Isaac /*@
21320cf1dd8SToby Isaac   PetscDualSpaceView - Views a PetscDualSpace
21420cf1dd8SToby Isaac 
215d083f849SBarry Smith   Collective on sp
21620cf1dd8SToby Isaac 
21720cf1dd8SToby Isaac   Input Parameter:
21820cf1dd8SToby Isaac + sp - the PetscDualSpace object to view
21920cf1dd8SToby Isaac - v  - the viewer
22020cf1dd8SToby Isaac 
221*a4ce7ad1SMatthew G. Knepley   Level: beginner
22220cf1dd8SToby Isaac 
22320cf1dd8SToby Isaac .seealso PetscDualSpaceDestroy()
22420cf1dd8SToby Isaac @*/
22520cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceView(PetscDualSpace sp, PetscViewer v)
22620cf1dd8SToby Isaac {
227d9bac1caSLisandro Dalcin   PetscBool      iascii;
22820cf1dd8SToby Isaac   PetscErrorCode ierr;
22920cf1dd8SToby Isaac 
23020cf1dd8SToby Isaac   PetscFunctionBegin;
23120cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
232d9bac1caSLisandro Dalcin   if (v) PetscValidHeaderSpecific(v, PETSC_VIEWER_CLASSID, 2);
23320cf1dd8SToby Isaac   if (!v) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) sp), &v);CHKERRQ(ierr);}
234d9bac1caSLisandro Dalcin   ierr = PetscObjectTypeCompare((PetscObject) v, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
235221d6281SMatthew G. Knepley   if (iascii) {ierr = PetscDualSpaceView_ASCII(sp, v);CHKERRQ(ierr);}
23620cf1dd8SToby Isaac   PetscFunctionReturn(0);
23720cf1dd8SToby Isaac }
23820cf1dd8SToby Isaac 
23920cf1dd8SToby Isaac /*@
24020cf1dd8SToby Isaac   PetscDualSpaceSetFromOptions - sets parameters in a PetscDualSpace from the options database
24120cf1dd8SToby Isaac 
242d083f849SBarry Smith   Collective on sp
24320cf1dd8SToby Isaac 
24420cf1dd8SToby Isaac   Input Parameter:
24520cf1dd8SToby Isaac . sp - the PetscDualSpace object to set options for
24620cf1dd8SToby Isaac 
24720cf1dd8SToby Isaac   Options Database:
2487be5e748SToby Isaac . -petscspace_degree the approximation order of the space
24920cf1dd8SToby Isaac 
250*a4ce7ad1SMatthew G. Knepley   Level: intermediate
25120cf1dd8SToby Isaac 
25220cf1dd8SToby Isaac .seealso PetscDualSpaceView()
25320cf1dd8SToby Isaac @*/
25420cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetFromOptions(PetscDualSpace sp)
25520cf1dd8SToby Isaac {
256063ee4adSMatthew G. Knepley   PetscDualSpaceReferenceCell refCell = PETSCDUALSPACE_REFCELL_SIMPLEX;
257063ee4adSMatthew G. Knepley   PetscInt                    refDim  = 0;
258063ee4adSMatthew G. Knepley   PetscBool                   flg;
25920cf1dd8SToby Isaac   const char                 *defaultType;
26020cf1dd8SToby Isaac   char                        name[256];
26120cf1dd8SToby Isaac   PetscErrorCode              ierr;
26220cf1dd8SToby Isaac 
26320cf1dd8SToby Isaac   PetscFunctionBegin;
26420cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
26520cf1dd8SToby Isaac   if (!((PetscObject) sp)->type_name) {
26620cf1dd8SToby Isaac     defaultType = PETSCDUALSPACELAGRANGE;
26720cf1dd8SToby Isaac   } else {
26820cf1dd8SToby Isaac     defaultType = ((PetscObject) sp)->type_name;
26920cf1dd8SToby Isaac   }
27020cf1dd8SToby Isaac   if (!PetscSpaceRegisterAllCalled) {ierr = PetscSpaceRegisterAll();CHKERRQ(ierr);}
27120cf1dd8SToby Isaac 
27220cf1dd8SToby Isaac   ierr = PetscObjectOptionsBegin((PetscObject) sp);CHKERRQ(ierr);
27320cf1dd8SToby Isaac   ierr = PetscOptionsFList("-petscdualspace_type", "Dual space", "PetscDualSpaceSetType", PetscDualSpaceList, defaultType, name, 256, &flg);CHKERRQ(ierr);
27420cf1dd8SToby Isaac   if (flg) {
27520cf1dd8SToby Isaac     ierr = PetscDualSpaceSetType(sp, name);CHKERRQ(ierr);
27620cf1dd8SToby Isaac   } else if (!((PetscObject) sp)->type_name) {
27720cf1dd8SToby Isaac     ierr = PetscDualSpaceSetType(sp, defaultType);CHKERRQ(ierr);
27820cf1dd8SToby Isaac   }
2795a856986SBarry Smith   ierr = PetscOptionsBoundedInt("-petscdualspace_degree", "The approximation order", "PetscDualSpaceSetOrder", sp->order, &sp->order, NULL,0);CHKERRQ(ierr);
2805a856986SBarry Smith   ierr = PetscOptionsBoundedInt("-petscdualspace_components", "The number of components", "PetscDualSpaceSetNumComponents", sp->Nc, &sp->Nc, NULL,1);CHKERRQ(ierr);
28120cf1dd8SToby Isaac   if (sp->ops->setfromoptions) {
28220cf1dd8SToby Isaac     ierr = (*sp->ops->setfromoptions)(PetscOptionsObject,sp);CHKERRQ(ierr);
28320cf1dd8SToby Isaac   }
2845a856986SBarry Smith   ierr = PetscOptionsBoundedInt("-petscdualspace_refdim", "The spatial dimension of the reference cell", "PetscDualSpaceSetReferenceCell", refDim, &refDim, NULL,0);CHKERRQ(ierr);
285063ee4adSMatthew G. Knepley   ierr = PetscOptionsEnum("-petscdualspace_refcell", "Reference cell", "PetscDualSpaceSetReferenceCell", PetscDualSpaceReferenceCells, (PetscEnum) refCell, (PetscEnum *) &refCell, &flg);CHKERRQ(ierr);
286063ee4adSMatthew G. Knepley   if (flg) {
287063ee4adSMatthew G. Knepley     DM K;
288063ee4adSMatthew G. Knepley 
289063ee4adSMatthew G. Knepley     if (!refDim) SETERRQ(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_INCOMP, "Reference cell specified without a dimension. Use -petscdualspace_refdim.");
290063ee4adSMatthew G. Knepley     ierr = PetscDualSpaceCreateReferenceCell(sp, refDim, refCell == PETSCDUALSPACE_REFCELL_SIMPLEX ? PETSC_TRUE : PETSC_FALSE, &K);CHKERRQ(ierr);
291063ee4adSMatthew G. Knepley     ierr = PetscDualSpaceSetDM(sp, K);CHKERRQ(ierr);
292063ee4adSMatthew G. Knepley     ierr = DMDestroy(&K);CHKERRQ(ierr);
293063ee4adSMatthew G. Knepley   }
294063ee4adSMatthew G. Knepley 
29520cf1dd8SToby Isaac   /* process any options handlers added with PetscObjectAddOptionsHandler() */
29620cf1dd8SToby Isaac   ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) sp);CHKERRQ(ierr);
29720cf1dd8SToby Isaac   ierr = PetscOptionsEnd();CHKERRQ(ierr);
298063ee4adSMatthew G. Knepley   sp->setfromoptionscalled = PETSC_TRUE;
29920cf1dd8SToby Isaac   PetscFunctionReturn(0);
30020cf1dd8SToby Isaac }
30120cf1dd8SToby Isaac 
30220cf1dd8SToby Isaac /*@
30320cf1dd8SToby Isaac   PetscDualSpaceSetUp - Construct a basis for the PetscDualSpace
30420cf1dd8SToby Isaac 
305d083f849SBarry Smith   Collective on sp
30620cf1dd8SToby Isaac 
30720cf1dd8SToby Isaac   Input Parameter:
30820cf1dd8SToby Isaac . sp - the PetscDualSpace object to setup
30920cf1dd8SToby Isaac 
310*a4ce7ad1SMatthew G. Knepley   Level: intermediate
31120cf1dd8SToby Isaac 
31220cf1dd8SToby Isaac .seealso PetscDualSpaceView(), PetscDualSpaceDestroy()
31320cf1dd8SToby Isaac @*/
31420cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetUp(PetscDualSpace sp)
31520cf1dd8SToby Isaac {
31620cf1dd8SToby Isaac   PetscErrorCode ierr;
31720cf1dd8SToby Isaac 
31820cf1dd8SToby Isaac   PetscFunctionBegin;
31920cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
32020cf1dd8SToby Isaac   if (sp->setupcalled) PetscFunctionReturn(0);
32120cf1dd8SToby Isaac   sp->setupcalled = PETSC_TRUE;
32220cf1dd8SToby Isaac   if (sp->ops->setup) {ierr = (*sp->ops->setup)(sp);CHKERRQ(ierr);}
323063ee4adSMatthew G. Knepley   if (sp->setfromoptionscalled) {ierr = PetscDualSpaceViewFromOptions(sp, NULL, "-petscdualspace_view");CHKERRQ(ierr);}
32420cf1dd8SToby Isaac   PetscFunctionReturn(0);
32520cf1dd8SToby Isaac }
32620cf1dd8SToby Isaac 
32720cf1dd8SToby Isaac /*@
32820cf1dd8SToby Isaac   PetscDualSpaceDestroy - Destroys a PetscDualSpace object
32920cf1dd8SToby Isaac 
330d083f849SBarry Smith   Collective on sp
33120cf1dd8SToby Isaac 
33220cf1dd8SToby Isaac   Input Parameter:
33320cf1dd8SToby Isaac . sp - the PetscDualSpace object to destroy
33420cf1dd8SToby Isaac 
335*a4ce7ad1SMatthew G. Knepley   Level: beginner
33620cf1dd8SToby Isaac 
33720cf1dd8SToby Isaac .seealso PetscDualSpaceView()
33820cf1dd8SToby Isaac @*/
33920cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceDestroy(PetscDualSpace *sp)
34020cf1dd8SToby Isaac {
34120cf1dd8SToby Isaac   PetscInt       dim, f;
34220cf1dd8SToby Isaac   PetscErrorCode ierr;
34320cf1dd8SToby Isaac 
34420cf1dd8SToby Isaac   PetscFunctionBegin;
34520cf1dd8SToby Isaac   if (!*sp) PetscFunctionReturn(0);
34620cf1dd8SToby Isaac   PetscValidHeaderSpecific((*sp), PETSCDUALSPACE_CLASSID, 1);
34720cf1dd8SToby Isaac 
34820cf1dd8SToby Isaac   if (--((PetscObject)(*sp))->refct > 0) {*sp = 0; PetscFunctionReturn(0);}
34920cf1dd8SToby Isaac   ((PetscObject) (*sp))->refct = 0;
35020cf1dd8SToby Isaac 
35120cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDimension(*sp, &dim);CHKERRQ(ierr);
35220cf1dd8SToby Isaac   for (f = 0; f < dim; ++f) {
35320cf1dd8SToby Isaac     ierr = PetscQuadratureDestroy(&(*sp)->functional[f]);CHKERRQ(ierr);
35420cf1dd8SToby Isaac   }
35520cf1dd8SToby Isaac   ierr = PetscFree((*sp)->functional);CHKERRQ(ierr);
35620cf1dd8SToby Isaac   ierr = PetscQuadratureDestroy(&(*sp)->allPoints);CHKERRQ(ierr);
35720cf1dd8SToby Isaac   ierr = DMDestroy(&(*sp)->dm);CHKERRQ(ierr);
35820cf1dd8SToby Isaac 
35920cf1dd8SToby Isaac   if ((*sp)->ops->destroy) {ierr = (*(*sp)->ops->destroy)(*sp);CHKERRQ(ierr);}
36020cf1dd8SToby Isaac   ierr = PetscHeaderDestroy(sp);CHKERRQ(ierr);
36120cf1dd8SToby Isaac   PetscFunctionReturn(0);
36220cf1dd8SToby Isaac }
36320cf1dd8SToby Isaac 
36420cf1dd8SToby Isaac /*@
36520cf1dd8SToby Isaac   PetscDualSpaceCreate - Creates an empty PetscDualSpace object. The type can then be set with PetscDualSpaceSetType().
36620cf1dd8SToby Isaac 
367d083f849SBarry Smith   Collective
36820cf1dd8SToby Isaac 
36920cf1dd8SToby Isaac   Input Parameter:
37020cf1dd8SToby Isaac . comm - The communicator for the PetscDualSpace object
37120cf1dd8SToby Isaac 
37220cf1dd8SToby Isaac   Output Parameter:
37320cf1dd8SToby Isaac . sp - The PetscDualSpace object
37420cf1dd8SToby Isaac 
37520cf1dd8SToby Isaac   Level: beginner
37620cf1dd8SToby Isaac 
37720cf1dd8SToby Isaac .seealso: PetscDualSpaceSetType(), PETSCDUALSPACELAGRANGE
37820cf1dd8SToby Isaac @*/
37920cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceCreate(MPI_Comm comm, PetscDualSpace *sp)
38020cf1dd8SToby Isaac {
38120cf1dd8SToby Isaac   PetscDualSpace s;
38220cf1dd8SToby Isaac   PetscErrorCode ierr;
38320cf1dd8SToby Isaac 
38420cf1dd8SToby Isaac   PetscFunctionBegin;
38520cf1dd8SToby Isaac   PetscValidPointer(sp, 2);
38620cf1dd8SToby Isaac   ierr = PetscCitationsRegister(FECitation,&FEcite);CHKERRQ(ierr);
38720cf1dd8SToby Isaac   *sp  = NULL;
38820cf1dd8SToby Isaac   ierr = PetscFEInitializePackage();CHKERRQ(ierr);
38920cf1dd8SToby Isaac 
39020cf1dd8SToby Isaac   ierr = PetscHeaderCreate(s, PETSCDUALSPACE_CLASSID, "PetscDualSpace", "Dual Space", "PetscDualSpace", comm, PetscDualSpaceDestroy, PetscDualSpaceView);CHKERRQ(ierr);
39120cf1dd8SToby Isaac 
39220cf1dd8SToby Isaac   s->order = 0;
39320cf1dd8SToby Isaac   s->Nc    = 1;
3944bee2e38SMatthew G. Knepley   s->k     = 0;
39520cf1dd8SToby Isaac   s->setupcalled = PETSC_FALSE;
39620cf1dd8SToby Isaac 
39720cf1dd8SToby Isaac   *sp = s;
39820cf1dd8SToby Isaac   PetscFunctionReturn(0);
39920cf1dd8SToby Isaac }
40020cf1dd8SToby Isaac 
40120cf1dd8SToby Isaac /*@
40220cf1dd8SToby Isaac   PetscDualSpaceDuplicate - Creates a duplicate PetscDualSpace object, however it is not setup.
40320cf1dd8SToby Isaac 
404d083f849SBarry Smith   Collective on sp
40520cf1dd8SToby Isaac 
40620cf1dd8SToby Isaac   Input Parameter:
40720cf1dd8SToby Isaac . sp - The original PetscDualSpace
40820cf1dd8SToby Isaac 
40920cf1dd8SToby Isaac   Output Parameter:
41020cf1dd8SToby Isaac . spNew - The duplicate PetscDualSpace
41120cf1dd8SToby Isaac 
41220cf1dd8SToby Isaac   Level: beginner
41320cf1dd8SToby Isaac 
41420cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate(), PetscDualSpaceSetType()
41520cf1dd8SToby Isaac @*/
41620cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceDuplicate(PetscDualSpace sp, PetscDualSpace *spNew)
41720cf1dd8SToby Isaac {
41820cf1dd8SToby Isaac   PetscErrorCode ierr;
41920cf1dd8SToby Isaac 
42020cf1dd8SToby Isaac   PetscFunctionBegin;
42120cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
42220cf1dd8SToby Isaac   PetscValidPointer(spNew, 2);
42320cf1dd8SToby Isaac   ierr = (*sp->ops->duplicate)(sp, spNew);CHKERRQ(ierr);
42420cf1dd8SToby Isaac   PetscFunctionReturn(0);
42520cf1dd8SToby Isaac }
42620cf1dd8SToby Isaac 
42720cf1dd8SToby Isaac /*@
42820cf1dd8SToby Isaac   PetscDualSpaceGetDM - Get the DM representing the reference cell
42920cf1dd8SToby Isaac 
43020cf1dd8SToby Isaac   Not collective
43120cf1dd8SToby Isaac 
43220cf1dd8SToby Isaac   Input Parameter:
43320cf1dd8SToby Isaac . sp - The PetscDualSpace
43420cf1dd8SToby Isaac 
43520cf1dd8SToby Isaac   Output Parameter:
43620cf1dd8SToby Isaac . dm - The reference cell
43720cf1dd8SToby Isaac 
43820cf1dd8SToby Isaac   Level: intermediate
43920cf1dd8SToby Isaac 
44020cf1dd8SToby Isaac .seealso: PetscDualSpaceSetDM(), PetscDualSpaceCreate()
44120cf1dd8SToby Isaac @*/
44220cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetDM(PetscDualSpace sp, DM *dm)
44320cf1dd8SToby Isaac {
44420cf1dd8SToby Isaac   PetscFunctionBegin;
44520cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
44620cf1dd8SToby Isaac   PetscValidPointer(dm, 2);
44720cf1dd8SToby Isaac   *dm = sp->dm;
44820cf1dd8SToby Isaac   PetscFunctionReturn(0);
44920cf1dd8SToby Isaac }
45020cf1dd8SToby Isaac 
45120cf1dd8SToby Isaac /*@
45220cf1dd8SToby Isaac   PetscDualSpaceSetDM - Get the DM representing the reference cell
45320cf1dd8SToby Isaac 
45420cf1dd8SToby Isaac   Not collective
45520cf1dd8SToby Isaac 
45620cf1dd8SToby Isaac   Input Parameters:
45720cf1dd8SToby Isaac + sp - The PetscDualSpace
45820cf1dd8SToby Isaac - dm - The reference cell
45920cf1dd8SToby Isaac 
46020cf1dd8SToby Isaac   Level: intermediate
46120cf1dd8SToby Isaac 
46220cf1dd8SToby Isaac .seealso: PetscDualSpaceGetDM(), PetscDualSpaceCreate()
46320cf1dd8SToby Isaac @*/
46420cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetDM(PetscDualSpace sp, DM dm)
46520cf1dd8SToby Isaac {
46620cf1dd8SToby Isaac   PetscErrorCode ierr;
46720cf1dd8SToby Isaac 
46820cf1dd8SToby Isaac   PetscFunctionBegin;
46920cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
47020cf1dd8SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
47120cf1dd8SToby Isaac   ierr = DMDestroy(&sp->dm);CHKERRQ(ierr);
47220cf1dd8SToby Isaac   ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);
47320cf1dd8SToby Isaac   sp->dm = dm;
47420cf1dd8SToby Isaac   PetscFunctionReturn(0);
47520cf1dd8SToby Isaac }
47620cf1dd8SToby Isaac 
47720cf1dd8SToby Isaac /*@
47820cf1dd8SToby Isaac   PetscDualSpaceGetOrder - Get the order of the dual space
47920cf1dd8SToby Isaac 
48020cf1dd8SToby Isaac   Not collective
48120cf1dd8SToby Isaac 
48220cf1dd8SToby Isaac   Input Parameter:
48320cf1dd8SToby Isaac . sp - The PetscDualSpace
48420cf1dd8SToby Isaac 
48520cf1dd8SToby Isaac   Output Parameter:
48620cf1dd8SToby Isaac . order - The order
48720cf1dd8SToby Isaac 
48820cf1dd8SToby Isaac   Level: intermediate
48920cf1dd8SToby Isaac 
49020cf1dd8SToby Isaac .seealso: PetscDualSpaceSetOrder(), PetscDualSpaceCreate()
49120cf1dd8SToby Isaac @*/
49220cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetOrder(PetscDualSpace sp, PetscInt *order)
49320cf1dd8SToby Isaac {
49420cf1dd8SToby Isaac   PetscFunctionBegin;
49520cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
49620cf1dd8SToby Isaac   PetscValidPointer(order, 2);
49720cf1dd8SToby Isaac   *order = sp->order;
49820cf1dd8SToby Isaac   PetscFunctionReturn(0);
49920cf1dd8SToby Isaac }
50020cf1dd8SToby Isaac 
50120cf1dd8SToby Isaac /*@
50220cf1dd8SToby Isaac   PetscDualSpaceSetOrder - Set the order of the dual space
50320cf1dd8SToby Isaac 
50420cf1dd8SToby Isaac   Not collective
50520cf1dd8SToby Isaac 
50620cf1dd8SToby Isaac   Input Parameters:
50720cf1dd8SToby Isaac + sp - The PetscDualSpace
50820cf1dd8SToby Isaac - order - The order
50920cf1dd8SToby Isaac 
51020cf1dd8SToby Isaac   Level: intermediate
51120cf1dd8SToby Isaac 
51220cf1dd8SToby Isaac .seealso: PetscDualSpaceGetOrder(), PetscDualSpaceCreate()
51320cf1dd8SToby Isaac @*/
51420cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetOrder(PetscDualSpace sp, PetscInt order)
51520cf1dd8SToby Isaac {
51620cf1dd8SToby Isaac   PetscFunctionBegin;
51720cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
51820cf1dd8SToby Isaac   sp->order = order;
51920cf1dd8SToby Isaac   PetscFunctionReturn(0);
52020cf1dd8SToby Isaac }
52120cf1dd8SToby Isaac 
52220cf1dd8SToby Isaac /*@
52320cf1dd8SToby Isaac   PetscDualSpaceGetNumComponents - Return the number of components for this space
52420cf1dd8SToby Isaac 
52520cf1dd8SToby Isaac   Input Parameter:
52620cf1dd8SToby Isaac . sp - The PetscDualSpace
52720cf1dd8SToby Isaac 
52820cf1dd8SToby Isaac   Output Parameter:
52920cf1dd8SToby Isaac . Nc - The number of components
53020cf1dd8SToby Isaac 
53120cf1dd8SToby Isaac   Note: A vector space, for example, will have d components, where d is the spatial dimension
53220cf1dd8SToby Isaac 
53320cf1dd8SToby Isaac   Level: intermediate
53420cf1dd8SToby Isaac 
53520cf1dd8SToby Isaac .seealso: PetscDualSpaceSetNumComponents(), PetscDualSpaceGetDimension(), PetscDualSpaceCreate(), PetscDualSpace
53620cf1dd8SToby Isaac @*/
53720cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetNumComponents(PetscDualSpace sp, PetscInt *Nc)
53820cf1dd8SToby Isaac {
53920cf1dd8SToby Isaac   PetscFunctionBegin;
54020cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
54120cf1dd8SToby Isaac   PetscValidPointer(Nc, 2);
54220cf1dd8SToby Isaac   *Nc = sp->Nc;
54320cf1dd8SToby Isaac   PetscFunctionReturn(0);
54420cf1dd8SToby Isaac }
54520cf1dd8SToby Isaac 
54620cf1dd8SToby Isaac /*@
54720cf1dd8SToby Isaac   PetscDualSpaceSetNumComponents - Set the number of components for this space
54820cf1dd8SToby Isaac 
54920cf1dd8SToby Isaac   Input Parameters:
55020cf1dd8SToby Isaac + sp - The PetscDualSpace
55120cf1dd8SToby Isaac - order - The number of components
55220cf1dd8SToby Isaac 
55320cf1dd8SToby Isaac   Level: intermediate
55420cf1dd8SToby Isaac 
55520cf1dd8SToby Isaac .seealso: PetscDualSpaceGetNumComponents(), PetscDualSpaceCreate(), PetscDualSpace
55620cf1dd8SToby Isaac @*/
55720cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetNumComponents(PetscDualSpace sp, PetscInt Nc)
55820cf1dd8SToby Isaac {
55920cf1dd8SToby Isaac   PetscFunctionBegin;
56020cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
56120cf1dd8SToby Isaac   sp->Nc = Nc;
56220cf1dd8SToby Isaac   PetscFunctionReturn(0);
56320cf1dd8SToby Isaac }
56420cf1dd8SToby Isaac 
56520cf1dd8SToby Isaac /*@
56620cf1dd8SToby Isaac   PetscDualSpaceGetFunctional - Get the i-th basis functional in the dual space
56720cf1dd8SToby Isaac 
56820cf1dd8SToby Isaac   Not collective
56920cf1dd8SToby Isaac 
57020cf1dd8SToby Isaac   Input Parameters:
57120cf1dd8SToby Isaac + sp - The PetscDualSpace
57220cf1dd8SToby Isaac - i  - The basis number
57320cf1dd8SToby Isaac 
57420cf1dd8SToby Isaac   Output Parameter:
57520cf1dd8SToby Isaac . functional - The basis functional
57620cf1dd8SToby Isaac 
57720cf1dd8SToby Isaac   Level: intermediate
57820cf1dd8SToby Isaac 
57920cf1dd8SToby Isaac .seealso: PetscDualSpaceGetDimension(), PetscDualSpaceCreate()
58020cf1dd8SToby Isaac @*/
58120cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetFunctional(PetscDualSpace sp, PetscInt i, PetscQuadrature *functional)
58220cf1dd8SToby Isaac {
58320cf1dd8SToby Isaac   PetscInt       dim;
58420cf1dd8SToby Isaac   PetscErrorCode ierr;
58520cf1dd8SToby Isaac 
58620cf1dd8SToby Isaac   PetscFunctionBegin;
58720cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
58820cf1dd8SToby Isaac   PetscValidPointer(functional, 3);
58920cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDimension(sp, &dim);CHKERRQ(ierr);
59020cf1dd8SToby Isaac   if ((i < 0) || (i >= dim)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Functional index %d must be in [0, %d)", i, dim);
59120cf1dd8SToby Isaac   *functional = sp->functional[i];
59220cf1dd8SToby Isaac   PetscFunctionReturn(0);
59320cf1dd8SToby Isaac }
59420cf1dd8SToby Isaac 
59520cf1dd8SToby Isaac /*@
59620cf1dd8SToby Isaac   PetscDualSpaceGetDimension - Get the dimension of the dual space, i.e. the number of basis functionals
59720cf1dd8SToby Isaac 
59820cf1dd8SToby Isaac   Not collective
59920cf1dd8SToby Isaac 
60020cf1dd8SToby Isaac   Input Parameter:
60120cf1dd8SToby Isaac . sp - The PetscDualSpace
60220cf1dd8SToby Isaac 
60320cf1dd8SToby Isaac   Output Parameter:
60420cf1dd8SToby Isaac . dim - The dimension
60520cf1dd8SToby Isaac 
60620cf1dd8SToby Isaac   Level: intermediate
60720cf1dd8SToby Isaac 
60820cf1dd8SToby Isaac .seealso: PetscDualSpaceGetFunctional(), PetscDualSpaceCreate()
60920cf1dd8SToby Isaac @*/
61020cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetDimension(PetscDualSpace sp, PetscInt *dim)
61120cf1dd8SToby Isaac {
61220cf1dd8SToby Isaac   PetscErrorCode ierr;
61320cf1dd8SToby Isaac 
61420cf1dd8SToby Isaac   PetscFunctionBegin;
61520cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
61620cf1dd8SToby Isaac   PetscValidPointer(dim, 2);
61720cf1dd8SToby Isaac   *dim = 0;
61820cf1dd8SToby Isaac   if (sp->ops->getdimension) {ierr = (*sp->ops->getdimension)(sp, dim);CHKERRQ(ierr);}
61920cf1dd8SToby Isaac   PetscFunctionReturn(0);
62020cf1dd8SToby Isaac }
62120cf1dd8SToby Isaac 
62220cf1dd8SToby Isaac /*@C
62320cf1dd8SToby Isaac   PetscDualSpaceGetNumDof - Get the number of degrees of freedom for each spatial (topological) dimension
62420cf1dd8SToby Isaac 
62520cf1dd8SToby Isaac   Not collective
62620cf1dd8SToby Isaac 
62720cf1dd8SToby Isaac   Input Parameter:
62820cf1dd8SToby Isaac . sp - The PetscDualSpace
62920cf1dd8SToby Isaac 
63020cf1dd8SToby Isaac   Output Parameter:
63120cf1dd8SToby Isaac . numDof - An array of length dim+1 which holds the number of dofs for each dimension
63220cf1dd8SToby Isaac 
63320cf1dd8SToby Isaac   Level: intermediate
63420cf1dd8SToby Isaac 
63520cf1dd8SToby Isaac .seealso: PetscDualSpaceGetFunctional(), PetscDualSpaceCreate()
63620cf1dd8SToby Isaac @*/
63720cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetNumDof(PetscDualSpace sp, const PetscInt **numDof)
63820cf1dd8SToby Isaac {
63920cf1dd8SToby Isaac   PetscErrorCode ierr;
64020cf1dd8SToby Isaac 
64120cf1dd8SToby Isaac   PetscFunctionBegin;
64220cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
64320cf1dd8SToby Isaac   PetscValidPointer(numDof, 2);
64420cf1dd8SToby Isaac   ierr = (*sp->ops->getnumdof)(sp, numDof);CHKERRQ(ierr);
64520cf1dd8SToby Isaac   if (!*numDof) SETERRQ(PetscObjectComm((PetscObject) sp), PETSC_ERR_LIB, "Empty numDof[] returned from dual space implementation");
64620cf1dd8SToby Isaac   PetscFunctionReturn(0);
64720cf1dd8SToby Isaac }
64820cf1dd8SToby Isaac 
649*a4ce7ad1SMatthew G. Knepley /*@
650*a4ce7ad1SMatthew G. Knepley   PetscDualSpaceCreateSection - Create a PetscSection over the reference cell with the layout from this space
651*a4ce7ad1SMatthew G. Knepley 
652*a4ce7ad1SMatthew G. Knepley   Collective on sp
653*a4ce7ad1SMatthew G. Knepley 
654*a4ce7ad1SMatthew G. Knepley   Input Parameters:
655*a4ce7ad1SMatthew G. Knepley + sp      - The PetscDualSpace
656*a4ce7ad1SMatthew G. Knepley 
657*a4ce7ad1SMatthew G. Knepley   Output Parameter:
658*a4ce7ad1SMatthew G. Knepley . section - The section
659*a4ce7ad1SMatthew G. Knepley 
660*a4ce7ad1SMatthew G. Knepley   Level: advanced
661*a4ce7ad1SMatthew G. Knepley 
662*a4ce7ad1SMatthew G. Knepley .seealso: PetscDualSpaceCreate(), DMPLEX
663*a4ce7ad1SMatthew G. Knepley @*/
66420cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceCreateSection(PetscDualSpace sp, PetscSection *section)
66520cf1dd8SToby Isaac {
66620cf1dd8SToby Isaac   DM             dm;
66720cf1dd8SToby Isaac   PetscInt       pStart, pEnd, depth, h, offset;
66820cf1dd8SToby Isaac   const PetscInt *numDof;
66920cf1dd8SToby Isaac   PetscErrorCode ierr;
67020cf1dd8SToby Isaac 
67120cf1dd8SToby Isaac   PetscFunctionBegin;
67220cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDM(sp,&dm);CHKERRQ(ierr);
67320cf1dd8SToby Isaac   ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr);
67420cf1dd8SToby Isaac   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)sp),section);CHKERRQ(ierr);
67520cf1dd8SToby Isaac   ierr = PetscSectionSetChart(*section,pStart,pEnd);CHKERRQ(ierr);
67620cf1dd8SToby Isaac   ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr);
67720cf1dd8SToby Isaac   ierr = PetscDualSpaceGetNumDof(sp,&numDof);CHKERRQ(ierr);
67820cf1dd8SToby Isaac   for (h = 0; h <= depth; h++) {
67920cf1dd8SToby Isaac     PetscInt hStart, hEnd, p, dof;
68020cf1dd8SToby Isaac 
68120cf1dd8SToby Isaac     ierr = DMPlexGetHeightStratum(dm,h,&hStart,&hEnd);CHKERRQ(ierr);
68220cf1dd8SToby Isaac     dof = numDof[depth - h];
68320cf1dd8SToby Isaac     for (p = hStart; p < hEnd; p++) {
68420cf1dd8SToby Isaac       ierr = PetscSectionSetDof(*section,p,dof);CHKERRQ(ierr);
68520cf1dd8SToby Isaac     }
68620cf1dd8SToby Isaac   }
68720cf1dd8SToby Isaac   ierr = PetscSectionSetUp(*section);CHKERRQ(ierr);
68820cf1dd8SToby Isaac   for (h = 0, offset = 0; h <= depth; h++) {
68920cf1dd8SToby Isaac     PetscInt hStart, hEnd, p, dof;
69020cf1dd8SToby Isaac 
69120cf1dd8SToby Isaac     ierr = DMPlexGetHeightStratum(dm,h,&hStart,&hEnd);CHKERRQ(ierr);
69220cf1dd8SToby Isaac     dof = numDof[depth - h];
69320cf1dd8SToby Isaac     for (p = hStart; p < hEnd; p++) {
69420cf1dd8SToby Isaac       ierr = PetscSectionGetDof(*section,p,&dof);CHKERRQ(ierr);
69520cf1dd8SToby Isaac       ierr = PetscSectionSetOffset(*section,p,offset);CHKERRQ(ierr);
69620cf1dd8SToby Isaac       offset += dof;
69720cf1dd8SToby Isaac     }
69820cf1dd8SToby Isaac   }
69920cf1dd8SToby Isaac   PetscFunctionReturn(0);
70020cf1dd8SToby Isaac }
70120cf1dd8SToby Isaac 
70220cf1dd8SToby Isaac /*@
70320cf1dd8SToby Isaac   PetscDualSpaceCreateReferenceCell - Create a DMPLEX with the appropriate FEM reference cell
70420cf1dd8SToby Isaac 
705d083f849SBarry Smith   Collective on sp
70620cf1dd8SToby Isaac 
70720cf1dd8SToby Isaac   Input Parameters:
70820cf1dd8SToby Isaac + sp      - The PetscDualSpace
70920cf1dd8SToby Isaac . dim     - The spatial dimension
71020cf1dd8SToby Isaac - simplex - Flag for simplex, otherwise use a tensor-product cell
71120cf1dd8SToby Isaac 
71220cf1dd8SToby Isaac   Output Parameter:
71320cf1dd8SToby Isaac . refdm - The reference cell
71420cf1dd8SToby Isaac 
715*a4ce7ad1SMatthew G. Knepley   Level: intermediate
71620cf1dd8SToby Isaac 
71720cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate(), DMPLEX
71820cf1dd8SToby Isaac @*/
71920cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceCreateReferenceCell(PetscDualSpace sp, PetscInt dim, PetscBool simplex, DM *refdm)
72020cf1dd8SToby Isaac {
72120cf1dd8SToby Isaac   PetscErrorCode ierr;
72220cf1dd8SToby Isaac 
72320cf1dd8SToby Isaac   PetscFunctionBeginUser;
72420cf1dd8SToby Isaac   ierr = DMPlexCreateReferenceCell(PetscObjectComm((PetscObject) sp), dim, simplex, refdm);CHKERRQ(ierr);
72520cf1dd8SToby Isaac   PetscFunctionReturn(0);
72620cf1dd8SToby Isaac }
72720cf1dd8SToby Isaac 
72820cf1dd8SToby Isaac /*@C
72920cf1dd8SToby Isaac   PetscDualSpaceApply - Apply a functional from the dual space basis to an input function
73020cf1dd8SToby Isaac 
73120cf1dd8SToby Isaac   Input Parameters:
73220cf1dd8SToby Isaac + sp      - The PetscDualSpace object
73320cf1dd8SToby Isaac . f       - The basis functional index
73420cf1dd8SToby Isaac . time    - The time
73520cf1dd8SToby 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)
73620cf1dd8SToby Isaac . numComp - The number of components for the function
73720cf1dd8SToby Isaac . func    - The input function
73820cf1dd8SToby Isaac - ctx     - A context for the function
73920cf1dd8SToby Isaac 
74020cf1dd8SToby Isaac   Output Parameter:
74120cf1dd8SToby Isaac . value   - numComp output values
74220cf1dd8SToby Isaac 
74320cf1dd8SToby Isaac   Note: The calling sequence for the callback func is given by:
74420cf1dd8SToby Isaac 
74520cf1dd8SToby Isaac $ func(PetscInt dim, PetscReal time, const PetscReal x[],
74620cf1dd8SToby Isaac $      PetscInt numComponents, PetscScalar values[], void *ctx)
74720cf1dd8SToby Isaac 
748*a4ce7ad1SMatthew G. Knepley   Level: beginner
74920cf1dd8SToby Isaac 
75020cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate()
75120cf1dd8SToby Isaac @*/
75220cf1dd8SToby 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)
75320cf1dd8SToby Isaac {
75420cf1dd8SToby Isaac   PetscErrorCode ierr;
75520cf1dd8SToby Isaac 
75620cf1dd8SToby Isaac   PetscFunctionBegin;
75720cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
75820cf1dd8SToby Isaac   PetscValidPointer(cgeom, 4);
75920cf1dd8SToby Isaac   PetscValidPointer(value, 8);
76020cf1dd8SToby Isaac   ierr = (*sp->ops->apply)(sp, f, time, cgeom, numComp, func, ctx, value);CHKERRQ(ierr);
76120cf1dd8SToby Isaac   PetscFunctionReturn(0);
76220cf1dd8SToby Isaac }
76320cf1dd8SToby Isaac 
76420cf1dd8SToby Isaac /*@C
76520cf1dd8SToby Isaac   PetscDualSpaceApplyAll - Apply all functionals from the dual space basis to the result of an evaluation at the points returned by PetscDualSpaceGetAllPoints()
76620cf1dd8SToby Isaac 
76720cf1dd8SToby Isaac   Input Parameters:
76820cf1dd8SToby Isaac + sp        - The PetscDualSpace object
76920cf1dd8SToby Isaac - pointEval - Evaluation at the points returned by PetscDualSpaceGetAllPoints()
77020cf1dd8SToby Isaac 
77120cf1dd8SToby Isaac   Output Parameter:
77220cf1dd8SToby Isaac . spValue   - The values of all dual space functionals
77320cf1dd8SToby Isaac 
774*a4ce7ad1SMatthew G. Knepley   Level: beginner
77520cf1dd8SToby Isaac 
77620cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate()
77720cf1dd8SToby Isaac @*/
77820cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceApplyAll(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue)
77920cf1dd8SToby Isaac {
78020cf1dd8SToby Isaac   PetscErrorCode ierr;
78120cf1dd8SToby Isaac 
78220cf1dd8SToby Isaac   PetscFunctionBegin;
78320cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
78420cf1dd8SToby Isaac   ierr = (*sp->ops->applyall)(sp, pointEval, spValue);CHKERRQ(ierr);
78520cf1dd8SToby Isaac   PetscFunctionReturn(0);
78620cf1dd8SToby Isaac }
78720cf1dd8SToby Isaac 
78820cf1dd8SToby Isaac /*@C
78920cf1dd8SToby Isaac   PetscDualSpaceApplyDefault - Apply a functional from the dual space basis to an input function by assuming a point evaluation functional.
79020cf1dd8SToby Isaac 
79120cf1dd8SToby Isaac   Input Parameters:
79220cf1dd8SToby Isaac + sp    - The PetscDualSpace object
79320cf1dd8SToby Isaac . f     - The basis functional index
79420cf1dd8SToby Isaac . time  - The time
79520cf1dd8SToby Isaac . cgeom - A context with geometric information for this cell, we use v0 (the initial vertex) and J (the Jacobian)
79620cf1dd8SToby Isaac . Nc    - The number of components for the function
79720cf1dd8SToby Isaac . func  - The input function
79820cf1dd8SToby Isaac - ctx   - A context for the function
79920cf1dd8SToby Isaac 
80020cf1dd8SToby Isaac   Output Parameter:
80120cf1dd8SToby Isaac . value   - The output value
80220cf1dd8SToby Isaac 
80320cf1dd8SToby Isaac   Note: The calling sequence for the callback func is given by:
80420cf1dd8SToby Isaac 
80520cf1dd8SToby Isaac $ func(PetscInt dim, PetscReal time, const PetscReal x[],
80620cf1dd8SToby Isaac $      PetscInt numComponents, PetscScalar values[], void *ctx)
80720cf1dd8SToby Isaac 
80820cf1dd8SToby Isaac and the idea is to evaluate the functional as an integral
80920cf1dd8SToby Isaac 
81020cf1dd8SToby Isaac $ n(f) = int dx n(x) . f(x)
81120cf1dd8SToby Isaac 
81220cf1dd8SToby Isaac where both n and f have Nc components.
81320cf1dd8SToby Isaac 
814*a4ce7ad1SMatthew G. Knepley   Level: beginner
81520cf1dd8SToby Isaac 
81620cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate()
81720cf1dd8SToby Isaac @*/
81820cf1dd8SToby 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)
81920cf1dd8SToby Isaac {
82020cf1dd8SToby Isaac   DM               dm;
82120cf1dd8SToby Isaac   PetscQuadrature  n;
82220cf1dd8SToby Isaac   const PetscReal *points, *weights;
82320cf1dd8SToby Isaac   PetscReal        x[3];
82420cf1dd8SToby Isaac   PetscScalar     *val;
82520cf1dd8SToby Isaac   PetscInt         dim, dE, qNc, c, Nq, q;
82620cf1dd8SToby Isaac   PetscBool        isAffine;
82720cf1dd8SToby Isaac   PetscErrorCode   ierr;
82820cf1dd8SToby Isaac 
82920cf1dd8SToby Isaac   PetscFunctionBegin;
83020cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
83120cf1dd8SToby Isaac   PetscValidPointer(value, 5);
83220cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr);
83320cf1dd8SToby Isaac   ierr = PetscDualSpaceGetFunctional(sp, f, &n);CHKERRQ(ierr);
83420cf1dd8SToby Isaac   ierr = PetscQuadratureGetData(n, &dim, &qNc, &Nq, &points, &weights);CHKERRQ(ierr);
83520cf1dd8SToby 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);
83620cf1dd8SToby Isaac   if (qNc != Nc) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature components %D != function components %D", qNc, Nc);
83720cf1dd8SToby Isaac   ierr = DMGetWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr);
83820cf1dd8SToby Isaac   *value = 0.0;
83920cf1dd8SToby Isaac   isAffine = cgeom->isAffine;
84020cf1dd8SToby Isaac   dE = cgeom->dimEmbed;
84120cf1dd8SToby Isaac   for (q = 0; q < Nq; ++q) {
84220cf1dd8SToby Isaac     if (isAffine) {
84320cf1dd8SToby Isaac       CoordinatesRefToReal(dE, cgeom->dim, cgeom->xi, cgeom->v, cgeom->J, &points[q*dim], x);
84420cf1dd8SToby Isaac       ierr = (*func)(dE, time, x, Nc, val, ctx);CHKERRQ(ierr);
84520cf1dd8SToby Isaac     } else {
84620cf1dd8SToby Isaac       ierr = (*func)(dE, time, &cgeom->v[dE*q], Nc, val, ctx);CHKERRQ(ierr);
84720cf1dd8SToby Isaac     }
84820cf1dd8SToby Isaac     for (c = 0; c < Nc; ++c) {
84920cf1dd8SToby Isaac       *value += val[c]*weights[q*Nc+c];
85020cf1dd8SToby Isaac     }
85120cf1dd8SToby Isaac   }
85220cf1dd8SToby Isaac   ierr = DMRestoreWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr);
85320cf1dd8SToby Isaac   PetscFunctionReturn(0);
85420cf1dd8SToby Isaac }
85520cf1dd8SToby Isaac 
85620cf1dd8SToby Isaac /*@C
85720cf1dd8SToby Isaac   PetscDualSpaceApplyAllDefault - Apply all functionals from the dual space basis to the result of an evaluation at the points returned by PetscDualSpaceGetAllPoints()
85820cf1dd8SToby Isaac 
85920cf1dd8SToby Isaac   Input Parameters:
86020cf1dd8SToby Isaac + sp        - The PetscDualSpace object
86120cf1dd8SToby Isaac - pointEval - Evaluation at the points returned by PetscDualSpaceGetAllPoints()
86220cf1dd8SToby Isaac 
86320cf1dd8SToby Isaac   Output Parameter:
86420cf1dd8SToby Isaac . spValue   - The values of all dual space functionals
86520cf1dd8SToby Isaac 
866*a4ce7ad1SMatthew G. Knepley   Level: beginner
86720cf1dd8SToby Isaac 
86820cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate()
86920cf1dd8SToby Isaac @*/
87020cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceApplyAllDefault(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue)
87120cf1dd8SToby Isaac {
87220cf1dd8SToby Isaac   PetscQuadrature  n;
87320cf1dd8SToby Isaac   const PetscReal *points, *weights;
87420cf1dd8SToby Isaac   PetscInt         qNc, c, Nq, q, f, spdim, Nc;
87520cf1dd8SToby Isaac   PetscInt         offset;
87620cf1dd8SToby Isaac   PetscErrorCode   ierr;
87720cf1dd8SToby Isaac 
87820cf1dd8SToby Isaac   PetscFunctionBegin;
87920cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
88020cf1dd8SToby Isaac   PetscValidScalarPointer(pointEval, 2);
88120cf1dd8SToby Isaac   PetscValidScalarPointer(spValue, 5);
88220cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDimension(sp, &spdim);CHKERRQ(ierr);
88320cf1dd8SToby Isaac   ierr = PetscDualSpaceGetNumComponents(sp, &Nc);CHKERRQ(ierr);
88420cf1dd8SToby Isaac   for (f = 0, offset = 0; f < spdim; f++) {
88520cf1dd8SToby Isaac     ierr = PetscDualSpaceGetFunctional(sp, f, &n);CHKERRQ(ierr);
88620cf1dd8SToby Isaac     ierr = PetscQuadratureGetData(n, NULL, &qNc, &Nq, &points, &weights);CHKERRQ(ierr);
88720cf1dd8SToby Isaac     if (qNc != Nc) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature components %D != function components %D", qNc, Nc);
88820cf1dd8SToby Isaac     spValue[f] = 0.0;
88920cf1dd8SToby Isaac     for (q = 0; q < Nq; ++q) {
89020cf1dd8SToby Isaac       for (c = 0; c < Nc; ++c) {
89120cf1dd8SToby Isaac         spValue[f] += pointEval[offset++]*weights[q*Nc+c];
89220cf1dd8SToby Isaac       }
89320cf1dd8SToby Isaac     }
89420cf1dd8SToby Isaac   }
89520cf1dd8SToby Isaac   PetscFunctionReturn(0);
89620cf1dd8SToby Isaac }
89720cf1dd8SToby Isaac 
898*a4ce7ad1SMatthew G. Knepley /*@
899*a4ce7ad1SMatthew G. Knepley   PetscDualSpaceGetAllPoints - Get all quadrature points from this space
900*a4ce7ad1SMatthew G. Knepley 
901*a4ce7ad1SMatthew G. Knepley   Input Parameter:
902*a4ce7ad1SMatthew G. Knepley . sp - The dualspace
903*a4ce7ad1SMatthew G. Knepley 
904*a4ce7ad1SMatthew G. Knepley   Output Parameter:
905*a4ce7ad1SMatthew G. Knepley . allPoints - A PetscQuadrature object containing all evaluation points
906*a4ce7ad1SMatthew G. Knepley 
907*a4ce7ad1SMatthew G. Knepley   Level: advanced
908*a4ce7ad1SMatthew G. Knepley 
909*a4ce7ad1SMatthew G. Knepley .seealso: PetscDualSpaceCreate()
910*a4ce7ad1SMatthew G. Knepley @*/
91120cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetAllPoints(PetscDualSpace sp, PetscQuadrature *allPoints)
91220cf1dd8SToby Isaac {
91320cf1dd8SToby Isaac   PetscErrorCode ierr;
91420cf1dd8SToby Isaac 
91520cf1dd8SToby Isaac   PetscFunctionBegin;
91620cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
91720cf1dd8SToby Isaac   PetscValidPointer(allPoints,2);
91820cf1dd8SToby Isaac   if (!sp->allPoints && sp->ops->createallpoints) {
91920cf1dd8SToby Isaac     ierr = (*sp->ops->createallpoints)(sp,&sp->allPoints);CHKERRQ(ierr);
92020cf1dd8SToby Isaac   }
92120cf1dd8SToby Isaac   *allPoints = sp->allPoints;
92220cf1dd8SToby Isaac   PetscFunctionReturn(0);
92320cf1dd8SToby Isaac }
92420cf1dd8SToby Isaac 
925*a4ce7ad1SMatthew G. Knepley /*@
926*a4ce7ad1SMatthew G. Knepley   PetscDualSpaceCreateAllPointsDefault - Create all evaluation points by examining functionals
927*a4ce7ad1SMatthew G. Knepley 
928*a4ce7ad1SMatthew G. Knepley   Input Parameter:
929*a4ce7ad1SMatthew G. Knepley . sp - The dualspace
930*a4ce7ad1SMatthew G. Knepley 
931*a4ce7ad1SMatthew G. Knepley   Output Parameter:
932*a4ce7ad1SMatthew G. Knepley . allPoints - A PetscQuadrature object containing all evaluation points
933*a4ce7ad1SMatthew G. Knepley 
934*a4ce7ad1SMatthew G. Knepley   Level: advanced
935*a4ce7ad1SMatthew G. Knepley 
936*a4ce7ad1SMatthew G. Knepley .seealso: PetscDualSpaceCreate()
937*a4ce7ad1SMatthew G. Knepley @*/
93820cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceCreateAllPointsDefault(PetscDualSpace sp, PetscQuadrature *allPoints)
93920cf1dd8SToby Isaac {
94020cf1dd8SToby Isaac   PetscInt        spdim;
94120cf1dd8SToby Isaac   PetscInt        numPoints, offset;
94220cf1dd8SToby Isaac   PetscReal       *points;
94320cf1dd8SToby Isaac   PetscInt        f, dim;
94420cf1dd8SToby Isaac   PetscQuadrature q;
94520cf1dd8SToby Isaac   PetscErrorCode  ierr;
94620cf1dd8SToby Isaac 
94720cf1dd8SToby Isaac   PetscFunctionBegin;
94820cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDimension(sp,&spdim);CHKERRQ(ierr);
94920cf1dd8SToby Isaac   if (!spdim) {
95020cf1dd8SToby Isaac     ierr = PetscQuadratureCreate(PETSC_COMM_SELF,allPoints);CHKERRQ(ierr);
95120cf1dd8SToby Isaac     ierr = PetscQuadratureSetData(*allPoints,0,0,0,NULL,NULL);CHKERRQ(ierr);
95220cf1dd8SToby Isaac   }
95320cf1dd8SToby Isaac   ierr = PetscDualSpaceGetFunctional(sp,0,&q);CHKERRQ(ierr);
95420cf1dd8SToby Isaac   ierr = PetscQuadratureGetData(q,&dim,NULL,&numPoints,NULL,NULL);CHKERRQ(ierr);
95520cf1dd8SToby Isaac   for (f = 1; f < spdim; f++) {
95620cf1dd8SToby Isaac     PetscInt Np;
95720cf1dd8SToby Isaac 
95820cf1dd8SToby Isaac     ierr = PetscDualSpaceGetFunctional(sp,f,&q);CHKERRQ(ierr);
95920cf1dd8SToby Isaac     ierr = PetscQuadratureGetData(q,NULL,NULL,&Np,NULL,NULL);CHKERRQ(ierr);
96020cf1dd8SToby Isaac     numPoints += Np;
96120cf1dd8SToby Isaac   }
96220cf1dd8SToby Isaac   ierr = PetscMalloc1(dim*numPoints,&points);CHKERRQ(ierr);
96320cf1dd8SToby Isaac   for (f = 0, offset = 0; f < spdim; f++) {
96420cf1dd8SToby Isaac     const PetscReal *p;
96520cf1dd8SToby Isaac     PetscInt        Np, i;
96620cf1dd8SToby Isaac 
96720cf1dd8SToby Isaac     ierr = PetscDualSpaceGetFunctional(sp,f,&q);CHKERRQ(ierr);
96820cf1dd8SToby Isaac     ierr = PetscQuadratureGetData(q,NULL,NULL,&Np,&p,NULL);CHKERRQ(ierr);
96920cf1dd8SToby Isaac     for (i = 0; i < Np * dim; i++) {
97020cf1dd8SToby Isaac       points[offset + i] = p[i];
97120cf1dd8SToby Isaac     }
97220cf1dd8SToby Isaac     offset += Np * dim;
97320cf1dd8SToby Isaac   }
97420cf1dd8SToby Isaac   ierr = PetscQuadratureCreate(PETSC_COMM_SELF,allPoints);CHKERRQ(ierr);
97520cf1dd8SToby Isaac   ierr = PetscQuadratureSetData(*allPoints,dim,0,numPoints,points,NULL);CHKERRQ(ierr);
97620cf1dd8SToby Isaac   PetscFunctionReturn(0);
97720cf1dd8SToby Isaac }
97820cf1dd8SToby Isaac 
97920cf1dd8SToby Isaac /*@C
98020cf1dd8SToby Isaac   PetscDualSpaceApplyFVM - Apply a functional from the dual space basis to an input function by assuming a point evaluation functional at the cell centroid.
98120cf1dd8SToby Isaac 
98220cf1dd8SToby Isaac   Input Parameters:
98320cf1dd8SToby Isaac + sp    - The PetscDualSpace object
98420cf1dd8SToby Isaac . f     - The basis functional index
98520cf1dd8SToby Isaac . time  - The time
98620cf1dd8SToby Isaac . cgeom - A context with geometric information for this cell, we currently just use the centroid
98720cf1dd8SToby Isaac . Nc    - The number of components for the function
98820cf1dd8SToby Isaac . func  - The input function
98920cf1dd8SToby Isaac - ctx   - A context for the function
99020cf1dd8SToby Isaac 
99120cf1dd8SToby Isaac   Output Parameter:
99220cf1dd8SToby Isaac . value - The output value (scalar)
99320cf1dd8SToby Isaac 
99420cf1dd8SToby Isaac   Note: The calling sequence for the callback func is given by:
99520cf1dd8SToby Isaac 
99620cf1dd8SToby Isaac $ func(PetscInt dim, PetscReal time, const PetscReal x[],
99720cf1dd8SToby Isaac $      PetscInt numComponents, PetscScalar values[], void *ctx)
99820cf1dd8SToby Isaac 
99920cf1dd8SToby Isaac and the idea is to evaluate the functional as an integral
100020cf1dd8SToby Isaac 
100120cf1dd8SToby Isaac $ n(f) = int dx n(x) . f(x)
100220cf1dd8SToby Isaac 
100320cf1dd8SToby Isaac where both n and f have Nc components.
100420cf1dd8SToby Isaac 
1005*a4ce7ad1SMatthew G. Knepley   Level: beginner
100620cf1dd8SToby Isaac 
100720cf1dd8SToby Isaac .seealso: PetscDualSpaceCreate()
100820cf1dd8SToby Isaac @*/
100920cf1dd8SToby 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)
101020cf1dd8SToby Isaac {
101120cf1dd8SToby Isaac   DM               dm;
101220cf1dd8SToby Isaac   PetscQuadrature  n;
101320cf1dd8SToby Isaac   const PetscReal *points, *weights;
101420cf1dd8SToby Isaac   PetscScalar     *val;
101520cf1dd8SToby Isaac   PetscInt         dimEmbed, qNc, c, Nq, q;
101620cf1dd8SToby Isaac   PetscErrorCode   ierr;
101720cf1dd8SToby Isaac 
101820cf1dd8SToby Isaac   PetscFunctionBegin;
101920cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
102020cf1dd8SToby Isaac   PetscValidPointer(value, 5);
102120cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr);
102220cf1dd8SToby Isaac   ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr);
102320cf1dd8SToby Isaac   ierr = PetscDualSpaceGetFunctional(sp, f, &n);CHKERRQ(ierr);
102420cf1dd8SToby Isaac   ierr = PetscQuadratureGetData(n, NULL, &qNc, &Nq, &points, &weights);CHKERRQ(ierr);
102520cf1dd8SToby Isaac   if (qNc != Nc) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature components %D != function components %D", qNc, Nc);
102620cf1dd8SToby Isaac   ierr = DMGetWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr);
102720cf1dd8SToby Isaac   *value = 0.;
102820cf1dd8SToby Isaac   for (q = 0; q < Nq; ++q) {
102920cf1dd8SToby Isaac     ierr = (*func)(dimEmbed, time, cgeom->centroid, Nc, val, ctx);CHKERRQ(ierr);
103020cf1dd8SToby Isaac     for (c = 0; c < Nc; ++c) {
103120cf1dd8SToby Isaac       *value += val[c]*weights[q*Nc+c];
103220cf1dd8SToby Isaac     }
103320cf1dd8SToby Isaac   }
103420cf1dd8SToby Isaac   ierr = DMRestoreWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr);
103520cf1dd8SToby Isaac   PetscFunctionReturn(0);
103620cf1dd8SToby Isaac }
103720cf1dd8SToby Isaac 
103820cf1dd8SToby Isaac /*@
103920cf1dd8SToby Isaac   PetscDualSpaceGetHeightSubspace - Get the subset of the dual space basis that is supported on a mesh point of a
104020cf1dd8SToby Isaac   given height.  This assumes that the reference cell is symmetric over points of this height.
104120cf1dd8SToby Isaac 
104220cf1dd8SToby Isaac   If the dual space is not defined on mesh points of the given height (e.g. if the space is discontinuous and
104320cf1dd8SToby Isaac   pointwise values are not defined on the element boundaries), or if the implementation of PetscDualSpace does not
104420cf1dd8SToby Isaac   support extracting subspaces, then NULL is returned.
104520cf1dd8SToby Isaac 
104620cf1dd8SToby Isaac   This does not increment the reference count on the returned dual space, and the user should not destroy it.
104720cf1dd8SToby Isaac 
104820cf1dd8SToby Isaac   Not collective
104920cf1dd8SToby Isaac 
105020cf1dd8SToby Isaac   Input Parameters:
105120cf1dd8SToby Isaac + sp - the PetscDualSpace object
105220cf1dd8SToby Isaac - height - the height of the mesh point for which the subspace is desired
105320cf1dd8SToby Isaac 
105420cf1dd8SToby Isaac   Output Parameter:
105520cf1dd8SToby Isaac . subsp - the subspace.  Note that the functionals in the subspace are with respect to the intrinsic geometry of the
105620cf1dd8SToby Isaac   point, which will be of lesser dimension if height > 0.
105720cf1dd8SToby Isaac 
105820cf1dd8SToby Isaac   Level: advanced
105920cf1dd8SToby Isaac 
106020cf1dd8SToby Isaac .seealso: PetscSpaceGetHeightSubspace(), PetscDualSpace
106120cf1dd8SToby Isaac @*/
106220cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetHeightSubspace(PetscDualSpace sp, PetscInt height, PetscDualSpace *subsp)
106320cf1dd8SToby Isaac {
106420cf1dd8SToby Isaac   PetscErrorCode ierr;
106520cf1dd8SToby Isaac 
106620cf1dd8SToby Isaac   PetscFunctionBegin;
106720cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
106820cf1dd8SToby Isaac   PetscValidPointer(subsp, 3);
106920cf1dd8SToby Isaac   *subsp = NULL;
1070aa545514SMatthew G. Knepley   if (sp->ops->getheightsubspace) {ierr = (*sp->ops->getheightsubspace)(sp, height, subsp);CHKERRQ(ierr);}
107120cf1dd8SToby Isaac   PetscFunctionReturn(0);
107220cf1dd8SToby Isaac }
107320cf1dd8SToby Isaac 
107420cf1dd8SToby Isaac /*@
107520cf1dd8SToby Isaac   PetscDualSpaceGetPointSubspace - Get the subset of the dual space basis that is supported on a particular mesh point.
107620cf1dd8SToby Isaac 
107720cf1dd8SToby Isaac   If the dual space is not defined on the mesh point (e.g. if the space is discontinuous and pointwise values are not
107820cf1dd8SToby Isaac   defined on the element boundaries), or if the implementation of PetscDualSpace does not support extracting
107920cf1dd8SToby Isaac   subspaces, then NULL is returned.
108020cf1dd8SToby Isaac 
108120cf1dd8SToby Isaac   This does not increment the reference count on the returned dual space, and the user should not destroy it.
108220cf1dd8SToby Isaac 
108320cf1dd8SToby Isaac   Not collective
108420cf1dd8SToby Isaac 
108520cf1dd8SToby Isaac   Input Parameters:
108620cf1dd8SToby Isaac + sp - the PetscDualSpace object
108720cf1dd8SToby Isaac - point - the point (in the dual space's DM) for which the subspace is desired
108820cf1dd8SToby Isaac 
108920cf1dd8SToby Isaac   Output Parameters:
109020cf1dd8SToby Isaac   bdsp - the subspace.  Note that the functionals in the subspace are with respect to the intrinsic geometry of the
109120cf1dd8SToby Isaac   point, which will be of lesser dimension if height > 0.
109220cf1dd8SToby Isaac 
109320cf1dd8SToby Isaac   Level: advanced
109420cf1dd8SToby Isaac 
109520cf1dd8SToby Isaac .seealso: PetscDualSpace
109620cf1dd8SToby Isaac @*/
109720cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetPointSubspace(PetscDualSpace sp, PetscInt point, PetscDualSpace *bdsp)
109820cf1dd8SToby Isaac {
109920cf1dd8SToby Isaac   PetscErrorCode ierr;
110020cf1dd8SToby Isaac 
110120cf1dd8SToby Isaac   PetscFunctionBegin;
110220cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
110320cf1dd8SToby Isaac   PetscValidPointer(bdsp,2);
110420cf1dd8SToby Isaac   *bdsp = NULL;
110520cf1dd8SToby Isaac   if (sp->ops->getpointsubspace) {
110620cf1dd8SToby Isaac     ierr = (*sp->ops->getpointsubspace)(sp,point,bdsp);CHKERRQ(ierr);
110720cf1dd8SToby Isaac   } else if (sp->ops->getheightsubspace) {
110820cf1dd8SToby Isaac     DM       dm;
110920cf1dd8SToby Isaac     DMLabel  label;
111020cf1dd8SToby Isaac     PetscInt dim, depth, height;
111120cf1dd8SToby Isaac 
111220cf1dd8SToby Isaac     ierr = PetscDualSpaceGetDM(sp,&dm);CHKERRQ(ierr);
111320cf1dd8SToby Isaac     ierr = DMPlexGetDepth(dm,&dim);CHKERRQ(ierr);
111420cf1dd8SToby Isaac     ierr = DMPlexGetDepthLabel(dm,&label);CHKERRQ(ierr);
111520cf1dd8SToby Isaac     ierr = DMLabelGetValue(label,point,&depth);CHKERRQ(ierr);
111620cf1dd8SToby Isaac     height = dim - depth;
111720cf1dd8SToby Isaac     ierr = (*sp->ops->getheightsubspace)(sp,height,bdsp);CHKERRQ(ierr);
111820cf1dd8SToby Isaac   }
111920cf1dd8SToby Isaac   PetscFunctionReturn(0);
112020cf1dd8SToby Isaac }
112120cf1dd8SToby Isaac 
11226f905325SMatthew G. Knepley /*@C
11236f905325SMatthew G. Knepley   PetscDualSpaceGetSymmetries - Returns a description of the symmetries of this basis
11246f905325SMatthew G. Knepley 
11256f905325SMatthew G. Knepley   Not collective
11266f905325SMatthew G. Knepley 
11276f905325SMatthew G. Knepley   Input Parameter:
11286f905325SMatthew G. Knepley . sp - the PetscDualSpace object
11296f905325SMatthew G. Knepley 
11306f905325SMatthew G. Knepley   Output Parameters:
11316f905325SMatthew G. Knepley + perms - Permutations of the local degrees of freedom, parameterized by the point orientation
11326f905325SMatthew G. Knepley - flips - Sign reversal of the local degrees of freedom, parameterized by the point orientation
11336f905325SMatthew G. Knepley 
11346f905325SMatthew G. Knepley   Note: The permutation and flip arrays are organized in the following way
11356f905325SMatthew G. Knepley $ perms[p][ornt][dof # on point] = new local dof #
11366f905325SMatthew G. Knepley $ flips[p][ornt][dof # on point] = reversal or not
11376f905325SMatthew G. Knepley 
11386f905325SMatthew G. Knepley   Level: developer
11396f905325SMatthew G. Knepley 
11406f905325SMatthew G. Knepley .seealso: PetscDualSpaceSetSymmetries()
11416f905325SMatthew G. Knepley @*/
11426f905325SMatthew G. Knepley PetscErrorCode PetscDualSpaceGetSymmetries(PetscDualSpace sp, const PetscInt ****perms, const PetscScalar ****flips)
11436f905325SMatthew G. Knepley {
11446f905325SMatthew G. Knepley   PetscErrorCode ierr;
11456f905325SMatthew G. Knepley 
11466f905325SMatthew G. Knepley   PetscFunctionBegin;
11476f905325SMatthew G. Knepley   PetscValidHeaderSpecific(sp,PETSCDUALSPACE_CLASSID,1);
11486f905325SMatthew G. Knepley   if (perms) {PetscValidPointer(perms,2); *perms = NULL;}
11496f905325SMatthew G. Knepley   if (flips) {PetscValidPointer(flips,3); *flips = NULL;}
11506f905325SMatthew G. Knepley   if (sp->ops->getsymmetries) {ierr = (sp->ops->getsymmetries)(sp,perms,flips);CHKERRQ(ierr);}
11516f905325SMatthew G. Knepley   PetscFunctionReturn(0);
11526f905325SMatthew G. Knepley }
11534bee2e38SMatthew G. Knepley 
11544bee2e38SMatthew G. Knepley /*@
11554bee2e38SMatthew G. Knepley   PetscDualSpaceGetDeRahm - Get the k-simplex associated with the functionals in this dual space
11564bee2e38SMatthew G. Knepley 
11574bee2e38SMatthew G. Knepley   Input Parameter:
11584bee2e38SMatthew G. Knepley . dsp - The PetscDualSpace
11594bee2e38SMatthew G. Knepley 
11604bee2e38SMatthew G. Knepley   Output Parameter:
11614bee2e38SMatthew G. Knepley . k   - The simplex dimension
11624bee2e38SMatthew G. Knepley 
1163*a4ce7ad1SMatthew G. Knepley   Level: developer
11644bee2e38SMatthew G. Knepley 
11654bee2e38SMatthew G. Knepley   Note: Currently supported values are
11664bee2e38SMatthew G. Knepley $ 0: These are H_1 methods that only transform coordinates
11674bee2e38SMatthew G. Knepley $ 1: These are Hcurl methods that transform functions using the covariant Piola transform (COVARIANT_PIOLA_TRANSFORM)
11684bee2e38SMatthew G. Knepley $ 2: These are the same as 1
11694bee2e38SMatthew G. Knepley $ 3: These are Hdiv methods that transform functions using the contravariant Piola transform (CONTRAVARIANT_PIOLA_TRANSFORM)
11704bee2e38SMatthew G. Knepley 
11714bee2e38SMatthew G. Knepley .seealso: PetscDualSpacePullback(), PetscDualSpacePushforward(), PetscDualSpaceTransform(), PetscDualSpaceTransformType
11724bee2e38SMatthew G. Knepley @*/
11734bee2e38SMatthew G. Knepley PetscErrorCode PetscDualSpaceGetDeRahm(PetscDualSpace dsp, PetscInt *k)
11744bee2e38SMatthew G. Knepley {
11754bee2e38SMatthew G. Knepley   PetscFunctionBeginHot;
11764bee2e38SMatthew G. Knepley   PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1);
11774bee2e38SMatthew G. Knepley   PetscValidPointer(k, 2);
11784bee2e38SMatthew G. Knepley   *k = dsp->k;
11794bee2e38SMatthew G. Knepley   PetscFunctionReturn(0);
11804bee2e38SMatthew G. Knepley }
11814bee2e38SMatthew G. Knepley 
11824bee2e38SMatthew G. Knepley /*@C
11834bee2e38SMatthew G. Knepley   PetscDualSpaceTransform - Transform the function values
11844bee2e38SMatthew G. Knepley 
11854bee2e38SMatthew G. Knepley   Input Parameters:
11864bee2e38SMatthew G. Knepley + dsp       - The PetscDualSpace
11874bee2e38SMatthew G. Knepley . trans     - The type of transform
11884bee2e38SMatthew G. Knepley . isInverse - Flag to invert the transform
11894bee2e38SMatthew G. Knepley . fegeom    - The cell geometry
11904bee2e38SMatthew G. Knepley . Nv        - The number of function samples
11914bee2e38SMatthew G. Knepley . Nc        - The number of function components
11924bee2e38SMatthew G. Knepley - vals      - The function values
11934bee2e38SMatthew G. Knepley 
11944bee2e38SMatthew G. Knepley   Output Parameter:
11954bee2e38SMatthew G. Knepley . vals      - The transformed function values
11964bee2e38SMatthew G. Knepley 
1197*a4ce7ad1SMatthew G. Knepley   Level: intermediate
11984bee2e38SMatthew G. Knepley 
1199625e0966SMatthew G. Knepley .seealso: PetscDualSpaceTransformGradient(), PetscDualSpacePullback(), PetscDualSpacePushforward(), PetscDualSpaceTransformType
12004bee2e38SMatthew G. Knepley @*/
12014bee2e38SMatthew G. Knepley PetscErrorCode PetscDualSpaceTransform(PetscDualSpace dsp, PetscDualSpaceTransformType trans, PetscBool isInverse, PetscFEGeom *fegeom, PetscInt Nv, PetscInt Nc, PetscScalar vals[])
12024bee2e38SMatthew G. Knepley {
12034bee2e38SMatthew G. Knepley   PetscInt dim, v, c;
12044bee2e38SMatthew G. Knepley 
12054bee2e38SMatthew G. Knepley   PetscFunctionBeginHot;
12064bee2e38SMatthew G. Knepley   PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1);
12074bee2e38SMatthew G. Knepley   PetscValidPointer(fegeom, 4);
12084bee2e38SMatthew G. Knepley   PetscValidPointer(vals, 7);
12092ae266adSMatthew G. Knepley   dim = dsp->dm->dim;
12104bee2e38SMatthew G. Knepley   /* Assume its a vector, otherwise assume its a bunch of scalars */
12114bee2e38SMatthew G. Knepley   if (Nc == 1 || Nc != dim) PetscFunctionReturn(0);
12124bee2e38SMatthew G. Knepley   switch (trans) {
12134bee2e38SMatthew G. Knepley     case IDENTITY_TRANSFORM: break;
12144bee2e38SMatthew G. Knepley     case COVARIANT_PIOLA_TRANSFORM: /* Covariant Piola mapping $\sigma^*(F) = J^{-T} F \circ \phi^{-1)$ */
12154bee2e38SMatthew G. Knepley     if (isInverse) {
12164bee2e38SMatthew G. Knepley       for (v = 0; v < Nv; ++v) {
12174bee2e38SMatthew G. Knepley         switch (dim)
12184bee2e38SMatthew G. Knepley         {
12196142fa51SMatthew G. Knepley           case 2: DMPlex_MultTranspose2DReal_Internal(fegeom->J, 1, &vals[v*Nc], &vals[v*Nc]);break;
12206142fa51SMatthew G. Knepley           case 3: DMPlex_MultTranspose3DReal_Internal(fegeom->J, 1, &vals[v*Nc], &vals[v*Nc]);break;
12214bee2e38SMatthew G. Knepley           default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim);
12224bee2e38SMatthew G. Knepley         }
12234bee2e38SMatthew G. Knepley       }
12244bee2e38SMatthew G. Knepley     } else {
12254bee2e38SMatthew G. Knepley       for (v = 0; v < Nv; ++v) {
12264bee2e38SMatthew G. Knepley         switch (dim)
12274bee2e38SMatthew G. Knepley         {
12286142fa51SMatthew G. Knepley           case 2: DMPlex_MultTranspose2DReal_Internal(fegeom->invJ, 1, &vals[v*Nc], &vals[v*Nc]);break;
12296142fa51SMatthew G. Knepley           case 3: DMPlex_MultTranspose3DReal_Internal(fegeom->invJ, 1, &vals[v*Nc], &vals[v*Nc]);break;
12304bee2e38SMatthew G. Knepley           default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim);
12314bee2e38SMatthew G. Knepley         }
12324bee2e38SMatthew G. Knepley       }
12334bee2e38SMatthew G. Knepley     }
12344bee2e38SMatthew G. Knepley     break;
12354bee2e38SMatthew G. Knepley     case CONTRAVARIANT_PIOLA_TRANSFORM: /* Contravariant Piola mapping $\sigma^*(F) = \frac{1}{|\det J|} J F \circ \phi^{-1}$ */
12364bee2e38SMatthew G. Knepley     if (isInverse) {
12374bee2e38SMatthew G. Knepley       for (v = 0; v < Nv; ++v) {
12384bee2e38SMatthew G. Knepley         switch (dim)
12394bee2e38SMatthew G. Knepley         {
12406142fa51SMatthew G. Knepley           case 2: DMPlex_Mult2DReal_Internal(fegeom->invJ, 1, &vals[v*Nc], &vals[v*Nc]);break;
12416142fa51SMatthew G. Knepley           case 3: DMPlex_Mult3DReal_Internal(fegeom->invJ, 1, &vals[v*Nc], &vals[v*Nc]);break;
12424bee2e38SMatthew G. Knepley           default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim);
12434bee2e38SMatthew G. Knepley         }
12444bee2e38SMatthew G. Knepley         for (c = 0; c < Nc; ++c) vals[v*Nc+c] *= fegeom->detJ[0];
12454bee2e38SMatthew G. Knepley       }
12464bee2e38SMatthew G. Knepley     } else {
12474bee2e38SMatthew G. Knepley       for (v = 0; v < Nv; ++v) {
12484bee2e38SMatthew G. Knepley         switch (dim)
12494bee2e38SMatthew G. Knepley         {
12506142fa51SMatthew G. Knepley           case 2: DMPlex_Mult2DReal_Internal(fegeom->J, 1, &vals[v*Nc], &vals[v*Nc]);break;
12516142fa51SMatthew G. Knepley           case 3: DMPlex_Mult3DReal_Internal(fegeom->J, 1, &vals[v*Nc], &vals[v*Nc]);break;
12524bee2e38SMatthew G. Knepley           default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim);
12534bee2e38SMatthew G. Knepley         }
12544bee2e38SMatthew G. Knepley         for (c = 0; c < Nc; ++c) vals[v*Nc+c] /= fegeom->detJ[0];
12554bee2e38SMatthew G. Knepley       }
12564bee2e38SMatthew G. Knepley     }
12574bee2e38SMatthew G. Knepley     break;
12584bee2e38SMatthew G. Knepley   }
12594bee2e38SMatthew G. Knepley   PetscFunctionReturn(0);
12604bee2e38SMatthew G. Knepley }
12614bee2e38SMatthew G. Knepley /*@C
12624bee2e38SMatthew G. Knepley   PetscDualSpaceTransformGradient - Transform the function gradient values
12634bee2e38SMatthew G. Knepley 
12644bee2e38SMatthew G. Knepley   Input Parameters:
12654bee2e38SMatthew G. Knepley + dsp       - The PetscDualSpace
12664bee2e38SMatthew G. Knepley . trans     - The type of transform
12674bee2e38SMatthew G. Knepley . isInverse - Flag to invert the transform
12684bee2e38SMatthew G. Knepley . fegeom    - The cell geometry
12694bee2e38SMatthew G. Knepley . Nv        - The number of function gradient samples
12704bee2e38SMatthew G. Knepley . Nc        - The number of function components
12714bee2e38SMatthew G. Knepley - vals      - The function gradient values
12724bee2e38SMatthew G. Knepley 
12734bee2e38SMatthew G. Knepley   Output Parameter:
12744bee2e38SMatthew G. Knepley . vals      - The transformed function values
12754bee2e38SMatthew G. Knepley 
1276*a4ce7ad1SMatthew G. Knepley   Level: intermediate
12774bee2e38SMatthew G. Knepley 
1278625e0966SMatthew G. Knepley .seealso: PetscDualSpaceTransform(), PetscDualSpacePullback(), PetscDualSpacePushforward(), PetscDualSpaceTransformType
12794bee2e38SMatthew G. Knepley @*/
12804bee2e38SMatthew G. Knepley PetscErrorCode PetscDualSpaceTransformGradient(PetscDualSpace dsp, PetscDualSpaceTransformType trans, PetscBool isInverse, PetscFEGeom *fegeom, PetscInt Nv, PetscInt Nc, PetscScalar vals[])
12814bee2e38SMatthew G. Knepley {
12824bee2e38SMatthew G. Knepley   PetscInt dim, v, c, d;
12834bee2e38SMatthew G. Knepley 
12844bee2e38SMatthew G. Knepley   PetscFunctionBeginHot;
12854bee2e38SMatthew G. Knepley   PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1);
12864bee2e38SMatthew G. Knepley   PetscValidPointer(fegeom, 4);
12874bee2e38SMatthew G. Knepley   PetscValidPointer(vals, 7);
12882ae266adSMatthew G. Knepley   dim = dsp->dm->dim;
12894bee2e38SMatthew G. Knepley   /* Transform gradient */
12904bee2e38SMatthew G. Knepley   for (v = 0; v < Nv; ++v) {
12914bee2e38SMatthew G. Knepley     for (c = 0; c < Nc; ++c) {
12924bee2e38SMatthew G. Knepley       switch (dim)
12934bee2e38SMatthew G. Knepley       {
1294d70dd48eSMatthew G. Knepley         case 1: vals[(v*Nc+c)*dim] *= fegeom->invJ[0];
12956142fa51SMatthew G. Knepley         case 2: DMPlex_MultTranspose2DReal_Internal(fegeom->invJ, 1, &vals[(v*Nc+c)*dim], &vals[(v*Nc+c)*dim]);break;
12966142fa51SMatthew G. Knepley         case 3: DMPlex_MultTranspose3DReal_Internal(fegeom->invJ, 1, &vals[(v*Nc+c)*dim], &vals[(v*Nc+c)*dim]);break;
12974bee2e38SMatthew G. Knepley         default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim);
12984bee2e38SMatthew G. Knepley       }
12994bee2e38SMatthew G. Knepley     }
13004bee2e38SMatthew G. Knepley   }
13014bee2e38SMatthew G. Knepley   /* Assume its a vector, otherwise assume its a bunch of scalars */
13024bee2e38SMatthew G. Knepley   if (Nc == 1 || Nc != dim) PetscFunctionReturn(0);
13034bee2e38SMatthew G. Knepley   switch (trans) {
13044bee2e38SMatthew G. Knepley     case IDENTITY_TRANSFORM: break;
13054bee2e38SMatthew G. Knepley     case COVARIANT_PIOLA_TRANSFORM: /* Covariant Piola mapping $\sigma^*(F) = J^{-T} F \circ \phi^{-1)$ */
13064bee2e38SMatthew G. Knepley     if (isInverse) {
13074bee2e38SMatthew G. Knepley       for (v = 0; v < Nv; ++v) {
13084bee2e38SMatthew G. Knepley         for (d = 0; d < dim; ++d) {
13094bee2e38SMatthew G. Knepley           switch (dim)
13104bee2e38SMatthew G. Knepley           {
13116142fa51SMatthew G. Knepley             case 2: DMPlex_MultTranspose2DReal_Internal(fegeom->J, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
13126142fa51SMatthew G. Knepley             case 3: DMPlex_MultTranspose3DReal_Internal(fegeom->J, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
13134bee2e38SMatthew G. Knepley             default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim);
13144bee2e38SMatthew G. Knepley           }
13154bee2e38SMatthew G. Knepley         }
13164bee2e38SMatthew G. Knepley       }
13174bee2e38SMatthew G. Knepley     } else {
13184bee2e38SMatthew G. Knepley       for (v = 0; v < Nv; ++v) {
13194bee2e38SMatthew G. Knepley         for (d = 0; d < dim; ++d) {
13204bee2e38SMatthew G. Knepley           switch (dim)
13214bee2e38SMatthew G. Knepley           {
13226142fa51SMatthew G. Knepley             case 2: DMPlex_MultTranspose2DReal_Internal(fegeom->invJ, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
13236142fa51SMatthew G. Knepley             case 3: DMPlex_MultTranspose3DReal_Internal(fegeom->invJ, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
13244bee2e38SMatthew G. Knepley             default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim);
13254bee2e38SMatthew G. Knepley           }
13264bee2e38SMatthew G. Knepley         }
13274bee2e38SMatthew G. Knepley       }
13284bee2e38SMatthew G. Knepley     }
13294bee2e38SMatthew G. Knepley     break;
13304bee2e38SMatthew G. Knepley     case CONTRAVARIANT_PIOLA_TRANSFORM: /* Contravariant Piola mapping $\sigma^*(F) = \frac{1}{|\det J|} J F \circ \phi^{-1}$ */
13314bee2e38SMatthew G. Knepley     if (isInverse) {
13324bee2e38SMatthew G. Knepley       for (v = 0; v < Nv; ++v) {
13334bee2e38SMatthew G. Knepley         for (d = 0; d < dim; ++d) {
13344bee2e38SMatthew G. Knepley           switch (dim)
13354bee2e38SMatthew G. Knepley           {
13366142fa51SMatthew G. Knepley             case 2: DMPlex_Mult2DReal_Internal(fegeom->invJ, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
13376142fa51SMatthew G. Knepley             case 3: DMPlex_Mult3DReal_Internal(fegeom->invJ, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
13384bee2e38SMatthew G. Knepley             default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim);
13394bee2e38SMatthew G. Knepley           }
13404bee2e38SMatthew G. Knepley           for (c = 0; c < Nc; ++c) vals[(v*Nc+c)*dim+d] *= fegeom->detJ[0];
13414bee2e38SMatthew G. Knepley         }
13424bee2e38SMatthew G. Knepley       }
13434bee2e38SMatthew G. Knepley     } else {
13444bee2e38SMatthew G. Knepley       for (v = 0; v < Nv; ++v) {
13454bee2e38SMatthew G. Knepley         for (d = 0; d < dim; ++d) {
13464bee2e38SMatthew G. Knepley           switch (dim)
13474bee2e38SMatthew G. Knepley           {
13486142fa51SMatthew G. Knepley             case 2: DMPlex_Mult2DReal_Internal(fegeom->J, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
13496142fa51SMatthew G. Knepley             case 3: DMPlex_Mult3DReal_Internal(fegeom->J, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break;
13504bee2e38SMatthew G. Knepley             default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim);
13514bee2e38SMatthew G. Knepley           }
13524bee2e38SMatthew G. Knepley           for (c = 0; c < Nc; ++c) vals[(v*Nc+c)*dim+d] /= fegeom->detJ[0];
13534bee2e38SMatthew G. Knepley         }
13544bee2e38SMatthew G. Knepley       }
13554bee2e38SMatthew G. Knepley     }
13564bee2e38SMatthew G. Knepley     break;
13574bee2e38SMatthew G. Knepley   }
13584bee2e38SMatthew G. Knepley   PetscFunctionReturn(0);
13594bee2e38SMatthew G. Knepley }
13604bee2e38SMatthew G. Knepley 
13614bee2e38SMatthew G. Knepley /*@C
13624bee2e38SMatthew 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.
13634bee2e38SMatthew G. Knepley 
13644bee2e38SMatthew G. Knepley   Input Parameters:
13654bee2e38SMatthew G. Knepley + dsp        - The PetscDualSpace
13664bee2e38SMatthew G. Knepley . fegeom     - The geometry for this cell
13674bee2e38SMatthew G. Knepley . Nq         - The number of function samples
13684bee2e38SMatthew G. Knepley . Nc         - The number of function components
13694bee2e38SMatthew G. Knepley - pointEval  - The function values
13704bee2e38SMatthew G. Knepley 
13714bee2e38SMatthew G. Knepley   Output Parameter:
13724bee2e38SMatthew G. Knepley . pointEval  - The transformed function values
13734bee2e38SMatthew G. Knepley 
13744bee2e38SMatthew G. Knepley   Level: advanced
13754bee2e38SMatthew G. Knepley 
13764bee2e38SMatthew 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.
13774bee2e38SMatthew G. Knepley 
13784bee2e38SMatthew G. Knepley .seealso: PetscDualSpacePushforward(), PetscDualSpaceTransform(), PetscDualSpaceGetDeRahm()
13794bee2e38SMatthew G. Knepley @*/
13804bee2e38SMatthew G. Knepley PetscErrorCode PetscDualSpacePullback(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, PetscScalar pointEval[])
13814bee2e38SMatthew G. Knepley {
13824bee2e38SMatthew G. Knepley   PetscDualSpaceTransformType trans;
13834bee2e38SMatthew G. Knepley   PetscErrorCode              ierr;
13844bee2e38SMatthew G. Knepley 
13854bee2e38SMatthew G. Knepley   PetscFunctionBeginHot;
13864bee2e38SMatthew G. Knepley   PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1);
13874bee2e38SMatthew G. Knepley   PetscValidPointer(fegeom, 2);
13884bee2e38SMatthew G. Knepley   PetscValidPointer(pointEval, 5);
13894bee2e38SMatthew G. Knepley   /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k.
13904bee2e38SMatthew G. Knepley      This determines their transformation properties. */
13912ae266adSMatthew G. Knepley   switch (dsp->k)
13924bee2e38SMatthew G. Knepley   {
13934bee2e38SMatthew G. Knepley     case 0: /* H^1 point evaluations */
13944bee2e38SMatthew G. Knepley     trans = IDENTITY_TRANSFORM;break;
13954bee2e38SMatthew G. Knepley     case 1: /* Hcurl preserves tangential edge traces  */
13964bee2e38SMatthew G. Knepley     case 2:
13974bee2e38SMatthew G. Knepley     trans = COVARIANT_PIOLA_TRANSFORM;break;
13984bee2e38SMatthew G. Knepley     case 3: /* Hdiv preserve normal traces */
13994bee2e38SMatthew G. Knepley     trans = CONTRAVARIANT_PIOLA_TRANSFORM;break;
14002ae266adSMatthew G. Knepley     default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %D for transformation", dsp->k);
14014bee2e38SMatthew G. Knepley   }
14024bee2e38SMatthew G. Knepley   ierr = PetscDualSpaceTransform(dsp, trans, PETSC_TRUE, fegeom, Nq, Nc, pointEval);CHKERRQ(ierr);
14034bee2e38SMatthew G. Knepley   PetscFunctionReturn(0);
14044bee2e38SMatthew G. Knepley }
14054bee2e38SMatthew G. Knepley 
14064bee2e38SMatthew G. Knepley /*@C
14074bee2e38SMatthew 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.
14084bee2e38SMatthew G. Knepley 
14094bee2e38SMatthew G. Knepley   Input Parameters:
14104bee2e38SMatthew G. Knepley + dsp        - The PetscDualSpace
14114bee2e38SMatthew G. Knepley . fegeom     - The geometry for this cell
14124bee2e38SMatthew G. Knepley . Nq         - The number of function samples
14134bee2e38SMatthew G. Knepley . Nc         - The number of function components
14144bee2e38SMatthew G. Knepley - pointEval  - The function values
14154bee2e38SMatthew G. Knepley 
14164bee2e38SMatthew G. Knepley   Output Parameter:
14174bee2e38SMatthew G. Knepley . pointEval  - The transformed function values
14184bee2e38SMatthew G. Knepley 
14194bee2e38SMatthew G. Knepley   Level: advanced
14204bee2e38SMatthew G. Knepley 
14214bee2e38SMatthew 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.
14224bee2e38SMatthew G. Knepley 
14234bee2e38SMatthew G. Knepley .seealso: PetscDualSpacePullback(), PetscDualSpaceTransform(), PetscDualSpaceGetDeRahm()
14244bee2e38SMatthew G. Knepley @*/
14254bee2e38SMatthew G. Knepley PetscErrorCode PetscDualSpacePushforward(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, PetscScalar pointEval[])
14264bee2e38SMatthew G. Knepley {
14274bee2e38SMatthew G. Knepley   PetscDualSpaceTransformType trans;
14284bee2e38SMatthew G. Knepley   PetscErrorCode              ierr;
14294bee2e38SMatthew G. Knepley 
14304bee2e38SMatthew G. Knepley   PetscFunctionBeginHot;
14314bee2e38SMatthew G. Knepley   PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1);
14324bee2e38SMatthew G. Knepley   PetscValidPointer(fegeom, 2);
14334bee2e38SMatthew G. Knepley   PetscValidPointer(pointEval, 5);
14344bee2e38SMatthew G. Knepley   /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k.
14354bee2e38SMatthew G. Knepley      This determines their transformation properties. */
14362ae266adSMatthew G. Knepley   switch (dsp->k)
14374bee2e38SMatthew G. Knepley   {
14384bee2e38SMatthew G. Knepley     case 0: /* H^1 point evaluations */
14394bee2e38SMatthew G. Knepley     trans = IDENTITY_TRANSFORM;break;
14404bee2e38SMatthew G. Knepley     case 1: /* Hcurl preserves tangential edge traces  */
14414bee2e38SMatthew G. Knepley     case 2:
14424bee2e38SMatthew G. Knepley     trans = COVARIANT_PIOLA_TRANSFORM;break;
14434bee2e38SMatthew G. Knepley     case 3: /* Hdiv preserve normal traces */
14444bee2e38SMatthew G. Knepley     trans = CONTRAVARIANT_PIOLA_TRANSFORM;break;
14452ae266adSMatthew G. Knepley     default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %D for transformation", dsp->k);
14464bee2e38SMatthew G. Knepley   }
14474bee2e38SMatthew G. Knepley   ierr = PetscDualSpaceTransform(dsp, trans, PETSC_FALSE, fegeom, Nq, Nc, pointEval);CHKERRQ(ierr);
14484bee2e38SMatthew G. Knepley   PetscFunctionReturn(0);
14494bee2e38SMatthew G. Knepley }
14504bee2e38SMatthew G. Knepley 
14514bee2e38SMatthew G. Knepley /*@C
14524bee2e38SMatthew 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.
14534bee2e38SMatthew G. Knepley 
14544bee2e38SMatthew G. Knepley   Input Parameters:
14554bee2e38SMatthew G. Knepley + dsp        - The PetscDualSpace
14564bee2e38SMatthew G. Knepley . fegeom     - The geometry for this cell
14574bee2e38SMatthew G. Knepley . Nq         - The number of function gradient samples
14584bee2e38SMatthew G. Knepley . Nc         - The number of function components
14594bee2e38SMatthew G. Knepley - pointEval  - The function gradient values
14604bee2e38SMatthew G. Knepley 
14614bee2e38SMatthew G. Knepley   Output Parameter:
14624bee2e38SMatthew G. Knepley . pointEval  - The transformed function gradient values
14634bee2e38SMatthew G. Knepley 
14644bee2e38SMatthew G. Knepley   Level: advanced
14654bee2e38SMatthew G. Knepley 
14664bee2e38SMatthew 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.
14674bee2e38SMatthew G. Knepley 
14684bee2e38SMatthew G. Knepley .seealso: PetscDualSpacePushforward(), PPetscDualSpacePullback(), PetscDualSpaceTransform(), PetscDualSpaceGetDeRahm()
1469dc0529c6SBarry Smith @*/
1470dc0529c6SBarry Smith PetscErrorCode PetscDualSpacePushforwardGradient(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, PetscScalar pointEval[])
14714bee2e38SMatthew G. Knepley {
14724bee2e38SMatthew G. Knepley   PetscDualSpaceTransformType trans;
14734bee2e38SMatthew G. Knepley   PetscErrorCode              ierr;
14744bee2e38SMatthew G. Knepley 
14754bee2e38SMatthew G. Knepley   PetscFunctionBeginHot;
14764bee2e38SMatthew G. Knepley   PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1);
14774bee2e38SMatthew G. Knepley   PetscValidPointer(fegeom, 2);
14784bee2e38SMatthew G. Knepley   PetscValidPointer(pointEval, 5);
14794bee2e38SMatthew G. Knepley   /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k.
14804bee2e38SMatthew G. Knepley      This determines their transformation properties. */
14812ae266adSMatthew G. Knepley   switch (dsp->k)
14824bee2e38SMatthew G. Knepley   {
14834bee2e38SMatthew G. Knepley     case 0: /* H^1 point evaluations */
14844bee2e38SMatthew G. Knepley     trans = IDENTITY_TRANSFORM;break;
14854bee2e38SMatthew G. Knepley     case 1: /* Hcurl preserves tangential edge traces  */
14864bee2e38SMatthew G. Knepley     case 2:
14874bee2e38SMatthew G. Knepley     trans = COVARIANT_PIOLA_TRANSFORM;break;
14884bee2e38SMatthew G. Knepley     case 3: /* Hdiv preserve normal traces */
14894bee2e38SMatthew G. Knepley     trans = CONTRAVARIANT_PIOLA_TRANSFORM;break;
14902ae266adSMatthew G. Knepley     default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %D for transformation", dsp->k);
14914bee2e38SMatthew G. Knepley   }
14924bee2e38SMatthew G. Knepley   ierr = PetscDualSpaceTransformGradient(dsp, trans, PETSC_FALSE, fegeom, Nq, Nc, pointEval);CHKERRQ(ierr);
14934bee2e38SMatthew G. Knepley   PetscFunctionReturn(0);
14944bee2e38SMatthew G. Knepley }
1495