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