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