xref: /petsc/src/dm/dt/interface/dt.c (revision 707aa5c505b0ca8de3125920f23b39da005b42e7)
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;
22bfa639d9SMatthew G. Knepley   ierr = PetscFree(q->quadPoints);CHKERRQ(ierr);
23bfa639d9SMatthew G. Knepley   ierr = PetscFree(q->quadWeights);CHKERRQ(ierr);
24bfa639d9SMatthew G. Knepley   PetscFunctionReturn(0);
25bfa639d9SMatthew G. Knepley }
26bfa639d9SMatthew G. Knepley 
27bfa639d9SMatthew G. Knepley #undef __FUNCT__
2837045ce4SJed Brown #define __FUNCT__ "PetscDTLegendreEval"
2937045ce4SJed Brown /*@
3037045ce4SJed Brown    PetscDTLegendreEval - evaluate Legendre polynomial at points
3137045ce4SJed Brown 
3237045ce4SJed Brown    Not Collective
3337045ce4SJed Brown 
3437045ce4SJed Brown    Input Arguments:
3537045ce4SJed Brown +  npoints - number of spatial points to evaluate at
3637045ce4SJed Brown .  points - array of locations to evaluate at
3737045ce4SJed Brown .  ndegree - number of basis degrees to evaluate
3837045ce4SJed Brown -  degrees - sorted array of degrees to evaluate
3937045ce4SJed Brown 
4037045ce4SJed Brown    Output Arguments:
410298fd71SBarry Smith +  B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL)
420298fd71SBarry Smith .  D - row-oriented derivative evaluation matrix (or NULL)
430298fd71SBarry Smith -  D2 - row-oriented second derivative evaluation matrix (or NULL)
4437045ce4SJed Brown 
4537045ce4SJed Brown    Level: intermediate
4637045ce4SJed Brown 
4737045ce4SJed Brown .seealso: PetscDTGaussQuadrature()
4837045ce4SJed Brown @*/
4937045ce4SJed Brown PetscErrorCode PetscDTLegendreEval(PetscInt npoints,const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2)
5037045ce4SJed Brown {
5137045ce4SJed Brown   PetscInt i,maxdegree;
5237045ce4SJed Brown 
5337045ce4SJed Brown   PetscFunctionBegin;
5437045ce4SJed Brown   if (!npoints || !ndegree) PetscFunctionReturn(0);
5537045ce4SJed Brown   maxdegree = degrees[ndegree-1];
5637045ce4SJed Brown   for (i=0; i<npoints; i++) {
5737045ce4SJed Brown     PetscReal pm1,pm2,pd1,pd2,pdd1,pdd2,x;
5837045ce4SJed Brown     PetscInt  j,k;
5937045ce4SJed Brown     x    = points[i];
6037045ce4SJed Brown     pm2  = 0;
6137045ce4SJed Brown     pm1  = 1;
6237045ce4SJed Brown     pd2  = 0;
6337045ce4SJed Brown     pd1  = 0;
6437045ce4SJed Brown     pdd2 = 0;
6537045ce4SJed Brown     pdd1 = 0;
6637045ce4SJed Brown     k    = 0;
6737045ce4SJed Brown     if (degrees[k] == 0) {
6837045ce4SJed Brown       if (B) B[i*ndegree+k] = pm1;
6937045ce4SJed Brown       if (D) D[i*ndegree+k] = pd1;
7037045ce4SJed Brown       if (D2) D2[i*ndegree+k] = pdd1;
7137045ce4SJed Brown       k++;
7237045ce4SJed Brown     }
7337045ce4SJed Brown     for (j=1; j<=maxdegree; j++,k++) {
7437045ce4SJed Brown       PetscReal p,d,dd;
7537045ce4SJed Brown       p    = ((2*j-1)*x*pm1 - (j-1)*pm2)/j;
7637045ce4SJed Brown       d    = pd2 + (2*j-1)*pm1;
7737045ce4SJed Brown       dd   = pdd2 + (2*j-1)*pd1;
7837045ce4SJed Brown       pm2  = pm1;
7937045ce4SJed Brown       pm1  = p;
8037045ce4SJed Brown       pd2  = pd1;
8137045ce4SJed Brown       pd1  = d;
8237045ce4SJed Brown       pdd2 = pdd1;
8337045ce4SJed Brown       pdd1 = dd;
8437045ce4SJed Brown       if (degrees[k] == j) {
8537045ce4SJed Brown         if (B) B[i*ndegree+k] = p;
8637045ce4SJed Brown         if (D) D[i*ndegree+k] = d;
8737045ce4SJed Brown         if (D2) D2[i*ndegree+k] = dd;
8837045ce4SJed Brown       }
8937045ce4SJed Brown     }
9037045ce4SJed Brown   }
9137045ce4SJed Brown   PetscFunctionReturn(0);
9237045ce4SJed Brown }
9337045ce4SJed Brown 
9437045ce4SJed Brown #undef __FUNCT__
9537045ce4SJed Brown #define __FUNCT__ "PetscDTGaussQuadrature"
9637045ce4SJed Brown /*@
9737045ce4SJed Brown    PetscDTGaussQuadrature - create Gauss quadrature
9837045ce4SJed Brown 
9937045ce4SJed Brown    Not Collective
10037045ce4SJed Brown 
10137045ce4SJed Brown    Input Arguments:
10237045ce4SJed Brown +  npoints - number of points
10337045ce4SJed Brown .  a - left end of interval (often-1)
10437045ce4SJed Brown -  b - right end of interval (often +1)
10537045ce4SJed Brown 
10637045ce4SJed Brown    Output Arguments:
10737045ce4SJed Brown +  x - quadrature points
10837045ce4SJed Brown -  w - quadrature weights
10937045ce4SJed Brown 
11037045ce4SJed Brown    Level: intermediate
11137045ce4SJed Brown 
11237045ce4SJed Brown    References:
11337045ce4SJed Brown    Golub and Welsch, Calculation of Quadrature Rules, Math. Comp. 23(106), 221--230, 1969.
11437045ce4SJed Brown 
11537045ce4SJed Brown .seealso: PetscDTLegendreEval()
11637045ce4SJed Brown @*/
11737045ce4SJed Brown PetscErrorCode PetscDTGaussQuadrature(PetscInt npoints,PetscReal a,PetscReal b,PetscReal *x,PetscReal *w)
11837045ce4SJed Brown {
11937045ce4SJed Brown   PetscErrorCode ierr;
12037045ce4SJed Brown   PetscInt       i;
12137045ce4SJed Brown   PetscReal      *work;
12237045ce4SJed Brown   PetscScalar    *Z;
12337045ce4SJed Brown   PetscBLASInt   N,LDZ,info;
12437045ce4SJed Brown 
12537045ce4SJed Brown   PetscFunctionBegin;
12637045ce4SJed Brown   /* Set up the Golub-Welsch system */
12737045ce4SJed Brown   for (i=0; i<npoints; i++) {
12837045ce4SJed Brown     x[i] = 0;                   /* diagonal is 0 */
12937045ce4SJed Brown     if (i) w[i-1] = 0.5 / PetscSqrtReal(1 - 1./PetscSqr(2*i));
13037045ce4SJed Brown   }
13137045ce4SJed Brown   ierr = PetscRealView(npoints-1,w,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
13237045ce4SJed Brown   ierr = PetscMalloc2(npoints*npoints,PetscScalar,&Z,PetscMax(1,2*npoints-2),PetscReal,&work);CHKERRQ(ierr);
133c5df96a5SBarry Smith   ierr = PetscBLASIntCast(npoints,&N);CHKERRQ(ierr);
13437045ce4SJed Brown   LDZ  = N;
13537045ce4SJed Brown   ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
1368b83055fSJed Brown   PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("I",&N,x,w,Z,&LDZ,work,&info));
13737045ce4SJed Brown   ierr = PetscFPTrapPop();CHKERRQ(ierr);
1381c3d6f74SJed Brown   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEQR error");
13937045ce4SJed Brown 
14037045ce4SJed Brown   for (i=0; i<(npoints+1)/2; i++) {
14137045ce4SJed Brown     PetscReal y = 0.5 * (-x[i] + x[npoints-i-1]); /* enforces symmetry */
14237045ce4SJed Brown     x[i]           = (a+b)/2 - y*(b-a)/2;
14337045ce4SJed Brown     x[npoints-i-1] = (a+b)/2 + y*(b-a)/2;
1440d644c17SKarl Rupp 
14537045ce4SJed Brown     w[i] = w[npoints-1-i] = (b-a)*PetscSqr(0.5*PetscAbsScalar(Z[i*npoints] + Z[(npoints-i-1)*npoints]));
14637045ce4SJed Brown   }
14737045ce4SJed Brown   ierr = PetscFree2(Z,work);CHKERRQ(ierr);
14837045ce4SJed Brown   PetscFunctionReturn(0);
14937045ce4SJed Brown }
150194825f6SJed Brown 
151194825f6SJed Brown #undef __FUNCT__
152494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTFactorial_Internal"
153494e7359SMatthew G. Knepley /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x.
154494e7359SMatthew G. Knepley    Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */
155494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorial_Internal(PetscInt n, PetscReal *factorial)
156494e7359SMatthew G. Knepley {
157494e7359SMatthew G. Knepley   PetscReal f = 1.0;
158494e7359SMatthew G. Knepley   PetscInt  i;
159494e7359SMatthew G. Knepley 
160494e7359SMatthew G. Knepley   PetscFunctionBegin;
161494e7359SMatthew G. Knepley   for (i = 1; i < n+1; ++i) f *= i;
162494e7359SMatthew G. Knepley   *factorial = f;
163494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
164494e7359SMatthew G. Knepley }
165494e7359SMatthew G. Knepley 
166494e7359SMatthew G. Knepley #undef __FUNCT__
167494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTComputeJacobi"
168494e7359SMatthew G. Knepley /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x.
169494e7359SMatthew G. Knepley    Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */
170494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobi(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P)
171494e7359SMatthew G. Knepley {
172494e7359SMatthew G. Knepley   PetscReal apb, pn1, pn2;
173494e7359SMatthew G. Knepley   PetscInt  k;
174494e7359SMatthew G. Knepley 
175494e7359SMatthew G. Knepley   PetscFunctionBegin;
176494e7359SMatthew G. Knepley   if (!n) {*P = 1.0; PetscFunctionReturn(0);}
177494e7359SMatthew G. Knepley   if (n == 1) {*P = 0.5 * (a - b + (a + b + 2.0) * x); PetscFunctionReturn(0);}
178494e7359SMatthew G. Knepley   apb = a + b;
179494e7359SMatthew G. Knepley   pn2 = 1.0;
180494e7359SMatthew G. Knepley   pn1 = 0.5 * (a - b + (apb + 2.0) * x);
181494e7359SMatthew G. Knepley   *P  = 0.0;
182494e7359SMatthew G. Knepley   for (k = 2; k < n+1; ++k) {
183494e7359SMatthew G. Knepley     PetscReal a1 = 2.0 * k * (k + apb) * (2.0*k + apb - 2.0);
184494e7359SMatthew G. Knepley     PetscReal a2 = (2.0 * k + apb - 1.0) * (a*a - b*b);
185494e7359SMatthew G. Knepley     PetscReal a3 = (2.0 * k + apb - 2.0) * (2.0 * k + apb - 1.0) * (2.0 * k + apb);
186494e7359SMatthew G. Knepley     PetscReal a4 = 2.0 * (k + a - 1.0) * (k + b - 1.0) * (2.0 * k + apb);
187494e7359SMatthew G. Knepley 
188494e7359SMatthew G. Knepley     a2  = a2 / a1;
189494e7359SMatthew G. Knepley     a3  = a3 / a1;
190494e7359SMatthew G. Knepley     a4  = a4 / a1;
191494e7359SMatthew G. Knepley     *P  = (a2 + a3 * x) * pn1 - a4 * pn2;
192494e7359SMatthew G. Knepley     pn2 = pn1;
193494e7359SMatthew G. Knepley     pn1 = *P;
194494e7359SMatthew G. Knepley   }
195494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
196494e7359SMatthew G. Knepley }
197494e7359SMatthew G. Knepley 
198494e7359SMatthew G. Knepley #undef __FUNCT__
199494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTComputeJacobiDerivative"
200494e7359SMatthew G. Knepley /* Evaluates the first derivative of P_{n}^{a,b} at a point x. */
201494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobiDerivative(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P)
202494e7359SMatthew G. Knepley {
203494e7359SMatthew G. Knepley   PetscReal      nP;
204494e7359SMatthew G. Knepley   PetscErrorCode ierr;
205494e7359SMatthew G. Knepley 
206494e7359SMatthew G. Knepley   PetscFunctionBegin;
207494e7359SMatthew G. Knepley   if (!n) {*P = 0.0; PetscFunctionReturn(0);}
208494e7359SMatthew G. Knepley   ierr = PetscDTComputeJacobi(a+1, b+1, n-1, x, &nP);CHKERRQ(ierr);
209494e7359SMatthew G. Knepley   *P   = 0.5 * (a + b + n + 1) * nP;
210494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
211494e7359SMatthew G. Knepley }
212494e7359SMatthew G. Knepley 
213494e7359SMatthew G. Knepley #undef __FUNCT__
214494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTMapSquareToTriangle_Internal"
215494e7359SMatthew G. Knepley /* Maps from [-1,1]^2 to the (-1,1) reference triangle */
216494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTMapSquareToTriangle_Internal(PetscReal x, PetscReal y, PetscReal *xi, PetscReal *eta)
217494e7359SMatthew G. Knepley {
218494e7359SMatthew G. Knepley   PetscFunctionBegin;
219494e7359SMatthew G. Knepley   *xi  = 0.5 * (1.0 + x) * (1.0 - y) - 1.0;
220494e7359SMatthew G. Knepley   *eta = y;
221494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
222494e7359SMatthew G. Knepley }
223494e7359SMatthew G. Knepley 
224494e7359SMatthew G. Knepley #undef __FUNCT__
225494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTMapCubeToTetrahedron_Internal"
226494e7359SMatthew G. Knepley /* Maps from [-1,1]^2 to the (-1,1) reference triangle */
227494e7359SMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode PetscDTMapCubeToTetrahedron_Internal(PetscReal x, PetscReal y, PetscReal z, PetscReal *xi, PetscReal *eta, PetscReal *zeta)
228494e7359SMatthew G. Knepley {
229494e7359SMatthew G. Knepley   PetscFunctionBegin;
230494e7359SMatthew G. Knepley   *xi   = 0.25 * (1.0 + x) * (1.0 - y) * (1.0 - z) - 1.0;
231494e7359SMatthew G. Knepley   *eta  = 0.5  * (1.0 + y) * (1.0 - z) - 1.0;
232494e7359SMatthew G. Knepley   *zeta = z;
233494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
234494e7359SMatthew G. Knepley }
235494e7359SMatthew G. Knepley 
236494e7359SMatthew G. Knepley #undef __FUNCT__
237494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTGaussJacobiQuadrature1D_Internal"
238494e7359SMatthew G. Knepley static PetscErrorCode PetscDTGaussJacobiQuadrature1D_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal *x, PetscReal *w)
239494e7359SMatthew G. Knepley {
240494e7359SMatthew G. Knepley   PetscInt       maxIter = 100;
241494e7359SMatthew G. Knepley   PetscReal      eps     = 1.0e-8;
242a8291ba1SSatish Balay   PetscReal      a1, a2, a3, a4, a5, a6;
243494e7359SMatthew G. Knepley   PetscInt       k;
244494e7359SMatthew G. Knepley   PetscErrorCode ierr;
245494e7359SMatthew G. Knepley 
246494e7359SMatthew G. Knepley   PetscFunctionBegin;
247a8291ba1SSatish Balay 
248a8291ba1SSatish Balay   a1      = pow(2, a+b+1);
249a8291ba1SSatish Balay #if defined(PETSC_HAVE_TGAMMA)
250a8291ba1SSatish Balay   a2      = tgamma(a + npoints + 1);
251a8291ba1SSatish Balay   a3      = tgamma(b + npoints + 1);
252a8291ba1SSatish Balay   a4      = tgamma(a + b + npoints + 1);
253a8291ba1SSatish Balay #else
254a8291ba1SSatish Balay   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"tgamma() - math routine is unavailable.");
255a8291ba1SSatish Balay #endif
256a8291ba1SSatish Balay 
257494e7359SMatthew G. Knepley   ierr = PetscDTFactorial_Internal(npoints, &a5);CHKERRQ(ierr);
258494e7359SMatthew G. Knepley   a6   = a1 * a2 * a3 / a4 / a5;
259494e7359SMatthew G. Knepley   /* Computes the m roots of P_{m}^{a,b} on [-1,1] by Newton's method with Chebyshev points as initial guesses.
260494e7359SMatthew G. Knepley    Algorithm implemented from the pseudocode given by Karniadakis and Sherwin and Python in FIAT */
261494e7359SMatthew G. Knepley   for (k = 0; k < npoints; ++k) {
2627f1c68b3SMatthew G. Knepley     PetscReal r = -cos((2.0*k + 1.0) * PETSC_PI / (2.0 * npoints)), dP;
263494e7359SMatthew G. Knepley     PetscInt  j;
264494e7359SMatthew G. Knepley 
265494e7359SMatthew G. Knepley     if (k > 0) r = 0.5 * (r + x[k-1]);
266494e7359SMatthew G. Knepley     for (j = 0; j < maxIter; ++j) {
267494e7359SMatthew G. Knepley       PetscReal s = 0.0, delta, f, fp;
268494e7359SMatthew G. Knepley       PetscInt  i;
269494e7359SMatthew G. Knepley 
270494e7359SMatthew G. Knepley       for (i = 0; i < k; ++i) s = s + 1.0 / (r - x[i]);
271494e7359SMatthew G. Knepley       ierr = PetscDTComputeJacobi(a, b, npoints, r, &f);CHKERRQ(ierr);
272494e7359SMatthew G. Knepley       ierr = PetscDTComputeJacobiDerivative(a, b, npoints, r, &fp);CHKERRQ(ierr);
273494e7359SMatthew G. Knepley       delta = f / (fp - f * s);
274494e7359SMatthew G. Knepley       r     = r - delta;
275494e7359SMatthew G. Knepley       if (fabs(delta) < eps) break;
276494e7359SMatthew G. Knepley     }
277494e7359SMatthew G. Knepley     x[k] = r;
278494e7359SMatthew G. Knepley     ierr = PetscDTComputeJacobiDerivative(a, b, npoints, x[k], &dP);CHKERRQ(ierr);
279494e7359SMatthew G. Knepley     w[k] = a6 / (1.0 - PetscSqr(x[k])) / PetscSqr(dP);
280494e7359SMatthew G. Knepley   }
281494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
282494e7359SMatthew G. Knepley }
283494e7359SMatthew G. Knepley 
284494e7359SMatthew G. Knepley #undef __FUNCT__
285494e7359SMatthew G. Knepley #define __FUNCT__ "PetscDTGaussJacobiQuadrature"
286fd9d31fbSMatthew G. Knepley /*@C
287494e7359SMatthew G. Knepley   PetscDTGaussJacobiQuadrature - create Gauss-Jacobi quadrature for a simplex
288494e7359SMatthew G. Knepley 
289494e7359SMatthew G. Knepley   Not Collective
290494e7359SMatthew G. Knepley 
291494e7359SMatthew G. Knepley   Input Arguments:
292494e7359SMatthew G. Knepley + dim - The simplex dimension
293552aa4f7SMatthew G. Knepley . order - The quadrature order
294494e7359SMatthew G. Knepley . a - left end of interval (often-1)
295494e7359SMatthew G. Knepley - b - right end of interval (often +1)
296494e7359SMatthew G. Knepley 
297494e7359SMatthew G. Knepley   Output Arguments:
298552aa4f7SMatthew G. Knepley . q - A PetscQuadrature object
299494e7359SMatthew G. Knepley 
300494e7359SMatthew G. Knepley   Level: intermediate
301494e7359SMatthew G. Knepley 
302494e7359SMatthew G. Knepley   References:
303494e7359SMatthew G. Knepley   Karniadakis and Sherwin.
304494e7359SMatthew G. Knepley   FIAT
305494e7359SMatthew G. Knepley 
306494e7359SMatthew G. Knepley .seealso: PetscDTGaussQuadrature()
307494e7359SMatthew G. Knepley @*/
308552aa4f7SMatthew G. Knepley PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt dim, PetscInt order, PetscReal a, PetscReal b, PetscQuadrature *q)
309494e7359SMatthew G. Knepley {
310552aa4f7SMatthew G. Knepley   PetscInt       npoints = dim > 1 ? dim > 2 ? order*PetscSqr(order) : PetscSqr(order) : order;
311494e7359SMatthew G. Knepley   PetscReal     *px, *wx, *py, *wy, *pz, *wz, *x, *w;
312494e7359SMatthew G. Knepley   PetscInt       i, j, k;
313494e7359SMatthew G. Knepley   PetscErrorCode ierr;
314494e7359SMatthew G. Knepley 
315494e7359SMatthew G. Knepley   PetscFunctionBegin;
316494e7359SMatthew G. Knepley   if ((a != -1.0) || (b != 1.0)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must use default internal right now");
317552aa4f7SMatthew G. Knepley   ierr = PetscMalloc(npoints*dim * sizeof(PetscReal), &x);CHKERRQ(ierr);
318552aa4f7SMatthew G. Knepley   ierr = PetscMalloc(npoints     * sizeof(PetscReal), &w);CHKERRQ(ierr);
319494e7359SMatthew G. Knepley   switch (dim) {
320*707aa5c5SMatthew G. Knepley   case 0:
321*707aa5c5SMatthew G. Knepley     ierr = PetscFree(x);CHKERRQ(ierr);
322*707aa5c5SMatthew G. Knepley     ierr = PetscFree(w);CHKERRQ(ierr);
323*707aa5c5SMatthew G. Knepley     ierr = PetscMalloc(1 * sizeof(PetscReal), &x);CHKERRQ(ierr);
324*707aa5c5SMatthew G. Knepley     ierr = PetscMalloc(1 * sizeof(PetscReal), &w);CHKERRQ(ierr);
325*707aa5c5SMatthew G. Knepley     x[0] = 0.0;
326*707aa5c5SMatthew G. Knepley     w[0] = 1.0;
327*707aa5c5SMatthew G. Knepley     break;
328494e7359SMatthew G. Knepley   case 1:
329552aa4f7SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, x, w);CHKERRQ(ierr);
330494e7359SMatthew G. Knepley     break;
331494e7359SMatthew G. Knepley   case 2:
332552aa4f7SMatthew G. Knepley     ierr = PetscMalloc4(order,PetscReal,&px,order,PetscReal,&wx,order,PetscReal,&py,order,PetscReal,&wy);CHKERRQ(ierr);
333552aa4f7SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, px, wx);CHKERRQ(ierr);
334552aa4f7SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 1.0, 0.0, py, wy);CHKERRQ(ierr);
335552aa4f7SMatthew G. Knepley     for (i = 0; i < order; ++i) {
336552aa4f7SMatthew G. Knepley       for (j = 0; j < order; ++j) {
337552aa4f7SMatthew G. Knepley         ierr = PetscDTMapSquareToTriangle_Internal(px[i], py[j], &x[(i*order+j)*2+0], &x[(i*order+j)*2+1]);CHKERRQ(ierr);
338552aa4f7SMatthew G. Knepley         w[i*order+j] = 0.5 * wx[i] * wy[j];
339494e7359SMatthew G. Knepley       }
340494e7359SMatthew G. Knepley     }
341494e7359SMatthew G. Knepley     ierr = PetscFree4(px,wx,py,wy);CHKERRQ(ierr);
342494e7359SMatthew G. Knepley     break;
343494e7359SMatthew G. Knepley   case 3:
344552aa4f7SMatthew G. Knepley     ierr = PetscMalloc6(order,PetscReal,&px,order,PetscReal,&wx,order,PetscReal,&py,order,PetscReal,&wy,order,PetscReal,&pz,order,PetscReal,&wz);CHKERRQ(ierr);
345552aa4f7SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 0.0, 0.0, px, wx);CHKERRQ(ierr);
346552aa4f7SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 1.0, 0.0, py, wy);CHKERRQ(ierr);
347552aa4f7SMatthew G. Knepley     ierr = PetscDTGaussJacobiQuadrature1D_Internal(order, 2.0, 0.0, pz, wz);CHKERRQ(ierr);
348552aa4f7SMatthew G. Knepley     for (i = 0; i < order; ++i) {
349552aa4f7SMatthew G. Knepley       for (j = 0; j < order; ++j) {
350552aa4f7SMatthew G. Knepley         for (k = 0; k < order; ++k) {
351552aa4f7SMatthew 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);
352552aa4f7SMatthew G. Knepley           w[(i*order+j)*order+k] = 0.125 * wx[i] * wy[j] * wz[k];
353494e7359SMatthew G. Knepley         }
354494e7359SMatthew G. Knepley       }
355494e7359SMatthew G. Knepley     }
356494e7359SMatthew G. Knepley     ierr = PetscFree6(px,wx,py,wy,pz,wz);CHKERRQ(ierr);
357494e7359SMatthew G. Knepley     break;
358494e7359SMatthew G. Knepley   default:
359494e7359SMatthew G. Knepley     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim);
360494e7359SMatthew G. Knepley   }
361552aa4f7SMatthew G. Knepley   q->numQuadPoints = npoints;
362552aa4f7SMatthew G. Knepley   q->quadPoints    = x;
363552aa4f7SMatthew G. Knepley   q->quadWeights   = w;
364494e7359SMatthew G. Knepley   PetscFunctionReturn(0);
365494e7359SMatthew G. Knepley }
366494e7359SMatthew G. Knepley 
367494e7359SMatthew G. Knepley #undef __FUNCT__
368194825f6SJed Brown #define __FUNCT__ "PetscDTPseudoInverseQR"
369194825f6SJed Brown /* Overwrites A. Can only handle full-rank problems with m>=n
370194825f6SJed Brown  * A in column-major format
371194825f6SJed Brown  * Ainv in row-major format
372194825f6SJed Brown  * tau has length m
373194825f6SJed Brown  * worksize must be >= max(1,n)
374194825f6SJed Brown  */
375194825f6SJed Brown static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m,PetscInt mstride,PetscInt n,PetscReal *A_in,PetscReal *Ainv_out,PetscScalar *tau,PetscInt worksize,PetscScalar *work)
376194825f6SJed Brown {
377194825f6SJed Brown   PetscErrorCode ierr;
378194825f6SJed Brown   PetscBLASInt M,N,K,lda,ldb,ldwork,info;
379194825f6SJed Brown   PetscScalar *A,*Ainv,*R,*Q,Alpha;
380194825f6SJed Brown 
381194825f6SJed Brown   PetscFunctionBegin;
382194825f6SJed Brown #if defined(PETSC_USE_COMPLEX)
383194825f6SJed Brown   {
384194825f6SJed Brown     PetscInt i,j;
385194825f6SJed Brown     ierr = PetscMalloc2(m*n,PetscScalar,&A,m*n,PetscScalar,&Ainv);CHKERRQ(ierr);
386194825f6SJed Brown     for (j=0; j<n; j++) {
387194825f6SJed Brown       for (i=0; i<m; i++) A[i+m*j] = A_in[i+mstride*j];
388194825f6SJed Brown     }
389194825f6SJed Brown     mstride = m;
390194825f6SJed Brown   }
391194825f6SJed Brown #else
392194825f6SJed Brown   A = A_in;
393194825f6SJed Brown   Ainv = Ainv_out;
394194825f6SJed Brown #endif
395194825f6SJed Brown 
396194825f6SJed Brown   ierr = PetscBLASIntCast(m,&M);CHKERRQ(ierr);
397194825f6SJed Brown   ierr = PetscBLASIntCast(n,&N);CHKERRQ(ierr);
398194825f6SJed Brown   ierr = PetscBLASIntCast(mstride,&lda);CHKERRQ(ierr);
399194825f6SJed Brown   ierr = PetscBLASIntCast(worksize,&ldwork);CHKERRQ(ierr);
400194825f6SJed Brown   ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
401194825f6SJed Brown   LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info);
402194825f6SJed Brown   ierr = PetscFPTrapPop();CHKERRQ(ierr);
403194825f6SJed Brown   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error");
404194825f6SJed Brown   R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */
405194825f6SJed Brown 
406194825f6SJed Brown   /* Extract an explicit representation of Q */
407194825f6SJed Brown   Q = Ainv;
408194825f6SJed Brown   ierr = PetscMemcpy(Q,A,mstride*n*sizeof(PetscScalar));CHKERRQ(ierr);
409194825f6SJed Brown   K = N;                        /* full rank */
410194825f6SJed Brown   LAPACKungqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info);
411194825f6SJed Brown   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error");
412194825f6SJed Brown 
413194825f6SJed Brown   /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */
414194825f6SJed Brown   Alpha = 1.0;
415194825f6SJed Brown   ldb = lda;
416194825f6SJed Brown   BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb);
417194825f6SJed Brown   /* Ainv is Q, overwritten with inverse */
418194825f6SJed Brown 
419194825f6SJed Brown #if defined(PETSC_USE_COMPLEX)
420194825f6SJed Brown   {
421194825f6SJed Brown     PetscInt i;
422194825f6SJed Brown     for (i=0; i<m*n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]);
423194825f6SJed Brown     ierr = PetscFree2(A,Ainv);CHKERRQ(ierr);
424194825f6SJed Brown   }
425194825f6SJed Brown #endif
426194825f6SJed Brown   PetscFunctionReturn(0);
427194825f6SJed Brown }
428194825f6SJed Brown 
429194825f6SJed Brown #undef __FUNCT__
430194825f6SJed Brown #define __FUNCT__ "PetscDTLegendreIntegrate"
431194825f6SJed Brown /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */
432194825f6SJed Brown static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval,const PetscReal *x,PetscInt ndegree,const PetscInt *degrees,PetscBool Transpose,PetscReal *B)
433194825f6SJed Brown {
434194825f6SJed Brown   PetscErrorCode ierr;
435194825f6SJed Brown   PetscReal *Bv;
436194825f6SJed Brown   PetscInt i,j;
437194825f6SJed Brown 
438194825f6SJed Brown   PetscFunctionBegin;
439194825f6SJed Brown   ierr = PetscMalloc((ninterval+1)*ndegree*sizeof(PetscReal),&Bv);CHKERRQ(ierr);
440194825f6SJed Brown   /* Point evaluation of L_p on all the source vertices */
441194825f6SJed Brown   ierr = PetscDTLegendreEval(ninterval+1,x,ndegree,degrees,Bv,NULL,NULL);CHKERRQ(ierr);
442194825f6SJed Brown   /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */
443194825f6SJed Brown   for (i=0; i<ninterval; i++) {
444194825f6SJed Brown     for (j=0; j<ndegree; j++) {
445194825f6SJed Brown       if (Transpose) B[i+ninterval*j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j];
446194825f6SJed Brown       else           B[i*ndegree+j]   = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j];
447194825f6SJed Brown     }
448194825f6SJed Brown   }
449194825f6SJed Brown   ierr = PetscFree(Bv);CHKERRQ(ierr);
450194825f6SJed Brown   PetscFunctionReturn(0);
451194825f6SJed Brown }
452194825f6SJed Brown 
453194825f6SJed Brown #undef __FUNCT__
454194825f6SJed Brown #define __FUNCT__ "PetscDTReconstructPoly"
455194825f6SJed Brown /*@
456194825f6SJed Brown    PetscDTReconstructPoly - create matrix representing polynomial reconstruction using cell intervals and evaluation at target intervals
457194825f6SJed Brown 
458194825f6SJed Brown    Not Collective
459194825f6SJed Brown 
460194825f6SJed Brown    Input Arguments:
461194825f6SJed Brown +  degree - degree of reconstruction polynomial
462194825f6SJed Brown .  nsource - number of source intervals
463194825f6SJed Brown .  sourcex - sorted coordinates of source cell boundaries (length nsource+1)
464194825f6SJed Brown .  ntarget - number of target intervals
465194825f6SJed Brown -  targetx - sorted coordinates of target cell boundaries (length ntarget+1)
466194825f6SJed Brown 
467194825f6SJed Brown    Output Arguments:
468194825f6SJed Brown .  R - reconstruction matrix, utarget = sum_s R[t*nsource+s] * usource[s]
469194825f6SJed Brown 
470194825f6SJed Brown    Level: advanced
471194825f6SJed Brown 
472194825f6SJed Brown .seealso: PetscDTLegendreEval()
473194825f6SJed Brown @*/
474194825f6SJed Brown PetscErrorCode PetscDTReconstructPoly(PetscInt degree,PetscInt nsource,const PetscReal *sourcex,PetscInt ntarget,const PetscReal *targetx,PetscReal *R)
475194825f6SJed Brown {
476194825f6SJed Brown   PetscErrorCode ierr;
477194825f6SJed Brown   PetscInt i,j,k,*bdegrees,worksize;
478194825f6SJed Brown   PetscReal xmin,xmax,center,hscale,*sourcey,*targety,*Bsource,*Bsinv,*Btarget;
479194825f6SJed Brown   PetscScalar *tau,*work;
480194825f6SJed Brown 
481194825f6SJed Brown   PetscFunctionBegin;
482194825f6SJed Brown   PetscValidRealPointer(sourcex,3);
483194825f6SJed Brown   PetscValidRealPointer(targetx,5);
484194825f6SJed Brown   PetscValidRealPointer(R,6);
485194825f6SJed 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);
486194825f6SJed Brown #if defined(PETSC_USE_DEBUG)
487194825f6SJed Brown   for (i=0; i<nsource; i++) {
488194825f6SJed 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]);
489194825f6SJed Brown   }
490194825f6SJed Brown   for (i=0; i<ntarget; i++) {
491194825f6SJed 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]);
492194825f6SJed Brown   }
493194825f6SJed Brown #endif
494194825f6SJed Brown   xmin = PetscMin(sourcex[0],targetx[0]);
495194825f6SJed Brown   xmax = PetscMax(sourcex[nsource],targetx[ntarget]);
496194825f6SJed Brown   center = (xmin + xmax)/2;
497194825f6SJed Brown   hscale = (xmax - xmin)/2;
498194825f6SJed Brown   worksize = nsource;
499194825f6SJed Brown   ierr = PetscMalloc4(degree+1,PetscInt,&bdegrees,nsource+1,PetscReal,&sourcey,nsource*(degree+1),PetscReal,&Bsource,worksize,PetscScalar,&work);CHKERRQ(ierr);
50082772646SJed Brown   ierr = PetscMalloc4(nsource,PetscScalar,&tau,nsource*(degree+1),PetscReal,&Bsinv,ntarget+1,PetscReal,&targety,ntarget*(degree+1),PetscReal,&Btarget);CHKERRQ(ierr);
501194825f6SJed Brown   for (i=0; i<=nsource; i++) sourcey[i] = (sourcex[i]-center)/hscale;
502194825f6SJed Brown   for (i=0; i<=degree; i++) bdegrees[i] = i+1;
503194825f6SJed Brown   ierr = PetscDTLegendreIntegrate(nsource,sourcey,degree+1,bdegrees,PETSC_TRUE,Bsource);CHKERRQ(ierr);
504194825f6SJed Brown   ierr = PetscDTPseudoInverseQR(nsource,nsource,degree+1,Bsource,Bsinv,tau,nsource,work);CHKERRQ(ierr);
505194825f6SJed Brown   for (i=0; i<=ntarget; i++) targety[i] = (targetx[i]-center)/hscale;
506194825f6SJed Brown   ierr = PetscDTLegendreIntegrate(ntarget,targety,degree+1,bdegrees,PETSC_FALSE,Btarget);CHKERRQ(ierr);
507194825f6SJed Brown   for (i=0; i<ntarget; i++) {
508194825f6SJed Brown     PetscReal rowsum = 0;
509194825f6SJed Brown     for (j=0; j<nsource; j++) {
510194825f6SJed Brown       PetscReal sum = 0;
511194825f6SJed Brown       for (k=0; k<degree+1; k++) {
512194825f6SJed Brown         sum += Btarget[i*(degree+1)+k] * Bsinv[k*nsource+j];
513194825f6SJed Brown       }
514194825f6SJed Brown       R[i*nsource+j] = sum;
515194825f6SJed Brown       rowsum += sum;
516194825f6SJed Brown     }
517194825f6SJed Brown     for (j=0; j<nsource; j++) R[i*nsource+j] /= rowsum; /* normalize each row */
518194825f6SJed Brown   }
519194825f6SJed Brown   ierr = PetscFree4(bdegrees,sourcey,Bsource,work);CHKERRQ(ierr);
520194825f6SJed Brown   ierr = PetscFree4(tau,Bsinv,targety,Btarget);CHKERRQ(ierr);
521194825f6SJed Brown   PetscFunctionReturn(0);
522194825f6SJed Brown }
523