Actual source code: dt.c

petsc-3.9.4 2018-09-11
Report Typos and Errors
  1: /* Discretization tools */

  3:  #include <petscdt.h>
  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>

 11: #if defined(PETSC_HAVE_MPFR)
 12: #include <mpfr.h>
 13: #endif

 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";

 25: /*@
 26:   PetscQuadratureCreate - Create a PetscQuadrature object

 28:   Collective on MPI_Comm

 30:   Input Parameter:
 31: . comm - The communicator for the PetscQuadrature object

 33:   Output Parameter:
 34: . q  - The PetscQuadrature object

 36:   Level: beginner

 38: .keywords: PetscQuadrature, quadrature, create
 39: .seealso: PetscQuadratureDestroy(), PetscQuadratureGetData()
 40: @*/
 41: PetscErrorCode PetscQuadratureCreate(MPI_Comm comm, PetscQuadrature *q)
 42: {

 47:   PetscSysInitializePackage();
 48:   PetscHeaderCreate(*q,PETSC_OBJECT_CLASSID,"PetscQuadrature","Quadrature","DT",comm,PetscQuadratureDestroy,PetscQuadratureView);
 49:   (*q)->dim       = -1;
 50:   (*q)->Nc        =  1;
 51:   (*q)->order     = -1;
 52:   (*q)->numPoints = 0;
 53:   (*q)->points    = NULL;
 54:   (*q)->weights   = NULL;
 55:   return(0);
 56: }

 58: /*@
 59:   PetscQuadratureDuplicate - Create a deep copy of the PetscQuadrature object

 61:   Collective on PetscQuadrature

 63:   Input Parameter:
 64: . q  - The PetscQuadrature object

 66:   Output Parameter:
 67: . r  - The new PetscQuadrature object

 69:   Level: beginner

 71: .keywords: PetscQuadrature, quadrature, clone
 72: .seealso: PetscQuadratureCreate(), PetscQuadratureDestroy(), PetscQuadratureGetData()
 73: @*/
 74: PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature q, PetscQuadrature *r)
 75: {
 76:   PetscInt         order, dim, Nc, Nq;
 77:   const PetscReal *points, *weights;
 78:   PetscReal       *p, *w;
 79:   PetscErrorCode   ierr;

 83:   PetscQuadratureCreate(PetscObjectComm((PetscObject) q), r);
 84:   PetscQuadratureGetOrder(q, &order);
 85:   PetscQuadratureSetOrder(*r, order);
 86:   PetscQuadratureGetData(q, &dim, &Nc, &Nq, &points, &weights);
 87:   PetscMalloc1(Nq*dim, &p);
 88:   PetscMalloc1(Nq*Nc, &w);
 89:   PetscMemcpy(p, points, Nq*dim * sizeof(PetscReal));
 90:   PetscMemcpy(w, weights, Nc * Nq * sizeof(PetscReal));
 91:   PetscQuadratureSetData(*r, dim, Nc, Nq, p, w);
 92:   return(0);
 93: }

 95: /*@
 96:   PetscQuadratureDestroy - Destroys a PetscQuadrature object

 98:   Collective on PetscQuadrature

100:   Input Parameter:
101: . q  - The PetscQuadrature object

103:   Level: beginner

105: .keywords: PetscQuadrature, quadrature, destroy
106: .seealso: PetscQuadratureCreate(), PetscQuadratureGetData()
107: @*/
108: PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *q)
109: {

113:   if (!*q) return(0);
115:   if (--((PetscObject)(*q))->refct > 0) {
116:     *q = NULL;
117:     return(0);
118:   }
119:   PetscFree((*q)->points);
120:   PetscFree((*q)->weights);
121:   PetscHeaderDestroy(q);
122:   return(0);
123: }

125: /*@
126:   PetscQuadratureGetOrder - Return the order of the method

128:   Not collective

130:   Input Parameter:
131: . q - The PetscQuadrature object

133:   Output Parameter:
134: . order - The order of the quadrature, i.e. the highest degree polynomial that is exactly integrated

136:   Level: intermediate

138: .seealso: PetscQuadratureSetOrder(), PetscQuadratureGetData(), PetscQuadratureSetData()
139: @*/
140: PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature q, PetscInt *order)
141: {
145:   *order = q->order;
146:   return(0);
147: }

149: /*@
150:   PetscQuadratureSetOrder - Return the order of the method

152:   Not collective

154:   Input Parameters:
155: + q - The PetscQuadrature object
156: - order - The order of the quadrature, i.e. the highest degree polynomial that is exactly integrated

158:   Level: intermediate

160: .seealso: PetscQuadratureGetOrder(), PetscQuadratureGetData(), PetscQuadratureSetData()
161: @*/
162: PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature q, PetscInt order)
163: {
166:   q->order = order;
167:   return(0);
168: }

170: /*@
171:   PetscQuadratureGetNumComponents - Return the number of components for functions to be integrated

173:   Not collective

175:   Input Parameter:
176: . q - The PetscQuadrature object

178:   Output Parameter:
179: . Nc - The number of components

181:   Note: We are performing an integral int f(x) . w(x) dx, where both f and w (the weight) have Nc components.

183:   Level: intermediate

185: .seealso: PetscQuadratureSetNumComponents(), PetscQuadratureGetData(), PetscQuadratureSetData()
186: @*/
187: PetscErrorCode PetscQuadratureGetNumComponents(PetscQuadrature q, PetscInt *Nc)
188: {
192:   *Nc = q->Nc;
193:   return(0);
194: }

196: /*@
197:   PetscQuadratureSetNumComponents - Return the number of components for functions to be integrated

199:   Not collective

201:   Input Parameters:
202: + q  - The PetscQuadrature object
203: - Nc - The number of components

205:   Note: We are performing an integral int f(x) . w(x) dx, where both f and w (the weight) have Nc components.

207:   Level: intermediate

209: .seealso: PetscQuadratureGetNumComponents(), PetscQuadratureGetData(), PetscQuadratureSetData()
210: @*/
211: PetscErrorCode PetscQuadratureSetNumComponents(PetscQuadrature q, PetscInt Nc)
212: {
215:   q->Nc = Nc;
216:   return(0);
217: }

