xref: /petsc/src/dm/dt/space/impls/poly/spacepoly.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
120cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
220cf1dd8SToby Isaac 
329b5c115SMatthew G. Knepley static PetscErrorCode PetscSpaceSetFromOptions_Polynomial(PetscOptionItems *PetscOptionsObject,PetscSpace sp)
420cf1dd8SToby Isaac {
520cf1dd8SToby Isaac   PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
620cf1dd8SToby Isaac   PetscErrorCode   ierr;
720cf1dd8SToby Isaac 
820cf1dd8SToby Isaac   PetscFunctionBegin;
920cf1dd8SToby Isaac   ierr = PetscOptionsHead(PetscOptionsObject,"PetscSpace polynomial options");CHKERRQ(ierr);
1020cf1dd8SToby Isaac   ierr = PetscOptionsBool("-petscspace_poly_tensor", "Use the tensor product polynomials", "PetscSpacePolynomialSetTensor", poly->tensor, &poly->tensor, NULL);CHKERRQ(ierr);
1120cf1dd8SToby Isaac   ierr = PetscOptionsTail();CHKERRQ(ierr);
1220cf1dd8SToby Isaac   PetscFunctionReturn(0);
1320cf1dd8SToby Isaac }
1420cf1dd8SToby Isaac 
15d9bac1caSLisandro Dalcin static PetscErrorCode PetscSpacePolynomialView_Ascii(PetscSpace sp, PetscViewer v)
1620cf1dd8SToby Isaac {
1720cf1dd8SToby Isaac   PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
1820cf1dd8SToby Isaac   PetscErrorCode   ierr;
1920cf1dd8SToby Isaac 
2020cf1dd8SToby Isaac   PetscFunctionBegin;
21f1436e55SToby Isaac   ierr = PetscViewerASCIIPrintf(v, "%s space of degree %D\n", poly->tensor ? "Tensor polynomial" : "Polynomial", sp->degree);CHKERRQ(ierr);
2220cf1dd8SToby Isaac   PetscFunctionReturn(0);
2320cf1dd8SToby Isaac }
2420cf1dd8SToby Isaac 
2529b5c115SMatthew G. Knepley static PetscErrorCode PetscSpaceView_Polynomial(PetscSpace sp, PetscViewer viewer)
2620cf1dd8SToby Isaac {
2720cf1dd8SToby Isaac   PetscBool      iascii;
2820cf1dd8SToby Isaac   PetscErrorCode ierr;
2920cf1dd8SToby Isaac 
3020cf1dd8SToby Isaac   PetscFunctionBegin;
3120cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
3220cf1dd8SToby Isaac   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
3320cf1dd8SToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
3420cf1dd8SToby Isaac   if (iascii) {ierr = PetscSpacePolynomialView_Ascii(sp, viewer);CHKERRQ(ierr);}
3520cf1dd8SToby Isaac   PetscFunctionReturn(0);
3620cf1dd8SToby Isaac }
3720cf1dd8SToby Isaac 
3829b5c115SMatthew G. Knepley static PetscErrorCode PetscSpaceDestroy_Polynomial(PetscSpace sp)
3920cf1dd8SToby Isaac {
4020cf1dd8SToby Isaac   PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
4120cf1dd8SToby Isaac   PetscErrorCode   ierr;
4220cf1dd8SToby Isaac 
4320cf1dd8SToby Isaac   PetscFunctionBegin;
4420cf1dd8SToby Isaac   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePolynomialGetTensor_C", NULL);CHKERRQ(ierr);
4520cf1dd8SToby Isaac   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePolynomialSetTensor_C", NULL);CHKERRQ(ierr);
4620cf1dd8SToby Isaac   if (poly->subspaces) {
4720cf1dd8SToby Isaac     PetscInt d;
4820cf1dd8SToby Isaac 
4920cf1dd8SToby Isaac     for (d = 0; d < sp->Nv; ++d) {
5020cf1dd8SToby Isaac       ierr = PetscSpaceDestroy(&poly->subspaces[d]);CHKERRQ(ierr);
5120cf1dd8SToby Isaac     }
5220cf1dd8SToby Isaac   }
5320cf1dd8SToby Isaac   ierr = PetscFree(poly->subspaces);CHKERRQ(ierr);
5420cf1dd8SToby Isaac   ierr = PetscFree(poly);CHKERRQ(ierr);
5520cf1dd8SToby Isaac   PetscFunctionReturn(0);
5620cf1dd8SToby Isaac }
5720cf1dd8SToby Isaac 
58f1436e55SToby Isaac static PetscErrorCode PetscSpaceSetUp_Polynomial(PetscSpace sp)
59f1436e55SToby Isaac {
60f1436e55SToby Isaac   PetscSpace_Poly *poly    = (PetscSpace_Poly *) sp->data;
61f1436e55SToby Isaac   PetscErrorCode   ierr;
62f1436e55SToby Isaac 
63f1436e55SToby Isaac   PetscFunctionBegin;
64f1436e55SToby Isaac   if (poly->setupCalled) PetscFunctionReturn(0);
65f1436e55SToby Isaac   if (sp->Nv <= 1) {
66f1436e55SToby Isaac     poly->tensor = PETSC_FALSE;
67f1436e55SToby Isaac   }
68f1436e55SToby Isaac   if (sp->Nc != 1) {
69f1436e55SToby Isaac     PetscInt    Nc = sp->Nc;
70f1436e55SToby Isaac     PetscBool   tensor = poly->tensor;
71f1436e55SToby Isaac     PetscInt    Nv = sp->Nv;
72f1436e55SToby Isaac     PetscInt    degree = sp->degree;
73417c287bSToby Isaac     const char *prefix;
74417c287bSToby Isaac     const char *name;
75417c287bSToby Isaac     char        subname[PETSC_MAX_PATH_LEN];
76f1436e55SToby Isaac     PetscSpace  subsp;
77f1436e55SToby Isaac 
78f1436e55SToby Isaac     ierr = PetscSpaceSetType(sp, PETSCSPACESUM);CHKERRQ(ierr);
79f1436e55SToby Isaac     ierr = PetscSpaceSumSetNumSubspaces(sp, Nc);CHKERRQ(ierr);
80f1436e55SToby Isaac     ierr = PetscSpaceCreate(PetscObjectComm((PetscObject)sp), &subsp);CHKERRQ(ierr);
81417c287bSToby Isaac     ierr = PetscObjectGetOptionsPrefix((PetscObject)sp, &prefix);CHKERRQ(ierr);
82417c287bSToby Isaac     ierr = PetscObjectSetOptionsPrefix((PetscObject)subsp, prefix);CHKERRQ(ierr);
83417c287bSToby Isaac     ierr = PetscObjectAppendOptionsPrefix((PetscObject)subsp, "sumcomp_");CHKERRQ(ierr);
84417c287bSToby Isaac     if (((PetscObject)sp)->name) {
85417c287bSToby Isaac       ierr = PetscObjectGetName((PetscObject)sp, &name);CHKERRQ(ierr);
86417c287bSToby Isaac       ierr = PetscSNPrintf(subname, PETSC_MAX_PATH_LEN-1, "%s sum component", name);CHKERRQ(ierr);
87417c287bSToby Isaac       ierr = PetscObjectSetName((PetscObject)subsp, subname);CHKERRQ(ierr);
88417c287bSToby Isaac     } else {
89417c287bSToby Isaac       ierr = PetscObjectSetName((PetscObject)subsp, "sum component");CHKERRQ(ierr);
90417c287bSToby Isaac     }
91f1436e55SToby Isaac     ierr = PetscSpaceSetType(subsp, PETSCSPACEPOLYNOMIAL);CHKERRQ(ierr);
92f1436e55SToby Isaac     ierr = PetscSpaceSetDegree(subsp, degree, PETSC_DETERMINE);CHKERRQ(ierr);
93f1436e55SToby Isaac     ierr = PetscSpaceSetNumComponents(subsp, 1);CHKERRQ(ierr);
94f1436e55SToby Isaac     ierr = PetscSpaceSetNumVariables(subsp, Nv);CHKERRQ(ierr);
95f1436e55SToby Isaac     ierr = PetscSpacePolynomialSetTensor(subsp, tensor);CHKERRQ(ierr);
96f1436e55SToby Isaac     ierr = PetscSpaceSetUp(subsp);CHKERRQ(ierr);
97f1436e55SToby Isaac     for (PetscInt i = 0; i < Nc; i++) {
98f1436e55SToby Isaac       ierr = PetscSpaceSumSetSubspace(sp, i, subsp);CHKERRQ(ierr);
99f1436e55SToby Isaac     }
100f1436e55SToby Isaac     ierr = PetscSpaceDestroy(&subsp);CHKERRQ(ierr);
101f1436e55SToby Isaac     ierr = PetscSpaceSetUp(sp);CHKERRQ(ierr);
102f1436e55SToby Isaac     PetscFunctionReturn(0);
103f1436e55SToby Isaac   }
104f1436e55SToby Isaac   if (poly->tensor) {
105f1436e55SToby Isaac     sp->maxDegree = PETSC_DETERMINE;
106f1436e55SToby Isaac     ierr = PetscSpaceSetType(sp, PETSCSPACETENSOR);CHKERRQ(ierr);
107f1436e55SToby Isaac     ierr = PetscSpaceSetUp(sp);CHKERRQ(ierr);
108f1436e55SToby Isaac     PetscFunctionReturn(0);
109f1436e55SToby Isaac   }
110*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(sp->degree < 0,PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_OUTOFRANGE, "Negative degree %D invalid", sp->degree);
111f1436e55SToby Isaac   sp->maxDegree = sp->degree;
112f1436e55SToby Isaac   poly->setupCalled = PETSC_TRUE;
113f1436e55SToby Isaac   PetscFunctionReturn(0);
114f1436e55SToby Isaac }
115f1436e55SToby Isaac 
11629b5c115SMatthew G. Knepley static PetscErrorCode PetscSpaceGetDimension_Polynomial(PetscSpace sp, PetscInt *dim)
11720cf1dd8SToby Isaac {
11820cf1dd8SToby Isaac   PetscInt         deg  = sp->degree;
119f1436e55SToby Isaac   PetscInt         n    = sp->Nv;
12020cf1dd8SToby Isaac   PetscErrorCode   ierr;
12120cf1dd8SToby Isaac 
12220cf1dd8SToby Isaac   PetscFunctionBegin;
123f1436e55SToby Isaac   ierr = PetscDTBinomialInt(n + deg, n, dim);CHKERRQ(ierr);
124f1436e55SToby Isaac   *dim *= sp->Nc;
12520cf1dd8SToby Isaac   PetscFunctionReturn(0);
12620cf1dd8SToby Isaac }
12720cf1dd8SToby Isaac 
128f1436e55SToby Isaac static PetscErrorCode CoordinateBasis(PetscInt dim, PetscInt npoints, const PetscReal points[], PetscInt jet, PetscInt Njet, PetscReal pScalar[])
12920cf1dd8SToby Isaac {
13020cf1dd8SToby Isaac   PetscErrorCode ierr;
13120cf1dd8SToby Isaac 
13220cf1dd8SToby Isaac   PetscFunctionBegin;
133f1436e55SToby Isaac   ierr = PetscArrayzero(pScalar, (1 + dim) * Njet * npoints);CHKERRQ(ierr);
134f1436e55SToby Isaac   for (PetscInt b = 0; b < 1 + dim; b++) {
135f1436e55SToby Isaac     for (PetscInt j = 0; j < PetscMin(1 + dim, Njet); j++) {
136f1436e55SToby Isaac       if (j == 0) {
137f1436e55SToby Isaac         if (b == 0) {
138f1436e55SToby Isaac           for (PetscInt pt = 0; pt < npoints; pt++) {
139f1436e55SToby Isaac             pScalar[b * Njet * npoints + j * npoints + pt] = 1.;
140f1436e55SToby Isaac           }
14120cf1dd8SToby Isaac         } else {
142f1436e55SToby Isaac           for (PetscInt pt = 0; pt < npoints; pt++) {
143f1436e55SToby Isaac             pScalar[b * Njet * npoints + j * npoints + pt] = points[pt * dim + (b-1)];
144f1436e55SToby Isaac           }
145f1436e55SToby Isaac         }
146f1436e55SToby Isaac       } else if (j == b) {
147f1436e55SToby Isaac         for (PetscInt pt = 0; pt < npoints; pt++) {
148f1436e55SToby Isaac           pScalar[b * Njet * npoints + j * npoints + pt] = 1.;
149f1436e55SToby Isaac         }
150f1436e55SToby Isaac       }
15120cf1dd8SToby Isaac     }
15220cf1dd8SToby Isaac   }
15320cf1dd8SToby Isaac   PetscFunctionReturn(0);
15420cf1dd8SToby Isaac }
15520cf1dd8SToby Isaac 
15629b5c115SMatthew G. Knepley static PetscErrorCode PetscSpaceEvaluate_Polynomial(PetscSpace sp, PetscInt npoints, const PetscReal points[], PetscReal B[], PetscReal D[], PetscReal H[])
15720cf1dd8SToby Isaac {
15820cf1dd8SToby Isaac   PetscSpace_Poly *poly    = (PetscSpace_Poly *) sp->data;
15920cf1dd8SToby Isaac   DM               dm      = sp->dm;
16020cf1dd8SToby Isaac   PetscInt         dim     = sp->Nv;
161f1436e55SToby Isaac   PetscInt         Nb, jet, Njet;
162f1436e55SToby Isaac   PetscReal       *pScalar;
16320cf1dd8SToby Isaac   PetscErrorCode   ierr;
16420cf1dd8SToby Isaac 
16520cf1dd8SToby Isaac   PetscFunctionBegin;
166f1436e55SToby Isaac   if (!poly->setupCalled) {
167f1436e55SToby Isaac     ierr = PetscSpaceSetUp(sp);CHKERRQ(ierr);
168f1436e55SToby Isaac     ierr = PetscSpaceEvaluate(sp, npoints, points, B, D, H);CHKERRQ(ierr);
169f1436e55SToby Isaac     PetscFunctionReturn(0);
17020cf1dd8SToby Isaac   }
171*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(poly->tensor || sp->Nc != 1,PETSC_COMM_SELF, PETSC_ERR_PLIB, "tensor and multicomponent spaces should have been converted");
172f1436e55SToby Isaac   ierr = PetscDTBinomialInt(dim + sp->degree, dim, &Nb);CHKERRQ(ierr);
173f1436e55SToby Isaac   if (H) {
174f1436e55SToby Isaac     jet = 2;
175f1436e55SToby Isaac   } else if (D) {
176f1436e55SToby Isaac     jet = 1;
177f1436e55SToby Isaac   } else {
178f1436e55SToby Isaac     jet = 0;
17920cf1dd8SToby Isaac   }
180f1436e55SToby Isaac   ierr = PetscDTBinomialInt(dim + jet, dim, &Njet);CHKERRQ(ierr);
181f1436e55SToby Isaac   ierr = DMGetWorkArray(dm, Nb * Njet * npoints, MPIU_REAL, &pScalar);CHKERRQ(ierr);
182f1436e55SToby Isaac   // Why are we handling the case degree == 1 specially?  Because we don't want numerical noise when we evaluate hat
183f1436e55SToby Isaac   // functions at the vertices of a simplex, which happens when we invert the Vandermonde matrix of the PKD basis.
184f1436e55SToby Isaac   // We don't make any promise about which basis is used.
185f1436e55SToby Isaac   if (sp->degree == 1) {
186f1436e55SToby Isaac     ierr = CoordinateBasis(dim, npoints, points, jet, Njet, pScalar);CHKERRQ(ierr);
187f1436e55SToby Isaac   } else {
188f1436e55SToby Isaac     ierr = PetscDTPKDEvalJet(dim, npoints, points, sp->degree, jet, pScalar);CHKERRQ(ierr);
18920cf1dd8SToby Isaac   }
19020cf1dd8SToby Isaac   if (B) {
191f1436e55SToby Isaac     PetscInt p_strl = Nb;
192f1436e55SToby Isaac     PetscInt b_strl = 1;
1933596293dSMatthew G. Knepley 
194f1436e55SToby Isaac     PetscInt b_strr = Njet*npoints;
195f1436e55SToby Isaac     PetscInt p_strr = 1;
196f1436e55SToby Isaac 
197f1436e55SToby Isaac     ierr = PetscArrayzero(B, npoints * Nb);CHKERRQ(ierr);
198f1436e55SToby Isaac     for (PetscInt b = 0; b < Nb; b++) {
199f1436e55SToby Isaac       for (PetscInt p = 0; p < npoints; p++) {
200f1436e55SToby Isaac         B[p*p_strl + b*b_strl] = pScalar[b*b_strr + p*p_strr];
2013596293dSMatthew G. Knepley       }
2023596293dSMatthew G. Knepley     }
20320cf1dd8SToby Isaac   }
20420cf1dd8SToby Isaac   if (D) {
205f1436e55SToby Isaac     PetscInt p_strl = dim*Nb;
206f1436e55SToby Isaac     PetscInt b_strl = dim;
207f1436e55SToby Isaac     PetscInt d_strl = 1;
208f1436e55SToby Isaac 
209f1436e55SToby Isaac     PetscInt b_strr = Njet*npoints;
210f1436e55SToby Isaac     PetscInt d_strr = npoints;
211f1436e55SToby Isaac     PetscInt p_strr = 1;
212f1436e55SToby Isaac 
213f1436e55SToby Isaac     ierr = PetscArrayzero(D, npoints * Nb * dim);CHKERRQ(ierr);
214f1436e55SToby Isaac     for (PetscInt d = 0; d < dim; d++) {
215f1436e55SToby Isaac       for (PetscInt b = 0; b < Nb; b++) {
216f1436e55SToby Isaac         for (PetscInt p = 0; p < npoints; p++) {
217f1436e55SToby Isaac           D[p*p_strl + b*b_strl + d*d_strl] = pScalar[b*b_strr + (1+d)*d_strr + p*p_strr];
21820cf1dd8SToby Isaac         }
21920cf1dd8SToby Isaac       }
22020cf1dd8SToby Isaac     }
22120cf1dd8SToby Isaac   }
22220cf1dd8SToby Isaac   if (H) {
223f1436e55SToby Isaac     PetscInt p_strl  = dim*dim*Nb;
224f1436e55SToby Isaac     PetscInt b_strl  = dim*dim;
225f1436e55SToby Isaac     PetscInt d1_strl = dim;
226f1436e55SToby Isaac     PetscInt d2_strl = 1;
227f1436e55SToby Isaac 
228f1436e55SToby Isaac     PetscInt b_strr = Njet*npoints;
229f1436e55SToby Isaac     PetscInt j_strr = npoints;
230f1436e55SToby Isaac     PetscInt p_strr = 1;
231f1436e55SToby Isaac 
232f1436e55SToby Isaac     PetscInt *derivs;
233f1436e55SToby Isaac     ierr = PetscCalloc1(dim, &derivs);CHKERRQ(ierr);
234f1436e55SToby Isaac     ierr = PetscArrayzero(H, npoints * Nb * dim * dim);CHKERRQ(ierr);
235f1436e55SToby Isaac     for (PetscInt d1 = 0; d1 < dim; d1++) {
236f1436e55SToby Isaac       for (PetscInt d2 = 0; d2 < dim; d2++) {
237f1436e55SToby Isaac         PetscInt j;
238f1436e55SToby Isaac         derivs[d1]++;
239f1436e55SToby Isaac         derivs[d2]++;
240f1436e55SToby Isaac         ierr = PetscDTGradedOrderToIndex(dim, derivs, &j);CHKERRQ(ierr);
241f1436e55SToby Isaac         derivs[d1]--;
242f1436e55SToby Isaac         derivs[d2]--;
243f1436e55SToby Isaac         for (PetscInt b = 0; b < Nb; b++) {
244f1436e55SToby Isaac           for (PetscInt p = 0; p < npoints; p++) {
245f1436e55SToby Isaac             H[p*p_strl + b*b_strl + d1*d1_strl + d2*d2_strl] = pScalar[b*b_strr + j*j_strr + p*p_strr];
24620cf1dd8SToby Isaac           }
24720cf1dd8SToby Isaac         }
24820cf1dd8SToby Isaac       }
24920cf1dd8SToby Isaac     }
250f1436e55SToby Isaac     ierr = PetscFree(derivs);CHKERRQ(ierr);
25120cf1dd8SToby Isaac   }
252f1436e55SToby Isaac   ierr = DMRestoreWorkArray(dm, Nb * Njet * npoints, MPIU_REAL, &pScalar);CHKERRQ(ierr);
25320cf1dd8SToby Isaac   PetscFunctionReturn(0);
25420cf1dd8SToby Isaac }
25520cf1dd8SToby Isaac 
25620cf1dd8SToby Isaac /*@
25720cf1dd8SToby Isaac   PetscSpacePolynomialSetTensor - Set whether a function space is a space of tensor polynomials (the space is spanned
258f1436e55SToby Isaac   by polynomials whose degree in each variable is bounded by the given order), as opposed to polynomials (the space is
25920cf1dd8SToby Isaac   spanned by polynomials whose total degree---summing over all variables---is bounded by the given order).
26020cf1dd8SToby Isaac 
26120cf1dd8SToby Isaac   Input Parameters:
26220cf1dd8SToby Isaac + sp     - the function space object
26320cf1dd8SToby Isaac - tensor - PETSC_TRUE for a tensor polynomial space, PETSC_FALSE for a polynomial space
26420cf1dd8SToby Isaac 
2654ab77754SMatthew G. Knepley   Options Database:
2664ab77754SMatthew G. Knepley . -petscspace_poly_tensor <bool> - Whether to use tensor product polynomials in higher dimension
2674ab77754SMatthew G. Knepley 
26829b5c115SMatthew G. Knepley   Level: intermediate
26920cf1dd8SToby Isaac 
27020cf1dd8SToby Isaac .seealso: PetscSpacePolynomialGetTensor(), PetscSpaceSetDegree(), PetscSpaceSetNumVariables()
27120cf1dd8SToby Isaac @*/
27220cf1dd8SToby Isaac PetscErrorCode PetscSpacePolynomialSetTensor(PetscSpace sp, PetscBool tensor)
27320cf1dd8SToby Isaac {
27420cf1dd8SToby Isaac   PetscErrorCode ierr;
27520cf1dd8SToby Isaac 
27620cf1dd8SToby Isaac   PetscFunctionBegin;
27720cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
27820cf1dd8SToby Isaac   ierr = PetscTryMethod(sp,"PetscSpacePolynomialSetTensor_C",(PetscSpace,PetscBool),(sp,tensor));CHKERRQ(ierr);
27920cf1dd8SToby Isaac   PetscFunctionReturn(0);
28020cf1dd8SToby Isaac }
28120cf1dd8SToby Isaac 
28220cf1dd8SToby Isaac /*@
28320cf1dd8SToby Isaac   PetscSpacePolynomialGetTensor - Get whether a function space is a space of tensor polynomials (the space is spanned
28420cf1dd8SToby Isaac   by polynomials whose degree in each variabl is bounded by the given order), as opposed to polynomials (the space is
28520cf1dd8SToby Isaac   spanned by polynomials whose total degree---summing over all variables---is bounded by the given order).
28620cf1dd8SToby Isaac 
28720cf1dd8SToby Isaac   Input Parameters:
28820cf1dd8SToby Isaac . sp     - the function space object
28920cf1dd8SToby Isaac 
29020cf1dd8SToby Isaac   Output Parameters:
29120cf1dd8SToby Isaac . tensor - PETSC_TRUE for a tensor polynomial space, PETSC_FALSE for a polynomial space
29220cf1dd8SToby Isaac 
29329b5c115SMatthew G. Knepley   Level: intermediate
29420cf1dd8SToby Isaac 
29520cf1dd8SToby Isaac .seealso: PetscSpacePolynomialSetTensor(), PetscSpaceSetDegree(), PetscSpaceSetNumVariables()
29620cf1dd8SToby Isaac @*/
29720cf1dd8SToby Isaac PetscErrorCode PetscSpacePolynomialGetTensor(PetscSpace sp, PetscBool *tensor)
29820cf1dd8SToby Isaac {
29920cf1dd8SToby Isaac   PetscErrorCode ierr;
30020cf1dd8SToby Isaac 
30120cf1dd8SToby Isaac   PetscFunctionBegin;
30220cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
30320cf1dd8SToby Isaac   PetscValidPointer(tensor, 2);
30420cf1dd8SToby Isaac   ierr = PetscTryMethod(sp,"PetscSpacePolynomialGetTensor_C",(PetscSpace,PetscBool*),(sp,tensor));CHKERRQ(ierr);
30520cf1dd8SToby Isaac   PetscFunctionReturn(0);
30620cf1dd8SToby Isaac }
30720cf1dd8SToby Isaac 
30820cf1dd8SToby Isaac static PetscErrorCode PetscSpacePolynomialSetTensor_Polynomial(PetscSpace sp, PetscBool tensor)
30920cf1dd8SToby Isaac {
31020cf1dd8SToby Isaac   PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
31120cf1dd8SToby Isaac 
31220cf1dd8SToby Isaac   PetscFunctionBegin;
31320cf1dd8SToby Isaac   poly->tensor = tensor;
31420cf1dd8SToby Isaac   PetscFunctionReturn(0);
31520cf1dd8SToby Isaac }
31620cf1dd8SToby Isaac 
31720cf1dd8SToby Isaac static PetscErrorCode PetscSpacePolynomialGetTensor_Polynomial(PetscSpace sp, PetscBool *tensor)
31820cf1dd8SToby Isaac {
31920cf1dd8SToby Isaac   PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
32020cf1dd8SToby Isaac 
32120cf1dd8SToby Isaac   PetscFunctionBegin;
32220cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
32320cf1dd8SToby Isaac   PetscValidPointer(tensor, 2);
32420cf1dd8SToby Isaac   *tensor = poly->tensor;
32520cf1dd8SToby Isaac   PetscFunctionReturn(0);
32620cf1dd8SToby Isaac }
32720cf1dd8SToby Isaac 
32820cf1dd8SToby Isaac static PetscErrorCode PetscSpaceGetHeightSubspace_Polynomial(PetscSpace sp, PetscInt height, PetscSpace *subsp)
32920cf1dd8SToby Isaac {
33020cf1dd8SToby Isaac   PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
33120cf1dd8SToby Isaac   PetscInt         Nc, dim, order;
33220cf1dd8SToby Isaac   PetscBool        tensor;
33320cf1dd8SToby Isaac   PetscErrorCode   ierr;
33420cf1dd8SToby Isaac 
33520cf1dd8SToby Isaac   PetscFunctionBegin;
33620cf1dd8SToby Isaac   ierr = PetscSpaceGetNumComponents(sp, &Nc);CHKERRQ(ierr);
33720cf1dd8SToby Isaac   ierr = PetscSpaceGetNumVariables(sp, &dim);CHKERRQ(ierr);
33820cf1dd8SToby Isaac   ierr = PetscSpaceGetDegree(sp, &order, NULL);CHKERRQ(ierr);
33920cf1dd8SToby Isaac   ierr = PetscSpacePolynomialGetTensor(sp, &tensor);CHKERRQ(ierr);
340*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(height > dim || height < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Asked for space at height %D for dimension %D space", height, dim);
34120cf1dd8SToby Isaac   if (!poly->subspaces) {ierr = PetscCalloc1(dim, &poly->subspaces);CHKERRQ(ierr);}
34220cf1dd8SToby Isaac   if (height <= dim) {
34320cf1dd8SToby Isaac     if (!poly->subspaces[height-1]) {
34420cf1dd8SToby Isaac       PetscSpace  sub;
3453f6b16c7SMatthew G. Knepley       const char *name;
34620cf1dd8SToby Isaac 
34720cf1dd8SToby Isaac       ierr = PetscSpaceCreate(PetscObjectComm((PetscObject) sp), &sub);CHKERRQ(ierr);
3483f6b16c7SMatthew G. Knepley       ierr = PetscObjectGetName((PetscObject) sp,  &name);CHKERRQ(ierr);
3493f6b16c7SMatthew G. Knepley       ierr = PetscObjectSetName((PetscObject) sub,  name);CHKERRQ(ierr);
35020cf1dd8SToby Isaac       ierr = PetscSpaceSetType(sub, PETSCSPACEPOLYNOMIAL);CHKERRQ(ierr);
351d39dd5f5SToby Isaac       ierr = PetscSpaceSetNumComponents(sub, Nc);CHKERRQ(ierr);
352d39dd5f5SToby Isaac       ierr = PetscSpaceSetDegree(sub, order, PETSC_DETERMINE);CHKERRQ(ierr);
35320cf1dd8SToby Isaac       ierr = PetscSpaceSetNumVariables(sub, dim-height);CHKERRQ(ierr);
35420cf1dd8SToby Isaac       ierr = PetscSpacePolynomialSetTensor(sub, tensor);CHKERRQ(ierr);
35520cf1dd8SToby Isaac       ierr = PetscSpaceSetUp(sub);CHKERRQ(ierr);
35620cf1dd8SToby Isaac       poly->subspaces[height-1] = sub;
35720cf1dd8SToby Isaac     }
35820cf1dd8SToby Isaac     *subsp = poly->subspaces[height-1];
35920cf1dd8SToby Isaac   } else {
36020cf1dd8SToby Isaac     *subsp = NULL;
36120cf1dd8SToby Isaac   }
36220cf1dd8SToby Isaac   PetscFunctionReturn(0);
36320cf1dd8SToby Isaac }
36420cf1dd8SToby Isaac 
36529b5c115SMatthew G. Knepley static PetscErrorCode PetscSpaceInitialize_Polynomial(PetscSpace sp)
36620cf1dd8SToby Isaac {
36720cf1dd8SToby Isaac   PetscErrorCode ierr;
36820cf1dd8SToby Isaac 
36920cf1dd8SToby Isaac   PetscFunctionBegin;
37020cf1dd8SToby Isaac   sp->ops->setfromoptions    = PetscSpaceSetFromOptions_Polynomial;
37120cf1dd8SToby Isaac   sp->ops->setup             = PetscSpaceSetUp_Polynomial;
37220cf1dd8SToby Isaac   sp->ops->view              = PetscSpaceView_Polynomial;
37320cf1dd8SToby Isaac   sp->ops->destroy           = PetscSpaceDestroy_Polynomial;
37420cf1dd8SToby Isaac   sp->ops->getdimension      = PetscSpaceGetDimension_Polynomial;
37520cf1dd8SToby Isaac   sp->ops->evaluate          = PetscSpaceEvaluate_Polynomial;
37620cf1dd8SToby Isaac   sp->ops->getheightsubspace = PetscSpaceGetHeightSubspace_Polynomial;
37720cf1dd8SToby Isaac   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePolynomialGetTensor_C", PetscSpacePolynomialGetTensor_Polynomial);CHKERRQ(ierr);
37820cf1dd8SToby Isaac   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePolynomialSetTensor_C", PetscSpacePolynomialSetTensor_Polynomial);CHKERRQ(ierr);
37920cf1dd8SToby Isaac   PetscFunctionReturn(0);
38020cf1dd8SToby Isaac }
38120cf1dd8SToby Isaac 
38220cf1dd8SToby Isaac /*MC
38320cf1dd8SToby Isaac   PETSCSPACEPOLYNOMIAL = "poly" - A PetscSpace object that encapsulates a polynomial space, e.g. P1 is the space of
38420cf1dd8SToby Isaac   linear polynomials. The space is replicated for each component.
38520cf1dd8SToby Isaac 
38620cf1dd8SToby Isaac   Level: intermediate
38720cf1dd8SToby Isaac 
38820cf1dd8SToby Isaac .seealso: PetscSpaceType, PetscSpaceCreate(), PetscSpaceSetType()
38920cf1dd8SToby Isaac M*/
39020cf1dd8SToby Isaac 
39120cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode PetscSpaceCreate_Polynomial(PetscSpace sp)
39220cf1dd8SToby Isaac {
39320cf1dd8SToby Isaac   PetscSpace_Poly *poly;
39420cf1dd8SToby Isaac   PetscErrorCode   ierr;
39520cf1dd8SToby Isaac 
39620cf1dd8SToby Isaac   PetscFunctionBegin;
39720cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
39820cf1dd8SToby Isaac   ierr     = PetscNewLog(sp,&poly);CHKERRQ(ierr);
39920cf1dd8SToby Isaac   sp->data = poly;
40020cf1dd8SToby Isaac 
40120cf1dd8SToby Isaac   poly->tensor    = PETSC_FALSE;
40220cf1dd8SToby Isaac   poly->subspaces = NULL;
40320cf1dd8SToby Isaac 
40420cf1dd8SToby Isaac   ierr = PetscSpaceInitialize_Polynomial(sp);CHKERRQ(ierr);
40520cf1dd8SToby Isaac   PetscFunctionReturn(0);
40620cf1dd8SToby Isaac }
40720cf1dd8SToby Isaac 
408