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