219: /*@C
220:   PetscQuadratureGetData - Returns the data defining the quadrature

222:   Not collective

224:   Input Parameter:
225: . q  - The PetscQuadrature object

227:   Output Parameters:
228: + dim - The spatial dimension
229: , Nc - The number of components
230: . npoints - The number of quadrature points
231: . points - The coordinates of each quadrature point
232: - weights - The weight of each quadrature point

234:   Level: intermediate

236:   Fortran Notes: From Fortran you must call PetscQuadratureRestoreData() when you are done with the data

238: .keywords: PetscQuadrature, quadrature
239: .seealso: PetscQuadratureCreate(), PetscQuadratureSetData()
240: @*/
241: PetscErrorCode PetscQuadratureGetData(PetscQuadrature q, PetscInt *dim, PetscInt *Nc, PetscInt *npoints, const PetscReal *points[], const PetscReal *weights[])
242: {
245:   if (dim) {
247:     *dim = q->dim;
248:   }
249:   if (Nc) {
251:     *Nc = q->Nc;
252:   }
253:   if (npoints) {
255:     *npoints = q->numPoints;
256:   }
257:   if (points) {
259:     *points = q->points;
260:   }
261:   if (weights) {
263:     *weights = q->weights;
264:   }
265:   return(0);
266: }

268: /*@C
269:   PetscQuadratureSetData - Sets the data defining the quadrature

271:   Not collective

273:   Input Parameters:
274: + q  - The PetscQuadrature object
275: . dim - The spatial dimension
276: , Nc - The number of components
277: . npoints - The number of quadrature points
278: . points - The coordinates of each quadrature point
279: - weights - The weight of each quadrature point

281:   Note: This routine owns the references to points and weights, so they msut be allocated using PetscMalloc() and the user should not free them.

283:   Level: intermediate

285: .keywords: PetscQuadrature, quadrature
286: .seealso: PetscQuadratureCreate(), PetscQuadratureGetData()
287: @*/
288: PetscErrorCode PetscQuadratureSetData(PetscQuadrature q, PetscInt dim, PetscInt Nc, PetscInt npoints, const PetscReal points[], const PetscReal weights[])
289: {
292:   if (dim >= 0)     q->dim       = dim;
293:   if (Nc >= 0)      q->Nc        = Nc;
294:   if (npoints >= 0) q->numPoints = npoints;
295:   if (points) {
297:     q->points = points;
298:   }
299:   if (weights) {
301:     q->weights = weights;
302:   }
303:   return(0);
304: }

306: /*@C
307:   PetscQuadratureView - Views a PetscQuadrature object

309:   Collective on PetscQuadrature

311:   Input Parameters:
312: + q  - The PetscQuadrature object
313: - viewer - The PetscViewer object

315:   Level: beginner

317: .keywords: PetscQuadrature, quadrature, view
318: .seealso: PetscQuadratureCreate(), PetscQuadratureGetData()
319: @*/
320: PetscErrorCode PetscQuadratureView(PetscQuadrature quad, PetscViewer viewer)
321: {
322:   PetscInt       q, d, c;

326:   PetscObjectPrintClassNamePrefixType((PetscObject)quad,viewer);
327:   if (quad->Nc > 1) {PetscViewerASCIIPrintf(viewer, "Quadrature on %D points with %D components\n  (", quad->numPoints, quad->Nc);}
328:   else              {PetscViewerASCIIPrintf(viewer, "Quadrature on %D points\n  (", quad->numPoints);}
329:   for (q = 0; q < quad->numPoints; ++q) {
330:     for (d = 0; d < quad->dim; ++d) {
331:       if (d) PetscViewerASCIIPrintf(viewer, ", ");
332:       PetscViewerASCIIPrintf(viewer, "%g\n", (double)quad->points[q*quad->dim+d]);
333:     }
334:     if (quad->Nc > 1) {
335:       PetscViewerASCIIPrintf(viewer, ") (");
336:       for (c = 0; c < quad->Nc; ++c) {
337:         if (c) PetscViewerASCIIPrintf(viewer, ", ");
338:         PetscViewerASCIIPrintf(viewer, "%g", (double)quad->weights[q*quad->Nc+c]);
339:       }
340:       PetscViewerASCIIPrintf(viewer, ")\n");
341:     } else {
342:       PetscViewerASCIIPrintf(viewer, ") %g\n", (double)quad->weights[q]);
343:     }
344:   }
345:   return(0);
346: }

348: /*@C
349:   PetscQuadratureExpandComposite - Return a quadrature over the composite element, which has the original quadrature in each subelement

351:   Not collective

353:   Input Parameter:
354: + q - The original PetscQuadrature
355: . numSubelements - The number of subelements the original element is divided into
356: . v0 - An array of the initial points for each subelement
357: - jac - An array of the Jacobian mappings from the reference to each subelement

359:   Output Parameters:
360: . dim - The dimension

362:   Note: Together v0 and jac define an affine mapping from the original reference element to each subelement

364:  Not available from Fortran

366:   Level: intermediate

368: .seealso: PetscFECreate(), PetscSpaceGetDimension(), PetscDualSpaceGetDimension()
369: @*/
370: PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature q, PetscInt numSubelements, const PetscReal v0[], const PetscReal jac[], PetscQuadrature *qref)
371: {
372:   const PetscReal *points,    *weights;
373:   PetscReal       *pointsRef, *weightsRef;
374:   PetscInt         dim, Nc, order, npoints, npointsRef, c, p, cp, d, e;
375:   PetscErrorCode   ierr;

382:   PetscQuadratureCreate(PETSC_COMM_SELF, qref);
383:   PetscQuadratureGetOrder(q, &order);
384:   PetscQuadratureGetData(q, &dim, &Nc, &npoints, &points, &weights);
385:   npointsRef = npoints*numSubelements;
386:   PetscMalloc1(npointsRef*dim,&pointsRef);
387:   PetscMalloc1(npointsRef*Nc, &weightsRef);
388:   for (c = 0; c < numSubelements; ++c) {
389:     for (p = 0; p < npoints; ++p) {
390:       for (d = 0; d < dim; ++d) {
391:         pointsRef[(c*npoints + p)*dim+d] = v0[c*dim+d];
392:         for (e = 0; e < dim; ++e) {
393:           pointsRef[(c*npoints + p)*dim+d] += jac[(c*dim + d)*dim+e]*(points[p*dim+e] + 1.0);
394:         }
395:       }
396:       /* Could also use detJ here */
397:       for (cp = 0; cp < Nc; ++cp) weightsRef[(c*npoints+p)*Nc+cp] = weights[p*Nc+cp]/numSubelements;
398:     }
399:   }
400:   PetscQuadratureSetOrder(*qref, order);
401:   PetscQuadratureSetData(*qref, dim, Nc, npointsRef, pointsRef, weightsRef);
402:   return(0);
403: }

