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