xref: /petsc/src/dm/dt/dualspace/impls/simple/dspacesimple.c (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
120cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
220cf1dd8SToby Isaac #include <petscdmplex.h>
320cf1dd8SToby Isaac 
49371c9d4SSatish Balay static PetscErrorCode PetscDualSpaceSetUp_Simple(PetscDualSpace sp) {
520cf1dd8SToby Isaac   PetscDualSpace_Simple *s  = (PetscDualSpace_Simple *)sp->data;
620cf1dd8SToby Isaac   DM                     dm = sp->dm;
7b4457527SToby Isaac   PetscInt               dim, pStart, pEnd;
8b4457527SToby Isaac   PetscSection           section;
920cf1dd8SToby Isaac 
1020cf1dd8SToby Isaac   PetscFunctionBegin;
119566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
129566063dSJacob Faibussowitsch   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
139566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(PETSC_COMM_SELF, &section));
149566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(section, pStart, pEnd));
159566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetDof(section, pStart, s->dim));
169566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(section));
17b4457527SToby Isaac   sp->pointSection = section;
1820cf1dd8SToby Isaac   PetscFunctionReturn(0);
1920cf1dd8SToby Isaac }
2020cf1dd8SToby Isaac 
219371c9d4SSatish Balay static PetscErrorCode PetscDualSpaceDestroy_Simple(PetscDualSpace sp) {
2220cf1dd8SToby Isaac   PetscDualSpace_Simple *s = (PetscDualSpace_Simple *)sp->data;
2320cf1dd8SToby Isaac 
2420cf1dd8SToby Isaac   PetscFunctionBegin;
259566063dSJacob Faibussowitsch   PetscCall(PetscFree(s->numDof));
269566063dSJacob Faibussowitsch   PetscCall(PetscFree(s));
279566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSimpleSetDimension_C", NULL));
289566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSimpleSetFunctional_C", NULL));
2920cf1dd8SToby Isaac   PetscFunctionReturn(0);
3020cf1dd8SToby Isaac }
3120cf1dd8SToby Isaac 
329371c9d4SSatish Balay static PetscErrorCode PetscDualSpaceDuplicate_Simple(PetscDualSpace sp, PetscDualSpace spNew) {
33b4457527SToby Isaac   PetscInt dim, d;
3420cf1dd8SToby Isaac 
3520cf1dd8SToby Isaac   PetscFunctionBegin;
369566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDimension(sp, &dim));
379566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSimpleSetDimension(spNew, dim));
3820cf1dd8SToby Isaac   for (d = 0; d < dim; ++d) {
3920cf1dd8SToby Isaac     PetscQuadrature q;
4020cf1dd8SToby Isaac 
419566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetFunctional(sp, d, &q));
429566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceSimpleSetFunctional(spNew, d, q));
4320cf1dd8SToby Isaac   }
4420cf1dd8SToby Isaac   PetscFunctionReturn(0);
4520cf1dd8SToby Isaac }
4620cf1dd8SToby Isaac 
479371c9d4SSatish Balay static PetscErrorCode PetscDualSpaceSetFromOptions_Simple(PetscDualSpace sp, PetscOptionItems *PetscOptionsObject) {
4820cf1dd8SToby Isaac   PetscFunctionBegin;
4920cf1dd8SToby Isaac   PetscFunctionReturn(0);
5020cf1dd8SToby Isaac }
5120cf1dd8SToby Isaac 
529371c9d4SSatish Balay static PetscErrorCode PetscDualSpaceSimpleSetDimension_Simple(PetscDualSpace sp, const PetscInt dim) {
5320cf1dd8SToby Isaac   PetscDualSpace_Simple *s = (PetscDualSpace_Simple *)sp->data;
5420cf1dd8SToby Isaac   DM                     dm;
5520cf1dd8SToby Isaac   PetscInt               spatialDim, f;
5620cf1dd8SToby Isaac 
5720cf1dd8SToby Isaac   PetscFunctionBegin;
589566063dSJacob Faibussowitsch   for (f = 0; f < s->dim; ++f) PetscCall(PetscQuadratureDestroy(&sp->functional[f]));
599566063dSJacob Faibussowitsch   PetscCall(PetscFree(sp->functional));
6020cf1dd8SToby Isaac   s->dim = dim;
619566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(s->dim, &sp->functional));
629566063dSJacob Faibussowitsch   PetscCall(PetscFree(s->numDof));
639566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDM(sp, &dm));
649566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &spatialDim));
659566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(spatialDim + 1, &s->numDof));
6620cf1dd8SToby Isaac   s->numDof[spatialDim] = dim;
6720cf1dd8SToby Isaac   PetscFunctionReturn(0);
6820cf1dd8SToby Isaac }
6920cf1dd8SToby Isaac 
709371c9d4SSatish Balay static PetscErrorCode PetscDualSpaceSimpleSetFunctional_Simple(PetscDualSpace sp, PetscInt f, PetscQuadrature q) {
7120cf1dd8SToby Isaac   PetscDualSpace_Simple *s = (PetscDualSpace_Simple *)sp->data;
7220cf1dd8SToby Isaac   PetscReal             *weights;
7320cf1dd8SToby Isaac   PetscInt               Nc, c, Nq, p;
7420cf1dd8SToby Isaac 
7520cf1dd8SToby Isaac   PetscFunctionBegin;
7663a3b9bcSJacob Faibussowitsch   PetscCheck(!(f < 0) && !(f >= s->dim), PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_OUTOFRANGE, "Basis index %" PetscInt_FMT " not in [0, %" PetscInt_FMT ")", f, s->dim);
779566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureDuplicate(q, &sp->functional[f]));
7820cf1dd8SToby Isaac   /* Reweight so that it has unit volume: Do we want to do this for Nc > 1? */
799566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(sp->functional[f], NULL, &Nc, &Nq, NULL, (const PetscReal **)&weights));
8020cf1dd8SToby Isaac   for (c = 0; c < Nc; ++c) {
8120cf1dd8SToby Isaac     PetscReal vol = 0.0;
8220cf1dd8SToby Isaac 
8320cf1dd8SToby Isaac     for (p = 0; p < Nq; ++p) vol += weights[p * Nc + c];
8420cf1dd8SToby Isaac     for (p = 0; p < Nq; ++p) weights[p * Nc + c] /= (vol == 0.0 ? 1.0 : vol);
8520cf1dd8SToby Isaac   }
8620cf1dd8SToby Isaac   PetscFunctionReturn(0);
8720cf1dd8SToby Isaac }
8820cf1dd8SToby Isaac 
8920cf1dd8SToby Isaac /*@
9020cf1dd8SToby Isaac   PetscDualSpaceSimpleSetDimension - Set the number of functionals in the dual space basis
9120cf1dd8SToby Isaac 
92d083f849SBarry Smith   Logically Collective on sp
9320cf1dd8SToby Isaac 
9420cf1dd8SToby Isaac   Input Parameters:
9520cf1dd8SToby Isaac + sp  - the PetscDualSpace
9620cf1dd8SToby Isaac - dim - the basis dimension
9720cf1dd8SToby Isaac 
9820cf1dd8SToby Isaac   Level: intermediate
9920cf1dd8SToby Isaac 
100db781477SPatrick Sanan .seealso: `PetscDualSpaceSimpleSetFunctional()`
10120cf1dd8SToby Isaac @*/
1029371c9d4SSatish Balay PetscErrorCode PetscDualSpaceSimpleSetDimension(PetscDualSpace sp, PetscInt dim) {
10320cf1dd8SToby Isaac   PetscFunctionBegin;
10420cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
10520cf1dd8SToby Isaac   PetscValidLogicalCollectiveInt(sp, dim, 2);
10628b400f6SJacob Faibussowitsch   PetscCheck(!sp->setupcalled, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change dimension after dual space is set up");
107cac4c232SBarry Smith   PetscTryMethod(sp, "PetscDualSpaceSimpleSetDimension_C", (PetscDualSpace, PetscInt), (sp, dim));
10820cf1dd8SToby Isaac   PetscFunctionReturn(0);
10920cf1dd8SToby Isaac }
11020cf1dd8SToby Isaac 
11120cf1dd8SToby Isaac /*@
11220cf1dd8SToby Isaac   PetscDualSpaceSimpleSetFunctional - Set the given basis element for this dual space
11320cf1dd8SToby Isaac 
11420cf1dd8SToby Isaac   Not Collective
11520cf1dd8SToby Isaac 
11620cf1dd8SToby Isaac   Input Parameters:
11720cf1dd8SToby Isaac + sp  - the PetscDualSpace
11820cf1dd8SToby Isaac . f - the basis index
11920cf1dd8SToby Isaac - q - the basis functional
12020cf1dd8SToby Isaac 
12120cf1dd8SToby Isaac   Level: intermediate
12220cf1dd8SToby Isaac 
12320cf1dd8SToby Isaac   Note: The quadrature will be reweighted so that it has unit volume.
12420cf1dd8SToby Isaac 
125db781477SPatrick Sanan .seealso: `PetscDualSpaceSimpleSetDimension()`
12620cf1dd8SToby Isaac @*/
1279371c9d4SSatish Balay PetscErrorCode PetscDualSpaceSimpleSetFunctional(PetscDualSpace sp, PetscInt func, PetscQuadrature q) {
12820cf1dd8SToby Isaac   PetscFunctionBegin;
12920cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
130cac4c232SBarry Smith   PetscTryMethod(sp, "PetscDualSpaceSimpleSetFunctional_C", (PetscDualSpace, PetscInt, PetscQuadrature), (sp, func, q));
13120cf1dd8SToby Isaac   PetscFunctionReturn(0);
13220cf1dd8SToby Isaac }
13320cf1dd8SToby Isaac 
1349371c9d4SSatish Balay static PetscErrorCode PetscDualSpaceInitialize_Simple(PetscDualSpace sp) {
13520cf1dd8SToby Isaac   PetscFunctionBegin;
13620cf1dd8SToby Isaac   sp->ops->setfromoptions       = PetscDualSpaceSetFromOptions_Simple;
13720cf1dd8SToby Isaac   sp->ops->setup                = PetscDualSpaceSetUp_Simple;
13820cf1dd8SToby Isaac   sp->ops->view                 = NULL;
13920cf1dd8SToby Isaac   sp->ops->destroy              = PetscDualSpaceDestroy_Simple;
14020cf1dd8SToby Isaac   sp->ops->duplicate            = PetscDualSpaceDuplicate_Simple;
141b4457527SToby Isaac   sp->ops->createheightsubspace = NULL;
142b4457527SToby Isaac   sp->ops->createpointsubspace  = NULL;
14320cf1dd8SToby Isaac   sp->ops->getsymmetries        = NULL;
14420cf1dd8SToby Isaac   sp->ops->apply                = PetscDualSpaceApplyDefault;
14520cf1dd8SToby Isaac   sp->ops->applyall             = PetscDualSpaceApplyAllDefault;
146b4457527SToby Isaac   sp->ops->applyint             = PetscDualSpaceApplyInteriorDefault;
147b4457527SToby Isaac   sp->ops->createalldata        = PetscDualSpaceCreateAllDataDefault;
148b4457527SToby Isaac   sp->ops->createintdata        = PetscDualSpaceCreateInteriorDataDefault;
14920cf1dd8SToby Isaac   PetscFunctionReturn(0);
15020cf1dd8SToby Isaac }
15120cf1dd8SToby Isaac 
15220cf1dd8SToby Isaac /*MC
15320cf1dd8SToby Isaac   PETSCDUALSPACESIMPLE = "simple" - A PetscDualSpace object that encapsulates a dual space of arbitrary functionals
15420cf1dd8SToby Isaac 
15520cf1dd8SToby Isaac   Level: intermediate
15620cf1dd8SToby Isaac 
157db781477SPatrick Sanan .seealso: `PetscDualSpaceType`, `PetscDualSpaceCreate()`, `PetscDualSpaceSetType()`
15820cf1dd8SToby Isaac M*/
15920cf1dd8SToby Isaac 
1609371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PetscDualSpaceCreate_Simple(PetscDualSpace sp) {
16120cf1dd8SToby Isaac   PetscDualSpace_Simple *s;
16220cf1dd8SToby Isaac 
16320cf1dd8SToby Isaac   PetscFunctionBegin;
16420cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
165*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&s));
16620cf1dd8SToby Isaac   sp->data = s;
16720cf1dd8SToby Isaac 
16820cf1dd8SToby Isaac   s->dim    = 0;
16920cf1dd8SToby Isaac   s->numDof = NULL;
17020cf1dd8SToby Isaac 
1719566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceInitialize_Simple(sp));
1729566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSimpleSetDimension_C", PetscDualSpaceSimpleSetDimension_Simple));
1739566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSimpleSetFunctional_C", PetscDualSpaceSimpleSetFunctional_Simple));
17420cf1dd8SToby Isaac   PetscFunctionReturn(0);
17520cf1dd8SToby Isaac }
176