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