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