405: /*@
406:    PetscDTLegendreEval - evaluate Legendre polynomial at points

408:    Not Collective

410:    Input Arguments:
411: +  npoints - number of spatial points to evaluate at
412: .  points - array of locations to evaluate at
413: .  ndegree - number of basis degrees to evaluate
414: -  degrees - sorted array of degrees to evaluate

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

421:    Level: intermediate

423: .seealso: PetscDTGaussQuadrature()
424: @*/
425: PetscErrorCode PetscDTLegendreEval(PetscInt npoints,const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2)
426: {
427:   PetscInt i,maxdegree;

430:   if (!npoints || !ndegree) return(0);
431:   maxdegree = degrees[ndegree-1];
432:   for (i=0; i<npoints; i++) {
433:     PetscReal pm1,pm2,pd1,pd2,pdd1,pdd2,x;
434:     PetscInt  j,k;
435:     x    = points[i];
436:     pm2  = 0;
437:     pm1  = 1;
438:     pd2  = 0;
439:     pd1  = 0;
440:     pdd2 = 0;
441:     pdd1 = 0;
442:     k    = 0;
443:     if (degrees[k] == 0) {
444:       if (B) B[i*ndegree+k] = pm1;
445:       if (D) D[i*ndegree+k] = pd1;
446:       if (D2) D2[i*ndegree+k] = pdd1;
447:       k++;
448:     }
449:     for (j=1; j<=maxdegree; j++,k++) {
450:       PetscReal p,d,dd;
451:       p    = ((2*j-1)*x*pm1 - (j-1)*pm2)/j;
452:       d    = pd2 + (2*j-1)*pm1;
453:       dd   = pdd2 + (2*j-1)*pd1;
454:       pm2  = pm1;
455:       pm1  = p;
456:       pd2  = pd1;
457:       pd1  = d;
458:       pdd2 = pdd1;
459:       pdd1 = dd;
460:       if (degrees[k] == j) {
461:         if (B) B[i*ndegree+k] = p;
462:         if (D) D[i*ndegree+k] = d;
463:         if (D2) D2[i*ndegree+k] = dd;
464:       }
465:     }
466:   }
467:   return(0);
468: }

470: /*@
471:    PetscDTGaussQuadrature - create Gauss quadrature

473:    Not Collective

475:    Input Arguments:
476: +  npoints - number of points
477: .  a - left end of interval (often-1)
478: -  b - right end of interval (often +1)

480:    Output Arguments:
481: +  x - quadrature points
482: -  w - quadrature weights

484:    Level: intermediate

486:    References:
487: .   1. - Golub and Welsch, Calculation of Quadrature Rules, Math. Comp. 23(106), 1969.

489: .seealso: PetscDTLegendreEval()
490: @*/
491: PetscErrorCode PetscDTGaussQuadrature(PetscInt npoints,PetscReal a,PetscReal b,PetscReal *x,PetscReal *w)
492: {
494:   PetscInt       i;
495:   PetscReal      *work;
496:   PetscScalar    *Z;
497:   PetscBLASInt   N,LDZ,info;

500:   PetscCitationsRegister(GaussCitation, &GaussCite);
501:   /* Set up the Golub-Welsch system */
502:   for (i=0; i<npoints; i++) {
503:     x[i] = 0;                   /* diagonal is 0 */
504:     if (i) w[i-1] = 0.5 / PetscSqrtReal(1 - 1./PetscSqr(2*i));
505:   }
506:   PetscMalloc2(npoints*npoints,&Z,PetscMax(1,2*npoints-2),&work);
507:   PetscBLASIntCast(npoints,&N);
508:   LDZ  = N;
509:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);
510:   PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("I",&N,x,w,Z,&LDZ,work,&info));
511:   PetscFPTrapPop();
512:   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEQR error");

514:   for (i=0; i<(npoints+1)/2; i++) {
515:     PetscReal y = 0.5 * (-x[i] + x[npoints-i-1]); /* enforces symmetry */
516:     x[i]           = (a+b)/2 - y*(b-a)/2;
517:     if (x[i] == -0.0) x[i] = 0.0;
518:     x[npoints-i-1] = (a+b)/2 + y*(b-a)/2;

520:     w[i] = w[npoints-1-i] = 0.5*(b-a)*(PetscSqr(PetscAbsScalar(Z[i*npoints])) + PetscSqr(PetscAbsScalar(Z[(npoints-i-1)*npoints])));
521:   }
522:   PetscFree2(Z,work);
523:   return(0);
524: }

