xref: /petsc/src/dm/dt/dualspace/impls/simple/dspacesimple.c (revision dce8aeba1c9b69b19f651c53d8a6b674bd7e9cbd)
120cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
220cf1dd8SToby Isaac #include <petscdmplex.h>
320cf1dd8SToby Isaac 
4d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDualSpaceSetUp_Simple(PetscDualSpace sp)
5d71ae5a4SJacob Faibussowitsch {
620cf1dd8SToby Isaac   PetscDualSpace_Simple *s  = (PetscDualSpace_Simple *)sp->data;
720cf1dd8SToby Isaac   DM                     dm = sp->dm;
8b4457527SToby Isaac   PetscInt               dim, pStart, pEnd;
9b4457527SToby Isaac   PetscSection           section;
1020cf1dd8SToby Isaac 
1120cf1dd8SToby Isaac   PetscFunctionBegin;
129566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
139566063dSJacob Faibussowitsch   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
149566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(PETSC_COMM_SELF, &section));
159566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(section, pStart, pEnd));
169566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetDof(section, pStart, s->dim));
179566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(section));
18b4457527SToby Isaac   sp->pointSection = section;
1920cf1dd8SToby Isaac   PetscFunctionReturn(0);
2020cf1dd8SToby Isaac }
2120cf1dd8SToby Isaac 
22d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDualSpaceDestroy_Simple(PetscDualSpace sp)
23d71ae5a4SJacob Faibussowitsch {
2420cf1dd8SToby Isaac   PetscDualSpace_Simple *s = (PetscDualSpace_Simple *)sp->data;
2520cf1dd8SToby Isaac 
2620cf1dd8SToby Isaac   PetscFunctionBegin;
279566063dSJacob Faibussowitsch   PetscCall(PetscFree(s->numDof));
289566063dSJacob Faibussowitsch   PetscCall(PetscFree(s));
299566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSimpleSetDimension_C", NULL));
309566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSimpleSetFunctional_C", NULL));
3120cf1dd8SToby Isaac   PetscFunctionReturn(0);
3220cf1dd8SToby Isaac }
3320cf1dd8SToby Isaac 
34d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDualSpaceDuplicate_Simple(PetscDualSpace sp, PetscDualSpace spNew)
35d71ae5a4SJacob Faibussowitsch {
36b4457527SToby Isaac   PetscInt dim, d;
3720cf1dd8SToby Isaac 
3820cf1dd8SToby Isaac   PetscFunctionBegin;
399566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDimension(sp, &dim));
409566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSimpleSetDimension(spNew, dim));
4120cf1dd8SToby Isaac   for (d = 0; d < dim; ++d) {
4220cf1dd8SToby Isaac     PetscQuadrature q;
4320cf1dd8SToby Isaac 
449566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetFunctional(sp, d, &q));
459566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceSimpleSetFunctional(spNew, d, q));
4620cf1dd8SToby Isaac   }
4720cf1dd8SToby Isaac   PetscFunctionReturn(0);
4820cf1dd8SToby Isaac }
4920cf1dd8SToby Isaac 
50d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDualSpaceSetFromOptions_Simple(PetscDualSpace sp, PetscOptionItems *PetscOptionsObject)
51d71ae5a4SJacob Faibussowitsch {
5220cf1dd8SToby Isaac   PetscFunctionBegin;
5320cf1dd8SToby Isaac   PetscFunctionReturn(0);
5420cf1dd8SToby Isaac }
5520cf1dd8SToby Isaac 
56d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDualSpaceSimpleSetDimension_Simple(PetscDualSpace sp, const PetscInt dim)
57d71ae5a4SJacob Faibussowitsch {
5820cf1dd8SToby Isaac   PetscDualSpace_Simple *s = (PetscDualSpace_Simple *)sp->data;
5920cf1dd8SToby Isaac   DM                     dm;
6020cf1dd8SToby Isaac   PetscInt               spatialDim, f;
6120cf1dd8SToby Isaac 
6220cf1dd8SToby Isaac   PetscFunctionBegin;
639566063dSJacob Faibussowitsch   for (f = 0; f < s->dim; ++f) PetscCall(PetscQuadratureDestroy(&sp->functional[f]));
649566063dSJacob Faibussowitsch   PetscCall(PetscFree(sp->functional));
6520cf1dd8SToby Isaac   s->dim = dim;
669566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(s->dim, &sp->functional));
679566063dSJacob Faibussowitsch   PetscCall(PetscFree(s->numDof));
689566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDM(sp, &dm));
699566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &spatialDim));
709566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(spatialDim + 1, &s->numDof));
7120cf1dd8SToby Isaac   s->numDof[spatialDim] = dim;
7220cf1dd8SToby Isaac   PetscFunctionReturn(0);
7320cf1dd8SToby Isaac }
7420cf1dd8SToby Isaac 
75d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDualSpaceSimpleSetFunctional_Simple(PetscDualSpace sp, PetscInt f, PetscQuadrature q)
76d71ae5a4SJacob Faibussowitsch {
7720cf1dd8SToby Isaac   PetscDualSpace_Simple *s = (PetscDualSpace_Simple *)sp->data;
7820cf1dd8SToby Isaac   PetscReal             *weights;
7920cf1dd8SToby Isaac   PetscInt               Nc, c, Nq, p;
8020cf1dd8SToby Isaac 
8120cf1dd8SToby Isaac   PetscFunctionBegin;
8263a3b9bcSJacob 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);
839566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureDuplicate(q, &sp->functional[f]));
8420cf1dd8SToby Isaac   /* Reweight so that it has unit volume: Do we want to do this for Nc > 1? */
859566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(sp->functional[f], NULL, &Nc, &Nq, NULL, (const PetscReal **)&weights));
8620cf1dd8SToby Isaac   for (c = 0; c < Nc; ++c) {
8720cf1dd8SToby Isaac     PetscReal vol = 0.0;
8820cf1dd8SToby Isaac 
8920cf1dd8SToby Isaac     for (p = 0; p < Nq; ++p) vol += weights[p * Nc + c];
9020cf1dd8SToby Isaac     for (p = 0; p < Nq; ++p) weights[p * Nc + c] /= (vol == 0.0 ? 1.0 : vol);
9120cf1dd8SToby Isaac   }
9220cf1dd8SToby Isaac   PetscFunctionReturn(0);
9320cf1dd8SToby Isaac }
9420cf1dd8SToby Isaac 
9520cf1dd8SToby Isaac /*@
9620cf1dd8SToby Isaac   PetscDualSpaceSimpleSetDimension - Set the number of functionals in the dual space basis
9720cf1dd8SToby Isaac 
98d083f849SBarry Smith   Logically Collective on sp
9920cf1dd8SToby Isaac 
10020cf1dd8SToby Isaac   Input Parameters:
101*dce8aebaSBarry Smith + sp  - the `PetscDualSpace`
10220cf1dd8SToby Isaac - dim - the basis dimension
10320cf1dd8SToby Isaac 
10420cf1dd8SToby Isaac   Level: intermediate
10520cf1dd8SToby Isaac 
106*dce8aebaSBarry Smith .seealso: `PetscDualSpace`, `PetscDualSpaceSimpleSetFunctional()`
10720cf1dd8SToby Isaac @*/
108d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDualSpaceSimpleSetDimension(PetscDualSpace sp, PetscInt dim)
109d71ae5a4SJacob Faibussowitsch {
11020cf1dd8SToby Isaac   PetscFunctionBegin;
11120cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
11220cf1dd8SToby Isaac   PetscValidLogicalCollectiveInt(sp, dim, 2);
11328b400f6SJacob Faibussowitsch   PetscCheck(!sp->setupcalled, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change dimension after dual space is set up");
114cac4c232SBarry Smith   PetscTryMethod(sp, "PetscDualSpaceSimpleSetDimension_C", (PetscDualSpace, PetscInt), (sp, dim));
11520cf1dd8SToby Isaac   PetscFunctionReturn(0);
11620cf1dd8SToby Isaac }
11720cf1dd8SToby Isaac 
11820cf1dd8SToby Isaac /*@
11920cf1dd8SToby Isaac   PetscDualSpaceSimpleSetFunctional - Set the given basis element for this dual space
12020cf1dd8SToby Isaac 
12120cf1dd8SToby Isaac   Not Collective
12220cf1dd8SToby Isaac 
12320cf1dd8SToby Isaac   Input Parameters:
124*dce8aebaSBarry Smith + sp  - the `PetscDualSpace`
12520cf1dd8SToby Isaac . f - the basis index
12620cf1dd8SToby Isaac - q - the basis functional
12720cf1dd8SToby Isaac 
12820cf1dd8SToby Isaac   Level: intermediate
12920cf1dd8SToby Isaac 
130*dce8aebaSBarry Smith   Note:
131*dce8aebaSBarry Smith   The quadrature will be reweighted so that it has unit volume.
13220cf1dd8SToby Isaac 
133*dce8aebaSBarry Smith .seealso: `PetscDualSpace`, `PetscDualSpaceSimpleSetDimension()`
13420cf1dd8SToby Isaac @*/
135d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDualSpaceSimpleSetFunctional(PetscDualSpace sp, PetscInt func, PetscQuadrature q)
136d71ae5a4SJacob Faibussowitsch {
13720cf1dd8SToby Isaac   PetscFunctionBegin;
13820cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
139cac4c232SBarry Smith   PetscTryMethod(sp, "PetscDualSpaceSimpleSetFunctional_C", (PetscDualSpace, PetscInt, PetscQuadrature), (sp, func, q));
14020cf1dd8SToby Isaac   PetscFunctionReturn(0);
14120cf1dd8SToby Isaac }
14220cf1dd8SToby Isaac 
143d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDualSpaceInitialize_Simple(PetscDualSpace sp)
144d71ae5a4SJacob Faibussowitsch {
14520cf1dd8SToby Isaac   PetscFunctionBegin;
14620cf1dd8SToby Isaac   sp->ops->setfromoptions       = PetscDualSpaceSetFromOptions_Simple;
14720cf1dd8SToby Isaac   sp->ops->setup                = PetscDualSpaceSetUp_Simple;
14820cf1dd8SToby Isaac   sp->ops->view                 = NULL;
14920cf1dd8SToby Isaac   sp->ops->destroy              = PetscDualSpaceDestroy_Simple;
15020cf1dd8SToby Isaac   sp->ops->duplicate            = PetscDualSpaceDuplicate_Simple;
151b4457527SToby Isaac   sp->ops->createheightsubspace = NULL;
152b4457527SToby Isaac   sp->ops->createpointsubspace  = NULL;
15320cf1dd8SToby Isaac   sp->ops->getsymmetries        = NULL;
15420cf1dd8SToby Isaac   sp->ops->apply                = PetscDualSpaceApplyDefault;
15520cf1dd8SToby Isaac   sp->ops->applyall             = PetscDualSpaceApplyAllDefault;
156b4457527SToby Isaac   sp->ops->applyint             = PetscDualSpaceApplyInteriorDefault;
157b4457527SToby Isaac   sp->ops->createalldata        = PetscDualSpaceCreateAllDataDefault;
158b4457527SToby Isaac   sp->ops->createintdata        = PetscDualSpaceCreateInteriorDataDefault;
15920cf1dd8SToby Isaac   PetscFunctionReturn(0);
16020cf1dd8SToby Isaac }
16120cf1dd8SToby Isaac 
16220cf1dd8SToby Isaac /*MC
163*dce8aebaSBarry Smith   PETSCDUALSPACESIMPLE = "simple" - A `PetscDualSpaceType` that encapsulates a dual space of arbitrary functionals
16420cf1dd8SToby Isaac 
16520cf1dd8SToby Isaac   Level: intermediate
16620cf1dd8SToby Isaac 
167*dce8aebaSBarry Smith .seealso: `PetscDualSpace`, `PetscDualSpaceType`, `PetscDualSpaceCreate()`, `PetscDualSpaceSetType()`
16820cf1dd8SToby Isaac M*/
16920cf1dd8SToby Isaac 
170d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscDualSpaceCreate_Simple(PetscDualSpace sp)
171d71ae5a4SJacob Faibussowitsch {
17220cf1dd8SToby Isaac   PetscDualSpace_Simple *s;
17320cf1dd8SToby Isaac 
17420cf1dd8SToby Isaac   PetscFunctionBegin;
17520cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1764dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&s));
17720cf1dd8SToby Isaac   sp->data = s;
17820cf1dd8SToby Isaac 
17920cf1dd8SToby Isaac   s->dim    = 0;
18020cf1dd8SToby Isaac   s->numDof = NULL;
18120cf1dd8SToby Isaac 
1829566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceInitialize_Simple(sp));
1839566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSimpleSetDimension_C", PetscDualSpaceSimpleSetDimension_Simple));
1849566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSimpleSetFunctional_C", PetscDualSpaceSimpleSetFunctional_Simple));
18520cf1dd8SToby Isaac   PetscFunctionReturn(0);
18620cf1dd8SToby Isaac }
187