xref: /petsc/src/dm/dt/space/impls/poly/spacepoly.c (revision ebead697dbf761eb322f829370bbe90b3bd93fa3)
1 #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
2 
3 static PetscErrorCode PetscSpaceSetFromOptions_Polynomial(PetscOptionItems *PetscOptionsObject,PetscSpace sp)
4 {
5   PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
6 
7   PetscFunctionBegin;
8   PetscOptionsHeadBegin(PetscOptionsObject,"PetscSpace polynomial options");
9   PetscCall(PetscOptionsBool("-petscspace_poly_tensor", "Use the tensor product polynomials", "PetscSpacePolynomialSetTensor", poly->tensor, &poly->tensor, NULL));
10   PetscOptionsHeadEnd();
11   PetscFunctionReturn(0);
12 }
13 
14 static PetscErrorCode PetscSpacePolynomialView_Ascii(PetscSpace sp, PetscViewer v)
15 {
16   PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
17 
18   PetscFunctionBegin;
19   PetscCall(PetscViewerASCIIPrintf(v, "%s space of degree %" PetscInt_FMT "\n", poly->tensor ? "Tensor polynomial" : "Polynomial", sp->degree));
20   PetscFunctionReturn(0);
21 }
22 
23 static PetscErrorCode PetscSpaceView_Polynomial(PetscSpace sp, PetscViewer viewer)
24 {
25   PetscBool      iascii;
26 
27   PetscFunctionBegin;
28   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
29   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
30   PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii));
31   if (iascii) PetscCall(PetscSpacePolynomialView_Ascii(sp, viewer));
32   PetscFunctionReturn(0);
33 }
34 
35 static PetscErrorCode PetscSpaceDestroy_Polynomial(PetscSpace sp)
36 {
37   PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
38 
39   PetscFunctionBegin;
40   PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePolynomialGetTensor_C", NULL));
41   PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePolynomialSetTensor_C", NULL));
42   if (poly->subspaces) {
43     PetscInt d;
44 
45     for (d = 0; d < sp->Nv; ++d) {
46       PetscCall(PetscSpaceDestroy(&poly->subspaces[d]));
47     }
48   }
49   PetscCall(PetscFree(poly->subspaces));
50   PetscCall(PetscFree(poly));
51   PetscFunctionReturn(0);
52 }
53 
54 static PetscErrorCode PetscSpaceSetUp_Polynomial(PetscSpace sp)
55 {
56   PetscSpace_Poly *poly    = (PetscSpace_Poly *) sp->data;
57 
58   PetscFunctionBegin;
59   if (poly->setupCalled) PetscFunctionReturn(0);
60   if (sp->Nv <= 1) {
61     poly->tensor = PETSC_FALSE;
62   }
63   if (sp->Nc != 1) {
64     PetscInt    Nc = sp->Nc;
65     PetscBool   tensor = poly->tensor;
66     PetscInt    Nv = sp->Nv;
67     PetscInt    degree = sp->degree;
68     const char *prefix;
69     const char *name;
70     char        subname[PETSC_MAX_PATH_LEN];
71     PetscSpace  subsp;
72 
73     PetscCall(PetscSpaceSetType(sp, PETSCSPACESUM));
74     PetscCall(PetscSpaceSumSetNumSubspaces(sp, Nc));
75     PetscCall(PetscSpaceCreate(PetscObjectComm((PetscObject)sp), &subsp));
76     PetscCall(PetscObjectGetOptionsPrefix((PetscObject)sp, &prefix));
77     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)subsp, prefix));
78     PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)subsp, "sumcomp_"));
79     if (((PetscObject)sp)->name) {
80       PetscCall(PetscObjectGetName((PetscObject)sp, &name));
81       PetscCall(PetscSNPrintf(subname, PETSC_MAX_PATH_LEN-1, "%s sum component", name));
82       PetscCall(PetscObjectSetName((PetscObject)subsp, subname));
83     } else PetscCall(PetscObjectSetName((PetscObject)subsp, "sum component"));
84     PetscCall(PetscSpaceSetType(subsp, PETSCSPACEPOLYNOMIAL));
85     PetscCall(PetscSpaceSetDegree(subsp, degree, PETSC_DETERMINE));
86     PetscCall(PetscSpaceSetNumComponents(subsp, 1));
87     PetscCall(PetscSpaceSetNumVariables(subsp, Nv));
88     PetscCall(PetscSpacePolynomialSetTensor(subsp, tensor));
89     PetscCall(PetscSpaceSetUp(subsp));
90     for (PetscInt i = 0; i < Nc; i++) {
91       PetscCall(PetscSpaceSumSetSubspace(sp, i, subsp));
92     }
93     PetscCall(PetscSpaceDestroy(&subsp));
94     PetscCall(PetscSpaceSetUp(sp));
95     PetscFunctionReturn(0);
96   }
97   if (poly->tensor) {
98     sp->maxDegree = PETSC_DETERMINE;
99     PetscCall(PetscSpaceSetType(sp, PETSCSPACETENSOR));
100     PetscCall(PetscSpaceSetUp(sp));
101     PetscFunctionReturn(0);
102   }
103   PetscCheck(sp->degree >= 0,PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_OUTOFRANGE, "Negative degree %" PetscInt_FMT " invalid", sp->degree);
104   sp->maxDegree = sp->degree;
105   poly->setupCalled = PETSC_TRUE;
106   PetscFunctionReturn(0);
107 }
108 
109 static PetscErrorCode PetscSpaceGetDimension_Polynomial(PetscSpace sp, PetscInt *dim)
110 {
111   PetscInt         deg  = sp->degree;
112   PetscInt         n    = sp->Nv;
113 
114   PetscFunctionBegin;
115   PetscCall(PetscDTBinomialInt(n + deg, n, dim));
116   *dim *= sp->Nc;
117   PetscFunctionReturn(0);
118 }
119 
120 static PetscErrorCode CoordinateBasis(PetscInt dim, PetscInt npoints, const PetscReal points[], PetscInt jet, PetscInt Njet, PetscReal pScalar[])
121 {
122   PetscFunctionBegin;
123   PetscCall(PetscArrayzero(pScalar, (1 + dim) * Njet * npoints));
124   for (PetscInt b = 0; b < 1 + dim; b++) {
125     for (PetscInt j = 0; j < PetscMin(1 + dim, Njet); j++) {
126       if (j == 0) {
127         if (b == 0) {
128           for (PetscInt pt = 0; pt < npoints; pt++) {
129             pScalar[b * Njet * npoints + j * npoints + pt] = 1.;
130           }
131         } else {
132           for (PetscInt pt = 0; pt < npoints; pt++) {
133             pScalar[b * Njet * npoints + j * npoints + pt] = points[pt * dim + (b-1)];
134           }
135         }
136       } else if (j == b) {
137         for (PetscInt pt = 0; pt < npoints; pt++) {
138           pScalar[b * Njet * npoints + j * npoints + pt] = 1.;
139         }
140       }
141     }
142   }
143   PetscFunctionReturn(0);
144 }
145 
146 static PetscErrorCode PetscSpaceEvaluate_Polynomial(PetscSpace sp, PetscInt npoints, const PetscReal points[], PetscReal B[], PetscReal D[], PetscReal H[])
147 {
148   PetscSpace_Poly *poly    = (PetscSpace_Poly *) sp->data;
149   DM               dm      = sp->dm;
150   PetscInt         dim     = sp->Nv;
151   PetscInt         Nb, jet, Njet;
152   PetscReal       *pScalar;
153 
154   PetscFunctionBegin;
155   if (!poly->setupCalled) {
156     PetscCall(PetscSpaceSetUp(sp));
157     PetscCall(PetscSpaceEvaluate(sp, npoints, points, B, D, H));
158     PetscFunctionReturn(0);
159   }
160   PetscCheck(!poly->tensor && sp->Nc == 1,PETSC_COMM_SELF, PETSC_ERR_PLIB, "tensor and multicomponent spaces should have been converted");
161   PetscCall(PetscDTBinomialInt(dim + sp->degree, dim, &Nb));
162   if (H) {
163     jet = 2;
164   } else if (D) {
165     jet = 1;
166   } else {
167     jet = 0;
168   }
169   PetscCall(PetscDTBinomialInt(dim + jet, dim, &Njet));
170   PetscCall(DMGetWorkArray(dm, Nb * Njet * npoints, MPIU_REAL, &pScalar));
171   // Why are we handling the case degree == 1 specially?  Because we don't want numerical noise when we evaluate hat
172   // functions at the vertices of a simplex, which happens when we invert the Vandermonde matrix of the PKD basis.
173   // We don't make any promise about which basis is used.
174   if (sp->degree == 1) {
175     PetscCall(CoordinateBasis(dim, npoints, points, jet, Njet, pScalar));
176   } else {
177     PetscCall(PetscDTPKDEvalJet(dim, npoints, points, sp->degree, jet, pScalar));
178   }
179   if (B) {
180     PetscInt p_strl = Nb;
181     PetscInt b_strl = 1;
182 
183     PetscInt b_strr = Njet*npoints;
184     PetscInt p_strr = 1;
185 
186     PetscCall(PetscArrayzero(B, npoints * Nb));
187     for (PetscInt b = 0; b < Nb; b++) {
188       for (PetscInt p = 0; p < npoints; p++) {
189         B[p*p_strl + b*b_strl] = pScalar[b*b_strr + p*p_strr];
190       }
191     }
192   }
193   if (D) {
194     PetscInt p_strl = dim*Nb;
195     PetscInt b_strl = dim;
196     PetscInt d_strl = 1;
197 
198     PetscInt b_strr = Njet*npoints;
199     PetscInt d_strr = npoints;
200     PetscInt p_strr = 1;
201 
202     PetscCall(PetscArrayzero(D, npoints * Nb * dim));
203     for (PetscInt d = 0; d < dim; d++) {
204       for (PetscInt b = 0; b < Nb; b++) {
205         for (PetscInt p = 0; p < npoints; p++) {
206           D[p*p_strl + b*b_strl + d*d_strl] = pScalar[b*b_strr + (1+d)*d_strr + p*p_strr];
207         }
208       }
209     }
210   }
211   if (H) {
212     PetscInt p_strl  = dim*dim*Nb;
213     PetscInt b_strl  = dim*dim;
214     PetscInt d1_strl = dim;
215     PetscInt d2_strl = 1;
216 
217     PetscInt b_strr = Njet*npoints;
218     PetscInt j_strr = npoints;
219     PetscInt p_strr = 1;
220 
221     PetscInt *derivs;
222     PetscCall(PetscCalloc1(dim, &derivs));
223     PetscCall(PetscArrayzero(H, npoints * Nb * dim * dim));
224     for (PetscInt d1 = 0; d1 < dim; d1++) {
225       for (PetscInt d2 = 0; d2 < dim; d2++) {
226         PetscInt j;
227         derivs[d1]++;
228         derivs[d2]++;
229         PetscCall(PetscDTGradedOrderToIndex(dim, derivs, &j));
230         derivs[d1]--;
231         derivs[d2]--;
232         for (PetscInt b = 0; b < Nb; b++) {
233           for (PetscInt p = 0; p < npoints; p++) {
234             H[p*p_strl + b*b_strl + d1*d1_strl + d2*d2_strl] = pScalar[b*b_strr + j*j_strr + p*p_strr];
235           }
236         }
237       }
238     }
239     PetscCall(PetscFree(derivs));
240   }
241   PetscCall(DMRestoreWorkArray(dm, Nb * Njet * npoints, MPIU_REAL, &pScalar));
242   PetscFunctionReturn(0);
243 }
244 
245 /*@
246   PetscSpacePolynomialSetTensor - Set whether a function space is a space of tensor polynomials (the space is spanned
247   by polynomials whose degree in each variable is bounded by the given order), as opposed to polynomials (the space is
248   spanned by polynomials whose total degree---summing over all variables---is bounded by the given order).
249 
250   Input Parameters:
251 + sp     - the function space object
252 - tensor - PETSC_TRUE for a tensor polynomial space, PETSC_FALSE for a polynomial space
253 
254   Options Database:
255 . -petscspace_poly_tensor <bool> - Whether to use tensor product polynomials in higher dimension
256 
257   Level: intermediate
258 
259 .seealso: `PetscSpacePolynomialGetTensor()`, `PetscSpaceSetDegree()`, `PetscSpaceSetNumVariables()`
260 @*/
261 PetscErrorCode PetscSpacePolynomialSetTensor(PetscSpace sp, PetscBool tensor)
262 {
263   PetscFunctionBegin;
264   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
265   PetscTryMethod(sp,"PetscSpacePolynomialSetTensor_C",(PetscSpace,PetscBool),(sp,tensor));
266   PetscFunctionReturn(0);
267 }
268 
269 /*@
270   PetscSpacePolynomialGetTensor - Get whether a function space is a space of tensor polynomials (the space is spanned
271   by polynomials whose degree in each variabl is bounded by the given order), as opposed to polynomials (the space is
272   spanned by polynomials whose total degree---summing over all variables---is bounded by the given order).
273 
274   Input Parameters:
275 . sp     - the function space object
276 
277   Output Parameters:
278 . tensor - PETSC_TRUE for a tensor polynomial space, PETSC_FALSE for a polynomial space
279 
280   Level: intermediate
281 
282 .seealso: `PetscSpacePolynomialSetTensor()`, `PetscSpaceSetDegree()`, `PetscSpaceSetNumVariables()`
283 @*/
284 PetscErrorCode PetscSpacePolynomialGetTensor(PetscSpace sp, PetscBool *tensor)
285 {
286   PetscFunctionBegin;
287   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
288   PetscValidBoolPointer(tensor, 2);
289   PetscTryMethod(sp,"PetscSpacePolynomialGetTensor_C",(PetscSpace,PetscBool*),(sp,tensor));
290   PetscFunctionReturn(0);
291 }
292 
293 static PetscErrorCode PetscSpacePolynomialSetTensor_Polynomial(PetscSpace sp, PetscBool tensor)
294 {
295   PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
296 
297   PetscFunctionBegin;
298   poly->tensor = tensor;
299   PetscFunctionReturn(0);
300 }
301 
302 static PetscErrorCode PetscSpacePolynomialGetTensor_Polynomial(PetscSpace sp, PetscBool *tensor)
303 {
304   PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
305 
306   PetscFunctionBegin;
307   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
308   PetscValidBoolPointer(tensor, 2);
309   *tensor = poly->tensor;
310   PetscFunctionReturn(0);
311 }
312 
313 static PetscErrorCode PetscSpaceGetHeightSubspace_Polynomial(PetscSpace sp, PetscInt height, PetscSpace *subsp)
314 {
315   PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
316   PetscInt         Nc, dim, order;
317   PetscBool        tensor;
318 
319   PetscFunctionBegin;
320   PetscCall(PetscSpaceGetNumComponents(sp, &Nc));
321   PetscCall(PetscSpaceGetNumVariables(sp, &dim));
322   PetscCall(PetscSpaceGetDegree(sp, &order, NULL));
323   PetscCall(PetscSpacePolynomialGetTensor(sp, &tensor));
324   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);
325   if (!poly->subspaces) PetscCall(PetscCalloc1(dim, &poly->subspaces));
326   if (height <= dim) {
327     if (!poly->subspaces[height-1]) {
328       PetscSpace  sub;
329       const char *name;
330 
331       PetscCall(PetscSpaceCreate(PetscObjectComm((PetscObject) sp), &sub));
332       PetscCall(PetscObjectGetName((PetscObject) sp,  &name));
333       PetscCall(PetscObjectSetName((PetscObject) sub,  name));
334       PetscCall(PetscSpaceSetType(sub, PETSCSPACEPOLYNOMIAL));
335       PetscCall(PetscSpaceSetNumComponents(sub, Nc));
336       PetscCall(PetscSpaceSetDegree(sub, order, PETSC_DETERMINE));
337       PetscCall(PetscSpaceSetNumVariables(sub, dim-height));
338       PetscCall(PetscSpacePolynomialSetTensor(sub, tensor));
339       PetscCall(PetscSpaceSetUp(sub));
340       poly->subspaces[height-1] = sub;
341     }
342     *subsp = poly->subspaces[height-1];
343   } else {
344     *subsp = NULL;
345   }
346   PetscFunctionReturn(0);
347 }
348 
349 static PetscErrorCode PetscSpaceInitialize_Polynomial(PetscSpace sp)
350 {
351   PetscFunctionBegin;
352   sp->ops->setfromoptions    = PetscSpaceSetFromOptions_Polynomial;
353   sp->ops->setup             = PetscSpaceSetUp_Polynomial;
354   sp->ops->view              = PetscSpaceView_Polynomial;
355   sp->ops->destroy           = PetscSpaceDestroy_Polynomial;
356   sp->ops->getdimension      = PetscSpaceGetDimension_Polynomial;
357   sp->ops->evaluate          = PetscSpaceEvaluate_Polynomial;
358   sp->ops->getheightsubspace = PetscSpaceGetHeightSubspace_Polynomial;
359   PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePolynomialGetTensor_C", PetscSpacePolynomialGetTensor_Polynomial));
360   PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePolynomialSetTensor_C", PetscSpacePolynomialSetTensor_Polynomial));
361   PetscFunctionReturn(0);
362 }
363 
364 /*MC
365   PETSCSPACEPOLYNOMIAL = "poly" - A PetscSpace object that encapsulates a polynomial space, e.g. P1 is the space of
366   linear polynomials. The space is replicated for each component.
367 
368   Level: intermediate
369 
370 .seealso: `PetscSpaceType`, `PetscSpaceCreate()`, `PetscSpaceSetType()`
371 M*/
372 
373 PETSC_EXTERN PetscErrorCode PetscSpaceCreate_Polynomial(PetscSpace sp)
374 {
375   PetscSpace_Poly *poly;
376 
377   PetscFunctionBegin;
378   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
379   PetscCall(PetscNewLog(sp,&poly));
380   sp->data = poly;
381 
382   poly->tensor    = PETSC_FALSE;
383   poly->subspaces = NULL;
384 
385   PetscCall(PetscSpaceInitialize_Polynomial(sp));
386   PetscFunctionReturn(0);
387 }
388