526: /*@
527:   PetscDTGaussTensorQuadrature - creates a tensor-product Gauss quadrature

529:   Not Collective

531:   Input Arguments:
532: + dim     - The spatial dimension
533: . Nc      - The number of components
534: . npoints - number of points in one dimension
535: . a       - left end of interval (often-1)
536: - b       - right end of interval (often +1)

538:   Output Argument:
539: . q - A PetscQuadrature object

541:   Level: intermediate

543: .seealso: PetscDTGaussQuadrature(), PetscDTLegendreEval()
544: @*/
545: PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt dim, PetscInt Nc, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q)
546: {
547:   PetscInt       totpoints = dim > 1 ? dim > 2 ? npoints*PetscSqr(npoints) : PetscSqr(npoints) : npoints, i, j, k, c;
548:   PetscReal     *x, *w, *xw, *ww;

552:   PetscMalloc1(totpoints*dim,&x);
553:   PetscMalloc1(totpoints*Nc,&w);
554:   /* Set up the Golub-Welsch system */
555:   switch (dim) {
556:   case 0:
557:     PetscFree(x);
558:     PetscFree(w);
559:     PetscMalloc1(1, &x);
560:     PetscMalloc1(Nc, &w);
561:     x[0] = 0.0;
562:     for (c = 0; c < Nc; ++c) w[c] = 1.0;
563:     break;
564:   case 1:
565:     PetscMalloc1(npoints,&ww);
566:     PetscDTGaussQuadrature(npoints, a, b, x, ww);
567:     for (i = 0; i < npoints; ++i) for (c = 0; c < Nc; ++c) w[i*Nc+c] = ww[i];
568:     PetscFree(ww);
569:     break;
570:   case 2:
571:     PetscMalloc2(npoints,&xw,npoints,&ww);
572:     PetscDTGaussQuadrature(npoints, a, b, xw, ww);
573:     for (i = 0; i < npoints; ++i) {
574:       for (j = 0; j < npoints; ++j) {
575:         x[(i*npoints+j)*dim+0] = xw[i];
576:         x[(i*npoints+j)*dim+1] = xw[j];
577:         for (c = 0; c < Nc; ++c) w[(i*npoints+j)*Nc+c] = ww[i] * ww[j];
578:       }
579:     }
580:     PetscFree2(xw,ww);
581:     break;
582:   case 3:
583:     PetscMalloc2(npoints,&xw,npoints,&ww);
584:     PetscDTGaussQuadrature(npoints, a, b, xw, ww);
585:     for (i = 0; i < npoints; ++i) {
586:       for (j = 0; j < npoints; ++j) {
587:         for (k = 0; k < npoints; ++k) {
588:           x[((i*npoints+j)*npoints+k)*dim+0] = xw[i];
589:           x[((i*npoints+j)*npoints+k)*dim+1] = xw[j];
590:           x[((i*npoints+j)*npoints+k)*dim+2] = xw[k];
591:           for (c = 0; c < Nc; ++c) w[((i*npoints+j)*npoints+k)*Nc+c] = ww[i] * ww[j] * ww[k];
592:         }
593:       }
594:     }
595:     PetscFree2(xw,ww);
596:     break;
597:   default:
598:     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim);
599:   }
600:   PetscQuadratureCreate(PETSC_COMM_SELF, q);
601:   PetscQuadratureSetOrder(*q, npoints-1);
602:   PetscQuadratureSetData(*q, dim, Nc, totpoints, x, w);
603:   return(0);
604: }

606: /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x.
607:    Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */
608: PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorial_Internal(PetscInt n, PetscReal *factorial)
609: {
610:   PetscReal f = 1.0;
611:   PetscInt  i;

614:   for (i = 1; i < n+1; ++i) f *= i;
615:   *factorial = f;
616:   return(0);
617: }

619: /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x.
620:    Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */
621: PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobi(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P)
622: {
623:   PetscReal apb, pn1, pn2;
624:   PetscInt  k;

627:   if (!n) {*P = 1.0; return(0);}
628:   if (n == 1) {*P = 0.5 * (a - b + (a + b + 2.0) * x); return(0);}
629:   apb = a + b;
630:   pn2 = 1.0;
631:   pn1 = 0.5 * (a - b + (apb + 2.0) * x);
632:   *P  = 0.0;
633:   for (k = 2; k < n+1; ++k) {
634:     PetscReal a1 = 2.0 * k * (k + apb) * (2.0*k + apb - 2.0);
635:     PetscReal a2 = (2.0 * k + apb - 1.0) * (a*a - b*b);
636:     PetscReal a3 = (2.0 * k + apb - 2.0) * (2.0 * k + apb - 1.0) * (2.0 * k + apb);
637:     PetscReal a4 = 2.0 * (k + a - 1.0) * (k + b - 1.0) * (2.0 * k + apb);

639:     a2  = a2 / a1;
640:     a3  = a3 / a1;
641:     a4  = a4 / a1;
642:     *P  = (a2 + a3 * x) * pn1 - a4 * pn2;
643:     pn2 = pn1;
644:     pn1 = *P;
645:   }
646:   return(0);
647: }

649: /* Evaluates the first derivative of P_{n}^{a,b} at a point x. */
650: PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobiDerivative(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P)
651: {
652:   PetscReal      nP;

656:   if (!n) {*P = 0.0; return(0);}
657:   PetscDTComputeJacobi(a+1, b+1, n-1, x, &nP);
658:   *P   = 0.5 * (a + b + n + 1) * nP;
659:   return(0);
660: }

662: /* Maps from [-1,1]^2 to the (-1,1) reference triangle */
663: PETSC_STATIC_INLINE PetscErrorCode PetscDTMapSquareToTriangle_Internal(PetscReal x, PetscReal y, PetscReal *xi, PetscReal *eta)
664: {
666:   *xi  = 0.5 * (1.0 + x) * (1.0 - y) - 1.0;
667:   *eta = y;
668:   return(0);
669: }

671: /* Maps from [-1,1]^2 to the (-1,1) reference triangle */
672: PETSC_STATIC_INLINE PetscErrorCode PetscDTMapCubeToTetrahedron_Internal(PetscReal x, PetscReal y, PetscReal z, PetscReal *xi, PetscReal *eta, PetscReal *zeta)
673: {
675:   *xi   = 0.25 * (1.0 + x) * (1.0 - y) * (1.0 - z) - 1.0;
676:   *eta  = 0.5  * (1.0 + y) * (1.0 - z) - 1.0;
677:   *zeta = z;
678:   return(0);
679: }

681: static PetscErrorCode PetscDTGaussJacobiQuadrature1D_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal *x, PetscReal *w)
682: {
683:   PetscInt       maxIter = 100;
684:   PetscReal      eps     = 1.0e-8;
685:   PetscReal      a1, a2, a3, a4, a5, a6;
686:   PetscInt       k;


691:   a1      = PetscPowReal(2.0, a+b+1);
692: #if defined(PETSC_HAVE_TGAMMA)
693:   a2      = PetscTGamma(a + npoints + 1);
694:   a3      = PetscTGamma(b + npoints + 1);
695:   a4      = PetscTGamma(a + b + npoints + 1);
696: #else
697:   {
698:     PetscInt ia, ib;

700:     ia = (PetscInt) a;
701:     ib = (PetscInt) b;
702:     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 */
703:       PetscDTFactorial_Internal(ia + npoints, &a2);
704:       PetscDTFactorial_Internal(ib + npoints, &a3);
705:       PetscDTFactorial_Internal(ia + ib + npoints, &a4);
706:     } else {
707:       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"tgamma() - math routine is unavailable.");
708:     }
709:   }
710: #endif

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

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

