xref: /petsc/src/dm/dt/dualspace/impls/sum/dualspacesum.c (revision 2dce792e531186164765a9583d36d03ffc15e9ea)
1*2dce792eSToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
2*2dce792eSToby Isaac 
3*2dce792eSToby Isaac /*@
4*2dce792eSToby Isaac   PetscDualSpaceSumGetNumSubspaces - Get the number of spaces in the sum space
5*2dce792eSToby Isaac 
6*2dce792eSToby Isaac   Input Parameter:
7*2dce792eSToby Isaac . sp - the dual space object
8*2dce792eSToby Isaac 
9*2dce792eSToby Isaac   Output Parameter:
10*2dce792eSToby Isaac . numSumSpaces - the number of spaces
11*2dce792eSToby Isaac 
12*2dce792eSToby Isaac   Level: intermediate
13*2dce792eSToby Isaac 
14*2dce792eSToby Isaac   Note:
15*2dce792eSToby Isaac   The name NumSubspaces is slightly misleading because it is actually getting the number of defining spaces of the sum, not a number of Subspaces of it
16*2dce792eSToby Isaac 
17*2dce792eSToby Isaac .seealso: `PETSCDUALSPACESUM`, `PetscDualSpace`, `PetscDualSpaceSumSetNumSubspaces()`
18*2dce792eSToby Isaac @*/
19*2dce792eSToby Isaac PetscErrorCode PetscDualSpaceSumGetNumSubspaces(PetscDualSpace sp, PetscInt *numSumSpaces)
20*2dce792eSToby Isaac {
21*2dce792eSToby Isaac   PetscFunctionBegin;
22*2dce792eSToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
23*2dce792eSToby Isaac   PetscAssertPointer(numSumSpaces, 2);
24*2dce792eSToby Isaac   PetscTryMethod(sp, "PetscDualSpaceSumGetNumSubspaces_C", (PetscDualSpace, PetscInt *), (sp, numSumSpaces));
25*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
26*2dce792eSToby Isaac }
27*2dce792eSToby Isaac 
28*2dce792eSToby Isaac /*@
29*2dce792eSToby Isaac   PetscDualSpaceSumSetNumSubspaces - Set the number of spaces in the sum space
30*2dce792eSToby Isaac 
31*2dce792eSToby Isaac   Input Parameters:
32*2dce792eSToby Isaac + sp           - the dual space object
33*2dce792eSToby Isaac - numSumSpaces - the number of spaces
34*2dce792eSToby Isaac 
35*2dce792eSToby Isaac   Level: intermediate
36*2dce792eSToby Isaac 
37*2dce792eSToby Isaac   Note:
38*2dce792eSToby Isaac   The name NumSubspaces is slightly misleading because it is actually setting the number of defining spaces of the sum, not a number of Subspaces of it
39*2dce792eSToby Isaac 
40*2dce792eSToby Isaac .seealso: `PETSCDUALSPACESUM`, `PetscDualSpace`, `PetscDualSpaceSumGetNumSubspaces()`
41*2dce792eSToby Isaac @*/
42*2dce792eSToby Isaac PetscErrorCode PetscDualSpaceSumSetNumSubspaces(PetscDualSpace sp, PetscInt numSumSpaces)
43*2dce792eSToby Isaac {
44*2dce792eSToby Isaac   PetscFunctionBegin;
45*2dce792eSToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
46*2dce792eSToby Isaac   PetscTryMethod(sp, "PetscDualSpaceSumSetNumSubspaces_C", (PetscDualSpace, PetscInt), (sp, numSumSpaces));
47*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
48*2dce792eSToby Isaac }
49*2dce792eSToby Isaac 
50*2dce792eSToby Isaac /*@
51*2dce792eSToby Isaac   PetscDualSpaceSumGetConcatenate - Get the concatenate flag for this space.
52*2dce792eSToby Isaac 
53*2dce792eSToby Isaac   Input Parameter:
54*2dce792eSToby Isaac . sp - the dual space object
55*2dce792eSToby Isaac 
56*2dce792eSToby Isaac   Output Parameter:
57*2dce792eSToby Isaac . concatenate - flag indicating whether subspaces are concatenated.
58*2dce792eSToby Isaac 
59*2dce792eSToby Isaac   Level: intermediate
60*2dce792eSToby Isaac 
61*2dce792eSToby Isaac   Note:
62*2dce792eSToby Isaac   A concatenated sum space will have the number of components equal to the sum of the number of
63*2dce792eSToby Isaac   components of all subspaces. A non-concatenated, or direct sum space will have the same
64*2dce792eSToby Isaac   number of components as its subspaces.
65*2dce792eSToby Isaac 
66*2dce792eSToby Isaac .seealso: `PETSCDUALSPACESUM`, `PetscDualSpace`, `PetscDualSpaceSumSetConcatenate()`
67*2dce792eSToby Isaac @*/
68*2dce792eSToby Isaac PetscErrorCode PetscDualSpaceSumGetConcatenate(PetscDualSpace sp, PetscBool *concatenate)
69*2dce792eSToby Isaac {
70*2dce792eSToby Isaac   PetscFunctionBegin;
71*2dce792eSToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
72*2dce792eSToby Isaac   PetscTryMethod(sp, "PetscDualSpaceSumGetConcatenate_C", (PetscDualSpace, PetscBool *), (sp, concatenate));
73*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
74*2dce792eSToby Isaac }
75*2dce792eSToby Isaac 
76*2dce792eSToby Isaac /*@
77*2dce792eSToby Isaac   PetscDualSpaceSumSetConcatenate - Sets the concatenate flag for this space.
78*2dce792eSToby Isaac 
79*2dce792eSToby Isaac   Input Parameters:
80*2dce792eSToby Isaac + sp          - the dual space object
81*2dce792eSToby Isaac - concatenate - are subspaces concatenated components (true) or direct summands (false)
82*2dce792eSToby Isaac 
83*2dce792eSToby Isaac   Level: intermediate
84*2dce792eSToby Isaac 
85*2dce792eSToby Isaac   Notes:
86*2dce792eSToby Isaac   A concatenated sum space will have the number of components equal to the sum of the number of
87*2dce792eSToby Isaac   components of all subspaces. A non-concatenated, or direct sum space will have the same
88*2dce792eSToby Isaac   number of components as its subspaces .
89*2dce792eSToby Isaac 
90*2dce792eSToby Isaac .seealso: `PETSCDUALSPACESUM`, `PetscDualSpace`, `PetscDualSpaceSumGetConcatenate()`
91*2dce792eSToby Isaac @*/
92*2dce792eSToby Isaac PetscErrorCode PetscDualSpaceSumSetConcatenate(PetscDualSpace sp, PetscBool concatenate)
93*2dce792eSToby Isaac {
94*2dce792eSToby Isaac   PetscFunctionBegin;
95*2dce792eSToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
96*2dce792eSToby Isaac   PetscTryMethod(sp, "PetscDualSpaceSumSetConcatenate_C", (PetscDualSpace, PetscBool), (sp, concatenate));
97*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
98*2dce792eSToby Isaac }
99*2dce792eSToby Isaac 
100*2dce792eSToby Isaac /*@
101*2dce792eSToby Isaac   PetscDualSpaceSumGetSubspace - Get a space in the sum space
102*2dce792eSToby Isaac 
103*2dce792eSToby Isaac   Input Parameters:
104*2dce792eSToby Isaac + sp - the dual space object
105*2dce792eSToby Isaac - s  - The space number
106*2dce792eSToby Isaac 
107*2dce792eSToby Isaac   Output Parameter:
108*2dce792eSToby Isaac . subsp - the `PetscDualSpace`
109*2dce792eSToby Isaac 
110*2dce792eSToby Isaac   Level: intermediate
111*2dce792eSToby Isaac 
112*2dce792eSToby Isaac   Note:
113*2dce792eSToby Isaac   The name GetSubspace is slightly misleading because it is actually getting one of the defining spaces of the sum, not a Subspace of it
114*2dce792eSToby Isaac 
115*2dce792eSToby Isaac .seealso: `PETSCDUALSPACESUM`, `PetscDualSpace`, `PetscDualSpaceSumSetSubspace()`
116*2dce792eSToby Isaac @*/
117*2dce792eSToby Isaac PetscErrorCode PetscDualSpaceSumGetSubspace(PetscDualSpace sp, PetscInt s, PetscDualSpace *subsp)
118*2dce792eSToby Isaac {
119*2dce792eSToby Isaac   PetscFunctionBegin;
120*2dce792eSToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
121*2dce792eSToby Isaac   PetscAssertPointer(subsp, 3);
122*2dce792eSToby Isaac   PetscTryMethod(sp, "PetscDualSpaceSumGetSubspace_C", (PetscDualSpace, PetscInt, PetscDualSpace *), (sp, s, subsp));
123*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
124*2dce792eSToby Isaac }
125*2dce792eSToby Isaac 
126*2dce792eSToby Isaac /*@
127*2dce792eSToby Isaac   PetscDualSpaceSumSetSubspace - Set a space in the sum space
128*2dce792eSToby Isaac 
129*2dce792eSToby Isaac   Input Parameters:
130*2dce792eSToby Isaac + sp    - the dual space object
131*2dce792eSToby Isaac . s     - The space number
132*2dce792eSToby Isaac - subsp - the number of spaces
133*2dce792eSToby Isaac 
134*2dce792eSToby Isaac   Level: intermediate
135*2dce792eSToby Isaac 
136*2dce792eSToby Isaac   Note:
137*2dce792eSToby Isaac   The name SetSubspace is slightly misleading because it is actually setting one of the defining spaces of the sum, not a Subspace of it
138*2dce792eSToby Isaac 
139*2dce792eSToby Isaac .seealso: `PETSCDUALSPACESUM`, `PetscDualSpace`, `PetscDualSpaceSumGetSubspace()`
140*2dce792eSToby Isaac @*/
141*2dce792eSToby Isaac PetscErrorCode PetscDualSpaceSumSetSubspace(PetscDualSpace sp, PetscInt s, PetscDualSpace subsp)
142*2dce792eSToby Isaac {
143*2dce792eSToby Isaac   PetscFunctionBegin;
144*2dce792eSToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
145*2dce792eSToby Isaac   if (subsp) PetscValidHeaderSpecific(subsp, PETSCDUALSPACE_CLASSID, 3);
146*2dce792eSToby Isaac   PetscTryMethod(sp, "PetscDualSpaceSumSetSubspace_C", (PetscDualSpace, PetscInt, PetscDualSpace), (sp, s, subsp));
147*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
148*2dce792eSToby Isaac }
149*2dce792eSToby Isaac 
150*2dce792eSToby Isaac static PetscErrorCode PetscDualSpaceSumGetNumSubspaces_Sum(PetscDualSpace space, PetscInt *numSumSpaces)
151*2dce792eSToby Isaac {
152*2dce792eSToby Isaac   PetscDualSpace_Sum *sum = (PetscDualSpace_Sum *)space->data;
153*2dce792eSToby Isaac 
154*2dce792eSToby Isaac   PetscFunctionBegin;
155*2dce792eSToby Isaac   *numSumSpaces = sum->numSumSpaces;
156*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
157*2dce792eSToby Isaac }
158*2dce792eSToby Isaac 
159*2dce792eSToby Isaac static PetscErrorCode PetscDualSpaceSumSetNumSubspaces_Sum(PetscDualSpace space, PetscInt numSumSpaces)
160*2dce792eSToby Isaac {
161*2dce792eSToby Isaac   PetscDualSpace_Sum *sum = (PetscDualSpace_Sum *)space->data;
162*2dce792eSToby Isaac   PetscInt            Ns  = sum->numSumSpaces;
163*2dce792eSToby Isaac 
164*2dce792eSToby Isaac   PetscFunctionBegin;
165*2dce792eSToby Isaac   PetscCheck(!sum->setupCalled, PetscObjectComm((PetscObject)space), PETSC_ERR_ARG_WRONGSTATE, "Cannot change number of subspaces after setup called");
166*2dce792eSToby Isaac   if (numSumSpaces == Ns) PetscFunctionReturn(PETSC_SUCCESS);
167*2dce792eSToby Isaac   for (PetscInt s = 0; s < Ns; ++s) PetscCall(PetscDualSpaceDestroy(&sum->sumspaces[s]));
168*2dce792eSToby Isaac   PetscCall(PetscFree(sum->sumspaces));
169*2dce792eSToby Isaac 
170*2dce792eSToby Isaac   Ns = sum->numSumSpaces = numSumSpaces;
171*2dce792eSToby Isaac   PetscCall(PetscCalloc1(Ns, &sum->sumspaces));
172*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
173*2dce792eSToby Isaac }
174*2dce792eSToby Isaac 
175*2dce792eSToby Isaac static PetscErrorCode PetscDualSpaceSumGetConcatenate_Sum(PetscDualSpace sp, PetscBool *concatenate)
176*2dce792eSToby Isaac {
177*2dce792eSToby Isaac   PetscDualSpace_Sum *sum = (PetscDualSpace_Sum *)sp->data;
178*2dce792eSToby Isaac 
179*2dce792eSToby Isaac   PetscFunctionBegin;
180*2dce792eSToby Isaac   *concatenate = sum->concatenate;
181*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
182*2dce792eSToby Isaac }
183*2dce792eSToby Isaac 
184*2dce792eSToby Isaac static PetscErrorCode PetscDualSpaceSumSetConcatenate_Sum(PetscDualSpace sp, PetscBool concatenate)
185*2dce792eSToby Isaac {
186*2dce792eSToby Isaac   PetscDualSpace_Sum *sum = (PetscDualSpace_Sum *)sp->data;
187*2dce792eSToby Isaac 
188*2dce792eSToby Isaac   PetscFunctionBegin;
189*2dce792eSToby Isaac   PetscCheck(!sum->setupCalled, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change space concatenation after setup called.");
190*2dce792eSToby Isaac 
191*2dce792eSToby Isaac   sum->concatenate = concatenate;
192*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
193*2dce792eSToby Isaac }
194*2dce792eSToby Isaac 
195*2dce792eSToby Isaac static PetscErrorCode PetscDualSpaceSumGetSubspace_Sum(PetscDualSpace space, PetscInt s, PetscDualSpace *subspace)
196*2dce792eSToby Isaac {
197*2dce792eSToby Isaac   PetscDualSpace_Sum *sum = (PetscDualSpace_Sum *)space->data;
198*2dce792eSToby Isaac   PetscInt            Ns  = sum->numSumSpaces;
199*2dce792eSToby Isaac 
200*2dce792eSToby Isaac   PetscFunctionBegin;
201*2dce792eSToby Isaac   PetscCheck(Ns >= 0, PetscObjectComm((PetscObject)space), PETSC_ERR_ARG_WRONGSTATE, "Must call PetscDualSpaceSumSetNumSubspaces() first");
202*2dce792eSToby Isaac   PetscCheck(s >= 0 && s < Ns, PetscObjectComm((PetscObject)space), PETSC_ERR_ARG_OUTOFRANGE, "Invalid subspace number %" PetscInt_FMT, s);
203*2dce792eSToby Isaac 
204*2dce792eSToby Isaac   *subspace = sum->sumspaces[s];
205*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
206*2dce792eSToby Isaac }
207*2dce792eSToby Isaac 
208*2dce792eSToby Isaac static PetscErrorCode PetscDualSpaceSumSetSubspace_Sum(PetscDualSpace space, PetscInt s, PetscDualSpace subspace)
209*2dce792eSToby Isaac {
210*2dce792eSToby Isaac   PetscDualSpace_Sum *sum = (PetscDualSpace_Sum *)space->data;
211*2dce792eSToby Isaac   PetscInt            Ns  = sum->numSumSpaces;
212*2dce792eSToby Isaac 
213*2dce792eSToby Isaac   PetscFunctionBegin;
214*2dce792eSToby Isaac   PetscCheck(!sum->setupCalled, PetscObjectComm((PetscObject)space), PETSC_ERR_ARG_WRONGSTATE, "Cannot change subspace after setup called");
215*2dce792eSToby Isaac   PetscCheck(Ns >= 0, PetscObjectComm((PetscObject)space), PETSC_ERR_ARG_WRONGSTATE, "Must call PetscDualSpaceSumSetNumSubspaces() first");
216*2dce792eSToby Isaac   PetscCheck(s >= 0 && s < Ns, PetscObjectComm((PetscObject)space), PETSC_ERR_ARG_OUTOFRANGE, "Invalid subspace number %" PetscInt_FMT, s);
217*2dce792eSToby Isaac 
218*2dce792eSToby Isaac   PetscCall(PetscObjectReference((PetscObject)subspace));
219*2dce792eSToby Isaac   PetscCall(PetscDualSpaceDestroy(&sum->sumspaces[s]));
220*2dce792eSToby Isaac   sum->sumspaces[s] = subspace;
221*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
222*2dce792eSToby Isaac }
223*2dce792eSToby Isaac 
224*2dce792eSToby Isaac static PetscErrorCode PetscDualSpaceDuplicate_Sum(PetscDualSpace sp, PetscDualSpace spNew)
225*2dce792eSToby Isaac {
226*2dce792eSToby Isaac   PetscInt       num_subspaces, Nc;
227*2dce792eSToby Isaac   PetscBool      concatenate, interleave_basis, interleave_components;
228*2dce792eSToby Isaac   PetscDualSpace subsp_first     = NULL;
229*2dce792eSToby Isaac   PetscDualSpace subsp_dup_first = NULL;
230*2dce792eSToby Isaac   DM             K;
231*2dce792eSToby Isaac 
232*2dce792eSToby Isaac   PetscFunctionBegin;
233*2dce792eSToby Isaac   PetscCall(PetscDualSpaceSumGetNumSubspaces(sp, &num_subspaces));
234*2dce792eSToby Isaac   PetscCall(PetscDualSpaceSumSetNumSubspaces(spNew, num_subspaces));
235*2dce792eSToby Isaac   PetscCall(PetscDualSpaceSumGetConcatenate(sp, &concatenate));
236*2dce792eSToby Isaac   PetscCall(PetscDualSpaceSumSetConcatenate(spNew, concatenate));
237*2dce792eSToby Isaac   PetscCall(PetscDualSpaceSumGetInterleave(sp, &interleave_basis, &interleave_components));
238*2dce792eSToby Isaac   PetscCall(PetscDualSpaceSumSetInterleave(spNew, interleave_basis, interleave_components));
239*2dce792eSToby Isaac   PetscCall(PetscDualSpaceGetDM(sp, &K));
240*2dce792eSToby Isaac   PetscCall(PetscDualSpaceSetDM(spNew, K));
241*2dce792eSToby Isaac   PetscCall(PetscDualSpaceGetNumComponents(sp, &Nc));
242*2dce792eSToby Isaac   PetscCall(PetscDualSpaceSetNumComponents(spNew, Nc));
243*2dce792eSToby Isaac   for (PetscInt s = 0; s < num_subspaces; s++) {
244*2dce792eSToby Isaac     PetscDualSpace subsp, subspNew;
245*2dce792eSToby Isaac 
246*2dce792eSToby Isaac     PetscCall(PetscDualSpaceSumGetSubspace(sp, s, &subsp));
247*2dce792eSToby Isaac     if (s == 0) {
248*2dce792eSToby Isaac       subsp_first = subsp;
249*2dce792eSToby Isaac       PetscCall(PetscDualSpaceDuplicate(subsp, &subsp_dup_first));
250*2dce792eSToby Isaac       subspNew = subsp_dup_first;
251*2dce792eSToby Isaac     } else if (subsp == subsp_first) {
252*2dce792eSToby Isaac       PetscCall(PetscObjectReference((PetscObject)subsp_dup_first));
253*2dce792eSToby Isaac       subspNew = subsp_dup_first;
254*2dce792eSToby Isaac     } else {
255*2dce792eSToby Isaac       PetscCall(PetscDualSpaceDuplicate(subsp, &subspNew));
256*2dce792eSToby Isaac     }
257*2dce792eSToby Isaac     PetscCall(PetscDualSpaceSumSetSubspace(spNew, s, subspNew));
258*2dce792eSToby Isaac     PetscCall(PetscDualSpaceDestroy(&subspNew));
259*2dce792eSToby Isaac   }
260*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
261*2dce792eSToby Isaac }
262*2dce792eSToby Isaac 
263*2dce792eSToby Isaac static PetscErrorCode PetscDualSpaceSumCreateQuadrature(PetscInt Ns, PetscInt cdim, PetscBool uniform_points, PetscQuadrature subquads[], PetscQuadrature *fullquad)
264*2dce792eSToby Isaac {
265*2dce792eSToby Isaac   PetscReal *points;
266*2dce792eSToby Isaac   PetscInt   Npoints;
267*2dce792eSToby Isaac 
268*2dce792eSToby Isaac   PetscFunctionBegin;
269*2dce792eSToby Isaac   if (uniform_points) {
270*2dce792eSToby Isaac     PetscCall(PetscObjectReference((PetscObject)subquads[0]));
271*2dce792eSToby Isaac     *fullquad = subquads[0];
272*2dce792eSToby Isaac     PetscFunctionReturn(PETSC_SUCCESS);
273*2dce792eSToby Isaac   }
274*2dce792eSToby Isaac   Npoints = 0;
275*2dce792eSToby Isaac   for (PetscInt s = 0; s < Ns; s++) {
276*2dce792eSToby Isaac     PetscInt subNpoints;
277*2dce792eSToby Isaac 
278*2dce792eSToby Isaac     if (!subquads[s]) continue;
279*2dce792eSToby Isaac     PetscCall(PetscQuadratureGetData(subquads[s], NULL, NULL, &subNpoints, NULL, NULL));
280*2dce792eSToby Isaac     Npoints += subNpoints;
281*2dce792eSToby Isaac   }
282*2dce792eSToby Isaac   *fullquad = NULL;
283*2dce792eSToby Isaac   if (!Npoints) PetscFunctionReturn(PETSC_SUCCESS);
284*2dce792eSToby Isaac   PetscCall(PetscMalloc1(Npoints * cdim, &points));
285*2dce792eSToby Isaac   for (PetscInt s = 0, offset = 0; s < Ns; s++) {
286*2dce792eSToby Isaac     PetscInt         subNpoints;
287*2dce792eSToby Isaac     const PetscReal *subpoints;
288*2dce792eSToby Isaac 
289*2dce792eSToby Isaac     PetscCall(PetscQuadratureGetData(subquads[s], NULL, NULL, &subNpoints, &subpoints, NULL));
290*2dce792eSToby Isaac     PetscCall(PetscArraycpy(&points[offset], subpoints, cdim * subNpoints));
291*2dce792eSToby Isaac     offset += cdim * subNpoints;
292*2dce792eSToby Isaac   }
293*2dce792eSToby Isaac   PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, fullquad));
294*2dce792eSToby Isaac   PetscCall(PetscQuadratureSetData(*fullquad, cdim, 0, Npoints, points, NULL));
295*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
296*2dce792eSToby Isaac }
297*2dce792eSToby Isaac 
298*2dce792eSToby Isaac static PetscErrorCode PetscDualSpaceSumCreateMatrix(PetscDualSpace sp, Mat submats[], ISLocalToGlobalMapping map_rows[], ISLocalToGlobalMapping map_cols[], PetscQuadrature fullquad, Mat *fullmat)
299*2dce792eSToby Isaac {
300*2dce792eSToby Isaac   Mat          mat;
301*2dce792eSToby Isaac   PetscInt    *i = NULL, *j = NULL;
302*2dce792eSToby Isaac   PetscScalar *v = NULL;
303*2dce792eSToby Isaac   PetscInt     npoints;
304*2dce792eSToby Isaac   PetscInt     nrows, ncols, nnz;
305*2dce792eSToby Isaac   PetscInt     Nc, num_subspaces;
306*2dce792eSToby Isaac 
307*2dce792eSToby Isaac   PetscFunctionBegin;
308*2dce792eSToby Isaac   PetscCall(PetscDualSpaceGetNumComponents(sp, &Nc));
309*2dce792eSToby Isaac   PetscCall(PetscDualSpaceSumGetNumSubspaces(sp, &num_subspaces));
310*2dce792eSToby Isaac   nrows = 0;
311*2dce792eSToby Isaac   ncols = 0;
312*2dce792eSToby Isaac   nnz   = 0;
313*2dce792eSToby Isaac   PetscCall(PetscQuadratureGetData(fullquad, NULL, NULL, &npoints, NULL, NULL));
314*2dce792eSToby Isaac   ncols = Nc * npoints;
315*2dce792eSToby Isaac   for (PetscInt s = 0; s < num_subspaces; s++) {
316*2dce792eSToby Isaac     // Get the COO data for each matrix, map the is and js, and append to growing COO data
317*2dce792eSToby Isaac     PetscInt        sNrows, sNcols;
318*2dce792eSToby Isaac     Mat             smat;
319*2dce792eSToby Isaac     const PetscInt *si;
320*2dce792eSToby Isaac     const PetscInt *sj;
321*2dce792eSToby Isaac     PetscScalar    *sv;
322*2dce792eSToby Isaac     PetscMemType    memtype;
323*2dce792eSToby Isaac     PetscInt        snz;
324*2dce792eSToby Isaac     PetscInt        snz_actual;
325*2dce792eSToby Isaac     PetscInt       *cooi;
326*2dce792eSToby Isaac     PetscInt       *cooj;
327*2dce792eSToby Isaac     PetscScalar    *coov;
328*2dce792eSToby Isaac     PetscScalar    *v_new;
329*2dce792eSToby Isaac     PetscInt       *i_new;
330*2dce792eSToby Isaac     PetscInt       *j_new;
331*2dce792eSToby Isaac 
332*2dce792eSToby Isaac     if (!submats[s]) continue;
333*2dce792eSToby Isaac     PetscCall(MatGetSize(submats[s], &sNrows, &sNcols));
334*2dce792eSToby Isaac     nrows += sNrows;
335*2dce792eSToby Isaac     PetscCall(MatConvert(submats[s], MATSEQAIJ, MAT_INITIAL_MATRIX, &smat));
336*2dce792eSToby Isaac     PetscCall(MatBindToCPU(smat, PETSC_TRUE));
337*2dce792eSToby Isaac     PetscCall(MatSeqAIJGetCSRAndMemType(smat, &si, &sj, &sv, &memtype));
338*2dce792eSToby Isaac     PetscCheck(memtype == PETSC_MEMTYPE_HOST, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not convert matrix to host memory");
339*2dce792eSToby Isaac     snz = si[sNrows];
340*2dce792eSToby Isaac 
341*2dce792eSToby Isaac     PetscCall(PetscMalloc1(nnz + snz, &v_new));
342*2dce792eSToby Isaac     PetscCall(PetscArraycpy(v_new, v, nnz));
343*2dce792eSToby Isaac     PetscCall(PetscFree(v));
344*2dce792eSToby Isaac     v = v_new;
345*2dce792eSToby Isaac 
346*2dce792eSToby Isaac     PetscCall(PetscMalloc1(nnz + snz, &i_new));
347*2dce792eSToby Isaac     PetscCall(PetscArraycpy(i_new, i, nnz));
348*2dce792eSToby Isaac     PetscCall(PetscFree(i));
349*2dce792eSToby Isaac     i = i_new;
350*2dce792eSToby Isaac 
351*2dce792eSToby Isaac     PetscCall(PetscMalloc1(nnz + snz, &j_new));
352*2dce792eSToby Isaac     PetscCall(PetscArraycpy(j_new, j, nnz));
353*2dce792eSToby Isaac     PetscCall(PetscFree(j));
354*2dce792eSToby Isaac     j = j_new;
355*2dce792eSToby Isaac 
356*2dce792eSToby Isaac     PetscCall(PetscMalloc2(snz, &cooi, snz, &cooj));
357*2dce792eSToby Isaac 
358*2dce792eSToby Isaac     coov = &v[nnz];
359*2dce792eSToby Isaac 
360*2dce792eSToby Isaac     snz_actual = 0;
361*2dce792eSToby Isaac     for (PetscInt row = 0; row < sNrows; row++) {
362*2dce792eSToby Isaac       for (PetscInt k = si[row]; k < si[row + 1]; k++, snz_actual++) {
363*2dce792eSToby Isaac         cooi[snz_actual] = row;
364*2dce792eSToby Isaac         cooj[snz_actual] = sj[k];
365*2dce792eSToby Isaac         coov[snz_actual] = sv[k];
366*2dce792eSToby Isaac       }
367*2dce792eSToby Isaac     }
368*2dce792eSToby Isaac     PetscCall(MatDestroy(&smat));
369*2dce792eSToby Isaac 
370*2dce792eSToby Isaac     if (snz_actual > 0) {
371*2dce792eSToby Isaac       PetscCall(ISLocalToGlobalMappingApply(map_rows[s], snz_actual, cooi, &i[nnz]));
372*2dce792eSToby Isaac       PetscCall(ISLocalToGlobalMappingApply(map_cols[s], snz_actual, cooj, &j[nnz]));
373*2dce792eSToby Isaac       nnz += snz_actual;
374*2dce792eSToby Isaac     }
375*2dce792eSToby Isaac     PetscCall(PetscFree2(cooi, cooj));
376*2dce792eSToby Isaac   }
377*2dce792eSToby Isaac   PetscCall(MatCreate(PETSC_COMM_SELF, fullmat));
378*2dce792eSToby Isaac   mat = *fullmat;
379*2dce792eSToby Isaac   PetscCall(MatSetSizes(mat, nrows, ncols, nrows, ncols));
380*2dce792eSToby Isaac   PetscCall(MatSetType(mat, MATSEQAIJ));
381*2dce792eSToby Isaac   PetscCall(MatSetPreallocationCOO(mat, nnz, i, j));
382*2dce792eSToby Isaac   PetscCall(MatSetValuesCOO(mat, v, INSERT_VALUES));
383*2dce792eSToby Isaac   PetscCall(PetscFree(i));
384*2dce792eSToby Isaac   PetscCall(PetscFree(j));
385*2dce792eSToby Isaac   PetscCall(PetscFree(v));
386*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
387*2dce792eSToby Isaac }
388*2dce792eSToby Isaac 
389*2dce792eSToby Isaac static PetscErrorCode PetscDualSpaceSumCreateMappings(PetscDualSpace sp, PetscBool interior, PetscBool uniform_points, ISLocalToGlobalMapping map_row[], ISLocalToGlobalMapping map_col[])
390*2dce792eSToby Isaac {
391*2dce792eSToby Isaac   PetscSection section;
392*2dce792eSToby Isaac   PetscInt     pStart, pEnd, Ns, Nc;
393*2dce792eSToby Isaac   PetscInt     num_points = 0;
394*2dce792eSToby Isaac   PetscBool    interleave_basis, interleave_components, concatenate;
395*2dce792eSToby Isaac   PetscInt    *roffset;
396*2dce792eSToby Isaac 
397*2dce792eSToby Isaac   PetscFunctionBegin;
398*2dce792eSToby Isaac   if (interior) {
399*2dce792eSToby Isaac     PetscCall(PetscDualSpaceGetInteriorSection(sp, &section));
400*2dce792eSToby Isaac   } else {
401*2dce792eSToby Isaac     PetscCall(PetscDualSpaceGetSection(sp, &section));
402*2dce792eSToby Isaac   }
403*2dce792eSToby Isaac   PetscCall(PetscDualSpaceGetNumComponents(sp, &Nc));
404*2dce792eSToby Isaac   PetscCall(PetscDualSpaceSumGetNumSubspaces(sp, &Ns));
405*2dce792eSToby Isaac   PetscCall(PetscDualSpaceSumGetInterleave(sp, &interleave_basis, &interleave_components));
406*2dce792eSToby Isaac   PetscCall(PetscDualSpaceSumGetConcatenate(sp, &concatenate));
407*2dce792eSToby Isaac   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
408*2dce792eSToby Isaac   for (PetscInt p = pStart; p < pEnd; p++) {
409*2dce792eSToby Isaac     PetscInt dof;
410*2dce792eSToby Isaac 
411*2dce792eSToby Isaac     PetscCall(PetscSectionGetDof(section, p, &dof));
412*2dce792eSToby Isaac     num_points += (dof > 0);
413*2dce792eSToby Isaac   }
414*2dce792eSToby Isaac   PetscCall(PetscCalloc1(pEnd - pStart, &roffset));
415*2dce792eSToby Isaac   for (PetscInt s = 0, coffset = 0; s < Ns; s++) {
416*2dce792eSToby Isaac     PetscDualSpace  subsp;
417*2dce792eSToby Isaac     PetscSection    subsection;
418*2dce792eSToby Isaac     IS              is_row, is_col;
419*2dce792eSToby Isaac     PetscInt        subNb;
420*2dce792eSToby Isaac     PetscInt        sNc, sNpoints, sNcols;
421*2dce792eSToby Isaac     PetscQuadrature quad;
422*2dce792eSToby Isaac 
423*2dce792eSToby Isaac     PetscCall(PetscDualSpaceSumGetSubspace(sp, s, &subsp));
424*2dce792eSToby Isaac     PetscCall(PetscDualSpaceGetNumComponents(subsp, &sNc));
425*2dce792eSToby Isaac     if (interior) {
426*2dce792eSToby Isaac       PetscCall(PetscDualSpaceGetInteriorSection(subsp, &subsection));
427*2dce792eSToby Isaac       PetscCall(PetscDualSpaceGetInteriorData(subsp, &quad, NULL));
428*2dce792eSToby Isaac     } else {
429*2dce792eSToby Isaac       PetscCall(PetscDualSpaceGetSection(subsp, &subsection));
430*2dce792eSToby Isaac       PetscCall(PetscDualSpaceGetAllData(subsp, &quad, NULL));
431*2dce792eSToby Isaac     }
432*2dce792eSToby Isaac     PetscCall(PetscSectionGetStorageSize(subsection, &subNb));
433*2dce792eSToby Isaac     if (num_points == 1) {
434*2dce792eSToby Isaac       PetscInt rstride;
435*2dce792eSToby Isaac 
436*2dce792eSToby Isaac       rstride = interleave_basis ? Ns : 1;
437*2dce792eSToby Isaac       PetscCall(ISCreateStride(PETSC_COMM_SELF, subNb, roffset[0], rstride, &is_row));
438*2dce792eSToby Isaac       roffset[0] += interleave_basis ? 1 : subNb;
439*2dce792eSToby Isaac     } else {
440*2dce792eSToby Isaac       PetscInt *rows;
441*2dce792eSToby Isaac 
442*2dce792eSToby Isaac       PetscCall(PetscMalloc1(subNb, &rows));
443*2dce792eSToby Isaac       for (PetscInt p = pStart; p < pEnd; p++) {
444*2dce792eSToby Isaac         PetscInt subdof, dof, off, suboff, stride;
445*2dce792eSToby Isaac 
446*2dce792eSToby Isaac         PetscCall(PetscSectionGetOffset(subsection, p, &suboff));
447*2dce792eSToby Isaac         PetscCall(PetscSectionGetDof(subsection, p, &subdof));
448*2dce792eSToby Isaac         PetscCall(PetscSectionGetOffset(section, p, &off));
449*2dce792eSToby Isaac         PetscCall(PetscSectionGetDof(section, p, &dof));
450*2dce792eSToby Isaac         PetscCheck(subdof * Ns == dof || !interleave_basis, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Basis cannot be interleaved\n");
451*2dce792eSToby Isaac         stride = interleave_basis ? Ns : 1;
452*2dce792eSToby Isaac         for (PetscInt k = 0; k < subdof; k++) { rows[suboff + k] = off + roffset[p - pStart] + k * stride; }
453*2dce792eSToby Isaac         roffset[p - pStart] += interleave_basis ? 1 : subdof;
454*2dce792eSToby Isaac       }
455*2dce792eSToby Isaac       PetscCall(ISCreateGeneral(PETSC_COMM_SELF, subNb, rows, PETSC_OWN_POINTER, &is_row));
456*2dce792eSToby Isaac     }
457*2dce792eSToby Isaac 
458*2dce792eSToby Isaac     sNpoints = 0;
459*2dce792eSToby Isaac     if (quad) PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &sNpoints, NULL, NULL));
460*2dce792eSToby Isaac     sNcols = sNpoints * sNc;
461*2dce792eSToby Isaac 
462*2dce792eSToby Isaac     if (!concatenate) {
463*2dce792eSToby Isaac       PetscCall(ISCreateStride(PETSC_COMM_SELF, sNcols, coffset, 1, &is_col));
464*2dce792eSToby Isaac       coffset += uniform_points ? 0 : sNcols;
465*2dce792eSToby Isaac     } else {
466*2dce792eSToby Isaac       if (uniform_points && interleave_components) {
467*2dce792eSToby Isaac         PetscCall(ISCreateStride(PETSC_COMM_SELF, sNcols, coffset, Ns, &is_col));
468*2dce792eSToby Isaac         coffset += 1;
469*2dce792eSToby Isaac       } else {
470*2dce792eSToby Isaac         PetscInt *cols;
471*2dce792eSToby Isaac 
472*2dce792eSToby Isaac         PetscCall(PetscMalloc1(sNcols, &cols));
473*2dce792eSToby Isaac         for (PetscInt p = 0, r = 0; p < sNpoints; p++) {
474*2dce792eSToby Isaac           for (PetscInt c = 0; c < sNc; c++, r++) { cols[r] = coffset + p * Nc + c; }
475*2dce792eSToby Isaac         }
476*2dce792eSToby Isaac         coffset += uniform_points ? sNc : Nc * sNpoints + sNc;
477*2dce792eSToby Isaac         PetscCall(ISCreateGeneral(PETSC_COMM_SELF, sNcols, cols, PETSC_OWN_POINTER, &is_col));
478*2dce792eSToby Isaac       }
479*2dce792eSToby Isaac     }
480*2dce792eSToby Isaac     PetscCall(ISLocalToGlobalMappingCreateIS(is_row, &map_row[s]));
481*2dce792eSToby Isaac     PetscCall(ISLocalToGlobalMappingCreateIS(is_col, &map_col[s]));
482*2dce792eSToby Isaac     PetscCall(ISDestroy(&is_row));
483*2dce792eSToby Isaac     PetscCall(ISDestroy(&is_col));
484*2dce792eSToby Isaac   }
485*2dce792eSToby Isaac   PetscCall(PetscFree(roffset));
486*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
487*2dce792eSToby Isaac }
488*2dce792eSToby Isaac 
489*2dce792eSToby Isaac static PetscErrorCode PetscDualSpaceCreateFacetSubspace_Sum(PetscDualSpace sp, PetscInt f, PetscDualSpace *bdsp)
490*2dce792eSToby Isaac {
491*2dce792eSToby Isaac   PetscInt       k, Nc, Nk, Nknew, Ncnew, Ns;
492*2dce792eSToby Isaac   PetscInt       dim, pointDim = -1;
493*2dce792eSToby Isaac   PetscInt       depth, Ncopies;
494*2dce792eSToby Isaac   PetscBool      any;
495*2dce792eSToby Isaac   DM             dm, K;
496*2dce792eSToby Isaac   DMPolytopeType ct;
497*2dce792eSToby Isaac 
498*2dce792eSToby Isaac   PetscFunctionBegin;
499*2dce792eSToby Isaac   PetscCall(PetscDualSpaceSumGetNumSubspaces(sp, &Ns));
500*2dce792eSToby Isaac   any = PETSC_FALSE;
501*2dce792eSToby Isaac   for (PetscInt s = 0; s < Ns; s++) {
502*2dce792eSToby Isaac     PetscDualSpace subsp, subpointsp;
503*2dce792eSToby Isaac 
504*2dce792eSToby Isaac     PetscCall(PetscDualSpaceSumGetSubspace(sp, s, &subsp));
505*2dce792eSToby Isaac     PetscCall(PetscDualSpaceGetPointSubspace(subsp, f, &subpointsp));
506*2dce792eSToby Isaac     if (subpointsp) any = PETSC_TRUE;
507*2dce792eSToby Isaac   }
508*2dce792eSToby Isaac   if (!any) {
509*2dce792eSToby Isaac     *bdsp = NULL;
510*2dce792eSToby Isaac     PetscFunctionReturn(PETSC_SUCCESS);
511*2dce792eSToby Isaac   }
512*2dce792eSToby Isaac   PetscCall(PetscDualSpaceGetDM(sp, &dm));
513*2dce792eSToby Isaac   PetscCall(DMGetDimension(dm, &dim));
514*2dce792eSToby Isaac   PetscCall(DMPlexGetDepth(dm, &depth));
515*2dce792eSToby Isaac   PetscCall(PetscDualSpaceDuplicate(sp, bdsp));
516*2dce792eSToby Isaac   PetscCheck((depth == dim) || (depth == 1), PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unsupported interpolation state of reference element");
517*2dce792eSToby Isaac   pointDim = (depth == dim) ? (dim - 1) : 0;
518*2dce792eSToby Isaac   PetscCall(DMPlexGetCellType(dm, f, &ct));
519*2dce792eSToby Isaac   PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, ct, &K));
520*2dce792eSToby Isaac   PetscCall(PetscDualSpaceSetDM(*bdsp, K));
521*2dce792eSToby Isaac   PetscCall(PetscDualSpaceGetFormDegree(sp, &k));
522*2dce792eSToby Isaac   PetscCall(PetscDualSpaceGetNumComponents(sp, &Nc));
523*2dce792eSToby Isaac   PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(k), &Nk));
524*2dce792eSToby Isaac   Ncopies = Nc / Nk;
525*2dce792eSToby Isaac   PetscCall(PetscDTBinomialInt(pointDim, PetscAbsInt(k), &Nknew));
526*2dce792eSToby Isaac   Ncnew = Nknew * Ncopies;
527*2dce792eSToby Isaac   PetscCall(PetscDualSpaceSetNumComponents(*bdsp, Ncnew));
528*2dce792eSToby Isaac   for (PetscInt s = 0; s < Ns; s++) {
529*2dce792eSToby Isaac     PetscDualSpace subsp, subpointsp;
530*2dce792eSToby Isaac 
531*2dce792eSToby Isaac     PetscCall(PetscDualSpaceSumGetSubspace(sp, s, &subsp));
532*2dce792eSToby Isaac     PetscCall(PetscDualSpaceGetPointSubspace(subsp, f, &subpointsp));
533*2dce792eSToby Isaac     if (subpointsp) {
534*2dce792eSToby Isaac       PetscCall(PetscObjectReference((PetscObject)subpointsp));
535*2dce792eSToby Isaac     } else {
536*2dce792eSToby Isaac       // make an empty dual space
537*2dce792eSToby Isaac       PetscInt subNc, subNcopies, subpointNc;
538*2dce792eSToby Isaac 
539*2dce792eSToby Isaac       PetscCall(PetscDualSpaceCreate(PetscObjectComm((PetscObject)subsp), &subpointsp));
540*2dce792eSToby Isaac       PetscCall(PetscDualSpaceGetNumComponents(subsp, &subNc));
541*2dce792eSToby Isaac       subNcopies = subNc / Nk;
542*2dce792eSToby Isaac       subpointNc = subNcopies * Nknew;
543*2dce792eSToby Isaac       PetscCall(PetscDualSpaceSetType(subpointsp, PETSCDUALSPACESIMPLE));
544*2dce792eSToby Isaac       PetscCall(PetscDualSpaceSimpleSetDimension(subpointsp, 0));
545*2dce792eSToby Isaac       PetscCall(PetscDualSpaceSetFormDegree(subpointsp, k));
546*2dce792eSToby Isaac       PetscCall(PetscDualSpaceSetNumComponents(subpointsp, subpointNc));
547*2dce792eSToby Isaac     }
548*2dce792eSToby Isaac     // K should be equal to subpointsp->dm, but we want it to be exactly the same
549*2dce792eSToby Isaac     PetscCall(PetscObjectReference((PetscObject)K));
550*2dce792eSToby Isaac     PetscCall(DMDestroy(&subpointsp->dm));
551*2dce792eSToby Isaac     subpointsp->dm = K;
552*2dce792eSToby Isaac     PetscCall(PetscDualSpaceSetUp(subpointsp));
553*2dce792eSToby Isaac     PetscCall(PetscDualSpaceSumSetSubspace(*bdsp, s, subpointsp));
554*2dce792eSToby Isaac     PetscCall(PetscDualSpaceDestroy(&subpointsp));
555*2dce792eSToby Isaac   }
556*2dce792eSToby Isaac   PetscCall(DMDestroy(&K));
557*2dce792eSToby Isaac   PetscCall(PetscDualSpaceSetUp(*bdsp));
558*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
559*2dce792eSToby Isaac }
560*2dce792eSToby Isaac 
561*2dce792eSToby Isaac static PetscErrorCode PetscDualSpaceSumIsUniform(PetscDualSpace sp, PetscBool *is_uniform)
562*2dce792eSToby Isaac {
563*2dce792eSToby Isaac   PetscDualSpace_Sum *sum     = (PetscDualSpace_Sum *)sp->data;
564*2dce792eSToby Isaac   PetscBool           uniform = PETSC_TRUE;
565*2dce792eSToby Isaac 
566*2dce792eSToby Isaac   PetscFunctionBegin;
567*2dce792eSToby Isaac   for (PetscInt s = 1; s < sum->numSumSpaces; s++) {
568*2dce792eSToby Isaac     if (sum->sumspaces[s] != sum->sumspaces[0]) {
569*2dce792eSToby Isaac       uniform = PETSC_FALSE;
570*2dce792eSToby Isaac       break;
571*2dce792eSToby Isaac     }
572*2dce792eSToby Isaac   }
573*2dce792eSToby Isaac   *is_uniform = uniform;
574*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
575*2dce792eSToby Isaac }
576*2dce792eSToby Isaac 
577*2dce792eSToby Isaac static PetscErrorCode PetscDualSpaceGetSymmetries_Sum(PetscDualSpace sp, const PetscInt ****perms, const PetscScalar ****flips)
578*2dce792eSToby Isaac {
579*2dce792eSToby Isaac   PetscDualSpace_Sum *sum = (PetscDualSpace_Sum *)sp->data;
580*2dce792eSToby Isaac 
581*2dce792eSToby Isaac   PetscFunctionBegin;
582*2dce792eSToby Isaac   if (!sum->symComputed) {
583*2dce792eSToby Isaac     PetscInt       Ns;
584*2dce792eSToby Isaac     PetscBool      any_perms = PETSC_FALSE;
585*2dce792eSToby Isaac     PetscBool      any_flips = PETSC_FALSE;
586*2dce792eSToby Isaac     PetscInt    ***symperms  = NULL;
587*2dce792eSToby Isaac     PetscScalar ***symflips  = NULL;
588*2dce792eSToby Isaac 
589*2dce792eSToby Isaac     sum->symComputed = PETSC_TRUE;
590*2dce792eSToby Isaac     PetscCall(PetscDualSpaceSumGetNumSubspaces(sp, &Ns));
591*2dce792eSToby Isaac     for (PetscInt s = 0; s < Ns; s++) {
592*2dce792eSToby Isaac       PetscDualSpace       subsp;
593*2dce792eSToby Isaac       const PetscInt    ***sub_perms;
594*2dce792eSToby Isaac       const PetscScalar ***sub_flips;
595*2dce792eSToby Isaac 
596*2dce792eSToby Isaac       PetscCall(PetscDualSpaceSumGetSubspace(sp, s, &subsp));
597*2dce792eSToby Isaac       PetscCall(PetscDualSpaceGetSymmetries(subsp, &sub_perms, &sub_flips));
598*2dce792eSToby Isaac       if (sub_perms) any_perms = PETSC_TRUE;
599*2dce792eSToby Isaac       if (sub_flips) any_flips = PETSC_TRUE;
600*2dce792eSToby Isaac     }
601*2dce792eSToby Isaac     if (any_perms || any_flips) {
602*2dce792eSToby Isaac       DM       K;
603*2dce792eSToby Isaac       PetscInt pStart, pEnd, numPoints;
604*2dce792eSToby Isaac       PetscInt spintdim;
605*2dce792eSToby Isaac 
606*2dce792eSToby Isaac       PetscCall(PetscDualSpaceGetDM(sp, &K));
607*2dce792eSToby Isaac       PetscCall(DMPlexGetChart(K, &pStart, &pEnd));
608*2dce792eSToby Isaac       numPoints = pEnd - pStart;
609*2dce792eSToby Isaac       PetscCall(PetscCalloc1(numPoints, &symperms));
610*2dce792eSToby Isaac       PetscCall(PetscCalloc1(numPoints, &symflips));
611*2dce792eSToby Isaac       PetscCall(PetscDualSpaceGetBoundarySymmetries_Internal(sp, symperms, symflips));
612*2dce792eSToby Isaac       // get interior symmetries
613*2dce792eSToby Isaac       PetscCall(PetscDualSpaceGetInteriorDimension(sp, &spintdim));
614*2dce792eSToby Isaac       if (spintdim) {
615*2dce792eSToby Isaac         PetscInt       groupSize;
616*2dce792eSToby Isaac         PetscInt     **cellPerms;
617*2dce792eSToby Isaac         PetscScalar  **cellFlips;
618*2dce792eSToby Isaac         DMPolytopeType ct;
619*2dce792eSToby Isaac 
620*2dce792eSToby Isaac         PetscCall(DMPlexGetCellType(K, 0, &ct));
621*2dce792eSToby Isaac         groupSize       = DMPolytopeTypeGetNumArrangments(ct);
622*2dce792eSToby Isaac         sum->numSelfSym = groupSize;
623*2dce792eSToby Isaac         sum->selfSymOff = groupSize / 2;
624*2dce792eSToby Isaac         PetscCall(PetscCalloc1(groupSize, &cellPerms));
625*2dce792eSToby Isaac         PetscCall(PetscCalloc1(groupSize, &cellFlips));
626*2dce792eSToby Isaac         symperms[0] = &cellPerms[groupSize / 2];
627*2dce792eSToby Isaac         symflips[0] = &cellFlips[groupSize / 2];
628*2dce792eSToby Isaac         for (PetscInt o = -groupSize / 2; o < groupSize / 2; o++) {
629*2dce792eSToby Isaac           PetscBool any_o_perms = PETSC_FALSE;
630*2dce792eSToby Isaac           PetscBool any_o_flips = PETSC_FALSE;
631*2dce792eSToby Isaac 
632*2dce792eSToby Isaac           for (PetscInt s = 0; s < Ns; s++) {
633*2dce792eSToby Isaac             PetscDualSpace       subsp;
634*2dce792eSToby Isaac             const PetscInt    ***sub_perms;
635*2dce792eSToby Isaac             const PetscScalar ***sub_flips;
636*2dce792eSToby Isaac 
637*2dce792eSToby Isaac             PetscCall(PetscDualSpaceSumGetSubspace(sp, s, &subsp));
638*2dce792eSToby Isaac             PetscCall(PetscDualSpaceGetSymmetries(subsp, &sub_perms, &sub_flips));
639*2dce792eSToby Isaac             if (sub_perms && sub_perms[0] && sub_perms[0][o]) any_o_perms = PETSC_TRUE;
640*2dce792eSToby Isaac             if (sub_flips && sub_flips[0] && sub_flips[0][o]) any_o_flips = PETSC_TRUE;
641*2dce792eSToby Isaac           }
642*2dce792eSToby Isaac           if (any_o_perms) {
643*2dce792eSToby Isaac             PetscInt *o_perm;
644*2dce792eSToby Isaac 
645*2dce792eSToby Isaac             PetscCall(PetscMalloc1(spintdim, &o_perm));
646*2dce792eSToby Isaac             for (PetscInt i = 0; i < spintdim; i++) o_perm[i] = i;
647*2dce792eSToby Isaac             for (PetscInt s = 0; s < Ns; s++) {
648*2dce792eSToby Isaac               PetscDualSpace       subsp;
649*2dce792eSToby Isaac               const PetscInt    ***sub_perms;
650*2dce792eSToby Isaac               const PetscScalar ***sub_flips;
651*2dce792eSToby Isaac 
652*2dce792eSToby Isaac               PetscCall(PetscDualSpaceSumGetSubspace(sp, s, &subsp));
653*2dce792eSToby Isaac               PetscCall(PetscDualSpaceGetSymmetries(subsp, &sub_perms, &sub_flips));
654*2dce792eSToby Isaac               if (sub_perms && sub_perms[0] && sub_perms[0][o]) {
655*2dce792eSToby Isaac                 PetscInt  subspdim;
656*2dce792eSToby Isaac                 PetscInt *range, *domain;
657*2dce792eSToby Isaac                 PetscInt *range_mapped, *domain_mapped;
658*2dce792eSToby Isaac 
659*2dce792eSToby Isaac                 PetscCall(PetscDualSpaceGetInteriorDimension(subsp, &subspdim));
660*2dce792eSToby Isaac                 PetscCall(PetscMalloc4(subspdim, &range, subspdim, &range_mapped, subspdim, &domain, subspdim, &domain_mapped));
661*2dce792eSToby Isaac                 for (PetscInt i = 0; i < subspdim; i++) domain[i] = i;
662*2dce792eSToby Isaac                 PetscCall(PetscArraycpy(range, sub_perms[0][o], subspdim));
663*2dce792eSToby Isaac                 PetscCall(ISLocalToGlobalMappingApply(sum->int_rows[s], subspdim, domain, domain_mapped));
664*2dce792eSToby Isaac                 PetscCall(ISLocalToGlobalMappingApply(sum->int_rows[s], subspdim, range, range_mapped));
665*2dce792eSToby Isaac                 for (PetscInt i = 0; i < subspdim; i++) o_perm[domain_mapped[i]] = range_mapped[i];
666*2dce792eSToby Isaac                 PetscCall(PetscFree4(range, range_mapped, domain, domain_mapped));
667*2dce792eSToby Isaac               }
668*2dce792eSToby Isaac             }
669*2dce792eSToby Isaac             symperms[0][o] = o_perm;
670*2dce792eSToby Isaac           }
671*2dce792eSToby Isaac           if (any_o_flips) {
672*2dce792eSToby Isaac             PetscScalar *o_flip;
673*2dce792eSToby Isaac 
674*2dce792eSToby Isaac             PetscCall(PetscMalloc1(spintdim, &o_flip));
675*2dce792eSToby Isaac             for (PetscInt i = 0; i < spintdim; i++) o_flip[i] = 1.0;
676*2dce792eSToby Isaac             for (PetscInt s = 0; s < Ns; s++) {
677*2dce792eSToby Isaac               PetscDualSpace       subsp;
678*2dce792eSToby Isaac               const PetscInt    ***sub_perms;
679*2dce792eSToby Isaac               const PetscScalar ***sub_flips;
680*2dce792eSToby Isaac 
681*2dce792eSToby Isaac               PetscCall(PetscDualSpaceSumGetSubspace(sp, s, &subsp));
682*2dce792eSToby Isaac               PetscCall(PetscDualSpaceGetSymmetries(subsp, &sub_perms, &sub_flips));
683*2dce792eSToby Isaac               if (sub_perms && sub_perms[0] && sub_perms[0][o]) {
684*2dce792eSToby Isaac                 PetscInt  subspdim;
685*2dce792eSToby Isaac                 PetscInt *domain;
686*2dce792eSToby Isaac                 PetscInt *domain_mapped;
687*2dce792eSToby Isaac 
688*2dce792eSToby Isaac                 PetscCall(PetscDualSpaceGetInteriorDimension(subsp, &subspdim));
689*2dce792eSToby Isaac                 PetscCall(PetscMalloc2(subspdim, &domain, subspdim, &domain_mapped));
690*2dce792eSToby Isaac                 for (PetscInt i = 0; i < subspdim; i++) domain[i] = i;
691*2dce792eSToby Isaac                 PetscCall(ISLocalToGlobalMappingApply(sum->int_rows[s], subspdim, domain, domain_mapped));
692*2dce792eSToby Isaac                 for (PetscInt i = 0; i < subspdim; i++) o_flip[domain_mapped[i]] = sub_perms[0][o][i];
693*2dce792eSToby Isaac                 PetscCall(PetscFree2(domain, domain_mapped));
694*2dce792eSToby Isaac               }
695*2dce792eSToby Isaac             }
696*2dce792eSToby Isaac             symflips[0][o] = o_flip;
697*2dce792eSToby Isaac           }
698*2dce792eSToby Isaac         }
699*2dce792eSToby Isaac         {
700*2dce792eSToby Isaac           PetscBool any_perms = PETSC_FALSE;
701*2dce792eSToby Isaac           PetscBool any_flips = PETSC_FALSE;
702*2dce792eSToby Isaac           for (PetscInt o = -groupSize / 2; o < groupSize / 2; o++) {
703*2dce792eSToby Isaac             if (symperms[0][o]) any_perms = PETSC_TRUE;
704*2dce792eSToby Isaac             if (symflips[0][o]) any_flips = PETSC_TRUE;
705*2dce792eSToby Isaac           }
706*2dce792eSToby Isaac           if (!any_perms) {
707*2dce792eSToby Isaac             PetscCall(PetscFree(cellPerms));
708*2dce792eSToby Isaac             symperms[0] = NULL;
709*2dce792eSToby Isaac           }
710*2dce792eSToby Isaac           if (!any_flips) {
711*2dce792eSToby Isaac             PetscCall(PetscFree(cellFlips));
712*2dce792eSToby Isaac             symflips[0] = NULL;
713*2dce792eSToby Isaac           }
714*2dce792eSToby Isaac         }
715*2dce792eSToby Isaac       }
716*2dce792eSToby Isaac       if (!any_perms) {
717*2dce792eSToby Isaac         PetscCall(PetscFree(symperms));
718*2dce792eSToby Isaac         symperms = NULL;
719*2dce792eSToby Isaac       }
720*2dce792eSToby Isaac       if (!any_flips) {
721*2dce792eSToby Isaac         PetscCall(PetscFree(symflips));
722*2dce792eSToby Isaac         symflips = NULL;
723*2dce792eSToby Isaac       }
724*2dce792eSToby Isaac     }
725*2dce792eSToby Isaac     sum->symperms = symperms;
726*2dce792eSToby Isaac     sum->symflips = symflips;
727*2dce792eSToby Isaac   }
728*2dce792eSToby Isaac   if (perms) *perms = (const PetscInt ***)sum->symperms;
729*2dce792eSToby Isaac   if (flips) *flips = (const PetscScalar ***)sum->symflips;
730*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
731*2dce792eSToby Isaac }
732*2dce792eSToby Isaac 
733*2dce792eSToby Isaac static PetscErrorCode PetscDualSpaceSetUp_Sum(PetscDualSpace sp)
734*2dce792eSToby Isaac {
735*2dce792eSToby Isaac   PetscDualSpace_Sum *sum         = (PetscDualSpace_Sum *)sp->data;
736*2dce792eSToby Isaac   PetscBool           concatenate = PETSC_TRUE;
737*2dce792eSToby Isaac   PetscBool           uniform;
738*2dce792eSToby Isaac   PetscInt            Ns, Nc, i, sum_Nc = 0;
739*2dce792eSToby Isaac   PetscInt            minNc, maxNc;
740*2dce792eSToby Isaac   PetscInt            minForm, maxForm, cdim, depth;
741*2dce792eSToby Isaac   DM                  K;
742*2dce792eSToby Isaac   PetscQuadrature    *all_quads = NULL;
743*2dce792eSToby Isaac   PetscQuadrature    *int_quads = NULL;
744*2dce792eSToby Isaac   Mat                *all_mats  = NULL;
745*2dce792eSToby Isaac   Mat                *int_mats  = NULL;
746*2dce792eSToby Isaac 
747*2dce792eSToby Isaac   PetscFunctionBegin;
748*2dce792eSToby Isaac   if (sum->setupCalled) PetscFunctionReturn(PETSC_SUCCESS);
749*2dce792eSToby Isaac   sum->setupCalled = PETSC_TRUE;
750*2dce792eSToby Isaac 
751*2dce792eSToby Isaac   PetscCall(PetscDualSpaceSumGetNumSubspaces(sp, &Ns));
752*2dce792eSToby Isaac   PetscCheck(Ns >= 0, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_OUTOFRANGE, "Cannot have %" PetscInt_FMT " subspaces", Ns);
753*2dce792eSToby Isaac 
754*2dce792eSToby Isaac   // step 1: make sure they share a DM
755*2dce792eSToby Isaac   PetscCall(PetscDualSpaceGetDM(sp, &K));
756*2dce792eSToby Isaac   PetscCall(DMGetCoordinateDim(K, &cdim));
757*2dce792eSToby Isaac   PetscCall(DMPlexGetDepth(K, &depth));
758*2dce792eSToby Isaac   PetscCall(PetscDualSpaceSumIsUniform(sp, &sp->uniform));
759*2dce792eSToby Isaac   uniform = sp->uniform;
760*2dce792eSToby Isaac   {
761*2dce792eSToby Isaac     for (PetscInt s = 0; s < Ns; s++) {
762*2dce792eSToby Isaac       PetscDualSpace subsp;
763*2dce792eSToby Isaac       DM             sub_K;
764*2dce792eSToby Isaac 
765*2dce792eSToby Isaac       PetscCall(PetscDualSpaceSumGetSubspace(sp, s, &subsp));
766*2dce792eSToby Isaac       PetscCall(PetscDualSpaceSetUp(subsp));
767*2dce792eSToby Isaac       PetscCall(PetscDualSpaceGetDM(subsp, &sub_K));
768*2dce792eSToby Isaac       if (s == 0 && K == NULL) {
769*2dce792eSToby Isaac         PetscCall(PetscDualSpaceSetDM(sp, sub_K));
770*2dce792eSToby Isaac         K = sub_K;
771*2dce792eSToby Isaac       }
772*2dce792eSToby Isaac       PetscCheck(sub_K == K, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Subspace %d does not have the same DM as the sum space", (int)s);
773*2dce792eSToby Isaac     }
774*2dce792eSToby Isaac   }
775*2dce792eSToby Isaac 
776*2dce792eSToby Isaac   // step 2: count components
777*2dce792eSToby Isaac   PetscCall(PetscDualSpaceGetNumComponents(sp, &Nc));
778*2dce792eSToby Isaac   PetscCall(PetscDualSpaceSumGetConcatenate(sp, &concatenate));
779*2dce792eSToby Isaac   minNc   = Nc;
780*2dce792eSToby Isaac   maxNc   = Nc;
781*2dce792eSToby Isaac   minForm = sp->k == PETSC_FORM_DEGREE_UNDEFINED ? PETSC_INT_MAX : sp->k;
782*2dce792eSToby Isaac   maxForm = sp->k == PETSC_FORM_DEGREE_UNDEFINED ? PETSC_INT_MIN : sp->k;
783*2dce792eSToby Isaac   for (i = 0; i < Ns; ++i) {
784*2dce792eSToby Isaac     PetscInt       sNc, formDegree;
785*2dce792eSToby Isaac     PetscDualSpace si;
786*2dce792eSToby Isaac 
787*2dce792eSToby Isaac     PetscCall(PetscDualSpaceSumGetSubspace(sp, i, &si));
788*2dce792eSToby Isaac     PetscCall(PetscDualSpaceSetUp(si));
789*2dce792eSToby Isaac     PetscCall(PetscDualSpaceGetNumComponents(si, &sNc));
790*2dce792eSToby Isaac     if (sNc != Nc) PetscCheck(concatenate, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Subspace as a different number of components but space does not concatenate components");
791*2dce792eSToby Isaac     minNc = PetscMin(minNc, sNc);
792*2dce792eSToby Isaac     maxNc = PetscMax(maxNc, sNc);
793*2dce792eSToby Isaac     sum_Nc += sNc;
794*2dce792eSToby Isaac     PetscCall(PetscDualSpaceGetFormDegree(si, &formDegree));
795*2dce792eSToby Isaac     minForm = PetscMin(minForm, formDegree);
796*2dce792eSToby Isaac     maxForm = PetscMax(maxForm, formDegree);
797*2dce792eSToby Isaac   }
798*2dce792eSToby Isaac   sp->k = (minForm != maxForm) ? PETSC_FORM_DEGREE_UNDEFINED : minForm;
799*2dce792eSToby Isaac 
800*2dce792eSToby Isaac   if (concatenate) PetscCheck(sum_Nc == Nc, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_OUTOFRANGE, "Total number of subspace components (%" PetscInt_FMT ") does not match number of target space components (%" PetscInt_FMT ").", sum_Nc, Nc);
801*2dce792eSToby Isaac   else PetscCheck(minNc == Nc && maxNc == Nc, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_OUTOFRANGE, "Subspaces must have same number of components as the target space.");
802*2dce792eSToby Isaac 
803*2dce792eSToby Isaac   /* A PetscDualSpace should have a fixed number of components, but
804*2dce792eSToby Isaac      if the spaces we are combining have different form degrees, they will not
805*2dce792eSToby Isaac      have the same number of components on subcomponents of the boundary,
806*2dce792eSToby Isaac      so we do not try to create boundary dual spaces in this case */
807*2dce792eSToby Isaac   if (sp->k != PETSC_FORM_DEGREE_UNDEFINED && depth > 0) {
808*2dce792eSToby Isaac     PetscInt  pStart, pEnd;
809*2dce792eSToby Isaac     PetscInt *pStratStart, *pStratEnd;
810*2dce792eSToby Isaac 
811*2dce792eSToby Isaac     PetscCall(DMPlexGetDepth(K, &depth));
812*2dce792eSToby Isaac     PetscCall(DMPlexGetChart(K, &pStart, &pEnd));
813*2dce792eSToby Isaac     PetscCall(PetscCalloc1(pEnd, &(sp->pointSpaces)));
814*2dce792eSToby Isaac     PetscCall(PetscMalloc2(depth + 1, &pStratStart, depth + 1, &pStratEnd));
815*2dce792eSToby Isaac     for (PetscInt d = 0; d <= depth; ++d) PetscCall(DMPlexGetDepthStratum(K, d, &pStratStart[d], &pStratEnd[d]));
816*2dce792eSToby Isaac 
817*2dce792eSToby Isaac     for (PetscInt p = pStratStart[depth - 1]; p < pStratEnd[depth - 1]; p++) { /* calculate the facet dual spaces */
818*2dce792eSToby Isaac       PetscReal      v0[3];
819*2dce792eSToby Isaac       DMPolytopeType ptype;
820*2dce792eSToby Isaac       PetscReal      J[9], detJ;
821*2dce792eSToby Isaac       PetscInt       q;
822*2dce792eSToby Isaac 
823*2dce792eSToby Isaac       PetscCall(DMPlexComputeCellGeometryAffineFEM(K, p, v0, J, NULL, &detJ));
824*2dce792eSToby Isaac       PetscCall(DMPlexGetCellType(K, p, &ptype));
825*2dce792eSToby Isaac 
826*2dce792eSToby Isaac       /* compare to previous facets: if computed, reference that dualspace */
827*2dce792eSToby Isaac       for (q = pStratStart[depth - 1]; q < p; q++) {
828*2dce792eSToby Isaac         DMPolytopeType qtype;
829*2dce792eSToby Isaac 
830*2dce792eSToby Isaac         PetscCall(DMPlexGetCellType(K, q, &qtype));
831*2dce792eSToby Isaac         if (qtype == ptype) break;
832*2dce792eSToby Isaac       }
833*2dce792eSToby Isaac       if (q < p) { /* this facet has the same dual space as that one */
834*2dce792eSToby Isaac         PetscCall(PetscObjectReference((PetscObject)sp->pointSpaces[q]));
835*2dce792eSToby Isaac         sp->pointSpaces[p] = sp->pointSpaces[q];
836*2dce792eSToby Isaac         continue;
837*2dce792eSToby Isaac       }
838*2dce792eSToby Isaac       /* if not, recursively compute this dual space */
839*2dce792eSToby Isaac       PetscCall(PetscDualSpaceCreateFacetSubspace_Sum(sp, p, &sp->pointSpaces[p]));
840*2dce792eSToby Isaac     }
841*2dce792eSToby Isaac     for (PetscInt h = 2; h <= depth; h++) { /* get the higher subspaces from the facet subspaces */
842*2dce792eSToby Isaac       PetscInt hd = depth - h;
843*2dce792eSToby Isaac 
844*2dce792eSToby Isaac       for (PetscInt p = pStratStart[hd]; p < pStratEnd[hd]; p++) {
845*2dce792eSToby Isaac         PetscInt        suppSize;
846*2dce792eSToby Isaac         const PetscInt *supp;
847*2dce792eSToby Isaac         DM              qdm;
848*2dce792eSToby Isaac         PetscDualSpace  qsp, psp;
849*2dce792eSToby Isaac         PetscInt        c, coneSize, q;
850*2dce792eSToby Isaac         const PetscInt *cone;
851*2dce792eSToby Isaac         const PetscInt *refCone;
852*2dce792eSToby Isaac 
853*2dce792eSToby Isaac         PetscCall(DMPlexGetSupportSize(K, p, &suppSize));
854*2dce792eSToby Isaac         PetscCall(DMPlexGetSupport(K, p, &supp));
855*2dce792eSToby Isaac         q   = supp[0];
856*2dce792eSToby Isaac         qsp = sp->pointSpaces[q];
857*2dce792eSToby Isaac         PetscCall(DMPlexGetConeSize(K, q, &coneSize));
858*2dce792eSToby Isaac         PetscCall(DMPlexGetCone(K, q, &cone));
859*2dce792eSToby Isaac         for (c = 0; c < coneSize; c++)
860*2dce792eSToby Isaac           if (cone[c] == p) break;
861*2dce792eSToby Isaac         PetscCheck(c != coneSize, PetscObjectComm((PetscObject)K), PETSC_ERR_PLIB, "cone/support mismatch");
862*2dce792eSToby Isaac         if (!qsp) {
863*2dce792eSToby Isaac           sp->pointSpaces[p] = NULL;
864*2dce792eSToby Isaac           continue;
865*2dce792eSToby Isaac         }
866*2dce792eSToby Isaac         PetscCall(PetscDualSpaceGetDM(qsp, &qdm));
867*2dce792eSToby Isaac         PetscCall(DMPlexGetCone(qdm, 0, &refCone));
868*2dce792eSToby Isaac         /* get the equivalent dual space from the support dual space */
869*2dce792eSToby Isaac         PetscCall(PetscDualSpaceGetPointSubspace(qsp, refCone[c], &psp));
870*2dce792eSToby Isaac         PetscCall(PetscObjectReference((PetscObject)psp));
871*2dce792eSToby Isaac         sp->pointSpaces[p] = psp;
872*2dce792eSToby Isaac       }
873*2dce792eSToby Isaac     }
874*2dce792eSToby Isaac     PetscCall(PetscFree2(pStratStart, pStratEnd));
875*2dce792eSToby Isaac   }
876*2dce792eSToby Isaac 
877*2dce792eSToby Isaac   sum->uniform = uniform;
878*2dce792eSToby Isaac   PetscCall(PetscCalloc1(Ns, &sum->all_rows));
879*2dce792eSToby Isaac   PetscCall(PetscCalloc1(Ns, &sum->all_cols));
880*2dce792eSToby Isaac   PetscCall(PetscCalloc1(Ns, &sum->int_rows));
881*2dce792eSToby Isaac   PetscCall(PetscCalloc1(Ns, &sum->int_cols));
882*2dce792eSToby Isaac   PetscCall(PetscMalloc4(Ns, &all_quads, Ns, &all_mats, Ns, &int_quads, Ns, &int_mats));
883*2dce792eSToby Isaac   {
884*2dce792eSToby Isaac     // test for uniform all points and uniform interior points
885*2dce792eSToby Isaac     PetscBool       uniform_all         = PETSC_TRUE;
886*2dce792eSToby Isaac     PetscBool       uniform_interior    = PETSC_TRUE;
887*2dce792eSToby Isaac     PetscQuadrature quad_all_first      = NULL;
888*2dce792eSToby Isaac     PetscQuadrature quad_interior_first = NULL;
889*2dce792eSToby Isaac     PetscInt        pStart, pEnd;
890*2dce792eSToby Isaac 
891*2dce792eSToby Isaac     PetscCall(DMPlexGetChart(K, &pStart, &pEnd));
892*2dce792eSToby Isaac     PetscCall(PetscDualSpaceSectionCreate_Internal(sp, &sp->pointSection));
893*2dce792eSToby Isaac 
894*2dce792eSToby Isaac     for (PetscInt p = pStart; p < pEnd; p++) {
895*2dce792eSToby Isaac       PetscInt full_dof = 0;
896*2dce792eSToby Isaac 
897*2dce792eSToby Isaac       for (PetscInt s = 0; s < Ns; s++) {
898*2dce792eSToby Isaac         PetscDualSpace subsp;
899*2dce792eSToby Isaac         PetscSection   subsection;
900*2dce792eSToby Isaac         PetscInt       dof;
901*2dce792eSToby Isaac 
902*2dce792eSToby Isaac         PetscCall(PetscDualSpaceSumGetSubspace(sp, s, &subsp));
903*2dce792eSToby Isaac         PetscCall(PetscDualSpaceGetSection(subsp, &subsection));
904*2dce792eSToby Isaac         PetscCall(PetscSectionGetDof(subsection, p, &dof));
905*2dce792eSToby Isaac         full_dof += dof;
906*2dce792eSToby Isaac       }
907*2dce792eSToby Isaac       PetscCall(PetscSectionSetDof(sp->pointSection, p, full_dof));
908*2dce792eSToby Isaac     }
909*2dce792eSToby Isaac     PetscCall(PetscDualSpaceSectionSetUp_Internal(sp, sp->pointSection));
910*2dce792eSToby Isaac 
911*2dce792eSToby Isaac     for (PetscInt s = 0; s < Ns; s++) {
912*2dce792eSToby Isaac       PetscDualSpace  subsp;
913*2dce792eSToby Isaac       PetscQuadrature subquad_all;
914*2dce792eSToby Isaac       PetscQuadrature subquad_interior;
915*2dce792eSToby Isaac       Mat             submat_all;
916*2dce792eSToby Isaac       Mat             submat_interior;
917*2dce792eSToby Isaac 
918*2dce792eSToby Isaac       PetscCall(PetscDualSpaceSumGetSubspace(sp, s, &subsp));
919*2dce792eSToby Isaac       PetscCall(PetscDualSpaceGetAllData(subsp, &subquad_all, &submat_all));
920*2dce792eSToby Isaac       PetscCall(PetscDualSpaceGetInteriorData(subsp, &subquad_interior, &submat_interior));
921*2dce792eSToby Isaac       if (!s) {
922*2dce792eSToby Isaac         quad_all_first      = subquad_all;
923*2dce792eSToby Isaac         quad_interior_first = subquad_interior;
924*2dce792eSToby Isaac       } else {
925*2dce792eSToby Isaac         if (subquad_all != quad_all_first) uniform_all = PETSC_FALSE;
926*2dce792eSToby Isaac         if (subquad_interior != quad_interior_first) uniform_interior = PETSC_FALSE;
927*2dce792eSToby Isaac       }
928*2dce792eSToby Isaac       all_quads[s] = subquad_all;
929*2dce792eSToby Isaac       int_quads[s] = subquad_interior;
930*2dce792eSToby Isaac       all_mats[s]  = submat_all;
931*2dce792eSToby Isaac       int_mats[s]  = submat_interior;
932*2dce792eSToby Isaac     }
933*2dce792eSToby Isaac     sum->uniform_all_points      = uniform_all;
934*2dce792eSToby Isaac     sum->uniform_interior_points = uniform_interior;
935*2dce792eSToby Isaac 
936*2dce792eSToby Isaac     PetscCall(PetscDualSpaceSumCreateMappings(sp, PETSC_TRUE, uniform_interior, sum->int_rows, sum->int_cols));
937*2dce792eSToby Isaac     PetscCall(PetscDualSpaceSumCreateQuadrature(Ns, cdim, uniform_interior, int_quads, &sp->intNodes));
938*2dce792eSToby Isaac     if (sp->intNodes) PetscCall(PetscDualSpaceSumCreateMatrix(sp, int_mats, sum->int_rows, sum->int_cols, sp->intNodes, &sp->intMat));
939*2dce792eSToby Isaac 
940*2dce792eSToby Isaac     PetscCall(PetscDualSpaceSumCreateMappings(sp, PETSC_FALSE, uniform_all, sum->all_rows, sum->all_cols));
941*2dce792eSToby Isaac     PetscCall(PetscDualSpaceSumCreateQuadrature(Ns, cdim, uniform_all, all_quads, &sp->allNodes));
942*2dce792eSToby Isaac     if (sp->allNodes) PetscCall(PetscDualSpaceSumCreateMatrix(sp, all_mats, sum->all_rows, sum->all_cols, sp->allNodes, &sp->allMat));
943*2dce792eSToby Isaac   }
944*2dce792eSToby Isaac   PetscCall(PetscFree4(all_quads, all_mats, int_quads, int_mats));
945*2dce792eSToby Isaac   PetscCall(PetscDualSpaceComputeFunctionalsFromAllData(sp));
946*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
947*2dce792eSToby Isaac }
948*2dce792eSToby Isaac 
949*2dce792eSToby Isaac static PetscErrorCode PetscDualSpaceSumView_Ascii(PetscDualSpace sp, PetscViewer v)
950*2dce792eSToby Isaac {
951*2dce792eSToby Isaac   PetscDualSpace_Sum *sum         = (PetscDualSpace_Sum *)sp->data;
952*2dce792eSToby Isaac   PetscBool           concatenate = sum->concatenate;
953*2dce792eSToby Isaac   PetscInt            i, Ns = sum->numSumSpaces;
954*2dce792eSToby Isaac 
955*2dce792eSToby Isaac   PetscFunctionBegin;
956*2dce792eSToby Isaac   if (concatenate) PetscCall(PetscViewerASCIIPrintf(v, "Sum dual space of %" PetscInt_FMT " concatenated subspaces%s\n", Ns, sum->uniform ? " (all identical)" : ""));
957*2dce792eSToby Isaac   else PetscCall(PetscViewerASCIIPrintf(v, "Sum dual space of %" PetscInt_FMT " subspaces%s\n", Ns, sum->uniform ? " (all identical)" : ""));
958*2dce792eSToby Isaac   for (i = 0; i < (sum->uniform ? (Ns > 0 ? 1 : 0) : Ns); ++i) {
959*2dce792eSToby Isaac     PetscCall(PetscViewerASCIIPushTab(v));
960*2dce792eSToby Isaac     PetscCall(PetscDualSpaceView(sum->sumspaces[i], v));
961*2dce792eSToby Isaac     PetscCall(PetscViewerASCIIPopTab(v));
962*2dce792eSToby Isaac   }
963*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
964*2dce792eSToby Isaac }
965*2dce792eSToby Isaac 
966*2dce792eSToby Isaac static PetscErrorCode PetscDualSpaceView_Sum(PetscDualSpace sp, PetscViewer viewer)
967*2dce792eSToby Isaac {
968*2dce792eSToby Isaac   PetscBool iascii;
969*2dce792eSToby Isaac 
970*2dce792eSToby Isaac   PetscFunctionBegin;
971*2dce792eSToby Isaac   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
972*2dce792eSToby Isaac   if (iascii) PetscCall(PetscDualSpaceSumView_Ascii(sp, viewer));
973*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
974*2dce792eSToby Isaac }
975*2dce792eSToby Isaac 
976*2dce792eSToby Isaac static PetscErrorCode PetscDualSpaceDestroy_Sum(PetscDualSpace sp)
977*2dce792eSToby Isaac {
978*2dce792eSToby Isaac   PetscDualSpace_Sum *sum = (PetscDualSpace_Sum *)sp->data;
979*2dce792eSToby Isaac   PetscInt            i, Ns = sum->numSumSpaces;
980*2dce792eSToby Isaac 
981*2dce792eSToby Isaac   PetscFunctionBegin;
982*2dce792eSToby Isaac   if (sum->symperms) {
983*2dce792eSToby Isaac     PetscInt **selfSyms = sum->symperms[0];
984*2dce792eSToby Isaac 
985*2dce792eSToby Isaac     if (selfSyms) {
986*2dce792eSToby Isaac       PetscInt i, **allocated = &selfSyms[-sum->selfSymOff];
987*2dce792eSToby Isaac 
988*2dce792eSToby Isaac       for (i = 0; i < sum->numSelfSym; i++) PetscCall(PetscFree(allocated[i]));
989*2dce792eSToby Isaac       PetscCall(PetscFree(allocated));
990*2dce792eSToby Isaac     }
991*2dce792eSToby Isaac     PetscCall(PetscFree(sum->symperms));
992*2dce792eSToby Isaac   }
993*2dce792eSToby Isaac   if (sum->symflips) {
994*2dce792eSToby Isaac     PetscScalar **selfSyms = sum->symflips[0];
995*2dce792eSToby Isaac 
996*2dce792eSToby Isaac     if (selfSyms) {
997*2dce792eSToby Isaac       PetscInt      i;
998*2dce792eSToby Isaac       PetscScalar **allocated = &selfSyms[-sum->selfSymOff];
999*2dce792eSToby Isaac 
1000*2dce792eSToby Isaac       for (i = 0; i < sum->numSelfSym; i++) PetscCall(PetscFree(allocated[i]));
1001*2dce792eSToby Isaac       PetscCall(PetscFree(allocated));
1002*2dce792eSToby Isaac     }
1003*2dce792eSToby Isaac     PetscCall(PetscFree(sum->symflips));
1004*2dce792eSToby Isaac   }
1005*2dce792eSToby Isaac   for (i = 0; i < Ns; ++i) {
1006*2dce792eSToby Isaac     PetscCall(PetscDualSpaceDestroy(&sum->sumspaces[i]));
1007*2dce792eSToby Isaac     if (sum->all_rows) PetscCall(ISLocalToGlobalMappingDestroy(&sum->all_rows[i]));
1008*2dce792eSToby Isaac     if (sum->all_cols) PetscCall(ISLocalToGlobalMappingDestroy(&sum->all_cols[i]));
1009*2dce792eSToby Isaac     if (sum->int_rows) PetscCall(ISLocalToGlobalMappingDestroy(&sum->int_rows[i]));
1010*2dce792eSToby Isaac     if (sum->int_cols) PetscCall(ISLocalToGlobalMappingDestroy(&sum->int_cols[i]));
1011*2dce792eSToby Isaac   }
1012*2dce792eSToby Isaac   PetscCall(PetscFree(sum->sumspaces));
1013*2dce792eSToby Isaac   PetscCall(PetscFree(sum->all_rows));
1014*2dce792eSToby Isaac   PetscCall(PetscFree(sum->all_cols));
1015*2dce792eSToby Isaac   PetscCall(PetscFree(sum->int_rows));
1016*2dce792eSToby Isaac   PetscCall(PetscFree(sum->int_cols));
1017*2dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumSetSubspace_C", NULL));
1018*2dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumGetSubspace_C", NULL));
1019*2dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumSetNumSubspaces_C", NULL));
1020*2dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumGetNumSubspaces_C", NULL));
1021*2dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumGetConcatenate_C", NULL));
1022*2dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumSetConcatenate_C", NULL));
1023*2dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumGetInterleave_C", NULL));
1024*2dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumSetInterleave_C", NULL));
1025*2dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetContinuity_C", NULL));
1026*2dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetContinuity_C", NULL));
1027*2dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetMomentOrder_C", NULL));
1028*2dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetMomentOrder_C", NULL));
1029*2dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetNodeType_C", NULL));
1030*2dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetNodeType_C", NULL));
1031*2dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetTensor_C", NULL));
1032*2dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetTensor_C", NULL));
1033*2dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetTrimmed_C", NULL));
1034*2dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetTrimmed_C", NULL));
1035*2dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetUseMoments_C", NULL));
1036*2dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetUseMoments_C", NULL));
1037*2dce792eSToby Isaac   PetscCall(PetscFree(sum));
1038*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1039*2dce792eSToby Isaac }
1040*2dce792eSToby Isaac 
1041*2dce792eSToby Isaac /*@
1042*2dce792eSToby Isaac   PetscDualSpaceSumSetInterleave - Set whether the basis functions and components of a uniform sum are interleaved
1043*2dce792eSToby Isaac 
1044*2dce792eSToby Isaac   Logically collective
1045*2dce792eSToby Isaac 
1046*2dce792eSToby Isaac   Input Parameters:
1047*2dce792eSToby Isaac + sp                    - a `PetscDualSpace` of type `PETSCDUALSPACESUM`
1048*2dce792eSToby Isaac . interleave_basis      - if `PETSC_TRUE`, the basis vectors of the subspaces are interleaved
1049*2dce792eSToby Isaac - interleave_components - if `PETSC_TRUE` and the space concatenates components (`PetscDualSpaceSumGetConcatenate()`),
1050*2dce792eSToby Isaac                           interleave the concatenated components
1051*2dce792eSToby Isaac 
1052*2dce792eSToby Isaac   Level: developer
1053*2dce792eSToby Isaac 
1054*2dce792eSToby Isaac .seealso: `PetscDualSpace`, `PETSCDUALSPACESUM`, `PETSCFEVECTOR`, `PetscDualSpaceSumGetInterleave()`
1055*2dce792eSToby Isaac @*/
1056*2dce792eSToby Isaac PetscErrorCode PetscDualSpaceSumSetInterleave(PetscDualSpace sp, PetscBool interleave_basis, PetscBool interleave_components)
1057*2dce792eSToby Isaac {
1058*2dce792eSToby Isaac   PetscFunctionBegin;
1059*2dce792eSToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1060*2dce792eSToby Isaac   PetscTryMethod(sp, "PetscDualSpaceSumSetInterleave_C", (PetscDualSpace, PetscBool, PetscBool), (sp, interleave_basis, interleave_components));
1061*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1062*2dce792eSToby Isaac }
1063*2dce792eSToby Isaac 
1064*2dce792eSToby Isaac static PetscErrorCode PetscDualSpaceSumSetInterleave_Sum(PetscDualSpace sp, PetscBool interleave_basis, PetscBool interleave_components)
1065*2dce792eSToby Isaac {
1066*2dce792eSToby Isaac   PetscDualSpace_Sum *sum = (PetscDualSpace_Sum *)sp->data;
1067*2dce792eSToby Isaac   PetscFunctionBegin;
1068*2dce792eSToby Isaac   sum->interleave_basis      = interleave_basis;
1069*2dce792eSToby Isaac   sum->interleave_components = interleave_components;
1070*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1071*2dce792eSToby Isaac }
1072*2dce792eSToby Isaac 
1073*2dce792eSToby Isaac /*@
1074*2dce792eSToby Isaac   PetscDualSpaceSumGetInterleave - Get whether the basis functions and components of a uniform sum are interleaved
1075*2dce792eSToby Isaac 
1076*2dce792eSToby Isaac   Logically collective
1077*2dce792eSToby Isaac 
1078*2dce792eSToby Isaac   Input Parameter:
1079*2dce792eSToby Isaac . sp - a `PetscDualSpace` of type `PETSCDUALSPACESUM`
1080*2dce792eSToby Isaac 
1081*2dce792eSToby Isaac   Output Parameters:
1082*2dce792eSToby Isaac + interleave_basis      - if `PETSC_TRUE`, the basis vectors of the subspaces are interleaved
1083*2dce792eSToby Isaac - interleave_components - if `PETSC_TRUE` and the space concatenates components (`PetscDualSpaceSumGetConcatenate()`),
1084*2dce792eSToby Isaac                           interleave the concatenated components
1085*2dce792eSToby Isaac 
1086*2dce792eSToby Isaac   Level: developer
1087*2dce792eSToby Isaac 
1088*2dce792eSToby Isaac .seealso: `PetscDualSpace`, `PETSCDUALSPACESUM`, `PETSCFEVECTOR`, `PetscDualSpaceSumSetInterleave()`
1089*2dce792eSToby Isaac @*/
1090*2dce792eSToby Isaac PetscErrorCode PetscDualSpaceSumGetInterleave(PetscDualSpace sp, PetscBool *interleave_basis, PetscBool *interleave_components)
1091*2dce792eSToby Isaac {
1092*2dce792eSToby Isaac   PetscFunctionBegin;
1093*2dce792eSToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1094*2dce792eSToby Isaac   if (interleave_basis) PetscAssertPointer(interleave_basis, 2);
1095*2dce792eSToby Isaac   if (interleave_components) PetscAssertPointer(interleave_components, 3);
1096*2dce792eSToby Isaac   PetscTryMethod(sp, "PetscDualSpaceSumGetInterleave_C", (PetscDualSpace, PetscBool *, PetscBool *), (sp, interleave_basis, interleave_components));
1097*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1098*2dce792eSToby Isaac }
1099*2dce792eSToby Isaac 
1100*2dce792eSToby Isaac static PetscErrorCode PetscDualSpaceSumGetInterleave_Sum(PetscDualSpace sp, PetscBool *interleave_basis, PetscBool *interleave_components)
1101*2dce792eSToby Isaac {
1102*2dce792eSToby Isaac   PetscDualSpace_Sum *sum = (PetscDualSpace_Sum *)sp->data;
1103*2dce792eSToby Isaac   PetscFunctionBegin;
1104*2dce792eSToby Isaac   if (interleave_basis) *interleave_basis = sum->interleave_basis;
1105*2dce792eSToby Isaac   if (interleave_components) *interleave_components = sum->interleave_components;
1106*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1107*2dce792eSToby Isaac }
1108*2dce792eSToby Isaac 
1109*2dce792eSToby Isaac #define PetscDualSpaceSumPassthrough(sp, func, ...) \
1110*2dce792eSToby Isaac   do { \
1111*2dce792eSToby Isaac     PetscDualSpace_Sum *sum = (PetscDualSpace_Sum *)sp->data; \
1112*2dce792eSToby Isaac     PetscBool           is_uniform; \
1113*2dce792eSToby Isaac     PetscCall(PetscDualSpaceSumIsUniform(sp, &is_uniform)); \
1114*2dce792eSToby Isaac     if (is_uniform && sum->numSumSpaces > 0) { \
1115*2dce792eSToby Isaac       PetscDualSpace subsp; \
1116*2dce792eSToby Isaac       PetscCall(PetscDualSpaceSumGetSubspace(sp, 0, &subsp)); \
1117*2dce792eSToby Isaac       PetscCall(func(subsp, __VA_ARGS__)); \
1118*2dce792eSToby Isaac     } \
1119*2dce792eSToby Isaac   } while (0)
1120*2dce792eSToby Isaac 
1121*2dce792eSToby Isaac static PetscErrorCode PetscDualSpaceLagrangeGetContinuity_Sum(PetscDualSpace sp, PetscBool *value)
1122*2dce792eSToby Isaac {
1123*2dce792eSToby Isaac   PetscFunctionBegin;
1124*2dce792eSToby Isaac   PetscDualSpaceSumPassthrough(sp, PetscDualSpaceLagrangeGetContinuity, value);
1125*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1126*2dce792eSToby Isaac }
1127*2dce792eSToby Isaac 
1128*2dce792eSToby Isaac static PetscErrorCode PetscDualSpaceLagrangeSetContinuity_Sum(PetscDualSpace sp, PetscBool value)
1129*2dce792eSToby Isaac {
1130*2dce792eSToby Isaac   PetscFunctionBegin;
1131*2dce792eSToby Isaac   PetscDualSpaceSumPassthrough(sp, PetscDualSpaceLagrangeSetContinuity, value);
1132*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1133*2dce792eSToby Isaac }
1134*2dce792eSToby Isaac 
1135*2dce792eSToby Isaac static PetscErrorCode PetscDualSpaceLagrangeGetTensor_Sum(PetscDualSpace sp, PetscBool *value)
1136*2dce792eSToby Isaac {
1137*2dce792eSToby Isaac   PetscFunctionBegin;
1138*2dce792eSToby Isaac   PetscDualSpaceSumPassthrough(sp, PetscDualSpaceLagrangeGetTensor, value);
1139*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1140*2dce792eSToby Isaac }
1141*2dce792eSToby Isaac 
1142*2dce792eSToby Isaac static PetscErrorCode PetscDualSpaceLagrangeSetTensor_Sum(PetscDualSpace sp, PetscBool value)
1143*2dce792eSToby Isaac {
1144*2dce792eSToby Isaac   PetscFunctionBegin;
1145*2dce792eSToby Isaac   PetscDualSpaceSumPassthrough(sp, PetscDualSpaceLagrangeSetTensor, value);
1146*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1147*2dce792eSToby Isaac }
1148*2dce792eSToby Isaac 
1149*2dce792eSToby Isaac static PetscErrorCode PetscDualSpaceLagrangeGetTrimmed_Sum(PetscDualSpace sp, PetscBool *value)
1150*2dce792eSToby Isaac {
1151*2dce792eSToby Isaac   PetscFunctionBegin;
1152*2dce792eSToby Isaac   PetscDualSpaceSumPassthrough(sp, PetscDualSpaceLagrangeGetTrimmed, value);
1153*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1154*2dce792eSToby Isaac }
1155*2dce792eSToby Isaac 
1156*2dce792eSToby Isaac static PetscErrorCode PetscDualSpaceLagrangeSetTrimmed_Sum(PetscDualSpace sp, PetscBool value)
1157*2dce792eSToby Isaac {
1158*2dce792eSToby Isaac   PetscFunctionBegin;
1159*2dce792eSToby Isaac   PetscDualSpaceSumPassthrough(sp, PetscDualSpaceLagrangeSetTrimmed, value);
1160*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1161*2dce792eSToby Isaac }
1162*2dce792eSToby Isaac 
1163*2dce792eSToby Isaac static PetscErrorCode PetscDualSpaceLagrangeGetUseMoments_Sum(PetscDualSpace sp, PetscBool *value)
1164*2dce792eSToby Isaac {
1165*2dce792eSToby Isaac   PetscFunctionBegin;
1166*2dce792eSToby Isaac   PetscDualSpaceSumPassthrough(sp, PetscDualSpaceLagrangeGetUseMoments, value);
1167*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1168*2dce792eSToby Isaac }
1169*2dce792eSToby Isaac 
1170*2dce792eSToby Isaac static PetscErrorCode PetscDualSpaceLagrangeSetUseMoments_Sum(PetscDualSpace sp, PetscBool value)
1171*2dce792eSToby Isaac {
1172*2dce792eSToby Isaac   PetscFunctionBegin;
1173*2dce792eSToby Isaac   PetscDualSpaceSumPassthrough(sp, PetscDualSpaceLagrangeSetUseMoments, value);
1174*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1175*2dce792eSToby Isaac }
1176*2dce792eSToby Isaac 
1177*2dce792eSToby Isaac static PetscErrorCode PetscDualSpaceLagrangeGetMomentOrder_Sum(PetscDualSpace sp, PetscInt *value)
1178*2dce792eSToby Isaac {
1179*2dce792eSToby Isaac   PetscFunctionBegin;
1180*2dce792eSToby Isaac   PetscDualSpaceSumPassthrough(sp, PetscDualSpaceLagrangeGetMomentOrder, value);
1181*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1182*2dce792eSToby Isaac }
1183*2dce792eSToby Isaac 
1184*2dce792eSToby Isaac static PetscErrorCode PetscDualSpaceLagrangeSetMomentOrder_Sum(PetscDualSpace sp, PetscInt value)
1185*2dce792eSToby Isaac {
1186*2dce792eSToby Isaac   PetscFunctionBegin;
1187*2dce792eSToby Isaac   PetscDualSpaceSumPassthrough(sp, PetscDualSpaceLagrangeSetMomentOrder, value);
1188*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1189*2dce792eSToby Isaac }
1190*2dce792eSToby Isaac 
1191*2dce792eSToby Isaac static PetscErrorCode PetscDualSpaceLagrangeGetNodeType_Sum(PetscDualSpace sp, PetscDTNodeType *node_type, PetscBool *include_endpoints, PetscReal *exponent)
1192*2dce792eSToby Isaac {
1193*2dce792eSToby Isaac   PetscFunctionBegin;
1194*2dce792eSToby Isaac   PetscDualSpaceSumPassthrough(sp, PetscDualSpaceLagrangeGetNodeType, node_type, include_endpoints, exponent);
1195*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1196*2dce792eSToby Isaac }
1197*2dce792eSToby Isaac 
1198*2dce792eSToby Isaac static PetscErrorCode PetscDualSpaceLagrangeSetNodeType_Sum(PetscDualSpace sp, PetscDTNodeType node_type, PetscBool include_endpoints, PetscReal exponent)
1199*2dce792eSToby Isaac {
1200*2dce792eSToby Isaac   PetscFunctionBegin;
1201*2dce792eSToby Isaac   PetscDualSpaceSumPassthrough(sp, PetscDualSpaceLagrangeSetNodeType, node_type, include_endpoints, exponent);
1202*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1203*2dce792eSToby Isaac }
1204*2dce792eSToby Isaac 
1205*2dce792eSToby Isaac static PetscErrorCode PetscDualSpaceInitialize_Sum(PetscDualSpace sp)
1206*2dce792eSToby Isaac {
1207*2dce792eSToby Isaac   PetscFunctionBegin;
1208*2dce792eSToby Isaac   sp->ops->destroy              = PetscDualSpaceDestroy_Sum;
1209*2dce792eSToby Isaac   sp->ops->view                 = PetscDualSpaceView_Sum;
1210*2dce792eSToby Isaac   sp->ops->setfromoptions       = NULL;
1211*2dce792eSToby Isaac   sp->ops->duplicate            = PetscDualSpaceDuplicate_Sum;
1212*2dce792eSToby Isaac   sp->ops->setup                = PetscDualSpaceSetUp_Sum;
1213*2dce792eSToby Isaac   sp->ops->createheightsubspace = NULL;
1214*2dce792eSToby Isaac   sp->ops->createpointsubspace  = NULL;
1215*2dce792eSToby Isaac   sp->ops->getsymmetries        = PetscDualSpaceGetSymmetries_Sum;
1216*2dce792eSToby Isaac   sp->ops->apply                = PetscDualSpaceApplyDefault;
1217*2dce792eSToby Isaac   sp->ops->applyall             = PetscDualSpaceApplyAllDefault;
1218*2dce792eSToby Isaac   sp->ops->applyint             = PetscDualSpaceApplyInteriorDefault;
1219*2dce792eSToby Isaac   sp->ops->createalldata        = PetscDualSpaceCreateAllDataDefault;
1220*2dce792eSToby Isaac   sp->ops->createintdata        = PetscDualSpaceCreateInteriorDataDefault;
1221*2dce792eSToby Isaac 
1222*2dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumGetNumSubspaces_C", PetscDualSpaceSumGetNumSubspaces_Sum));
1223*2dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumSetNumSubspaces_C", PetscDualSpaceSumSetNumSubspaces_Sum));
1224*2dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumGetSubspace_C", PetscDualSpaceSumGetSubspace_Sum));
1225*2dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumSetSubspace_C", PetscDualSpaceSumSetSubspace_Sum));
1226*2dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumGetConcatenate_C", PetscDualSpaceSumGetConcatenate_Sum));
1227*2dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumSetConcatenate_C", PetscDualSpaceSumSetConcatenate_Sum));
1228*2dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumGetInterleave_C", PetscDualSpaceSumGetInterleave_Sum));
1229*2dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumSetInterleave_C", PetscDualSpaceSumSetInterleave_Sum));
1230*2dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetContinuity_C", PetscDualSpaceLagrangeGetContinuity_Sum));
1231*2dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetContinuity_C", PetscDualSpaceLagrangeSetContinuity_Sum));
1232*2dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetMomentOrder_C", PetscDualSpaceLagrangeGetMomentOrder_Sum));
1233*2dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetMomentOrder_C", PetscDualSpaceLagrangeSetMomentOrder_Sum));
1234*2dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetNodeType_C", PetscDualSpaceLagrangeGetNodeType_Sum));
1235*2dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetNodeType_C", PetscDualSpaceLagrangeSetNodeType_Sum));
1236*2dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetTensor_C", PetscDualSpaceLagrangeGetTensor_Sum));
1237*2dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetTensor_C", PetscDualSpaceLagrangeSetTensor_Sum));
1238*2dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetTrimmed_C", PetscDualSpaceLagrangeGetTrimmed_Sum));
1239*2dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetTrimmed_C", PetscDualSpaceLagrangeSetTrimmed_Sum));
1240*2dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetUseMoments_C", PetscDualSpaceLagrangeGetUseMoments_Sum));
1241*2dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetUseMoments_C", PetscDualSpaceLagrangeSetUseMoments_Sum));
1242*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1243*2dce792eSToby Isaac }
1244*2dce792eSToby Isaac 
1245*2dce792eSToby Isaac /*MC
1246*2dce792eSToby Isaac   PETSCDUALSPACESUM = "sum" - A `PetscDualSpace` object that encapsulates a sum of subspaces.
1247*2dce792eSToby Isaac 
1248*2dce792eSToby Isaac   Level: intermediate
1249*2dce792eSToby Isaac 
1250*2dce792eSToby Isaac   Note:
1251*2dce792eSToby Isaac   That sum can either be direct or a concatenation. For example if A and B are spaces each with 2 components,
1252*2dce792eSToby Isaac   the direct sum of A and B will also have 2 components while the concatenated sum will have 4 components. In both cases A and B must be defined over the
1253*2dce792eSToby Isaac   same reference element.
1254*2dce792eSToby Isaac 
1255*2dce792eSToby Isaac .seealso: `PetscDualSpace`, `PetscDualSpaceType`, `PetscDualSpaceCreate()`, `PetscDualSpaceSetType()`, `PetscDualSpaceSumGetNumSubspaces()`, `PetscDualSpaceSumSetNumSubspaces()`,
1256*2dce792eSToby Isaac           `PetscDualSpaceSumGetConcatenate()`, `PetscDualSpaceSumSetConcatenate()`, `PetscDualSpaceSumSetInterleave()`, `PetscDualSpaceSumGetInterleave()`
1257*2dce792eSToby Isaac M*/
1258*2dce792eSToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceCreate_Sum(PetscDualSpace sp)
1259*2dce792eSToby Isaac {
1260*2dce792eSToby Isaac   PetscDualSpace_Sum *sum;
1261*2dce792eSToby Isaac 
1262*2dce792eSToby Isaac   PetscFunctionBegin;
1263*2dce792eSToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1264*2dce792eSToby Isaac   PetscCall(PetscNew(&sum));
1265*2dce792eSToby Isaac   sum->numSumSpaces = PETSC_DEFAULT;
1266*2dce792eSToby Isaac   sp->data          = sum;
1267*2dce792eSToby Isaac   sp->k             = PETSC_FORM_DEGREE_UNDEFINED;
1268*2dce792eSToby Isaac   PetscCall(PetscDualSpaceInitialize_Sum(sp));
1269*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1270*2dce792eSToby Isaac }
1271*2dce792eSToby Isaac 
1272*2dce792eSToby Isaac /*@
1273*2dce792eSToby Isaac   PetscDualSpaceCreateSum - Create a finite element dual basis that is the sum of other dual bases
1274*2dce792eSToby Isaac 
1275*2dce792eSToby Isaac   Collective
1276*2dce792eSToby Isaac 
1277*2dce792eSToby Isaac   Input Parameters:
1278*2dce792eSToby Isaac + numSubspaces - the number of spaces that will be added together
1279*2dce792eSToby Isaac . subspaces    - an array of length `numSubspaces` of spaces
1280*2dce792eSToby Isaac - concatenate  - if `PETSC_FALSE`, the sum-space has the same components as the individual dual spaces (`PetscDualSpaceGetNumComponents()`); if `PETSC_TRUE`, the individual components are concatenated to create a dual space with more components
1281*2dce792eSToby Isaac 
1282*2dce792eSToby Isaac   Output Parameter:
1283*2dce792eSToby Isaac . sumSpace - a `PetscDualSpace` of type `PETSCDUALSPACESUM`
1284*2dce792eSToby Isaac 
1285*2dce792eSToby Isaac   Level: advanced
1286*2dce792eSToby Isaac 
1287*2dce792eSToby Isaac .seealso: `PetscDualSpace`, `PETSCDUALSPACESUM`, `PETSCSPACESUM`
1288*2dce792eSToby Isaac @*/
1289*2dce792eSToby Isaac PetscErrorCode PetscDualSpaceCreateSum(PetscInt numSubspaces, const PetscDualSpace subspaces[], PetscBool concatenate, PetscDualSpace *sumSpace)
1290*2dce792eSToby Isaac {
1291*2dce792eSToby Isaac   PetscInt i, Nc = 0;
1292*2dce792eSToby Isaac 
1293*2dce792eSToby Isaac   PetscFunctionBegin;
1294*2dce792eSToby Isaac   PetscCall(PetscDualSpaceCreate(PetscObjectComm((PetscObject)subspaces[0]), sumSpace));
1295*2dce792eSToby Isaac   PetscCall(PetscDualSpaceSetType(*sumSpace, PETSCDUALSPACESUM));
1296*2dce792eSToby Isaac   PetscCall(PetscDualSpaceSumSetNumSubspaces(*sumSpace, numSubspaces));
1297*2dce792eSToby Isaac   PetscCall(PetscDualSpaceSumSetConcatenate(*sumSpace, concatenate));
1298*2dce792eSToby Isaac   for (i = 0; i < numSubspaces; ++i) {
1299*2dce792eSToby Isaac     PetscInt sNc;
1300*2dce792eSToby Isaac 
1301*2dce792eSToby Isaac     PetscCall(PetscDualSpaceSumSetSubspace(*sumSpace, i, subspaces[i]));
1302*2dce792eSToby Isaac     PetscCall(PetscDualSpaceGetNumComponents(subspaces[i], &sNc));
1303*2dce792eSToby Isaac     if (concatenate) Nc += sNc;
1304*2dce792eSToby Isaac     else Nc = sNc;
1305*2dce792eSToby Isaac 
1306*2dce792eSToby Isaac     if (i == 0) {
1307*2dce792eSToby Isaac       DM dm;
1308*2dce792eSToby Isaac 
1309*2dce792eSToby Isaac       PetscCall(PetscDualSpaceGetDM(subspaces[i], &dm));
1310*2dce792eSToby Isaac       PetscCall(PetscDualSpaceSetDM(*sumSpace, dm));
1311*2dce792eSToby Isaac     }
1312*2dce792eSToby Isaac   }
1313*2dce792eSToby Isaac   PetscCall(PetscDualSpaceSetNumComponents(*sumSpace, Nc));
1314*2dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1315*2dce792eSToby Isaac }
1316