120cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/ 220cf1dd8SToby Isaac #include <petscdmplex.h> 320cf1dd8SToby Isaac 4*a4ce7ad1SMatthew G. Knepley static PetscErrorCode PetscDualSpaceSetUp_Simple(PetscDualSpace sp) 520cf1dd8SToby Isaac { 620cf1dd8SToby Isaac PetscDualSpace_Simple *s = (PetscDualSpace_Simple *) sp->data; 720cf1dd8SToby Isaac DM dm = sp->dm; 820cf1dd8SToby Isaac PetscInt dim; 920cf1dd8SToby Isaac PetscErrorCode ierr; 1020cf1dd8SToby Isaac 1120cf1dd8SToby Isaac PetscFunctionBegin; 1220cf1dd8SToby Isaac ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1320cf1dd8SToby Isaac ierr = PetscCalloc1(dim+1, &s->numDof);CHKERRQ(ierr); 1420cf1dd8SToby Isaac PetscFunctionReturn(0); 1520cf1dd8SToby Isaac } 1620cf1dd8SToby Isaac 17*a4ce7ad1SMatthew G. Knepley static PetscErrorCode PetscDualSpaceDestroy_Simple(PetscDualSpace sp) 1820cf1dd8SToby Isaac { 1920cf1dd8SToby Isaac PetscDualSpace_Simple *s = (PetscDualSpace_Simple *) sp->data; 2020cf1dd8SToby Isaac PetscErrorCode ierr; 2120cf1dd8SToby Isaac 2220cf1dd8SToby Isaac PetscFunctionBegin; 2320cf1dd8SToby Isaac ierr = PetscFree(s->numDof);CHKERRQ(ierr); 2420cf1dd8SToby Isaac ierr = PetscFree(s);CHKERRQ(ierr); 2520cf1dd8SToby Isaac ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceSimpleSetDimension_C", NULL);CHKERRQ(ierr); 2620cf1dd8SToby Isaac ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceSimpleSetFunctional_C", NULL);CHKERRQ(ierr); 2720cf1dd8SToby Isaac PetscFunctionReturn(0); 2820cf1dd8SToby Isaac } 2920cf1dd8SToby Isaac 30*a4ce7ad1SMatthew G. Knepley static PetscErrorCode PetscDualSpaceDuplicate_Simple(PetscDualSpace sp, PetscDualSpace *spNew) 3120cf1dd8SToby Isaac { 3220cf1dd8SToby Isaac PetscInt dim, d, Nc; 3320cf1dd8SToby Isaac PetscErrorCode ierr; 3420cf1dd8SToby Isaac 3520cf1dd8SToby Isaac PetscFunctionBegin; 3620cf1dd8SToby Isaac ierr = PetscDualSpaceCreate(PetscObjectComm((PetscObject) sp), spNew);CHKERRQ(ierr); 3720cf1dd8SToby Isaac ierr = PetscDualSpaceSetType(*spNew, PETSCDUALSPACESIMPLE);CHKERRQ(ierr); 3820cf1dd8SToby Isaac ierr = PetscDualSpaceGetNumComponents(sp, &Nc);CHKERRQ(ierr); 3920cf1dd8SToby Isaac ierr = PetscDualSpaceSetNumComponents(sp, Nc);CHKERRQ(ierr); 4020cf1dd8SToby Isaac ierr = PetscDualSpaceGetDimension(sp, &dim);CHKERRQ(ierr); 4120cf1dd8SToby Isaac ierr = PetscDualSpaceSimpleSetDimension(*spNew, dim);CHKERRQ(ierr); 4220cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 4320cf1dd8SToby Isaac PetscQuadrature q; 4420cf1dd8SToby Isaac 4520cf1dd8SToby Isaac ierr = PetscDualSpaceGetFunctional(sp, d, &q);CHKERRQ(ierr); 4620cf1dd8SToby Isaac ierr = PetscDualSpaceSimpleSetFunctional(*spNew, d, q);CHKERRQ(ierr); 4720cf1dd8SToby Isaac } 4820cf1dd8SToby Isaac PetscFunctionReturn(0); 4920cf1dd8SToby Isaac } 5020cf1dd8SToby Isaac 51*a4ce7ad1SMatthew G. Knepley static PetscErrorCode PetscDualSpaceSetFromOptions_Simple(PetscOptionItems *PetscOptionsObject,PetscDualSpace sp) 5220cf1dd8SToby Isaac { 5320cf1dd8SToby Isaac PetscFunctionBegin; 5420cf1dd8SToby Isaac PetscFunctionReturn(0); 5520cf1dd8SToby Isaac } 5620cf1dd8SToby Isaac 57*a4ce7ad1SMatthew G. Knepley static PetscErrorCode PetscDualSpaceGetDimension_Simple(PetscDualSpace sp, PetscInt *dim) 5820cf1dd8SToby Isaac { 5920cf1dd8SToby Isaac PetscDualSpace_Simple *s = (PetscDualSpace_Simple *) sp->data; 6020cf1dd8SToby Isaac 6120cf1dd8SToby Isaac PetscFunctionBegin; 6220cf1dd8SToby Isaac *dim = s->dim; 6320cf1dd8SToby Isaac PetscFunctionReturn(0); 6420cf1dd8SToby Isaac } 6520cf1dd8SToby Isaac 66*a4ce7ad1SMatthew G. Knepley static PetscErrorCode PetscDualSpaceSimpleSetDimension_Simple(PetscDualSpace sp, const PetscInt dim) 6720cf1dd8SToby Isaac { 6820cf1dd8SToby Isaac PetscDualSpace_Simple *s = (PetscDualSpace_Simple *) sp->data; 6920cf1dd8SToby Isaac DM dm; 7020cf1dd8SToby Isaac PetscInt spatialDim, f; 7120cf1dd8SToby Isaac PetscErrorCode ierr; 7220cf1dd8SToby Isaac 7320cf1dd8SToby Isaac PetscFunctionBegin; 7420cf1dd8SToby Isaac for (f = 0; f < s->dim; ++f) {ierr = PetscQuadratureDestroy(&sp->functional[f]);CHKERRQ(ierr);} 7520cf1dd8SToby Isaac ierr = PetscFree(sp->functional);CHKERRQ(ierr); 7620cf1dd8SToby Isaac s->dim = dim; 7720cf1dd8SToby Isaac ierr = PetscCalloc1(s->dim, &sp->functional);CHKERRQ(ierr); 7820cf1dd8SToby Isaac ierr = PetscFree(s->numDof);CHKERRQ(ierr); 7920cf1dd8SToby Isaac ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr); 8020cf1dd8SToby Isaac ierr = DMGetCoordinateDim(dm, &spatialDim);CHKERRQ(ierr); 8120cf1dd8SToby Isaac ierr = PetscCalloc1(spatialDim+1, &s->numDof);CHKERRQ(ierr); 8220cf1dd8SToby Isaac s->numDof[spatialDim] = dim; 8320cf1dd8SToby Isaac PetscFunctionReturn(0); 8420cf1dd8SToby Isaac } 8520cf1dd8SToby Isaac 86*a4ce7ad1SMatthew G. Knepley static PetscErrorCode PetscDualSpaceGetNumDof_Simple(PetscDualSpace sp, const PetscInt **numDof) 8720cf1dd8SToby Isaac { 8820cf1dd8SToby Isaac PetscDualSpace_Simple *s = (PetscDualSpace_Simple *) sp->data; 8920cf1dd8SToby Isaac 9020cf1dd8SToby Isaac PetscFunctionBegin; 9120cf1dd8SToby Isaac *numDof = s->numDof; 9220cf1dd8SToby Isaac PetscFunctionReturn(0); 9320cf1dd8SToby Isaac } 9420cf1dd8SToby Isaac 95*a4ce7ad1SMatthew G. Knepley static PetscErrorCode PetscDualSpaceSimpleSetFunctional_Simple(PetscDualSpace sp, PetscInt f, PetscQuadrature q) 9620cf1dd8SToby Isaac { 9720cf1dd8SToby Isaac PetscDualSpace_Simple *s = (PetscDualSpace_Simple *) sp->data; 9820cf1dd8SToby Isaac PetscReal *weights; 9920cf1dd8SToby Isaac PetscInt Nc, c, Nq, p; 10020cf1dd8SToby Isaac PetscErrorCode ierr; 10120cf1dd8SToby Isaac 10220cf1dd8SToby Isaac PetscFunctionBegin; 10320cf1dd8SToby Isaac if ((f < 0) || (f >= s->dim)) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_OUTOFRANGE, "Basis index %d not in [0, %d)", f, s->dim); 10420cf1dd8SToby Isaac ierr = PetscQuadratureDuplicate(q, &sp->functional[f]);CHKERRQ(ierr); 10520cf1dd8SToby Isaac /* Reweight so that it has unit volume: Do we want to do this for Nc > 1? */ 10620cf1dd8SToby Isaac ierr = PetscQuadratureGetData(sp->functional[f], NULL, &Nc, &Nq, NULL, (const PetscReal **) &weights);CHKERRQ(ierr); 10720cf1dd8SToby Isaac for (c = 0; c < Nc; ++c) { 10820cf1dd8SToby Isaac PetscReal vol = 0.0; 10920cf1dd8SToby Isaac 11020cf1dd8SToby Isaac for (p = 0; p < Nq; ++p) vol += weights[p*Nc+c]; 11120cf1dd8SToby Isaac for (p = 0; p < Nq; ++p) weights[p*Nc+c] /= (vol == 0.0 ? 1.0 : vol); 11220cf1dd8SToby Isaac } 11320cf1dd8SToby Isaac PetscFunctionReturn(0); 11420cf1dd8SToby Isaac } 11520cf1dd8SToby Isaac 11620cf1dd8SToby Isaac /*@ 11720cf1dd8SToby Isaac PetscDualSpaceSimpleSetDimension - Set the number of functionals in the dual space basis 11820cf1dd8SToby Isaac 119d083f849SBarry Smith Logically Collective on sp 12020cf1dd8SToby Isaac 12120cf1dd8SToby Isaac Input Parameters: 12220cf1dd8SToby Isaac + sp - the PetscDualSpace 12320cf1dd8SToby Isaac - dim - the basis dimension 12420cf1dd8SToby Isaac 12520cf1dd8SToby Isaac Level: intermediate 12620cf1dd8SToby Isaac 12720cf1dd8SToby Isaac .seealso: PetscDualSpaceSimpleSetFunctional() 12820cf1dd8SToby Isaac @*/ 12920cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSimpleSetDimension(PetscDualSpace sp, PetscInt dim) 13020cf1dd8SToby Isaac { 13120cf1dd8SToby Isaac PetscErrorCode ierr; 13220cf1dd8SToby Isaac 13320cf1dd8SToby Isaac PetscFunctionBegin; 13420cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 13520cf1dd8SToby Isaac PetscValidLogicalCollectiveInt(sp, dim, 2); 13620cf1dd8SToby Isaac ierr = PetscTryMethod(sp, "PetscDualSpaceSimpleSetDimension_C", (PetscDualSpace,PetscInt),(sp,dim));CHKERRQ(ierr); 13720cf1dd8SToby Isaac PetscFunctionReturn(0); 13820cf1dd8SToby Isaac } 13920cf1dd8SToby Isaac 14020cf1dd8SToby Isaac /*@ 14120cf1dd8SToby Isaac PetscDualSpaceSimpleSetFunctional - Set the given basis element for this dual space 14220cf1dd8SToby Isaac 14320cf1dd8SToby Isaac Not Collective 14420cf1dd8SToby Isaac 14520cf1dd8SToby Isaac Input Parameters: 14620cf1dd8SToby Isaac + sp - the PetscDualSpace 14720cf1dd8SToby Isaac . f - the basis index 14820cf1dd8SToby Isaac - q - the basis functional 14920cf1dd8SToby Isaac 15020cf1dd8SToby Isaac Level: intermediate 15120cf1dd8SToby Isaac 15220cf1dd8SToby Isaac Note: The quadrature will be reweighted so that it has unit volume. 15320cf1dd8SToby Isaac 15420cf1dd8SToby Isaac .seealso: PetscDualSpaceSimpleSetDimension() 15520cf1dd8SToby Isaac @*/ 15620cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSimpleSetFunctional(PetscDualSpace sp, PetscInt func, PetscQuadrature q) 15720cf1dd8SToby Isaac { 15820cf1dd8SToby Isaac PetscErrorCode ierr; 15920cf1dd8SToby Isaac 16020cf1dd8SToby Isaac PetscFunctionBegin; 16120cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 16220cf1dd8SToby Isaac ierr = PetscTryMethod(sp, "PetscDualSpaceSimpleSetFunctional_C", (PetscDualSpace,PetscInt,PetscQuadrature),(sp,func,q));CHKERRQ(ierr); 16320cf1dd8SToby Isaac PetscFunctionReturn(0); 16420cf1dd8SToby Isaac } 16520cf1dd8SToby Isaac 166*a4ce7ad1SMatthew G. Knepley static PetscErrorCode PetscDualSpaceInitialize_Simple(PetscDualSpace sp) 16720cf1dd8SToby Isaac { 16820cf1dd8SToby Isaac PetscFunctionBegin; 16920cf1dd8SToby Isaac sp->ops->setfromoptions = PetscDualSpaceSetFromOptions_Simple; 17020cf1dd8SToby Isaac sp->ops->setup = PetscDualSpaceSetUp_Simple; 17120cf1dd8SToby Isaac sp->ops->view = NULL; 17220cf1dd8SToby Isaac sp->ops->destroy = PetscDualSpaceDestroy_Simple; 17320cf1dd8SToby Isaac sp->ops->duplicate = PetscDualSpaceDuplicate_Simple; 17420cf1dd8SToby Isaac sp->ops->getdimension = PetscDualSpaceGetDimension_Simple; 17520cf1dd8SToby Isaac sp->ops->getnumdof = PetscDualSpaceGetNumDof_Simple; 17620cf1dd8SToby Isaac sp->ops->getheightsubspace = NULL; 17720cf1dd8SToby Isaac sp->ops->getsymmetries = NULL; 17820cf1dd8SToby Isaac sp->ops->apply = PetscDualSpaceApplyDefault; 17920cf1dd8SToby Isaac sp->ops->applyall = PetscDualSpaceApplyAllDefault; 18020cf1dd8SToby Isaac sp->ops->createallpoints = PetscDualSpaceCreateAllPointsDefault; 18120cf1dd8SToby Isaac PetscFunctionReturn(0); 18220cf1dd8SToby Isaac } 18320cf1dd8SToby Isaac 18420cf1dd8SToby Isaac /*MC 18520cf1dd8SToby Isaac PETSCDUALSPACESIMPLE = "simple" - A PetscDualSpace object that encapsulates a dual space of arbitrary functionals 18620cf1dd8SToby Isaac 18720cf1dd8SToby Isaac Level: intermediate 18820cf1dd8SToby Isaac 18920cf1dd8SToby Isaac .seealso: PetscDualSpaceType, PetscDualSpaceCreate(), PetscDualSpaceSetType() 19020cf1dd8SToby Isaac M*/ 19120cf1dd8SToby Isaac 19220cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceCreate_Simple(PetscDualSpace sp) 19320cf1dd8SToby Isaac { 19420cf1dd8SToby Isaac PetscDualSpace_Simple *s; 19520cf1dd8SToby Isaac PetscErrorCode ierr; 19620cf1dd8SToby Isaac 19720cf1dd8SToby Isaac PetscFunctionBegin; 19820cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 19920cf1dd8SToby Isaac ierr = PetscNewLog(sp,&s);CHKERRQ(ierr); 20020cf1dd8SToby Isaac sp->data = s; 20120cf1dd8SToby Isaac 20220cf1dd8SToby Isaac s->dim = 0; 20320cf1dd8SToby Isaac s->numDof = NULL; 20420cf1dd8SToby Isaac 20520cf1dd8SToby Isaac ierr = PetscDualSpaceInitialize_Simple(sp);CHKERRQ(ierr); 20620cf1dd8SToby Isaac ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceSimpleSetDimension_C", PetscDualSpaceSimpleSetDimension_Simple);CHKERRQ(ierr); 20720cf1dd8SToby Isaac ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceSimpleSetFunctional_C", PetscDualSpaceSimpleSetFunctional_Simple);CHKERRQ(ierr); 20820cf1dd8SToby Isaac PetscFunctionReturn(0); 20920cf1dd8SToby Isaac } 210