xref: /petsc/src/dm/dt/interface/dt.c (revision 46e7e8de0ebd9938610bcaa9bbaa8e3209b7e39b)
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__ "PetscQuadratureDestroy"
17 PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *q)
18 {
19   PetscErrorCode ierr;
20 
21   PetscFunctionBegin;
22   ierr = PetscFree(q->quadPoints);CHKERRQ(ierr);
23   ierr = PetscFree(q->quadWeights);CHKERRQ(ierr);
24   PetscFunctionReturn(0);
25 }
26 
27 #undef __FUNCT__
28 #define __FUNCT__ "PetscDTLegendreEval"
29 /*@
30    PetscDTLegendreEval - evaluate Legendre polynomial at points
31 
32    Not Collective
33 
34    Input Arguments:
35 +  npoints - number of spatial points to evaluate at
36 .  points - array of locations to evaluate at
37 .  ndegree - number of basis degrees to evaluate
38 -  degrees - sorted array of degrees to evaluate
39 
40    Output Arguments:
41 +  B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL)
42 .  D - row-oriented derivative evaluation matrix (or NULL)
43 -  D2 - row-oriented second derivative evaluation matrix (or NULL)
44 
45    Level: intermediate
46 
47 .seealso: PetscDTGaussQuadrature()
48 @*/
49 PetscErrorCode PetscDTLegendreEval(PetscInt npoints,const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2)
50 {
51   PetscInt i,maxdegree;
52 
53   PetscFunctionBegin;
54   if (!npoints || !ndegree) PetscFunctionReturn(0);
55   maxdegree = degrees[ndegree-1];
56   for (i=0; i<npoints; i++) {
57     PetscReal pm1,pm2,pd1,pd2,pdd1,pdd2,x;
58     PetscInt  j,k;
59     x    = points[i];
60     pm2  = 0;
61     pm1  = 1;
62     pd2  = 0;
63     pd1  = 0;
64     pdd2 = 0;
65     pdd1 = 0;
66     k    = 0;
67     if (degrees[k] == 0) {
68       if (B) B[i*ndegree+k] = pm1;
69       if (D) D[i*ndegree+k] = pd1;
70       if (D2) D2[i*ndegree+k] = pdd1;
71       k++;
72     }
73     for (j=1; j<=maxdegree; j++,k++) {
74       PetscReal p,d,dd;
75       p    = ((2*j-1)*x*pm1 - (j-1)*pm2)/j;
76       d    = pd2 + (2*j-1)*pm1;
77       dd   = pdd2 + (2*j-1)*pd1;
78       pm2  = pm1;
79       pm1  = p;
80       pd2  = pd1;
81       pd1  = d;
82       pdd2 = pdd1;
83       pdd1 = dd;
84       if (degrees[k] == j) {
85         if (B) B[i*ndegree+k] = p;
86         if (D) D[i*ndegree+k] = d;
87         if (D2) D2[i*ndegree+k] = dd;
88       }
89     }
90   }
91   PetscFunctionReturn(0);
92 }
93 
94 #undef __FUNCT__
95 #define __FUNCT__ "PetscDTGaussQuadrature"
96 /*@
97    PetscDTGaussQuadrature - create Gauss quadrature
98 
99    Not Collective
100 
101    Input Arguments:
102 +  npoints - number of points
103 .  a - left end of interval (often-1)
104 -  b - right end of interval (often +1)
105 
106    Output Arguments:
107 +  x - quadrature points
108 -  w - quadrature weights
109 
110    Level: intermediate
111 
112    References:
113    Golub and Welsch, Calculation of Quadrature Rules, Math. Comp. 23(106), 221--230, 1969.
114 
115 .seealso: PetscDTLegendreEval()
116 @*/
117 PetscErrorCode PetscDTGaussQuadrature(PetscInt npoints,PetscReal a,PetscReal b,PetscReal *x,PetscReal *w)
118 {
119   PetscErrorCode ierr;
120   PetscInt       i;
121   PetscReal      *work;
122   PetscScalar    *Z;
123   PetscBLASInt   N,LDZ,info;
124 
125   PetscFunctionBegin;
126   /* Set up the Golub-Welsch system */
127   for (i=0; i<npoints; i++) {
128     x[i] = 0;                   /* diagonal is 0 */
129     if (i) w[i-1] = 0.5 / PetscSqrtReal(1 - 1./PetscSqr(2*i));
130   }
131   ierr = PetscRealView(npoints-1,w,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
132   ierr = PetscMalloc2(npoints*npoints,PetscScalar,&Z,PetscMax(1,2*npoints-2),PetscReal,&work);CHKERRQ(ierr);
133   ierr = PetscBLASIntCast(npoints,&N);CHKERRQ(ierr);
134   LDZ  = N;
135   ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
136   PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("I",&N,x,w,Z,&LDZ,work,&info));
137   ierr = PetscFPTrapPop();CHKERRQ(ierr);
138   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEQR error");
139 
140   for (i=0; i<(npoints+1)/2; i++) {
141     PetscReal y = 0.5 * (-x[i] + x[npoints-i-1]); /* enforces symmetry */
142     x[i]           = (a+b)/2 - y*(b-a)/2;
143     x[npoints-i-1] = (a+b)/2 + y*(b-a)/2;
144 
145     w[i] = w[npoints-1-i] = (b-a)*PetscSqr(0.5*PetscAbsScalar(Z[i*npoints] + Z[(npoints-i-1)*npoints]));
146   }
147   ierr = PetscFree2(Z,work);CHKERRQ(ierr);
148   PetscFunctionReturn(0);
149 }
150 
151 #undef __FUNCT__
152 #define __FUNCT__ "PetscDTFactorial_Internal"
153 /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x.
154    Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */
155 PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorial_Internal(PetscInt n, PetscReal *factorial)
156 {
157   PetscReal f = 1.0;
158   PetscInt  i;
159 
160   PetscFunctionBegin;
161   for (i = 1; i < n+1; ++i) f *= i;
162   *factorial = f;
163   PetscFunctionReturn(0);
164 }
165 
166 #undef __FUNCT__
167 #define __FUNCT__ "PetscDTComputeJacobi"
168 /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x.
169    Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */
170 PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobi(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P)
171 {
172   PetscReal apb, pn1, pn2;
173   PetscInt  k;
174 
175   PetscFunctionBegin;
176   if (!n) {*P = 1.0; PetscFunctionReturn(0);}
177   if (n == 1) {*P = 0.5 * (a - b + (a + b + 2.0) * x); PetscFunctionReturn(0);}
178   apb = a + b;
179   pn2 = 1.0;
180   pn1 = 0.5 * (a - b + (apb + 2.0) * x);
181   *P  = 0.0;
182   for (k = 2; k < n+1; ++k) {
183     PetscReal a1 = 2.0 * k * (k + apb) * (2.0*k + apb - 2.0);
184     PetscReal a2 = (2.0 * k + apb - 1.0) * (a*a - b*b);
185     PetscReal a3 = (2.0 * k + apb - 2.0) * (2.0 * k + apb - 1.0) * (2.0 * k + apb);
186     PetscReal a4 = 2.0 * (k + a - 1.0) * (k + b - 1.0) * (2.0 * k + apb);
187 
188     a2  = a2 / a1;
189     a3  = a3 / a1;
190     a4  = a4 / a1;
191     *P  = (a2 + a3 * x) * pn1 - a4 * pn2;
192     pn2 = pn1;
193     pn1 = *P;
194   }
195   PetscFunctionReturn(0);
196 }
197 
198 #undef __FUNCT__
199 #define __FUNCT__ "PetscDTComputeJacobiDerivative"
200 /* Evaluates the first derivative of P_{n}^{a,b} at a point x. */
201 PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobiDerivative(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P)
202 {
203   PetscReal      nP;
204   PetscErrorCode ierr;
205 
206   PetscFunctionBegin;
207   if (!n) {*P = 0.0; PetscFunctionReturn(0);}
208   ierr = PetscDTComputeJacobi(a+1, b+1, n-1, x, &nP);CHKERRQ(ierr);
209   *P   = 0.5 * (a + b + n + 1) * nP;
210   PetscFunctionReturn(0);
211 }
212 
213 #undef __FUNCT__
214 #define __FUNCT__ "PetscDTMapSquareToTriangle_Internal"
215 /* Maps from [-1,1]^2 to the (-1,1) reference triangle */
216 PETSC_STATIC_INLINE PetscErrorCode PetscDTMapSquareToTriangle_Internal(PetscReal x, PetscReal y, PetscReal *xi, PetscReal *eta)
217 {
218   PetscFunctionBegin;
219   *xi  = 0.5 * (1.0 + x) * (1.0 - y) - 1.0;
220   *eta = y;
221   PetscFunctionReturn(0);
222 }
223 
224 #undef __FUNCT__
225 #define __FUNCT__ "PetscDTMapCubeToTetrahedron_Internal"
226 /* Maps from [-1,1]^2 to the (-1,1) reference triangle */
227 PETSC_STATIC_INLINE PetscErrorCode PetscDTMapCubeToTetrahedron_Internal(PetscReal x, PetscReal y, PetscReal z, PetscReal *xi, PetscReal *eta, PetscReal *zeta)
228 {
229   PetscFunctionBegin;
230   *xi   = 0.25 * (1.0 + x) * (1.0 - y) * (1.0 - z) - 1.0;
231   *eta  = 0.5  * (1.0 + y) * (1.0 - z) - 1.0;
232   *zeta = z;
233   PetscFunctionReturn(0);
234 }
235 
236 #undef __FUNCT__
237 #define __FUNCT__ "PetscDTGaussJacobiQuadrature1D_Internal"
238 static PetscErrorCode PetscDTGaussJacobiQuadrature1D_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal *x, PetscReal *w)
239 {
240   PetscInt       maxIter = 100;
241   PetscReal      eps     = 1.0e-8;
242   PetscReal      a1, a2, a3, a4, a5, a6;
243   PetscInt       k;
244   PetscErrorCode ierr;
245 
246   PetscFunctionBegin;
247 
248   a1      = pow(2, a+b+1);
249 #if defined(PETSC_HAVE_TGAMMA)
250   a2      = tgamma(a + npoints + 1);
251   a3      = tgamma(b + npoints + 1);
252   a4      = tgamma(a + b + npoints + 1);
253 #else
254   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"tgamma() - math routine is unavailable.");
255 #endif
256 
257   ierr = PetscDTFactorial_Internal(npoints, &a5);CHKERRQ(ierr);
258   a6   = a1 * a2 * a3 / a4 / a5;
259   /* Computes the m roots of P_{m}^{a,b} on [-1,1] by Newton's method with Chebyshev points as initial guesses.
260    Algorithm implemented from the pseudocode given by Karniadakis and Sherwin and Python in FIAT */
261   for (k = 0; k < npoints; ++k) {
262     PetscReal r = -cos((2.0*k + 1.0) * PETSC_PI / (2.0 * npoints)), dP;
263     PetscInt  j;
264 
265     if (k > 0) r = 0.5 * (r + x[k-1]);
266     for (j = 0; j < maxIter; ++j) {
267       PetscReal s = 0.0, delta, f, fp;
268       PetscInt  i;
269 
270       for (i = 0; i < k; ++i) s = s + 1.0 / (r - x[i]);
271       ierr = PetscDTComputeJacobi(a, b, npoints, r, &f);CHKERRQ(ierr);
272       ierr = PetscDTComputeJacobiDerivative(a, b, npoints, r, &fp);CHKERRQ(ierr);
273       delta = f / (fp - f * s);
274       r     = r - delta;
275       if (fabs(delta) < eps) break;
276     }
277     x[k] = r;
278     ierr = PetscDTComputeJacobiDerivative(a, b, npoints, x[k], &dP);CHKERRQ(ierr);
279     w[k] = a6 / (1.0 - PetscSqr(x[k])) / PetscSqr(dP);
280   }
281   PetscFunctionReturn(0);
282 }
283 
284 #undef __FUNCT__
285 #define __FUNCT__ "PetscDTGaussJacobiQuadrature"
286 /*@C
287   PetscDTGaussJacobiQuadrature - create Gauss-Jacobi quadrature for a simplex
288 
289   Not Collective
290 
291   Input Arguments:
292 + dim - The simplex dimension
293 . npoints - number of points
294 . a - left end of interval (often-1)
295 - b - right end of interval (often +1)
296 
297   Output Arguments:
298 + points - quadrature points
299 - weights - quadrature weights
300 
301   Level: intermediate
302 
303   References:
304   Karniadakis and Sherwin.
305   FIAT
306 
307 .seealso: PetscDTGaussQuadrature()
308 @*/
309 PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt dim, PetscInt npoints, PetscReal a, PetscReal b, PetscReal *points[], PetscReal *weights[])
310 {
311   PetscReal     *px, *wx, *py, *wy, *pz, *wz, *x, *w;
312   PetscInt       i, j, k;
313   PetscErrorCode ierr;
314 
315   PetscFunctionBegin;
316   if ((a != -1.0) || (b != 1.0)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must use default internal right now");
317   switch (dim) {
318   case 1:
319     ierr = PetscMalloc(npoints * sizeof(PetscReal), &x);CHKERRQ(ierr);
320     ierr = PetscMalloc(npoints * sizeof(PetscReal), &w);CHKERRQ(ierr);
321     ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, x, w);CHKERRQ(ierr);
322     break;
323   case 2:
324     ierr = PetscMalloc(npoints*npoints*2 * sizeof(PetscReal), &x);CHKERRQ(ierr);
325     ierr = PetscMalloc(npoints*npoints   * sizeof(PetscReal), &w);CHKERRQ(ierr);
326     ierr = PetscMalloc4(npoints,PetscReal,&px,npoints,PetscReal,&wx,npoints,PetscReal,&py,npoints,PetscReal,&wy);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     for (i = 0; i < npoints; ++i) {
330       for (j = 0; j < npoints; ++j) {
331         ierr = PetscDTMapSquareToTriangle_Internal(px[i], py[j], &x[(i*npoints+j)*2+0], &x[(i*npoints+j)*2+1]);CHKERRQ(ierr);
332         w[i*npoints+j] = 0.5 * wx[i] * wy[j];
333       }
334     }
335     ierr = PetscFree4(px,wx,py,wy);CHKERRQ(ierr);
336     break;
337   case 3:
338     ierr = PetscMalloc(npoints*npoints*3 * sizeof(PetscReal), &x);CHKERRQ(ierr);
339     ierr = PetscMalloc(npoints*npoints   * sizeof(PetscReal), &w);CHKERRQ(ierr);
340     ierr = PetscMalloc6(npoints,PetscReal,&px,npoints,PetscReal,&wx,npoints,PetscReal,&py,npoints,PetscReal,&wy,npoints,PetscReal,&pz,npoints,PetscReal,&wz);CHKERRQ(ierr);
341     ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, px, wx);CHKERRQ(ierr);
342     ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 1.0, 0.0, py, wy);CHKERRQ(ierr);
343     ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 2.0, 0.0, pz, wz);CHKERRQ(ierr);
344     for (i = 0; i < npoints; ++i) {
345       for (j = 0; j < npoints; ++j) {
346         for (k = 0; k < npoints; ++k) {
347           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);
348           w[(i*npoints+j)*npoints+k] = 0.125 * wx[i] * wy[j] * wz[k];
349         }
350       }
351     }
352     ierr = PetscFree6(px,wx,py,wy,pz,wz);CHKERRQ(ierr);
353     break;
354   default:
355     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim);
356   }
357   if (points)  *points  = x;
358   if (weights) *weights = w;
359   PetscFunctionReturn(0);
360 }
361 
362 #undef __FUNCT__
363 #define __FUNCT__ "PetscDTPseudoInverseQR"
364 /* Overwrites A. Can only handle full-rank problems with m>=n
365  * A in column-major format
366  * Ainv in row-major format
367  * tau has length m
368  * worksize must be >= max(1,n)
369  */
370 static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m,PetscInt mstride,PetscInt n,PetscReal *A_in,PetscReal *Ainv_out,PetscScalar *tau,PetscInt worksize,PetscScalar *work)
371 {
372   PetscErrorCode ierr;
373   PetscBLASInt M,N,K,lda,ldb,ldwork,info;
374   PetscScalar *A,*Ainv,*R,*Q,Alpha;
375 
376   PetscFunctionBegin;
377 #if defined(PETSC_USE_COMPLEX)
378   {
379     PetscInt i,j;
380     ierr = PetscMalloc2(m*n,PetscScalar,&A,m*n,PetscScalar,&Ainv);CHKERRQ(ierr);
381     for (j=0; j<n; j++) {
382       for (i=0; i<m; i++) A[i+m*j] = A_in[i+mstride*j];
383     }
384     mstride = m;
385   }
386 #else
387   A = A_in;
388   Ainv = Ainv_out;
389 #endif
390 
391   ierr = PetscBLASIntCast(m,&M);CHKERRQ(ierr);
392   ierr = PetscBLASIntCast(n,&N);CHKERRQ(ierr);
393   ierr = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr);
394   ierr = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr);
395   ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
396   LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info);
397   ierr = PetscFPTrapPop();CHKERRQ(ierr);
398   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error");
399   R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */
400 
401   /* Extract an explicit representation of Q */
402   Q = Ainv;
403   ierr = PetscMemcpy(Q,A,mstride*n*sizeof(PetscScalar));CHKERRQ(ierr);
404   K = N;                        /* full rank */
405   LAPACKungqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info);
406   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error");
407 
408   /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */
409   Alpha = 1.0;
410   ldb = lda;
411   BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb);
412   /* Ainv is Q, overwritten with inverse */
413 
414 #if defined(PETSC_USE_COMPLEX)
415   {
416     PetscInt i;
417     for (i=0; i<m*n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]);
418     ierr = PetscFree2(A,Ainv);CHKERRQ(ierr);
419   }
420 #endif
421   PetscFunctionReturn(0);
422 }
423 
424 #undef __FUNCT__
425 #define __FUNCT__ "PetscDTLegendreIntegrate"
426 /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */
427 static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval,const PetscReal *x,PetscInt ndegree,const PetscInt *degrees,PetscBool Transpose,PetscReal *B)
428 {
429   PetscErrorCode ierr;
430   PetscReal *Bv;
431   PetscInt i,j;
432 
433   PetscFunctionBegin;
434   ierr = PetscMalloc((ninterval+1)*ndegree*sizeof(PetscReal),&Bv);CHKERRQ(ierr);
435   /* Point evaluation of L_p on all the source vertices */
436   ierr = PetscDTLegendreEval(ninterval+1,x,ndegree,degrees,Bv,NULL,NULL);CHKERRQ(ierr);
437   /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */
438   for (i=0; i<ninterval; i++) {
439     for (j=0; j<ndegree; j++) {
440       if (Transpose) B[i+ninterval*j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j];
441       else           B[i*ndegree+j]   = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j];
442     }
443   }
444   ierr = PetscFree(Bv);CHKERRQ(ierr);
445   PetscFunctionReturn(0);
446 }
447 
448 #undef __FUNCT__
449 #define __FUNCT__ "PetscDTReconstructPoly"
450 /*@
451    PetscDTReconstructPoly - create matrix representing polynomial reconstruction using cell intervals and evaluation at target intervals
452 
453    Not Collective
454 
455    Input Arguments:
456 +  degree - degree of reconstruction polynomial
457 .  nsource - number of source intervals
458 .  sourcex - sorted coordinates of source cell boundaries (length nsource+1)
459 .  ntarget - number of target intervals
460 -  targetx - sorted coordinates of target cell boundaries (length ntarget+1)
461 
462    Output Arguments:
463 .  R - reconstruction matrix, utarget = sum_s R[t*nsource+s] * usource[s]
464 
465    Level: advanced
466 
467 .seealso: PetscDTLegendreEval()
468 @*/
469 PetscErrorCode PetscDTReconstructPoly(PetscInt degree,PetscInt nsource,const PetscReal *sourcex,PetscInt ntarget,const PetscReal *targetx,PetscReal *R)
470 {
471   PetscErrorCode ierr;
472   PetscInt i,j,k,*bdegrees,worksize;
473   PetscReal xmin,xmax,center,hscale,*sourcey,*targety,*Bsource,*Bsinv,*Btarget;
474   PetscScalar *tau,*work;
475 
476   PetscFunctionBegin;
477   PetscValidRealPointer(sourcex,3);
478   PetscValidRealPointer(targetx,5);
479   PetscValidRealPointer(R,6);
480   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);
481 #if defined(PETSC_USE_DEBUG)
482   for (i=0; i<nsource; i++) {
483     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]);
484   }
485   for (i=0; i<ntarget; i++) {
486     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]);
487   }
488 #endif
489   xmin = PetscMin(sourcex[0],targetx[0]);
490   xmax = PetscMax(sourcex[nsource],targetx[ntarget]);
491   center = (xmin + xmax)/2;
492   hscale = (xmax - xmin)/2;
493   worksize = nsource;
494   ierr = PetscMalloc4(degree+1,PetscInt,&bdegrees,nsource+1,PetscReal,&sourcey,nsource*(degree+1),PetscReal,&Bsource,worksize,PetscScalar,&work);CHKERRQ(ierr);
495   ierr = PetscMalloc4(nsource,PetscScalar,&tau,nsource*(degree+1),PetscReal,&Bsinv,ntarget+1,PetscReal,&targety,ntarget*(degree+1),PetscReal,&Btarget);CHKERRQ(ierr);
496   for (i=0; i<=nsource; i++) sourcey[i] = (sourcex[i]-center)/hscale;
497   for (i=0; i<=degree; i++) bdegrees[i] = i+1;
498   ierr = PetscDTLegendreIntegrate(nsource,sourcey,degree+1,bdegrees,PETSC_TRUE,Bsource);CHKERRQ(ierr);
499   ierr = PetscDTPseudoInverseQR(nsource,nsource,degree+1,Bsource,Bsinv,tau,nsource,work);CHKERRQ(ierr);
500   for (i=0; i<=ntarget; i++) targety[i] = (targetx[i]-center)/hscale;
501   ierr = PetscDTLegendreIntegrate(ntarget,targety,degree+1,bdegrees,PETSC_FALSE,Btarget);CHKERRQ(ierr);
502   for (i=0; i<ntarget; i++) {
503     PetscReal rowsum = 0;
504     for (j=0; j<nsource; j++) {
505       PetscReal sum = 0;
506       for (k=0; k<degree+1; k++) {
507         sum += Btarget[i*(degree+1)+k] * Bsinv[k*nsource+j];
508       }
509       R[i*nsource+j] = sum;
510       rowsum += sum;
511     }
512     for (j=0; j<nsource; j++) R[i*nsource+j] /= rowsum; /* normalize each row */
513   }
514   ierr = PetscFree4(bdegrees,sourcey,Bsource,work);CHKERRQ(ierr);
515   ierr = PetscFree4(tau,Bsinv,targety,Btarget);CHKERRQ(ierr);
516   PetscFunctionReturn(0);
517 }
518 
519 /* Basis Jet Tabulation
520 
521 We would like to tabulate the nodal basis functions and derivatives at a set of points, usually quadrature points. We
522 follow here the derviation in http://www.math.ttu.edu/~kirby/papers/fiat-toms-2004.pdf. The nodal basis $\psi_i$ can
523 be expressed in terms of a prime basis $\phi_i$ which can be stably evaluated. In PETSc, we will use the Legendre basis
524 as a prime basis.
525 
526   \psi_i = \sum_k \alpha_{ki} \phi_k
527 
528 Our nodal basis is defined in terms of the dual basis $n_j$
529 
530   n_j \cdot \psi_i = \delta_{ji}
531 
532 and we may act on the first equation to obtain
533 
534   n_j \cdot \psi_i = \sum_k \alpha_{ki} n_j \cdot \phi_k
535        \delta_{ji} = \sum_k \alpha_{ki} V_{jk}
536                  I = V \alpha
537 
538 so the coefficients of the nodal basis in the prime basis are
539 
540    \alpha = V^{-1}
541 
542 We will define the dual basis vectors $n_j$ using a quadrature rule.
543 
544 Right now, we will just use the polynomial spaces P^k. I know some elements use the space of symmetric polynomials
545 (I think Nedelec), but we will neglect this for now. Constraints in the space, e.g. Arnold-Winther elements, can
546 be implemented exactly as in FIAT using functionals $L_j$.
547 
548 I will have to count the degrees correctly for the Legendre product when we are on simplices.
549 
550 We will have three objects:
551  - Space, P: this just need point evaluation I think
552  - 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
553  - FEM: This keeps {P, P', Q}
554 */
555 #include <petsc-private/petscfeimpl.h>
556 
557 PetscInt PETSCSPACE_CLASSID = 0;
558 
559 PetscFunctionList PetscSpaceList              = NULL;
560 PetscBool         PetscSpaceRegisterAllCalled = PETSC_FALSE;
561 
562 #undef __FUNCT__
563 #define __FUNCT__ "PetscSpaceRegister"
564 /*@C
565   PetscSpaceRegister - Adds a new PetscSpace implementation
566 
567   Not Collective
568 
569   Input Parameters:
570 + name        - The name of a new user-defined creation routine
571 - create_func - The creation routine itself
572 
573   Notes:
574   PetscSpaceRegister() may be called multiple times to add several user-defined PetscSpaces
575 
576   Sample usage:
577 .vb
578     PetscSpaceRegister("my_space", MyPetscSpaceCreate);
579 .ve
580 
581   Then, your PetscSpace type can be chosen with the procedural interface via
582 .vb
583     PetscSpaceCreate(MPI_Comm, PetscSpace *);
584     PetscSpaceSetType(PetscSpace, "my_space");
585 .ve
586    or at runtime via the option
587 .vb
588     -petscspace_type my_space
589 .ve
590 
591   Level: advanced
592 
593 .keywords: PetscSpace, register
594 .seealso: PetscSpaceRegisterAll(), PetscSpaceRegisterDestroy()
595 
596 @*/
597 PetscErrorCode PetscSpaceRegister(const char sname[], PetscErrorCode (*function)(PetscSpace))
598 {
599   PetscErrorCode ierr;
600 
601   PetscFunctionBegin;
602   ierr = PetscFunctionListAdd(&PetscSpaceList, sname, function);CHKERRQ(ierr);
603   PetscFunctionReturn(0);
604 }
605 
606 #undef __FUNCT__
607 #define __FUNCT__ "PetscSpaceSetType"
608 /*@C
609   PetscSpaceSetType - Builds a particular PetscSpace
610 
611   Collective on PetscSpace
612 
613   Input Parameters:
614 + sp   - The PetscSpace object
615 - name - The kind of space
616 
617   Options Database Key:
618 . -petscspace_type <type> - Sets the PetscSpace type; use -help for a list of available types
619 
620   Level: intermediate
621 
622 .keywords: PetscSpace, set, type
623 .seealso: PetscSpaceGetType(), PetscSpaceCreate()
624 @*/
625 PetscErrorCode PetscSpaceSetType(PetscSpace sp, PetscSpaceType name)
626 {
627   PetscErrorCode (*r)(PetscSpace);
628   PetscBool      match;
629   PetscErrorCode ierr;
630 
631   PetscFunctionBegin;
632   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
633   ierr = PetscObjectTypeCompare((PetscObject) sp, name, &match);CHKERRQ(ierr);
634   if (match) PetscFunctionReturn(0);
635 
636   if (!PetscSpaceRegisterAllCalled) {ierr = PetscSpaceRegisterAll();CHKERRQ(ierr);}
637   ierr = PetscFunctionListFind(PetscSpaceList, name, &r);CHKERRQ(ierr);
638   if (!r) SETERRQ1(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscSpace type: %s", name);
639 
640   if (sp->ops->destroy) {
641     ierr             = (*sp->ops->destroy)(sp);CHKERRQ(ierr);
642     sp->ops->destroy = NULL;
643   }
644   ierr = (*r)(sp);CHKERRQ(ierr);
645   ierr = PetscObjectChangeTypeName((PetscObject) sp, name);CHKERRQ(ierr);
646   PetscFunctionReturn(0);
647 }
648 
649 #undef __FUNCT__
650 #define __FUNCT__ "PetscSpaceGetType"
651 /*@C
652   PetscSpaceGetType - Gets the PetscSpace type name (as a string) from the object.
653 
654   Not Collective
655 
656   Input Parameter:
657 . dm  - The PetscSpace
658 
659   Output Parameter:
660 . name - The PetscSpace type name
661 
662   Level: intermediate
663 
664 .keywords: PetscSpace, get, type, name
665 .seealso: PetscSpaceSetType(), PetscSpaceCreate()
666 @*/
667 PetscErrorCode PetscSpaceGetType(PetscSpace sp, PetscSpaceType *name)
668 {
669   PetscErrorCode ierr;
670 
671   PetscFunctionBegin;
672   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
673   PetscValidCharPointer(name, 2);
674   if (!PetscSpaceRegisterAllCalled) {
675     ierr = PetscSpaceRegisterAll();CHKERRQ(ierr);
676   }
677   *name = ((PetscObject) sp)->type_name;
678   PetscFunctionReturn(0);
679 }
680 
681 #undef __FUNCT__
682 #define __FUNCT__ "PetscSpaceView"
683 /*@C
684   PetscSpaceView - Views a PetscSpace
685 
686   Collective on PetscSpace
687 
688   Input Parameter:
689 + sp - the PetscSpace object to view
690 - v  - the viewer
691 
692   Level: developer
693 
694 .seealso PetscSpaceDestroy()
695 @*/
696 PetscErrorCode PetscSpaceView(PetscSpace sp, PetscViewer v)
697 {
698   PetscErrorCode ierr;
699 
700   PetscFunctionBegin;
701   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
702   if (!v) {
703     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) sp), &v);CHKERRQ(ierr);
704   }
705   if (sp->ops->view) {
706     ierr = (*sp->ops->view)(sp, v);CHKERRQ(ierr);
707   }
708   PetscFunctionReturn(0);
709 }
710 
711 #undef __FUNCT__
712 #define __FUNCT__ "PetscSpaceViewFromOptions"
713 /*
714   PetscSpaceViewFromOptions - Processes command line options to determine if/how a PetscSpace is to be viewed.
715 
716   Collective on PetscSpace
717 
718   Input Parameters:
719 + sp   - the PetscSpace
720 . prefix - prefix to use for viewing, or NULL to use prefix of 'rnd'
721 - optionname - option to activate viewing
722 
723   Level: intermediate
724 
725 .keywords: PetscSpace, view, options, database
726 .seealso: VecViewFromOptions(), MatViewFromOptions()
727 */
728 PetscErrorCode PetscSpaceViewFromOptions(PetscSpace sp, const char prefix[], const char optionname[])
729 {
730   PetscViewer       viewer;
731   PetscViewerFormat format;
732   PetscBool         flg;
733   PetscErrorCode    ierr;
734 
735   PetscFunctionBegin;
736   if (prefix) {
737     ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject) sp), prefix, optionname, &viewer, &format, &flg);CHKERRQ(ierr);
738   } else {
739     ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject) sp), ((PetscObject) sp)->prefix, optionname, &viewer, &format, &flg);CHKERRQ(ierr);
740   }
741   if (flg) {
742     ierr = PetscViewerPushFormat(viewer, format);CHKERRQ(ierr);
743     ierr = PetscSpaceView(sp, viewer);CHKERRQ(ierr);
744     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
745     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
746   }
747   PetscFunctionReturn(0);
748 }
749 
750 #undef __FUNCT__
751 #define __FUNCT__ "PetscSpaceSetFromOptions"
752 /*@
753   PetscSpaceSetFromOptions - sets parameters in a PetscSpace from the options database
754 
755   Collective on PetscSpace
756 
757   Input Parameter:
758 . sp - the PetscSpace object to set options for
759 
760   Options Database:
761 . -petscspace_order the approximation order of the space
762 
763   Level: developer
764 
765 .seealso PetscSpaceView()
766 @*/
767 PetscErrorCode PetscSpaceSetFromOptions(PetscSpace sp)
768 {
769   const char    *defaultType;
770   char           typename[256];
771   PetscBool      flg;
772   PetscErrorCode ierr;
773 
774   PetscFunctionBegin;
775   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
776   if (!((PetscObject) sp)->type_name) {
777     defaultType = PETSCSPACEPOLYNOMIAL;
778   } else {
779     defaultType = ((PetscObject) sp)->type_name;
780   }
781   if (!PetscSpaceRegisterAllCalled) {ierr = PetscSpaceRegisterAll();CHKERRQ(ierr);}
782 
783   ierr = PetscObjectOptionsBegin((PetscObject) sp);CHKERRQ(ierr);
784   ierr = PetscOptionsList("-petscspace_type", "Linear space", "PetscSpaceSetType", PetscSpaceList, defaultType, typename, 256, &flg);CHKERRQ(ierr);
785   if (flg) {
786     ierr = PetscSpaceSetType(sp, typename);CHKERRQ(ierr);
787   } else if (!((PetscObject) sp)->type_name) {
788     ierr = PetscSpaceSetType(sp, defaultType);CHKERRQ(ierr);
789   }
790   ierr = PetscOptionsInt("-petscspace_order", "The approximation order", "PetscSpaceSetOrder", sp->order, &sp->order, NULL);CHKERRQ(ierr);
791   if (sp->ops->setfromoptions) {
792     ierr = (*sp->ops->setfromoptions)(sp);CHKERRQ(ierr);
793   }
794   /* process any options handlers added with PetscObjectAddOptionsHandler() */
795   ierr = PetscObjectProcessOptionsHandlers((PetscObject) sp);CHKERRQ(ierr);
796   ierr = PetscOptionsEnd();CHKERRQ(ierr);
797   ierr = PetscSpaceViewFromOptions(sp, NULL, "-petscspace_view");CHKERRQ(ierr);
798   PetscFunctionReturn(0);
799 }
800 
801 #undef __FUNCT__
802 #define __FUNCT__ "PetscSpaceSetUp"
803 /*@C
804   PetscSpaceSetUp - Construct data structures for the PetscSpace
805 
806   Collective on PetscSpace
807 
808   Input Parameter:
809 . sp - the PetscSpace object to setup
810 
811   Level: developer
812 
813 .seealso PetscSpaceView(), PetscSpaceDestroy()
814 @*/
815 PetscErrorCode PetscSpaceSetUp(PetscSpace sp)
816 {
817   PetscErrorCode ierr;
818 
819   PetscFunctionBegin;
820   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
821   if (sp->ops->setup) {ierr = (*sp->ops->setup)(sp);CHKERRQ(ierr);}
822   PetscFunctionReturn(0);
823 }
824 
825 #undef __FUNCT__
826 #define __FUNCT__ "PetscSpaceDestroy"
827 /*@
828   PetscSpaceDestroy - Destroys a PetscSpace object
829 
830   Collective on PetscSpace
831 
832   Input Parameter:
833 . sp - the PetscSpace object to destroy
834 
835   Level: developer
836 
837 .seealso PetscSpaceView()
838 @*/
839 PetscErrorCode PetscSpaceDestroy(PetscSpace *sp)
840 {
841   PetscErrorCode ierr;
842 
843   PetscFunctionBegin;
844   if (!*sp) PetscFunctionReturn(0);
845   PetscValidHeaderSpecific((*sp), PETSCSPACE_CLASSID, 1);
846 
847   if (--((PetscObject)(*sp))->refct > 0) {*sp = 0; PetscFunctionReturn(0);}
848   ((PetscObject) (*sp))->refct = 0;
849   /* if memory was published with AMS then destroy it */
850   ierr = PetscObjectAMSViewOff((PetscObject) *sp);CHKERRQ(ierr);
851 
852   ierr = DMDestroy(&(*sp)->dm);CHKERRQ(ierr);
853 
854   ierr = (*(*sp)->ops->destroy)(*sp);CHKERRQ(ierr);
855   ierr = PetscHeaderDestroy(sp);CHKERRQ(ierr);
856   PetscFunctionReturn(0);
857 }
858 
859 #undef __FUNCT__
860 #define __FUNCT__ "PetscSpaceCreate"
861 /*@
862   PetscSpaceCreate - Creates an empty PetscSpace object. The type can then be set with PetscSpaceSetType().
863 
864   Collective on MPI_Comm
865 
866   Input Parameter:
867 . comm - The communicator for the PetscSpace object
868 
869   Output Parameter:
870 . sp - The PetscSpace object
871 
872   Level: beginner
873 
874 .seealso: PetscSpaceSetType(), PETSCSPACEPOLYNOMIAL
875 @*/
876 PetscErrorCode PetscSpaceCreate(MPI_Comm comm, PetscSpace *sp)
877 {
878   PetscSpace     s;
879   PetscErrorCode ierr;
880 
881   PetscFunctionBegin;
882   PetscValidPointer(sp, 2);
883   *sp = NULL;
884 #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
885   ierr = PetscFEInitializePackage();CHKERRQ(ierr);
886 #endif
887 
888   ierr = PetscHeaderCreate(s, _p_PetscSpace, struct _PetscSpaceOps, PETSCSPACE_CLASSID, "PetscSpace", "Linear Space", "PetscSpace", comm, PetscSpaceDestroy, PetscSpaceView);CHKERRQ(ierr);
889   ierr = PetscMemzero(s->ops, sizeof(struct _PetscSpaceOps));CHKERRQ(ierr);
890 
891   s->order = 0;
892   ierr = DMShellCreate(comm, &s->dm);CHKERRQ(ierr);
893 
894   *sp = s;
895   PetscFunctionReturn(0);
896 }
897 
898 #undef __FUNCT__
899 #define __FUNCT__ "PetscSpaceGetDimension"
900 /* Dimension of the space, i.e. number of basis vectors */
901 PetscErrorCode PetscSpaceGetDimension(PetscSpace sp, PetscInt *dim)
902 {
903   PetscErrorCode ierr;
904 
905   PetscFunctionBegin;
906   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
907   PetscValidPointer(dim, 2);
908   *dim = 0;
909   if (sp->ops->getdimension) {ierr = (*sp->ops->getdimension)(sp, dim);CHKERRQ(ierr);}
910   PetscFunctionReturn(0);
911 }
912 
913 #undef __FUNCT__
914 #define __FUNCT__ "PetscSpaceGetOrder"
915 PetscErrorCode PetscSpaceGetOrder(PetscSpace sp, PetscInt *order)
916 {
917   PetscFunctionBegin;
918   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
919   PetscValidPointer(order, 2);
920   *order = sp->order;
921   PetscFunctionReturn(0);
922 }
923 
924 #undef __FUNCT__
925 #define __FUNCT__ "PetscSpaceSetOrder"
926 PetscErrorCode PetscSpaceSetOrder(PetscSpace sp, PetscInt order)
927 {
928   PetscFunctionBegin;
929   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
930   sp->order = order;
931   PetscFunctionReturn(0);
932 }
933 
934 #undef __FUNCT__
935 #define __FUNCT__ "PetscSpaceEvaluate"
936 PetscErrorCode PetscSpaceEvaluate(PetscSpace sp, PetscInt npoints, const PetscReal points[], PetscReal B[], PetscReal D[], PetscReal H[])
937 {
938   PetscErrorCode ierr;
939 
940   PetscFunctionBegin;
941   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
942   PetscValidPointer(points, 3);
943   if (B) PetscValidPointer(B, 4);
944   if (D) PetscValidPointer(D, 5);
945   if (H) PetscValidPointer(H, 6);
946   if (sp->ops->evaluate) {ierr = (*sp->ops->evaluate)(sp, npoints, points, B, D, H);CHKERRQ(ierr);}
947   PetscFunctionReturn(0);
948 }
949 
950 #undef __FUNCT__
951 #define __FUNCT__ "PetscSpaceSetFromOptions_Polynomial"
952 PetscErrorCode PetscSpaceSetFromOptions_Polynomial(PetscSpace sp)
953 {
954   PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
955   PetscErrorCode   ierr;
956 
957   PetscFunctionBegin;
958   ierr = PetscObjectOptionsBegin((PetscObject) sp);CHKERRQ(ierr);
959   ierr = PetscOptionsInt("-petscspace_num_variables", "The number of different variables, e.g. x and y", "PetscSpacePolynomialSetNumVariables", poly->numVariables, &poly->numVariables, NULL);CHKERRQ(ierr);
960   ierr = PetscOptionsBool("-petscspace_poly_sym", "Use only symmetric polynomials", "PetscSpacePolynomialSetSymmetric", poly->symmetric, &poly->symmetric, NULL);CHKERRQ(ierr);
961   ierr = PetscOptionsEnd();CHKERRQ(ierr);
962   PetscFunctionReturn(0);
963 }
964 
965 #undef __FUNCT__
966 #define __FUNCT__ "PetscSpacePolynomialView_Ascii"
967 PetscErrorCode PetscSpacePolynomialView_Ascii(PetscSpace sp, PetscViewer viewer)
968 {
969   PetscSpace_Poly  *poly = (PetscSpace_Poly *) sp->data;
970   PetscViewerFormat format;
971   PetscErrorCode    ierr;
972 
973   PetscFunctionBegin;
974   ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
975   if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
976     ierr = PetscViewerASCIIPrintf(viewer, "Polynomial space in %d variables of order %d", poly->numVariables, sp->order);CHKERRQ(ierr);
977   } else {
978     ierr = PetscViewerASCIIPrintf(viewer, "Polynomial space in %d variables of order %d", poly->numVariables, sp->order);CHKERRQ(ierr);
979   }
980   PetscFunctionReturn(0);
981 }
982 
983 #undef __FUNCT__
984 #define __FUNCT__ "PetscSpaceView_Polynomial"
985 PetscErrorCode PetscSpaceView_Polynomial(PetscSpace sp, PetscViewer viewer)
986 {
987   PetscBool      iascii;
988   PetscErrorCode ierr;
989 
990   PetscFunctionBegin;
991   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
992   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
993   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
994   if (iascii) {ierr = PetscSpacePolynomialView_Ascii(sp, viewer);CHKERRQ(ierr);}
995   PetscFunctionReturn(0);
996 }
997 
998 #undef __FUNCT__
999 #define __FUNCT__ "PetscSpaceSetUp_Polynomial"
1000 PetscErrorCode PetscSpaceSetUp_Polynomial(PetscSpace sp)
1001 {
1002   PetscSpace_Poly *poly    = (PetscSpace_Poly *) sp->data;
1003   PetscInt         ndegree = sp->order+1;
1004   PetscInt         deg;
1005   PetscErrorCode   ierr;
1006 
1007   PetscFunctionBegin;
1008   ierr = PetscMalloc(ndegree * sizeof(PetscInt), &poly->degrees);CHKERRQ(ierr);
1009   for (deg = 0; deg < ndegree; ++deg) poly->degrees[deg] = deg;
1010   PetscFunctionReturn(0);
1011 }
1012 
1013 #undef __FUNCT__
1014 #define __FUNCT__ "PetscSpaceDestroy_Polynomial"
1015 PetscErrorCode PetscSpaceDestroy_Polynomial(PetscSpace sp)
1016 {
1017   PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
1018   PetscErrorCode   ierr;
1019 
1020   PetscFunctionBegin;
1021   ierr = PetscFree(poly->degrees);CHKERRQ(ierr);
1022   ierr = PetscFree(poly);CHKERRQ(ierr);
1023   PetscFunctionReturn(0);
1024 }
1025 
1026 #undef __FUNCT__
1027 #define __FUNCT__ "PetscSpaceGetDimension_Polynomial"
1028 PetscErrorCode PetscSpaceGetDimension_Polynomial(PetscSpace sp, PetscInt *dim)
1029 {
1030   PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
1031   PetscInt         deg  = sp->order;
1032   PetscInt         n    = poly->numVariables, i;
1033   PetscReal        D    = 1.0;
1034 
1035   PetscFunctionBegin;
1036   for (i = 1; i <= n; ++i) {
1037     D *= ((PetscReal) (deg+i))/i;
1038   }
1039   *dim = (PetscInt) (D + 0.5);
1040   PetscFunctionReturn(0);
1041 }
1042 
1043 #undef __FUNCT__
1044 #define __FUNCT__ "LatticePoint_Internal"
1045 /*
1046   LatticePoint_Internal - Returns all tuples of size 'len' with nonnegative integers that sum up to 'sum'.
1047 
1048   Input Parameters:
1049 + len - The length of the tuple
1050 . sum - The sum of all entries in the tuple
1051 - ind - The current multi-index of the tuple, initialized to the 0 tuple
1052 
1053   Output Parameter:
1054 + ind - The multi-index of the tuple, -1 indicates the iteration has terminated
1055 . tup - A tuple of len integers addig to sum
1056 
1057   Level: developer
1058 
1059 .seealso:
1060 */
1061 static PetscErrorCode LatticePoint_Internal(PetscInt len, PetscInt sum, PetscInt ind[], PetscInt tup[])
1062 {
1063   PetscInt       i;
1064   PetscErrorCode ierr;
1065 
1066   PetscFunctionBegin;
1067   if (len == 1) {
1068     ind[0] = -1;
1069     tup[0] = sum;
1070   } else if (sum == 0) {
1071     for (i = 0; i < len; ++i) {ind[0] = -1; tup[i] = 0;}
1072   } else {
1073     tup[0] = sum - ind[0];
1074     ierr = LatticePoint_Internal(len-1, ind[0], &ind[1], &tup[1]);CHKERRQ(ierr);
1075     if (ind[1] < 0) {
1076       if (ind[0] == sum) {ind[0] = -1;}
1077       else               {ind[1] = 0; ++ind[0];}
1078     }
1079   }
1080   PetscFunctionReturn(0);
1081 }
1082 
1083 #undef __FUNCT__
1084 #define __FUNCT__ "PetscSpaceEvaluate_Polynomial"
1085 PetscErrorCode PetscSpaceEvaluate_Polynomial(PetscSpace sp, PetscInt npoints, const PetscReal points[], PetscReal B[], PetscReal D[], PetscReal H[])
1086 {
1087   PetscSpace_Poly *poly    = (PetscSpace_Poly *) sp->data;
1088   DM               dm      = sp->dm;
1089   PetscInt         ndegree = sp->order+1;
1090   PetscInt        *degrees = poly->degrees;
1091   PetscInt         dim     = poly->numVariables;
1092   PetscReal       *lpoints, *tmp, *LB, *LD, *LH;
1093   PetscInt        *ind, *tup;
1094   PetscInt         pdim, d, i, p, deg, o;
1095   PetscErrorCode   ierr;
1096 
1097   PetscFunctionBegin;
1098   ierr = PetscSpaceGetDimension(sp, &pdim);CHKERRQ(ierr);
1099   ierr = DMGetWorkArray(dm, npoints, PETSC_REAL, &lpoints);CHKERRQ(ierr);
1100   ierr = DMGetWorkArray(dm, npoints*ndegree*3, PETSC_REAL, &tmp);CHKERRQ(ierr);
1101   if (B) {ierr = DMGetWorkArray(dm, npoints*dim*ndegree, PETSC_REAL, &LB);CHKERRQ(ierr);}
1102   if (D) {ierr = DMGetWorkArray(dm, npoints*dim*ndegree, PETSC_REAL, &LD);CHKERRQ(ierr);}
1103   if (H) {ierr = DMGetWorkArray(dm, npoints*dim*ndegree, PETSC_REAL, &LH);CHKERRQ(ierr);}
1104   for (d = 0; d < dim; ++d) {
1105     for (p = 0; p < npoints; ++p) {
1106       lpoints[p] = points[p*dim+d];
1107     }
1108     ierr = PetscDTLegendreEval(npoints, lpoints, ndegree, degrees, tmp, &tmp[1*npoints*ndegree], &tmp[2*npoints*ndegree]);CHKERRQ(ierr);
1109     /* LB, LD, LH (ndegree * dim x npoints) */
1110     for (deg = 0; deg < ndegree; ++deg) {
1111       for (p = 0; p < npoints; ++p) {
1112         if (B) LB[(deg*dim + d)*npoints + p] = tmp[(0*npoints + p)*ndegree+deg];
1113         if (D) LD[(deg*dim + d)*npoints + p] = tmp[(1*npoints + p)*ndegree+deg];
1114         if (H) LH[(deg*dim + d)*npoints + p] = tmp[(2*npoints + p)*ndegree+deg];
1115       }
1116     }
1117   }
1118   /* Multiply by A (pdim x ndegree * dim) */
1119   ierr = PetscMalloc2(dim,PetscInt,&ind,dim,PetscInt,&tup);CHKERRQ(ierr);
1120   if (B) {
1121     /* B (npoints x pdim) */
1122     i = 0;
1123     for (o = 0; o <= sp->order; ++o) {
1124       ierr = PetscMemzero(ind, dim * sizeof(PetscInt));CHKERRQ(ierr);
1125       while (ind[0] >= 0) {
1126         ierr = LatticePoint_Internal(dim, o, ind, tup);CHKERRQ(ierr);
1127         for (p = 0; p < npoints; ++p) {
1128           B[p*pdim + i] = 1.0;
1129           for (d = 0; d < dim; ++d) {
1130             B[p*pdim + i] *= LB[(tup[d]*dim + d)*npoints + p];
1131           }
1132         }
1133         ++i;
1134       }
1135     }
1136   }
1137   ierr = PetscFree2(ind,tup);CHKERRQ(ierr);
1138   if (D) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to code first derivatives");
1139   if (H) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to code second derivatives");
1140   if (B) {ierr = DMRestoreWorkArray(dm, npoints*dim*ndegree, PETSC_REAL, &LB);CHKERRQ(ierr);}
1141   if (D) {ierr = DMRestoreWorkArray(dm, npoints*dim*ndegree, PETSC_REAL, &LD);CHKERRQ(ierr);}
1142   if (H) {ierr = DMRestoreWorkArray(dm, npoints*dim*ndegree, PETSC_REAL, &LH);CHKERRQ(ierr);}
1143   ierr = DMRestoreWorkArray(dm, npoints*ndegree*3, PETSC_REAL, &tmp);CHKERRQ(ierr);
1144   ierr = DMRestoreWorkArray(dm, npoints, PETSC_REAL, &lpoints);CHKERRQ(ierr);
1145   PetscFunctionReturn(0);
1146 }
1147 
1148 #undef __FUNCT__
1149 #define __FUNCT__ "PetscSpaceInitialize_Polynomial"
1150 PetscErrorCode PetscSpaceInitialize_Polynomial(PetscSpace sp)
1151 {
1152   PetscFunctionBegin;
1153   sp->ops->setfromoptions = PetscSpaceSetFromOptions_Polynomial;
1154   sp->ops->setup          = PetscSpaceSetUp_Polynomial;
1155   sp->ops->view           = PetscSpaceView_Polynomial;
1156   sp->ops->destroy        = PetscSpaceDestroy_Polynomial;
1157   sp->ops->getdimension   = PetscSpaceGetDimension_Polynomial;
1158   sp->ops->evaluate       = PetscSpaceEvaluate_Polynomial;
1159   PetscFunctionReturn(0);
1160 }
1161 
1162 /*MC
1163   PETSCSPACEPOLYNOMIAL = "poly" - A PetscSpace object that encapsulates a polynomial space, e.g. P1 is the space of linear polynomials.
1164 
1165   Level: intermediate
1166 
1167 .seealso: PetscSpaceType, PetscSpaceCreate(), PetscSpaceSetType()
1168 M*/
1169 
1170 #undef __FUNCT__
1171 #define __FUNCT__ "PetscSpaceCreate_Polynomial"
1172 PETSC_EXTERN PetscErrorCode PetscSpaceCreate_Polynomial(PetscSpace sp)
1173 {
1174   PetscSpace_Poly *poly;
1175   PetscErrorCode   ierr;
1176 
1177   PetscFunctionBegin;
1178   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
1179   ierr     = PetscNewLog(sp, PetscSpace_Poly, &poly);CHKERRQ(ierr);
1180   sp->data = poly;
1181 
1182   poly->numVariables = 0;
1183   poly->symmetric    = PETSC_FALSE;
1184   poly->degrees      = 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 = PetscDualSpaceGetDimension(*sp, &dim);CHKERRQ(ierr);
1538   for (f = 0; f < dim; ++f) {
1539     ierr = PetscQuadratureDestroy(&(*sp)->functional[f]);CHKERRQ(ierr);
1540   }
1541   ierr = PetscFree((*sp)->functional);CHKERRQ(ierr);
1542   ierr = DMDestroy(&(*sp)->dm);CHKERRQ(ierr);
1543 
1544   if ((*sp)->ops->destroy) {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   PetscReal     *qpoints, *qweights;
1674   PetscInt       depth, dim, pdim, *pStart, *pEnd, coneSize, d, n, f = 0;
1675   PetscErrorCode ierr;
1676 
1677   PetscFunctionBegin;
1678   ierr = PetscDualSpaceGetDimension(sp, &pdim);CHKERRQ(ierr);
1679   ierr = PetscMalloc(pdim * sizeof(PetscQuadrature), &sp->functional);CHKERRQ(ierr);
1680   /* Classify element type */
1681   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
1682   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
1683   ierr = PetscMalloc2(depth+1,PetscInt,&pStart,depth+1,PetscInt,&pEnd);CHKERRQ(ierr);
1684   for (d = 0; d < depth; ++d) {
1685     ierr = DMPlexGetDepthStratum(dm, d, &pStart[d], &pEnd[d]);CHKERRQ(ierr);
1686   }
1687   ierr = DMPlexGetConeSize(dm, pStart[depth], &coneSize);CHKERRQ(ierr);
1688   ierr = DMPlexGetCoordinateSection(dm, &csection);CHKERRQ(ierr);
1689   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
1690   if (coneSize == dim+1) {
1691     PetscInt *closure = NULL, closureSize, c;
1692 
1693     /* Simplex */
1694     ierr = DMPlexGetTransitiveClosure(dm, pStart[depth], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1695     for (c = 0; c < closureSize*2; c += 2) {
1696       const PetscInt p = closure[c];
1697 
1698       if ((p >= pStart[0]) && (p < pEnd[0])) {
1699         /* Vertices */
1700         const PetscScalar *coords;
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         PetscScalar *coords;
1720         PetscInt     k;
1721 
1722         if (order < 2) continue;
1723         coords = NULL;
1724         ierr = DMPlexVecGetClosure(dm, csection, coordinates, p, &n, &coords);CHKERRQ(ierr);
1725         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);
1726         for (k = 1; k < order; ++k) {
1727           sp->functional[f].numQuadPoints = 1;
1728           ierr = PetscMalloc(sp->functional[f].numQuadPoints*dim * sizeof(PetscReal), &qpoints);CHKERRQ(ierr);
1729           ierr = PetscMalloc(sp->functional[f].numQuadPoints     * sizeof(PetscReal), &qweights);CHKERRQ(ierr);
1730           for (d = 0; d < dim; ++d) {qpoints[d] = k*(coords[1*dim+d] - coords[0*dim+d])/order + coords[0*dim+d];}
1731           qweights[0] = 1.0;
1732           sp->functional[f].quadPoints  = qpoints;
1733           sp->functional[f].quadWeights = qweights;
1734           ++f;
1735         }
1736         ierr = DMPlexVecRestoreClosure(dm, csection, coordinates, p, &n, &coords);CHKERRQ(ierr);
1737       } else if ((p >= pStart[depth-1]) && (p < pEnd[depth-1])) {
1738         /* Faces */
1739 
1740         if (order < 3) continue;
1741         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to implement faces");
1742       } else if ((p >= pStart[depth]) && (p < pEnd[depth])) {
1743         /* Cells */
1744 
1745         if ((order > 0) && (order < 3)) continue;
1746         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to implement cells");
1747       }
1748     }
1749     ierr = DMPlexRestoreTransitiveClosure(dm, pStart[depth], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1750   } else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot handle cells with cone size %d", coneSize);
1751   ierr = PetscFree2(pStart,pEnd);CHKERRQ(ierr);
1752   if (f != pdim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of dual basis vector %d not equal to dimension %d", f, pdim);
1753   PetscFunctionReturn(0);
1754 }
1755 
1756 #undef __FUNCT__
1757 #define __FUNCT__ "PetscDualSpaceDestroy_Lagrange"
1758 PetscErrorCode PetscDualSpaceDestroy_Lagrange(PetscDualSpace sp)
1759 {
1760   PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data;
1761   PetscErrorCode      ierr;
1762 
1763   PetscFunctionBegin;
1764   ierr = PetscFree(lag);CHKERRQ(ierr);
1765   PetscFunctionReturn(0);
1766 }
1767 
1768 #undef __FUNCT__
1769 #define __FUNCT__ "PetscDualSpaceGetDimension_Lagrange"
1770 PetscErrorCode PetscDualSpaceGetDimension_Lagrange(PetscDualSpace sp, PetscInt *dim)
1771 {
1772   PetscInt            deg = sp->order;
1773   PetscReal           D   = 1.0;
1774   PetscInt            n, i;
1775   PetscErrorCode      ierr;
1776 
1777   PetscFunctionBegin;
1778   /* TODO: Assumes simplices */
1779   ierr = DMPlexGetDimension(sp->dm, &n);CHKERRQ(ierr);
1780   for (i = 1; i <= n; ++i) {
1781     D *= ((PetscReal) (deg+i))/i;
1782   }
1783   *dim = (PetscInt) (D + 0.5);
1784   PetscFunctionReturn(0);
1785 }
1786 
1787 #undef __FUNCT__
1788 #define __FUNCT__ "PetscDualSpaceInitialize_Lagrange"
1789 PetscErrorCode PetscDualSpaceInitialize_Lagrange(PetscDualSpace sp)
1790 {
1791   PetscFunctionBegin;
1792   sp->ops->setfromoptions = NULL;
1793   sp->ops->setup          = PetscDualSpaceSetUp_Lagrange;
1794   sp->ops->view           = NULL;
1795   sp->ops->destroy        = PetscDualSpaceDestroy_Lagrange;
1796   sp->ops->getdimension   = PetscDualSpaceGetDimension_Lagrange;
1797   PetscFunctionReturn(0);
1798 }
1799 
1800 /*MC
1801   PETSCDUALSPACELAGRANGE = "lagrange" - A PetscDualSpace object that encapsulates a dual space of pointwise evaluation functionals
1802 
1803   Level: intermediate
1804 
1805 .seealso: PetscDualSpaceType, PetscDualSpaceCreate(), PetscDualSpaceSetType()
1806 M*/
1807 
1808 #undef __FUNCT__
1809 #define __FUNCT__ "PetscDualSpaceCreate_Lagrange"
1810 PETSC_EXTERN PetscErrorCode PetscDualSpaceCreate_Lagrange(PetscDualSpace sp)
1811 {
1812   PetscDualSpace_Lag *lag;
1813   PetscErrorCode      ierr;
1814 
1815   PetscFunctionBegin;
1816   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1817   ierr     = PetscNewLog(sp, PetscDualSpace_Lag, &lag);CHKERRQ(ierr);
1818   sp->data = lag;
1819 
1820   /* lag->n = 0; */
1821 
1822   ierr = PetscDualSpaceInitialize_Lagrange(sp);CHKERRQ(ierr);
1823   PetscFunctionReturn(0);
1824 }
1825 
1826 
1827 PetscInt PETSCFE_CLASSID = 0;
1828 
1829 #undef __FUNCT__
1830 #define __FUNCT__ "PetscFEView"
1831 /*@C
1832   PetscFEView - Views a PetscFE
1833 
1834   Collective on PetscFE
1835 
1836   Input Parameter:
1837 + sp - the PetscFE object to view
1838 - v  - the viewer
1839 
1840   Level: developer
1841 
1842 .seealso PetscFEDestroy()
1843 @*/
1844 PetscErrorCode PetscFEView(PetscFE fem, PetscViewer v)
1845 {
1846   PetscErrorCode ierr;
1847 
1848   PetscFunctionBegin;
1849   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
1850   if (!v) {
1851     ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) fem), &v);CHKERRQ(ierr);
1852   }
1853   if (fem->ops->view) {
1854     ierr = (*fem->ops->view)(fem, v);CHKERRQ(ierr);
1855   }
1856   PetscFunctionReturn(0);
1857 }
1858 
1859 #undef __FUNCT__
1860 #define __FUNCT__ "PetscFEDestroy"
1861 /*@
1862   PetscFEDestroy - Destroys a PetscFE object
1863 
1864   Collective on PetscFE
1865 
1866   Input Parameter:
1867 . fem - the PetscFE object to destroy
1868 
1869   Level: developer
1870 
1871 .seealso PetscFEView()
1872 @*/
1873 PetscErrorCode PetscFEDestroy(PetscFE *fem)
1874 {
1875   PetscErrorCode ierr;
1876 
1877   PetscFunctionBegin;
1878   if (!*fem) PetscFunctionReturn(0);
1879   PetscValidHeaderSpecific((*fem), PETSCFE_CLASSID, 1);
1880 
1881   if (--((PetscObject)(*fem))->refct > 0) {*fem = 0; PetscFunctionReturn(0);}
1882   ((PetscObject) (*fem))->refct = 0;
1883   /* if memory was published with AMS then destroy it */
1884   ierr = PetscObjectAMSViewOff((PetscObject) *fem);CHKERRQ(ierr);
1885 
1886   ierr = PetscSpaceDestroy(&(*fem)->basisSpace);CHKERRQ(ierr);
1887   ierr = PetscDualSpaceDestroy(&(*fem)->dualSpace);CHKERRQ(ierr);
1888   ierr = PetscQuadratureDestroy(&(*fem)->quadrature);CHKERRQ(ierr);
1889 
1890   if ((*fem)->ops->destroy) {ierr = (*(*fem)->ops->destroy)(*fem);CHKERRQ(ierr);}
1891   ierr = PetscHeaderDestroy(fem);CHKERRQ(ierr);
1892   PetscFunctionReturn(0);
1893 }
1894 
1895 #undef __FUNCT__
1896 #define __FUNCT__ "PetscFECreate"
1897 /*@
1898   PetscFECreate - Creates an empty PetscFE object. The type can then be set with PetscFESetType().
1899 
1900   Collective on MPI_Comm
1901 
1902   Input Parameter:
1903 . comm - The communicator for the PetscFE object
1904 
1905   Output Parameter:
1906 . fem - The PetscFE object
1907 
1908   Level: beginner
1909 
1910 .seealso: PetscFESetType(), PETSCFEGALERKIN
1911 @*/
1912 PetscErrorCode PetscFECreate(MPI_Comm comm, PetscFE *fem)
1913 {
1914   PetscFE        f;
1915   PetscErrorCode ierr;
1916 
1917   PetscFunctionBegin;
1918   PetscValidPointer(fem, 2);
1919   *fem = NULL;
1920 #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
1921   ierr = PetscFEInitializePackage();CHKERRQ(ierr);
1922 #endif
1923 
1924   ierr = PetscHeaderCreate(f, _p_PetscFE, struct _PetscFEOps, PETSCFE_CLASSID, "PetscFE", "Finite Element", "PetscFE", comm, PetscFEDestroy, PetscFEView);CHKERRQ(ierr);
1925   ierr = PetscMemzero(f->ops, sizeof(struct _PetscFEOps));CHKERRQ(ierr);
1926 
1927   f->basisSpace    = NULL;
1928   f->dualSpace     = NULL;
1929   f->numComponents = 1;
1930   ierr = PetscMemzero(&f->quadrature, sizeof(PetscQuadrature));CHKERRQ(ierr);
1931 
1932   *fem = f;
1933   PetscFunctionReturn(0);
1934 }
1935 
1936 #undef __FUNCT__
1937 #define __FUNCT__ "PetscFEGetDimension"
1938 PetscErrorCode PetscFEGetDimension(PetscFE fem, PetscInt *dim)
1939 {
1940   PetscErrorCode ierr;
1941 
1942   PetscFunctionBegin;
1943   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
1944   PetscValidPointer(dim, 2);
1945   ierr = PetscSpaceGetDimension(fem->basisSpace, dim);CHKERRQ(ierr);
1946   PetscFunctionReturn(0);
1947 }
1948 
1949 #undef __FUNCT__
1950 #define __FUNCT__ "PetscFESetNumComponents"
1951 PetscErrorCode PetscFESetNumComponents(PetscFE fem, PetscInt comp)
1952 {
1953   PetscFunctionBegin;
1954   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
1955   fem->numComponents = comp;
1956   PetscFunctionReturn(0);
1957 }
1958 
1959 #undef __FUNCT__
1960 #define __FUNCT__ "PetscFEGetNumComponents"
1961 PetscErrorCode PetscFEGetNumComponents(PetscFE fem, PetscInt *comp)
1962 {
1963   PetscFunctionBegin;
1964   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
1965   PetscValidPointer(comp, 2);
1966   *comp = fem->numComponents;
1967   PetscFunctionReturn(0);
1968 }
1969 
1970 #undef __FUNCT__
1971 #define __FUNCT__ "PetscFEGetBasisSpace"
1972 PetscErrorCode PetscFEGetBasisSpace(PetscFE fem, PetscSpace *sp)
1973 {
1974   PetscFunctionBegin;
1975   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
1976   PetscValidPointer(sp, 2);
1977   *sp = fem->basisSpace;
1978   PetscFunctionReturn(0);
1979 }
1980 
1981 #undef __FUNCT__
1982 #define __FUNCT__ "PetscFESetBasisSpace"
1983 PetscErrorCode PetscFESetBasisSpace(PetscFE fem, PetscSpace sp)
1984 {
1985   PetscErrorCode ierr;
1986 
1987   PetscFunctionBegin;
1988   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
1989   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 2);
1990   ierr = PetscSpaceDestroy(&fem->basisSpace);CHKERRQ(ierr);
1991   fem->basisSpace = sp;
1992   ierr = PetscObjectReference((PetscObject) fem->basisSpace);CHKERRQ(ierr);
1993   PetscFunctionReturn(0);
1994 }
1995 
1996 #undef __FUNCT__
1997 #define __FUNCT__ "PetscFEGetDualSpace"
1998 PetscErrorCode PetscFEGetDualSpace(PetscFE fem, PetscDualSpace *sp)
1999 {
2000   PetscFunctionBegin;
2001   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
2002   PetscValidPointer(sp, 2);
2003   *sp = fem->dualSpace;
2004   PetscFunctionReturn(0);
2005 }
2006 
2007 #undef __FUNCT__
2008 #define __FUNCT__ "PetscFESetDualSpace"
2009 PetscErrorCode PetscFESetDualSpace(PetscFE fem, PetscDualSpace sp)
2010 {
2011   PetscErrorCode ierr;
2012 
2013   PetscFunctionBegin;
2014   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
2015   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 2);
2016   ierr = PetscDualSpaceDestroy(&fem->dualSpace);CHKERRQ(ierr);
2017   fem->dualSpace = sp;
2018   ierr = PetscObjectReference((PetscObject) fem->dualSpace);CHKERRQ(ierr);
2019   PetscFunctionReturn(0);
2020 }
2021 
2022 #undef __FUNCT__
2023 #define __FUNCT__ "PetscFEGetQuadrature"
2024 PetscErrorCode PetscFEGetQuadrature(PetscFE fem, PetscQuadrature *q)
2025 {
2026   PetscFunctionBegin;
2027   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
2028   PetscValidPointer(q, 2);
2029   *q = fem->quadrature;
2030   PetscFunctionReturn(0);
2031 }
2032 
2033 #undef __FUNCT__
2034 #define __FUNCT__ "PetscFESetQuadrature"
2035 PetscErrorCode PetscFESetQuadrature(PetscFE fem, PetscQuadrature q)
2036 {
2037   PetscErrorCode ierr;
2038 
2039   PetscFunctionBegin;
2040   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
2041   ierr = PetscQuadratureDestroy(&fem->quadrature);CHKERRQ(ierr);
2042   fem->quadrature = q;
2043   PetscFunctionReturn(0);
2044 }
2045 
2046 #undef __FUNCT__
2047 #define __FUNCT__ "PetscFEGetTabulation"
2048 PetscErrorCode PetscFEGetTabulation(PetscFE fem, PetscReal **B, PetscReal **D, PetscReal **H)
2049 {
2050   DM              dm;
2051   PetscInt        pdim; /* Dimension of FE space P */
2052   PetscInt        dim;  /* Spatial dimension */
2053   PetscInt        comp; /* Field components */
2054   PetscInt        npoints = fem->quadrature.numQuadPoints;
2055   PetscReal      *points  = fem->quadrature.quadPoints;
2056   PetscReal      *tmpB, *invV;
2057   PetscInt        p, j, k;
2058   PetscErrorCode  ierr;
2059 
2060   PetscFunctionBegin;
2061   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
2062   PetscValidPointer(points, 3);
2063   if (B) PetscValidPointer(B, 4);
2064   if (D) PetscValidPointer(D, 5);
2065   if (H) PetscValidPointer(H, 6);
2066   ierr = PetscDualSpaceGetDM(fem->dualSpace, &dm);CHKERRQ(ierr);
2067 
2068   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
2069   ierr = PetscSpaceGetDimension(fem->basisSpace, &pdim);CHKERRQ(ierr);
2070   ierr = PetscFEGetNumComponents(fem, &comp);CHKERRQ(ierr);
2071   /* 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); */
2072 
2073   if (B) {
2074     ierr = DMGetWorkArray(dm, npoints*pdim*comp, PETSC_REAL, B);CHKERRQ(ierr);
2075     ierr = DMGetWorkArray(dm, npoints*pdim, PETSC_REAL, &tmpB);CHKERRQ(ierr);
2076   }
2077   if (D) {ierr = DMGetWorkArray(dm, npoints*pdim*dim, PETSC_REAL, D);CHKERRQ(ierr);}
2078   if (H) {ierr = DMGetWorkArray(dm, npoints*pdim*dim*dim, PETSC_REAL, H);CHKERRQ(ierr);}
2079   ierr = PetscSpaceEvaluate(fem->basisSpace, npoints, points, B ? tmpB : NULL, D ? *D : NULL, H ? *H : NULL);CHKERRQ(ierr);
2080 
2081   ierr = DMGetWorkArray(dm, pdim*pdim, PETSC_REAL, &invV);CHKERRQ(ierr);
2082   for (j = 0; j < pdim; ++j) {
2083     PetscReal      *Bf;
2084     PetscQuadrature f;
2085     PetscInt        q;
2086 
2087     ierr = PetscDualSpaceGetFunctional(fem->dualSpace, j, &f);
2088     ierr = DMGetWorkArray(dm, f.numQuadPoints*pdim, PETSC_REAL, &Bf);CHKERRQ(ierr);
2089     ierr = PetscSpaceEvaluate(fem->basisSpace, f.numQuadPoints, f.quadPoints, Bf, NULL, NULL);CHKERRQ(ierr);
2090     for (k = 0; k < pdim; ++k) {
2091       /* n_j \cdot \phi_k */
2092       invV[j*pdim+k] = 0.0;
2093       for (q = 0; q < f.numQuadPoints; ++q) {
2094         invV[j*pdim+k] += Bf[q*pdim+k]*f.quadWeights[q];
2095       }
2096     }
2097     ierr = DMRestoreWorkArray(dm, f.numQuadPoints*pdim, PETSC_REAL, &Bf);CHKERRQ(ierr);
2098   }
2099   {
2100     PetscReal    *work;
2101     PetscBLASInt *pivots;
2102     PetscBLASInt  n = pdim, info;
2103 
2104     ierr = DMGetWorkArray(dm, pdim, PETSC_INT, &pivots);CHKERRQ(ierr);
2105     ierr = DMGetWorkArray(dm, pdim, PETSC_REAL, &work);CHKERRQ(ierr);
2106     PetscStackCallBLAS("LAPACKgetrf", LAPACKgetrf_(&n, &n, invV, &n, pivots, &info));
2107     PetscStackCallBLAS("LAPACKgetri", LAPACKgetri_(&n, invV, &n, pivots, work, &n, &info));
2108     ierr = DMRestoreWorkArray(dm, pdim, PETSC_INT, &pivots);CHKERRQ(ierr);
2109     ierr = DMRestoreWorkArray(dm, pdim, PETSC_REAL, &work);CHKERRQ(ierr);
2110   }
2111   for (p = 0; p < npoints; ++p) {
2112     if (B) {
2113       /* Multiply by V^{-1} (pdim x pdim) */
2114       for (j = 0; j < pdim; ++j) {
2115         const PetscInt i = (p*pdim + j)*comp;
2116         PetscInt       c;
2117 
2118         (*B)[i] = 0.0;
2119         for (k = 0; k < pdim; ++k) {
2120           (*B)[i] += invV[k*pdim+j] * tmpB[p*pdim + k];
2121         }
2122         for (c = 1; c < comp; ++c) {
2123           (*B)[i+c] = (*B)[i];
2124         }
2125       }
2126     }
2127   }
2128   ierr = DMRestoreWorkArray(dm, pdim*pdim, PETSC_REAL, &invV);CHKERRQ(ierr);
2129   if (B) {
2130     ierr = DMRestoreWorkArray(dm, npoints*pdim, PETSC_REAL, &tmpB);CHKERRQ(ierr);
2131   }
2132   PetscFunctionReturn(0);
2133 }
2134 
2135 #undef __FUNCT__
2136 #define __FUNCT__ "PetscFERestoreTabulation"
2137 PetscErrorCode PetscFERestoreTabulation(PetscFE fem, PetscReal **B, PetscReal **D, PetscReal **D2)
2138 {
2139   DM             dm;
2140   PetscErrorCode ierr;
2141 
2142   PetscFunctionBegin;
2143   PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1);
2144   ierr = PetscDualSpaceGetDM(fem->dualSpace, &dm);CHKERRQ(ierr);
2145   if (B)  {ierr = DMRestoreWorkArray(dm, 0, PETSC_REAL, B);CHKERRQ(ierr);}
2146   if (D)  {ierr = DMRestoreWorkArray(dm, 0, PETSC_REAL, D);CHKERRQ(ierr);}
2147   if (D2) {ierr = DMRestoreWorkArray(dm, 0, PETSC_REAL, D2);CHKERRQ(ierr);}
2148   PetscFunctionReturn(0);
2149 }
2150