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