xref: /petsc/src/dm/dt/interface/dt.c (revision 785e854f82a3c614b452fca2cf5ad4f2afe8bdde)
137045ce4SJed Brown /* Discretization tools */
237045ce4SJed Brown 
3a6fc04d9SSatish Balay #include <petscconf.h>
4a6fc04d9SSatish Balay #if defined(PETSC_HAVE_MATHIMF_H)
5a6fc04d9SSatish Balay #include <mathimf.h>           /* this needs to be included before math.h */
6a6fc04d9SSatish Balay #endif
7a6fc04d9SSatish Balay 
80c35b76eSJed Brown #include <petscdt.h>            /*I "petscdt.h" I*/
937045ce4SJed Brown #include <petscblaslapack.h>
10194825f6SJed Brown #include <petsc-private/petscimpl.h>
11665c2dedSJed Brown #include <petscviewer.h>
1259804f93SMatthew G. Knepley #include <petscdmplex.h>
1359804f93SMatthew G. Knepley #include <petscdmshell.h>
1437045ce4SJed Brown 
1537045ce4SJed Brown #undef __FUNCT__
16bfa639d9SMatthew G. Knepley #define __FUNCT__ "PetscQuadratureDestroy"
17bfa639d9SMatthew G. Knepley PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *q)
18bfa639d9SMatthew G. Knepley {
19bfa639d9SMatthew G. Knepley   PetscErrorCode ierr;
20bfa639d9SMatthew G. Knepley 
21bfa639d9SMatthew G. Knepley   PetscFunctionBegin;
22f9fd7fdbSMatthew G. Knepley   ierr = PetscFree(q->points);CHKERRQ(ierr);
23f9fd7fdbSMatthew G. Knepley   ierr = PetscFree(q->weights);CHKERRQ(ierr);
24f9fd7fdbSMatthew G. Knepley   PetscFunctionReturn(0);
25f9fd7fdbSMatthew G. Knepley }
26f9fd7fdbSMatthew G. Knepley 
27f9fd7fdbSMatthew G. Knepley #undef __FUNCT__
28f9fd7fdbSMatthew G. Knepley #define __FUNCT__ "PetscQuadratureView"
29f9fd7fdbSMatthew G. Knepley PetscErrorCode PetscQuadratureView(PetscQuadrature quad, PetscViewer viewer)
30f9fd7fdbSMatthew G. Knepley {
31f9fd7fdbSMatthew G. Knepley   PetscInt       q, d;
32f9fd7fdbSMatthew G. Knepley   PetscErrorCode ierr;
33f9fd7fdbSMatthew G. Knepley 
34f9fd7fdbSMatthew G. Knepley   PetscFunctionBegin;
35f9fd7fdbSMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "Quadrature on %d points\n  (", quad.numPoints);CHKERRQ(ierr);
36f9fd7fdbSMatthew G. Knepley   for (q = 0; q < quad.numPoints; ++q) {
37f9fd7fdbSMatthew G. Knepley     for (d = 0; d < quad.dim; ++d) {
38f9fd7fdbSMatthew G. Knepley       if (d) ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr);
39f9fd7fdbSMatthew G. Knepley       ierr = PetscViewerASCIIPrintf(viewer, "%g\n", quad.points[q*quad.dim+d]);CHKERRQ(ierr);
40f9fd7fdbSMatthew G. Knepley     }
41f9fd7fdbSMatthew G. Knepley     ierr = PetscViewerASCIIPrintf(viewer, ") %g\n", quad.weights[q]);CHKERRQ(ierr);
42f9fd7fdbSMatthew G. Knepley   }
43bfa639d9SMatthew G. Knepley   PetscFunctionReturn(0);
44bfa639d9SMatthew G. Knepley }
45bfa639d9SMatthew G. Knepley 
46bfa639d9SMatthew G. Knepley #undef __FUNCT__
4737045ce4SJed Brown #define __FUNCT__ "PetscDTLegendreEval"
4837045ce4SJed Brown /*@
4937045ce4SJed Brown    PetscDTLegendreEval - evaluate Legendre polynomial at points
5037045ce4SJed Brown 
5137045ce4SJed Brown    Not Collective
5237045ce4SJed Brown 
5337045ce4SJed Brown    Input Arguments:
5437045ce4SJed Brown +  npoints - number of spatial points to evaluate at
5537045ce4SJed Brown .  points - array of locations to evaluate at
5637045ce4SJed Brown .  ndegree - number of basis degrees to evaluate
5737045ce4SJed Brown -  degrees - sorted array of degrees to evaluate
5837045ce4SJed Brown 
5937045ce4SJed Brown    Output Arguments:
600298fd71SBarry Smith +  B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL)
610298fd71SBarry Smith .  D - row-oriented derivative evaluation matrix (or NULL)
620298fd71SBarry Smith -  D2 - row-oriented second derivative evaluation matrix (or NULL)
6337045ce4SJed Brown 
6437045ce4SJed Brown    Level: intermediate
6537045ce4SJed Brown 
6637045ce4SJed Brown .seealso: PetscDTGaussQuadrature()
6737045ce4SJed Brown @*/
6837045ce4SJed Brown PetscErrorCode PetscDTLegendreEval(PetscInt npoints,const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2)
6937045ce4SJed Brown {
7037045ce4SJed Brown   PetscInt i,maxdegree;
7137045ce4SJed Brown 
7237045ce4SJed Brown   PetscFunctionBegin;
7337045ce4SJed Brown   if (!npoints || !ndegree) PetscFunctionReturn(0);
7437045ce4SJed Brown   maxdegree = degrees[ndegree-1];
7537045ce4SJed Brown   for (i=0; i<npoints; i++) {
7637045ce4SJed Brown     PetscReal pm1,pm2,pd1,pd2,pdd1,pdd2,x;
7737045ce4SJed Brown     PetscInt  j,k;
7837045ce4SJed Brown     x    = points[i];
7937045ce4SJed Brown     pm2  = 0;
8037045ce4SJed Brown     pm1  = 1;
8137045ce4SJed Brown     pd2  = 0;
8237045ce4SJed Brown     pd1  = 0;
8337045ce4SJed Brown     pdd2 = 0;
8437045ce4SJed Brown     pdd1 = 0;
8537045ce4SJed Brown     k    = 0;
8637045ce4SJed Brown     if (degrees[k] == 0) {
8737045ce4SJed Brown       if (B) B[i*ndegree+k] = pm1;
8837045ce4SJed Brown       if (D) D[i*ndegree+k] = pd1;
8937045ce4SJed Brown       if (D2) D2[i*ndegree+k] = pdd1;
9037045ce4SJed Brown       k++;
9137045ce4SJed Brown     }
9237045ce4SJed Brown     for (j=1; j<=maxdegree; j++,k++) {
9337045ce4SJed Brown       PetscReal p,d,dd;
9437045ce4SJed Brown       p    = ((2*j-1)*x*pm1 - (j-1)*pm2)/j;
9537045ce4SJed Brown       d    = pd2 + (2*j-1)*pm1;
9637045ce4SJed Brown       dd   = pdd2 + (2*j-1)*pd1;
9737045ce4SJed Brown       pm2  = pm1;
9837045ce4SJed Brown       pm1  = p;
9937045ce4SJed Brown       pd2  = pd1;
10037045ce4SJed Brown       pd1  = d;
10137045ce4SJed Brown       pdd2 = pdd1;
10237045ce4SJed Brown       pdd1 = dd;
10337045ce4SJed Brown       if (degrees[k] == j) {
10437045ce4SJed Brown         if (B) B[i*ndegree+k] = p;
10537045ce4SJed Brown         if (D) D[i*ndegree+k] = d;
10637045ce4SJed Brown         if (D2) D2[i*ndegree+k] = dd;
10737045ce4SJed Brown       }
10837045ce4SJed Brown     }
10937045ce4SJed Brown   }
11037045ce4SJed Brown   PetscFunctionReturn(0);
11137045ce4SJed Brown }
11237045ce4SJed Brown 
11337045ce4SJed Brown #undef __FUNCT__
11437045ce4SJed Brown #define __FUNCT__ "PetscDTGaussQuadrature"
11537045ce4SJed Brown /*@
11637045ce4SJed Brown    PetscDTGaussQuadrature - create Gauss quadrature
11737045ce4SJed Brown 
11837045ce4SJed Brown    Not Collective
11937045ce4SJed Brown 
12037045ce4SJed Brown    Input Arguments:
12137045ce4SJed Brown +  npoints - number of points
12237045ce4SJed Brown .  a - left end of interval (often-1)
12337045ce4SJed Brown -  b - right end of interval (often +1)
12437045ce4SJed Brown 
12537045ce4SJed Brown    Output Arguments:
12637045ce4SJed Brown +  x - quadrature points
12737045ce4SJed Brown -  w - quadrature weights
12837045ce4SJed Brown 
12937045ce4SJed Brown    Level: intermediate
13037045ce4SJed Brown 
13137045ce4SJed Brown    References:
13237045ce4SJed Brown    Golub and Welsch, Calculation of Quadrature Rules, Math. Comp. 23(106), 221--230, 1969.
13337045ce4SJed Brown 
13437045ce4SJed Brown .seealso: PetscDTLegendreEval()
13537045ce4SJed Brown @*/
13637045ce4SJed Brown PetscErrorCode PetscDTGaussQuadrature(PetscInt npoints,PetscReal a,PetscReal b,PetscReal *x,PetscReal *w)
13737045ce4SJed Brown {
13837045ce4SJed Brown   PetscErrorCode ierr;
13937045ce4SJed Brown   PetscInt       i;
14037045ce4SJed Brown   PetscReal      *work;
14137045ce4SJed Brown   PetscScalar    *Z;
14237045ce4SJed Brown   PetscBLASInt   N,LDZ,info;
14337045ce4SJed Brown 
14437045ce4SJed Brown   PetscFunctionBegin;
14537045ce4SJed Brown   /* Set up the Golub-Welsch system */
14637045ce4SJed Brown   for (i=0; i<npoints; i++) {
14737045ce4SJed Brown     x[i] = 0;                   /* diagonal is 0 */
14837045ce4SJed Brown     if (i) w[i-1] = 0.5 / PetscSqrtReal(1 - 1./PetscSqr(2*i));
14937045ce4SJed Brown   }
15037045ce4SJed Brown   ierr = PetscRealView(npoints-1,w,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
151dcca6d9dSJed Brown   ierr = PetscMalloc2(npoints*npoints,&Z,PetscMax(1,2*npoints-2),&work);CHKERRQ(ierr);
152c5df96a5SBarry Smith   ierr = PetscBLASIntCast(npoints,&N);CHKERRQ(ierr);
15337045ce4SJed Brown   LDZ  = N;
15437045ce4SJed Brown   ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1558b83055fSJed Brown   PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("I",&N,x,w,Z,&LDZ,work,&info));
15637045ce4SJed Brown   ierr = PetscFPTrapPop();CHKERRQ(ierr);
1571c3d6f74SJed Brown   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEQR error");
15837045ce4SJed Brown 
15937045ce4SJed Brown   for (i=0; i<(npoints+1)/2; i++) {
16037045ce4SJed Brown     PetscReal y = 0.5 * (-x[i] + x[npoints-i-1]); /* enforces symmetry */
16137045ce4SJed Brown     x[i]           = (a+b)/2 - y*(b-a)/2;
16237045ce4SJed Brown     x[npoints-i-1] = (a+b)/2 + y*(b-a)/2;
1630d644c17SKarl Rupp 
16437045ce4SJed Brown     w[i] = w[npoints-1-i] = (b-a)*PetscSqr(0.5*PetscAbsScalar(Z[i*npoints] + Z[(npoints-i-1)*npoints]));
16537045ce4SJed Brown   }
16637045ce4SJed Brown   ierr = PetscFree2(Z,work);CHKERRQ(ierr);
16737045ce4SJed Brown   PetscFunctionReturn(0);
16837045ce4SJed Brown }
169194825f6SJed Brown 
170194825f6SJed Brown #undef __FUNCT__
171494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTFactorial_Internal"
172494e7359SMatthew G. Knepley /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x.
173494e7359SMatthew G. Knepley    Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */
174494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorial_Internal(PetscInt n, PetscReal *factorial)
175494e7359SMatthew G. Knepley {
176494e7359SMatthew G. Knepley   PetscReal f = 1.0;
177494e7359SMatthew G. Knepley   PetscInt  i;
178494e7359SMatthew G. Knepley 
179494e7359SMatthew G. Knepley   PetscFunctionBegin;
180494e7359SMatthew G. Knepley   for (i = 1; i < n+1; ++i) f *= i;
181494e7359SMatthew G. Knepley   *factorial = f;
182494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
183494e7359SMatthew G. Knepley }
184494e7359SMatthew G. Knepley 
185494e7359SMatthew G. Knepley #undef __FUNCT__
186494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTComputeJacobi"
187494e7359SMatthew G. Knepley /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x.
188494e7359SMatthew G. Knepley    Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */
189494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobi(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P)
190494e7359SMatthew G. Knepley {
191494e7359SMatthew G. Knepley   PetscReal apb, pn1, pn2;
192494e7359SMatthew G. Knepley   PetscInt  k;
193494e7359SMatthew G. Knepley 
194494e7359SMatthew G. Knepley   PetscFunctionBegin;
195494e7359SMatthew G. Knepley   if (!n) {*P = 1.0; PetscFunctionReturn(0);}
196494e7359SMatthew G. Knepley   if (n == 1) {*P = 0.5 * (a - b + (a + b + 2.0) * x); PetscFunctionReturn(0);}
197494e7359SMatthew G. Knepley   apb = a + b;
198494e7359SMatthew G. Knepley   pn2 = 1.0;
199494e7359SMatthew G. Knepley   pn1 = 0.5 * (a - b + (apb + 2.0) * x);
200494e7359SMatthew G. Knepley   *P  = 0.0;
201494e7359SMatthew G. Knepley   for (k = 2; k < n+1; ++k) {
202494e7359SMatthew G. Knepley     PetscReal a1 = 2.0 * k * (k + apb) * (2.0*k + apb - 2.0);
203494e7359SMatthew G. Knepley     PetscReal a2 = (2.0 * k + apb - 1.0) * (a*a - b*b);
204494e7359SMatthew G. Knepley     PetscReal a3 = (2.0 * k + apb - 2.0) * (2.0 * k + apb - 1.0) * (2.0 * k + apb);
205494e7359SMatthew G. Knepley     PetscReal a4 = 2.0 * (k + a - 1.0) * (k + b - 1.0) * (2.0 * k + apb);
206494e7359SMatthew G. Knepley 
207494e7359SMatthew G. Knepley     a2  = a2 / a1;
208494e7359SMatthew G. Knepley     a3  = a3 / a1;
209494e7359SMatthew G. Knepley     a4  = a4 / a1;
210494e7359SMatthew G. Knepley     *P  = (a2 + a3 * x) * pn1 - a4 * pn2;
211494e7359SMatthew G. Knepley     pn2 = pn1;
212494e7359SMatthew G. Knepley     pn1 = *P;
213494e7359SMatthew G. Knepley   }
214494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
215494e7359SMatthew G. Knepley }
216494e7359SMatthew G. Knepley 
217494e7359SMatthew G. Knepley #undef __FUNCT__
218494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTComputeJacobiDerivative"
219494e7359SMatthew G. Knepley /* Evaluates the first derivative of P_{n}^{a,b} at a point x. */
220494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobiDerivative(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P)
221494e7359SMatthew G. Knepley {
222494e7359SMatthew G. Knepley   PetscReal      nP;
223494e7359SMatthew G. Knepley   PetscErrorCode ierr;
224494e7359SMatthew G. Knepley 
225494e7359SMatthew G. Knepley   PetscFunctionBegin;
226494e7359SMatthew G. Knepley   if (!n) {*P = 0.0; PetscFunctionReturn(0);}
227494e7359SMatthew G. Knepley   ierr = PetscDTComputeJacobi(a+1, b+1, n-1, x, &nP);CHKERRQ(ierr);
228494e7359SMatthew G. Knepley   *P   = 0.5 * (a + b + n + 1) * nP;
229494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
230494e7359SMatthew G. Knepley }
231494e7359SMatthew G. Knepley 
232494e7359SMatthew G. Knepley #undef __FUNCT__
233494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTMapSquareToTriangle_Internal"
234494e7359SMatthew G. Knepley /* Maps from [-1,1]^2 to the (-1,1) reference triangle */
235494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTMapSquareToTriangle_Internal(PetscReal x, PetscReal y, PetscReal *xi, PetscReal *eta)
236494e7359SMatthew G. Knepley {
237494e7359SMatthew G. Knepley   PetscFunctionBegin;
238494e7359SMatthew G. Knepley   *xi  = 0.5 * (1.0 + x) * (1.0 - y) - 1.0;
239494e7359SMatthew G. Knepley   *eta = y;
240494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
241494e7359SMatthew G. Knepley }
242494e7359SMatthew G. Knepley 
243494e7359SMatthew G. Knepley #undef __FUNCT__
244494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTMapCubeToTetrahedron_Internal"
245494e7359SMatthew G. Knepley /* Maps from [-1,1]^2 to the (-1,1) reference triangle */
246494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTMapCubeToTetrahedron_Internal(PetscReal x, PetscReal y, PetscReal z, PetscReal *xi, PetscReal *eta, PetscReal *zeta)
247494e7359SMatthew G. Knepley {
248494e7359SMatthew G. Knepley   PetscFunctionBegin;
249494e7359SMatthew G. Knepley   *xi   = 0.25 * (1.0 + x) * (1.0 - y) * (1.0 - z) - 1.0;
250494e7359SMatthew G. Knepley   *eta  = 0.5  * (1.0 + y) * (1.0 - z) - 1.0;
251494e7359SMatthew G. Knepley   *zeta = z;
252494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
253494e7359SMatthew G. Knepley }
254494e7359SMatthew G. Knepley 
255494e7359SMatthew G. Knepley #undef __FUNCT__
256494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTGaussJacobiQuadrature1D_Internal"
257494e7359SMatthew G. Knepley static PetscErrorCode PetscDTGaussJacobiQuadrature1D_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal *x, PetscReal *w)
258494e7359SMatthew G. Knepley {
259494e7359SMatthew G. Knepley   PetscInt       maxIter = 100;
260494e7359SMatthew G. Knepley   PetscReal      eps     = 1.0e-8;
261a8291ba1SSatish Balay   PetscReal      a1, a2, a3, a4, a5, a6;
262494e7359SMatthew G. Knepley   PetscInt       k;
263494e7359SMatthew G. Knepley   PetscErrorCode ierr;
264494e7359SMatthew G. Knepley 
265494e7359SMatthew G. Knepley   PetscFunctionBegin;
266a8291ba1SSatish Balay 
267a8291ba1SSatish Balay   a1      = pow(2, a+b+1);
268a8291ba1SSatish Balay #if defined(PETSC_HAVE_TGAMMA)
269a8291ba1SSatish Balay   a2      = tgamma(a + npoints + 1);
270a8291ba1SSatish Balay   a3      = tgamma(b + npoints + 1);
271a8291ba1SSatish Balay   a4      = tgamma(a + b + npoints + 1);
272a8291ba1SSatish Balay #else
273a8291ba1SSatish Balay   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"tgamma() - math routine is unavailable.");
274a8291ba1SSatish Balay #endif
275a8291ba1SSatish Balay 
276494e7359SMatthew G. Knepley   ierr = PetscDTFactorial_Internal(npoints, &a5);CHKERRQ(ierr);
277494e7359SMatthew G. Knepley   a6   = a1 * a2 * a3 / a4 / a5;
278494e7359SMatthew G. Knepley   /* Computes the m roots of P_{m}^{a,b} on [-1,1] by Newton's method with Chebyshev points as initial guesses.
279494e7359SMatthew G. Knepley    Algorithm implemented from the pseudocode given by Karniadakis and Sherwin and Python in FIAT */
280494e7359SMatthew G. Knepley   for (k = 0; k < npoints; ++k) {
2817f1c68b3SMatthew G. Knepley     PetscReal r = -cos((2.0*k + 1.0) * PETSC_PI / (2.0 * npoints)), dP;
282494e7359SMatthew G. Knepley     PetscInt  j;
283494e7359SMatthew G. Knepley 
284494e7359SMatthew G. Knepley     if (k > 0) r = 0.5 * (r + x[k-1]);
285494e7359SMatthew G. Knepley     for (j = 0; j < maxIter; ++j) {
286494e7359SMatthew G. Knepley       PetscReal s = 0.0, delta, f, fp;
287494e7359SMatthew G. Knepley       PetscInt  i;
288494e7359SMatthew G. Knepley 
289494e7359SMatthew G. Knepley       for (i = 0; i < k; ++i) s = s + 1.0 / (r - x[i]);
290494e7359SMatthew G. Knepley       ierr = PetscDTComputeJacobi(a, b, npoints, r, &f);CHKERRQ(ierr);
291494e7359SMatthew G. Knepley       ierr = PetscDTComputeJacobiDerivative(a, b, npoints, r, &fp);CHKERRQ(ierr);
292494e7359SMatthew G. Knepley       delta = f / (fp - f * s);
293494e7359SMatthew G. Knepley       r     = r - delta;
294494e7359SMatthew G. Knepley       if (fabs(delta) < eps) break;
295494e7359SMatthew G. Knepley     }
296494e7359SMatthew G. Knepley     x[k] = r;
297494e7359SMatthew G. Knepley     ierr = PetscDTComputeJacobiDerivative(a, b, npoints, x[k], &dP);CHKERRQ(ierr);
298494e7359SMatthew G. Knepley     w[k] = a6 / (1.0 - PetscSqr(x[k])) / PetscSqr(dP);
299494e7359SMatthew G. Knepley   }
300494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
301494e7359SMatthew G. Knepley }
302494e7359SMatthew G. Knepley 
303494e7359SMatthew G. Knepley #undef __FUNCT__
304494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTGaussJacobiQuadrature"
305fd9d31fbSMatthew G. Knepley /*@C
306494e7359SMatthew G. Knepley   PetscDTGaussJacobiQuadrature - create Gauss-Jacobi quadrature for a simplex
307494e7359SMatthew G. Knepley 
308494e7359SMatthew G. Knepley   Not Collective
309494e7359SMatthew G. Knepley 
310494e7359SMatthew G. Knepley   Input Arguments:
311494e7359SMatthew G. Knepley + dim - The simplex dimension
312552aa4f7SMatthew G. Knepley . order - The quadrature order
313494e7359SMatthew G. Knepley . a - left end of interval (often-1)
314494e7359SMatthew G. Knepley - b - right end of interval (often +1)
315494e7359SMatthew G. Knepley 
316494e7359SMatthew G. Knepley   Output Arguments:
317552aa4f7SMatthew G. Knepley . q - A PetscQuadrature object
318494e7359SMatthew G. Knepley 
319494e7359SMatthew G. Knepley   Level: intermediate
320494e7359SMatthew G. Knepley 
321494e7359SMatthew G. Knepley   References:
322494e7359SMatthew G. Knepley   Karniadakis and Sherwin.
323494e7359SMatthew G. Knepley   FIAT
324494e7359SMatthew G. Knepley 
325494e7359SMatthew G. Knepley .seealso: PetscDTGaussQuadrature()
326494e7359SMatthew G. Knepley @*/
327552aa4f7SMatthew G. Knepley PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt dim, PetscInt order, PetscReal a, PetscReal b, PetscQuadrature *q)
328494e7359SMatthew G. Knepley {
329552aa4f7SMatthew G. Knepley   PetscInt       npoints = dim > 1 ? dim > 2 ? order*PetscSqr(order) : PetscSqr(order) : order;
330494e7359SMatthew G. Knepley   PetscReal     *px, *wx, *py, *wy, *pz, *wz, *x, *w;
331494e7359SMatthew G. Knepley   PetscInt       i, j, k;
332494e7359SMatthew G. Knepley   PetscErrorCode ierr;
333494e7359SMatthew G. Knepley 
334494e7359SMatthew G. Knepley   PetscFunctionBegin;
335494e7359SMatthew G. Knepley   if ((a != -1.0) || (b != 1.0)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must use default internal right now");
336*785e854fSJed Brown   ierr = PetscMalloc1(npoints*dim, &x);CHKERRQ(ierr);
337*785e854fSJed Brown   ierr = PetscMalloc1(npoints, &w);CHKERRQ(ierr);
338494e7359SMatthew G. Knepley   switch (dim) {
339707aa5c5SMatthew G. Knepley   case 0:
340707aa5c5SMatthew G. Knepley     ierr = PetscFree(x);CHKERRQ(ierr);
341707aa5c5SMatthew G. Knepley     ierr = PetscFree(w);CHKERRQ(ierr);
342*785e854fSJed Brown     ierr = PetscMalloc1(1, &x);CHKERRQ(ierr);
343*785e854fSJed Brown     ierr = PetscMalloc1(1, &w);CHKERRQ(ierr);
344707aa5c5SMatthew G. Knepley     x[0] = 0.0;
345707aa5c5SMatthew G. Knepley     w[0] = 1.0;
346707aa5c5SMatthew G. Knepley     break;
347494e7359SMatthew G. Knepley   case 1:
348552aa4f7SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, x, w);CHKERRQ(ierr);
349494e7359SMatthew G. Knepley     break;
350494e7359SMatthew G. Knepley   case 2:
351dcca6d9dSJed Brown     ierr = PetscMalloc4(order,&px,order,&wx,order,&py,order,&wy);CHKERRQ(ierr);
352552aa4f7SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, px, wx);CHKERRQ(ierr);
353552aa4f7SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 1.0, 0.0, py, wy);CHKERRQ(ierr);
354552aa4f7SMatthew G. Knepley     for (i = 0; i < order; ++i) {
355552aa4f7SMatthew G. Knepley       for (j = 0; j < order; ++j) {
356552aa4f7SMatthew G. Knepley         ierr = PetscDTMapSquareToTriangle_Internal(px[i], py[j], &x[(i*order+j)*2+0], &x[(i*order+j)*2+1]);CHKERRQ(ierr);
357552aa4f7SMatthew G. Knepley         w[i*order+j] = 0.5 * wx[i] * wy[j];
358494e7359SMatthew G. Knepley       }
359494e7359SMatthew G. Knepley     }
360494e7359SMatthew G. Knepley     ierr = PetscFree4(px,wx,py,wy);CHKERRQ(ierr);
361494e7359SMatthew G. Knepley     break;
362494e7359SMatthew G. Knepley   case 3:
363dcca6d9dSJed Brown     ierr = PetscMalloc6(order,&px,order,&wx,order,&py,order,&wy,order,&pz,order,&wz);CHKERRQ(ierr);
364552aa4f7SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, px, wx);CHKERRQ(ierr);
365552aa4f7SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 1.0, 0.0, py, wy);CHKERRQ(ierr);
366552aa4f7SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 2.0, 0.0, pz, wz);CHKERRQ(ierr);
367552aa4f7SMatthew G. Knepley     for (i = 0; i < order; ++i) {
368552aa4f7SMatthew G. Knepley       for (j = 0; j < order; ++j) {
369552aa4f7SMatthew G. Knepley         for (k = 0; k < order; ++k) {
370552aa4f7SMatthew G. Knepley           ierr = PetscDTMapCubeToTetrahedron_Internal(px[i], py[j], pz[k], &x[((i*order+j)*order+k)*3+0], &x[((i*order+j)*order+k)*3+1], &x[((i*order+j)*order+k)*3+2]);CHKERRQ(ierr);
371552aa4f7SMatthew G. Knepley           w[(i*order+j)*order+k] = 0.125 * wx[i] * wy[j] * wz[k];
372494e7359SMatthew G. Knepley         }
373494e7359SMatthew G. Knepley       }
374494e7359SMatthew G. Knepley     }
375494e7359SMatthew G. Knepley     ierr = PetscFree6(px,wx,py,wy,pz,wz);CHKERRQ(ierr);
376494e7359SMatthew G. Knepley     break;
377494e7359SMatthew G. Knepley   default:
378494e7359SMatthew G. Knepley     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim);
379494e7359SMatthew G. Knepley   }
380f9fd7fdbSMatthew G. Knepley   q->dim       = dim;
381f9fd7fdbSMatthew G. Knepley   q->numPoints = npoints;
382f9fd7fdbSMatthew G. Knepley   q->points    = x;
383f9fd7fdbSMatthew G. Knepley   q->weights   = w;
384494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
385494e7359SMatthew G. Knepley }
386494e7359SMatthew G. Knepley 
387494e7359SMatthew G. Knepley #undef __FUNCT__
388194825f6SJed Brown #define __FUNCT__ "PetscDTPseudoInverseQR"
389194825f6SJed Brown /* Overwrites A. Can only handle full-rank problems with m>=n
390194825f6SJed Brown  * A in column-major format
391194825f6SJed Brown  * Ainv in row-major format
392194825f6SJed Brown  * tau has length m
393194825f6SJed Brown  * worksize must be >= max(1,n)
394194825f6SJed Brown  */
395194825f6SJed Brown static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m,PetscInt mstride,PetscInt n,PetscReal *A_in,PetscReal *Ainv_out,PetscScalar *tau,PetscInt worksize,PetscScalar *work)
396194825f6SJed Brown {
397194825f6SJed Brown   PetscErrorCode ierr;
398194825f6SJed Brown   PetscBLASInt M,N,K,lda,ldb,ldwork,info;
399194825f6SJed Brown   PetscScalar *A,*Ainv,*R,*Q,Alpha;
400194825f6SJed Brown 
401194825f6SJed Brown   PetscFunctionBegin;
402194825f6SJed Brown #if defined(PETSC_USE_COMPLEX)
403194825f6SJed Brown   {
404194825f6SJed Brown     PetscInt i,j;
405dcca6d9dSJed Brown     ierr = PetscMalloc2(m*n,&A,m*n,&Ainv);CHKERRQ(ierr);
406194825f6SJed Brown     for (j=0; j<n; j++) {
407194825f6SJed Brown       for (i=0; i<m; i++) A[i+m*j] = A_in[i+mstride*j];
408194825f6SJed Brown     }
409194825f6SJed Brown     mstride = m;
410194825f6SJed Brown   }
411194825f6SJed Brown #else
412194825f6SJed Brown   A = A_in;
413194825f6SJed Brown   Ainv = Ainv_out;
414194825f6SJed Brown #endif
415194825f6SJed Brown 
416194825f6SJed Brown   ierr = PetscBLASIntCast(m,&M);CHKERRQ(ierr);
417194825f6SJed Brown   ierr = PetscBLASIntCast(n,&N);CHKERRQ(ierr);
418194825f6SJed Brown   ierr = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr);
419194825f6SJed Brown   ierr = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr);
420194825f6SJed Brown   ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
421194825f6SJed Brown   LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info);
422194825f6SJed Brown   ierr = PetscFPTrapPop();CHKERRQ(ierr);
423194825f6SJed Brown   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error");
424194825f6SJed Brown   R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */
425194825f6SJed Brown 
426194825f6SJed Brown   /* Extract an explicit representation of Q */
427194825f6SJed Brown   Q = Ainv;
428194825f6SJed Brown   ierr = PetscMemcpy(Q,A,mstride*n*sizeof(PetscScalar));CHKERRQ(ierr);
429194825f6SJed Brown   K = N;                        /* full rank */
430194825f6SJed Brown   LAPACKungqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info);
431194825f6SJed Brown   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error");
432194825f6SJed Brown 
433194825f6SJed Brown   /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */
434194825f6SJed Brown   Alpha = 1.0;
435194825f6SJed Brown   ldb = lda;
436194825f6SJed Brown   BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb);
437194825f6SJed Brown   /* Ainv is Q, overwritten with inverse */
438194825f6SJed Brown 
439194825f6SJed Brown #if defined(PETSC_USE_COMPLEX)
440194825f6SJed Brown   {
441194825f6SJed Brown     PetscInt i;
442194825f6SJed Brown     for (i=0; i<m*n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]);
443194825f6SJed Brown     ierr = PetscFree2(A,Ainv);CHKERRQ(ierr);
444194825f6SJed Brown   }
445194825f6SJed Brown #endif
446194825f6SJed Brown   PetscFunctionReturn(0);
447194825f6SJed Brown }
448194825f6SJed Brown 
449194825f6SJed Brown #undef __FUNCT__
450194825f6SJed Brown #define __FUNCT__ "PetscDTLegendreIntegrate"
451194825f6SJed Brown /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */
452194825f6SJed Brown static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval,const PetscReal *x,PetscInt ndegree,const PetscInt *degrees,PetscBool Transpose,PetscReal *B)
453194825f6SJed Brown {
454194825f6SJed Brown   PetscErrorCode ierr;
455194825f6SJed Brown   PetscReal *Bv;
456194825f6SJed Brown   PetscInt i,j;
457194825f6SJed Brown 
458194825f6SJed Brown   PetscFunctionBegin;
459*785e854fSJed Brown   ierr = PetscMalloc1((ninterval+1)*ndegree,&Bv);CHKERRQ(ierr);
460194825f6SJed Brown   /* Point evaluation of L_p on all the source vertices */
461194825f6SJed Brown   ierr = PetscDTLegendreEval(ninterval+1,x,ndegree,degrees,Bv,NULL,NULL);CHKERRQ(ierr);
462194825f6SJed Brown   /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */
463194825f6SJed Brown   for (i=0; i<ninterval; i++) {
464194825f6SJed Brown     for (j=0; j<ndegree; j++) {
465194825f6SJed Brown       if (Transpose) B[i+ninterval*j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j];
466194825f6SJed Brown       else           B[i*ndegree+j]   = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j];
467194825f6SJed Brown     }
468194825f6SJed Brown   }
469194825f6SJed Brown   ierr = PetscFree(Bv);CHKERRQ(ierr);
470194825f6SJed Brown   PetscFunctionReturn(0);
471194825f6SJed Brown }
472194825f6SJed Brown 
473194825f6SJed Brown #undef __FUNCT__
474194825f6SJed Brown #define __FUNCT__ "PetscDTReconstructPoly"
475194825f6SJed Brown /*@
476194825f6SJed Brown    PetscDTReconstructPoly - create matrix representing polynomial reconstruction using cell intervals and evaluation at target intervals
477194825f6SJed Brown 
478194825f6SJed Brown    Not Collective
479194825f6SJed Brown 
480194825f6SJed Brown    Input Arguments:
481194825f6SJed Brown +  degree - degree of reconstruction polynomial
482194825f6SJed Brown .  nsource - number of source intervals
483194825f6SJed Brown .  sourcex - sorted coordinates of source cell boundaries (length nsource+1)
484194825f6SJed Brown .  ntarget - number of target intervals
485194825f6SJed Brown -  targetx - sorted coordinates of target cell boundaries (length ntarget+1)
486194825f6SJed Brown 
487194825f6SJed Brown    Output Arguments:
488194825f6SJed Brown .  R - reconstruction matrix, utarget = sum_s R[t*nsource+s] * usource[s]
489194825f6SJed Brown 
490194825f6SJed Brown    Level: advanced
491194825f6SJed Brown 
492194825f6SJed Brown .seealso: PetscDTLegendreEval()
493194825f6SJed Brown @*/
494194825f6SJed Brown PetscErrorCode PetscDTReconstructPoly(PetscInt degree,PetscInt nsource,const PetscReal *sourcex,PetscInt ntarget,const PetscReal *targetx,PetscReal *R)
495194825f6SJed Brown {
496194825f6SJed Brown   PetscErrorCode ierr;
497194825f6SJed Brown   PetscInt i,j,k,*bdegrees,worksize;
498194825f6SJed Brown   PetscReal xmin,xmax,center,hscale,*sourcey,*targety,*Bsource,*Bsinv,*Btarget;
499194825f6SJed Brown   PetscScalar *tau,*work;
500194825f6SJed Brown 
501194825f6SJed Brown   PetscFunctionBegin;
502194825f6SJed Brown   PetscValidRealPointer(sourcex,3);
503194825f6SJed Brown   PetscValidRealPointer(targetx,5);
504194825f6SJed Brown   PetscValidRealPointer(R,6);
505194825f6SJed Brown   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);
506194825f6SJed Brown #if defined(PETSC_USE_DEBUG)
507194825f6SJed Brown   for (i=0; i<nsource; i++) {
508194825f6SJed Brown     if (sourcex[i] >= sourcex[i+1]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Source interval %D has negative orientation (%G,%G)",i,sourcex[i],sourcex[i+1]);
509194825f6SJed Brown   }
510194825f6SJed Brown   for (i=0; i<ntarget; i++) {
511194825f6SJed Brown     if (targetx[i] >= targetx[i+1]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"Target interval %D has negative orientation (%G,%G)",i,targetx[i],targetx[i+1]);
512194825f6SJed Brown   }
513194825f6SJed Brown #endif
514194825f6SJed Brown   xmin = PetscMin(sourcex[0],targetx[0]);
515194825f6SJed Brown   xmax = PetscMax(sourcex[nsource],targetx[ntarget]);
516194825f6SJed Brown   center = (xmin + xmax)/2;
517194825f6SJed Brown   hscale = (xmax - xmin)/2;
518194825f6SJed Brown   worksize = nsource;
519dcca6d9dSJed Brown   ierr = PetscMalloc4(degree+1,&bdegrees,nsource+1,&sourcey,nsource*(degree+1),&Bsource,worksize,&work);CHKERRQ(ierr);
520dcca6d9dSJed Brown   ierr = PetscMalloc4(nsource,&tau,nsource*(degree+1),&Bsinv,ntarget+1,&targety,ntarget*(degree+1),&Btarget);CHKERRQ(ierr);
521194825f6SJed Brown   for (i=0; i<=nsource; i++) sourcey[i] = (sourcex[i]-center)/hscale;
522194825f6SJed Brown   for (i=0; i<=degree; i++) bdegrees[i] = i+1;
523194825f6SJed Brown   ierr = PetscDTLegendreIntegrate(nsource,sourcey,degree+1,bdegrees,PETSC_TRUE,Bsource);CHKERRQ(ierr);
524194825f6SJed Brown   ierr = PetscDTPseudoInverseQR(nsource,nsource,degree+1,Bsource,Bsinv,tau,nsource,work);CHKERRQ(ierr);
525194825f6SJed Brown   for (i=0; i<=ntarget; i++) targety[i] = (targetx[i]-center)/hscale;
526194825f6SJed Brown   ierr = PetscDTLegendreIntegrate(ntarget,targety,degree+1,bdegrees,PETSC_FALSE,Btarget);CHKERRQ(ierr);
527194825f6SJed Brown   for (i=0; i<ntarget; i++) {
528194825f6SJed Brown     PetscReal rowsum = 0;
529194825f6SJed Brown     for (j=0; j<nsource; j++) {
530194825f6SJed Brown       PetscReal sum = 0;
531194825f6SJed Brown       for (k=0; k<degree+1; k++) {
532194825f6SJed Brown         sum += Btarget[i*(degree+1)+k] * Bsinv[k*nsource+j];
533194825f6SJed Brown       }
534194825f6SJed Brown       R[i*nsource+j] = sum;
535194825f6SJed Brown       rowsum += sum;
536194825f6SJed Brown     }
537194825f6SJed Brown     for (j=0; j<nsource; j++) R[i*nsource+j] /= rowsum; /* normalize each row */
538194825f6SJed Brown   }
539194825f6SJed Brown   ierr = PetscFree4(bdegrees,sourcey,Bsource,work);CHKERRQ(ierr);
540194825f6SJed Brown   ierr = PetscFree4(tau,Bsinv,targety,Btarget);CHKERRQ(ierr);
541194825f6SJed Brown   PetscFunctionReturn(0);
542194825f6SJed Brown }
543