xref: /petsc/src/dm/dt/space/impls/ptrimmed/spaceptrimmed.c (revision 9d0fdff2e8b5ebdd0e57bf4dfe77e45fec5df614)
1 #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
2 
3 static PetscErrorCode PetscSpaceSetFromOptions_Ptrimmed(PetscSpace sp, PetscOptionItems *PetscOptionsObject) {
4   PetscSpace_Ptrimmed *pt = (PetscSpace_Ptrimmed *)sp->data;
5 
6   PetscFunctionBegin;
7   PetscOptionsHeadBegin(PetscOptionsObject, "PetscSpace polynomial options");
8   PetscCall(PetscOptionsInt("-petscspace_ptrimmed_form_degree", "form degree of trimmed space", "PetscSpacePTrimmedSetFormDegree", pt->formDegree, &(pt->formDegree), NULL));
9   PetscOptionsHeadEnd();
10   PetscFunctionReturn(0);
11 }
12 
13 static PetscErrorCode PetscSpacePTrimmedView_Ascii(PetscSpace sp, PetscViewer v) {
14   PetscSpace_Ptrimmed *pt = (PetscSpace_Ptrimmed *)sp->data;
15   PetscInt             f, tdegree;
16 
17   PetscFunctionBegin;
18   f       = pt->formDegree;
19   tdegree = f == 0 ? sp->degree : sp->degree + 1;
20   PetscCall(PetscViewerASCIIPrintf(v, "Trimmed polynomials %" PetscInt_FMT "%s-forms of degree %" PetscInt_FMT " (P-%" PetscInt_FMT "/\\%" PetscInt_FMT ")\n", PetscAbsInt(f), f < 0 ? "*" : "", sp->degree, tdegree, PetscAbsInt(f)));
21   PetscFunctionReturn(0);
22 }
23 
24 static PetscErrorCode PetscSpaceView_Ptrimmed(PetscSpace sp, PetscViewer viewer) {
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(PetscSpacePTrimmedView_Ascii(sp, viewer));
32   PetscFunctionReturn(0);
33 }
34 
35 static PetscErrorCode PetscSpaceDestroy_Ptrimmed(PetscSpace sp) {
36   PetscSpace_Ptrimmed *pt = (PetscSpace_Ptrimmed *)sp->data;
37 
38   PetscFunctionBegin;
39   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscSpacePTrimmedGetFormDegree_C", NULL));
40   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscSpacePTrimmedSetFormDegree_C", NULL));
41   if (pt->subspaces) {
42     PetscInt d;
43 
44     for (d = 0; d < sp->Nv; ++d) PetscCall(PetscSpaceDestroy(&pt->subspaces[d]));
45   }
46   PetscCall(PetscFree(pt->subspaces));
47   PetscCall(PetscFree(pt));
48   PetscFunctionReturn(0);
49 }
50 
51 static PetscErrorCode PetscSpaceSetUp_Ptrimmed(PetscSpace sp) {
52   PetscSpace_Ptrimmed *pt = (PetscSpace_Ptrimmed *)sp->data;
53   PetscInt             Nf;
54 
55   PetscFunctionBegin;
56   if (pt->setupCalled) PetscFunctionReturn(0);
57   PetscCheck(pt->formDegree >= -sp->Nv && pt->formDegree <= sp->Nv, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_OUTOFRANGE, "Form degree %" PetscInt_FMT " not in valid range [%" PetscInt_FMT ",%" PetscInt_FMT "]", pt->formDegree, sp->Nv, sp->Nv);
58   PetscCall(PetscDTBinomialInt(sp->Nv, PetscAbsInt(pt->formDegree), &Nf));
59   if (sp->Nc == PETSC_DETERMINE) sp->Nc = Nf;
60   PetscCheck(sp->Nc % Nf == 0, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_INCOMP, "Number of components %" PetscInt_FMT " is not a multiple of form dimension %" PetscInt_FMT, sp->Nc, Nf);
61   if (sp->Nc != Nf) {
62     PetscSpace  subsp;
63     PetscInt    nCopies = sp->Nc / Nf;
64     PetscInt    Nv, deg, maxDeg;
65     PetscInt    formDegree = pt->formDegree;
66     const char *prefix;
67     const char *name;
68     char        subname[PETSC_MAX_PATH_LEN];
69 
70     PetscCall(PetscSpaceSetType(sp, PETSCSPACESUM));
71     PetscCall(PetscSpaceSumSetConcatenate(sp, PETSC_TRUE));
72     PetscCall(PetscSpaceSumSetNumSubspaces(sp, nCopies));
73     PetscCall(PetscSpaceCreate(PetscObjectComm((PetscObject)sp), &subsp));
74     PetscCall(PetscObjectGetOptionsPrefix((PetscObject)sp, &prefix));
75     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)subsp, prefix));
76     PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)subsp, "sumcomp_"));
77     if (((PetscObject)sp)->name) {
78       PetscCall(PetscObjectGetName((PetscObject)sp, &name));
79       PetscCall(PetscSNPrintf(subname, PETSC_MAX_PATH_LEN - 1, "%s sum component", name));
80       PetscCall(PetscObjectSetName((PetscObject)subsp, subname));
81     } else PetscCall(PetscObjectSetName((PetscObject)subsp, "sum component"));
82     PetscCall(PetscSpaceSetType(subsp, PETSCSPACEPTRIMMED));
83     PetscCall(PetscSpaceGetNumVariables(sp, &Nv));
84     PetscCall(PetscSpaceSetNumVariables(subsp, Nv));
85     PetscCall(PetscSpaceSetNumComponents(subsp, Nf));
86     PetscCall(PetscSpaceGetDegree(sp, &deg, &maxDeg));
87     PetscCall(PetscSpaceSetDegree(subsp, deg, maxDeg));
88     PetscCall(PetscSpacePTrimmedSetFormDegree(subsp, formDegree));
89     PetscCall(PetscSpaceSetUp(subsp));
90     for (PetscInt i = 0; i < nCopies; i++) PetscCall(PetscSpaceSumSetSubspace(sp, i, subsp));
91     PetscCall(PetscSpaceDestroy(&subsp));
92     PetscCall(PetscSpaceSetUp(sp));
93     PetscFunctionReturn(0);
94   }
95   if (sp->degree == PETSC_DEFAULT) {
96     sp->degree = 0;
97   } else if (sp->degree < 0) {
98     SETERRQ(PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_OUTOFRANGE, "Invalid negative degree %" PetscInt_FMT, sp->degree);
99   }
100   sp->maxDegree = (pt->formDegree == 0 || PetscAbsInt(pt->formDegree) == sp->Nv) ? sp->degree : sp->degree + 1;
101   if (pt->formDegree == 0 || PetscAbsInt(pt->formDegree) == sp->Nv) {
102     // Convert to regular polynomial space
103     PetscCall(PetscSpaceSetType(sp, PETSCSPACEPOLYNOMIAL));
104     PetscCall(PetscSpaceSetUp(sp));
105     PetscFunctionReturn(0);
106   }
107   pt->setupCalled = PETSC_TRUE;
108   PetscFunctionReturn(0);
109 }
110 
111 static PetscErrorCode PetscSpaceGetDimension_Ptrimmed(PetscSpace sp, PetscInt *dim) {
112   PetscSpace_Ptrimmed *pt = (PetscSpace_Ptrimmed *)sp->data;
113   PetscInt             f;
114   PetscInt             Nf;
115 
116   PetscFunctionBegin;
117   f = pt->formDegree;
118   // For PetscSpace, degree refers to the largest complete polynomial degree contained in the space which
119   // is equal to the index of a P trimmed space only for 0-forms: otherwise, the index is degree + 1
120   PetscCall(PetscDTPTrimmedSize(sp->Nv, f == 0 ? sp->degree : sp->degree + 1, pt->formDegree, dim));
121   PetscCall(PetscDTBinomialInt(sp->Nv, PetscAbsInt(pt->formDegree), &Nf));
122   *dim *= (sp->Nc / Nf);
123   PetscFunctionReturn(0);
124 }
125 
126 /*
127   p in [0, npoints), i in [0, pdim), c in [0, Nc)
128 
129   B[p][i][c] = B[p][i_scalar][c][c]
130 */
131 static PetscErrorCode PetscSpaceEvaluate_Ptrimmed(PetscSpace sp, PetscInt npoints, const PetscReal points[], PetscReal B[], PetscReal D[], PetscReal H[]) {
132   PetscSpace_Ptrimmed *pt = (PetscSpace_Ptrimmed *)sp->data;
133   DM                   dm = sp->dm;
134   PetscInt             jet, degree, Nf, Ncopies, Njet;
135   PetscInt             Nc = sp->Nc;
136   PetscInt             f;
137   PetscInt             dim = sp->Nv;
138   PetscReal           *eval;
139   PetscInt             Nb;
140 
141   PetscFunctionBegin;
142   if (!pt->setupCalled) {
143     PetscCall(PetscSpaceSetUp(sp));
144     PetscCall(PetscSpaceEvaluate(sp, npoints, points, B, D, H));
145     PetscFunctionReturn(0);
146   }
147   if (H) {
148     jet = 2;
149   } else if (D) {
150     jet = 1;
151   } else {
152     jet = 0;
153   }
154   f      = pt->formDegree;
155   degree = f == 0 ? sp->degree : sp->degree + 1;
156   PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(f), &Nf));
157   Ncopies = Nc / Nf;
158   PetscCheck(Ncopies == 1, PetscObjectComm((PetscObject)sp), PETSC_ERR_PLIB, "Multicopy spaces should have been converted to PETSCSPACESUM");
159   PetscCall(PetscDTBinomialInt(dim + jet, dim, &Njet));
160   PetscCall(PetscDTPTrimmedSize(dim, degree, f, &Nb));
161   PetscCall(DMGetWorkArray(dm, Nb * Nf * Njet * npoints, MPIU_REAL, &eval));
162   PetscCall(PetscDTPTrimmedEvalJet(dim, npoints, points, degree, f, jet, eval));
163   if (B) {
164     PetscInt p_strl = Nf * Nb;
165     PetscInt b_strl = Nf;
166     PetscInt v_strl = 1;
167 
168     PetscInt b_strr = Nf * Njet * npoints;
169     PetscInt v_strr = Njet * npoints;
170     PetscInt p_strr = 1;
171 
172     for (PetscInt v = 0; v < Nf; v++) {
173       for (PetscInt b = 0; b < Nb; b++) {
174         for (PetscInt p = 0; p < npoints; p++) B[p * p_strl + b * b_strl + v * v_strl] = eval[b * b_strr + v * v_strr + p * p_strr];
175       }
176     }
177   }
178   if (D) {
179     PetscInt p_strl = dim * Nf * Nb;
180     PetscInt b_strl = dim * Nf;
181     PetscInt v_strl = dim;
182     PetscInt d_strl = 1;
183 
184     PetscInt b_strr = Nf * Njet * npoints;
185     PetscInt v_strr = Njet * npoints;
186     PetscInt d_strr = npoints;
187     PetscInt p_strr = 1;
188 
189     for (PetscInt v = 0; v < Nf; v++) {
190       for (PetscInt d = 0; d < dim; d++) {
191         for (PetscInt b = 0; b < Nb; b++) {
192           for (PetscInt p = 0; p < npoints; p++) D[p * p_strl + b * b_strl + v * v_strl + d * d_strl] = eval[b * b_strr + v * v_strr + (1 + d) * d_strr + p * p_strr];
193         }
194       }
195     }
196   }
197   if (H) {
198     PetscInt p_strl  = dim * dim * Nf * Nb;
199     PetscInt b_strl  = dim * dim * Nf;
200     PetscInt v_strl  = dim * dim;
201     PetscInt d1_strl = dim;
202     PetscInt d2_strl = 1;
203 
204     PetscInt b_strr = Nf * Njet * npoints;
205     PetscInt v_strr = Njet * npoints;
206     PetscInt j_strr = npoints;
207     PetscInt p_strr = 1;
208 
209     PetscInt *derivs;
210     PetscCall(PetscCalloc1(dim, &derivs));
211     for (PetscInt d1 = 0; d1 < dim; d1++) {
212       for (PetscInt d2 = 0; d2 < dim; d2++) {
213         PetscInt j;
214         derivs[d1]++;
215         derivs[d2]++;
216         PetscCall(PetscDTGradedOrderToIndex(dim, derivs, &j));
217         derivs[d1]--;
218         derivs[d2]--;
219         for (PetscInt v = 0; v < Nf; v++) {
220           for (PetscInt b = 0; b < Nb; b++) {
221             for (PetscInt p = 0; p < npoints; p++) H[p * p_strl + b * b_strl + v * v_strl + d1 * d1_strl + d2 * d2_strl] = eval[b * b_strr + v * v_strr + j * j_strr + p * p_strr];
222           }
223         }
224       }
225     }
226     PetscCall(PetscFree(derivs));
227   }
228   PetscCall(DMRestoreWorkArray(dm, Nb * Nf * Njet * npoints, MPIU_REAL, &eval));
229   PetscFunctionReturn(0);
230 }
231 
232 /*@
233   PetscSpacePTrimmedSetFormDegree - Set the form degree of the trimmed polynomials.
234 
235   Input Parameters:
236 + sp         - the function space object
237 - formDegree - the form degree
238 
239   Options Database:
240 . -petscspace_ptrimmed_form_degree <int> - The trimmed polynomial form degree
241 
242   Level: intermediate
243 
244 .seealso: `PetscDTAltV`, `PetscDTPTrimmedEvalJet()`, `PetscSpacePTrimmedGetFormDegree()`
245 @*/
246 PetscErrorCode PetscSpacePTrimmedSetFormDegree(PetscSpace sp, PetscInt formDegree) {
247   PetscFunctionBegin;
248   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
249   PetscTryMethod(sp, "PetscSpacePTrimmedSetFormDegree_C", (PetscSpace, PetscInt), (sp, formDegree));
250   PetscFunctionReturn(0);
251 }
252 
253 /*@
254   PetscSpacePTrimmedGetFormDegree - Get the form degree of the trimmed polynomials.
255 
256   Input Parameters:
257 . sp     - the function space object
258 
259   Output Parameters:
260 . formDegee - the form degree
261 
262   Level: intermediate
263 
264 .seealso: `PetscDTAltV`, `PetscDTPTrimmedEvalJet()`, `PetscSpacePTrimmedSetFormDegree()`
265 @*/
266 PetscErrorCode PetscSpacePTrimmedGetFormDegree(PetscSpace sp, PetscInt *formDegree) {
267   PetscFunctionBegin;
268   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
269   PetscValidIntPointer(formDegree, 2);
270   PetscTryMethod(sp, "PetscSpacePTrimmedGetFormDegree_C", (PetscSpace, PetscInt *), (sp, formDegree));
271   PetscFunctionReturn(0);
272 }
273 
274 static PetscErrorCode PetscSpacePTrimmedSetFormDegree_Ptrimmed(PetscSpace sp, PetscInt formDegree) {
275   PetscSpace_Ptrimmed *pt = (PetscSpace_Ptrimmed *)sp->data;
276 
277   PetscFunctionBegin;
278   pt->formDegree = formDegree;
279   PetscFunctionReturn(0);
280 }
281 
282 static PetscErrorCode PetscSpacePTrimmedGetFormDegree_Ptrimmed(PetscSpace sp, PetscInt *formDegree) {
283   PetscSpace_Ptrimmed *pt = (PetscSpace_Ptrimmed *)sp->data;
284 
285   PetscFunctionBegin;
286   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
287   PetscValidIntPointer(formDegree, 2);
288   *formDegree = pt->formDegree;
289   PetscFunctionReturn(0);
290 }
291 
292 static PetscErrorCode PetscSpaceGetHeightSubspace_Ptrimmed(PetscSpace sp, PetscInt height, PetscSpace *subsp) {
293   PetscSpace_Ptrimmed *pt = (PetscSpace_Ptrimmed *)sp->data;
294   PetscInt             dim;
295 
296   PetscFunctionBegin;
297   PetscCall(PetscSpaceGetNumVariables(sp, &dim));
298   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);
299   if (!pt->subspaces) PetscCall(PetscCalloc1(dim, &(pt->subspaces)));
300   if ((dim - height) <= PetscAbsInt(pt->formDegree)) {
301     if (!pt->subspaces[height - 1]) {
302       PetscInt    Nc, degree, Nf, Ncopies, Nfsub;
303       PetscSpace  sub;
304       const char *name;
305 
306       PetscCall(PetscSpaceGetNumComponents(sp, &Nc));
307       PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(pt->formDegree), &Nf));
308       PetscCall(PetscDTBinomialInt((dim - height), PetscAbsInt(pt->formDegree), &Nfsub));
309       Ncopies = Nf / Nc;
310       PetscCall(PetscSpaceGetDegree(sp, &degree, NULL));
311 
312       PetscCall(PetscSpaceCreate(PetscObjectComm((PetscObject)sp), &sub));
313       PetscCall(PetscObjectGetName((PetscObject)sp, &name));
314       PetscCall(PetscObjectSetName((PetscObject)sub, name));
315       PetscCall(PetscSpaceSetType(sub, PETSCSPACEPTRIMMED));
316       PetscCall(PetscSpaceSetNumComponents(sub, Nfsub * Ncopies));
317       PetscCall(PetscSpaceSetDegree(sub, degree, PETSC_DETERMINE));
318       PetscCall(PetscSpaceSetNumVariables(sub, dim - height));
319       PetscCall(PetscSpacePTrimmedSetFormDegree(sub, pt->formDegree));
320       PetscCall(PetscSpaceSetUp(sub));
321       pt->subspaces[height - 1] = sub;
322     }
323     *subsp = pt->subspaces[height - 1];
324   } else {
325     *subsp = NULL;
326   }
327   PetscFunctionReturn(0);
328 }
329 
330 static PetscErrorCode PetscSpaceInitialize_Ptrimmed(PetscSpace sp) {
331   PetscFunctionBegin;
332   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscSpacePTrimmedGetFormDegree_C", PetscSpacePTrimmedGetFormDegree_Ptrimmed));
333   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscSpacePTrimmedSetFormDegree_C", PetscSpacePTrimmedSetFormDegree_Ptrimmed));
334   sp->ops->setfromoptions    = PetscSpaceSetFromOptions_Ptrimmed;
335   sp->ops->setup             = PetscSpaceSetUp_Ptrimmed;
336   sp->ops->view              = PetscSpaceView_Ptrimmed;
337   sp->ops->destroy           = PetscSpaceDestroy_Ptrimmed;
338   sp->ops->getdimension      = PetscSpaceGetDimension_Ptrimmed;
339   sp->ops->evaluate          = PetscSpaceEvaluate_Ptrimmed;
340   sp->ops->getheightsubspace = PetscSpaceGetHeightSubspace_Ptrimmed;
341   PetscFunctionReturn(0);
342 }
343 
344 /*MC
345   PETSCSPACEPTRIMMED = "ptrimmed" - A PetscSpace object that encapsulates a trimmed polynomial space.
346 
347   Level: intermediate
348 
349 .seealso: `PetscSpaceType`, `PetscSpaceCreate()`, `PetscSpaceSetType()`, `PetscDTPTrimmedEvalJet()`
350 M*/
351 
352 PETSC_EXTERN PetscErrorCode PetscSpaceCreate_Ptrimmed(PetscSpace sp) {
353   PetscSpace_Ptrimmed *pt;
354 
355   PetscFunctionBegin;
356   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
357   PetscCall(PetscNew(&pt));
358   sp->data = pt;
359 
360   pt->subspaces = NULL;
361   sp->Nc        = PETSC_DETERMINE;
362 
363   PetscCall(PetscSpaceInitialize_Ptrimmed(sp));
364   PetscFunctionReturn(0);
365 }
366