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