725:       for (i = 0; i < k; ++i) s = s + 1.0 / (r - x[i]);
726:       PetscDTComputeJacobi(a, b, npoints, r, &f);
727:       PetscDTComputeJacobiDerivative(a, b, npoints, r, &fp);
728:       delta = f / (fp - f * s);
729:       r     = r - delta;
730:       if (PetscAbsReal(delta) < eps) break;
731:     }
732:     x[k] = r;
733:     PetscDTComputeJacobiDerivative(a, b, npoints, x[k], &dP);
734:     w[k] = a6 / (1.0 - PetscSqr(x[k])) / PetscSqr(dP);
735:   }
736:   return(0);
737: }

739: /*@
740:   PetscDTGaussJacobiQuadrature - create Gauss-Jacobi quadrature for a simplex

742:   Not Collective

744:   Input Arguments:
745: + dim     - The simplex dimension
746: . Nc      - The number of components
747: . npoints - The number of points in one dimension
748: . a       - left end of interval (often-1)
749: - b       - right end of interval (often +1)

751:   Output Argument:
752: . q - A PetscQuadrature object

754:   Level: intermediate

756:   References:
757: .  1. - Karniadakis and Sherwin.  FIAT

759: .seealso: PetscDTGaussTensorQuadrature(), PetscDTGaussQuadrature()
760: @*/
761: PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt dim, PetscInt Nc, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q)
762: {
763:   PetscInt       totpoints = dim > 1 ? dim > 2 ? npoints*PetscSqr(npoints) : PetscSqr(npoints) : npoints;
764:   PetscReal     *px, *wx, *py, *wy, *pz, *wz, *x, *w;
765:   PetscInt       i, j, k, c;

769:   if ((a != -1.0) || (b != 1.0)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must use default internal right now");
770:   PetscMalloc1(totpoints*dim, &x);
771:   PetscMalloc1(totpoints*Nc, &w);
772:   switch (dim) {
773:   case 0:
774:     PetscFree(x);
775:     PetscFree(w);
776:     PetscMalloc1(1, &x);
777:     PetscMalloc1(Nc, &w);
778:     x[0] = 0.0;
779:     for (c = 0; c < Nc; ++c) w[c] = 1.0;
780:     break;
781:   case 1:
782:     PetscMalloc1(npoints,&wx);
783:     PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, x, wx);
784:     for (i = 0; i < npoints; ++i) for (c = 0; c < Nc; ++c) w[i*Nc+c] = wx[i];
785:     PetscFree(wx);
786:     break;
787:   case 2:
788:     PetscMalloc4(npoints,&px,npoints,&wx,npoints,&py,npoints,&wy);
789:     PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, px, wx);
790:     PetscDTGaussJacobiQuadrature1D_Internal(npoints, 1.0, 0.0, py, wy);
791:     for (i = 0; i < npoints; ++i) {
792:       for (j = 0; j < npoints; ++j) {
793:         PetscDTMapSquareToTriangle_Internal(px[i], py[j], &x[(i*npoints+j)*2+0], &x[(i*npoints+j)*2+1]);
794:         for (c = 0; c < Nc; ++c) w[(i*npoints+j)*Nc+c] = 0.5 * wx[i] * wy[j];
795:       }
796:     }
797:     PetscFree4(px,wx,py,wy);
798:     break;
799:   case 3:
800:     PetscMalloc6(npoints,&px,npoints,&wx,npoints,&py,npoints,&wy,npoints,&pz,npoints,&wz);
801:     PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, px, wx);
802:     PetscDTGaussJacobiQuadrature1D_Internal(npoints, 1.0, 0.0, py, wy);
803:     PetscDTGaussJacobiQuadrature1D_Internal(npoints, 2.0, 0.0, pz, wz);
804:     for (i = 0; i < npoints; ++i) {
805:       for (j = 0; j < npoints; ++j) {
806:         for (k = 0; k < npoints; ++k) {
807:           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]);
808:           for (c = 0; c < Nc; ++c) w[((i*npoints+j)*npoints+k)*Nc+c] = 0.125 * wx[i] * wy[j] * wz[k];
809:         }
810:       }
811:     }
812:     PetscFree6(px,wx,py,wy,pz,wz);
813:     break;
814:   default:
815:     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim);
816:   }
817:   PetscQuadratureCreate(PETSC_COMM_SELF, q);
818:   PetscQuadratureSetOrder(*q, npoints-1);
819:   PetscQuadratureSetData(*q, dim, Nc, totpoints, x, w);
820:   return(0);
821: }

823: /*@
824:   PetscDTTanhSinhTensorQuadrature - create tanh-sinh quadrature for a tensor product cell

826:   Not Collective

828:   Input Arguments:
829: + dim   - The cell dimension
830: . level - The number of points in one dimension, 2^l
831: . a     - left end of interval (often-1)
832: - b     - right end of interval (often +1)

834:   Output Argument:
835: . q - A PetscQuadrature object

837:   Level: intermediate

839: .seealso: PetscDTGaussTensorQuadrature()
840: @*/
841: PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt dim, PetscInt level, PetscReal a, PetscReal b, PetscQuadrature *q)
842: {
843:   const PetscInt  p     = 16;                        /* Digits of precision in the evaluation */
844:   const PetscReal alpha = (b-a)/2.;                  /* Half-width of the integration interval */
845:   const PetscReal beta  = (b+a)/2.;                  /* Center of the integration interval */
846:   const PetscReal h     = PetscPowReal(2.0, -level); /* Step size, length between x_k */
847:   PetscReal       xk;                                /* Quadrature point x_k on reference domain [-1, 1] */
848:   PetscReal       wk    = 0.5*PETSC_PI;              /* Quadrature weight at x_k */
849:   PetscReal      *x, *w;
850:   PetscInt        K, k, npoints;
851:   PetscErrorCode  ierr;

854:   if (dim > 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Dimension %d not yet implemented", dim);
855:   if (!level) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a number of significant digits");
856:   /* Find K such that the weights are < 32 digits of precision */
857:   for (K = 1; PetscAbsReal(PetscLog10Real(wk)) < 2*p; ++K) {
858:     wk = 0.5*h*PETSC_PI*PetscCoshReal(K*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(K*h)));
859:   }
860:   PetscQuadratureCreate(PETSC_COMM_SELF, q);
861:   PetscQuadratureSetOrder(*q, 2*K+1);
862:   npoints = 2*K-1;
863:   PetscMalloc1(npoints*dim, &x);
864:   PetscMalloc1(npoints, &w);
865:   /* Center term */
866:   x[0] = beta;
867:   w[0] = 0.5*alpha*PETSC_PI;
868:   for (k = 1; k < K; ++k) {
869:     wk = 0.5*alpha*h*PETSC_PI*PetscCoshReal(k*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h)));
870:     xk = PetscTanhReal(0.5*PETSC_PI*PetscSinhReal(k*h));
871:     x[2*k-1] = -alpha*xk+beta;
872:     w[2*k-1] = wk;
873:     x[2*k+0] =  alpha*xk+beta;
874:     w[2*k+0] = wk;
875:   }
876:   PetscQuadratureSetData(*q, dim, 1, npoints, x, w);
877:   return(0);
878: }

