Actual source code: dt.c

petsc-3.4.1 2013-06-10
  1: /* Discretization tools */

  3: #include <petscdt.h>            /*I "petscdt.h" I*/
  4: #include <petscblaslapack.h>
  5: #include <petsc-private/petscimpl.h>
  6: #include <petscviewer.h>

 10: /*@
 11:    PetscDTLegendreEval - evaluate Legendre polynomial at points

 13:    Not Collective

 15:    Input Arguments:
 16: +  npoints - number of spatial points to evaluate at
 17: .  points - array of locations to evaluate at
 18: .  ndegree - number of basis degrees to evaluate
 19: -  degrees - sorted array of degrees to evaluate

 21:    Output Arguments:
 22: +  B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL)
 23: .  D - row-oriented derivative evaluation matrix (or NULL)
 24: -  D2 - row-oriented second derivative evaluation matrix (or NULL)

 26:    Level: intermediate

 28: .seealso: PetscDTGaussQuadrature()
 29: @*/
 30: PetscErrorCode PetscDTLegendreEval(PetscInt npoints,const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2)
 31: {
 32:   PetscInt i,maxdegree;

 35:   if (!npoints || !ndegree) return(0);
 36:   maxdegree = degrees[ndegree-1];
 37:   for (i=0; i<npoints; i++) {
 38:     PetscReal pm1,pm2,pd1,pd2,pdd1,pdd2,x;
 39:     PetscInt  j,k;
 40:     x    = points[i];
 41:     pm2  = 0;
 42:     pm1  = 1;
 43:     pd2  = 0;
 44:     pd1  = 0;
 45:     pdd2 = 0;
 46:     pdd1 = 0;
 47:     k    = 0;
 48:     if (degrees[k] == 0) {
 49:       if (B) B[i*ndegree+k] = pm1;
 50:       if (D) D[i*ndegree+k] = pd1;
 51:       if (D2) D2[i*ndegree+k] = pdd1;
 52:       k++;
 53:     }
 54:     for (j=1; j<=maxdegree; j++,k++) {
 55:       PetscReal p,d,dd;
 56:       p    = ((2*j-1)*x*pm1 - (j-1)*pm2)/j;
 57:       d    = pd2 + (2*j-1)*pm1;
 58:       dd   = pdd2 + (2*j-1)*pd1;
 59:       pm2  = pm1;
 60:       pm1  = p;
 61:       pd2  = pd1;
 62:       pd1  = d;
 63:       pdd2 = pdd1;
 64:       pdd1 = dd;
 65:       if (degrees[k] == j) {
 66:         if (B) B[i*ndegree+k] = p;
 67:         if (D) D[i*ndegree+k] = d;
 68:         if (D2) D2[i*ndegree+k] = dd;
 69:       }
 70:     }
 71:   }
 72:   return(0);
 73: }

 77: /*@
 78:    PetscDTGaussQuadrature - create Gauss quadrature

 80:    Not Collective

 82:    Input Arguments:
 83: +  npoints - number of points
 84: .  a - left end of interval (often-1)
 85: -  b - right end of interval (often +1)

 87:    Output Arguments:
 88: +  x - quadrature points
 89: -  w - quadrature weights

 91:    Level: intermediate

 93:    References:
 94:    Golub and Welsch, Calculation of Quadrature Rules, Math. Comp. 23(106), 221--230, 1969.

 96: .seealso: PetscDTLegendreEval()
 97: @*/
 98: PetscErrorCode PetscDTGaussQuadrature(PetscInt npoints,PetscReal a,PetscReal b,PetscReal *x,PetscReal *w)
 99: {
101:   PetscInt       i;
102:   PetscReal      *work;
103:   PetscScalar    *Z;
104:   PetscBLASInt   N,LDZ,info;

107:   /* Set up the Golub-Welsch system */
108:   for (i=0; i<npoints; i++) {
109:     x[i] = 0;                   /* diagonal is 0 */
110:     if (i) w[i-1] = 0.5 / PetscSqrtReal(1 - 1./PetscSqr(2*i));
111:   }
112:   PetscRealView(npoints-1,w,PETSC_VIEWER_STDOUT_SELF);
113:   PetscMalloc2(npoints*npoints,PetscScalar,&Z,PetscMax(1,2*npoints-2),PetscReal,&work);
114:   PetscBLASIntCast(npoints,&N);
115:   LDZ  = N;
116:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);
117:   PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("I",&N,x,w,Z,&LDZ,work,&info));
118:   PetscFPTrapPop();
119:   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEQR error");

