xref: /petsc/src/dm/dt/interface/dt.c (revision 83b1d2c3652270b8f524ddc204331e3bf2caac59)
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 static PetscBool GaussCite       = PETSC_FALSE;
17 const char       GaussCitation[] = "@article{GolubWelsch1969,\n"
18                                    "  author  = {Golub and Welsch},\n"
19                                    "  title   = {Calculation of Quadrature Rules},\n"
20                                    "  journal = {Math. Comp.},\n"
21                                    "  volume  = {23},\n"
22                                    "  number  = {106},\n"
23                                    "  pages   = {221--230},\n"
24                                    "  year    = {1969}\n}\n";
25 
26 #undef __FUNCT__
27 #define __FUNCT__ "PetscQuadratureCreate"
28 /*@
29   PetscQuadratureCreate - Create a PetscQuadrature object
30 
31   Collective on MPI_Comm
32 
33   Input Parameter:
34 . comm - The comm
35 
36   Output Parameter:
37 . q - The PetscQuadrature object
38 
39   Level: beginner
40 
41 .seealso: PetscQuadratureDestroy(), PetscQuadratureGetData(), PetscQuadratureSetData()
42 @*/
43 PetscErrorCode PetscQuadratureCreate(MPI_Comm comm, PetscQuadrature *q)
44 {
45   PetscErrorCode ierr;
46 
47   PetscFunctionBegin;
48   PetscValidPointer(q, 2);
49   ierr = DMInitializePackage();CHKERRQ(ierr);
50   ierr = PetscHeaderCreate(*q,_p_PetscQuadrature,int,PETSC_OBJECT_CLASSID,"PetscQuadrature","Quadrature","DT",comm,PetscQuadratureDestroy,PetscQuadratureView);CHKERRQ(ierr);
51   (*q)->dim       = -1;
52   (*q)->order     = -1;
53   (*q)->numPoints = 0;
54   (*q)->points    = NULL;
55   (*q)->weights   = NULL;
56   PetscFunctionReturn(0);
57 }
58 
59 #undef __FUNCT__
60 #define __FUNCT__ "PetscQuadratureDestroy"
61 /*@
62   PetscQuadratureDestroy - Destroy a PetscQuadrature object
63 
64   Collective on PetscQuadrature
65 
66   Input Parameter:
67 . q - The PetscQuadrature object
68 
69   Level: beginner
70 
71 .seealso: PetscQuadratureCreate(), PetscQuadratureGetData(), PetscQuadratureSetData()
72 @*/
73 PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *q)
74 {
75   PetscErrorCode ierr;
76 
77   PetscFunctionBegin;
78   if (!*q) PetscFunctionReturn(0);
79   PetscValidHeaderSpecific((*q),PETSC_OBJECT_CLASSID,1);
80   if (--((PetscObject)(*q))->refct > 0) {
81     *q = NULL;
82     PetscFunctionReturn(0);
83   }
84   ierr = PetscFree((*q)->points);CHKERRQ(ierr);
85   ierr = PetscFree((*q)->weights);CHKERRQ(ierr);
86   ierr = PetscHeaderDestroy(q);CHKERRQ(ierr);
87   PetscFunctionReturn(0);
88 }
89 
90 #undef __FUNCT__
91 #define __FUNCT__ "PetscQuadratureGetOrder"
92 /*@
93   PetscQuadratureGetOrder - Return the quadrature information
94 
95   Not collective
96 
97   Input Parameter:
98 . q - The PetscQuadrature object
99 
100   Output Parameter:
101 . order - The order of the quadrature, i.e. the highest degree polynomial that is exactly integrated
102 
103   Output Parameter:
104 
105   Level: intermediate
106 
107 .seealso: PetscQuadratureSetOrder(), PetscQuadratureGetData(), PetscQuadratureSetData()
108 @*/
109 PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature q, PetscInt *order)
110 {
111   PetscFunctionBegin;
112   PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1);
113   PetscValidPointer(order, 2);
114   *order = q->order;
115   PetscFunctionReturn(0);
116 }
117 
118 #undef __FUNCT__
119 #define __FUNCT__ "PetscQuadratureSetOrder"
120 /*@
121   PetscQuadratureSetOrder - Return the quadrature information
122 
123   Not collective
124 
125   Input Parameters:
126 + q - The PetscQuadrature object
127 - order - The order of the quadrature, i.e. the highest degree polynomial that is exactly integrated
128 
129   Level: intermediate
130 
131 .seealso: PetscQuadratureGetOrder(), PetscQuadratureGetData(), PetscQuadratureSetData()
132 @*/
133 PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature q, PetscInt order)
134 {
135   PetscFunctionBegin;
136   PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1);
137   q->order = order;
138   PetscFunctionReturn(0);
139 }
140 
141 #undef __FUNCT__
142 #define __FUNCT__ "PetscQuadratureGetData"
143 /*@C
144   PetscQuadratureGetData - Return the quadrature information
145 
146   Not collective
147 
148   Input Parameter:
149 . q - The PetscQuadrature object
150 
151   Output Parameters:
152 + dim     - The spatial dimension
153 . npoints - The number of quadrature points
154 . points  - The coordinates of the quadrature points
155 - weights - The quadrature weights
156 
157   Level: intermediate
158 
159 .seealso: PetscQuadratureSetData(), PetscQuadratureGetOrder(), PetscQuadratureSetOrder()
160 @*/
161 PetscErrorCode PetscQuadratureGetData(PetscQuadrature q, PetscInt *dim, PetscInt *npoints, const PetscReal *points[], const PetscReal *weights[])
162 {
163   PetscFunctionBegin;
164   PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1);
165   if (dim) {
166     PetscValidPointer(dim, 2);
167     *dim = q->dim;
168   }
169   if (npoints) {
170     PetscValidPointer(npoints, 3);
171     *npoints = q->numPoints;
172   }
173   if (points) {
174     PetscValidPointer(points, 4);
175     *points = q->points;
176   }
177   if (weights) {
178     PetscValidPointer(weights, 5);
179     *weights = q->weights;
180   }
181   PetscFunctionReturn(0);
182 }
183 
184 #undef __FUNCT__
185 #define __FUNCT__ "PetscQuadratureSetData"
186 /*@C
187   PetscQuadratureSetData - Set the quadrature information
188 
189   Not collective
190 
191   Input Parameters:
192 + q       - The PetscQuadrature object
193 . dim     - The spatial dimension
194 . npoints - The number of quadrature points
195 . points  - The coordinates of the quadrature points
196 - weights - The quadrature weights
197 
198   Level: intermediate
199 
200 .seealso: PetscQuadratureGetData(), PetscQuadratureGetOrder(), PetscQuadratureSetOrder()
201 @*/
202 PetscErrorCode PetscQuadratureSetData(PetscQuadrature q, PetscInt dim, PetscInt npoints, const PetscReal points[], const PetscReal weights[])
203 {
204   PetscFunctionBegin;
205   PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1);
206   if (dim >= 0)     q->dim       = dim;
207   if (npoints >= 0) q->numPoints = npoints;
208   if (points) {
209     PetscValidPointer(points, 4);
210     q->points = points;
211   }
212   if (weights) {
213     PetscValidPointer(weights, 5);
214     q->weights = weights;
215   }
216   PetscFunctionReturn(0);
217 }
218 
219 #undef __FUNCT__
220 #define __FUNCT__ "PetscQuadratureView"
221 /*@
222   PetscQuadratureView - View a PetscQuadrature object
223 
224   Collective on MPI_Comm
225 
226   Input Parameters:
227 + q - The PetscQuadrature object
228 - viewer - The PetscViewer object
229 
230   Level: beginner
231 
232 .seealso: PetscQuadratureCreate(), PetscQuadratureGetData(), PetscQuadratureSetData()
233 @*/
234 PetscErrorCode PetscQuadratureView(PetscQuadrature quad, PetscViewer viewer)
235 {
236   PetscInt       q, d;
237   PetscErrorCode ierr;
238 
239   PetscFunctionBegin;
240   ierr = PetscObjectPrintClassNamePrefixType((PetscObject)quad,viewer);CHKERRQ(ierr);
241   ierr = PetscViewerASCIIPrintf(viewer, "Quadrature on %d points\n  (", quad->numPoints);CHKERRQ(ierr);
242   for (q = 0; q < quad->numPoints; ++q) {
243     for (d = 0; d < quad->dim; ++d) {
244       if (d) ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr);
245       ierr = PetscViewerASCIIPrintf(viewer, "%g\n", (double)quad->points[q*quad->dim+d]);CHKERRQ(ierr);
246     }
247     ierr = PetscViewerASCIIPrintf(viewer, ") %g\n", (double)quad->weights[q]);CHKERRQ(ierr);
248   }
249   PetscFunctionReturn(0);
250 }
251 
252 #undef __FUNCT__
253 #define __FUNCT__ "PetscDTLegendreEval"
254 /*@
255    PetscDTLegendreEval - evaluate Legendre polynomial at points
256 
257    Not Collective
258 
259    Input Arguments:
260 +  npoints - number of spatial points to evaluate at
261 .  points - array of locations to evaluate at
262 .  ndegree - number of basis degrees to evaluate
263 -  degrees - sorted array of degrees to evaluate
264 
265    Output Arguments:
266 +  B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL)
267 .  D - row-oriented derivative evaluation matrix (or NULL)
268 -  D2 - row-oriented second derivative evaluation matrix (or NULL)
269 
270    Level: intermediate
271 
272 .seealso: PetscDTGaussQuadrature()
273 @*/
274 PetscErrorCode PetscDTLegendreEval(PetscInt npoints,const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2)
275 {
276   PetscInt i,maxdegree;
277 
278   PetscFunctionBegin;
279   if (!npoints || !ndegree) PetscFunctionReturn(0);
280   maxdegree = degrees[ndegree-1];
281   for (i=0; i<npoints; i++) {
282     PetscReal pm1,pm2,pd1,pd2,pdd1,pdd2,x;
283     PetscInt  j,k;
284     x    = points[i];
285     pm2  = 0;
286     pm1  = 1;
287     pd2  = 0;
288     pd1  = 0;
289     pdd2 = 0;
290     pdd1 = 0;
291     k    = 0;
292     if (degrees[k] == 0) {
293       if (B) B[i*ndegree+k] = pm1;
294       if (D) D[i*ndegree+k] = pd1;
295       if (D2) D2[i*ndegree+k] = pdd1;
296       k++;
297     }
298     for (j=1; j<=maxdegree; j++,k++) {
299       PetscReal p,d,dd;
300       p    = ((2*j-1)*x*pm1 - (j-1)*pm2)/j;
301       d    = pd2 + (2*j-1)*pm1;
302       dd   = pdd2 + (2*j-1)*pd1;
303       pm2  = pm1;
304       pm1  = p;
305       pd2  = pd1;
306       pd1  = d;
307       pdd2 = pdd1;
308       pdd1 = dd;
309       if (degrees[k] == j) {
310         if (B) B[i*ndegree+k] = p;
311         if (D) D[i*ndegree+k] = d;
312         if (D2) D2[i*ndegree+k] = dd;
313       }
314     }
315   }
316   PetscFunctionReturn(0);
317 }
318 
319 #undef __FUNCT__
320 #define __FUNCT__ "PetscDTGaussQuadrature"
321 /*@
322    PetscDTGaussQuadrature - create Gauss quadrature
323 
324    Not Collective
325 
326    Input Arguments:
327 +  npoints - number of points
328 .  a - left end of interval (often-1)
329 -  b - right end of interval (often +1)
330 
331    Output Arguments:
332 +  x - quadrature points
333 -  w - quadrature weights
334 
335    Level: intermediate
336 
337    References:
338    Golub and Welsch, Calculation of Quadrature Rules, Math. Comp. 23(106), 221--230, 1969.
339 
340 .seealso: PetscDTLegendreEval()
341 @*/
342 PetscErrorCode PetscDTGaussQuadrature(PetscInt npoints,PetscReal a,PetscReal b,PetscReal *x,PetscReal *w)
343 {
344   PetscErrorCode ierr;
345   PetscInt       i;
346   PetscReal      *work;
347   PetscScalar    *Z;
348   PetscBLASInt   N,LDZ,info;
349 
350   PetscFunctionBegin;
351   ierr = PetscCitationsRegister(GaussCitation, &GaussCite);CHKERRQ(ierr);
352   /* Set up the Golub-Welsch system */
353   for (i=0; i<npoints; i++) {
354     x[i] = 0;                   /* diagonal is 0 */
355     if (i) w[i-1] = 0.5 / PetscSqrtReal(1 - 1./PetscSqr(2*i));
356   }
357   ierr = PetscMalloc2(npoints*npoints,&Z,PetscMax(1,2*npoints-2),&work);CHKERRQ(ierr);
358   ierr = PetscBLASIntCast(npoints,&N);CHKERRQ(ierr);
359   LDZ  = N;
360   ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
361   PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("I",&N,x,w,Z,&LDZ,work,&info));
362   ierr = PetscFPTrapPop();CHKERRQ(ierr);
363   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEQR error");
364 
365   for (i=0; i<(npoints+1)/2; i++) {
366     PetscReal y = 0.5 * (-x[i] + x[npoints-i-1]); /* enforces symmetry */
367     x[i]           = (a+b)/2 - y*(b-a)/2;
368     x[npoints-i-1] = (a+b)/2 + y*(b-a)/2;
369 
370     w[i] = w[npoints-1-i] = 0.5*(b-a)*(PetscSqr(PetscAbsScalar(Z[i*npoints])) + PetscSqr(PetscAbsScalar(Z[(npoints-i-1)*npoints])));
371   }
372   ierr = PetscFree2(Z,work);CHKERRQ(ierr);
373   PetscFunctionReturn(0);
374 }
375 
376 #undef __FUNCT__
377 #define __FUNCT__ "PetscDTGaussTensorQuadrature"
378 /*@
379   PetscDTGaussTensorQuadrature - creates a tensor-product Gauss quadrature
380 
381   Not Collective
382 
383   Input Arguments:
384 + dim     - The spatial dimension
385 . npoints - number of points in one dimension
386 . a       - left end of interval (often-1)
387 - b       - right end of interval (often +1)
388 
389   Output Argument:
390 . q - A PetscQuadrature object
391 
392   Level: intermediate
393 
394 .seealso: PetscDTGaussQuadrature(), PetscDTLegendreEval()
395 @*/
396 PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt dim, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q)
397 {
398   PetscInt       totpoints = dim > 1 ? dim > 2 ? npoints*PetscSqr(npoints) : PetscSqr(npoints) : npoints, i, j, k;
399   PetscReal     *x, *w, *xw, *ww;
400   PetscErrorCode ierr;
401 
402   PetscFunctionBegin;
403   ierr = PetscMalloc1(totpoints*dim,&x);CHKERRQ(ierr);
404   ierr = PetscMalloc1(totpoints,&w);CHKERRQ(ierr);
405   /* Set up the Golub-Welsch system */
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 = PetscDTGaussQuadrature(npoints, a, b, x, w);CHKERRQ(ierr);
417     break;
418   case 2:
419     ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr);
420     ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr);
421     for (i = 0; i < npoints; ++i) {
422       for (j = 0; j < npoints; ++j) {
423         x[(i*npoints+j)*dim+0] = xw[i];
424         x[(i*npoints+j)*dim+1] = xw[j];
425         w[i*npoints+j]         = ww[i] * ww[j];
426       }
427     }
428     ierr = PetscFree2(xw,ww);CHKERRQ(ierr);
429     break;
430   case 3:
431     ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr);
432     ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr);
433     for (i = 0; i < npoints; ++i) {
434       for (j = 0; j < npoints; ++j) {
435         for (k = 0; k < npoints; ++k) {
436           x[((i*npoints+j)*npoints+k)*dim+0] = xw[i];
437           x[((i*npoints+j)*npoints+k)*dim+1] = xw[j];
438           x[((i*npoints+j)*npoints+k)*dim+2] = xw[k];
439           w[(i*npoints+j)*npoints+k]         = ww[i] * ww[j] * ww[k];
440         }
441       }
442     }
443     ierr = PetscFree2(xw,ww);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 = PetscQuadratureSetOrder(*q, npoints);CHKERRQ(ierr);
450   ierr = PetscQuadratureSetData(*q, dim, totpoints, x, w);CHKERRQ(ierr);
451   PetscFunctionReturn(0);
452 }
453 
454 #undef __FUNCT__
455 #define __FUNCT__ "PetscDTFactorial_Internal"
456 /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x.
457    Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */
458 PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorial_Internal(PetscInt n, PetscReal *factorial)
459 {
460   PetscReal f = 1.0;
461   PetscInt  i;
462 
463   PetscFunctionBegin;
464   for (i = 1; i < n+1; ++i) f *= i;
465   *factorial = f;
466   PetscFunctionReturn(0);
467 }
468 
469 #undef __FUNCT__
470 #define __FUNCT__ "PetscDTComputeJacobi"
471 /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x.
472    Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */
473 PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobi(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P)
474 {
475   PetscReal apb, pn1, pn2;
476   PetscInt  k;
477 
478   PetscFunctionBegin;
479   if (!n) {*P = 1.0; PetscFunctionReturn(0);}
480   if (n == 1) {*P = 0.5 * (a - b + (a + b + 2.0) * x); PetscFunctionReturn(0);}
481   apb = a + b;
482   pn2 = 1.0;
483   pn1 = 0.5 * (a - b + (apb + 2.0) * x);
484   *P  = 0.0;
485   for (k = 2; k < n+1; ++k) {
486     PetscReal a1 = 2.0 * k * (k + apb) * (2.0*k + apb - 2.0);
487     PetscReal a2 = (2.0 * k + apb - 1.0) * (a*a - b*b);
488     PetscReal a3 = (2.0 * k + apb - 2.0) * (2.0 * k + apb - 1.0) * (2.0 * k + apb);
489     PetscReal a4 = 2.0 * (k + a - 1.0) * (k + b - 1.0) * (2.0 * k + apb);
490 
491     a2  = a2 / a1;
492     a3  = a3 / a1;
493     a4  = a4 / a1;
494     *P  = (a2 + a3 * x) * pn1 - a4 * pn2;
495     pn2 = pn1;
496     pn1 = *P;
497   }
498   PetscFunctionReturn(0);
499 }
500 
501 #undef __FUNCT__
502 #define __FUNCT__ "PetscDTComputeJacobiDerivative"
503 /* Evaluates the first derivative of P_{n}^{a,b} at a point x. */
504 PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobiDerivative(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P)
505 {
506   PetscReal      nP;
507   PetscErrorCode ierr;
508 
509   PetscFunctionBegin;
510   if (!n) {*P = 0.0; PetscFunctionReturn(0);}
511   ierr = PetscDTComputeJacobi(a+1, b+1, n-1, x, &nP);CHKERRQ(ierr);
512   *P   = 0.5 * (a + b + n + 1) * nP;
513   PetscFunctionReturn(0);
514 }
515 
516 #undef __FUNCT__
517 #define __FUNCT__ "PetscDTMapSquareToTriangle_Internal"
518 /* Maps from [-1,1]^2 to the (-1,1) reference triangle */
519 PETSC_STATIC_INLINE PetscErrorCode PetscDTMapSquareToTriangle_Internal(PetscReal x, PetscReal y, PetscReal *xi, PetscReal *eta)
520 {
521   PetscFunctionBegin;
522   *xi  = 0.5 * (1.0 + x) * (1.0 - y) - 1.0;
523   *eta = y;
524   PetscFunctionReturn(0);
525 }
526 
527 #undef __FUNCT__
528 #define __FUNCT__ "PetscDTMapCubeToTetrahedron_Internal"
529 /* Maps from [-1,1]^2 to the (-1,1) reference triangle */
530 PETSC_STATIC_INLINE PetscErrorCode PetscDTMapCubeToTetrahedron_Internal(PetscReal x, PetscReal y, PetscReal z, PetscReal *xi, PetscReal *eta, PetscReal *zeta)
531 {
532   PetscFunctionBegin;
533   *xi   = 0.25 * (1.0 + x) * (1.0 - y) * (1.0 - z) - 1.0;
534   *eta  = 0.5  * (1.0 + y) * (1.0 - z) - 1.0;
535   *zeta = z;
536   PetscFunctionReturn(0);
537 }
538 
539 #undef __FUNCT__
540 #define __FUNCT__ "PetscDTGaussJacobiQuadrature1D_Internal"
541 static PetscErrorCode PetscDTGaussJacobiQuadrature1D_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal *x, PetscReal *w)
542 {
543   PetscInt       maxIter = 100;
544   PetscReal      eps     = 1.0e-8;
545   PetscReal      a1, a2, a3, a4, a5, a6;
546   PetscInt       k;
547   PetscErrorCode ierr;
548 
549   PetscFunctionBegin;
550 
551   a1      = PetscPowReal(2.0, a+b+1);
552 #if defined(PETSC_HAVE_TGAMMA)
553   a2      = PetscTGamma(a + npoints + 1);
554   a3      = PetscTGamma(b + npoints + 1);
555   a4      = PetscTGamma(a + b + npoints + 1);
556 #else
557   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"tgamma() - math routine is unavailable.");
558 #endif
559 
560   ierr = PetscDTFactorial_Internal(npoints, &a5);CHKERRQ(ierr);
561   a6   = a1 * a2 * a3 / a4 / a5;
562   /* Computes the m roots of P_{m}^{a,b} on [-1,1] by Newton's method with Chebyshev points as initial guesses.
563    Algorithm implemented from the pseudocode given by Karniadakis and Sherwin and Python in FIAT */
564   for (k = 0; k < npoints; ++k) {
565     PetscReal r = -PetscCosReal((2.0*k + 1.0) * PETSC_PI / (2.0 * npoints)), dP;
566     PetscInt  j;
567 
568     if (k > 0) r = 0.5 * (r + x[k-1]);
569     for (j = 0; j < maxIter; ++j) {
570       PetscReal s = 0.0, delta, f, fp;
571       PetscInt  i;
572 
573       for (i = 0; i < k; ++i) s = s + 1.0 / (r - x[i]);
574       ierr = PetscDTComputeJacobi(a, b, npoints, r, &f);CHKERRQ(ierr);
575       ierr = PetscDTComputeJacobiDerivative(a, b, npoints, r, &fp);CHKERRQ(ierr);
576       delta = f / (fp - f * s);
577       r     = r - delta;
578       if (PetscAbsReal(delta) < eps) break;
579     }
580     x[k] = r;
581     ierr = PetscDTComputeJacobiDerivative(a, b, npoints, x[k], &dP);CHKERRQ(ierr);
582     w[k] = a6 / (1.0 - PetscSqr(x[k])) / PetscSqr(dP);
583   }
584   PetscFunctionReturn(0);
585 }
586 
587 #undef __FUNCT__
588 #define __FUNCT__ "PetscDTGaussJacobiQuadrature"
589 /*@C
590   PetscDTGaussJacobiQuadrature - create Gauss-Jacobi quadrature for a simplex
591 
592   Not Collective
593 
594   Input Arguments:
595 + dim   - The simplex dimension
596 . order - The number of points in one dimension
597 . a     - left end of interval (often-1)
598 - b     - right end of interval (often +1)
599 
600   Output Argument:
601 . q - A PetscQuadrature object
602 
603   Level: intermediate
604 
605   References:
606   Karniadakis and Sherwin.
607   FIAT
608 
609 .seealso: PetscDTGaussTensorQuadrature(), PetscDTGaussQuadrature()
610 @*/
611 PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt dim, PetscInt order, PetscReal a, PetscReal b, PetscQuadrature *q)
612 {
613   PetscInt       npoints = dim > 1 ? dim > 2 ? order*PetscSqr(order) : PetscSqr(order) : order;
614   PetscReal     *px, *wx, *py, *wy, *pz, *wz, *x, *w;
615   PetscInt       i, j, k;
616   PetscErrorCode ierr;
617 
618   PetscFunctionBegin;
619   if ((a != -1.0) || (b != 1.0)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must use default internal right now");
620   ierr = PetscMalloc1(npoints*dim, &x);CHKERRQ(ierr);
621   ierr = PetscMalloc1(npoints, &w);CHKERRQ(ierr);
622   switch (dim) {
623   case 0:
624     ierr = PetscFree(x);CHKERRQ(ierr);
625     ierr = PetscFree(w);CHKERRQ(ierr);
626     ierr = PetscMalloc1(1, &x);CHKERRQ(ierr);
627     ierr = PetscMalloc1(1, &w);CHKERRQ(ierr);
628     x[0] = 0.0;
629     w[0] = 1.0;
630     break;
631   case 1:
632     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, x, w);CHKERRQ(ierr);
633     break;
634   case 2:
635     ierr = PetscMalloc4(order,&px,order,&wx,order,&py,order,&wy);CHKERRQ(ierr);
636     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, px, wx);CHKERRQ(ierr);
637     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 1.0, 0.0, py, wy);CHKERRQ(ierr);
638     for (i = 0; i < order; ++i) {
639       for (j = 0; j < order; ++j) {
640         ierr = PetscDTMapSquareToTriangle_Internal(px[i], py[j], &x[(i*order+j)*2+0], &x[(i*order+j)*2+1]);CHKERRQ(ierr);
641         w[i*order+j] = 0.5 * wx[i] * wy[j];
642       }
643     }
644     ierr = PetscFree4(px,wx,py,wy);CHKERRQ(ierr);
645     break;
646   case 3:
647     ierr = PetscMalloc6(order,&px,order,&wx,order,&py,order,&wy,order,&pz,order,&wz);CHKERRQ(ierr);
648     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, px, wx);CHKERRQ(ierr);
649     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 1.0, 0.0, py, wy);CHKERRQ(ierr);
650     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 2.0, 0.0, pz, wz);CHKERRQ(ierr);
651     for (i = 0; i < order; ++i) {
652       for (j = 0; j < order; ++j) {
653         for (k = 0; k < order; ++k) {
654           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);
655           w[(i*order+j)*order+k] = 0.125 * wx[i] * wy[j] * wz[k];
656         }
657       }
658     }
659     ierr = PetscFree6(px,wx,py,wy,pz,wz);CHKERRQ(ierr);
660     break;
661   default:
662     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim);
663   }
664   ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr);
665   ierr = PetscQuadratureSetOrder(*q, order);CHKERRQ(ierr);
666   ierr = PetscQuadratureSetData(*q, dim, npoints, x, w);CHKERRQ(ierr);
667   PetscFunctionReturn(0);
668 }
669 
670 #undef __FUNCT__
671 #define __FUNCT__ "PetscDTPseudoInverseQR"
672 /* Overwrites A. Can only handle full-rank problems with m>=n
673  * A in column-major format
674  * Ainv in row-major format
675  * tau has length m
676  * worksize must be >= max(1,n)
677  */
678 static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m,PetscInt mstride,PetscInt n,PetscReal *A_in,PetscReal *Ainv_out,PetscScalar *tau,PetscInt worksize,PetscScalar *work)
679 {
680   PetscErrorCode ierr;
681   PetscBLASInt   M,N,K,lda,ldb,ldwork,info;
682   PetscScalar    *A,*Ainv,*R,*Q,Alpha;
683 
684   PetscFunctionBegin;
685 #if defined(PETSC_USE_COMPLEX)
686   {
687     PetscInt i,j;
688     ierr = PetscMalloc2(m*n,&A,m*n,&Ainv);CHKERRQ(ierr);
689     for (j=0; j<n; j++) {
690       for (i=0; i<m; i++) A[i+m*j] = A_in[i+mstride*j];
691     }
692     mstride = m;
693   }
694 #else
695   A = A_in;
696   Ainv = Ainv_out;
697 #endif
698 
699   ierr = PetscBLASIntCast(m,&M);CHKERRQ(ierr);
700   ierr = PetscBLASIntCast(n,&N);CHKERRQ(ierr);
701   ierr = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr);
702   ierr = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr);
703   ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
704   PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info));
705   ierr = PetscFPTrapPop();CHKERRQ(ierr);
706   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error");
707   R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */
708 
709   /* Extract an explicit representation of Q */
710   Q = Ainv;
711   ierr = PetscMemcpy(Q,A,mstride*n*sizeof(PetscScalar));CHKERRQ(ierr);
712   K = N;                        /* full rank */
713   PetscStackCallBLAS("LAPACKungqr",LAPACKungqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info));
714   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error");
715 
716   /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */
717   Alpha = 1.0;
718   ldb = lda;
719   PetscStackCallBLAS("BLAStrsm",BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb));
720   /* Ainv is Q, overwritten with inverse */
721 
722 #if defined(PETSC_USE_COMPLEX)
723   {
724     PetscInt i;
725     for (i=0; i<m*n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]);
726     ierr = PetscFree2(A,Ainv);CHKERRQ(ierr);
727   }
728 #endif
729   PetscFunctionReturn(0);
730 }
731 
732 #undef __FUNCT__
733 #define __FUNCT__ "PetscDTLegendreIntegrate"
734 /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */
735 static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval,const PetscReal *x,PetscInt ndegree,const PetscInt *degrees,PetscBool Transpose,PetscReal *B)
736 {
737   PetscErrorCode ierr;
738   PetscReal      *Bv;
739   PetscInt       i,j;
740 
741   PetscFunctionBegin;
742   ierr = PetscMalloc1((ninterval+1)*ndegree,&Bv);CHKERRQ(ierr);
743   /* Point evaluation of L_p on all the source vertices */
744   ierr = PetscDTLegendreEval(ninterval+1,x,ndegree,degrees,Bv,NULL,NULL);CHKERRQ(ierr);
745   /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */
746   for (i=0; i<ninterval; i++) {
747     for (j=0; j<ndegree; j++) {
748       if (Transpose) B[i+ninterval*j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j];
749       else           B[i*ndegree+j]   = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j];
750     }
751   }
752   ierr = PetscFree(Bv);CHKERRQ(ierr);
753   PetscFunctionReturn(0);
754 }
755 
756 #undef __FUNCT__
757 #define __FUNCT__ "PetscDTReconstructPoly"
758 /*@
759    PetscDTReconstructPoly - create matrix representing polynomial reconstruction using cell intervals and evaluation at target intervals
760 
761    Not Collective
762 
763    Input Arguments:
764 +  degree - degree of reconstruction polynomial
765 .  nsource - number of source intervals
766 .  sourcex - sorted coordinates of source cell boundaries (length nsource+1)
767 .  ntarget - number of target intervals
768 -  targetx - sorted coordinates of target cell boundaries (length ntarget+1)
769 
770    Output Arguments:
771 .  R - reconstruction matrix, utarget = sum_s R[t*nsource+s] * usource[s]
772 
773    Level: advanced
774 
775 .seealso: PetscDTLegendreEval()
776 @*/
777 PetscErrorCode PetscDTReconstructPoly(PetscInt degree,PetscInt nsource,const PetscReal *sourcex,PetscInt ntarget,const PetscReal *targetx,PetscReal *R)
778 {
779   PetscErrorCode ierr;
780   PetscInt       i,j,k,*bdegrees,worksize;
781   PetscReal      xmin,xmax,center,hscale,*sourcey,*targety,*Bsource,*Bsinv,*Btarget;
782   PetscScalar    *tau,*work;
783 
784   PetscFunctionBegin;
785   PetscValidRealPointer(sourcex,3);
786   PetscValidRealPointer(targetx,5);
787   PetscValidRealPointer(R,6);
788   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);
789 #if defined(PETSC_USE_DEBUG)
790   for (i=0; i<nsource; i++) {
791     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]);
792   }
793   for (i=0; i<ntarget; i++) {
794     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]);
795   }
796 #endif
797   xmin = PetscMin(sourcex[0],targetx[0]);
798   xmax = PetscMax(sourcex[nsource],targetx[ntarget]);
799   center = (xmin + xmax)/2;
800   hscale = (xmax - xmin)/2;
801   worksize = nsource;
802   ierr = PetscMalloc4(degree+1,&bdegrees,nsource+1,&sourcey,nsource*(degree+1),&Bsource,worksize,&work);CHKERRQ(ierr);
803   ierr = PetscMalloc4(nsource,&tau,nsource*(degree+1),&Bsinv,ntarget+1,&targety,ntarget*(degree+1),&Btarget);CHKERRQ(ierr);
804   for (i=0; i<=nsource; i++) sourcey[i] = (sourcex[i]-center)/hscale;
805   for (i=0; i<=degree; i++) bdegrees[i] = i+1;
806   ierr = PetscDTLegendreIntegrate(nsource,sourcey,degree+1,bdegrees,PETSC_TRUE,Bsource);CHKERRQ(ierr);
807   ierr = PetscDTPseudoInverseQR(nsource,nsource,degree+1,Bsource,Bsinv,tau,nsource,work);CHKERRQ(ierr);
808   for (i=0; i<=ntarget; i++) targety[i] = (targetx[i]-center)/hscale;
809   ierr = PetscDTLegendreIntegrate(ntarget,targety,degree+1,bdegrees,PETSC_FALSE,Btarget);CHKERRQ(ierr);
810   for (i=0; i<ntarget; i++) {
811     PetscReal rowsum = 0;
812     for (j=0; j<nsource; j++) {
813       PetscReal sum = 0;
814       for (k=0; k<degree+1; k++) {
815         sum += Btarget[i*(degree+1)+k] * Bsinv[k*nsource+j];
816       }
817       R[i*nsource+j] = sum;
818       rowsum += sum;
819     }
820     for (j=0; j<nsource; j++) R[i*nsource+j] /= rowsum; /* normalize each row */
821   }
822   ierr = PetscFree4(bdegrees,sourcey,Bsource,work);CHKERRQ(ierr);
823   ierr = PetscFree4(tau,Bsinv,targety,Btarget);CHKERRQ(ierr);
824   PetscFunctionReturn(0);
825 }
826