xref: /petsc/src/dm/dt/interface/dt.c (revision 7fdeb8b9be692273ce99bf256278a44482a052ca)
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*/
9 #include <petscblaslapack.h>
10 #include <petsc-private/petscimpl.h>
11 #include <petsc-private/dtimpl.h>
12 #include <petscviewer.h>
13 #include <petscdmplex.h>
14 #include <petscdmshell.h>
15 
16 #undef __FUNCT__
17 #define __FUNCT__ "PetscQuadratureCreate"
18 PetscErrorCode PetscQuadratureCreate(MPI_Comm comm, PetscQuadrature *q)
19 {
20   PetscErrorCode ierr;
21 
22   PetscFunctionBegin;
23   PetscValidPointer(q, 2);
24   ierr = DMInitializePackage();CHKERRQ(ierr);
25   ierr = PetscHeaderCreate(*q,_p_PetscQuadrature,int,PETSC_OBJECT_CLASSID,"PetscQuadrature","Quadrature","DT",comm,PetscQuadratureDestroy,PetscQuadratureView);CHKERRQ(ierr);
26   (*q)->dim       = -1;
27   (*q)->numPoints = 0;
28   (*q)->points    = NULL;
29   (*q)->weights   = NULL;
30   PetscFunctionReturn(0);
31 }
32 
33 #undef __FUNCT__
34 #define __FUNCT__ "PetscQuadratureDestroy"
35 PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *q)
36 {
37   PetscErrorCode ierr;
38 
39   PetscFunctionBegin;
40   if (!*q) PetscFunctionReturn(0);
41   PetscValidHeaderSpecific((*q),PETSC_OBJECT_CLASSID,1);
42   if (--((PetscObject)(*q))->refct > 0) {
43     *q = NULL;
44     PetscFunctionReturn(0);
45   }
46   ierr = PetscFree((*q)->points);CHKERRQ(ierr);
47   ierr = PetscFree((*q)->weights);CHKERRQ(ierr);
48   ierr = PetscHeaderDestroy(q);CHKERRQ(ierr);
49   PetscFunctionReturn(0);
50 }
51 
52 #undef __FUNCT__
53 #define __FUNCT__ "PetscQuadratureGetData"
54 PetscErrorCode PetscQuadratureGetData(PetscQuadrature q, PetscInt *dim, PetscInt *npoints, const PetscReal *points[], const PetscReal *weights[])
55 {
56   PetscFunctionBegin;
57   PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1);
58   if (dim) {
59     PetscValidPointer(dim, 2);
60     *dim = q->dim;
61   }
62   if (npoints) {
63     PetscValidPointer(npoints, 3);
64     *npoints = q->numPoints;
65   }
66   if (points) {
67     PetscValidPointer(points, 4);
68     *points = q->points;
69   }
70   if (weights) {
71     PetscValidPointer(weights, 5);
72     *weights = q->weights;
73   }
74   PetscFunctionReturn(0);
75 }
76 
77 #undef __FUNCT__
78 #define __FUNCT__ "PetscQuadratureSetData"
79 PetscErrorCode PetscQuadratureSetData(PetscQuadrature q, PetscInt dim, PetscInt npoints, const PetscReal points[], const PetscReal weights[])
80 {
81   PetscFunctionBegin;
82   PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1);
83   if (dim >= 0)     q->dim       = dim;
84   if (npoints >= 0) q->numPoints = npoints;
85   if (points) {
86     PetscValidPointer(points, 4);
87     q->points = points;
88   }
89   if (weights) {
90     PetscValidPointer(weights, 5);
91     q->weights = weights;
92   }
93   PetscFunctionReturn(0);
94 }
95 
96 #undef __FUNCT__
97 #define __FUNCT__ "PetscQuadratureView"
98 PetscErrorCode PetscQuadratureView(PetscQuadrature quad, PetscViewer viewer)
99 {
100   PetscInt       q, d;
101   PetscErrorCode ierr;
102 
103   PetscFunctionBegin;
104   ierr = PetscViewerASCIIPrintf(viewer, "Quadrature on %d points\n  (", quad->numPoints);CHKERRQ(ierr);
105   for (q = 0; q < quad->numPoints; ++q) {
106     for (d = 0; d < quad->dim; ++d) {
107       if (d) ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr);
108       ierr = PetscViewerASCIIPrintf(viewer, "%g\n", (double)quad->points[q*quad->dim+d]);CHKERRQ(ierr);
109     }
110     ierr = PetscViewerASCIIPrintf(viewer, ") %g\n", (double)quad->weights[q]);CHKERRQ(ierr);
111   }
112   PetscFunctionReturn(0);
113 }
114 
115 #undef __FUNCT__
116 #define __FUNCT__ "PetscDTLegendreEval"
117 /*@
118    PetscDTLegendreEval - evaluate Legendre polynomial at points
119 
120    Not Collective
121 
122    Input Arguments:
123 +  npoints - number of spatial points to evaluate at
124 .  points - array of locations to evaluate at
125 .  ndegree - number of basis degrees to evaluate
126 -  degrees - sorted array of degrees to evaluate
127 
128    Output Arguments:
129 +  B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL)
130 .  D - row-oriented derivative evaluation matrix (or NULL)
131 -  D2 - row-oriented second derivative evaluation matrix (or NULL)
132 
133    Level: intermediate
134 
135 .seealso: PetscDTGaussQuadrature()
136 @*/
137 PetscErrorCode PetscDTLegendreEval(PetscInt npoints,const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2)
138 {
139   PetscInt i,maxdegree;
140 
141   PetscFunctionBegin;
142   if (!npoints || !ndegree) PetscFunctionReturn(0);
143   maxdegree = degrees[ndegree-1];
144   for (i=0; i<npoints; i++) {
145     PetscReal pm1,pm2,pd1,pd2,pdd1,pdd2,x;
146     PetscInt  j,k;
147     x    = points[i];
148     pm2  = 0;
149     pm1  = 1;
150     pd2  = 0;
151     pd1  = 0;
152     pdd2 = 0;
153     pdd1 = 0;
154     k    = 0;
155     if (degrees[k] == 0) {
156       if (B) B[i*ndegree+k] = pm1;
157       if (D) D[i*ndegree+k] = pd1;
158       if (D2) D2[i*ndegree+k] = pdd1;
159       k++;
160     }
161     for (j=1; j<=maxdegree; j++,k++) {
162       PetscReal p,d,dd;
163       p    = ((2*j-1)*x*pm1 - (j-1)*pm2)/j;
164       d    = pd2 + (2*j-1)*pm1;
165       dd   = pdd2 + (2*j-1)*pd1;
166       pm2  = pm1;
167       pm1  = p;
168       pd2  = pd1;
169       pd1  = d;
170       pdd2 = pdd1;
171       pdd1 = dd;
172       if (degrees[k] == j) {
173         if (B) B[i*ndegree+k] = p;
174         if (D) D[i*ndegree+k] = d;
175         if (D2) D2[i*ndegree+k] = dd;
176       }
177     }
178   }
179   PetscFunctionReturn(0);
180 }
181 
182 #undef __FUNCT__
183 #define __FUNCT__ "PetscDTGaussQuadrature"
184 /*@
185    PetscDTGaussQuadrature - create Gauss quadrature
186 
187    Not Collective
188 
189    Input Arguments:
190 +  npoints - number of points
191 .  a - left end of interval (often-1)
192 -  b - right end of interval (often +1)
193 
194    Output Arguments:
195 +  x - quadrature points
196 -  w - quadrature weights
197 
198    Level: intermediate
199 
200    References:
201    Golub and Welsch, Calculation of Quadrature Rules, Math. Comp. 23(106), 221--230, 1969.
202 
203 .seealso: PetscDTLegendreEval()
204 @*/
205 PetscErrorCode PetscDTGaussQuadrature(PetscInt npoints,PetscReal a,PetscReal b,PetscReal *x,PetscReal *w)
206 {
207   PetscErrorCode ierr;
208   PetscInt       i;
209   PetscReal      *work;
210   PetscScalar    *Z;
211   PetscBLASInt   N,LDZ,info;
212 
213   PetscFunctionBegin;
214   /* Set up the Golub-Welsch system */
215   for (i=0; i<npoints; i++) {
216     x[i] = 0;                   /* diagonal is 0 */
217     if (i) w[i-1] = 0.5 / PetscSqrtReal(1 - 1./PetscSqr(2*i));
218   }
219   ierr = PetscMalloc2(npoints*npoints,&Z,PetscMax(1,2*npoints-2),&work);CHKERRQ(ierr);
220   ierr = PetscBLASIntCast(npoints,&N);CHKERRQ(ierr);
221   LDZ  = N;
222   ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
223   PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("I",&N,x,w,Z,&LDZ,work,&info));
224   ierr = PetscFPTrapPop();CHKERRQ(ierr);
225   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEQR error");
226 
227   for (i=0; i<(npoints+1)/2; i++) {
228     PetscReal y = 0.5 * (-x[i] + x[npoints-i-1]); /* enforces symmetry */
229     x[i]           = (a+b)/2 - y*(b-a)/2;
230     x[npoints-i-1] = (a+b)/2 + y*(b-a)/2;
231 
232     w[i] = w[npoints-1-i] = (b-a)*PetscSqr(0.5*PetscAbsScalar(Z[i*npoints] + Z[(npoints-i-1)*npoints]));
233   }
234   ierr = PetscFree2(Z,work);CHKERRQ(ierr);
235   PetscFunctionReturn(0);
236 }
237 
238 #undef __FUNCT__
239 #define __FUNCT__ "PetscDTFactorial_Internal"
240 /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x.
241    Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */
242 PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorial_Internal(PetscInt n, PetscReal *factorial)
243 {
244   PetscReal f = 1.0;
245   PetscInt  i;
246 
247   PetscFunctionBegin;
248   for (i = 1; i < n+1; ++i) f *= i;
249   *factorial = f;
250   PetscFunctionReturn(0);
251 }
252 
253 #undef __FUNCT__
254 #define __FUNCT__ "PetscDTComputeJacobi"
255 /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x.
256    Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */
257 PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobi(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P)
258 {
259   PetscReal apb, pn1, pn2;
260   PetscInt  k;
261 
262   PetscFunctionBegin;
263   if (!n) {*P = 1.0; PetscFunctionReturn(0);}
264   if (n == 1) {*P = 0.5 * (a - b + (a + b + 2.0) * x); PetscFunctionReturn(0);}
265   apb = a + b;
266   pn2 = 1.0;
267   pn1 = 0.5 * (a - b + (apb + 2.0) * x);
268   *P  = 0.0;
269   for (k = 2; k < n+1; ++k) {
270     PetscReal a1 = 2.0 * k * (k + apb) * (2.0*k + apb - 2.0);
271     PetscReal a2 = (2.0 * k + apb - 1.0) * (a*a - b*b);
272     PetscReal a3 = (2.0 * k + apb - 2.0) * (2.0 * k + apb - 1.0) * (2.0 * k + apb);
273     PetscReal a4 = 2.0 * (k + a - 1.0) * (k + b - 1.0) * (2.0 * k + apb);
274 
275     a2  = a2 / a1;
276     a3  = a3 / a1;
277     a4  = a4 / a1;
278     *P  = (a2 + a3 * x) * pn1 - a4 * pn2;
279     pn2 = pn1;
280     pn1 = *P;
281   }
282   PetscFunctionReturn(0);
283 }
284 
285 #undef __FUNCT__
286 #define __FUNCT__ "PetscDTComputeJacobiDerivative"
287 /* Evaluates the first derivative of P_{n}^{a,b} at a point x. */
288 PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobiDerivative(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P)
289 {
290   PetscReal      nP;
291   PetscErrorCode ierr;
292 
293   PetscFunctionBegin;
294   if (!n) {*P = 0.0; PetscFunctionReturn(0);}
295   ierr = PetscDTComputeJacobi(a+1, b+1, n-1, x, &nP);CHKERRQ(ierr);
296   *P   = 0.5 * (a + b + n + 1) * nP;
297   PetscFunctionReturn(0);
298 }
299 
300 #undef __FUNCT__
301 #define __FUNCT__ "PetscDTMapSquareToTriangle_Internal"
302 /* Maps from [-1,1]^2 to the (-1,1) reference triangle */
303 PETSC_STATIC_INLINE PetscErrorCode PetscDTMapSquareToTriangle_Internal(PetscReal x, PetscReal y, PetscReal *xi, PetscReal *eta)
304 {
305   PetscFunctionBegin;
306   *xi  = 0.5 * (1.0 + x) * (1.0 - y) - 1.0;
307   *eta = y;
308   PetscFunctionReturn(0);
309 }
310 
311 #undef __FUNCT__
312 #define __FUNCT__ "PetscDTMapCubeToTetrahedron_Internal"
313 /* Maps from [-1,1]^2 to the (-1,1) reference triangle */
314 PETSC_STATIC_INLINE PetscErrorCode PetscDTMapCubeToTetrahedron_Internal(PetscReal x, PetscReal y, PetscReal z, PetscReal *xi, PetscReal *eta, PetscReal *zeta)
315 {
316   PetscFunctionBegin;
317   *xi   = 0.25 * (1.0 + x) * (1.0 - y) * (1.0 - z) - 1.0;
318   *eta  = 0.5  * (1.0 + y) * (1.0 - z) - 1.0;
319   *zeta = z;
320   PetscFunctionReturn(0);
321 }
322 
323 #undef __FUNCT__
324 #define __FUNCT__ "PetscDTGaussJacobiQuadrature1D_Internal"
325 static PetscErrorCode PetscDTGaussJacobiQuadrature1D_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal *x, PetscReal *w)
326 {
327   PetscInt       maxIter = 100;
328   PetscReal      eps     = 1.0e-8;
329   PetscReal      a1, a2, a3, a4, a5, a6;
330   PetscInt       k;
331   PetscErrorCode ierr;
332 
333   PetscFunctionBegin;
334 
335   a1      = PetscPowReal(2.0, a+b+1);
336 #if defined(PETSC_HAVE_TGAMMA)
337   a2      = PetscTGamma(a + npoints + 1);
338   a3      = PetscTGamma(b + npoints + 1);
339   a4      = PetscTGamma(a + b + npoints + 1);
340 #else
341   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"tgamma() - math routine is unavailable.");
342 #endif
343 
344   ierr = PetscDTFactorial_Internal(npoints, &a5);CHKERRQ(ierr);
345   a6   = a1 * a2 * a3 / a4 / a5;
346   /* Computes the m roots of P_{m}^{a,b} on [-1,1] by Newton's method with Chebyshev points as initial guesses.
347    Algorithm implemented from the pseudocode given by Karniadakis and Sherwin and Python in FIAT */
348   for (k = 0; k < npoints; ++k) {
349     PetscReal r = -PetscCosReal((2.0*k + 1.0) * PETSC_PI / (2.0 * npoints)), dP;
350     PetscInt  j;
351 
352     if (k > 0) r = 0.5 * (r + x[k-1]);
353     for (j = 0; j < maxIter; ++j) {
354       PetscReal s = 0.0, delta, f, fp;
355       PetscInt  i;
356 
357       for (i = 0; i < k; ++i) s = s + 1.0 / (r - x[i]);
358       ierr = PetscDTComputeJacobi(a, b, npoints, r, &f);CHKERRQ(ierr);
359       ierr = PetscDTComputeJacobiDerivative(a, b, npoints, r, &fp);CHKERRQ(ierr);
360       delta = f / (fp - f * s);
361       r     = r - delta;
362       if (PetscAbs(delta) < eps) break;
363     }
364     x[k] = r;
365     ierr = PetscDTComputeJacobiDerivative(a, b, npoints, x[k], &dP);CHKERRQ(ierr);
366     w[k] = a6 / (1.0 - PetscSqr(x[k])) / PetscSqr(dP);
367   }
368   PetscFunctionReturn(0);
369 }
370 
371 #undef __FUNCT__
372 #define __FUNCT__ "PetscDTGaussJacobiQuadrature"
373 /*@C
374   PetscDTGaussJacobiQuadrature - create Gauss-Jacobi quadrature for a simplex
375 
376   Not Collective
377 
378   Input Arguments:
379 + dim - The simplex dimension
380 . order - The quadrature order
381 . a - left end of interval (often-1)
382 - b - right end of interval (often +1)
383 
384   Output Arguments:
385 . q - A PetscQuadrature object
386 
387   Level: intermediate
388 
389   References:
390   Karniadakis and Sherwin.
391   FIAT
392 
393 .seealso: PetscDTGaussQuadrature()
394 @*/
395 PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt dim, PetscInt order, PetscReal a, PetscReal b, PetscQuadrature *q)
396 {
397   PetscInt       npoints = dim > 1 ? dim > 2 ? order*PetscSqr(order) : PetscSqr(order) : order;
398   PetscReal     *px, *wx, *py, *wy, *pz, *wz, *x, *w;
399   PetscInt       i, j, k;
400   PetscErrorCode ierr;
401 
402   PetscFunctionBegin;
403   if ((a != -1.0) || (b != 1.0)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must use default internal right now");
404   ierr = PetscMalloc1(npoints*dim, &x);CHKERRQ(ierr);
405   ierr = PetscMalloc1(npoints, &w);CHKERRQ(ierr);
406   switch (dim) {
407   case 0:
408     ierr = PetscFree(x);CHKERRQ(ierr);
409     ierr = PetscFree(w);CHKERRQ(ierr);
410     ierr = PetscMalloc1(1, &x);CHKERRQ(ierr);
411     ierr = PetscMalloc1(1, &w);CHKERRQ(ierr);
412     x[0] = 0.0;
413     w[0] = 1.0;
414     break;
415   case 1:
416     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, x, w);CHKERRQ(ierr);
417     break;
418   case 2:
419     ierr = PetscMalloc4(order,&px,order,&wx,order,&py,order,&wy);CHKERRQ(ierr);
420     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, px, wx);CHKERRQ(ierr);
421     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 1.0, 0.0, py, wy);CHKERRQ(ierr);
422     for (i = 0; i < order; ++i) {
423       for (j = 0; j < order; ++j) {
424         ierr = PetscDTMapSquareToTriangle_Internal(px[i], py[j], &x[(i*order+j)*2+0], &x[(i*order+j)*2+1]);CHKERRQ(ierr);
425         w[i*order+j] = 0.5 * wx[i] * wy[j];
426       }
427     }
428     ierr = PetscFree4(px,wx,py,wy);CHKERRQ(ierr);
429     break;
430   case 3:
431     ierr = PetscMalloc6(order,&px,order,&wx,order,&py,order,&wy,order,&pz,order,&wz);CHKERRQ(ierr);
432     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, px, wx);CHKERRQ(ierr);
433     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 1.0, 0.0, py, wy);CHKERRQ(ierr);
434     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 2.0, 0.0, pz, wz);CHKERRQ(ierr);
435     for (i = 0; i < order; ++i) {
436       for (j = 0; j < order; ++j) {
437         for (k = 0; k < order; ++k) {
438           ierr = PetscDTMapCubeToTetrahedron_Internal(px[i], py[j], pz[k], &x[((i*order+j)*order+k)*3+0], &x[((i*order+j)*order+k)*3+1], &x[((i*order+j)*order+k)*3+2]);CHKERRQ(ierr);
439           w[(i*order+j)*order+k] = 0.125 * wx[i] * wy[j] * wz[k];
440         }
441       }
442     }
443     ierr = PetscFree6(px,wx,py,wy,pz,wz);CHKERRQ(ierr);
444     break;
445   default:
446     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim);
447   }
448   ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr);
449   ierr = PetscQuadratureSetData(*q, dim, npoints, x, w);CHKERRQ(ierr);
450   PetscFunctionReturn(0);
451 }
452 
453 #undef __FUNCT__
454 #define __FUNCT__ "PetscDTPseudoInverseQR"
455 /* Overwrites A. Can only handle full-rank problems with m>=n
456  * A in column-major format
457  * Ainv in row-major format
458  * tau has length m
459  * worksize must be >= max(1,n)
460  */
461 static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m,PetscInt mstride,PetscInt n,PetscReal *A_in,PetscReal *Ainv_out,PetscScalar *tau,PetscInt worksize,PetscScalar *work)
462 {
463   PetscErrorCode ierr;
464   PetscBLASInt   M,N,K,lda,ldb,ldwork,info;
465   PetscScalar    *A,*Ainv,*R,*Q,Alpha;
466 
467   PetscFunctionBegin;
468 #if defined(PETSC_USE_COMPLEX)
469   {
470     PetscInt i,j;
471     ierr = PetscMalloc2(m*n,&A,m*n,&Ainv);CHKERRQ(ierr);
472     for (j=0; j<n; j++) {
473       for (i=0; i<m; i++) A[i+m*j] = A_in[i+mstride*j];
474     }
475     mstride = m;
476   }
477 #else
478   A = A_in;
479   Ainv = Ainv_out;
480 #endif
481 
482   ierr = PetscBLASIntCast(m,&M);CHKERRQ(ierr);
483   ierr = PetscBLASIntCast(n,&N);CHKERRQ(ierr);
484   ierr = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr);
485   ierr = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr);
486   ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
487   PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info));
488   ierr = PetscFPTrapPop();CHKERRQ(ierr);
489   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error");
490   R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */
491 
492   /* Extract an explicit representation of Q */
493   Q = Ainv;
494   ierr = PetscMemcpy(Q,A,mstride*n*sizeof(PetscScalar));CHKERRQ(ierr);
495   K = N;                        /* full rank */
496   PetscStackCallBLAS("LAPACKungqr",LAPACKungqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info));
497   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error");
498 
499   /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */
500   Alpha = 1.0;
501   ldb = lda;
502   PetscStackCallBLAS("BLAStrsm",BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb));
503   /* Ainv is Q, overwritten with inverse */
504 
505 #if defined(PETSC_USE_COMPLEX)
506   {
507     PetscInt i;
508     for (i=0; i<m*n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]);
509     ierr = PetscFree2(A,Ainv);CHKERRQ(ierr);
510   }
511 #endif
512   PetscFunctionReturn(0);
513 }
514 
515 #undef __FUNCT__
516 #define __FUNCT__ "PetscDTLegendreIntegrate"
517 /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */
518 static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval,const PetscReal *x,PetscInt ndegree,const PetscInt *degrees,PetscBool Transpose,PetscReal *B)
519 {
520   PetscErrorCode ierr;
521   PetscReal      *Bv;
522   PetscInt       i,j;
523 
524   PetscFunctionBegin;
525   ierr = PetscMalloc1((ninterval+1)*ndegree,&Bv);CHKERRQ(ierr);
526   /* Point evaluation of L_p on all the source vertices */
527   ierr = PetscDTLegendreEval(ninterval+1,x,ndegree,degrees,Bv,NULL,NULL);CHKERRQ(ierr);
528   /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */
529   for (i=0; i<ninterval; i++) {
530     for (j=0; j<ndegree; j++) {
531       if (Transpose) B[i+ninterval*j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j];
532       else           B[i*ndegree+j]   = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j];
533     }
534   }
535   ierr = PetscFree(Bv);CHKERRQ(ierr);
536   PetscFunctionReturn(0);
537 }
538 
539 #undef __FUNCT__
540 #define __FUNCT__ "PetscDTReconstructPoly"
541 /*@
542    PetscDTReconstructPoly - create matrix representing polynomial reconstruction using cell intervals and evaluation at target intervals
543 
544    Not Collective
545 
546    Input Arguments:
547 +  degree - degree of reconstruction polynomial
548 .  nsource - number of source intervals
549 .  sourcex - sorted coordinates of source cell boundaries (length nsource+1)
550 .  ntarget - number of target intervals
551 -  targetx - sorted coordinates of target cell boundaries (length ntarget+1)
552 
553    Output Arguments:
554 .  R - reconstruction matrix, utarget = sum_s R[t*nsource+s] * usource[s]
555 
556    Level: advanced
557 
558 .seealso: PetscDTLegendreEval()
559 @*/
560 PetscErrorCode PetscDTReconstructPoly(PetscInt degree,PetscInt nsource,const PetscReal *sourcex,PetscInt ntarget,const PetscReal *targetx,PetscReal *R)
561 {
562   PetscErrorCode ierr;
563   PetscInt       i,j,k,*bdegrees,worksize;
564   PetscReal      xmin,xmax,center,hscale,*sourcey,*targety,*Bsource,*Bsinv,*Btarget;
565   PetscScalar    *tau,*work;
566 
567   PetscFunctionBegin;
568   PetscValidRealPointer(sourcex,3);
569   PetscValidRealPointer(targetx,5);
570   PetscValidRealPointer(R,6);
571   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);
572 #if defined(PETSC_USE_DEBUG)
573   for (i=0; i<nsource; i++) {
574     if (sourcex[i] >= sourcex[i+1]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Source interval %D has negative orientation (%g,%g)",i,(double)sourcex[i],(double)sourcex[i+1]);
575   }
576   for (i=0; i<ntarget; i++) {
577     if (targetx[i] >= targetx[i+1]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Target interval %D has negative orientation (%g,%g)",i,(double)targetx[i],(double)targetx[i+1]);
578   }
579 #endif
580   xmin = PetscMin(sourcex[0],targetx[0]);
581   xmax = PetscMax(sourcex[nsource],targetx[ntarget]);
582   center = (xmin + xmax)/2;
583   hscale = (xmax - xmin)/2;
584   worksize = nsource;
585   ierr = PetscMalloc4(degree+1,&bdegrees,nsource+1,&sourcey,nsource*(degree+1),&Bsource,worksize,&work);CHKERRQ(ierr);
586   ierr = PetscMalloc4(nsource,&tau,nsource*(degree+1),&Bsinv,ntarget+1,&targety,ntarget*(degree+1),&Btarget);CHKERRQ(ierr);
587   for (i=0; i<=nsource; i++) sourcey[i] = (sourcex[i]-center)/hscale;
588   for (i=0; i<=degree; i++) bdegrees[i] = i+1;
589   ierr = PetscDTLegendreIntegrate(nsource,sourcey,degree+1,bdegrees,PETSC_TRUE,Bsource);CHKERRQ(ierr);
590   ierr = PetscDTPseudoInverseQR(nsource,nsource,degree+1,Bsource,Bsinv,tau,nsource,work);CHKERRQ(ierr);
591   for (i=0; i<=ntarget; i++) targety[i] = (targetx[i]-center)/hscale;
592   ierr = PetscDTLegendreIntegrate(ntarget,targety,degree+1,bdegrees,PETSC_FALSE,Btarget);CHKERRQ(ierr);
593   for (i=0; i<ntarget; i++) {
594     PetscReal rowsum = 0;
595     for (j=0; j<nsource; j++) {
596       PetscReal sum = 0;
597       for (k=0; k<degree+1; k++) {
598         sum += Btarget[i*(degree+1)+k] * Bsinv[k*nsource+j];
599       }
600       R[i*nsource+j] = sum;
601       rowsum += sum;
602     }
603     for (j=0; j<nsource; j++) R[i*nsource+j] /= rowsum; /* normalize each row */
604   }
605   ierr = PetscFree4(bdegrees,sourcey,Bsource,work);CHKERRQ(ierr);
606   ierr = PetscFree4(tau,Bsinv,targety,Btarget);CHKERRQ(ierr);
607   PetscFunctionReturn(0);
608 }
609