121:   for (i=0; i<(npoints+1)/2; i++) {
122:     PetscReal y = 0.5 * (-x[i] + x[npoints-i-1]); /* enforces symmetry */
123:     x[i]           = (a+b)/2 - y*(b-a)/2;
124:     x[npoints-i-1] = (a+b)/2 + y*(b-a)/2;

126:     w[i] = w[npoints-1-i] = (b-a)*PetscSqr(0.5*PetscAbsScalar(Z[i*npoints] + Z[(npoints-i-1)*npoints]));
127:   }
128:   PetscFree2(Z,work);
129:   return(0);
130: }

134: /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x.
135:    Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */
136: PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorial_Internal(PetscInt n, PetscReal *factorial)
137: {
138:   PetscReal f = 1.0;
139:   PetscInt  i;

142:   for (i = 1; i < n+1; ++i) f *= i;
143:   *factorial = f;
144:   return(0);
145: }

149: /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x.
150:    Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */
151: PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobi(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P)
152: {
153:   PetscReal apb, pn1, pn2;
154:   PetscInt  k;

157:   if (!n) {*P = 1.0; return(0);}
158:   if (n == 1) {*P = 0.5 * (a - b + (a + b + 2.0) * x); return(0);}
159:   apb = a + b;
160:   pn2 = 1.0;
161:   pn1 = 0.5 * (a - b + (apb + 2.0) * x);
162:   *P  = 0.0;
163:   for (k = 2; k < n+1; ++k) {
164:     PetscReal a1 = 2.0 * k * (k + apb) * (2.0*k + apb - 2.0);
165:     PetscReal a2 = (2.0 * k + apb - 1.0) * (a*a - b*b);
166:     PetscReal a3 = (2.0 * k + apb - 2.0) * (2.0 * k + apb - 1.0) * (2.0 * k + apb);
167:     PetscReal a4 = 2.0 * (k + a - 1.0) * (k + b - 1.0) * (2.0 * k + apb);

169:     a2  = a2 / a1;
170:     a3  = a3 / a1;
171:     a4  = a4 / a1;
172:     *P  = (a2 + a3 * x) * pn1 - a4 * pn2;
173:     pn2 = pn1;
174:     pn1 = *P;
175:   }
176:   return(0);
177: }

181: /* Evaluates the first derivative of P_{n}^{a,b} at a point x. */
182: PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobiDerivative(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P)
183: {
184:   PetscReal      nP;

188:   if (!n) {*P = 0.0; return(0);}
189:   PetscDTComputeJacobi(a+1, b+1, n-1, x, &nP);
190:   *P   = 0.5 * (a + b + n + 1) * nP;
191:   return(0);
192: }

196: /* Maps from [-1,1]^2 to the (-1,1) reference triangle */
197: PETSC_STATIC_INLINE PetscErrorCode PetscDTMapSquareToTriangle_Internal(PetscReal x, PetscReal y, PetscReal *xi, PetscReal *eta)
198: {
200:   *xi  = 0.5 * (1.0 + x) * (1.0 - y) - 1.0;
201:   *eta = y;
202:   return(0);
203: }

207: /* Maps from [-1,1]^2 to the (-1,1) reference triangle */
208: PETSC_STATIC_INLINE PetscErrorCode PetscDTMapCubeToTetrahedron_Internal(PetscReal x, PetscReal y, PetscReal z, PetscReal *xi, PetscReal *eta, PetscReal *zeta)
209: {
211:   *xi   = 0.25 * (1.0 + x) * (1.0 - y) * (1.0 - z) - 1.0;
212:   *eta  = 0.5  * (1.0 + y) * (1.0 - z) - 1.0;
213:   *zeta = z;
214:   return(0);
215: }

