xref: /petsc/src/dm/dt/space/impls/poly/spacepoly.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
120cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
220cf1dd8SToby Isaac 
39371c9d4SSatish Balay static PetscErrorCode PetscSpaceSetFromOptions_Polynomial(PetscSpace sp, PetscOptionItems *PetscOptionsObject) {
420cf1dd8SToby Isaac   PetscSpace_Poly *poly = (PetscSpace_Poly *)sp->data;
520cf1dd8SToby Isaac 
620cf1dd8SToby Isaac   PetscFunctionBegin;
7d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "PetscSpace polynomial options");
89566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-petscspace_poly_tensor", "Use the tensor product polynomials", "PetscSpacePolynomialSetTensor", poly->tensor, &poly->tensor, NULL));
9d0609cedSBarry Smith   PetscOptionsHeadEnd();
1020cf1dd8SToby Isaac   PetscFunctionReturn(0);
1120cf1dd8SToby Isaac }
1220cf1dd8SToby Isaac 
139371c9d4SSatish Balay static PetscErrorCode PetscSpacePolynomialView_Ascii(PetscSpace sp, PetscViewer v) {
1420cf1dd8SToby Isaac   PetscSpace_Poly *poly = (PetscSpace_Poly *)sp->data;
1520cf1dd8SToby Isaac 
1620cf1dd8SToby Isaac   PetscFunctionBegin;
1763a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(v, "%s space of degree %" PetscInt_FMT "\n", poly->tensor ? "Tensor polynomial" : "Polynomial", sp->degree));
1820cf1dd8SToby Isaac   PetscFunctionReturn(0);
1920cf1dd8SToby Isaac }
2020cf1dd8SToby Isaac 
219371c9d4SSatish Balay static PetscErrorCode PetscSpaceView_Polynomial(PetscSpace sp, PetscViewer viewer) {
2220cf1dd8SToby Isaac   PetscBool iascii;
2320cf1dd8SToby Isaac 
2420cf1dd8SToby Isaac   PetscFunctionBegin;
2520cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
2620cf1dd8SToby Isaac   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
279566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
289566063dSJacob Faibussowitsch   if (iascii) PetscCall(PetscSpacePolynomialView_Ascii(sp, viewer));
2920cf1dd8SToby Isaac   PetscFunctionReturn(0);
3020cf1dd8SToby Isaac }
3120cf1dd8SToby Isaac 
329371c9d4SSatish Balay static PetscErrorCode PetscSpaceDestroy_Polynomial(PetscSpace sp) {
3320cf1dd8SToby Isaac   PetscSpace_Poly *poly = (PetscSpace_Poly *)sp->data;
3420cf1dd8SToby Isaac 
3520cf1dd8SToby Isaac   PetscFunctionBegin;
369566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscSpacePolynomialGetTensor_C", NULL));
379566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscSpacePolynomialSetTensor_C", NULL));
3820cf1dd8SToby Isaac   if (poly->subspaces) {
3920cf1dd8SToby Isaac     PetscInt d;
4020cf1dd8SToby Isaac 
41*48a46eb9SPierre Jolivet     for (d = 0; d < sp->Nv; ++d) PetscCall(PetscSpaceDestroy(&poly->subspaces[d]));
4220cf1dd8SToby Isaac   }
439566063dSJacob Faibussowitsch   PetscCall(PetscFree(poly->subspaces));
449566063dSJacob Faibussowitsch   PetscCall(PetscFree(poly));
4520cf1dd8SToby Isaac   PetscFunctionReturn(0);
4620cf1dd8SToby Isaac }
4720cf1dd8SToby Isaac 
489371c9d4SSatish Balay static PetscErrorCode PetscSpaceSetUp_Polynomial(PetscSpace sp) {
49f1436e55SToby Isaac   PetscSpace_Poly *poly = (PetscSpace_Poly *)sp->data;
50f1436e55SToby Isaac 
51f1436e55SToby Isaac   PetscFunctionBegin;
52f1436e55SToby Isaac   if (poly->setupCalled) PetscFunctionReturn(0);
539371c9d4SSatish Balay   if (sp->Nv <= 1) { poly->tensor = PETSC_FALSE; }
54f1436e55SToby Isaac   if (sp->Nc != 1) {
55f1436e55SToby Isaac     PetscInt    Nc     = sp->Nc;
56f1436e55SToby Isaac     PetscBool   tensor = poly->tensor;
57f1436e55SToby Isaac     PetscInt    Nv     = sp->Nv;
58f1436e55SToby Isaac     PetscInt    degree = sp->degree;
59417c287bSToby Isaac     const char *prefix;
60417c287bSToby Isaac     const char *name;
61417c287bSToby Isaac     char        subname[PETSC_MAX_PATH_LEN];
62f1436e55SToby Isaac     PetscSpace  subsp;
63f1436e55SToby Isaac 
649566063dSJacob Faibussowitsch     PetscCall(PetscSpaceSetType(sp, PETSCSPACESUM));
659566063dSJacob Faibussowitsch     PetscCall(PetscSpaceSumSetNumSubspaces(sp, Nc));
669566063dSJacob Faibussowitsch     PetscCall(PetscSpaceCreate(PetscObjectComm((PetscObject)sp), &subsp));
679566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetOptionsPrefix((PetscObject)sp, &prefix));
689566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)subsp, prefix));
699566063dSJacob Faibussowitsch     PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)subsp, "sumcomp_"));
70417c287bSToby Isaac     if (((PetscObject)sp)->name) {
719566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetName((PetscObject)sp, &name));
729566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(subname, PETSC_MAX_PATH_LEN - 1, "%s sum component", name));
739566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)subsp, subname));
741baa6e33SBarry Smith     } else PetscCall(PetscObjectSetName((PetscObject)subsp, "sum component"));
759566063dSJacob Faibussowitsch     PetscCall(PetscSpaceSetType(subsp, PETSCSPACEPOLYNOMIAL));
769566063dSJacob Faibussowitsch     PetscCall(PetscSpaceSetDegree(subsp, degree, PETSC_DETERMINE));
779566063dSJacob Faibussowitsch     PetscCall(PetscSpaceSetNumComponents(subsp, 1));
789566063dSJacob Faibussowitsch     PetscCall(PetscSpaceSetNumVariables(subsp, Nv));
799566063dSJacob Faibussowitsch     PetscCall(PetscSpacePolynomialSetTensor(subsp, tensor));
809566063dSJacob Faibussowitsch     PetscCall(PetscSpaceSetUp(subsp));
81*48a46eb9SPierre Jolivet     for (PetscInt i = 0; i < Nc; i++) PetscCall(PetscSpaceSumSetSubspace(sp, i, subsp));
829566063dSJacob Faibussowitsch     PetscCall(PetscSpaceDestroy(&subsp));
839566063dSJacob Faibussowitsch     PetscCall(PetscSpaceSetUp(sp));
84f1436e55SToby Isaac     PetscFunctionReturn(0);
85f1436e55SToby Isaac   }
86f1436e55SToby Isaac   if (poly->tensor) {
87f1436e55SToby Isaac     sp->maxDegree = PETSC_DETERMINE;
889566063dSJacob Faibussowitsch     PetscCall(PetscSpaceSetType(sp, PETSCSPACETENSOR));
899566063dSJacob Faibussowitsch     PetscCall(PetscSpaceSetUp(sp));
90f1436e55SToby Isaac     PetscFunctionReturn(0);
91f1436e55SToby Isaac   }
9263a3b9bcSJacob Faibussowitsch   PetscCheck(sp->degree >= 0, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_OUTOFRANGE, "Negative degree %" PetscInt_FMT " invalid", sp->degree);
93f1436e55SToby Isaac   sp->maxDegree     = sp->degree;
94f1436e55SToby Isaac   poly->setupCalled = PETSC_TRUE;
95f1436e55SToby Isaac   PetscFunctionReturn(0);
96f1436e55SToby Isaac }
97f1436e55SToby Isaac 
989371c9d4SSatish Balay static PetscErrorCode PetscSpaceGetDimension_Polynomial(PetscSpace sp, PetscInt *dim) {
9920cf1dd8SToby Isaac   PetscInt deg = sp->degree;
100f1436e55SToby Isaac   PetscInt n   = sp->Nv;
10120cf1dd8SToby Isaac 
10220cf1dd8SToby Isaac   PetscFunctionBegin;
1039566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(n + deg, n, dim));
104f1436e55SToby Isaac   *dim *= sp->Nc;
10520cf1dd8SToby Isaac   PetscFunctionReturn(0);
10620cf1dd8SToby Isaac }
10720cf1dd8SToby Isaac 
1089371c9d4SSatish Balay static PetscErrorCode CoordinateBasis(PetscInt dim, PetscInt npoints, const PetscReal points[], PetscInt jet, PetscInt Njet, PetscReal pScalar[]) {
10920cf1dd8SToby Isaac   PetscFunctionBegin;
1109566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(pScalar, (1 + dim) * Njet * npoints));
111f1436e55SToby Isaac   for (PetscInt b = 0; b < 1 + dim; b++) {
112f1436e55SToby Isaac     for (PetscInt j = 0; j < PetscMin(1 + dim, Njet); j++) {
113f1436e55SToby Isaac       if (j == 0) {
114f1436e55SToby Isaac         if (b == 0) {
1159371c9d4SSatish Balay           for (PetscInt pt = 0; pt < npoints; pt++) { pScalar[b * Njet * npoints + j * npoints + pt] = 1.; }
11620cf1dd8SToby Isaac         } else {
1179371c9d4SSatish Balay           for (PetscInt pt = 0; pt < npoints; pt++) { pScalar[b * Njet * npoints + j * npoints + pt] = points[pt * dim + (b - 1)]; }
118f1436e55SToby Isaac         }
119f1436e55SToby Isaac       } else if (j == b) {
1209371c9d4SSatish Balay         for (PetscInt pt = 0; pt < npoints; pt++) { pScalar[b * Njet * npoints + j * npoints + pt] = 1.; }
121f1436e55SToby Isaac       }
12220cf1dd8SToby Isaac     }
12320cf1dd8SToby Isaac   }
12420cf1dd8SToby Isaac   PetscFunctionReturn(0);
12520cf1dd8SToby Isaac }
12620cf1dd8SToby Isaac 
1279371c9d4SSatish Balay static PetscErrorCode PetscSpaceEvaluate_Polynomial(PetscSpace sp, PetscInt npoints, const PetscReal points[], PetscReal B[], PetscReal D[], PetscReal H[]) {
12820cf1dd8SToby Isaac   PetscSpace_Poly *poly = (PetscSpace_Poly *)sp->data;
12920cf1dd8SToby Isaac   DM               dm   = sp->dm;
13020cf1dd8SToby Isaac   PetscInt         dim  = sp->Nv;
131f1436e55SToby Isaac   PetscInt         Nb, jet, Njet;
132f1436e55SToby Isaac   PetscReal       *pScalar;
13320cf1dd8SToby Isaac 
13420cf1dd8SToby Isaac   PetscFunctionBegin;
135f1436e55SToby Isaac   if (!poly->setupCalled) {
1369566063dSJacob Faibussowitsch     PetscCall(PetscSpaceSetUp(sp));
1379566063dSJacob Faibussowitsch     PetscCall(PetscSpaceEvaluate(sp, npoints, points, B, D, H));
138f1436e55SToby Isaac     PetscFunctionReturn(0);
13920cf1dd8SToby Isaac   }
1401dca8a05SBarry Smith   PetscCheck(!poly->tensor && sp->Nc == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "tensor and multicomponent spaces should have been converted");
1419566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(dim + sp->degree, dim, &Nb));
142f1436e55SToby Isaac   if (H) {
143f1436e55SToby Isaac     jet = 2;
144f1436e55SToby Isaac   } else if (D) {
145f1436e55SToby Isaac     jet = 1;
146f1436e55SToby Isaac   } else {
147f1436e55SToby Isaac     jet = 0;
14820cf1dd8SToby Isaac   }
1499566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(dim + jet, dim, &Njet));
1509566063dSJacob Faibussowitsch   PetscCall(DMGetWorkArray(dm, Nb * Njet * npoints, MPIU_REAL, &pScalar));
151f1436e55SToby Isaac   // Why are we handling the case degree == 1 specially?  Because we don't want numerical noise when we evaluate hat
152f1436e55SToby Isaac   // functions at the vertices of a simplex, which happens when we invert the Vandermonde matrix of the PKD basis.
153f1436e55SToby Isaac   // We don't make any promise about which basis is used.
154f1436e55SToby Isaac   if (sp->degree == 1) {
1559566063dSJacob Faibussowitsch     PetscCall(CoordinateBasis(dim, npoints, points, jet, Njet, pScalar));
156f1436e55SToby Isaac   } else {
1579566063dSJacob Faibussowitsch     PetscCall(PetscDTPKDEvalJet(dim, npoints, points, sp->degree, jet, pScalar));
15820cf1dd8SToby Isaac   }
15920cf1dd8SToby Isaac   if (B) {
160f1436e55SToby Isaac     PetscInt p_strl = Nb;
161f1436e55SToby Isaac     PetscInt b_strl = 1;
1623596293dSMatthew G. Knepley 
163f1436e55SToby Isaac     PetscInt b_strr = Njet * npoints;
164f1436e55SToby Isaac     PetscInt p_strr = 1;
165f1436e55SToby Isaac 
1669566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(B, npoints * Nb));
167f1436e55SToby Isaac     for (PetscInt b = 0; b < Nb; b++) {
1689371c9d4SSatish Balay       for (PetscInt p = 0; p < npoints; p++) { B[p * p_strl + b * b_strl] = pScalar[b * b_strr + p * p_strr]; }
1693596293dSMatthew G. Knepley     }
17020cf1dd8SToby Isaac   }
17120cf1dd8SToby Isaac   if (D) {
172f1436e55SToby Isaac     PetscInt p_strl = dim * Nb;
173f1436e55SToby Isaac     PetscInt b_strl = dim;
174f1436e55SToby Isaac     PetscInt d_strl = 1;
175f1436e55SToby Isaac 
176f1436e55SToby Isaac     PetscInt b_strr = Njet * npoints;
177f1436e55SToby Isaac     PetscInt d_strr = npoints;
178f1436e55SToby Isaac     PetscInt p_strr = 1;
179f1436e55SToby Isaac 
1809566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(D, npoints * Nb * dim));
181f1436e55SToby Isaac     for (PetscInt d = 0; d < dim; d++) {
182f1436e55SToby Isaac       for (PetscInt b = 0; b < Nb; b++) {
1839371c9d4SSatish Balay         for (PetscInt p = 0; p < npoints; p++) { D[p * p_strl + b * b_strl + d * d_strl] = pScalar[b * b_strr + (1 + d) * d_strr + p * p_strr]; }
18420cf1dd8SToby Isaac       }
18520cf1dd8SToby Isaac     }
18620cf1dd8SToby Isaac   }
18720cf1dd8SToby Isaac   if (H) {
188f1436e55SToby Isaac     PetscInt p_strl  = dim * dim * Nb;
189f1436e55SToby Isaac     PetscInt b_strl  = dim * dim;
190f1436e55SToby Isaac     PetscInt d1_strl = dim;
191f1436e55SToby Isaac     PetscInt d2_strl = 1;
192f1436e55SToby Isaac 
193f1436e55SToby Isaac     PetscInt b_strr = Njet * npoints;
194f1436e55SToby Isaac     PetscInt j_strr = npoints;
195f1436e55SToby Isaac     PetscInt p_strr = 1;
196f1436e55SToby Isaac 
197f1436e55SToby Isaac     PetscInt *derivs;
1989566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(dim, &derivs));
1999566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(H, npoints * Nb * dim * dim));
200f1436e55SToby Isaac     for (PetscInt d1 = 0; d1 < dim; d1++) {
201f1436e55SToby Isaac       for (PetscInt d2 = 0; d2 < dim; d2++) {
202f1436e55SToby Isaac         PetscInt j;
203f1436e55SToby Isaac         derivs[d1]++;
204f1436e55SToby Isaac         derivs[d2]++;
2059566063dSJacob Faibussowitsch         PetscCall(PetscDTGradedOrderToIndex(dim, derivs, &j));
206f1436e55SToby Isaac         derivs[d1]--;
207f1436e55SToby Isaac         derivs[d2]--;
208f1436e55SToby Isaac         for (PetscInt b = 0; b < Nb; b++) {
2099371c9d4SSatish Balay           for (PetscInt p = 0; p < npoints; p++) { H[p * p_strl + b * b_strl + d1 * d1_strl + d2 * d2_strl] = pScalar[b * b_strr + j * j_strr + p * p_strr]; }
21020cf1dd8SToby Isaac         }
21120cf1dd8SToby Isaac       }
21220cf1dd8SToby Isaac     }
2139566063dSJacob Faibussowitsch     PetscCall(PetscFree(derivs));
21420cf1dd8SToby Isaac   }
2159566063dSJacob Faibussowitsch   PetscCall(DMRestoreWorkArray(dm, Nb * Njet * npoints, MPIU_REAL, &pScalar));
21620cf1dd8SToby Isaac   PetscFunctionReturn(0);
21720cf1dd8SToby Isaac }
21820cf1dd8SToby Isaac 
21920cf1dd8SToby Isaac /*@
22020cf1dd8SToby Isaac   PetscSpacePolynomialSetTensor - Set whether a function space is a space of tensor polynomials (the space is spanned
221f1436e55SToby Isaac   by polynomials whose degree in each variable is bounded by the given order), as opposed to polynomials (the space is
22220cf1dd8SToby Isaac   spanned by polynomials whose total degree---summing over all variables---is bounded by the given order).
22320cf1dd8SToby Isaac 
22420cf1dd8SToby Isaac   Input Parameters:
22520cf1dd8SToby Isaac + sp     - the function space object
22620cf1dd8SToby Isaac - tensor - PETSC_TRUE for a tensor polynomial space, PETSC_FALSE for a polynomial space
22720cf1dd8SToby Isaac 
2284ab77754SMatthew G. Knepley   Options Database:
2294ab77754SMatthew G. Knepley . -petscspace_poly_tensor <bool> - Whether to use tensor product polynomials in higher dimension
2304ab77754SMatthew G. Knepley 
23129b5c115SMatthew G. Knepley   Level: intermediate
23220cf1dd8SToby Isaac 
233db781477SPatrick Sanan .seealso: `PetscSpacePolynomialGetTensor()`, `PetscSpaceSetDegree()`, `PetscSpaceSetNumVariables()`
23420cf1dd8SToby Isaac @*/
2359371c9d4SSatish Balay PetscErrorCode PetscSpacePolynomialSetTensor(PetscSpace sp, PetscBool tensor) {
23620cf1dd8SToby Isaac   PetscFunctionBegin;
23720cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
238cac4c232SBarry Smith   PetscTryMethod(sp, "PetscSpacePolynomialSetTensor_C", (PetscSpace, PetscBool), (sp, tensor));
23920cf1dd8SToby Isaac   PetscFunctionReturn(0);
24020cf1dd8SToby Isaac }
24120cf1dd8SToby Isaac 
24220cf1dd8SToby Isaac /*@
24320cf1dd8SToby Isaac   PetscSpacePolynomialGetTensor - Get whether a function space is a space of tensor polynomials (the space is spanned
24420cf1dd8SToby Isaac   by polynomials whose degree in each variabl is bounded by the given order), as opposed to polynomials (the space is
24520cf1dd8SToby Isaac   spanned by polynomials whose total degree---summing over all variables---is bounded by the given order).
24620cf1dd8SToby Isaac 
24720cf1dd8SToby Isaac   Input Parameters:
24820cf1dd8SToby Isaac . sp     - the function space object
24920cf1dd8SToby Isaac 
25020cf1dd8SToby Isaac   Output Parameters:
25120cf1dd8SToby Isaac . tensor - PETSC_TRUE for a tensor polynomial space, PETSC_FALSE for a polynomial space
25220cf1dd8SToby Isaac 
25329b5c115SMatthew G. Knepley   Level: intermediate
25420cf1dd8SToby Isaac 
255db781477SPatrick Sanan .seealso: `PetscSpacePolynomialSetTensor()`, `PetscSpaceSetDegree()`, `PetscSpaceSetNumVariables()`
25620cf1dd8SToby Isaac @*/
2579371c9d4SSatish Balay PetscErrorCode PetscSpacePolynomialGetTensor(PetscSpace sp, PetscBool *tensor) {
25820cf1dd8SToby Isaac   PetscFunctionBegin;
25920cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
260dadcf809SJacob Faibussowitsch   PetscValidBoolPointer(tensor, 2);
261cac4c232SBarry Smith   PetscTryMethod(sp, "PetscSpacePolynomialGetTensor_C", (PetscSpace, PetscBool *), (sp, tensor));
26220cf1dd8SToby Isaac   PetscFunctionReturn(0);
26320cf1dd8SToby Isaac }
26420cf1dd8SToby Isaac 
2659371c9d4SSatish Balay static PetscErrorCode PetscSpacePolynomialSetTensor_Polynomial(PetscSpace sp, PetscBool tensor) {
26620cf1dd8SToby Isaac   PetscSpace_Poly *poly = (PetscSpace_Poly *)sp->data;
26720cf1dd8SToby Isaac 
26820cf1dd8SToby Isaac   PetscFunctionBegin;
26920cf1dd8SToby Isaac   poly->tensor = tensor;
27020cf1dd8SToby Isaac   PetscFunctionReturn(0);
27120cf1dd8SToby Isaac }
27220cf1dd8SToby Isaac 
2739371c9d4SSatish Balay static PetscErrorCode PetscSpacePolynomialGetTensor_Polynomial(PetscSpace sp, PetscBool *tensor) {
27420cf1dd8SToby Isaac   PetscSpace_Poly *poly = (PetscSpace_Poly *)sp->data;
27520cf1dd8SToby Isaac 
27620cf1dd8SToby Isaac   PetscFunctionBegin;
27720cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
278dadcf809SJacob Faibussowitsch   PetscValidBoolPointer(tensor, 2);
27920cf1dd8SToby Isaac   *tensor = poly->tensor;
28020cf1dd8SToby Isaac   PetscFunctionReturn(0);
28120cf1dd8SToby Isaac }
28220cf1dd8SToby Isaac 
2839371c9d4SSatish Balay static PetscErrorCode PetscSpaceGetHeightSubspace_Polynomial(PetscSpace sp, PetscInt height, PetscSpace *subsp) {
28420cf1dd8SToby Isaac   PetscSpace_Poly *poly = (PetscSpace_Poly *)sp->data;
28520cf1dd8SToby Isaac   PetscInt         Nc, dim, order;
28620cf1dd8SToby Isaac   PetscBool        tensor;
28720cf1dd8SToby Isaac 
28820cf1dd8SToby Isaac   PetscFunctionBegin;
2899566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetNumComponents(sp, &Nc));
2909566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetNumVariables(sp, &dim));
2919566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetDegree(sp, &order, NULL));
2929566063dSJacob Faibussowitsch   PetscCall(PetscSpacePolynomialGetTensor(sp, &tensor));
2931dca8a05SBarry Smith   PetscCheck(height <= dim && height >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Asked for space at height %" PetscInt_FMT " for dimension %" PetscInt_FMT " space", height, dim);
2949566063dSJacob Faibussowitsch   if (!poly->subspaces) PetscCall(PetscCalloc1(dim, &poly->subspaces));
29520cf1dd8SToby Isaac   if (height <= dim) {
29620cf1dd8SToby Isaac     if (!poly->subspaces[height - 1]) {
29720cf1dd8SToby Isaac       PetscSpace  sub;
2983f6b16c7SMatthew G. Knepley       const char *name;
29920cf1dd8SToby Isaac 
3009566063dSJacob Faibussowitsch       PetscCall(PetscSpaceCreate(PetscObjectComm((PetscObject)sp), &sub));
3019566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetName((PetscObject)sp, &name));
3029566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)sub, name));
3039566063dSJacob Faibussowitsch       PetscCall(PetscSpaceSetType(sub, PETSCSPACEPOLYNOMIAL));
3049566063dSJacob Faibussowitsch       PetscCall(PetscSpaceSetNumComponents(sub, Nc));
3059566063dSJacob Faibussowitsch       PetscCall(PetscSpaceSetDegree(sub, order, PETSC_DETERMINE));
3069566063dSJacob Faibussowitsch       PetscCall(PetscSpaceSetNumVariables(sub, dim - height));
3079566063dSJacob Faibussowitsch       PetscCall(PetscSpacePolynomialSetTensor(sub, tensor));
3089566063dSJacob Faibussowitsch       PetscCall(PetscSpaceSetUp(sub));
30920cf1dd8SToby Isaac       poly->subspaces[height - 1] = sub;
31020cf1dd8SToby Isaac     }
31120cf1dd8SToby Isaac     *subsp = poly->subspaces[height - 1];
31220cf1dd8SToby Isaac   } else {
31320cf1dd8SToby Isaac     *subsp = NULL;
31420cf1dd8SToby Isaac   }
31520cf1dd8SToby Isaac   PetscFunctionReturn(0);
31620cf1dd8SToby Isaac }
31720cf1dd8SToby Isaac 
3189371c9d4SSatish Balay static PetscErrorCode PetscSpaceInitialize_Polynomial(PetscSpace sp) {
31920cf1dd8SToby Isaac   PetscFunctionBegin;
32020cf1dd8SToby Isaac   sp->ops->setfromoptions    = PetscSpaceSetFromOptions_Polynomial;
32120cf1dd8SToby Isaac   sp->ops->setup             = PetscSpaceSetUp_Polynomial;
32220cf1dd8SToby Isaac   sp->ops->view              = PetscSpaceView_Polynomial;
32320cf1dd8SToby Isaac   sp->ops->destroy           = PetscSpaceDestroy_Polynomial;
32420cf1dd8SToby Isaac   sp->ops->getdimension      = PetscSpaceGetDimension_Polynomial;
32520cf1dd8SToby Isaac   sp->ops->evaluate          = PetscSpaceEvaluate_Polynomial;
32620cf1dd8SToby Isaac   sp->ops->getheightsubspace = PetscSpaceGetHeightSubspace_Polynomial;
3279566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscSpacePolynomialGetTensor_C", PetscSpacePolynomialGetTensor_Polynomial));
3289566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscSpacePolynomialSetTensor_C", PetscSpacePolynomialSetTensor_Polynomial));
32920cf1dd8SToby Isaac   PetscFunctionReturn(0);
33020cf1dd8SToby Isaac }
33120cf1dd8SToby Isaac 
33220cf1dd8SToby Isaac /*MC
33320cf1dd8SToby Isaac   PETSCSPACEPOLYNOMIAL = "poly" - A PetscSpace object that encapsulates a polynomial space, e.g. P1 is the space of
33420cf1dd8SToby Isaac   linear polynomials. The space is replicated for each component.
33520cf1dd8SToby Isaac 
33620cf1dd8SToby Isaac   Level: intermediate
33720cf1dd8SToby Isaac 
338db781477SPatrick Sanan .seealso: `PetscSpaceType`, `PetscSpaceCreate()`, `PetscSpaceSetType()`
33920cf1dd8SToby Isaac M*/
34020cf1dd8SToby Isaac 
3419371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PetscSpaceCreate_Polynomial(PetscSpace sp) {
34220cf1dd8SToby Isaac   PetscSpace_Poly *poly;
34320cf1dd8SToby Isaac 
34420cf1dd8SToby Isaac   PetscFunctionBegin;
34520cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
3469566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(sp, &poly));
34720cf1dd8SToby Isaac   sp->data = poly;
34820cf1dd8SToby Isaac 
34920cf1dd8SToby Isaac   poly->tensor    = PETSC_FALSE;
35020cf1dd8SToby Isaac   poly->subspaces = NULL;
35120cf1dd8SToby Isaac 
3529566063dSJacob Faibussowitsch   PetscCall(PetscSpaceInitialize_Polynomial(sp));
35320cf1dd8SToby Isaac   PetscFunctionReturn(0);
35420cf1dd8SToby Isaac }
355