120cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/ 220cf1dd8SToby Isaac #include <petscdmplex.h> 320cf1dd8SToby Isaac 4a4ce7ad1SMatthew 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; 8*b4457527SToby Isaac PetscInt dim, pStart, pEnd; 9*b4457527SToby Isaac PetscSection section; 1020cf1dd8SToby Isaac PetscErrorCode ierr; 1120cf1dd8SToby Isaac 1220cf1dd8SToby Isaac PetscFunctionBegin; 1320cf1dd8SToby Isaac ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 14*b4457527SToby Isaac ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 15*b4457527SToby Isaac ierr = PetscSectionCreate(PETSC_COMM_SELF, §ion);CHKERRQ(ierr); 16*b4457527SToby Isaac ierr = PetscSectionSetChart(section, pStart, pEnd);CHKERRQ(ierr); 17*b4457527SToby Isaac ierr = PetscSectionSetDof(section, pStart, s->dim);CHKERRQ(ierr); 18*b4457527SToby Isaac ierr = PetscSectionSetUp(section);CHKERRQ(ierr); 19*b4457527SToby Isaac sp->pointSection = section; 2020cf1dd8SToby Isaac PetscFunctionReturn(0); 2120cf1dd8SToby Isaac } 2220cf1dd8SToby Isaac 23a4ce7ad1SMatthew G. Knepley static PetscErrorCode PetscDualSpaceDestroy_Simple(PetscDualSpace sp) 2420cf1dd8SToby Isaac { 2520cf1dd8SToby Isaac PetscDualSpace_Simple *s = (PetscDualSpace_Simple *) sp->data; 2620cf1dd8SToby Isaac PetscErrorCode ierr; 2720cf1dd8SToby Isaac 2820cf1dd8SToby Isaac PetscFunctionBegin; 2920cf1dd8SToby Isaac ierr = PetscFree(s->numDof);CHKERRQ(ierr); 3020cf1dd8SToby Isaac ierr = PetscFree(s);CHKERRQ(ierr); 3120cf1dd8SToby Isaac ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceSimpleSetDimension_C", NULL);CHKERRQ(ierr); 3220cf1dd8SToby Isaac ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceSimpleSetFunctional_C", NULL);CHKERRQ(ierr); 3320cf1dd8SToby Isaac PetscFunctionReturn(0); 3420cf1dd8SToby Isaac } 3520cf1dd8SToby Isaac 36*b4457527SToby Isaac static PetscErrorCode PetscDualSpaceDuplicate_Simple(PetscDualSpace sp, PetscDualSpace spNew) 3720cf1dd8SToby Isaac { 38*b4457527SToby Isaac PetscInt dim, d; 3920cf1dd8SToby Isaac PetscErrorCode ierr; 4020cf1dd8SToby Isaac 4120cf1dd8SToby Isaac PetscFunctionBegin; 4220cf1dd8SToby Isaac ierr = PetscDualSpaceGetDimension(sp, &dim);CHKERRQ(ierr); 43*b4457527SToby Isaac ierr = PetscDualSpaceSimpleSetDimension(spNew, dim);CHKERRQ(ierr); 4420cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 4520cf1dd8SToby Isaac PetscQuadrature q; 4620cf1dd8SToby Isaac 4720cf1dd8SToby Isaac ierr = PetscDualSpaceGetFunctional(sp, d, &q);CHKERRQ(ierr); 48*b4457527SToby Isaac ierr = PetscDualSpaceSimpleSetFunctional(spNew, d, q);CHKERRQ(ierr); 4920cf1dd8SToby Isaac } 5020cf1dd8SToby Isaac PetscFunctionReturn(0); 5120cf1dd8SToby Isaac } 5220cf1dd8SToby Isaac 53a4ce7ad1SMatthew G. Knepley static PetscErrorCode PetscDualSpaceSetFromOptions_Simple(PetscOptionItems *PetscOptionsObject,PetscDualSpace sp) 5420cf1dd8SToby Isaac { 5520cf1dd8SToby Isaac PetscFunctionBegin; 5620cf1dd8SToby Isaac PetscFunctionReturn(0); 5720cf1dd8SToby Isaac } 5820cf1dd8SToby Isaac 59a4ce7ad1SMatthew G. Knepley static PetscErrorCode PetscDualSpaceSimpleSetDimension_Simple(PetscDualSpace sp, const PetscInt dim) 6020cf1dd8SToby Isaac { 6120cf1dd8SToby Isaac PetscDualSpace_Simple *s = (PetscDualSpace_Simple *) sp->data; 6220cf1dd8SToby Isaac DM dm; 6320cf1dd8SToby Isaac PetscInt spatialDim, f; 6420cf1dd8SToby Isaac PetscErrorCode ierr; 6520cf1dd8SToby Isaac 6620cf1dd8SToby Isaac PetscFunctionBegin; 6720cf1dd8SToby Isaac for (f = 0; f < s->dim; ++f) {ierr = PetscQuadratureDestroy(&sp->functional[f]);CHKERRQ(ierr);} 6820cf1dd8SToby Isaac ierr = PetscFree(sp->functional);CHKERRQ(ierr); 6920cf1dd8SToby Isaac s->dim = dim; 7020cf1dd8SToby Isaac ierr = PetscCalloc1(s->dim, &sp->functional);CHKERRQ(ierr); 7120cf1dd8SToby Isaac ierr = PetscFree(s->numDof);CHKERRQ(ierr); 7220cf1dd8SToby Isaac ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr); 7320cf1dd8SToby Isaac ierr = DMGetCoordinateDim(dm, &spatialDim);CHKERRQ(ierr); 7420cf1dd8SToby Isaac ierr = PetscCalloc1(spatialDim+1, &s->numDof);CHKERRQ(ierr); 7520cf1dd8SToby Isaac s->numDof[spatialDim] = dim; 7620cf1dd8SToby Isaac PetscFunctionReturn(0); 7720cf1dd8SToby Isaac } 7820cf1dd8SToby Isaac 79a4ce7ad1SMatthew G. Knepley static PetscErrorCode PetscDualSpaceSimpleSetFunctional_Simple(PetscDualSpace sp, PetscInt f, PetscQuadrature q) 8020cf1dd8SToby Isaac { 8120cf1dd8SToby Isaac PetscDualSpace_Simple *s = (PetscDualSpace_Simple *) sp->data; 8220cf1dd8SToby Isaac PetscReal *weights; 8320cf1dd8SToby Isaac PetscInt Nc, c, Nq, p; 8420cf1dd8SToby Isaac PetscErrorCode ierr; 8520cf1dd8SToby Isaac 8620cf1dd8SToby Isaac PetscFunctionBegin; 8720cf1dd8SToby 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); 8820cf1dd8SToby Isaac ierr = PetscQuadratureDuplicate(q, &sp->functional[f]);CHKERRQ(ierr); 8920cf1dd8SToby Isaac /* Reweight so that it has unit volume: Do we want to do this for Nc > 1? */ 9020cf1dd8SToby Isaac ierr = PetscQuadratureGetData(sp->functional[f], NULL, &Nc, &Nq, NULL, (const PetscReal **) &weights);CHKERRQ(ierr); 9120cf1dd8SToby Isaac for (c = 0; c < Nc; ++c) { 9220cf1dd8SToby Isaac PetscReal vol = 0.0; 9320cf1dd8SToby Isaac 9420cf1dd8SToby Isaac for (p = 0; p < Nq; ++p) vol += weights[p*Nc+c]; 9520cf1dd8SToby Isaac for (p = 0; p < Nq; ++p) weights[p*Nc+c] /= (vol == 0.0 ? 1.0 : vol); 9620cf1dd8SToby Isaac } 9720cf1dd8SToby Isaac PetscFunctionReturn(0); 9820cf1dd8SToby Isaac } 9920cf1dd8SToby Isaac 10020cf1dd8SToby Isaac /*@ 10120cf1dd8SToby Isaac PetscDualSpaceSimpleSetDimension - Set the number of functionals in the dual space basis 10220cf1dd8SToby Isaac 103d083f849SBarry Smith Logically Collective on sp 10420cf1dd8SToby Isaac 10520cf1dd8SToby Isaac Input Parameters: 10620cf1dd8SToby Isaac + sp - the PetscDualSpace 10720cf1dd8SToby Isaac - dim - the basis dimension 10820cf1dd8SToby Isaac 10920cf1dd8SToby Isaac Level: intermediate 11020cf1dd8SToby Isaac 11120cf1dd8SToby Isaac .seealso: PetscDualSpaceSimpleSetFunctional() 11220cf1dd8SToby Isaac @*/ 11320cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSimpleSetDimension(PetscDualSpace sp, PetscInt dim) 11420cf1dd8SToby Isaac { 11520cf1dd8SToby Isaac PetscErrorCode ierr; 11620cf1dd8SToby Isaac 11720cf1dd8SToby Isaac PetscFunctionBegin; 11820cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 11920cf1dd8SToby Isaac PetscValidLogicalCollectiveInt(sp, dim, 2); 120*b4457527SToby Isaac if (sp->setupcalled) SETERRQ(PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change dimension after dual space is set up"); 12120cf1dd8SToby Isaac ierr = PetscTryMethod(sp, "PetscDualSpaceSimpleSetDimension_C", (PetscDualSpace,PetscInt),(sp,dim));CHKERRQ(ierr); 12220cf1dd8SToby Isaac PetscFunctionReturn(0); 12320cf1dd8SToby Isaac } 12420cf1dd8SToby Isaac 12520cf1dd8SToby Isaac /*@ 12620cf1dd8SToby Isaac PetscDualSpaceSimpleSetFunctional - Set the given basis element for this dual space 12720cf1dd8SToby Isaac 12820cf1dd8SToby Isaac Not Collective 12920cf1dd8SToby Isaac 13020cf1dd8SToby Isaac Input Parameters: 13120cf1dd8SToby Isaac + sp - the PetscDualSpace 13220cf1dd8SToby Isaac . f - the basis index 13320cf1dd8SToby Isaac - q - the basis functional 13420cf1dd8SToby Isaac 13520cf1dd8SToby Isaac Level: intermediate 13620cf1dd8SToby Isaac 13720cf1dd8SToby Isaac Note: The quadrature will be reweighted so that it has unit volume. 13820cf1dd8SToby Isaac 13920cf1dd8SToby Isaac .seealso: PetscDualSpaceSimpleSetDimension() 14020cf1dd8SToby Isaac @*/ 14120cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSimpleSetFunctional(PetscDualSpace sp, PetscInt func, PetscQuadrature q) 14220cf1dd8SToby Isaac { 14320cf1dd8SToby Isaac PetscErrorCode ierr; 14420cf1dd8SToby Isaac 14520cf1dd8SToby Isaac PetscFunctionBegin; 14620cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 14720cf1dd8SToby Isaac ierr = PetscTryMethod(sp, "PetscDualSpaceSimpleSetFunctional_C", (PetscDualSpace,PetscInt,PetscQuadrature),(sp,func,q));CHKERRQ(ierr); 14820cf1dd8SToby Isaac PetscFunctionReturn(0); 14920cf1dd8SToby Isaac } 15020cf1dd8SToby Isaac 151a4ce7ad1SMatthew G. Knepley static PetscErrorCode PetscDualSpaceInitialize_Simple(PetscDualSpace sp) 15220cf1dd8SToby Isaac { 15320cf1dd8SToby Isaac PetscFunctionBegin; 15420cf1dd8SToby Isaac sp->ops->setfromoptions = PetscDualSpaceSetFromOptions_Simple; 15520cf1dd8SToby Isaac sp->ops->setup = PetscDualSpaceSetUp_Simple; 15620cf1dd8SToby Isaac sp->ops->view = NULL; 15720cf1dd8SToby Isaac sp->ops->destroy = PetscDualSpaceDestroy_Simple; 15820cf1dd8SToby Isaac sp->ops->duplicate = PetscDualSpaceDuplicate_Simple; 159*b4457527SToby Isaac sp->ops->createheightsubspace = NULL; 160*b4457527SToby Isaac sp->ops->createpointsubspace = NULL; 16120cf1dd8SToby Isaac sp->ops->getsymmetries = NULL; 16220cf1dd8SToby Isaac sp->ops->apply = PetscDualSpaceApplyDefault; 16320cf1dd8SToby Isaac sp->ops->applyall = PetscDualSpaceApplyAllDefault; 164*b4457527SToby Isaac sp->ops->applyint = PetscDualSpaceApplyInteriorDefault; 165*b4457527SToby Isaac sp->ops->createalldata = PetscDualSpaceCreateAllDataDefault; 166*b4457527SToby Isaac sp->ops->createintdata = PetscDualSpaceCreateInteriorDataDefault; 16720cf1dd8SToby Isaac PetscFunctionReturn(0); 16820cf1dd8SToby Isaac } 16920cf1dd8SToby Isaac 17020cf1dd8SToby Isaac /*MC 17120cf1dd8SToby Isaac PETSCDUALSPACESIMPLE = "simple" - A PetscDualSpace object that encapsulates a dual space of arbitrary functionals 17220cf1dd8SToby Isaac 17320cf1dd8SToby Isaac Level: intermediate 17420cf1dd8SToby Isaac 17520cf1dd8SToby Isaac .seealso: PetscDualSpaceType, PetscDualSpaceCreate(), PetscDualSpaceSetType() 17620cf1dd8SToby Isaac M*/ 17720cf1dd8SToby Isaac 17820cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceCreate_Simple(PetscDualSpace sp) 17920cf1dd8SToby Isaac { 18020cf1dd8SToby Isaac PetscDualSpace_Simple *s; 18120cf1dd8SToby Isaac PetscErrorCode ierr; 18220cf1dd8SToby Isaac 18320cf1dd8SToby Isaac PetscFunctionBegin; 18420cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 18520cf1dd8SToby Isaac ierr = PetscNewLog(sp,&s);CHKERRQ(ierr); 18620cf1dd8SToby Isaac sp->data = s; 18720cf1dd8SToby Isaac 18820cf1dd8SToby Isaac s->dim = 0; 18920cf1dd8SToby Isaac s->numDof = NULL; 19020cf1dd8SToby Isaac 19120cf1dd8SToby Isaac ierr = PetscDualSpaceInitialize_Simple(sp);CHKERRQ(ierr); 19220cf1dd8SToby Isaac ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceSimpleSetDimension_C", PetscDualSpaceSimpleSetDimension_Simple);CHKERRQ(ierr); 19320cf1dd8SToby Isaac ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceSimpleSetFunctional_C", PetscDualSpaceSimpleSetFunctional_Simple);CHKERRQ(ierr); 19420cf1dd8SToby Isaac PetscFunctionReturn(0); 19520cf1dd8SToby Isaac } 196