219: static PetscErrorCode PetscDTGaussJacobiQuadrature1D_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal *x, PetscReal *w)
220: {
221:   PetscInt       maxIter = 100;
222:   PetscReal      eps     = 1.0e-8;
223:   PetscReal      a1, a2, a3, a4, a5, a6;
224:   PetscInt       k;


229:   a1      = pow(2, a+b+1);
230: #if defined(PETSC_HAVE_TGAMMA)
231:   a2      = tgamma(a + npoints + 1);
232:   a3      = tgamma(b + npoints + 1);
233:   a4      = tgamma(a + b + npoints + 1);
234: #else
235:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"tgamma() - math routine is unavailable.");
236: #endif

238:   PetscDTFactorial_Internal(npoints, &a5);
239:   a6   = a1 * a2 * a3 / a4 / a5;
240:   /* Computes the m roots of P_{m}^{a,b} on [-1,1] by Newton's method with Chebyshev points as initial guesses.
241:    Algorithm implemented from the pseudocode given by Karniadakis and Sherwin and Python in FIAT */
242:   for (k = 0; k < npoints; ++k) {
243:     PetscReal r = -cos((2.0*k + 1.0) * PETSC_PI / (2.0 * npoints)), dP;
244:     PetscInt  j;

246:     if (k > 0) r = 0.5 * (r + x[k-1]);
247:     for (j = 0; j < maxIter; ++j) {
248:       PetscReal s = 0.0, delta, f, fp;
249:       PetscInt  i;

251:       for (i = 0; i < k; ++i) s = s + 1.0 / (r - x[i]);
252:       PetscDTComputeJacobi(a, b, npoints, r, &f);
253:       PetscDTComputeJacobiDerivative(a, b, npoints, r, &fp);
254:       delta = f / (fp - f * s);
255:       r     = r - delta;
256:       if (fabs(delta) < eps) break;
257:     }
258:     x[k] = r;
259:     PetscDTComputeJacobiDerivative(a, b, npoints, x[k], &dP);
260:     w[k] = a6 / (1.0 - PetscSqr(x[k])) / PetscSqr(dP);
261:   }
262:   return(0);
263: }

267: /*@C
268:   PetscDTGaussJacobiQuadrature - create Gauss-Jacobi quadrature for a simplex

270:   Not Collective

272:   Input Arguments:
273: + dim - The simplex dimension
274: . npoints - number of points
275: . a - left end of interval (often-1)
276: - b - right end of interval (often +1)

278:   Output Arguments:
279: + points - quadrature points
280: - weights - quadrature weights

282:   Level: intermediate

284:   References:
285:   Karniadakis and Sherwin.
286:   FIAT

288: .seealso: PetscDTGaussQuadrature()
289: @*/
290: PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt dim, PetscInt npoints, PetscReal a, PetscReal b, PetscReal *points[], PetscReal *weights[])
291: {
292:   PetscReal     *px, *wx, *py, *wy, *pz, *wz, *x, *w;
293:   PetscInt       i, j, k;

297:   if ((a != -1.0) || (b != 1.0)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must use default internal right now");
298:   switch (dim) {
299:   case 1:
300:     PetscMalloc(npoints * sizeof(PetscReal), &x);
301:     PetscMalloc(npoints * sizeof(PetscReal), &w);
302:     PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, x, w);
303:     break;
304:   case 2:
305:     PetscMalloc(npoints*npoints*2 * sizeof(PetscReal), &x);
306:     PetscMalloc(npoints*npoints   * sizeof(PetscReal), &w);
307:     PetscMalloc4(npoints,PetscReal,&px,npoints,PetscReal,&wx,npoints,PetscReal,&py,npoints,PetscReal,&wy);
308:     PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, px, wx);
309:     PetscDTGaussJacobiQuadrature1D_Internal(npoints, 1.0, 0.0, py, wy);
310:     for (i = 0; i < npoints; ++i) {
311:       for (j = 0; j < npoints; ++j) {
312:         PetscDTMapSquareToTriangle_Internal(px[i], py[j], &x[(i*npoints+j)*2+0], &x[(i*npoints+j)*2+1]);
313:         w[i*npoints+j] = 0.5 * wx[i] * wy[j];
314:       }
315:     }
316:     PetscFree4(px,wx,py,wy);
317:     break;
318:   case 3:
319:     PetscMalloc(npoints*npoints*3 * sizeof(PetscReal), &x);
320:     PetscMalloc(npoints*npoints   * sizeof(PetscReal), &w);
321:     PetscMalloc6(npoints,PetscReal,&px,npoints,PetscReal,&wx,npoints,PetscReal,&py,npoints,PetscReal,&wy,npoints,PetscReal,&pz,npoints,PetscReal,&wz);
322:     PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, px, wx);
323:     PetscDTGaussJacobiQuadrature1D_Internal(npoints, 1.0, 0.0, py, wy);
324:     PetscDTGaussJacobiQuadrature1D_Internal(npoints, 2.0, 0.0, pz, wz);
325:     for (i = 0; i < npoints; ++i) {
326:       for (j = 0; j < npoints; ++j) {
327:         for (k = 0; k < npoints; ++k) {
328:           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]);
329:           w[(i*npoints+j)*npoints+k] = 0.125 * wx[i] * wy[j] * wz[k];
330:         }
331:       }
332:     }
333:     PetscFree6(px,wx,py,wy,pz,wz);
334:     break;
335:   default:
336:     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim);
337:   }
338:   if (points)  *points  = x;
339:   if (weights) *weights = w;
340:   return(0);
341: }

