xref: /petsc/src/dm/dt/dualspace/impls/simple/dspacesimple.c (revision 20cf1dd8cd656e27510885a71ea4be028bf48813)
1*20cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
2*20cf1dd8SToby Isaac #include <petscdmplex.h>
3*20cf1dd8SToby Isaac 
4*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetUp_Simple(PetscDualSpace sp)
5*20cf1dd8SToby Isaac {
6*20cf1dd8SToby Isaac   PetscDualSpace_Simple *s  = (PetscDualSpace_Simple *) sp->data;
7*20cf1dd8SToby Isaac   DM                     dm = sp->dm;
8*20cf1dd8SToby Isaac   PetscInt               dim;
9*20cf1dd8SToby Isaac   PetscErrorCode         ierr;
10*20cf1dd8SToby Isaac 
11*20cf1dd8SToby Isaac   PetscFunctionBegin;
12*20cf1dd8SToby Isaac   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
13*20cf1dd8SToby Isaac   ierr = PetscCalloc1(dim+1, &s->numDof);CHKERRQ(ierr);
14*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
15*20cf1dd8SToby Isaac }
16*20cf1dd8SToby Isaac 
17*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceDestroy_Simple(PetscDualSpace sp)
18*20cf1dd8SToby Isaac {
19*20cf1dd8SToby Isaac   PetscDualSpace_Simple *s = (PetscDualSpace_Simple *) sp->data;
20*20cf1dd8SToby Isaac   PetscErrorCode         ierr;
21*20cf1dd8SToby Isaac 
22*20cf1dd8SToby Isaac   PetscFunctionBegin;
23*20cf1dd8SToby Isaac   ierr = PetscFree(s->numDof);CHKERRQ(ierr);
24*20cf1dd8SToby Isaac   ierr = PetscFree(s);CHKERRQ(ierr);
25*20cf1dd8SToby Isaac   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceSimpleSetDimension_C", NULL);CHKERRQ(ierr);
26*20cf1dd8SToby Isaac   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceSimpleSetFunctional_C", NULL);CHKERRQ(ierr);
27*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
28*20cf1dd8SToby Isaac }
29*20cf1dd8SToby Isaac 
30*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceDuplicate_Simple(PetscDualSpace sp, PetscDualSpace *spNew)
31*20cf1dd8SToby Isaac {
32*20cf1dd8SToby Isaac   PetscInt       dim, d, Nc;
33*20cf1dd8SToby Isaac   PetscErrorCode ierr;
34*20cf1dd8SToby Isaac 
35*20cf1dd8SToby Isaac   PetscFunctionBegin;
36*20cf1dd8SToby Isaac   ierr = PetscDualSpaceCreate(PetscObjectComm((PetscObject) sp), spNew);CHKERRQ(ierr);
37*20cf1dd8SToby Isaac   ierr = PetscDualSpaceSetType(*spNew, PETSCDUALSPACESIMPLE);CHKERRQ(ierr);
38*20cf1dd8SToby Isaac   ierr = PetscDualSpaceGetNumComponents(sp, &Nc);CHKERRQ(ierr);
39*20cf1dd8SToby Isaac   ierr = PetscDualSpaceSetNumComponents(sp, Nc);CHKERRQ(ierr);
40*20cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDimension(sp, &dim);CHKERRQ(ierr);
41*20cf1dd8SToby Isaac   ierr = PetscDualSpaceSimpleSetDimension(*spNew, dim);CHKERRQ(ierr);
42*20cf1dd8SToby Isaac   for (d = 0; d < dim; ++d) {
43*20cf1dd8SToby Isaac     PetscQuadrature q;
44*20cf1dd8SToby Isaac 
45*20cf1dd8SToby Isaac     ierr = PetscDualSpaceGetFunctional(sp, d, &q);CHKERRQ(ierr);
46*20cf1dd8SToby Isaac     ierr = PetscDualSpaceSimpleSetFunctional(*spNew, d, q);CHKERRQ(ierr);
47*20cf1dd8SToby Isaac   }
48*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
49*20cf1dd8SToby Isaac }
50*20cf1dd8SToby Isaac 
51*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetFromOptions_Simple(PetscOptionItems *PetscOptionsObject,PetscDualSpace sp)
52*20cf1dd8SToby Isaac {
53*20cf1dd8SToby Isaac   PetscFunctionBegin;
54*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
55*20cf1dd8SToby Isaac }
56*20cf1dd8SToby Isaac 
57*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetDimension_Simple(PetscDualSpace sp, PetscInt *dim)
58*20cf1dd8SToby Isaac {
59*20cf1dd8SToby Isaac   PetscDualSpace_Simple *s = (PetscDualSpace_Simple *) sp->data;
60*20cf1dd8SToby Isaac 
61*20cf1dd8SToby Isaac   PetscFunctionBegin;
62*20cf1dd8SToby Isaac   *dim = s->dim;
63*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
64*20cf1dd8SToby Isaac }
65*20cf1dd8SToby Isaac 
66*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSimpleSetDimension_Simple(PetscDualSpace sp, const PetscInt dim)
67*20cf1dd8SToby Isaac {
68*20cf1dd8SToby Isaac   PetscDualSpace_Simple *s = (PetscDualSpace_Simple *) sp->data;
69*20cf1dd8SToby Isaac   DM                     dm;
70*20cf1dd8SToby Isaac   PetscInt               spatialDim, f;
71*20cf1dd8SToby Isaac   PetscErrorCode         ierr;
72*20cf1dd8SToby Isaac 
73*20cf1dd8SToby Isaac   PetscFunctionBegin;
74*20cf1dd8SToby Isaac   for (f = 0; f < s->dim; ++f) {ierr = PetscQuadratureDestroy(&sp->functional[f]);CHKERRQ(ierr);}
75*20cf1dd8SToby Isaac   ierr = PetscFree(sp->functional);CHKERRQ(ierr);
76*20cf1dd8SToby Isaac   s->dim = dim;
77*20cf1dd8SToby Isaac   ierr = PetscCalloc1(s->dim, &sp->functional);CHKERRQ(ierr);
78*20cf1dd8SToby Isaac   ierr = PetscFree(s->numDof);CHKERRQ(ierr);
79*20cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr);
80*20cf1dd8SToby Isaac   ierr = DMGetCoordinateDim(dm, &spatialDim);CHKERRQ(ierr);
81*20cf1dd8SToby Isaac   ierr = PetscCalloc1(spatialDim+1, &s->numDof);CHKERRQ(ierr);
82*20cf1dd8SToby Isaac   s->numDof[spatialDim] = dim;
83*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
84*20cf1dd8SToby Isaac }
85*20cf1dd8SToby Isaac 
86*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetNumDof_Simple(PetscDualSpace sp, const PetscInt **numDof)
87*20cf1dd8SToby Isaac {
88*20cf1dd8SToby Isaac   PetscDualSpace_Simple *s = (PetscDualSpace_Simple *) sp->data;
89*20cf1dd8SToby Isaac 
90*20cf1dd8SToby Isaac   PetscFunctionBegin;
91*20cf1dd8SToby Isaac   *numDof = s->numDof;
92*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
93*20cf1dd8SToby Isaac }
94*20cf1dd8SToby Isaac 
95*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSimpleSetFunctional_Simple(PetscDualSpace sp, PetscInt f, PetscQuadrature q)
96*20cf1dd8SToby Isaac {
97*20cf1dd8SToby Isaac   PetscDualSpace_Simple *s = (PetscDualSpace_Simple *) sp->data;
98*20cf1dd8SToby Isaac   PetscReal             *weights;
99*20cf1dd8SToby Isaac   PetscInt               Nc, c, Nq, p;
100*20cf1dd8SToby Isaac   PetscErrorCode         ierr;
101*20cf1dd8SToby Isaac 
102*20cf1dd8SToby Isaac   PetscFunctionBegin;
103*20cf1dd8SToby 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);
104*20cf1dd8SToby Isaac   ierr = PetscQuadratureDuplicate(q, &sp->functional[f]);CHKERRQ(ierr);
105*20cf1dd8SToby Isaac   /* Reweight so that it has unit volume: Do we want to do this for Nc > 1? */
106*20cf1dd8SToby Isaac   ierr = PetscQuadratureGetData(sp->functional[f], NULL, &Nc, &Nq, NULL, (const PetscReal **) &weights);CHKERRQ(ierr);
107*20cf1dd8SToby Isaac   for (c = 0; c < Nc; ++c) {
108*20cf1dd8SToby Isaac     PetscReal vol = 0.0;
109*20cf1dd8SToby Isaac 
110*20cf1dd8SToby Isaac     for (p = 0; p < Nq; ++p) vol += weights[p*Nc+c];
111*20cf1dd8SToby Isaac     for (p = 0; p < Nq; ++p) weights[p*Nc+c] /= (vol == 0.0 ? 1.0 : vol);
112*20cf1dd8SToby Isaac   }
113*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
114*20cf1dd8SToby Isaac }
115*20cf1dd8SToby Isaac 
116*20cf1dd8SToby Isaac /*@
117*20cf1dd8SToby Isaac   PetscDualSpaceSimpleSetDimension - Set the number of functionals in the dual space basis
118*20cf1dd8SToby Isaac 
119*20cf1dd8SToby Isaac   Logically Collective on PetscDualSpace
120*20cf1dd8SToby Isaac 
121*20cf1dd8SToby Isaac   Input Parameters:
122*20cf1dd8SToby Isaac + sp  - the PetscDualSpace
123*20cf1dd8SToby Isaac - dim - the basis dimension
124*20cf1dd8SToby Isaac 
125*20cf1dd8SToby Isaac   Level: intermediate
126*20cf1dd8SToby Isaac 
127*20cf1dd8SToby Isaac .keywords: PetscDualSpace, dimension
128*20cf1dd8SToby Isaac .seealso: PetscDualSpaceSimpleSetFunctional()
129*20cf1dd8SToby Isaac @*/
130*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSimpleSetDimension(PetscDualSpace sp, PetscInt dim)
131*20cf1dd8SToby Isaac {
132*20cf1dd8SToby Isaac   PetscErrorCode ierr;
133*20cf1dd8SToby Isaac 
134*20cf1dd8SToby Isaac   PetscFunctionBegin;
135*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
136*20cf1dd8SToby Isaac   PetscValidLogicalCollectiveInt(sp, dim, 2);
137*20cf1dd8SToby Isaac   ierr = PetscTryMethod(sp, "PetscDualSpaceSimpleSetDimension_C", (PetscDualSpace,PetscInt),(sp,dim));CHKERRQ(ierr);
138*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
139*20cf1dd8SToby Isaac }
140*20cf1dd8SToby Isaac 
141*20cf1dd8SToby Isaac /*@
142*20cf1dd8SToby Isaac   PetscDualSpaceSimpleSetFunctional - Set the given basis element for this dual space
143*20cf1dd8SToby Isaac 
144*20cf1dd8SToby Isaac   Not Collective
145*20cf1dd8SToby Isaac 
146*20cf1dd8SToby Isaac   Input Parameters:
147*20cf1dd8SToby Isaac + sp  - the PetscDualSpace
148*20cf1dd8SToby Isaac . f - the basis index
149*20cf1dd8SToby Isaac - q - the basis functional
150*20cf1dd8SToby Isaac 
151*20cf1dd8SToby Isaac   Level: intermediate
152*20cf1dd8SToby Isaac 
153*20cf1dd8SToby Isaac   Note: The quadrature will be reweighted so that it has unit volume.
154*20cf1dd8SToby Isaac 
155*20cf1dd8SToby Isaac .keywords: PetscDualSpace, functional
156*20cf1dd8SToby Isaac .seealso: PetscDualSpaceSimpleSetDimension()
157*20cf1dd8SToby Isaac @*/
158*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSimpleSetFunctional(PetscDualSpace sp, PetscInt func, PetscQuadrature q)
159*20cf1dd8SToby Isaac {
160*20cf1dd8SToby Isaac   PetscErrorCode ierr;
161*20cf1dd8SToby Isaac 
162*20cf1dd8SToby Isaac   PetscFunctionBegin;
163*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
164*20cf1dd8SToby Isaac   ierr = PetscTryMethod(sp, "PetscDualSpaceSimpleSetFunctional_C", (PetscDualSpace,PetscInt,PetscQuadrature),(sp,func,q));CHKERRQ(ierr);
165*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
166*20cf1dd8SToby Isaac }
167*20cf1dd8SToby Isaac 
168*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceInitialize_Simple(PetscDualSpace sp)
169*20cf1dd8SToby Isaac {
170*20cf1dd8SToby Isaac   PetscFunctionBegin;
171*20cf1dd8SToby Isaac   sp->ops->setfromoptions    = PetscDualSpaceSetFromOptions_Simple;
172*20cf1dd8SToby Isaac   sp->ops->setup             = PetscDualSpaceSetUp_Simple;
173*20cf1dd8SToby Isaac   sp->ops->view              = NULL;
174*20cf1dd8SToby Isaac   sp->ops->destroy           = PetscDualSpaceDestroy_Simple;
175*20cf1dd8SToby Isaac   sp->ops->duplicate         = PetscDualSpaceDuplicate_Simple;
176*20cf1dd8SToby Isaac   sp->ops->getdimension      = PetscDualSpaceGetDimension_Simple;
177*20cf1dd8SToby Isaac   sp->ops->getnumdof         = PetscDualSpaceGetNumDof_Simple;
178*20cf1dd8SToby Isaac   sp->ops->getheightsubspace = NULL;
179*20cf1dd8SToby Isaac   sp->ops->getsymmetries     = NULL;
180*20cf1dd8SToby Isaac   sp->ops->apply             = PetscDualSpaceApplyDefault;
181*20cf1dd8SToby Isaac   sp->ops->applyall          = PetscDualSpaceApplyAllDefault;
182*20cf1dd8SToby Isaac   sp->ops->createallpoints   = PetscDualSpaceCreateAllPointsDefault;
183*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
184*20cf1dd8SToby Isaac }
185*20cf1dd8SToby Isaac 
186*20cf1dd8SToby Isaac /*MC
187*20cf1dd8SToby Isaac   PETSCDUALSPACESIMPLE = "simple" - A PetscDualSpace object that encapsulates a dual space of arbitrary functionals
188*20cf1dd8SToby Isaac 
189*20cf1dd8SToby Isaac   Level: intermediate
190*20cf1dd8SToby Isaac 
191*20cf1dd8SToby Isaac .seealso: PetscDualSpaceType, PetscDualSpaceCreate(), PetscDualSpaceSetType()
192*20cf1dd8SToby Isaac M*/
193*20cf1dd8SToby Isaac 
194*20cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceCreate_Simple(PetscDualSpace sp)
195*20cf1dd8SToby Isaac {
196*20cf1dd8SToby Isaac   PetscDualSpace_Simple *s;
197*20cf1dd8SToby Isaac   PetscErrorCode         ierr;
198*20cf1dd8SToby Isaac 
199*20cf1dd8SToby Isaac   PetscFunctionBegin;
200*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
201*20cf1dd8SToby Isaac   ierr     = PetscNewLog(sp,&s);CHKERRQ(ierr);
202*20cf1dd8SToby Isaac   sp->data = s;
203*20cf1dd8SToby Isaac 
204*20cf1dd8SToby Isaac   s->dim    = 0;
205*20cf1dd8SToby Isaac   s->numDof = NULL;
206*20cf1dd8SToby Isaac 
207*20cf1dd8SToby Isaac   ierr = PetscDualSpaceInitialize_Simple(sp);CHKERRQ(ierr);
208*20cf1dd8SToby Isaac   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceSimpleSetDimension_C", PetscDualSpaceSimpleSetDimension_Simple);CHKERRQ(ierr);
209*20cf1dd8SToby Isaac   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceSimpleSetFunctional_C", PetscDualSpaceSimpleSetFunctional_Simple);CHKERRQ(ierr);
210*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
211*20cf1dd8SToby Isaac }
212*20cf1dd8SToby Isaac 
213*20cf1dd8SToby Isaac 
214*20cf1dd8SToby Isaac 
215