880: PetscErrorCode PetscDTTanhSinhIntegrate(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol)
881: {
882:   const PetscInt  p     = 16;        /* Digits of precision in the evaluation */
883:   const PetscReal alpha = (b-a)/2.;  /* Half-width of the integration interval */
884:   const PetscReal beta  = (b+a)/2.;  /* Center of the integration interval */
885:   PetscReal       h     = 1.0;       /* Step size, length between x_k */
886:   PetscInt        l     = 0;         /* Level of refinement, h = 2^{-l} */
887:   PetscReal       osum  = 0.0;       /* Integral on last level */
888:   PetscReal       psum  = 0.0;       /* Integral on the level before the last level */
889:   PetscReal       sum;               /* Integral on current level */
890:   PetscReal       yk;                /* Quadrature point 1 - x_k on reference domain [-1, 1] */
891:   PetscReal       lx, rx;            /* Quadrature points to the left and right of 0 on the real domain [a, b] */
892:   PetscReal       wk;                /* Quadrature weight at x_k */
893:   PetscReal       lval, rval;        /* Terms in the quadature sum to the left and right of 0 */
894:   PetscInt        d;                 /* Digits of precision in the integral */

897:   if (digits <= 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a positive number of significant digits");
898:   /* Center term */
899:   func(beta, &lval);
900:   sum = 0.5*alpha*PETSC_PI*lval;
901:   /* */
902:   do {
903:     PetscReal lterm, rterm, maxTerm = 0.0, d1, d2, d3, d4;
904:     PetscInt  k = 1;

906:     ++l;
907:     /* PetscPrintf(PETSC_COMM_SELF, "LEVEL %D sum: %15.15f\n", l, sum); */
908:     /* At each level of refinement, h --> h/2 and sum --> sum/2 */
909:     psum = osum;
910:     osum = sum;
911:     h   *= 0.5;
912:     sum *= 0.5;
913:     do {
914:       wk = 0.5*h*PETSC_PI*PetscCoshReal(k*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h)));
915:       yk = 1.0/(PetscExpReal(0.5*PETSC_PI*PetscSinhReal(k*h)) * PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h)));
916:       lx = -alpha*(1.0 - yk)+beta;
917:       rx =  alpha*(1.0 - yk)+beta;
918:       func(lx, &lval);
919:       func(rx, &rval);
920:       lterm   = alpha*wk*lval;
921:       maxTerm = PetscMax(PetscAbsReal(lterm), maxTerm);
922:       sum    += lterm;
923:       rterm   = alpha*wk*rval;
924:       maxTerm = PetscMax(PetscAbsReal(rterm), maxTerm);
925:       sum    += rterm;
926:       ++k;
927:       /* Only need to evaluate every other point on refined levels */
928:       if (l != 1) ++k;
929:     } while (PetscAbsReal(PetscLog10Real(wk)) < p); /* Only need to evaluate sum until weights are < 32 digits of precision */

931:     d1 = PetscLog10Real(PetscAbsReal(sum - osum));
932:     d2 = PetscLog10Real(PetscAbsReal(sum - psum));
933:     d3 = PetscLog10Real(maxTerm) - p;
934:     if (PetscMax(PetscAbsReal(lterm), PetscAbsReal(rterm)) == 0.0) d4 = 0.0;
935:     else d4 = PetscLog10Real(PetscMax(PetscAbsReal(lterm), PetscAbsReal(rterm)));
936:     d  = PetscAbsInt(PetscMin(0, PetscMax(PetscMax(PetscMax(PetscSqr(d1)/d2, 2*d1), d3), d4)));
937:   } while (d < digits && l < 12);
938:   *sol = sum;

940:   return(0);
941: }

