xref: /petsc/src/dm/dt/space/impls/poly/spacepoly.c (revision f1436e559c5ffc6e8099008682a09bacce1cc973)
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;
21*f1436e55SToby 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 
58*f1436e55SToby Isaac static PetscErrorCode PetscSpaceSetUp_Polynomial(PetscSpace sp)
59*f1436e55SToby Isaac {
60*f1436e55SToby Isaac   PetscSpace_Poly *poly    = (PetscSpace_Poly *) sp->data;
61*f1436e55SToby Isaac   PetscErrorCode   ierr;
62*f1436e55SToby Isaac 
63*f1436e55SToby Isaac   PetscFunctionBegin;
64*f1436e55SToby Isaac   if (poly->setupCalled) PetscFunctionReturn(0);
65*f1436e55SToby Isaac   if (sp->Nv <= 1) {
66*f1436e55SToby Isaac     poly->tensor = PETSC_FALSE;
67*f1436e55SToby Isaac   }
68*f1436e55SToby Isaac   if (sp->Nc != 1) {
69*f1436e55SToby Isaac     PetscInt  Nc = sp->Nc;
70*f1436e55SToby Isaac     PetscBool tensor = poly->tensor;
71*f1436e55SToby Isaac     PetscInt  Nv = sp->Nv;
72*f1436e55SToby Isaac     PetscInt  degree = sp->degree;
73*f1436e55SToby Isaac     PetscSpace subsp;
74*f1436e55SToby Isaac 
75*f1436e55SToby Isaac     ierr = PetscSpaceSetType(sp, PETSCSPACESUM);CHKERRQ(ierr);
76*f1436e55SToby Isaac     ierr = PetscSpaceSumSetNumSubspaces(sp, Nc);CHKERRQ(ierr);
77*f1436e55SToby Isaac     ierr = PetscSpaceCreate(PetscObjectComm((PetscObject)sp), &subsp);CHKERRQ(ierr);
78*f1436e55SToby Isaac     ierr = PetscSpaceSetType(subsp, PETSCSPACEPOLYNOMIAL);CHKERRQ(ierr);
79*f1436e55SToby Isaac     ierr = PetscSpaceSetDegree(subsp, degree, PETSC_DETERMINE);CHKERRQ(ierr);
80*f1436e55SToby Isaac     ierr = PetscSpaceSetNumComponents(subsp, 1);CHKERRQ(ierr);
81*f1436e55SToby Isaac     ierr = PetscSpaceSetNumVariables(subsp, Nv);CHKERRQ(ierr);
82*f1436e55SToby Isaac     ierr = PetscSpacePolynomialSetTensor(subsp, tensor);CHKERRQ(ierr);
83*f1436e55SToby Isaac     ierr = PetscSpaceSetUp(subsp);CHKERRQ(ierr);
84*f1436e55SToby Isaac     for (PetscInt i = 0; i < Nc; i++) {
85*f1436e55SToby Isaac       ierr = PetscSpaceSumSetSubspace(sp, i, subsp);CHKERRQ(ierr);
86*f1436e55SToby Isaac     }
87*f1436e55SToby Isaac     ierr = PetscSpaceDestroy(&subsp);CHKERRQ(ierr);
88*f1436e55SToby Isaac     ierr = PetscSpaceSetUp(sp);CHKERRQ(ierr);
89*f1436e55SToby Isaac     PetscFunctionReturn(0);
90*f1436e55SToby Isaac   }
91*f1436e55SToby Isaac   if (poly->tensor) {
92*f1436e55SToby Isaac     sp->maxDegree = PETSC_DETERMINE;
93*f1436e55SToby Isaac     ierr = PetscSpaceSetType(sp, PETSCSPACETENSOR);CHKERRQ(ierr);
94*f1436e55SToby Isaac     ierr = PetscSpaceSetUp(sp);CHKERRQ(ierr);
95*f1436e55SToby Isaac     PetscFunctionReturn(0);
96*f1436e55SToby Isaac   }
97*f1436e55SToby Isaac   if (sp->degree < 0) SETERRQ1(PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_OUTOFRANGE, "Negative degree %D invalid\n", sp->degree);
98*f1436e55SToby Isaac   sp->maxDegree = sp->degree;
99*f1436e55SToby Isaac   poly->setupCalled = PETSC_TRUE;
100*f1436e55SToby Isaac   PetscFunctionReturn(0);
101*f1436e55SToby Isaac }
102*f1436e55SToby Isaac 
10329b5c115SMatthew G. Knepley static PetscErrorCode PetscSpaceGetDimension_Polynomial(PetscSpace sp, PetscInt *dim)
10420cf1dd8SToby Isaac {
10520cf1dd8SToby Isaac   PetscInt         deg  = sp->degree;
106*f1436e55SToby Isaac   PetscInt         n    = sp->Nv;
10720cf1dd8SToby Isaac   PetscErrorCode   ierr;
10820cf1dd8SToby Isaac 
10920cf1dd8SToby Isaac   PetscFunctionBegin;
110*f1436e55SToby Isaac   ierr = PetscDTBinomialInt(n + deg, n, dim);CHKERRQ(ierr);
111*f1436e55SToby Isaac   *dim *= sp->Nc;
11220cf1dd8SToby Isaac   PetscFunctionReturn(0);
11320cf1dd8SToby Isaac }
11420cf1dd8SToby Isaac 
115*f1436e55SToby Isaac static PetscErrorCode CoordinateBasis(PetscInt dim, PetscInt npoints, const PetscReal points[], PetscInt jet, PetscInt Njet, PetscReal pScalar[])
11620cf1dd8SToby Isaac {
11720cf1dd8SToby Isaac   PetscErrorCode ierr;
11820cf1dd8SToby Isaac 
11920cf1dd8SToby Isaac   PetscFunctionBegin;
120*f1436e55SToby Isaac   ierr = PetscArrayzero(pScalar, (1 + dim) * Njet * npoints);CHKERRQ(ierr);
121*f1436e55SToby Isaac   for (PetscInt b = 0; b < 1 + dim; b++) {
122*f1436e55SToby Isaac     for (PetscInt j = 0; j < PetscMin(1 + dim, Njet); j++) {
123*f1436e55SToby Isaac       if (j == 0) {
124*f1436e55SToby Isaac         if (b == 0) {
125*f1436e55SToby Isaac           for (PetscInt pt = 0; pt < npoints; pt++) {
126*f1436e55SToby Isaac             pScalar[b * Njet * npoints + j * npoints + pt] = 1.;
127*f1436e55SToby Isaac           }
12820cf1dd8SToby Isaac         } else {
129*f1436e55SToby Isaac           for (PetscInt pt = 0; pt < npoints; pt++) {
130*f1436e55SToby Isaac             pScalar[b * Njet * npoints + j * npoints + pt] = points[pt * dim + (b-1)];
131*f1436e55SToby Isaac           }
132*f1436e55SToby Isaac         }
133*f1436e55SToby Isaac       } else if (j == b) {
134*f1436e55SToby Isaac         for (PetscInt pt = 0; pt < npoints; pt++) {
135*f1436e55SToby Isaac           pScalar[b * Njet * npoints + j * npoints + pt] = 1.;
136*f1436e55SToby Isaac         }
137*f1436e55SToby Isaac       }
13820cf1dd8SToby Isaac     }
13920cf1dd8SToby Isaac   }
14020cf1dd8SToby Isaac   PetscFunctionReturn(0);
14120cf1dd8SToby Isaac }
14220cf1dd8SToby Isaac 
14329b5c115SMatthew G. Knepley static PetscErrorCode PetscSpaceEvaluate_Polynomial(PetscSpace sp, PetscInt npoints, const PetscReal points[], PetscReal B[], PetscReal D[], PetscReal H[])
14420cf1dd8SToby Isaac {
14520cf1dd8SToby Isaac   PetscSpace_Poly *poly    = (PetscSpace_Poly *) sp->data;
14620cf1dd8SToby Isaac   DM               dm      = sp->dm;
14720cf1dd8SToby Isaac   PetscInt         dim     = sp->Nv;
148*f1436e55SToby Isaac   PetscInt         Nb, jet, Njet;
149*f1436e55SToby Isaac   PetscReal       *pScalar;
15020cf1dd8SToby Isaac   PetscErrorCode   ierr;
15120cf1dd8SToby Isaac 
15220cf1dd8SToby Isaac   PetscFunctionBegin;
153*f1436e55SToby Isaac   if (!poly->setupCalled) {
154*f1436e55SToby Isaac     ierr = PetscSpaceSetUp(sp);CHKERRQ(ierr);
155*f1436e55SToby Isaac     ierr = PetscSpaceEvaluate(sp, npoints, points, B, D, H);CHKERRQ(ierr);
156*f1436e55SToby Isaac     PetscFunctionReturn(0);
15720cf1dd8SToby Isaac   }
158*f1436e55SToby Isaac   if (poly->tensor || sp->Nc != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "tensor and multicomponent spaces should have been converted");
159*f1436e55SToby Isaac   ierr = PetscDTBinomialInt(dim + sp->degree, dim, &Nb);CHKERRQ(ierr);
160*f1436e55SToby Isaac   if (H) {
161*f1436e55SToby Isaac     jet = 2;
162*f1436e55SToby Isaac   } else if (D) {
163*f1436e55SToby Isaac     jet = 1;
164*f1436e55SToby Isaac   } else {
165*f1436e55SToby Isaac     jet = 0;
16620cf1dd8SToby Isaac   }
167*f1436e55SToby Isaac   ierr = PetscDTBinomialInt(dim + jet, dim, &Njet);CHKERRQ(ierr);
168*f1436e55SToby Isaac   ierr = DMGetWorkArray(dm, Nb * Njet * npoints, MPIU_REAL, &pScalar);CHKERRQ(ierr);
169*f1436e55SToby Isaac   // Why are we handling the case degree == 1 specially?  Because we don't want numerical noise when we evaluate hat
170*f1436e55SToby Isaac   // functions at the vertices of a simplex, which happens when we invert the Vandermonde matrix of the PKD basis.
171*f1436e55SToby Isaac   // We don't make any promise about which basis is used.
172*f1436e55SToby Isaac   if (sp->degree == 1) {
173*f1436e55SToby Isaac     ierr = CoordinateBasis(dim, npoints, points, jet, Njet, pScalar);CHKERRQ(ierr);
174*f1436e55SToby Isaac   } else {
175*f1436e55SToby Isaac     ierr = PetscDTPKDEvalJet(dim, npoints, points, sp->degree, jet, pScalar);CHKERRQ(ierr);
17620cf1dd8SToby Isaac   }
17720cf1dd8SToby Isaac   if (B) {
178*f1436e55SToby Isaac     PetscInt p_strl = Nb;
179*f1436e55SToby Isaac     PetscInt b_strl = 1;
1803596293dSMatthew G. Knepley 
181*f1436e55SToby Isaac     PetscInt b_strr = Njet*npoints;
182*f1436e55SToby Isaac     PetscInt p_strr = 1;
183*f1436e55SToby Isaac 
184*f1436e55SToby Isaac     ierr = PetscArrayzero(B, npoints * Nb);CHKERRQ(ierr);
185*f1436e55SToby Isaac     for (PetscInt b = 0; b < Nb; b++) {
186*f1436e55SToby Isaac       for (PetscInt p = 0; p < npoints; p++) {
187*f1436e55SToby Isaac         B[p*p_strl + b*b_strl] = pScalar[b*b_strr + p*p_strr];
1883596293dSMatthew G. Knepley       }
1893596293dSMatthew G. Knepley     }
19020cf1dd8SToby Isaac   }
19120cf1dd8SToby Isaac   if (D) {
192*f1436e55SToby Isaac     PetscInt p_strl = dim*Nb;
193*f1436e55SToby Isaac     PetscInt b_strl = dim;
194*f1436e55SToby Isaac     PetscInt d_strl = 1;
195*f1436e55SToby Isaac 
196*f1436e55SToby Isaac     PetscInt b_strr = Njet*npoints;
197*f1436e55SToby Isaac     PetscInt d_strr = npoints;
198*f1436e55SToby Isaac     PetscInt p_strr = 1;
199*f1436e55SToby Isaac 
200*f1436e55SToby Isaac     ierr = PetscArrayzero(D, npoints * Nb * dim);CHKERRQ(ierr);
201*f1436e55SToby Isaac     for (PetscInt d = 0; d < dim; d++) {
202*f1436e55SToby Isaac       for (PetscInt b = 0; b < Nb; b++) {
203*f1436e55SToby Isaac         for (PetscInt p = 0; p < npoints; p++) {
204*f1436e55SToby Isaac           D[p*p_strl + b*b_strl + d*d_strl] = pScalar[b*b_strr + (1+d)*d_strr + p*p_strr];
20520cf1dd8SToby Isaac         }
20620cf1dd8SToby Isaac       }
20720cf1dd8SToby Isaac     }
20820cf1dd8SToby Isaac   }
20920cf1dd8SToby Isaac   if (H) {
210*f1436e55SToby Isaac     PetscInt p_strl  = dim*dim*Nb;
211*f1436e55SToby Isaac     PetscInt b_strl  = dim*dim;
212*f1436e55SToby Isaac     PetscInt d1_strl = dim;
213*f1436e55SToby Isaac     PetscInt d2_strl = 1;
214*f1436e55SToby Isaac 
215*f1436e55SToby Isaac     PetscInt b_strr = Njet*npoints;
216*f1436e55SToby Isaac     PetscInt j_strr = npoints;
217*f1436e55SToby Isaac     PetscInt p_strr = 1;
218*f1436e55SToby Isaac 
219*f1436e55SToby Isaac     PetscInt *derivs;
220*f1436e55SToby Isaac     ierr = PetscCalloc1(dim, &derivs);CHKERRQ(ierr);
221*f1436e55SToby Isaac     ierr = PetscArrayzero(H, npoints * Nb * dim * dim);CHKERRQ(ierr);
222*f1436e55SToby Isaac     for (PetscInt d1 = 0; d1 < dim; d1++) {
223*f1436e55SToby Isaac       for (PetscInt d2 = 0; d2 < dim; d2++) {
224*f1436e55SToby Isaac         PetscInt j;
225*f1436e55SToby Isaac         derivs[d1]++;
226*f1436e55SToby Isaac         derivs[d2]++;
227*f1436e55SToby Isaac         ierr = PetscDTGradedOrderToIndex(dim, derivs, &j);CHKERRQ(ierr);
228*f1436e55SToby Isaac         derivs[d1]--;
229*f1436e55SToby Isaac         derivs[d2]--;
230*f1436e55SToby Isaac         for (PetscInt b = 0; b < Nb; b++) {
231*f1436e55SToby Isaac           for (PetscInt p = 0; p < npoints; p++) {
232*f1436e55SToby Isaac             H[p*p_strl + b*b_strl + d1*d1_strl + d2*d2_strl] = pScalar[b*b_strr + j*j_strr + p*p_strr];
23320cf1dd8SToby Isaac           }
23420cf1dd8SToby Isaac         }
23520cf1dd8SToby Isaac       }
23620cf1dd8SToby Isaac     }
237*f1436e55SToby Isaac     ierr = PetscFree(derivs);CHKERRQ(ierr);
23820cf1dd8SToby Isaac   }
239*f1436e55SToby Isaac   ierr = DMRestoreWorkArray(dm, Nb * Njet * npoints, MPIU_REAL, &pScalar);CHKERRQ(ierr);
24020cf1dd8SToby Isaac   PetscFunctionReturn(0);
24120cf1dd8SToby Isaac }
24220cf1dd8SToby Isaac 
24320cf1dd8SToby Isaac /*@
24420cf1dd8SToby Isaac   PetscSpacePolynomialSetTensor - Set whether a function space is a space of tensor polynomials (the space is spanned
245*f1436e55SToby Isaac   by polynomials whose degree in each variable is bounded by the given order), as opposed to polynomials (the space is
24620cf1dd8SToby Isaac   spanned by polynomials whose total degree---summing over all variables---is bounded by the given order).
24720cf1dd8SToby Isaac 
24820cf1dd8SToby Isaac   Input Parameters:
24920cf1dd8SToby Isaac + sp     - the function space object
25020cf1dd8SToby Isaac - tensor - PETSC_TRUE for a tensor polynomial space, PETSC_FALSE for a polynomial space
25120cf1dd8SToby Isaac 
2524ab77754SMatthew G. Knepley   Options Database:
2534ab77754SMatthew G. Knepley . -petscspace_poly_tensor <bool> - Whether to use tensor product polynomials in higher dimension
2544ab77754SMatthew G. Knepley 
25529b5c115SMatthew G. Knepley   Level: intermediate
25620cf1dd8SToby Isaac 
25720cf1dd8SToby Isaac .seealso: PetscSpacePolynomialGetTensor(), PetscSpaceSetDegree(), PetscSpaceSetNumVariables()
25820cf1dd8SToby Isaac @*/
25920cf1dd8SToby Isaac PetscErrorCode PetscSpacePolynomialSetTensor(PetscSpace sp, PetscBool tensor)
26020cf1dd8SToby Isaac {
26120cf1dd8SToby Isaac   PetscErrorCode ierr;
26220cf1dd8SToby Isaac 
26320cf1dd8SToby Isaac   PetscFunctionBegin;
26420cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
26520cf1dd8SToby Isaac   ierr = PetscTryMethod(sp,"PetscSpacePolynomialSetTensor_C",(PetscSpace,PetscBool),(sp,tensor));CHKERRQ(ierr);
26620cf1dd8SToby Isaac   PetscFunctionReturn(0);
26720cf1dd8SToby Isaac }
26820cf1dd8SToby Isaac 
26920cf1dd8SToby Isaac /*@
27020cf1dd8SToby Isaac   PetscSpacePolynomialGetTensor - Get whether a function space is a space of tensor polynomials (the space is spanned
27120cf1dd8SToby Isaac   by polynomials whose degree in each variabl is bounded by the given order), as opposed to polynomials (the space is
27220cf1dd8SToby Isaac   spanned by polynomials whose total degree---summing over all variables---is bounded by the given order).
27320cf1dd8SToby Isaac 
27420cf1dd8SToby Isaac   Input Parameters:
27520cf1dd8SToby Isaac . sp     - the function space object
27620cf1dd8SToby Isaac 
27720cf1dd8SToby Isaac   Output Parameters:
27820cf1dd8SToby Isaac . tensor - PETSC_TRUE for a tensor polynomial space, PETSC_FALSE for a polynomial space
27920cf1dd8SToby Isaac 
28029b5c115SMatthew G. Knepley   Level: intermediate
28120cf1dd8SToby Isaac 
28220cf1dd8SToby Isaac .seealso: PetscSpacePolynomialSetTensor(), PetscSpaceSetDegree(), PetscSpaceSetNumVariables()
28320cf1dd8SToby Isaac @*/
28420cf1dd8SToby Isaac PetscErrorCode PetscSpacePolynomialGetTensor(PetscSpace sp, PetscBool *tensor)
28520cf1dd8SToby Isaac {
28620cf1dd8SToby Isaac   PetscErrorCode ierr;
28720cf1dd8SToby Isaac 
28820cf1dd8SToby Isaac   PetscFunctionBegin;
28920cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
29020cf1dd8SToby Isaac   PetscValidPointer(tensor, 2);
29120cf1dd8SToby Isaac   ierr = PetscTryMethod(sp,"PetscSpacePolynomialGetTensor_C",(PetscSpace,PetscBool*),(sp,tensor));CHKERRQ(ierr);
29220cf1dd8SToby Isaac   PetscFunctionReturn(0);
29320cf1dd8SToby Isaac }
29420cf1dd8SToby Isaac 
29520cf1dd8SToby Isaac static PetscErrorCode PetscSpacePolynomialSetTensor_Polynomial(PetscSpace sp, PetscBool tensor)
29620cf1dd8SToby Isaac {
29720cf1dd8SToby Isaac   PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
29820cf1dd8SToby Isaac 
29920cf1dd8SToby Isaac   PetscFunctionBegin;
30020cf1dd8SToby Isaac   poly->tensor = tensor;
30120cf1dd8SToby Isaac   PetscFunctionReturn(0);
30220cf1dd8SToby Isaac }
30320cf1dd8SToby Isaac 
30420cf1dd8SToby Isaac static PetscErrorCode PetscSpacePolynomialGetTensor_Polynomial(PetscSpace sp, PetscBool *tensor)
30520cf1dd8SToby Isaac {
30620cf1dd8SToby Isaac   PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
30720cf1dd8SToby Isaac 
30820cf1dd8SToby Isaac   PetscFunctionBegin;
30920cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
31020cf1dd8SToby Isaac   PetscValidPointer(tensor, 2);
31120cf1dd8SToby Isaac   *tensor = poly->tensor;
31220cf1dd8SToby Isaac   PetscFunctionReturn(0);
31320cf1dd8SToby Isaac }
31420cf1dd8SToby Isaac 
31520cf1dd8SToby Isaac static PetscErrorCode PetscSpaceGetHeightSubspace_Polynomial(PetscSpace sp, PetscInt height, PetscSpace *subsp)
31620cf1dd8SToby Isaac {
31720cf1dd8SToby Isaac   PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
31820cf1dd8SToby Isaac   PetscInt         Nc, dim, order;
31920cf1dd8SToby Isaac   PetscBool        tensor;
32020cf1dd8SToby Isaac   PetscErrorCode   ierr;
32120cf1dd8SToby Isaac 
32220cf1dd8SToby Isaac   PetscFunctionBegin;
32320cf1dd8SToby Isaac   ierr = PetscSpaceGetNumComponents(sp, &Nc);CHKERRQ(ierr);
32420cf1dd8SToby Isaac   ierr = PetscSpaceGetNumVariables(sp, &dim);CHKERRQ(ierr);
32520cf1dd8SToby Isaac   ierr = PetscSpaceGetDegree(sp, &order, NULL);CHKERRQ(ierr);
32620cf1dd8SToby Isaac   ierr = PetscSpacePolynomialGetTensor(sp, &tensor);CHKERRQ(ierr);
3272da392ccSBarry Smith   if (height > dim || height < 0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Asked for space at height %D for dimension %D space", height, dim);
32820cf1dd8SToby Isaac   if (!poly->subspaces) {ierr = PetscCalloc1(dim, &poly->subspaces);CHKERRQ(ierr);}
32920cf1dd8SToby Isaac   if (height <= dim) {
33020cf1dd8SToby Isaac     if (!poly->subspaces[height-1]) {
33120cf1dd8SToby Isaac       PetscSpace  sub;
3323f6b16c7SMatthew G. Knepley       const char *name;
33320cf1dd8SToby Isaac 
33420cf1dd8SToby Isaac       ierr = PetscSpaceCreate(PetscObjectComm((PetscObject) sp), &sub);CHKERRQ(ierr);
3353f6b16c7SMatthew G. Knepley       ierr = PetscObjectGetName((PetscObject) sp,  &name);CHKERRQ(ierr);
3363f6b16c7SMatthew G. Knepley       ierr = PetscObjectSetName((PetscObject) sub,  name);CHKERRQ(ierr);
33720cf1dd8SToby Isaac       ierr = PetscSpaceSetType(sub, PETSCSPACEPOLYNOMIAL);CHKERRQ(ierr);
338d39dd5f5SToby Isaac       ierr = PetscSpaceSetNumComponents(sub, Nc);CHKERRQ(ierr);
339d39dd5f5SToby Isaac       ierr = PetscSpaceSetDegree(sub, order, PETSC_DETERMINE);CHKERRQ(ierr);
34020cf1dd8SToby Isaac       ierr = PetscSpaceSetNumVariables(sub, dim-height);CHKERRQ(ierr);
34120cf1dd8SToby Isaac       ierr = PetscSpacePolynomialSetTensor(sub, tensor);CHKERRQ(ierr);
34220cf1dd8SToby Isaac       ierr = PetscSpaceSetUp(sub);CHKERRQ(ierr);
34320cf1dd8SToby Isaac       poly->subspaces[height-1] = sub;
34420cf1dd8SToby Isaac     }
34520cf1dd8SToby Isaac     *subsp = poly->subspaces[height-1];
34620cf1dd8SToby Isaac   } else {
34720cf1dd8SToby Isaac     *subsp = NULL;
34820cf1dd8SToby Isaac   }
34920cf1dd8SToby Isaac   PetscFunctionReturn(0);
35020cf1dd8SToby Isaac }
35120cf1dd8SToby Isaac 
35229b5c115SMatthew G. Knepley static PetscErrorCode PetscSpaceInitialize_Polynomial(PetscSpace sp)
35320cf1dd8SToby Isaac {
35420cf1dd8SToby Isaac   PetscErrorCode ierr;
35520cf1dd8SToby Isaac 
35620cf1dd8SToby Isaac   PetscFunctionBegin;
35720cf1dd8SToby Isaac   sp->ops->setfromoptions    = PetscSpaceSetFromOptions_Polynomial;
35820cf1dd8SToby Isaac   sp->ops->setup             = PetscSpaceSetUp_Polynomial;
35920cf1dd8SToby Isaac   sp->ops->view              = PetscSpaceView_Polynomial;
36020cf1dd8SToby Isaac   sp->ops->destroy           = PetscSpaceDestroy_Polynomial;
36120cf1dd8SToby Isaac   sp->ops->getdimension      = PetscSpaceGetDimension_Polynomial;
36220cf1dd8SToby Isaac   sp->ops->evaluate          = PetscSpaceEvaluate_Polynomial;
36320cf1dd8SToby Isaac   sp->ops->getheightsubspace = PetscSpaceGetHeightSubspace_Polynomial;
36420cf1dd8SToby Isaac   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePolynomialGetTensor_C", PetscSpacePolynomialGetTensor_Polynomial);CHKERRQ(ierr);
36520cf1dd8SToby Isaac   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePolynomialSetTensor_C", PetscSpacePolynomialSetTensor_Polynomial);CHKERRQ(ierr);
36620cf1dd8SToby Isaac   PetscFunctionReturn(0);
36720cf1dd8SToby Isaac }
36820cf1dd8SToby Isaac 
36920cf1dd8SToby Isaac /*MC
37020cf1dd8SToby Isaac   PETSCSPACEPOLYNOMIAL = "poly" - A PetscSpace object that encapsulates a polynomial space, e.g. P1 is the space of
37120cf1dd8SToby Isaac   linear polynomials. The space is replicated for each component.
37220cf1dd8SToby Isaac 
37320cf1dd8SToby Isaac   Level: intermediate
37420cf1dd8SToby Isaac 
37520cf1dd8SToby Isaac .seealso: PetscSpaceType, PetscSpaceCreate(), PetscSpaceSetType()
37620cf1dd8SToby Isaac M*/
37720cf1dd8SToby Isaac 
37820cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode PetscSpaceCreate_Polynomial(PetscSpace sp)
37920cf1dd8SToby Isaac {
38020cf1dd8SToby Isaac   PetscSpace_Poly *poly;
38120cf1dd8SToby Isaac   PetscErrorCode   ierr;
38220cf1dd8SToby Isaac 
38320cf1dd8SToby Isaac   PetscFunctionBegin;
38420cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
38520cf1dd8SToby Isaac   ierr     = PetscNewLog(sp,&poly);CHKERRQ(ierr);
38620cf1dd8SToby Isaac   sp->data = poly;
38720cf1dd8SToby Isaac 
38820cf1dd8SToby Isaac   poly->tensor    = PETSC_FALSE;
38920cf1dd8SToby Isaac   poly->subspaces = NULL;
39020cf1dd8SToby Isaac 
39120cf1dd8SToby Isaac   ierr = PetscSpaceInitialize_Polynomial(sp);CHKERRQ(ierr);
39220cf1dd8SToby Isaac   PetscFunctionReturn(0);
39320cf1dd8SToby Isaac }
39420cf1dd8SToby Isaac 
395