xref: /petsc/src/dm/dt/interface/dt.c (revision eb8bc6c8b0a47cf79b83038d572b2ecbe73f8b61)
1 /* Discretization tools */
2 
3 #include <petscconf.h>
4 #if defined(PETSC_HAVE_MATHIMF_H)
5 #include <mathimf.h>           /* this needs to be included before math.h */
6 #endif
7 
8 #include <petscdt.h>            /*I "petscdt.h" I*/ /*I "petscfe.h" I*/
9 #include <petscblaslapack.h>
10 #include <petsc-private/petscimpl.h>
11 #include <petscviewer.h>
12 #include <petscdmplex.h>
13 #include <petscdmshell.h>
14 
15 #undef __FUNCT__
16 #define __FUNCT__ "PetscDTLegendreEval"
17 /*@
18    PetscDTLegendreEval - evaluate Legendre polynomial at points
19 
20    Not Collective
21 
22    Input Arguments:
23 +  npoints - number of spatial points to evaluate at
24 .  points - array of locations to evaluate at
25 .  ndegree - number of basis degrees to evaluate
26 -  degrees - sorted array of degrees to evaluate
27 
28    Output Arguments:
29 +  B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL)
30 .  D - row-oriented derivative evaluation matrix (or NULL)
31 -  D2 - row-oriented second derivative evaluation matrix (or NULL)
32 
33    Level: intermediate
34 
35 .seealso: PetscDTGaussQuadrature()
36 @*/
37 PetscErrorCode PetscDTLegendreEval(PetscInt npoints,const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2)
38 {
39   PetscInt i,maxdegree;
40 
41   PetscFunctionBegin;
42   if (!npoints || !ndegree) PetscFunctionReturn(0);
43   maxdegree = degrees[ndegree-1];
44   for (i=0; i<npoints; i++) {
45     PetscReal pm1,pm2,pd1,pd2,pdd1,pdd2,x;
46     PetscInt  j,k;
47     x    = points[i];
48     pm2  = 0;
49     pm1  = 1;
50     pd2  = 0;
51     pd1  = 0;
52     pdd2 = 0;
53     pdd1 = 0;
54     k    = 0;
55     if (degrees[k] == 0) {
56       if (B) B[i*ndegree+k] = pm1;
57       if (D) D[i*ndegree+k] = pd1;
58       if (D2) D2[i*ndegree+k] = pdd1;
59       k++;
60     }
61     for (j=1; j<=maxdegree; j++,k++) {
62       PetscReal p,d,dd;
63       p    = ((2*j-1)*x*pm1 - (j-1)*pm2)/j;
64       d    = pd2 + (2*j-1)*pm1;
65       dd   = pdd2 + (2*j-1)*pd1;
66       pm2  = pm1;
67       pm1  = p;
68       pd2  = pd1;
69       pd1  = d;
70       pdd2 = pdd1;
71       pdd1 = dd;
72       if (degrees[k] == j) {
73         if (B) B[i*ndegree+k] = p;
74         if (D) D[i*ndegree+k] = d;
75         if (D2) D2[i*ndegree+k] = dd;
76       }
77     }
78   }
79   PetscFunctionReturn(0);
80 }
81 
82 #undef __FUNCT__
83 #define __FUNCT__ "PetscDTGaussQuadrature"
84 /*@
85    PetscDTGaussQuadrature - create Gauss quadrature
86 
87    Not Collective
88 
89    Input Arguments:
90 +  npoints - number of points
91 .  a - left end of interval (often-1)
92 -  b - right end of interval (often +1)
93 
94    Output Arguments:
95 +  x - quadrature points
96 -  w - quadrature weights
97 
98    Level: intermediate
99 
100    References:
101    Golub and Welsch, Calculation of Quadrature Rules, Math. Comp. 23(106), 221--230, 1969.
102 
103 .seealso: PetscDTLegendreEval()
104 @*/
105 PetscErrorCode PetscDTGaussQuadrature(PetscInt npoints,PetscReal a,PetscReal b,PetscReal *x,PetscReal *w)
106 {
107   PetscErrorCode ierr;
108   PetscInt       i;
109   PetscReal      *work;
110   PetscScalar    *Z;
111   PetscBLASInt   N,LDZ,info;
112 
113   PetscFunctionBegin;
114   /* Set up the Golub-Welsch system */
115   for (i=0; i<npoints; i++) {
116     x[i] = 0;                   /* diagonal is 0 */
117     if (i) w[i-1] = 0.5 / PetscSqrtReal(1 - 1./PetscSqr(2*i));
118   }
119   ierr = PetscRealView(npoints-1,w,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
120   ierr = PetscMalloc2(npoints*npoints,PetscScalar,&Z,PetscMax(1,2*npoints-2),PetscReal,&work);CHKERRQ(ierr);
121   ierr = PetscBLASIntCast(npoints,&N);CHKERRQ(ierr);
122   LDZ  = N;
123   ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
124   PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("I",&N,x,w,Z,&LDZ,work,&info));
125   ierr = PetscFPTrapPop();CHKERRQ(ierr);
126   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEQR error");
127 
128   for (i=0; i<(npoints+1)/2; i++) {
129     PetscReal y = 0.5 * (-x[i] + x[npoints-i-1]); /* enforces symmetry */
130     x[i]           = (a+b)/2 - y*(b-a)/2;
131     x[npoints-i-1] = (a+b)/2 + y*(b-a)/2;
132 
133     w[i] = w[npoints-1-i] = (b-a)*PetscSqr(0.5*PetscAbsScalar(Z[i*npoints] + Z[(npoints-i-1)*npoints]));
134   }
135   ierr = PetscFree2(Z,work);CHKERRQ(ierr);
136   PetscFunctionReturn(0);
137 }
138 
139 #undef __FUNCT__
140 #define __FUNCT__ "PetscDTFactorial_Internal"
141 /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x.
142    Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */
143 PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorial_Internal(PetscInt n, PetscReal *factorial)
144 {
145   PetscReal f = 1.0;
146   PetscInt  i;
147 
148   PetscFunctionBegin;
149   for (i = 1; i < n+1; ++i) f *= i;
150   *factorial = f;
151   PetscFunctionReturn(0);
152 }
153 
154 #undef __FUNCT__
155 #define __FUNCT__ "PetscDTComputeJacobi"
156 /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x.
157    Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */
158 PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobi(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P)
159 {
160   PetscReal apb, pn1, pn2;
161   PetscInt  k;
162 
163   PetscFunctionBegin;
164   if (!n) {*P = 1.0; PetscFunctionReturn(0);}
165   if (n == 1) {*P = 0.5 * (a - b + (a + b + 2.0) * x); PetscFunctionReturn(0);}
166   apb = a + b;
167   pn2 = 1.0;
168   pn1 = 0.5 * (a - b + (apb + 2.0) * x);
169   *P  = 0.0;
170   for (k = 2; k < n+1; ++k) {
171     PetscReal a1 = 2.0 * k * (k + apb) * (2.0*k + apb - 2.0);
172     PetscReal a2 = (2.0 * k + apb - 1.0) * (a*a - b*b);
173     PetscReal a3 = (2.0 * k + apb - 2.0) * (2.0 * k + apb - 1.0) * (2.0 * k + apb);
174     PetscReal a4 = 2.0 * (k + a - 1.0) * (k + b - 1.0) * (2.0 * k + apb);
175 
176     a2  = a2 / a1;
177     a3  = a3 / a1;
178     a4  = a4 / a1;
179     *P  = (a2 + a3 * x) * pn1 - a4 * pn2;
180     pn2 = pn1;
181     pn1 = *P;
182   }
183   PetscFunctionReturn(0);
184 }
185 
186 #undef __FUNCT__
187 #define __FUNCT__ "PetscDTComputeJacobiDerivative"
188 /* Evaluates the first derivative of P_{n}^{a,b} at a point x. */
189 PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobiDerivative(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P)
190 {
191   PetscReal      nP;
192   PetscErrorCode ierr;
193 
194   PetscFunctionBegin;
195   if (!n) {*P = 0.0; PetscFunctionReturn(0);}
196   ierr = PetscDTComputeJacobi(a+1, b+1, n-1, x, &nP);CHKERRQ(ierr);
197   *P   = 0.5 * (a + b + n + 1) * nP;
198   PetscFunctionReturn(0);
199 }
200 
201 #undef __FUNCT__
202 #define __FUNCT__ "PetscDTMapSquareToTriangle_Internal"
203 /* Maps from [-1,1]^2 to the (-1,1) reference triangle */
204 PETSC_STATIC_INLINE PetscErrorCode PetscDTMapSquareToTriangle_Internal(PetscReal x, PetscReal y, PetscReal *xi, PetscReal *eta)
205 {
206   PetscFunctionBegin;
207   *xi  = 0.5 * (1.0 + x) * (1.0 - y) - 1.0;
208   *eta = y;
209   PetscFunctionReturn(0);
210 }
211 
212 #undef __FUNCT__
213 #define __FUNCT__ "PetscDTMapCubeToTetrahedron_Internal"
214 /* Maps from [-1,1]^2 to the (-1,1) reference triangle */
215 PETSC_STATIC_INLINE PetscErrorCode PetscDTMapCubeToTetrahedron_Internal(PetscReal x, PetscReal y, PetscReal z, PetscReal *xi, PetscReal *eta, PetscReal *zeta)
216 {
217   PetscFunctionBegin;
218   *xi   = 0.25 * (1.0 + x) * (1.0 - y) * (1.0 - z) - 1.0;
219   *eta  = 0.5  * (1.0 + y) * (1.0 - z) - 1.0;
220   *zeta = z;
221   PetscFunctionReturn(0);
222 }
223 
224 #undef __FUNCT__
225 #define __FUNCT__ "PetscDTGaussJacobiQuadrature1D_Internal"
226 static PetscErrorCode PetscDTGaussJacobiQuadrature1D_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal *x, PetscReal *w)
227 {
228   PetscInt       maxIter = 100;
229   PetscReal      eps     = 1.0e-8;
230   PetscReal      a1, a2, a3, a4, a5, a6;
231   PetscInt       k;
232   PetscErrorCode ierr;
233 
234   PetscFunctionBegin;
235 
236   a1      = pow(2, a+b+1);
237 #if defined(PETSC_HAVE_TGAMMA)
238   a2      = tgamma(a + npoints + 1);
239   a3      = tgamma(b + npoints + 1);
240   a4      = tgamma(a + b + npoints + 1);
241 #else
242   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"tgamma() - math routine is unavailable.");
243 #endif
244 
245   ierr = PetscDTFactorial_Internal(npoints, &a5);CHKERRQ(ierr);
246   a6   = a1 * a2 * a3 / a4 / a5;
247   /* Computes the m roots of P_{m}^{a,b} on [-1,1] by Newton's method with Chebyshev points as initial guesses.
248    Algorithm implemented from the pseudocode given by Karniadakis and Sherwin and Python in FIAT */
249   for (k = 0; k < npoints; ++k) {
250     PetscReal r = -cos((2.0*k + 1.0) * PETSC_PI / (2.0 * npoints)), dP;
251     PetscInt  j;
252 
253     if (k > 0) r = 0.5 * (r + x[k-1]);
254     for (j = 0; j < maxIter; ++j) {
255       PetscReal s = 0.0, delta, f, fp;
256       PetscInt  i;
257 
258       for (i = 0; i < k; ++i) s = s + 1.0 / (r - x[i]);
259       ierr = PetscDTComputeJacobi(a, b, npoints, r, &f);CHKERRQ(ierr);
260       ierr = PetscDTComputeJacobiDerivative(a, b, npoints, r, &fp);CHKERRQ(ierr);
261       delta = f / (fp - f * s);
262       r     = r - delta;
263       if (fabs(delta) < eps) break;
264     }
265     x[k] = r;
266     ierr = PetscDTComputeJacobiDerivative(a, b, npoints, x[k], &dP);CHKERRQ(ierr);
267     w[k] = a6 / (1.0 - PetscSqr(x[k])) / PetscSqr(dP);
268   }
269   PetscFunctionReturn(0);
270 }
271 
272 #undef __FUNCT__
273 #define __FUNCT__ "PetscDTGaussJacobiQuadrature"
274 /*@C
275   PetscDTGaussJacobiQuadrature - create Gauss-Jacobi quadrature for a simplex
276 
277   Not Collective
278 
279   Input Arguments:
280 + dim - The simplex dimension
281 . npoints - number of points
282 . a - left end of interval (often-1)
283 - b - right end of interval (often +1)
284 
285   Output Arguments:
286 + points - quadrature points
287 - weights - quadrature weights
288 
289   Level: intermediate
290 
291   References:
292   Karniadakis and Sherwin.
293   FIAT
294 
295 .seealso: PetscDTGaussQuadrature()
296 @*/
297 PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt dim, PetscInt npoints, PetscReal a, PetscReal b, PetscReal *points[], PetscReal *weights[])
298 {
299   PetscReal     *px, *wx, *py, *wy, *pz, *wz, *x, *w;
300   PetscInt       i, j, k;
301   PetscErrorCode ierr;
302 
303   PetscFunctionBegin;
304   if ((a != -1.0) || (b != 1.0)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must use default internal right now");
305   switch (dim) {
306   case 1:
307     ierr = PetscMalloc(npoints * sizeof(PetscReal), &x);CHKERRQ(ierr);
308     ierr = PetscMalloc(npoints * sizeof(PetscReal), &w);CHKERRQ(ierr);
309     ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, x, w);CHKERRQ(ierr);
310     break;
311   case 2:
312     ierr = PetscMalloc(npoints*npoints*2 * sizeof(PetscReal), &x);CHKERRQ(ierr);
313     ierr = PetscMalloc(npoints*npoints   * sizeof(PetscReal), &w);CHKERRQ(ierr);
314     ierr = PetscMalloc4(npoints,PetscReal,&px,npoints,PetscReal,&wx,npoints,PetscReal,&py,npoints,PetscReal,&wy);CHKERRQ(ierr);
315     ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, px, wx);CHKERRQ(ierr);
316     ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 1.0, 0.0, py, wy);CHKERRQ(ierr);
317     for (i = 0; i < npoints; ++i) {
318       for (j = 0; j < npoints; ++j) {
319         ierr = PetscDTMapSquareToTriangle_Internal(px[i], py[j], &x[(i*npoints+j)*2+0], &x[(i*npoints+j)*2+1]);CHKERRQ(ierr);
320         w[i*npoints+j] = 0.5 * wx[i] * wy[j];
321       }
322     }
323     ierr = PetscFree4(px,wx,py,wy);CHKERRQ(ierr);
324     break;
325   case 3:
326     ierr = PetscMalloc(npoints*npoints*3 * sizeof(PetscReal), &x);CHKERRQ(ierr);
327     ierr = PetscMalloc(npoints*npoints   * sizeof(PetscReal), &w);CHKERRQ(ierr);
328     ierr = PetscMalloc6(npoints,PetscReal,&px,npoints,PetscReal,&wx,npoints,PetscReal,&py,npoints,PetscReal,&wy,npoints,PetscReal,&pz,npoints,PetscReal,&wz);CHKERRQ(ierr);
329     ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, px, wx);CHKERRQ(ierr);
330     ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 1.0, 0.0, py, wy);CHKERRQ(ierr);
331     ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 2.0, 0.0, pz, wz);CHKERRQ(ierr);
332     for (i = 0; i < npoints; ++i) {
333       for (j = 0; j < npoints; ++j) {
334         for (k = 0; k < npoints; ++k) {
335           ierr = PetscDTMapCubeToTetrahedron_Internal(px[i], py[j], pz[k], &x[((i*npoints+j)*npoints+k)*3+0], &x[((i*npoints+j)*npoints+k)*3+1], &x[((i*npoints+j)*npoints+k)*3+2]);CHKERRQ(ierr);
336           w[(i*npoints+j)*npoints+k] = 0.125 * wx[i] * wy[j] * wz[k];
337         }
338       }
339     }
340     ierr = PetscFree6(px,wx,py,wy,pz,wz);CHKERRQ(ierr);
341     break;
342   default:
343     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim);
344   }
345   if (points)  *points  = x;
346   if (weights) *weights = w;
347   PetscFunctionReturn(0);
348 }
349 
350 #undef __FUNCT__
351 #define __FUNCT__ "PetscDTPseudoInverseQR"
352 /* Overwrites A. Can only handle full-rank problems with m>=n
353  * A in column-major format
354  * Ainv in row-major format
355  * tau has length m
356  * worksize must be >= max(1,n)
357  */
358 static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m,PetscInt mstride,PetscInt n,PetscReal *A_in,PetscReal *Ainv_out,PetscScalar *tau,PetscInt worksize,PetscScalar *work)
359 {
360   PetscErrorCode ierr;
361   PetscBLASInt M,N,K,lda,ldb,ldwork,info;
362   PetscScalar *A,*Ainv,*R,*Q,Alpha;
363 
364   PetscFunctionBegin;
365 #if defined(PETSC_USE_COMPLEX)
366   {
367     PetscInt i,j;
368     ierr = PetscMalloc2(m*n,PetscScalar,&A,m*n,PetscScalar,&Ainv);CHKERRQ(ierr);
369     for (j=0; j<n; j++) {
370       for (i=0; i<m; i++) A[i+m*j] = A_in[i+mstride*j];
371     }
372     mstride = m;
373   }
374 #else
375   A = A_in;
376   Ainv = Ainv_out;
377 #endif
378 
379   ierr = PetscBLASIntCast(m,&M);CHKERRQ(ierr);
380   ierr = PetscBLASIntCast(n,&N);CHKERRQ(ierr);
381   ierr = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr);
382   ierr = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr);
383   ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
384   LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info);
385   ierr = PetscFPTrapPop();CHKERRQ(ierr);
386   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error");
387   R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */
388 
389   /* Extract an explicit representation of Q */
390   Q = Ainv;
391   ierr = PetscMemcpy(Q,A,mstride*n*sizeof(PetscScalar));CHKERRQ(ierr);
392   K = N;                        /* full rank */
393   LAPACKungqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info);
394   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error");
395 
396   /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */
397   Alpha = 1.0;
398   ldb = lda;
399   BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb);
400   /* Ainv is Q, overwritten with inverse */
401 
402 #if defined(PETSC_USE_COMPLEX)
403   {
404     PetscInt i;
405     for (i=0; i<m*n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]);
406     ierr = PetscFree2(A,Ainv);CHKERRQ(ierr);
407   }
408 #endif
409   PetscFunctionReturn(0);
410 }
411 
412 #undef __FUNCT__
413 #define __FUNCT__ "PetscDTLegendreIntegrate"
414 /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */
415 static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval,const PetscReal *x,PetscInt ndegree,const PetscInt *degrees,PetscBool Transpose,PetscReal *B)
416 {
417   PetscErrorCode ierr;
418   PetscReal *Bv;
419   PetscInt i,j;
420 
421   PetscFunctionBegin;
422   ierr = PetscMalloc((ninterval+1)*ndegree*sizeof(PetscReal),&Bv);CHKERRQ(ierr);
423   /* Point evaluation of L_p on all the source vertices */
424   ierr = PetscDTLegendreEval(ninterval+1,x,ndegree,degrees,Bv,NULL,NULL);CHKERRQ(ierr);
425   /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */
426   for (i=0; i<ninterval; i++) {
427     for (j=0; j<ndegree; j++) {
428       if (Transpose) B[i+ninterval*j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j];
429       else           B[i*ndegree+j]   = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j];
430     }
431   }
432   ierr = PetscFree(Bv);CHKERRQ(ierr);
433   PetscFunctionReturn(0);
434 }
435 
436 #undef __FUNCT__
437 #define __FUNCT__ "PetscDTReconstructPoly"
438 /*@
439    PetscDTReconstructPoly - create matrix representing polynomial reconstruction using cell intervals and evaluation at target intervals
440 
441    Not Collective
442 
443    Input Arguments:
444 +  degree - degree of reconstruction polynomial
445 .  nsource - number of source intervals
446 .  sourcex - sorted coordinates of source cell boundaries (length nsource+1)
447 .  ntarget - number of target intervals
448 -  targetx - sorted coordinates of target cell boundaries (length ntarget+1)
449 
450    Output Arguments:
451 .  R - reconstruction matrix, utarget = sum_s R[t*nsource+s] * usource[s]
452 
453    Level: advanced
454 
455 .seealso: PetscDTLegendreEval()
456 @*/
457 PetscErrorCode PetscDTReconstructPoly(PetscInt degree,PetscInt nsource,const PetscReal *sourcex,PetscInt ntarget,const PetscReal *targetx,PetscReal *R)
458 {
459   PetscErrorCode ierr;
460   PetscInt i,j,k,*bdegrees,worksize;
461   PetscReal xmin,xmax,center,hscale,*sourcey,*targety,*Bsource,*Bsinv,*Btarget;
462   PetscScalar *tau,*work;
463 
464   PetscFunctionBegin;
465   PetscValidRealPointer(sourcex,3);
466   PetscValidRealPointer(targetx,5);
467   PetscValidRealPointer(R,6);
468   if (degree >= nsource) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Reconstruction degree %D must be less than number of source intervals %D",degree,nsource);
469 #if defined(PETSC_USE_DEBUG)
470   for (i=0; i<nsource; i++) {
471     if (sourcex[i] >= sourcex[i+1]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Source interval %D has negative orientation (%G,%G)",i,sourcex[i],sourcex[i+1]);
472   }
473   for (i=0; i<ntarget; i++) {
474     if (targetx[i] >= targetx[i+1]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Target interval %D has negative orientation (%G,%G)",i,targetx[i],targetx[i+1]);
475   }
476 #endif
477   xmin = PetscMin(sourcex[0],targetx[0]);
478   xmax = PetscMax(sourcex[nsource],targetx[ntarget]);
479   center = (xmin + xmax)/2;
480   hscale = (xmax - xmin)/2;
481   worksize = nsource;
482   ierr = PetscMalloc4(degree+1,PetscInt,&bdegrees,nsource+1,PetscReal,&sourcey,nsource*(degree+1),PetscReal,&Bsource,worksize,PetscScalar,&work);CHKERRQ(ierr);
483   ierr = PetscMalloc4(nsource,PetscScalar,&tau,nsource*(degree+1),PetscReal,&Bsinv,ntarget+1,PetscReal,&targety,ntarget*(degree+1),PetscReal,&Btarget);CHKERRQ(ierr);
484   for (i=0; i<=nsource; i++) sourcey[i] = (sourcex[i]-center)/hscale;
485   for (i=0; i<=degree; i++) bdegrees[i] = i+1;
486   ierr = PetscDTLegendreIntegrate(nsource,sourcey,degree+1,bdegrees,PETSC_TRUE,Bsource);CHKERRQ(ierr);
487   ierr = PetscDTPseudoInverseQR(nsource,nsource,degree+1,Bsource,Bsinv,tau,nsource,work);CHKERRQ(ierr);
488   for (i=0; i<=ntarget; i++) targety[i] = (targetx[i]-center)/hscale;
489   ierr = PetscDTLegendreIntegrate(ntarget,targety,degree+1,bdegrees,PETSC_FALSE,Btarget);CHKERRQ(ierr);
490   for (i=0; i<ntarget; i++) {
491     PetscReal rowsum = 0;
492     for (j=0; j<nsource; j++) {
493       PetscReal sum = 0;
494       for (k=0; k<degree+1; k++) {
495         sum += Btarget[i*(degree+1)+k] * Bsinv[k*nsource+j];
496       }
497       R[i*nsource+j] = sum;
498       rowsum += sum;
499     }
500     for (j=0; j<nsource; j++) R[i*nsource+j] /= rowsum; /* normalize each row */
501   }
502   ierr = PetscFree4(bdegrees,sourcey,Bsource,work);CHKERRQ(ierr);
503   ierr = PetscFree4(tau,Bsinv,targety,Btarget);CHKERRQ(ierr);
504   PetscFunctionReturn(0);
505 }
506 
507 /* Basis Jet Tabulation
508 
509 We would like to tabulate the nodal basis functions and derivatives at a set of points, usually quadrature points. We
510 follow here the derviation in http://www.math.ttu.edu/~kirby/papers/fiat-toms-2004.pdf. The nodal basis $\psi_i$ can
511 be expressed in terms of a prime basis $\phi_i$ which can be stably evaluated. In PETSc, we will use the Legendre basis
512 as a prime basis.
513 
514   \psi_i = \sum_k \alpha_{ki} \phi_k
515 
516 Our nodal basis is defined in terms of the dual basis $n_j$
517 
518   n_j \cdot \psi_i = \delta_{ji}
519 
520 and we may act on the first equation to obtain
521 
522   n_j \cdot \psi_i = \sum_k \alpha_{ki} n_j \cdot \phi_k
523        \delta_{ji} = \sum_k \alpha_{ki} V_{jk}
524                  I = V \alpha
525 
526 so the coefficients of the nodal basis in the prime basis are
527 
528    \alpha = V^{-1}
529 
530 We will define the dual basis vectors $n_j$ using a quadrature rule.
531 
532 Right now, we will just use the polynomial spaces P^k. I know some elements use the space of symmetric polynomials
533 (I think Nedelec), but we will neglect this for now. Constraints in the space, e.g. Arnold-Winther elements, can
534 be implemented exactly as in FIAT using functionals $L_j$.
535 
536 I will have to count the degrees correctly for the Legendre product when we are on simplices.
537 
538 We will have three objects:
539  - Space, P: this just need point evaluation I think
540  - Dual Space, P'+K: This looks like a set of functionals that can act on members of P, each n is defined by a Q
541  - FEM: This keeps {P, P', Q}
542 */
543 #include <petsc-private/petscfeimpl.h>
544 
545 PetscInt PETSCSPACE_CLASSID = 0;
546 
547 PetscFunctionList PetscSpaceList              = NULL;
548 PetscBool         PetscSpaceRegisterAllCalled = PETSC_FALSE;
549 
550 #undef __FUNCT__
551 #define __FUNCT__ "PetscSpaceRegister"
552 /*@C
553   PetscSpaceRegister - Adds a new PetscSpace implementation
554 
555   Not Collective
556 
557   Input Parameters:
558 + name        - The name of a new user-defined creation routine
559 - create_func - The creation routine itself
560 
561   Notes:
562   PetscSpaceRegister() may be called multiple times to add several user-defined PetscSpaces
563 
564   Sample usage:
565 .vb
566     PetscSpaceRegister("my_space", MyPetscSpaceCreate);
567 .ve
568 
569   Then, your PetscSpace type can be chosen with the procedural interface via
570 .vb
571     PetscSpaceCreate(MPI_Comm, PetscSpace *);
572     PetscSpaceSetType(PetscSpace, "my_space");
573 .ve
574    or at runtime via the option
575 .vb
576     -petscspace_type my_space
577 .ve
578 
579   Level: advanced
580 
581 .keywords: PetscSpace, register
582 .seealso: PetscSpaceRegisterAll(), PetscSpaceRegisterDestroy()
583 
584 @*/
585 PetscErrorCode PetscSpaceRegister(const char sname[], PetscErrorCode (*function)(PetscSpace))
586 {
587   PetscErrorCode ierr;
588 
589   PetscFunctionBegin;
590   ierr = PetscFunctionListAdd(&PetscSpaceList, sname, function);CHKERRQ(ierr);
591   PetscFunctionReturn(0);
592 }
593 
594 #undef __FUNCT__
595 #define __FUNCT__ "PetscSpaceSetType"
596 /*@C
597   PetscSpaceSetType - Builds a particular PetscSpace
598 
599   Collective on PetscSpace
600 
601   Input Parameters:
602 + sp   - The PetscSpace object
603 - name - The kind of space
604 
605   Options Database Key:
606 . -petscspace_type <type> - Sets the PetscSpace type; use -help for a list of available types
607 
608   Level: intermediate
609 
610 .keywords: PetscSpace, set, type
611 .seealso: PetscSpaceGetType(), PetscSpaceCreate()
612 @*/
613 PetscErrorCode PetscSpaceSetType(PetscSpace sp, PetscSpaceType name)
614 {
615   PetscErrorCode (*r)(PetscSpace);
616   PetscBool      match;
617   PetscErrorCode ierr;
618 
619   PetscFunctionBegin;
620   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
621   ierr = PetscObjectTypeCompare((PetscObject) sp, name, &match);CHKERRQ(ierr);
622   if (match) PetscFunctionReturn(0);
623 
624   if (!PetscSpaceRegisterAllCalled) {ierr = PetscSpaceRegisterAll();CHKERRQ(ierr);}
625   ierr = PetscFunctionListFind(PetscSpaceList, name, &r);CHKERRQ(ierr);
626   if (!r) SETERRQ1(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscSpace type: %s", name);
627 
628   if (sp->ops->destroy) {
629     ierr             = (*sp->ops->destroy)(sp);CHKERRQ(ierr);
630     sp->ops->destroy = NULL;
631   }
632   ierr = (*r)(sp);CHKERRQ(ierr);
633   ierr = PetscObjectChangeTypeName((PetscObject) sp, name);CHKERRQ(ierr);
634   PetscFunctionReturn(0);
635 }
636 
637 #undef __FUNCT__
638 #define __FUNCT__ "PetscSpaceGetType"
639 /*@C
640   PetscSpaceGetType - Gets the PetscSpace type name (as a string) from the object.
641 
642   Not Collective
643 
644   Input Parameter:
645 . dm  - The PetscSpace
646 
647   Output Parameter:
648 . name - The PetscSpace type name
649 
650   Level: intermediate
651 
652 .keywords: PetscSpace, get, type, name
653 .seealso: PetscSpaceSetType(), PetscSpaceCreate()
654 @*/
655 PetscErrorCode PetscSpaceGetType(PetscSpace sp, PetscSpaceType *name)
656 {
657   PetscErrorCode ierr;
658 
659   PetscFunctionBegin;
660   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
661   PetscValidCharPointer(name, 2);
662   if (!PetscSpaceRegisterAllCalled) {
663     ierr = PetscSpaceRegisterAll();CHKERRQ(ierr);
664   }
665   *name = ((PetscObject) sp)->type_name;
666   PetscFunctionReturn(0);
667 }
668 
669 #undef __FUNCT__
670 #define __FUNCT__ "PetscSpaceView"
671 /*@C
672   PetscSpaceView - Views a PetscSpace
673 
674   Collective on PetscSpace
675 
676   Input Parameter:
677 + sp - the PetscSpace object to view
678 - v  - the viewer
679 
680   Level: developer
681 
682 .seealso PetscSpaceDestroy()
683 @*/
684 PetscErrorCode PetscSpaceView(PetscSpace sp, PetscViewer v)
685 {
686   PetscErrorCode ierr;
687 
688   PetscFunctionBegin;
689   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
690   if (!v) {
691     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) sp), &v);CHKERRQ(ierr);
692   }
693   if (sp->ops->view) {
694     ierr = (*sp->ops->view)(sp, v);CHKERRQ(ierr);
695   }
696   PetscFunctionReturn(0);
697 }
698 
699 #undef __FUNCT__
700 #define __FUNCT__ "PetscSpaceViewFromOptions"
701 /*
702   PetscSpaceViewFromOptions - Processes command line options to determine if/how a PetscSpace is to be viewed.
703 
704   Collective on PetscSpace
705 
706   Input Parameters:
707 + sp   - the PetscSpace
708 . prefix - prefix to use for viewing, or NULL to use prefix of 'rnd'
709 - optionname - option to activate viewing
710 
711   Level: intermediate
712 
713 .keywords: PetscSpace, view, options, database
714 .seealso: VecViewFromOptions(), MatViewFromOptions()
715 */
716 PetscErrorCode PetscSpaceViewFromOptions(PetscSpace sp, const char prefix[], const char optionname[])
717 {
718   PetscViewer       viewer;
719   PetscViewerFormat format;
720   PetscBool         flg;
721   PetscErrorCode    ierr;
722 
723   PetscFunctionBegin;
724   if (prefix) {
725     ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject) sp), prefix, optionname, &viewer, &format, &flg);CHKERRQ(ierr);
726   } else {
727     ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject) sp), ((PetscObject) sp)->prefix, optionname, &viewer, &format, &flg);CHKERRQ(ierr);
728   }
729   if (flg) {
730     ierr = PetscViewerPushFormat(viewer, format);CHKERRQ(ierr);
731     ierr = PetscSpaceView(sp, viewer);CHKERRQ(ierr);
732     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
733     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
734   }
735   PetscFunctionReturn(0);
736 }
737 
738 #undef __FUNCT__
739 #define __FUNCT__ "PetscSpaceSetFromOptions"
740 /*@
741   PetscSpaceSetFromOptions - sets parameters in a PetscSpace from the options database
742 
743   Collective on PetscSpace
744 
745   Input Parameter:
746 . sp - the PetscSpace object to set options for
747 
748   Options Database:
749 . -petscspace_order the approximation order of the space
750 
751   Level: developer
752 
753 .seealso PetscSpaceView()
754 @*/
755 PetscErrorCode PetscSpaceSetFromOptions(PetscSpace sp)
756 {
757   const char    *defaultType;
758   char           typename[256];
759   PetscBool      flg;
760   PetscErrorCode ierr;
761 
762   PetscFunctionBegin;
763   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
764   if (!((PetscObject) sp)->type_name) {
765     defaultType = PETSCSPACEPOLYNOMIAL;
766   } else {
767     defaultType = ((PetscObject) sp)->type_name;
768   }
769   if (!PetscSpaceRegisterAllCalled) {ierr = PetscSpaceRegisterAll();CHKERRQ(ierr);}
770 
771   ierr = PetscObjectOptionsBegin((PetscObject) sp);CHKERRQ(ierr);
772   ierr = PetscOptionsList("-petscspace_type", "Linear space", "PetscSpaceSetType", PetscSpaceList, defaultType, typename, 256, &flg);CHKERRQ(ierr);
773   if (flg) {
774     ierr = PetscSpaceSetType(sp, typename);CHKERRQ(ierr);
775   } else if (!((PetscObject) sp)->type_name) {
776     ierr = PetscSpaceSetType(sp, defaultType);CHKERRQ(ierr);
777   }
778   ierr = PetscOptionsInt("-petscspace_order", "The approximation order", "PetscSpaceSetOrder", sp->order, &sp->order, NULL);CHKERRQ(ierr);
779   if (sp->ops->setfromoptions) {
780     ierr = (*sp->ops->setfromoptions)(sp);CHKERRQ(ierr);
781   }
782   /* process any options handlers added with PetscObjectAddOptionsHandler() */
783   ierr = PetscObjectProcessOptionsHandlers((PetscObject) sp);CHKERRQ(ierr);
784   ierr = PetscOptionsEnd();CHKERRQ(ierr);
785   ierr = PetscSpaceViewFromOptions(sp, NULL, "-petscspace_view");CHKERRQ(ierr);
786   PetscFunctionReturn(0);
787 }
788 
789 #undef __FUNCT__
790 #define __FUNCT__ "PetscSpaceSetUp"
791 /*@C
792   PetscSpaceSetUp - Construct data structures for the PetscSpace
793 
794   Collective on PetscSpace
795 
796   Input Parameter:
797 . sp - the PetscSpace object to setup
798 
799   Level: developer
800 
801 .seealso PetscSpaceView(), PetscSpaceDestroy()
802 @*/
803 PetscErrorCode PetscSpaceSetUp(PetscSpace sp)
804 {
805   PetscErrorCode ierr;
806 
807   PetscFunctionBegin;
808   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
809   if (sp->ops->setup) {ierr = (*sp->ops->setup)(sp);CHKERRQ(ierr);}
810   PetscFunctionReturn(0);
811 }
812 
813 #undef __FUNCT__
814 #define __FUNCT__ "PetscSpaceDestroy"
815 /*@
816   PetscSpaceDestroy - Destroys a PetscSpace object
817 
818   Collective on PetscSpace
819 
820   Input Parameter:
821 . sp - the PetscSpace object to destroy
822 
823   Level: developer
824 
825 .seealso PetscSpaceView()
826 @*/
827 PetscErrorCode PetscSpaceDestroy(PetscSpace *sp)
828 {
829   PetscErrorCode ierr;
830 
831   PetscFunctionBegin;
832   if (!*sp) PetscFunctionReturn(0);
833   PetscValidHeaderSpecific((*sp), PETSCSPACE_CLASSID, 1);
834 
835   if (--((PetscObject)(*sp))->refct > 0) {*sp = 0; PetscFunctionReturn(0);}
836   ((PetscObject) (*sp))->refct = 0;
837   /* if memory was published with AMS then destroy it */
838   ierr = PetscObjectAMSViewOff((PetscObject) *sp);CHKERRQ(ierr);
839 
840   ierr = DMDestroy(&(*sp)->dm);CHKERRQ(ierr);
841 
842   ierr = (*(*sp)->ops->destroy)(*sp);CHKERRQ(ierr);
843   ierr = PetscHeaderDestroy(sp);CHKERRQ(ierr);
844   PetscFunctionReturn(0);
845 }
846 
847 #undef __FUNCT__
848 #define __FUNCT__ "PetscSpaceCreate"
849 /*@
850   PetscSpaceCreate - Creates an empty PetscSpace object. The type can then be set with PetscSpaceSetType().
851 
852   Collective on MPI_Comm
853 
854   Input Parameter:
855 . comm - The communicator for the PetscSpace object
856 
857   Output Parameter:
858 . sp - The PetscSpace object
859 
860   Level: beginner
861 
862 .seealso: PetscSpaceSetType(), PETSCSPACEPOLYNOMIAL
863 @*/
864 PetscErrorCode PetscSpaceCreate(MPI_Comm comm, PetscSpace *sp)
865 {
866   PetscSpace     s;
867   PetscErrorCode ierr;
868 
869   PetscFunctionBegin;
870   PetscValidPointer(sp, 2);
871   *sp = NULL;
872 #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
873   ierr = PetscFEInitializePackage();CHKERRQ(ierr);
874 #endif
875 
876   ierr = PetscHeaderCreate(s, _p_PetscSpace, struct _PetscSpaceOps, PETSCSPACE_CLASSID, "PetscSpace", "Linear Space", "PetscSpace", comm, PetscSpaceDestroy, PetscSpaceView);CHKERRQ(ierr);
877   ierr = PetscMemzero(s->ops, sizeof(struct _PetscSpaceOps));CHKERRQ(ierr);
878 
879   s->order = 0;
880   ierr = DMShellCreate(comm, &s->dm);CHKERRQ(ierr);
881 
882   *sp = s;
883   PetscFunctionReturn(0);
884 }
885 
886 #undef __FUNCT__
887 #define __FUNCT__ "PetscSpaceGetDimension"
888 /* Dimension of the space, i.e. number of basis vectors */
889 PetscErrorCode PetscSpaceGetDimension(PetscSpace sp, PetscInt *dim)
890 {
891   PetscErrorCode ierr;
892 
893   PetscFunctionBegin;
894   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
895   PetscValidPointer(dim, 2);
896   *dim = 0;
897   if (sp->ops->getdimension) {ierr = (*sp->ops->getdimension)(sp, dim);CHKERRQ(ierr);}
898   PetscFunctionReturn(0);
899 }
900 
901 #undef __FUNCT__
902 #define __FUNCT__ "PetscSpaceGetOrder"
903 PetscErrorCode PetscSpaceGetOrder(PetscSpace sp, PetscInt *order)
904 {
905   PetscFunctionBegin;
906   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
907   PetscValidPointer(order, 2);
908   *order = sp->order;
909   PetscFunctionReturn(0);
910 }
911 
912 #undef __FUNCT__
913 #define __FUNCT__ "PetscSpaceSetOrder"
914 PetscErrorCode PetscSpaceSetOrder(PetscSpace sp, PetscInt order)
915 {
916   PetscFunctionBegin;
917   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
918   sp->order = order;
919   PetscFunctionReturn(0);
920 }
921 
922 #undef __FUNCT__
923 #define __FUNCT__ "PetscSpaceEvaluate"
924 PetscErrorCode PetscSpaceEvaluate(PetscSpace sp, PetscInt npoints, const PetscReal points[], PetscReal B[], PetscReal D[], PetscReal H[])
925 {
926   PetscErrorCode ierr;
927 
928   PetscFunctionBegin;
929   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
930   PetscValidPointer(points, 3);
931   if (B) PetscValidPointer(B, 4);
932   if (D) PetscValidPointer(D, 5);
933   if (H) PetscValidPointer(H, 6);
934   if (sp->ops->evaluate) {ierr = (*sp->ops->evaluate)(sp, npoints, points, B, D, H);CHKERRQ(ierr);}
935   PetscFunctionReturn(0);
936 }
937 
938 #undef __FUNCT__
939 #define __FUNCT__ "PetscSpaceSetFromOptions_Polynomial"
940 PetscErrorCode PetscSpaceSetFromOptions_Polynomial(PetscSpace sp)
941 {
942   PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
943   PetscErrorCode   ierr;
944 
945   PetscFunctionBegin;
946   ierr = PetscObjectOptionsBegin((PetscObject) sp);CHKERRQ(ierr);
947   ierr = PetscOptionsInt("-petscspace_num_variables", "The number of different variables, e.g. x and y", "PetscSpacePolynomialSetNumVariables", poly->numVariables, &poly->numVariables, NULL);CHKERRQ(ierr);
948   ierr = PetscOptionsBool("-petscspace_poly_sym", "Use only symmetric polynomials", "PetscSpacePolynomialSetSymmetric", poly->symmetric, &poly->symmetric, NULL);CHKERRQ(ierr);
949   ierr = PetscOptionsEnd();CHKERRQ(ierr);
950   PetscFunctionReturn(0);
951 }
952 
953 #undef __FUNCT__
954 #define __FUNCT__ "PetscSpacePolynomialView_Ascii"
955 PetscErrorCode PetscSpacePolynomialView_Ascii(PetscSpace sp, PetscViewer viewer)
956 {
957   PetscSpace_Poly  *poly = (PetscSpace_Poly *) sp->data;
958   PetscViewerFormat format;
959   PetscErrorCode    ierr;
960 
961   PetscFunctionBegin;
962   ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
963   if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
964     ierr = PetscViewerASCIIPrintf(viewer, "Polynomial space in %d variables of order %d", poly->numVariables, sp->order);CHKERRQ(ierr);
965   } else {
966     ierr = PetscViewerASCIIPrintf(viewer, "Polynomial space in %d variables of order %d", poly->numVariables, sp->order);CHKERRQ(ierr);
967   }
968   PetscFunctionReturn(0);
969 }
970 
971 #undef __FUNCT__
972 #define __FUNCT__ "PetscSpaceView_Polynomial"
973 PetscErrorCode PetscSpaceView_Polynomial(PetscSpace sp, PetscViewer viewer)
974 {
975   PetscBool      iascii;
976   PetscErrorCode ierr;
977 
978   PetscFunctionBegin;
979   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
980   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
981   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
982   if (iascii) {ierr = PetscSpacePolynomialView_Ascii(sp, viewer);CHKERRQ(ierr);}
983   PetscFunctionReturn(0);
984 }
985 
986 #undef __FUNCT__
987 #define __FUNCT__ "PetscSpaceSetUp_Polynomial"
988 PetscErrorCode PetscSpaceSetUp_Polynomial(PetscSpace sp)
989 {
990   PetscSpace_Poly *poly    = (PetscSpace_Poly *) sp->data;
991   PetscInt         ndegree = sp->order+1;
992   PetscInt         dim     = poly->numVariables;
993   PetscInt         pdim, deg;
994   PetscErrorCode   ierr;
995 
996   PetscFunctionBegin;
997   ierr = PetscMalloc(ndegree * sizeof(PetscInt), &poly->degrees);CHKERRQ(ierr);
998   for (deg = 0; deg < ndegree; ++deg) poly->degrees[deg] = deg;
999   ierr = PetscSpaceGetDimension(sp, &pdim);CHKERRQ(ierr);
1000   /* Create A (pdim x ndegree * dim) */
1001   ierr = PetscMalloc(pdim*ndegree*dim * sizeof(PetscReal), &poly->A);CHKERRQ(ierr);
1002   ierr = PetscMemzero(poly->A, pdim*ndegree*dim * sizeof(PetscReal));CHKERRQ(ierr);
1003   /* Hardcode P_1: Here we need a way to iterate through the basis */
1004   poly->A[(0*ndegree + 0)*dim + 0] = 1.0; /* 1 */
1005   poly->A[(1*ndegree + 1)*dim + 0] = 1.0; /* x */
1006   poly->A[(2*ndegree + 1)*dim + 1] = 1.0; /* y */
1007   PetscFunctionReturn(0);
1008 }
1009 
1010 #undef __FUNCT__
1011 #define __FUNCT__ "PetscSpaceDestroy_Polynomial"
1012 PetscErrorCode PetscSpaceDestroy_Polynomial(PetscSpace sp)
1013 {
1014   PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
1015   PetscErrorCode   ierr;
1016 
1017   PetscFunctionBegin;
1018   ierr = PetscFree(poly->A);CHKERRQ(ierr);
1019   ierr = PetscFree(poly->degrees);CHKERRQ(ierr);
1020   ierr = PetscFree(poly);CHKERRQ(ierr);
1021   PetscFunctionReturn(0);
1022 }
1023 
1024 #undef __FUNCT__
1025 #define __FUNCT__ "PetscSpaceGetDimension_Polynomial"
1026 PetscErrorCode PetscSpaceGetDimension_Polynomial(PetscSpace sp, PetscInt *dim)
1027 {
1028   PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
1029   PetscInt         deg  = sp->order;
1030   PetscInt         n    = poly->numVariables, i;
1031   PetscReal        D    = 1.0;
1032 
1033   PetscFunctionBegin;
1034   for (i = 1; i <= n; ++i) {
1035     D *= ((PetscReal) (deg+i))/i;
1036   }
1037   *dim = (PetscInt) (D + 0.5);
1038   PetscFunctionReturn(0);
1039 }
1040 
1041 #undef __FUNCT__
1042 #define __FUNCT__ "LatticePoint_Internal"
1043 /*
1044   LatticePoint_Internal - Returns all tuples of size 'len' with nonnegative integers that sum up to 'sum'.
1045 
1046   Input Parameters:
1047 + len - The length of the tuple
1048 . sum - The sum of all entries in the tuple
1049 - ind - The current multi-index of the tuple, initialized to the 0 tuple
1050 
1051   Output Parameter:
1052 + ind - The multi-index of the tuple, -1 indicates the iteration has terminated
1053 . tup - A tuple of len integers addig to sum
1054 
1055   Level: developer
1056 
1057 .seealso:
1058 */
1059 static PetscErrorCode LatticePoint_Internal(PetscInt len, PetscInt sum, PetscInt ind[], PetscInt tup[])
1060 {
1061   PetscInt       i;
1062   PetscErrorCode ierr;
1063 
1064   PetscFunctionBegin;
1065   if (len == 1) {
1066     ind[0] = -1;
1067     tup[0] = sum;
1068   } else if (sum == 0) {
1069     for (i = 0; i < len; ++i) {ind[0] = -1; tup[i] = 0;}
1070   } else {
1071     tup[0] = sum - ind[0];
1072     ierr = LatticePoint_Internal(len-1, ind[0], &ind[1], &tup[1]);CHKERRQ(ierr);
1073     if (ind[1] < 0) {
1074       if (ind[0] == sum) {ind[0] = -1;}
1075       else               {ind[1] = 0; ++ind[0];}
1076     }
1077   }
1078   PetscFunctionReturn(0);
1079 }
1080 
1081 #undef __FUNCT__
1082 #define __FUNCT__ "PetscSpaceEvaluate_Polynomial"
1083 PetscErrorCode PetscSpaceEvaluate_Polynomial(PetscSpace sp, PetscInt npoints, const PetscReal points[], PetscReal B[], PetscReal D[], PetscReal H[])
1084 {
1085   PetscSpace_Poly *poly    = (PetscSpace_Poly *) sp->data;
1086   DM               dm      = sp->dm;
1087   PetscInt         ndegree = sp->order+1;
1088   PetscInt        *degrees = poly->degrees;
1089   PetscInt         dim     = poly->numVariables;
1090   PetscReal       *A       = poly->A;
1091   PetscReal       *lpoints, *tmp, *LB, *LD, *LH;
1092   PetscInt        *ind, *tup;
1093   PetscInt         pdim, d, i, p, deg, o;
1094   PetscErrorCode   ierr;
1095 
1096   PetscFunctionBegin;
1097   ierr = PetscSpaceGetDimension(sp, &pdim);CHKERRQ(ierr);
1098   ierr = DMGetWorkArray(dm, npoints, PETSC_REAL, &lpoints);CHKERRQ(ierr);
1099   ierr = DMGetWorkArray(dm, npoints*ndegree*3, PETSC_REAL, &tmp);CHKERRQ(ierr);
1100   if (B) {ierr = DMGetWorkArray(dm, npoints*dim*ndegree, PETSC_REAL, &LB);CHKERRQ(ierr);}
1101   if (D) {ierr = DMGetWorkArray(dm, npoints*dim*ndegree, PETSC_REAL, &LD);CHKERRQ(ierr);}
1102   if (H) {ierr = DMGetWorkArray(dm, npoints*dim*ndegree, PETSC_REAL, &LH);CHKERRQ(ierr);}
1103   for (d = 0; d < dim; ++d) {
1104     for (p = 0; p < npoints; ++p) {
1105       lpoints[p] = points[p*dim+d];
1106     }
1107     ierr = PetscDTLegendreEval(npoints, lpoints, ndegree, degrees, tmp, &tmp[1*npoints*ndegree], &tmp[2*npoints*ndegree]);CHKERRQ(ierr);
1108     /* LB, LD, LH (ndegree * dim x npoints) */
1109     for (deg = 0; deg < ndegree; ++deg) {
1110       for (p = 0; p < npoints; ++p) {
1111         if (B) LB[(deg*dim + d)*npoints + p] = tmp[(0*npoints + p)*ndegree+deg];
1112         if (D) LD[(deg*dim + d)*npoints + p] = tmp[(1*npoints + p)*ndegree+deg];
1113         if (H) LH[(deg*dim + d)*npoints + p] = tmp[(2*npoints + p)*ndegree+deg];
1114       }
1115     }
1116   }
1117   /* Multiply by A (pdim x ndegree * dim) */
1118   ierr = PetscMalloc2(dim,PetscInt,&ind,dim,PetscInt,&tup);CHKERRQ(ierr);
1119   if (B) {
1120     /* B (npoints x pdim) */
1121     i = 0;
1122     for (o = 0; o <= sp->order; ++o) {
1123       ierr = PetscMemzero(ind, dim * sizeof(PetscInt));CHKERRQ(ierr);
1124       while (ind[0] >= 0) {
1125         ierr = LatticePoint_Internal(dim, o, ind, tup);CHKERRQ(ierr);
1126         for (p = 0; p < npoints; ++p) {
1127           B[p*pdim + i] = 1.0;
1128           for (d = 0; d < dim; ++d) {
1129             B[p*pdim + i] *= LB[(tup[d]*dim + d)*npoints + p];
1130           }
1131         }
1132         ++i;
1133       }
1134     }
1135   }
1136   ierr = PetscFree2(ind,tup);CHKERRQ(ierr);
1137   if (D) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to code first derivatives");
1138   if (H) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to code second derivatives");
1139   if (B) {ierr = DMRestoreWorkArray(dm, npoints*dim*ndegree, PETSC_REAL, &LB);CHKERRQ(ierr);}
1140   if (D) {ierr = DMRestoreWorkArray(dm, npoints*dim*ndegree, PETSC_REAL, &LD);CHKERRQ(ierr);}
1141   if (H) {ierr = DMRestoreWorkArray(dm, npoints*dim*ndegree, PETSC_REAL, &LH);CHKERRQ(ierr);}
1142   ierr = DMRestoreWorkArray(dm, npoints*ndegree*3, PETSC_REAL, &tmp);CHKERRQ(ierr);
1143   ierr = DMRestoreWorkArray(dm, npoints, PETSC_REAL, &lpoints);CHKERRQ(ierr);
1144   PetscFunctionReturn(0);
1145 }
1146 
1147 #undef __FUNCT__
1148 #define __FUNCT__ "PetscSpaceInitialize_Polynomial"
1149 PetscErrorCode PetscSpaceInitialize_Polynomial(PetscSpace sp)
1150 {
1151   PetscFunctionBegin;
1152   sp->ops->setfromoptions = PetscSpaceSetFromOptions_Polynomial;
1153   sp->ops->setup          = PetscSpaceSetUp_Polynomial;
1154   sp->ops->view           = PetscSpaceView_Polynomial;
1155   sp->ops->destroy        = PetscSpaceDestroy_Polynomial;
1156   sp->ops->getdimension   = PetscSpaceGetDimension_Polynomial;
1157   sp->ops->evaluate       = PetscSpaceEvaluate_Polynomial;
1158   PetscFunctionReturn(0);
1159 }
1160 
1161 /*MC
1162   PETSCSPACEPOLYNOMIAL = "poly" - A PetscSpace object that encapsulates a polynomial space, e.g. P1 is the space of linear polynomials.
1163 
1164   Level: intermediate
1165 
1166 .seealso: PetscSpaceType, PetscSpaceCreate(), PetscSpaceSetType()
1167 M*/
1168 
1169 #undef __FUNCT__
1170 #define __FUNCT__ "PetscSpaceCreate_Polynomial"
1171 PETSC_EXTERN PetscErrorCode PetscSpaceCreate_Polynomial(PetscSpace sp)
1172 {
1173   PetscSpace_Poly *poly;
1174   PetscErrorCode   ierr;
1175 
1176   PetscFunctionBegin;
1177   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
1178   ierr     = PetscNewLog(sp, PetscSpace_Poly, &poly);CHKERRQ(ierr);
1179   sp->data = poly;
1180 
1181   poly->numVariables = 0;
1182   poly->symmetric    = PETSC_FALSE;
1183   poly->degrees      = NULL;
1184   poly->A            = NULL;
1185 
1186   ierr = PetscSpaceInitialize_Polynomial(sp);CHKERRQ(ierr);
1187   PetscFunctionReturn(0);
1188 }
1189 
1190 #undef __FUNCT__
1191 #define __FUNCT__ "PetscSpacePolynomialSetSymmetric"
1192 PetscErrorCode PetscSpacePolynomialSetSymmetric(PetscSpace sp, PetscBool sym)
1193 {
1194   PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
1195 
1196   PetscFunctionBegin;
1197   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
1198   poly->symmetric = sym;
1199   PetscFunctionReturn(0);
1200 }
1201 
1202 #undef __FUNCT__
1203 #define __FUNCT__ "PetscSpacePolynomialGetSymmetric"
1204 PetscErrorCode PetscSpacePolynomialGetSymmetric(PetscSpace sp, PetscBool *sym)
1205 {
1206   PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
1207 
1208   PetscFunctionBegin;
1209   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
1210   PetscValidPointer(sym, 2);
1211   *sym = poly->symmetric;
1212   PetscFunctionReturn(0);
1213 }
1214 
1215 #undef __FUNCT__
1216 #define __FUNCT__ "PetscSpacePolynomialSetNumVariables"
1217 PetscErrorCode PetscSpacePolynomialSetNumVariables(PetscSpace sp, PetscInt n)
1218 {
1219   PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
1220 
1221   PetscFunctionBegin;
1222   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
1223   poly->numVariables = n;
1224   PetscFunctionReturn(0);
1225 }
1226 
1227 #undef __FUNCT__
1228 #define __FUNCT__ "PetscSpacePolynomialGetNumVariables"
1229 PetscErrorCode PetscSpacePolynomialGetNumVariables(PetscSpace sp, PetscInt *n)
1230 {
1231   PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
1232 
1233   PetscFunctionBegin;
1234   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
1235   PetscValidPointer(n, 2);
1236   *n = poly->numVariables;
1237   PetscFunctionReturn(0);
1238 }
1239 
1240 
1241 PetscInt PETSCDUALSPACE_CLASSID = 0;
1242 
1243 PetscFunctionList PetscDualSpaceList              = NULL;
1244 PetscBool         PetscDualSpaceRegisterAllCalled = PETSC_FALSE;
1245 
1246 #undef __FUNCT__
1247 #define __FUNCT__ "PetscDualSpaceRegister"
1248 /*@C
1249   PetscDualSpaceRegister - Adds a new PetscDualSpace implementation
1250 
1251   Not Collective
1252 
1253   Input Parameters:
1254 + name        - The name of a new user-defined creation routine
1255 - create_func - The creation routine itself
1256 
1257   Notes:
1258   PetscDualSpaceRegister() may be called multiple times to add several user-defined PetscDualSpaces
1259 
1260   Sample usage:
1261 .vb
1262     PetscDualSpaceRegister("my_space", MyPetscDualSpaceCreate);
1263 .ve
1264 
1265   Then, your PetscDualSpace type can be chosen with the procedural interface via
1266 .vb
1267     PetscDualSpaceCreate(MPI_Comm, PetscDualSpace *);
1268     PetscDualSpaceSetType(PetscDualSpace, "my_dual_space");
1269 .ve
1270    or at runtime via the option
1271 .vb
1272     -petscdualspace_type my_dual_space
1273 .ve
1274 
1275   Level: advanced
1276 
1277 .keywords: PetscDualSpace, register
1278 .seealso: PetscDualSpaceRegisterAll(), PetscDualSpaceRegisterDestroy()
1279 
1280 @*/
1281 PetscErrorCode PetscDualSpaceRegister(const char sname[], PetscErrorCode (*function)(PetscDualSpace))
1282 {
1283   PetscErrorCode ierr;
1284 
1285   PetscFunctionBegin;
1286   ierr = PetscFunctionListAdd(&PetscDualSpaceList, sname, function);CHKERRQ(ierr);
1287   PetscFunctionReturn(0);
1288 }
1289 
1290 #undef __FUNCT__
1291 #define __FUNCT__ "PetscDualSpaceSetType"
1292 /*@C
1293   PetscDualSpaceSetType - Builds a particular PetscDualSpace
1294 
1295   Collective on PetscDualSpace
1296 
1297   Input Parameters:
1298 + sp   - The PetscDualSpace object
1299 - name - The kind of space
1300 
1301   Options Database Key:
1302 . -petscdualspace_type <type> - Sets the PetscDualSpace type; use -help for a list of available types
1303 
1304   Level: intermediate
1305 
1306 .keywords: PetscDualSpace, set, type
1307 .seealso: PetscDualSpaceGetType(), PetscDualSpaceCreate()
1308 @*/
1309 PetscErrorCode PetscDualSpaceSetType(PetscDualSpace sp, PetscDualSpaceType name)
1310 {
1311   PetscErrorCode (*r)(PetscDualSpace);
1312   PetscBool      match;
1313   PetscErrorCode ierr;
1314 
1315   PetscFunctionBegin;
1316   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1317   ierr = PetscObjectTypeCompare((PetscObject) sp, name, &match);CHKERRQ(ierr);
1318   if (match) PetscFunctionReturn(0);
1319 
1320   if (!PetscDualSpaceRegisterAllCalled) {ierr = PetscDualSpaceRegisterAll();CHKERRQ(ierr);}
1321   ierr = PetscFunctionListFind(PetscDualSpaceList, name, &r);CHKERRQ(ierr);
1322   if (!r) SETERRQ1(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscDualSpace type: %s", name);
1323 
1324   if (sp->ops->destroy) {
1325     ierr             = (*sp->ops->destroy)(sp);CHKERRQ(ierr);
1326     sp->ops->destroy = NULL;
1327   }
1328   ierr = (*r)(sp);CHKERRQ(ierr);
1329   ierr = PetscObjectChangeTypeName((PetscObject) sp, name);CHKERRQ(ierr);
1330   PetscFunctionReturn(0);
1331 }
1332 
1333 #undef __FUNCT__
1334 #define __FUNCT__ "PetscDualSpaceGetType"
1335 /*@C
1336   PetscDualSpaceGetType - Gets the PetscDualSpace type name (as a string) from the object.
1337 
1338   Not Collective
1339 
1340   Input Parameter:
1341 . dm  - The PetscDualSpace
1342 
1343   Output Parameter:
1344 . name - The PetscDualSpace type name
1345 
1346   Level: intermediate
1347 
1348 .keywords: PetscDualSpace, get, type, name
1349 .seealso: PetscDualSpaceSetType(), PetscDualSpaceCreate()
1350 @*/
1351 PetscErrorCode PetscDualSpaceGetType(PetscDualSpace sp, PetscDualSpaceType *name)
1352 {
1353   PetscErrorCode ierr;
1354 
1355   PetscFunctionBegin;
1356   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1357   PetscValidCharPointer(name, 2);
1358   if (!PetscDualSpaceRegisterAllCalled) {
1359     ierr = PetscDualSpaceRegisterAll();CHKERRQ(ierr);
1360   }
1361   *name = ((PetscObject) sp)->type_name;
1362   PetscFunctionReturn(0);
1363 }
1364 
1365 #undef __FUNCT__
1366 #define __FUNCT__ "PetscDualSpaceView"
1367 /*@C
1368   PetscDualSpaceView - Views a PetscDualSpace
1369 
1370   Collective on PetscDualSpace
1371 
1372   Input Parameter:
1373 + sp - the PetscDualSpace object to view
1374 - v  - the viewer
1375 
1376   Level: developer
1377 
1378 .seealso PetscDualSpaceDestroy()
1379 @*/
1380 PetscErrorCode PetscDualSpaceView(PetscDualSpace sp, PetscViewer v)
1381 {
1382   PetscErrorCode ierr;
1383 
1384   PetscFunctionBegin;
1385   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1386   if (!v) {
1387     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) sp), &v);CHKERRQ(ierr);
1388   }
1389   if (sp->ops->view) {
1390     ierr = (*sp->ops->view)(sp, v);CHKERRQ(ierr);
1391   }
1392   PetscFunctionReturn(0);
1393 }
1394 
1395 #undef __FUNCT__
1396 #define __FUNCT__ "PetscDualSpaceViewFromOptions"
1397 /*
1398   PetscDualSpaceViewFromOptions - Processes command line options to determine if/how a PetscDualSpace is to be viewed.
1399 
1400   Collective on PetscDualSpace
1401 
1402   Input Parameters:
1403 + sp   - the PetscDualSpace
1404 . prefix - prefix to use for viewing, or NULL to use prefix of 'rnd'
1405 - optionname - option to activate viewing
1406 
1407   Level: intermediate
1408 
1409 .keywords: PetscDualSpace, view, options, database
1410 .seealso: VecViewFromOptions(), MatViewFromOptions()
1411 */
1412 PetscErrorCode PetscDualSpaceViewFromOptions(PetscDualSpace sp, const char prefix[], const char optionname[])
1413 {
1414   PetscViewer       viewer;
1415   PetscViewerFormat format;
1416   PetscBool         flg;
1417   PetscErrorCode    ierr;
1418 
1419   PetscFunctionBegin;
1420   if (prefix) {
1421     ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject) sp), prefix, optionname, &viewer, &format, &flg);CHKERRQ(ierr);
1422   } else {
1423     ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject) sp), ((PetscObject) sp)->prefix, optionname, &viewer, &format, &flg);CHKERRQ(ierr);
1424   }
1425   if (flg) {
1426     ierr = PetscViewerPushFormat(viewer, format);CHKERRQ(ierr);
1427     ierr = PetscDualSpaceView(sp, viewer);CHKERRQ(ierr);
1428     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
1429     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1430   }
1431   PetscFunctionReturn(0);
1432 }
1433 
1434 #undef __FUNCT__
1435 #define __FUNCT__ "PetscDualSpaceSetFromOptions"
1436 /*@
1437   PetscDualSpaceSetFromOptions - sets parameters in a PetscDualSpace from the options database
1438 
1439   Collective on PetscDualSpace
1440 
1441   Input Parameter:
1442 . sp - the PetscDualSpace object to set options for
1443 
1444   Options Database:
1445 . -petscspace_order the approximation order of the space
1446 
1447   Level: developer
1448 
1449 .seealso PetscDualSpaceView()
1450 @*/
1451 PetscErrorCode PetscDualSpaceSetFromOptions(PetscDualSpace sp)
1452 {
1453   const char    *defaultType;
1454   char           typename[256];
1455   PetscBool      flg;
1456   PetscErrorCode ierr;
1457 
1458   PetscFunctionBegin;
1459   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1460   if (!((PetscObject) sp)->type_name) {
1461     defaultType = PETSCDUALSPACELAGRANGE;
1462   } else {
1463     defaultType = ((PetscObject) sp)->type_name;
1464   }
1465   if (!PetscSpaceRegisterAllCalled) {ierr = PetscSpaceRegisterAll();CHKERRQ(ierr);}
1466 
1467   ierr = PetscObjectOptionsBegin((PetscObject) sp);CHKERRQ(ierr);
1468   ierr = PetscOptionsList("-petscdualspace_type", "Dual space", "PetscDualSpaceSetType", PetscDualSpaceList, defaultType, typename, 256, &flg);CHKERRQ(ierr);
1469   if (flg) {
1470     ierr = PetscDualSpaceSetType(sp, typename);CHKERRQ(ierr);
1471   } else if (!((PetscObject) sp)->type_name) {
1472     ierr = PetscDualSpaceSetType(sp, defaultType);CHKERRQ(ierr);
1473   }
1474   ierr = PetscOptionsInt("-petscdualspace_order", "The approximation order", "PetscDualSpaceSetOrder", sp->order, &sp->order, NULL);CHKERRQ(ierr);
1475   if (sp->ops->setfromoptions) {
1476     ierr = (*sp->ops->setfromoptions)(sp);CHKERRQ(ierr);
1477   }
1478   /* process any options handlers added with PetscObjectAddOptionsHandler() */
1479   ierr = PetscObjectProcessOptionsHandlers((PetscObject) sp);CHKERRQ(ierr);
1480   ierr = PetscOptionsEnd();CHKERRQ(ierr);
1481   ierr = PetscDualSpaceViewFromOptions(sp, NULL, "-petscdualspace_view");CHKERRQ(ierr);
1482   PetscFunctionReturn(0);
1483 }
1484 
1485 #undef __FUNCT__
1486 #define __FUNCT__ "PetscDualSpaceSetUp"
1487 /*@C
1488   PetscDualSpaceSetUp - Construct a basis for the PetscDualSpace
1489 
1490   Collective on PetscDualSpace
1491 
1492   Input Parameter:
1493 . sp - the PetscDualSpace object to setup
1494 
1495   Level: developer
1496 
1497 .seealso PetscDualSpaceView(), PetscDualSpaceDestroy()
1498 @*/
1499 PetscErrorCode PetscDualSpaceSetUp(PetscDualSpace sp)
1500 {
1501   PetscErrorCode ierr;
1502 
1503   PetscFunctionBegin;
1504   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1505   if (sp->ops->setup) {ierr = (*sp->ops->setup)(sp);CHKERRQ(ierr);}
1506   PetscFunctionReturn(0);
1507 }
1508 
1509 #undef __FUNCT__
1510 #define __FUNCT__ "PetscDualSpaceDestroy"
1511 /*@
1512   PetscDualSpaceDestroy - Destroys a PetscDualSpace object
1513 
1514   Collective on PetscDualSpace
1515 
1516   Input Parameter:
1517 . sp - the PetscDualSpace object to destroy
1518 
1519   Level: developer
1520 
1521 .seealso PetscDualSpaceView()
1522 @*/
1523 PetscErrorCode PetscDualSpaceDestroy(PetscDualSpace *sp)
1524 {
1525   PetscInt       dim, f;
1526   PetscErrorCode ierr;
1527 
1528   PetscFunctionBegin;
1529   if (!*sp) PetscFunctionReturn(0);
1530   PetscValidHeaderSpecific((*sp), PETSCDUALSPACE_CLASSID, 1);
1531 
1532   if (--((PetscObject)(*sp))->refct > 0) {*sp = 0; PetscFunctionReturn(0);}
1533   ((PetscObject) (*sp))->refct = 0;
1534   /* if memory was published with AMS then destroy it */
1535   ierr = PetscObjectAMSViewOff((PetscObject) *sp);CHKERRQ(ierr);
1536 
1537   ierr = DMDestroy(&(*sp)->dm);CHKERRQ(ierr);
1538   ierr = PetscDualSpaceGetDimension(*sp, &dim);CHKERRQ(ierr);
1539   for (f = 0; f < dim; ++f) {
1540     /* ierr = PetscQuadratureDestroy((*sp)->functional[f]);CHKERRQ(ierr); */
1541   }
1542   ierr = PetscFree((*sp)->functional);CHKERRQ(ierr);
1543 
1544   ierr = (*(*sp)->ops->destroy)(*sp);CHKERRQ(ierr);
1545   ierr = PetscHeaderDestroy(sp);CHKERRQ(ierr);
1546   PetscFunctionReturn(0);
1547 }
1548 
1549 #undef __FUNCT__
1550 #define __FUNCT__ "PetscDualSpaceCreate"
1551 /*@
1552   PetscDualSpaceCreate - Creates an empty PetscDualSpace object. The type can then be set with PetscDualSpaceSetType().
1553 
1554   Collective on MPI_Comm
1555 
1556   Input Parameter:
1557 . comm - The communicator for the PetscDualSpace object
1558 
1559   Output Parameter:
1560 . sp - The PetscDualSpace object
1561 
1562   Level: beginner
1563 
1564 .seealso: PetscDualSpaceSetType(), PETSCDUALSPACELAGRANGE
1565 @*/
1566 PetscErrorCode PetscDualSpaceCreate(MPI_Comm comm, PetscDualSpace *sp)
1567 {
1568   PetscDualSpace s;
1569   PetscErrorCode ierr;
1570 
1571   PetscFunctionBegin;
1572   PetscValidPointer(sp, 2);
1573   *sp = NULL;
1574 #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
1575   ierr = PetscFEInitializePackage();CHKERRQ(ierr);
1576 #endif
1577 
1578   ierr = PetscHeaderCreate(s, _p_PetscDualSpace, struct _PetscDualSpaceOps, PETSCDUALSPACE_CLASSID, "PetscDualSpace", "Dual Space", "PetscDualSpace", comm, PetscDualSpaceDestroy, PetscDualSpaceView);CHKERRQ(ierr);
1579   ierr = PetscMemzero(s->ops, sizeof(struct _PetscDualSpaceOps));CHKERRQ(ierr);
1580 
1581   s->order = 0;
1582 
1583   *sp = s;
1584   PetscFunctionReturn(0);
1585 }
1586 
1587 #undef __FUNCT__
1588 #define __FUNCT__ "PetscDualSpaceGetDM"
1589 PetscErrorCode PetscDualSpaceGetDM(PetscDualSpace sp, DM *dm)
1590 {
1591   PetscFunctionBegin;
1592   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1593   PetscValidPointer(dm, 2);
1594   *dm = sp->dm;
1595   PetscFunctionReturn(0);
1596 }
1597 
1598 #undef __FUNCT__
1599 #define __FUNCT__ "PetscDualSpaceSetDM"
1600 PetscErrorCode PetscDualSpaceSetDM(PetscDualSpace sp, DM dm)
1601 {
1602   PetscErrorCode ierr;
1603 
1604   PetscFunctionBegin;
1605   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1606   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
1607   ierr = DMDestroy(&sp->dm);CHKERRQ(ierr);
1608   ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr);
1609   sp->dm = dm;
1610   PetscFunctionReturn(0);
1611 }
1612 
1613 #undef __FUNCT__
1614 #define __FUNCT__ "PetscDualSpaceGetOrder"
1615 PetscErrorCode PetscDualSpaceGetOrder(PetscDualSpace sp, PetscInt *order)
1616 {
1617   PetscFunctionBegin;
1618   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1619   PetscValidPointer(order, 2);
1620   *order = sp->order;
1621   PetscFunctionReturn(0);
1622 }
1623 
1624 #undef __FUNCT__
1625 #define __FUNCT__ "PetscDualSpaceSetOrder"
1626 PetscErrorCode PetscDualSpaceSetOrder(PetscDualSpace sp, PetscInt order)
1627 {
1628   PetscFunctionBegin;
1629   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1630   sp->order = order;
1631   PetscFunctionReturn(0);
1632 }
1633 
1634 #undef __FUNCT__
1635 #define __FUNCT__ "PetscDualSpaceGetFunctional"
1636 PetscErrorCode PetscDualSpaceGetFunctional(PetscDualSpace sp, PetscInt i, PetscQuadrature *functional)
1637 {
1638   PetscInt       dim;
1639   PetscErrorCode ierr;
1640 
1641   PetscFunctionBegin;
1642   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1643   PetscValidPointer(functional, 3);
1644   ierr = PetscDualSpaceGetDimension(sp, &dim);CHKERRQ(ierr);
1645   if ((i < 0) || (i >= dim)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Functional index %d must be in [0, %d)", i, dim);
1646   *functional = sp->functional[i];
1647   PetscFunctionReturn(0);
1648 }
1649 
1650 #undef __FUNCT__
1651 #define __FUNCT__ "PetscDualSpaceGetDimension"
1652 /* Dimension of the space, i.e. number of basis vectors */
1653 PetscErrorCode PetscDualSpaceGetDimension(PetscDualSpace sp, PetscInt *dim)
1654 {
1655   PetscErrorCode ierr;
1656 
1657   PetscFunctionBegin;
1658   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1659   PetscValidPointer(dim, 2);
1660   *dim = 0;
1661   if (sp->ops->getdimension) {ierr = (*sp->ops->getdimension)(sp, dim);CHKERRQ(ierr);}
1662   PetscFunctionReturn(0);
1663 }
1664 
1665 #undef __FUNCT__
1666 #define __FUNCT__ "PetscDualSpaceSetUp_Lagrange"
1667 PetscErrorCode PetscDualSpaceSetUp_Lagrange(PetscDualSpace sp)
1668 {
1669   DM                 dm    = sp->dm;
1670   PetscInt           order = sp->order;
1671   PetscSection       csection;
1672   Vec                coordinates;
1673   const PetscScalar *coords;
1674   PetscReal         *qpoints, *qweights;
1675   PetscInt           depth, dim, pdim, *pStart, *pEnd, coneSize, d, n, f = 0;
1676   PetscErrorCode     ierr;
1677 
1678   PetscFunctionBegin;
1679   ierr = PetscDualSpaceGetDimension(sp, &pdim);CHKERRQ(ierr);
1680   ierr = PetscMalloc(pdim * sizeof(PetscQuadrature), &sp->functional);CHKERRQ(ierr);
1681   /* Classify element type */
1682   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
1683   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
1684   ierr = PetscMalloc2(depth+1,PetscInt,&pStart,depth+1,PetscInt,&pEnd);CHKERRQ(ierr);
1685   for (d = 0; d < depth; ++d) {
1686     ierr = DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);CHKERRQ(ierr);
1687   }
1688   ierr = DMPlexGetConeSize(dm, pStart[depth], &coneSize);CHKERRQ(ierr);
1689   ierr = DMPlexGetCoordinateSection(dm, &csection);CHKERRQ(ierr);
1690   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
1691   if (coneSize == dim+1) {
1692     PetscInt *closure = NULL, closureSize, c;
1693 
1694     /* Simplex */
1695     ierr = DMPlexGetTransitiveClosure(dm, pStart[depth], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1696     for (c = 0; c < closureSize*2; c += 2) {
1697       const PetscInt p = closure[c];
1698 
1699       if ((p >= pStart[0]) && (p < pEnd[0])) {
1700         /* Vertices */
1701         PetscInt dof, off, d;
1702 
1703         if (order < 1) continue;
1704         sp->functional[f].numQuadPoints = 1;
1705         ierr = PetscMalloc(sp->functional[f].numQuadPoints*dim * sizeof(PetscReal), &qpoints);CHKERRQ(ierr);
1706         ierr = PetscMalloc(sp->functional[f].numQuadPoints     * sizeof(PetscReal), &qweights);CHKERRQ(ierr);
1707         ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr);
1708         ierr = PetscSectionGetDof(csection, p, &dof);CHKERRQ(ierr);
1709         ierr = PetscSectionGetOffset(csection, p, &off);CHKERRQ(ierr);
1710         if (dof != dim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of coordinates %d does not match spatial dimension %d", dof, dim);
1711         for (d = 0; d < dof; ++d) {qpoints[d] = coords[off+d];}
1712         qweights[0] = 1.0;
1713         sp->functional[f].quadPoints  = qpoints;
1714         sp->functional[f].quadWeights = qweights;
1715         ++f;
1716         ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr);
1717       } else if ((p >= pStart[1]) && (p < pEnd[1])) {
1718         /* Edges */
1719         PetscInt k;
1720 
1721         if (order < 2) continue;
1722         coords = NULL;
1723         ierr = DMPlexVecGetClosure(dm, csection, coordinates, p, &n, &coords);CHKERRQ(ierr);
1724         if (n != dim*2) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Point %d has %d coordinate values instead of %d", p, n, dim*2);
1725         for (k = 1; k < order; ++k) {
1726           sp->functional[f].numQuadPoints = 1;
1727           ierr = PetscMalloc(sp->functional[f].numQuadPoints*dim * sizeof(PetscReal), &qpoints);CHKERRQ(ierr);
1728           ierr = PetscMalloc(sp->functional[f].numQuadPoints     * sizeof(PetscReal), &qweights);CHKERRQ(ierr);
1729           for (d = 0; d < dim; ++d) {qpoints[d] = k*(coords[1*dim+d] - coords[0*dim+d])/order + coords[0*dim+d];}
1730           qweights[0] = 1.0;
1731           sp->functional[f].quadPoints  = qpoints;
1732           sp->functional[f].quadWeights = qweights;
1733           ++f;
1734         }
1735         ierr = DMPlexVecRestoreClosure(dm, csection, coordinates, p, &n, &coords);CHKERRQ(ierr);
1736       } else if ((p >= pStart[depth-1]) && (p < pEnd[depth-1])) {
1737         /* Faces */
1738 
1739         if (order < 3) continue;
1740         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to implement faces");
1741       } else if ((p >= pStart[depth]) && (p < pEnd[depth])) {
1742         /* Cells */
1743 
1744         if ((order > 0) && (order < 3)) continue;
1745         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to implement cells");
1746       }
1747     }
1748     ierr = DMPlexRestoreTransitiveClosure(dm, pStart[depth], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1749   } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot handle cells with cone size %d", coneSize);
1750   ierr = PetscFree2(pStart,pEnd);CHKERRQ(ierr);
1751   if (f != pdim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of dual basis vector %d not equal to dimension %d", f, pdim);
1752   PetscFunctionReturn(0);
1753 }
1754 
1755 #undef __FUNCT__
1756 #define __FUNCT__ "PetscDualSpaceGetDimension_Lagrange"
1757 PetscErrorCode PetscDualSpaceGetDimension_Lagrange(PetscDualSpace sp, PetscInt *dim)
1758 {
1759   PetscInt            deg = sp->order;
1760   PetscReal           D   = 1.0;
1761   PetscInt            n, i;
1762   PetscErrorCode      ierr;
1763 
1764   PetscFunctionBegin;
1765   /* TODO: Assumes simplices */
1766   ierr = DMPlexGetDimension(sp->dm, &n);CHKERRQ(ierr);
1767   for (i = 1; i <= n; ++i) {
1768     D *= ((PetscReal) (deg+i))/i;
1769   }
1770   *dim = (PetscInt) (D + 0.5);
1771   PetscFunctionReturn(0);
1772 }
1773 
1774 #undef __FUNCT__
1775 #define __FUNCT__ "PetscDualSpaceInitialize_Lagrange"
1776 PetscErrorCode PetscDualSpaceInitialize_Lagrange(PetscDualSpace sp)
1777 {
1778   PetscFunctionBegin;
1779   sp->ops->setfromoptions = NULL;
1780   sp->ops->setup          = PetscDualSpaceSetUp_Lagrange;
1781   sp->ops->view           = NULL;
1782   sp->ops->destroy        = NULL;
1783   sp->ops->getdimension   = PetscDualSpaceGetDimension_Lagrange;
1784   PetscFunctionReturn(0);
1785 }
1786 
1787 /*MC
1788   PETSCDUALSPACELAGRANGE = "lagrange" - A PetscDualSpace object that encapsulates a dual space of pointwise evaluation functionals
1789 
1790   Level: intermediate
1791 
1792 .seealso: PetscDualSpaceType, PetscDualSpaceCreate(), PetscDualSpaceSetType()
1793 M*/
1794 
1795 #undef __FUNCT__
1796 #define __FUNCT__ "PetscDualSpaceCreate_Lagrange"
1797 PETSC_EXTERN PetscErrorCode PetscDualSpaceCreate_Lagrange(PetscDualSpace sp)
1798 {
1799   PetscDualSpace_Lag *lag;
1800   PetscErrorCode      ierr;
1801 
1802   PetscFunctionBegin;
1803   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1804   ierr     = PetscNewLog(sp, PetscDualSpace_Lag, &lag);CHKERRQ(ierr);
1805   sp->data = lag;
1806 
1807   /* lag->n = 0; */
1808 
1809   ierr = PetscDualSpaceInitialize_Lagrange(sp);CHKERRQ(ierr);
1810   PetscFunctionReturn(0);
1811 }
1812 
1813 
1814 PetscInt PETSCFE_CLASSID = 0;
1815 
1816 #undef __FUNCT__
1817 #define __FUNCT__ "PetscFEView"
1818 /*@C
1819   PetscFEView - Views a PetscFE
1820 
1821   Collective on PetscFE
1822 
1823   Input Parameter:
1824 + sp - the PetscFE object to view
1825 - v  - the viewer
1826 
1827   Level: developer
1828 
1829 .seealso PetscFEDestroy()
1830 @*/
1831 PetscErrorCode PetscFEView(PetscFE fem, PetscViewer v)
1832 {
1833   PetscErrorCode ierr;
1834 
1835   PetscFunctionBegin;
1836   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
1837   if (!v) {
1838     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) fem), &v);CHKERRQ(ierr);
1839   }
1840   if (fem->ops->view) {
1841     ierr = (*fem->ops->view)(fem, v);CHKERRQ(ierr);
1842   }
1843   PetscFunctionReturn(0);
1844 }
1845 
1846 #undef __FUNCT__
1847 #define __FUNCT__ "PetscFEDestroy"
1848 /*@
1849   PetscFEDestroy - Destroys a PetscFE object
1850 
1851   Collective on PetscFE
1852 
1853   Input Parameter:
1854 . fem - the PetscFE object to destroy
1855 
1856   Level: developer
1857 
1858 .seealso PetscFEView()
1859 @*/
1860 PetscErrorCode PetscFEDestroy(PetscFE *fem)
1861 {
1862   PetscErrorCode ierr;
1863 
1864   PetscFunctionBegin;
1865   if (!*fem) PetscFunctionReturn(0);
1866   PetscValidHeaderSpecific((*fem), PETSCFE_CLASSID, 1);
1867 
1868   if (--((PetscObject)(*fem))->refct > 0) {*fem = 0; PetscFunctionReturn(0);}
1869   ((PetscObject) (*fem))->refct = 0;
1870   /* if memory was published with AMS then destroy it */
1871   ierr = PetscObjectAMSViewOff((PetscObject) *fem);CHKERRQ(ierr);
1872 
1873   ierr = (*(*fem)->ops->destroy)(*fem);CHKERRQ(ierr);
1874   ierr = PetscHeaderDestroy(fem);CHKERRQ(ierr);
1875   PetscFunctionReturn(0);
1876 }
1877 
1878 #undef __FUNCT__
1879 #define __FUNCT__ "PetscFECreate"
1880 /*@
1881   PetscFECreate - Creates an empty PetscFE object. The type can then be set with PetscFESetType().
1882 
1883   Collective on MPI_Comm
1884 
1885   Input Parameter:
1886 . comm - The communicator for the PetscFE object
1887 
1888   Output Parameter:
1889 . fem - The PetscFE object
1890 
1891   Level: beginner
1892 
1893 .seealso: PetscFESetType(), PETSCFEGALERKIN
1894 @*/
1895 PetscErrorCode PetscFECreate(MPI_Comm comm, PetscFE *fem)
1896 {
1897   PetscFE        f;
1898   PetscErrorCode ierr;
1899 
1900   PetscFunctionBegin;
1901   PetscValidPointer(fem, 2);
1902   *fem = NULL;
1903 #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
1904   ierr = PetscFEInitializePackage();CHKERRQ(ierr);
1905 #endif
1906 
1907   ierr = PetscHeaderCreate(f, _p_PetscFE, struct _PetscFEOps, PETSCFE_CLASSID, "PetscFE", "Finite Element", "PetscFE", comm, PetscFEDestroy, PetscFEView);CHKERRQ(ierr);
1908   ierr = PetscMemzero(f->ops, sizeof(struct _PetscFEOps));CHKERRQ(ierr);
1909 
1910   f->basisSpace = NULL;
1911   f->dualSpace  = NULL;
1912 
1913   *fem = f;
1914   PetscFunctionReturn(0);
1915 }
1916 
1917 #undef __FUNCT__
1918 #define __FUNCT__ "PetscFEGetDimension"
1919 PetscErrorCode PetscFEGetDimension(PetscFE fem, PetscInt *dim)
1920 {
1921   PetscErrorCode ierr;
1922 
1923   PetscFunctionBegin;
1924   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
1925   PetscValidPointer(dim, 2);
1926   ierr = PetscSpaceGetDimension(fem->basisSpace, dim);CHKERRQ(ierr);
1927   PetscFunctionReturn(0);
1928 }
1929 
1930 #undef __FUNCT__
1931 #define __FUNCT__ "PetscFEGetNumComponents"
1932 PetscErrorCode PetscFEGetNumComponents(PetscFE fem, PetscInt *comp)
1933 {
1934   PetscFunctionBegin;
1935   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
1936   PetscValidPointer(comp, 2);
1937   *comp = 1;
1938   PetscFunctionReturn(0);
1939 }
1940 
1941 #undef __FUNCT__
1942 #define __FUNCT__ "PetscFEGetBasisSpace"
1943 PetscErrorCode PetscFEGetBasisSpace(PetscFE fem, PetscSpace *sp)
1944 {
1945   PetscFunctionBegin;
1946   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
1947   PetscValidPointer(sp, 2);
1948   *sp = fem->basisSpace;
1949   PetscFunctionReturn(0);
1950 }
1951 
1952 #undef __FUNCT__
1953 #define __FUNCT__ "PetscFESetBasisSpace"
1954 PetscErrorCode PetscFESetBasisSpace(PetscFE fem, PetscSpace sp)
1955 {
1956   PetscFunctionBegin;
1957   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
1958   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 2);
1959   fem->basisSpace = sp;
1960   PetscFunctionReturn(0);
1961 }
1962 
1963 #undef __FUNCT__
1964 #define __FUNCT__ "PetscFEGetDualSpace"
1965 PetscErrorCode PetscFEGetDualSpace(PetscFE fem, PetscDualSpace *sp)
1966 {
1967   PetscFunctionBegin;
1968   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
1969   PetscValidPointer(sp, 2);
1970   *sp = fem->dualSpace;
1971   PetscFunctionReturn(0);
1972 }
1973 
1974 #undef __FUNCT__
1975 #define __FUNCT__ "PetscFESetDualSpace"
1976 PetscErrorCode PetscFESetDualSpace(PetscFE fem, PetscDualSpace sp)
1977 {
1978   PetscFunctionBegin;
1979   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
1980   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 2);
1981   fem->dualSpace = sp;
1982   PetscFunctionReturn(0);
1983 }
1984 
1985 #undef __FUNCT__
1986 #define __FUNCT__ "PetscFEGetTabulation"
1987 PetscErrorCode PetscFEGetTabulation(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscReal **B, PetscReal **D, PetscReal **H)
1988 {
1989   DM             dm;
1990   PetscReal     *tmpB, *invV;
1991   PetscInt       pdim; /* Dimension of FE space P */
1992   PetscInt       dim;  /* Spatial dimension */
1993   PetscInt       p, j, k;
1994   PetscErrorCode ierr;
1995 
1996   PetscFunctionBegin;
1997   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
1998   PetscValidPointer(points, 3);
1999   if (B) PetscValidPointer(B, 4);
2000   if (D) PetscValidPointer(D, 5);
2001   if (H) PetscValidPointer(H, 6);
2002   ierr = PetscDualSpaceGetDM(fem->dualSpace, &dm);CHKERRQ(ierr);
2003 
2004   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
2005   ierr = PetscSpaceGetDimension(fem->basisSpace, &pdim);CHKERRQ(ierr);
2006   /* if (nvalues%dim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The number of coordinate values %d must be divisible by the spatial dimension %d", nvalues, dim); */
2007 
2008   if (B) {
2009     ierr = DMGetWorkArray(dm, npoints*pdim, PETSC_REAL, B);CHKERRQ(ierr);
2010     ierr = DMGetWorkArray(dm, npoints*pdim, PETSC_REAL, &tmpB);CHKERRQ(ierr);
2011   }
2012   if (D) {ierr = DMGetWorkArray(dm, npoints*pdim*dim, PETSC_REAL, D);CHKERRQ(ierr);}
2013   if (H) {ierr = DMGetWorkArray(dm, npoints*pdim*dim*dim, PETSC_REAL, H);CHKERRQ(ierr);}
2014   ierr = PetscSpaceEvaluate(fem->basisSpace, npoints, points, B ? tmpB : NULL, D ? *D : NULL, H ? *H : NULL);CHKERRQ(ierr);
2015 
2016   ierr = DMGetWorkArray(dm, pdim*pdim, PETSC_REAL, &invV);CHKERRQ(ierr);
2017   for (j = 0; j < pdim; ++j) {
2018     PetscReal      *Bf;
2019     PetscQuadrature f;
2020     PetscInt        q;
2021 
2022     ierr = PetscDualSpaceGetFunctional(fem->dualSpace, j, &f);
2023     ierr = DMGetWorkArray(dm, f.numQuadPoints*pdim, PETSC_REAL, &Bf);CHKERRQ(ierr);
2024     ierr = PetscSpaceEvaluate(fem->basisSpace, f.numQuadPoints, f.quadPoints, Bf, NULL, NULL);CHKERRQ(ierr);
2025     for (k = 0; k < pdim; ++k) {
2026       /* n_j \cdot \phi_k */
2027       invV[j*pdim+k] = 0.0;
2028       for (q = 0; q < f.numQuadPoints; ++q) {
2029         invV[j*pdim+k] += Bf[q*pdim+k]*f.quadWeights[q];
2030       }
2031     }
2032     ierr = DMRestoreWorkArray(dm, f.numQuadPoints*pdim, PETSC_REAL, &Bf);CHKERRQ(ierr);
2033   }
2034   {
2035     PetscReal    *work;
2036     PetscBLASInt *pivots;
2037     PetscBLASInt  n = pdim, info;
2038 
2039     ierr = DMGetWorkArray(dm, pdim, PETSC_INT, &pivots);CHKERRQ(ierr);
2040     ierr = DMGetWorkArray(dm, pdim, PETSC_REAL, &work);CHKERRQ(ierr);
2041     PetscStackCallBLAS("LAPACKgetrf", LAPACKgetrf_(&n, &n, invV, &n, pivots, &info));
2042     PetscStackCallBLAS("LAPACKgetri", LAPACKgetri_(&n, invV, &n, pivots, work, &n, &info));
2043     ierr = DMRestoreWorkArray(dm, pdim, PETSC_INT, &pivots);CHKERRQ(ierr);
2044     ierr = DMRestoreWorkArray(dm, pdim, PETSC_REAL, &work);CHKERRQ(ierr);
2045   }
2046   for (p = 0; p < npoints; ++p) {
2047     if (B) {
2048       /* Multiply by V^{-1} (pdim x pdim) */
2049       for (j = 0; j < pdim; ++j) {
2050         const PetscInt i = p*pdim + j;
2051 
2052         (*B)[i] = 0.0;
2053         for (k = 0; k < pdim; ++k) {
2054           (*B)[i] += invV[k*pdim+j] * tmpB[p*pdim + k];
2055         }
2056       }
2057     }
2058   }
2059   ierr = DMRestoreWorkArray(dm, pdim*pdim, PETSC_REAL, &invV);CHKERRQ(ierr);
2060   if (B) {
2061     ierr = DMRestoreWorkArray(dm, npoints*pdim, PETSC_REAL, &tmpB);CHKERRQ(ierr);
2062   }
2063   PetscFunctionReturn(0);
2064 }
2065 
2066 #undef __FUNCT__
2067 #define __FUNCT__ "PetscFERestoreTabulation"
2068 PetscErrorCode PetscFERestoreTabulation(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscReal **B, PetscReal **D, PetscReal **D2)
2069 {
2070   DM             dm;
2071   PetscErrorCode ierr;
2072 
2073   PetscFunctionBegin;
2074   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
2075   ierr = PetscDualSpaceGetDM(fem->dualSpace, &dm);CHKERRQ(ierr);
2076   if (B)  {ierr = DMRestoreWorkArray(dm, 0, PETSC_REAL, B);CHKERRQ(ierr);}
2077   if (D)  {ierr = DMRestoreWorkArray(dm, 0, PETSC_REAL, D);CHKERRQ(ierr);}
2078   if (D2) {ierr = DMRestoreWorkArray(dm, 0, PETSC_REAL, D2);CHKERRQ(ierr);}
2079   PetscFunctionReturn(0);
2080 }
2081