943: #ifdef PETSC_HAVE_MPFR
944: PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol)
945: {
946:   const PetscInt  safetyFactor = 2;  /* Calculate abcissa until 2*p digits */
947:   PetscInt        l            = 0;  /* Level of refinement, h = 2^{-l} */
948:   mpfr_t          alpha;             /* Half-width of the integration interval */
949:   mpfr_t          beta;              /* Center of the integration interval */
950:   mpfr_t          h;                 /* Step size, length between x_k */
951:   mpfr_t          osum;              /* Integral on last level */
952:   mpfr_t          psum;              /* Integral on the level before the last level */
953:   mpfr_t          sum;               /* Integral on current level */
954:   mpfr_t          yk;                /* Quadrature point 1 - x_k on reference domain [-1, 1] */
955:   mpfr_t          lx, rx;            /* Quadrature points to the left and right of 0 on the real domain [a, b] */
956:   mpfr_t          wk;                /* Quadrature weight at x_k */
957:   PetscReal       lval, rval;        /* Terms in the quadature sum to the left and right of 0 */
958:   PetscInt        d;                 /* Digits of precision in the integral */
959:   mpfr_t          pi2, kh, msinh, mcosh, maxTerm, curTerm, tmp;

962:   if (digits <= 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a positive number of significant digits");
963:   /* Create high precision storage */
964:   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);
965:   /* Initialization */
966:   mpfr_set_d(alpha, 0.5*(b-a), MPFR_RNDN);
967:   mpfr_set_d(beta,  0.5*(b+a), MPFR_RNDN);
968:   mpfr_set_d(osum,  0.0,       MPFR_RNDN);
969:   mpfr_set_d(psum,  0.0,       MPFR_RNDN);
970:   mpfr_set_d(h,     1.0,       MPFR_RNDN);
971:   mpfr_const_pi(pi2, MPFR_RNDN);
972:   mpfr_mul_d(pi2, pi2, 0.5, MPFR_RNDN);
973:   /* Center term */
974:   func(0.5*(b+a), &lval);
975:   mpfr_set(sum, pi2, MPFR_RNDN);
976:   mpfr_mul(sum, sum, alpha, MPFR_RNDN);
977:   mpfr_mul_d(sum, sum, lval, MPFR_RNDN);
978:   /* */
979:   do {
980:     PetscReal d1, d2, d3, d4;
981:     PetscInt  k = 1;

983:     ++l;
984:     mpfr_set_d(maxTerm, 0.0, MPFR_RNDN);
985:     /* PetscPrintf(PETSC_COMM_SELF, "LEVEL %D sum: %15.15f\n", l, sum); */
986:     /* At each level of refinement, h --> h/2 and sum --> sum/2 */
987:     mpfr_set(psum, osum, MPFR_RNDN);
988:     mpfr_set(osum,  sum, MPFR_RNDN);
989:     mpfr_mul_d(h,   h,   0.5, MPFR_RNDN);
990:     mpfr_mul_d(sum, sum, 0.5, MPFR_RNDN);
991:     do {
992:       mpfr_set_si(kh, k, MPFR_RNDN);
993:       mpfr_mul(kh, kh, h, MPFR_RNDN);
994:       /* Weight */
995:       mpfr_set(wk, h, MPFR_RNDN);
996:       mpfr_sinh_cosh(msinh, mcosh, kh, MPFR_RNDN);
997:       mpfr_mul(msinh, msinh, pi2, MPFR_RNDN);
998:       mpfr_mul(mcosh, mcosh, pi2, MPFR_RNDN);
999:       mpfr_cosh(tmp, msinh, MPFR_RNDN);
1000:       mpfr_sqr(tmp, tmp, MPFR_RNDN);
1001:       mpfr_mul(wk, wk, mcosh, MPFR_RNDN);
1002:       mpfr_div(wk, wk, tmp, MPFR_RNDN);
1003:       /* Abscissa */
1004:       mpfr_set_d(yk, 1.0, MPFR_RNDZ);
1005:       mpfr_cosh(tmp, msinh, MPFR_RNDN);
1006:       mpfr_div(yk, yk, tmp, MPFR_RNDZ);
1007:       mpfr_exp(tmp, msinh, MPFR_RNDN);
1008:       mpfr_div(yk, yk, tmp, MPFR_RNDZ);
1009:       /* Quadrature points */
1010:       mpfr_sub_d(lx, yk, 1.0, MPFR_RNDZ);
1011:       mpfr_mul(lx, lx, alpha, MPFR_RNDU);
1012:       mpfr_add(lx, lx, beta, MPFR_RNDU);
1013:       mpfr_d_sub(rx, 1.0, yk, MPFR_RNDZ);
1014:       mpfr_mul(rx, rx, alpha, MPFR_RNDD);
1015:       mpfr_add(rx, rx, beta, MPFR_RNDD);
1016:       /* Evaluation */
1017:       func(mpfr_get_d(lx, MPFR_RNDU), &lval);
1018:       func(mpfr_get_d(rx, MPFR_RNDD), &rval);
1019:       /* Update */
1020:       mpfr_mul(tmp, wk, alpha, MPFR_RNDN);
1021:       mpfr_mul_d(tmp, tmp, lval, MPFR_RNDN);
1022:       mpfr_add(sum, sum, tmp, MPFR_RNDN);
1023:       mpfr_abs(tmp, tmp, MPFR_RNDN);
1024:       mpfr_max(maxTerm, maxTerm, tmp, MPFR_RNDN);
1025:       mpfr_set(curTerm, tmp, MPFR_RNDN);
1026:       mpfr_mul(tmp, wk, alpha, MPFR_RNDN);
1027:       mpfr_mul_d(tmp, tmp, rval, MPFR_RNDN);
1028:       mpfr_add(sum, sum, tmp, MPFR_RNDN);
1029:       mpfr_abs(tmp, tmp, MPFR_RNDN);
1030:       mpfr_max(maxTerm, maxTerm, tmp, MPFR_RNDN);
1031:       mpfr_max(curTerm, curTerm, tmp, MPFR_RNDN);
1032:       ++k;
1033:       /* Only need to evaluate every other point on refined levels */
1034:       if (l != 1) ++k;
1035:       mpfr_log10(tmp, wk, MPFR_RNDN);
1036:       mpfr_abs(tmp, tmp, MPFR_RNDN);
1037:     } while (mpfr_get_d(tmp, MPFR_RNDN) < safetyFactor*digits); /* Only need to evaluate sum until weights are < 32 digits of precision */
1038:     mpfr_sub(tmp, sum, osum, MPFR_RNDN);
1039:     mpfr_abs(tmp, tmp, MPFR_RNDN);
1040:     mpfr_log10(tmp, tmp, MPFR_RNDN);
1041:     d1 = mpfr_get_d(tmp, MPFR_RNDN);
1042:     mpfr_sub(tmp, sum, psum, MPFR_RNDN);
1043:     mpfr_abs(tmp, tmp, MPFR_RNDN);
1044:     mpfr_log10(tmp, tmp, MPFR_RNDN);
1045:     d2 = mpfr_get_d(tmp, MPFR_RNDN);
1046:     mpfr_log10(tmp, maxTerm, MPFR_RNDN);
1047:     d3 = mpfr_get_d(tmp, MPFR_RNDN) - digits;
1048:     mpfr_log10(tmp, curTerm, MPFR_RNDN);
1049:     d4 = mpfr_get_d(tmp, MPFR_RNDN);
1050:     d  = PetscAbsInt(PetscMin(0, PetscMax(PetscMax(PetscMax(PetscSqr(d1)/d2, 2*d1), d3), d4)));
1051:   } while (d < digits && l < 8);
1052:   *sol = mpfr_get_d(sum, MPFR_RNDN);
1053:   /* Cleanup */
1054:   mpfr_clears(alpha, beta, h, sum, osum, psum, yk, wk, lx, rx, tmp, maxTerm, curTerm, pi2, kh, msinh, mcosh, NULL);
1055:   return(0);
1056: }
1057: #else

1059: PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol)
1060: {
1061:   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "This method will not work without MPFR. Reconfigure using --download-mpfr --download-gmp");
1062: }
1063: #endif

