xref: /petsc/src/dm/dt/interface/dt.c (revision 2cd228614d2c80ef874f383cf4c4b629b8f30ed9)
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 static PetscBool GaussCite       = PETSC_FALSE;
16 const char       GaussCitation[] = "@article{GolubWelsch1969,\n"
17                                    "  author  = {Golub and Welsch},\n"
18                                    "  title   = {Calculation of Quadrature Rules},\n"
19                                    "  journal = {Math. Comp.},\n"
20                                    "  volume  = {23},\n"
21                                    "  number  = {106},\n"
22                                    "  pages   = {221--230},\n"
23                                    "  year    = {1969}\n}\n";
24 
25 
26 PetscClassId PETSCQUADRATURE_CLASSID = 0;
27 
28 /*@
29   PetscQuadratureCreate - Create a PetscQuadrature object
30 
31   Collective
32 
33   Input Parameter:
34 . comm - The communicator for the PetscQuadrature object
35 
36   Output Parameter:
37 . q  - The PetscQuadrature object
38 
39   Level: beginner
40 
41 .seealso: PetscQuadratureDestroy(), PetscQuadratureGetData()
42 @*/
43 PetscErrorCode PetscQuadratureCreate(MPI_Comm comm, PetscQuadrature *q)
44 {
45   PetscErrorCode ierr;
46 
47   PetscFunctionBegin;
48   PetscValidPointer(q, 2);
49   ierr = DMInitializePackage();CHKERRQ(ierr);
50   ierr = PetscHeaderCreate(*q,PETSCQUADRATURE_CLASSID,"PetscQuadrature","Quadrature","DT",comm,PetscQuadratureDestroy,PetscQuadratureView);CHKERRQ(ierr);
51   (*q)->dim       = -1;
52   (*q)->Nc        =  1;
53   (*q)->order     = -1;
54   (*q)->numPoints = 0;
55   (*q)->points    = NULL;
56   (*q)->weights   = NULL;
57   PetscFunctionReturn(0);
58 }
59 
60 /*@
61   PetscQuadratureDuplicate - Create a deep copy of the PetscQuadrature object
62 
63   Collective on q
64 
65   Input Parameter:
66 . q  - The PetscQuadrature object
67 
68   Output Parameter:
69 . r  - The new PetscQuadrature object
70 
71   Level: beginner
72 
73 .seealso: PetscQuadratureCreate(), PetscQuadratureDestroy(), PetscQuadratureGetData()
74 @*/
75 PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature q, PetscQuadrature *r)
76 {
77   PetscInt         order, dim, Nc, Nq;
78   const PetscReal *points, *weights;
79   PetscReal       *p, *w;
80   PetscErrorCode   ierr;
81 
82   PetscFunctionBegin;
83   PetscValidPointer(q, 2);
84   ierr = PetscQuadratureCreate(PetscObjectComm((PetscObject) q), r);CHKERRQ(ierr);
85   ierr = PetscQuadratureGetOrder(q, &order);CHKERRQ(ierr);
86   ierr = PetscQuadratureSetOrder(*r, order);CHKERRQ(ierr);
87   ierr = PetscQuadratureGetData(q, &dim, &Nc, &Nq, &points, &weights);CHKERRQ(ierr);
88   ierr = PetscMalloc1(Nq*dim, &p);CHKERRQ(ierr);
89   ierr = PetscMalloc1(Nq*Nc, &w);CHKERRQ(ierr);
90   ierr = PetscArraycpy(p, points, Nq*dim);CHKERRQ(ierr);
91   ierr = PetscArraycpy(w, weights, Nc * Nq);CHKERRQ(ierr);
92   ierr = PetscQuadratureSetData(*r, dim, Nc, Nq, p, w);CHKERRQ(ierr);
93   PetscFunctionReturn(0);
94 }
95 
96 /*@
97   PetscQuadratureDestroy - Destroys a PetscQuadrature object
98 
99   Collective on q
100 
101   Input Parameter:
102 . q  - The PetscQuadrature object
103 
104   Level: beginner
105 
106 .seealso: PetscQuadratureCreate(), PetscQuadratureGetData()
107 @*/
108 PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *q)
109 {
110   PetscErrorCode ierr;
111 
112   PetscFunctionBegin;
113   if (!*q) PetscFunctionReturn(0);
114   PetscValidHeaderSpecific((*q),PETSCQUADRATURE_CLASSID,1);
115   if (--((PetscObject)(*q))->refct > 0) {
116     *q = NULL;
117     PetscFunctionReturn(0);
118   }
119   ierr = PetscFree((*q)->points);CHKERRQ(ierr);
120   ierr = PetscFree((*q)->weights);CHKERRQ(ierr);
121   ierr = PetscHeaderDestroy(q);CHKERRQ(ierr);
122   PetscFunctionReturn(0);
123 }
124 
125 /*@
126   PetscQuadratureGetOrder - Return the order of the method
127 
128   Not collective
129 
130   Input Parameter:
131 . q - The PetscQuadrature object
132 
133   Output Parameter:
134 . order - The order of the quadrature, i.e. the highest degree polynomial that is exactly integrated
135 
136   Level: intermediate
137 
138 .seealso: PetscQuadratureSetOrder(), PetscQuadratureGetData(), PetscQuadratureSetData()
139 @*/
140 PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature q, PetscInt *order)
141 {
142   PetscFunctionBegin;
143   PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1);
144   PetscValidPointer(order, 2);
145   *order = q->order;
146   PetscFunctionReturn(0);
147 }
148 
149 /*@
150   PetscQuadratureSetOrder - Return the order of the method
151 
152   Not collective
153 
154   Input Parameters:
155 + q - The PetscQuadrature object
156 - order - The order of the quadrature, i.e. the highest degree polynomial that is exactly integrated
157 
158   Level: intermediate
159 
160 .seealso: PetscQuadratureGetOrder(), PetscQuadratureGetData(), PetscQuadratureSetData()
161 @*/
162 PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature q, PetscInt order)
163 {
164   PetscFunctionBegin;
165   PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1);
166   q->order = order;
167   PetscFunctionReturn(0);
168 }
169 
170 /*@
171   PetscQuadratureGetNumComponents - Return the number of components for functions to be integrated
172 
173   Not collective
174 
175   Input Parameter:
176 . q - The PetscQuadrature object
177 
178   Output Parameter:
179 . Nc - The number of components
180 
181   Note: We are performing an integral int f(x) . w(x) dx, where both f and w (the weight) have Nc components.
182 
183   Level: intermediate
184 
185 .seealso: PetscQuadratureSetNumComponents(), PetscQuadratureGetData(), PetscQuadratureSetData()
186 @*/
187 PetscErrorCode PetscQuadratureGetNumComponents(PetscQuadrature q, PetscInt *Nc)
188 {
189   PetscFunctionBegin;
190   PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1);
191   PetscValidPointer(Nc, 2);
192   *Nc = q->Nc;
193   PetscFunctionReturn(0);
194 }
195 
196 /*@
197   PetscQuadratureSetNumComponents - Return the number of components for functions to be integrated
198 
199   Not collective
200 
201   Input Parameters:
202 + q  - The PetscQuadrature object
203 - Nc - The number of components
204 
205   Note: We are performing an integral int f(x) . w(x) dx, where both f and w (the weight) have Nc components.
206 
207   Level: intermediate
208 
209 .seealso: PetscQuadratureGetNumComponents(), PetscQuadratureGetData(), PetscQuadratureSetData()
210 @*/
211 PetscErrorCode PetscQuadratureSetNumComponents(PetscQuadrature q, PetscInt Nc)
212 {
213   PetscFunctionBegin;
214   PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1);
215   q->Nc = Nc;
216   PetscFunctionReturn(0);
217 }
218 
219 /*@C
220   PetscQuadratureGetData - Returns the data defining the quadrature
221 
222   Not collective
223 
224   Input Parameter:
225 . q  - The PetscQuadrature object
226 
227   Output Parameters:
228 + dim - The spatial dimension
229 . Nc - The number of components
230 . npoints - The number of quadrature points
231 . points - The coordinates of each quadrature point
232 - weights - The weight of each quadrature point
233 
234   Level: intermediate
235 
236   Fortran Notes:
237     From Fortran you must call PetscQuadratureRestoreData() when you are done with the data
238 
239 .seealso: PetscQuadratureCreate(), PetscQuadratureSetData()
240 @*/
241 PetscErrorCode PetscQuadratureGetData(PetscQuadrature q, PetscInt *dim, PetscInt *Nc, PetscInt *npoints, const PetscReal *points[], const PetscReal *weights[])
242 {
243   PetscFunctionBegin;
244   PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1);
245   if (dim) {
246     PetscValidPointer(dim, 2);
247     *dim = q->dim;
248   }
249   if (Nc) {
250     PetscValidPointer(Nc, 3);
251     *Nc = q->Nc;
252   }
253   if (npoints) {
254     PetscValidPointer(npoints, 4);
255     *npoints = q->numPoints;
256   }
257   if (points) {
258     PetscValidPointer(points, 5);
259     *points = q->points;
260   }
261   if (weights) {
262     PetscValidPointer(weights, 6);
263     *weights = q->weights;
264   }
265   PetscFunctionReturn(0);
266 }
267 
268 /*@C
269   PetscQuadratureSetData - Sets the data defining the quadrature
270 
271   Not collective
272 
273   Input Parameters:
274 + q  - The PetscQuadrature object
275 . dim - The spatial dimension
276 . Nc - The number of components
277 . npoints - The number of quadrature points
278 . points - The coordinates of each quadrature point
279 - weights - The weight of each quadrature point
280 
281   Note: This routine owns the references to points and weights, so they must be allocated using PetscMalloc() and the user should not free them.
282 
283   Level: intermediate
284 
285 .seealso: PetscQuadratureCreate(), PetscQuadratureGetData()
286 @*/
287 PetscErrorCode PetscQuadratureSetData(PetscQuadrature q, PetscInt dim, PetscInt Nc, PetscInt npoints, const PetscReal points[], const PetscReal weights[])
288 {
289   PetscFunctionBegin;
290   PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1);
291   if (dim >= 0)     q->dim       = dim;
292   if (Nc >= 0)      q->Nc        = Nc;
293   if (npoints >= 0) q->numPoints = npoints;
294   if (points) {
295     PetscValidPointer(points, 4);
296     q->points = points;
297   }
298   if (weights) {
299     PetscValidPointer(weights, 5);
300     q->weights = weights;
301   }
302   PetscFunctionReturn(0);
303 }
304 
305 static PetscErrorCode PetscQuadratureView_Ascii(PetscQuadrature quad, PetscViewer v)
306 {
307   PetscInt          q, d, c;
308   PetscViewerFormat format;
309   PetscErrorCode    ierr;
310 
311   PetscFunctionBegin;
312   if (quad->Nc > 1) {ierr = PetscViewerASCIIPrintf(v, "Quadrature of order %D on %D points (dim %D) with %D components\n", quad->order, quad->numPoints, quad->dim, quad->Nc);CHKERRQ(ierr);}
313   else              {ierr = PetscViewerASCIIPrintf(v, "Quadrature of order %D on %D points (dim %D)\n", quad->order, quad->numPoints, quad->dim);CHKERRQ(ierr);}
314   ierr = PetscViewerGetFormat(v, &format);CHKERRQ(ierr);
315   if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) PetscFunctionReturn(0);
316   for (q = 0; q < quad->numPoints; ++q) {
317     ierr = PetscViewerASCIIPrintf(v, "p%D (", q);CHKERRQ(ierr);
318     ierr = PetscViewerASCIIUseTabs(v, PETSC_FALSE);CHKERRQ(ierr);
319     for (d = 0; d < quad->dim; ++d) {
320       if (d) {ierr = PetscViewerASCIIPrintf(v, ", ");CHKERRQ(ierr);}
321       ierr = PetscViewerASCIIPrintf(v, "%+g", (double)quad->points[q*quad->dim+d]);CHKERRQ(ierr);
322     }
323     ierr = PetscViewerASCIIPrintf(v, ") ");CHKERRQ(ierr);
324     if (quad->Nc > 1) {ierr = PetscViewerASCIIPrintf(v, "w%D (", q);CHKERRQ(ierr);}
325     for (c = 0; c < quad->Nc; ++c) {
326       if (c) {ierr = PetscViewerASCIIPrintf(v, ", ");CHKERRQ(ierr);}
327       ierr = PetscViewerASCIIPrintf(v, "%+g", (double)quad->weights[q*quad->Nc+c]);CHKERRQ(ierr);
328     }
329     if (quad->Nc > 1) {ierr = PetscViewerASCIIPrintf(v, ")");CHKERRQ(ierr);}
330     ierr = PetscViewerASCIIPrintf(v, "\n");CHKERRQ(ierr);
331     ierr = PetscViewerASCIIUseTabs(v, PETSC_TRUE);CHKERRQ(ierr);
332   }
333   PetscFunctionReturn(0);
334 }
335 
336 /*@C
337   PetscQuadratureView - Views a PetscQuadrature object
338 
339   Collective on quad
340 
341   Input Parameters:
342 + quad  - The PetscQuadrature object
343 - viewer - The PetscViewer object
344 
345   Level: beginner
346 
347 .seealso: PetscQuadratureCreate(), PetscQuadratureGetData()
348 @*/
349 PetscErrorCode PetscQuadratureView(PetscQuadrature quad, PetscViewer viewer)
350 {
351   PetscBool      iascii;
352   PetscErrorCode ierr;
353 
354   PetscFunctionBegin;
355   PetscValidHeader(quad, 1);
356   if (viewer) PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
357   if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) quad), &viewer);CHKERRQ(ierr);}
358   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
359   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
360   if (iascii) {ierr = PetscQuadratureView_Ascii(quad, viewer);CHKERRQ(ierr);}
361   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
362   PetscFunctionReturn(0);
363 }
364 
365 /*@C
366   PetscQuadratureExpandComposite - Return a quadrature over the composite element, which has the original quadrature in each subelement
367 
368   Not collective
369 
370   Input Parameter:
371 + q - The original PetscQuadrature
372 . numSubelements - The number of subelements the original element is divided into
373 . v0 - An array of the initial points for each subelement
374 - jac - An array of the Jacobian mappings from the reference to each subelement
375 
376   Output Parameters:
377 . dim - The dimension
378 
379   Note: Together v0 and jac define an affine mapping from the original reference element to each subelement
380 
381  Not available from Fortran
382 
383   Level: intermediate
384 
385 .seealso: PetscFECreate(), PetscSpaceGetDimension(), PetscDualSpaceGetDimension()
386 @*/
387 PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature q, PetscInt numSubelements, const PetscReal v0[], const PetscReal jac[], PetscQuadrature *qref)
388 {
389   const PetscReal *points,    *weights;
390   PetscReal       *pointsRef, *weightsRef;
391   PetscInt         dim, Nc, order, npoints, npointsRef, c, p, cp, d, e;
392   PetscErrorCode   ierr;
393 
394   PetscFunctionBegin;
395   PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 1);
396   PetscValidPointer(v0, 3);
397   PetscValidPointer(jac, 4);
398   PetscValidPointer(qref, 5);
399   ierr = PetscQuadratureCreate(PETSC_COMM_SELF, qref);CHKERRQ(ierr);
400   ierr = PetscQuadratureGetOrder(q, &order);CHKERRQ(ierr);
401   ierr = PetscQuadratureGetData(q, &dim, &Nc, &npoints, &points, &weights);CHKERRQ(ierr);
402   npointsRef = npoints*numSubelements;
403   ierr = PetscMalloc1(npointsRef*dim,&pointsRef);CHKERRQ(ierr);
404   ierr = PetscMalloc1(npointsRef*Nc, &weightsRef);CHKERRQ(ierr);
405   for (c = 0; c < numSubelements; ++c) {
406     for (p = 0; p < npoints; ++p) {
407       for (d = 0; d < dim; ++d) {
408         pointsRef[(c*npoints + p)*dim+d] = v0[c*dim+d];
409         for (e = 0; e < dim; ++e) {
410           pointsRef[(c*npoints + p)*dim+d] += jac[(c*dim + d)*dim+e]*(points[p*dim+e] + 1.0);
411         }
412       }
413       /* Could also use detJ here */
414       for (cp = 0; cp < Nc; ++cp) weightsRef[(c*npoints+p)*Nc+cp] = weights[p*Nc+cp]/numSubelements;
415     }
416   }
417   ierr = PetscQuadratureSetOrder(*qref, order);CHKERRQ(ierr);
418   ierr = PetscQuadratureSetData(*qref, dim, Nc, npointsRef, pointsRef, weightsRef);CHKERRQ(ierr);
419   PetscFunctionReturn(0);
420 }
421 
422 /*@
423    PetscDTLegendreEval - evaluate Legendre polynomial at points
424 
425    Not Collective
426 
427    Input Arguments:
428 +  npoints - number of spatial points to evaluate at
429 .  points - array of locations to evaluate at
430 .  ndegree - number of basis degrees to evaluate
431 -  degrees - sorted array of degrees to evaluate
432 
433    Output Arguments:
434 +  B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL)
435 .  D - row-oriented derivative evaluation matrix (or NULL)
436 -  D2 - row-oriented second derivative evaluation matrix (or NULL)
437 
438    Level: intermediate
439 
440 .seealso: PetscDTGaussQuadrature()
441 @*/
442 PetscErrorCode PetscDTLegendreEval(PetscInt npoints,const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2)
443 {
444   PetscInt i,maxdegree;
445 
446   PetscFunctionBegin;
447   if (!npoints || !ndegree) PetscFunctionReturn(0);
448   maxdegree = degrees[ndegree-1];
449   for (i=0; i<npoints; i++) {
450     PetscReal pm1,pm2,pd1,pd2,pdd1,pdd2,x;
451     PetscInt  j,k;
452     x    = points[i];
453     pm2  = 0;
454     pm1  = 1;
455     pd2  = 0;
456     pd1  = 0;
457     pdd2 = 0;
458     pdd1 = 0;
459     k    = 0;
460     if (degrees[k] == 0) {
461       if (B) B[i*ndegree+k] = pm1;
462       if (D) D[i*ndegree+k] = pd1;
463       if (D2) D2[i*ndegree+k] = pdd1;
464       k++;
465     }
466     for (j=1; j<=maxdegree; j++,k++) {
467       PetscReal p,d,dd;
468       p    = ((2*j-1)*x*pm1 - (j-1)*pm2)/j;
469       d    = pd2 + (2*j-1)*pm1;
470       dd   = pdd2 + (2*j-1)*pd1;
471       pm2  = pm1;
472       pm1  = p;
473       pd2  = pd1;
474       pd1  = d;
475       pdd2 = pdd1;
476       pdd1 = dd;
477       if (degrees[k] == j) {
478         if (B) B[i*ndegree+k] = p;
479         if (D) D[i*ndegree+k] = d;
480         if (D2) D2[i*ndegree+k] = dd;
481       }
482     }
483   }
484   PetscFunctionReturn(0);
485 }
486 
487 /*@
488    PetscDTGaussQuadrature - create Gauss quadrature
489 
490    Not Collective
491 
492    Input Arguments:
493 +  npoints - number of points
494 .  a - left end of interval (often-1)
495 -  b - right end of interval (often +1)
496 
497    Output Arguments:
498 +  x - quadrature points
499 -  w - quadrature weights
500 
501    Level: intermediate
502 
503    References:
504 .   1. - Golub and Welsch, Calculation of Quadrature Rules, Math. Comp. 23(106), 1969.
505 
506 .seealso: PetscDTLegendreEval()
507 @*/
508 PetscErrorCode PetscDTGaussQuadrature(PetscInt npoints,PetscReal a,PetscReal b,PetscReal *x,PetscReal *w)
509 {
510   PetscErrorCode ierr;
511   PetscInt       i;
512   PetscReal      *work;
513   PetscScalar    *Z;
514   PetscBLASInt   N,LDZ,info;
515 
516   PetscFunctionBegin;
517   ierr = PetscCitationsRegister(GaussCitation, &GaussCite);CHKERRQ(ierr);
518   /* Set up the Golub-Welsch system */
519   for (i=0; i<npoints; i++) {
520     x[i] = 0;                   /* diagonal is 0 */
521     if (i) w[i-1] = 0.5 / PetscSqrtReal(1 - 1./PetscSqr(2*i));
522   }
523   ierr = PetscMalloc2(npoints*npoints,&Z,PetscMax(1,2*npoints-2),&work);CHKERRQ(ierr);
524   ierr = PetscBLASIntCast(npoints,&N);CHKERRQ(ierr);
525   LDZ  = N;
526   ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
527   PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("I",&N,x,w,Z,&LDZ,work,&info));
528   ierr = PetscFPTrapPop();CHKERRQ(ierr);
529   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEQR error");
530 
531   for (i=0; i<(npoints+1)/2; i++) {
532     PetscReal y = 0.5 * (-x[i] + x[npoints-i-1]); /* enforces symmetry */
533     x[i]           = (a+b)/2 - y*(b-a)/2;
534     if (x[i] == -0.0) x[i] = 0.0;
535     x[npoints-i-1] = (a+b)/2 + y*(b-a)/2;
536 
537     w[i] = w[npoints-1-i] = 0.5*(b-a)*(PetscSqr(PetscAbsScalar(Z[i*npoints])) + PetscSqr(PetscAbsScalar(Z[(npoints-i-1)*npoints])));
538   }
539   ierr = PetscFree2(Z,work);CHKERRQ(ierr);
540   PetscFunctionReturn(0);
541 }
542 
543 static void qAndLEvaluation(PetscInt n, PetscReal x, PetscReal *q, PetscReal *qp, PetscReal *Ln)
544 /*
545   Compute the polynomial q(x) = L_{N+1}(x) - L_{n-1}(x) and its derivative in
546   addition to L_N(x) as these are needed for computing the GLL points via Newton's method.
547   Reference: "Implementing Spectral Methods for Partial Differential Equations: Algorithms
548   for Scientists and Engineers" by David A. Kopriva.
549 */
550 {
551   PetscInt k;
552 
553   PetscReal Lnp;
554   PetscReal Lnp1, Lnp1p;
555   PetscReal Lnm1, Lnm1p;
556   PetscReal Lnm2, Lnm2p;
557 
558   Lnm1  = 1.0;
559   *Ln   = x;
560   Lnm1p = 0.0;
561   Lnp   = 1.0;
562 
563   for (k=2; k<=n; ++k) {
564     Lnm2  = Lnm1;
565     Lnm1  = *Ln;
566     Lnm2p = Lnm1p;
567     Lnm1p = Lnp;
568     *Ln   = (2.*((PetscReal)k)-1.)/(1.0*((PetscReal)k))*x*Lnm1 - (((PetscReal)k)-1.)/((PetscReal)k)*Lnm2;
569     Lnp   = Lnm2p + (2.0*((PetscReal)k)-1.)*Lnm1;
570   }
571   k     = n+1;
572   Lnp1  = (2.*((PetscReal)k)-1.)/(((PetscReal)k))*x*(*Ln) - (((PetscReal)k)-1.)/((PetscReal)k)*Lnm1;
573   Lnp1p = Lnm1p + (2.0*((PetscReal)k)-1.)*(*Ln);
574   *q    = Lnp1 - Lnm1;
575   *qp   = Lnp1p - Lnm1p;
576 }
577 
578 /*@C
579    PetscDTGaussLobattoLegendreQuadrature - creates a set of the locations and weights of the Gauss-Lobatto-Legendre
580                       nodes of a given size on the domain [-1,1]
581 
582    Not Collective
583 
584    Input Parameter:
585 +  n - number of grid nodes
586 -  type - PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA or PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON
587 
588    Output Arguments:
589 +  x - quadrature points
590 -  w - quadrature weights
591 
592    Notes:
593     For n > 30  the Newton approach computes duplicate (incorrect) values for some nodes because the initial guess is apparently not
594           close enough to the desired solution
595 
596    These are useful for implementing spectral methods based on Gauss-Lobatto-Legendre (GLL) nodes
597 
598    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
599 
600    Level: intermediate
601 
602 .seealso: PetscDTGaussQuadrature()
603 
604 @*/
605 PetscErrorCode PetscDTGaussLobattoLegendreQuadrature(PetscInt npoints,PetscGaussLobattoLegendreCreateType type,PetscReal *x,PetscReal *w)
606 {
607   PetscErrorCode ierr;
608 
609   PetscFunctionBegin;
610   if (npoints < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must provide at least 2 grid points per element");
611 
612   if (type == PETSCGAUSSLOBATTOLEGENDRE_VIA_LINEAR_ALGEBRA) {
613     PetscReal      *M,si;
614     PetscBLASInt   bn,lierr;
615     PetscReal      x0,z0,z1,z2;
616     PetscInt       i,p = npoints - 1,nn;
617 
618     x[0]   =-1.0;
619     x[npoints-1] = 1.0;
620     if (npoints-2 > 0){
621       ierr = PetscMalloc1(npoints-1,&M);CHKERRQ(ierr);
622       for (i=0; i<npoints-2; i++) {
623         si  = ((PetscReal)i)+1.0;
624         M[i]=0.5*PetscSqrtReal(si*(si+2.0)/((si+0.5)*(si+1.5)));
625       }
626       ierr = PetscBLASIntCast(npoints-2,&bn);CHKERRQ(ierr);
627       ierr = PetscArrayzero(&x[1],bn);CHKERRQ(ierr);
628       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
629       x0=0;
630       PetscStackCallBLAS("LAPACKsteqr",LAPACKREALsteqr_("N",&bn,&x[1],M,&x0,&bn,M,&lierr));
631       if (lierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in STERF Lapack routine %d",(int)lierr);
632       ierr = PetscFPTrapPop();CHKERRQ(ierr);
633       ierr = PetscFree(M);CHKERRQ(ierr);
634     }
635     if ((npoints-1)%2==0) {
636       x[(npoints-1)/2]   = 0.0; /* hard wire to exactly 0.0 since linear algebra produces nonzero */
637     }
638 
639     w[0] = w[p] = 2.0/(((PetscReal)(p))*(((PetscReal)p)+1.0));
640     z2 = -1.;                      /* Dummy value to avoid -Wmaybe-initialized */
641     for (i=1; i<p; i++) {
642       x0  = x[i];
643       z0 = 1.0;
644       z1 = x0;
645       for (nn=1; nn<p; nn++) {
646         z2 = x0*z1*(2.0*((PetscReal)nn)+1.0)/(((PetscReal)nn)+1.0)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.0));
647         z0 = z1;
648         z1 = z2;
649       }
650       w[i]=2.0/(((PetscReal)p)*(((PetscReal)p)+1.0)*z2*z2);
651     }
652   } else {
653     PetscInt  j,m;
654     PetscReal z1,z,q,qp,Ln;
655     PetscReal *pt;
656     ierr = PetscMalloc1(npoints,&pt);CHKERRQ(ierr);
657 
658     if (npoints > 30) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"PETSCGAUSSLOBATTOLEGENDRE_VIA_NEWTON produces incorrect answers for n > 30");
659     x[0]     = -1.0;
660     x[npoints-1]   = 1.0;
661     w[0]   = w[npoints-1] = 2./(((PetscReal)npoints)*(((PetscReal)npoints)-1.0));
662     m  = (npoints-1)/2; /* The roots are symmetric, so we only find half of them. */
663     for (j=1; j<=m; j++) { /* Loop over the desired roots. */
664       z = -1.0*PetscCosReal((PETSC_PI*((PetscReal)j)+0.25)/(((PetscReal)npoints)-1.0))-(3.0/(8.0*(((PetscReal)npoints)-1.0)*PETSC_PI))*(1.0/(((PetscReal)j)+0.25));
665       /* Starting with the above approximation to the ith root, we enter */
666       /* the main loop of refinement by Newton's method.                 */
667       do {
668         qAndLEvaluation(npoints-1,z,&q,&qp,&Ln);
669         z1 = z;
670         z  = z1-q/qp; /* Newton's method. */
671       } while (PetscAbs(z-z1) > 10.*PETSC_MACHINE_EPSILON);
672       qAndLEvaluation(npoints-1,z,&q,&qp,&Ln);
673 
674       x[j]       = z;
675       x[npoints-1-j]   = -z;      /* and put in its symmetric counterpart.   */
676       w[j]     = 2.0/(((PetscReal)npoints)*(((PetscReal)npoints)-1.)*Ln*Ln);  /* Compute the weight */
677       w[npoints-1-j] = w[j];                 /* and its symmetric counterpart. */
678       pt[j]=qp;
679     }
680 
681     if ((npoints-1)%2==0) {
682       qAndLEvaluation(npoints-1,0.0,&q,&qp,&Ln);
683       x[(npoints-1)/2]   = 0.0;
684       w[(npoints-1)/2] = 2.0/(((PetscReal)npoints)*(((PetscReal)npoints)-1.)*Ln*Ln);
685     }
686     ierr = PetscFree(pt);CHKERRQ(ierr);
687   }
688   PetscFunctionReturn(0);
689 }
690 
691 /*@
692   PetscDTGaussTensorQuadrature - creates a tensor-product Gauss quadrature
693 
694   Not Collective
695 
696   Input Arguments:
697 + dim     - The spatial dimension
698 . Nc      - The number of components
699 . npoints - number of points in one dimension
700 . a       - left end of interval (often-1)
701 - b       - right end of interval (often +1)
702 
703   Output Argument:
704 . q - A PetscQuadrature object
705 
706   Level: intermediate
707 
708 .seealso: PetscDTGaussQuadrature(), PetscDTLegendreEval()
709 @*/
710 PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt dim, PetscInt Nc, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q)
711 {
712   PetscInt       totpoints = dim > 1 ? dim > 2 ? npoints*PetscSqr(npoints) : PetscSqr(npoints) : npoints, i, j, k, c;
713   PetscReal     *x, *w, *xw, *ww;
714   PetscErrorCode ierr;
715 
716   PetscFunctionBegin;
717   ierr = PetscMalloc1(totpoints*dim,&x);CHKERRQ(ierr);
718   ierr = PetscMalloc1(totpoints*Nc,&w);CHKERRQ(ierr);
719   /* Set up the Golub-Welsch system */
720   switch (dim) {
721   case 0:
722     ierr = PetscFree(x);CHKERRQ(ierr);
723     ierr = PetscFree(w);CHKERRQ(ierr);
724     ierr = PetscMalloc1(1, &x);CHKERRQ(ierr);
725     ierr = PetscMalloc1(Nc, &w);CHKERRQ(ierr);
726     x[0] = 0.0;
727     for (c = 0; c < Nc; ++c) w[c] = 1.0;
728     break;
729   case 1:
730     ierr = PetscMalloc1(npoints,&ww);CHKERRQ(ierr);
731     ierr = PetscDTGaussQuadrature(npoints, a, b, x, ww);CHKERRQ(ierr);
732     for (i = 0; i < npoints; ++i) for (c = 0; c < Nc; ++c) w[i*Nc+c] = ww[i];
733     ierr = PetscFree(ww);CHKERRQ(ierr);
734     break;
735   case 2:
736     ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr);
737     ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr);
738     for (i = 0; i < npoints; ++i) {
739       for (j = 0; j < npoints; ++j) {
740         x[(i*npoints+j)*dim+0] = xw[i];
741         x[(i*npoints+j)*dim+1] = xw[j];
742         for (c = 0; c < Nc; ++c) w[(i*npoints+j)*Nc+c] = ww[i] * ww[j];
743       }
744     }
745     ierr = PetscFree2(xw,ww);CHKERRQ(ierr);
746     break;
747   case 3:
748     ierr = PetscMalloc2(npoints,&xw,npoints,&ww);CHKERRQ(ierr);
749     ierr = PetscDTGaussQuadrature(npoints, a, b, xw, ww);CHKERRQ(ierr);
750     for (i = 0; i < npoints; ++i) {
751       for (j = 0; j < npoints; ++j) {
752         for (k = 0; k < npoints; ++k) {
753           x[((i*npoints+j)*npoints+k)*dim+0] = xw[i];
754           x[((i*npoints+j)*npoints+k)*dim+1] = xw[j];
755           x[((i*npoints+j)*npoints+k)*dim+2] = xw[k];
756           for (c = 0; c < Nc; ++c) w[((i*npoints+j)*npoints+k)*Nc+c] = ww[i] * ww[j] * ww[k];
757         }
758       }
759     }
760     ierr = PetscFree2(xw,ww);CHKERRQ(ierr);
761     break;
762   default:
763     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim);
764   }
765   ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr);
766   ierr = PetscQuadratureSetOrder(*q, 2*npoints-1);CHKERRQ(ierr);
767   ierr = PetscQuadratureSetData(*q, dim, Nc, totpoints, x, w);CHKERRQ(ierr);
768   ierr = PetscObjectChangeTypeName((PetscObject)*q,"GaussTensor");CHKERRQ(ierr);
769   PetscFunctionReturn(0);
770 }
771 
772 /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x.
773    Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */
774 PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorial_Internal(PetscInt n, PetscReal *factorial)
775 {
776   PetscReal f = 1.0;
777   PetscInt  i;
778 
779   PetscFunctionBegin;
780   for (i = 1; i < n+1; ++i) f *= i;
781   *factorial = f;
782   PetscFunctionReturn(0);
783 }
784 
785 /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x.
786    Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */
787 PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobi(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P)
788 {
789   PetscReal apb, pn1, pn2;
790   PetscInt  k;
791 
792   PetscFunctionBegin;
793   if (!n) {*P = 1.0; PetscFunctionReturn(0);}
794   if (n == 1) {*P = 0.5 * (a - b + (a + b + 2.0) * x); PetscFunctionReturn(0);}
795   apb = a + b;
796   pn2 = 1.0;
797   pn1 = 0.5 * (a - b + (apb + 2.0) * x);
798   *P  = 0.0;
799   for (k = 2; k < n+1; ++k) {
800     PetscReal a1 = 2.0 * k * (k + apb) * (2.0*k + apb - 2.0);
801     PetscReal a2 = (2.0 * k + apb - 1.0) * (a*a - b*b);
802     PetscReal a3 = (2.0 * k + apb - 2.0) * (2.0 * k + apb - 1.0) * (2.0 * k + apb);
803     PetscReal a4 = 2.0 * (k + a - 1.0) * (k + b - 1.0) * (2.0 * k + apb);
804 
805     a2  = a2 / a1;
806     a3  = a3 / a1;
807     a4  = a4 / a1;
808     *P  = (a2 + a3 * x) * pn1 - a4 * pn2;
809     pn2 = pn1;
810     pn1 = *P;
811   }
812   PetscFunctionReturn(0);
813 }
814 
815 /* Evaluates the first derivative of P_{n}^{a,b} at a point x. */
816 PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobiDerivative(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P)
817 {
818   PetscReal      nP;
819   PetscErrorCode ierr;
820 
821   PetscFunctionBegin;
822   if (!n) {*P = 0.0; PetscFunctionReturn(0);}
823   ierr = PetscDTComputeJacobi(a+1, b+1, n-1, x, &nP);CHKERRQ(ierr);
824   *P   = 0.5 * (a + b + n + 1) * nP;
825   PetscFunctionReturn(0);
826 }
827 
828 /* Maps from [-1,1]^2 to the (-1,1) reference triangle */
829 PETSC_STATIC_INLINE PetscErrorCode PetscDTMapSquareToTriangle_Internal(PetscReal x, PetscReal y, PetscReal *xi, PetscReal *eta)
830 {
831   PetscFunctionBegin;
832   *xi  = 0.5 * (1.0 + x) * (1.0 - y) - 1.0;
833   *eta = y;
834   PetscFunctionReturn(0);
835 }
836 
837 /* Maps from [-1,1]^2 to the (-1,1) reference triangle */
838 PETSC_STATIC_INLINE PetscErrorCode PetscDTMapCubeToTetrahedron_Internal(PetscReal x, PetscReal y, PetscReal z, PetscReal *xi, PetscReal *eta, PetscReal *zeta)
839 {
840   PetscFunctionBegin;
841   *xi   = 0.25 * (1.0 + x) * (1.0 - y) * (1.0 - z) - 1.0;
842   *eta  = 0.5  * (1.0 + y) * (1.0 - z) - 1.0;
843   *zeta = z;
844   PetscFunctionReturn(0);
845 }
846 
847 static PetscErrorCode PetscDTGaussJacobiQuadrature1D_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal *x, PetscReal *w)
848 {
849   PetscInt       maxIter = 100;
850   PetscReal      eps     = 1.0e-8;
851   PetscReal      a1, a2, a3, a4, a5, a6;
852   PetscInt       k;
853   PetscErrorCode ierr;
854 
855   PetscFunctionBegin;
856 
857   a1      = PetscPowReal(2.0, a+b+1);
858 #if defined(PETSC_HAVE_TGAMMA)
859   a2      = PetscTGamma(a + npoints + 1);
860   a3      = PetscTGamma(b + npoints + 1);
861   a4      = PetscTGamma(a + b + npoints + 1);
862 #else
863   {
864     PetscInt ia, ib;
865 
866     ia = (PetscInt) a;
867     ib = (PetscInt) b;
868     if (ia == a && ib == b && ia + npoints + 1 > 0 && ib + npoints + 1 > 0 && ia + ib + npoints + 1 > 0) { /* All gamma(x) terms are (x-1)! terms */
869       ierr = PetscDTFactorial_Internal(ia + npoints, &a2);CHKERRQ(ierr);
870       ierr = PetscDTFactorial_Internal(ib + npoints, &a3);CHKERRQ(ierr);
871       ierr = PetscDTFactorial_Internal(ia + ib + npoints, &a4);CHKERRQ(ierr);
872     } else {
873       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"tgamma() - math routine is unavailable.");
874     }
875   }
876 #endif
877 
878   ierr = PetscDTFactorial_Internal(npoints, &a5);CHKERRQ(ierr);
879   a6   = a1 * a2 * a3 / a4 / a5;
880   /* Computes the m roots of P_{m}^{a,b} on [-1,1] by Newton's method with Chebyshev points as initial guesses.
881    Algorithm implemented from the pseudocode given by Karniadakis and Sherwin and Python in FIAT */
882   for (k = 0; k < npoints; ++k) {
883     PetscReal r = -PetscCosReal((2.0*k + 1.0) * PETSC_PI / (2.0 * npoints)), dP;
884     PetscInt  j;
885 
886     if (k > 0) r = 0.5 * (r + x[k-1]);
887     for (j = 0; j < maxIter; ++j) {
888       PetscReal s = 0.0, delta, f, fp;
889       PetscInt  i;
890 
891       for (i = 0; i < k; ++i) s = s + 1.0 / (r - x[i]);
892       ierr = PetscDTComputeJacobi(a, b, npoints, r, &f);CHKERRQ(ierr);
893       ierr = PetscDTComputeJacobiDerivative(a, b, npoints, r, &fp);CHKERRQ(ierr);
894       delta = f / (fp - f * s);
895       r     = r - delta;
896       if (PetscAbsReal(delta) < eps) break;
897     }
898     x[k] = r;
899     ierr = PetscDTComputeJacobiDerivative(a, b, npoints, x[k], &dP);CHKERRQ(ierr);
900     w[k] = a6 / (1.0 - PetscSqr(x[k])) / PetscSqr(dP);
901   }
902   PetscFunctionReturn(0);
903 }
904 
905 /*@
906   PetscDTGaussJacobiQuadrature - create Gauss-Jacobi quadrature for a simplex
907 
908   Not Collective
909 
910   Input Arguments:
911 + dim     - The simplex dimension
912 . Nc      - The number of components
913 . npoints - The number of points in one dimension
914 . a       - left end of interval (often-1)
915 - b       - right end of interval (often +1)
916 
917   Output Argument:
918 . q - A PetscQuadrature object
919 
920   Level: intermediate
921 
922   References:
923 .  1. - Karniadakis and Sherwin.  FIAT
924 
925 .seealso: PetscDTGaussTensorQuadrature(), PetscDTGaussQuadrature()
926 @*/
927 PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt dim, PetscInt Nc, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q)
928 {
929   PetscInt       totpoints = dim > 1 ? dim > 2 ? npoints*PetscSqr(npoints) : PetscSqr(npoints) : npoints;
930   PetscReal     *px, *wx, *py, *wy, *pz, *wz, *x, *w;
931   PetscInt       i, j, k, c;
932   PetscErrorCode ierr;
933 
934   PetscFunctionBegin;
935   if ((a != -1.0) || (b != 1.0)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must use default internal right now");
936   ierr = PetscMalloc1(totpoints*dim, &x);CHKERRQ(ierr);
937   ierr = PetscMalloc1(totpoints*Nc, &w);CHKERRQ(ierr);
938   switch (dim) {
939   case 0:
940     ierr = PetscFree(x);CHKERRQ(ierr);
941     ierr = PetscFree(w);CHKERRQ(ierr);
942     ierr = PetscMalloc1(1, &x);CHKERRQ(ierr);
943     ierr = PetscMalloc1(Nc, &w);CHKERRQ(ierr);
944     x[0] = 0.0;
945     for (c = 0; c < Nc; ++c) w[c] = 1.0;
946     break;
947   case 1:
948     ierr = PetscMalloc1(npoints,&wx);CHKERRQ(ierr);
949     ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, x, wx);CHKERRQ(ierr);
950     for (i = 0; i < npoints; ++i) for (c = 0; c < Nc; ++c) w[i*Nc+c] = wx[i];
951     ierr = PetscFree(wx);CHKERRQ(ierr);
952     break;
953   case 2:
954     ierr = PetscMalloc4(npoints,&px,npoints,&wx,npoints,&py,npoints,&wy);CHKERRQ(ierr);
955     ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, px, wx);CHKERRQ(ierr);
956     ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 1.0, 0.0, py, wy);CHKERRQ(ierr);
957     for (i = 0; i < npoints; ++i) {
958       for (j = 0; j < npoints; ++j) {
959         ierr = PetscDTMapSquareToTriangle_Internal(px[i], py[j], &x[(i*npoints+j)*2+0], &x[(i*npoints+j)*2+1]);CHKERRQ(ierr);
960         for (c = 0; c < Nc; ++c) w[(i*npoints+j)*Nc+c] = 0.5 * wx[i] * wy[j];
961       }
962     }
963     ierr = PetscFree4(px,wx,py,wy);CHKERRQ(ierr);
964     break;
965   case 3:
966     ierr = PetscMalloc6(npoints,&px,npoints,&wx,npoints,&py,npoints,&wy,npoints,&pz,npoints,&wz);CHKERRQ(ierr);
967     ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, px, wx);CHKERRQ(ierr);
968     ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 1.0, 0.0, py, wy);CHKERRQ(ierr);
969     ierr = PetscDTGaussJacobiQuadrature1D_Internal(npoints, 2.0, 0.0, pz, wz);CHKERRQ(ierr);
970     for (i = 0; i < npoints; ++i) {
971       for (j = 0; j < npoints; ++j) {
972         for (k = 0; k < npoints; ++k) {
973           ierr = PetscDTMapCubeToTetrahedron_Internal(px[i], py[j], pz[k], &x[((i*npoints+j)*npoints+k)*3+0], &x[((i*npoints+j)*npoints+k)*3+1], &x[((i*npoints+j)*npoints+k)*3+2]);CHKERRQ(ierr);
974           for (c = 0; c < Nc; ++c) w[((i*npoints+j)*npoints+k)*Nc+c] = 0.125 * wx[i] * wy[j] * wz[k];
975         }
976       }
977     }
978     ierr = PetscFree6(px,wx,py,wy,pz,wz);CHKERRQ(ierr);
979     break;
980   default:
981     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim);
982   }
983   ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr);
984   ierr = PetscQuadratureSetOrder(*q, 2*npoints-1);CHKERRQ(ierr);
985   ierr = PetscQuadratureSetData(*q, dim, Nc, totpoints, x, w);CHKERRQ(ierr);
986   ierr = PetscObjectChangeTypeName((PetscObject)*q,"GaussJacobi");CHKERRQ(ierr);
987   PetscFunctionReturn(0);
988 }
989 
990 /*@
991   PetscDTTanhSinhTensorQuadrature - create tanh-sinh quadrature for a tensor product cell
992 
993   Not Collective
994 
995   Input Arguments:
996 + dim   - The cell dimension
997 . level - The number of points in one dimension, 2^l
998 . a     - left end of interval (often-1)
999 - b     - right end of interval (often +1)
1000 
1001   Output Argument:
1002 . q - A PetscQuadrature object
1003 
1004   Level: intermediate
1005 
1006 .seealso: PetscDTGaussTensorQuadrature()
1007 @*/
1008 PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt dim, PetscInt level, PetscReal a, PetscReal b, PetscQuadrature *q)
1009 {
1010   const PetscInt  p     = 16;                        /* Digits of precision in the evaluation */
1011   const PetscReal alpha = (b-a)/2.;                  /* Half-width of the integration interval */
1012   const PetscReal beta  = (b+a)/2.;                  /* Center of the integration interval */
1013   const PetscReal h     = PetscPowReal(2.0, -level); /* Step size, length between x_k */
1014   PetscReal       xk;                                /* Quadrature point x_k on reference domain [-1, 1] */
1015   PetscReal       wk    = 0.5*PETSC_PI;              /* Quadrature weight at x_k */
1016   PetscReal      *x, *w;
1017   PetscInt        K, k, npoints;
1018   PetscErrorCode  ierr;
1019 
1020   PetscFunctionBegin;
1021   if (dim > 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Dimension %d not yet implemented", dim);
1022   if (!level) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a number of significant digits");
1023   /* Find K such that the weights are < 32 digits of precision */
1024   for (K = 1; PetscAbsReal(PetscLog10Real(wk)) < 2*p; ++K) {
1025     wk = 0.5*h*PETSC_PI*PetscCoshReal(K*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(K*h)));
1026   }
1027   ierr = PetscQuadratureCreate(PETSC_COMM_SELF, q);CHKERRQ(ierr);
1028   ierr = PetscQuadratureSetOrder(*q, 2*K+1);CHKERRQ(ierr);
1029   npoints = 2*K-1;
1030   ierr = PetscMalloc1(npoints*dim, &x);CHKERRQ(ierr);
1031   ierr = PetscMalloc1(npoints, &w);CHKERRQ(ierr);
1032   /* Center term */
1033   x[0] = beta;
1034   w[0] = 0.5*alpha*PETSC_PI;
1035   for (k = 1; k < K; ++k) {
1036     wk = 0.5*alpha*h*PETSC_PI*PetscCoshReal(k*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h)));
1037     xk = PetscTanhReal(0.5*PETSC_PI*PetscSinhReal(k*h));
1038     x[2*k-1] = -alpha*xk+beta;
1039     w[2*k-1] = wk;
1040     x[2*k+0] =  alpha*xk+beta;
1041     w[2*k+0] = wk;
1042   }
1043   ierr = PetscQuadratureSetData(*q, dim, 1, npoints, x, w);CHKERRQ(ierr);
1044   PetscFunctionReturn(0);
1045 }
1046 
1047 PetscErrorCode PetscDTTanhSinhIntegrate(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol)
1048 {
1049   const PetscInt  p     = 16;        /* Digits of precision in the evaluation */
1050   const PetscReal alpha = (b-a)/2.;  /* Half-width of the integration interval */
1051   const PetscReal beta  = (b+a)/2.;  /* Center of the integration interval */
1052   PetscReal       h     = 1.0;       /* Step size, length between x_k */
1053   PetscInt        l     = 0;         /* Level of refinement, h = 2^{-l} */
1054   PetscReal       osum  = 0.0;       /* Integral on last level */
1055   PetscReal       psum  = 0.0;       /* Integral on the level before the last level */
1056   PetscReal       sum;               /* Integral on current level */
1057   PetscReal       yk;                /* Quadrature point 1 - x_k on reference domain [-1, 1] */
1058   PetscReal       lx, rx;            /* Quadrature points to the left and right of 0 on the real domain [a, b] */
1059   PetscReal       wk;                /* Quadrature weight at x_k */
1060   PetscReal       lval, rval;        /* Terms in the quadature sum to the left and right of 0 */
1061   PetscInt        d;                 /* Digits of precision in the integral */
1062 
1063   PetscFunctionBegin;
1064   if (digits <= 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a positive number of significant digits");
1065   /* Center term */
1066   func(beta, &lval);
1067   sum = 0.5*alpha*PETSC_PI*lval;
1068   /* */
1069   do {
1070     PetscReal lterm, rterm, maxTerm = 0.0, d1, d2, d3, d4;
1071     PetscInt  k = 1;
1072 
1073     ++l;
1074     /* PetscPrintf(PETSC_COMM_SELF, "LEVEL %D sum: %15.15f\n", l, sum); */
1075     /* At each level of refinement, h --> h/2 and sum --> sum/2 */
1076     psum = osum;
1077     osum = sum;
1078     h   *= 0.5;
1079     sum *= 0.5;
1080     do {
1081       wk = 0.5*h*PETSC_PI*PetscCoshReal(k*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h)));
1082       yk = 1.0/(PetscExpReal(0.5*PETSC_PI*PetscSinhReal(k*h)) * PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h)));
1083       lx = -alpha*(1.0 - yk)+beta;
1084       rx =  alpha*(1.0 - yk)+beta;
1085       func(lx, &lval);
1086       func(rx, &rval);
1087       lterm   = alpha*wk*lval;
1088       maxTerm = PetscMax(PetscAbsReal(lterm), maxTerm);
1089       sum    += lterm;
1090       rterm   = alpha*wk*rval;
1091       maxTerm = PetscMax(PetscAbsReal(rterm), maxTerm);
1092       sum    += rterm;
1093       ++k;
1094       /* Only need to evaluate every other point on refined levels */
1095       if (l != 1) ++k;
1096     } while (PetscAbsReal(PetscLog10Real(wk)) < p); /* Only need to evaluate sum until weights are < 32 digits of precision */
1097 
1098     d1 = PetscLog10Real(PetscAbsReal(sum - osum));
1099     d2 = PetscLog10Real(PetscAbsReal(sum - psum));
1100     d3 = PetscLog10Real(maxTerm) - p;
1101     if (PetscMax(PetscAbsReal(lterm), PetscAbsReal(rterm)) == 0.0) d4 = 0.0;
1102     else d4 = PetscLog10Real(PetscMax(PetscAbsReal(lterm), PetscAbsReal(rterm)));
1103     d  = PetscAbsInt(PetscMin(0, PetscMax(PetscMax(PetscMax(PetscSqr(d1)/d2, 2*d1), d3), d4)));
1104   } while (d < digits && l < 12);
1105   *sol = sum;
1106 
1107   PetscFunctionReturn(0);
1108 }
1109 
1110 #if defined(PETSC_HAVE_MPFR)
1111 PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol)
1112 {
1113   const PetscInt  safetyFactor = 2;  /* Calculate abcissa until 2*p digits */
1114   PetscInt        l            = 0;  /* Level of refinement, h = 2^{-l} */
1115   mpfr_t          alpha;             /* Half-width of the integration interval */
1116   mpfr_t          beta;              /* Center of the integration interval */
1117   mpfr_t          h;                 /* Step size, length between x_k */
1118   mpfr_t          osum;              /* Integral on last level */
1119   mpfr_t          psum;              /* Integral on the level before the last level */
1120   mpfr_t          sum;               /* Integral on current level */
1121   mpfr_t          yk;                /* Quadrature point 1 - x_k on reference domain [-1, 1] */
1122   mpfr_t          lx, rx;            /* Quadrature points to the left and right of 0 on the real domain [a, b] */
1123   mpfr_t          wk;                /* Quadrature weight at x_k */
1124   PetscReal       lval, rval;        /* Terms in the quadature sum to the left and right of 0 */
1125   PetscInt        d;                 /* Digits of precision in the integral */
1126   mpfr_t          pi2, kh, msinh, mcosh, maxTerm, curTerm, tmp;
1127 
1128   PetscFunctionBegin;
1129   if (digits <= 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a positive number of significant digits");
1130   /* Create high precision storage */
1131   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);
1132   /* Initialization */
1133   mpfr_set_d(alpha, 0.5*(b-a), MPFR_RNDN);
1134   mpfr_set_d(beta,  0.5*(b+a), MPFR_RNDN);
1135   mpfr_set_d(osum,  0.0,       MPFR_RNDN);
1136   mpfr_set_d(psum,  0.0,       MPFR_RNDN);
1137   mpfr_set_d(h,     1.0,       MPFR_RNDN);
1138   mpfr_const_pi(pi2, MPFR_RNDN);
1139   mpfr_mul_d(pi2, pi2, 0.5, MPFR_RNDN);
1140   /* Center term */
1141   func(0.5*(b+a), &lval);
1142   mpfr_set(sum, pi2, MPFR_RNDN);
1143   mpfr_mul(sum, sum, alpha, MPFR_RNDN);
1144   mpfr_mul_d(sum, sum, lval, MPFR_RNDN);
1145   /* */
1146   do {
1147     PetscReal d1, d2, d3, d4;
1148     PetscInt  k = 1;
1149 
1150     ++l;
1151     mpfr_set_d(maxTerm, 0.0, MPFR_RNDN);
1152     /* PetscPrintf(PETSC_COMM_SELF, "LEVEL %D sum: %15.15f\n", l, sum); */
1153     /* At each level of refinement, h --> h/2 and sum --> sum/2 */
1154     mpfr_set(psum, osum, MPFR_RNDN);
1155     mpfr_set(osum,  sum, MPFR_RNDN);
1156     mpfr_mul_d(h,   h,   0.5, MPFR_RNDN);
1157     mpfr_mul_d(sum, sum, 0.5, MPFR_RNDN);
1158     do {
1159       mpfr_set_si(kh, k, MPFR_RNDN);
1160       mpfr_mul(kh, kh, h, MPFR_RNDN);
1161       /* Weight */
1162       mpfr_set(wk, h, MPFR_RNDN);
1163       mpfr_sinh_cosh(msinh, mcosh, kh, MPFR_RNDN);
1164       mpfr_mul(msinh, msinh, pi2, MPFR_RNDN);
1165       mpfr_mul(mcosh, mcosh, pi2, MPFR_RNDN);
1166       mpfr_cosh(tmp, msinh, MPFR_RNDN);
1167       mpfr_sqr(tmp, tmp, MPFR_RNDN);
1168       mpfr_mul(wk, wk, mcosh, MPFR_RNDN);
1169       mpfr_div(wk, wk, tmp, MPFR_RNDN);
1170       /* Abscissa */
1171       mpfr_set_d(yk, 1.0, MPFR_RNDZ);
1172       mpfr_cosh(tmp, msinh, MPFR_RNDN);
1173       mpfr_div(yk, yk, tmp, MPFR_RNDZ);
1174       mpfr_exp(tmp, msinh, MPFR_RNDN);
1175       mpfr_div(yk, yk, tmp, MPFR_RNDZ);
1176       /* Quadrature points */
1177       mpfr_sub_d(lx, yk, 1.0, MPFR_RNDZ);
1178       mpfr_mul(lx, lx, alpha, MPFR_RNDU);
1179       mpfr_add(lx, lx, beta, MPFR_RNDU);
1180       mpfr_d_sub(rx, 1.0, yk, MPFR_RNDZ);
1181       mpfr_mul(rx, rx, alpha, MPFR_RNDD);
1182       mpfr_add(rx, rx, beta, MPFR_RNDD);
1183       /* Evaluation */
1184       func(mpfr_get_d(lx, MPFR_RNDU), &lval);
1185       func(mpfr_get_d(rx, MPFR_RNDD), &rval);
1186       /* Update */
1187       mpfr_mul(tmp, wk, alpha, MPFR_RNDN);
1188       mpfr_mul_d(tmp, tmp, lval, MPFR_RNDN);
1189       mpfr_add(sum, sum, tmp, MPFR_RNDN);
1190       mpfr_abs(tmp, tmp, MPFR_RNDN);
1191       mpfr_max(maxTerm, maxTerm, tmp, MPFR_RNDN);
1192       mpfr_set(curTerm, tmp, MPFR_RNDN);
1193       mpfr_mul(tmp, wk, alpha, MPFR_RNDN);
1194       mpfr_mul_d(tmp, tmp, rval, MPFR_RNDN);
1195       mpfr_add(sum, sum, tmp, MPFR_RNDN);
1196       mpfr_abs(tmp, tmp, MPFR_RNDN);
1197       mpfr_max(maxTerm, maxTerm, tmp, MPFR_RNDN);
1198       mpfr_max(curTerm, curTerm, tmp, MPFR_RNDN);
1199       ++k;
1200       /* Only need to evaluate every other point on refined levels */
1201       if (l != 1) ++k;
1202       mpfr_log10(tmp, wk, MPFR_RNDN);
1203       mpfr_abs(tmp, tmp, MPFR_RNDN);
1204     } while (mpfr_get_d(tmp, MPFR_RNDN) < safetyFactor*digits); /* Only need to evaluate sum until weights are < 32 digits of precision */
1205     mpfr_sub(tmp, sum, osum, MPFR_RNDN);
1206     mpfr_abs(tmp, tmp, MPFR_RNDN);
1207     mpfr_log10(tmp, tmp, MPFR_RNDN);
1208     d1 = mpfr_get_d(tmp, MPFR_RNDN);
1209     mpfr_sub(tmp, sum, psum, MPFR_RNDN);
1210     mpfr_abs(tmp, tmp, MPFR_RNDN);
1211     mpfr_log10(tmp, tmp, MPFR_RNDN);
1212     d2 = mpfr_get_d(tmp, MPFR_RNDN);
1213     mpfr_log10(tmp, maxTerm, MPFR_RNDN);
1214     d3 = mpfr_get_d(tmp, MPFR_RNDN) - digits;
1215     mpfr_log10(tmp, curTerm, MPFR_RNDN);
1216     d4 = mpfr_get_d(tmp, MPFR_RNDN);
1217     d  = PetscAbsInt(PetscMin(0, PetscMax(PetscMax(PetscMax(PetscSqr(d1)/d2, 2*d1), d3), d4)));
1218   } while (d < digits && l < 8);
1219   *sol = mpfr_get_d(sum, MPFR_RNDN);
1220   /* Cleanup */
1221   mpfr_clears(alpha, beta, h, sum, osum, psum, yk, wk, lx, rx, tmp, maxTerm, curTerm, pi2, kh, msinh, mcosh, NULL);
1222   PetscFunctionReturn(0);
1223 }
1224 #else
1225 
1226 PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol)
1227 {
1228   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "This method will not work without MPFR. Reconfigure using --download-mpfr --download-gmp");
1229 }
1230 #endif
1231 
1232 /* Overwrites A. Can only handle full-rank problems with m>=n
1233  * A in column-major format
1234  * Ainv in row-major format
1235  * tau has length m
1236  * worksize must be >= max(1,n)
1237  */
1238 static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m,PetscInt mstride,PetscInt n,PetscReal *A_in,PetscReal *Ainv_out,PetscScalar *tau,PetscInt worksize,PetscScalar *work)
1239 {
1240   PetscErrorCode ierr;
1241   PetscBLASInt   M,N,K,lda,ldb,ldwork,info;
1242   PetscScalar    *A,*Ainv,*R,*Q,Alpha;
1243 
1244   PetscFunctionBegin;
1245 #if defined(PETSC_USE_COMPLEX)
1246   {
1247     PetscInt i,j;
1248     ierr = PetscMalloc2(m*n,&A,m*n,&Ainv);CHKERRQ(ierr);
1249     for (j=0; j<n; j++) {
1250       for (i=0; i<m; i++) A[i+m*j] = A_in[i+mstride*j];
1251     }
1252     mstride = m;
1253   }
1254 #else
1255   A = A_in;
1256   Ainv = Ainv_out;
1257 #endif
1258 
1259   ierr = PetscBLASIntCast(m,&M);CHKERRQ(ierr);
1260   ierr = PetscBLASIntCast(n,&N);CHKERRQ(ierr);
1261   ierr = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr);
1262   ierr = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr);
1263   ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1264   PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info));
1265   ierr = PetscFPTrapPop();CHKERRQ(ierr);
1266   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error");
1267   R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */
1268 
1269   /* Extract an explicit representation of Q */
1270   Q = Ainv;
1271   ierr = PetscArraycpy(Q,A,mstride*n);CHKERRQ(ierr);
1272   K = N;                        /* full rank */
1273   PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info));
1274   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error");
1275 
1276   /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */
1277   Alpha = 1.0;
1278   ldb = lda;
1279   PetscStackCallBLAS("BLAStrsm",BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb));
1280   /* Ainv is Q, overwritten with inverse */
1281 
1282 #if defined(PETSC_USE_COMPLEX)
1283   {
1284     PetscInt i;
1285     for (i=0; i<m*n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]);
1286     ierr = PetscFree2(A,Ainv);CHKERRQ(ierr);
1287   }
1288 #endif
1289   PetscFunctionReturn(0);
1290 }
1291 
1292 /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */
1293 static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval,const PetscReal *x,PetscInt ndegree,const PetscInt *degrees,PetscBool Transpose,PetscReal *B)
1294 {
1295   PetscErrorCode ierr;
1296   PetscReal      *Bv;
1297   PetscInt       i,j;
1298 
1299   PetscFunctionBegin;
1300   ierr = PetscMalloc1((ninterval+1)*ndegree,&Bv);CHKERRQ(ierr);
1301   /* Point evaluation of L_p on all the source vertices */
1302   ierr = PetscDTLegendreEval(ninterval+1,x,ndegree,degrees,Bv,NULL,NULL);CHKERRQ(ierr);
1303   /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */
1304   for (i=0; i<ninterval; i++) {
1305     for (j=0; j<ndegree; j++) {
1306       if (Transpose) B[i+ninterval*j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j];
1307       else           B[i*ndegree+j]   = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j];
1308     }
1309   }
1310   ierr = PetscFree(Bv);CHKERRQ(ierr);
1311   PetscFunctionReturn(0);
1312 }
1313 
1314 /*@
1315    PetscDTReconstructPoly - create matrix representing polynomial reconstruction using cell intervals and evaluation at target intervals
1316 
1317    Not Collective
1318 
1319    Input Arguments:
1320 +  degree - degree of reconstruction polynomial
1321 .  nsource - number of source intervals
1322 .  sourcex - sorted coordinates of source cell boundaries (length nsource+1)
1323 .  ntarget - number of target intervals
1324 -  targetx - sorted coordinates of target cell boundaries (length ntarget+1)
1325 
1326    Output Arguments:
1327 .  R - reconstruction matrix, utarget = sum_s R[t*nsource+s] * usource[s]
1328 
1329    Level: advanced
1330 
1331 .seealso: PetscDTLegendreEval()
1332 @*/
1333 PetscErrorCode PetscDTReconstructPoly(PetscInt degree,PetscInt nsource,const PetscReal *sourcex,PetscInt ntarget,const PetscReal *targetx,PetscReal *R)
1334 {
1335   PetscErrorCode ierr;
1336   PetscInt       i,j,k,*bdegrees,worksize;
1337   PetscReal      xmin,xmax,center,hscale,*sourcey,*targety,*Bsource,*Bsinv,*Btarget;
1338   PetscScalar    *tau,*work;
1339 
1340   PetscFunctionBegin;
1341   PetscValidRealPointer(sourcex,3);
1342   PetscValidRealPointer(targetx,5);
1343   PetscValidRealPointer(R,6);
1344   if (degree >= nsource) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Reconstruction degree %D must be less than number of source intervals %D",degree,nsource);
1345 #if defined(PETSC_USE_DEBUG)
1346   for (i=0; i<nsource; i++) {
1347     if (sourcex[i] >= sourcex[i+1]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Source interval %D has negative orientation (%g,%g)",i,(double)sourcex[i],(double)sourcex[i+1]);
1348   }
1349   for (i=0; i<ntarget; i++) {
1350     if (targetx[i] >= targetx[i+1]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Target interval %D has negative orientation (%g,%g)",i,(double)targetx[i],(double)targetx[i+1]);
1351   }
1352 #endif
1353   xmin = PetscMin(sourcex[0],targetx[0]);
1354   xmax = PetscMax(sourcex[nsource],targetx[ntarget]);
1355   center = (xmin + xmax)/2;
1356   hscale = (xmax - xmin)/2;
1357   worksize = nsource;
1358   ierr = PetscMalloc4(degree+1,&bdegrees,nsource+1,&sourcey,nsource*(degree+1),&Bsource,worksize,&work);CHKERRQ(ierr);
1359   ierr = PetscMalloc4(nsource,&tau,nsource*(degree+1),&Bsinv,ntarget+1,&targety,ntarget*(degree+1),&Btarget);CHKERRQ(ierr);
1360   for (i=0; i<=nsource; i++) sourcey[i] = (sourcex[i]-center)/hscale;
1361   for (i=0; i<=degree; i++) bdegrees[i] = i+1;
1362   ierr = PetscDTLegendreIntegrate(nsource,sourcey,degree+1,bdegrees,PETSC_TRUE,Bsource);CHKERRQ(ierr);
1363   ierr = PetscDTPseudoInverseQR(nsource,nsource,degree+1,Bsource,Bsinv,tau,nsource,work);CHKERRQ(ierr);
1364   for (i=0; i<=ntarget; i++) targety[i] = (targetx[i]-center)/hscale;
1365   ierr = PetscDTLegendreIntegrate(ntarget,targety,degree+1,bdegrees,PETSC_FALSE,Btarget);CHKERRQ(ierr);
1366   for (i=0; i<ntarget; i++) {
1367     PetscReal rowsum = 0;
1368     for (j=0; j<nsource; j++) {
1369       PetscReal sum = 0;
1370       for (k=0; k<degree+1; k++) {
1371         sum += Btarget[i*(degree+1)+k] * Bsinv[k*nsource+j];
1372       }
1373       R[i*nsource+j] = sum;
1374       rowsum += sum;
1375     }
1376     for (j=0; j<nsource; j++) R[i*nsource+j] /= rowsum; /* normalize each row */
1377   }
1378   ierr = PetscFree4(bdegrees,sourcey,Bsource,work);CHKERRQ(ierr);
1379   ierr = PetscFree4(tau,Bsinv,targety,Btarget);CHKERRQ(ierr);
1380   PetscFunctionReturn(0);
1381 }
1382 
1383 /*@C
1384    PetscGaussLobattoLegendreIntegrate - Compute the L2 integral of a function on the GLL points
1385 
1386    Not Collective
1387 
1388    Input Parameter:
1389 +  n - the number of GLL nodes
1390 .  nodes - the GLL nodes
1391 .  weights - the GLL weights
1392 .  f - the function values at the nodes
1393 
1394    Output Parameter:
1395 .  in - the value of the integral
1396 
1397    Level: beginner
1398 
1399 .seealso: PetscDTGaussLobattoLegendreQuadrature()
1400 
1401 @*/
1402 PetscErrorCode PetscGaussLobattoLegendreIntegrate(PetscInt n,PetscReal *nodes,PetscReal *weights,const PetscReal *f,PetscReal *in)
1403 {
1404   PetscInt          i;
1405 
1406   PetscFunctionBegin;
1407   *in = 0.;
1408   for (i=0; i<n; i++) {
1409     *in += f[i]*f[i]*weights[i];
1410   }
1411   PetscFunctionReturn(0);
1412 }
1413 
1414 /*@C
1415    PetscGaussLobattoLegendreElementLaplacianCreate - computes the Laplacian for a single 1d GLL element
1416 
1417    Not Collective
1418 
1419    Input Parameter:
1420 +  n - the number of GLL nodes
1421 .  nodes - the GLL nodes
1422 .  weights - the GLL weights
1423 
1424    Output Parameter:
1425 .  A - the stiffness element
1426 
1427    Level: beginner
1428 
1429    Notes:
1430     Destroy this with PetscGaussLobattoLegendreElementLaplacianDestroy()
1431 
1432    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)
1433 
1434 .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianDestroy()
1435 
1436 @*/
1437 PetscErrorCode PetscGaussLobattoLegendreElementLaplacianCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA)
1438 {
1439   PetscReal        **A;
1440   PetscErrorCode  ierr;
1441   const PetscReal  *gllnodes = nodes;
1442   const PetscInt   p = n-1;
1443   PetscReal        z0,z1,z2 = -1,x,Lpj,Lpr;
1444   PetscInt         i,j,nn,r;
1445 
1446   PetscFunctionBegin;
1447   ierr = PetscMalloc1(n,&A);CHKERRQ(ierr);
1448   ierr = PetscMalloc1(n*n,&A[0]);CHKERRQ(ierr);
1449   for (i=1; i<n; i++) A[i] = A[i-1]+n;
1450 
1451   for (j=1; j<p; j++) {
1452     x  = gllnodes[j];
1453     z0 = 1.;
1454     z1 = x;
1455     for (nn=1; nn<p; nn++) {
1456       z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.));
1457       z0 = z1;
1458       z1 = z2;
1459     }
1460     Lpj=z2;
1461     for (r=1; r<p; r++) {
1462       if (r == j) {
1463         A[j][j]=2./(3.*(1.-gllnodes[j]*gllnodes[j])*Lpj*Lpj);
1464       } else {
1465         x  = gllnodes[r];
1466         z0 = 1.;
1467         z1 = x;
1468         for (nn=1; nn<p; nn++) {
1469           z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.));
1470           z0 = z1;
1471           z1 = z2;
1472         }
1473         Lpr     = z2;
1474         A[r][j] = 4./(((PetscReal)p)*(((PetscReal)p)+1.)*Lpj*Lpr*(gllnodes[j]-gllnodes[r])*(gllnodes[j]-gllnodes[r]));
1475       }
1476     }
1477   }
1478   for (j=1; j<p+1; j++) {
1479     x  = gllnodes[j];
1480     z0 = 1.;
1481     z1 = x;
1482     for (nn=1; nn<p; nn++) {
1483       z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.));
1484       z0 = z1;
1485       z1 = z2;
1486     }
1487     Lpj     = z2;
1488     A[j][0] = 4.*PetscPowRealInt(-1.,p)/(((PetscReal)p)*(((PetscReal)p)+1.)*Lpj*(1.+gllnodes[j])*(1.+gllnodes[j]));
1489     A[0][j] = A[j][0];
1490   }
1491   for (j=0; j<p; j++) {
1492     x  = gllnodes[j];
1493     z0 = 1.;
1494     z1 = x;
1495     for (nn=1; nn<p; nn++) {
1496       z2 = x*z1*(2.*((PetscReal)nn)+1.)/(((PetscReal)nn)+1.)-z0*(((PetscReal)nn)/(((PetscReal)nn)+1.));
1497       z0 = z1;
1498       z1 = z2;
1499     }
1500     Lpj=z2;
1501 
1502     A[p][j] = 4./(((PetscReal)p)*(((PetscReal)p)+1.)*Lpj*(1.-gllnodes[j])*(1.-gllnodes[j]));
1503     A[j][p] = A[p][j];
1504   }
1505   A[0][0]=0.5+(((PetscReal)p)*(((PetscReal)p)+1.)-2.)/6.;
1506   A[p][p]=A[0][0];
1507   *AA = A;
1508   PetscFunctionReturn(0);
1509 }
1510 
1511 /*@C
1512    PetscGaussLobattoLegendreElementLaplacianDestroy - frees the Laplacian for a single 1d GLL element
1513 
1514    Not Collective
1515 
1516    Input Parameter:
1517 +  n - the number of GLL nodes
1518 .  nodes - the GLL nodes
1519 .  weights - the GLL weightss
1520 -  A - the stiffness element
1521 
1522    Level: beginner
1523 
1524 .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianCreate()
1525 
1526 @*/
1527 PetscErrorCode PetscGaussLobattoLegendreElementLaplacianDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA)
1528 {
1529   PetscErrorCode ierr;
1530 
1531   PetscFunctionBegin;
1532   ierr = PetscFree((*AA)[0]);CHKERRQ(ierr);
1533   ierr = PetscFree(*AA);CHKERRQ(ierr);
1534   *AA  = NULL;
1535   PetscFunctionReturn(0);
1536 }
1537 
1538 /*@C
1539    PetscGaussLobattoLegendreElementGradientCreate - computes the gradient for a single 1d GLL element
1540 
1541    Not Collective
1542 
1543    Input Parameter:
1544 +  n - the number of GLL nodes
1545 .  nodes - the GLL nodes
1546 .  weights - the GLL weights
1547 
1548    Output Parameter:
1549 .  AA - the stiffness element
1550 -  AAT - the transpose of AA (pass in NULL if you do not need this array)
1551 
1552    Level: beginner
1553 
1554    Notes:
1555     Destroy this with PetscGaussLobattoLegendreElementGradientDestroy()
1556 
1557    You can access entries in these arrays with AA[i][j] but in memory it is stored in contiguous memory, row oriented
1558 
1559 .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianDestroy()
1560 
1561 @*/
1562 PetscErrorCode PetscGaussLobattoLegendreElementGradientCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA,PetscReal ***AAT)
1563 {
1564   PetscReal        **A, **AT = NULL;
1565   PetscErrorCode  ierr;
1566   const PetscReal  *gllnodes = nodes;
1567   const PetscInt   p = n-1;
1568   PetscReal        q,qp,Li, Lj,d0;
1569   PetscInt         i,j;
1570 
1571   PetscFunctionBegin;
1572   ierr = PetscMalloc1(n,&A);CHKERRQ(ierr);
1573   ierr = PetscMalloc1(n*n,&A[0]);CHKERRQ(ierr);
1574   for (i=1; i<n; i++) A[i] = A[i-1]+n;
1575 
1576   if (AAT) {
1577     ierr = PetscMalloc1(n,&AT);CHKERRQ(ierr);
1578     ierr = PetscMalloc1(n*n,&AT[0]);CHKERRQ(ierr);
1579     for (i=1; i<n; i++) AT[i] = AT[i-1]+n;
1580   }
1581 
1582   if (n==1) {A[0][0] = 0.;}
1583   d0 = (PetscReal)p*((PetscReal)p+1.)/4.;
1584   for  (i=0; i<n; i++) {
1585     for  (j=0; j<n; j++) {
1586       A[i][j] = 0.;
1587       qAndLEvaluation(p,gllnodes[i],&q,&qp,&Li);
1588       qAndLEvaluation(p,gllnodes[j],&q,&qp,&Lj);
1589       if (i!=j)             A[i][j] = Li/(Lj*(gllnodes[i]-gllnodes[j]));
1590       if ((j==i) && (i==0)) A[i][j] = -d0;
1591       if (j==i && i==p)     A[i][j] = d0;
1592       if (AT) AT[j][i] = A[i][j];
1593     }
1594   }
1595   if (AAT) *AAT = AT;
1596   *AA  = A;
1597   PetscFunctionReturn(0);
1598 }
1599 
1600 /*@C
1601    PetscGaussLobattoLegendreElementGradientDestroy - frees the gradient for a single 1d GLL element obtained with PetscGaussLobattoLegendreElementGradientCreate()
1602 
1603    Not Collective
1604 
1605    Input Parameter:
1606 +  n - the number of GLL nodes
1607 .  nodes - the GLL nodes
1608 .  weights - the GLL weights
1609 .  AA - the stiffness element
1610 -  AAT - the transpose of the element
1611 
1612    Level: beginner
1613 
1614 .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianCreate(), PetscGaussLobattoLegendreElementAdvectionCreate()
1615 
1616 @*/
1617 PetscErrorCode PetscGaussLobattoLegendreElementGradientDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA,PetscReal ***AAT)
1618 {
1619   PetscErrorCode ierr;
1620 
1621   PetscFunctionBegin;
1622   ierr = PetscFree((*AA)[0]);CHKERRQ(ierr);
1623   ierr = PetscFree(*AA);CHKERRQ(ierr);
1624   *AA  = NULL;
1625   if (*AAT) {
1626     ierr = PetscFree((*AAT)[0]);CHKERRQ(ierr);
1627     ierr = PetscFree(*AAT);CHKERRQ(ierr);
1628     *AAT  = NULL;
1629   }
1630   PetscFunctionReturn(0);
1631 }
1632 
1633 /*@C
1634    PetscGaussLobattoLegendreElementAdvectionCreate - computes the advection operator for a single 1d GLL element
1635 
1636    Not Collective
1637 
1638    Input Parameter:
1639 +  n - the number of GLL nodes
1640 .  nodes - the GLL nodes
1641 .  weights - the GLL weightss
1642 
1643    Output Parameter:
1644 .  AA - the stiffness element
1645 
1646    Level: beginner
1647 
1648    Notes:
1649     Destroy this with PetscGaussLobattoLegendreElementAdvectionDestroy()
1650 
1651    This is the same as the Gradient operator multiplied by the diagonal mass matrix
1652 
1653    You can access entries in this array with AA[i][j] but in memory it is stored in contiguous memory, row oriented
1654 
1655 .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementLaplacianCreate(), PetscGaussLobattoLegendreElementAdvectionDestroy()
1656 
1657 @*/
1658 PetscErrorCode PetscGaussLobattoLegendreElementAdvectionCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA)
1659 {
1660   PetscReal       **D;
1661   PetscErrorCode  ierr;
1662   const PetscReal  *gllweights = weights;
1663   const PetscInt   glln = n;
1664   PetscInt         i,j;
1665 
1666   PetscFunctionBegin;
1667   ierr = PetscGaussLobattoLegendreElementGradientCreate(n,nodes,weights,&D,NULL);CHKERRQ(ierr);
1668   for (i=0; i<glln; i++){
1669     for (j=0; j<glln; j++) {
1670       D[i][j] = gllweights[i]*D[i][j];
1671     }
1672   }
1673   *AA = D;
1674   PetscFunctionReturn(0);
1675 }
1676 
1677 /*@C
1678    PetscGaussLobattoLegendreElementAdvectionDestroy - frees the advection stiffness for a single 1d GLL element
1679 
1680    Not Collective
1681 
1682    Input Parameter:
1683 +  n - the number of GLL nodes
1684 .  nodes - the GLL nodes
1685 .  weights - the GLL weights
1686 -  A - advection
1687 
1688    Level: beginner
1689 
1690 .seealso: PetscDTGaussLobattoLegendreQuadrature(), PetscGaussLobattoLegendreElementAdvectionCreate()
1691 
1692 @*/
1693 PetscErrorCode PetscGaussLobattoLegendreElementAdvectionDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA)
1694 {
1695   PetscErrorCode ierr;
1696 
1697   PetscFunctionBegin;
1698   ierr = PetscFree((*AA)[0]);CHKERRQ(ierr);
1699   ierr = PetscFree(*AA);CHKERRQ(ierr);
1700   *AA  = NULL;
1701   PetscFunctionReturn(0);
1702 }
1703 
1704 PetscErrorCode PetscGaussLobattoLegendreElementMassCreate(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA)
1705 {
1706   PetscReal        **A;
1707   PetscErrorCode  ierr;
1708   const PetscReal  *gllweights = weights;
1709   const PetscInt   glln = n;
1710   PetscInt         i,j;
1711 
1712   PetscFunctionBegin;
1713   ierr = PetscMalloc1(glln,&A);CHKERRQ(ierr);
1714   ierr = PetscMalloc1(glln*glln,&A[0]);CHKERRQ(ierr);
1715   for (i=1; i<glln; i++) A[i] = A[i-1]+glln;
1716   if (glln==1) {A[0][0] = 0.;}
1717   for  (i=0; i<glln; i++) {
1718     for  (j=0; j<glln; j++) {
1719       A[i][j] = 0.;
1720       if (j==i)     A[i][j] = gllweights[i];
1721     }
1722   }
1723   *AA  = A;
1724   PetscFunctionReturn(0);
1725 }
1726 
1727 PetscErrorCode PetscGaussLobattoLegendreElementMassDestroy(PetscInt n,PetscReal *nodes,PetscReal *weights,PetscReal ***AA)
1728 {
1729   PetscErrorCode ierr;
1730 
1731   PetscFunctionBegin;
1732   ierr = PetscFree((*AA)[0]);CHKERRQ(ierr);
1733   ierr = PetscFree(*AA);CHKERRQ(ierr);
1734   *AA  = NULL;
1735   PetscFunctionReturn(0);
1736 }
1737