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