xref: /petsc/src/dm/dt/interface/dt.c (revision 61bb57ed3979a08b29eaff0911e2dad60879110d)
1 /* Discretization tools */
2 
3 #include <petscdt.h>            /*I "petscdt.h" I*/
4 #include <petscblaslapack.h>
5 #include <petsc/private/petscimpl.h>
6 #include <petsc/private/dtimpl.h>
7 #include <petscviewer.h>
8 #include <petscdmplex.h>
9 #include <petscdmshell.h>
10 
11 #if defined(PETSC_HAVE_MPFR)
12 #include <mpfr.h>
13 #endif
14 
15 static PetscBool GolubWelschCite       = PETSC_FALSE;
16 const char       GolubWelschCitation[] = "@article{GolubWelsch1969,\n"
17                                          "  author  = {Golub and Welsch},\n"
18                                          "  title   = {Calculation of Quadrature Rules},\n"
19                                          "  journal = {Math. Comp.},\n"
20                                          "  volume  = {23},\n"
21                                          "  number  = {106},\n"
22                                          "  pages   = {221--230},\n"
23                                          "  year    = {1969}\n}\n";
24 
25 /* Numerical tests in src/dm/dt/tests/ex1.c show that when computing the nodes and weights of Gauss-Jacobi
26    quadrature rules:
27 
28    - in double precision, Newton's method and Golub & Welsch both work for moderate degrees (< 100),
29    - in single precision, Newton's method starts producing incorrect roots around n = 15, but
30      the weights from Golub & Welsch become a problem before then: they produces errors
31      in computing the Jacobi-polynomial Gram matrix around n = 6.
32 
33    So we default to Newton's method (required fewer dependencies) */
34 PetscBool PetscDTGaussQuadratureNewton_Internal = PETSC_TRUE;
35 
36 PetscClassId PETSCQUADRATURE_CLASSID = 0;
37 
38 /*@
39   PetscQuadratureCreate - Create a PetscQuadrature object
40 
41   Collective
42 
43   Input Parameter:
44 . comm - The communicator for the PetscQuadrature object
45 
46   Output Parameter:
47 . q  - The PetscQuadrature object
48 
49   Level: beginner
50 
51 .seealso: PetscQuadratureDestroy(), PetscQuadratureGetData()
52 @*/
53 PetscErrorCode PetscQuadratureCreate(MPI_Comm comm, PetscQuadrature *q)
54 {
55   PetscErrorCode ierr;
56 
57   PetscFunctionBegin;
58   PetscValidPointer(q, 2);
59   ierr = DMInitializePackage();CHKERRQ(ierr);
60   ierr = PetscHeaderCreate(*q,PETSCQUADRATURE_CLASSID,"PetscQuadrature","Quadrature","DT",comm,PetscQuadratureDestroy,PetscQuadratureView);CHKERRQ(ierr);
61   (*q)->dim       = -1;
62   (*q)->Nc        =  1;
63   (*q)->order     = -1;
64   (*q)->numPoints = 0;
65   (*q)->points    = NULL;
66   (*q)->weights   = NULL;
67   PetscFunctionReturn(0);
68 }
69 
70 /*@
71   PetscQuadratureDuplicate - Create a deep copy of the PetscQuadrature object
72 
73   Collective on q
74 
75   Input Parameter:
76 . q  - The PetscQuadrature object
77 
78   Output Parameter:
79 . r  - The new PetscQuadrature object
80 
81   Level: beginner
82 
83 .seealso: PetscQuadratureCreate(), PetscQuadratureDestroy(), PetscQuadratureGetData()
84 @*/
85 PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature q, PetscQuadrature *r)
86 {
87   PetscInt         order, dim, Nc, Nq;
88   const PetscReal *points, *weights;
89   PetscReal       *p, *w;
90   PetscErrorCode   ierr;
91 
92   PetscFunctionBegin;
93   PetscValidPointer(q, 2);
94   ierr = PetscQuadratureCreate(PetscObjectComm((PetscObject) q), r);CHKERRQ(ierr);
95   ierr = PetscQuadratureGetOrder(q, &order);CHKERRQ(ierr);
96   ierr = PetscQuadratureSetOrder(*r, order);CHKERRQ(ierr);
97   ierr = PetscQuadratureGetData(q, &dim, &Nc, &Nq, &points, &weights);CHKERRQ(ierr);
98   ierr = PetscMalloc1(Nq*dim, &p);CHKERRQ(ierr);
99   ierr = PetscMalloc1(Nq*Nc, &w);CHKERRQ(ierr);
100   ierr = PetscArraycpy(p, points, Nq*dim);CHKERRQ(ierr);
101   ierr = PetscArraycpy(w, weights, Nc * Nq);CHKERRQ(ierr);
102   ierr = PetscQuadratureSetData(*r, dim, Nc, Nq, p, w);CHKERRQ(ierr);
103   PetscFunctionReturn(0);
104 }
105 
106 /*@
107   PetscQuadratureDestroy - Destroys a PetscQuadrature object
108 
109   Collective on q
110 
111   Input Parameter:
112 . q  - The PetscQuadrature object
113 
114   Level: beginner
115 
116 .seealso: PetscQuadratureCreate(), PetscQuadratureGetData()
117 @*/
118 PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *q)
119 {
120   PetscErrorCode ierr;
121 
122   PetscFunctionBegin;
123   if (!*q) PetscFunctionReturn(0);
124   PetscValidHeaderSpecific((*q),PETSCQUADRATURE_CLASSID,1);
125   if (--((PetscObject)(*q))->refct > 0) {
126     *q = NULL;
127     PetscFunctionReturn(0);
128   }
129   ierr = PetscFree((*q)->points);CHKERRQ(ierr);
130   ierr = PetscFree((*q)->weights);CHKERRQ(ierr);
131   ierr = PetscHeaderDestroy(q);CHKERRQ(ierr);
132   PetscFunctionReturn(0);
133 }
134 
135 /*@
136   PetscQuadratureGetOrder - Return the order of the method
137 
138   Not collective
139 
140   Input Parameter:
141 . q - The PetscQuadrature object
142 
143   Output Parameter:
144 . order - The order of the quadrature, i.e. the highest degree polynomial that is exactly integrated
145 
146   Level: intermediate
147 
148 .seealso: PetscQuadratureSetOrder(), PetscQuadratureGetData(), PetscQuadratureSetData()
149 @*/
150 PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature q, PetscInt *order)
151 {
152   PetscFunctionBegin;
153   PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1);
154   PetscValidPointer(order, 2);
155   *order = q->order;
156   PetscFunctionReturn(0);
157 }
158 
159 /*@
160   PetscQuadratureSetOrder - Return the order of the method
161 
162   Not collective
163 
164   Input Parameters:
165 + q - The PetscQuadrature object
166 - order - The order of the quadrature, i.e. the highest degree polynomial that is exactly integrated
167 
168   Level: intermediate
169 
170 .seealso: PetscQuadratureGetOrder(), PetscQuadratureGetData(), PetscQuadratureSetData()
171 @*/
172 PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature q, PetscInt order)
173 {
174   PetscFunctionBegin;
175   PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1);
176   q->order = order;
177   PetscFunctionReturn(0);
178 }
179 
180 /*@
181   PetscQuadratureGetNumComponents - Return the number of components for functions to be integrated
182 
183   Not collective
184 
185   Input Parameter:
186 . q - The PetscQuadrature object
187 
188   Output Parameter:
189 . Nc - The number of components
190 
191   Note: We are performing an integral int f(x) . w(x) dx, where both f and w (the weight) have Nc components.
192 
193   Level: intermediate
194 
195 .seealso: PetscQuadratureSetNumComponents(), PetscQuadratureGetData(), PetscQuadratureSetData()
196 @*/
197 PetscErrorCode PetscQuadratureGetNumComponents(PetscQuadrature q, PetscInt *Nc)
198 {
199   PetscFunctionBegin;
200   PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1);
201   PetscValidPointer(Nc, 2);
202   *Nc = q->Nc;
203   PetscFunctionReturn(0);
204 }
205 
206 /*@
207   PetscQuadratureSetNumComponents - Return the number of components for functions to be integrated
208 
209   Not collective
210 
211   Input Parameters:
212 + q  - The PetscQuadrature object
213 - Nc - The number of components
214 
215   Note: We are performing an integral int f(x) . w(x) dx, where both f and w (the weight) have Nc components.
216 
217   Level: intermediate
218 
219 .seealso: PetscQuadratureGetNumComponents(), PetscQuadratureGetData(), PetscQuadratureSetData()
220 @*/
221 PetscErrorCode PetscQuadratureSetNumComponents(PetscQuadrature q, PetscInt Nc)
222 {
223   PetscFunctionBegin;
224   PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1);
225   q->Nc = Nc;
226   PetscFunctionReturn(0);
227 }
228 
229 /*@C
230   PetscQuadratureGetData - Returns the data defining the quadrature
231 
232   Not collective
233 
234   Input Parameter:
235 . q  - The PetscQuadrature object
236 
237   Output Parameters:
238 + dim - The spatial dimension
239 . Nc - The number of components
240 . npoints - The number of quadrature points
241 . points - The coordinates of each quadrature point
242 - weights - The weight of each quadrature point
243 
244   Level: intermediate
245 
246   Fortran Notes:
247     From Fortran you must call PetscQuadratureRestoreData() when you are done with the data
248 
249 .seealso: PetscQuadratureCreate(), PetscQuadratureSetData()
250 @*/
251 PetscErrorCode PetscQuadratureGetData(PetscQuadrature q, PetscInt *dim, PetscInt *Nc, PetscInt *npoints, const PetscReal *points[], const PetscReal *weights[])
252 {
253   PetscFunctionBegin;
254   PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1);
255   if (dim) {
256     PetscValidPointer(dim, 2);
257     *dim = q->dim;
258   }
259   if (Nc) {
260     PetscValidPointer(Nc, 3);
261     *Nc = q->Nc;
262   }
263   if (npoints) {
264     PetscValidPointer(npoints, 4);
265     *npoints = q->numPoints;
266   }
267   if (points) {
268     PetscValidPointer(points, 5);
269     *points = q->points;
270   }
271   if (weights) {
272     PetscValidPointer(weights, 6);
273     *weights = q->weights;
274   }
275   PetscFunctionReturn(0);
276 }
277 
278 static PetscErrorCode PetscDTJacobianInverse_Internal(PetscInt m, PetscInt n, const PetscReal J[], PetscReal Jinv[])
279 {
280   PetscScalar    *Js, *Jinvs;
281   PetscInt       i, j, k;
282   PetscBLASInt   bm, bn, info;
283   PetscErrorCode ierr;
284 
285   PetscFunctionBegin;
286   ierr = PetscBLASIntCast(m, &bm);CHKERRQ(ierr);
287   ierr = PetscBLASIntCast(n, &bn);CHKERRQ(ierr);
288 #if defined(PETSC_USE_COMPLEX)
289   ierr = PetscMalloc2(m*n, &Js, m*n, &Jinvs);CHKERRQ(ierr);
290   for (i = 0; i < m*n; i++) Js[i] = J[i];
291 #else
292   Js = (PetscReal *) J;
293   Jinvs = Jinv;
294 #endif
295   if (m == n) {
296     PetscBLASInt *pivots;
297     PetscScalar *W;
298 
299     ierr = PetscMalloc2(m, &pivots, m, &W);CHKERRQ(ierr);
300 
301     ierr = PetscArraycpy(Jinvs, Js, m * m);CHKERRQ(ierr);
302     PetscStackCallBLAS("LAPACKgetrf", LAPACKgetrf_(&bm, &bm, Jinvs, &bm, pivots, &info));
303     if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetrf %D",(PetscInt)info);
304     PetscStackCallBLAS("LAPACKgetri", LAPACKgetri_(&bm, Jinvs, &bm, pivots, W, &bm, &info));
305     if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetri %D",(PetscInt)info);
306     ierr = PetscFree2(pivots, W);CHKERRQ(ierr);
307   } else if (m < n) {
308     PetscScalar *JJT;
309     PetscBLASInt *pivots;
310     PetscScalar *W;
311 
312     ierr = PetscMalloc1(m*m, &JJT);CHKERRQ(ierr);
313     ierr = PetscMalloc2(m, &pivots, m, &W);CHKERRQ(ierr);
314     for (i = 0; i < m; i++) {
315       for (j = 0; j < m; j++) {
316         PetscScalar val = 0.;
317 
318         for (k = 0; k < n; k++) val += Js[i * n + k] * Js[j * n + k];
319         JJT[i * m + j] = val;
320       }
321     }
322 
323     PetscStackCallBLAS("LAPACKgetrf", LAPACKgetrf_(&bm, &bm, JJT, &bm, pivots, &info));
324     if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetrf %D",(PetscInt)info);
325     PetscStackCallBLAS("LAPACKgetri", LAPACKgetri_(&bm, JJT, &bm, pivots, W, &bm, &info));
326     if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetri %D",(PetscInt)info);
327     for (i = 0; i < n; i++) {
328       for (j = 0; j < m; j++) {
329         PetscScalar val = 0.;
330 
331         for (k = 0; k < m; k++) val += Js[k * n + i] * JJT[k * m + j];
332         Jinvs[i * m + j] = val;
333       }
334     }
335     ierr = PetscFree2(pivots, W);CHKERRQ(ierr);
336     ierr = PetscFree(JJT);CHKERRQ(ierr);
337   } else {
338     PetscScalar *JTJ;
339     PetscBLASInt *pivots;
340     PetscScalar *W;
341 
342     ierr = PetscMalloc1(n*n, &JTJ);CHKERRQ(ierr);
343     ierr = PetscMalloc2(n, &pivots, n, &W);CHKERRQ(ierr);
344     for (i = 0; i < n; i++) {
345       for (j = 0; j < n; j++) {
346         PetscScalar val = 0.;
347 
348         for (k = 0; k < m; k++) val += Js[k * n + i] * Js[k * n + j];
349         JTJ[i * n + j] = val;
350       }
351     }
352 
353     PetscStackCallBLAS("LAPACKgetrf", LAPACKgetrf_(&bn, &bn, JTJ, &bm, pivots, &info));
354     if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetrf %D",(PetscInt)info);
355     PetscStackCallBLAS("LAPACKgetri", LAPACKgetri_(&bn, JTJ, &bn, pivots, W, &bn, &info));
356     if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error returned from LAPACKgetri %D",(PetscInt)info);
357     for (i = 0; i < n; i++) {
358       for (j = 0; j < m; j++) {
359         PetscScalar val = 0.;
360 
361         for (k = 0; k < n; k++) val += JTJ[i * n + k] * Js[j * n + k];
362         Jinvs[i * m + j] = val;
363       }
364     }
365     ierr = PetscFree2(pivots, W);CHKERRQ(ierr);
366     ierr = PetscFree(JTJ);CHKERRQ(ierr);
367   }
368 #if defined(PETSC_USE_COMPLEX)
369   for (i = 0; i < m*n; i++) Jinv[i] = PetscRealPart(Jinvs[i]);
370   ierr = PetscFree2(Js, Jinvs);CHKERRQ(ierr);
371 #endif
372   PetscFunctionReturn(0);
373 }
374 
375 /*@
376    PetscQuadraturePushForward - Push forward a quadrature functional under an affine transformation.
377 
378    Collecive on PetscQuadrature
379 
380    Input Arguments:
381 +  q - the quadrature functional
382 .  imageDim - the dimension of the image of the transformation
383 .  origin - a point in the original space
384 .  originImage - the image of the origin under the transformation
385 .  J - the Jacobian of the image: an [imageDim x dim] matrix in row major order
386 -  formDegree - transform the quadrature weights as k-forms of this form degree (if the number of components is a multiple of (dim choose formDegree), it is assumed that they represent multiple k-forms) [see PetscDTAltVPullback() for interpretation of formDegree]
387 
388    Output Arguments:
389 .  Jinvstarq - a quadrature rule where each point is the image of a point in the original quadrature rule, and where the k-form weights have been pulled-back by the pseudoinverse of J to the k-form weights in the image space.
390 
391    Note: the new quadrature rule will have a different number of components if spaces have different dimensions.  For example, pushing a 2-form forward from a two dimensional space to a three dimensional space changes the number of components from 1 to 3.
392 
393    Level: intermediate
394 
395 .seealso: PetscDTAltVPullback(), PetscDTAltVPullbackMatrix()
396 @*/
397 PetscErrorCode PetscQuadraturePushForward(PetscQuadrature q, PetscInt imageDim, const PetscReal origin[], const PetscReal originImage[], const PetscReal J[], PetscInt formDegree, PetscQuadrature *Jinvstarq)
398 {
399   PetscInt         dim, Nc, imageNc, formSize, Ncopies, imageFormSize, Npoints, pt, i, j, c;
400   const PetscReal *points;
401   const PetscReal *weights;
402   PetscReal       *imagePoints, *imageWeights;
403   PetscReal       *Jinv;
404   PetscReal       *Jinvstar;
405   PetscErrorCode   ierr;
406 
407   PetscFunctionBegin;
408   PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 1);
409   if (imageDim < PetscAbsInt(formDegree)) SETERRQ2(PetscObjectComm((PetscObject)q), PETSC_ERR_ARG_INCOMP, "Cannot represent a %D-form in %D dimensions", PetscAbsInt(formDegree), imageDim);
410   ierr = PetscQuadratureGetData(q, &dim, &Nc, &Npoints, &points, &weights);CHKERRQ(ierr);
411   ierr = PetscDTBinomialInt(dim, PetscAbsInt(formDegree), &formSize);CHKERRQ(ierr);
412   if (Nc % formSize) SETERRQ2(PetscObjectComm((PetscObject)q), PETSC_ERR_ARG_INCOMP, "Number of components %D is not a multiple of formSize %D\n", Nc, formSize);
413   Ncopies = Nc / formSize;
414   ierr = PetscDTBinomialInt(imageDim, PetscAbsInt(formDegree), &imageFormSize);CHKERRQ(ierr);
415   imageNc = Ncopies * imageFormSize;
416   ierr = PetscMalloc1(Npoints * imageDim, &imagePoints);CHKERRQ(ierr);
417   ierr = PetscMalloc1(Npoints * imageNc, &imageWeights);CHKERRQ(ierr);
418   ierr = PetscMalloc2(imageDim * dim, &Jinv, formSize * imageFormSize, &Jinvstar);CHKERRQ(ierr);
419   ierr = PetscDTJacobianInverse_Internal(dim, imageDim, J, Jinv);CHKERRQ(ierr);
420   ierr = PetscDTAltVPullbackMatrix(imageDim, dim, Jinv, formDegree, Jinvstar);CHKERRQ(ierr);
421   for (pt = 0; pt < Npoints; pt++) {
422     const PetscReal *point = &points[pt * dim];
423     PetscReal       *imagePoint = &imagePoints[pt * imageDim];
424 
425     for (i = 0; i < imageDim; i++) {
426       PetscReal val = originImage[i];
427 
428       for (j = 0; j < dim; j++) val += J[i * dim + j] * (point[j] - origin[j]);
429       imagePoint[i] = val;
430     }
431     for (c = 0; c < Ncopies; c++) {
432       const PetscReal *form = &weights[pt * Nc + c * formSize];
433       PetscReal       *imageForm = &imageWeights[pt * imageNc + c * imageFormSize];
434 
435       for (i = 0; i < imageFormSize; i++) {
436         PetscReal val = 0.;
437 
438         for (j = 0; j < formSize; j++) val += Jinvstar[i * formSize + j] * form[j];
439         imageForm[i] = val;
440       }
441     }
442   }
443   ierr = PetscQuadratureCreate(PetscObjectComm((PetscObject)q), Jinvstarq);CHKERRQ(ierr);
444   ierr = PetscQuadratureSetData(*Jinvstarq, imageDim, imageNc, Npoints, imagePoints, imageWeights);CHKERRQ(ierr);
445   ierr = PetscFree2(Jinv, Jinvstar);CHKERRQ(ierr);
446   PetscFunctionReturn(0);
447 }
448 
449 /*@C
450   PetscQuadratureSetData - Sets the data defining the quadrature
451 
452   Not collective
453 
454   Input Parameters:
455 + q  - The PetscQuadrature object
456 . dim - The spatial dimension
457 . Nc - The number of components
458 . npoints - The number of quadrature points
459 . points - The coordinates of each quadrature point
460 - weights - The weight of each quadrature point
461 
462   Note: This routine owns the references to points and weights, so they must be allocated using PetscMalloc() and the user should not free them.
463 
464   Level: intermediate
465 
466 .seealso: PetscQuadratureCreate(), PetscQuadratureGetData()
467 @*/
468 PetscErrorCode PetscQuadratureSetData(PetscQuadrature q, PetscInt dim, PetscInt Nc, PetscInt npoints, const PetscReal points[], const PetscReal weights[])
469 {
470   PetscFunctionBegin;
471   PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1);
472   if (dim >= 0)     q->dim       = dim;
473   if (Nc >= 0)      q->Nc        = Nc;
474   if (npoints >= 0) q->numPoints = npoints;
475   if (points) {
476     PetscValidPointer(points, 4);
477     q->points = points;
478   }
479   if (weights) {
480     PetscValidPointer(weights, 5);
481     q->weights = weights;
482   }
483   PetscFunctionReturn(0);
484 }
485 
486 static PetscErrorCode PetscQuadratureView_Ascii(PetscQuadrature quad, PetscViewer v)
487 {
488   PetscInt          q, d, c;
489   PetscViewerFormat format;
490   PetscErrorCode    ierr;
491 
492   PetscFunctionBegin;
493   if (quad->Nc > 1) {ierr = PetscViewerASCIIPrintf(v, "Quadrature of order %D on %D points (dim %D) with %D components\n", quad->order, quad->numPoints, quad->dim, quad->Nc);CHKERRQ(ierr);}
494   else              {ierr = PetscViewerASCIIPrintf(v, "Quadrature of order %D on %D points (dim %D)\n", quad->order, quad->numPoints, quad->dim);CHKERRQ(ierr);}
495   ierr = PetscViewerGetFormat(v, &format);CHKERRQ(ierr);
496   if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) PetscFunctionReturn(0);
497   for (q = 0; q < quad->numPoints; ++q) {
498     ierr = PetscViewerASCIIPrintf(v, "p%D (", q);CHKERRQ(ierr);
499     ierr = PetscViewerASCIIUseTabs(v, PETSC_FALSE);CHKERRQ(ierr);
500     for (d = 0; d < quad->dim; ++d) {
501       if (d) {ierr = PetscViewerASCIIPrintf(v, ", ");CHKERRQ(ierr);}
502       ierr = PetscViewerASCIIPrintf(v, "%+g", (double)quad->points[q*quad->dim+d]);CHKERRQ(ierr);
503     }
504     ierr = PetscViewerASCIIPrintf(v, ") ");CHKERRQ(ierr);
505     if (quad->Nc > 1) {ierr = PetscViewerASCIIPrintf(v, "w%D (", q);CHKERRQ(ierr);}
506     for (c = 0; c < quad->Nc; ++c) {
507       if (c) {ierr = PetscViewerASCIIPrintf(v, ", ");CHKERRQ(ierr);}
508       ierr = PetscViewerASCIIPrintf(v, "%+g", (double)quad->weights[q*quad->Nc+c]);CHKERRQ(ierr);
509     }
510     if (quad->Nc > 1) {ierr = PetscViewerASCIIPrintf(v, ")");CHKERRQ(ierr);}
511     ierr = PetscViewerASCIIPrintf(v, "\n");CHKERRQ(ierr);
512     ierr = PetscViewerASCIIUseTabs(v, PETSC_TRUE);CHKERRQ(ierr);
513   }
514   PetscFunctionReturn(0);
515 }
516 
517 /*@C
518   PetscQuadratureView - Views a PetscQuadrature object
519 
520   Collective on quad
521 
522   Input Parameters:
523 + quad  - The PetscQuadrature object
524 - viewer - The PetscViewer object
525 
526   Level: beginner
527 
528 .seealso: PetscQuadratureCreate(), PetscQuadratureGetData()
529 @*/
530 PetscErrorCode PetscQuadratureView(PetscQuadrature quad, PetscViewer viewer)
531 {
532   PetscBool      iascii;
533   PetscErrorCode ierr;
534 
535   PetscFunctionBegin;
536   PetscValidHeader(quad, 1);
537   if (viewer) PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
538   if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) quad), &viewer);CHKERRQ(ierr);}
539   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
540   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
541   if (iascii) {ierr = PetscQuadratureView_Ascii(quad, viewer);CHKERRQ(ierr);}
542   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
543   PetscFunctionReturn(0);
544 }
545 
546 /*@C
547   PetscQuadratureExpandComposite - Return a quadrature over the composite element, which has the original quadrature in each subelement
548 
549   Not collective
550 
551   Input Parameter:
552 + q - The original PetscQuadrature
553 . numSubelements - The number of subelements the original element is divided into
554 . v0 - An array of the initial points for each subelement
555 - jac - An array of the Jacobian mappings from the reference to each subelement
556 
557   Output Parameters:
558 . dim - The dimension
559 
560   Note: Together v0 and jac define an affine mapping from the original reference element to each subelement
561 
562  Not available from Fortran
563 
564   Level: intermediate
565 
566 .seealso: PetscFECreate(), PetscSpaceGetDimension(), PetscDualSpaceGetDimension()
567 @*/
568 PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature q, PetscInt numSubelements, const PetscReal v0[], const PetscReal jac[], PetscQuadrature *qref)
569 {
570   const PetscReal *points,    *weights;
571   PetscReal       *pointsRef, *weightsRef;
572   PetscInt         dim, Nc, order, npoints, npointsRef, c, p, cp, d, e;
573   PetscErrorCode   ierr;
574 
575   PetscFunctionBegin;
576   PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1);
577   PetscValidPointer(v0, 3);
578   PetscValidPointer(jac, 4);
579   PetscValidPointer(qref, 5);
580   ierr = PetscQuadratureCreate(PETSC_COMM_SELF, qref);CHKERRQ(ierr);
581   ierr = PetscQuadratureGetOrder(q, &order);CHKERRQ(ierr);
582   ierr = PetscQuadratureGetData(q, &dim, &Nc, &npoints, &points, &weights);CHKERRQ(ierr);
583   npointsRef = npoints*numSubelements;
584   ierr = PetscMalloc1(npointsRef*dim,&pointsRef);CHKERRQ(ierr);
585   ierr = PetscMalloc1(npointsRef*Nc, &weightsRef);CHKERRQ(ierr);
586   for (c = 0; c < numSubelements; ++c) {
587     for (p = 0; p < npoints; ++p) {
588       for (d = 0; d < dim; ++d) {
589         pointsRef[(c*npoints + p)*dim+d] = v0[c*dim+d];
590         for (e = 0; e < dim; ++e) {
591           pointsRef[(c*npoints + p)*dim+d] += jac[(c*dim + d)*dim+e]*(points[p*dim+e] + 1.0);
592         }
593       }
594       /* Could also use detJ here */
595       for (cp = 0; cp < Nc; ++cp) weightsRef[(c*npoints+p)*Nc+cp] = weights[p*Nc+cp]/numSubelements;
596     }
597   }
598   ierr = PetscQuadratureSetOrder(*qref, order);CHKERRQ(ierr);
599   ierr = PetscQuadratureSetData(*qref, dim, Nc, npointsRef, pointsRef, weightsRef);CHKERRQ(ierr);
600   PetscFunctionReturn(0);
601 }
602 
603 /* Compute the coefficients for the Jacobi polynomial recurrence,
604  *
605  * J^{a,b}_n(x) = (cnm1 + cnm1x * x) * J^{a,b}_{n-1}(x) - cnm2 * J^{a,b}_{n-2}(x).
606  */
607 #define PetscDTJacobiRecurrence_Internal(n,a,b,cnm1,cnm1x,cnm2) \
608 do {                                                            \
609   PetscReal _a = (a);                                           \
610   PetscReal _b = (b);                                           \
611   PetscReal _n = (n);                                           \
612   if (n == 1) {                                                 \
613     (cnm1) = (_a-_b) * 0.5;                                     \
614     (cnm1x) = (_a+_b+2.)*0.5;                                   \
615     (cnm2) = 0.;                                                \
616   } else {                                                      \
617     PetscReal _2n = _n+_n;                                      \
618     PetscReal _d = (_2n*(_n+_a+_b)*(_2n+_a+_b-2));              \
619     PetscReal _n1 = (_2n+_a+_b-1.)*(_a*_a-_b*_b);               \
620     PetscReal _n1x = (_2n+_a+_b-1.)*(_2n+_a+_b)*(_2n+_a+_b-2);  \
621     PetscReal _n2 = 2.*((_n+_a-1.)*(_n+_b-1.)*(_2n+_a+_b));     \
622     (cnm1) = _n1 / _d;                                          \
623     (cnm1x) = _n1x / _d;                                        \
624     (cnm2) = _n2 / _d;                                          \
625   }                                                             \
626 } while (0)
627 
628 static PetscErrorCode PetscDTJacobiEval_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscInt k, const PetscReal *points, PetscInt ndegree, const PetscInt *degrees, PetscReal *p)
629 {
630   PetscReal ak, bk;
631   PetscReal abk1;
632   PetscInt i,l,maxdegree;
633 
634   PetscFunctionBegin;
635   maxdegree = degrees[ndegree-1] - k;
636   ak = a + k;
637   bk = b + k;
638   abk1 = a + b + k + 1.;
639   if (maxdegree < 0) {
640     for (i = 0; i < npoints; i++) for (l = 0; l < ndegree; l++) p[i*ndegree+l] = 0.;
641     PetscFunctionReturn(0);
642   }
643   for (i=0; i<npoints; i++) {
644     PetscReal pm1,pm2,x;
645     PetscReal cnm1, cnm1x, cnm2;
646     PetscInt  j,m;
647 
648     x    = points[i];
649     pm2  = 1.;
650     PetscDTJacobiRecurrence_Internal(1,ak,bk,cnm1,cnm1x,cnm2);
651     pm1 = (cnm1 + cnm1x*x);
652     l    = 0;
653     while (l < ndegree && degrees[l] - k < 0) {
654       p[l++] = 0.;
655     }
656     while (l < ndegree && degrees[l] - k == 0) {
657       p[l] = pm2;
658       for (m = 0; m < k; m++) p[l] *= (abk1 + m) * 0.5;
659       l++;
660     }
661     while (l < ndegree && degrees[l] - k == 1) {
662       p[l] = pm1;
663       for (m = 0; m < k; m++) p[l] *= (abk1 + 1 + m) * 0.5;
664       l++;
665     }
666     for (j=2; j<=maxdegree; j++) {
667       PetscReal pp;
668 
669       PetscDTJacobiRecurrence_Internal(j,ak,bk,cnm1,cnm1x,cnm2);
670       pp   = (cnm1 + cnm1x*x)*pm1 - cnm2*pm2;
671       pm2  = pm1;
672       pm1  = pp;
673       while (l < ndegree && degrees[l] - k == j) {
674         p[l] = pp;
675         for (m = 0; m < k; m++) p[l] *= (abk1 + j + m) * 0.5;
676         l++;
677       }
678     }
679     p += ndegree;
680   }
681   PetscFunctionReturn(0);
682 }
683 
684 /*@
685    PetscDTJacobiEval - evaluate Jacobi polynomials for the weight function $(1.+x)^{\alpha} (1.-x)^{\beta}$
686                        at points
687 
688    Not Collective
689 
690    Input Arguments:
691 +  npoints - number of spatial points to evaluate at
692 .  alpha - the left exponent > -1
693 .  beta - the right exponent > -1
694 .  points - array of locations to evaluate at
695 .  ndegree - number of basis degrees to evaluate
696 -  degrees - sorted array of degrees to evaluate
697 
698    Output Arguments:
699 +  B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL)
700 .  D - row-oriented derivative evaluation matrix (or NULL)
701 -  D2 - row-oriented second derivative evaluation matrix (or NULL)
702 
703    Level: intermediate
704 
705 .seealso: PetscDTGaussQuadrature()
706 @*/
707 PetscErrorCode PetscDTJacobiEval(PetscInt npoints,PetscReal alpha, PetscReal beta, const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2)
708 {
709   PetscErrorCode ierr;
710 
711   PetscFunctionBegin;
712   if (alpha <= -1.) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"alpha must be > -1.");
713   if (beta <= -1.) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"beta must be > -1.");
714   if (!npoints || !ndegree) PetscFunctionReturn(0);
715   if (B)  {ierr = PetscDTJacobiEval_Internal(npoints, alpha, beta, 0, points, ndegree, degrees, B);CHKERRQ(ierr);}
716   if (D)  {ierr = PetscDTJacobiEval_Internal(npoints, alpha, beta, 1, points, ndegree, degrees, D);CHKERRQ(ierr);}
717   if (D2) {ierr = PetscDTJacobiEval_Internal(npoints, alpha, beta, 2, points, ndegree, degrees, D2);CHKERRQ(ierr);}
718   PetscFunctionReturn(0);
719 }
720 
721 /*@
722    PetscDTLegendreEval - evaluate Legendre polynomials at points
723 
724    Not Collective
725 
726    Input Arguments:
727 +  npoints - number of spatial points to evaluate at
728 .  points - array of locations to evaluate at
729 .  ndegree - number of basis degrees to evaluate
730 -  degrees - sorted array of degrees to evaluate
731 
732    Output Arguments:
733 +  B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL)
734 .  D - row-oriented derivative evaluation matrix (or NULL)
735 -  D2 - row-oriented second derivative evaluation matrix (or NULL)
736 
737    Level: intermediate
738 
739 .seealso: PetscDTGaussQuadrature()
740 @*/
741 PetscErrorCode PetscDTLegendreEval(PetscInt npoints,const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2)
742 {
743   PetscErrorCode ierr;
744 
745   PetscFunctionBegin;
746   ierr = PetscDTJacobiEval(npoints, 0., 0., points, ndegree, degrees, B, D, D2);CHKERRQ(ierr);
747   PetscFunctionReturn(0);
748 }
749 
750 /* solve the symmetric tridiagonal eigenvalue system, writing the eigenvalues into eigs and the eigenvectors into V
751  * with lds n; diag and subdiag are overwritten */
752 static PetscErrorCode PetscDTSymmetricTridiagonalEigensolve(PetscInt n, PetscReal diag[], PetscReal subdiag[],
753                                                             PetscReal eigs[], PetscScalar V[])
754 {
755   char jobz = 'V'; /* eigenvalues and eigenvectors */
756   char range = 'A'; /* all eigenvalues will be found */
757   PetscReal VL = 0.; /* ignored because range is 'A' */
758   PetscReal VU = 0.; /* ignored because range is 'A' */
759   PetscBLASInt IL = 0; /* ignored because range is 'A' */
760   PetscBLASInt IU = 0; /* ignored because range is 'A' */
761   PetscReal abstol = 0.; /* unused */
762   PetscBLASInt bn, bm, ldz; /* bm will equal bn on exit */
763   PetscBLASInt *isuppz;
764   PetscBLASInt lwork, liwork;
765   PetscReal workquery;
766   PetscBLASInt  iworkquery;
767   PetscBLASInt *iwork;
768   PetscBLASInt info;
769   PetscReal *work = NULL;
770   PetscErrorCode ierr;
771 
772   PetscFunctionBegin;
773 #if !defined(PETSCDTGAUSSIANQUADRATURE_EIG)
774   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP_SYS, "A LAPACK symmetric tridiagonal eigensolver could not be found");
775 #endif
776   ierr = PetscBLASIntCast(n, &bn);CHKERRQ(ierr);
777   ierr = PetscBLASIntCast(n, &ldz);CHKERRQ(ierr);
778 #if !defined(PETSC_MISSING_LAPACK_STEGR)
779   ierr = PetscMalloc1(2 * n, &isuppz);CHKERRQ(ierr);
780   lwork = -1;
781   liwork = -1;
782   PetscStackCallBLAS("LAPACKstegr",LAPACKstegr_(&jobz,&range,&bn,diag,subdiag,&VL,&VU,&IL,&IU,&abstol,&bm,eigs,V,&ldz,isuppz,&workquery,&lwork,&iworkquery,&liwork,&info));
783   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEGR error");
784   lwork = (PetscBLASInt) workquery;
785   liwork = (PetscBLASInt) iworkquery;
786   ierr = PetscMalloc2(lwork, &work, liwork, &iwork);CHKERRQ(ierr);
787   ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
788   PetscStackCallBLAS("LAPACKstegr",LAPACKstegr_(&jobz,&range,&bn,diag,subdiag,&VL,&VU,&IL,&IU,&abstol,&bm,eigs,V,&ldz,isuppz,work,&lwork,iwork,&liwork,&info));
789   ierr = PetscFPTrapPop();CHKERRQ(ierr);
790   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEGR error");
791   ierr = PetscFree2(work, iwork);CHKERRQ(ierr);
792   ierr = PetscFree(isuppz);CHKERRQ(ierr);
793 #elif !defined(PETSC_MISSING_LAPACK_STEQR)
794   jobz = 'I'; /* Compute eigenvalues and eigenvectors of the
795                  tridiagonal matrix.  Z is initialized to the identity
796                  matrix. */
797   ierr = PetscMalloc1(PetscMax(1,2*n-2),&work);CHKERRQ(ierr);
798   PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("I",&bn,diag,subdiag,V,&ldz,work,&info));
799   ierr = PetscFPTrapPop();CHKERRQ(ierr);
800   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEQR error");
801   ierr = PetscFree(work);CHKERRQ(ierr);
802   ierr = PetscArraycpy(eigs,diag,n);CHKERRQ(ierr);
803 #endif
804   PetscFunctionReturn(0);
805 }
806 
807 /* Formula for the weights at the endpoints (-1 and 1) of Gauss-Lobatto-Jacobi
808  * quadrature rules on the interval [-1, 1] */
809 static PetscErrorCode PetscDTGaussLobattoJacobiEndweights_Internal(PetscInt n, PetscReal alpha, PetscReal beta, PetscReal *leftw, PetscReal *rightw)
810 {
811   PetscReal twoab1;
812   PetscInt  m = n - 2;
813   PetscReal a = alpha + 1.;
814   PetscReal b = beta + 1.;
815   PetscReal gra, grb;
816 
817   PetscFunctionBegin;
818   twoab1 = PetscPowReal(2., a + b - 1.);
819 #if defined(PETSC_HAVE_LGAMMA)
820   grb = PetscExpReal(2. * PetscLGamma(b+1.) + PetscLGamma(m+1.) + PetscLGamma(m+a+1.) -
821                      (PetscLGamma(m+b+1) + PetscLGamma(m+a+b+1.)));
822   gra = PetscExpReal(2. * PetscLGamma(a+1.) + PetscLGamma(m+1.) + PetscLGamma(m+b+1.) -
823                      (PetscLGamma(m+a+1) + PetscLGamma(m+a+b+1.)));
824 #else
825   {
826     PetscInt alphai = (PetscInt) alpha;
827     PetscInt betai = (PetscInt) beta;
828     PetscErrorCode ierr;
829 
830     if ((PetscReal) alphai == alpha && (PetscReal) betai == beta) {
831       PetscReal binom1, binom2;
832 
833       ierr = PetscDTBinomial(m+b, b, &binom1);CHKERRQ(ierr);
834       ierr = PetscDTBinomial(m+a+b, b, &binom2);CHKERRQ(ierr);
835       grb = 1./ (binom1 * binom2);
836       ierr = PetscDTBinomial(m+a, a, &binom1);CHKERRQ(ierr);
837       ierr = PetscDTBinomial(m+a+b, a, &binom2);CHKERRQ(ierr);
838       gra = 1./ (binom1 * binom2);
839     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"lgamma() - math routine is unavailable.");
840   }
841 #endif
842   *leftw = twoab1 * grb / b;
843   *rightw = twoab1 * gra / a;
844   PetscFunctionReturn(0);
845 }
846 
847 /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x.
848    Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */
849 PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobi(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P)
850 {
851   PetscReal pn1, pn2;
852   PetscReal cnm1, cnm1x, cnm2;
853   PetscInt  k;
854 
855   PetscFunctionBegin;
856   if (!n) {*P = 1.0; PetscFunctionReturn(0);}
857   PetscDTJacobiRecurrence_Internal(1,a,b,cnm1,cnm1x,cnm2);
858   pn2 = 1.;
859   pn1 = cnm1 + cnm1x*x;
860   if (n == 1) {*P = pn1; PetscFunctionReturn(0);}
861   *P  = 0.0;
862   for (k = 2; k < n+1; ++k) {
863     PetscDTJacobiRecurrence_Internal(k,a,b,cnm1,cnm1x,cnm2);
864 
865     *P  = (cnm1 + cnm1x*x)*pn1 - cnm2*pn2;
866     pn2 = pn1;
867     pn1 = *P;
868   }
869   PetscFunctionReturn(0);
870 }
871 
872 /* Evaluates the first derivative of P_{n}^{a,b} at a point x. */
873 PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobiDerivative(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscInt k, PetscReal *P)
874 {
875   PetscReal      nP;
876   PetscInt       i;
877   PetscErrorCode ierr;
878 
879   PetscFunctionBegin;
880   if (k > n) {*P = 0.0; PetscFunctionReturn(0);}
881   ierr = PetscDTComputeJacobi(a+k, b+k, n-k, x, &nP);CHKERRQ(ierr);
882   for (i = 0; i < k; i++) nP *= (a + b + n + 1. + i) * 0.5;
883   *P = nP;
884   PetscFunctionReturn(0);
885 }
886 
887 /* Maps from [-1,1]^2 to the (-1,1) reference triangle */
888 PETSC_STATIC_INLINE PetscErrorCode PetscDTMapSquareToTriangle_Internal(PetscReal x, PetscReal y, PetscReal *xi, PetscReal *eta)
889 {
890   PetscFunctionBegin;
891   *xi  = 0.5 * (1.0 + x) * (1.0 - y) - 1.0;
892   *eta = y;
893   PetscFunctionReturn(0);
894 }
895 
896 static PetscErrorCode PetscDTGaussJacobiQuadrature_Newton_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal x[], PetscReal w[])
897 {
898   PetscInt       maxIter = 100;
899   PetscReal      eps     = PetscExpReal(0.75 * PetscLogReal(PETSC_MACHINE_EPSILON));
900   PetscReal      a1, a6, gf;
901   PetscInt       k;
902   PetscErrorCode ierr;
903 
904   PetscFunctionBegin;
905 
906   a1 = PetscPowReal(2.0, a+b+1);
907 #if defined(PETSC_HAVE_LGAMMA)
908   {
909     PetscReal a2, a3, a4, a5;
910     a2 = PetscLGamma(a + npoints + 1);
911     a3 = PetscLGamma(b + npoints + 1);
912     a4 = PetscLGamma(a + b + npoints + 1);
913     a5 = PetscLGamma(npoints + 1);
914     gf = PetscExpReal(a2 + a3 - (a4 + a5));
915   }
916 #else
917   {
918     PetscInt ia, ib;
919 
920     ia = (PetscInt) a;
921     ib = (PetscInt) b;
922     gf = 1.;
923     if (ia == a && ia >= 0) { /* compute ratio of rising factorals wrt a */
924       for (k = 0; k < ia; k++) gf *= (npoints + 1. + k) / (npoints + b + 1. + k);
925     } else if (b == b && ib >= 0) { /* compute ratio of rising factorials wrt b */
926       for (k = 0; k < ib; k++) gf *= (npoints + 1. + k) / (npoints + a + 1. + k);
927     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"lgamma() - math routine is unavailable.");
928   }
929 #endif
930 
931   a6   = a1 * gf;
932   /* Computes the m roots of P_{m}^{a,b} on [-1,1] by Newton's method with Chebyshev points as initial guesses.
933    Algorithm implemented from the pseudocode given by Karniadakis and Sherwin and Python in FIAT */
934   for (k = 0; k < npoints; ++k) {
935     PetscReal r = PetscCosReal(PETSC_PI * (1. - (4.*k + 3. + 2.*b) / (4.*npoints + 2.*(a + b + 1.)))), dP;
936     PetscInt  j;
937 
938     if (k > 0) r = 0.5 * (r + x[k-1]);
939     for (j = 0; j < maxIter; ++j) {
940       PetscReal s = 0.0, delta, f, fp;
941       PetscInt  i;
942 
943       for (i = 0; i < k; ++i) s = s + 1.0 / (r - x[i]);
944       ierr = PetscDTComputeJacobi(a, b, npoints, r, &f);CHKERRQ(ierr);
945       ierr = PetscDTComputeJacobiDerivative(a, b, npoints, r, 1, &fp);CHKERRQ(ierr);
946       delta = f / (fp - f * s);
947       r     = r - delta;
948       if (PetscAbsReal(delta) < eps) break;
949     }
950     x[k] = r;
951     ierr = PetscDTComputeJacobiDerivative(a, b, npoints, x[k], 1, &dP);CHKERRQ(ierr);
952     w[k] = a6 / (1.0 - PetscSqr(x[k])) / PetscSqr(dP);
953   }
954   PetscFunctionReturn(0);
955 }
956 
957 /* Compute the diagonals of the Jacobi matrix used in Golub & Welsch algorithms for Gauss-Jacobi
958  * quadrature weight calculations on [-1,1] for exponents (1. + x)^a (1.-x)^b */
959 static PetscErrorCode PetscDTJacobiMatrix_Internal(PetscInt nPoints, PetscReal a, PetscReal b, PetscReal *d, PetscReal *s)
960 {
961   PetscInt       i;
962 
963   PetscFunctionBegin;
964   for (i = 0; i < nPoints; i++) {
965     PetscReal A, B, C;
966 
967     PetscDTJacobiRecurrence_Internal(i+1,a,b,A,B,C);
968     d[i] = -A / B;
969     if (i) s[i-1] *= C / B;
970     if (i < nPoints - 1) s[i] = 1. / B;
971   }
972   PetscFunctionReturn(0);
973 }
974 
975 static PetscErrorCode PetscDTGaussJacobiQuadrature_GolubWelsch_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal *x, PetscReal *w)
976 {
977   PetscReal mu0;
978   PetscReal ga, gb, gab;
979   PetscInt i;
980   PetscErrorCode ierr;
981 
982   PetscFunctionBegin;
983   ierr = PetscCitationsRegister(GolubWelschCitation, &GolubWelschCite);CHKERRQ(ierr);
984 
985 #if defined(PETSC_HAVE_TGAMMA)
986   ga  = PetscTGamma(a + 1);
987   gb  = PetscTGamma(b + 1);
988   gab = PetscTGamma(a + b + 2);
989 #else
990   {
991     PetscInt ia, ib;
992 
993     ia = (PetscInt) a;
994     ib = (PetscInt) b;
995     if (ia == a && ib == b && ia + 1 > 0 && ib + 1 > 0 && ia + ib + 2 > 0) { /* All gamma(x) terms are (x-1)! terms */
996       ierr = PetscDTFactorial(ia, &ga);CHKERRQ(ierr);
997       ierr = PetscDTFactorial(ib, &gb);CHKERRQ(ierr);
998       ierr = PetscDTFactorial(ia + ib + 1, &gb);CHKERRQ(ierr);
999     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"tgamma() - math routine is unavailable.");
1000   }
1001 #endif
1002   mu0 = PetscPowReal(2.,a + b + 1.) * ga * gb / gab;
1003 
1004 #if defined(PETSCDTGAUSSIANQUADRATURE_EIG)
1005   {
1006     PetscReal *diag, *subdiag;
1007     PetscScalar *V;
1008 
1009     ierr = PetscMalloc2(npoints, &diag, npoints, &subdiag);CHKERRQ(ierr);
1010     ierr = PetscMalloc1(npoints*npoints, &V);CHKERRQ(ierr);
1011     ierr = PetscDTJacobiMatrix_Internal(npoints, a, b, diag, subdiag);CHKERRQ(ierr);
1012     for (i = 0; i < npoints - 1; i++) subdiag[i] = PetscSqrtReal(subdiag[i]);
1013     ierr = PetscDTSymmetricTridiagonalEigensolve(npoints, diag, subdiag, x, V);CHKERRQ(ierr);
1014     for (i = 0; i < npoints; i++) w[i] = PetscSqr(PetscRealPart(V[i * npoints])) * mu0;
1015     ierr = PetscFree(V);CHKERRQ(ierr);
1016     ierr = PetscFree2(diag, subdiag);CHKERRQ(ierr);
1017   }
1018 #else
1019   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP_SYS, "A LAPACK symmetric tridiagonal eigensolver could not be found");
1020 #endif
1021   { /* As of March 2, 2020, The Sun Performance Library breaks the LAPACK contract for xstegr and xsteqr: the
1022        eigenvalues are not guaranteed to be in ascending order.  So we heave a passive aggressive sigh and check that
1023        the eigenvalues are sorted */
1024     PetscBool sorted;
1025 
1026     ierr = PetscSortedReal(npoints, x, &sorted);CHKERRQ(ierr);
1027     if (!sorted) {
1028       PetscInt *order, i;
1029       PetscReal *tmp;
1030 
1031       ierr = PetscMalloc2(npoints, &order, npoints, &tmp);CHKERRQ(ierr);
1032       for (i = 0; i < npoints; i++) order[i] = i;
1033       ierr = PetscSortRealWithPermutation(npoints, x, order);CHKERRQ(ierr);
1034       ierr = PetscArraycpy(tmp, x, npoints);CHKERRQ(ierr);
1035       for (i = 0; i < npoints; i++) x[i] = tmp[order[i]];
1036       ierr = PetscArraycpy(tmp, w, npoints);CHKERRQ(ierr);
1037       for (i = 0; i < npoints; i++) w[i] = tmp[order[i]];
1038       ierr = PetscFree2(order, tmp);CHKERRQ(ierr);
1039     }
1040   }
1041   PetscFunctionReturn(0);
1042 }
1043 
1044 static PetscErrorCode PetscDTGaussJacobiQuadrature_Internal(PetscInt npoints,PetscReal alpha, PetscReal beta, PetscReal x[], PetscReal w[], PetscBool newton)
1045 {
1046   PetscErrorCode ierr;
1047 
1048   PetscFunctionBegin;
1049   if (npoints < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of points must be positive");
1050   /* If asking for a 1D Lobatto point, just return the non-Lobatto 1D point */
1051   if (alpha <= -1.) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"alpha must be > -1.");
1052   if (beta <= -1.) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"beta must be > -1.");
1053 
1054   if (newton) {
1055     ierr = PetscDTGaussJacobiQuadrature_Newton_Internal(npoints, alpha, beta, x, w);CHKERRQ(ierr);
1056   } else {
1057     ierr = PetscDTGaussJacobiQuadrature_GolubWelsch_Internal(npoints, alpha, beta, x, w);CHKERRQ(ierr);
1058   }
1059   if (alpha == beta) { /* symmetrize */
1060     PetscInt i;
1061     for (i = 0; i < (npoints + 1) / 2; i++) {
1062       PetscInt  j  = npoints - 1 - i;
1063       PetscReal xi = x[i];
1064       PetscReal xj = x[j];
1065       PetscReal wi = w[i];
1066       PetscReal wj = w[j];
1067 
1068       x[i] = (xi - xj) / 2.;
1069       x[j] = (xj - xi) / 2.;
1070       w[i] = w[j] = (wi + wj) / 2.;
1071     }
1072   }
1073   PetscFunctionReturn(0);
1074 }
1075 
1076 /*@
1077   PetscDTGaussJacobiQuadrature - quadrature for the interval [a, b] with the weight function
1078   $(x-a)^\alpha (x-b)^\beta$.
1079 
1080   Not collective
1081 
1082   Input Parameters:
1083 + npoints - the number of points in the quadrature rule
1084 . a - the left endpoint of the interval
1085 . b - the right endpoint of the interval
1086 . alpha - the left exponent
1087 - beta - the right exponent
1088 
1089   Output Parameters:
1090 + x - array of length npoints, the locations of the quadrature points
1091 - w - array of length npoints, the weights of the quadrature points
1092 
1093   Level: intermediate
1094 
1095   Note: this quadrature rule is exact for polynomials up to degree 2*npoints - 1.
1096 @*/
1097 PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt npoints,PetscReal a, PetscReal b, PetscReal alpha, PetscReal beta, PetscReal x[], PetscReal w[])
1098 {
1099   PetscInt       i;
1100   PetscErrorCode ierr;
1101 
1102   PetscFunctionBegin;
1103   ierr = PetscDTGaussJacobiQuadrature_Internal(npoints, alpha, beta, x, w, PetscDTGaussQuadratureNewton_Internal);CHKERRQ(ierr);
1104   if (a != -1. || b != 1.) { /* shift */
1105     for (i = 0; i < npoints; i++) {
1106       x[i] = (x[i] + 1.) * ((b - a) / 2.) + a;
1107       w[i] *= (b - a) / 2.;
1108     }
1109   }
1110   PetscFunctionReturn(0);
1111 }
1112 
1113 static PetscErrorCode PetscDTGaussLobattoJacobiQuadrature_Internal(PetscInt npoints,PetscReal alpha, PetscReal beta, PetscReal x[], PetscReal w[], PetscBool newton)
1114 {
1115   PetscInt       i;
1116   PetscErrorCode ierr;
1117 
1118   PetscFunctionBegin;
1119   if (npoints < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of points must be positive");
1120   /* If asking for a 1D Lobatto point, just return the non-Lobatto 1D point */
1121   if (alpha <= -1.) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"alpha must be > -1.");
1122   if (beta <= -1.) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"beta must be > -1.");
1123 
1124   x[0] = -1.;
1125   x[npoints-1] = 1.;
1126   if (npoints > 2) {
1127     ierr = PetscDTGaussJacobiQuadrature_Internal(npoints-2, alpha+1., beta+1., &x[1], &w[1], newton);CHKERRQ(ierr);
1128   }
1129   for (i = 1; i < npoints - 1; i++) {
1130     w[i] /= (1. - x[i]*x[i]);
1131   }
1132   ierr = PetscDTGaussLobattoJacobiEndweights_Internal(npoints, alpha, beta, &w[0], &w[npoints-1]);CHKERRQ(ierr);
1133   PetscFunctionReturn(0);
1134 }
1135 
1136 /*@
1137   PetscDTGaussLobattoJacobiQuadrature - quadrature for the interval [a, b] with the weight function
1138   $(x-a)^\alpha (x-b)^\beta$, with endpoints a and b included as quadrature points.
1139 
1140   Not collective
1141 
1142   Input Parameters:
1143 + npoints - the number of points in the quadrature rule
1144 . a - the left endpoint of the interval
1145 . b - the right endpoint of the interval
1146 . alpha - the left exponent
1147 - beta - the right exponent
1148 
1149   Output Parameters:
1150 + x - array of length npoints, the locations of the quadrature points
1151 - w - array of length npoints, the weights of the quadrature points
1152 
1153   Level: intermediate
1154 
1155   Note: this quadrature rule is exact for polynomials up to degree 2*npoints - 3.
1156 @*/
1157 PetscErrorCode PetscDTGaussLobattoJacobiQuadrature(PetscInt npoints,PetscReal a, PetscReal b, PetscReal alpha, PetscReal beta, PetscReal x[], PetscReal w[])
1158 {
1159   PetscInt       i;
1160   PetscErrorCode ierr;
1161 
1162   PetscFunctionBegin;
1163   ierr = PetscDTGaussLobattoJacobiQuadrature_Internal(npoints, alpha, beta, x, w, PetscDTGaussQuadratureNewton_Internal);CHKERRQ(ierr);
1164   if (a != -1. || b != 1.) { /* shift */
1165     for (i = 0; i < npoints; i++) {
1166       x[i] = (x[i] + 1.) * ((b - a) / 2.) + a;
1167       w[i] *= (b - a) / 2.;
1168     }
1169   }
1170   PetscFunctionReturn(0);
1171 }
1172 
1173 /*@
1174    PetscDTGaussQuadrature - create Gauss-Legendre quadrature
1175 
1176    Not Collective
1177 
1178    Input Arguments:
1179 +  npoints - number of points
1180 .  a - left end of interval (often-1)
1181 -  b - right end of interval (often +1)
1182 
1183    Output Arguments:
1184 +  x - quadrature points
1185 -  w - quadrature weights
1186 
1187    Level: intermediate
1188 
1189    References:
1190 .   1. - Golub and Welsch, Calculation of Quadrature Rules, Math. Comp. 23(106), 1969.
1191 
1192 .seealso: PetscDTLegendreEval()
1193 @*/
1194 PetscErrorCode PetscDTGaussQuadrature(PetscInt npoints,PetscReal a,PetscReal b,PetscReal *x,PetscReal *w)
1195 {
1196   PetscInt       i;
1197   PetscErrorCode ierr;
1198 
1199   PetscFunctionBegin;
1200   ierr = PetscDTGaussJacobiQuadrature_Internal(npoints, 0., 0., x, w, PetscDTGaussQuadratureNewton_Internal);CHKERRQ(ierr);
1201   if (a != -1. || b != 1.) { /* shift */
1202     for (i = 0; i < npoints; i++) {
1203       x[i] = (x[i] + 1.) * ((b - a) / 2.) + a;
1204       w[i] *= (b - a) / 2.;
1205     }
1206   }
1207   PetscFunctionReturn(0);
1208 }
1209 
1210 /*@C
1211    PetscDTGaussLobattoLegendreQuadrature - creates a set of the locations and weights of the Gauss-Lobatto-Legendre
1212                       nodes of a given size on the domain [-1,1]
1213 
1214    Not Collective
1215 
1216    Input Parameter:
1217 +  n - number of grid nodes
1218 -  type - PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA or PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON
1219 
1220    Output Arguments:
1221 +  x - quadrature points
1222 -  w - quadrature weights
1223 
1224    Notes:
1225     For n > 30  the Newton approach computes duplicate (incorrect) values for some nodes because the initial guess is apparently not
1226           close enough to the desired solution
1227 
1228    These are useful for implementing spectral methods based on Gauss-Lobatto-Legendre (GLL) nodes
1229 
1230    See  https://epubs.siam.org/doi/abs/10.1137/110855442  https://epubs.siam.org/doi/abs/10.1137/120889873 for better ways to compute GLL nodes
1231 
1232    Level: intermediate
1233 
1234 .seealso: PetscDTGaussQuadrature()
1235 
1236 @*/
1237 PetscErrorCode PetscDTGaussLobattoLegendreQuadrature(PetscInt npoints,PetscGaussLobattoLegendreCreateType type,PetscReal *x,PetscReal *w)
1238 {
1239   PetscBool      newton;
1240   PetscErrorCode ierr;
1241 
1242   PetscFunctionBegin;
1243   if (npoints < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must provide at least 2 grid points per element");
1244   newton = (PetscBool) (type == PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON);
1245   ierr = PetscDTGaussLobattoJacobiQuadrature_Internal(npoints, 0., 0., x, w, newton);CHKERRQ(ierr);
1246   PetscFunctionReturn(0);
1247 }
1248 
1249 /*@
1250   PetscDTGaussTensorQuadrature - creates a tensor-product Gauss quadrature
1251 
1252   Not Collective
1253 
1254   Input Arguments:
1255 + dim     - The spatial dimension
1256 . Nc      - The number of components
1257 . npoints - number of points in one dimension
1258 . a       - left end of interval (often-1)
1259 - b       - right end of interval (often +1)
1260 
1261   Output Argument:
1262 . q - A PetscQuadrature object
1263 
1264   Level: intermediate
1265 
1266 .seealso: PetscDTGaussQuadrature(), PetscDTLegendreEval()
1267 @*/
1268 PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt dim, PetscInt Nc, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q)
1269 {
1270   PetscInt       totpoints = dim > 1 ? dim > 2 ? npoints*PetscSqr(npoints) : PetscSqr(npoints) : npoints, i, j, k, c;
1271   PetscReal     *x, *w, *xw, *ww;
1272   PetscErrorCode ierr;
1273 
1274   PetscFunctionBegin;
1275   ierr = PetscMalloc1(totpoints*dim,&x);CHKERRQ(ierr);
1276   ierr = PetscMalloc1(totpoints*Nc,&w);CHKERRQ(ierr);
1277   /* Set up the Golub-Welsch system */
1278   switch (dim) {
1279   case 0:
1280     ierr = PetscFree(x);CHKERRQ(ierr);
1281     ierr = PetscFree(w);CHKERRQ(ierr);
1282     ierr = PetscMalloc1(1, &x);CHKERRQ(ierr);
1283     ierr = PetscMalloc1(Nc, &w);CHKERRQ(ierr);
1284     x[0] = 0.0;
1285     for (c = 0; c < Nc; ++c) w[c] = 1.0;
1286     break;
1287   case 1:
1288     ierr = PetscMalloc1(npoints,&ww);CHKERRQ(ierr);
1289     ierr = PetscDTGaussQuadrature(npoints, a, b, x, ww);CHKERRQ(ierr);
1290     for (i = 0; i < npoints; ++i) for (c = 0; c < Nc; ++c) w[i*Nc+c] = ww[i];
1291     ierr = PetscFree(ww);CHKERRQ(ierr);
1292     break;
1293   case 2:
1294     ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr);
1295     ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr);
1296     for (i = 0; i < npoints; ++i) {
1297       for (j = 0; j < npoints; ++j) {
1298         x[(i*npoints+j)*dim+0] = xw[i];
1299         x[(i*npoints+j)*dim+1] = xw[j];
1300         for (c = 0; c < Nc; ++c) w[(i*npoints+j)*Nc+c] = ww[i] * ww[j];
1301       }
1302     }
1303     ierr = PetscFree2(xw,ww);CHKERRQ(ierr);
1304     break;
1305   case 3:
1306     ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr);
1307     ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr);
1308     for (i = 0; i < npoints; ++i) {
1309       for (j = 0; j < npoints; ++j) {
1310         for (k = 0; k < npoints; ++k) {
1311           x[((i*npoints+j)*npoints+k)*dim+0] = xw[i];
1312           x[((i*npoints+j)*npoints+k)*dim+1] = xw[j];
1313           x[((i*npoints+j)*npoints+k)*dim+2] = xw[k];
1314           for (c = 0; c < Nc; ++c) w[((i*npoints+j)*npoints+k)*Nc+c] = ww[i] * ww[j] * ww[k];
1315         }
1316       }
1317     }
1318     ierr = PetscFree2(xw,ww);CHKERRQ(ierr);
1319     break;
1320   default:
1321     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim);
1322   }
1323   ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr);
1324   ierr = PetscQuadratureSetOrder(*q, 2*npoints-1);CHKERRQ(ierr);
1325   ierr = PetscQuadratureSetData(*q, dim, Nc, totpoints, x, w);CHKERRQ(ierr);
1326   ierr = PetscObjectChangeTypeName((PetscObject)*q,"GaussTensor");CHKERRQ(ierr);
1327   PetscFunctionReturn(0);
1328 }
1329 
1330 /* Maps from [-1,1]^2 to the (-1,1) reference triangle */
1331 PETSC_STATIC_INLINE PetscErrorCode PetscDTMapCubeToTetrahedron_Internal(PetscReal x, PetscReal y, PetscReal z, PetscReal *xi, PetscReal *eta, PetscReal *zeta)
1332 {
1333   PetscFunctionBegin;
1334   *xi   = 0.25 * (1.0 + x) * (1.0 - y) * (1.0 - z) - 1.0;
1335   *eta  = 0.5  * (1.0 + y) * (1.0 - z) - 1.0;
1336   *zeta = z;
1337   PetscFunctionReturn(0);
1338 }
1339 
1340 
1341 /*@
1342   PetscDTStroudConicalQuadrature - create Stroud conical quadrature for a simplex
1343 
1344   Not Collective
1345 
1346   Input Arguments:
1347 + dim     - The simplex dimension
1348 . Nc      - The number of components
1349 . npoints - The number of points in one dimension
1350 . a       - left end of interval (often-1)
1351 - b       - right end of interval (often +1)
1352 
1353   Output Argument:
1354 . q - A PetscQuadrature object
1355 
1356   Level: intermediate
1357 
1358   References:
1359 .  1. - Karniadakis and Sherwin.  FIAT
1360 
1361   Note: For dim == 1, this is Gauss-Legendre quadrature
1362 
1363 .seealso: PetscDTGaussTensorQuadrature(), PetscDTGaussQuadrature()
1364 @*/
1365 PetscErrorCode PetscDTStroudConicalQuadrature(PetscInt dim, PetscInt Nc, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q)
1366 {
1367   PetscInt       totpoints = dim > 1 ? dim > 2 ? npoints*PetscSqr(npoints) : PetscSqr(npoints) : npoints;
1368   PetscReal     *px, *wx, *py, *wy, *pz, *wz, *x, *w;
1369   PetscInt       i, j, k, c; PetscErrorCode ierr;
1370 
1371   PetscFunctionBegin;
1372   if ((a != -1.0) || (b != 1.0)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must use default internal right now");
1373   ierr = PetscMalloc1(totpoints*dim, &x);CHKERRQ(ierr);
1374   ierr = PetscMalloc1(totpoints*Nc, &w);CHKERRQ(ierr);
1375   switch (dim) {
1376   case 0:
1377     ierr = PetscFree(x);CHKERRQ(ierr);
1378     ierr = PetscFree(w);CHKERRQ(ierr);
1379     ierr = PetscMalloc1(1, &x);CHKERRQ(ierr);
1380     ierr = PetscMalloc1(Nc, &w);CHKERRQ(ierr);
1381     x[0] = 0.0;
1382     for (c = 0; c < Nc; ++c) w[c] = 1.0;
1383     break;
1384   case 1:
1385     ierr = PetscMalloc1(npoints,&wx);CHKERRQ(ierr);
1386     ierr = PetscDTGaussJacobiQuadrature(npoints, -1., 1., 0.0, 0.0, x, wx);CHKERRQ(ierr);
1387     for (i = 0; i < npoints; ++i) for (c = 0; c < Nc; ++c) w[i*Nc+c] = wx[i];
1388     ierr = PetscFree(wx);CHKERRQ(ierr);
1389     break;
1390   case 2:
1391     ierr = PetscMalloc4(npoints,&px,npoints,&wx,npoints,&py,npoints,&wy);CHKERRQ(ierr);
1392     ierr = PetscDTGaussJacobiQuadrature(npoints, -1., 1., 0.0, 0.0, px, wx);CHKERRQ(ierr);
1393     ierr = PetscDTGaussJacobiQuadrature(npoints, -1., 1., 1.0, 0.0, py, wy);CHKERRQ(ierr);
1394     for (i = 0; i < npoints; ++i) {
1395       for (j = 0; j < npoints; ++j) {
1396         ierr = PetscDTMapSquareToTriangle_Internal(px[i], py[j], &x[(i*npoints+j)*2+0], &x[(i*npoints+j)*2+1]);CHKERRQ(ierr);
1397         for (c = 0; c < Nc; ++c) w[(i*npoints+j)*Nc+c] = 0.5 * wx[i] * wy[j];
1398       }
1399     }
1400     ierr = PetscFree4(px,wx,py,wy);CHKERRQ(ierr);
1401     break;
1402   case 3:
1403     ierr = PetscMalloc6(npoints,&px,npoints,&wx,npoints,&py,npoints,&wy,npoints,&pz,npoints,&wz);CHKERRQ(ierr);
1404     ierr = PetscDTGaussJacobiQuadrature(npoints, -1., 1., 0.0, 0.0, px, wx);CHKERRQ(ierr);
1405     ierr = PetscDTGaussJacobiQuadrature(npoints, -1., 1., 1.0, 0.0, py, wy);CHKERRQ(ierr);
1406     ierr = PetscDTGaussJacobiQuadrature(npoints, -1., 1., 2.0, 0.0, pz, wz);CHKERRQ(ierr);
1407     for (i = 0; i < npoints; ++i) {
1408       for (j = 0; j < npoints; ++j) {
1409         for (k = 0; k < npoints; ++k) {
1410           ierr = PetscDTMapCubeToTetrahedron_Internal(px[i], py[j], pz[k], &x[((i*npoints+j)*npoints+k)*3+0], &x[((i*npoints+j)*npoints+k)*3+1], &x[((i*npoints+j)*npoints+k)*3+2]);CHKERRQ(ierr);
1411           for (c = 0; c < Nc; ++c) w[((i*npoints+j)*npoints+k)*Nc+c] = 0.125 * wx[i] * wy[j] * wz[k];
1412         }
1413       }
1414     }
1415     ierr = PetscFree6(px,wx,py,wy,pz,wz);CHKERRQ(ierr);
1416     break;
1417   default:
1418     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim);
1419   }
1420   ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr);
1421   ierr = PetscQuadratureSetOrder(*q, 2*npoints-1);CHKERRQ(ierr);
1422   ierr = PetscQuadratureSetData(*q, dim, Nc, totpoints, x, w);CHKERRQ(ierr);
1423   ierr = PetscObjectChangeTypeName((PetscObject)*q,"GaussJacobi");CHKERRQ(ierr);
1424   PetscFunctionReturn(0);
1425 }
1426 
1427 /*@
1428   PetscDTTanhSinhTensorQuadrature - create tanh-sinh quadrature for a tensor product cell
1429 
1430   Not Collective
1431 
1432   Input Arguments:
1433 + dim   - The cell dimension
1434 . level - The number of points in one dimension, 2^l
1435 . a     - left end of interval (often-1)
1436 - b     - right end of interval (often +1)
1437 
1438   Output Argument:
1439 . q - A PetscQuadrature object
1440 
1441   Level: intermediate
1442 
1443 .seealso: PetscDTGaussTensorQuadrature()
1444 @*/
1445 PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt dim, PetscInt level, PetscReal a, PetscReal b, PetscQuadrature *q)
1446 {
1447   const PetscInt  p     = 16;                        /* Digits of precision in the evaluation */
1448   const PetscReal alpha = (b-a)/2.;                  /* Half-width of the integration interval */
1449   const PetscReal beta  = (b+a)/2.;                  /* Center of the integration interval */
1450   const PetscReal h     = PetscPowReal(2.0, -level); /* Step size, length between x_k */
1451   PetscReal       xk;                                /* Quadrature point x_k on reference domain [-1, 1] */
1452   PetscReal       wk    = 0.5*PETSC_PI;              /* Quadrature weight at x_k */
1453   PetscReal      *x, *w;
1454   PetscInt        K, k, npoints;
1455   PetscErrorCode  ierr;
1456 
1457   PetscFunctionBegin;
1458   if (dim > 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Dimension %d not yet implemented", dim);
1459   if (!level) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a number of significant digits");
1460   /* Find K such that the weights are < 32 digits of precision */
1461   for (K = 1; PetscAbsReal(PetscLog10Real(wk)) < 2*p; ++K) {
1462     wk = 0.5*h*PETSC_PI*PetscCoshReal(K*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(K*h)));
1463   }
1464   ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr);
1465   ierr = PetscQuadratureSetOrder(*q, 2*K+1);CHKERRQ(ierr);
1466   npoints = 2*K-1;
1467   ierr = PetscMalloc1(npoints*dim, &x);CHKERRQ(ierr);
1468   ierr = PetscMalloc1(npoints, &w);CHKERRQ(ierr);
1469   /* Center term */
1470   x[0] = beta;
1471   w[0] = 0.5*alpha*PETSC_PI;
1472   for (k = 1; k < K; ++k) {
1473     wk = 0.5*alpha*h*PETSC_PI*PetscCoshReal(k*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h)));
1474     xk = PetscTanhReal(0.5*PETSC_PI*PetscSinhReal(k*h));
1475     x[2*k-1] = -alpha*xk+beta;
1476     w[2*k-1] = wk;
1477     x[2*k+0] =  alpha*xk+beta;
1478     w[2*k+0] = wk;
1479   }
1480   ierr = PetscQuadratureSetData(*q, dim, 1, npoints, x, w);CHKERRQ(ierr);
1481   PetscFunctionReturn(0);
1482 }
1483 
1484 PetscErrorCode PetscDTTanhSinhIntegrate(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol)
1485 {
1486   const PetscInt  p     = 16;        /* Digits of precision in the evaluation */
1487   const PetscReal alpha = (b-a)/2.;  /* Half-width of the integration interval */
1488   const PetscReal beta  = (b+a)/2.;  /* Center of the integration interval */
1489   PetscReal       h     = 1.0;       /* Step size, length between x_k */
1490   PetscInt        l     = 0;         /* Level of refinement, h = 2^{-l} */
1491   PetscReal       osum  = 0.0;       /* Integral on last level */
1492   PetscReal       psum  = 0.0;       /* Integral on the level before the last level */
1493   PetscReal       sum;               /* Integral on current level */
1494   PetscReal       yk;                /* Quadrature point 1 - x_k on reference domain [-1, 1] */
1495   PetscReal       lx, rx;            /* Quadrature points to the left and right of 0 on the real domain [a, b] */
1496   PetscReal       wk;                /* Quadrature weight at x_k */
1497   PetscReal       lval, rval;        /* Terms in the quadature sum to the left and right of 0 */
1498   PetscInt        d;                 /* Digits of precision in the integral */
1499 
1500   PetscFunctionBegin;
1501   if (digits <= 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a positive number of significant digits");
1502   /* Center term */
1503   func(beta, &lval);
1504   sum = 0.5*alpha*PETSC_PI*lval;
1505   /* */
1506   do {
1507     PetscReal lterm, rterm, maxTerm = 0.0, d1, d2, d3, d4;
1508     PetscInt  k = 1;
1509 
1510     ++l;
1511     /* PetscPrintf(PETSC_COMM_SELF, "LEVEL %D sum: %15.15f\n", l, sum); */
1512     /* At each level of refinement, h --> h/2 and sum --> sum/2 */
1513     psum = osum;
1514     osum = sum;
1515     h   *= 0.5;
1516     sum *= 0.5;
1517     do {
1518       wk = 0.5*h*PETSC_PI*PetscCoshReal(k*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h)));
1519       yk = 1.0/(PetscExpReal(0.5*PETSC_PI*PetscSinhReal(k*h)) * PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h)));
1520       lx = -alpha*(1.0 - yk)+beta;
1521       rx =  alpha*(1.0 - yk)+beta;
1522       func(lx, &lval);
1523       func(rx, &rval);
1524       lterm   = alpha*wk*lval;
1525       maxTerm = PetscMax(PetscAbsReal(lterm), maxTerm);
1526       sum    += lterm;
1527       rterm   = alpha*wk*rval;
1528       maxTerm = PetscMax(PetscAbsReal(rterm), maxTerm);
1529       sum    += rterm;
1530       ++k;
1531       /* Only need to evaluate every other point on refined levels */
1532       if (l != 1) ++k;
1533     } while (PetscAbsReal(PetscLog10Real(wk)) < p); /* Only need to evaluate sum until weights are < 32 digits of precision */
1534 
1535     d1 = PetscLog10Real(PetscAbsReal(sum - osum));
1536     d2 = PetscLog10Real(PetscAbsReal(sum - psum));
1537     d3 = PetscLog10Real(maxTerm) - p;
1538     if (PetscMax(PetscAbsReal(lterm), PetscAbsReal(rterm)) == 0.0) d4 = 0.0;
1539     else d4 = PetscLog10Real(PetscMax(PetscAbsReal(lterm), PetscAbsReal(rterm)));
1540     d  = PetscAbsInt(PetscMin(0, PetscMax(PetscMax(PetscMax(PetscSqr(d1)/d2, 2*d1), d3), d4)));
1541   } while (d < digits && l < 12);
1542   *sol = sum;
1543 
1544   PetscFunctionReturn(0);
1545 }
1546 
1547 #if defined(PETSC_HAVE_MPFR)
1548 PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol)
1549 {
1550   const PetscInt  safetyFactor = 2;  /* Calculate abcissa until 2*p digits */
1551   PetscInt        l            = 0;  /* Level of refinement, h = 2^{-l} */
1552   mpfr_t          alpha;             /* Half-width of the integration interval */
1553   mpfr_t          beta;              /* Center of the integration interval */
1554   mpfr_t          h;                 /* Step size, length between x_k */
1555   mpfr_t          osum;              /* Integral on last level */
1556   mpfr_t          psum;              /* Integral on the level before the last level */
1557   mpfr_t          sum;               /* Integral on current level */
1558   mpfr_t          yk;                /* Quadrature point 1 - x_k on reference domain [-1, 1] */
1559   mpfr_t          lx, rx;            /* Quadrature points to the left and right of 0 on the real domain [a, b] */
1560   mpfr_t          wk;                /* Quadrature weight at x_k */
1561   PetscReal       lval, rval;        /* Terms in the quadature sum to the left and right of 0 */
1562   PetscInt        d;                 /* Digits of precision in the integral */
1563   mpfr_t          pi2, kh, msinh, mcosh, maxTerm, curTerm, tmp;
1564 
1565   PetscFunctionBegin;
1566   if (digits <= 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a positive number of significant digits");
1567   /* Create high precision storage */
1568   mpfr_inits2(PetscCeilReal(safetyFactor*digits*PetscLogReal(10.)/PetscLogReal(2.)), alpha, beta, h, sum, osum, psum, yk, wk, lx, rx, tmp, maxTerm, curTerm, pi2, kh, msinh, mcosh, NULL);
1569   /* Initialization */
1570   mpfr_set_d(alpha, 0.5*(b-a), MPFR_RNDN);
1571   mpfr_set_d(beta,  0.5*(b+a), MPFR_RNDN);
1572   mpfr_set_d(osum,  0.0,       MPFR_RNDN);
1573   mpfr_set_d(psum,  0.0,       MPFR_RNDN);
1574   mpfr_set_d(h,     1.0,       MPFR_RNDN);
1575   mpfr_const_pi(pi2, MPFR_RNDN);
1576   mpfr_mul_d(pi2, pi2, 0.5, MPFR_RNDN);
1577   /* Center term */
1578   func(0.5*(b+a), &lval);
1579   mpfr_set(sum, pi2, MPFR_RNDN);
1580   mpfr_mul(sum, sum, alpha, MPFR_RNDN);
1581   mpfr_mul_d(sum, sum, lval, MPFR_RNDN);
1582   /* */
1583   do {
1584     PetscReal d1, d2, d3, d4;
1585     PetscInt  k = 1;
1586 
1587     ++l;
1588     mpfr_set_d(maxTerm, 0.0, MPFR_RNDN);
1589     /* PetscPrintf(PETSC_COMM_SELF, "LEVEL %D sum: %15.15f\n", l, sum); */
1590     /* At each level of refinement, h --> h/2 and sum --> sum/2 */
1591     mpfr_set(psum, osum, MPFR_RNDN);
1592     mpfr_set(osum,  sum, MPFR_RNDN);
1593     mpfr_mul_d(h,   h,   0.5, MPFR_RNDN);
1594     mpfr_mul_d(sum, sum, 0.5, MPFR_RNDN);
1595     do {
1596       mpfr_set_si(kh, k, MPFR_RNDN);
1597       mpfr_mul(kh, kh, h, MPFR_RNDN);
1598       /* Weight */
1599       mpfr_set(wk, h, MPFR_RNDN);
1600       mpfr_sinh_cosh(msinh, mcosh, kh, MPFR_RNDN);
1601       mpfr_mul(msinh, msinh, pi2, MPFR_RNDN);
1602       mpfr_mul(mcosh, mcosh, pi2, MPFR_RNDN);
1603       mpfr_cosh(tmp, msinh, MPFR_RNDN);
1604       mpfr_sqr(tmp, tmp, MPFR_RNDN);
1605       mpfr_mul(wk, wk, mcosh, MPFR_RNDN);
1606       mpfr_div(wk, wk, tmp, MPFR_RNDN);
1607       /* Abscissa */
1608       mpfr_set_d(yk, 1.0, MPFR_RNDZ);
1609       mpfr_cosh(tmp, msinh, MPFR_RNDN);
1610       mpfr_div(yk, yk, tmp, MPFR_RNDZ);
1611       mpfr_exp(tmp, msinh, MPFR_RNDN);
1612       mpfr_div(yk, yk, tmp, MPFR_RNDZ);
1613       /* Quadrature points */
1614       mpfr_sub_d(lx, yk, 1.0, MPFR_RNDZ);
1615       mpfr_mul(lx, lx, alpha, MPFR_RNDU);
1616       mpfr_add(lx, lx, beta, MPFR_RNDU);
1617       mpfr_d_sub(rx, 1.0, yk, MPFR_RNDZ);
1618       mpfr_mul(rx, rx, alpha, MPFR_RNDD);
1619       mpfr_add(rx, rx, beta, MPFR_RNDD);
1620       /* Evaluation */
1621       func(mpfr_get_d(lx, MPFR_RNDU), &lval);
1622       func(mpfr_get_d(rx, MPFR_RNDD), &rval);
1623       /* Update */
1624       mpfr_mul(tmp, wk, alpha, MPFR_RNDN);
1625       mpfr_mul_d(tmp, tmp, lval, MPFR_RNDN);
1626       mpfr_add(sum, sum, tmp, MPFR_RNDN);
1627       mpfr_abs(tmp, tmp, MPFR_RNDN);
1628       mpfr_max(maxTerm, maxTerm, tmp, MPFR_RNDN);
1629       mpfr_set(curTerm, tmp, MPFR_RNDN);
1630       mpfr_mul(tmp, wk, alpha, MPFR_RNDN);
1631       mpfr_mul_d(tmp, tmp, rval, MPFR_RNDN);
1632       mpfr_add(sum, sum, tmp, MPFR_RNDN);
1633       mpfr_abs(tmp, tmp, MPFR_RNDN);
1634       mpfr_max(maxTerm, maxTerm, tmp, MPFR_RNDN);
1635       mpfr_max(curTerm, curTerm, tmp, MPFR_RNDN);
1636       ++k;
1637       /* Only need to evaluate every other point on refined levels */
1638       if (l != 1) ++k;
1639       mpfr_log10(tmp, wk, MPFR_RNDN);
1640       mpfr_abs(tmp, tmp, MPFR_RNDN);
1641     } while (mpfr_get_d(tmp, MPFR_RNDN) < safetyFactor*digits); /* Only need to evaluate sum until weights are < 32 digits of precision */
1642     mpfr_sub(tmp, sum, osum, MPFR_RNDN);
1643     mpfr_abs(tmp, tmp, MPFR_RNDN);
1644     mpfr_log10(tmp, tmp, MPFR_RNDN);
1645     d1 = mpfr_get_d(tmp, MPFR_RNDN);
1646     mpfr_sub(tmp, sum, psum, MPFR_RNDN);
1647     mpfr_abs(tmp, tmp, MPFR_RNDN);
1648     mpfr_log10(tmp, tmp, MPFR_RNDN);
1649     d2 = mpfr_get_d(tmp, MPFR_RNDN);
1650     mpfr_log10(tmp, maxTerm, MPFR_RNDN);
1651     d3 = mpfr_get_d(tmp, MPFR_RNDN) - digits;
1652     mpfr_log10(tmp, curTerm, MPFR_RNDN);
1653     d4 = mpfr_get_d(tmp, MPFR_RNDN);
1654     d  = PetscAbsInt(PetscMin(0, PetscMax(PetscMax(PetscMax(PetscSqr(d1)/d2, 2*d1), d3), d4)));
1655   } while (d < digits && l < 8);
1656   *sol = mpfr_get_d(sum, MPFR_RNDN);
1657   /* Cleanup */
1658   mpfr_clears(alpha, beta, h, sum, osum, psum, yk, wk, lx, rx, tmp, maxTerm, curTerm, pi2, kh, msinh, mcosh, NULL);
1659   PetscFunctionReturn(0);
1660 }
1661 #else
1662 
1663 PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol)
1664 {
1665   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "This method will not work without MPFR. Reconfigure using --download-mpfr --download-gmp");
1666 }
1667 #endif
1668 
1669 /* Overwrites A. Can only handle full-rank problems with m>=n
1670  * A in column-major format
1671  * Ainv in row-major format
1672  * tau has length m
1673  * worksize must be >= max(1,n)
1674  */
1675 static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m,PetscInt mstride,PetscInt n,PetscReal *A_in,PetscReal *Ainv_out,PetscScalar *tau,PetscInt worksize,PetscScalar *work)
1676 {
1677   PetscErrorCode ierr;
1678   PetscBLASInt   M,N,K,lda,ldb,ldwork,info;
1679   PetscScalar    *A,*Ainv,*R,*Q,Alpha;
1680 
1681   PetscFunctionBegin;
1682 #if defined(PETSC_USE_COMPLEX)
1683   {
1684     PetscInt i,j;
1685     ierr = PetscMalloc2(m*n,&A,m*n,&Ainv);CHKERRQ(ierr);
1686     for (j=0; j<n; j++) {
1687       for (i=0; i<m; i++) A[i+m*j] = A_in[i+mstride*j];
1688     }
1689     mstride = m;
1690   }
1691 #else
1692   A = A_in;
1693   Ainv = Ainv_out;
1694 #endif
1695 
1696   ierr = PetscBLASIntCast(m,&M);CHKERRQ(ierr);
1697   ierr = PetscBLASIntCast(n,&N);CHKERRQ(ierr);
1698   ierr = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr);
1699   ierr = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr);
1700   ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1701   PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info));
1702   ierr = PetscFPTrapPop();CHKERRQ(ierr);
1703   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error");
1704   R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */
1705 
1706   /* Extract an explicit representation of Q */
1707   Q = Ainv;
1708   ierr = PetscArraycpy(Q,A,mstride*n);CHKERRQ(ierr);
1709   K = N;                        /* full rank */
1710   PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info));
1711   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error");
1712 
1713   /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */
1714   Alpha = 1.0;
1715   ldb = lda;
1716   PetscStackCallBLAS("BLAStrsm",BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb));
1717   /* Ainv is Q, overwritten with inverse */
1718 
1719 #if defined(PETSC_USE_COMPLEX)
1720   {
1721     PetscInt i;
1722     for (i=0; i<m*n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]);
1723     ierr = PetscFree2(A,Ainv);CHKERRQ(ierr);
1724   }
1725 #endif
1726   PetscFunctionReturn(0);
1727 }
1728 
1729 /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */
1730 static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval,const PetscReal *x,PetscInt ndegree,const PetscInt *degrees,PetscBool Transpose,PetscReal *B)
1731 {
1732   PetscErrorCode ierr;
1733   PetscReal      *Bv;
1734   PetscInt       i,j;
1735 
1736   PetscFunctionBegin;
1737   ierr = PetscMalloc1((ninterval+1)*ndegree,&Bv);CHKERRQ(ierr);
1738   /* Point evaluation of L_p on all the source vertices */
1739   ierr = PetscDTLegendreEval(ninterval+1,x,ndegree,degrees,Bv,NULL,NULL);CHKERRQ(ierr);
1740   /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */
1741   for (i=0; i<ninterval; i++) {
1742     for (j=0; j<ndegree; j++) {
1743       if (Transpose) B[i+ninterval*j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j];
1744       else           B[i*ndegree+j]   = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j];
1745     }
1746   }
1747   ierr = PetscFree(Bv);CHKERRQ(ierr);
1748   PetscFunctionReturn(0);
1749 }
1750 
1751 /*@
1752    PetscDTReconstructPoly - create matrix representing polynomial reconstruction using cell intervals and evaluation at target intervals
1753 
1754    Not Collective
1755 
1756    Input Arguments:
1757 +  degree - degree of reconstruction polynomial
1758 .  nsource - number of source intervals
1759 .  sourcex - sorted coordinates of source cell boundaries (length nsource+1)
1760 .  ntarget - number of target intervals
1761 -  targetx - sorted coordinates of target cell boundaries (length ntarget+1)
1762 
1763    Output Arguments:
1764 .  R - reconstruction matrix, utarget = sum_s R[t*nsource+s] * usource[s]
1765 
1766    Level: advanced
1767 
1768 .seealso: PetscDTLegendreEval()
1769 @*/
1770 PetscErrorCode PetscDTReconstructPoly(PetscInt degree,PetscInt nsource,const PetscReal *sourcex,PetscInt ntarget,const PetscReal *targetx,PetscReal *R)
1771 {
1772   PetscErrorCode ierr;
1773   PetscInt       i,j,k,*bdegrees,worksize;
1774   PetscReal      xmin,xmax,center,hscale,*sourcey,*targety,*Bsource,*Bsinv,*Btarget;
1775   PetscScalar    *tau,*work;
1776 
1777   PetscFunctionBegin;
1778   PetscValidRealPointer(sourcex,3);
1779   PetscValidRealPointer(targetx,5);
1780   PetscValidRealPointer(R,6);
1781   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);
1782 #if defined(PETSC_USE_DEBUG)
1783   for (i=0; i<nsource; i++) {
1784     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]);
1785   }
1786   for (i=0; i<ntarget; i++) {
1787     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]);
1788   }
1789 #endif
1790   xmin = PetscMin(sourcex[0],targetx[0]);
1791   xmax = PetscMax(sourcex[nsource],targetx[ntarget]);
1792   center = (xmin + xmax)/2;
1793   hscale = (xmax - xmin)/2;
1794   worksize = nsource;
1795   ierr = PetscMalloc4(degree+1,&bdegrees,nsource+1,&sourcey,nsource*(degree+1),&Bsource,worksize,&work);CHKERRQ(ierr);
1796   ierr = PetscMalloc4(nsource,&tau,nsource*(degree+1),&Bsinv,ntarget+1,&targety,ntarget*(degree+1),&Btarget);CHKERRQ(ierr);
1797   for (i=0; i<=nsource; i++) sourcey[i] = (sourcex[i]-center)/hscale;
1798   for (i=0; i<=degree; i++) bdegrees[i] = i+1;
1799   ierr = PetscDTLegendreIntegrate(nsource,sourcey,degree+1,bdegrees,PETSC_TRUE,Bsource);CHKERRQ(ierr);
1800   ierr = PetscDTPseudoInverseQR(nsource,nsource,degree+1,Bsource,Bsinv,tau,nsource,work);CHKERRQ(ierr);
1801   for (i=0; i<=ntarget; i++) targety[i] = (targetx[i]-center)/hscale;
1802   ierr = PetscDTLegendreIntegrate(ntarget,targety,degree+1,bdegrees,PETSC_FALSE,Btarget);CHKERRQ(ierr);
1803   for (i=0; i<ntarget; i++) {
1804     PetscReal rowsum = 0;
1805     for (j=0; j<nsource; j++) {
1806       PetscReal sum = 0;
1807       for (k=0; k<degree+1; k++) {
1808         sum += Btarget[i*(degree+1)+k] * Bsinv[k*nsource+j];
1809       }
1810       R[i*nsource+j] = sum;
1811       rowsum += sum;
1812     }
1813     for (j=0; j<nsource; j++) R[i*nsource+j] /= rowsum; /* normalize each row */
1814   }
1815   ierr = PetscFree4(bdegrees,sourcey,Bsource,work);CHKERRQ(ierr);
1816   ierr = PetscFree4(tau,Bsinv,targety,Btarget);CHKERRQ(ierr);
1817   PetscFunctionReturn(0);
1818 }
1819 
1820 /*@C
1821    PetscGaussLobattoLegendreIntegrate - Compute the L2 integral of a function on the GLL points
1822 
1823    Not Collective
1824 
1825    Input Parameter:
1826 +  n - the number of GLL nodes
1827 .  nodes - the GLL nodes
1828 .  weights - the GLL weights
1829 .  f - the function values at the nodes
1830 
1831    Output Parameter:
1832 .  in - the value of the integral
1833 
1834    Level: beginner
1835 
1836 .seealso: PetscDTGaussLobattoLegendreQuadrature()
1837 
1838 @*/
1839 PetscErrorCode PetscGaussLobattoLegendreIntegrate(PetscInt n,PetscReal *nodes,PetscReal *weights,const PetscReal *f,PetscReal *in)
1840 {
1841   PetscInt          i;
1842 
1843   PetscFunctionBegin;
1844   *in = 0.;
1845   for (i=0; i<n; i++) {
1846     *in += f[i]*f[i]*weights[i];
1847   }
1848   PetscFunctionReturn(0);
1849 }
1850 
1851 /*@C
1852    PetscGaussLobattoLegendreElementLaplacianCreate - computes the Laplacian for a single 1d GLL element
1853 
1854    Not Collective
1855 
1856    Input Parameter:
1857 +  n - the number of GLL nodes
1858 .  nodes - the GLL nodes
1859 .  weights - the GLL weights
1860 
1861    Output Parameter:
1862 .  A - the stiffness element
1863 
1864    Level: beginner
1865 
1866    Notes:
1867     Destroy this with PetscGaussLobattoLegendreElementLaplacianDestroy()
1868 
1869    You can access entries in this array with AA[i][j] but in memory it is stored in contiguous memory, row oriented (the array is symmetric)
1870 
1871 .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianDestroy()
1872 
1873 @*/
1874 PetscErrorCode PetscGaussLobattoLegendreElementLaplacianCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA)
1875 {
1876   PetscReal        **A;
1877   PetscErrorCode  ierr;
1878   const PetscReal  *gllnodes = nodes;
1879   const PetscInt   p = n-1;
1880   PetscReal        z0,z1,z2 = -1,x,Lpj,Lpr;
1881   PetscInt         i,j,nn,r;
1882 
1883   PetscFunctionBegin;
1884   ierr = PetscMalloc1(n,&A);CHKERRQ(ierr);
1885   ierr = PetscMalloc1(n*n,&A[0]);CHKERRQ(ierr);
1886   for (i=1; i<n; i++) A[i] = A[i-1]+n;
1887 
1888   for (j=1; j<p; j++) {
1889     x  = gllnodes[j];
1890     z0 = 1.;
1891     z1 = x;
1892     for (nn=1; nn<p; nn++) {
1893       z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.));
1894       z0 = z1;
1895       z1 = z2;
1896     }
1897     Lpj=z2;
1898     for (r=1; r<p; r++) {
1899       if (r == j) {
1900         A[j][j]=2./(3.*(1.-gllnodes[j]*gllnodes[j])*Lpj*Lpj);
1901       } else {
1902         x  = gllnodes[r];
1903         z0 = 1.;
1904         z1 = x;
1905         for (nn=1; nn<p; nn++) {
1906           z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.));
1907           z0 = z1;
1908           z1 = z2;
1909         }
1910         Lpr     = z2;
1911         A[r][j] = 4./(((PetscReal)p)*(((PetscReal)p)+1.)*Lpj*Lpr*(gllnodes[j]-gllnodes[r])*(gllnodes[j]-gllnodes[r]));
1912       }
1913     }
1914   }
1915   for (j=1; j<p+1; j++) {
1916     x  = gllnodes[j];
1917     z0 = 1.;
1918     z1 = x;
1919     for (nn=1; nn<p; nn++) {
1920       z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.));
1921       z0 = z1;
1922       z1 = z2;
1923     }
1924     Lpj     = z2;
1925     A[j][0] = 4.*PetscPowRealInt(-1.,p)/(((PetscReal)p)*(((PetscReal)p)+1.)*Lpj*(1.+gllnodes[j])*(1.+gllnodes[j]));
1926     A[0][j] = A[j][0];
1927   }
1928   for (j=0; j<p; j++) {
1929     x  = gllnodes[j];
1930     z0 = 1.;
1931     z1 = x;
1932     for (nn=1; nn<p; nn++) {
1933       z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.));
1934       z0 = z1;
1935       z1 = z2;
1936     }
1937     Lpj=z2;
1938 
1939     A[p][j] = 4./(((PetscReal)p)*(((PetscReal)p)+1.)*Lpj*(1.-gllnodes[j])*(1.-gllnodes[j]));
1940     A[j][p] = A[p][j];
1941   }
1942   A[0][0]=0.5+(((PetscReal)p)*(((PetscReal)p)+1.)-2.)/6.;
1943   A[p][p]=A[0][0];
1944   *AA = A;
1945   PetscFunctionReturn(0);
1946 }
1947 
1948 /*@C
1949    PetscGaussLobattoLegendreElementLaplacianDestroy - frees the Laplacian for a single 1d GLL element
1950 
1951    Not Collective
1952 
1953    Input Parameter:
1954 +  n - the number of GLL nodes
1955 .  nodes - the GLL nodes
1956 .  weights - the GLL weightss
1957 -  A - the stiffness element
1958 
1959    Level: beginner
1960 
1961 .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianCreate()
1962 
1963 @*/
1964 PetscErrorCode PetscGaussLobattoLegendreElementLaplacianDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA)
1965 {
1966   PetscErrorCode ierr;
1967 
1968   PetscFunctionBegin;
1969   ierr = PetscFree((*AA)[0]);CHKERRQ(ierr);
1970   ierr = PetscFree(*AA);CHKERRQ(ierr);
1971   *AA  = NULL;
1972   PetscFunctionReturn(0);
1973 }
1974 
1975 /*@C
1976    PetscGaussLobattoLegendreElementGradientCreate - computes the gradient for a single 1d GLL element
1977 
1978    Not Collective
1979 
1980    Input Parameter:
1981 +  n - the number of GLL nodes
1982 .  nodes - the GLL nodes
1983 .  weights - the GLL weights
1984 
1985    Output Parameter:
1986 .  AA - the stiffness element
1987 -  AAT - the transpose of AA (pass in NULL if you do not need this array)
1988 
1989    Level: beginner
1990 
1991    Notes:
1992     Destroy this with PetscGaussLobattoLegendreElementGradientDestroy()
1993 
1994    You can access entries in these arrays with AA[i][j] but in memory it is stored in contiguous memory, row oriented
1995 
1996 .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianDestroy()
1997 
1998 @*/
1999 PetscErrorCode PetscGaussLobattoLegendreElementGradientCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA,PetscReal ***AAT)
2000 {
2001   PetscReal        **A, **AT = NULL;
2002   PetscErrorCode  ierr;
2003   const PetscReal  *gllnodes = nodes;
2004   const PetscInt   p = n-1;
2005   PetscReal        Li, Lj,d0;
2006   PetscInt         i,j;
2007 
2008   PetscFunctionBegin;
2009   ierr = PetscMalloc1(n,&A);CHKERRQ(ierr);
2010   ierr = PetscMalloc1(n*n,&A[0]);CHKERRQ(ierr);
2011   for (i=1; i<n; i++) A[i] = A[i-1]+n;
2012 
2013   if (AAT) {
2014     ierr = PetscMalloc1(n,&AT);CHKERRQ(ierr);
2015     ierr = PetscMalloc1(n*n,&AT[0]);CHKERRQ(ierr);
2016     for (i=1; i<n; i++) AT[i] = AT[i-1]+n;
2017   }
2018 
2019   if (n==1) {A[0][0] = 0.;}
2020   d0 = (PetscReal)p*((PetscReal)p+1.)/4.;
2021   for  (i=0; i<n; i++) {
2022     for  (j=0; j<n; j++) {
2023       A[i][j] = 0.;
2024       ierr = PetscDTComputeJacobi(0., 0., p, gllnodes[i], &Li);CHKERRQ(ierr);
2025       ierr = PetscDTComputeJacobi(0., 0., p, gllnodes[j], &Lj);CHKERRQ(ierr);
2026       if (i!=j)             A[i][j] = Li/(Lj*(gllnodes[i]-gllnodes[j]));
2027       if ((j==i) && (i==0)) A[i][j] = -d0;
2028       if (j==i && i==p)     A[i][j] = d0;
2029       if (AT) AT[j][i] = A[i][j];
2030     }
2031   }
2032   if (AAT) *AAT = AT;
2033   *AA  = A;
2034   PetscFunctionReturn(0);
2035 }
2036 
2037 /*@C
2038    PetscGaussLobattoLegendreElementGradientDestroy - frees the gradient for a single 1d GLL element obtained with PetscGaussLobattoLegendreElementGradientCreate()
2039 
2040    Not Collective
2041 
2042    Input Parameter:
2043 +  n - the number of GLL nodes
2044 .  nodes - the GLL nodes
2045 .  weights - the GLL weights
2046 .  AA - the stiffness element
2047 -  AAT - the transpose of the element
2048 
2049    Level: beginner
2050 
2051 .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianCreate(), PetscGaussLobattoLegendreElementAdvectionCreate()
2052 
2053 @*/
2054 PetscErrorCode PetscGaussLobattoLegendreElementGradientDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA,PetscReal ***AAT)
2055 {
2056   PetscErrorCode ierr;
2057 
2058   PetscFunctionBegin;
2059   ierr = PetscFree((*AA)[0]);CHKERRQ(ierr);
2060   ierr = PetscFree(*AA);CHKERRQ(ierr);
2061   *AA  = NULL;
2062   if (*AAT) {
2063     ierr = PetscFree((*AAT)[0]);CHKERRQ(ierr);
2064     ierr = PetscFree(*AAT);CHKERRQ(ierr);
2065     *AAT  = NULL;
2066   }
2067   PetscFunctionReturn(0);
2068 }
2069 
2070 /*@C
2071    PetscGaussLobattoLegendreElementAdvectionCreate - computes the advection operator for a single 1d GLL element
2072 
2073    Not Collective
2074 
2075    Input Parameter:
2076 +  n - the number of GLL nodes
2077 .  nodes - the GLL nodes
2078 .  weights - the GLL weightss
2079 
2080    Output Parameter:
2081 .  AA - the stiffness element
2082 
2083    Level: beginner
2084 
2085    Notes:
2086     Destroy this with PetscGaussLobattoLegendreElementAdvectionDestroy()
2087 
2088    This is the same as the Gradient operator multiplied by the diagonal mass matrix
2089 
2090    You can access entries in this array with AA[i][j] but in memory it is stored in contiguous memory, row oriented
2091 
2092 .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianCreate(), PetscGaussLobattoLegendreElementAdvectionDestroy()
2093 
2094 @*/
2095 PetscErrorCode PetscGaussLobattoLegendreElementAdvectionCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA)
2096 {
2097   PetscReal       **D;
2098   PetscErrorCode  ierr;
2099   const PetscReal  *gllweights = weights;
2100   const PetscInt   glln = n;
2101   PetscInt         i,j;
2102 
2103   PetscFunctionBegin;
2104   ierr = PetscGaussLobattoLegendreElementGradientCreate(n,nodes,weights,&D,NULL);CHKERRQ(ierr);
2105   for (i=0; i<glln; i++){
2106     for (j=0; j<glln; j++) {
2107       D[i][j] = gllweights[i]*D[i][j];
2108     }
2109   }
2110   *AA = D;
2111   PetscFunctionReturn(0);
2112 }
2113 
2114 /*@C
2115    PetscGaussLobattoLegendreElementAdvectionDestroy - frees the advection stiffness for a single 1d GLL element
2116 
2117    Not Collective
2118 
2119    Input Parameter:
2120 +  n - the number of GLL nodes
2121 .  nodes - the GLL nodes
2122 .  weights - the GLL weights
2123 -  A - advection
2124 
2125    Level: beginner
2126 
2127 .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementAdvectionCreate()
2128 
2129 @*/
2130 PetscErrorCode PetscGaussLobattoLegendreElementAdvectionDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA)
2131 {
2132   PetscErrorCode ierr;
2133 
2134   PetscFunctionBegin;
2135   ierr = PetscFree((*AA)[0]);CHKERRQ(ierr);
2136   ierr = PetscFree(*AA);CHKERRQ(ierr);
2137   *AA  = NULL;
2138   PetscFunctionReturn(0);
2139 }
2140 
2141 PetscErrorCode PetscGaussLobattoLegendreElementMassCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA)
2142 {
2143   PetscReal        **A;
2144   PetscErrorCode  ierr;
2145   const PetscReal  *gllweights = weights;
2146   const PetscInt   glln = n;
2147   PetscInt         i,j;
2148 
2149   PetscFunctionBegin;
2150   ierr = PetscMalloc1(glln,&A);CHKERRQ(ierr);
2151   ierr = PetscMalloc1(glln*glln,&A[0]);CHKERRQ(ierr);
2152   for (i=1; i<glln; i++) A[i] = A[i-1]+glln;
2153   if (glln==1) {A[0][0] = 0.;}
2154   for  (i=0; i<glln; i++) {
2155     for  (j=0; j<glln; j++) {
2156       A[i][j] = 0.;
2157       if (j==i)     A[i][j] = gllweights[i];
2158     }
2159   }
2160   *AA  = A;
2161   PetscFunctionReturn(0);
2162 }
2163 
2164 PetscErrorCode PetscGaussLobattoLegendreElementMassDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA)
2165 {
2166   PetscErrorCode ierr;
2167 
2168   PetscFunctionBegin;
2169   ierr = PetscFree((*AA)[0]);CHKERRQ(ierr);
2170   ierr = PetscFree(*AA);CHKERRQ(ierr);
2171   *AA  = NULL;
2172   PetscFunctionReturn(0);
2173 }
2174