xref: /petsc/src/dm/dt/interface/dt.c (revision decb47aa5ac6684705207c1cbd56e8b9ae5b704e)
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", quad->points[q*quad->dim+d]);CHKERRQ(ierr);
109     }
110     ierr = PetscViewerASCIIPrintf(viewer, ") %g\n", 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 = PetscRealView(npoints-1,w,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
220   ierr = PetscMalloc2(npoints*npoints,&Z,PetscMax(1,2*npoints-2),&work);CHKERRQ(ierr);
221   ierr = PetscBLASIntCast(npoints,&N);CHKERRQ(ierr);
222   LDZ  = N;
223   ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
224   PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("I",&N,x,w,Z,&LDZ,work,&info));
225   ierr = PetscFPTrapPop();CHKERRQ(ierr);
226   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEQR error");
227 
228   for (i=0; i<(npoints+1)/2; i++) {
229     PetscReal y = 0.5 * (-x[i] + x[npoints-i-1]); /* enforces symmetry */
230     x[i]           = (a+b)/2 - y*(b-a)/2;
231     x[npoints-i-1] = (a+b)/2 + y*(b-a)/2;
232 
233     w[i] = w[npoints-1-i] = (b-a)*PetscSqr(0.5*PetscAbsScalar(Z[i*npoints] + Z[(npoints-i-1)*npoints]));
234   }
235   ierr = PetscFree2(Z,work);CHKERRQ(ierr);
236   PetscFunctionReturn(0);
237 }
238 
239 #undef __FUNCT__
240 #define __FUNCT__ "PetscDTFactorial_Internal"
241 /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x.
242    Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */
243 PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorial_Internal(PetscInt n, PetscReal *factorial)
244 {
245   PetscReal f = 1.0;
246   PetscInt  i;
247 
248   PetscFunctionBegin;
249   for (i = 1; i < n+1; ++i) f *= i;
250   *factorial = f;
251   PetscFunctionReturn(0);
252 }
253 
254 #undef __FUNCT__
255 #define __FUNCT__ "PetscDTComputeJacobi"
256 /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x.
257    Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */
258 PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobi(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P)
259 {
260   PetscReal apb, pn1, pn2;
261   PetscInt  k;
262 
263   PetscFunctionBegin;
264   if (!n) {*P = 1.0; PetscFunctionReturn(0);}
265   if (n == 1) {*P = 0.5 * (a - b + (a + b + 2.0) * x); PetscFunctionReturn(0);}
266   apb = a + b;
267   pn2 = 1.0;
268   pn1 = 0.5 * (a - b + (apb + 2.0) * x);
269   *P  = 0.0;
270   for (k = 2; k < n+1; ++k) {
271     PetscReal a1 = 2.0 * k * (k + apb) * (2.0*k + apb - 2.0);
272     PetscReal a2 = (2.0 * k + apb - 1.0) * (a*a - b*b);
273     PetscReal a3 = (2.0 * k + apb - 2.0) * (2.0 * k + apb - 1.0) * (2.0 * k + apb);
274     PetscReal a4 = 2.0 * (k + a - 1.0) * (k + b - 1.0) * (2.0 * k + apb);
275 
276     a2  = a2 / a1;
277     a3  = a3 / a1;
278     a4  = a4 / a1;
279     *P  = (a2 + a3 * x) * pn1 - a4 * pn2;
280     pn2 = pn1;
281     pn1 = *P;
282   }
283   PetscFunctionReturn(0);
284 }
285 
286 #undef __FUNCT__
287 #define __FUNCT__ "PetscDTComputeJacobiDerivative"
288 /* Evaluates the first derivative of P_{n}^{a,b} at a point x. */
289 PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobiDerivative(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P)
290 {
291   PetscReal      nP;
292   PetscErrorCode ierr;
293 
294   PetscFunctionBegin;
295   if (!n) {*P = 0.0; PetscFunctionReturn(0);}
296   ierr = PetscDTComputeJacobi(a+1, b+1, n-1, x, &nP);CHKERRQ(ierr);
297   *P   = 0.5 * (a + b + n + 1) * nP;
298   PetscFunctionReturn(0);
299 }
300 
301 #undef __FUNCT__
302 #define __FUNCT__ "PetscDTMapSquareToTriangle_Internal"
303 /* Maps from [-1,1]^2 to the (-1,1) reference triangle */
304 PETSC_STATIC_INLINE PetscErrorCode PetscDTMapSquareToTriangle_Internal(PetscReal x, PetscReal y, PetscReal *xi, PetscReal *eta)
305 {
306   PetscFunctionBegin;
307   *xi  = 0.5 * (1.0 + x) * (1.0 - y) - 1.0;
308   *eta = y;
309   PetscFunctionReturn(0);
310 }
311 
312 #undef __FUNCT__
313 #define __FUNCT__ "PetscDTMapCubeToTetrahedron_Internal"
314 /* Maps from [-1,1]^2 to the (-1,1) reference triangle */
315 PETSC_STATIC_INLINE PetscErrorCode PetscDTMapCubeToTetrahedron_Internal(PetscReal x, PetscReal y, PetscReal z, PetscReal *xi, PetscReal *eta, PetscReal *zeta)
316 {
317   PetscFunctionBegin;
318   *xi   = 0.25 * (1.0 + x) * (1.0 - y) * (1.0 - z) - 1.0;
319   *eta  = 0.5  * (1.0 + y) * (1.0 - z) - 1.0;
320   *zeta = z;
321   PetscFunctionReturn(0);
322 }
323 
324 #undef __FUNCT__
325 #define __FUNCT__ "PetscDTGaussJacobiQuadrature1D_Internal"
326 static PetscErrorCode PetscDTGaussJacobiQuadrature1D_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal *x, PetscReal *w)
327 {
328   PetscInt       maxIter = 100;
329   PetscReal      eps     = 1.0e-8;
330   PetscReal      a1, a2, a3, a4, a5, a6;
331   PetscInt       k;
332   PetscErrorCode ierr;
333 
334   PetscFunctionBegin;
335 
336   a1      = PetscPowReal(2.0, a+b+1);
337 #if defined(PETSC_HAVE_TGAMMA)
338   a2      = tgamma(a + npoints + 1);
339   a3      = tgamma(b + npoints + 1);
340   a4      = tgamma(a + b + npoints + 1);
341 #else
342   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"tgamma() - math routine is unavailable.");
343 #endif
344 
345   ierr = PetscDTFactorial_Internal(npoints, &a5);CHKERRQ(ierr);
346   a6   = a1 * a2 * a3 / a4 / a5;
347   /* Computes the m roots of P_{m}^{a,b} on [-1,1] by Newton's method with Chebyshev points as initial guesses.
348    Algorithm implemented from the pseudocode given by Karniadakis and Sherwin and Python in FIAT */
349   for (k = 0; k < npoints; ++k) {
350     PetscReal r = -PetscCosReal((2.0*k + 1.0) * PETSC_PI / (2.0 * npoints)), dP;
351     PetscInt  j;
352 
353     if (k > 0) r = 0.5 * (r + x[k-1]);
354     for (j = 0; j < maxIter; ++j) {
355       PetscReal s = 0.0, delta, f, fp;
356       PetscInt  i;
357 
358       for (i = 0; i < k; ++i) s = s + 1.0 / (r - x[i]);
359       ierr = PetscDTComputeJacobi(a, b, npoints, r, &f);CHKERRQ(ierr);
360       ierr = PetscDTComputeJacobiDerivative(a, b, npoints, r, &fp);CHKERRQ(ierr);
361       delta = f / (fp - f * s);
362       r     = r - delta;
363       if (fabs(delta) < eps) break;
364     }
365     x[k] = r;
366     ierr = PetscDTComputeJacobiDerivative(a, b, npoints, x[k], &dP);CHKERRQ(ierr);
367     w[k] = a6 / (1.0 - PetscSqr(x[k])) / PetscSqr(dP);
368   }
369   PetscFunctionReturn(0);
370 }
371 
372 #undef __FUNCT__
373 #define __FUNCT__ "PetscDTGaussJacobiQuadrature"
374 /*@C
375   PetscDTGaussJacobiQuadrature - create Gauss-Jacobi quadrature for a simplex
376 
377   Not Collective
378 
379   Input Arguments:
380 + dim - The simplex dimension
381 . order - The quadrature order
382 . a - left end of interval (often-1)
383 - b - right end of interval (often +1)
384 
385   Output Arguments:
386 . q - A PetscQuadrature object
387 
388   Level: intermediate
389 
390   References:
391   Karniadakis and Sherwin.
392   FIAT
393 
394 .seealso: PetscDTGaussQuadrature()
395 @*/
396 PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt dim, PetscInt order, PetscReal a, PetscReal b, PetscQuadrature *q)
397 {
398   PetscInt       npoints = dim > 1 ? dim > 2 ? order*PetscSqr(order) : PetscSqr(order) : order;
399   PetscReal     *px, *wx, *py, *wy, *pz, *wz, *x, *w;
400   PetscInt       i, j, k;
401   PetscErrorCode ierr;
402 
403   PetscFunctionBegin;
404   if ((a != -1.0) || (b != 1.0)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must use default internal right now");
405   ierr = PetscMalloc1(npoints*dim, &x);CHKERRQ(ierr);
406   ierr = PetscMalloc1(npoints, &w);CHKERRQ(ierr);
407   switch (dim) {
408   case 0:
409     ierr = PetscFree(x);CHKERRQ(ierr);
410     ierr = PetscFree(w);CHKERRQ(ierr);
411     ierr = PetscMalloc1(1, &x);CHKERRQ(ierr);
412     ierr = PetscMalloc1(1, &w);CHKERRQ(ierr);
413     x[0] = 0.0;
414     w[0] = 1.0;
415     break;
416   case 1:
417     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, x, w);CHKERRQ(ierr);
418     break;
419   case 2:
420     ierr = PetscMalloc4(order,&px,order,&wx,order,&py,order,&wy);CHKERRQ(ierr);
421     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, px, wx);CHKERRQ(ierr);
422     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 1.0, 0.0, py, wy);CHKERRQ(ierr);
423     for (i = 0; i < order; ++i) {
424       for (j = 0; j < order; ++j) {
425         ierr = PetscDTMapSquareToTriangle_Internal(px[i], py[j], &x[(i*order+j)*2+0], &x[(i*order+j)*2+1]);CHKERRQ(ierr);
426         w[i*order+j] = 0.5 * wx[i] * wy[j];
427       }
428     }
429     ierr = PetscFree4(px,wx,py,wy);CHKERRQ(ierr);
430     break;
431   case 3:
432     ierr = PetscMalloc6(order,&px,order,&wx,order,&py,order,&wy,order,&pz,order,&wz);CHKERRQ(ierr);
433     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, px, wx);CHKERRQ(ierr);
434     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 1.0, 0.0, py, wy);CHKERRQ(ierr);
435     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 2.0, 0.0, pz, wz);CHKERRQ(ierr);
436     for (i = 0; i < order; ++i) {
437       for (j = 0; j < order; ++j) {
438         for (k = 0; k < order; ++k) {
439           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);
440           w[(i*order+j)*order+k] = 0.125 * wx[i] * wy[j] * wz[k];
441         }
442       }
443     }
444     ierr = PetscFree6(px,wx,py,wy,pz,wz);CHKERRQ(ierr);
445     break;
446   default:
447     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim);
448   }
449   ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr);
450   ierr = PetscQuadratureSetData(*q, dim, npoints, x, w);CHKERRQ(ierr);
451   PetscFunctionReturn(0);
452 }
453 
454 #undef __FUNCT__
455 #define __FUNCT__ "PetscDTPseudoInverseQR"
456 /* Overwrites A. Can only handle full-rank problems with m>=n
457  * A in column-major format
458  * Ainv in row-major format
459  * tau has length m
460  * worksize must be >= max(1,n)
461  */
462 static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m,PetscInt mstride,PetscInt n,PetscReal *A_in,PetscReal *Ainv_out,PetscScalar *tau,PetscInt worksize,PetscScalar *work)
463 {
464   PetscErrorCode ierr;
465   PetscBLASInt M,N,K,lda,ldb,ldwork,info;
466   PetscScalar *A,*Ainv,*R,*Q,Alpha;
467 
468   PetscFunctionBegin;
469 #if defined(PETSC_USE_COMPLEX)
470   {
471     PetscInt i,j;
472     ierr = PetscMalloc2(m*n,&A,m*n,&Ainv);CHKERRQ(ierr);
473     for (j=0; j<n; j++) {
474       for (i=0; i<m; i++) A[i+m*j] = A_in[i+mstride*j];
475     }
476     mstride = m;
477   }
478 #else
479   A = A_in;
480   Ainv = Ainv_out;
481 #endif
482 
483   ierr = PetscBLASIntCast(m,&M);CHKERRQ(ierr);
484   ierr = PetscBLASIntCast(n,&N);CHKERRQ(ierr);
485   ierr = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr);
486   ierr = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr);
487   ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
488   LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info);
489   ierr = PetscFPTrapPop();CHKERRQ(ierr);
490   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error");
491   R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */
492 
493   /* Extract an explicit representation of Q */
494   Q = Ainv;
495   ierr = PetscMemcpy(Q,A,mstride*n*sizeof(PetscScalar));CHKERRQ(ierr);
496   K = N;                        /* full rank */
497   LAPACKungqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info);
498   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error");
499 
500   /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */
501   Alpha = 1.0;
502   ldb = lda;
503   BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb);
504   /* Ainv is Q, overwritten with inverse */
505 
506 #if defined(PETSC_USE_COMPLEX)
507   {
508     PetscInt i;
509     for (i=0; i<m*n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]);
510     ierr = PetscFree2(A,Ainv);CHKERRQ(ierr);
511   }
512 #endif
513   PetscFunctionReturn(0);
514 }
515 
516 #undef __FUNCT__
517 #define __FUNCT__ "PetscDTLegendreIntegrate"
518 /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */
519 static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval,const PetscReal *x,PetscInt ndegree,const PetscInt *degrees,PetscBool Transpose,PetscReal *B)
520 {
521   PetscErrorCode ierr;
522   PetscReal *Bv;
523   PetscInt i,j;
524 
525   PetscFunctionBegin;
526   ierr = PetscMalloc1((ninterval+1)*ndegree,&Bv);CHKERRQ(ierr);
527   /* Point evaluation of L_p on all the source vertices */
528   ierr = PetscDTLegendreEval(ninterval+1,x,ndegree,degrees,Bv,NULL,NULL);CHKERRQ(ierr);
529   /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */
530   for (i=0; i<ninterval; i++) {
531     for (j=0; j<ndegree; j++) {
532       if (Transpose) B[i+ninterval*j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j];
533       else           B[i*ndegree+j]   = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j];
534     }
535   }
536   ierr = PetscFree(Bv);CHKERRQ(ierr);
537   PetscFunctionReturn(0);
538 }
539 
540 #undef __FUNCT__
541 #define __FUNCT__ "PetscDTReconstructPoly"
542 /*@
543    PetscDTReconstructPoly - create matrix representing polynomial reconstruction using cell intervals and evaluation at target intervals
544 
545    Not Collective
546 
547    Input Arguments:
548 +  degree - degree of reconstruction polynomial
549 .  nsource - number of source intervals
550 .  sourcex - sorted coordinates of source cell boundaries (length nsource+1)
551 .  ntarget - number of target intervals
552 -  targetx - sorted coordinates of target cell boundaries (length ntarget+1)
553 
554    Output Arguments:
555 .  R - reconstruction matrix, utarget = sum_s R[t*nsource+s] * usource[s]
556 
557    Level: advanced
558 
559 .seealso: PetscDTLegendreEval()
560 @*/
561 PetscErrorCode PetscDTReconstructPoly(PetscInt degree,PetscInt nsource,const PetscReal *sourcex,PetscInt ntarget,const PetscReal *targetx,PetscReal *R)
562 {
563   PetscErrorCode ierr;
564   PetscInt i,j,k,*bdegrees,worksize;
565   PetscReal xmin,xmax,center,hscale,*sourcey,*targety,*Bsource,*Bsinv,*Btarget;
566   PetscScalar *tau,*work;
567 
568   PetscFunctionBegin;
569   PetscValidRealPointer(sourcex,3);
570   PetscValidRealPointer(targetx,5);
571   PetscValidRealPointer(R,6);
572   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);
573 #if defined(PETSC_USE_DEBUG)
574   for (i=0; i<nsource; i++) {
575     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]);
576   }
577   for (i=0; i<ntarget; i++) {
578     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]);
579   }
580 #endif
581   xmin = PetscMin(sourcex[0],targetx[0]);
582   xmax = PetscMax(sourcex[nsource],targetx[ntarget]);
583   center = (xmin + xmax)/2;
584   hscale = (xmax - xmin)/2;
585   worksize = nsource;
586   ierr = PetscMalloc4(degree+1,&bdegrees,nsource+1,&sourcey,nsource*(degree+1),&Bsource,worksize,&work);CHKERRQ(ierr);
587   ierr = PetscMalloc4(nsource,&tau,nsource*(degree+1),&Bsinv,ntarget+1,&targety,ntarget*(degree+1),&Btarget);CHKERRQ(ierr);
588   for (i=0; i<=nsource; i++) sourcey[i] = (sourcex[i]-center)/hscale;
589   for (i=0; i<=degree; i++) bdegrees[i] = i+1;
590   ierr = PetscDTLegendreIntegrate(nsource,sourcey,degree+1,bdegrees,PETSC_TRUE,Bsource);CHKERRQ(ierr);
591   ierr = PetscDTPseudoInverseQR(nsource,nsource,degree+1,Bsource,Bsinv,tau,nsource,work);CHKERRQ(ierr);
592   for (i=0; i<=ntarget; i++) targety[i] = (targetx[i]-center)/hscale;
593   ierr = PetscDTLegendreIntegrate(ntarget,targety,degree+1,bdegrees,PETSC_FALSE,Btarget);CHKERRQ(ierr);
594   for (i=0; i<ntarget; i++) {
595     PetscReal rowsum = 0;
596     for (j=0; j<nsource; j++) {
597       PetscReal sum = 0;
598       for (k=0; k<degree+1; k++) {
599         sum += Btarget[i*(degree+1)+k] * Bsinv[k*nsource+j];
600       }
601       R[i*nsource+j] = sum;
602       rowsum += sum;
603     }
604     for (j=0; j<nsource; j++) R[i*nsource+j] /= rowsum; /* normalize each row */
605   }
606   ierr = PetscFree4(bdegrees,sourcey,Bsource,work);CHKERRQ(ierr);
607   ierr = PetscFree4(tau,Bsinv,targety,Btarget);CHKERRQ(ierr);
608   PetscFunctionReturn(0);
609 }
610