xref: /petsc/src/dm/dt/dualspace/impls/lagrange/dspacelagrange.c (revision 20cf1dd8cd656e27510885a71ea4be028bf48813)
1*20cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
2*20cf1dd8SToby Isaac #include <petscdmplex.h>
3*20cf1dd8SToby Isaac 
4*20cf1dd8SToby Isaac static PetscErrorCode PetscDualSpaceLagrangeGetTensor_Lagrange(PetscDualSpace sp, PetscBool *tensor)
5*20cf1dd8SToby Isaac {
6*20cf1dd8SToby Isaac   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data;
7*20cf1dd8SToby Isaac 
8*20cf1dd8SToby Isaac   PetscFunctionBegin;
9*20cf1dd8SToby Isaac   *tensor = lag->tensorSpace;
10*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
11*20cf1dd8SToby Isaac }
12*20cf1dd8SToby Isaac 
13*20cf1dd8SToby Isaac static PetscErrorCode PetscDualSpaceLagrangeSetTensor_Lagrange(PetscDualSpace sp, PetscBool tensor)
14*20cf1dd8SToby Isaac {
15*20cf1dd8SToby Isaac   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data;
16*20cf1dd8SToby Isaac 
17*20cf1dd8SToby Isaac   PetscFunctionBegin;
18*20cf1dd8SToby Isaac   lag->tensorSpace = tensor;
19*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
20*20cf1dd8SToby Isaac }
21*20cf1dd8SToby Isaac 
22*20cf1dd8SToby Isaac /*
23*20cf1dd8SToby Isaac   LatticePointLexicographic_Internal - Returns all tuples of size 'len' with nonnegative integers that sum up to at most 'max'.
24*20cf1dd8SToby Isaac                                        Ordering is lexicographic with lowest index as least significant in ordering.
25*20cf1dd8SToby Isaac                                        e.g. for len == 2 and max == 2, this will return, in order, {0,0}, {1,0}, {2,0}, {0,1}, {1,1}, {2,0}.
26*20cf1dd8SToby Isaac 
27*20cf1dd8SToby Isaac   Input Parameters:
28*20cf1dd8SToby Isaac + len - The length of the tuple
29*20cf1dd8SToby Isaac . max - The maximum sum
30*20cf1dd8SToby Isaac - tup - A tuple of length len+1: tup[len] > 0 indicates a stopping condition
31*20cf1dd8SToby Isaac 
32*20cf1dd8SToby Isaac   Output Parameter:
33*20cf1dd8SToby Isaac . tup - A tuple of len integers whos sum is at most 'max'
34*20cf1dd8SToby Isaac */
35*20cf1dd8SToby Isaac static PetscErrorCode LatticePointLexicographic_Internal(PetscInt len, PetscInt max, PetscInt tup[])
36*20cf1dd8SToby Isaac {
37*20cf1dd8SToby Isaac   PetscFunctionBegin;
38*20cf1dd8SToby Isaac   while (len--) {
39*20cf1dd8SToby Isaac     max -= tup[len];
40*20cf1dd8SToby Isaac     if (!max) {
41*20cf1dd8SToby Isaac       tup[len] = 0;
42*20cf1dd8SToby Isaac       break;
43*20cf1dd8SToby Isaac     }
44*20cf1dd8SToby Isaac   }
45*20cf1dd8SToby Isaac   tup[++len]++;
46*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
47*20cf1dd8SToby Isaac }
48*20cf1dd8SToby Isaac 
49*20cf1dd8SToby Isaac /*
50*20cf1dd8SToby Isaac   TensorPointLexicographic_Internal - Returns all tuples of size 'len' with nonnegative integers that are all less than or equal to 'max'.
51*20cf1dd8SToby Isaac                                       Ordering is lexicographic with lowest index as least significant in ordering.
52*20cf1dd8SToby Isaac                                       e.g. for len == 2 and max == 2, this will return, in order, {0,0}, {1,0}, {2,0}, {0,1}, {1,1}, {2,1}, {0,2}, {1,2}, {2,2}.
53*20cf1dd8SToby Isaac 
54*20cf1dd8SToby Isaac   Input Parameters:
55*20cf1dd8SToby Isaac + len - The length of the tuple
56*20cf1dd8SToby Isaac . max - The maximum value
57*20cf1dd8SToby Isaac - tup - A tuple of length len+1: tup[len] > 0 indicates a stopping condition
58*20cf1dd8SToby Isaac 
59*20cf1dd8SToby Isaac   Output Parameter:
60*20cf1dd8SToby Isaac . tup - A tuple of len integers whos sum is at most 'max'
61*20cf1dd8SToby Isaac */
62*20cf1dd8SToby Isaac static PetscErrorCode TensorPointLexicographic_Internal(PetscInt len, PetscInt max, PetscInt tup[])
63*20cf1dd8SToby Isaac {
64*20cf1dd8SToby Isaac   PetscInt       i;
65*20cf1dd8SToby Isaac 
66*20cf1dd8SToby Isaac   PetscFunctionBegin;
67*20cf1dd8SToby Isaac   for (i = 0; i < len; i++) {
68*20cf1dd8SToby Isaac     if (tup[i] < max) {
69*20cf1dd8SToby Isaac       break;
70*20cf1dd8SToby Isaac     } else {
71*20cf1dd8SToby Isaac       tup[i] = 0;
72*20cf1dd8SToby Isaac     }
73*20cf1dd8SToby Isaac   }
74*20cf1dd8SToby Isaac   tup[i]++;
75*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
76*20cf1dd8SToby Isaac }
77*20cf1dd8SToby Isaac 
78*20cf1dd8SToby Isaac 
79*20cf1dd8SToby Isaac #define BaryIndex(perEdge,a,b,c) (((b)*(2*perEdge+1-(b)))/2)+(c)
80*20cf1dd8SToby Isaac 
81*20cf1dd8SToby Isaac #define CartIndex(perEdge,a,b) (perEdge*(a)+b)
82*20cf1dd8SToby Isaac 
83*20cf1dd8SToby Isaac static PetscErrorCode PetscDualSpaceGetSymmetries_Lagrange(PetscDualSpace sp, const PetscInt ****perms, const PetscScalar ****flips)
84*20cf1dd8SToby Isaac {
85*20cf1dd8SToby Isaac 
86*20cf1dd8SToby Isaac   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data;
87*20cf1dd8SToby Isaac   PetscInt           dim, order, p, Nc;
88*20cf1dd8SToby Isaac   PetscErrorCode     ierr;
89*20cf1dd8SToby Isaac 
90*20cf1dd8SToby Isaac   PetscFunctionBegin;
91*20cf1dd8SToby Isaac   ierr = PetscDualSpaceGetOrder(sp,&order);CHKERRQ(ierr);
92*20cf1dd8SToby Isaac   ierr = PetscDualSpaceGetNumComponents(sp,&Nc);CHKERRQ(ierr);
93*20cf1dd8SToby Isaac   ierr = DMGetDimension(sp->dm,&dim);CHKERRQ(ierr);
94*20cf1dd8SToby Isaac   if (!dim || !lag->continuous || order < 3) PetscFunctionReturn(0);
95*20cf1dd8SToby Isaac   if (dim > 3) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Lagrange symmetries not implemented for dim = %D > 3",dim);
96*20cf1dd8SToby Isaac   if (!lag->symmetries) { /* store symmetries */
97*20cf1dd8SToby Isaac     PetscDualSpace hsp;
98*20cf1dd8SToby Isaac     DM             K;
99*20cf1dd8SToby Isaac     PetscInt       numPoints = 1, d;
100*20cf1dd8SToby Isaac     PetscInt       numFaces;
101*20cf1dd8SToby Isaac     PetscInt       ***symmetries;
102*20cf1dd8SToby Isaac     const PetscInt ***hsymmetries;
103*20cf1dd8SToby Isaac 
104*20cf1dd8SToby Isaac     if (lag->simplexCell) {
105*20cf1dd8SToby Isaac       numFaces = 1 + dim;
106*20cf1dd8SToby Isaac       for (d = 0; d < dim; d++) numPoints = numPoints * 2 + 1;
107*20cf1dd8SToby Isaac     }
108*20cf1dd8SToby Isaac     else {
109*20cf1dd8SToby Isaac       numPoints = PetscPowInt(3,dim);
110*20cf1dd8SToby Isaac       numFaces  = 2 * dim;
111*20cf1dd8SToby Isaac     }
112*20cf1dd8SToby Isaac     ierr = PetscCalloc1(numPoints,&symmetries);CHKERRQ(ierr);
113*20cf1dd8SToby Isaac     if (0 < dim && dim < 3) { /* compute self symmetries */
114*20cf1dd8SToby Isaac       PetscInt **cellSymmetries;
115*20cf1dd8SToby Isaac 
116*20cf1dd8SToby Isaac       lag->numSelfSym = 2 * numFaces;
117*20cf1dd8SToby Isaac       lag->selfSymOff = numFaces;
118*20cf1dd8SToby Isaac       ierr = PetscCalloc1(2*numFaces,&cellSymmetries);CHKERRQ(ierr);
119*20cf1dd8SToby Isaac       /* we want to be able to index symmetries directly with the orientations, which range from [-numFaces,numFaces) */
120*20cf1dd8SToby Isaac       symmetries[0] = &cellSymmetries[numFaces];
121*20cf1dd8SToby Isaac       if (dim == 1) {
122*20cf1dd8SToby Isaac         PetscInt dofPerEdge = order - 1;
123*20cf1dd8SToby Isaac 
124*20cf1dd8SToby Isaac         if (dofPerEdge > 1) {
125*20cf1dd8SToby Isaac           PetscInt i, j, *reverse;
126*20cf1dd8SToby Isaac 
127*20cf1dd8SToby Isaac           ierr = PetscMalloc1(dofPerEdge*Nc,&reverse);CHKERRQ(ierr);
128*20cf1dd8SToby Isaac           for (i = 0; i < dofPerEdge; i++) {
129*20cf1dd8SToby Isaac             for (j = 0; j < Nc; j++) {
130*20cf1dd8SToby Isaac               reverse[i*Nc + j] = Nc * (dofPerEdge - 1 - i) + j;
131*20cf1dd8SToby Isaac             }
132*20cf1dd8SToby Isaac           }
133*20cf1dd8SToby Isaac           symmetries[0][-2] = reverse;
134*20cf1dd8SToby Isaac 
135*20cf1dd8SToby Isaac           /* yes, this is redundant, but it makes it easier to cleanup if I don't have to worry about what not to free */
136*20cf1dd8SToby Isaac           ierr = PetscMalloc1(dofPerEdge*Nc,&reverse);CHKERRQ(ierr);
137*20cf1dd8SToby Isaac           for (i = 0; i < dofPerEdge; i++) {
138*20cf1dd8SToby Isaac             for (j = 0; j < Nc; j++) {
139*20cf1dd8SToby Isaac               reverse[i*Nc + j] = Nc * (dofPerEdge - 1 - i) + j;
140*20cf1dd8SToby Isaac             }
141*20cf1dd8SToby Isaac           }
142*20cf1dd8SToby Isaac           symmetries[0][1] = reverse;
143*20cf1dd8SToby Isaac         }
144*20cf1dd8SToby Isaac       } else {
145*20cf1dd8SToby Isaac         PetscInt dofPerEdge = lag->simplexCell ? (order - 2) : (order - 1), s;
146*20cf1dd8SToby Isaac         PetscInt dofPerFace;
147*20cf1dd8SToby Isaac 
148*20cf1dd8SToby Isaac         if (dofPerEdge > 1) {
149*20cf1dd8SToby Isaac           for (s = -numFaces; s < numFaces; s++) {
150*20cf1dd8SToby Isaac             PetscInt *sym, i, j, k, l;
151*20cf1dd8SToby Isaac 
152*20cf1dd8SToby Isaac             if (!s) continue;
153*20cf1dd8SToby Isaac             if (lag->simplexCell) {
154*20cf1dd8SToby Isaac               dofPerFace = (dofPerEdge * (dofPerEdge + 1))/2;
155*20cf1dd8SToby Isaac               ierr = PetscMalloc1(Nc*dofPerFace,&sym);CHKERRQ(ierr);
156*20cf1dd8SToby Isaac               for (j = 0, l = 0; j < dofPerEdge; j++) {
157*20cf1dd8SToby Isaac                 for (k = 0; k < dofPerEdge - j; k++, l++) {
158*20cf1dd8SToby Isaac                   i = dofPerEdge - 1 - j - k;
159*20cf1dd8SToby Isaac                   switch (s) {
160*20cf1dd8SToby Isaac                   case -3:
161*20cf1dd8SToby Isaac                     sym[Nc*l] = BaryIndex(dofPerEdge,i,k,j);
162*20cf1dd8SToby Isaac                     break;
163*20cf1dd8SToby Isaac                   case -2:
164*20cf1dd8SToby Isaac                     sym[Nc*l] = BaryIndex(dofPerEdge,j,i,k);
165*20cf1dd8SToby Isaac                     break;
166*20cf1dd8SToby Isaac                   case -1:
167*20cf1dd8SToby Isaac                     sym[Nc*l] = BaryIndex(dofPerEdge,k,j,i);
168*20cf1dd8SToby Isaac                     break;
169*20cf1dd8SToby Isaac                   case 1:
170*20cf1dd8SToby Isaac                     sym[Nc*l] = BaryIndex(dofPerEdge,k,i,j);
171*20cf1dd8SToby Isaac                     break;
172*20cf1dd8SToby Isaac                   case 2:
173*20cf1dd8SToby Isaac                     sym[Nc*l] = BaryIndex(dofPerEdge,j,k,i);
174*20cf1dd8SToby Isaac                     break;
175*20cf1dd8SToby Isaac                   }
176*20cf1dd8SToby Isaac                 }
177*20cf1dd8SToby Isaac               }
178*20cf1dd8SToby Isaac             } else {
179*20cf1dd8SToby Isaac               dofPerFace = dofPerEdge * dofPerEdge;
180*20cf1dd8SToby Isaac               ierr = PetscMalloc1(Nc*dofPerFace,&sym);CHKERRQ(ierr);
181*20cf1dd8SToby Isaac               for (j = 0, l = 0; j < dofPerEdge; j++) {
182*20cf1dd8SToby Isaac                 for (k = 0; k < dofPerEdge; k++, l++) {
183*20cf1dd8SToby Isaac                   switch (s) {
184*20cf1dd8SToby Isaac                   case -4:
185*20cf1dd8SToby Isaac                     sym[Nc*l] = CartIndex(dofPerEdge,k,j);
186*20cf1dd8SToby Isaac                     break;
187*20cf1dd8SToby Isaac                   case -3:
188*20cf1dd8SToby Isaac                     sym[Nc*l] = CartIndex(dofPerEdge,(dofPerEdge - 1 - j),k);
189*20cf1dd8SToby Isaac                     break;
190*20cf1dd8SToby Isaac                   case -2:
191*20cf1dd8SToby Isaac                     sym[Nc*l] = CartIndex(dofPerEdge,(dofPerEdge - 1 - k),(dofPerEdge - 1 - j));
192*20cf1dd8SToby Isaac                     break;
193*20cf1dd8SToby Isaac                   case -1:
194*20cf1dd8SToby Isaac                     sym[Nc*l] = CartIndex(dofPerEdge,j,(dofPerEdge - 1 - k));
195*20cf1dd8SToby Isaac                     break;
196*20cf1dd8SToby Isaac                   case 1:
197*20cf1dd8SToby Isaac                     sym[Nc*l] = CartIndex(dofPerEdge,(dofPerEdge - 1 - k),j);
198*20cf1dd8SToby Isaac                     break;
199*20cf1dd8SToby Isaac                   case 2:
200*20cf1dd8SToby Isaac                     sym[Nc*l] = CartIndex(dofPerEdge,(dofPerEdge - 1 - j),(dofPerEdge - 1 - k));
201*20cf1dd8SToby Isaac                     break;
202*20cf1dd8SToby Isaac                   case 3:
203*20cf1dd8SToby Isaac                     sym[Nc*l] = CartIndex(dofPerEdge,k,(dofPerEdge - 1 - j));
204*20cf1dd8SToby Isaac                     break;
205*20cf1dd8SToby Isaac                   }
206*20cf1dd8SToby Isaac                 }
207*20cf1dd8SToby Isaac               }
208*20cf1dd8SToby Isaac             }
209*20cf1dd8SToby Isaac             for (i = 0; i < dofPerFace; i++) {
210*20cf1dd8SToby Isaac               sym[Nc*i] *= Nc;
211*20cf1dd8SToby Isaac               for (j = 1; j < Nc; j++) {
212*20cf1dd8SToby Isaac                 sym[Nc*i+j] = sym[Nc*i] + j;
213*20cf1dd8SToby Isaac               }
214*20cf1dd8SToby Isaac             }
215*20cf1dd8SToby Isaac             symmetries[0][s] = sym;
216*20cf1dd8SToby Isaac           }
217*20cf1dd8SToby Isaac         }
218*20cf1dd8SToby Isaac       }
219*20cf1dd8SToby Isaac     }
220*20cf1dd8SToby Isaac     ierr = PetscDualSpaceGetHeightSubspace(sp,1,&hsp);CHKERRQ(ierr);
221*20cf1dd8SToby Isaac     ierr = PetscDualSpaceGetSymmetries(hsp,&hsymmetries,NULL);CHKERRQ(ierr);
222*20cf1dd8SToby Isaac     if (hsymmetries) {
223*20cf1dd8SToby Isaac       PetscBool      *seen;
224*20cf1dd8SToby Isaac       const PetscInt *cone;
225*20cf1dd8SToby Isaac       PetscInt       KclosureSize, *Kclosure = NULL;
226*20cf1dd8SToby Isaac 
227*20cf1dd8SToby Isaac       ierr = PetscDualSpaceGetDM(sp,&K);CHKERRQ(ierr);
228*20cf1dd8SToby Isaac       ierr = PetscCalloc1(numPoints,&seen);CHKERRQ(ierr);
229*20cf1dd8SToby Isaac       ierr = DMPlexGetCone(K,0,&cone);CHKERRQ(ierr);
230*20cf1dd8SToby Isaac       ierr = DMPlexGetTransitiveClosure(K,0,PETSC_TRUE,&KclosureSize,&Kclosure);CHKERRQ(ierr);
231*20cf1dd8SToby Isaac       for (p = 0; p < numFaces; p++) {
232*20cf1dd8SToby Isaac         PetscInt closureSize, *closure = NULL, q;
233*20cf1dd8SToby Isaac 
234*20cf1dd8SToby Isaac         ierr = DMPlexGetTransitiveClosure(K,cone[p],PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
235*20cf1dd8SToby Isaac         for (q = 0; q < closureSize; q++) {
236*20cf1dd8SToby Isaac           PetscInt point = closure[2*q], r;
237*20cf1dd8SToby Isaac 
238*20cf1dd8SToby Isaac           if(!seen[point]) {
239*20cf1dd8SToby Isaac             for (r = 0; r < KclosureSize; r++) {
240*20cf1dd8SToby Isaac               if (Kclosure[2 * r] == point) break;
241*20cf1dd8SToby Isaac             }
242*20cf1dd8SToby Isaac             seen[point] = PETSC_TRUE;
243*20cf1dd8SToby Isaac             symmetries[r] = (PetscInt **) hsymmetries[q];
244*20cf1dd8SToby Isaac           }
245*20cf1dd8SToby Isaac         }
246*20cf1dd8SToby Isaac         ierr = DMPlexRestoreTransitiveClosure(K,cone[p],PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
247*20cf1dd8SToby Isaac       }
248*20cf1dd8SToby Isaac       ierr = DMPlexRestoreTransitiveClosure(K,0,PETSC_TRUE,&KclosureSize,&Kclosure);CHKERRQ(ierr);
249*20cf1dd8SToby Isaac       ierr = PetscFree(seen);CHKERRQ(ierr);
250*20cf1dd8SToby Isaac     }
251*20cf1dd8SToby Isaac     lag->symmetries = symmetries;
252*20cf1dd8SToby Isaac   }
253*20cf1dd8SToby Isaac   if (perms) *perms = (const PetscInt ***) lag->symmetries;
254*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
255*20cf1dd8SToby Isaac }
256*20cf1dd8SToby Isaac 
257*20cf1dd8SToby Isaac /*@C
258*20cf1dd8SToby Isaac   PetscDualSpaceGetSymmetries - Returns a description of the symmetries of this basis
259*20cf1dd8SToby Isaac 
260*20cf1dd8SToby Isaac   Not collective
261*20cf1dd8SToby Isaac 
262*20cf1dd8SToby Isaac   Input Parameter:
263*20cf1dd8SToby Isaac . sp - the PetscDualSpace object
264*20cf1dd8SToby Isaac 
265*20cf1dd8SToby Isaac   Output Parameters:
266*20cf1dd8SToby Isaac + perms - Permutations of the local degrees of freedom, parameterized by the point orientation
267*20cf1dd8SToby Isaac - flips - Sign reversal of the local degrees of freedom, parameterized by the point orientation
268*20cf1dd8SToby Isaac 
269*20cf1dd8SToby Isaac   Note: The permutation and flip arrays are organized in the following way
270*20cf1dd8SToby Isaac $ perms[p][ornt][dof # on point] = new local dof #
271*20cf1dd8SToby Isaac $ flips[p][ornt][dof # on point] = reversal or not
272*20cf1dd8SToby Isaac 
273*20cf1dd8SToby Isaac   Level: developer
274*20cf1dd8SToby Isaac 
275*20cf1dd8SToby Isaac .seealso: PetscDualSpaceSetSymmetries()
276*20cf1dd8SToby Isaac @*/
277*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetSymmetries(PetscDualSpace sp, const PetscInt ****perms, const PetscScalar ****flips)
278*20cf1dd8SToby Isaac {
279*20cf1dd8SToby Isaac   PetscErrorCode ierr;
280*20cf1dd8SToby Isaac 
281*20cf1dd8SToby Isaac   PetscFunctionBegin;
282*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp,PETSCDUALSPACE_CLASSID,1);
283*20cf1dd8SToby Isaac   if (perms) {
284*20cf1dd8SToby Isaac     PetscValidPointer(perms,2);
285*20cf1dd8SToby Isaac     *perms = NULL;
286*20cf1dd8SToby Isaac   }
287*20cf1dd8SToby Isaac   if (flips) {
288*20cf1dd8SToby Isaac     PetscValidPointer(flips,3);
289*20cf1dd8SToby Isaac     *flips = NULL;
290*20cf1dd8SToby Isaac   }
291*20cf1dd8SToby Isaac   if (sp->ops->getsymmetries) {
292*20cf1dd8SToby Isaac     ierr = (sp->ops->getsymmetries)(sp,perms,flips);CHKERRQ(ierr);
293*20cf1dd8SToby Isaac   }
294*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
295*20cf1dd8SToby Isaac }
296*20cf1dd8SToby Isaac 
297*20cf1dd8SToby Isaac static PetscErrorCode PetscDualSpaceLagrangeView_Ascii(PetscDualSpace sp, PetscViewer viewer)
298*20cf1dd8SToby Isaac {
299*20cf1dd8SToby Isaac   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data;
300*20cf1dd8SToby Isaac   PetscErrorCode      ierr;
301*20cf1dd8SToby Isaac 
302*20cf1dd8SToby Isaac   PetscFunctionBegin;
303*20cf1dd8SToby Isaac   ierr = PetscViewerASCIIPrintf(viewer, "%s %sLagrange dual space of order %D", lag->continuous ? "Continuous" : "Discontinuous", lag->tensorSpace ? "Tensor " : "", sp->order, sp->Nc);CHKERRQ(ierr);
304*20cf1dd8SToby Isaac   if (sp->Nc > 1) {ierr = PetscViewerASCIIPrintf(viewer, " with %D components", sp->Nc);CHKERRQ(ierr);}
305*20cf1dd8SToby Isaac   ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr);
306*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
307*20cf1dd8SToby Isaac }
308*20cf1dd8SToby Isaac 
309*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceView_Lagrange(PetscDualSpace sp, PetscViewer viewer)
310*20cf1dd8SToby Isaac {
311*20cf1dd8SToby Isaac   PetscBool      iascii;
312*20cf1dd8SToby Isaac   PetscErrorCode ierr;
313*20cf1dd8SToby Isaac 
314*20cf1dd8SToby Isaac   PetscFunctionBegin;
315*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
316*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
317*20cf1dd8SToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
318*20cf1dd8SToby Isaac   if (iascii) {ierr = PetscDualSpaceLagrangeView_Ascii(sp, viewer);CHKERRQ(ierr);}
319*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
320*20cf1dd8SToby Isaac }
321*20cf1dd8SToby Isaac 
322*20cf1dd8SToby Isaac static PetscErrorCode PetscDualSpaceGetDimension_SingleCell_Lagrange(PetscDualSpace sp, PetscInt order, PetscInt *dim)
323*20cf1dd8SToby Isaac {
324*20cf1dd8SToby Isaac   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data;
325*20cf1dd8SToby Isaac   PetscReal           D   = 1.0;
326*20cf1dd8SToby Isaac   PetscInt            n, i;
327*20cf1dd8SToby Isaac   PetscErrorCode      ierr;
328*20cf1dd8SToby Isaac 
329*20cf1dd8SToby Isaac   PetscFunctionBegin;
330*20cf1dd8SToby Isaac   *dim = -1;                    /* Ensure that the compiler knows *dim is set. */
331*20cf1dd8SToby Isaac   ierr = DMGetDimension(sp->dm, &n);CHKERRQ(ierr);
332*20cf1dd8SToby Isaac   if (!lag->tensorSpace) {
333*20cf1dd8SToby Isaac     for (i = 1; i <= n; ++i) {
334*20cf1dd8SToby Isaac       D *= ((PetscReal) (order+i))/i;
335*20cf1dd8SToby Isaac     }
336*20cf1dd8SToby Isaac     *dim = (PetscInt) (D + 0.5);
337*20cf1dd8SToby Isaac   } else {
338*20cf1dd8SToby Isaac     *dim = 1;
339*20cf1dd8SToby Isaac     for (i = 0; i < n; ++i) *dim *= (order+1);
340*20cf1dd8SToby Isaac   }
341*20cf1dd8SToby Isaac   *dim *= sp->Nc;
342*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
343*20cf1dd8SToby Isaac }
344*20cf1dd8SToby Isaac 
345*20cf1dd8SToby Isaac static PetscErrorCode PetscDualSpaceCreateHeightSubspace_Lagrange(PetscDualSpace sp, PetscInt height, PetscDualSpace *bdsp)
346*20cf1dd8SToby Isaac {
347*20cf1dd8SToby Isaac   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data;
348*20cf1dd8SToby Isaac   PetscBool          continuous, tensor;
349*20cf1dd8SToby Isaac   PetscInt           order;
350*20cf1dd8SToby Isaac   PetscErrorCode     ierr;
351*20cf1dd8SToby Isaac 
352*20cf1dd8SToby Isaac   PetscFunctionBegin;
353*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
354*20cf1dd8SToby Isaac   PetscValidPointer(bdsp,2);
355*20cf1dd8SToby Isaac   ierr = PetscDualSpaceLagrangeGetContinuity(sp,&continuous);CHKERRQ(ierr);
356*20cf1dd8SToby Isaac   ierr = PetscDualSpaceGetOrder(sp,&order);CHKERRQ(ierr);
357*20cf1dd8SToby Isaac   if (height == 0) {
358*20cf1dd8SToby Isaac     ierr = PetscObjectReference((PetscObject)sp);CHKERRQ(ierr);
359*20cf1dd8SToby Isaac     *bdsp = sp;
360*20cf1dd8SToby Isaac   } else if (continuous == PETSC_FALSE || !order) {
361*20cf1dd8SToby Isaac     *bdsp = NULL;
362*20cf1dd8SToby Isaac   } else {
363*20cf1dd8SToby Isaac     DM dm, K;
364*20cf1dd8SToby Isaac     PetscInt dim;
365*20cf1dd8SToby Isaac 
366*20cf1dd8SToby Isaac     ierr = PetscDualSpaceGetDM(sp,&dm);CHKERRQ(ierr);
367*20cf1dd8SToby Isaac     ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
368*20cf1dd8SToby Isaac     if (height > dim || height < 0) {SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Asked for dual space at height %d for dimension %d reference element\n",height,dim);}
369*20cf1dd8SToby Isaac     ierr = PetscDualSpaceDuplicate(sp,bdsp);CHKERRQ(ierr);
370*20cf1dd8SToby Isaac     ierr = PetscDualSpaceCreateReferenceCell(*bdsp, dim-height, lag->simplexCell, &K);CHKERRQ(ierr);
371*20cf1dd8SToby Isaac     ierr = PetscDualSpaceSetDM(*bdsp, K);CHKERRQ(ierr);
372*20cf1dd8SToby Isaac     ierr = DMDestroy(&K);CHKERRQ(ierr);
373*20cf1dd8SToby Isaac     ierr = PetscDualSpaceLagrangeGetTensor(sp,&tensor);CHKERRQ(ierr);
374*20cf1dd8SToby Isaac     ierr = PetscDualSpaceLagrangeSetTensor(*bdsp,tensor);CHKERRQ(ierr);
375*20cf1dd8SToby Isaac     ierr = PetscDualSpaceSetUp(*bdsp);CHKERRQ(ierr);
376*20cf1dd8SToby Isaac   }
377*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
378*20cf1dd8SToby Isaac }
379*20cf1dd8SToby Isaac 
380*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetUp_Lagrange(PetscDualSpace sp)
381*20cf1dd8SToby Isaac {
382*20cf1dd8SToby Isaac   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data;
383*20cf1dd8SToby Isaac   DM                  dm    = sp->dm;
384*20cf1dd8SToby Isaac   PetscInt            order = sp->order;
385*20cf1dd8SToby Isaac   PetscInt            Nc    = sp->Nc;
386*20cf1dd8SToby Isaac   PetscBool           continuous;
387*20cf1dd8SToby Isaac   PetscSection        csection;
388*20cf1dd8SToby Isaac   Vec                 coordinates;
389*20cf1dd8SToby Isaac   PetscReal          *qpoints, *qweights;
390*20cf1dd8SToby Isaac   PetscInt            depth, dim, pdimMax, pStart, pEnd, p, *pStratStart, *pStratEnd, coneSize, d, f = 0, c;
391*20cf1dd8SToby Isaac   PetscBool           simplex, tensorSpace;
392*20cf1dd8SToby Isaac   PetscErrorCode      ierr;
393*20cf1dd8SToby Isaac 
394*20cf1dd8SToby Isaac   PetscFunctionBegin;
395*20cf1dd8SToby Isaac   /* Classify element type */
396*20cf1dd8SToby Isaac   if (!order) lag->continuous = PETSC_FALSE;
397*20cf1dd8SToby Isaac   continuous = lag->continuous;
398*20cf1dd8SToby Isaac   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
399*20cf1dd8SToby Isaac   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
400*20cf1dd8SToby Isaac   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
401*20cf1dd8SToby Isaac   ierr = PetscCalloc1(dim+1, &lag->numDof);CHKERRQ(ierr);
402*20cf1dd8SToby Isaac   ierr = PetscMalloc2(depth+1,&pStratStart,depth+1,&pStratEnd);CHKERRQ(ierr);
403*20cf1dd8SToby Isaac   for (d = 0; d <= depth; ++d) {ierr = DMPlexGetDepthStratum(dm, d, &pStratStart[d], &pStratEnd[d]);CHKERRQ(ierr);}
404*20cf1dd8SToby Isaac   ierr = DMPlexGetConeSize(dm, pStratStart[depth], &coneSize);CHKERRQ(ierr);
405*20cf1dd8SToby Isaac   ierr = DMGetCoordinateSection(dm, &csection);CHKERRQ(ierr);
406*20cf1dd8SToby Isaac   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
407*20cf1dd8SToby Isaac   if (depth == 1) {
408*20cf1dd8SToby Isaac     if      (coneSize == dim+1)    simplex = PETSC_TRUE;
409*20cf1dd8SToby Isaac     else if (coneSize == 1 << dim) simplex = PETSC_FALSE;
410*20cf1dd8SToby Isaac     else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only support simplices and tensor product cells");
411*20cf1dd8SToby Isaac   } else if (depth == dim) {
412*20cf1dd8SToby Isaac     if      (coneSize == dim+1)   simplex = PETSC_TRUE;
413*20cf1dd8SToby Isaac     else if (coneSize == 2 * dim) simplex = PETSC_FALSE;
414*20cf1dd8SToby Isaac     else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only support simplices and tensor product cells");
415*20cf1dd8SToby Isaac   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only support cell-vertex meshes or interpolated meshes");
416*20cf1dd8SToby Isaac   lag->simplexCell = simplex;
417*20cf1dd8SToby Isaac   if (dim > 1 && continuous && lag->simplexCell == lag->tensorSpace) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "Mismatching simplex/tensor cells and spaces only allowed for discontinuous elements");
418*20cf1dd8SToby Isaac   tensorSpace    = lag->tensorSpace;
419*20cf1dd8SToby Isaac   lag->height    = 0;
420*20cf1dd8SToby Isaac   lag->subspaces = NULL;
421*20cf1dd8SToby Isaac   if (continuous && sp->order > 0 && dim > 0) {
422*20cf1dd8SToby Isaac     PetscInt i;
423*20cf1dd8SToby Isaac 
424*20cf1dd8SToby Isaac     lag->height = dim;
425*20cf1dd8SToby Isaac     ierr = PetscMalloc1(dim,&lag->subspaces);CHKERRQ(ierr);
426*20cf1dd8SToby Isaac     ierr = PetscDualSpaceCreateHeightSubspace_Lagrange(sp,1,&lag->subspaces[0]);CHKERRQ(ierr);
427*20cf1dd8SToby Isaac     ierr = PetscDualSpaceSetUp(lag->subspaces[0]);CHKERRQ(ierr);
428*20cf1dd8SToby Isaac     for (i = 1; i < dim; i++) {
429*20cf1dd8SToby Isaac       ierr = PetscDualSpaceGetHeightSubspace(lag->subspaces[i-1],1,&lag->subspaces[i]);CHKERRQ(ierr);
430*20cf1dd8SToby Isaac       ierr = PetscObjectReference((PetscObject)(lag->subspaces[i]));CHKERRQ(ierr);
431*20cf1dd8SToby Isaac     }
432*20cf1dd8SToby Isaac   }
433*20cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDimension_SingleCell_Lagrange(sp, sp->order, &pdimMax);CHKERRQ(ierr);
434*20cf1dd8SToby Isaac   pdimMax *= (pStratEnd[depth] - pStratStart[depth]);
435*20cf1dd8SToby Isaac   ierr = PetscMalloc1(pdimMax, &sp->functional);CHKERRQ(ierr);
436*20cf1dd8SToby Isaac   if (!dim) {
437*20cf1dd8SToby Isaac     for (c = 0; c < Nc; ++c) {
438*20cf1dd8SToby Isaac       ierr = PetscQuadratureCreate(PETSC_COMM_SELF, &sp->functional[f]);CHKERRQ(ierr);
439*20cf1dd8SToby Isaac       ierr = PetscCalloc1(Nc, &qweights);CHKERRQ(ierr);
440*20cf1dd8SToby Isaac       ierr = PetscQuadratureSetOrder(sp->functional[f], 0);CHKERRQ(ierr);
441*20cf1dd8SToby Isaac       ierr = PetscQuadratureSetData(sp->functional[f], 0, Nc, 1, NULL, qweights);CHKERRQ(ierr);
442*20cf1dd8SToby Isaac       qweights[c] = 1.0;
443*20cf1dd8SToby Isaac       ++f;
444*20cf1dd8SToby Isaac       lag->numDof[0]++;
445*20cf1dd8SToby Isaac     }
446*20cf1dd8SToby Isaac   } else {
447*20cf1dd8SToby Isaac     PetscInt     *tup;
448*20cf1dd8SToby Isaac     PetscReal    *v0, *hv0, *J, *invJ, detJ, hdetJ;
449*20cf1dd8SToby Isaac     PetscSection section;
450*20cf1dd8SToby Isaac 
451*20cf1dd8SToby Isaac     ierr = PetscSectionCreate(PETSC_COMM_SELF,&section);CHKERRQ(ierr);
452*20cf1dd8SToby Isaac     ierr = PetscSectionSetChart(section,pStart,pEnd);CHKERRQ(ierr);
453*20cf1dd8SToby Isaac     ierr = PetscCalloc5(dim+1,&tup,dim,&v0,dim,&hv0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
454*20cf1dd8SToby Isaac     for (p = pStart; p < pEnd; p++) {
455*20cf1dd8SToby Isaac       PetscInt       pointDim, d, nFunc = 0;
456*20cf1dd8SToby Isaac       PetscDualSpace hsp;
457*20cf1dd8SToby Isaac 
458*20cf1dd8SToby Isaac       ierr = DMPlexComputeCellGeometryFEM(dm, p, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
459*20cf1dd8SToby Isaac       for (d = 0; d < depth; d++) {if (p >= pStratStart[d] && p < pStratEnd[d]) break;}
460*20cf1dd8SToby Isaac       pointDim = (depth == 1 && d == 1) ? dim : d;
461*20cf1dd8SToby Isaac       hsp = ((pointDim < dim) && lag->subspaces) ? lag->subspaces[dim - pointDim - 1] : NULL;
462*20cf1dd8SToby Isaac       if (hsp) {
463*20cf1dd8SToby Isaac         PetscDualSpace_Lag *hlag = (PetscDualSpace_Lag *) hsp->data;
464*20cf1dd8SToby Isaac         DM                 hdm;
465*20cf1dd8SToby Isaac 
466*20cf1dd8SToby Isaac         ierr = PetscDualSpaceGetDM(hsp,&hdm);CHKERRQ(ierr);
467*20cf1dd8SToby Isaac         ierr = DMPlexComputeCellGeometryFEM(hdm, 0, NULL, hv0, NULL, NULL, &hdetJ);CHKERRQ(ierr);
468*20cf1dd8SToby Isaac         nFunc = lag->numDof[pointDim] = hlag->numDof[pointDim];
469*20cf1dd8SToby Isaac       }
470*20cf1dd8SToby Isaac       if (pointDim == dim) {
471*20cf1dd8SToby Isaac         /* Cells, create for self */
472*20cf1dd8SToby Isaac         PetscInt     orderEff = continuous ? (!tensorSpace ? order-1-dim : order-2) : order;
473*20cf1dd8SToby Isaac         PetscReal    denom    = continuous ? order : (!tensorSpace ? order+1+dim : order+2);
474*20cf1dd8SToby Isaac         PetscReal    numer    = (!simplex || !tensorSpace) ? 2. : (2./dim);
475*20cf1dd8SToby Isaac         PetscReal    dx = numer/denom;
476*20cf1dd8SToby Isaac         PetscInt     cdim, d, d2;
477*20cf1dd8SToby Isaac 
478*20cf1dd8SToby Isaac         if (orderEff < 0) continue;
479*20cf1dd8SToby Isaac         ierr = PetscDualSpaceGetDimension_SingleCell_Lagrange(sp, orderEff, &cdim);CHKERRQ(ierr);
480*20cf1dd8SToby Isaac         ierr = PetscMemzero(tup,(dim+1)*sizeof(PetscInt));CHKERRQ(ierr);
481*20cf1dd8SToby Isaac         if (!tensorSpace) {
482*20cf1dd8SToby Isaac           while (!tup[dim]) {
483*20cf1dd8SToby Isaac             for (c = 0; c < Nc; ++c) {
484*20cf1dd8SToby Isaac               ierr = PetscQuadratureCreate(PETSC_COMM_SELF, &sp->functional[f]);CHKERRQ(ierr);
485*20cf1dd8SToby Isaac               ierr = PetscMalloc1(dim, &qpoints);CHKERRQ(ierr);
486*20cf1dd8SToby Isaac               ierr = PetscCalloc1(Nc,  &qweights);CHKERRQ(ierr);
487*20cf1dd8SToby Isaac               ierr = PetscQuadratureSetOrder(sp->functional[f], 0);CHKERRQ(ierr);
488*20cf1dd8SToby Isaac               ierr = PetscQuadratureSetData(sp->functional[f], dim, Nc, 1, qpoints, qweights);CHKERRQ(ierr);
489*20cf1dd8SToby Isaac               for (d = 0; d < dim; ++d) {
490*20cf1dd8SToby Isaac                 qpoints[d] = v0[d];
491*20cf1dd8SToby Isaac                 for (d2 = 0; d2 < dim; ++d2) qpoints[d] += J[d*dim+d2]*((tup[d2]+1)*dx);
492*20cf1dd8SToby Isaac               }
493*20cf1dd8SToby Isaac               qweights[c] = 1.0;
494*20cf1dd8SToby Isaac               ++f;
495*20cf1dd8SToby Isaac             }
496*20cf1dd8SToby Isaac             ierr = LatticePointLexicographic_Internal(dim, orderEff, tup);CHKERRQ(ierr);
497*20cf1dd8SToby Isaac           }
498*20cf1dd8SToby Isaac         } else {
499*20cf1dd8SToby Isaac           while (!tup[dim]) {
500*20cf1dd8SToby Isaac             for (c = 0; c < Nc; ++c) {
501*20cf1dd8SToby Isaac               ierr = PetscQuadratureCreate(PETSC_COMM_SELF, &sp->functional[f]);CHKERRQ(ierr);
502*20cf1dd8SToby Isaac               ierr = PetscMalloc1(dim, &qpoints);CHKERRQ(ierr);
503*20cf1dd8SToby Isaac               ierr = PetscCalloc1(Nc,  &qweights);CHKERRQ(ierr);
504*20cf1dd8SToby Isaac               ierr = PetscQuadratureSetOrder(sp->functional[f], 0);CHKERRQ(ierr);
505*20cf1dd8SToby Isaac               ierr = PetscQuadratureSetData(sp->functional[f], dim, Nc, 1, qpoints, qweights);CHKERRQ(ierr);
506*20cf1dd8SToby Isaac               for (d = 0; d < dim; ++d) {
507*20cf1dd8SToby Isaac                 qpoints[d] = v0[d];
508*20cf1dd8SToby Isaac                 for (d2 = 0; d2 < dim; ++d2) qpoints[d] += J[d*dim+d2]*((tup[d2]+1)*dx);
509*20cf1dd8SToby Isaac               }
510*20cf1dd8SToby Isaac               qweights[c] = 1.0;
511*20cf1dd8SToby Isaac               ++f;
512*20cf1dd8SToby Isaac             }
513*20cf1dd8SToby Isaac             ierr = TensorPointLexicographic_Internal(dim, orderEff, tup);CHKERRQ(ierr);
514*20cf1dd8SToby Isaac           }
515*20cf1dd8SToby Isaac         }
516*20cf1dd8SToby Isaac         lag->numDof[dim] = cdim;
517*20cf1dd8SToby Isaac       } else { /* transform functionals from subspaces */
518*20cf1dd8SToby Isaac         PetscInt q;
519*20cf1dd8SToby Isaac 
520*20cf1dd8SToby Isaac         for (q = 0; q < nFunc; q++, f++) {
521*20cf1dd8SToby Isaac           PetscQuadrature fn;
522*20cf1dd8SToby Isaac           PetscInt        fdim, Nc, c, nPoints, i;
523*20cf1dd8SToby Isaac           const PetscReal *points;
524*20cf1dd8SToby Isaac           const PetscReal *weights;
525*20cf1dd8SToby Isaac           PetscReal       *qpoints;
526*20cf1dd8SToby Isaac           PetscReal       *qweights;
527*20cf1dd8SToby Isaac 
528*20cf1dd8SToby Isaac           ierr = PetscDualSpaceGetFunctional(hsp, q, &fn);CHKERRQ(ierr);
529*20cf1dd8SToby Isaac           ierr = PetscQuadratureGetData(fn,&fdim,&Nc,&nPoints,&points,&weights);CHKERRQ(ierr);
530*20cf1dd8SToby Isaac           if (fdim != pointDim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Expected height dual space dim %D, got %D",pointDim,fdim);
531*20cf1dd8SToby Isaac           ierr = PetscMalloc1(nPoints * dim, &qpoints);CHKERRQ(ierr);
532*20cf1dd8SToby Isaac           ierr = PetscCalloc1(nPoints * Nc,  &qweights);CHKERRQ(ierr);
533*20cf1dd8SToby Isaac           for (i = 0; i < nPoints; i++) {
534*20cf1dd8SToby Isaac             PetscInt  j, k;
535*20cf1dd8SToby Isaac             PetscReal *qp = &qpoints[i * dim];
536*20cf1dd8SToby Isaac 
537*20cf1dd8SToby Isaac             for (c = 0; c < Nc; ++c) qweights[i*Nc+c] = weights[i*Nc+c];
538*20cf1dd8SToby Isaac             for (j = 0; j < dim; ++j) qp[j] = v0[j];
539*20cf1dd8SToby Isaac             for (j = 0; j < dim; ++j) {
540*20cf1dd8SToby Isaac               for (k = 0; k < pointDim; k++) qp[j] += J[dim * j + k] * (points[pointDim * i + k] - hv0[k]);
541*20cf1dd8SToby Isaac             }
542*20cf1dd8SToby Isaac           }
543*20cf1dd8SToby Isaac           ierr = PetscQuadratureCreate(PETSC_COMM_SELF, &sp->functional[f]);CHKERRQ(ierr);
544*20cf1dd8SToby Isaac           ierr = PetscQuadratureSetOrder(sp->functional[f],0);CHKERRQ(ierr);
545*20cf1dd8SToby Isaac           ierr = PetscQuadratureSetData(sp->functional[f],dim,Nc,nPoints,qpoints,qweights);CHKERRQ(ierr);
546*20cf1dd8SToby Isaac         }
547*20cf1dd8SToby Isaac       }
548*20cf1dd8SToby Isaac       ierr = PetscSectionSetDof(section,p,lag->numDof[pointDim]);CHKERRQ(ierr);
549*20cf1dd8SToby Isaac     }
550*20cf1dd8SToby Isaac     ierr = PetscFree5(tup,v0,hv0,J,invJ);CHKERRQ(ierr);
551*20cf1dd8SToby Isaac     ierr = PetscSectionSetUp(section);CHKERRQ(ierr);
552*20cf1dd8SToby Isaac     { /* reorder to closure order */
553*20cf1dd8SToby Isaac       PetscInt *key, count;
554*20cf1dd8SToby Isaac       PetscQuadrature *reorder = NULL;
555*20cf1dd8SToby Isaac 
556*20cf1dd8SToby Isaac       ierr = PetscCalloc1(f,&key);CHKERRQ(ierr);
557*20cf1dd8SToby Isaac       ierr = PetscMalloc1(f*sp->Nc,&reorder);CHKERRQ(ierr);
558*20cf1dd8SToby Isaac 
559*20cf1dd8SToby Isaac       for (p = pStratStart[depth], count = 0; p < pStratEnd[depth]; p++) {
560*20cf1dd8SToby Isaac         PetscInt *closure = NULL, closureSize, c;
561*20cf1dd8SToby Isaac 
562*20cf1dd8SToby Isaac         ierr = DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
563*20cf1dd8SToby Isaac         for (c = 0; c < closureSize; c++) {
564*20cf1dd8SToby Isaac           PetscInt point = closure[2 * c], dof, off, i;
565*20cf1dd8SToby Isaac 
566*20cf1dd8SToby Isaac           ierr = PetscSectionGetDof(section,point,&dof);CHKERRQ(ierr);
567*20cf1dd8SToby Isaac           ierr = PetscSectionGetOffset(section,point,&off);CHKERRQ(ierr);
568*20cf1dd8SToby Isaac           for (i = 0; i < dof; i++) {
569*20cf1dd8SToby Isaac             PetscInt fi = i + off;
570*20cf1dd8SToby Isaac             if (!key[fi]) {
571*20cf1dd8SToby Isaac               key[fi] = 1;
572*20cf1dd8SToby Isaac               reorder[count++] = sp->functional[fi];
573*20cf1dd8SToby Isaac             }
574*20cf1dd8SToby Isaac           }
575*20cf1dd8SToby Isaac         }
576*20cf1dd8SToby Isaac         ierr = DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
577*20cf1dd8SToby Isaac       }
578*20cf1dd8SToby Isaac       ierr = PetscFree(sp->functional);CHKERRQ(ierr);
579*20cf1dd8SToby Isaac       sp->functional = reorder;
580*20cf1dd8SToby Isaac       ierr = PetscFree(key);CHKERRQ(ierr);
581*20cf1dd8SToby Isaac     }
582*20cf1dd8SToby Isaac     ierr = PetscSectionDestroy(&section);CHKERRQ(ierr);
583*20cf1dd8SToby Isaac   }
584*20cf1dd8SToby Isaac   if (pStratEnd[depth] == 1 && f != pdimMax) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of dual basis vectors %d not equal to dimension %d", f, pdimMax);
585*20cf1dd8SToby Isaac   ierr = PetscFree2(pStratStart, pStratEnd);CHKERRQ(ierr);
586*20cf1dd8SToby Isaac   if (f > pdimMax) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of dual basis vectors %d is greater than dimension %d", f, pdimMax);
587*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
588*20cf1dd8SToby Isaac }
589*20cf1dd8SToby Isaac 
590*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceDestroy_Lagrange(PetscDualSpace sp)
591*20cf1dd8SToby Isaac {
592*20cf1dd8SToby Isaac   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data;
593*20cf1dd8SToby Isaac   PetscInt            i;
594*20cf1dd8SToby Isaac   PetscErrorCode      ierr;
595*20cf1dd8SToby Isaac 
596*20cf1dd8SToby Isaac   PetscFunctionBegin;
597*20cf1dd8SToby Isaac   if (lag->symmetries) {
598*20cf1dd8SToby Isaac     PetscInt **selfSyms = lag->symmetries[0];
599*20cf1dd8SToby Isaac 
600*20cf1dd8SToby Isaac     if (selfSyms) {
601*20cf1dd8SToby Isaac       PetscInt i, **allocated = &selfSyms[-lag->selfSymOff];
602*20cf1dd8SToby Isaac 
603*20cf1dd8SToby Isaac       for (i = 0; i < lag->numSelfSym; i++) {
604*20cf1dd8SToby Isaac         ierr = PetscFree(allocated[i]);CHKERRQ(ierr);
605*20cf1dd8SToby Isaac       }
606*20cf1dd8SToby Isaac       ierr = PetscFree(allocated);CHKERRQ(ierr);
607*20cf1dd8SToby Isaac     }
608*20cf1dd8SToby Isaac     ierr = PetscFree(lag->symmetries);CHKERRQ(ierr);
609*20cf1dd8SToby Isaac   }
610*20cf1dd8SToby Isaac   for (i = 0; i < lag->height; i++) {
611*20cf1dd8SToby Isaac     ierr = PetscDualSpaceDestroy(&lag->subspaces[i]);CHKERRQ(ierr);
612*20cf1dd8SToby Isaac   }
613*20cf1dd8SToby Isaac   ierr = PetscFree(lag->subspaces);CHKERRQ(ierr);
614*20cf1dd8SToby Isaac   ierr = PetscFree(lag->numDof);CHKERRQ(ierr);
615*20cf1dd8SToby Isaac   ierr = PetscFree(lag);CHKERRQ(ierr);
616*20cf1dd8SToby Isaac   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetContinuity_C", NULL);CHKERRQ(ierr);
617*20cf1dd8SToby Isaac   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetContinuity_C", NULL);CHKERRQ(ierr);
618*20cf1dd8SToby Isaac   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetTensor_C", NULL);CHKERRQ(ierr);
619*20cf1dd8SToby Isaac   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetTensor_C", NULL);CHKERRQ(ierr);
620*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
621*20cf1dd8SToby Isaac }
622*20cf1dd8SToby Isaac 
623*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceDuplicate_Lagrange(PetscDualSpace sp, PetscDualSpace *spNew)
624*20cf1dd8SToby Isaac {
625*20cf1dd8SToby Isaac   PetscInt       order, Nc;
626*20cf1dd8SToby Isaac   PetscBool      cont, tensor;
627*20cf1dd8SToby Isaac   PetscErrorCode ierr;
628*20cf1dd8SToby Isaac 
629*20cf1dd8SToby Isaac   PetscFunctionBegin;
630*20cf1dd8SToby Isaac   ierr = PetscDualSpaceCreate(PetscObjectComm((PetscObject) sp), spNew);CHKERRQ(ierr);
631*20cf1dd8SToby Isaac   ierr = PetscDualSpaceSetType(*spNew, PETSCDUALSPACELAGRANGE);CHKERRQ(ierr);
632*20cf1dd8SToby Isaac   ierr = PetscDualSpaceGetOrder(sp, &order);CHKERRQ(ierr);
633*20cf1dd8SToby Isaac   ierr = PetscDualSpaceSetOrder(*spNew, order);CHKERRQ(ierr);
634*20cf1dd8SToby Isaac   ierr = PetscDualSpaceGetNumComponents(sp, &Nc);CHKERRQ(ierr);
635*20cf1dd8SToby Isaac   ierr = PetscDualSpaceSetNumComponents(*spNew, Nc);CHKERRQ(ierr);
636*20cf1dd8SToby Isaac   ierr = PetscDualSpaceLagrangeGetContinuity(sp, &cont);CHKERRQ(ierr);
637*20cf1dd8SToby Isaac   ierr = PetscDualSpaceLagrangeSetContinuity(*spNew, cont);CHKERRQ(ierr);
638*20cf1dd8SToby Isaac   ierr = PetscDualSpaceLagrangeGetTensor(sp, &tensor);CHKERRQ(ierr);
639*20cf1dd8SToby Isaac   ierr = PetscDualSpaceLagrangeSetTensor(*spNew, tensor);CHKERRQ(ierr);
640*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
641*20cf1dd8SToby Isaac }
642*20cf1dd8SToby Isaac 
643*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceSetFromOptions_Lagrange(PetscOptionItems *PetscOptionsObject,PetscDualSpace sp)
644*20cf1dd8SToby Isaac {
645*20cf1dd8SToby Isaac   PetscBool      continuous, tensor, flg;
646*20cf1dd8SToby Isaac   PetscErrorCode ierr;
647*20cf1dd8SToby Isaac 
648*20cf1dd8SToby Isaac   PetscFunctionBegin;
649*20cf1dd8SToby Isaac   ierr = PetscDualSpaceLagrangeGetContinuity(sp, &continuous);CHKERRQ(ierr);
650*20cf1dd8SToby Isaac   ierr = PetscDualSpaceLagrangeGetTensor(sp, &tensor);CHKERRQ(ierr);
651*20cf1dd8SToby Isaac   ierr = PetscOptionsHead(PetscOptionsObject,"PetscDualSpace Lagrange Options");CHKERRQ(ierr);
652*20cf1dd8SToby Isaac   ierr = PetscOptionsBool("-petscdualspace_lagrange_continuity", "Flag for continuous element", "PetscDualSpaceLagrangeSetContinuity", continuous, &continuous, &flg);CHKERRQ(ierr);
653*20cf1dd8SToby Isaac   if (flg) {ierr = PetscDualSpaceLagrangeSetContinuity(sp, continuous);CHKERRQ(ierr);}
654*20cf1dd8SToby Isaac   ierr = PetscOptionsBool("-petscdualspace_lagrange_tensor", "Flag for tensor dual space", "PetscDualSpaceLagrangeSetContinuity", tensor, &tensor, &flg);CHKERRQ(ierr);
655*20cf1dd8SToby Isaac   if (flg) {ierr = PetscDualSpaceLagrangeSetTensor(sp, tensor);CHKERRQ(ierr);}
656*20cf1dd8SToby Isaac   ierr = PetscOptionsTail();CHKERRQ(ierr);
657*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
658*20cf1dd8SToby Isaac }
659*20cf1dd8SToby Isaac 
660*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetDimension_Lagrange(PetscDualSpace sp, PetscInt *dim)
661*20cf1dd8SToby Isaac {
662*20cf1dd8SToby Isaac   DM              K;
663*20cf1dd8SToby Isaac   const PetscInt *numDof;
664*20cf1dd8SToby Isaac   PetscInt        spatialDim, Nc, size = 0, d;
665*20cf1dd8SToby Isaac   PetscErrorCode  ierr;
666*20cf1dd8SToby Isaac 
667*20cf1dd8SToby Isaac   PetscFunctionBegin;
668*20cf1dd8SToby Isaac   ierr = PetscDualSpaceGetDM(sp, &K);CHKERRQ(ierr);
669*20cf1dd8SToby Isaac   ierr = PetscDualSpaceGetNumDof(sp, &numDof);CHKERRQ(ierr);
670*20cf1dd8SToby Isaac   ierr = DMGetDimension(K, &spatialDim);CHKERRQ(ierr);
671*20cf1dd8SToby Isaac   ierr = DMPlexGetHeightStratum(K, 0, NULL, &Nc);CHKERRQ(ierr);
672*20cf1dd8SToby Isaac   if (Nc == 1) {ierr = PetscDualSpaceGetDimension_SingleCell_Lagrange(sp, sp->order, dim);CHKERRQ(ierr); PetscFunctionReturn(0);}
673*20cf1dd8SToby Isaac   for (d = 0; d <= spatialDim; ++d) {
674*20cf1dd8SToby Isaac     PetscInt pStart, pEnd;
675*20cf1dd8SToby Isaac 
676*20cf1dd8SToby Isaac     ierr = DMPlexGetDepthStratum(K, d, &pStart, &pEnd);CHKERRQ(ierr);
677*20cf1dd8SToby Isaac     size += (pEnd-pStart)*numDof[d];
678*20cf1dd8SToby Isaac   }
679*20cf1dd8SToby Isaac   *dim = size;
680*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
681*20cf1dd8SToby Isaac }
682*20cf1dd8SToby Isaac 
683*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetNumDof_Lagrange(PetscDualSpace sp, const PetscInt **numDof)
684*20cf1dd8SToby Isaac {
685*20cf1dd8SToby Isaac   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data;
686*20cf1dd8SToby Isaac 
687*20cf1dd8SToby Isaac   PetscFunctionBegin;
688*20cf1dd8SToby Isaac   *numDof = lag->numDof;
689*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
690*20cf1dd8SToby Isaac }
691*20cf1dd8SToby Isaac 
692*20cf1dd8SToby Isaac static PetscErrorCode PetscDualSpaceLagrangeGetContinuity_Lagrange(PetscDualSpace sp, PetscBool *continuous)
693*20cf1dd8SToby Isaac {
694*20cf1dd8SToby Isaac   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data;
695*20cf1dd8SToby Isaac 
696*20cf1dd8SToby Isaac   PetscFunctionBegin;
697*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
698*20cf1dd8SToby Isaac   PetscValidPointer(continuous, 2);
699*20cf1dd8SToby Isaac   *continuous = lag->continuous;
700*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
701*20cf1dd8SToby Isaac }
702*20cf1dd8SToby Isaac 
703*20cf1dd8SToby Isaac static PetscErrorCode PetscDualSpaceLagrangeSetContinuity_Lagrange(PetscDualSpace sp, PetscBool continuous)
704*20cf1dd8SToby Isaac {
705*20cf1dd8SToby Isaac   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data;
706*20cf1dd8SToby Isaac 
707*20cf1dd8SToby Isaac   PetscFunctionBegin;
708*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
709*20cf1dd8SToby Isaac   lag->continuous = continuous;
710*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
711*20cf1dd8SToby Isaac }
712*20cf1dd8SToby Isaac 
713*20cf1dd8SToby Isaac /*@
714*20cf1dd8SToby Isaac   PetscDualSpaceLagrangeGetContinuity - Retrieves the flag for element continuity
715*20cf1dd8SToby Isaac 
716*20cf1dd8SToby Isaac   Not Collective
717*20cf1dd8SToby Isaac 
718*20cf1dd8SToby Isaac   Input Parameter:
719*20cf1dd8SToby Isaac . sp         - the PetscDualSpace
720*20cf1dd8SToby Isaac 
721*20cf1dd8SToby Isaac   Output Parameter:
722*20cf1dd8SToby Isaac . continuous - flag for element continuity
723*20cf1dd8SToby Isaac 
724*20cf1dd8SToby Isaac   Level: intermediate
725*20cf1dd8SToby Isaac 
726*20cf1dd8SToby Isaac .keywords: PetscDualSpace, Lagrange, continuous, discontinuous
727*20cf1dd8SToby Isaac .seealso: PetscDualSpaceLagrangeSetContinuity()
728*20cf1dd8SToby Isaac @*/
729*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceLagrangeGetContinuity(PetscDualSpace sp, PetscBool *continuous)
730*20cf1dd8SToby Isaac {
731*20cf1dd8SToby Isaac   PetscErrorCode ierr;
732*20cf1dd8SToby Isaac 
733*20cf1dd8SToby Isaac   PetscFunctionBegin;
734*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
735*20cf1dd8SToby Isaac   PetscValidPointer(continuous, 2);
736*20cf1dd8SToby Isaac   ierr = PetscTryMethod(sp, "PetscDualSpaceLagrangeGetContinuity_C", (PetscDualSpace,PetscBool*),(sp,continuous));CHKERRQ(ierr);
737*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
738*20cf1dd8SToby Isaac }
739*20cf1dd8SToby Isaac 
740*20cf1dd8SToby Isaac /*@
741*20cf1dd8SToby Isaac   PetscDualSpaceLagrangeSetContinuity - Indicate whether the element is continuous
742*20cf1dd8SToby Isaac 
743*20cf1dd8SToby Isaac   Logically Collective on PetscDualSpace
744*20cf1dd8SToby Isaac 
745*20cf1dd8SToby Isaac   Input Parameters:
746*20cf1dd8SToby Isaac + sp         - the PetscDualSpace
747*20cf1dd8SToby Isaac - continuous - flag for element continuity
748*20cf1dd8SToby Isaac 
749*20cf1dd8SToby Isaac   Options Database:
750*20cf1dd8SToby Isaac . -petscdualspace_lagrange_continuity <bool>
751*20cf1dd8SToby Isaac 
752*20cf1dd8SToby Isaac   Level: intermediate
753*20cf1dd8SToby Isaac 
754*20cf1dd8SToby Isaac .keywords: PetscDualSpace, Lagrange, continuous, discontinuous
755*20cf1dd8SToby Isaac .seealso: PetscDualSpaceLagrangeGetContinuity()
756*20cf1dd8SToby Isaac @*/
757*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceLagrangeSetContinuity(PetscDualSpace sp, PetscBool continuous)
758*20cf1dd8SToby Isaac {
759*20cf1dd8SToby Isaac   PetscErrorCode ierr;
760*20cf1dd8SToby Isaac 
761*20cf1dd8SToby Isaac   PetscFunctionBegin;
762*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
763*20cf1dd8SToby Isaac   PetscValidLogicalCollectiveBool(sp, continuous, 2);
764*20cf1dd8SToby Isaac   ierr = PetscTryMethod(sp, "PetscDualSpaceLagrangeSetContinuity_C", (PetscDualSpace,PetscBool),(sp,continuous));CHKERRQ(ierr);
765*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
766*20cf1dd8SToby Isaac }
767*20cf1dd8SToby Isaac 
768*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceGetHeightSubspace_Lagrange(PetscDualSpace sp, PetscInt height, PetscDualSpace *bdsp)
769*20cf1dd8SToby Isaac {
770*20cf1dd8SToby Isaac   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data;
771*20cf1dd8SToby Isaac   PetscErrorCode     ierr;
772*20cf1dd8SToby Isaac 
773*20cf1dd8SToby Isaac   PetscFunctionBegin;
774*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
775*20cf1dd8SToby Isaac   PetscValidPointer(bdsp,2);
776*20cf1dd8SToby Isaac   if (height == 0) {
777*20cf1dd8SToby Isaac     *bdsp = sp;
778*20cf1dd8SToby Isaac   }
779*20cf1dd8SToby Isaac   else {
780*20cf1dd8SToby Isaac     DM dm;
781*20cf1dd8SToby Isaac     PetscInt dim;
782*20cf1dd8SToby Isaac 
783*20cf1dd8SToby Isaac     ierr = PetscDualSpaceGetDM(sp,&dm);CHKERRQ(ierr);
784*20cf1dd8SToby Isaac     ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
785*20cf1dd8SToby Isaac     if (height > dim || height < 0) {SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Asked for dual space at height %d for dimension %d reference element\n",height,dim);}
786*20cf1dd8SToby Isaac     if (height <= lag->height) {
787*20cf1dd8SToby Isaac       *bdsp = lag->subspaces[height-1];
788*20cf1dd8SToby Isaac     }
789*20cf1dd8SToby Isaac     else {
790*20cf1dd8SToby Isaac       *bdsp = NULL;
791*20cf1dd8SToby Isaac     }
792*20cf1dd8SToby Isaac   }
793*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
794*20cf1dd8SToby Isaac }
795*20cf1dd8SToby Isaac 
796*20cf1dd8SToby Isaac PetscErrorCode PetscDualSpaceInitialize_Lagrange(PetscDualSpace sp)
797*20cf1dd8SToby Isaac {
798*20cf1dd8SToby Isaac   PetscFunctionBegin;
799*20cf1dd8SToby Isaac   sp->ops->setfromoptions    = PetscDualSpaceSetFromOptions_Lagrange;
800*20cf1dd8SToby Isaac   sp->ops->setup             = PetscDualSpaceSetUp_Lagrange;
801*20cf1dd8SToby Isaac   sp->ops->view              = PetscDualSpaceView_Lagrange;
802*20cf1dd8SToby Isaac   sp->ops->destroy           = PetscDualSpaceDestroy_Lagrange;
803*20cf1dd8SToby Isaac   sp->ops->duplicate         = PetscDualSpaceDuplicate_Lagrange;
804*20cf1dd8SToby Isaac   sp->ops->getdimension      = PetscDualSpaceGetDimension_Lagrange;
805*20cf1dd8SToby Isaac   sp->ops->getnumdof         = PetscDualSpaceGetNumDof_Lagrange;
806*20cf1dd8SToby Isaac   sp->ops->getheightsubspace = PetscDualSpaceGetHeightSubspace_Lagrange;
807*20cf1dd8SToby Isaac   sp->ops->getsymmetries     = PetscDualSpaceGetSymmetries_Lagrange;
808*20cf1dd8SToby Isaac   sp->ops->apply             = PetscDualSpaceApplyDefault;
809*20cf1dd8SToby Isaac   sp->ops->applyall          = PetscDualSpaceApplyAllDefault;
810*20cf1dd8SToby Isaac   sp->ops->createallpoints   = PetscDualSpaceCreateAllPointsDefault;
811*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
812*20cf1dd8SToby Isaac }
813*20cf1dd8SToby Isaac 
814*20cf1dd8SToby Isaac /*MC
815*20cf1dd8SToby Isaac   PETSCDUALSPACELAGRANGE = "lagrange" - A PetscDualSpace object that encapsulates a dual space of pointwise evaluation functionals
816*20cf1dd8SToby Isaac 
817*20cf1dd8SToby Isaac   Level: intermediate
818*20cf1dd8SToby Isaac 
819*20cf1dd8SToby Isaac .seealso: PetscDualSpaceType, PetscDualSpaceCreate(), PetscDualSpaceSetType()
820*20cf1dd8SToby Isaac M*/
821*20cf1dd8SToby Isaac 
822*20cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceCreate_Lagrange(PetscDualSpace sp)
823*20cf1dd8SToby Isaac {
824*20cf1dd8SToby Isaac   PetscDualSpace_Lag *lag;
825*20cf1dd8SToby Isaac   PetscErrorCode      ierr;
826*20cf1dd8SToby Isaac 
827*20cf1dd8SToby Isaac   PetscFunctionBegin;
828*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
829*20cf1dd8SToby Isaac   ierr     = PetscNewLog(sp,&lag);CHKERRQ(ierr);
830*20cf1dd8SToby Isaac   sp->data = lag;
831*20cf1dd8SToby Isaac 
832*20cf1dd8SToby Isaac   lag->numDof      = NULL;
833*20cf1dd8SToby Isaac   lag->simplexCell = PETSC_TRUE;
834*20cf1dd8SToby Isaac   lag->tensorSpace = PETSC_FALSE;
835*20cf1dd8SToby Isaac   lag->continuous  = PETSC_TRUE;
836*20cf1dd8SToby Isaac 
837*20cf1dd8SToby Isaac   ierr = PetscDualSpaceInitialize_Lagrange(sp);CHKERRQ(ierr);
838*20cf1dd8SToby Isaac   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetContinuity_C", PetscDualSpaceLagrangeGetContinuity_Lagrange);CHKERRQ(ierr);
839*20cf1dd8SToby Isaac   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetContinuity_C", PetscDualSpaceLagrangeSetContinuity_Lagrange);CHKERRQ(ierr);
840*20cf1dd8SToby Isaac   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetTensor_C", PetscDualSpaceLagrangeGetTensor_Lagrange);CHKERRQ(ierr);
841*20cf1dd8SToby Isaac   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetTensor_C", PetscDualSpaceLagrangeSetTensor_Lagrange);CHKERRQ(ierr);
842*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
843*20cf1dd8SToby Isaac }
844*20cf1dd8SToby Isaac 
845