xref: /petsc/src/dm/dt/space/impls/poly/spacepoly.c (revision cac4c232dc4f93991e342196e27ef7b0655dac7b)
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 
720cf1dd8SToby Isaac   PetscFunctionBegin;
89566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHead(PetscOptionsObject,"PetscSpace polynomial options"));
99566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-petscspace_poly_tensor", "Use the tensor product polynomials", "PetscSpacePolynomialSetTensor", poly->tensor, &poly->tensor, NULL));
109566063dSJacob Faibussowitsch   PetscCall(PetscOptionsTail());
1120cf1dd8SToby Isaac   PetscFunctionReturn(0);
1220cf1dd8SToby Isaac }
1320cf1dd8SToby Isaac 
14d9bac1caSLisandro Dalcin static PetscErrorCode PetscSpacePolynomialView_Ascii(PetscSpace sp, PetscViewer v)
1520cf1dd8SToby Isaac {
1620cf1dd8SToby Isaac   PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
1720cf1dd8SToby Isaac 
1820cf1dd8SToby Isaac   PetscFunctionBegin;
199566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(v, "%s space of degree %D\n", poly->tensor ? "Tensor polynomial" : "Polynomial", sp->degree));
2020cf1dd8SToby Isaac   PetscFunctionReturn(0);
2120cf1dd8SToby Isaac }
2220cf1dd8SToby Isaac 
2329b5c115SMatthew G. Knepley static PetscErrorCode PetscSpaceView_Polynomial(PetscSpace sp, PetscViewer viewer)
2420cf1dd8SToby Isaac {
2520cf1dd8SToby Isaac   PetscBool      iascii;
2620cf1dd8SToby Isaac 
2720cf1dd8SToby Isaac   PetscFunctionBegin;
2820cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
2920cf1dd8SToby Isaac   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
309566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii));
319566063dSJacob Faibussowitsch   if (iascii) PetscCall(PetscSpacePolynomialView_Ascii(sp, viewer));
3220cf1dd8SToby Isaac   PetscFunctionReturn(0);
3320cf1dd8SToby Isaac }
3420cf1dd8SToby Isaac 
3529b5c115SMatthew G. Knepley static PetscErrorCode PetscSpaceDestroy_Polynomial(PetscSpace sp)
3620cf1dd8SToby Isaac {
3720cf1dd8SToby Isaac   PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
3820cf1dd8SToby Isaac 
3920cf1dd8SToby Isaac   PetscFunctionBegin;
409566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePolynomialGetTensor_C", NULL));
419566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePolynomialSetTensor_C", NULL));
4220cf1dd8SToby Isaac   if (poly->subspaces) {
4320cf1dd8SToby Isaac     PetscInt d;
4420cf1dd8SToby Isaac 
4520cf1dd8SToby Isaac     for (d = 0; d < sp->Nv; ++d) {
469566063dSJacob Faibussowitsch       PetscCall(PetscSpaceDestroy(&poly->subspaces[d]));
4720cf1dd8SToby Isaac     }
4820cf1dd8SToby Isaac   }
499566063dSJacob Faibussowitsch   PetscCall(PetscFree(poly->subspaces));
509566063dSJacob Faibussowitsch   PetscCall(PetscFree(poly));
5120cf1dd8SToby Isaac   PetscFunctionReturn(0);
5220cf1dd8SToby Isaac }
5320cf1dd8SToby Isaac 
54f1436e55SToby Isaac static PetscErrorCode PetscSpaceSetUp_Polynomial(PetscSpace sp)
55f1436e55SToby Isaac {
56f1436e55SToby Isaac   PetscSpace_Poly *poly    = (PetscSpace_Poly *) sp->data;
57f1436e55SToby Isaac 
58f1436e55SToby Isaac   PetscFunctionBegin;
59f1436e55SToby Isaac   if (poly->setupCalled) PetscFunctionReturn(0);
60f1436e55SToby Isaac   if (sp->Nv <= 1) {
61f1436e55SToby Isaac     poly->tensor = PETSC_FALSE;
62f1436e55SToby Isaac   }
63f1436e55SToby Isaac   if (sp->Nc != 1) {
64f1436e55SToby Isaac     PetscInt    Nc = sp->Nc;
65f1436e55SToby Isaac     PetscBool   tensor = poly->tensor;
66f1436e55SToby Isaac     PetscInt    Nv = sp->Nv;
67f1436e55SToby Isaac     PetscInt    degree = sp->degree;
68417c287bSToby Isaac     const char *prefix;
69417c287bSToby Isaac     const char *name;
70417c287bSToby Isaac     char        subname[PETSC_MAX_PATH_LEN];
71f1436e55SToby Isaac     PetscSpace  subsp;
72f1436e55SToby Isaac 
739566063dSJacob Faibussowitsch     PetscCall(PetscSpaceSetType(sp, PETSCSPACESUM));
749566063dSJacob Faibussowitsch     PetscCall(PetscSpaceSumSetNumSubspaces(sp, Nc));
759566063dSJacob Faibussowitsch     PetscCall(PetscSpaceCreate(PetscObjectComm((PetscObject)sp), &subsp));
769566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetOptionsPrefix((PetscObject)sp, &prefix));
779566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)subsp, prefix));
789566063dSJacob Faibussowitsch     PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)subsp, "sumcomp_"));
79417c287bSToby Isaac     if (((PetscObject)sp)->name) {
809566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetName((PetscObject)sp, &name));
819566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(subname, PETSC_MAX_PATH_LEN-1, "%s sum component", name));
829566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)subsp, subname));
83417c287bSToby Isaac     } else {
849566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)subsp, "sum component"));
85417c287bSToby Isaac     }
869566063dSJacob Faibussowitsch     PetscCall(PetscSpaceSetType(subsp, PETSCSPACEPOLYNOMIAL));
879566063dSJacob Faibussowitsch     PetscCall(PetscSpaceSetDegree(subsp, degree, PETSC_DETERMINE));
889566063dSJacob Faibussowitsch     PetscCall(PetscSpaceSetNumComponents(subsp, 1));
899566063dSJacob Faibussowitsch     PetscCall(PetscSpaceSetNumVariables(subsp, Nv));
909566063dSJacob Faibussowitsch     PetscCall(PetscSpacePolynomialSetTensor(subsp, tensor));
919566063dSJacob Faibussowitsch     PetscCall(PetscSpaceSetUp(subsp));
92f1436e55SToby Isaac     for (PetscInt i = 0; i < Nc; i++) {
939566063dSJacob Faibussowitsch       PetscCall(PetscSpaceSumSetSubspace(sp, i, subsp));
94f1436e55SToby Isaac     }
959566063dSJacob Faibussowitsch     PetscCall(PetscSpaceDestroy(&subsp));
969566063dSJacob Faibussowitsch     PetscCall(PetscSpaceSetUp(sp));
97f1436e55SToby Isaac     PetscFunctionReturn(0);
98f1436e55SToby Isaac   }
99f1436e55SToby Isaac   if (poly->tensor) {
100f1436e55SToby Isaac     sp->maxDegree = PETSC_DETERMINE;
1019566063dSJacob Faibussowitsch     PetscCall(PetscSpaceSetType(sp, PETSCSPACETENSOR));
1029566063dSJacob Faibussowitsch     PetscCall(PetscSpaceSetUp(sp));
103f1436e55SToby Isaac     PetscFunctionReturn(0);
104f1436e55SToby Isaac   }
1052c71b3e2SJacob Faibussowitsch   PetscCheckFalse(sp->degree < 0,PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_OUTOFRANGE, "Negative degree %D invalid", sp->degree);
106f1436e55SToby Isaac   sp->maxDegree = sp->degree;
107f1436e55SToby Isaac   poly->setupCalled = PETSC_TRUE;
108f1436e55SToby Isaac   PetscFunctionReturn(0);
109f1436e55SToby Isaac }
110f1436e55SToby Isaac 
11129b5c115SMatthew G. Knepley static PetscErrorCode PetscSpaceGetDimension_Polynomial(PetscSpace sp, PetscInt *dim)
11220cf1dd8SToby Isaac {
11320cf1dd8SToby Isaac   PetscInt         deg  = sp->degree;
114f1436e55SToby Isaac   PetscInt         n    = sp->Nv;
11520cf1dd8SToby Isaac 
11620cf1dd8SToby Isaac   PetscFunctionBegin;
1179566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(n + deg, n, dim));
118f1436e55SToby Isaac   *dim *= sp->Nc;
11920cf1dd8SToby Isaac   PetscFunctionReturn(0);
12020cf1dd8SToby Isaac }
12120cf1dd8SToby Isaac 
122f1436e55SToby Isaac static PetscErrorCode CoordinateBasis(PetscInt dim, PetscInt npoints, const PetscReal points[], PetscInt jet, PetscInt Njet, PetscReal pScalar[])
12320cf1dd8SToby Isaac {
12420cf1dd8SToby Isaac   PetscFunctionBegin;
1259566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(pScalar, (1 + dim) * Njet * npoints));
126f1436e55SToby Isaac   for (PetscInt b = 0; b < 1 + dim; b++) {
127f1436e55SToby Isaac     for (PetscInt j = 0; j < PetscMin(1 + dim, Njet); j++) {
128f1436e55SToby Isaac       if (j == 0) {
129f1436e55SToby Isaac         if (b == 0) {
130f1436e55SToby Isaac           for (PetscInt pt = 0; pt < npoints; pt++) {
131f1436e55SToby Isaac             pScalar[b * Njet * npoints + j * npoints + pt] = 1.;
132f1436e55SToby Isaac           }
13320cf1dd8SToby Isaac         } else {
134f1436e55SToby Isaac           for (PetscInt pt = 0; pt < npoints; pt++) {
135f1436e55SToby Isaac             pScalar[b * Njet * npoints + j * npoints + pt] = points[pt * dim + (b-1)];
136f1436e55SToby Isaac           }
137f1436e55SToby Isaac         }
138f1436e55SToby Isaac       } else if (j == b) {
139f1436e55SToby Isaac         for (PetscInt pt = 0; pt < npoints; pt++) {
140f1436e55SToby Isaac           pScalar[b * Njet * npoints + j * npoints + pt] = 1.;
141f1436e55SToby Isaac         }
142f1436e55SToby Isaac       }
14320cf1dd8SToby Isaac     }
14420cf1dd8SToby Isaac   }
14520cf1dd8SToby Isaac   PetscFunctionReturn(0);
14620cf1dd8SToby Isaac }
14720cf1dd8SToby Isaac 
14829b5c115SMatthew G. Knepley static PetscErrorCode PetscSpaceEvaluate_Polynomial(PetscSpace sp, PetscInt npoints, const PetscReal points[], PetscReal B[], PetscReal D[], PetscReal H[])
14920cf1dd8SToby Isaac {
15020cf1dd8SToby Isaac   PetscSpace_Poly *poly    = (PetscSpace_Poly *) sp->data;
15120cf1dd8SToby Isaac   DM               dm      = sp->dm;
15220cf1dd8SToby Isaac   PetscInt         dim     = sp->Nv;
153f1436e55SToby Isaac   PetscInt         Nb, jet, Njet;
154f1436e55SToby Isaac   PetscReal       *pScalar;
15520cf1dd8SToby Isaac 
15620cf1dd8SToby Isaac   PetscFunctionBegin;
157f1436e55SToby Isaac   if (!poly->setupCalled) {
1589566063dSJacob Faibussowitsch     PetscCall(PetscSpaceSetUp(sp));
1599566063dSJacob Faibussowitsch     PetscCall(PetscSpaceEvaluate(sp, npoints, points, B, D, H));
160f1436e55SToby Isaac     PetscFunctionReturn(0);
16120cf1dd8SToby Isaac   }
1622c71b3e2SJacob Faibussowitsch   PetscCheckFalse(poly->tensor || sp->Nc != 1,PETSC_COMM_SELF, PETSC_ERR_PLIB, "tensor and multicomponent spaces should have been converted");
1639566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(dim + sp->degree, dim, &Nb));
164f1436e55SToby Isaac   if (H) {
165f1436e55SToby Isaac     jet = 2;
166f1436e55SToby Isaac   } else if (D) {
167f1436e55SToby Isaac     jet = 1;
168f1436e55SToby Isaac   } else {
169f1436e55SToby Isaac     jet = 0;
17020cf1dd8SToby Isaac   }
1719566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(dim + jet, dim, &Njet));
1729566063dSJacob Faibussowitsch   PetscCall(DMGetWorkArray(dm, Nb * Njet * npoints, MPIU_REAL, &pScalar));
173f1436e55SToby Isaac   // Why are we handling the case degree == 1 specially?  Because we don't want numerical noise when we evaluate hat
174f1436e55SToby Isaac   // functions at the vertices of a simplex, which happens when we invert the Vandermonde matrix of the PKD basis.
175f1436e55SToby Isaac   // We don't make any promise about which basis is used.
176f1436e55SToby Isaac   if (sp->degree == 1) {
1779566063dSJacob Faibussowitsch     PetscCall(CoordinateBasis(dim, npoints, points, jet, Njet, pScalar));
178f1436e55SToby Isaac   } else {
1799566063dSJacob Faibussowitsch     PetscCall(PetscDTPKDEvalJet(dim, npoints, points, sp->degree, jet, pScalar));
18020cf1dd8SToby Isaac   }
18120cf1dd8SToby Isaac   if (B) {
182f1436e55SToby Isaac     PetscInt p_strl = Nb;
183f1436e55SToby Isaac     PetscInt b_strl = 1;
1843596293dSMatthew G. Knepley 
185f1436e55SToby Isaac     PetscInt b_strr = Njet*npoints;
186f1436e55SToby Isaac     PetscInt p_strr = 1;
187f1436e55SToby Isaac 
1889566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(B, npoints * Nb));
189f1436e55SToby Isaac     for (PetscInt b = 0; b < Nb; b++) {
190f1436e55SToby Isaac       for (PetscInt p = 0; p < npoints; p++) {
191f1436e55SToby Isaac         B[p*p_strl + b*b_strl] = pScalar[b*b_strr + p*p_strr];
1923596293dSMatthew G. Knepley       }
1933596293dSMatthew G. Knepley     }
19420cf1dd8SToby Isaac   }
19520cf1dd8SToby Isaac   if (D) {
196f1436e55SToby Isaac     PetscInt p_strl = dim*Nb;
197f1436e55SToby Isaac     PetscInt b_strl = dim;
198f1436e55SToby Isaac     PetscInt d_strl = 1;
199f1436e55SToby Isaac 
200f1436e55SToby Isaac     PetscInt b_strr = Njet*npoints;
201f1436e55SToby Isaac     PetscInt d_strr = npoints;
202f1436e55SToby Isaac     PetscInt p_strr = 1;
203f1436e55SToby Isaac 
2049566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(D, npoints * Nb * dim));
205f1436e55SToby Isaac     for (PetscInt d = 0; d < dim; d++) {
206f1436e55SToby Isaac       for (PetscInt b = 0; b < Nb; b++) {
207f1436e55SToby Isaac         for (PetscInt p = 0; p < npoints; p++) {
208f1436e55SToby Isaac           D[p*p_strl + b*b_strl + d*d_strl] = pScalar[b*b_strr + (1+d)*d_strr + p*p_strr];
20920cf1dd8SToby Isaac         }
21020cf1dd8SToby Isaac       }
21120cf1dd8SToby Isaac     }
21220cf1dd8SToby Isaac   }
21320cf1dd8SToby Isaac   if (H) {
214f1436e55SToby Isaac     PetscInt p_strl  = dim*dim*Nb;
215f1436e55SToby Isaac     PetscInt b_strl  = dim*dim;
216f1436e55SToby Isaac     PetscInt d1_strl = dim;
217f1436e55SToby Isaac     PetscInt d2_strl = 1;
218f1436e55SToby Isaac 
219f1436e55SToby Isaac     PetscInt b_strr = Njet*npoints;
220f1436e55SToby Isaac     PetscInt j_strr = npoints;
221f1436e55SToby Isaac     PetscInt p_strr = 1;
222f1436e55SToby Isaac 
223f1436e55SToby Isaac     PetscInt *derivs;
2249566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(dim, &derivs));
2259566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(H, npoints * Nb * dim * dim));
226f1436e55SToby Isaac     for (PetscInt d1 = 0; d1 < dim; d1++) {
227f1436e55SToby Isaac       for (PetscInt d2 = 0; d2 < dim; d2++) {
228f1436e55SToby Isaac         PetscInt j;
229f1436e55SToby Isaac         derivs[d1]++;
230f1436e55SToby Isaac         derivs[d2]++;
2319566063dSJacob Faibussowitsch         PetscCall(PetscDTGradedOrderToIndex(dim, derivs, &j));
232f1436e55SToby Isaac         derivs[d1]--;
233f1436e55SToby Isaac         derivs[d2]--;
234f1436e55SToby Isaac         for (PetscInt b = 0; b < Nb; b++) {
235f1436e55SToby Isaac           for (PetscInt p = 0; p < npoints; p++) {
236f1436e55SToby Isaac             H[p*p_strl + b*b_strl + d1*d1_strl + d2*d2_strl] = pScalar[b*b_strr + j*j_strr + p*p_strr];
23720cf1dd8SToby Isaac           }
23820cf1dd8SToby Isaac         }
23920cf1dd8SToby Isaac       }
24020cf1dd8SToby Isaac     }
2419566063dSJacob Faibussowitsch     PetscCall(PetscFree(derivs));
24220cf1dd8SToby Isaac   }
2439566063dSJacob Faibussowitsch   PetscCall(DMRestoreWorkArray(dm, Nb * Njet * npoints, MPIU_REAL, &pScalar));
24420cf1dd8SToby Isaac   PetscFunctionReturn(0);
24520cf1dd8SToby Isaac }
24620cf1dd8SToby Isaac 
24720cf1dd8SToby Isaac /*@
24820cf1dd8SToby Isaac   PetscSpacePolynomialSetTensor - Set whether a function space is a space of tensor polynomials (the space is spanned
249f1436e55SToby Isaac   by polynomials whose degree in each variable is bounded by the given order), as opposed to polynomials (the space is
25020cf1dd8SToby Isaac   spanned by polynomials whose total degree---summing over all variables---is bounded by the given order).
25120cf1dd8SToby Isaac 
25220cf1dd8SToby Isaac   Input Parameters:
25320cf1dd8SToby Isaac + sp     - the function space object
25420cf1dd8SToby Isaac - tensor - PETSC_TRUE for a tensor polynomial space, PETSC_FALSE for a polynomial space
25520cf1dd8SToby Isaac 
2564ab77754SMatthew G. Knepley   Options Database:
2574ab77754SMatthew G. Knepley . -petscspace_poly_tensor <bool> - Whether to use tensor product polynomials in higher dimension
2584ab77754SMatthew G. Knepley 
25929b5c115SMatthew G. Knepley   Level: intermediate
26020cf1dd8SToby Isaac 
26120cf1dd8SToby Isaac .seealso: PetscSpacePolynomialGetTensor(), PetscSpaceSetDegree(), PetscSpaceSetNumVariables()
26220cf1dd8SToby Isaac @*/
26320cf1dd8SToby Isaac PetscErrorCode PetscSpacePolynomialSetTensor(PetscSpace sp, PetscBool tensor)
26420cf1dd8SToby Isaac {
26520cf1dd8SToby Isaac   PetscFunctionBegin;
26620cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
267*cac4c232SBarry Smith   PetscTryMethod(sp,"PetscSpacePolynomialSetTensor_C",(PetscSpace,PetscBool),(sp,tensor));
26820cf1dd8SToby Isaac   PetscFunctionReturn(0);
26920cf1dd8SToby Isaac }
27020cf1dd8SToby Isaac 
27120cf1dd8SToby Isaac /*@
27220cf1dd8SToby Isaac   PetscSpacePolynomialGetTensor - Get whether a function space is a space of tensor polynomials (the space is spanned
27320cf1dd8SToby Isaac   by polynomials whose degree in each variabl is bounded by the given order), as opposed to polynomials (the space is
27420cf1dd8SToby Isaac   spanned by polynomials whose total degree---summing over all variables---is bounded by the given order).
27520cf1dd8SToby Isaac 
27620cf1dd8SToby Isaac   Input Parameters:
27720cf1dd8SToby Isaac . sp     - the function space object
27820cf1dd8SToby Isaac 
27920cf1dd8SToby Isaac   Output Parameters:
28020cf1dd8SToby Isaac . tensor - PETSC_TRUE for a tensor polynomial space, PETSC_FALSE for a polynomial space
28120cf1dd8SToby Isaac 
28229b5c115SMatthew G. Knepley   Level: intermediate
28320cf1dd8SToby Isaac 
28420cf1dd8SToby Isaac .seealso: PetscSpacePolynomialSetTensor(), PetscSpaceSetDegree(), PetscSpaceSetNumVariables()
28520cf1dd8SToby Isaac @*/
28620cf1dd8SToby Isaac PetscErrorCode PetscSpacePolynomialGetTensor(PetscSpace sp, PetscBool *tensor)
28720cf1dd8SToby Isaac {
28820cf1dd8SToby Isaac   PetscFunctionBegin;
28920cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
290dadcf809SJacob Faibussowitsch   PetscValidBoolPointer(tensor, 2);
291*cac4c232SBarry Smith   PetscTryMethod(sp,"PetscSpacePolynomialGetTensor_C",(PetscSpace,PetscBool*),(sp,tensor));
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);
310dadcf809SJacob Faibussowitsch   PetscValidBoolPointer(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 
32120cf1dd8SToby Isaac   PetscFunctionBegin;
3229566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetNumComponents(sp, &Nc));
3239566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetNumVariables(sp, &dim));
3249566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetDegree(sp, &order, NULL));
3259566063dSJacob Faibussowitsch   PetscCall(PetscSpacePolynomialGetTensor(sp, &tensor));
3262c71b3e2SJacob 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);
3279566063dSJacob Faibussowitsch   if (!poly->subspaces) PetscCall(PetscCalloc1(dim, &poly->subspaces));
32820cf1dd8SToby Isaac   if (height <= dim) {
32920cf1dd8SToby Isaac     if (!poly->subspaces[height-1]) {
33020cf1dd8SToby Isaac       PetscSpace  sub;
3313f6b16c7SMatthew G. Knepley       const char *name;
33220cf1dd8SToby Isaac 
3339566063dSJacob Faibussowitsch       PetscCall(PetscSpaceCreate(PetscObjectComm((PetscObject) sp), &sub));
3349566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetName((PetscObject) sp,  &name));
3359566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject) sub,  name));
3369566063dSJacob Faibussowitsch       PetscCall(PetscSpaceSetType(sub, PETSCSPACEPOLYNOMIAL));
3379566063dSJacob Faibussowitsch       PetscCall(PetscSpaceSetNumComponents(sub, Nc));
3389566063dSJacob Faibussowitsch       PetscCall(PetscSpaceSetDegree(sub, order, PETSC_DETERMINE));
3399566063dSJacob Faibussowitsch       PetscCall(PetscSpaceSetNumVariables(sub, dim-height));
3409566063dSJacob Faibussowitsch       PetscCall(PetscSpacePolynomialSetTensor(sub, tensor));
3419566063dSJacob Faibussowitsch       PetscCall(PetscSpaceSetUp(sub));
34220cf1dd8SToby Isaac       poly->subspaces[height-1] = sub;
34320cf1dd8SToby Isaac     }
34420cf1dd8SToby Isaac     *subsp = poly->subspaces[height-1];
34520cf1dd8SToby Isaac   } else {
34620cf1dd8SToby Isaac     *subsp = NULL;
34720cf1dd8SToby Isaac   }
34820cf1dd8SToby Isaac   PetscFunctionReturn(0);
34920cf1dd8SToby Isaac }
35020cf1dd8SToby Isaac 
35129b5c115SMatthew G. Knepley static PetscErrorCode PetscSpaceInitialize_Polynomial(PetscSpace sp)
35220cf1dd8SToby Isaac {
35320cf1dd8SToby Isaac   PetscFunctionBegin;
35420cf1dd8SToby Isaac   sp->ops->setfromoptions    = PetscSpaceSetFromOptions_Polynomial;
35520cf1dd8SToby Isaac   sp->ops->setup             = PetscSpaceSetUp_Polynomial;
35620cf1dd8SToby Isaac   sp->ops->view              = PetscSpaceView_Polynomial;
35720cf1dd8SToby Isaac   sp->ops->destroy           = PetscSpaceDestroy_Polynomial;
35820cf1dd8SToby Isaac   sp->ops->getdimension      = PetscSpaceGetDimension_Polynomial;
35920cf1dd8SToby Isaac   sp->ops->evaluate          = PetscSpaceEvaluate_Polynomial;
36020cf1dd8SToby Isaac   sp->ops->getheightsubspace = PetscSpaceGetHeightSubspace_Polynomial;
3619566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePolynomialGetTensor_C", PetscSpacePolynomialGetTensor_Polynomial));
3629566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePolynomialSetTensor_C", PetscSpacePolynomialSetTensor_Polynomial));
36320cf1dd8SToby Isaac   PetscFunctionReturn(0);
36420cf1dd8SToby Isaac }
36520cf1dd8SToby Isaac 
36620cf1dd8SToby Isaac /*MC
36720cf1dd8SToby Isaac   PETSCSPACEPOLYNOMIAL = "poly" - A PetscSpace object that encapsulates a polynomial space, e.g. P1 is the space of
36820cf1dd8SToby Isaac   linear polynomials. The space is replicated for each component.
36920cf1dd8SToby Isaac 
37020cf1dd8SToby Isaac   Level: intermediate
37120cf1dd8SToby Isaac 
37220cf1dd8SToby Isaac .seealso: PetscSpaceType, PetscSpaceCreate(), PetscSpaceSetType()
37320cf1dd8SToby Isaac M*/
37420cf1dd8SToby Isaac 
37520cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode PetscSpaceCreate_Polynomial(PetscSpace sp)
37620cf1dd8SToby Isaac {
37720cf1dd8SToby Isaac   PetscSpace_Poly *poly;
37820cf1dd8SToby Isaac 
37920cf1dd8SToby Isaac   PetscFunctionBegin;
38020cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
3819566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(sp,&poly));
38220cf1dd8SToby Isaac   sp->data = poly;
38320cf1dd8SToby Isaac 
38420cf1dd8SToby Isaac   poly->tensor    = PETSC_FALSE;
38520cf1dd8SToby Isaac   poly->subspaces = NULL;
38620cf1dd8SToby Isaac 
3879566063dSJacob Faibussowitsch   PetscCall(PetscSpaceInitialize_Polynomial(sp));
38820cf1dd8SToby Isaac   PetscFunctionReturn(0);
38920cf1dd8SToby Isaac }
390