345: /* Overwrites A. Can only handle full-rank problems with m>=n
346:  * A in column-major format
347:  * Ainv in row-major format
348:  * tau has length m
349:  * worksize must be >= max(1,n)
350:  */
351: static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m,PetscInt mstride,PetscInt n,PetscReal *A_in,PetscReal *Ainv_out,PetscScalar *tau,PetscInt worksize,PetscScalar *work)
352: {
354:   PetscBLASInt M,N,K,lda,ldb,ldwork,info;
355:   PetscScalar *A,*Ainv,*R,*Q,Alpha;

358: #if defined(PETSC_USE_COMPLEX)
359:   {
360:     PetscInt i,j;
361:     PetscMalloc2(m*n,PetscScalar,&A,m*n,PetscScalar,&Ainv);
362:     for (j=0; j<n; j++) {
363:       for (i=0; i<m; i++) A[i+m*j] = A_in[i+mstride*j];
364:     }
365:     mstride = m;
366:   }
367: #else
368:   A = A_in;
369:   Ainv = Ainv_out;
370: #endif

372:   PetscBLASIntCast(m,&M);
373:   PetscBLASIntCast(n,&N);
374:   PetscBLASIntCast(mstride,&lda);
375:   PetscBLASIntCast(worksize,&ldwork);
376:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);
377:   LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info);
378:   PetscFPTrapPop();
379:   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error");
380:   R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */

382:   /* Extract an explicit representation of Q */
383:   Q = Ainv;
384:   PetscMemcpy(Q,A,mstride*n*sizeof(PetscScalar));
385:   K = N;                        /* full rank */
386:   LAPACKungqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info);
387:   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error");

389:   /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */
390:   Alpha = 1.0;
391:   ldb = lda;
392:   BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb);
393:   /* Ainv is Q, overwritten with inverse */

395: #if defined(PETSC_USE_COMPLEX)
396:   {
397:     PetscInt i;
398:     for (i=0; i<m*n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]);
399:     PetscFree2(A,Ainv);
400:   }
401: #endif
402:   return(0);
403: }

