xref: /petsc/src/dm/dt/dualspace/impls/simple/dspacesimple.c (revision 813083eb7473ec21f4ff50a178709600a8f71c4c)
120cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
220cf1dd8SToby Isaac #include <petscdmplex.h>
320cf1dd8SToby Isaac 
PetscDualSpaceSetUp_Simple(PetscDualSpace sp)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;
193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2020cf1dd8SToby Isaac }
2120cf1dd8SToby Isaac 
PetscDualSpaceDestroy_Simple(PetscDualSpace sp)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));
313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3220cf1dd8SToby Isaac }
3320cf1dd8SToby Isaac 
PetscDualSpaceDuplicate_Simple(PetscDualSpace sp,PetscDualSpace spNew)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   }
473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4820cf1dd8SToby Isaac }
4920cf1dd8SToby Isaac 
PetscDualSpaceSimpleSetDimension_Simple(PetscDualSpace sp,const PetscInt dim)50d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDualSpaceSimpleSetDimension_Simple(PetscDualSpace sp, const PetscInt dim)
51d71ae5a4SJacob Faibussowitsch {
5220cf1dd8SToby Isaac   PetscDualSpace_Simple *s = (PetscDualSpace_Simple *)sp->data;
5320cf1dd8SToby Isaac   DM                     dm;
5420cf1dd8SToby Isaac   PetscInt               spatialDim, f;
5520cf1dd8SToby Isaac 
5620cf1dd8SToby Isaac   PetscFunctionBegin;
579566063dSJacob Faibussowitsch   for (f = 0; f < s->dim; ++f) PetscCall(PetscQuadratureDestroy(&sp->functional[f]));
589566063dSJacob Faibussowitsch   PetscCall(PetscFree(sp->functional));
5920cf1dd8SToby Isaac   s->dim = dim;
609566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(s->dim, &sp->functional));
619566063dSJacob Faibussowitsch   PetscCall(PetscFree(s->numDof));
629566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDM(sp, &dm));
639566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &spatialDim));
649566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(spatialDim + 1, &s->numDof));
6520cf1dd8SToby Isaac   s->numDof[spatialDim] = dim;
663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6720cf1dd8SToby Isaac }
6820cf1dd8SToby Isaac 
PetscDualSpaceSimpleSetFunctional_Simple(PetscDualSpace sp,PetscInt f,PetscQuadrature q)69d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDualSpaceSimpleSetFunctional_Simple(PetscDualSpace sp, PetscInt f, PetscQuadrature q)
70d71ae5a4SJacob Faibussowitsch {
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   }
863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8720cf1dd8SToby Isaac }
8820cf1dd8SToby Isaac 
8920cf1dd8SToby Isaac /*@
9020cf1dd8SToby Isaac   PetscDualSpaceSimpleSetDimension - Set the number of functionals in the dual space basis
9120cf1dd8SToby Isaac 
9220f4b53cSBarry Smith   Logically Collective
9320cf1dd8SToby Isaac 
9420cf1dd8SToby Isaac   Input Parameters:
95dce8aebaSBarry Smith + sp  - the `PetscDualSpace`
9620cf1dd8SToby Isaac - dim - the basis dimension
9720cf1dd8SToby Isaac 
9820cf1dd8SToby Isaac   Level: intermediate
9920cf1dd8SToby Isaac 
100*b24fb147SBarry Smith .seealso: `PETSCDUALSPACESIMPLE`, `PetscDualSpace`, `PetscDualSpaceSimpleSetFunctional()`
10120cf1dd8SToby Isaac @*/
PetscDualSpaceSimpleSetDimension(PetscDualSpace sp,PetscInt dim)102d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDualSpaceSimpleSetDimension(PetscDualSpace sp, PetscInt dim)
103d71ae5a4SJacob Faibussowitsch {
10420cf1dd8SToby Isaac   PetscFunctionBegin;
10520cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
10620cf1dd8SToby Isaac   PetscValidLogicalCollectiveInt(sp, dim, 2);
10728b400f6SJacob Faibussowitsch   PetscCheck(!sp->setupcalled, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change dimension after dual space is set up");
108cac4c232SBarry Smith   PetscTryMethod(sp, "PetscDualSpaceSimpleSetDimension_C", (PetscDualSpace, PetscInt), (sp, dim));
1093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11020cf1dd8SToby Isaac }
11120cf1dd8SToby Isaac 
11220cf1dd8SToby Isaac /*@
113*b24fb147SBarry Smith   PetscDualSpaceSimpleSetFunctional - Set the given basis functional for this dual space
11420cf1dd8SToby Isaac 
11520cf1dd8SToby Isaac   Not Collective
11620cf1dd8SToby Isaac 
11720cf1dd8SToby Isaac   Input Parameters:
118dce8aebaSBarry Smith + sp   - the `PetscDualSpace`
119*b24fb147SBarry Smith . func - the basis index
12020cf1dd8SToby Isaac - q    - the basis functional
12120cf1dd8SToby Isaac 
12220cf1dd8SToby Isaac   Level: intermediate
12320cf1dd8SToby Isaac 
124dce8aebaSBarry Smith   Note:
125dce8aebaSBarry Smith   The quadrature will be reweighted so that it has unit volume.
12620cf1dd8SToby Isaac 
127*b24fb147SBarry Smith .seealso: `PETSCDUALSPACESIMPLE`, `PetscDualSpace`, `PetscDualSpaceSimpleSetDimension()`
12820cf1dd8SToby Isaac @*/
PetscDualSpaceSimpleSetFunctional(PetscDualSpace sp,PetscInt func,PetscQuadrature q)129d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDualSpaceSimpleSetFunctional(PetscDualSpace sp, PetscInt func, PetscQuadrature q)
130d71ae5a4SJacob Faibussowitsch {
13120cf1dd8SToby Isaac   PetscFunctionBegin;
13220cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
133cac4c232SBarry Smith   PetscTryMethod(sp, "PetscDualSpaceSimpleSetFunctional_C", (PetscDualSpace, PetscInt, PetscQuadrature), (sp, func, q));
1343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13520cf1dd8SToby Isaac }
13620cf1dd8SToby Isaac 
PetscDualSpaceInitialize_Simple(PetscDualSpace sp)137d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscDualSpaceInitialize_Simple(PetscDualSpace sp)
138d71ae5a4SJacob Faibussowitsch {
13920cf1dd8SToby Isaac   PetscFunctionBegin;
14020cf1dd8SToby Isaac   sp->ops->setup                = PetscDualSpaceSetUp_Simple;
14120cf1dd8SToby Isaac   sp->ops->view                 = NULL;
14220cf1dd8SToby Isaac   sp->ops->destroy              = PetscDualSpaceDestroy_Simple;
14320cf1dd8SToby Isaac   sp->ops->duplicate            = PetscDualSpaceDuplicate_Simple;
144b4457527SToby Isaac   sp->ops->createheightsubspace = NULL;
145b4457527SToby Isaac   sp->ops->createpointsubspace  = NULL;
14620cf1dd8SToby Isaac   sp->ops->getsymmetries        = NULL;
14720cf1dd8SToby Isaac   sp->ops->apply                = PetscDualSpaceApplyDefault;
14820cf1dd8SToby Isaac   sp->ops->applyall             = PetscDualSpaceApplyAllDefault;
149b4457527SToby Isaac   sp->ops->applyint             = PetscDualSpaceApplyInteriorDefault;
150b4457527SToby Isaac   sp->ops->createalldata        = PetscDualSpaceCreateAllDataDefault;
151b4457527SToby Isaac   sp->ops->createintdata        = PetscDualSpaceCreateInteriorDataDefault;
1523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15320cf1dd8SToby Isaac }
15420cf1dd8SToby Isaac 
15520cf1dd8SToby Isaac /*MC
156*b24fb147SBarry Smith   PETSCDUALSPACESIMPLE = "simple" - A `PetscDualSpaceType` that encapsulates a dual space of functionals provided with `PetscDualSpaceSimpleSetFunctional()`
15720cf1dd8SToby Isaac 
15820cf1dd8SToby Isaac   Level: intermediate
15920cf1dd8SToby Isaac 
160*b24fb147SBarry Smith   Developer Note:
161*b24fb147SBarry Smith   It is not clear this has a good name
162*b24fb147SBarry Smith 
163*b24fb147SBarry Smith .seealso: `PetscDualSpace`, `PetscDualSpaceSimpleSetFunctional()`, `PetscDualSpaceType`, `PetscDualSpaceCreate()`, `PetscDualSpaceSetType()`
16420cf1dd8SToby Isaac M*/
16520cf1dd8SToby Isaac 
PetscDualSpaceCreate_Simple(PetscDualSpace sp)166d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscDualSpaceCreate_Simple(PetscDualSpace sp)
167d71ae5a4SJacob Faibussowitsch {
16820cf1dd8SToby Isaac   PetscDualSpace_Simple *s;
16920cf1dd8SToby Isaac 
17020cf1dd8SToby Isaac   PetscFunctionBegin;
17120cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1724dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&s));
17320cf1dd8SToby Isaac   sp->data = s;
17420cf1dd8SToby Isaac 
17520cf1dd8SToby Isaac   s->dim    = 0;
17620cf1dd8SToby Isaac   s->numDof = NULL;
17720cf1dd8SToby Isaac 
1789566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceInitialize_Simple(sp));
1799566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSimpleSetDimension_C", PetscDualSpaceSimpleSetDimension_Simple));
1809566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSimpleSetFunctional_C", PetscDualSpaceSimpleSetFunctional_Simple));
1813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
18220cf1dd8SToby Isaac }
183