1065: /* Overwrites A. Can only handle full-rank problems with m>=n
1066:  * A in column-major format
1067:  * Ainv in row-major format
1068:  * tau has length m
1069:  * worksize must be >= max(1,n)
1070:  */
1071: static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m,PetscInt mstride,PetscInt n,PetscReal *A_in,PetscReal *Ainv_out,PetscScalar *tau,PetscInt worksize,PetscScalar *work)
1072: {
1074:   PetscBLASInt   M,N,K,lda,ldb,ldwork,info;
1075:   PetscScalar    *A,*Ainv,*R,*Q,Alpha;

1078: #if defined(PETSC_USE_COMPLEX)
1079:   {
1080:     PetscInt i,j;
1081:     PetscMalloc2(m*n,&A,m*n,&Ainv);
1082:     for (j=0; j<n; j++) {
1083:       for (i=0; i<m; i++) A[i+m*j] = A_in[i+mstride*j];
1084:     }
1085:     mstride = m;
1086:   }
1087: #else
1088:   A = A_in;
1089:   Ainv = Ainv_out;
1090: #endif

1092:   PetscBLASIntCast(m,&M);
1093:   PetscBLASIntCast(n,&N);
1094:   PetscBLASIntCast(mstride,&lda);
1095:   PetscBLASIntCast(worksize,&ldwork);
1096:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);
1097:   PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info));
1098:   PetscFPTrapPop();
1099:   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error");
1100:   R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */

1102:   /* Extract an explicit representation of Q */
1103:   Q = Ainv;
1104:   PetscMemcpy(Q,A,mstride*n*sizeof(PetscScalar));
1105:   K = N;                        /* full rank */
1106:   PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info));
1107:   if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error");

1109:   /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */
1110:   Alpha = 1.0;
1111:   ldb = lda;
1112:   PetscStackCallBLAS("BLAStrsm",BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb));
1113:   /* Ainv is Q, overwritten with inverse */

1115: #if defined(PETSC_USE_COMPLEX)
1116:   {
1117:     PetscInt i;
1118:     for (i=0; i<m*n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]);
1119:     PetscFree2(A,Ainv);
1120:   }
1121: #endif
1122:   return(0);
1123: }

1125: /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */
1126: static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval,const PetscReal *x,PetscInt ndegree,const PetscInt *degrees,PetscBool Transpose,PetscReal *B)
1127: {
1129:   PetscReal      *Bv;
1130:   PetscInt       i,j;

1133:   PetscMalloc1((ninterval+1)*ndegree,&Bv);
1134:   /* Point evaluation of L_p on all the source vertices */
1135:   PetscDTLegendreEval(ninterval+1,x,ndegree,degrees,Bv,NULL,NULL);
1136:   /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */
1137:   for (i=0; i<ninterval; i++) {
1138:     for (j=0; j<ndegree; j++) {
1139:       if (Transpose) B[i+ninterval*j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j];
1140:       else           B[i*ndegree+j]   = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j];
1141:     }
1142:   }
1143:   PetscFree(Bv);
1144:   return(0);
1145: }

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

1150:    Not Collective

1152:    Input Arguments:
1153: +  degree - degree of reconstruction polynomial
1154: .  nsource - number of source intervals
1155: .  sourcex - sorted coordinates of source cell boundaries (length nsource+1)
1156: .  ntarget - number of target intervals
1157: -  targetx - sorted coordinates of target cell boundaries (length ntarget+1)

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

1162:    Level: advanced

1164: .seealso: PetscDTLegendreEval()
1165: @*/
1166: PetscErrorCode PetscDTReconstructPoly(PetscInt degree,PetscInt nsource,const PetscReal *sourcex,PetscInt ntarget,const PetscReal *targetx,PetscReal *R)
1167: {
1169:   PetscInt       i,j,k,*bdegrees,worksize;
1170:   PetscReal      xmin,xmax,center,hscale,*sourcey,*targety,*Bsource,*Bsinv,*Btarget;
1171:   PetscScalar    *tau,*work;

1177:   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);
1178: #if defined(PETSC_USE_DEBUG)
1179:   for (i=0; i<nsource; i++) {
1180:     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]);
1181:   }
1182:   for (i=0; i<ntarget; i++) {
1183:     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]);
1184:   }
1185: #endif
1186:   xmin = PetscMin(sourcex[0],targetx[0]);
1187:   xmax = PetscMax(sourcex[nsource],targetx[ntarget]);
1188:   center = (xmin + xmax)/2;
1189:   hscale = (xmax - xmin)/2;
1190:   worksize = nsource;
1191:   PetscMalloc4(degree+1,&bdegrees,nsource+1,&sourcey,nsource*(degree+1),&Bsource,worksize,&work);
1192:   PetscMalloc4(nsource,&tau,nsource*(degree+1),&Bsinv,ntarget+1,&targety,ntarget*(degree+1),&Btarget);
1193:   for (i=0; i<=nsource; i++) sourcey[i] = (sourcex[i]-center)/hscale;
1194:   for (i=0; i<=degree; i++) bdegrees[i] = i+1;
1195:   PetscDTLegendreIntegrate(nsource,sourcey,degree+1,bdegrees,PETSC_TRUE,Bsource);
1196:   PetscDTPseudoInverseQR(nsource,nsource,degree+1,Bsource,Bsinv,tau,nsource,work);
1197:   for (i=0; i<=ntarget; i++) targety[i] = (targetx[i]-center)/hscale;
1198:   PetscDTLegendreIntegrate(ntarget,targety,degree+1,bdegrees,PETSC_FALSE,Btarget);
1199:   for (i=0; i<ntarget; i++) {
1200:     PetscReal rowsum = 0;
1201:     for (j=0; j<nsource; j++) {
1202:       PetscReal sum = 0;
1203:       for (k=0; k<degree+1; k++) {
1204:         sum += Btarget[i*(degree+1)+k] * Bsinv[k*nsource+j];
1205:       }
1206:       R[i*nsource+j] = sum;
1207:       rowsum += sum;
1208:     }
1209:     for (j=0; j<nsource; j++) R[i*nsource+j] /= rowsum; /* normalize each row */
1210:   }
1211:   PetscFree4(bdegrees,sourcey,Bsource,work);
1212:   PetscFree4(tau,Bsinv,targety,Btarget);
1213:   return(0);
1214: }