120cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
220cf1dd8SToby Isaac
PetscSpaceSetFromOptions_Polynomial(PetscSpace sp,PetscOptionItems PetscOptionsObject)3ce78bad3SBarry Smith 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
PetscSpacePolynomialView_Ascii(PetscSpace sp,PetscViewer v)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
PetscSpaceView_Polynomial(PetscSpace sp,PetscViewer viewer)23d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSpaceView_Polynomial(PetscSpace sp, PetscViewer viewer)
24d71ae5a4SJacob Faibussowitsch {
25*9f196a02SMartin Diehl PetscBool isascii;
2620cf1dd8SToby Isaac
2720cf1dd8SToby Isaac PetscFunctionBegin;
2820cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
2920cf1dd8SToby Isaac PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
30*9f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
31*9f196a02SMartin Diehl if (isascii) PetscCall(PetscSpacePolynomialView_Ascii(sp, viewer));
323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3320cf1dd8SToby Isaac }
3420cf1dd8SToby Isaac
PetscSpaceDestroy_Polynomial(PetscSpace sp)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
PetscSpaceSetUp_Polynomial(PetscSpace sp)52d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSpaceSetUp_Polynomial(PetscSpace sp)
53d71ae5a4SJacob Faibussowitsch {
54f1436e55SToby Isaac PetscSpace_Poly *poly = (PetscSpace_Poly *)sp->data;
55f1436e55SToby Isaac
56f1436e55SToby Isaac PetscFunctionBegin;
57371d2eb7SMartin Diehl 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));
712dce792eSToby Isaac PetscCall(PetscSpaceSumSetInterleave(sp, PETSC_TRUE, PETSC_FALSE));
729566063dSJacob Faibussowitsch PetscCall(PetscSpaceCreate(PetscObjectComm((PetscObject)sp), &subsp));
739566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)sp, &prefix));
749566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)subsp, prefix));
759566063dSJacob Faibussowitsch PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)subsp, "sumcomp_"));
76417c287bSToby Isaac if (((PetscObject)sp)->name) {
779566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)sp, &name));
789566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(subname, PETSC_MAX_PATH_LEN - 1, "%s sum component", name));
799566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)subsp, subname));
801baa6e33SBarry Smith } else PetscCall(PetscObjectSetName((PetscObject)subsp, "sum component"));
819566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetType(subsp, PETSCSPACEPOLYNOMIAL));
829566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetDegree(subsp, degree, PETSC_DETERMINE));
839566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetNumComponents(subsp, 1));
849566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetNumVariables(subsp, Nv));
859566063dSJacob Faibussowitsch PetscCall(PetscSpacePolynomialSetTensor(subsp, tensor));
869566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetUp(subsp));
8748a46eb9SPierre Jolivet for (PetscInt i = 0; i < Nc; i++) PetscCall(PetscSpaceSumSetSubspace(sp, i, subsp));
889566063dSJacob Faibussowitsch PetscCall(PetscSpaceDestroy(&subsp));
899566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetUp(sp));
903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
91f1436e55SToby Isaac }
92f1436e55SToby Isaac if (poly->tensor) {
93f1436e55SToby Isaac sp->maxDegree = PETSC_DETERMINE;
949566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetType(sp, PETSCSPACETENSOR));
959566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetUp(sp));
963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
97f1436e55SToby Isaac }
9863a3b9bcSJacob Faibussowitsch PetscCheck(sp->degree >= 0, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_OUTOFRANGE, "Negative degree %" PetscInt_FMT " invalid", sp->degree);
99f1436e55SToby Isaac sp->maxDegree = sp->degree;
100371d2eb7SMartin Diehl poly->setupcalled = PETSC_TRUE;
1013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
102f1436e55SToby Isaac }
103f1436e55SToby Isaac
PetscSpaceGetDimension_Polynomial(PetscSpace sp,PetscInt * dim)104d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSpaceGetDimension_Polynomial(PetscSpace sp, PetscInt *dim)
105d71ae5a4SJacob Faibussowitsch {
10620cf1dd8SToby Isaac PetscInt deg = sp->degree;
107f1436e55SToby Isaac PetscInt n = sp->Nv;
10820cf1dd8SToby Isaac
10920cf1dd8SToby Isaac PetscFunctionBegin;
1109566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(n + deg, n, dim));
111f1436e55SToby Isaac *dim *= sp->Nc;
1123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
11320cf1dd8SToby Isaac }
11420cf1dd8SToby Isaac
CoordinateBasis(PetscInt dim,PetscInt npoints,const PetscReal points[],PetscInt jet,PetscInt Njet,PetscReal pScalar[])115d71ae5a4SJacob Faibussowitsch static PetscErrorCode CoordinateBasis(PetscInt dim, PetscInt npoints, const PetscReal points[], PetscInt jet, PetscInt Njet, PetscReal pScalar[])
116d71ae5a4SJacob Faibussowitsch {
11720cf1dd8SToby Isaac PetscFunctionBegin;
1189566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(pScalar, (1 + dim) * Njet * npoints));
119f1436e55SToby Isaac for (PetscInt b = 0; b < 1 + dim; b++) {
120f1436e55SToby Isaac for (PetscInt j = 0; j < PetscMin(1 + dim, Njet); j++) {
121f1436e55SToby Isaac if (j == 0) {
122f1436e55SToby Isaac if (b == 0) {
123ad540459SPierre Jolivet for (PetscInt pt = 0; pt < npoints; pt++) pScalar[b * Njet * npoints + j * npoints + pt] = 1.;
12420cf1dd8SToby Isaac } else {
125ad540459SPierre Jolivet for (PetscInt pt = 0; pt < npoints; pt++) pScalar[b * Njet * npoints + j * npoints + pt] = points[pt * dim + (b - 1)];
126f1436e55SToby Isaac }
127f1436e55SToby Isaac } else if (j == b) {
128ad540459SPierre Jolivet for (PetscInt pt = 0; pt < npoints; pt++) pScalar[b * Njet * npoints + j * npoints + pt] = 1.;
129f1436e55SToby Isaac }
13020cf1dd8SToby Isaac }
13120cf1dd8SToby Isaac }
1323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
13320cf1dd8SToby Isaac }
13420cf1dd8SToby Isaac
PetscSpaceEvaluate_Polynomial(PetscSpace sp,PetscInt npoints,const PetscReal points[],PetscReal B[],PetscReal D[],PetscReal H[])135d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSpaceEvaluate_Polynomial(PetscSpace sp, PetscInt npoints, const PetscReal points[], PetscReal B[], PetscReal D[], PetscReal H[])
136d71ae5a4SJacob Faibussowitsch {
13720cf1dd8SToby Isaac PetscSpace_Poly *poly = (PetscSpace_Poly *)sp->data;
13820cf1dd8SToby Isaac DM dm = sp->dm;
13920cf1dd8SToby Isaac PetscInt dim = sp->Nv;
140f1436e55SToby Isaac PetscInt Nb, jet, Njet;
141f1436e55SToby Isaac PetscReal *pScalar;
14220cf1dd8SToby Isaac
14320cf1dd8SToby Isaac PetscFunctionBegin;
144371d2eb7SMartin Diehl if (!poly->setupcalled) {
1459566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetUp(sp));
1469566063dSJacob Faibussowitsch PetscCall(PetscSpaceEvaluate(sp, npoints, points, B, D, H));
1473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
14820cf1dd8SToby Isaac }
1491dca8a05SBarry Smith PetscCheck(!poly->tensor && sp->Nc == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "tensor and multicomponent spaces should have been converted");
1509566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(dim + sp->degree, dim, &Nb));
151f1436e55SToby Isaac if (H) {
152f1436e55SToby Isaac jet = 2;
153f1436e55SToby Isaac } else if (D) {
154f1436e55SToby Isaac jet = 1;
155f1436e55SToby Isaac } else {
156f1436e55SToby Isaac jet = 0;
15720cf1dd8SToby Isaac }
1589566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(dim + jet, dim, &Njet));
1599566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(dm, Nb * Njet * npoints, MPIU_REAL, &pScalar));
160f1436e55SToby Isaac // Why are we handling the case degree == 1 specially? Because we don't want numerical noise when we evaluate hat
161f1436e55SToby Isaac // functions at the vertices of a simplex, which happens when we invert the Vandermonde matrix of the PKD basis.
162f1436e55SToby Isaac // We don't make any promise about which basis is used.
163f1436e55SToby Isaac if (sp->degree == 1) {
1649566063dSJacob Faibussowitsch PetscCall(CoordinateBasis(dim, npoints, points, jet, Njet, pScalar));
165f1436e55SToby Isaac } else {
1669566063dSJacob Faibussowitsch PetscCall(PetscDTPKDEvalJet(dim, npoints, points, sp->degree, jet, pScalar));
16720cf1dd8SToby Isaac }
16820cf1dd8SToby Isaac if (B) {
169f1436e55SToby Isaac PetscInt p_strl = Nb;
170f1436e55SToby Isaac PetscInt b_strl = 1;
1713596293dSMatthew G. Knepley
172f1436e55SToby Isaac PetscInt b_strr = Njet * npoints;
173f1436e55SToby Isaac PetscInt p_strr = 1;
174f1436e55SToby Isaac
1759566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(B, npoints * Nb));
176f1436e55SToby Isaac for (PetscInt b = 0; b < Nb; b++) {
177ad540459SPierre Jolivet for (PetscInt p = 0; p < npoints; p++) B[p * p_strl + b * b_strl] = pScalar[b * b_strr + p * p_strr];
1783596293dSMatthew G. Knepley }
17920cf1dd8SToby Isaac }
18020cf1dd8SToby Isaac if (D) {
181f1436e55SToby Isaac PetscInt p_strl = dim * Nb;
182f1436e55SToby Isaac PetscInt b_strl = dim;
183f1436e55SToby Isaac PetscInt d_strl = 1;
184f1436e55SToby Isaac
185f1436e55SToby Isaac PetscInt b_strr = Njet * npoints;
186f1436e55SToby Isaac PetscInt d_strr = npoints;
187f1436e55SToby Isaac PetscInt p_strr = 1;
188f1436e55SToby Isaac
1899566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(D, npoints * Nb * dim));
190f1436e55SToby Isaac for (PetscInt d = 0; d < dim; d++) {
191f1436e55SToby Isaac for (PetscInt b = 0; b < Nb; b++) {
192ad540459SPierre 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];
19320cf1dd8SToby Isaac }
19420cf1dd8SToby Isaac }
19520cf1dd8SToby Isaac }
19620cf1dd8SToby Isaac if (H) {
197f1436e55SToby Isaac PetscInt p_strl = dim * dim * Nb;
198f1436e55SToby Isaac PetscInt b_strl = dim * dim;
199f1436e55SToby Isaac PetscInt d1_strl = dim;
200f1436e55SToby Isaac PetscInt d2_strl = 1;
201f1436e55SToby Isaac
202f1436e55SToby Isaac PetscInt b_strr = Njet * npoints;
203f1436e55SToby Isaac PetscInt j_strr = npoints;
204f1436e55SToby Isaac PetscInt p_strr = 1;
205f1436e55SToby Isaac
206f1436e55SToby Isaac PetscInt *derivs;
2079566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(dim, &derivs));
2089566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(H, npoints * Nb * dim * dim));
209f1436e55SToby Isaac for (PetscInt d1 = 0; d1 < dim; d1++) {
210f1436e55SToby Isaac for (PetscInt d2 = 0; d2 < dim; d2++) {
211f1436e55SToby Isaac PetscInt j;
212f1436e55SToby Isaac derivs[d1]++;
213f1436e55SToby Isaac derivs[d2]++;
2149566063dSJacob Faibussowitsch PetscCall(PetscDTGradedOrderToIndex(dim, derivs, &j));
215f1436e55SToby Isaac derivs[d1]--;
216f1436e55SToby Isaac derivs[d2]--;
217f1436e55SToby Isaac for (PetscInt b = 0; b < Nb; b++) {
218ad540459SPierre 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];
21920cf1dd8SToby Isaac }
22020cf1dd8SToby Isaac }
22120cf1dd8SToby Isaac }
2229566063dSJacob Faibussowitsch PetscCall(PetscFree(derivs));
22320cf1dd8SToby Isaac }
2249566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(dm, Nb * Njet * npoints, MPIU_REAL, &pScalar));
2253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
22620cf1dd8SToby Isaac }
22720cf1dd8SToby Isaac
22820cf1dd8SToby Isaac /*@
229a4e35b19SJacob Faibussowitsch PetscSpacePolynomialSetTensor - Set whether a function space is a space of tensor polynomials.
23020cf1dd8SToby Isaac
23120cf1dd8SToby Isaac Input Parameters:
23220cf1dd8SToby Isaac + sp - the function space object
233dce8aebaSBarry Smith - tensor - `PETSC_TRUE` for a tensor polynomial space, `PETSC_FALSE` for a polynomial space
23420cf1dd8SToby Isaac
235dce8aebaSBarry Smith Options Database Key:
2364ab77754SMatthew G. Knepley . -petscspace_poly_tensor <bool> - Whether to use tensor product polynomials in higher dimension
2374ab77754SMatthew G. Knepley
23829b5c115SMatthew G. Knepley Level: intermediate
23920cf1dd8SToby Isaac
240a4e35b19SJacob Faibussowitsch Notes:
241a4e35b19SJacob Faibussowitsch It is a tensor space if it is spanned by polynomials whose degree in each variable is
242a4e35b19SJacob Faibussowitsch bounded by the given order, as opposed to the space spanned by polynomials
243a4e35b19SJacob Faibussowitsch whose total degree---summing over all variables---is bounded by the given order.
244a4e35b19SJacob Faibussowitsch
245dce8aebaSBarry Smith .seealso: `PetscSpace`, `PetscSpacePolynomialGetTensor()`, `PetscSpaceSetDegree()`, `PetscSpaceSetNumVariables()`
24620cf1dd8SToby Isaac @*/
PetscSpacePolynomialSetTensor(PetscSpace sp,PetscBool tensor)247d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSpacePolynomialSetTensor(PetscSpace sp, PetscBool tensor)
248d71ae5a4SJacob Faibussowitsch {
24920cf1dd8SToby Isaac PetscFunctionBegin;
25020cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
251cac4c232SBarry Smith PetscTryMethod(sp, "PetscSpacePolynomialSetTensor_C", (PetscSpace, PetscBool), (sp, tensor));
2523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
25320cf1dd8SToby Isaac }
25420cf1dd8SToby Isaac
25520cf1dd8SToby Isaac /*@
256a4e35b19SJacob Faibussowitsch PetscSpacePolynomialGetTensor - Get whether a function space is a space of tensor
257a4e35b19SJacob Faibussowitsch polynomials.
25820cf1dd8SToby Isaac
2592fe279fdSBarry Smith Input Parameter:
26020cf1dd8SToby Isaac . sp - the function space object
26120cf1dd8SToby Isaac
2622fe279fdSBarry Smith Output Parameter:
263dce8aebaSBarry Smith . tensor - `PETSC_TRUE` for a tensor polynomial space, `PETSC_FALSE` for a polynomial space
26420cf1dd8SToby Isaac
26529b5c115SMatthew G. Knepley Level: intermediate
26620cf1dd8SToby Isaac
267a4e35b19SJacob Faibussowitsch Notes:
268a4e35b19SJacob Faibussowitsch The space is a tensor space if it is spanned by polynomials whose degree in each variable is
269a4e35b19SJacob Faibussowitsch bounded by the given order, as opposed to the space spanned by polynomials
270a4e35b19SJacob Faibussowitsch whose total degree---summing over all variables---is bounded by the given order.
271a4e35b19SJacob Faibussowitsch
272dce8aebaSBarry Smith .seealso: `PetscSpace`, `PetscSpacePolynomialSetTensor()`, `PetscSpaceSetDegree()`, `PetscSpaceSetNumVariables()`
27320cf1dd8SToby Isaac @*/
PetscSpacePolynomialGetTensor(PetscSpace sp,PetscBool * tensor)274d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSpacePolynomialGetTensor(PetscSpace sp, PetscBool *tensor)
275d71ae5a4SJacob Faibussowitsch {
27620cf1dd8SToby Isaac PetscFunctionBegin;
27720cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
2784f572ea9SToby Isaac PetscAssertPointer(tensor, 2);
279cac4c232SBarry Smith PetscTryMethod(sp, "PetscSpacePolynomialGetTensor_C", (PetscSpace, PetscBool *), (sp, tensor));
2803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
28120cf1dd8SToby Isaac }
28220cf1dd8SToby Isaac
PetscSpacePolynomialSetTensor_Polynomial(PetscSpace sp,PetscBool tensor)283d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSpacePolynomialSetTensor_Polynomial(PetscSpace sp, PetscBool tensor)
284d71ae5a4SJacob Faibussowitsch {
28520cf1dd8SToby Isaac PetscSpace_Poly *poly = (PetscSpace_Poly *)sp->data;
28620cf1dd8SToby Isaac
28720cf1dd8SToby Isaac PetscFunctionBegin;
28820cf1dd8SToby Isaac poly->tensor = tensor;
2893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
29020cf1dd8SToby Isaac }
29120cf1dd8SToby Isaac
PetscSpacePolynomialGetTensor_Polynomial(PetscSpace sp,PetscBool * tensor)292d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSpacePolynomialGetTensor_Polynomial(PetscSpace sp, PetscBool *tensor)
293d71ae5a4SJacob Faibussowitsch {
29420cf1dd8SToby Isaac PetscSpace_Poly *poly = (PetscSpace_Poly *)sp->data;
29520cf1dd8SToby Isaac
29620cf1dd8SToby Isaac PetscFunctionBegin;
29720cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
2984f572ea9SToby Isaac PetscAssertPointer(tensor, 2);
29920cf1dd8SToby Isaac *tensor = poly->tensor;
3003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
30120cf1dd8SToby Isaac }
30220cf1dd8SToby Isaac
PetscSpaceGetHeightSubspace_Polynomial(PetscSpace sp,PetscInt height,PetscSpace * subsp)303d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSpaceGetHeightSubspace_Polynomial(PetscSpace sp, PetscInt height, PetscSpace *subsp)
304d71ae5a4SJacob Faibussowitsch {
30520cf1dd8SToby Isaac PetscSpace_Poly *poly = (PetscSpace_Poly *)sp->data;
30620cf1dd8SToby Isaac PetscInt Nc, dim, order;
30720cf1dd8SToby Isaac PetscBool tensor;
30820cf1dd8SToby Isaac
30920cf1dd8SToby Isaac PetscFunctionBegin;
3109566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetNumComponents(sp, &Nc));
3119566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetNumVariables(sp, &dim));
3129566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetDegree(sp, &order, NULL));
3139566063dSJacob Faibussowitsch PetscCall(PetscSpacePolynomialGetTensor(sp, &tensor));
3141dca8a05SBarry 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);
3159566063dSJacob Faibussowitsch if (!poly->subspaces) PetscCall(PetscCalloc1(dim, &poly->subspaces));
31620cf1dd8SToby Isaac if (height <= dim) {
31720cf1dd8SToby Isaac if (!poly->subspaces[height - 1]) {
31820cf1dd8SToby Isaac PetscSpace sub;
3193f6b16c7SMatthew G. Knepley const char *name;
32020cf1dd8SToby Isaac
3219566063dSJacob Faibussowitsch PetscCall(PetscSpaceCreate(PetscObjectComm((PetscObject)sp), &sub));
3229566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)sp, &name));
3239566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)sub, name));
3249566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetType(sub, PETSCSPACEPOLYNOMIAL));
3259566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetNumComponents(sub, Nc));
3269566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetDegree(sub, order, PETSC_DETERMINE));
3279566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetNumVariables(sub, dim - height));
3289566063dSJacob Faibussowitsch PetscCall(PetscSpacePolynomialSetTensor(sub, tensor));
3299566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetUp(sub));
33020cf1dd8SToby Isaac poly->subspaces[height - 1] = sub;
33120cf1dd8SToby Isaac }
33220cf1dd8SToby Isaac *subsp = poly->subspaces[height - 1];
33320cf1dd8SToby Isaac } else {
33420cf1dd8SToby Isaac *subsp = NULL;
33520cf1dd8SToby Isaac }
3363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
33720cf1dd8SToby Isaac }
33820cf1dd8SToby Isaac
PetscSpaceInitialize_Polynomial(PetscSpace sp)339d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSpaceInitialize_Polynomial(PetscSpace sp)
340d71ae5a4SJacob Faibussowitsch {
34120cf1dd8SToby Isaac PetscFunctionBegin;
34220cf1dd8SToby Isaac sp->ops->setfromoptions = PetscSpaceSetFromOptions_Polynomial;
34320cf1dd8SToby Isaac sp->ops->setup = PetscSpaceSetUp_Polynomial;
34420cf1dd8SToby Isaac sp->ops->view = PetscSpaceView_Polynomial;
34520cf1dd8SToby Isaac sp->ops->destroy = PetscSpaceDestroy_Polynomial;
34620cf1dd8SToby Isaac sp->ops->getdimension = PetscSpaceGetDimension_Polynomial;
34720cf1dd8SToby Isaac sp->ops->evaluate = PetscSpaceEvaluate_Polynomial;
34820cf1dd8SToby Isaac sp->ops->getheightsubspace = PetscSpaceGetHeightSubspace_Polynomial;
3499566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscSpacePolynomialGetTensor_C", PetscSpacePolynomialGetTensor_Polynomial));
3509566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscSpacePolynomialSetTensor_C", PetscSpacePolynomialSetTensor_Polynomial));
3513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
35220cf1dd8SToby Isaac }
35320cf1dd8SToby Isaac
35420cf1dd8SToby Isaac /*MC
355dce8aebaSBarry Smith PETSCSPACEPOLYNOMIAL = "poly" - A `PetscSpace` object that encapsulates a polynomial space, e.g. P1 is the space of
35620cf1dd8SToby Isaac linear polynomials. The space is replicated for each component.
35720cf1dd8SToby Isaac
35820cf1dd8SToby Isaac Level: intermediate
35920cf1dd8SToby Isaac
360dce8aebaSBarry Smith .seealso: `PetscSpace`, `PetscSpaceType`, `PetscSpaceCreate()`, `PetscSpaceSetType()`
36120cf1dd8SToby Isaac M*/
36220cf1dd8SToby Isaac
PetscSpaceCreate_Polynomial(PetscSpace sp)363d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscSpaceCreate_Polynomial(PetscSpace sp)
364d71ae5a4SJacob Faibussowitsch {
36520cf1dd8SToby Isaac PetscSpace_Poly *poly;
36620cf1dd8SToby Isaac
36720cf1dd8SToby Isaac PetscFunctionBegin;
36820cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
3694dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&poly));
37020cf1dd8SToby Isaac sp->data = poly;
37120cf1dd8SToby Isaac
37220cf1dd8SToby Isaac poly->tensor = PETSC_FALSE;
37320cf1dd8SToby Isaac poly->subspaces = NULL;
37420cf1dd8SToby Isaac
3759566063dSJacob Faibussowitsch PetscCall(PetscSpaceInitialize_Polynomial(sp));
3763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
37720cf1dd8SToby Isaac }
378