407: /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */
408: static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval,const PetscReal *x,PetscInt ndegree,const PetscInt *degrees,PetscBool Transpose,PetscReal *B)
409: {
411:   PetscReal *Bv;
412:   PetscInt i,j;

415:   PetscMalloc((ninterval+1)*ndegree*sizeof(PetscReal),&Bv);
416:   /* Point evaluation of L_p on all the source vertices */
417:   PetscDTLegendreEval(ninterval+1,x,ndegree,degrees,Bv,NULL,NULL);
418:   /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */
419:   for (i=0; i<ninterval; i++) {
420:     for (j=0; j<ndegree; j++) {
421:       if (Transpose) B[i+ninterval*j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j];
422:       else           B[i*ndegree+j]   = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j];
423:     }
424:   }
425:   PetscFree(Bv);
426:   return(0);
427: }

431: /*@
432:    PetscDTReconstructPoly - create matrix representing polynomial reconstruction using cell intervals and evaluation at target intervals

434:    Not Collective

436:    Input Arguments:
437: +  degree - degree of reconstruction polynomial
438: .  nsource - number of source intervals
439: .  sourcex - sorted coordinates of source cell boundaries (length nsource+1)
440: .  ntarget - number of target intervals
441: -  targetx - sorted coordinates of target cell boundaries (length ntarget+1)

443:    Output Arguments:
444: .  R - reconstruction matrix, utarget = sum_s R[t*nsource+s] * usource[s]

446:    Level: advanced

448: .seealso: PetscDTLegendreEval()
449: @*/
450: PetscErrorCode PetscDTReconstructPoly(PetscInt degree,PetscInt nsource,const PetscReal *sourcex,PetscInt ntarget,const PetscReal *targetx,PetscReal *R)
451: {
453:   PetscInt i,j,k,*bdegrees,worksize;
454:   PetscReal xmin,xmax,center,hscale,*sourcey,*targety,*Bsource,*Bsinv,*Btarget;
455:   PetscScalar *tau,*work;

461:   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);
462: #if defined(PETSC_USE_DEBUG)
463:   for (i=0; i<nsource; i++) {
464:     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]);
465:   }
466:   for (i=0; i<ntarget; i++) {
467:     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]);
468:   }
469: #endif
470:   xmin = PetscMin(sourcex[0],targetx[0]);
471:   xmax = PetscMax(sourcex[nsource],targetx[ntarget]);
472:   center = (xmin + xmax)/2;
473:   hscale = (xmax - xmin)/2;
474:   worksize = nsource;
475:   PetscMalloc4(degree+1,PetscInt,&bdegrees,nsource+1,PetscReal,&sourcey,nsource*(degree+1),PetscReal,&Bsource,worksize,PetscScalar,&work);
476:   PetscMalloc4(nsource,PetscScalar,&tau,nsource*(degree+1),PetscReal,&Bsinv,ntarget+1,PetscReal,&targety,ntarget*(degree+1),PetscReal,&Btarget);
477:   for (i=0; i<=nsource; i++) sourcey[i] = (sourcex[i]-center)/hscale;
478:   for (i=0; i<=degree; i++) bdegrees[i] = i+1;
479:   PetscDTLegendreIntegrate(nsource,sourcey,degree+1,bdegrees,PETSC_TRUE,Bsource);
480:   PetscDTPseudoInverseQR(nsource,nsource,degree+1,Bsource,Bsinv,tau,nsource,work);
481:   for (i=0; i<=ntarget; i++) targety[i] = (targetx[i]-center)/hscale;
482:   PetscDTLegendreIntegrate(ntarget,targety,degree+1,bdegrees,PETSC_FALSE,Btarget);
483:   for (i=0; i<ntarget; i++) {
484:     PetscReal rowsum = 0;
485:     for (j=0; j<nsource; j++) {
486:       PetscReal sum = 0;
487:       for (k=0; k<degree+1; k++) {
488:         sum += Btarget[i*(degree+1)+k] * Bsinv[k*nsource+j];
489:       }
490:       R[i*nsource+j] = sum;
491:       rowsum += sum;
492:     }
493:     for (j=0; j<nsource; j++) R[i*nsource+j] /= rowsum; /* normalize each row */
494:   }
495:   PetscFree4(bdegrees,sourcey,Bsource,work);
496:   PetscFree4(tau,Bsinv,targety,Btarget);
497:   return(0);
498: }