Actual source code: dt.c
petsc-3.10.5 2019-03-28
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:
237: From Fortran you must call PetscQuadratureRestoreData() when you are done with the data
239: .keywords: PetscQuadrature, quadrature
240: .seealso: PetscQuadratureCreate(), PetscQuadratureSetData()
241: @*/
242: PetscErrorCode PetscQuadratureGetData(PetscQuadrature q, PetscInt *dim, PetscInt *Nc, PetscInt *npoints, const PetscReal *points[], const PetscReal *weights[])
243: {
246: if (dim) {
248: *dim = q->dim;
249: }
250: if (Nc) {
252: *Nc = q->Nc;
253: }
254: if (npoints) {
256: *npoints = q->numPoints;
257: }
258: if (points) {
260: *points = q->points;
261: }
262: if (weights) {
264: *weights = q->weights;
265: }
266: return(0);
267: }
269: /*@C
270: PetscQuadratureSetData - Sets the data defining the quadrature
272: Not collective
274: Input Parameters:
275: + q - The PetscQuadrature object
276: . dim - The spatial dimension
277: . Nc - The number of components
278: . npoints - The number of quadrature points
279: . points - The coordinates of each quadrature point
280: - weights - The weight of each quadrature point
282: Note: This routine owns the references to points and weights, so they must be allocated using PetscMalloc() and the user should not free them.
284: Level: intermediate
286: .keywords: PetscQuadrature, quadrature
287: .seealso: PetscQuadratureCreate(), PetscQuadratureGetData()
288: @*/
289: PetscErrorCode PetscQuadratureSetData(PetscQuadrature q, PetscInt dim, PetscInt Nc, PetscInt npoints, const PetscReal points[], const PetscReal weights[])
290: {
293: if (dim >= 0) q->dim = dim;
294: if (Nc >= 0) q->Nc = Nc;
295: if (npoints >= 0) q->numPoints = npoints;
296: if (points) {
298: q->points = points;
299: }
300: if (weights) {
302: q->weights = weights;
303: }
304: return(0);
305: }
307: /*@C
308: PetscQuadratureView - Views a PetscQuadrature object
310: Collective on PetscQuadrature
312: Input Parameters:
313: + q - The PetscQuadrature object
314: - viewer - The PetscViewer object
316: Level: beginner
318: .keywords: PetscQuadrature, quadrature, view
319: .seealso: PetscQuadratureCreate(), PetscQuadratureGetData()
320: @*/
321: PetscErrorCode PetscQuadratureView(PetscQuadrature quad, PetscViewer viewer)
322: {
323: PetscInt q, d, c;
327: PetscObjectPrintClassNamePrefixType((PetscObject)quad,viewer);
328: if (quad->Nc > 1) {PetscViewerASCIIPrintf(viewer, "Quadrature on %D points with %D components\n (", quad->numPoints, quad->Nc);}
329: else {PetscViewerASCIIPrintf(viewer, "Quadrature on %D points\n (", quad->numPoints);}
330: for (q = 0; q < quad->numPoints; ++q) {
331: for (d = 0; d < quad->dim; ++d) {
332: if (d) PetscViewerASCIIPrintf(viewer, ", ");
333: PetscViewerASCIIPrintf(viewer, "%g\n", (double)quad->points[q*quad->dim+d]);
334: }
335: if (quad->Nc > 1) {
336: PetscViewerASCIIPrintf(viewer, ") (");
337: for (c = 0; c < quad->Nc; ++c) {
338: if (c) PetscViewerASCIIPrintf(viewer, ", ");
339: PetscViewerASCIIPrintf(viewer, "%g", (double)quad->weights[q*quad->Nc+c]);
340: }
341: PetscViewerASCIIPrintf(viewer, ")\n");
342: } else {
343: PetscViewerASCIIPrintf(viewer, ") %g\n", (double)quad->weights[q]);
344: }
345: }
346: return(0);
347: }
349: /*@C
350: PetscQuadratureExpandComposite - Return a quadrature over the composite element, which has the original quadrature in each subelement
352: Not collective
354: Input Parameter:
355: + q - The original PetscQuadrature
356: . numSubelements - The number of subelements the original element is divided into
357: . v0 - An array of the initial points for each subelement
358: - jac - An array of the Jacobian mappings from the reference to each subelement
360: Output Parameters:
361: . dim - The dimension
363: Note: Together v0 and jac define an affine mapping from the original reference element to each subelement
365: Not available from Fortran
367: Level: intermediate
369: .seealso: PetscFECreate(), PetscSpaceGetDimension(), PetscDualSpaceGetDimension()
370: @*/
371: PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature q, PetscInt numSubelements, const PetscReal v0[], const PetscReal jac[], PetscQuadrature *qref)
372: {
373: const PetscReal *points, *weights;
374: PetscReal *pointsRef, *weightsRef;
375: PetscInt dim, Nc, order, npoints, npointsRef, c, p, cp, d, e;
376: PetscErrorCode ierr;
383: PetscQuadratureCreate(PETSC_COMM_SELF, qref);
384: PetscQuadratureGetOrder(q, &order);
385: PetscQuadratureGetData(q, &dim, &Nc, &npoints, &points, &weights);
386: npointsRef = npoints*numSubelements;
387: PetscMalloc1(npointsRef*dim,&pointsRef);
388: PetscMalloc1(npointsRef*Nc, &weightsRef);
389: for (c = 0; c < numSubelements; ++c) {
390: for (p = 0; p < npoints; ++p) {
391: for (d = 0; d < dim; ++d) {
392: pointsRef[(c*npoints + p)*dim+d] = v0[c*dim+d];
393: for (e = 0; e < dim; ++e) {
394: pointsRef[(c*npoints + p)*dim+d] += jac[(c*dim + d)*dim+e]*(points[p*dim+e] + 1.0);
395: }
396: }
397: /* Could also use detJ here */
398: for (cp = 0; cp < Nc; ++cp) weightsRef[(c*npoints+p)*Nc+cp] = weights[p*Nc+cp]/numSubelements;
399: }
400: }
401: PetscQuadratureSetOrder(*qref, order);
402: PetscQuadratureSetData(*qref, dim, Nc, npointsRef, pointsRef, weightsRef);
403: return(0);
404: }
406: /*@
407: PetscDTLegendreEval - evaluate Legendre polynomial at points
409: Not Collective
411: Input Arguments:
412: + npoints - number of spatial points to evaluate at
413: . points - array of locations to evaluate at
414: . ndegree - number of basis degrees to evaluate
415: - degrees - sorted array of degrees to evaluate
417: Output Arguments:
418: + B - row-oriented basis evaluation matrix B[point*ndegree + degree] (dimension npoints*ndegrees, allocated by caller) (or NULL)
419: . D - row-oriented derivative evaluation matrix (or NULL)
420: - D2 - row-oriented second derivative evaluation matrix (or NULL)
422: Level: intermediate
424: .seealso: PetscDTGaussQuadrature()
425: @*/
426: PetscErrorCode PetscDTLegendreEval(PetscInt npoints,const PetscReal *points,PetscInt ndegree,const PetscInt *degrees,PetscReal *B,PetscReal *D,PetscReal *D2)
427: {
428: PetscInt i,maxdegree;
431: if (!npoints || !ndegree) return(0);
432: maxdegree = degrees[ndegree-1];
433: for (i=0; i<npoints; i++) {
434: PetscReal pm1,pm2,pd1,pd2,pdd1,pdd2,x;
435: PetscInt j,k;
436: x = points[i];
437: pm2 = 0;
438: pm1 = 1;
439: pd2 = 0;
440: pd1 = 0;
441: pdd2 = 0;
442: pdd1 = 0;
443: k = 0;
444: if (degrees[k] == 0) {
445: if (B) B[i*ndegree+k] = pm1;
446: if (D) D[i*ndegree+k] = pd1;
447: if (D2) D2[i*ndegree+k] = pdd1;
448: k++;
449: }
450: for (j=1; j<=maxdegree; j++,k++) {
451: PetscReal p,d,dd;
452: p = ((2*j-1)*x*pm1 - (j-1)*pm2)/j;
453: d = pd2 + (2*j-1)*pm1;
454: dd = pdd2 + (2*j-1)*pd1;
455: pm2 = pm1;
456: pm1 = p;
457: pd2 = pd1;
458: pd1 = d;
459: pdd2 = pdd1;
460: pdd1 = dd;
461: if (degrees[k] == j) {
462: if (B) B[i*ndegree+k] = p;
463: if (D) D[i*ndegree+k] = d;
464: if (D2) D2[i*ndegree+k] = dd;
465: }
466: }
467: }
468: return(0);
469: }
471: /*@
472: PetscDTGaussQuadrature - create Gauss quadrature
474: Not Collective
476: Input Arguments:
477: + npoints - number of points
478: . a - left end of interval (often-1)
479: - b - right end of interval (often +1)
481: Output Arguments:
482: + x - quadrature points
483: - w - quadrature weights
485: Level: intermediate
487: References:
488: . 1. - Golub and Welsch, Calculation of Quadrature Rules, Math. Comp. 23(106), 1969.
490: .seealso: PetscDTLegendreEval()
491: @*/
492: PetscErrorCode PetscDTGaussQuadrature(PetscInt npoints,PetscReal a,PetscReal b,PetscReal *x,PetscReal *w)
493: {
495: PetscInt i;
496: PetscReal *work;
497: PetscScalar *Z;
498: PetscBLASInt N,LDZ,info;
501: PetscCitationsRegister(GaussCitation, &GaussCite);
502: /* Set up the Golub-Welsch system */
503: for (i=0; i<npoints; i++) {
504: x[i] = 0; /* diagonal is 0 */
505: if (i) w[i-1] = 0.5 / PetscSqrtReal(1 - 1./PetscSqr(2*i));
506: }
507: PetscMalloc2(npoints*npoints,&Z,PetscMax(1,2*npoints-2),&work);
508: PetscBLASIntCast(npoints,&N);
509: LDZ = N;
510: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
511: PetscStackCallBLAS("LAPACKsteqr",LAPACKsteqr_("I",&N,x,w,Z,&LDZ,work,&info));
512: PetscFPTrapPop();
513: if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"xSTEQR error");
515: for (i=0; i<(npoints+1)/2; i++) {
516: PetscReal y = 0.5 * (-x[i] + x[npoints-i-1]); /* enforces symmetry */
517: x[i] = (a+b)/2 - y*(b-a)/2;
518: if (x[i] == -0.0) x[i] = 0.0;
519: x[npoints-i-1] = (a+b)/2 + y*(b-a)/2;
521: w[i] = w[npoints-1-i] = 0.5*(b-a)*(PetscSqr(PetscAbsScalar(Z[i*npoints])) + PetscSqr(PetscAbsScalar(Z[(npoints-i-1)*npoints])));
522: }
523: PetscFree2(Z,work);
524: return(0);
525: }
527: /*@
528: PetscDTGaussTensorQuadrature - creates a tensor-product Gauss quadrature
530: Not Collective
532: Input Arguments:
533: + dim - The spatial dimension
534: . Nc - The number of components
535: . npoints - number of points in one dimension
536: . a - left end of interval (often-1)
537: - b - right end of interval (often +1)
539: Output Argument:
540: . q - A PetscQuadrature object
542: Level: intermediate
544: .seealso: PetscDTGaussQuadrature(), PetscDTLegendreEval()
545: @*/
546: PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt dim, PetscInt Nc, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q)
547: {
548: PetscInt totpoints = dim > 1 ? dim > 2 ? npoints*PetscSqr(npoints) : PetscSqr(npoints) : npoints, i, j, k, c;
549: PetscReal *x, *w, *xw, *ww;
553: PetscMalloc1(totpoints*dim,&x);
554: PetscMalloc1(totpoints*Nc,&w);
555: /* Set up the Golub-Welsch system */
556: switch (dim) {
557: case 0:
558: PetscFree(x);
559: PetscFree(w);
560: PetscMalloc1(1, &x);
561: PetscMalloc1(Nc, &w);
562: x[0] = 0.0;
563: for (c = 0; c < Nc; ++c) w[c] = 1.0;
564: break;
565: case 1:
566: PetscMalloc1(npoints,&ww);
567: PetscDTGaussQuadrature(npoints, a, b, x, ww);
568: for (i = 0; i < npoints; ++i) for (c = 0; c < Nc; ++c) w[i*Nc+c] = ww[i];
569: PetscFree(ww);
570: break;
571: case 2:
572: PetscMalloc2(npoints,&xw,npoints,&ww);
573: PetscDTGaussQuadrature(npoints, a, b, xw, ww);
574: for (i = 0; i < npoints; ++i) {
575: for (j = 0; j < npoints; ++j) {
576: x[(i*npoints+j)*dim+0] = xw[i];
577: x[(i*npoints+j)*dim+1] = xw[j];
578: for (c = 0; c < Nc; ++c) w[(i*npoints+j)*Nc+c] = ww[i] * ww[j];
579: }
580: }
581: PetscFree2(xw,ww);
582: break;
583: case 3:
584: PetscMalloc2(npoints,&xw,npoints,&ww);
585: PetscDTGaussQuadrature(npoints, a, b, xw, ww);
586: for (i = 0; i < npoints; ++i) {
587: for (j = 0; j < npoints; ++j) {
588: for (k = 0; k < npoints; ++k) {
589: x[((i*npoints+j)*npoints+k)*dim+0] = xw[i];
590: x[((i*npoints+j)*npoints+k)*dim+1] = xw[j];
591: x[((i*npoints+j)*npoints+k)*dim+2] = xw[k];
592: for (c = 0; c < Nc; ++c) w[((i*npoints+j)*npoints+k)*Nc+c] = ww[i] * ww[j] * ww[k];
593: }
594: }
595: }
596: PetscFree2(xw,ww);
597: break;
598: default:
599: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim);
600: }
601: PetscQuadratureCreate(PETSC_COMM_SELF, q);
602: PetscQuadratureSetOrder(*q, 2*npoints-1);
603: PetscQuadratureSetData(*q, dim, Nc, totpoints, x, w);
604: return(0);
605: }
607: /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x.
608: Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */
609: PETSC_STATIC_INLINE PetscErrorCode PetscDTFactorial_Internal(PetscInt n, PetscReal *factorial)
610: {
611: PetscReal f = 1.0;
612: PetscInt i;
615: for (i = 1; i < n+1; ++i) f *= i;
616: *factorial = f;
617: return(0);
618: }
620: /* Evaluates the nth jacobi polynomial with weight parameters a,b at a point x.
621: Recurrence relations implemented from the pseudocode given in Karniadakis and Sherwin, Appendix B */
622: PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobi(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P)
623: {
624: PetscReal apb, pn1, pn2;
625: PetscInt k;
628: if (!n) {*P = 1.0; return(0);}
629: if (n == 1) {*P = 0.5 * (a - b + (a + b + 2.0) * x); return(0);}
630: apb = a + b;
631: pn2 = 1.0;
632: pn1 = 0.5 * (a - b + (apb + 2.0) * x);
633: *P = 0.0;
634: for (k = 2; k < n+1; ++k) {
635: PetscReal a1 = 2.0 * k * (k + apb) * (2.0*k + apb - 2.0);
636: PetscReal a2 = (2.0 * k + apb - 1.0) * (a*a - b*b);
637: PetscReal a3 = (2.0 * k + apb - 2.0) * (2.0 * k + apb - 1.0) * (2.0 * k + apb);
638: PetscReal a4 = 2.0 * (k + a - 1.0) * (k + b - 1.0) * (2.0 * k + apb);
640: a2 = a2 / a1;
641: a3 = a3 / a1;
642: a4 = a4 / a1;
643: *P = (a2 + a3 * x) * pn1 - a4 * pn2;
644: pn2 = pn1;
645: pn1 = *P;
646: }
647: return(0);
648: }
650: /* Evaluates the first derivative of P_{n}^{a,b} at a point x. */
651: PETSC_STATIC_INLINE PetscErrorCode PetscDTComputeJacobiDerivative(PetscReal a, PetscReal b, PetscInt n, PetscReal x, PetscReal *P)
652: {
653: PetscReal nP;
657: if (!n) {*P = 0.0; return(0);}
658: PetscDTComputeJacobi(a+1, b+1, n-1, x, &nP);
659: *P = 0.5 * (a + b + n + 1) * nP;
660: return(0);
661: }
663: /* Maps from [-1,1]^2 to the (-1,1) reference triangle */
664: PETSC_STATIC_INLINE PetscErrorCode PetscDTMapSquareToTriangle_Internal(PetscReal x, PetscReal y, PetscReal *xi, PetscReal *eta)
665: {
667: *xi = 0.5 * (1.0 + x) * (1.0 - y) - 1.0;
668: *eta = y;
669: return(0);
670: }
672: /* Maps from [-1,1]^2 to the (-1,1) reference triangle */
673: PETSC_STATIC_INLINE PetscErrorCode PetscDTMapCubeToTetrahedron_Internal(PetscReal x, PetscReal y, PetscReal z, PetscReal *xi, PetscReal *eta, PetscReal *zeta)
674: {
676: *xi = 0.25 * (1.0 + x) * (1.0 - y) * (1.0 - z) - 1.0;
677: *eta = 0.5 * (1.0 + y) * (1.0 - z) - 1.0;
678: *zeta = z;
679: return(0);
680: }
682: static PetscErrorCode PetscDTGaussJacobiQuadrature1D_Internal(PetscInt npoints, PetscReal a, PetscReal b, PetscReal *x, PetscReal *w)
683: {
684: PetscInt maxIter = 100;
685: PetscReal eps = 1.0e-8;
686: PetscReal a1, a2, a3, a4, a5, a6;
687: PetscInt k;
692: a1 = PetscPowReal(2.0, a+b+1);
693: #if defined(PETSC_HAVE_TGAMMA)
694: a2 = PetscTGamma(a + npoints + 1);
695: a3 = PetscTGamma(b + npoints + 1);
696: a4 = PetscTGamma(a + b + npoints + 1);
697: #else
698: {
699: PetscInt ia, ib;
701: ia = (PetscInt) a;
702: ib = (PetscInt) b;
703: 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 */
704: PetscDTFactorial_Internal(ia + npoints, &a2);
705: PetscDTFactorial_Internal(ib + npoints, &a3);
706: PetscDTFactorial_Internal(ia + ib + npoints, &a4);
707: } else {
708: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"tgamma() - math routine is unavailable.");
709: }
710: }
711: #endif
713: PetscDTFactorial_Internal(npoints, &a5);
714: a6 = a1 * a2 * a3 / a4 / a5;
715: /* Computes the m roots of P_{m}^{a,b} on [-1,1] by Newton's method with Chebyshev points as initial guesses.
716: Algorithm implemented from the pseudocode given by Karniadakis and Sherwin and Python in FIAT */
717: for (k = 0; k < npoints; ++k) {
718: PetscReal r = -PetscCosReal((2.0*k + 1.0) * PETSC_PI / (2.0 * npoints)), dP;
719: PetscInt j;
721: if (k > 0) r = 0.5 * (r + x[k-1]);
722: for (j = 0; j < maxIter; ++j) {
723: PetscReal s = 0.0, delta, f, fp;
724: PetscInt i;
726: for (i = 0; i < k; ++i) s = s + 1.0 / (r - x[i]);
727: PetscDTComputeJacobi(a, b, npoints, r, &f);
728: PetscDTComputeJacobiDerivative(a, b, npoints, r, &fp);
729: delta = f / (fp - f * s);
730: r = r - delta;
731: if (PetscAbsReal(delta) < eps) break;
732: }
733: x[k] = r;
734: PetscDTComputeJacobiDerivative(a, b, npoints, x[k], &dP);
735: w[k] = a6 / (1.0 - PetscSqr(x[k])) / PetscSqr(dP);
736: }
737: return(0);
738: }
740: /*@
741: PetscDTGaussJacobiQuadrature - create Gauss-Jacobi quadrature for a simplex
743: Not Collective
745: Input Arguments:
746: + dim - The simplex dimension
747: . Nc - The number of components
748: . npoints - The number of points in one dimension
749: . a - left end of interval (often-1)
750: - b - right end of interval (often +1)
752: Output Argument:
753: . q - A PetscQuadrature object
755: Level: intermediate
757: References:
758: . 1. - Karniadakis and Sherwin. FIAT
760: .seealso: PetscDTGaussTensorQuadrature(), PetscDTGaussQuadrature()
761: @*/
762: PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt dim, PetscInt Nc, PetscInt npoints, PetscReal a, PetscReal b, PetscQuadrature *q)
763: {
764: PetscInt totpoints = dim > 1 ? dim > 2 ? npoints*PetscSqr(npoints) : PetscSqr(npoints) : npoints;
765: PetscReal *px, *wx, *py, *wy, *pz, *wz, *x, *w;
766: PetscInt i, j, k, c;
770: if ((a != -1.0) || (b != 1.0)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must use default internal right now");
771: PetscMalloc1(totpoints*dim, &x);
772: PetscMalloc1(totpoints*Nc, &w);
773: switch (dim) {
774: case 0:
775: PetscFree(x);
776: PetscFree(w);
777: PetscMalloc1(1, &x);
778: PetscMalloc1(Nc, &w);
779: x[0] = 0.0;
780: for (c = 0; c < Nc; ++c) w[c] = 1.0;
781: break;
782: case 1:
783: PetscMalloc1(npoints,&wx);
784: PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, x, wx);
785: for (i = 0; i < npoints; ++i) for (c = 0; c < Nc; ++c) w[i*Nc+c] = wx[i];
786: PetscFree(wx);
787: break;
788: case 2:
789: PetscMalloc4(npoints,&px,npoints,&wx,npoints,&py,npoints,&wy);
790: PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, px, wx);
791: PetscDTGaussJacobiQuadrature1D_Internal(npoints, 1.0, 0.0, py, wy);
792: for (i = 0; i < npoints; ++i) {
793: for (j = 0; j < npoints; ++j) {
794: PetscDTMapSquareToTriangle_Internal(px[i], py[j], &x[(i*npoints+j)*2+0], &x[(i*npoints+j)*2+1]);
795: for (c = 0; c < Nc; ++c) w[(i*npoints+j)*Nc+c] = 0.5 * wx[i] * wy[j];
796: }
797: }
798: PetscFree4(px,wx,py,wy);
799: break;
800: case 3:
801: PetscMalloc6(npoints,&px,npoints,&wx,npoints,&py,npoints,&wy,npoints,&pz,npoints,&wz);
802: PetscDTGaussJacobiQuadrature1D_Internal(npoints, 0.0, 0.0, px, wx);
803: PetscDTGaussJacobiQuadrature1D_Internal(npoints, 1.0, 0.0, py, wy);
804: PetscDTGaussJacobiQuadrature1D_Internal(npoints, 2.0, 0.0, pz, wz);
805: for (i = 0; i < npoints; ++i) {
806: for (j = 0; j < npoints; ++j) {
807: for (k = 0; k < npoints; ++k) {
808: 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]);
809: for (c = 0; c < Nc; ++c) w[((i*npoints+j)*npoints+k)*Nc+c] = 0.125 * wx[i] * wy[j] * wz[k];
810: }
811: }
812: }
813: PetscFree6(px,wx,py,wy,pz,wz);
814: break;
815: default:
816: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot construct quadrature rule for dimension %d", dim);
817: }
818: PetscQuadratureCreate(PETSC_COMM_SELF, q);
819: PetscQuadratureSetOrder(*q, 2*npoints-1);
820: PetscQuadratureSetData(*q, dim, Nc, totpoints, x, w);
821: return(0);
822: }
824: /*@
825: PetscDTTanhSinhTensorQuadrature - create tanh-sinh quadrature for a tensor product cell
827: Not Collective
829: Input Arguments:
830: + dim - The cell dimension
831: . level - The number of points in one dimension, 2^l
832: . a - left end of interval (often-1)
833: - b - right end of interval (often +1)
835: Output Argument:
836: . q - A PetscQuadrature object
838: Level: intermediate
840: .seealso: PetscDTGaussTensorQuadrature()
841: @*/
842: PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt dim, PetscInt level, PetscReal a, PetscReal b, PetscQuadrature *q)
843: {
844: const PetscInt p = 16; /* Digits of precision in the evaluation */
845: const PetscReal alpha = (b-a)/2.; /* Half-width of the integration interval */
846: const PetscReal beta = (b+a)/2.; /* Center of the integration interval */
847: const PetscReal h = PetscPowReal(2.0, -level); /* Step size, length between x_k */
848: PetscReal xk; /* Quadrature point x_k on reference domain [-1, 1] */
849: PetscReal wk = 0.5*PETSC_PI; /* Quadrature weight at x_k */
850: PetscReal *x, *w;
851: PetscInt K, k, npoints;
852: PetscErrorCode ierr;
855: if (dim > 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Dimension %d not yet implemented", dim);
856: if (!level) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a number of significant digits");
857: /* Find K such that the weights are < 32 digits of precision */
858: for (K = 1; PetscAbsReal(PetscLog10Real(wk)) < 2*p; ++K) {
859: wk = 0.5*h*PETSC_PI*PetscCoshReal(K*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(K*h)));
860: }
861: PetscQuadratureCreate(PETSC_COMM_SELF, q);
862: PetscQuadratureSetOrder(*q, 2*K+1);
863: npoints = 2*K-1;
864: PetscMalloc1(npoints*dim, &x);
865: PetscMalloc1(npoints, &w);
866: /* Center term */
867: x[0] = beta;
868: w[0] = 0.5*alpha*PETSC_PI;
869: for (k = 1; k < K; ++k) {
870: wk = 0.5*alpha*h*PETSC_PI*PetscCoshReal(k*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h)));
871: xk = PetscTanhReal(0.5*PETSC_PI*PetscSinhReal(k*h));
872: x[2*k-1] = -alpha*xk+beta;
873: w[2*k-1] = wk;
874: x[2*k+0] = alpha*xk+beta;
875: w[2*k+0] = wk;
876: }
877: PetscQuadratureSetData(*q, dim, 1, npoints, x, w);
878: return(0);
879: }
881: PetscErrorCode PetscDTTanhSinhIntegrate(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol)
882: {
883: const PetscInt p = 16; /* Digits of precision in the evaluation */
884: const PetscReal alpha = (b-a)/2.; /* Half-width of the integration interval */
885: const PetscReal beta = (b+a)/2.; /* Center of the integration interval */
886: PetscReal h = 1.0; /* Step size, length between x_k */
887: PetscInt l = 0; /* Level of refinement, h = 2^{-l} */
888: PetscReal osum = 0.0; /* Integral on last level */
889: PetscReal psum = 0.0; /* Integral on the level before the last level */
890: PetscReal sum; /* Integral on current level */
891: PetscReal yk; /* Quadrature point 1 - x_k on reference domain [-1, 1] */
892: PetscReal lx, rx; /* Quadrature points to the left and right of 0 on the real domain [a, b] */
893: PetscReal wk; /* Quadrature weight at x_k */
894: PetscReal lval, rval; /* Terms in the quadature sum to the left and right of 0 */
895: PetscInt d; /* Digits of precision in the integral */
898: if (digits <= 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a positive number of significant digits");
899: /* Center term */
900: func(beta, &lval);
901: sum = 0.5*alpha*PETSC_PI*lval;
902: /* */
903: do {
904: PetscReal lterm, rterm, maxTerm = 0.0, d1, d2, d3, d4;
905: PetscInt k = 1;
907: ++l;
908: /* PetscPrintf(PETSC_COMM_SELF, "LEVEL %D sum: %15.15f\n", l, sum); */
909: /* At each level of refinement, h --> h/2 and sum --> sum/2 */
910: psum = osum;
911: osum = sum;
912: h *= 0.5;
913: sum *= 0.5;
914: do {
915: wk = 0.5*h*PETSC_PI*PetscCoshReal(k*h)/PetscSqr(PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h)));
916: yk = 1.0/(PetscExpReal(0.5*PETSC_PI*PetscSinhReal(k*h)) * PetscCoshReal(0.5*PETSC_PI*PetscSinhReal(k*h)));
917: lx = -alpha*(1.0 - yk)+beta;
918: rx = alpha*(1.0 - yk)+beta;
919: func(lx, &lval);
920: func(rx, &rval);
921: lterm = alpha*wk*lval;
922: maxTerm = PetscMax(PetscAbsReal(lterm), maxTerm);
923: sum += lterm;
924: rterm = alpha*wk*rval;
925: maxTerm = PetscMax(PetscAbsReal(rterm), maxTerm);
926: sum += rterm;
927: ++k;
928: /* Only need to evaluate every other point on refined levels */
929: if (l != 1) ++k;
930: } while (PetscAbsReal(PetscLog10Real(wk)) < p); /* Only need to evaluate sum until weights are < 32 digits of precision */
932: d1 = PetscLog10Real(PetscAbsReal(sum - osum));
933: d2 = PetscLog10Real(PetscAbsReal(sum - psum));
934: d3 = PetscLog10Real(maxTerm) - p;
935: if (PetscMax(PetscAbsReal(lterm), PetscAbsReal(rterm)) == 0.0) d4 = 0.0;
936: else d4 = PetscLog10Real(PetscMax(PetscAbsReal(lterm), PetscAbsReal(rterm)));
937: d = PetscAbsInt(PetscMin(0, PetscMax(PetscMax(PetscMax(PetscSqr(d1)/d2, 2*d1), d3), d4)));
938: } while (d < digits && l < 12);
939: *sol = sum;
941: return(0);
942: }
944: #if defined(PETSC_HAVE_MPFR)
945: PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol)
946: {
947: const PetscInt safetyFactor = 2; /* Calculate abcissa until 2*p digits */
948: PetscInt l = 0; /* Level of refinement, h = 2^{-l} */
949: mpfr_t alpha; /* Half-width of the integration interval */
950: mpfr_t beta; /* Center of the integration interval */
951: mpfr_t h; /* Step size, length between x_k */
952: mpfr_t osum; /* Integral on last level */
953: mpfr_t psum; /* Integral on the level before the last level */
954: mpfr_t sum; /* Integral on current level */
955: mpfr_t yk; /* Quadrature point 1 - x_k on reference domain [-1, 1] */
956: mpfr_t lx, rx; /* Quadrature points to the left and right of 0 on the real domain [a, b] */
957: mpfr_t wk; /* Quadrature weight at x_k */
958: PetscReal lval, rval; /* Terms in the quadature sum to the left and right of 0 */
959: PetscInt d; /* Digits of precision in the integral */
960: mpfr_t pi2, kh, msinh, mcosh, maxTerm, curTerm, tmp;
963: if (digits <= 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must give a positive number of significant digits");
964: /* Create high precision storage */
965: 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);
966: /* Initialization */
967: mpfr_set_d(alpha, 0.5*(b-a), MPFR_RNDN);
968: mpfr_set_d(beta, 0.5*(b+a), MPFR_RNDN);
969: mpfr_set_d(osum, 0.0, MPFR_RNDN);
970: mpfr_set_d(psum, 0.0, MPFR_RNDN);
971: mpfr_set_d(h, 1.0, MPFR_RNDN);
972: mpfr_const_pi(pi2, MPFR_RNDN);
973: mpfr_mul_d(pi2, pi2, 0.5, MPFR_RNDN);
974: /* Center term */
975: func(0.5*(b+a), &lval);
976: mpfr_set(sum, pi2, MPFR_RNDN);
977: mpfr_mul(sum, sum, alpha, MPFR_RNDN);
978: mpfr_mul_d(sum, sum, lval, MPFR_RNDN);
979: /* */
980: do {
981: PetscReal d1, d2, d3, d4;
982: PetscInt k = 1;
984: ++l;
985: mpfr_set_d(maxTerm, 0.0, MPFR_RNDN);
986: /* PetscPrintf(PETSC_COMM_SELF, "LEVEL %D sum: %15.15f\n", l, sum); */
987: /* At each level of refinement, h --> h/2 and sum --> sum/2 */
988: mpfr_set(psum, osum, MPFR_RNDN);
989: mpfr_set(osum, sum, MPFR_RNDN);
990: mpfr_mul_d(h, h, 0.5, MPFR_RNDN);
991: mpfr_mul_d(sum, sum, 0.5, MPFR_RNDN);
992: do {
993: mpfr_set_si(kh, k, MPFR_RNDN);
994: mpfr_mul(kh, kh, h, MPFR_RNDN);
995: /* Weight */
996: mpfr_set(wk, h, MPFR_RNDN);
997: mpfr_sinh_cosh(msinh, mcosh, kh, MPFR_RNDN);
998: mpfr_mul(msinh, msinh, pi2, MPFR_RNDN);
999: mpfr_mul(mcosh, mcosh, pi2, MPFR_RNDN);
1000: mpfr_cosh(tmp, msinh, MPFR_RNDN);
1001: mpfr_sqr(tmp, tmp, MPFR_RNDN);
1002: mpfr_mul(wk, wk, mcosh, MPFR_RNDN);
1003: mpfr_div(wk, wk, tmp, MPFR_RNDN);
1004: /* Abscissa */
1005: mpfr_set_d(yk, 1.0, MPFR_RNDZ);
1006: mpfr_cosh(tmp, msinh, MPFR_RNDN);
1007: mpfr_div(yk, yk, tmp, MPFR_RNDZ);
1008: mpfr_exp(tmp, msinh, MPFR_RNDN);
1009: mpfr_div(yk, yk, tmp, MPFR_RNDZ);
1010: /* Quadrature points */
1011: mpfr_sub_d(lx, yk, 1.0, MPFR_RNDZ);
1012: mpfr_mul(lx, lx, alpha, MPFR_RNDU);
1013: mpfr_add(lx, lx, beta, MPFR_RNDU);
1014: mpfr_d_sub(rx, 1.0, yk, MPFR_RNDZ);
1015: mpfr_mul(rx, rx, alpha, MPFR_RNDD);
1016: mpfr_add(rx, rx, beta, MPFR_RNDD);
1017: /* Evaluation */
1018: func(mpfr_get_d(lx, MPFR_RNDU), &lval);
1019: func(mpfr_get_d(rx, MPFR_RNDD), &rval);
1020: /* Update */
1021: mpfr_mul(tmp, wk, alpha, MPFR_RNDN);
1022: mpfr_mul_d(tmp, tmp, lval, MPFR_RNDN);
1023: mpfr_add(sum, sum, tmp, MPFR_RNDN);
1024: mpfr_abs(tmp, tmp, MPFR_RNDN);
1025: mpfr_max(maxTerm, maxTerm, tmp, MPFR_RNDN);
1026: mpfr_set(curTerm, tmp, MPFR_RNDN);
1027: mpfr_mul(tmp, wk, alpha, MPFR_RNDN);
1028: mpfr_mul_d(tmp, tmp, rval, MPFR_RNDN);
1029: mpfr_add(sum, sum, tmp, MPFR_RNDN);
1030: mpfr_abs(tmp, tmp, MPFR_RNDN);
1031: mpfr_max(maxTerm, maxTerm, tmp, MPFR_RNDN);
1032: mpfr_max(curTerm, curTerm, tmp, MPFR_RNDN);
1033: ++k;
1034: /* Only need to evaluate every other point on refined levels */
1035: if (l != 1) ++k;
1036: mpfr_log10(tmp, wk, MPFR_RNDN);
1037: mpfr_abs(tmp, tmp, MPFR_RNDN);
1038: } while (mpfr_get_d(tmp, MPFR_RNDN) < safetyFactor*digits); /* Only need to evaluate sum until weights are < 32 digits of precision */
1039: mpfr_sub(tmp, sum, osum, MPFR_RNDN);
1040: mpfr_abs(tmp, tmp, MPFR_RNDN);
1041: mpfr_log10(tmp, tmp, MPFR_RNDN);
1042: d1 = mpfr_get_d(tmp, MPFR_RNDN);
1043: mpfr_sub(tmp, sum, psum, MPFR_RNDN);
1044: mpfr_abs(tmp, tmp, MPFR_RNDN);
1045: mpfr_log10(tmp, tmp, MPFR_RNDN);
1046: d2 = mpfr_get_d(tmp, MPFR_RNDN);
1047: mpfr_log10(tmp, maxTerm, MPFR_RNDN);
1048: d3 = mpfr_get_d(tmp, MPFR_RNDN) - digits;
1049: mpfr_log10(tmp, curTerm, MPFR_RNDN);
1050: d4 = mpfr_get_d(tmp, MPFR_RNDN);
1051: d = PetscAbsInt(PetscMin(0, PetscMax(PetscMax(PetscMax(PetscSqr(d1)/d2, 2*d1), d3), d4)));
1052: } while (d < digits && l < 8);
1053: *sol = mpfr_get_d(sum, MPFR_RNDN);
1054: /* Cleanup */
1055: mpfr_clears(alpha, beta, h, sum, osum, psum, yk, wk, lx, rx, tmp, maxTerm, curTerm, pi2, kh, msinh, mcosh, NULL);
1056: return(0);
1057: }
1058: #else
1060: PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*func)(PetscReal, PetscReal *), PetscReal a, PetscReal b, PetscInt digits, PetscReal *sol)
1061: {
1062: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "This method will not work without MPFR. Reconfigure using --download-mpfr --download-gmp");
1063: }
1064: #endif
1066: /* Overwrites A. Can only handle full-rank problems with m>=n
1067: * A in column-major format
1068: * Ainv in row-major format
1069: * tau has length m
1070: * worksize must be >= max(1,n)
1071: */
1072: static PetscErrorCode PetscDTPseudoInverseQR(PetscInt m,PetscInt mstride,PetscInt n,PetscReal *A_in,PetscReal *Ainv_out,PetscScalar *tau,PetscInt worksize,PetscScalar *work)
1073: {
1075: PetscBLASInt M,N,K,lda,ldb,ldwork,info;
1076: PetscScalar *A,*Ainv,*R,*Q,Alpha;
1079: #if defined(PETSC_USE_COMPLEX)
1080: {
1081: PetscInt i,j;
1082: PetscMalloc2(m*n,&A,m*n,&Ainv);
1083: for (j=0; j<n; j++) {
1084: for (i=0; i<m; i++) A[i+m*j] = A_in[i+mstride*j];
1085: }
1086: mstride = m;
1087: }
1088: #else
1089: A = A_in;
1090: Ainv = Ainv_out;
1091: #endif
1093: PetscBLASIntCast(m,&M);
1094: PetscBLASIntCast(n,&N);
1095: PetscBLASIntCast(mstride,&lda);
1096: PetscBLASIntCast(worksize,&ldwork);
1097: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
1098: PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info));
1099: PetscFPTrapPop();
1100: if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error");
1101: R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */
1103: /* Extract an explicit representation of Q */
1104: Q = Ainv;
1105: PetscMemcpy(Q,A,mstride*n*sizeof(PetscScalar));
1106: K = N; /* full rank */
1107: PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info));
1108: if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error");
1110: /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */
1111: Alpha = 1.0;
1112: ldb = lda;
1113: PetscStackCallBLAS("BLAStrsm",BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb));
1114: /* Ainv is Q, overwritten with inverse */
1116: #if defined(PETSC_USE_COMPLEX)
1117: {
1118: PetscInt i;
1119: for (i=0; i<m*n; i++) Ainv_out[i] = PetscRealPart(Ainv[i]);
1120: PetscFree2(A,Ainv);
1121: }
1122: #endif
1123: return(0);
1124: }
1126: /* Computes integral of L_p' over intervals {(x0,x1),(x1,x2),...} */
1127: static PetscErrorCode PetscDTLegendreIntegrate(PetscInt ninterval,const PetscReal *x,PetscInt ndegree,const PetscInt *degrees,PetscBool Transpose,PetscReal *B)
1128: {
1130: PetscReal *Bv;
1131: PetscInt i,j;
1134: PetscMalloc1((ninterval+1)*ndegree,&Bv);
1135: /* Point evaluation of L_p on all the source vertices */
1136: PetscDTLegendreEval(ninterval+1,x,ndegree,degrees,Bv,NULL,NULL);
1137: /* Integral over each interval: \int_a^b L_p' = L_p(b)-L_p(a) */
1138: for (i=0; i<ninterval; i++) {
1139: for (j=0; j<ndegree; j++) {
1140: if (Transpose) B[i+ninterval*j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j];
1141: else B[i*ndegree+j] = Bv[(i+1)*ndegree+j] - Bv[i*ndegree+j];
1142: }
1143: }
1144: PetscFree(Bv);
1145: return(0);
1146: }
1148: /*@
1149: PetscDTReconstructPoly - create matrix representing polynomial reconstruction using cell intervals and evaluation at target intervals
1151: Not Collective
1153: Input Arguments:
1154: + degree - degree of reconstruction polynomial
1155: . nsource - number of source intervals
1156: . sourcex - sorted coordinates of source cell boundaries (length nsource+1)
1157: . ntarget - number of target intervals
1158: - targetx - sorted coordinates of target cell boundaries (length ntarget+1)
1160: Output Arguments:
1161: . R - reconstruction matrix, utarget = sum_s R[t*nsource+s] * usource[s]
1163: Level: advanced
1165: .seealso: PetscDTLegendreEval()
1166: @*/
1167: PetscErrorCode PetscDTReconstructPoly(PetscInt degree,PetscInt nsource,const PetscReal *sourcex,PetscInt ntarget,const PetscReal *targetx,PetscReal *R)
1168: {
1170: PetscInt i,j,k,*bdegrees,worksize;
1171: PetscReal xmin,xmax,center,hscale,*sourcey,*targety,*Bsource,*Bsinv,*Btarget;
1172: PetscScalar *tau,*work;
1178: 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);
1179: #if defined(PETSC_USE_DEBUG)
1180: for (i=0; i<nsource; i++) {
1181: 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]);
1182: }
1183: for (i=0; i<ntarget; i++) {
1184: 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]);
1185: }
1186: #endif
1187: xmin = PetscMin(sourcex[0],targetx[0]);
1188: xmax = PetscMax(sourcex[nsource],targetx[ntarget]);
1189: center = (xmin + xmax)/2;
1190: hscale = (xmax - xmin)/2;
1191: worksize = nsource;
1192: PetscMalloc4(degree+1,&bdegrees,nsource+1,&sourcey,nsource*(degree+1),&Bsource,worksize,&work);
1193: PetscMalloc4(nsource,&tau,nsource*(degree+1),&Bsinv,ntarget+1,&targety,ntarget*(degree+1),&Btarget);
1194: for (i=0; i<=nsource; i++) sourcey[i] = (sourcex[i]-center)/hscale;
1195: for (i=0; i<=degree; i++) bdegrees[i] = i+1;
1196: PetscDTLegendreIntegrate(nsource,sourcey,degree+1,bdegrees,PETSC_TRUE,Bsource);
1197: PetscDTPseudoInverseQR(nsource,nsource,degree+1,Bsource,Bsinv,tau,nsource,work);
1198: for (i=0; i<=ntarget; i++) targety[i] = (targetx[i]-center)/hscale;
1199: PetscDTLegendreIntegrate(ntarget,targety,degree+1,bdegrees,PETSC_FALSE,Btarget);
1200: for (i=0; i<ntarget; i++) {
1201: PetscReal rowsum = 0;
1202: for (j=0; j<nsource; j++) {
1203: PetscReal sum = 0;
1204: for (k=0; k<degree+1; k++) {
1205: sum += Btarget[i*(degree+1)+k] * Bsinv[k*nsource+j];
1206: }
1207: R[i*nsource+j] = sum;
1208: rowsum += sum;
1209: }
1210: for (j=0; j<nsource; j++) R[i*nsource+j] /= rowsum; /* normalize each row */
1211: }
1212: PetscFree4(bdegrees,sourcey,Bsource,work);
1213: PetscFree4(tau,Bsinv,targety,Btarget);
1214: return(0);
1215: }