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