xref: /petsc/src/dm/dt/space/impls/poly/spacepoly.c (revision a4e35b1925eceef64945ea472b84f2bf06a67b5e)
120cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
220cf1dd8SToby Isaac 
3d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSpaceSetFromOptions_Polynomial(PetscSpace sp, PetscOptionItems *PetscOptionsObject)
4d71ae5a4SJacob Faibussowitsch {
520cf1dd8SToby Isaac   PetscSpace_Poly *poly = (PetscSpace_Poly *)sp->data;
620cf1dd8SToby Isaac 
720cf1dd8SToby Isaac   PetscFunctionBegin;
8d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "PetscSpace polynomial options");
99566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-petscspace_poly_tensor", "Use the tensor product polynomials", "PetscSpacePolynomialSetTensor", poly->tensor, &poly->tensor, NULL));
10d0609cedSBarry Smith   PetscOptionsHeadEnd();
113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1220cf1dd8SToby Isaac }
1320cf1dd8SToby Isaac 
14d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSpacePolynomialView_Ascii(PetscSpace sp, PetscViewer v)
15d71ae5a4SJacob Faibussowitsch {
1620cf1dd8SToby Isaac   PetscSpace_Poly *poly = (PetscSpace_Poly *)sp->data;
1720cf1dd8SToby Isaac 
1820cf1dd8SToby Isaac   PetscFunctionBegin;
1963a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(v, "%s space of degree %" PetscInt_FMT "\n", poly->tensor ? "Tensor polynomial" : "Polynomial", sp->degree));
203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2120cf1dd8SToby Isaac }
2220cf1dd8SToby Isaac 
23d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSpaceView_Polynomial(PetscSpace sp, PetscViewer viewer)
24d71ae5a4SJacob Faibussowitsch {
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));
323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3320cf1dd8SToby Isaac }
3420cf1dd8SToby Isaac 
35d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSpaceDestroy_Polynomial(PetscSpace sp)
36d71ae5a4SJacob Faibussowitsch {
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 
4548a46eb9SPierre Jolivet     for (d = 0; d < sp->Nv; ++d) PetscCall(PetscSpaceDestroy(&poly->subspaces[d]));
4620cf1dd8SToby Isaac   }
479566063dSJacob Faibussowitsch   PetscCall(PetscFree(poly->subspaces));
489566063dSJacob Faibussowitsch   PetscCall(PetscFree(poly));
493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5020cf1dd8SToby Isaac }
5120cf1dd8SToby Isaac 
52d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSpaceSetUp_Polynomial(PetscSpace sp)
53d71ae5a4SJacob Faibussowitsch {
54f1436e55SToby Isaac   PetscSpace_Poly *poly = (PetscSpace_Poly *)sp->data;
55f1436e55SToby Isaac 
56f1436e55SToby Isaac   PetscFunctionBegin;
573ba16761SJacob Faibussowitsch   if (poly->setupCalled) PetscFunctionReturn(PETSC_SUCCESS);
58ad540459SPierre Jolivet   if (sp->Nv <= 1) poly->tensor = PETSC_FALSE;
59f1436e55SToby Isaac   if (sp->Nc != 1) {
60f1436e55SToby Isaac     PetscInt    Nc     = sp->Nc;
61f1436e55SToby Isaac     PetscBool   tensor = poly->tensor;
62f1436e55SToby Isaac     PetscInt    Nv     = sp->Nv;
63f1436e55SToby Isaac     PetscInt    degree = sp->degree;
64417c287bSToby Isaac     const char *prefix;
65417c287bSToby Isaac     const char *name;
66417c287bSToby Isaac     char        subname[PETSC_MAX_PATH_LEN];
67f1436e55SToby Isaac     PetscSpace  subsp;
68f1436e55SToby Isaac 
699566063dSJacob Faibussowitsch     PetscCall(PetscSpaceSetType(sp, PETSCSPACESUM));
709566063dSJacob Faibussowitsch     PetscCall(PetscSpaceSumSetNumSubspaces(sp, Nc));
719566063dSJacob Faibussowitsch     PetscCall(PetscSpaceCreate(PetscObjectComm((PetscObject)sp), &subsp));
729566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetOptionsPrefix((PetscObject)sp, &prefix));
739566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)subsp, prefix));
749566063dSJacob Faibussowitsch     PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)subsp, "sumcomp_"));
75417c287bSToby Isaac     if (((PetscObject)sp)->name) {
769566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetName((PetscObject)sp, &name));
779566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(subname, PETSC_MAX_PATH_LEN - 1, "%s sum component", name));
789566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)subsp, subname));
791baa6e33SBarry Smith     } else PetscCall(PetscObjectSetName((PetscObject)subsp, "sum component"));
809566063dSJacob Faibussowitsch     PetscCall(PetscSpaceSetType(subsp, PETSCSPACEPOLYNOMIAL));
819566063dSJacob Faibussowitsch     PetscCall(PetscSpaceSetDegree(subsp, degree, PETSC_DETERMINE));
829566063dSJacob Faibussowitsch     PetscCall(PetscSpaceSetNumComponents(subsp, 1));
839566063dSJacob Faibussowitsch     PetscCall(PetscSpaceSetNumVariables(subsp, Nv));
849566063dSJacob Faibussowitsch     PetscCall(PetscSpacePolynomialSetTensor(subsp, tensor));
859566063dSJacob Faibussowitsch     PetscCall(PetscSpaceSetUp(subsp));
8648a46eb9SPierre Jolivet     for (PetscInt i = 0; i < Nc; i++) PetscCall(PetscSpaceSumSetSubspace(sp, i, subsp));
879566063dSJacob Faibussowitsch     PetscCall(PetscSpaceDestroy(&subsp));
889566063dSJacob Faibussowitsch     PetscCall(PetscSpaceSetUp(sp));
893ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
90f1436e55SToby Isaac   }
91f1436e55SToby Isaac   if (poly->tensor) {
92f1436e55SToby Isaac     sp->maxDegree = PETSC_DETERMINE;
939566063dSJacob Faibussowitsch     PetscCall(PetscSpaceSetType(sp, PETSCSPACETENSOR));
949566063dSJacob Faibussowitsch     PetscCall(PetscSpaceSetUp(sp));
953ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
96f1436e55SToby Isaac   }
9763a3b9bcSJacob Faibussowitsch   PetscCheck(sp->degree >= 0, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_OUTOFRANGE, "Negative degree %" PetscInt_FMT " invalid", sp->degree);
98f1436e55SToby Isaac   sp->maxDegree     = sp->degree;
99f1436e55SToby Isaac   poly->setupCalled = PETSC_TRUE;
1003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
101f1436e55SToby Isaac }
102f1436e55SToby Isaac 
103d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSpaceGetDimension_Polynomial(PetscSpace sp, PetscInt *dim)
104d71ae5a4SJacob Faibussowitsch {
10520cf1dd8SToby Isaac   PetscInt deg = sp->degree;
106f1436e55SToby Isaac   PetscInt n   = sp->Nv;
10720cf1dd8SToby Isaac 
10820cf1dd8SToby Isaac   PetscFunctionBegin;
1099566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(n + deg, n, dim));
110f1436e55SToby Isaac   *dim *= sp->Nc;
1113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11220cf1dd8SToby Isaac }
11320cf1dd8SToby Isaac 
114d71ae5a4SJacob Faibussowitsch static PetscErrorCode CoordinateBasis(PetscInt dim, PetscInt npoints, const PetscReal points[], PetscInt jet, PetscInt Njet, PetscReal pScalar[])
115d71ae5a4SJacob Faibussowitsch {
11620cf1dd8SToby Isaac   PetscFunctionBegin;
1179566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(pScalar, (1 + dim) * Njet * npoints));
118f1436e55SToby Isaac   for (PetscInt b = 0; b < 1 + dim; b++) {
119f1436e55SToby Isaac     for (PetscInt j = 0; j < PetscMin(1 + dim, Njet); j++) {
120f1436e55SToby Isaac       if (j == 0) {
121f1436e55SToby Isaac         if (b == 0) {
122ad540459SPierre Jolivet           for (PetscInt pt = 0; pt < npoints; pt++) pScalar[b * Njet * npoints + j * npoints + pt] = 1.;
12320cf1dd8SToby Isaac         } else {
124ad540459SPierre Jolivet           for (PetscInt pt = 0; pt < npoints; pt++) pScalar[b * Njet * npoints + j * npoints + pt] = points[pt * dim + (b - 1)];
125f1436e55SToby Isaac         }
126f1436e55SToby Isaac       } else if (j == b) {
127ad540459SPierre Jolivet         for (PetscInt pt = 0; pt < npoints; pt++) pScalar[b * Njet * npoints + j * npoints + pt] = 1.;
128f1436e55SToby Isaac       }
12920cf1dd8SToby Isaac     }
13020cf1dd8SToby Isaac   }
1313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13220cf1dd8SToby Isaac }
13320cf1dd8SToby Isaac 
134d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSpaceEvaluate_Polynomial(PetscSpace sp, PetscInt npoints, const PetscReal points[], PetscReal B[], PetscReal D[], PetscReal H[])
135d71ae5a4SJacob Faibussowitsch {
13620cf1dd8SToby Isaac   PetscSpace_Poly *poly = (PetscSpace_Poly *)sp->data;
13720cf1dd8SToby Isaac   DM               dm   = sp->dm;
13820cf1dd8SToby Isaac   PetscInt         dim  = sp->Nv;
139f1436e55SToby Isaac   PetscInt         Nb, jet, Njet;
140f1436e55SToby Isaac   PetscReal       *pScalar;
14120cf1dd8SToby Isaac 
14220cf1dd8SToby Isaac   PetscFunctionBegin;
143f1436e55SToby Isaac   if (!poly->setupCalled) {
1449566063dSJacob Faibussowitsch     PetscCall(PetscSpaceSetUp(sp));
1459566063dSJacob Faibussowitsch     PetscCall(PetscSpaceEvaluate(sp, npoints, points, B, D, H));
1463ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
14720cf1dd8SToby Isaac   }
1481dca8a05SBarry Smith   PetscCheck(!poly->tensor && sp->Nc == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "tensor and multicomponent spaces should have been converted");
1499566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(dim + sp->degree, dim, &Nb));
150f1436e55SToby Isaac   if (H) {
151f1436e55SToby Isaac     jet = 2;
152f1436e55SToby Isaac   } else if (D) {
153f1436e55SToby Isaac     jet = 1;
154f1436e55SToby Isaac   } else {
155f1436e55SToby Isaac     jet = 0;
15620cf1dd8SToby Isaac   }
1579566063dSJacob Faibussowitsch   PetscCall(PetscDTBinomialInt(dim + jet, dim, &Njet));
1589566063dSJacob Faibussowitsch   PetscCall(DMGetWorkArray(dm, Nb * Njet * npoints, MPIU_REAL, &pScalar));
159f1436e55SToby Isaac   // Why are we handling the case degree == 1 specially?  Because we don't want numerical noise when we evaluate hat
160f1436e55SToby Isaac   // functions at the vertices of a simplex, which happens when we invert the Vandermonde matrix of the PKD basis.
161f1436e55SToby Isaac   // We don't make any promise about which basis is used.
162f1436e55SToby Isaac   if (sp->degree == 1) {
1639566063dSJacob Faibussowitsch     PetscCall(CoordinateBasis(dim, npoints, points, jet, Njet, pScalar));
164f1436e55SToby Isaac   } else {
1659566063dSJacob Faibussowitsch     PetscCall(PetscDTPKDEvalJet(dim, npoints, points, sp->degree, jet, pScalar));
16620cf1dd8SToby Isaac   }
16720cf1dd8SToby Isaac   if (B) {
168f1436e55SToby Isaac     PetscInt p_strl = Nb;
169f1436e55SToby Isaac     PetscInt b_strl = 1;
1703596293dSMatthew G. Knepley 
171f1436e55SToby Isaac     PetscInt b_strr = Njet * npoints;
172f1436e55SToby Isaac     PetscInt p_strr = 1;
173f1436e55SToby Isaac 
1749566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(B, npoints * Nb));
175f1436e55SToby Isaac     for (PetscInt b = 0; b < Nb; b++) {
176ad540459SPierre Jolivet       for (PetscInt p = 0; p < npoints; p++) B[p * p_strl + b * b_strl] = pScalar[b * b_strr + p * p_strr];
1773596293dSMatthew G. Knepley     }
17820cf1dd8SToby Isaac   }
17920cf1dd8SToby Isaac   if (D) {
180f1436e55SToby Isaac     PetscInt p_strl = dim * Nb;
181f1436e55SToby Isaac     PetscInt b_strl = dim;
182f1436e55SToby Isaac     PetscInt d_strl = 1;
183f1436e55SToby Isaac 
184f1436e55SToby Isaac     PetscInt b_strr = Njet * npoints;
185f1436e55SToby Isaac     PetscInt d_strr = npoints;
186f1436e55SToby Isaac     PetscInt p_strr = 1;
187f1436e55SToby Isaac 
1889566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(D, npoints * Nb * dim));
189f1436e55SToby Isaac     for (PetscInt d = 0; d < dim; d++) {
190f1436e55SToby Isaac       for (PetscInt b = 0; b < Nb; b++) {
191ad540459SPierre Jolivet         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];
19220cf1dd8SToby Isaac       }
19320cf1dd8SToby Isaac     }
19420cf1dd8SToby Isaac   }
19520cf1dd8SToby Isaac   if (H) {
196f1436e55SToby Isaac     PetscInt p_strl  = dim * dim * Nb;
197f1436e55SToby Isaac     PetscInt b_strl  = dim * dim;
198f1436e55SToby Isaac     PetscInt d1_strl = dim;
199f1436e55SToby Isaac     PetscInt d2_strl = 1;
200f1436e55SToby Isaac 
201f1436e55SToby Isaac     PetscInt b_strr = Njet * npoints;
202f1436e55SToby Isaac     PetscInt j_strr = npoints;
203f1436e55SToby Isaac     PetscInt p_strr = 1;
204f1436e55SToby Isaac 
205f1436e55SToby Isaac     PetscInt *derivs;
2069566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(dim, &derivs));
2079566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(H, npoints * Nb * dim * dim));
208f1436e55SToby Isaac     for (PetscInt d1 = 0; d1 < dim; d1++) {
209f1436e55SToby Isaac       for (PetscInt d2 = 0; d2 < dim; d2++) {
210f1436e55SToby Isaac         PetscInt j;
211f1436e55SToby Isaac         derivs[d1]++;
212f1436e55SToby Isaac         derivs[d2]++;
2139566063dSJacob Faibussowitsch         PetscCall(PetscDTGradedOrderToIndex(dim, derivs, &j));
214f1436e55SToby Isaac         derivs[d1]--;
215f1436e55SToby Isaac         derivs[d2]--;
216f1436e55SToby Isaac         for (PetscInt b = 0; b < Nb; b++) {
217ad540459SPierre Jolivet           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];
21820cf1dd8SToby Isaac         }
21920cf1dd8SToby Isaac       }
22020cf1dd8SToby Isaac     }
2219566063dSJacob Faibussowitsch     PetscCall(PetscFree(derivs));
22220cf1dd8SToby Isaac   }
2239566063dSJacob Faibussowitsch   PetscCall(DMRestoreWorkArray(dm, Nb * Njet * npoints, MPIU_REAL, &pScalar));
2243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
22520cf1dd8SToby Isaac }
22620cf1dd8SToby Isaac 
22720cf1dd8SToby Isaac /*@
228*a4e35b19SJacob Faibussowitsch   PetscSpacePolynomialSetTensor - Set whether a function space is a space of tensor polynomials.
22920cf1dd8SToby Isaac 
23020cf1dd8SToby Isaac   Input Parameters:
23120cf1dd8SToby Isaac + sp     - the function space object
232dce8aebaSBarry Smith - tensor - `PETSC_TRUE` for a tensor polynomial space, `PETSC_FALSE` for a polynomial space
23320cf1dd8SToby Isaac 
234dce8aebaSBarry Smith   Options Database Key:
2354ab77754SMatthew G. Knepley . -petscspace_poly_tensor <bool> - Whether to use tensor product polynomials in higher dimension
2364ab77754SMatthew G. Knepley 
23729b5c115SMatthew G. Knepley   Level: intermediate
23820cf1dd8SToby Isaac 
239*a4e35b19SJacob Faibussowitsch   Notes:
240*a4e35b19SJacob Faibussowitsch   It is a tensor space if it is spanned by polynomials whose degree in each variable is
241*a4e35b19SJacob Faibussowitsch   bounded by the given order, as opposed to the space spanned by polynomials
242*a4e35b19SJacob Faibussowitsch   whose total degree---summing over all variables---is bounded by the given order.
243*a4e35b19SJacob Faibussowitsch 
244dce8aebaSBarry Smith .seealso: `PetscSpace`, `PetscSpacePolynomialGetTensor()`, `PetscSpaceSetDegree()`, `PetscSpaceSetNumVariables()`
24520cf1dd8SToby Isaac @*/
246d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSpacePolynomialSetTensor(PetscSpace sp, PetscBool tensor)
247d71ae5a4SJacob Faibussowitsch {
24820cf1dd8SToby Isaac   PetscFunctionBegin;
24920cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
250cac4c232SBarry Smith   PetscTryMethod(sp, "PetscSpacePolynomialSetTensor_C", (PetscSpace, PetscBool), (sp, tensor));
2513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
25220cf1dd8SToby Isaac }
25320cf1dd8SToby Isaac 
25420cf1dd8SToby Isaac /*@
255*a4e35b19SJacob Faibussowitsch   PetscSpacePolynomialGetTensor - Get whether a function space is a space of tensor
256*a4e35b19SJacob Faibussowitsch   polynomials.
25720cf1dd8SToby Isaac 
2582fe279fdSBarry Smith   Input Parameter:
25920cf1dd8SToby Isaac . sp - the function space object
26020cf1dd8SToby Isaac 
2612fe279fdSBarry Smith   Output Parameter:
262dce8aebaSBarry Smith . tensor - `PETSC_TRUE` for a tensor polynomial space, `PETSC_FALSE` for a polynomial space
26320cf1dd8SToby Isaac 
26429b5c115SMatthew G. Knepley   Level: intermediate
26520cf1dd8SToby Isaac 
266*a4e35b19SJacob Faibussowitsch   Notes:
267*a4e35b19SJacob Faibussowitsch   The space is a tensor space if it is spanned by polynomials whose degree in each variable is
268*a4e35b19SJacob Faibussowitsch   bounded by the given order, as opposed to the space spanned by polynomials
269*a4e35b19SJacob Faibussowitsch   whose total degree---summing over all variables---is bounded by the given order.
270*a4e35b19SJacob Faibussowitsch 
271dce8aebaSBarry Smith .seealso: `PetscSpace`, `PetscSpacePolynomialSetTensor()`, `PetscSpaceSetDegree()`, `PetscSpaceSetNumVariables()`
27220cf1dd8SToby Isaac @*/
273d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSpacePolynomialGetTensor(PetscSpace sp, PetscBool *tensor)
274d71ae5a4SJacob Faibussowitsch {
27520cf1dd8SToby Isaac   PetscFunctionBegin;
27620cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
2774f572ea9SToby Isaac   PetscAssertPointer(tensor, 2);
278cac4c232SBarry Smith   PetscTryMethod(sp, "PetscSpacePolynomialGetTensor_C", (PetscSpace, PetscBool *), (sp, tensor));
2793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
28020cf1dd8SToby Isaac }
28120cf1dd8SToby Isaac 
282d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSpacePolynomialSetTensor_Polynomial(PetscSpace sp, PetscBool tensor)
283d71ae5a4SJacob Faibussowitsch {
28420cf1dd8SToby Isaac   PetscSpace_Poly *poly = (PetscSpace_Poly *)sp->data;
28520cf1dd8SToby Isaac 
28620cf1dd8SToby Isaac   PetscFunctionBegin;
28720cf1dd8SToby Isaac   poly->tensor = tensor;
2883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
28920cf1dd8SToby Isaac }
29020cf1dd8SToby Isaac 
291d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSpacePolynomialGetTensor_Polynomial(PetscSpace sp, PetscBool *tensor)
292d71ae5a4SJacob Faibussowitsch {
29320cf1dd8SToby Isaac   PetscSpace_Poly *poly = (PetscSpace_Poly *)sp->data;
29420cf1dd8SToby Isaac 
29520cf1dd8SToby Isaac   PetscFunctionBegin;
29620cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
2974f572ea9SToby Isaac   PetscAssertPointer(tensor, 2);
29820cf1dd8SToby Isaac   *tensor = poly->tensor;
2993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
30020cf1dd8SToby Isaac }
30120cf1dd8SToby Isaac 
302d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSpaceGetHeightSubspace_Polynomial(PetscSpace sp, PetscInt height, PetscSpace *subsp)
303d71ae5a4SJacob Faibussowitsch {
30420cf1dd8SToby Isaac   PetscSpace_Poly *poly = (PetscSpace_Poly *)sp->data;
30520cf1dd8SToby Isaac   PetscInt         Nc, dim, order;
30620cf1dd8SToby Isaac   PetscBool        tensor;
30720cf1dd8SToby Isaac 
30820cf1dd8SToby Isaac   PetscFunctionBegin;
3099566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetNumComponents(sp, &Nc));
3109566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetNumVariables(sp, &dim));
3119566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetDegree(sp, &order, NULL));
3129566063dSJacob Faibussowitsch   PetscCall(PetscSpacePolynomialGetTensor(sp, &tensor));
3131dca8a05SBarry 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);
3149566063dSJacob Faibussowitsch   if (!poly->subspaces) PetscCall(PetscCalloc1(dim, &poly->subspaces));
31520cf1dd8SToby Isaac   if (height <= dim) {
31620cf1dd8SToby Isaac     if (!poly->subspaces[height - 1]) {
31720cf1dd8SToby Isaac       PetscSpace  sub;
3183f6b16c7SMatthew G. Knepley       const char *name;
31920cf1dd8SToby Isaac 
3209566063dSJacob Faibussowitsch       PetscCall(PetscSpaceCreate(PetscObjectComm((PetscObject)sp), &sub));
3219566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetName((PetscObject)sp, &name));
3229566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)sub, name));
3239566063dSJacob Faibussowitsch       PetscCall(PetscSpaceSetType(sub, PETSCSPACEPOLYNOMIAL));
3249566063dSJacob Faibussowitsch       PetscCall(PetscSpaceSetNumComponents(sub, Nc));
3259566063dSJacob Faibussowitsch       PetscCall(PetscSpaceSetDegree(sub, order, PETSC_DETERMINE));
3269566063dSJacob Faibussowitsch       PetscCall(PetscSpaceSetNumVariables(sub, dim - height));
3279566063dSJacob Faibussowitsch       PetscCall(PetscSpacePolynomialSetTensor(sub, tensor));
3289566063dSJacob Faibussowitsch       PetscCall(PetscSpaceSetUp(sub));
32920cf1dd8SToby Isaac       poly->subspaces[height - 1] = sub;
33020cf1dd8SToby Isaac     }
33120cf1dd8SToby Isaac     *subsp = poly->subspaces[height - 1];
33220cf1dd8SToby Isaac   } else {
33320cf1dd8SToby Isaac     *subsp = NULL;
33420cf1dd8SToby Isaac   }
3353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
33620cf1dd8SToby Isaac }
33720cf1dd8SToby Isaac 
338d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSpaceInitialize_Polynomial(PetscSpace sp)
339d71ae5a4SJacob Faibussowitsch {
34020cf1dd8SToby Isaac   PetscFunctionBegin;
34120cf1dd8SToby Isaac   sp->ops->setfromoptions    = PetscSpaceSetFromOptions_Polynomial;
34220cf1dd8SToby Isaac   sp->ops->setup             = PetscSpaceSetUp_Polynomial;
34320cf1dd8SToby Isaac   sp->ops->view              = PetscSpaceView_Polynomial;
34420cf1dd8SToby Isaac   sp->ops->destroy           = PetscSpaceDestroy_Polynomial;
34520cf1dd8SToby Isaac   sp->ops->getdimension      = PetscSpaceGetDimension_Polynomial;
34620cf1dd8SToby Isaac   sp->ops->evaluate          = PetscSpaceEvaluate_Polynomial;
34720cf1dd8SToby Isaac   sp->ops->getheightsubspace = PetscSpaceGetHeightSubspace_Polynomial;
3489566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscSpacePolynomialGetTensor_C", PetscSpacePolynomialGetTensor_Polynomial));
3499566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscSpacePolynomialSetTensor_C", PetscSpacePolynomialSetTensor_Polynomial));
3503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
35120cf1dd8SToby Isaac }
35220cf1dd8SToby Isaac 
35320cf1dd8SToby Isaac /*MC
354dce8aebaSBarry Smith   PETSCSPACEPOLYNOMIAL = "poly" - A `PetscSpace` object that encapsulates a polynomial space, e.g. P1 is the space of
35520cf1dd8SToby Isaac   linear polynomials. The space is replicated for each component.
35620cf1dd8SToby Isaac 
35720cf1dd8SToby Isaac   Level: intermediate
35820cf1dd8SToby Isaac 
359dce8aebaSBarry Smith .seealso: `PetscSpace`, `PetscSpaceType`, `PetscSpaceCreate()`, `PetscSpaceSetType()`
36020cf1dd8SToby Isaac M*/
36120cf1dd8SToby Isaac 
362d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscSpaceCreate_Polynomial(PetscSpace sp)
363d71ae5a4SJacob Faibussowitsch {
36420cf1dd8SToby Isaac   PetscSpace_Poly *poly;
36520cf1dd8SToby Isaac 
36620cf1dd8SToby Isaac   PetscFunctionBegin;
36720cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
3684dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&poly));
36920cf1dd8SToby Isaac   sp->data = poly;
37020cf1dd8SToby Isaac 
37120cf1dd8SToby Isaac   poly->tensor    = PETSC_FALSE;
37220cf1dd8SToby Isaac   poly->subspaces = NULL;
37320cf1dd8SToby Isaac 
3749566063dSJacob Faibussowitsch   PetscCall(PetscSpaceInitialize_Polynomial(sp));
3753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
37620cf1dd8SToby Isaac }
377